Loading...
Searching...
No Matches
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-2022 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>
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
46namespace Foam
47{
50 (
53 dict
54 );
56 (
59 dict,
60 sphere
61 );
62}
63
64
65// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
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{
74 if (token(dict.lookup(name)).isNumber())
75 {
76 return vector::uniform(dict.get<scalar>(name));
77 }
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; }
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 * * * * * * * * * * * * * * //
151// Searching
152namespace Foam
154
155// Max iterations for root finding
156static constexpr int maxIters = 100;
158// Relative ellipse size within the root finding (1)
159static constexpr scalar tolCloseness = 1e-3;
160
161
162// Find root for distance to ellipse
163static scalar findRootEllipseDistance
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
207 return s;
208}
209
210
211// Find root for distance to ellipsoid
212static scalar findRootEllipsoidDistance
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
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);
368 return 0;
369}
370
371
372// Distance (squared) to an ellipsoid (3D)
373static scalar distanceToEllipsoid
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.hitPoint
581 (
582 origin_
583 + (magN < ROOTVSMALL ? vector(radii_[0],0,0) : (radii_[0]*n/magN))
584 );
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.point();
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
686void 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 }
774}
775
776
777// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
778
779Foam::searchableSphere::searchableSphere
780(
781 const IOobject& io,
782 const point& origin,
783 const scalar radius
785:
786 searchableSphere(io, origin, vector::uniform(radius))
787{}
788
789
790Foam::searchableSphere::searchableSphere
791(
792 const IOobject& io,
793 const point& origin,
794 const vector& radii
795)
796:
798 origin_(origin),
799 radii_(radii),
800 order_(getOrdering(radii_)) // NB: use () not {} for copy initialization
802 bounds().min() = (centre() - radii_);
803 bounds().max() = (centre() + radii_);
804}
805
806
807Foam::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)
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
856bool Foam::searchableSphere::overlaps(const boundBox& bb) const
857{
858 if (order_.shape == shapeType::SPHERE)
859 {
860 return bb.overlaps(origin_, sqr(radius()));
861 }
862
863 if (!bb.good())
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 }
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());
931 // Add a bit to make sure all points are tested inside
932 radiusSqr += Foam::sqr(SMALL);
933}
934
935
936void Foam::searchableSphere::findNearest
937(
938 const pointField& samples,
939 const scalarField& nearestDistSqr,
941) const
942{
943 info.resize(samples.size());
944
945 forAll(samples, i)
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,
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
1004void Foam::searchableSphere::findLineAll
1005(
1006 const pointField& start,
1007 const pointField& end,
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();
1044 }
1045 }
1046}
1047
1048
1050(
1051 const List<pointIndexHit>& info,
1052 labelList& region
1053) const
1055 region.resize(info.size());
1056 region = 0;
1057}
1058
1059
1061(
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].point() - 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].point());
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
1100(
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// ************************************************************************* //
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
void normalise()
Inplace normalise this field. Generally a no-op except for vector fields.
Definition Field.C:600
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition IOobject.H:1041
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
bool hit() const noexcept
Is there a hit?
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static constexpr direction nComponents
Number of components in this vector space.
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition VectorI.H:114
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
const point & max() const noexcept
Maximum describing the bounding box.
Definition boundBoxI.H:168
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition boundBoxI.H:439
const point & min() const noexcept
Minimum describing the bounding box.
Definition boundBoxI.H:162
bool good() const
Bounding box is non-inverted.
Definition boundBoxI.H:156
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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 bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
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.
scalar radius() const noexcept
The radius of the sphere, or major radius of the spheroid.
virtual tmp< pointField > points() const
Get the points that define the surface.
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.
point surfacePoint(const scalar theta, const scalar phi) const
A point on the sphere at given location.
@ 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:70
bool isNumber() const noexcept
Token is (signed/unsigned) integer type, FLOAT or DOUBLE.
Definition tokenI.H:992
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:66
@ INSIDE
A location inside the volume.
Definition volumeType.H:65
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
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))
#define InfoInFunction
Report an information message using Foam::Info.
const char * end
Definition SVGTools.H:223
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
static unsigned getOctant(const point &p)
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
dimensionedScalar y0(const dimensionedScalar &ds)
static vector getRadius(const word &name, const dictionary &dict)
dimensionedScalar sin(const dimensionedScalar &ds)
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition label.H:180
static scalar vectorMag(const scalar x, const scalar y)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static constexpr scalar tolCloseness
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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)
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
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:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
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)
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
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:50
dictionary dict
volScalarField & b
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
scalarField samples(nIntervals, Zero)