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-2015 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 
41 {
42  if (debug)
43  {
44  Info<< "meshToMesh0::calculateAddressing() : "
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  Info<< "meshToMesh0::calculateAddressing() : "
66  << "Setting up rescue" << endl;
67  }
68 
69  List<bool> boundaryCell(fromCells.size(), false);
70 
71  // set reference to boundary
72  const polyPatchList& patchesFrom = fromMesh_.boundaryMesh();
73 
74  forAll(patchesFrom, patchI)
75  {
76  // get reference to cells next to the boundary
77  const labelUList& bCells = patchesFrom[patchI].faceCells();
78 
79  forAll(bCells, faceI)
80  {
81  boundaryCell[bCells[faceI]] = true;
82  }
83  }
84 
85  treeBoundBox meshBb(fromPoints);
86 
87  scalar typDim = meshBb.avgDim()/(2.0*cbrt(scalar(fromCells.size())));
88 
89  treeBoundBox shiftedBb
90  (
91  meshBb.min(),
92  meshBb.max() + vector(typDim, typDim, typDim)
93  );
94 
95  if (debug)
96  {
97  Info<< "\nMesh" << endl;
98  Info<< " bounding box : " << meshBb << endl;
99  Info<< " bounding box (shifted) : " << shiftedBb << endl;
100  Info<< " typical dimension :" << shiftedBb.typDim() << endl;
101  }
102 
104  (
106  shiftedBb, // overall bounding box
107  8, // maxLevel
108  10, // leafsize
109  6.0 // duplicity
110  );
111 
112  if (debug)
113  {
114  oc.print(Pout, false, 0);
115  }
116 
118  (
121  fromMesh_,
122  boundaryCell,
123  oc
124  );
125 
127  {
128  const polyPatch& toPatch = toMesh_.boundaryMesh()[patchi];
129 
130  if (cuttingPatches_.found(toPatch.name()))
131  {
132  boundaryAddressing_[patchi].setSize(toPatch.size());
133 
135  (
137  toPatch.faceCentres(),
138  fromMesh_,
139  boundaryCell,
140  oc
141  );
142  }
143  else if
144  (
145  patchMap_.found(toPatch.name())
146  && fromMeshPatches_.found(patchMap_.find(toPatch.name())())
147  )
148  {
149  const polyPatch& fromPatch = fromMesh_.boundaryMesh()
150  [
151  fromMeshPatches_.find(patchMap_.find(toPatch.name())())()
152  ];
153 
154  if (fromPatch.empty())
155  {
157  << "Source patch " << fromPatch.name()
158  << " has no faces. Not performing mapping for it."
159  << endl;
161  }
162  else
163  {
164  treeBoundBox wallBb(fromPatch.localPoints());
165  scalar typDim =
166  wallBb.avgDim()/(2.0*sqrt(scalar(fromPatch.size())));
167 
168  treeBoundBox shiftedBb
169  (
170  wallBb.min(),
171  wallBb.max() + vector(typDim, typDim, typDim)
172  );
173 
174  // Note: allow more levels than in meshSearch. Assume patch
175  // is not as big as all boundary faces
177  (
178  treeDataFace(false, fromPatch),
179  shiftedBb, // overall search domain
180  12, // maxLevel
181  10, // leafsize
182  6.0 // duplicity
183  );
184 
185  const vectorField::subField centresToBoundary =
186  toPatch.faceCentres();
187 
188  boundaryAddressing_[patchi].setSize(toPatch.size());
189 
190  scalar distSqr = sqr(wallBb.mag());
191 
192  forAll(toPatch, toi)
193  {
195  (
196  centresToBoundary[toi],
197  distSqr
198  ).index();
199  }
200  }
201  }
202  }
203 
204  if (debug)
205  {
206  Info<< "meshToMesh0::calculateAddressing() : "
207  << "finished calculating mesh-to-mesh cell addressing" << endl;
208  }
209 }
210 
211 
213 (
214  labelList& cellAddressing_,
215  const pointField& points,
216  const fvMesh& fromMesh,
217  const List<bool>& boundaryCell,
219 ) const
220 {
221  // the implemented search method is a simple neighbour array search.
222  // It starts from a cell zero, searches its neighbours and finds one
223  // which is nearer to the target point than the current position.
224  // The location of the "current position" is reset to that cell and
225  // search through the neighbours continues. The search is finished
226  // when all the neighbours of the cell are farther from the target
227  // point than the current cell
228 
229  // set curCell label to zero (start)
230  label curCell = 0;
231 
232  // set reference to cell to cell addressing
233  const vectorField& centresFrom = fromMesh.cellCentres();
234  const labelListList& cc = fromMesh.cellCells();
235 
236  forAll(points, toI)
237  {
238  // pick up target position
239  const vector& p = points[toI];
240 
241  // set the sqr-distance
242  scalar distSqr = magSqr(p - centresFrom[curCell]);
243 
244  bool closer;
245 
246  do
247  {
248  closer = false;
249 
250  // set the current list of neighbouring cells
251  const labelList& neighbours = cc[curCell];
252 
253  forAll(neighbours, nI)
254  {
255  scalar curDistSqr =
256  magSqr(p - centresFrom[neighbours[nI]]);
257 
258  // search through all the neighbours.
259  // If the cell is closer, reset current cell and distance
260  if (curDistSqr < (1 - SMALL)*distSqr)
261  {
262  curCell = neighbours[nI];
263  distSqr = curDistSqr;
264  closer = true; // a closer neighbour has been found
265  }
266  }
267  } while (closer);
268 
269  cellAddressing_[toI] = -1;
270 
271  // Check point is actually in the nearest cell
272  if (fromMesh.pointInCell(p, curCell))
273  {
274  cellAddressing_[toI] = curCell;
275  }
276  else
277  {
278  // If curCell is a boundary cell then the point maybe either outside
279  // the domain or in an other region of the doamin, either way use
280  // the octree search to find it.
281  if (boundaryCell[curCell])
282  {
283  cellAddressing_[toI] = oc.findInside(p);
284  }
285  else
286  {
287  // If not on the boundary search the neighbours
288  bool found = false;
289 
290  // set the current list of neighbouring cells
291  const labelList& neighbours = cc[curCell];
292 
293  forAll(neighbours, nI)
294  {
295  // search through all the neighbours.
296  // If point is in neighbour reset current cell
297  if (fromMesh.pointInCell(p, neighbours[nI]))
298  {
299  cellAddressing_[toI] = neighbours[nI];
300  found = true;
301  break;
302  }
303  }
304 
305  if (!found)
306  {
307  // If still not found search the neighbour-neighbours
308 
309  // set the current list of neighbouring cells
310  const labelList& neighbours = cc[curCell];
311 
312  forAll(neighbours, nI)
313  {
314  // set the current list of neighbour-neighbouring cells
315  const labelList& nn = cc[neighbours[nI]];
316 
317  forAll(nn, nI)
318  {
319  // search through all the neighbours.
320  // If point is in neighbour reset current cell
321  if (fromMesh.pointInCell(p, nn[nI]))
322  {
323  cellAddressing_[toI] = nn[nI];
324  found = true;
325  break;
326  }
327  }
328  if (found) break;
329  }
330  }
331 
332  if (!found)
333  {
334  // Still not found so us the octree
335  cellAddressing_[toI] = oc.findInside(p);
336  }
337  }
338  }
339  }
340 }
341 
342 
343 // ************************************************************************* //
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::boundBox::mag
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
p
p
Definition: pEqn.H:62
Foam::meshToMesh0::cellAddressing_
labelList cellAddressing_
Cell addressing.
Definition: meshToMesh0.H:85
SubField.H
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
meshToMesh0.H
Foam::meshToMesh0::cuttingPatches_
HashTable< label > cuttingPatches_
toMesh patch labels which cut the from-mesh
Definition: meshToMesh0.H:82
Foam::meshToMesh0::fromMesh_
const fvMesh & fromMesh_
Definition: meshToMesh0.H:69
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::meshToMesh0::patchMap_
HashTable< word > patchMap_
Patch map.
Definition: meshToMesh0.H:79
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
indexedOctree.H
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::polyMesh::pointInCell
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1292
meshBb
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
treeDataFace.H
Foam::meshToMesh0::calcAddressing
void calcAddressing()
Definition: calculateMeshToMesh0Addressing.C:40
Foam::SubField
Pre-declare related SubField type.
Definition: Field.H:61
Foam::treeDataCell
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e....
Definition: treeDataCell.H:54
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::treeBoundBox::typDim
scalar typDim() const
Typical dimension length,height,width.
Definition: treeBoundBoxI.H:57
Foam::primitiveMesh::cellCells
const labelListList & cellCells() const
Definition: primitiveMeshCellCells.C:100
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataTriSurface.H:47
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::indexedOctree::findInside
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
Definition: indexedOctree.C:2828
Foam::indexedOctree::print
void print(prefixOSstream &, const bool printContents, const label) const
Print tree. Either print all indices (printContent = true) or.
Definition: indexedOctree.C:2984
Foam::meshToMesh0::boundaryAddressing_
labelListList boundaryAddressing_
Boundary addressing.
Definition: meshToMesh0.H:88
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::meshToMesh0::toMesh_
const fvMesh & toMesh_
Definition: meshToMesh0.H:70
treeDataCell.H
Foam::HashTable::find
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:211
Foam::polyMesh::CELL_TETS
@ CELL_TETS
Definition: polyMesh.H:107
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:308
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::treeDataFace
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
patchi
label patchi
Definition: getPatchFieldScalar.H:1
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::boundBox::avgDim
scalar avgDim() const
Average length/height/width dimension.
Definition: boundBoxI.H:114
Foam::meshToMesh0::cellAddresses
void cellAddresses(labelList &cells, const pointField &points, const fvMesh &fromMesh, const List< bool > &boundaryCell, const indexedOctree< treeDataCell > &oc) const
Definition: calculateMeshToMesh0Addressing.C:213
Foam::indexedOctree::findNearest
void findNearest(const label nodeI, const linePointRef &ln, treeBoundBox &tightest, label &nearestShapeI, point &linePoint, point &nearestPoint, const FindNearestOp &fnOp) const
Find nearest point to line.
Definition: indexedOctree.C:590
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:153
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::meshToMesh0::fromMeshPatches_
HashTable< label > fromMeshPatches_
fromMesh patch labels
Definition: meshToMesh0.H:73