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