Loading...
Searching...
No Matches
calculateMeshToMesh0Addressing.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) 2011-2020 OpenFOAM Foundation
9 Copyright (C) 2020 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
27Description
28 private member of meshToMesh0.
29 Calculates mesh to mesh addressing pattern (for each cell from one mesh
30 find the closest cell centre in the other mesh).
31
32\*---------------------------------------------------------------------------*/
33
34#include "meshToMesh0.H"
35
36#include "treeDataCell.H"
37#include "treeDataFace.H"
38
39// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40
41void Foam::meshToMesh0::calcAddressing()
42{
44 << "Calculating mesh-to-mesh cell addressing" << endl;
45
46 // set reference to cells
47 const cellList& fromCells = fromMesh_.cells();
48 const pointField& fromPoints = fromMesh_.points();
49
50 // In an attempt to preserve the efficiency of linear search and prevent
51 // failure, a RESCUE mechanism will be set up. Here, we shall mark all
52 // cells next to the solid boundaries. If such a cell is found as the
53 // closest, the relationship between the origin and cell will be examined.
54 // If the origin is outside the cell, a global n-squared search is
55 // triggered.
56
57 // SETTING UP RESCUE
58
59 // visit all boundaries and mark the cell next to the boundary.
60
61 DebugInFunction << "Setting up rescue" << endl;
62
63 List<bool> boundaryCell(fromCells.size(), false);
64
65 // set reference to boundary
66 const polyPatchList& patchesFrom = fromMesh_.boundaryMesh();
67
68 forAll(patchesFrom, patchi)
69 {
70 // get reference to cells next to the boundary
71 const labelUList& bCells = patchesFrom[patchi].faceCells();
72
73 forAll(bCells, facei)
74 {
75 boundaryCell[bCells[facei]] = true;
76 }
77 }
78
79 treeBoundBox meshBb(fromPoints);
80 treeBoundBox shiftedBb(meshBb);
81
82 scalar typDim = meshBb.avgDim()/(2.0*cbrt(scalar(fromCells.size())));
83 shiftedBb.max() += vector::uniform(typDim);
84
86 << "\nMesh" << nl
87 << " bounding box : " << meshBb << nl
88 << " bounding box (shifted) : " << shiftedBb << nl
89 << " typical dimension : " << shiftedBb.avgDim() << endl;
90
91 indexedOctree<treeDataCell> cellTree
92 (
93 treeDataCell(fromMesh_, polyMesh::CELL_TETS),
94 shiftedBb, // overall bounding box
95 8, // maxLevel
96 10, // leafsize
97 6.0 // duplicity
98 );
99
100 if (debug)
101 {
102 cellTree.print(Pout, false, 0);
103 }
104
105 cellAddresses
106 (
107 cellAddressing_,
108 toMesh_.cellCentres(),
109 fromMesh_,
110 boundaryCell,
111 cellTree
112 );
113
114 forAll(toMesh_.boundaryMesh(), patchi)
115 {
116 const polyPatch& toPatch = toMesh_.boundaryMesh()[patchi];
117
118 if (cuttingPatches_.found(toPatch.name()))
119 {
120 boundaryAddressing_[patchi].setSize(toPatch.size());
121
122 cellAddresses
123 (
124 boundaryAddressing_[patchi],
125 toPatch.faceCentres(),
126 fromMesh_,
127 boundaryCell,
128 cellTree
129 );
130 }
131 else if
132 (
133 patchMap_.found(toPatch.name())
134 && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
135 )
136 {
137 const polyPatch& fromPatch = fromMesh_.boundaryMesh()
138 [
139 fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
140 ];
141
142 if (fromPatch.empty())
143 {
145 << "Source patch " << fromPatch.name()
146 << " has no faces. Not performing mapping for it."
147 << endl;
148 boundaryAddressing_[patchi].setSize(toPatch.size());
149 boundaryAddressing_[patchi] = -1;
150 }
151 else
152 {
153 treeBoundBox wallBb(fromPatch.localPoints());
154 treeBoundBox shiftedBb(wallBb);
155
156 scalar typDim =
157 wallBb.avgDim()/(2.0*sqrt(scalar(fromPatch.size())));
158
159 shiftedBb.max() += vector::uniform(typDim);
160
161 // Note: allow more levels than in meshSearch. Assume patch
162 // is not as big as all boundary faces
163 indexedOctree<treeDataFace> faceTree
164 (
165 treeDataFace(fromPatch),
166 shiftedBb, // overall search domain
167 12, // maxLevel
168 10, // leafsize
169 6.0 // duplicity
170 );
171
172 const vectorField::subField centresToBoundary =
173 toPatch.faceCentres();
174
175 boundaryAddressing_[patchi].setSize(toPatch.size());
176
177 scalar distSqr = wallBb.magSqr();
178
179 forAll(toPatch, toi)
180 {
181 boundaryAddressing_[patchi][toi] = faceTree.findNearest
182 (
183 centresToBoundary[toi],
184 distSqr
185 ).index();
186 }
187 }
188 }
189 }
190
192 << "Finished calculating mesh-to-mesh cell addressing" << endl;
193}
194
195
196void Foam::meshToMesh0::cellAddresses
197(
198 labelList& cellAddressing_,
199 const pointField& points,
200 const fvMesh& fromMesh,
201 const List<bool>& boundaryCell,
202 const indexedOctree<treeDataCell>& cellTree
203) const
204{
205 // the implemented search method is a simple neighbour array search.
206 // It starts from a cell zero, searches its neighbours and finds one
207 // which is nearer to the target point than the current position.
208 // The location of the "current position" is reset to that cell and
209 // search through the neighbours continues. The search is finished
210 // when all the neighbours of the cell are farther from the target
211 // point than the current cell
212
213 // set curCell label to zero (start)
214 label curCell = 0;
215
216 // set reference to cell to cell addressing
217 const vectorField& centresFrom = fromMesh.cellCentres();
218 const labelListList& cc = fromMesh.cellCells();
219
220 forAll(points, toI)
221 {
222 // pick up target position
223 const vector& p = points[toI];
224
225 // set the sqr-distance
226 scalar distSqr = p.distSqr(centresFrom[curCell]);
227
228 bool closer;
229
230 do
231 {
232 closer = false;
233
234 // set the current list of neighbouring cells
235 const labelList& neighbours = cc[curCell];
236
237 forAll(neighbours, nI)
238 {
239 scalar curDistSqr = p.distSqr(centresFrom[neighbours[nI]]);
240
241 // search through all the neighbours.
242 // If the cell is closer, reset current cell and distance
243 if (curDistSqr < (1 - SMALL)*distSqr)
244 {
245 curCell = neighbours[nI];
246 distSqr = curDistSqr;
247 closer = true; // a closer neighbour has been found
248 }
249 }
250 } while (closer);
251
252 cellAddressing_[toI] = -1;
253
254 // Check point is actually in the nearest cell
255 if (fromMesh.pointInCell(p, curCell))
256 {
257 cellAddressing_[toI] = curCell;
258 }
259 else
260 {
261 // If curCell is a boundary cell then the point maybe either outside
262 // the domain or in an other region of the doamin, either way use
263 // the octree search to find it.
264 if (boundaryCell[curCell])
265 {
266 cellAddressing_[toI] = cellTree.findInside(p);
267
268 if (cellAddressing_[toI] != -1)
269 {
270 curCell = cellAddressing_[toI];
271 }
272 }
273 else
274 {
275 // If not on the boundary search the neighbours
276 bool found = false;
277
278 // set the current list of neighbouring cells
279 const labelList& neighbours = cc[curCell];
280
281 forAll(neighbours, nI)
282 {
283 // search through all the neighbours.
284 // If point is in neighbour reset current cell
285 if (fromMesh.pointInCell(p, neighbours[nI]))
286 {
287 cellAddressing_[toI] = neighbours[nI];
288 found = true;
289 break;
290 }
291 }
292
293 if (!found)
294 {
295 // If still not found search the neighbour-neighbours
296
297 // set the current list of neighbouring cells
298 const labelList& neighbours = cc[curCell];
299
300 forAll(neighbours, nI)
301 {
302 // set the current list of neighbour-neighbouring cells
303 const labelList& nn = cc[neighbours[nI]];
304
305 forAll(nn, nI)
306 {
307 // search through all the neighbours.
308 // If point is in neighbour reset current cell
309 if (fromMesh.pointInCell(p, nn[nI]))
310 {
311 cellAddressing_[toI] = nn[nI];
312 found = true;
313 break;
314 }
315 }
316 if (found) break;
317 }
318 }
319
320 if (!found)
321 {
322 // Still not found so use the octree
323 cellAddressing_[toI] = cellTree.findInside(p);
324
325 if (cellAddressing_[toI] != -1)
326 {
327 curCell = cellAddressing_[toI];
328 }
329 }
330 }
331 }
332 }
333}
334
335
336// ************************************************************************* //
bool found
SubField< vector > subField
Definition Field.H:183
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
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Non-pointer based hierarchical recursive searching.
volScalarField & p
const pointField & points
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugInFunction
Report an information message using Foam::Info.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar cbrt(const dimensionedScalar &ds)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299