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