faceMapper.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "faceMapper.H"
30 #include "polyMesh.H"
31 #include "mapPolyMesh.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::faceMapper::calcAddressing() const
36 {
37  if
38  (
39  directAddrPtr_
40  || interpAddrPtr_
41  || weightsPtr_
42  || insertedObjectsPtr_
43  )
44  {
46  << "Addressing already calculated."
47  << abort(FatalError);
48  }
49 
50  if (direct())
51  {
52  // Direct addressing, no weights
53 
54  // Restrict addressing list to contain only live faces
55  directAddrPtr_ = std::make_unique<labelList>
56  (
57  labelList::subList(mpm_.faceMap(), mapperLen_)
58  );
59  auto& directAddr = *directAddrPtr_;
60 
61  insertedObjectsPtr_ = std::make_unique<labelList>();
62  auto& inserted = *insertedObjectsPtr_;
63 
64  // The nInsertedObjects_ already counted in the constructor
65  if (nInsertedObjects_)
66  {
67  inserted.resize(nInsertedObjects_);
68 
69  label nInserted = 0;
70  forAll(directAddr, i)
71  {
72  if (directAddr[i] < 0)
73  {
74  // Found inserted
75  directAddr[i] = 0;
76  inserted[nInserted] = i;
77  ++nInserted;
78 
79  // TBD: check (nInsertedObjects_ < nInserted)?
80  #ifdef FULLDEBUG
81  if (nInsertedObjects_ < nInserted)
82  {
84  << "Unexpected insert of more than "
85  << nInsertedObjects_ << " items\n"
86  << abort(FatalError);
87  }
88  #endif
89  }
90  }
91  // TBD: check (nInserted < nInsertedObjects_)?
92  #ifdef FULLDEBUG
93  if (nInserted < nInsertedObjects_)
94  {
96  << "Found " << nInserted << " instead of "
97  << nInsertedObjects_ << " items to insert\n";
98  }
99  #endif
100  // The resize should be unnecessary
101  inserted.resize(nInserted);
102  }
103  }
104  else
105  {
106  // Interpolative addressing
107 
108  interpAddrPtr_ = std::make_unique<labelListList>(mapperLen_);
109  auto& addr = *interpAddrPtr_;
110 
111  weightsPtr_ = std::make_unique<scalarListList>(mapperLen_);
112  auto& wght = *weightsPtr_;
113 
114 
115  // Set the addressing and uniform weight
116  const auto setAddrWeights = [&]
117  (
118  const List<objectMap>& maps,
119  const char * const nameOfMap
120  )
121  {
122  for (const objectMap& map : maps)
123  {
124  // Get index, addressing
125  const label facei = map.index();
126  const labelList& mo = map.masterObjects();
127  if (mo.empty()) continue; // safety
128 
129  if (addr[facei].size())
130  {
132  << "Master face " << facei
133  << " already mapped, cannot apply "
134  << nameOfMap
135  << flatOutput(mo) << abort(FatalError);
136  }
137 
138  // Map from masters, uniform weights
139  addr[facei] = mo;
140  wght[facei] = scalarList(mo.size(), 1.0/mo.size());
141  }
142  };
143 
144 
145  setAddrWeights(mpm_.facesFromPointsMap(), "point faces");
146  setAddrWeights(mpm_.facesFromEdgesMap(), "edge faces");
147  setAddrWeights(mpm_.facesFromFacesMap(), "face faces");
148 
149 
150  // Do mapped faces.
151  // - may already have been set, so check if addressing still empty().
152 
153  {
154  const labelList& map = mpm_.faceMap();
155 
156  // NB: faceMap can be longer than nFaces()
157  for (label facei = 0; facei < mapperLen_; ++facei)
158  {
159  const label mappedi = map[facei];
160 
161  if (mappedi >= 0 && addr[facei].empty())
162  {
163  // Mapped from a single face
164  addr[facei].resize(1, mappedi);
165  wght[facei].resize(1, 1.0);
166  }
167  }
168  }
169 
170 
171  // Grab inserted faces (for them the size of addressing is still zero)
172 
173  insertedObjectsPtr_ = std::make_unique<labelList>();
174  auto& inserted = *insertedObjectsPtr_;
175 
176  // The nInsertedObjects_ already counted in the constructor
177  if (nInsertedObjects_)
178  {
179  inserted.resize(nInsertedObjects_);
180 
181  label nInserted = 0;
182  forAll(addr, i)
183  {
184  if (addr[i].empty())
185  {
186  // Mapped from dummy face 0
187  addr[i].resize(1, 0);
188  wght[i].resize(1, 1.0);
189 
190  inserted[nInserted] = i;
191  ++nInserted;
192 
193  // TBD: check (nInsertedObjects_ < nInserted)?
194  #ifdef FULLDEBUG
195  if (nInsertedObjects_ < nInserted)
196  {
198  << "Unexpected insert of more than "
199  << nInsertedObjects_ << " items\n"
200  << abort(FatalError);
201  }
202  #endif
203  }
204  }
205  // TBD: check (nInserted < nInsertedObjects_)?
206  #ifdef FULLDEBUG
207  if (nInserted < nInsertedObjects_)
208  {
210  << "Found " << nInserted << " instead of "
211  << nInsertedObjects_ << " items to insert\n";
212  }
213  #endif
214  // The resize should be unnecessary
215  inserted.resize(nInserted);
216  }
217  }
218 }
219 
220 
221 // void Foam::faceMapper::clearOut()
222 // {
223 // directAddrPtr_.reset(nullptr);
224 // interpAddrPtr_.reset(nullptr);
225 // weightsPtr_.reset(nullptr);
226 // insertedObjectsPtr_.reset(nullptr);
227 // }
228 
229 
230 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
231 
233 :
234  mpm_(mpm),
235  mapperLen_(mpm.mesh().nFaces()),
236  nInsertedObjects_(0),
237  direct_
238  (
239  // Mapping without interpolation?
240  mpm.facesFromPointsMap().empty()
241  && mpm.facesFromEdgesMap().empty()
242  && mpm.facesFromFacesMap().empty()
243  )
244 {
245  const auto& directMap = mpm_.faceMap();
246 
247  if (!mapperLen_)
248  {
249  // Empty mesh
250  direct_ = true;
251  nInsertedObjects_ = 0;
252  }
253  else if (direct_)
254  {
255  // Number of inserted faces (-ve values)
256  nInsertedObjects_ = std::count_if
257  (
258  directMap.cbegin(),
259  directMap.cbegin(mapperLen_),
260  [](label i) { return (i < 0); }
261  );
262  }
263  else
264  {
265  // Check if there are inserted faces with no owner
266  // (check all lists)
267 
268  bitSet unmapped(mapperLen_, true);
269 
270  unmapped.unset(directMap); // direct mapped
271 
272  for (const auto& map : mpm_.facesFromPointsMap())
273  {
274  if (!map.empty()) unmapped.unset(map.index());
275  }
276 
277  for (const auto& map : mpm_.facesFromEdgesMap())
278  {
279  if (!map.empty()) unmapped.unset(map.index());
280  }
281 
282  for (const auto& map : mpm_.facesFromFacesMap())
283  {
284  if (!map.empty()) unmapped.unset(map.index());
285  }
286 
287  nInsertedObjects_ = label(unmapped.count());
288  }
289 }
290 
291 
292 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
295 {}
296 
297 
298 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
300 Foam::label Foam::faceMapper::size() const
301 {
302  return mapperLen_;
303 }
304 
306 Foam::label Foam::faceMapper::sizeBeforeMapping() const
307 {
308  return mpm_.nOldFaces();
309 }
310 
313 {
314  return mpm_.nOldInternalFaces();
315 }
316 
317 
319 {
320  if (!direct())
321  {
323  << "Requested direct addressing for an interpolative mapper."
324  << abort(FatalError);
325  }
326 
327  if (!insertedObjects())
328  {
329  // No inserted faces. Re-use faceMap
330  return mpm_.faceMap();
331  }
332  else
333  {
334  if (!directAddrPtr_)
335  {
336  calcAddressing();
337  }
338 
339  return *directAddrPtr_;
340  }
341 }
342 
343 
345 {
346  if (direct())
347  {
349  << "Requested interpolative addressing for a direct mapper."
350  << abort(FatalError);
351  }
352 
353  if (!interpAddrPtr_)
354  {
355  calcAddressing();
356  }
357 
358  return *interpAddrPtr_;
359 }
360 
361 
363 {
364  if (direct())
365  {
367  << "Requested interpolative weights for a direct mapper."
368  << abort(FatalError);
369  }
370 
371  if (!weightsPtr_)
372  {
373  calcAddressing();
374  }
375 
376  return *weightsPtr_;
377 }
378 
379 
381 {
382  if (!insertedObjectsPtr_)
383  {
384  if (!nInsertedObjects_)
385  {
386  // No inserted objects will be created
387  return labelList::null();
388  }
389 
390  calcAddressing();
391  }
392 
393  return *insertedObjectsPtr_;
394 }
395 
398 {
399  return mpm_.flipFaceFlux();
400 }
401 
403 Foam::label Foam::faceMapper::nOldInternalFaces() const
404 {
405  return mpm_.nOldInternalFaces();
406 }
407 
410 {
411  return mpm_.oldPatchStarts();
412 }
413 
414 
416 {
417  return mpm_.oldPatchSizes();
418 }
419 
420 
421 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:420
label count_if(const ListType &input, const UnaryPredicate &pred, const label start=0)
Count the number of matching entries.
virtual label internalSizeBeforeMapping() const
Return number of internal faces before mapping.
Definition: faceMapper.C:305
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:153
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
Definition: faceMapper.C:373
SubList< label > subList
Declare type of subList.
Definition: List.H:144
const List< objectMap > & facesFromPointsMap() const noexcept
Faces inflated from points.
Definition: mapPolyMesh.H:511
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
virtual label sizeBeforeMapping() const
Return size of field before mapping.
Definition: faceMapper.C:299
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
faceMapper(const faceMapper &)=delete
No copy construct.
dynamicFvMesh & mesh
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faceMapper.C:337
virtual const labelHashSet & flipFaceFlux() const
Return flux flip map.
Definition: faceMapper.C:390
bitSet & unset(const bitSet &other)
Unset (subtract) the bits specified in the other bitset, which is a set difference corresponds to the...
Definition: bitSetI.H:536
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: faceMapper.C:355
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const List< objectMap > & facesFromEdgesMap() const noexcept
Faces inflated from edges.
Definition: mapPolyMesh.H:519
virtual bool direct() const
Is the mapping direct.
Definition: faceMapper.H:159
virtual ~faceMapper()
Destructor.
Definition: faceMapper.C:287
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
virtual label size() const
The mapper size.
Definition: faceMapper.C:293
virtual const labelList & oldPatchStarts() const
Return old patch starts.
Definition: faceMapper.C:402
const labelList & faceMap() const noexcept
Old face map.
Definition: mapPolyMesh.H:503
virtual label nOldInternalFaces() const
Return number of old internalFaces.
Definition: faceMapper.C:396
List< label > labelList
A List of labels.
Definition: List.H:62
const List< objectMap > & facesFromFacesMap() const noexcept
Faces originating from faces.
Definition: mapPolyMesh.H:527
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faceMapper.C:311
static const List< label > & null() noexcept
Return a null List (reference to a nullObject). Behaves like an empty List.
Definition: List.H:153
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
virtual const labelList & oldPatchSizes() const
Return old patch sizes.
Definition: faceMapper.C:408