calculateMeshToMesh0Addressing.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "meshToMesh0.H"
27 #include "treeDataCell.H"
28 #include "treeDataFace.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 void Foam::meshToMesh0::calcAddressing()
33 {
34  if (debug)
35  {
37  << "Calculating mesh-to-mesh cell addressing" << endl;
38  }
39 
40  // Set reference to cells
41  const cellList& fromCells = fromMesh_.cells();
42  const pointField& fromPoints = fromMesh_.points();
43 
44  // In an attempt to preserve the efficiency of linear search and prevent
45  // failure, a RESCUE mechanism will be set up. Here, we shall mark all
46  // cells next to the solid boundaries. If such a cell is found as the
47  // closest, the relationship between the origin and cell will be examined.
48  // If the origin is outside the cell, a global n-squared search is
49  // triggered.
50 
51  // SETTING UP RESCUE
52 
53  // Visit all boundaries and mark the cell next to the boundary.
54 
55  if (debug)
56  {
57  InfoInFunction << "Setting up rescue" << endl;
58  }
59 
60  List<bool> boundaryCell(fromCells.size(), false);
61 
62  // Set reference to boundary
63  const polyPatchList& patchesFrom = fromMesh_.boundaryMesh();
64 
65  forAll(patchesFrom, patchi)
66  {
67  // Get reference to cells next to the boundary
68  const labelUList& bCells = patchesFrom[patchi].faceCells();
69 
70  forAll(bCells, facei)
71  {
72  boundaryCell[bCells[facei]] = true;
73  }
74  }
75 
76  treeBoundBox meshBb(fromPoints);
77 
78  scalar typDim = meshBb.avgDim()/(2.0*cbrt(scalar(fromCells.size())));
79 
80  treeBoundBox shiftedBb
81  (
82  meshBb.min(),
83  meshBb.max() + vector(typDim, typDim, typDim)
84  );
85 
86  if (debug)
87  {
88  Info<< "\nMesh\n"
89  << " bounding box : " << meshBb << nl
90  << " bounding box (shifted) : " << shiftedBb << nl
91  << " typical dimension :" << shiftedBb.typDim() << endl;
92  }
93 
94  indexedOctree<treeDataCell> oc
95  (
96  treeDataCell(false, fromMesh_, polyMesh::CELL_TETS),
97  shiftedBb, // overall bounding box
98  8, // maxLevel
99  10, // leafsize
100  6.0 // duplicity
101  );
102 
103  if (debug)
104  {
105  oc.print(Pout, false, 0);
106  }
107 
108  cellAddresses
109  (
110  cellAddressing_,
111  toMesh_.cellCentres(),
112  fromMesh_,
113  boundaryCell,
114  oc
115  );
116 
117  forAll(toMesh_.boundaryMesh(), patchi)
118  {
119  const polyPatch& toPatch = toMesh_.boundaryMesh()[patchi];
120 
121  if (cuttingPatches_.found(toPatch.name()))
122  {
123  boundaryAddressing_[patchi].setSize(toPatch.size());
124 
125  cellAddresses
126  (
127  boundaryAddressing_[patchi],
128  toPatch.faceCentres(),
129  fromMesh_,
130  boundaryCell,
131  oc
132  );
133  }
134  else if
135  (
136  patchMap_.found(toPatch.name())
137  && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
138  )
139  {
140  const polyPatch& fromPatch = fromMesh_.boundaryMesh()
141  [
142  fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
143  ];
144 
145  if (fromPatch.empty())
146  {
148  << "Source patch " << fromPatch.name()
149  << " has no faces. Not performing mapping for it."
150  << endl;
151  boundaryAddressing_[patchi].setSize(toPatch.size());
152  boundaryAddressing_[patchi] = -1;
153  }
154  else
155  {
156  treeBoundBox wallBb(fromPatch.localPoints());
157  scalar typDim =
158  wallBb.avgDim()/(2.0*sqrt(scalar(fromPatch.size())));
159 
160  treeBoundBox shiftedBb
161  (
162  wallBb.min(),
163  wallBb.max() + vector(typDim, typDim, typDim)
164  );
165 
166  // Note: allow more levels than in meshSearch. Assume patch
167  // is not as big as all boundary faces
168  indexedOctree<treeDataFace> oc
169  (
170  treeDataFace(false, fromPatch),
171  shiftedBb, // overall search domain
172  12, // maxLevel
173  10, // leafsize
174  6.0 // duplicity
175  );
176 
177  const vectorField::subField centresToBoundary =
178  toPatch.faceCentres();
179 
180  boundaryAddressing_[patchi].setSize(toPatch.size());
181 
182  const scalar distSqr = sqr(wallBb.mag());
183 
184  forAll(toPatch, toi)
185  {
186  boundaryAddressing_[patchi][toi] = oc.findNearest
187  (
188  centresToBoundary[toi],
189  distSqr
190  ).index();
191  }
192  }
193  }
194  }
195 
196  if (debug)
197  {
199  << "Finished calculating mesh-to-mesh cell addressing" << endl;
200  }
201 }
202 
203 
204 void Foam::meshToMesh0::cellAddresses
205 (
206  labelList& cellAddressing_,
207  const pointField& points,
208  const fvMesh& fromMesh,
209  const List<bool>& boundaryCell,
210  const indexedOctree<treeDataCell>& oc
211 ) const
212 {
213  // The implemented search method is a simple neighbour array search.
214  // It starts from a cell zero, searches its neighbours and finds one
215  // which is nearer to the target point than the current position.
216  // The location of the "current position" is reset to that cell and
217  // search through the neighbours continues. The search is finished
218  // when all the neighbours of the cell are farther from the target
219  // point than the current cell
220 
221  // Set curCell label to zero (start)
222  label curCell = 0;
223 
224  // Set reference to cell to cell addressing
225  const vectorField& centresFrom = fromMesh.cellCentres();
226  const labelListList& cc = fromMesh.cellCells();
227 
228  forAll(points, toi)
229  {
230  // Pick up target position
231  const vector& p = points[toi];
232 
233  // Set the sqr-distance
234  scalar distSqr = magSqr(p - centresFrom[curCell]);
235 
236  bool closer;
237 
238  do
239  {
240  closer = false;
241 
242  // Set the current list of neighbouring cells
243  const labelList& neighbours = cc[curCell];
244 
245  forAll(neighbours, ni)
246  {
247  const scalar curDistSqr =
248  magSqr(p - centresFrom[neighbours[ni]]);
249 
250  // Search through all the neighbours.
251  // If the cell is closer, reset current cell and distance
252  if (curDistSqr < (1 - small)*distSqr)
253  {
254  curCell = neighbours[ni];
255  distSqr = curDistSqr;
256  closer = true; // A closer neighbour has been found
257  }
258  }
259  } while (closer);
260 
261  cellAddressing_[toi] = -1;
262 
263  // Check point is actually in the nearest cell
264  if (fromMesh.pointInCell(p, curCell))
265  {
266  cellAddressing_[toi] = curCell;
267  }
268  else
269  {
270  // If curCell is a boundary cell then the point maybe either outside
271  // the domain or in an other region of the doamin, either way use
272  // the octree search to find it.
273  if (boundaryCell[curCell])
274  {
275  cellAddressing_[toi] = oc.findInside(p);
276  if (cellAddressing_[toi] != -1)
277  {
278  curCell = cellAddressing_[toi];
279  }
280  }
281  else
282  {
283  // If not on the boundary search the neighbours
284  bool found = false;
285 
286  // Set the current list of neighbouring cells
287  const labelList& neighbours = cc[curCell];
288 
289  forAll(neighbours, ni)
290  {
291  // Search through all the neighbours.
292  // If point is in neighbour reset current cell
293  if (fromMesh.pointInCell(p, neighbours[ni]))
294  {
295  cellAddressing_[toi] = neighbours[ni];
296  found = true;
297  break;
298  }
299  }
300 
301  if (!found)
302  {
303  // If still not found search the neighbour-neighbours
304 
305  // Set the current list of neighbouring cells
306  const labelList& neighbours = cc[curCell];
307 
308  forAll(neighbours, ni)
309  {
310  // Set the current list of neighbour-neighbouring cells
311  const labelList& nn = cc[neighbours[ni]];
312 
313  forAll(nn, ni)
314  {
315  // Search through all the neighbours.
316  // If point is in neighbour reset current cell
317  if (fromMesh.pointInCell(p, nn[ni]))
318  {
319  cellAddressing_[toi] = nn[ni];
320  found = true;
321  break;
322  }
323  }
324  if (found) break;
325  }
326  }
327 
328  if (!found)
329  {
330  // Still not found so use the octree
331  cellAddressing_[toi] = oc.findInside(p);
332 
333  if (cellAddressing_[toi] != -1)
334  {
335  curCell = cellAddressing_[toi];
336  }
337  }
338  }
339  }
340  }
341 }
342 
343 
344 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:45
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const cellList & cells() const
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
UList< label > labelUList
Definition: UList.H:65
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
dimensionedScalar cbrt(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition: labelList.H:56
const vectorField & cellCentres() const
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(1e-3))
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:260
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
messageStream Info
SubField< vector > subField
Declare type of subField.
Definition: Field.H:100
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< cell > cellList
list of cells
Definition: cellList.H:42
#define InfoInFunction
Report an information message using Foam::Info.