Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
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)
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> queuedCells(10);
221 DynamicList<label> visitedCells(10);
222
223 queuedCells.push_back(cell2);
224
225 scalar d = GREAT;
226
227 while (!queuedCells.empty())
228 {
229 // Process as LIFO
230 const label currCelli = queuedCells.back();
231 queuedCells.pop_back();
232 visitedCells.push_back(currCelli);
233
234 scalar dTest = p1.distSqr(Cc2[currCelli]);
235
236 if (dTest < d)
237 {
238 cell2 = currCelli;
239 d = dTest;
240 appendNbrCells(cell2, mesh2, visitedCells, queuedCells);
241 }
242 }
243}
244
245
247(
248 label& startSeedI,
249 label& srcCelli,
250 label& tgtCelli,
251 boolList& mapFlag,
252 const labelList& srcCellIDs
253) const
254{
255 const labelList& srcNbr = src_.cellCells()[srcCelli];
256
257 srcCelli = -1;
258 forAll(srcNbr, i)
259 {
260 label celli = srcNbr[i];
261 if (mapFlag[celli])
262 {
263 srcCelli = celli;
264 return;
265 }
266 }
267
268 for (label i = startSeedI; i < srcCellIDs.size(); i++)
269 {
270 label celli = srcCellIDs[i];
271 if (mapFlag[celli])
272 {
273 startSeedI = i;
274 break;
275 }
276 }
277
278 (void)findInitialSeeds
279 (
280 srcCellIDs,
281 mapFlag,
282 startSeedI,
283 srcCelli,
284 tgtCelli
285 );
286}
287
288
290(
291 const label tgtCelli,
292 const List<DynamicList<label>>& tgtToSrc
293) const
294{
295 DynamicList<label> queuedCells(16);
296 DynamicList<label> visitedCells(16);
297
298 queuedCells.push_back(tgtCelli);
299
300 while (!queuedCells.empty())
301 {
302 // Process as LIFO
303 const label tgtI = queuedCells.back();
304 queuedCells.pop_back();
305
306 // search target tgtCelli neighbours for match with source cell
307
308 if (!visitedCells.found(tgtI))
309 {
310 visitedCells.push_back(tgtI);
311
312 if (tgtToSrc[tgtI].size())
313 {
314 return tgtToSrc[tgtI][0];
315 }
316 else
317 {
318 const labelList& nbrCells = tgt_.cellCells()[tgtI];
319
320 for (const label nbrCelli : nbrCells)
321 {
322 if (!visitedCells.found(nbrCelli))
323 {
324 queuedCells.push_back(nbrCelli);
325 }
326 }
327 }
328 }
329 }
330
331 // did not find any match - should not be possible to get here!
332 return -1;
333}
334
335
336// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
337
339(
340 const polyMesh& src,
341 const polyMesh& tgt
342)
344 meshToMeshMethod(src, tgt)
345{}
346
347
348// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
351{}
352
353
354// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
355
357(
358 labelListList& srcToTgtAddr,
359 scalarListList& srcToTgtWght,
360 pointListList& srcToTgtVec,
361 labelListList& tgtToSrcAddr,
362 scalarListList& tgtToSrcWght,
363 pointListList& tgtToSrcVec
364)
365{
366 bool ok = initialise
367 (
368 srcToTgtAddr,
369 srcToTgtWght,
370 tgtToSrcAddr,
371 tgtToSrcWght
372 );
373
374 if (!ok)
375 {
376 return;
377 }
378
379 // (potentially) participating source mesh cells
380 const labelList srcCellIDs(maskCells());
381
382 // list to keep track of whether src cell can be mapped
383 boolList mapFlag(src_.nCells(), false);
384 boolUIndList(mapFlag, srcCellIDs) = true;
385
386 // find initial point in tgt mesh
387 label srcSeedI = -1;
388 label tgtSeedI = -1;
389 label startSeedI = 0;
390
391 bool startWalk =
392 findInitialSeeds
393 (
394 srcCellIDs,
395 mapFlag,
396 startSeedI,
397 srcSeedI,
398 tgtSeedI
399 );
400
401 if (startWalk)
402 {
403 calculateAddressing
404 (
405 srcToTgtAddr,
406 srcToTgtWght,
407 tgtToSrcAddr,
408 tgtToSrcWght,
409 srcSeedI,
410 tgtSeedI,
411 srcCellIDs,
412 mapFlag,
413 startSeedI
414 );
415 }
416 else
417 {
418 // if meshes are collocated, after inflating the source mesh bounding
419 // box tgt mesh cells may be transferred, but may still not overlap
420 // with the source mesh
421 return;
422 }
423}
424
425
426// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void pop_back(label n=1)
Reduce size by 1 or more elements. Can be called on an empty list.
void append(const T &val)
Copy append an element to the end of this list.
void push_back(const T &val)
Copy append an element to the end of this list.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
label index() const noexcept
Return the hit index.
bool hit() const noexcept
Is there a hit?
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
T & back()
Access last element of the list, position [size()-1].
Definition UListI.H:253
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
scalar distSqr(const Vector< Cmpt > &v2) const
The L2-norm distance squared from another vector. The magSqr() of the difference.
Definition VectorI.H:95
Map nearest mesh-to-mesh interpolation class.
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.
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.
virtual void setNextNearestCells(label &startSeedI, label &srcCelli, label &tgtCelli, boolList &mapFlag, const labelList &srcCellIDs) const
Set the next cells for the marching front algorithm.
virtual label findMappedSrcCell(const label tgtCelli, const List< DynamicList< label > > &tgtToSrc) const
Find a source cell mapped to target cell tgtCelli.
mapNearestMethod(const mapNearestMethod &)=delete
No copy construct.
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.
virtual void calculate(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, pointListList &srcToTgtVec, labelListList &tgtToSrcAddr, scalarListList &tgtToSrcWght, pointListList &tgtToSrcVec)
Calculate addressing and weights and optionally offset vectors.
virtual ~mapNearestMethod()
Destructor.
Base class for mesh-to-mesh calculation methods.
const polyMesh & tgt_
Reference to the target mesh.
const polyMesh & src() const
Return const access to the source mesh.
const polyMesh & src_
Reference to the source mesh.
const polyMesh & tgt() const
Return const access to the target mesh.
labelList maskCells() const
Return src cell IDs for the overlap region.
scalar V_
Cell total volume in overlap region [m3].
virtual bool initialise(labelListList &srcToTgtAddr, scalarListList &srcToTgtWght, labelListList &tgtToTgtAddr, scalarListList &tgtToTgtWght) const
virtual void appendNbrCells(const label tgtCelli, const polyMesh &mesh, const labelUList &visitedTgtCells, DynamicList< label > &nbrTgtCellIDs) const
Append target cell neighbour cells to cellIDs list.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition polyMesh.C:924
const vectorField & cellCentres() const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
UIndirectList< bool > boolUIndList
UIndirectList of bools.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
List< pointList > pointListList
List of pointList.
Definition pointList.H:35
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299