cellVolumeWeightMethod.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "cellVolumeWeightMethod.H"
27 #include "faceAreaWeightAMI.H"
28 #include "indexedOctree.H"
29 #include "treeDataCell.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(cellVolumeWeightMethod, 0);
38  (
39  meshToMeshMethod,
40  cellVolumeWeightMethod,
41  components
42  );
43 }
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 (
49  const labelList& srcCellIDs,
50  const boolList& mapFlag,
51  const label startSeedI,
52  label& srcSeedI,
53  label& tgtSeedI
54 ) const
55 {
56  const cellList& srcCells = src_.cells();
57  const faceList& srcFaces = src_.faces();
58  const pointField& srcPts = src_.points();
59 
60  for (label i = startSeedI; i < srcCellIDs.size(); i++)
61  {
62  const label srcI = srcCellIDs[i];
63 
64  if (mapFlag[srcI])
65  {
66  const labelList tgtIDs
67  (
68  tgt_.cellTree().findBox
69  (
70  treeBoundBox(srcCells[srcI].bb(srcPts, srcFaces))
71  )
72  );
73 
74  forAll(tgtIDs, j)
75  {
76  const label tgtI = tgtIDs[j];
77 
78  if (intersect(srcI, tgtI))
79  {
80  srcSeedI = srcI;
81  tgtSeedI = tgtI;
82 
83  return true;
84  }
85  }
86  }
87  }
88 
89  if (debug)
90  {
91  Pout<< "could not find starting seed" << endl;
92  }
93 
94  return false;
95 }
96 
97 
99 (
100  labelListList& srcToTgtCellAddr,
101  scalarListList& srcToTgtCellWght,
102  labelListList& tgtToSrcCellAddr,
103  scalarListList& tgtToSrcCellWght,
104  const label srcSeedI,
105  const label tgtSeedI,
106  const labelList& srcCellIDs,
107  boolList& mapFlag,
108  label& startSeedI
109 )
110 {
111  label srcCelli = srcSeedI;
112  label tgtCelli = tgtSeedI;
113 
114  List<DynamicList<label>> srcToTgtAddr(src_.nCells());
115  List<DynamicList<scalar>> srcToTgtWght(src_.nCells());
116 
117  List<DynamicList<label>> tgtToSrcAddr(tgt_.nCells());
118  List<DynamicList<scalar>> tgtToSrcWght(tgt_.nCells());
119 
120  // list of tgt cell neighbour cells
121  DynamicList<label> nbrTgtCells(10);
122 
123  // list of tgt cells currently visited for srcCelli to avoid multiple hits
124  DynamicList<label> visitedTgtCells(10);
125 
126  // list to keep track of tgt cells used to seed src cells
127  labelList seedCells(src_.nCells(), -1);
128  seedCells[srcCelli] = tgtCelli;
129 
130  const scalarField& srcVol = src_.cellVolumes();
131 
132  do
133  {
134  nbrTgtCells.clear();
135  visitedTgtCells.clear();
136 
137  // append initial target cell and neighbours
138  nbrTgtCells.append(tgtCelli);
139  appendNbrCells(tgtCelli, tgt_, visitedTgtCells, nbrTgtCells);
140 
141  do
142  {
143  tgtCelli = nbrTgtCells.remove();
144  visitedTgtCells.append(tgtCelli);
145 
146  scalar vol = interVol(srcCelli, tgtCelli);
147 
148  // accumulate addressing and weights for valid intersection
149  if (vol/srcVol[srcCelli] > tolerance_)
150  {
151  // store src/tgt cell pair
152  srcToTgtAddr[srcCelli].append(tgtCelli);
153  srcToTgtWght[srcCelli].append(vol);
154 
155  tgtToSrcAddr[tgtCelli].append(srcCelli);
156  tgtToSrcWght[tgtCelli].append(vol);
157 
158  appendNbrCells(tgtCelli, tgt_, visitedTgtCells, nbrTgtCells);
159 
160  // accumulate intersection volume
161  V_ += vol;
162  }
163  }
164  while (!nbrTgtCells.empty());
165 
166  mapFlag[srcCelli] = false;
167 
168  // find new source seed cell
169  setNextCells
170  (
171  startSeedI,
172  srcCelli,
173  tgtCelli,
174  srcCellIDs,
175  mapFlag,
176  visitedTgtCells,
177  seedCells
178  );
179  }
180  while (srcCelli != -1);
181 
182  // transfer addressing into persistent storage
183  forAll(srcToTgtCellAddr, i)
184  {
185  srcToTgtCellAddr[i].transfer(srcToTgtAddr[i]);
186  srcToTgtCellWght[i].transfer(srcToTgtWght[i]);
187  }
188 
189  forAll(tgtToSrcCellAddr, i)
190  {
191  tgtToSrcCellAddr[i].transfer(tgtToSrcAddr[i]);
192  tgtToSrcCellWght[i].transfer(tgtToSrcWght[i]);
193  }
194 }
195 
196 
198 (
199  label& startSeedI,
200  label& srcCelli,
201  label& tgtCelli,
202  const labelList& srcCellIDs,
203  const boolList& mapFlag,
204  const DynamicList<label>& visitedCells,
205  labelList& seedCells
206 ) const
207 {
208  const labelList& srcNbrCells = src_.cellCells()[srcCelli];
209 
210  // set possible seeds for later use by querying all src cell neighbours
211  // with all visited target cells
212  bool valuesSet = false;
213  forAll(srcNbrCells, i)
214  {
215  label cellS = srcNbrCells[i];
216 
217  if (mapFlag[cellS] && seedCells[cellS] == -1)
218  {
219  forAll(visitedCells, j)
220  {
221  label cellT = visitedCells[j];
222 
223  if (intersect(cellS, cellT))
224  {
225  seedCells[cellS] = cellT;
226 
227  if (!valuesSet)
228  {
229  srcCelli = cellS;
230  tgtCelli = cellT;
231  valuesSet = true;
232  }
233  }
234  }
235  }
236  }
237 
238  // set next src and tgt cells if not set above
239  if (valuesSet)
240  {
241  return;
242  }
243  else
244  {
245  // try to use existing seed
246  bool foundNextSeed = false;
247  for (label i = startSeedI; i < srcCellIDs.size(); i++)
248  {
249  label cellS = srcCellIDs[i];
250 
251  if (mapFlag[cellS])
252  {
253  if (!foundNextSeed)
254  {
255  startSeedI = i;
256  foundNextSeed = true;
257  }
258 
259  if (seedCells[cellS] != -1)
260  {
261  srcCelli = cellS;
262  tgtCelli = seedCells[cellS];
263 
264  return;
265  }
266  }
267  }
268 
269  // perform new search to find match
270  if (debug)
271  {
272  Pout<< "Advancing front stalled: searching for new "
273  << "target cell" << endl;
274  }
275 
276  bool restart =
277  findInitialSeeds
278  (
279  srcCellIDs,
280  mapFlag,
281  startSeedI,
282  srcCelli,
283  tgtCelli
284  );
285 
286  if (restart)
287  {
288  // successfully found new starting seed-pair
289  return;
290  }
291  }
292 
293  // if we have got to here, there are no more src/tgt cell intersections
294  srcCelli = -1;
295  tgtCelli = -1;
296 }
297 
298 
299 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
300 
302 (
303  const polyMesh& src,
304  const polyMesh& tgt
305 )
306 :
307  meshToMeshMethod(src, tgt)
308 {}
309 
310 
311 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
312 
314 {}
315 
316 
317 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
318 
320 {
321  return faceAreaWeightAMI::typeName;
322 }
323 
324 
326 (
327  labelListList& srcToTgtAddr,
328  scalarListList& srcToTgtWght,
329  labelListList& tgtToSrcAddr,
330  scalarListList& tgtToSrcWght
331 )
332 {
333  bool ok = initialise
334  (
335  srcToTgtAddr,
336  srcToTgtWght,
337  tgtToSrcAddr,
338  tgtToSrcWght
339  );
340 
341  if (!ok)
342  {
343  return;
344  }
345 
346  // (potentially) participating source mesh cells
347  const labelList srcCellIDs(maskCells());
348 
349  // list to keep track of whether src cell can be mapped
350  boolList mapFlag(src_.nCells(), false);
351  UIndirectList<bool>(mapFlag, srcCellIDs) = true;
352 
353  // find initial point in tgt mesh
354  label srcSeedI = -1;
355  label tgtSeedI = -1;
356  label startSeedI = 0;
357 
358  bool startWalk =
359  findInitialSeeds
360  (
361  srcCellIDs,
362  mapFlag,
363  startSeedI,
364  srcSeedI,
365  tgtSeedI
366  );
367 
368  if (startWalk)
369  {
370  calculateAddressing
371  (
372  srcToTgtAddr,
373  srcToTgtWght,
374  tgtToSrcAddr,
375  tgtToSrcWght,
376  srcSeedI,
377  tgtSeedI,
378  srcCellIDs,
379  mapFlag,
380  startSeedI
381  );
382  }
383  else
384  {
385  // if meshes are collocated, after inflating the source mesh bounding
386  // box tgt mesh cells may be transferred, but may still not overlap
387  // with the source mesh
388  return;
389  }
390 }
391 
392 
393 // ************************************************************************* //
cellVolumeWeightMethod(const polyMesh &src, const polyMesh &tgt)
Construct from source and target meshes.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
virtual const word & AMImethod() const
Return the corresponding AMI method for patch interpolation.
virtual ~cellVolumeWeightMethod()
Destructor.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void calculateAddressing(labelListList &srcToTgtCellAddr, scalarListList &srcToTgtCellWght, labelListList &tgtToSrcCellAddr, scalarListList &tgtToSrcCellWght, const label srcSeedI, const label tgtSeedI, const labelList &srcCellIDs, boolList &mapFlag, label &startSeedI)
Calculate the mesh-to-mesh addressing and weights.
Macros for easy insertion into run-time selection tables.
virtual void calculate(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, labelListList &tgtToTgtAddr, scalarListList &tgtToTgtWght)
Calculate addressing and weights.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
A class for handling words, derived from string.
Definition: word.H:59
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
Base class for mesh-to-mesh calculation methods.
void setNextCells(label &startSeedI, label &srcCelli, label &tgtCelli, const labelList &srcCellIDs, const boolList &mapFlag, const DynamicList< label > &visitedCells, labelList &seedCells) const
Set the next cells in the advancing front algorithm.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
A List with indirect addressing.
Definition: fvMatrix.H:106
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:236
bool findInitialSeeds(const labelList &srcCellIDs, const boolList &mapFlag, const label startSeedI, label &srcSeedI, label &tgtSeedI) const
Find indices of overlapping cells in src and tgt meshes - returns.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
Namespace for OpenFOAM.