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