matchingCellsToCells.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 "matchingCellsToCells.H"
27 #include "indexedOctree.H"
28 #include "treeDataCell.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace cellsToCellss
36 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 bool Foam::cellsToCellss::matching::intersect
46 (
47  const polyMesh& srcMesh,
48  const polyMesh& tgtMesh,
49  const label srcCelli,
50  const label tgtCelli
51 ) const
52 {
53  return tgtMesh.pointInCell
54  (
55  srcMesh.cellCentres()[srcCelli],
56  tgtCelli,
58  );
59 }
60 
61 
62 bool Foam::cellsToCellss::matching::findInitialSeeds
63 (
64  const polyMesh& srcMesh,
65  const polyMesh& tgtMesh,
66  const labelList& srcCellIDs,
67  const boolList& mapFlag,
68  const label startSeedI,
69  label& srcSeedI,
70  label& tgtSeedI
71 ) const
72 {
73  const cellList& srcCells = srcMesh.cells();
74  const faceList& srcFaces = srcMesh.faces();
75  const pointField& srcPts = srcMesh.points();
76 
77  for (label i = startSeedI; i < srcCellIDs.size(); i++)
78  {
79  label srcI = srcCellIDs[i];
80 
81  if (mapFlag[srcI])
82  {
83  const point srcCtr(srcCells[srcI].centre(srcPts, srcFaces));
84  label tgtI = tgtMesh.cellTree().findInside(srcCtr);
85 
86  if (tgtI != -1 && intersect(srcMesh, tgtMesh, srcI, tgtI))
87  {
88  srcSeedI = srcI;
89  tgtSeedI = tgtI;
90 
91  return true;
92  }
93  }
94  }
95 
96  if (debug)
97  {
98  Pout<< "could not find starting seed" << endl;
99  }
100 
101  return false;
102 }
103 
104 
105 Foam::scalar Foam::cellsToCellss::matching::calculateAddressing
106 (
107  const polyMesh& srcMesh,
108  const polyMesh& tgtMesh,
109  labelListList& srcToTgtCellAddr,
110  scalarListList& srcToTgtCellWght,
111  labelListList& tgtToSrcCellAddr,
112  scalarListList& tgtToSrcCellWght,
113  const label srcSeedI,
114  const label tgtSeedI,
115  const labelList& srcCellIDs, // not used
116  boolList& mapFlag,
117  label& startSeedI
118 )
119 {
120  scalar V = 0;
121 
122  // store a list of src cells already mapped
123  labelList srcTgtSeed(srcMesh.nCells(), -1);
124 
125  List<DynamicList<label>> srcToTgt(srcMesh.nCells());
126  List<DynamicList<label>> tgtToSrc(tgtMesh.nCells());
127 
128  DynamicList<label> srcSeeds(10);
129 
130  const scalarField& srcVc = srcMesh.cellVolumes();
131 
132  label srcCelli = srcSeedI;
133  label tgtCelli = tgtSeedI;
134 
135  do
136  {
137  // store src/tgt cell pair
138  srcToTgt[srcCelli].append(tgtCelli);
139  tgtToSrc[tgtCelli].append(srcCelli);
140 
141  // mark source cell srcSeedI as matched
142  mapFlag[srcCelli] = false;
143 
144  // accumulate intersection volume
145  V += srcVc[srcCelli];
146 
147  // find new source seed cell
148  appendToDirectSeeds
149  (
150  srcMesh,
151  tgtMesh,
152  mapFlag,
153  srcTgtSeed,
154  srcSeeds,
155  srcCelli,
156  tgtCelli
157  );
158  }
159  while (srcCelli >= 0);
160 
161  // transfer addressing into persistent storage
162  forAll(srcToTgtCellAddr, i)
163  {
164  srcToTgtCellWght[i] = scalarList(srcToTgt[i].size(), scalar(1));
165  srcToTgtCellAddr[i].transfer(srcToTgt[i]);
166  }
167 
168  forAll(tgtToSrcCellAddr, i)
169  {
170  tgtToSrcCellWght[i] = scalarList(tgtToSrc[i].size(), scalar(1));
171  tgtToSrcCellAddr[i].transfer(tgtToSrc[i]);
172  }
173 
174  return V;
175 }
176 
177 
178 void Foam::cellsToCellss::matching::appendToDirectSeeds
179 (
180  const polyMesh& srcMesh,
181  const polyMesh& tgtMesh,
182  boolList& mapFlag,
183  labelList& srcTgtSeed,
184  DynamicList<label>& srcSeeds,
185  label& srcSeedI,
186  label& tgtSeedI
187 ) const
188 {
189  const labelList& srcNbr = srcMesh.cellCells()[srcSeedI];
190  const labelList& tgtNbr = tgtMesh.cellCells()[tgtSeedI];
191 
192  forAll(srcNbr, i)
193  {
194  label srcI = srcNbr[i];
195 
196  if (mapFlag[srcI] && (srcTgtSeed[srcI] == -1))
197  {
198  // source cell srcI not yet mapped
199 
200  // identify if target cell exists for source cell srcI
201  bool found = false;
202  forAll(tgtNbr, j)
203  {
204  label tgtI = tgtNbr[j];
205 
206  if (intersect(srcMesh, tgtMesh, srcI, tgtI))
207  {
208  // new match - append to lists
209  found = true;
210 
211  srcTgtSeed[srcI] = tgtI;
212  srcSeeds.append(srcI);
213 
214  break;
215  }
216  }
217 
218  if (!found)
219  {
220  // no match available for source cell srcI
221  mapFlag[srcI] = false;
222  }
223  }
224  }
225 
226  if (srcSeeds.size())
227  {
228  srcSeedI = srcSeeds.remove();
229  tgtSeedI = srcTgtSeed[srcSeedI];
230  }
231  else
232  {
233  srcSeedI = -1;
234  tgtSeedI = -1;
235  }
236 }
237 
238 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
239 
241 (
242  const polyMesh& srcMesh,
243  const polyMesh& tgtMesh
244 )
245 {
246  initialise(srcMesh, tgtMesh);
247 
248  // Determine (potentially) participating source mesh cells
249  const labelList srcCellIDs(maskCells(srcMesh, tgtMesh));
250 
251  // Initialise list to keep track of whether src cell can be mapped
252  boolList mapFlag(srcMesh.nCells(), false);
253  UIndirectList<bool>(mapFlag, srcCellIDs) = true;
254 
255  // Find initial point in tgt mesh
256  label srcSeedI = -1;
257  label tgtSeedI = -1;
258  label startSeedI = 0;
259  bool startWalk =
260  findInitialSeeds
261  (
262  srcMesh,
263  tgtMesh,
264  srcCellIDs,
265  mapFlag,
266  startSeedI,
267  srcSeedI,
268  tgtSeedI
269  );
270 
271  if (startWalk)
272  {
273  return
274  calculateAddressing
275  (
276  srcMesh,
277  tgtMesh,
278  srcLocalTgtCells_,
279  srcWeights_,
280  tgtLocalSrcCells_,
281  tgtWeights_,
282  srcSeedI,
283  tgtSeedI,
284  srcCellIDs,
285  mapFlag,
286  startSeedI
287  );
288  }
289  else
290  {
291  return 0;
292  }
293 }
294 
295 
297 (
298  const polyMesh& srcMesh,
299  labelListList& srcToTgtAddr,
300  scalarListList& srcToTgtWght
301 ) const
302 {
303  forAll(srcToTgtWght, srcCelli)
304  {
305  if (srcToTgtWght[srcCelli].size() > 1)
306  {
307  srcToTgtAddr[srcCelli].resize(1);
308  srcToTgtWght[srcCelli].resize(1);
309  }
310  }
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
315 
317 :
318  cellsToCells()
319 {}
320 
321 
322 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
323 
325 {}
326 
327 
328 // ************************************************************************* //
bool found
#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
Matching, one-to-one, 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
defineTypeNameAndDebug(intersection, 0)
addToRunTimeSelectionTable(cellsToCells, intersection, word)
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
List< cell > cellList
list of cells
Definition: cellList.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
vector point
Point is a vector.
Definition: point.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
List< face > faceList
Definition: faceListFwd.H:41