searchableCylinder.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018-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 "searchableCylinder.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(searchableCylinder, 0);
38  (
39  searchableSurface,
40  searchableCylinder,
41  dict
42  );
44  (
45  searchableSurface,
46  searchableCylinder,
47  dict,
48  cylinder
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
56 {
57  return tmp<pointField>::New(1, 0.5*(point1_ + point2_));
58 }
59 
60 
62 (
63  pointField& centres,
64  scalarField& radiusSqr
65 ) const
66 {
67  centres.setSize(1);
68  centres[0] = 0.5*(point1_ + point2_);
69 
70  radiusSqr.setSize(1);
71  radiusSqr[0] = Foam::magSqr(point1_-centres[0]) + Foam::sqr(radius_);
72 
73  // Add a bit to make sure all points are tested inside
74  radiusSqr += Foam::sqr(SMALL);
75 }
76 
77 
79 {
80  auto tpts = tmp<pointField>::New(2);
81  auto& pts = tpts.ref();
82 
83  pts[0] = point1_;
84  pts[1] = point2_;
85 
86  return tpts;
87 }
88 
89 
90 Foam::pointIndexHit Foam::searchableCylinder::findNearest
91 (
92  const point& sample,
93  const scalar nearestDistSqr
94 ) const
95 {
96  pointIndexHit info(false, sample, -1);
97 
98  vector v(sample - point1_);
99 
100  // Decompose sample-point1 into normal and parallel component
101  scalar parallel = (v & unitDir_);
102 
103  // Remove the parallel component and normalise
104  v -= parallel*unitDir_;
105  scalar magV = mag(v);
106 
107  if (magV < ROOTVSMALL)
108  {
109  v = Zero;
110  }
111  else
112  {
113  v /= magV;
114  }
115 
116  if (parallel <= 0)
117  {
118  // nearest is at point1 end cap. Clip to radius.
119  info.setPoint(point1_ + min(magV, radius_)*v);
120  }
121  else if (parallel >= magDir_)
122  {
123  // nearest is at point2 end cap. Clip to radius.
124  info.setPoint(point2_ + min(magV, radius_)*v);
125  }
126  else
127  {
128  // inbetween endcaps. Might either be nearer endcaps or cylinder wall
129 
130  // distance to endpoint: parallel or parallel-magDir
131  // distance to cylinder wall: magV-radius_
132 
133  // Nearest cylinder point
134  point cylPt;
135  if (magV < ROOTVSMALL)
136  {
137  // Point exactly on centre line. Take any point on wall.
138  vector e1 = point(1,0,0) ^ unitDir_;
139  scalar magE1 = mag(e1);
140  if (magE1 < SMALL)
141  {
142  e1 = point(0,1,0) ^ unitDir_;
143  magE1 = mag(e1);
144  }
145  e1 /= magE1;
146  cylPt = sample + radius_*e1;
147  }
148  else
149  {
150  cylPt = sample + (radius_-magV)*v;
151  }
152 
153  if (parallel < 0.5*magDir_)
154  {
155  // Project onto p1 endcap
156  point end1Pt = point1_ + min(magV, radius_)*v;
157 
158  if (magSqr(sample-cylPt) < magSqr(sample-end1Pt))
159  {
160  info.setPoint(cylPt);
161  }
162  else
163  {
164  info.setPoint(end1Pt);
165  }
166  }
167  else
168  {
169  // Project onto p2 endcap
170  point end2Pt = point2_ + min(magV, radius_)*v;
171 
172  if (magSqr(sample-cylPt) < magSqr(sample-end2Pt))
173  {
174  info.setPoint(cylPt);
175  }
176  else
177  {
178  info.setPoint(end2Pt);
179  }
180  }
181  }
182 
183  if (info.point().distSqr(sample) < nearestDistSqr)
184  {
185  info.setHit();
186  info.setIndex(0);
187  }
188 
189  return info;
190 }
191 
192 
193 Foam::scalar Foam::searchableCylinder::radius2(const point& pt) const
194 {
195  const vector x = (pt-point1_) ^ unitDir_;
196  return (x & x);
197 }
198 
199 
200 // From http://www.gamedev.net/community/forums/topic.asp?topic_id=467789 -
201 // intersection of cylinder with ray
202 void Foam::searchableCylinder::findLineAll
203 (
204  const point& start,
205  const point& end,
206  pointIndexHit& near,
207  pointIndexHit& far
208 ) const
209 {
210  near.setMiss();
211  far.setMiss();
212 
213  vector point1Start(start-point1_);
214  vector point2Start(start-point2_);
215  vector point1End(end-point1_);
216 
217  // Quick rejection of complete vector outside endcaps
218  scalar s1 = point1Start&unitDir_;
219  scalar s2 = point1End&unitDir_;
220 
221  if ((s1 < 0 && s2 < 0) || (s1 > magDir_ && s2 > magDir_))
222  {
223  return;
224  }
225 
226  // Line as P = start+t*V where V is unit vector and t=[0..mag(end-start)]
227  vector V(end-start);
228  scalar magV = mag(V);
229  if (magV < ROOTVSMALL)
230  {
231  return;
232  }
233  V /= magV;
234 
235 
236  // We now get the nearest intersections to start. This can either be
237  // the intersection with the end plane or with the cylinder side.
238 
239  // Get the two points (expressed in t) on the end planes. This is to
240  // clip any cylinder intersection against.
241  scalar tPoint1;
242  scalar tPoint2;
243 
244  // Maintain the two intersections with the endcaps
245  scalar tNear = VGREAT;
246  scalar tFar = VGREAT;
247 
248  {
249  scalar s = (V&unitDir_);
250  if (mag(s) > VSMALL)
251  {
252  tPoint1 = -s1/s;
253  tPoint2 = -(point2Start&unitDir_)/s;
254  if (tPoint2 < tPoint1)
255  {
256  std::swap(tPoint1, tPoint2);
257  }
258  if (tPoint1 > magV || tPoint2 < 0)
259  {
260  return;
261  }
262 
263  // See if the points on the endcaps are actually inside the cylinder
264  if (tPoint1 >= 0 && tPoint1 <= magV)
265  {
266  if (radius2(start+tPoint1*V) <= sqr(radius_))
267  {
268  tNear = tPoint1;
269  }
270  }
271  if (tPoint2 >= 0 && tPoint2 <= magV)
272  {
273  if (radius2(start+tPoint2*V) <= sqr(radius_))
274  {
275  // Check if already have a near hit from point1
276  if (tNear <= 1)
277  {
278  tFar = tPoint2;
279  }
280  else
281  {
282  tNear = tPoint2;
283  }
284  }
285  }
286  }
287  else
288  {
289  // Vector perpendicular to cylinder. Check for outside already done
290  // above so just set tpoint to allow all.
291  tPoint1 = -VGREAT;
292  tPoint2 = VGREAT;
293  }
294  }
295 
296 
297  const vector x = point1Start ^ unitDir_;
298  const vector y = V ^ unitDir_;
299  const scalar d = sqr(radius_);
300 
301  // Second order equation of the form a*t^2 + b*t + c
302  const scalar a = (y&y);
303  const scalar b = 2*(x&y);
304  const scalar c = (x&x)-d;
305 
306  const scalar disc = b*b-4*a*c;
307 
308  scalar t1 = -VGREAT;
309  scalar t2 = VGREAT;
310 
311  if (disc < 0)
312  {
313  // Fully outside
314  return;
315  }
316  else if (disc < ROOTVSMALL)
317  {
318  // Single solution
319  if (mag(a) > ROOTVSMALL)
320  {
321  t1 = -b/(2*a);
322 
323  //Pout<< "single solution t:" << t1
324  // << " for start:" << start << " end:" << end
325  // << " c:" << c << endl;
326 
327  if (t1 >= 0 && t1 <= magV && t1 >= tPoint1 && t1 <= tPoint2)
328  {
329  // valid. Insert sorted.
330  if (t1 < tNear)
331  {
332  tFar = tNear;
333  tNear = t1;
334  }
335  else if (t1 < tFar)
336  {
337  tFar = t1;
338  }
339  }
340  else
341  {
342  return;
343  }
344  }
345  else
346  {
347  // Aligned with axis. Check if outside radius
348  //Pout<< "small discriminant:" << disc
349  // << " for start:" << start << " end:" << end
350  // << " magV:" << magV
351  // << " c:" << c << endl;
352  if (c > 0)
353  {
354  return;
355  }
356  }
357  }
358  else
359  {
360  if (mag(a) > ROOTVSMALL)
361  {
362  scalar sqrtDisc = sqrt(disc);
363 
364  t1 = (-b - sqrtDisc)/(2*a);
365  t2 = (-b + sqrtDisc)/(2*a);
366  if (t2 < t1)
367  {
368  std::swap(t1, t2);
369  }
370 
371  if (t1 >= 0 && t1 <= magV && t1 >= tPoint1 && t1 <= tPoint2)
372  {
373  // valid. Insert sorted.
374  if (t1 < tNear)
375  {
376  tFar = tNear;
377  tNear = t1;
378  }
379  else if (t1 < tFar)
380  {
381  tFar = t1;
382  }
383  }
384  if (t2 >= 0 && t2 <= magV && t2 >= tPoint1 && t2 <= tPoint2)
385  {
386  // valid. Insert sorted.
387  if (t2 < tNear)
388  {
389  tFar = tNear;
390  tNear = t2;
391  }
392  else if (t2 < tFar)
393  {
394  tFar = t2;
395  }
396  }
397  //Pout<< "two solutions t1:" << t1 << " t2:" << t2
398  // << " for start:" << start << " end:" << end
399  // << " magV:" << magV
400  // << " c:" << c << endl;
401  }
402  else
403  {
404  // Aligned with axis. Check if outside radius
405  //Pout<< "large discriminant:" << disc
406  // << " small a:" << a
407  // << " for start:" << start << " end:" << end
408  // << " magV:" << magV
409  // << " c:" << c << endl;
410  if (c > 0)
411  {
412  return;
413  }
414  }
415  }
416 
417  // Check tNear, tFar
418  if (tNear >= 0 && tNear <= magV)
419  {
420  near.hitPoint(start+tNear*V);
421  near.setIndex(0);
422 
423  if (tFar <= magV)
424  {
425  far.hitPoint(start+tFar*V);
426  far.setIndex(0);
427  }
428  }
429  else if (tFar >= 0 && tFar <= magV)
430  {
431  near.hitPoint(start+tFar*V);
432  near.setIndex(0);
433  }
434 }
435 
436 
437 Foam::boundBox Foam::searchableCylinder::calcBounds() const
438 {
439 
440  // Adapted from
441  // http://www.gamedev.net/community/forums
442  // /topic.asp?topic_id=338522&forum_id=20&gforum_id=0
443 
444  // Let cylinder have end points A,B and radius r,
445 
446  // Bounds in direction X (same for Y and Z) can be found as:
447  // Let A.X<B.X (otherwise swap points)
448  // Good approximate lowest bound is A.X-r and highest is B.X+r (precise for
449  // capsule). At worst, in one direction it can be larger than needed by 2*r.
450 
451  // Accurate bounds for cylinder is
452  // A.X-kx*r, B.X+kx*r
453  // where
454  // kx=sqrt(((A.Y-B.Y)^2+(A.Z-B.Z)^2)/((A.X-B.X)^2+(A.Y-B.Y)^2+(A.Z-B.Z)^2))
455 
456  // similar thing for Y and Z
457  // (i.e.
458  // ky=sqrt(((A.X-B.X)^2+(A.Z-B.Z)^2)/((A.X-B.X)^2+(A.Y-B.Y)^2+(A.Z-B.Z)^2))
459  // kz=sqrt(((A.X-B.X)^2+(A.Y-B.Y)^2)/((A.X-B.X)^2+(A.Y-B.Y)^2+(A.Z-B.Z)^2))
460  // )
461 
462  // How derived: geometric reasoning. Bounds of cylinder is same as for 2
463  // circles centered on A and B. This sqrt thingy gives sine of angle between
464  // axis and direction, used to find projection of radius.
465 
466  vector kr
467  (
468  sqrt(sqr(unitDir_.y()) + sqr(unitDir_.z())),
469  sqrt(sqr(unitDir_.x()) + sqr(unitDir_.z())),
470  sqrt(sqr(unitDir_.x()) + sqr(unitDir_.y()))
471  );
472 
473  kr *= radius_;
474 
475  point min = point1_ - kr;
476  point max = point1_ + kr;
477 
478  min = ::Foam::min(min, point2_ - kr);
479  max = ::Foam::max(max, point2_ + kr);
480 
481  return boundBox(min, max);
482 }
483 
484 
485 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
486 
487 Foam::searchableCylinder::searchableCylinder
488 (
489  const IOobject& io,
490  const point& point1,
491  const point& point2,
492  const scalar radius
493 )
494 :
496  point1_(point1),
497  point2_(point2),
498  magDir_(mag(point2_-point1_)),
499  unitDir_((point2_-point1_)/magDir_),
500  radius_(radius)
501 {
502  bounds() = calcBounds();
503 }
504 
505 
506 Foam::searchableCylinder::searchableCylinder
507 (
508  const IOobject& io,
509  const dictionary& dict
510 )
511 :
513  point1_(dict.get<point>("point1")),
514  point2_(dict.get<point>("point2")),
515  magDir_(mag(point2_-point1_)),
516  unitDir_((point2_-point1_)/magDir_),
517  radius_(dict.get<scalar>("radius"))
518 {
519  bounds() = calcBounds();
520 }
521 
522 
523 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
524 
526 {
527  if (regions_.empty())
528  {
529  regions_.setSize(1);
530  regions_[0] = "region0";
531  }
532  return regions_;
533 }
534 
535 
536 void Foam::searchableCylinder::findNearest
537 (
538  const pointField& samples,
539  const scalarField& nearestDistSqr,
540  List<pointIndexHit>& info
541 ) const
542 {
543  info.setSize(samples.size());
544 
545  forAll(samples, i)
546  {
547  info[i] = findNearest(samples[i], nearestDistSqr[i]);
548  }
549 }
550 
551 
553 (
554  const pointField& start,
555  const pointField& end,
556  List<pointIndexHit>& info
557 ) const
558 {
559  info.setSize(start.size());
560 
561  forAll(start, i)
562  {
563  // Pick nearest intersection. If none intersected take second one.
565  findLineAll(start[i], end[i], info[i], b);
566  if (!info[i].hit() && b.hit())
567  {
568  info[i] = b;
569  }
570  }
571 }
572 
573 
575 (
576  const pointField& start,
577  const pointField& end,
578  List<pointIndexHit>& info
579 ) const
580 {
581  info.setSize(start.size());
582 
583  forAll(start, i)
584  {
585  // Discard far intersection
587  findLineAll(start[i], end[i], info[i], b);
588  if (!info[i].hit() && b.hit())
589  {
590  info[i] = b;
591  }
592  }
593 }
594 
595 
596 void Foam::searchableCylinder::findLineAll
597 (
598  const pointField& start,
599  const pointField& end,
600  List<List<pointIndexHit>>& info
601 ) const
602 {
603  info.setSize(start.size());
604 
605  forAll(start, i)
606  {
607  pointIndexHit near, far;
608  findLineAll(start[i], end[i], near, far);
609 
610  if (near.hit())
611  {
612  if (far.hit())
613  {
614  info[i].setSize(2);
615  info[i][0] = near;
616  info[i][1] = far;
617  }
618  else
619  {
620  info[i].setSize(1);
621  info[i][0] = near;
622  }
623  }
624  else
625  {
626  if (far.hit())
627  {
628  info[i].setSize(1);
629  info[i][0] = far;
630  }
631  else
632  {
633  info[i].clear();
634  }
635  }
636  }
637 }
638 
639 
641 (
642  const List<pointIndexHit>& info,
643  labelList& region
644 ) const
645 {
646  region.setSize(info.size());
647  region = 0;
648 }
649 
650 
652 (
653  const List<pointIndexHit>& info,
654  vectorField& normal
655 ) const
656 {
657  normal.setSize(info.size());
658  normal = Zero;
659 
660  forAll(info, i)
661  {
662  if (info[i].hit())
663  {
664  vector v(info[i].point() - point1_);
665 
666  // Decompose sample-point1 into normal and parallel component
667  const scalar parallel = (v & unitDir_);
668 
669  // Remove the parallel component and normalise
670  v -= parallel*unitDir_;
671  scalar magV = mag(v);
672 
673  if (parallel <= 0)
674  {
675  if ((magV-radius_) < mag(parallel))
676  {
677  // either above endcap (magV<radius) or outside but closer
678  normal[i] = -unitDir_;
679  }
680  else
681  {
682  normal[i] = v/magV;
683  }
684  }
685  else if (parallel <= 0.5*magDir_)
686  {
687  // See if endcap closer or sidewall
688  if (magV >= radius_ || (radius_-magV) < parallel)
689  {
690  normal[i] = v/magV;
691  }
692  else
693  {
694  // closer to endcap
695  normal[i] = -unitDir_;
696  }
697  }
698  else if (parallel <= magDir_)
699  {
700  // See if endcap closer or sidewall
701  if (magV >= radius_ || (radius_-magV) < (magDir_-parallel))
702  {
703  normal[i] = v/magV;
704  }
705  else
706  {
707  // closer to endcap
708  normal[i] = unitDir_;
709  }
710  }
711  else // beyond cylinder
712  {
713  if ((magV-radius_) < (parallel-magDir_))
714  {
715  // above endcap
716  normal[i] = unitDir_;
717  }
718  else
719  {
720  normal[i] = v/magV;
721  }
722  }
723  }
724  }
725 }
726 
727 
729 (
730  const pointField& points,
731  List<volumeType>& volType
732 ) const
733 {
734  volType.setSize(points.size());
735 
736  forAll(points, pointi)
737  {
738  const point& pt = points[pointi];
739 
740  volType[pointi] = volumeType::OUTSIDE;
741 
742  vector v(pt - point1_);
743 
744  // Decompose sample-point1 into normal and parallel component
745  const scalar parallel = (v & unitDir_);
746 
747  // Quick rejection. Left of point1 endcap, or right of point2 endcap
748  if (parallel < 0 || parallel > magDir_)
749  {
750  volType[pointi] = volumeType::OUTSIDE;
751  }
752  else
753  {
754  // Remove the parallel component
755  v -= parallel*unitDir_;
756 
757  volType[pointi] =
758  (
759  mag(v) <= radius_
762  );
763  }
764  }
765 }
766 
767 
768 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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.
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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.
virtual const wordList & regions() const
Names of regions.
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
scalar y
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:316
const pointField & points
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< pointField > points() const
Get the points that define the surface.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Vector< scalar > vector
Definition: vector.H:57
A location inside the volume.
Definition: volumeType.H:65
A location outside the volume.
Definition: volumeType.H:66
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
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)
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
vector point
Point is a vector.
Definition: point.H:37
const dimensionedScalar c
Speed of light in a vacuum.
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
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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)
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127