cellMapper.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 "cellMapper.H"
30 #include "polyMesh.H"
31 #include "mapPolyMesh.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::cellMapper::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  directAddrPtr_ = std::make_unique<labelList>
55  (
56  // No retired cells, so cellMap().size() == mapperLen_ anyhow
57  labelList::subList(mpm_.cellMap(), 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 celli = map.index();
126  const labelList& mo = map.masterObjects();
127  if (mo.empty()) continue; // safety
128 
129  if (addr[celli].size())
130  {
132  << "Master cell " << celli
133  << " already mapped, cannot apply "
134  << nameOfMap
135  << flatOutput(mo) << abort(FatalError);
136  }
137 
138  // Map from masters, uniform weights
139  addr[celli] = mo;
140  wght[celli] = scalarList(mo.size(), 1.0/mo.size());
141  }
142  };
143 
144 
145  setAddrWeights(mpm_.cellsFromPointsMap(), "point cells");
146  setAddrWeights(mpm_.cellsFromEdgesMap(), "edge cells");
147  setAddrWeights(mpm_.cellsFromFacesMap(), "face cells");
148 
149  // Volume conservative mapping if possible
150 
151  const List<objectMap>& cellsFromCells = mpm_.cellsFromCellsMap();
152  setAddrWeights(cellsFromCells, "cell cells");
153 
154  if (mpm_.hasOldCellVolumes())
155  {
156  // Volume weighted
157 
158  const scalarField& V = mpm_.oldCellVolumes();
159 
160  if (V.size() != sizeBeforeMapping())
161  {
163  << "cellVolumes size " << V.size()
164  << " != old number of cells " << sizeBeforeMapping()
165  << ". Are your cellVolumes already mapped?"
166  << " (new number of cells " << size() << ")"
167  << abort(FatalError);
168  }
169 
170  for (const auto& map : cellsFromCells)
171  {
172  // Get index, addressing
173  const label celli = map.index();
174  const labelList& mo = map.masterObjects();
175  if (mo.empty()) continue; // safety
176 
177  // wght[celli] is already sized and uniform weighted
178  auto& wght_cell = wght[celli];
179 
180  scalar sumV = 0;
181  forAll(mo, ci)
182  {
183  wght_cell[ci] = V[mo[ci]];
184  sumV += V[mo[ci]];
185  }
186  if (sumV > VSMALL)
187  {
188  for (auto& w : wght_cell)
189  {
190  w /= sumV;
191  }
192  }
193  else
194  {
195  // Exception: zero volume. Use uniform mapping
196  wght_cell = (1.0/mo.size());
197  }
198  }
199  }
200 
201 
202  // Do mapped cells.
203  // - may already have been set, so check if addressing still empty().
204 
205  {
206  const labelList& map = mpm_.cellMap();
207 
208  // The cellMap.size() == nCells() anyhow
209  for (label celli = 0; celli < mapperLen_; ++celli)
210  {
211  const label mappedi = map[celli];
212 
213  if (mappedi >= 0 && addr[celli].empty())
214  {
215  // Mapped from a single cell
216  addr[celli].resize(1, mappedi);
217  wght[celli].resize(1, 1.0);
218  }
219  }
220  }
221 
222 
223  // Grab inserted points (for them the size of addressing is still zero)
224 
225  insertedObjectsPtr_ = std::make_unique<labelList>();
226  auto& inserted = *insertedObjectsPtr_;
227 
228  // The nInsertedObjects_ already counted in the constructor
229  if (nInsertedObjects_)
230  {
231  inserted.resize(nInsertedObjects_);
232 
233  label nInserted = 0;
234  forAll(addr, i)
235  {
236  if (addr[i].empty())
237  {
238  // Mapped from dummy cell 0
239  addr[i].resize(1, 0);
240  wght[i].resize(1, 1.0);
241 
242  inserted[nInserted] = i;
243  ++nInserted;
244 
245  // TBD: check (nInsertedObjects_ < nInserted)?
246  #ifdef FULLDEBUG
247  if (nInsertedObjects_ < nInserted)
248  {
250  << "Unexpected insert of more than "
251  << nInsertedObjects_ << " items\n"
252  << abort(FatalError);
253  }
254  #endif
255  }
256  }
257  // TBD: check (nInserted < nInsertedObjects_)?
258  #ifdef FULLDEBUG
259  if (nInserted < nInsertedObjects_)
260  {
262  << "Found " << nInserted << " instead of "
263  << nInsertedObjects_ << " items to insert\n";
264  }
265  #endif
266  // The resize should be unnecessary
267  inserted.resize(nInserted);
268  }
269  }
270 }
271 
272 
273 // void Foam::cellMapper::clearOut()
274 // {
275 // directAddrPtr_.reset(nullptr);
276 // interpAddrPtr_.reset(nullptr);
277 // weightsPtr_.reset(nullptr);
278 // insertedObjectsPtr_.reset(nullptr);
279 // }
280 
281 
282 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
283 
285 :
286  mpm_(mpm),
287  mapperLen_(mpm.mesh().nCells()),
288  nInsertedObjects_(0),
289  direct_
290  (
291  // Mapping without interpolation?
292  mpm.cellsFromPointsMap().empty()
293  && mpm.cellsFromEdgesMap().empty()
294  && mpm.cellsFromFacesMap().empty()
295  && mpm.cellsFromCellsMap().empty()
296  )
297 {
298  const auto& directMap = mpm_.cellMap();
299 
300  if (!mapperLen_)
301  {
302  // Empty mesh
303  direct_ = true;
304  nInsertedObjects_ = 0;
305  }
306  else if (direct_)
307  {
308  // Number of inserted cells (-ve values)
309  nInsertedObjects_ = std::count_if
310  (
311  directMap.cbegin(),
312  directMap.cbegin(mapperLen_),
313  [](label i) { return (i < 0); }
314  );
315  }
316  else
317  {
318  // Check if there are inserted cells with no owner
319  // (check all lists)
320 
321  bitSet unmapped(mapperLen_, true);
322 
323  unmapped.unset(directMap); // direct mapped
324 
325  for (const auto& map : mpm_.cellsFromPointsMap())
326  {
327  if (!map.empty()) unmapped.unset(map.index());
328  }
329 
330  for (const auto& map : mpm_.cellsFromEdgesMap())
331  {
332  if (!map.empty()) unmapped.unset(map.index());
333  }
334 
335  for (const auto& map : mpm_.cellsFromFacesMap())
336  {
337  if (!map.empty()) unmapped.unset(map.index());
338  }
339 
340  for (const auto& map : mpm_.cellsFromCellsMap())
341  {
342  if (!map.empty()) unmapped.unset(map.index());
343  }
344 
345  nInsertedObjects_ = label(unmapped.count());
346  }
347 }
348 
349 
350 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
353 {}
354 
355 
356 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
357 
358 Foam::label Foam::cellMapper::size() const
359 {
360  // OR: return mapperLen_;
361  return mpm_.cellMap().size();
362 }
363 
365 Foam::label Foam::cellMapper::sizeBeforeMapping() const
366 {
367  return mpm_.nOldCells();
368 }
369 
370 
372 {
373  if (!direct())
374  {
376  << "Requested direct addressing for an interpolative mapper."
377  << abort(FatalError);
378  }
379 
380  if (!insertedObjects())
381  {
382  // No inserted cells. Re-use cellMap
383  return mpm_.cellMap();
384  }
385  else
386  {
387  if (!directAddrPtr_)
388  {
389  calcAddressing();
390  }
391 
392  return *directAddrPtr_;
393  }
394 }
395 
396 
398 {
399  if (direct())
400  {
402  << "Requested interpolative addressing for a direct mapper."
403  << abort(FatalError);
404  }
405 
406  if (!interpAddrPtr_)
407  {
408  calcAddressing();
409  }
410 
411  return *interpAddrPtr_;
412 }
413 
414 
416 {
417  if (direct())
418  {
420  << "Requested interpolative weights for a direct mapper."
421  << abort(FatalError);
422  }
423 
424  if (!weightsPtr_)
425  {
426  calcAddressing();
427  }
428 
429  return *weightsPtr_;
430 }
431 
432 
434 {
435  if (!insertedObjectsPtr_)
436  {
437  if (!nInsertedObjects_)
438  {
439  // No inserted objects will be created
440  return labelList::null();
441  }
442 
443  calcAddressing();
444  }
445 
446  return *insertedObjectsPtr_;
447 }
448 
449 
450 // ************************************************************************* //
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.
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...
virtual label sizeBeforeMapping() const
Return size before mapping.
Definition: cellMapper.C:358
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: cellMapper.C:408
bool hasOldCellVolumes() const noexcept
Definition: mapPolyMesh.H:800
cellMapper(const cellMapper &)=delete
No copy construct.
SubList< label > subList
Declare type of subList.
Definition: List.H:144
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
const List< objectMap > & cellsFromCellsMap() const noexcept
Cells originating from cells.
Definition: mapPolyMesh.H:569
const labelList & cellMap() const noexcept
Old cell map.
Definition: mapPolyMesh.H:537
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const List< objectMap > & cellsFromPointsMap() const noexcept
Cells inflated from points.
Definition: mapPolyMesh.H:545
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: cellMapper.C:390
errorManip< error > abort(error &err)
Definition: errorManip.H:139
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: cellMapper.C:364
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.
const scalarField & oldCellVolumes() const
Definition: mapPolyMesh.H:805
const List< objectMap > & cellsFromEdgesMap() const noexcept
Cells inflated from edges.
Definition: mapPolyMesh.H:553
const List< objectMap > & cellsFromFacesMap() const noexcept
Cells inflated from faces.
Definition: mapPolyMesh.H:561
List< label > labelList
A List of labels.
Definition: List.H:62
virtual label size() const
The mapper size.
Definition: cellMapper.C:351
virtual const labelList & insertedObjectLabels() const
Return list of inserted cells.
Definition: cellMapper.C:426
virtual ~cellMapper()
Destructor.
Definition: cellMapper.C:345
static const List< label > & null() noexcept
Return a null List (reference to a nullObject). Behaves like an empty List.
Definition: List.H:153
virtual bool direct() const
Is the mapping direct.
Definition: cellMapper.H:153
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225