cylinder_searchableSurface.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  namespace searchableSurfaces
34  {
36 
38  (
40  cylinder,
42  );
43 
45  (
47  cylinder,
48  dictionary,
49  searchableCylinder,
50  "searchableCylinder"
51  );
52  }
53 }
54 
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
60 {
61  tmp<pointField> tCtrs(new pointField(1, 0.5*(point1_ + point2_)));
62 
63  return tCtrs;
64 }
65 
66 
68 (
69  pointField& centres,
70  scalarField& radiusSqr
71 ) const
72 {
73  centres.setSize(1);
74  centres[0] = 0.5*(point1_ + point2_);
75 
76  radiusSqr.setSize(1);
77  radiusSqr[0] = Foam::magSqr(point1_-centres[0]) + Foam::sqr(radius_);
78 
79  // Add a bit to make sure all points are tested inside
80  radiusSqr += Foam::sqr(small);
81 }
82 
83 
85 {
86  tmp<pointField> tPts(new pointField(2));
87  pointField& pts = tPts.ref();
88 
89  pts[0] = point1_;
90  pts[1] = point2_;
91 
92  return tPts;
93 }
94 
95 
96 Foam::pointIndexHit Foam::searchableSurfaces::cylinder::findNearest
97 (
98  const point& sample,
99  const scalar nearestDistSqr
100 ) const
101 {
102  pointIndexHit info(false, sample, -1);
103 
104  vector v(sample - point1_);
105 
106  // Decompose sample-point1 into normal and parallel component
107  scalar parallel = (v & unitDir_);
108 
109  // Remove the parallel component and normalise
110  v -= parallel*unitDir_;
111  scalar magV = mag(v);
112 
113  if (magV < rootVSmall)
114  {
115  v = Zero;
116  }
117  else
118  {
119  v /= magV;
120  }
121 
122  if (parallel <= 0)
123  {
124  // nearest is at point1 end cap. Clip to radius.
125  info.setPoint(point1_ + min(magV, radius_)*v);
126  }
127  else if (parallel >= magDir_)
128  {
129  // nearest is at point2 end cap. Clip to radius.
130  info.setPoint(point2_ + min(magV, radius_)*v);
131  }
132  else
133  {
134  // in between endcaps. Might either be nearer endcaps or cylinder wall
135 
136  // distance to endpoint: parallel or parallel-magDir
137  // distance to cylinder wall: magV-radius_
138 
139  // Nearest cylinder point
140  point cylPt;
141  if (magV < rootVSmall)
142  {
143  // Point exactly on centre line. Take any point on wall.
144  vector e1 = point(1,0,0) ^ unitDir_;
145  scalar magE1 = mag(e1);
146  if (magE1 < small)
147  {
148  e1 = point(0,1,0) ^ unitDir_;
149  magE1 = mag(e1);
150  }
151  e1 /= magE1;
152  cylPt = sample + radius_*e1;
153  }
154  else
155  {
156  cylPt = sample + (radius_-magV)*v;
157  }
158 
159  if (parallel < 0.5*magDir_)
160  {
161  // Project onto p1 endcap
162  point end1Pt = point1_ + min(magV, radius_)*v;
163 
164  if (magSqr(sample-cylPt) < magSqr(sample-end1Pt))
165  {
166  info.setPoint(cylPt);
167  }
168  else
169  {
170  info.setPoint(end1Pt);
171  }
172  }
173  else
174  {
175  // Project onto p2 endcap
176  point end2Pt = point2_ + min(magV, radius_)*v;
177 
178  if (magSqr(sample-cylPt) < magSqr(sample-end2Pt))
179  {
180  info.setPoint(cylPt);
181  }
182  else
183  {
184  info.setPoint(end2Pt);
185  }
186  }
187  }
188 
189  if (magSqr(sample - info.rawPoint()) < nearestDistSqr)
190  {
191  info.setHit();
192  info.setIndex(0);
193  }
194 
195  return info;
196 }
197 
198 
199 Foam::scalar Foam::searchableSurfaces::cylinder::radius2(const point& pt) const
200 {
201  const vector x = (pt-point1_) ^ unitDir_;
202  return x&x;
203 }
204 
205 
206 // From http://www.gamedev.net/community/forums/topic.asp?topic_id=467789 -
207 // intersection of cylinder with ray
208 void Foam::searchableSurfaces::cylinder::findLineAll
209 (
210  const point& start,
211  const point& end,
212  pointIndexHit& near,
213  pointIndexHit& far
214 ) const
215 {
216  near.setMiss();
217  far.setMiss();
218 
219  vector point1Start(start-point1_);
220  vector point2Start(start-point2_);
221  vector point1End(end-point1_);
222 
223  // Quick rejection of complete vector outside endcaps
224  scalar s1 = point1Start&unitDir_;
225  scalar s2 = point1End&unitDir_;
226 
227  if ((s1 < 0 && s2 < 0) || (s1 > magDir_ && s2 > magDir_))
228  {
229  return;
230  }
231 
232  // Line as P = start+t*V where V is unit vector and t=[0..mag(end-start)]
233  vector V(end-start);
234  scalar magV = mag(V);
235  if (magV < rootVSmall)
236  {
237  return;
238  }
239  V /= magV;
240 
241 
242  // We now get the nearest intersections to start. This can either be
243  // the intersection with the end plane or with the cylinder side.
244 
245  // Get the two points (expressed in t) on the end planes. This is to
246  // clip any cylinder intersection against.
247  scalar tPoint1;
248  scalar tPoint2;
249 
250  // Maintain the two intersections with the endcaps
251  scalar tNear = vGreat;
252  scalar tFar = vGreat;
253 
254  {
255  scalar s = (V&unitDir_);
256  if (mag(s) > vSmall)
257  {
258  tPoint1 = -s1/s;
259  tPoint2 = -(point2Start&unitDir_)/s;
260  if (tPoint2 < tPoint1)
261  {
262  Swap(tPoint1, tPoint2);
263  }
264  if (tPoint1 > magV || tPoint2 < 0)
265  {
266  return;
267  }
268 
269  // See if the points on the endcaps are actually inside the cylinder
270  if (tPoint1 >= 0 && tPoint1 <= magV)
271  {
272  if (radius2(start+tPoint1*V) <= sqr(radius_))
273  {
274  tNear = tPoint1;
275  }
276  }
277  if (tPoint2 >= 0 && tPoint2 <= magV)
278  {
279  if (radius2(start+tPoint2*V) <= sqr(radius_))
280  {
281  // Check if already have a near hit from point1
282  if (tNear <= 1)
283  {
284  tFar = tPoint2;
285  }
286  else
287  {
288  tNear = tPoint2;
289  }
290  }
291  }
292  }
293  else
294  {
295  // Vector perpendicular to cylinder. Check for outside already done
296  // above so just set tpoint to allow all.
297  tPoint1 = -vGreat;
298  tPoint2 = vGreat;
299  }
300  }
301 
302 
303  const vector x = point1Start ^ unitDir_;
304  const vector y = V ^ unitDir_;
305  const scalar d = sqr(radius_);
306 
307  // Second order equation of the form a*t^2 + b*t + c
308  const scalar a = (y&y);
309  const scalar b = 2*(x&y);
310  const scalar c = (x&x)-d;
311 
312  const scalar disc = b*b-4*a*c;
313 
314  scalar t1 = -vGreat;
315  scalar t2 = vGreat;
316 
317  if (disc < 0)
318  {
319  // Fully outside
320  return;
321  }
322  else if (disc < rootVSmall)
323  {
324  // Single solution
325  if (mag(a) > rootVSmall)
326  {
327  t1 = -b/(2*a);
328 
329  // Pout<< "single solution t:" << t1
330  // << " for start:" << start << " end:" << end
331  // << " c:" << c << endl;
332 
333  if (t1 >= 0 && t1 <= magV && t1 >= tPoint1 && t1 <= tPoint2)
334  {
335  // valid. Insert sorted.
336  if (t1 < tNear)
337  {
338  tFar = tNear;
339  tNear = t1;
340  }
341  else if (t1 < tFar)
342  {
343  tFar = t1;
344  }
345  }
346  else
347  {
348  return;
349  }
350  }
351  else
352  {
353  // Aligned with axis. Check if outside radius
354  // Pout<< "small discriminant:" << disc
355  // << " for start:" << start << " end:" << end
356  // << " magV:" << magV
357  // << " c:" << c << endl;
358  if (c > 0)
359  {
360  return;
361  }
362  }
363  }
364  else
365  {
366  if (mag(a) > rootVSmall)
367  {
368  scalar sqrtDisc = sqrt(disc);
369 
370  t1 = (-b - sqrtDisc)/(2*a);
371  t2 = (-b + sqrtDisc)/(2*a);
372  if (t2 < t1)
373  {
374  Swap(t1, t2);
375  }
376 
377  if (t1 >= 0 && t1 <= magV && t1 >= tPoint1 && t1 <= tPoint2)
378  {
379  // valid. Insert sorted.
380  if (t1 < tNear)
381  {
382  tFar = tNear;
383  tNear = t1;
384  }
385  else if (t1 < tFar)
386  {
387  tFar = t1;
388  }
389  }
390  if (t2 >= 0 && t2 <= magV && t2 >= tPoint1 && t2 <= tPoint2)
391  {
392  // valid. Insert sorted.
393  if (t2 < tNear)
394  {
395  tFar = tNear;
396  tNear = t2;
397  }
398  else if (t2 < tFar)
399  {
400  tFar = t2;
401  }
402  }
403  // Pout<< "two solutions t1:" << t1 << " t2:" << t2
404  // << " for start:" << start << " end:" << end
405  // << " magV:" << magV
406  // << " c:" << c << endl;
407  }
408  else
409  {
410  // Aligned with axis. Check if outside radius
411  // Pout<< "large discriminant:" << disc
412  // << " small a:" << a
413  // << " for start:" << start << " end:" << end
414  // << " magV:" << magV
415  // << " c:" << c << endl;
416  if (c > 0)
417  {
418  return;
419  }
420  }
421  }
422 
423  // Check tNear, tFar
424  if (tNear >= 0 && tNear <= magV)
425  {
426  near.setPoint(start+tNear*V);
427  near.setHit();
428  near.setIndex(0);
429 
430  if (tFar <= magV)
431  {
432  far.setPoint(start+tFar*V);
433  far.setHit();
434  far.setIndex(0);
435  }
436  }
437  else if (tFar >= 0 && tFar <= magV)
438  {
439  near.setPoint(start+tFar*V);
440  near.setHit();
441  near.setIndex(0);
442  }
443 }
444 
445 
446 Foam::boundBox Foam::searchableSurfaces::cylinder::calcBounds() const
447 {
448 
449  // Adapted from
450  // http://www.gamedev.net/community/forums
451  // /topic.asp?topic_id=338522&forum_id=20&gforum_id=0
452 
453  // Let cylinder have end points A,B and radius r,
454 
455  // Bounds in direction X (same for Y and Z) can be found as:
456  // Let A.X<B.X (otherwise swap points)
457  // Good approximate lowest bound is A.X-r and highest is B.X+r (precise for
458  // capsule). At worst, in one direction it can be larger than needed by 2*r.
459 
460  // Accurate bounds for cylinder is
461  // A.X-kx*r, B.X+kx*r
462  // where
463  // kx=sqrt(((A.Y-B.Y)^2+(A.Z-B.Z)^2)/((A.X-B.X)^2+(A.Y-B.Y)^2+(A.Z-B.Z)^2))
464 
465  // similar thing for Y and Z
466  // (i.e.
467  // ky=sqrt(((A.X-B.X)^2+(A.Z-B.Z)^2)/((A.X-B.X)^2+(A.Y-B.Y)^2+(A.Z-B.Z)^2))
468  // kz=sqrt(((A.X-B.X)^2+(A.Y-B.Y)^2)/((A.X-B.X)^2+(A.Y-B.Y)^2+(A.Z-B.Z)^2))
469  // )
470 
471  // How derived: geometric reasoning. Bounds of cylinder is same as for 2
472  // circles centered on A and B. This sqrt thingy gives sine of angle between
473  // axis and direction, used to find projection of radius.
474 
475  vector kr
476  (
477  sqrt(sqr(unitDir_.y()) + sqr(unitDir_.z())),
478  sqrt(sqr(unitDir_.x()) + sqr(unitDir_.z())),
479  sqrt(sqr(unitDir_.x()) + sqr(unitDir_.y()))
480  );
481 
482  kr *= radius_;
483 
484  point min = point1_ - kr;
485  point max = point1_ + kr;
486 
487  min = ::Foam::min(min, point2_ - kr);
488  max = ::Foam::max(max, point2_ + kr);
489 
490  return boundBox(min, max);
491 }
492 
493 
494 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
495 
497 (
498  const IOobject& io,
499  const point& point1,
500  const point& point2,
501  const scalar radius
502 )
503 :
504  searchableSurface(io),
505  point1_(point1),
506  point2_(point2),
507  magDir_(mag(point2_-point1_)),
508  unitDir_((point2_-point1_)/magDir_),
509  radius_(radius)
510 {
511  bounds() = calcBounds();
512 }
513 
514 
516 (
517  const IOobject& io,
518  const dictionary& dict
519 )
520 :
521  searchableSurface(io),
522  point1_(dict.lookup<point>("point1", dimLength)),
523  point2_(dict.lookup<point>("point2", dimLength)),
524  magDir_(mag(point2_ - point1_)),
525  unitDir_((point2_ - point1_)/magDir_),
526  radius_(dict.lookup<scalar>("radius", dimLength))
527 {
528  bounds() = calcBounds();
529 }
530 
531 
532 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
533 
535 {}
536 
537 
538 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
539 
541 {
542  if (regions_.empty())
543  {
544  regions_.setSize(1);
545  regions_[0] = "region0";
546  }
547  return regions_;
548 }
549 
550 
551 void Foam::searchableSurfaces::cylinder::findNearest
552 (
553  const pointField& samples,
554  const scalarField& nearestDistSqr,
555  List<pointIndexHit>& info
556 ) const
557 {
558  info.setSize(samples.size());
559 
560  forAll(samples, i)
561  {
562  info[i] = findNearest(samples[i], nearestDistSqr[i]);
563  }
564 }
565 
566 
568 (
569  const pointField& start,
570  const pointField& end,
571  List<pointIndexHit>& info
572 ) const
573 {
574  info.setSize(start.size());
575 
576  forAll(start, i)
577  {
578  // Pick nearest intersection. If none intersected take second one.
580  findLineAll(start[i], end[i], info[i], b);
581  if (!info[i].hit() && b.hit())
582  {
583  info[i] = b;
584  }
585  }
586 }
587 
588 
590 (
591  const pointField& start,
592  const pointField& end,
593  List<pointIndexHit>& info
594 ) const
595 {
596  info.setSize(start.size());
597 
598  forAll(start, i)
599  {
600  // Discard far intersection
602  findLineAll(start[i], end[i], info[i], b);
603  if (!info[i].hit() && b.hit())
604  {
605  info[i] = b;
606  }
607  }
608 }
609 
610 
611 void Foam::searchableSurfaces::cylinder::findLineAll
612 (
613  const pointField& start,
614  const pointField& end,
616 ) const
617 {
618  info.setSize(start.size());
619 
620  forAll(start, i)
621  {
622  pointIndexHit near, far;
623  findLineAll(start[i], end[i], near, far);
624 
625  if (near.hit())
626  {
627  if (far.hit())
628  {
629  info[i].setSize(2);
630  info[i][0] = near;
631  info[i][1] = far;
632  }
633  else
634  {
635  info[i].setSize(1);
636  info[i][0] = near;
637  }
638  }
639  else
640  {
641  if (far.hit())
642  {
643  info[i].setSize(1);
644  info[i][0] = far;
645  }
646  else
647  {
648  info[i].clear();
649  }
650  }
651  }
652 }
653 
654 
656 (
657  const List<pointIndexHit>& info,
658  labelList& region
659 ) const
660 {
661  region.setSize(info.size());
662  region = 0;
663 }
664 
665 
667 (
668  const List<pointIndexHit>& info,
669  vectorField& normal
670 ) const
671 {
672  normal.setSize(info.size());
673  normal = Zero;
674 
675  forAll(info, i)
676  {
677  if (info[i].hit())
678  {
679  vector v(info[i].hitPoint() - point1_);
680 
681  // Decompose sample-point1 into normal and parallel component
682  scalar parallel = (v & unitDir_);
683 
684  // Remove the parallel component and normalise
685  v -= parallel*unitDir_;
686  scalar magV = mag(v);
687 
688  if (parallel <= 0)
689  {
690  if ((magV-radius_) < mag(parallel))
691  {
692  // either above endcap (magV<radius) or outside but closer
693  normal[i] = -unitDir_;
694  }
695  else
696  {
697  normal[i] = v/magV;
698  }
699  }
700  else if (parallel <= 0.5*magDir_)
701  {
702  // See if endcap closer or sidewall
703  if (magV >= radius_ || (radius_-magV) < parallel)
704  {
705  normal[i] = v/magV;
706  }
707  else
708  {
709  // closer to endcap
710  normal[i] = -unitDir_;
711  }
712  }
713  else if (parallel <= magDir_)
714  {
715  // See if endcap closer or sidewall
716  if (magV >= radius_ || (radius_-magV) < (magDir_-parallel))
717  {
718  normal[i] = v/magV;
719  }
720  else
721  {
722  // closer to endcap
723  normal[i] = unitDir_;
724  }
725  }
726  else // beyond cylinder
727  {
728  if ((magV-radius_) < (parallel-magDir_))
729  {
730  // above endcap
731  normal[i] = unitDir_;
732  }
733  else
734  {
735  normal[i] = v/magV;
736  }
737  }
738  }
739  }
740 }
741 
742 
744 (
745  const pointField& points,
746  List<volumeType>& volType
747 ) const
748 {
749  volType.setSize(points.size());
750  volType = volumeType::inside;
751 
752  forAll(points, pointi)
753  {
754  const point& pt = points[pointi];
755 
756  vector v(pt - point1_);
757 
758  // Decompose sample-point1 into normal and parallel component
759  scalar parallel = v & unitDir_;
760 
761  if (parallel < 0)
762  {
763  // left of point1 endcap
764  volType[pointi] = volumeType::outside;
765  }
766  else if (parallel > magDir_)
767  {
768  // right of point2 endcap
769  volType[pointi] = volumeType::outside;
770  }
771  else
772  {
773  // Remove the parallel component
774  v -= parallel*unitDir_;
775 
776  if (mag(v) > radius_)
777  {
778  volType[pointi] = volumeType::outside;
779  }
780  else
781  {
782  volType[pointi] = volumeType::inside;
783  }
784  }
785  }
786 }
787 
788 
789 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:54
bool hit() const
Is there a hit.
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
const boundBox & bounds() const
Return const reference to boundBox.
Surface geometry with a cylinder shape, which can be used with snappyHexMesh.
virtual void getVolumeType(const pointField &, List< volumeType > &) const
Determine type (inside/outside/mixed) for point. unknown if.
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 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.
cylinder(const IOobject &io, const point &, const point &, const scalar radius)
Construct from components.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
virtual tmp< pointField > points() const
Get the points that define the surface.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
const dimensionedScalar c
Speed of light in a vacuum.
addToRunTimeSelectionTable(searchableSurface, box, dictionary)
addBackwardCompatibleToRunTimeSelectionTable(searchableSurface, box, dictionary, searchableBox, "searchableBox")
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
const dimensionSet dimLength
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
vector point
Point is a vector.
Definition: point.H:41
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void Swap(T &a, T &b)
Definition: Swap.H:43
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict
scalarField samples(nIntervals, 0)