searchableSphere.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-2020 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 Description
28  The search for nearest point on an ellipse or ellipsoid follows the
29  description given by Geometric Tools (David Eberly), which also
30  include some pseudo code. The content is CC-BY 4.0
31 
32 https://www.geometrictools.com/Documentation/DistancePointEllipseEllipsoid.pdf
33 
34  In the search algorithm, symmetry is exploited and the searching is
35  confined to the first (+x,+y,+z) octant, and the radii are ordered
36  from largest to smallest.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "searchableSphere.H"
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47  defineTypeNameAndDebug(searchableSphere, 0);
49  (
50  searchableSurface,
51  searchableSphere,
52  dict
53  );
55  (
56  searchableSurface,
57  searchableSphere,
58  dict,
59  sphere
60  );
61 }
62 
63 
64 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
65 
66 // General handling
67 namespace Foam
68 {
69 
70 // Dictionary entry with single scalar or vector quantity
71 inline static vector getRadius(const word& name, const dictionary& dict)
72 {
73  if (token(dict.lookup(name)).isNumber())
74  {
75  return vector::uniform(dict.get<scalar>(name));
76  }
77 
78  return dict.get<vector>(name);
79 }
80 
81 
82 // Test point for negative components, return the sign-changes
83 inline static unsigned getOctant(const point& p)
84 {
85  unsigned octant = 0;
86 
87  if (p.x() < 0) { octant |= 1; }
88  if (p.y() < 0) { octant |= 2; }
89  if (p.z() < 0) { octant |= 4; }
90 
91  return octant;
92 }
93 
94 
95 // Apply sign-changes to point
96 inline static void applyOctant(point& p, unsigned octant)
97 {
98  if (octant & 1) { p.x() = -p.x(); }
99  if (octant & 2) { p.y() = -p.y(); }
100  if (octant & 4) { p.z() = -p.z(); }
101 }
102 
103 
104 // Vector magnitudes
105 
106 inline static scalar vectorMagSqr
107 (
108  const scalar x,
109  const scalar y
110 )
111 {
112  return (sqr(x) + sqr(y));
113 }
114 
115 inline static scalar vectorMagSqr
116 (
117  const scalar x,
118  const scalar y,
119  const scalar z
120 )
121 {
122  return (sqr(x) + sqr(y) + sqr(z));
123 }
124 
125 inline static scalar vectorMag
126 (
127  const scalar x,
128  const scalar y
129 )
130 {
131  return hypot(x, y);
132 }
133 
134 inline static scalar vectorMag
135 (
136  const scalar x,
137  const scalar y,
138  const scalar z
139 )
140 {
142 }
143 
144 
145 } // End namespace Foam
146 
147 
148 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
149 
150 // Searching
151 namespace Foam
152 {
153 
154 // Max iterations for root finding
155 static constexpr int maxIters = 100;
156 
157 // Relative ellipse size within the root finding (1)
158 static constexpr scalar tolCloseness = 1e-3;
159 
160 
161 // Find root for distance to ellipse
162 static scalar findRootEllipseDistance
163 (
164  const scalar r0,
165  const scalar z0,
166  const scalar z1,
167  scalar g
168 )
169 {
170  const scalar n0 = r0*z0;
171 
172  scalar s0 = z1 - 1;
173  scalar s1 = (g < 0 ? 0 : vectorMag(n0, z1) - 1);
174  scalar s = 0;
175 
176  int nIters = 0;
177  while (nIters++ < maxIters)
178  {
179  s = (s0 + s1) / 2;
180  if (equal(s, s0) || equal(s, s1))
181  {
182  break;
183  }
184 
185  g = sqr(n0/(s+r0)) + sqr(z1/(s+1)) - 1;
186 
187  if (mag(g) < tolCloseness)
188  {
189  break;
190  }
191  else if (g > 0)
192  {
193  s0 = s;
194  }
195  else // g < 0
196  {
197  s1 = s;
198  }
199  }
200 
201  #ifdef FULLDEBUG
203  << "Located root in " << nIters << " iterations" << endl;
204  #endif
205 
206  return s;
207 }
208 
209 
210 // Find root for distance to ellipsoid
211 static scalar findRootEllipsoidDistance
212 (
213  const scalar r0,
214  const scalar r1,
215  const scalar z0,
216  const scalar z1,
217  const scalar z2,
218  scalar g
219 )
220 {
221  const scalar n0 = r0*z0;
222  const scalar n1 = r1*z1;
223 
224  scalar s0 = z2 - 1;
225  scalar s1 = (g < 0 ? 0 : vectorMag(n0, n1, z2) - 1);
226  scalar s = 0;
227 
228  int nIters = 0;
229  while (nIters++ < maxIters)
230  {
231  s = (s0 + s1) / 2;
232  if (equal(s, s0) || equal(s, s1))
233  {
234  break;
235  }
236 
237  g = vectorMagSqr(n0/(s+r0), n1/(s+r1), z2/(s+1)) - 1;
238 
239  if (mag(g) < tolCloseness)
240  {
241  break;
242  }
243  else if (g > 0)
244  {
245  s0 = s;
246  }
247  else // g < 0
248  {
249  s1 = s;
250  }
251  }
252 
253  #ifdef FULLDEBUG
255  << "root at " << s << " found in " << nIters
256  << " iterations" << endl;
257  #endif
258 
259  return s;
260 }
261 
262 
263 // Distance (squared) to an ellipse (2D)
264 static scalar distanceToEllipse
265 (
266  // [in] Ellipse characteristics. e0 >= e1
267  const scalar e0, const scalar e1,
268  // [in] search point. y0 >= 0, y1 >= 0
269  const scalar y0, const scalar y1,
270  // [out] nearest point on ellipse
271  scalar& x0, scalar& x1
272 )
273 {
274  if (equal(y1, 0))
275  {
276  // On the y1 = 0 axis
277 
278  const scalar numer0 = e0*y0;
279  const scalar denom0 = sqr(e0) - sqr(e1);
280 
281  if (numer0 < denom0)
282  {
283  const scalar xde0 = numer0/denom0; // Is always < 1
284 
285  x0 = e0*xde0;
286  x1 = e1*sqrt(1 - sqr(xde0));
287 
288  return vectorMagSqr((x0-y0), x1);
289  }
290 
291  // Fallthrough
292  x0 = e0;
293  x1 = 0;
294 
295  return sqr(y0-e0);
296  }
297  else if (equal(y0, 0))
298  {
299  // On the y0 = 0 axis, in the y1 > 0 half
300  x0 = 0;
301  x1 = e1;
302  return sqr(y1-e1);
303  }
304  else
305  {
306  // In the y0, y1 > 0 quadrant
307 
308  const scalar z0 = y0 / e0;
309  const scalar z1 = y1 / e1;
310  scalar eval = sqr(z0) + sqr(z1);
311 
312  scalar g = eval - 1;
313 
314  if (mag(g) < tolCloseness)
315  {
316  x0 = y0;
317  x1 = y1;
318 
319  if (!equal(eval, 1))
320  {
321  // Very close, scale accordingly.
322  eval = sqrt(eval);
323  x0 /= eval;
324  x1 /= eval;
325  }
326 
327  return sqr(x0-y0) + sqr(x1-y1);
328  }
329 
330 
331  // General search.
332  // Uses root find to get tbar of F(t) on (-e1^2,+inf)
333 
334  // Ratio major/minor
335  const scalar r0 = sqr(e0 / e1);
336 
337  const scalar sbar =
338  findRootEllipseDistance(r0, z0, z1, g);
339 
340  x0 = r0 * y0 / (sbar + r0);
341  x1 = y1 / (sbar + 1);
342 
343  // Re-evaluate
344  eval = sqr(x0/e0) + sqr(x1/e1);
345 
346  if (!equal(eval, 1))
347  {
348  // Very close, scale accordingly.
349  //
350  // This is not exact - the point is projected at a
351  // slight angle, but we are only correcting for
352  // rounding in the first place.
353 
354  eval = sqrt(eval);
355  x0 /= eval;
356  x1 /= eval;
357  }
358 
359  return sqr(x0-y0) + sqr(x1-y1);
360  }
361 
362  // Code never reaches here
364  << "Programming/logic error" << nl
365  << exit(FatalError);
366 
367  return 0;
368 }
369 
370 
371 // Distance (squared) to an ellipsoid (3D)
372 static scalar distanceToEllipsoid
373 (
374  // [in] Ellipsoid characteristics. e0 >= e1 >= e2
375  const scalar e0, const scalar e1, const scalar e2,
376  // [in] search point. y0 >= 0, y1 >= 0, y2 >= 0
377  const scalar y0, const scalar y1, const scalar y2,
378  // [out] nearest point on ellipsoid
379  scalar& x0, scalar& x1, scalar& x2
380 )
381 {
382  if (equal(y2, 0))
383  {
384  // On the y2 = 0 plane. Can use 2D ellipse finding
385 
386  const scalar numer0 = e0*y0;
387  const scalar numer1 = e1*y1;
388  const scalar denom0 = sqr(e0) - sqr(e2);
389  const scalar denom1 = sqr(e1) - sqr(e2);
390 
391  if (numer0 < denom0 && numer1 < denom1)
392  {
393  const scalar xde0 = numer0/denom0; // Is always < 1
394  const scalar xde1 = numer1/denom1; // Is always < 1
395 
396  const scalar disc = (1 - sqr(xde0) - sqr(xde1));
397 
398  if (disc > 0)
399  {
400  x0 = e0*xde0;
401  x1 = e1*xde1;
402  x2 = e2*sqrt(disc);
403 
404  return vectorMagSqr((x0-y0), (x1-y1), x2);
405  }
406  }
407 
408  // Fallthrough - use 2D form
409  x2 = 0;
410  return distanceToEllipse(e0,e1, y0,y1, x0,x1);
411  }
412  else if (equal(y1, 0))
413  {
414  // On the y1 = 0 plane, in the y2 > 0 half
415  x1 = 0;
416  if (equal(y0, 0))
417  {
418  x0 = 0;
419  x2 = e2;
420  return sqr(y2-e2);
421  }
422  else // y0 > 0
423  {
424  return distanceToEllipse(e0,e2, y0,y2, x0,x2);
425  }
426  }
427  else if (equal(y0, 0))
428  {
429  // On the y1 = 0 plane, in the y1, y2 > 0 quadrant
430  x0 = 0;
431  return distanceToEllipse(e1,e2, y1,y2, x1,x2);
432  }
433  else
434  {
435  // In the y0, y1, y2 > 0 octant
436 
437  const scalar z0 = y0/e0;
438  const scalar z1 = y1/e1;
439  const scalar z2 = y2/e2;
440  scalar eval = vectorMagSqr(z0, z1, z2);
441 
442  scalar g = eval - 1;
443 
444  if (mag(g) < tolCloseness)
445  {
446  x0 = y0;
447  x1 = y1;
448  x2 = y2;
449 
450  if (equal(eval, 1))
451  {
452  // Exactly on the ellipsoid - we are done
453  return 0;
454  }
455 
456  // Very close, scale accordingly.
457  eval = sqrt(eval);
458  x0 /= eval;
459  x1 /= eval;
460  x2 /= eval;
461 
462  return vectorMagSqr((x0-y0), (x1-y1), (x2-y2));
463  }
464 
465 
466  // General search.
467  // Compute the unique root tbar of F(t) on (-e2^2,+inf)
468 
469  const scalar r0 = sqr(e0/e2);
470  const scalar r1 = sqr(e1/e2);
471 
472  const scalar sbar =
473  findRootEllipsoidDistance(r0,r1, z0,z1,z2, g);
474 
475  x0 = r0*y0/(sbar+r0);
476  x1 = r1*y1/(sbar+r1);
477  x2 = y2/(sbar+1);
478 
479  // Reevaluate
480  eval = vectorMagSqr((x0/e0), (x1/e1), (x2/e2));
481 
482  if (!equal(eval, 1))
483  {
484  // Not exactly on ellipsoid?
485  //
486  // Scale accordingly. This is not exact - the point
487  // is projected at a slight angle, but we are only
488  // correcting for rounding in the first place.
489 
490  eval = sqrt(eval);
491  x0 /= eval;
492  x1 /= eval;
493  x2 /= eval;
494  }
495 
496  return vectorMagSqr((x0-y0), (x1-y1), (x2-y2));
497  }
498 
499  // Code never reaches here
501  << "Programming/logic error" << nl
502  << exit(FatalError);
503 
504  return 0;
505 }
506 
507 } // End namespace Foam
508 
509 
510 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
511 
512 inline Foam::searchableSphere::componentOrder
513 Foam::searchableSphere::getOrdering(const vector& radii)
514 {
515  #ifdef FULLDEBUG
516  for (direction cmpt = 0; cmpt < vector::nComponents; ++cmpt)
517  {
518  if (radii[cmpt] <= 0)
519  {
521  << "Radii must be positive, non-zero: " << radii << endl
522  << exit(FatalError);
523  }
524  }
525  #endif
526 
527  std::array<uint8_t, 3> idx{0, 1, 2};
528 
529  // Reverse sort by magnitude (largest first...)
530  // Radii are positive (checked above, or just always true)
531  std::stable_sort
532  (
533  idx.begin(),
534  idx.end(),
535  [&](uint8_t a, uint8_t b){ return radii[a] > radii[b]; }
536  );
537 
538  componentOrder order{idx[0], idx[1], idx[2], shapeType::GENERAL};
539 
540  if (equal(radii[order.major], radii[order.minor]))
541  {
542  order.shape = shapeType::SPHERE;
543  }
544  else if (equal(radii[order.major], radii[order.mezzo]))
545  {
546  order.shape = shapeType::OBLATE;
547  }
548  else if (equal(radii[order.mezzo], radii[order.minor]))
549  {
550  order.shape = shapeType::PROLATE;
551  }
552 
553  return order;
554 }
555 
556 
557 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
558 
559 Foam::pointIndexHit Foam::searchableSphere::findNearest
560 (
561  const point& sample,
562  const scalar nearestDistSqr
563 ) const
564 {
565  pointIndexHit info(false, sample, -1);
566 
567  // Handle special cases first
568 
569  if (order_.shape == shapeType::SPHERE)
570  {
571  // Point relative to origin, simultaneously the normal on the sphere
572  const vector n(sample - origin_);
573  const scalar magN = mag(n);
574 
575  // It is a sphere, all radii are identical
576 
577  if (nearestDistSqr >= sqr(magN - radii_[0]))
578  {
579  info.setPoint
580  (
581  origin_
582  + (magN < ROOTVSMALL ? vector(radii_[0],0,0) : (radii_[0]*n/magN))
583  );
584  info.setHit();
585  info.setIndex(0);
586  }
587 
588  return info;
589  }
590 
591 
592  //
593  // Non-sphere
594  //
595 
596  // Local point relative to the origin
597  vector relPt(sample - origin_);
598 
599  // Detect -ve octants
600  const unsigned octant = getOctant(relPt);
601 
602  // Flip everything into positive octant.
603  // That is what the algorithm expects.
604  applyOctant(relPt, octant);
605 
606 
607  // TODO - quick reject for things that are too far away
608 
609  point& near = info.rawPoint();
610  scalar distSqr{0};
611 
612  if (order_.shape == shapeType::OBLATE)
613  {
614  // Oblate (major = mezzo > minor) - use 2D algorithm
615  // Distance from the minor axis to relPt
616  const scalar axisDist = hypot(relPt[order_.major], relPt[order_.mezzo]);
617 
618  // Distance from the minor axis to near
619  scalar nearAxis;
620 
621  distSqr = distanceToEllipse
622  (
623  radii_[order_.major], radii_[order_.minor],
624  axisDist, relPt[order_.minor],
625  nearAxis, near[order_.minor]
626  );
627 
628  // Now nearAxis is the ratio, by which their components have changed
629  nearAxis /= (axisDist + VSMALL);
630 
631  near[order_.major] = relPt[order_.major] * nearAxis;
632  near[order_.mezzo] = relPt[order_.mezzo] * nearAxis;
633  // near[order_.minor] = already calculated
634  }
635  else if (order_.shape == shapeType::PROLATE)
636  {
637  // Prolate (major > mezzo = minor) - use 2D algorithm
638  // Distance from the major axis to relPt
639  const scalar axisDist = hypot(relPt[order_.mezzo], relPt[order_.minor]);
640 
641  // Distance from the major axis to near
642  scalar nearAxis;
643 
644  distSqr = distanceToEllipse
645  (
646  radii_[order_.major], radii_[order_.minor],
647  relPt[order_.major], axisDist,
648  near[order_.major], nearAxis
649  );
650 
651  // Now nearAxis is the ratio, by which their components have changed
652  nearAxis /= (axisDist + VSMALL);
653 
654  // near[order_.major] = already calculated
655  near[order_.mezzo] = relPt[order_.mezzo] * nearAxis;
656  near[order_.minor] = relPt[order_.minor] * nearAxis;
657  }
658  else // General case
659  {
660  distSqr = distanceToEllipsoid
661  (
662  radii_[order_.major], radii_[order_.mezzo], radii_[order_.minor],
663  relPt[order_.major], relPt[order_.mezzo], relPt[order_.minor],
664  near[order_.major], near[order_.mezzo], near[order_.minor]
665  );
666  }
667 
668  // Flip everything back to original octant
669  applyOctant(near, octant);
670 
671  // From local to global
672  near += origin_;
673 
674 
675  // Accept/reject based on distance
676  if (distSqr <= nearestDistSqr)
677  {
678  info.setHit();
679  }
680 
681  return info;
682 }
683 
684 
685 // From Graphics Gems - intersection of sphere with ray
686 void Foam::searchableSphere::findLineAll
687 (
688  const point& start,
689  const point& end,
690  pointIndexHit& near,
691  pointIndexHit& far
692 ) const
693 {
694  near.setMiss();
695  far.setMiss();
696 
697  if (order_.shape == shapeType::SPHERE)
698  {
699  vector dir(end-start);
700  const scalar magSqrDir = magSqr(dir);
701 
702  if (magSqrDir > ROOTVSMALL)
703  {
704  dir /= Foam::sqrt(magSqrDir);
705 
706  const vector relStart(start - origin_);
707 
708  const scalar v = -(relStart & dir);
709 
710  const scalar disc = sqr(radius()) - (magSqr(relStart) - sqr(v));
711 
712  if (disc >= 0)
713  {
714  const scalar d = Foam::sqrt(disc);
715 
716  const scalar nearParam = v - d;
717  const scalar farParam = v + d;
718 
719  if (nearParam >= 0 && sqr(nearParam) <= magSqrDir)
720  {
721  near.hitPoint(start + nearParam*dir, 0);
722  }
723 
724  if (farParam >= 0 && sqr(farParam) <= magSqrDir)
725  {
726  far.hitPoint(start + farParam*dir, 0);
727  }
728  }
729  }
730 
731  return;
732  }
733 
734 
735  // General case
736 
737  // Similar to intersection of sphere with ray (Graphics Gems),
738  // but we scale x/y/z components according to radii
739  // to have a unit spheroid for the interactions.
740  // When finished, we unscale to get the real points
741 
742  // Note - can also be used for the spherical case
743 
744  const point relStart = scalePoint(start);
745 
746  vector dir(scalePoint(end) - relStart);
747  const scalar magSqrDir = magSqr(dir);
748 
749  if (magSqrDir > ROOTVSMALL)
750  {
751  dir /= Foam::sqrt(magSqrDir);
752 
753  const scalar v = -(relStart & dir);
754 
755  const scalar disc = scalar(1) - (magSqr(relStart) - sqr(v));
756 
757  if (disc >= 0)
758  {
759  const scalar d = Foam::sqrt(disc);
760 
761  const scalar nearParam = v - d;
762  const scalar farParam = v + d;
763 
764  if (nearParam >= 0 && sqr(nearParam) <= magSqrDir)
765  {
766  near.hitPoint(unscalePoint(relStart + nearParam*dir), 0);
767  }
768  if (farParam >= 0 && sqr(farParam) <= magSqrDir)
769  {
770  far.hitPoint(unscalePoint(relStart + farParam*dir), 0);
771  }
772  }
773  }
774 }
775 
776 
777 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
778 
779 Foam::searchableSphere::searchableSphere
780 (
781  const IOobject& io,
782  const point& origin,
783  const scalar radius
784 )
785 :
786  searchableSphere(io, origin, vector::uniform(radius))
787 {}
788 
789 
790 Foam::searchableSphere::searchableSphere
791 (
792  const IOobject& io,
793  const point& origin,
794  const vector& radii
795 )
796 :
797  searchableSurface(io),
798  origin_(origin),
799  radii_(radii),
800  order_(getOrdering(radii_)) // NB: use () not {} for copy initialization
801 {
802  bounds().min() = (centre() - radii_);
803  bounds().max() = (centre() + radii_);
804 }
805 
806 
807 Foam::searchableSphere::searchableSphere
808 (
809  const IOobject& io,
810  const dictionary& dict
811 )
812 :
814  (
815  io,
816  dict.getCompat<vector>("origin", {{"centre", -1806}}),
817  getRadius("radius", dict)
818  )
819 {}
820 
821 
822 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
823 
825 (
826  const scalar theta,
827  const scalar phi
828 ) const
829 {
830  return point
831  (
832  origin_.x() + radii_.x() * cos(theta)*sin(phi),
833  origin_.y() + radii_.y() * sin(theta)*sin(phi),
834  origin_.z() + radii_.z() * cos(phi)
835  );
836 }
837 
838 
840 (
841  const scalar theta,
842  const scalar phi
843 ) const
844 {
845  // Normal is (x0/r0^2, x1/r1^2, x2/r2^2)
846 
847  return vector
848  (
849  cos(theta)*sin(phi) / radii_.x(),
850  sin(theta)*sin(phi) / radii_.y(),
851  cos(phi) / radii_.z()
852  ).normalise();
853 }
854 
855 
857 {
858  if (order_.shape == shapeType::SPHERE)
859  {
860  return bb.overlaps(origin_, sqr(radius()));
861  }
862 
863  if (!bb.valid())
864  {
865  return false;
866  }
867 
868  // Code largely as per
869  // boundBox::overlaps(const point& centre, const scalar radiusSqr)
870  // but normalized for a unit size
871 
872  // Find out where centre is in relation to bb.
873  // Find nearest point on bb.
874 
875  // Note: no major advantage in treating sphere specially
876 
877  scalar distSqr = 0;
878  for (direction dir = 0; dir < vector::nComponents; ++dir)
879  {
880  const scalar d0 = bb.min()[dir] - origin_[dir];
881  const scalar d1 = bb.max()[dir] - origin_[dir];
882 
883  if ((d0 > 0) == (d1 > 0))
884  {
885  // Both min/max are on the same side of the origin
886  // ie, box does not span spheroid in this direction
887 
888  if (Foam::mag(d0) < Foam::mag(d1))
889  {
890  distSqr += Foam::sqr(d0/radii_[dir]);
891  }
892  else
893  {
894  distSqr += Foam::sqr(d1/radii_[dir]);
895  }
896 
897  if (distSqr > 1)
898  {
899  return false;
900  }
901  }
902  }
903 
904  return true;
905 }
906 
907 
909 {
910  if (regions_.empty())
911  {
912  regions_.resize(1);
913  regions_.first() = "region0";
914  }
915  return regions_;
916 }
917 
918 
920 (
921  pointField& centres,
922  scalarField& radiusSqr
923 ) const
924 {
925  centres.resize(1);
926  radiusSqr.resize(1);
927 
928  centres[0] = origin_;
929  radiusSqr[0] = Foam::sqr(radius());
930 
931  // Add a bit to make sure all points are tested inside
932  radiusSqr += Foam::sqr(SMALL);
933 }
934 
935 
936 void Foam::searchableSphere::findNearest
937 (
938  const pointField& samples,
939  const scalarField& nearestDistSqr,
940  List<pointIndexHit>& info
941 ) const
942 {
943  info.resize(samples.size());
944 
945  forAll(samples, i)
946  {
947  info[i] = findNearest(samples[i], nearestDistSqr[i]);
948  }
949 }
950 
951 
953 (
954  const pointField& start,
955  const pointField& end,
956  List<pointIndexHit>& info
957 ) const
958 {
959  info.resize(start.size());
960 
962 
963  forAll(start, i)
964  {
965  // Pick nearest intersection.
966  // If none intersected take second one.
967 
968  findLineAll(start[i], end[i], info[i], b);
969 
970  if (!info[i].hit() && b.hit())
971  {
972  info[i] = b;
973  }
974  }
975 }
976 
977 
979 (
980  const pointField& start,
981  const pointField& end,
982  List<pointIndexHit>& info
983 ) const
984 {
985  info.resize(start.size());
986 
988 
989  forAll(start, i)
990  {
991  // Pick nearest intersection.
992  // Discard far intersection
993 
994  findLineAll(start[i], end[i], info[i], b);
995 
996  if (!info[i].hit() && b.hit())
997  {
998  info[i] = b;
999  }
1000  }
1001 }
1002 
1003 
1004 void Foam::searchableSphere::findLineAll
1006  const pointField& start,
1007  const pointField& end,
1008  List<List<pointIndexHit>>& info
1009 ) const
1010 {
1011  info.resize(start.size());
1012 
1013  forAll(start, i)
1014  {
1015  pointIndexHit near, far;
1016 
1017  findLineAll(start[i], end[i], near, far);
1018 
1019  if (near.hit())
1020  {
1021  if (far.hit())
1022  {
1023  info[i].resize(2);
1024  info[i][0] = near;
1025  info[i][1] = far;
1026  }
1027  else
1028  {
1029  info[i].resize(1);
1030  info[i][0] = near;
1031  }
1032  }
1033  else
1034  {
1035  if (far.hit())
1036  {
1037  info[i].resize(1);
1038  info[i][0] = far;
1039  }
1040  else
1041  {
1042  info[i].clear();
1043  }
1044  }
1045  }
1046 }
1047 
1048 
1051  const List<pointIndexHit>& info,
1052  labelList& region
1053 ) const
1054 {
1055  region.resize(info.size());
1056  region = 0;
1057 }
1058 
1059 
1062  const List<pointIndexHit>& info,
1063  vectorField& normal
1064 ) const
1065 {
1066  normal.resize(info.size());
1067 
1068  forAll(info, i)
1069  {
1070  if (info[i].hit())
1071  {
1072  if (order_.shape == shapeType::SPHERE)
1073  {
1074  // Special case (sphere)
1075  normal[i] = normalised(info[i].hitPoint() - origin_);
1076  }
1077  else
1078  {
1079  // General case
1080  // Normal is (x0/r0^2, x1/r1^2, x2/r2^2)
1081 
1082  normal[i] = scalePoint(info[i].hitPoint());
1083 
1084  normal[i].x() /= radii_.x();
1085  normal[i].y() /= radii_.y();
1086  normal[i].z() /= radii_.z();
1087  normal[i].normalise();
1088  }
1089  }
1090  else
1091  {
1092  // Set to what?
1093  normal[i] = Zero;
1094  }
1095  }
1096 }
1097 
1098 
1101  const pointField& points,
1102  List<volumeType>& volType
1103 ) const
1104 {
1105  volType.resize(points.size());
1106 
1107  if (order_.shape == shapeType::SPHERE)
1108  {
1109  // Special case. Minor advantage in treating specially
1110 
1111  const scalar rad2 = sqr(radius());
1112 
1113  forAll(points, pointi)
1114  {
1115  const point& p = points[pointi];
1116 
1117  volType[pointi] =
1118  (
1119  (magSqr(p - origin_) <= rad2)
1121  );
1122  }
1123 
1124  return;
1125  }
1126 
1127  // General case - could also do component-wise (manually)
1128  // Evaluate: (x/r0)^2 + (y/r1)^2 + (z/r2)^2 - 1 = 0
1129  // [sphere]: x^2 + y^2 + z^2 - R^2 = 0
1130 
1131  forAll(points, pointi)
1132  {
1133  const point p = scalePoint(points[pointi]);
1134 
1135  volType[pointi] =
1136  (
1137  (magSqr(p) <= 1)
1139  );
1140  }
1141 }
1142 
1143 
1144 // ************************************************************************* //
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
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:325
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::uniform
static Vector< Cmpt > uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
s
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))
Definition: gmvOutputSpray.H:25
Foam::searchableSphere::findLine
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
Definition: searchableSphere.C:953
Foam::y1
dimensionedScalar y1(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:282
searchableSphere.H
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::searchableSphere::getRegion
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
Definition: searchableSphere.C:1050
Foam::maxIters
static constexpr int maxIters
Definition: searchableSphere.C:155
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::applyOctant
static void applyOctant(point &p, unsigned octant)
Definition: searchableSphere.C:96
Foam::searchableSphere::regions
virtual const wordList & regions() const
Names of regions.
Definition: searchableSphere.C:908
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::token
A token holds an item read from Istream.
Definition: token.H:68
Foam::vectorMagSqr
static scalar vectorMagSqr(const scalar x, const scalar y)
Definition: searchableSphere.C:107
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::dictionary::getCompat
T getCompat(const word &keyword, std::initializer_list< std::pair< const char *, int >> compat, enum keyType::option=keyType::REGEX) const
Definition: dictionaryTemplates.C:108
Foam::distanceToEllipse
static scalar distanceToEllipse(const scalar e0, const scalar e1, const scalar y0, const scalar y1, scalar &x0, scalar &x1)
Definition: searchableSphere.C:265
Foam::getRadius
static vector getRadius(const word &name, const dictionary &dict)
Definition: searchableSphere.C:71
Foam::searchableSphere::surfacePoint
point surfacePoint(const scalar theta, const scalar phi) const
A point on the sphere at given location.
Definition: searchableSphere.C:825
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::hypot
dimensionedScalar hypot(const dimensionedScalar &x, const dimensionedScalar &y)
Definition: dimensionedScalar.C:327
Foam::boundBox::overlaps
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:221
sphere
Specialization of rigidBody to construct a sphere given the mass and radius.
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::searchableSphere::radii
const vector & radii() const noexcept
The radii of the spheroid.
Definition: searchableSphere.H:249
Foam::getOctant
static unsigned getOctant(const point &p)
Definition: searchableSphere.C:83
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:52
Foam::PointIndexHit::hit
bool hit() const noexcept
Is there a hit?
Definition: PointIndexHit.H:130
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::searchableSphere::getNormal
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
Definition: searchableSphere.C:1061
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::vectorMag
static scalar vectorMag(const scalar x, const scalar y)
Definition: searchableSphere.C:126
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:69
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
Foam::searchableSphere::surfaceNormal
vector surfaceNormal(const scalar theta, const scalar phi) const
Surface normal on the sphere at given location.
Definition: searchableSphere.C:840
Foam::y0
dimensionedScalar y0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:281
samples
scalarField samples(nIntervals, Zero)
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::searchableSphere::boundingSpheres
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
Definition: searchableSphere.C:920
Foam::List::resize
void resize(const label newSize)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::searchableSphere::findLineAny
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
Definition: searchableSphere.C:979
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:424
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
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.
Foam::searchableSphere::getVolumeType
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
Definition: searchableSphere.C:1100
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
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::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:483
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Vector< scalar >
Foam::distanceToEllipsoid
static scalar distanceToEllipsoid(const scalar e0, const scalar e1, const scalar e2, const scalar y0, const scalar y1, const scalar y2, scalar &x0, scalar &x1, scalar &x2)
Definition: searchableSphere.C:373
Foam::findRootEllipseDistance
static scalar findRootEllipseDistance(const scalar r0, const scalar z0, const scalar z1, scalar g)
Definition: searchableSphere.C:163
Foam::List< word >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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::direction
uint8_t direction
Definition: direction.H:52
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::findRootEllipsoidDistance
static scalar findRootEllipsoidDistance(const scalar r0, const scalar r1, const scalar z0, const scalar z1, const scalar z2, scalar g)
Definition: searchableSphere.C:212
Foam::volumeType::INSIDE
A location inside the volume.
Definition: volumeType.H:68
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::boundBox::valid
bool valid() const
Bounding box is non-inverted.
Definition: boundBoxI.H:76
Foam::searchableSphere::overlaps
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
Definition: searchableSphere.C:856
Foam::equal
bool equal(const T &s1, const T &s2)
Compare two values for equality.
Definition: doubleFloat.H:46
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::volumeType::OUTSIDE
A location outside the volume.
Definition: volumeType.H:69
Foam::tolCloseness
static constexpr scalar tolCloseness
Definition: searchableSphere.C:158
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::searchableSphere
Searching on general spheroid.
Definition: searchableSphere.H:92