treeDataCell.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-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2022 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 \*---------------------------------------------------------------------------*/
28 
29 #include "treeDataCell.H"
30 #include "indexedOctree.H"
31 #include "polyMesh.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeName(treeDataCell);
38 }
39 
40 
41 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // Bound boxes corresponding to specified cells
47 template<class ElementIds>
49 (
50  const primitiveMesh& mesh,
51  const ElementIds& elemIds
52 )
53 {
54  treeBoundBoxList bbs(elemIds.size());
55 
56  if (mesh.hasCellPoints())
57  {
59  (
60  elemIds.cbegin(),
61  elemIds.cend(),
62  bbs.begin(),
63  [&](label celli)
64  {
65  return treeBoundBox(mesh.points(), mesh.cellPoints(celli));
66  }
67  );
68  }
69  else
70  {
72  (
73  elemIds.cbegin(),
74  elemIds.cend(),
75  bbs.begin(),
76  [&](label celli)
77  {
78  return treeBoundBox(mesh.cells()[celli].box(mesh));
79  }
80  );
81  }
82 
83  return bbs;
84 }
85 
86 } // End namespace Foam
87 
88 
89 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
90 
93 {
94  // All cells
95  return boxesImpl(mesh, labelRange(mesh.nCells()));
96 }
97 
98 
101 (
102  const primitiveMesh& mesh,
103  const labelUList& cellIds
104 )
105 {
106  return boxesImpl(mesh, cellIds);
107 }
108 
109 
112 (
113  const primitiveMesh& mesh,
114  const labelUList& cellIds
115 )
116 {
117  treeBoundBox bb;
118 
119  if (mesh.hasCellPoints())
120  {
121  for (const label celli : cellIds)
122  {
123  bb.add(mesh.points(), mesh.cellPoints(celli));
124  }
125  }
126  else
127  {
128  for (const label celli : cellIds)
129  {
130  bb.add(mesh.cells()[celli].box(mesh));
131  }
132  }
133 
134  return bb;
135 }
136 
137 
138 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
139 
140 inline Foam::treeBoundBox
141 Foam::treeDataCell::getBounds(const label index) const
142 {
143  return treeBoundBox(mesh_.cellBb(objectIndex(index)));
144 }
145 
146 
147 void Foam::treeDataCell::update()
148 {
149  if (cacheBb_)
150  {
151  if (useSubset_)
152  {
153  bbs_ = treeDataCell::boxes(mesh_, cellLabels_);
154  }
155  else
156  {
157  bbs_ = treeDataCell::boxes(mesh_);
158  }
159  }
160 }
161 
162 
163 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
164 
166 (
167  const bool cacheBb,
168  const polyMesh& mesh,
169  const polyMesh::cellDecomposition decompMode
170 )
171 :
172  mesh_(mesh),
173  cellLabels_(),
174  useSubset_(false),
175  cacheBb_(cacheBb),
176  decompMode_(decompMode)
177 {
178  update();
179 }
180 
181 
183 (
184  const bool cacheBb,
185  const polyMesh& mesh,
186  const labelUList& cellLabels,
187  const polyMesh::cellDecomposition decompMode
188 )
189 :
190  mesh_(mesh),
191  cellLabels_(cellLabels),
192  useSubset_(true),
193  cacheBb_(cacheBb),
194  decompMode_(decompMode)
195 {
196  update();
197 }
198 
199 
201 (
202  const bool cacheBb,
203  const polyMesh& mesh,
204  labelList&& cellLabels,
205  const polyMesh::cellDecomposition decompMode
206 )
207 :
208  mesh_(mesh),
209  cellLabels_(std::move(cellLabels)),
210  useSubset_(true),
211  cacheBb_(cacheBb),
212  decompMode_(decompMode)
213 {
214  update();
215 }
216 
217 
218 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
219 
221 Foam::treeDataCell::bounds(const labelUList& indices) const
222 {
223  if (useSubset_)
224  {
225  treeBoundBox bb;
226 
227  if (mesh_.hasCellPoints())
228  {
229  for (const label index : indices)
230  {
231  const label celli = cellLabels_[index];
232 
233  bb.add(mesh_.points(), mesh_.cellPoints(celli));
234  }
235  }
236  else
237  {
238  for (const label index : indices)
239  {
240  const label celli = cellLabels_[index];
241 
242  bb.add(mesh_.cells()[celli].box(mesh_));
243  }
244  }
245 
246  return bb;
247  }
248 
249  return treeDataCell::bounds(mesh_, indices);
250 }
251 
252 
254 {
255  if (useSubset_)
256  {
257  return tmp<pointField>::New(mesh_.cellCentres(), cellLabels_);
258  }
259 
260  return mesh_.cellCentres();
261 }
262 
263 
265 (
266  const label index,
267  const treeBoundBox& searchBox
268 ) const
269 {
270  return
271  (
272  cacheBb_
273  ? searchBox.overlaps(bbs_[index])
274  : searchBox.overlaps(getBounds(index))
275  );
276 }
277 
278 
280 (
281  const label index,
282  const point& sample
283 ) const
284 {
285  return mesh_.pointInCell(sample, objectIndex(index), decompMode_);
286 }
287 
288 
289 // * * * * * * * * * * * * * * * * Searching * * * * * * * * * * * * * * * * //
290 
292 (
294 )
295 :
296  tree_(tree)
297 {}
298 
299 
301 (
303 )
304 :
305  tree_(tree)
306 {}
307 
308 
310 (
311  const labelUList& indices,
312  const point& sample,
313 
314  scalar& nearestDistSqr,
315  label& minIndex,
316  point& nearestPoint
317 ) const
318 {
319  for (const label index : indices)
320  {
321  const point& pt = centre(index);
322 
323  const scalar distSqr = sample.distSqr(pt);
324 
325  if (distSqr < nearestDistSqr)
326  {
327  nearestDistSqr = distSqr;
328  minIndex = index;
329  nearestPoint = pt;
330  }
331  }
332 }
333 
334 
335 void Foam::treeDataCell::findNearestOp::operator()
336 (
337  const labelUList& indices,
338  const point& sample,
339 
340  scalar& nearestDistSqr,
341  label& minIndex,
342  point& nearestPoint
343 ) const
344 {
345  tree_.shapes().findNearest
346  (
347  indices,
348  sample,
349  nearestDistSqr,
350  minIndex,
351  nearestPoint
352  );
353 }
354 
355 
356 void Foam::treeDataCell::findNearestOp::operator()
357 (
358  const labelUList& indices,
359  const linePointRef& ln,
360 
361  treeBoundBox& tightest,
362  label& minIndex,
363  point& linePoint,
364  point& nearestPoint
365 ) const
366 {
368 }
369 
370 
371 bool Foam::treeDataCell::findIntersectOp::operator()
372 (
373  const label index,
374  const point& start,
375  const point& end,
376  point& intersectionPoint
377 ) const
378 {
379  const treeDataCell& shape = tree_.shapes();
380 
381  // Do quick rejection test
382  if (shape.cacheBb_)
383  {
384  const treeBoundBox& cellBb = shape.bbs_[index];
385 
386  if ((cellBb.posBits(start) & cellBb.posBits(end)) != 0)
387  {
388  // start and end in same block outside of cellBb.
389  return false;
390  }
391  }
392  else
393  {
394  const treeBoundBox cellBb = shape.getBounds(index);
395 
396  if ((cellBb.posBits(start) & cellBb.posBits(end)) != 0)
397  {
398  // start and end in same block outside of cellBb.
399  return false;
400  }
401  }
402 
403 
404  // Do intersection with all faces of cell
405  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
406 
407  // Disable picking up intersections behind us.
408  const scalar oldTol = intersection::setPlanarTol(0);
409 
410  const vector dir(end - start);
411  scalar minDistSqr = magSqr(dir);
412  bool hasMin = false;
413 
414  const label celli = shape.objectIndex(index);
415 
416  for (const label facei : shape.mesh().cells()[celli])
417  {
418  const face& f = shape.mesh().faces()[facei];
419 
420  pointHit inter = f.ray
421  (
422  start,
423  dir,
424  shape.mesh().points(),
426  );
427 
428  if (inter.hit() && sqr(inter.distance()) <= minDistSqr)
429  {
430  // Note: no extra test on whether intersection is in front of us
431  // since using half_ray AND zero tolerance. (note that tolerance
432  // is used to look behind us)
433  minDistSqr = sqr(inter.distance());
434  intersectionPoint = inter.point();
435  hasMin = true;
436  }
437  }
438 
439  // Restore picking tolerance
441 
442  return hasMin;
443 }
444 
445 
446 // ************************************************************************* //
A line primitive.
Definition: line.H:52
tmp< pointField > centres() const
Representative point cloud for contained shapes. One point per shape, corresponding to the cell centr...
Definition: treeDataCell.C:246
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
dimensionedSymmTensor sqr(const dimensionedVector &dv)
A range or interval of labels defined by a start and a size.
Definition: labelRange.H:52
const polyMesh & mesh() const noexcept
Reference to the supporting mesh.
Definition: treeDataCell.H:279
const cellList & cells() const
bool overlaps(const label index, const treeBoundBox &searchBox) const
Does (bb of) shape at index overlap searchBox.
Definition: treeDataCell.C:258
static treeBoundBoxList boxes(const primitiveMesh &mesh)
Calculate and return bounding boxes for all mesh cells.
Definition: treeDataCell.C:85
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:103
static treeBoundBoxList boxesImpl(const primitiveMesh &mesh, const ElementIds &elemIds)
Definition: treeDataFace.C:44
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
direction posBits(const point &pt) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:486
static scalar setPlanarTol(const scalar t)
Set the planar tolerance, returning the previous value.
Definition: intersection.H:101
defineTypeName(manifoldCellsMeshObject)
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:323
findIntersectOp(const indexedOctree< treeDataCell > &tree)
Definition: treeDataCell.C:294
dynamicFvMesh & mesh
Tree tree(triangles.begin(), triangles.end())
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
static treeBoundBox bounds(const primitiveMesh &mesh, const labelUList &cellIds)
Return bounding box of specified mesh cells.
Definition: treeDataCell.C:105
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
void findNearest(const labelUList &indices, const point &sample, scalar &nearestDistSqr, label &nearestIndex, point &nearestPoint) const
Calculates nearest (to sample) point in shape.
Definition: treeDataCell.C:303
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:1237
bool hasCellPoints() const noexcept
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
mesh update()
labelList f(nPoints)
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e...
Definition: treeDataCell.H:53
treeDataCell(const bool cacheBb, const polyMesh &mesh, polyMesh::cellDecomposition decompMode)
Construct from mesh, using all cells in mesh.
Definition: treeDataCell.C:159
findNearestOp(const indexedOctree< treeDataCell > &tree)
Definition: treeDataCell.C:285
vector point
Point is a vector.
Definition: point.H:37
Non-pointer based hierarchical recursive searching.
label objectIndex(const label index) const
Map from shape index to original (non-subset) cell label.
Definition: treeDataCell.H:318
label nCells() const noexcept
Number of mesh cells.
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: HashPtrTable.H:50
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const labelListList & cellPoints() const
bool contains(const label index, const point &sample) const
Does shape at index contain sample.
Definition: treeDataCell.C:273
PointHit< point > pointHit
A PointHit with a 3D point.
Definition: pointHit.H:43
Namespace for OpenFOAM.