mapNearestMethod.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-2021 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 "mapNearestMethod.H"
27 #include "mapNearestAMI.H"
28 #include "pointIndexHit.H"
29 #include "indexedOctree.H"
30 #include "treeDataCell.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(mapNearestMethod, 0);
38  addToRunTimeSelectionTable(meshToMeshMethod, mapNearestMethod, components);
39 }
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 (
45  const labelList& srcCellIDs,
46  const boolList& mapFlag,
47  const label startSeedI,
48  label& srcSeedI,
49  label& tgtSeedI
50 ) const
51 {
52  const vectorField& srcCcs = src_.cellCentres();
53 
54  for (label i = startSeedI; i < srcCellIDs.size(); i++)
55  {
56  label srcI = srcCellIDs[i];
57 
58  if (mapFlag[srcI])
59  {
60  const point& srcCc = srcCcs[srcI];
61  pointIndexHit hit =
62  tgt_.cellTree().findNearest(srcCc, great);
63 
64  if (hit.hit())
65  {
66  srcSeedI = srcI;
67  tgtSeedI = hit.index();
68 
69  return true;
70  }
71  else
72  {
74  << "Unable to find nearest target cell"
75  << " for source cell " << srcI
76  << " with centre " << srcCc
77  << abort(FatalError);
78  }
79 
80  break;
81  }
82  }
83 
84  if (debug)
85  {
86  Pout<< "could not find starting seed" << endl;
87  }
88 
89  return false;
90 }
91 
92 
94 (
95  labelListList& srcToTgtCellAddr,
96  scalarListList& srcToTgtCellWght,
97  labelListList& tgtToSrcCellAddr,
98  scalarListList& tgtToSrcCellWght,
99  const label srcSeedI,
100  const label tgtSeedI,
101  const labelList& srcCellIDs,
102  boolList& mapFlag,
103  label& startSeedI
104 )
105 {
106  List<DynamicList<label>> srcToTgt(src_.nCells());
107  List<DynamicList<label>> tgtToSrc(tgt_.nCells());
108 
109  const scalarField& srcVc = src_.cellVolumes();
110  const scalarField& tgtVc = tgt_.cellVolumes();
111 
112  {
113  label srcCelli = srcSeedI;
114  label tgtCelli = tgtSeedI;
115 
116  do
117  {
118  // find nearest tgt cell
119  findNearestCell(src_, tgt_, srcCelli, tgtCelli);
120 
121  // store src/tgt cell pair
122  srcToTgt[srcCelli].append(tgtCelli);
123  tgtToSrc[tgtCelli].append(srcCelli);
124 
125  // mark source cell srcCelli and tgtCelli as matched
126  mapFlag[srcCelli] = false;
127 
128  // accumulate intersection volume
129  V_ += srcVc[srcCelli];
130 
131  // find new source cell
132  setNextNearestCells
133  (
134  startSeedI,
135  srcCelli,
136  tgtCelli,
137  mapFlag,
138  srcCellIDs
139  );
140  }
141  while (srcCelli >= 0);
142  }
143 
144  // for the case of multiple source cells per target cell, select the
145  // nearest source cell only and discard the others
146  const vectorField& srcCc = src_.cellCentres();
147  const vectorField& tgtCc = tgt_.cellCentres();
148 
149  forAll(tgtToSrc, targetCelli)
150  {
151  if (tgtToSrc[targetCelli].size() > 1)
152  {
153  const vector& tgtC = tgtCc[targetCelli];
154 
155  DynamicList<label>& srcCells = tgtToSrc[targetCelli];
156 
157  label srcCelli = srcCells[0];
158  scalar d = magSqr(tgtC - srcCc[srcCelli]);
159 
160  for (label i = 1; i < srcCells.size(); i++)
161  {
162  label srcI = srcCells[i];
163  scalar dNew = magSqr(tgtC - srcCc[srcI]);
164  if (dNew < d)
165  {
166  d = dNew;
167  srcCelli = srcI;
168  }
169  }
170 
171  srcCells.clear();
172  srcCells.append(srcCelli);
173  }
174  }
175 
176  // If there are more target cells than source cells, some target cells
177  // might not yet be mapped
178  forAll(tgtToSrc, tgtCelli)
179  {
180  if (tgtToSrc[tgtCelli].empty())
181  {
182  label srcCelli = findMappedSrcCell(tgtCelli, tgtToSrc);
183 
184  findNearestCell(tgt_, src_, tgtCelli, srcCelli);
185 
186  tgtToSrc[tgtCelli].append(srcCelli);
187  }
188  }
189 
190  // transfer addressing into persistent storage
191  forAll(srcToTgtCellAddr, i)
192  {
193  srcToTgtCellWght[i] = scalarList(srcToTgt[i].size(), srcVc[i]);
194  srcToTgtCellAddr[i].transfer(srcToTgt[i]);
195  }
196 
197  forAll(tgtToSrcCellAddr, i)
198  {
199  tgtToSrcCellWght[i] = scalarList(tgtToSrc[i].size(), tgtVc[i]);
200  tgtToSrcCellAddr[i].transfer(tgtToSrc[i]);
201  }
202 }
203 
204 
206 (
207  const polyMesh& mesh1,
208  const polyMesh& mesh2,
209  const label cell1,
210  label& cell2
211 ) const
212 {
213  const vectorField& Cc1 = mesh1.cellCentres();
214  const vectorField& Cc2 = mesh2.cellCentres();
215 
216  const vector& p1 = Cc1[cell1];
217 
218  DynamicList<label> cells2(10);
219  cells2.append(cell2);
220 
221  DynamicList<label> visitedCells(10);
222 
223  scalar d = great;
224 
225  do
226  {
227  label c2 = cells2.remove();
228  visitedCells.append(c2);
229 
230  scalar dTest = magSqr(Cc2[c2] - p1);
231  if (dTest < d)
232  {
233  cell2 = c2;
234  d = dTest;
235  appendNbrCells(cell2, mesh2, visitedCells, cells2);
236  }
237 
238  } while (cells2.size() > 0);
239 }
240 
241 
243 (
244  label& startSeedI,
245  label& srcCelli,
246  label& tgtCelli,
247  boolList& mapFlag,
248  const labelList& srcCellIDs
249 ) const
250 {
251  const labelList& srcNbr = src_.cellCells()[srcCelli];
252 
253  srcCelli = -1;
254  forAll(srcNbr, i)
255  {
256  label celli = srcNbr[i];
257  if (mapFlag[celli])
258  {
259  srcCelli = celli;
260  return;
261  }
262  }
263 
264  for (label i = startSeedI; i < srcCellIDs.size(); i++)
265  {
266  label celli = srcCellIDs[i];
267  if (mapFlag[celli])
268  {
269  startSeedI = i;
270  break;
271  }
272  }
273 
274  (void)findInitialSeeds
275  (
276  srcCellIDs,
277  mapFlag,
278  startSeedI,
279  srcCelli,
280  tgtCelli
281  );
282 }
283 
284 
286 (
287  const label tgtCelli,
288  const List<DynamicList<label>>& tgtToSrc
289 ) const
290 {
291  DynamicList<label> testCells(10);
292  DynamicList<label> visitedCells(10);
293 
294  testCells.append(tgtCelli);
295 
296  do
297  {
298  // search target tgtCelli neighbours for match with source cell
299  label tgtI = testCells.remove();
300 
301  if (findIndex(visitedCells, tgtI) == -1)
302  {
303  visitedCells.append(tgtI);
304 
305  if (tgtToSrc[tgtI].size())
306  {
307  return tgtToSrc[tgtI][0];
308  }
309  else
310  {
311  const labelList& nbrCells = tgt_.cellCells()[tgtI];
312 
313  forAll(nbrCells, i)
314  {
315  if (findIndex(visitedCells, nbrCells[i]) == -1)
316  {
317  testCells.append(nbrCells[i]);
318  }
319  }
320  }
321  }
322  } while (testCells.size());
323 
324  // did not find any match - should not be possible to get here!
325  return -1;
326 }
327 
328 
329 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
330 
332 (
333  const polyMesh& src,
334  const polyMesh& tgt
335 )
336 :
337  meshToMeshMethod(src, tgt)
338 {}
339 
340 
341 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
342 
344 {}
345 
346 
347 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
348 
350 {
351  return mapNearestAMI::typeName;
352 }
353 
354 
356 (
357  labelListList& srcToTgtAddr,
358  scalarListList& srcToTgtWght,
359  labelListList& tgtToSrcAddr,
360  scalarListList& tgtToSrcWght
361 )
362 {
363  bool ok = initialise
364  (
365  srcToTgtAddr,
366  srcToTgtWght,
367  tgtToSrcAddr,
368  tgtToSrcWght
369  );
370 
371  if (!ok)
372  {
373  return;
374  }
375 
376  // (potentially) participating source mesh cells
377  const labelList srcCellIDs(maskCells());
378 
379  // list to keep track of whether src cell can be mapped
380  boolList mapFlag(src_.nCells(), false);
381  UIndirectList<bool>(mapFlag, srcCellIDs) = true;
382 
383  // find initial point in tgt mesh
384  label srcSeedI = -1;
385  label tgtSeedI = -1;
386  label startSeedI = 0;
387 
388  bool startWalk =
389  findInitialSeeds
390  (
391  srcCellIDs,
392  mapFlag,
393  startSeedI,
394  srcSeedI,
395  tgtSeedI
396  );
397 
398  if (startWalk)
399  {
400  calculateAddressing
401  (
402  srcToTgtAddr,
403  srcToTgtWght,
404  tgtToSrcAddr,
405  tgtToSrcWght,
406  srcSeedI,
407  tgtSeedI,
408  srcCellIDs,
409  mapFlag,
410  startSeedI
411  );
412  }
413  else
414  {
415  // if meshes are collocated, after inflating the source mesh bounding
416  // box tgt mesh cells may be transferred, but may still not overlap
417  // with the source mesh
418  return;
419  }
420 }
421 
422 
423 // ************************************************************************* //
virtual const word & AMImethod() const
Return the corresponding AMI method for patch interpolation.
virtual 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.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual void setNextNearestCells(label &startSeedI, label &srcCelli, label &tgtCelli, boolList &mapFlag, const labelList &srcCellIDs) const
Set the next cells for the marching front algorithm.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual void calculate(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, labelListList &tgtToTgtAddr, scalarListList &tgtToTgtWght)
Calculate addressing and weights.
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
virtual void findNearestCell(const polyMesh &mesh1, const polyMesh &mesh2, const label cell1, label &cell2) const
Find the nearest cell on mesh2 for cell1 on mesh1.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Macros for easy insertion into run-time selection tables.
const dimensionedScalar c2
Second radiation constant: default SI units: [m K].
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A class for handling words, derived from string.
Definition: word.H:59
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
bool hit() const
Is there a hit.
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensioned< scalar > magSqr(const dimensioned< Type > &)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
Base class for mesh-to-mesh calculation methods.
T remove()
Remove and return the top element.
Definition: DynamicListI.H:351
virtual ~mapNearestMethod()
Destructor.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
A List with indirect addressing.
Definition: fvMatrix.H:106
virtual 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.
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
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
label index() const
Return index.
virtual label findMappedSrcCell(const label tgtCelli, const List< DynamicList< label >> &tgtToSrc) const
Find a source cell mapped to target cell tgtCelli.
Namespace for OpenFOAM.
mapNearestMethod(const polyMesh &src, const polyMesh &tgt)
Construct from source and target meshes.