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 -------------------------------------------------------------------------------
11 License
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 Description
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 
41 void 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 
81  scalar typDim = meshBb.avgDim()/(2.0*cbrt(scalar(fromCells.size())));
82 
83  treeBoundBox shiftedBb
84  (
85  meshBb.min(),
86  meshBb.max() + vector(typDim, typDim, typDim)
87  );
88 
89  DebugInfo
90  << "\nMesh" << nl
91  << " bounding box : " << meshBb << nl
92  << " bounding box (shifted) : " << shiftedBb << nl
93  << " typical dimension : " << shiftedBb.typDim() << endl;
94 
95  indexedOctree<treeDataCell> oc
96  (
97  treeDataCell(false, fromMesh_, polyMesh::CELL_TETS),
98  shiftedBb, // overall bounding box
99  8, // maxLevel
100  10, // leafsize
101  6.0 // duplicity
102  );
103 
104  if (debug)
105  {
106  oc.print(Pout, false, 0);
107  }
108 
109  cellAddresses
110  (
111  cellAddressing_,
112  toMesh_.cellCentres(),
113  fromMesh_,
114  boundaryCell,
115  oc
116  );
117 
118  forAll(toMesh_.boundaryMesh(), patchi)
119  {
120  const polyPatch& toPatch = toMesh_.boundaryMesh()[patchi];
121 
122  if (cuttingPatches_.found(toPatch.name()))
123  {
124  boundaryAddressing_[patchi].setSize(toPatch.size());
125 
126  cellAddresses
127  (
128  boundaryAddressing_[patchi],
129  toPatch.faceCentres(),
130  fromMesh_,
131  boundaryCell,
132  oc
133  );
134  }
135  else if
136  (
137  patchMap_.found(toPatch.name())
138  && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
139  )
140  {
141  const polyPatch& fromPatch = fromMesh_.boundaryMesh()
142  [
143  fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
144  ];
145 
146  if (fromPatch.empty())
147  {
149  << "Source patch " << fromPatch.name()
150  << " has no faces. Not performing mapping for it."
151  << endl;
152  boundaryAddressing_[patchi].setSize(toPatch.size());
153  boundaryAddressing_[patchi] = -1;
154  }
155  else
156  {
157  treeBoundBox wallBb(fromPatch.localPoints());
158  scalar typDim =
159  wallBb.avgDim()/(2.0*sqrt(scalar(fromPatch.size())));
160 
161  treeBoundBox shiftedBb
162  (
163  wallBb.min(),
164  wallBb.max() + vector(typDim, typDim, typDim)
165  );
166 
167  // Note: allow more levels than in meshSearch. Assume patch
168  // is not as big as all boundary faces
169  indexedOctree<treeDataFace> oc
170  (
171  treeDataFace(false, fromPatch),
172  shiftedBb, // overall search domain
173  12, // maxLevel
174  10, // leafsize
175  6.0 // duplicity
176  );
177 
178  const vectorField::subField centresToBoundary =
179  toPatch.faceCentres();
180 
181  boundaryAddressing_[patchi].setSize(toPatch.size());
182 
183  scalar distSqr = sqr(wallBb.mag());
184 
185  forAll(toPatch, toi)
186  {
187  boundaryAddressing_[patchi][toi] = oc.findNearest
188  (
189  centresToBoundary[toi],
190  distSqr
191  ).index();
192  }
193  }
194  }
195  }
196 
198  << "Finished calculating mesh-to-mesh cell addressing" << endl;
199 }
200 
201 
202 void Foam::meshToMesh0::cellAddresses
203 (
204  labelList& cellAddressing_,
205  const pointField& points,
206  const fvMesh& fromMesh,
207  const List<bool>& boundaryCell,
208  const indexedOctree<treeDataCell>& oc
209 ) const
210 {
211  // the implemented search method is a simple neighbour array search.
212  // It starts from a cell zero, searches its neighbours and finds one
213  // which is nearer to the target point than the current position.
214  // The location of the "current position" is reset to that cell and
215  // search through the neighbours continues. The search is finished
216  // when all the neighbours of the cell are farther from the target
217  // point than the current cell
218 
219  // set curCell label to zero (start)
220  label curCell = 0;
221 
222  // set reference to cell to cell addressing
223  const vectorField& centresFrom = fromMesh.cellCentres();
224  const labelListList& cc = fromMesh.cellCells();
225 
226  forAll(points, toI)
227  {
228  // pick up target position
229  const vector& p = points[toI];
230 
231  // set the sqr-distance
232  scalar distSqr = magSqr(p - centresFrom[curCell]);
233 
234  bool closer;
235 
236  do
237  {
238  closer = false;
239 
240  // set the current list of neighbouring cells
241  const labelList& neighbours = cc[curCell];
242 
243  forAll(neighbours, nI)
244  {
245  scalar curDistSqr =
246  magSqr(p - centresFrom[neighbours[nI]]);
247 
248  // search through all the neighbours.
249  // If the cell is closer, reset current cell and distance
250  if (curDistSqr < (1 - SMALL)*distSqr)
251  {
252  curCell = neighbours[nI];
253  distSqr = curDistSqr;
254  closer = true; // a closer neighbour has been found
255  }
256  }
257  } while (closer);
258 
259  cellAddressing_[toI] = -1;
260 
261  // Check point is actually in the nearest cell
262  if (fromMesh.pointInCell(p, curCell))
263  {
264  cellAddressing_[toI] = curCell;
265  }
266  else
267  {
268  // If curCell is a boundary cell then the point maybe either outside
269  // the domain or in an other region of the doamin, either way use
270  // the octree search to find it.
271  if (boundaryCell[curCell])
272  {
273  cellAddressing_[toI] = oc.findInside(p);
274 
275  if (cellAddressing_[toI] != -1)
276  {
277  curCell = cellAddressing_[toI];
278  }
279  }
280  else
281  {
282  // If not on the boundary search the neighbours
283  bool found = false;
284 
285  // set the current list of neighbouring cells
286  const labelList& neighbours = cc[curCell];
287 
288  forAll(neighbours, nI)
289  {
290  // search through all the neighbours.
291  // If point is in neighbour reset current cell
292  if (fromMesh.pointInCell(p, neighbours[nI]))
293  {
294  cellAddressing_[toI] = neighbours[nI];
295  found = true;
296  break;
297  }
298  }
299 
300  if (!found)
301  {
302  // If still not found search the neighbour-neighbours
303 
304  // set the current list of neighbouring cells
305  const labelList& neighbours = cc[curCell];
306 
307  forAll(neighbours, nI)
308  {
309  // set the current list of neighbour-neighbouring cells
310  const labelList& nn = cc[neighbours[nI]];
311 
312  forAll(nn, nI)
313  {
314  // search through all the neighbours.
315  // If point is in neighbour reset current cell
316  if (fromMesh.pointInCell(p, nn[nI]))
317  {
318  cellAddressing_[toI] = nn[nI];
319  found = true;
320  break;
321  }
322  }
323  if (found) break;
324  }
325  }
326 
327  if (!found)
328  {
329  // Still not found so use the octree
330  cellAddressing_[toI] = oc.findInside(p);
331 
332  if (cellAddressing_[toI] != -1)
333  {
334  curCell = cellAddressing_[toI];
335  }
336  }
337  }
338  }
339  }
340 }
341 
342 
343 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
p
volScalarField & p
Definition: createFieldRefs.H:8
meshToMesh0.H
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
meshBb
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
treeDataFace.H
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::polyPatchList
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:47
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
treeDataCell.H
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::HashTable::find
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition: HashTableI.H:114
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::polyMesh::CELL_TETS
Definition: polyMesh.H:109
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::Field< vector >::subField
SubField< vector > subField
Declare type of subField.
Definition: Field.H:89
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Foam::HashTable::found
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328