searchableRotatedBox.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) 2014 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 "searchableRotatedBox.H"
28 #include "axesRotation.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(searchableRotatedBox, 0);
35  addToRunTimeSelectionTable(searchableSurface, searchableRotatedBox, dict);
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const IOobject& io,
44  const dictionary& dict
45 )
46 :
48  box_
49  (
50  IOobject
51  (
52  io.name() + "_box",
53  io.instance(),
54  io.local(),
55  io.db(),
56  io.readOpt(),
57  io.writeOpt(),
58  false //io.registerObject(),
59  ),
60  treeBoundBox(point::zero, dict.lookup("span"))
61  ),
62  transform_
63  (
64  "rotation",
65  dict.lookup("origin"),
66  dict.lookup("e3"),
67  dict.lookup("e1")
68  )
69 {
70  points_ = transform_.globalPosition(box_.points());
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75 
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
83 {
84  return box_.regions();
85 }
86 
87 
89 {
90  return transform_.globalPosition(box_.coordinates());
91 }
92 
93 
95 (
96  pointField& centres,
97  scalarField& radiusSqr
98 ) const
99 {
100  box_.boundingSpheres(centres, radiusSqr);
101  centres = transform_.globalPosition(centres);
102 }
103 
104 
106 {
107  return points_;
108 }
109 
110 
112 {
113  // (from treeDataPrimitivePatch.C)
114 
115  // 1. bounding box
116  if (!bb.overlaps(bounds()))
117  {
118  return false;
119  }
120 
121  // 2. Check if one or more face points inside
122  const faceList& fcs = treeBoundBox::faces;
123  forAll(fcs, faceI)
124  {
125  if (bb.containsAny(points_))
126  {
127  return true;
128  }
129  }
130 
131  // 3. Difficult case: all points are outside but connecting edges might
132  // go through cube.
133 
134  const treeBoundBox treeBb(bb);
135 
136  // 3a. my edges through bb faces
137  const edgeList& edges = treeBoundBox::edges;
138  forAll(edges, edgeI)
139  {
140  const edge& e = edges[edgeI];
141 
142  point inter;
143  if (treeBb.intersects(points_[e[0]], points_[e[1]], inter))
144  {
145  return true;
146  }
147  }
148 
149  // 3b. bb edges through my faces
150 
151  const pointField bbPoints(bb.points());
152 
153  forAll(fcs, faceI)
154  {
155  const face& f = fcs[faceI];
156  point fc = f.centre(points_);
157 
158  forAll(edges, edgeI)
159  {
160  const edge& e = edges[edgeI];
161 
162  pointHit inter = f.intersection
163  (
164  bbPoints[e[0]],
165  bbPoints[e[1]],
166  fc,
167  points_,
169  );
170 
171  if (inter.hit() && inter.distance() <= 1)
172  {
173  return true;
174  }
175  }
176  }
177 
178  return false;
179 }
180 
181 
183 (
184  const point& sample,
185  const scalar nearestDistSqr
186 ) const
187 {
188  pointIndexHit boxNearest
189  (
190  box_.findNearest
191  (
192  transform_.localPosition(sample),
193  nearestDistSqr
194  )
195  );
196 
197  boxNearest.rawPoint() = transform_.globalPosition(boxNearest.rawPoint());
198 
199  return boxNearest;
200 }
201 
202 
204 (
205  const linePointRef& ln,
206  treeBoundBox& tightest,
207  point& linePoint
208 ) const
209 {
211  return pointIndexHit();
212 }
213 
214 
216 (
217  const point& start,
218  const point& end
219 ) const
220 {
221  pointIndexHit boxHit
222  (
223  box_.findLine
224  (
225  transform_.localPosition(start),
226  transform_.localPosition(end)
227  )
228  );
229 
230  boxHit.rawPoint() = transform_.globalPosition(boxHit.rawPoint());
231 
232  return boxHit;
233 }
234 
235 
237 (
238  const point& start,
239  const point& end
240 ) const
241 {
242  return findLine(start, end);
243 }
244 
245 
247 (
248  const pointField& samples,
249  const scalarField& nearestDistSqr,
250  List<pointIndexHit>& info
251 ) const
252 {
253  info.setSize(samples.size());
254 
255  forAll(samples, i)
256  {
257  info[i] = findNearest(samples[i], nearestDistSqr[i]);
258  }
259 }
260 
261 
263 (
264  const pointField& start,
265  const pointField& end,
266  List<pointIndexHit>& info
267 ) const
268 {
269  info.setSize(start.size());
270 
271  forAll(start, i)
272  {
273  info[i] = findLine(start[i], end[i]);
274  }
275 }
276 
277 
279 (
280  const pointField& start,
281  const pointField& end,
282  List<pointIndexHit>& info
283 ) const
284 {
285  info.setSize(start.size());
286 
287  forAll(start, i)
288  {
289  info[i] = findLineAny(start[i], end[i]);
290  }
291 }
292 
293 
295 (
296  const pointField& start,
297  const pointField& end,
298  List<List<pointIndexHit> >& info
299 ) const
300 {
301  info.setSize(start.size());
302 
303  // Work array
305 
306  // Tolerances:
307  // To find all intersections we add a small vector to the last intersection
308  // This is chosen such that
309  // - it is significant (SMALL is smallest representative relative tolerance;
310  // we need something bigger since we're doing calculations)
311  // - if the start-end vector is zero we still progress
312  const vectorField dirVec(end-start);
313  const scalarField magSqrDirVec(magSqr(dirVec));
314  const vectorField smallVec
315  (
316  Foam::sqrt(SMALL)*dirVec
317  + vector(ROOTVSMALL,ROOTVSMALL,ROOTVSMALL)
318  );
319 
320  forAll(start, pointI)
321  {
322  // See if any intersection between pt and end
323  pointIndexHit inter = findLine(start[pointI], end[pointI]);
324 
325  if (inter.hit())
326  {
327  hits.clear();
328  hits.append(inter);
329 
330  point pt = inter.hitPoint() + smallVec[pointI];
331 
332  while (((pt-start[pointI])&dirVec[pointI]) <= magSqrDirVec[pointI])
333  {
334  // See if any intersection between pt and end
335  pointIndexHit inter = findLine(pt, end[pointI]);
336 
337  // Check for not hit or hit same face as before (can happen
338  // if vector along surface of face)
339  if
340  (
341  !inter.hit()
342  || (inter.index() == hits.last().index())
343  )
344  {
345  break;
346  }
347  hits.append(inter);
348 
349  pt = inter.hitPoint() + smallVec[pointI];
350  }
351 
352  info[pointI].transfer(hits);
353  }
354  else
355  {
356  info[pointI].clear();
357  }
358  }
359 }
360 
361 
363 (
364  const List<pointIndexHit>& info,
365  labelList& region
366 ) const
367 {
368  region.setSize(info.size());
369  region = 0;
370 }
371 
372 
374 (
375  const List<pointIndexHit>& info,
377 ) const
378 {
379  // searchableBox does not use hitPoints so no need to transform
380  box_.getNormal(info, normal);
381 
382  normal = transform_.globalVector(normal);
383 }
384 
385 
387 (
388  const pointField& points,
389  List<volumeType>& volType
390 ) const
391 {
392  box_.getVolumeType(transform_.localPosition(points), volType);
393 }
394 
395 
396 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:111
Foam::searchableRotatedBox::regions
virtual const wordList & regions() const
Names of regions.
Definition: searchableRotatedBox.C:82
Foam::PointHit::hit
bool hit() const
Is there a hit.
Definition: PointHit.H:120
Foam::intersection::HALF_RAY
@ HALF_RAY
Definition: intersection.H:72
Foam::searchableRotatedBox::points
virtual tmp< pointField > points() const
Get the points that define the surface.
Definition: searchableRotatedBox.C:105
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::PointHit
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
Definition: PointHit.H:51
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::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::searchableRotatedBox::findLineAny
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.
Definition: searchableRotatedBox.C:237
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
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::searchableRotatedBox::getRegion
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
Definition: searchableRotatedBox.C:363
Foam::searchableRotatedBox::getNormal
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
Definition: searchableRotatedBox.C:374
Foam::treeBoundBox::edges
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:157
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobject.H:350
Foam::IOobject::writeOpt
writeOption writeOpt() const
Definition: IOobject.H:327
Foam::searchableRotatedBox::coordinates
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
Definition: searchableRotatedBox.C:88
Foam::IOobject::db
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:239
samples
scalarField samples(nIntervals, 0)
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:365
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::PointIndexHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointIndexHit.H:143
Foam::PointHit::distance
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
Foam::searchableRotatedBox::overlaps
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
Definition: searchableRotatedBox.C:111
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::IOobject::readOpt
readOption readOpt() const
Definition: IOobject.H:317
Foam::treeBoundBox::intersects
bool intersects(const point &overallStart, const vector &overallVec, const point &start, const point &end, point &pt, direction &ptBits) const
Intersects segment; set point to intersection position and face,.
Definition: treeBoundBox.C:253
Foam::boundBox::overlaps
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:117
Foam::treeBoundBox::faces
static const faceList faces
Face to point addressing.
Definition: treeBoundBox.H:154
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:66
Foam::DynamicList::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:242
Foam::searchableRotatedBox::findNearest
pointIndexHit findNearest(const point &sample, const scalar nearestDistSqr) const
Calculate nearest point on surface.
Definition: searchableRotatedBox.C:183
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::boundBox::containsAny
bool containsAny(const UList< point > &) const
Contains any of the points? (inside or on edge)
Definition: boundBox.C:261
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
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::searchableRotatedBox::~searchableRotatedBox
virtual ~searchableRotatedBox()
Destructor.
Definition: searchableRotatedBox.C:76
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::searchableRotatedBox::findLine
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end.
Definition: searchableRotatedBox.C:216
Foam::searchableRotatedBox::boundingSpheres
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
Definition: searchableRotatedBox.C:95
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
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
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:379
Foam::line
A line primitive.
Definition: line.H:56
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::searchableRotatedBox::searchableRotatedBox
searchableRotatedBox(const searchableRotatedBox &)
Disallow default bitwise copy construct.
Foam::ln
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:854
Foam::searchableRotatedBox::getVolumeType
virtual void getVolumeType(const pointField &, List< volumeType > &) const
Determine type (inside/outside/mixed) for point. unknown if.
Definition: searchableRotatedBox.C:387
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
searchableRotatedBox.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::boundBox::points
tmp< pointField > points() const
Return corner points in an order corresponding to a 'hex' cell.
Definition: boundBox.C:149
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
axesRotation.H
Foam::IOobject::local
const fileName & local() const
Definition: IOobject.H:360
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::searchableRotatedBox::findLineAll
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &) const
Get all intersections in order from start to end.
Definition: searchableRotatedBox.C:295
normal
A normal distribution model.