triSurfaceSearch.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 | Copyright (C) 2015 OpenCFD Ltd.
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 "triSurfaceSearch.H"
27 #include "triSurface.H"
28 #include "PatchTools.H"
29 #include "volumeType.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
34 (
35  const pointIndexHit& currHit,
37  const vector& lineVec
38 ) const
39 {
40  // Classify the hit
41  label nearType = -1;
42  label nearLabel = -1;
43 
44  const labelledTri& f = surface()[currHit.index()];
45 
46  f.nearestPointClassify
47  (
48  currHit.hitPoint(),
49  surface().points(),
50  nearType,
51  nearLabel
52  );
53 
54  if (nearType == triPointRef::POINT)
55  {
56  // near point
57 
58  const label nearPointI = f[nearLabel];
59 
60  const labelList& pointFaces =
61  surface().pointFaces()[surface().meshPointMap()[nearPointI]];
62 
63  forAll(pointFaces, pI)
64  {
65  const label pointFaceI = pointFaces[pI];
66 
67  if (pointFaceI != currHit.index())
68  {
69  forAll(hits, hI)
70  {
71  const pointIndexHit& hit = hits[hI];
72 
73  if (hit.index() == pointFaceI)
74  {
75  return false;
76  }
77  }
78  }
79  }
80  }
81  else if (nearType == triPointRef::EDGE)
82  {
83  // near edge
84  // check if the other face of the edge is already hit
85 
86  const labelList& fEdges = surface().faceEdges()[currHit.index()];
87 
88  const label edgeI = fEdges[nearLabel];
89 
90  const labelList& edgeFaces = surface().edgeFaces()[edgeI];
91 
92  forAll(edgeFaces, fI)
93  {
94  const label edgeFaceI = edgeFaces[fI];
95 
96  if (edgeFaceI != currHit.index())
97  {
98  forAll(hits, hI)
99  {
100  const pointIndexHit& hit = hits[hI];
101 
102  if (hit.index() == edgeFaceI)
103  {
104  // Check normals
105  const vector currHitNormal =
106  surface().faceNormals()[currHit.index()];
107 
108  const vector existingHitNormal =
109  surface().faceNormals()[edgeFaceI];
110 
111  const label signCurrHit =
112  pos(currHitNormal & lineVec);
113 
114  const label signExistingHit =
115  pos(existingHitNormal & lineVec);
116 
117  if (signCurrHit == signExistingHit)
118  {
119  return false;
120  }
121  }
122  }
123  }
124  }
125  }
126 
127  return true;
128 }
129 
130 
131 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
132 
134 :
135  surface_(surface),
136  tolerance_(indexedOctree<treeDataTriSurface>::perturbTol()),
137  maxTreeDepth_(10),
138  treePtr_(NULL)
139 {}
140 
141 
143 (
144  const triSurface& surface,
145  const dictionary& dict
146 )
147 :
148  surface_(surface),
150  maxTreeDepth_(10),
151  treePtr_(NULL)
152 {
153  // Have optional non-standard search tolerance for gappy surfaces.
154  if (dict.readIfPresent("tolerance", tolerance_) && tolerance_ > 0)
155  {
156  Info<< " using intersection tolerance " << tolerance_ << endl;
157  }
158 
159  // Have optional non-standard tree-depth to limit storage.
160  if (dict.readIfPresent("maxTreeDepth", maxTreeDepth_) && maxTreeDepth_ > 0)
161  {
162  Info<< " using maximum tree depth " << maxTreeDepth_ << endl;
163  }
164 }
165 
166 
168 (
169  const triSurface& surface,
170  const scalar tolerance,
171  const label maxTreeDepth
172 )
173 :
174  surface_(surface),
175  tolerance_(tolerance),
176  maxTreeDepth_(maxTreeDepth),
177  treePtr_(NULL)
178 {}
179 
180 
181 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
182 
184 {
185  clearOut();
186 }
187 
188 
190 {
191  treePtr_.clear();
192 }
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
199 {
200  if (treePtr_.empty())
201  {
202  // Calculate bb without constructing local point numbering.
204 
205  if (surface().size())
206  {
207  label nPoints;
209 
210  if (nPoints != surface().points().size())
211  {
213  << "Surface does not have compact point numbering."
214  << " Of " << surface().points().size()
215  << " only " << nPoints
216  << " are used. This might give problems in some routines."
217  << endl;
218  }
219 
220  // Random number generator. Bit dodgy since not exactly random ;-)
221  Random rndGen(65431);
222 
223  // Slightly extended bb. Slightly off-centred just so on symmetric
224  // geometry there are less face/edge aligned items.
225  bb = bb.extend(rndGen, 1e-4);
226  bb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
227  bb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
228  }
229 
232 
233  treePtr_.reset
234  (
236  (
237  treeDataTriSurface(false, surface_, tolerance_),
238  bb,
239  maxTreeDepth_, // maxLevel
240  10, // leafsize
241  3.0 // duplicity
242  )
243  );
244 
246  }
247 
248  return treePtr_();
249 }
250 
251 
252 // Determine inside/outside for samples
254 (
255  const pointField& samples
256 ) const
257 {
258  boolList inside(samples.size());
259 
260  forAll(samples, sampleI)
261  {
262  const point& sample = samples[sampleI];
263 
264  if (!tree().bb().contains(sample))
265  {
266  inside[sampleI] = false;
267  }
268  else if (tree().getVolumeType(sample) == volumeType::INSIDE)
269  {
270  inside[sampleI] = true;
271  }
272  else
273  {
274  inside[sampleI] = false;
275  }
276  }
277  return inside;
278 }
279 
280 
282 (
283  const pointField& samples,
284  const scalarField& nearestDistSqr,
285  List<pointIndexHit>& info
286 ) const
287 {
290 
291  const indexedOctree<treeDataTriSurface>& octree = tree();
292 
293  const treeDataTriSurface::findNearestOp fOp(octree);
294 
295  info.setSize(samples.size());
296 
297  forAll(samples, i)
298  {
299  info[i] = octree.findNearest
300  (
301  samples[i],
302  nearestDistSqr[i],
303  fOp
304  );
305  }
306 
308 }
309 
310 
312 (
313  const point& pt,
314  const vector& span
315 )
316 const
317 {
318  const scalar nearestDistSqr = 0.25*magSqr(span);
319 
320  return tree().findNearest(pt, nearestDistSqr);
321 }
322 
323 
325 (
326  const pointField& start,
327  const pointField& end,
328  List<pointIndexHit>& info
329 ) const
330 {
331  const indexedOctree<treeDataTriSurface>& octree = tree();
332 
333  info.setSize(start.size());
334 
337 
338  forAll(start, i)
339  {
340  info[i] = octree.findLine(start[i], end[i]);
341  }
342 
344 }
345 
346 
348 (
349  const pointField& start,
350  const pointField& end,
351  List<pointIndexHit>& info
352 ) const
353 {
354  const indexedOctree<treeDataTriSurface>& octree = tree();
355 
356  info.setSize(start.size());
357 
360 
361  forAll(start, i)
362  {
363  info[i] = octree.findLineAny(start[i], end[i]);
364  }
365 
367 }
368 
369 
371 (
372  const pointField& start,
373  const pointField& end,
374  List<List<pointIndexHit> >& info
375 ) const
376 {
377  const indexedOctree<treeDataTriSurface>& octree = tree();
378 
379  info.setSize(start.size());
380 
383 
384  // Work array
386 
387  DynamicList<label> shapeMask;
388 
389  treeDataTriSurface::findAllIntersectOp allIntersectOp(octree, shapeMask);
390 
391  forAll(start, pointI)
392  {
393  hits.clear();
394  shapeMask.clear();
395 
396  while (true)
397  {
398  // See if any intersection between pt and end
399  pointIndexHit inter = octree.findLine
400  (
401  start[pointI],
402  end[pointI],
403  allIntersectOp
404  );
405 
406  if (inter.hit())
407  {
408  vector lineVec = end[pointI] - start[pointI];
409  lineVec /= mag(lineVec) + VSMALL;
410 
411  if
412  (
413  checkUniqueHit
414  (
415  inter,
416  hits,
417  lineVec
418  )
419  )
420  {
421  hits.append(inter);
422  }
423 
424  shapeMask.append(inter.index());
425  }
426  else
427  {
428  break;
429  }
430  }
431 
432  info[pointI].transfer(hits);
433  }
434 
436 }
437 
438 
439 // ************************************************************************* //
Foam::Random
Simple random number generator.
Definition: Random.H:49
Foam::volumeType::INSIDE
@ INSIDE
Definition: volumeType.H:63
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:111
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Foam::treeDataTriSurface
Encapsulates data for (indexedOc)tree searches on triSurface.
Definition: treeDataTriSurface.H:53
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::triSurfaceSearch::nearest
pointIndexHit nearest(const point &, const vector &span) const
Calculate nearest point on surface for single searchPoint. Returns.
Definition: triSurfaceSearch.C:312
Foam::PatchTools::calcBounds
static void calcBounds(const PrimitivePatch< Face, FaceList, PointField, PointType > &p, boundBox &bb, label &nPoints)
Definition: PatchToolsSearch.C:224
Foam::treeDataPrimitivePatch::findNearestOp
Definition: treeDataPrimitivePatch.H:92
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:105
Foam::List::transfer
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Foam::triSurfaceSearch::triSurfaceSearch
triSurfaceSearch(const triSurfaceSearch &)
Disallow default bitwise copy construct.
Foam::triSurfaceSearch::findNearest
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, List< pointIndexHit > &info) const
Definition: triSurfaceSearch.C:282
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::triSurfaceSearch::findLineAll
void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &info) const
Calculate all intersections from start to end.
Definition: triSurfaceSearch.C:371
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::indexedOctree::perturbTol
static scalar & perturbTol()
Get the perturbation tolerance.
Definition: indexedOctree.C:2550
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
samples
scalarField samples(nIntervals, 0)
volumeType.H
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
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::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
Foam::treeBoundBox::extend
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
Definition: treeBoundBoxI.H:317
Foam::indexedOctree::findLine
pointIndexHit findLine(const bool findAny, const point &treeStart, const point &treeEnd, const label startNodeI, const direction startOctantI, const FindIntersectOp &fiOp, const bool verbose=false) const
Find any or nearest intersection.
Foam::Info
messageStream Info
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:117
Foam::triSurfaceSearch::tree
const indexedOctree< treeDataTriSurface > & tree() const
Demand driven construction of the octree.
Definition: triSurfaceSearch.C:198
Foam::triSurfaceSearch::~triSurfaceSearch
~triSurfaceSearch()
Destructor.
Definition: triSurfaceSearch.C:183
Foam::DynamicList::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:242
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataTriSurface.H:47
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::triSurfaceSearch::clearOut
void clearOut()
Clear storage.
Definition: triSurfaceSearch.C:189
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::triSurfaceSearch::calcInside
boolList calcInside(const pointField &searchPoints) const
Calculate for each searchPoint inside/outside status.
Definition: triSurfaceSearch.C:254
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::treeDataTriSurface
treeDataPrimitivePatch< triSurface > treeDataTriSurface
Definition: treeDataTriSurface.H:47
Foam::triSurfaceSearch::findLine
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
Definition: triSurfaceSearch.C:325
f
labelList f(nPoints)
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::surface
Definition: surface.H:55
Foam::indexedOctree::findLineAny
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
triSurfaceSearch.H
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::triSurfaceSearch::findLineAny
void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
Definition: triSurfaceSearch.C:348
Foam::point
vector point
Point is a vector.
Definition: point.H:41
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::treeDataPrimitivePatch::findAllIntersectOp
Definition: treeDataPrimitivePatch.H:145
Foam::triSurfaceSearch::checkUniqueHit
bool checkUniqueHit(const pointIndexHit &currHit, const DynamicList< pointIndexHit, 1, 1 > &hits, const vector &lineVec) const
Check whether the current hit on the surface which lies on lineVec.
Definition: triSurfaceSearch.C:34
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
rndGen
cachedRandom rndGen(label(0), -1)
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190