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-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 "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 vector& span,
59  const coordSystem::cartesian& csys
60 )
61 :
63  box_
64  (
65  IOobject
66  (
67  io.name() + "_box",
68  io.instance(),
69  io.local(),
70  io.db(),
71  io.readOpt(),
72  io.writeOpt(),
73  false // never register
74  ),
75  treeBoundBox(Zero, span)
76  ),
77  transform_(csys.origin(), csys.e3(), csys.e1())
78 {
79  points_ = transform_.globalPosition(box_.points());
80 }
81 
82 
83 Foam::searchableRotatedBox::searchableRotatedBox
84 (
85  const IOobject& io,
86  const dictionary& dict
87 )
88 :
90  (
91  io,
92  dict.get<vector>("span"),
93  coordSystem::cartesian
94  (
95  dict.get<point>("origin"),
96  dict.get<vector>("e3"),
97  dict.get<vector>("e1")
98  )
99  )
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 {
107  return box_.regions();
108 }
109 
110 
112 {
113  return transform_.globalPosition(box_.coordinates());
114 }
115 
116 
118 (
119  pointField& centres,
120  scalarField& radiusSqr
121 ) const
122 {
123  box_.boundingSpheres(centres, radiusSqr);
124  centres = transform_.globalPosition(centres);
125 }
126 
129 {
130  return points_;
131 }
132 
133 
135 {
136  // (from treeDataPrimitivePatch.C)
137 
138  // 1. bounding box
139  if (!bb.overlaps(bounds()))
140  {
141  return false;
142  }
143 
144  // 2. Check if any corner points inside
145  if (bb.containsAny(points_))
146  {
147  return true;
148  }
149 
150  // 3. Difficult case: all points are outside but connecting edges might
151  // go through cube.
152 
153  const treeBoundBox treeBb(bb);
154 
155  // 3a. my edges through bb faces
156  for (const edge& e : treeBoundBox::edges)
157  {
158  point inter;
159  if (treeBb.intersects(points_[e[0]], points_[e[1]], inter))
160  {
161  return true;
162  }
163  }
164 
165  // 3b. bb edges through my faces
166 
167  const pointField bbPoints(bb.points());
168 
169  for (const face& f : treeBoundBox::faces)
170  {
171  const point fc = f.centre(points_);
172 
173  for (const edge& e : treeBoundBox::edges)
174  {
175  pointHit inter = f.intersection
176  (
177  bbPoints[e[0]],
178  bbPoints[e[1]],
179  fc,
180  points_,
182  );
183 
184  if (inter.hit() && inter.distance() <= 1)
185  {
186  return true;
187  }
188  }
189  }
190 
191  return false;
192 }
193 
194 
196 (
197  const point& sample,
198  const scalar nearestDistSqr
199 ) const
200 {
201  pointIndexHit boxNearest
202  (
203  box_.findNearest
204  (
205  transform_.localPosition(sample),
206  nearestDistSqr
207  )
208  );
209 
210  boxNearest.point() = transform_.globalPosition(boxNearest.point());
211 
212  return boxNearest;
213 }
214 
215 
217 (
218  const linePointRef& ln,
219  treeBoundBox& tightest,
220  point& linePoint
221 ) const
222 {
224  return pointIndexHit();
225 }
226 
227 
229 (
230  const point& start,
231  const point& end
232 ) const
233 {
234  pointIndexHit boxHit
235  (
236  box_.findLine
237  (
238  transform_.localPosition(start),
239  transform_.localPosition(end)
240  )
241  );
242 
243  boxHit.point() = transform_.globalPosition(boxHit.point());
244 
245  return boxHit;
246 }
247 
248 
250 (
251  const point& start,
252  const point& end
253 ) const
254 {
255  return findLine(start, end);
256 }
257 
258 
260 (
261  const pointField& samples,
262  const scalarField& nearestDistSqr,
263  List<pointIndexHit>& info
264 ) const
265 {
266  info.setSize(samples.size());
267 
268  forAll(samples, i)
269  {
270  info[i] = findNearest(samples[i], nearestDistSqr[i]);
271  }
272 }
273 
274 
276 (
277  const pointField& start,
278  const pointField& end,
279  List<pointIndexHit>& info
280 ) const
281 {
282  info.setSize(start.size());
283 
284  forAll(start, i)
285  {
286  info[i] = findLine(start[i], end[i]);
287  }
288 }
289 
290 
292 (
293  const pointField& start,
294  const pointField& end,
295  List<pointIndexHit>& info
296 ) const
297 {
298  info.setSize(start.size());
299 
300  forAll(start, i)
301  {
302  info[i] = findLineAny(start[i], end[i]);
303  }
304 }
305 
306 
308 (
309  const pointField& start,
310  const pointField& end,
311  List<List<pointIndexHit>>& info
312 ) const
313 {
314  info.setSize(start.size());
315 
316  // Work array
317  DynamicList<pointIndexHit> hits;
318 
319  // Tolerances:
320  // To find all intersections we add a small vector to the last intersection
321  // This is chosen such that
322  // - it is significant (SMALL is smallest representative relative tolerance;
323  // we need something bigger since we're doing calculations)
324  // - if the start-end vector is zero we still progress
325  const vectorField dirVec(end-start);
326  const scalarField magSqrDirVec(magSqr(dirVec));
327  const vectorField smallVec
328  (
329  Foam::sqrt(SMALL)*dirVec
330  + vector(ROOTVSMALL,ROOTVSMALL,ROOTVSMALL)
331  );
332 
333  forAll(start, pointI)
334  {
335  // See if any intersection between pt and end
336  pointIndexHit inter = findLine(start[pointI], end[pointI]);
337 
338  if (inter.hit())
339  {
340  hits.clear();
341  hits.append(inter);
342 
343  point pt = inter.point() + smallVec[pointI];
344 
345  while (((pt-start[pointI])&dirVec[pointI]) <= magSqrDirVec[pointI])
346  {
347  // See if any intersection between pt and end
348  pointIndexHit inter = findLine(pt, end[pointI]);
349 
350  // Check for not hit or hit same face as before (can happen
351  // if vector along surface of face)
352  if
353  (
354  !inter.hit()
355  || (inter.index() == hits.last().index())
356  )
357  {
358  break;
359  }
360  hits.append(inter);
361 
362  pt = inter.point() + smallVec[pointI];
363  }
364 
365  info[pointI].transfer(hits);
366  }
367  else
368  {
369  info[pointI].clear();
370  }
371  }
372 }
373 
374 
376 (
377  const List<pointIndexHit>& info,
378  labelList& region
379 ) const
380 {
381  region.setSize(info.size());
382  region = 0;
383 }
384 
385 
387 (
388  const List<pointIndexHit>& info,
389  vectorField& normal
390 ) const
391 {
392  // searchableBox does not use hitPoints so no need to transform
393  box_.getNormal(info, normal);
394 
395  normal = transform_.globalVector(normal);
396 }
397 
398 
400 (
401  const pointField& points,
402  List<volumeType>& volType
403 ) const
404 {
405  box_.getVolumeType(transform_.localPosition(points), volType);
406 }
407 
408 
409 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
A line primitive.
Definition: line.H:52
virtual const wordList & regions() const
Names of regions.
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
scalarField samples(nIntervals, Zero)
dimensionedScalar sqrt(const dimensionedScalar &ds)
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
static const faceList faces
Face to point addressing, using octant corner points.
Definition: treeBoundBox.H:174
Macros for easy insertion into run-time selection tables.
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
A Cartesian coordinate system.
Definition: cartesianCS.H:65
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point. unknown if.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Searching on a rotated box.
static const edgeList edges
Edge to point addressing, using octant corner points.
Definition: treeBoundBox.H:179
bool local
Definition: EEqn.H:20
Vector< scalar > vector
Definition: vector.H:57
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:1237
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
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
bool containsAny(const UList< point > &points) const
Contains any of the points? (inside or on edge)
Definition: boundBox.C:385
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.
pointIndexHit findNearest(const point &sample, const scalar nearestDistSqr) const
Calculate nearest point on surface.
vector point
Point is a vector.
Definition: point.H:37
tmp< pointField > points() const
Corner points in an order corresponding to a &#39;hex&#39; cell.
Definition: boundBox.H:374
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual tmp< pointField > points() const
Get the points that define the surface.
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:439
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
PointHit< point > pointHit
A PointHit with a 3D point.
Definition: pointHit.H:43
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127