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