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