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 | 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 \*---------------------------------------------------------------------------*/
25 
26 #include "treeDataCell.H"
27 #include "indexedOctree.H"
28 #include "polyMesh.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 defineTypeNameAndDebug(treeDataCell, 0);
35 }
36 
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 Foam::treeBoundBox Foam::treeDataCell::calcCellBb(const label celli) const
41 {
42  const cellList& cells = mesh_.cells();
43  const faceList& faces = mesh_.faces();
44  const pointField& points = mesh_.points();
45 
46  treeBoundBox cellBb
47  (
48  vector(GREAT, GREAT, GREAT),
49  vector(-GREAT, -GREAT, -GREAT)
50  );
51 
52  const cell& cFaces = cells[celli];
53 
54  forAll(cFaces, cFacei)
55  {
56  const face& f = faces[cFaces[cFacei]];
57 
58  forAll(f, fp)
59  {
60  const point& p = points[f[fp]];
61 
62  cellBb.min() = min(cellBb.min(), p);
63  cellBb.max() = max(cellBb.max(), p);
64  }
65  }
66  return cellBb;
67 }
68 
69 
70 void Foam::treeDataCell::update()
71 {
72  if (cacheBb_)
73  {
74  bbs_.setSize(cellLabels_.size());
75 
76  forAll(cellLabels_, i)
77  {
78  bbs_[i] = calcCellBb(cellLabels_[i]);
79  }
80  }
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
87 (
88  const bool cacheBb,
89  const polyMesh& mesh,
90  const labelUList& cellLabels,
91  const polyMesh::cellDecomposition decompMode
92 )
93 :
94  mesh_(mesh),
95  cellLabels_(cellLabels),
96  cacheBb_(cacheBb),
97  decompMode_(decompMode)
98 {
99  update();
100 }
101 
102 
104 (
105  const bool cacheBb,
106  const polyMesh& mesh,
107  const Xfer<labelList>& cellLabels,
108  const polyMesh::cellDecomposition decompMode
109 )
110 :
111  mesh_(mesh),
112  cellLabels_(cellLabels),
113  cacheBb_(cacheBb),
114  decompMode_(decompMode)
115 {
116  update();
117 }
118 
119 
121 (
122  const bool cacheBb,
123  const polyMesh& mesh,
124  const polyMesh::cellDecomposition decompMode
125 )
126 :
127  mesh_(mesh),
128  cellLabels_(identity(mesh_.nCells())),
129  cacheBb_(cacheBb),
130  decompMode_(decompMode)
131 {
132  update();
133 }
134 
135 
137 (
138  const indexedOctree<treeDataCell>& tree
139 )
140 :
141  tree_(tree)
142 {}
143 
144 
146 (
147  const indexedOctree<treeDataCell>& tree
148 )
149 :
150  tree_(tree)
151 {}
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
157 {
158  pointField cc(cellLabels_.size());
159 
160  forAll(cellLabels_, i)
161  {
162  cc[i] = mesh_.cellCentres()[cellLabels_[i]];
163  }
164 
165  return cc;
166 }
167 
168 
170 (
171  const label index,
172  const treeBoundBox& cubeBb
173 ) const
174 {
175  if (cacheBb_)
176  {
177  return cubeBb.overlaps(bbs_[index]);
178  }
179  else
180  {
181  return cubeBb.overlaps(calcCellBb(cellLabels_[index]));
182  }
183 }
184 
185 
187 (
188  const label index,
189  const point& sample
190 ) const
191 {
192  return mesh_.pointInCell(sample, cellLabels_[index], decompMode_);
193 }
194 
195 
196 void Foam::treeDataCell::findNearestOp::operator()
197 (
198  const labelUList& indices,
199  const point& sample,
200 
201  scalar& nearestDistSqr,
202  label& minIndex,
203  point& nearestPoint
204 ) const
205 {
206  const treeDataCell& shape = tree_.shapes();
207 
208  forAll(indices, i)
209  {
210  label index = indices[i];
211  label celli = shape.cellLabels()[index];
212  scalar distSqr = magSqr(sample - shape.mesh().cellCentres()[celli]);
213 
214  if (distSqr < nearestDistSqr)
215  {
216  nearestDistSqr = distSqr;
217  minIndex = index;
218  nearestPoint = shape.mesh().cellCentres()[celli];
219  }
220  }
221 }
222 
223 
224 void Foam::treeDataCell::findNearestOp::operator()
225 (
226  const labelUList& indices,
227  const linePointRef& ln,
228 
229  treeBoundBox& tightest,
230  label& minIndex,
231  point& linePoint,
232  point& nearestPoint
233 ) const
234 {
236 }
237 
238 
239 bool Foam::treeDataCell::findIntersectOp::operator()
240 (
241  const label index,
242  const point& start,
243  const point& end,
244  point& intersectionPoint
245 ) const
246 {
247  const treeDataCell& shape = tree_.shapes();
248 
249  // Do quick rejection test
250  if (shape.cacheBb_)
251  {
252  const treeBoundBox& cellBb = shape.bbs_[index];
253 
254  if ((cellBb.posBits(start) & cellBb.posBits(end)) != 0)
255  {
256  // start and end in same block outside of cellBb.
257  return false;
258  }
259  }
260  else
261  {
262  const treeBoundBox cellBb = shape.calcCellBb(shape.cellLabels_[index]);
263 
264  if ((cellBb.posBits(start) & cellBb.posBits(end)) != 0)
265  {
266  // start and end in same block outside of cellBb.
267  return false;
268  }
269  }
270 
271 
272  // Do intersection with all faces of cell
273  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
274 
275  // Disable picking up intersections behind us.
276  scalar oldTol = intersection::setPlanarTol(0.0);
277 
278  const cell& cFaces = shape.mesh_.cells()[shape.cellLabels_[index]];
279 
280  const vector dir(end - start);
281  scalar minDistSqr = magSqr(dir);
282  bool hasMin = false;
283 
284  forAll(cFaces, i)
285  {
286  const face& f = shape.mesh_.faces()[cFaces[i]];
287 
288  pointHit inter = f.ray
289  (
290  start,
291  dir,
292  shape.mesh_.points(),
294  );
295 
296  if (inter.hit() && sqr(inter.distance()) <= minDistSqr)
297  {
298  // Note: no extra test on whether intersection is in front of us
299  // since using half_ray AND zero tolerance. (note that tolerance
300  // is used to look behind us)
301  minDistSqr = sqr(inter.distance());
302  intersectionPoint = inter.hitPoint();
303  hasMin = true;
304  }
305  }
306 
307  // Restore picking tolerance
309 
310  return hasMin;
311 }
312 
313 
314 // ************************************************************************* //
bool overlaps(const label index, const treeBoundBox &sampleBb) const
Does (bb of) shape at index overlap bb.
Definition: treeDataCell.C:170
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A line primitive.
Definition: line.H:56
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const polyMesh & mesh() const
Definition: treeDataCell.H:175
List< face > faceList
Definition: faceListFwd.H:43
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
static const Form min
Definition: VectorSpace.H:113
const cellList & cells() const
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:98
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
static scalar setPlanarTol(const scalar t)
Set the planar tolerance, returning the previous value.
Definition: intersection.H:91
findIntersectOp(const indexedOctree< treeDataCell > &tree)
Definition: treeDataCell.C:146
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
direction posBits(const point &) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:464
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const vectorField & cellCentres() const
treeDataCell(const bool cacheBb, const polyMesh &, const labelUList &, const polyMesh::cellDecomposition decompMode)
Construct from mesh and subset of cells.
Definition: treeDataCell.C:87
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
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:720
dimensioned< scalar > magSqr(const dimensioned< Type > &)
bool contains(const label index, const point &sample) const
Does shape at index contain sample.
Definition: treeDataCell.C:187
defineTypeNameAndDebug(combustionModel, 0)
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e...
Definition: treeDataCell.H:54
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
pointHit ray(const point &p, const vector &n, const pointField &, const intersection::algorithm alg=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
Return potential intersection with face with a ray starting.
findNearestOp(const indexedOctree< treeDataCell > &tree)
Definition: treeDataCell.C:137
vector point
Point is a vector.
Definition: point.H:41
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:126
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
pointField shapePoints() const
Get representative point cloud for all shapes inside.
Definition: treeDataCell.C:156
bool hit() const
Is there a hit.
Definition: PointHit.H:120
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
const labelList & cellLabels() const
Definition: treeDataCell.H:170
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
List< cell > cellList
list of cells
Definition: cellList.H:42
Namespace for OpenFOAM.