Loading...
Searching...
No Matches
cellVolumeWeightMethod.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 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
30#include "indexedOctree.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
40 (
42 cellVolumeWeightMethod,
43 components
44 );
45}
46
47// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48
50(
51 const labelList& srcCellIDs,
52 const boolList& mapFlag,
53 const label startSeedI,
54 label& srcSeedI,
55 label& tgtSeedI
56) const
57{
58 const cellList& srcCells = src_.cells();
59 const faceList& srcFaces = src_.faces();
60 const pointField& srcPts = src_.points();
61
62 for (label i = startSeedI; i < srcCellIDs.size(); ++i)
63 {
64 const label srcI = srcCellIDs[i];
65
66 if (mapFlag[srcI])
67 {
68 const pointField pts(srcCells[srcI].points(srcFaces, srcPts));
69
70 for (const point& pt : pts)
71 {
72 const label tgtI = tgt_.cellTree().findInside(pt);
73
74 if (tgtI != -1 && intersect(srcI, tgtI))
75 {
76 srcSeedI = srcI;
77 tgtSeedI = tgtI;
78
79 return true;
80 }
81 }
82 }
83 }
84
85 if (debug)
86 {
87 Pout<< "could not find starting seed" << endl;
88 }
89
90 return false;
91}
92
93
95(
96 labelListList& srcToTgtCellAddr,
97 scalarListList& srcToTgtCellWght,
98 labelListList& tgtToSrcCellAddr,
99 scalarListList& tgtToSrcCellWght,
100 const label srcSeedI,
101 const label tgtSeedI,
102 const labelList& srcCellIDs,
103 boolList& mapFlag,
104 label& startSeedI
105)
106{
107 label srcCelli = srcSeedI;
108 label tgtCelli = tgtSeedI;
109
110 List<DynamicList<label>> srcToTgtAddr(src_.nCells());
111 List<DynamicList<scalar>> srcToTgtWght(src_.nCells());
112
113 List<DynamicList<label>> tgtToSrcAddr(tgt_.nCells());
114 List<DynamicList<scalar>> tgtToSrcWght(tgt_.nCells());
115
116 // List of tgt cell neighbour cells
117 DynamicList<label> queuedCells(10);
118
119 // List of tgt cells currently visited for srcCellI to avoid multiple hits
120 DynamicList<label> visitedCells(10);
121
122 // list to keep track of tgt cells used to seed src cells
123 labelList seedCells(src_.nCells(), -1);
124 seedCells[srcCelli] = tgtCelli;
125
126 const scalarField& srcVol = src_.cellVolumes();
127
128 do
129 {
130 queuedCells.clear();
131 visitedCells.clear();
132
133 // Initial target cell and neighbours
134 queuedCells.push_back(tgtCelli);
135 appendNbrCells(tgtCelli, tgt_, visitedCells, queuedCells);
136
137 while (!queuedCells.empty())
138 {
139 // Process new target cell as LIFO
140 tgtCelli = queuedCells.back();
141 queuedCells.pop_back();
142 visitedCells.push_back(tgtCelli);
143
144 scalar vol = interVol(srcCelli, tgtCelli);
145
146 // accumulate addressing and weights for valid intersection
147 if (vol/srcVol[srcCelli] > tolerance_)
148 {
149 // store src/tgt cell pair
150 srcToTgtAddr[srcCelli].append(tgtCelli);
151 srcToTgtWght[srcCelli].append(vol);
152
153 tgtToSrcAddr[tgtCelli].append(srcCelli);
154 tgtToSrcWght[tgtCelli].append(vol);
155
156 appendNbrCells(tgtCelli, tgt_, visitedCells, queuedCells);
157
158 // accumulate intersection volume
159 V_ += vol;
160 }
161 }
162
163 mapFlag[srcCelli] = false;
164
165 // find new source seed cell
166 setNextCells
167 (
168 startSeedI,
169 srcCelli,
170 tgtCelli,
171 srcCellIDs,
172 mapFlag,
173 visitedCells,
174 seedCells
175 );
176 }
177 while (srcCelli != -1);
178
179 // transfer addressing into persistent storage
180 forAll(srcToTgtCellAddr, i)
181 {
182 srcToTgtCellAddr[i].transfer(srcToTgtAddr[i]);
183 srcToTgtCellWght[i].transfer(srcToTgtWght[i]);
184 }
185
186 forAll(tgtToSrcCellAddr, i)
187 {
188 tgtToSrcCellAddr[i].transfer(tgtToSrcAddr[i]);
189 tgtToSrcCellWght[i].transfer(tgtToSrcWght[i]);
190 }
191
192
193 if (debug%2)
194 {
195 // At this point the overlaps are still in volume so we could
196 // get out the relative error
197 forAll(srcToTgtCellAddr, cellI)
198 {
199 scalar srcVol = src_.cellVolumes()[cellI];
200 scalar tgtVol = sum(srcToTgtCellWght[cellI]);
201
202 if (mag(srcVol) > ROOTVSMALL && mag((tgtVol-srcVol)/srcVol) > 1e-6)
203 {
205 << "At cell " << cellI << " cc:"
206 << src_.cellCentres()[cellI]
207 << " vol:" << srcVol
208 << " total overlap volume:" << tgtVol
209 << endl;
210 }
211 }
212
213 forAll(tgtToSrcCellAddr, cellI)
214 {
215 scalar tgtVol = tgt_.cellVolumes()[cellI];
216 scalar srcVol = sum(tgtToSrcCellWght[cellI]);
217
218 if (mag(tgtVol) > ROOTVSMALL && mag((srcVol-tgtVol)/tgtVol) > 1e-6)
219 {
221 << "At cell " << cellI << " cc:"
222 << tgt_.cellCentres()[cellI]
223 << " vol:" << tgtVol
224 << " total overlap volume:" << srcVol
226 }
227 }
228 }
229}
230
231
233(
234 label& startSeedI,
235 label& srcCelli,
236 label& tgtCelli,
237 const labelList& srcCellIDs,
238 const boolList& mapFlag,
239 const labelUList& visitedCells,
240 labelList& seedCells
241) const
242{
243 const labelList& srcNbrCells = src_.cellCells()[srcCelli];
244
245 // set possible seeds for later use by querying all src cell neighbours
246 // with all visited target cells
247 bool valuesSet = false;
248 forAll(srcNbrCells, i)
249 {
250 label cellS = srcNbrCells[i];
251
252 if (mapFlag[cellS] && seedCells[cellS] == -1)
253 {
254 forAll(visitedCells, j)
255 {
256 label cellT = visitedCells[j];
257
258 if (intersect(cellS, cellT))
259 {
260 seedCells[cellS] = cellT;
261
262 if (!valuesSet)
263 {
264 srcCelli = cellS;
265 tgtCelli = cellT;
266 valuesSet = true;
267 }
268 }
269 }
270 }
271 }
272
273 // set next src and tgt cells if not set above
274 if (valuesSet)
275 {
276 return;
277 }
278 else
279 {
280 // try to use existing seed
281 bool foundNextSeed = false;
282 for (label i = startSeedI; i < srcCellIDs.size(); i++)
283 {
284 label cellS = srcCellIDs[i];
285
286 if (mapFlag[cellS])
287 {
288 if (!foundNextSeed)
289 {
290 startSeedI = i;
291 foundNextSeed = true;
292 }
293
294 if (seedCells[cellS] != -1)
295 {
296 srcCelli = cellS;
297 tgtCelli = seedCells[cellS];
298
299 return;
300 }
301 }
302 }
303
304 // perform new search to find match
305 if (debug)
306 {
307 Pout<< "Advancing front stalled: searching for new "
308 << "target cell" << endl;
309 }
310
311 bool restart =
312 findInitialSeeds
313 (
314 srcCellIDs,
315 mapFlag,
316 startSeedI,
317 srcCelli,
318 tgtCelli
319 );
320
321 if (restart)
322 {
323 // successfully found new starting seed-pair
324 return;
325 }
326 }
327
328 // if we have got to here, there are no more src/tgt cell intersections
329 srcCelli = -1;
330 tgtCelli = -1;
331}
332
333
334// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
335
336Foam::cellVolumeWeightMethod::cellVolumeWeightMethod
337(
338 const polyMesh& src,
339 const polyMesh& tgt
340)
342 meshToMeshMethod(src, tgt)
343{}
344
345
346// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
349{}
350
351
352// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
353
355(
356 labelListList& srcToTgtAddr,
357 scalarListList& srcToTgtWght,
358 pointListList& srcToTgtVec,
359 labelListList& tgtToSrcAddr,
360 scalarListList& tgtToSrcWght,
361 pointListList& tgtToSrcVec
362)
363{
364 bool ok = initialise
365 (
366 srcToTgtAddr,
367 srcToTgtWght,
368 tgtToSrcAddr,
369 tgtToSrcWght
370 );
371
372 if (!ok)
373 {
374 return;
375 }
376
377 // (potentially) participating source mesh cells
378 const labelList srcCellIDs(maskCells());
379
380 // list to keep track of whether src cell can be mapped
381 boolList mapFlag(src_.nCells(), false);
382 boolUIndList(mapFlag, srcCellIDs) = true;
383
384 // find initial point in tgt mesh
385 label srcSeedI = -1;
386 label tgtSeedI = -1;
387 label startSeedI = 0;
388
389 bool startWalk =
390 findInitialSeeds
391 (
392 srcCellIDs,
393 mapFlag,
394 startSeedI,
395 srcSeedI,
396 tgtSeedI
397 );
398
399 if (startWalk)
400 {
401 calculateAddressing
402 (
403 srcToTgtAddr,
404 srcToTgtWght,
405 tgtToSrcAddr,
406 tgtToSrcWght,
407 srcSeedI,
408 tgtSeedI,
409 srcCellIDs,
410 mapFlag,
411 startSeedI
412 );
413 }
414 else
415 {
416 // if meshes are collocated, after inflating the source mesh bounding
417 // box tgt mesh cells may be transferred, but may still not overlap
418 // with the source mesh
419 return;
420 }
421}
422
423
424// ************************************************************************* //
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 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
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
Cell-volume-weighted mesh-to-mesh interpolation class.
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 ~cellVolumeWeightMethod()
Destructor.
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.
void setNextCells(label &startSeedI, label &srcCelli, label &tgtCelli, const labelList &srcCellIDs, const boolList &mapFlag, const labelUList &visitedCells, labelList &seedCells) const
Set the next cells in the advancing front algorithm.
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.
static scalar tolerance_
Tolerance used in volume overlap calculations.
virtual scalar interVol(const label srcCelli, const label tgtCelli) const
Return the intersection volume between two cells.
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 intersect(const label srcCelli, const label tgtCelli) const
Return the true if cells intersect.
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
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition polyMesh.C:924
const cellList & cells() const
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
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
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< pointList > pointListList
List of pointList.
Definition pointList.H:35
const pointField & pts
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299