nearestCellsToCells.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-2025 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 "nearestCellsToCells.H"
27 #include "pointIndexHit.H"
28 #include "indexedOctree.H"
29 #include "treeDataCell.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace cellsToCellss
37 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 bool Foam::cellsToCellss::nearest::findInitialSeeds
47 (
48  const polyMesh& srcMesh,
49  const polyMesh& tgtMesh,
50  const labelList& srcCellIDs,
51  const boolList& mapFlag,
52  const label startSeedI,
53  label& srcSeedI,
54  label& tgtSeedI
55 ) const
56 {
57  const vectorField& srcCcs = srcMesh.cellCentres();
58 
59  for (label i = startSeedI; i < srcCellIDs.size(); i++)
60  {
61  label srcI = srcCellIDs[i];
62 
63  if (mapFlag[srcI])
64  {
65  const point& srcCc = srcCcs[srcI];
66  pointIndexHit hit =
67  tgtMesh.cellTree().findNearest(srcCc, great);
68 
69  if (hit.hit())
70  {
71  srcSeedI = srcI;
72  tgtSeedI = hit.index();
73 
74  return true;
75  }
76  else
77  {
79  << "Unable to find nearest target cell"
80  << " for source cell " << srcI
81  << " with centre " << srcCc
82  << abort(FatalError);
83  }
84 
85  break;
86  }
87  }
88 
89  if (debug)
90  {
91  Pout<< "could not find starting seed" << endl;
92  }
93 
94  return false;
95 }
96 
97 
98 Foam::scalar Foam::cellsToCellss::nearest::calculateAddressing
99 (
100  const polyMesh& srcMesh,
101  const polyMesh& tgtMesh,
102  labelListList& srcToTgtCellAddr,
103  scalarListList& srcToTgtCellWght,
104  labelListList& tgtToSrcCellAddr,
105  scalarListList& tgtToSrcCellWght,
106  const label srcSeedI,
107  const label tgtSeedI,
108  const labelList& srcCellIDs,
109  boolList& mapFlag,
110  label& startSeedI
111 )
112 {
113  scalar V = 0;
114 
115  List<DynamicList<label>> srcToTgt(srcMesh.nCells());
116  List<DynamicList<label>> tgtToSrc(tgtMesh.nCells());
117 
118  const scalarField& srcVc = srcMesh.cellVolumes();
119 
120  {
121  label srcCelli = srcSeedI;
122  label tgtCelli = tgtSeedI;
123 
124  do
125  {
126  // find nearest tgt cell
127  findNearestCell(srcMesh, tgtMesh, srcCelli, tgtCelli);
128 
129  // store src/tgt cell pair
130  srcToTgt[srcCelli].append(tgtCelli);
131  tgtToSrc[tgtCelli].append(srcCelli);
132 
133  // mark source cell srcCelli and tgtCelli as matched
134  mapFlag[srcCelli] = false;
135 
136  // accumulate intersection volume
137  V += srcVc[srcCelli];
138 
139  // find new source cell
140  setNextNearestCells
141  (
142  srcMesh,
143  tgtMesh,
144  startSeedI,
145  srcCelli,
146  tgtCelli,
147  mapFlag,
148  srcCellIDs
149  );
150  }
151  while (srcCelli >= 0);
152  }
153 
154  // for the case of multiple source cells per target cell, select the
155  // nearest source cell only and discard the others
156  const vectorField& srcCc = srcMesh.cellCentres();
157  const vectorField& tgtCc = tgtMesh.cellCentres();
158 
159  forAll(tgtToSrc, targetCelli)
160  {
161  if (tgtToSrc[targetCelli].size() > 1)
162  {
163  const vector& tgtC = tgtCc[targetCelli];
164 
165  DynamicList<label>& srcCells = tgtToSrc[targetCelli];
166 
167  label srcCelli = srcCells[0];
168  scalar d = magSqr(tgtC - srcCc[srcCelli]);
169 
170  for (label i = 1; i < srcCells.size(); i++)
171  {
172  label srcI = srcCells[i];
173  scalar dNew = magSqr(tgtC - srcCc[srcI]);
174  if (dNew < d)
175  {
176  d = dNew;
177  srcCelli = srcI;
178  }
179  }
180 
181  srcCells.clear();
182  srcCells.append(srcCelli);
183  }
184  }
185 
186  // If there are more target cells than source cells, some target cells
187  // might not yet be mapped
188  forAll(tgtToSrc, tgtCelli)
189  {
190  if (tgtToSrc[tgtCelli].empty())
191  {
192  label srcCelli =
193  findMappedSrcCell(srcMesh, tgtMesh, tgtCelli, tgtToSrc);
194 
195  findNearestCell(tgtMesh, srcMesh, tgtCelli, srcCelli);
196 
197  tgtToSrc[tgtCelli].append(srcCelli);
198  }
199  }
200 
201  // transfer addressing into persistent storage
202  forAll(srcToTgtCellAddr, i)
203  {
204  srcToTgtCellWght[i] = scalarList(srcToTgt[i].size(), scalar(1));
205  srcToTgtCellAddr[i].transfer(srcToTgt[i]);
206  }
207 
208  forAll(tgtToSrcCellAddr, i)
209  {
210  tgtToSrcCellWght[i] = scalarList(tgtToSrc[i].size(), scalar(1));
211  tgtToSrcCellAddr[i].transfer(tgtToSrc[i]);
212  }
213 
214  return V;
215 }
216 
217 
218 void Foam::cellsToCellss::nearest::findNearestCell
219 (
220  const polyMesh& mesh1,
221  const polyMesh& mesh2,
222  const label cell1,
223  label& cell2
224 ) const
225 {
226  const vectorField& Cc1 = mesh1.cellCentres();
227  const vectorField& Cc2 = mesh2.cellCentres();
228 
229  const vector& p1 = Cc1[cell1];
230 
231  DynamicList<label> cells2(10);
232  cells2.append(cell2);
233 
234  DynamicList<label> visitedCells(10);
235 
236  scalar d = great;
237 
238  do
239  {
240  label c2 = cells2.remove();
241  visitedCells.append(c2);
242 
243  scalar dTest = magSqr(Cc2[c2] - p1);
244  if (dTest < d)
245  {
246  cell2 = c2;
247  d = dTest;
248  appendNbrCells(cell2, mesh2, visitedCells, cells2);
249  }
250 
251  } while (cells2.size() > 0);
252 }
253 
254 
255 void Foam::cellsToCellss::nearest::setNextNearestCells
256 (
257  const polyMesh& srcMesh,
258  const polyMesh& tgtMesh,
259  label& startSeedI,
260  label& srcCelli,
261  label& tgtCelli,
262  boolList& mapFlag,
263  const labelList& srcCellIDs
264 ) const
265 {
266  const labelList& srcNbr = srcMesh.cellCells()[srcCelli];
267 
268  srcCelli = -1;
269  forAll(srcNbr, i)
270  {
271  label celli = srcNbr[i];
272  if (mapFlag[celli])
273  {
274  srcCelli = celli;
275  return;
276  }
277  }
278 
279  for (label i = startSeedI; i < srcCellIDs.size(); i++)
280  {
281  label celli = srcCellIDs[i];
282  if (mapFlag[celli])
283  {
284  startSeedI = i;
285  break;
286  }
287  }
288 
289  findInitialSeeds
290  (
291  srcMesh,
292  tgtMesh,
293  srcCellIDs,
294  mapFlag,
295  startSeedI,
296  srcCelli,
297  tgtCelli
298  );
299 }
300 
301 
302 Foam::label Foam::cellsToCellss::nearest::findMappedSrcCell
303 (
304  const polyMesh& srcMesh,
305  const polyMesh& tgtMesh,
306  const label tgtCelli,
307  const List<DynamicList<label>>& tgtToSrc
308 ) const
309 {
310  DynamicList<label> testCells(10);
311  DynamicList<label> visitedCells(10);
312 
313  testCells.append(tgtCelli);
314 
315  do
316  {
317  // search target tgtCelli neighbours for match with source cell
318  label tgtI = testCells.remove();
319 
320  if (findIndex(visitedCells, tgtI) == -1)
321  {
322  visitedCells.append(tgtI);
323 
324  if (tgtToSrc[tgtI].size())
325  {
326  return tgtToSrc[tgtI][0];
327  }
328  else
329  {
330  const labelList& nbrCells = tgtMesh.cellCells()[tgtI];
331 
332  forAll(nbrCells, i)
333  {
334  if (findIndex(visitedCells, nbrCells[i]) == -1)
335  {
336  testCells.append(nbrCells[i]);
337  }
338  }
339  }
340  }
341  } while (testCells.size());
342 
343  // did not find any match - should not be possible to get here!
344  return -1;
345 }
346 
347 
348 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
349 
351 (
352  const polyMesh& srcMesh,
353  const polyMesh& tgtMesh
354 )
355 {
356  initialise(srcMesh, tgtMesh);
357 
358  // Determine (potentially) participating source mesh cells
359  const labelList srcCellIDs(maskCells(srcMesh, tgtMesh));
360 
361  // Initialise list to keep track of whether src cell can be mapped
362  boolList mapFlag(srcMesh.nCells(), false);
363  UIndirectList<bool>(mapFlag, srcCellIDs) = true;
364 
365  // Find initial point in tgt mesh
366  label srcSeedI = -1;
367  label tgtSeedI = -1;
368  label startSeedI = 0;
369  bool startWalk =
370  findInitialSeeds
371  (
372  srcMesh,
373  tgtMesh,
374  srcCellIDs,
375  mapFlag,
376  startSeedI,
377  srcSeedI,
378  tgtSeedI
379  );
380 
381  if (startWalk)
382  {
383  return
384  calculateAddressing
385  (
386  srcMesh,
387  tgtMesh,
388  srcLocalTgtCells_,
389  srcWeights_,
390  tgtLocalSrcCells_,
391  tgtWeights_,
392  srcSeedI,
393  tgtSeedI,
394  srcCellIDs,
395  mapFlag,
396  startSeedI
397  );
398  }
399  else
400  {
401  return 0;
402  }
403 }
404 
405 
407 (
408  const polyMesh& srcMesh,
409  labelListList& srcToTgtAddr,
410  scalarListList& srcToTgtWght
411 ) const
412 {
413  forAll(srcToTgtWght, srcCelli)
414  {
415  if (srcToTgtWght[srcCelli].size() > 1)
416  {
417  srcToTgtAddr[srcCelli].resize(1);
418  srcToTgtWght[srcCelli].resize(1);
419  }
420  }
421 }
422 
423 
424 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
425 
427 :
428  cellsToCells()
429 {}
430 
431 
432 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
433 
435 {}
436 
437 
438 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
A List with indirect addressing.
Definition: UIndirectList.H:60
Class to calculate interpolative addressing and weights between the cells of two overlapping meshes.
Definition: cellsToCells.H:56
Nearest cells-to-cells interpolation class.
virtual scalar calculate(const polyMesh &srcMesh, const polyMesh &tgtMesh)
Calculate the addressing and weights.
virtual void normalise(const polyMesh &mesh, labelListList &localOtherCells, scalarListList &weights) const
Normalise the weights for a given mesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
label nCells() const
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
defineTypeNameAndDebug(intersection, 0)
addToRunTimeSelectionTable(cellsToCells, intersection, word)
const dimensionedScalar c2
Second radiation constant: default SI units: [m K].
Namespace for OpenFOAM.
List< scalarList > scalarListList
Definition: scalarList.H:51
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)