Loading...
Searching...
No Matches
directMethod.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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "directMethod.H"
29#include "indexedOctree.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
37 defineTypeNameAndDebug(directMethod, 0);
38 addToRunTimeSelectionTable(meshToMeshMethod, directMethod, components);
39}
40
41// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42
44(
45 const label srcCelli,
46 const label tgtCelli
47) const
48{
49 return tgt_.pointInCell
50 (
51 src_.cellCentres()[srcCelli],
52 tgtCelli,
54 );
55}
56
57
59(
60 const labelList& srcCellIDs,
61 const boolList& mapFlag,
62 const label startSeedI,
63 label& srcSeedI,
64 label& tgtSeedI
65) const
66{
67 const cellList& srcCells = src_.cells();
68 const faceList& srcFaces = src_.faces();
69 const pointField& srcPts = src_.points();
70
71 for (label i = startSeedI; i < srcCellIDs.size(); i++)
72 {
73 label srcI = srcCellIDs[i];
74
75 if (mapFlag[srcI])
76 {
77 const point srcCtr(srcCells[srcI].centre(srcPts, srcFaces));
78 label tgtI = tgt_.cellTree().findInside(srcCtr);
79
80 if (tgtI != -1 && intersect(srcI, tgtI))
81 {
82 srcSeedI = srcI;
83 tgtSeedI = tgtI;
84
85 return true;
86 }
87 }
88 }
89
90 if (debug)
91 {
92 Pout<< "could not find starting seed" << endl;
93 }
94
95 return false;
96}
97
98
100(
101 labelListList& srcToTgtCellAddr,
102 scalarListList& srcToTgtCellWght,
103 labelListList& tgtToSrcCellAddr,
104 scalarListList& tgtToSrcCellWght,
105 const label srcSeedI,
106 const label tgtSeedI,
107 const labelList& srcCellIDs, // not used
108 boolList& mapFlag,
109 label& startSeedI
110)
111{
112 // store a list of src cells already mapped
113 labelList srcTgtSeed(src_.nCells(), -1);
114
115 List<DynamicList<label>> srcToTgt(src_.nCells());
116 List<DynamicList<label>> tgtToSrc(tgt_.nCells());
117
118 DynamicList<label> srcSeeds(10);
119
120 const scalarField& srcVc = src_.cellVolumes();
121 const scalarField& tgtVc = tgt_.cellVolumes();
122
123 label srcCelli = srcSeedI;
124 label tgtCelli = tgtSeedI;
125
126 do
127 {
128 // store src/tgt cell pair
129 srcToTgt[srcCelli].append(tgtCelli);
130 tgtToSrc[tgtCelli].append(srcCelli);
131
132 // mark source cell srcSeedI as matched
133 mapFlag[srcCelli] = false;
134
135 // accumulate intersection volume
136 V_ += srcVc[srcCelli];
137
138 // find new source seed cell
139 appendToDirectSeeds
140 (
141 mapFlag,
142 srcTgtSeed,
143 srcSeeds,
144 srcCelli,
145 tgtCelli
146 );
147 }
148 while (srcCelli >= 0);
149
150 // transfer addressing into persistent storage
151 forAll(srcToTgtCellAddr, i)
152 {
153 srcToTgtCellWght[i] = scalarList(srcToTgt[i].size(), srcVc[i]);
154 srcToTgtCellAddr[i].transfer(srcToTgt[i]);
155 }
156
157 forAll(tgtToSrcCellAddr, i)
159 tgtToSrcCellWght[i] = scalarList(tgtToSrc[i].size(), tgtVc[i]);
160 tgtToSrcCellAddr[i].transfer(tgtToSrc[i]);
161 }
162}
163
164
166(
167 boolList& mapFlag,
168 labelList& srcTgtSeed,
169 DynamicList<label>& srcSeeds,
170 label& srcSeedI,
171 label& tgtSeedI
172) const
173{
174 const labelList& srcNbr = src_.cellCells()[srcSeedI];
175 const labelList& tgtNbr = tgt_.cellCells()[tgtSeedI];
176
177 for (const label srcI : srcNbr)
178 {
179 if (mapFlag[srcI] && (srcTgtSeed[srcI] == -1))
180 {
181 // source cell srcI not yet mapped
182
183 // identify if target cell exists for source cell srcI
184 bool found = false;
185 for (const label tgtI : tgtNbr)
186 {
187 if (intersect(srcI, tgtI))
188 {
189 // new match - append to lists
190 found = true;
191
192 srcTgtSeed[srcI] = tgtI;
193 srcSeeds.push_back(srcI);
194
195 break;
196 }
197 }
198
199 if (!found)
200 {
201 // no match available for source cell srcI
202 mapFlag[srcI] = false;
203 }
204 }
205 }
206
207 if (!srcSeeds.empty())
208 {
209 srcSeedI = srcSeeds.back();
210 tgtSeedI = srcTgtSeed[srcSeedI];
211 srcSeeds.pop_back();
212 }
213 else
214 {
215 srcSeedI = -1;
216 tgtSeedI = -1;
217 }
218}
219
220// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
221
222Foam::directMethod::directMethod
223(
224 const polyMesh& src,
225 const polyMesh& tgt
226)
228 meshToMeshMethod(src, tgt)
229{}
230
231
232// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
235{}
236
237
238// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
239
241(
242 labelListList& srcToTgtAddr,
243 scalarListList& srcToTgtWght,
244 pointListList& srcToTgtVec,
245 labelListList& tgtToSrcAddr,
246 scalarListList& tgtToSrcWght,
247 pointListList& tgtToSrcVec
248)
249{
250 bool ok = initialise
251 (
252 srcToTgtAddr,
253 srcToTgtWght,
254 tgtToSrcAddr,
255 tgtToSrcWght
256 );
257
258 if (!ok)
259 {
260 return;
261 }
262
263 // (potentially) participating source mesh cells
264 const labelList srcCellIDs(maskCells());
265
266 // list to keep track of whether src cell can be mapped
267 boolList mapFlag(src_.nCells(), false);
268 boolUIndList(mapFlag, srcCellIDs) = true;
269
270 // find initial point in tgt mesh
271 label srcSeedI = -1;
272 label tgtSeedI = -1;
273 label startSeedI = 0;
274
275 bool startWalk =
276 findInitialSeeds
277 (
278 srcCellIDs,
279 mapFlag,
280 startSeedI,
281 srcSeedI,
282 tgtSeedI
283 );
284
285 if (startWalk)
286 {
287 calculateAddressing
288 (
289 srcToTgtAddr,
290 srcToTgtWght,
291 tgtToSrcAddr,
292 tgtToSrcWght,
293 srcSeedI,
294 tgtSeedI,
295 srcCellIDs,
296 mapFlag,
297 startSeedI
298 );
299 }
300 else
301 {
302 // if meshes are collocated, after inflating the source mesh bounding
303 // box tgt mesh cells may be transferred, but may still not overlap
304 // with the source mesh
305 return;
306 }
307}
308
309
310// ************************************************************************* //
bool found
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 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
Direct (one-to-one cell correspondence) 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 ~directMethod()
Destructor.
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 void appendToDirectSeeds(boolList &mapFlag, labelList &srcTgtSeed, DynamicList< label > &srcSeeds, label &srcSeedI, label &tgtSeedI) const
Append to list of src mesh seed indices.
virtual bool intersect(const label srcCelli, const label tgtCelli) const
Return the true if cells intersect.
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
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition polyMesh.C:1386
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
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
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.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
List< pointList > pointListList
List of pointList.
Definition pointList.H:35
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299