Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoParaboloid.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Mihaela Gheata 20/06/04
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TGeoParaboloid
13\ingroup Shapes_classes
14
15A paraboloid is defined by the revolution surface generated by a
16parabola and is bounded by two planes perpendicular to Z axis. The
17parabola equation is taken in the form: `z = a·r2 + b`, where:
18`r2 = x2 + y2`. Note the missing linear term (parabola symmetric with
19respect to Z axis).
20
21The coefficients a and b are computed from the input values which are
22the radii of the circular sections cut by the planes at `+/-dz`:
23
24 - `-dz = a*r2low + b`
25 - ` dz = a*r2high + b`
26
27~~~{.cpp}
28TGeoParaboloid(Double_t rlo,Double_t rhi,Double_t dz);
29~~~
30
31Begin_Macro
32{
33 new TGeoManager("parab", "paraboloid");
34 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
35 TGeoMedium *med = new TGeoMedium("MED",1,mat);
36 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
37 gGeoManager->SetTopVolume(top);
38 TGeoVolume *vol = gGeoManager->MakeParaboloid("PARAB",med,0, 40, 50);
39 TGeoParaboloid *par = (TGeoParaboloid*)vol->GetShape();
40 top->AddNode(vol,1);
41 gGeoManager->CloseGeometry();
42 gGeoManager->SetNsegments(80);
43 top->Draw();
44 TView *view = gPad->GetView();
45 view->ShowAxis();
46}
47End_Macro
48*/
49
50#include <iostream>
51#include "TGeoManager.h"
52#include "TGeoVolume.h"
53#include "TVirtualGeoPainter.h"
54#include "TGeoParaboloid.h"
55#include "TBuffer3D.h"
56#include "TBuffer3DTypes.h"
57#include "TMath.h"
58
60
61////////////////////////////////////////////////////////////////////////////////
62/// Dummy constructor
63
65{
66 fRlo = 0;
67 fRhi = 0;
68 fDz = 0;
69 fA = 0;
70 fB = 0;
72}
73
74////////////////////////////////////////////////////////////////////////////////
75/// Default constructor specifying X and Y semiaxis length
76
78 :TGeoBBox(0,0,0)
79{
80 fRlo = 0;
81 fRhi = 0;
82 fDz = 0;
83 fA = 0;
84 fB = 0;
86 SetParaboloidDimensions(rlo, rhi, dz);
88}
89
90////////////////////////////////////////////////////////////////////////////////
91/// Default constructor specifying X and Y semiaxis length
92
94 :TGeoBBox(name, 0, 0, 0)
95{
96 fRlo = 0;
97 fRhi = 0;
98 fDz = 0;
99 fA = 0;
100 fB = 0;
102 SetParaboloidDimensions(rlo, rhi, dz);
103 ComputeBBox();
104}
105
106////////////////////////////////////////////////////////////////////////////////
107/// Default constructor specifying minimum and maximum radius
108/// - param[0] = rlo
109/// - param[1] = rhi
110/// - param[2] = dz
111
113{
115 SetDimensions(param);
116 ComputeBBox();
117}
118
119////////////////////////////////////////////////////////////////////////////////
120/// destructor
121
123{
124}
125
126////////////////////////////////////////////////////////////////////////////////
127/// Computes capacity of the shape in [length^3]
128
130{
131 Double_t capacity = TMath::Pi()*fDz*(fRlo*fRlo+fRhi*fRhi);
132 return capacity;
133}
134
135////////////////////////////////////////////////////////////////////////////////
136/// compute bounding box of the tube
137
139{
141 fDY = fDX;
142 fDZ = fDz;
143}
144
145////////////////////////////////////////////////////////////////////////////////
146/// Compute normal to closest surface from POINT.
147
148void TGeoParaboloid::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
149{
150 norm[0] = norm[1] = 0.0;
151 if (TMath::Abs(point[2]) > fDz) {
152 norm[2] = TMath::Sign(1., dir[2]);
153 return;
154 }
155 Double_t safz = fDz-TMath::Abs(point[2]);
156 Double_t r = TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
157 Double_t safr = TMath::Abs(r-TMath::Sqrt((point[2]-fB)/fA));
158 if (safz<safr) {
159 norm[2] = TMath::Sign(1., dir[2]);
160 return;
161 }
162 Double_t talf = -2.*fA*r;
163 Double_t calf = 1./TMath::Sqrt(1.+talf*talf);
164 Double_t salf = talf * calf;
165 Double_t phi = TMath::ATan2(point[1], point[0]);
166
167 norm[0] = salf*TMath::Cos(phi);
168 norm[1] = salf*TMath::Sin(phi);
169 norm[2] = calf;
170 Double_t ndotd = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
171 if (ndotd < 0) {
172 norm[0] = -norm[0];
173 norm[1] = -norm[1];
174 norm[2] = -norm[2];
175 }
176}
177
178////////////////////////////////////////////////////////////////////////////////
179/// test if point is inside the elliptical tube
180
182{
183 if (TMath::Abs(point[2])>fDz) return kFALSE;
184 Double_t aa = fA*(point[2]-fB);
185 if (aa < 0) return kFALSE;
186 Double_t rsq = point[0]*point[0]+point[1]*point[1];
187 if (aa < fA*fA*rsq) return kFALSE;
188 return kTRUE;
189}
190
191////////////////////////////////////////////////////////////////////////////////
192/// compute closest distance from point px,py to each vertex
193
195{
197 const Int_t numPoints=n*(n+1)+2;
198 return ShapeDistancetoPrimitive(numPoints, px, py);
199}
200
201////////////////////////////////////////////////////////////////////////////////
202/// Compute distance from a point to the parabola given by:
203/// `z = a*rsq + b; rsq = x*x+y*y`
204
206{
207 Double_t rsq = point[0]*point[0]+point[1]*point[1];
208 Double_t a = fA * (dir[0]*dir[0] + dir[1]*dir[1]);
209 Double_t b = 2.*fA*(point[0]*dir[0]+point[1]*dir[1])-dir[2];
210 Double_t c = fA*rsq + fB - point[2];
211 Double_t dist = TGeoShape::Big();
213 if (TMath::Abs(b)<TGeoShape::Tolerance()) return dist; // big
214 dist = -c/b;
215 if (dist < 0) return TGeoShape::Big();
216 return dist; // OK
217 }
218 Double_t ainv = 1./a;
219 Double_t sum = - b*ainv;
220 Double_t prod = c*ainv;
221 Double_t delta = sum*sum - 4.*prod;
222 if (delta<0) return dist; // big
223 delta = TMath::Sqrt(delta);
224 Double_t sone = TMath::Sign(1.,ainv);
225 Int_t i = -1;
226 while (i<2) {
227 dist = 0.5*(sum+i*sone*delta);
228 i += 2;
229 if (dist<0) continue;
230 if (dist<1.E-8) {
231 Double_t talf = -2.*fA*TMath::Sqrt(rsq);
232 Double_t phi = TMath::ATan2(point[1], point[0]);
233 Double_t ndotd = talf*(TMath::Cos(phi)*dir[0]+TMath::Sin(phi)*dir[1])+dir[2];
234 if (!in) ndotd *= -1;
235 if (ndotd<0) return dist;
236 } else return dist;
237 }
238 return TGeoShape::Big();
239}
240
241////////////////////////////////////////////////////////////////////////////////
242/// compute distance from inside point to surface of the paraboloid
243
244Double_t TGeoParaboloid::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
245{
246 if (iact<3 && safe) {
247 // compute safe distance
248 *safe = Safety(point, kTRUE);
249 if (iact==0) return TGeoShape::Big();
250 if (iact==1 && step<*safe) return TGeoShape::Big();
251 }
252
254 if (dir[2]<0) {
255 dz = -(point[2]+fDz)/dir[2];
256 } else if (dir[2]>0) {
257 dz = (fDz-point[2])/dir[2];
258 }
259 Double_t dpara = DistToParaboloid(point, dir, kTRUE);
260 return TMath::Min(dz, dpara);
261}
262
263////////////////////////////////////////////////////////////////////////////////
264/// compute distance from outside point to surface of the paraboloid and safe distance
265
266Double_t TGeoParaboloid::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
267{
268 if (iact<3 && safe) {
269 // compute safe distance
270 *safe = Safety(point, kFALSE);
271 if (iact==0) return TGeoShape::Big();
272 if (iact==1 && step<*safe) return TGeoShape::Big();
273 }
274 Double_t xnew, ynew, znew;
275 if (point[2]<=-fDz) {
276 if (dir[2]<=0) return TGeoShape::Big();
277 Double_t snxt = -(fDz+point[2])/dir[2];
278 // find extrapolated X and Y
279 xnew = point[0]+snxt*dir[0];
280 ynew = point[1]+snxt*dir[1];
281 if ((xnew*xnew+ynew*ynew) <= fRlo*fRlo) return snxt;
282 } else if (point[2]>=fDz) {
283 if (dir[2]>=0) return TGeoShape::Big();
284 Double_t snxt = (fDz-point[2])/dir[2];
285 // find extrapolated X and Y
286 xnew = point[0]+snxt*dir[0];
287 ynew = point[1]+snxt*dir[1];
288 if ((xnew*xnew+ynew*ynew) <= fRhi*fRhi) return snxt;
289 }
290 Double_t snxt = DistToParaboloid(point, dir, kFALSE);
291 if (snxt > 1E20) return snxt;
292 znew = point[2]+snxt*dir[2];
293 if (TMath::Abs(znew) <= fDz) return snxt;
294 return TGeoShape::Big();
295}
296
297////////////////////////////////////////////////////////////////////////////////
298/// Divide the paraboloid along one axis.
299
300TGeoVolume *TGeoParaboloid::Divide(TGeoVolume * /*voldiv*/, const char * /*divname*/, Int_t /*iaxis*/, Int_t /*ndiv*/,
301 Double_t /*start*/, Double_t /*step*/)
302{
303 Error("Divide", "Paraboloid divisions not implemented");
304 return 0;
305}
306
307////////////////////////////////////////////////////////////////////////////////
308/// Fill vector param[4] with the bounding cylinder parameters. The order
309/// is the following : Rmin, Rmax, Phi1, Phi2
310
312{
313 param[0] = 0.; // Rmin
314 param[1] = fDX; // Rmax
315 param[1] *= param[1];
316 param[2] = 0.; // Phi1
317 param[3] = 360.; // Phi2
318}
319
320////////////////////////////////////////////////////////////////////////////////
321/// in case shape has some negative parameters, these has to be computed
322/// in order to fit the mother
323
325{
326 return 0;
327}
328
329////////////////////////////////////////////////////////////////////////////////
330/// print shape parameters
331
333{
334 printf("*** Shape %s: TGeoParaboloid ***\n", GetName());
335 printf(" rlo = %11.5f\n", fRlo);
336 printf(" rhi = %11.5f\n", fRhi);
337 printf(" dz = %11.5f\n", fDz);
338 printf(" Bounding box:\n");
340}
341
342////////////////////////////////////////////////////////////////////////////////
343/// Creates a TBuffer3D describing *this* shape.
344/// Coordinates are in local reference frame.
345
347{
349 Int_t nbPnts = n*(n+1)+2;
350 Int_t nbSegs = n*(2*n+3);
351 Int_t nbPols = n*(n+2);
352
354 nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 2*n*5 + n*n*6);
355
356 if (buff)
357 {
358 SetPoints(buff->fPnts);
359 SetSegsAndPols(*buff);
360 }
361
362 return buff;
363}
364
365////////////////////////////////////////////////////////////////////////////////
366/// Fill TBuffer3D structure for segments and polygons.
367
369{
370 Int_t indx, i, j;
372
374
375 Int_t nn1 = (n+1)*n+1;
376 indx = 0;
377 // Lower end-cap (n radial segments)
378 for (j=0; j<n; j++) {
379 buff.fSegs[indx++] = c+2;
380 buff.fSegs[indx++] = 0;
381 buff.fSegs[indx++] = j+1;
382 }
383 // Sectors (n)
384 for (i=0; i<n+1; i++) {
385 // lateral (circles) segments (n)
386 for (j=0; j<n; j++) {
387 buff.fSegs[indx++] = c;
388 buff.fSegs[indx++] = n*i+1+j;
389 buff.fSegs[indx++] = n*i+1+((j+1)%n);
390 }
391 if (i==n) break; // skip i=n for generators
392 // generator segments (n)
393 for (j=0; j<n; j++) {
394 buff.fSegs[indx++] = c;
395 buff.fSegs[indx++] = n*i+1+j;
396 buff.fSegs[indx++] = n*(i+1)+1+j;
397 }
398 }
399 // Upper end-cap
400 for (j=0; j<n; j++) {
401 buff.fSegs[indx++] = c+1;
402 buff.fSegs[indx++] = n*n+1+j;
403 buff.fSegs[indx++] = nn1;
404 }
405
406 indx = 0;
407
408 // lower end-cap (n polygons)
409 for (j=0; j<n; j++) {
410 buff.fPols[indx++] = c+2;
411 buff.fPols[indx++] = 3;
412 buff.fPols[indx++] = n+j;
413 buff.fPols[indx++] = (j+1)%n;
414 buff.fPols[indx++] = j;
415 }
416 // Sectors (n)
417 for (i=0; i<n; i++) {
418 // lateral faces (n)
419 for (j=0; j<n; j++) {
420 buff.fPols[indx++] = c;
421 buff.fPols[indx++] = 4;
422 buff.fPols[indx++] = (2*i+1)*n+j;
423 buff.fPols[indx++] = 2*(i+1)*n+j;
424 buff.fPols[indx++] = (2*i+3)*n+j;
425 buff.fPols[indx++] = 2*(i+1)*n+((j+1)%n);
426 }
427 }
428 // upper end-cap (n polygons)
429 for (j=0; j<n; j++) {
430 buff.fPols[indx++] = c+1;
431 buff.fPols[indx++] = 3;
432 buff.fPols[indx++] = 2*n*(n+1)+j;
433 buff.fPols[indx++] = 2*n*(n+1)+((j+1)%n);
434 buff.fPols[indx++] = (2*n+1)*n+j;
435 }
436}
437
438////////////////////////////////////////////////////////////////////////////////
439/// Computes the closest distance from given point to this shape.
440
442{
443 Double_t safz = fDz-TMath::Abs(point[2]);
444 if (!in) safz = -safz;
445 Double_t safr;
446 Double_t rsq = point[0]*point[0]+point[1]*point[1];
447 Double_t z0 = fA*rsq+fB;
448 Double_t r0sq = (point[2]-fB)/fA;
449 if (r0sq<0) {
450 if (in) return 0.;
451 return safz;
452 }
453 Double_t dr = TMath::Sqrt(rsq)-TMath::Sqrt(r0sq);
454 if (in) {
455 if (dr>-1.E-8) return 0.;
456 Double_t dz = TMath::Abs(point[2]-z0);
457 safr = -dr*dz/TMath::Sqrt(dr*dr+dz*dz);
458 } else {
459 if (dr<1.E-8) return safz;
460 Double_t talf = -2.*fA*TMath::Sqrt(r0sq);
461 Double_t salf = talf/TMath::Sqrt(1.+talf*talf);
462 safr = TMath::Abs(dr*salf);
463 }
464 if (in) return TMath::Min(safr,safz);
465 return TMath::Max(safr,safz);
466}
467
468////////////////////////////////////////////////////////////////////////////////
469/// Set paraboloid dimensions.
470
472{
473 if ((rlo<0) || (rhi<0) || (dz<=0) || TMath::Abs(rlo-rhi)<TGeoShape::Tolerance()) {
475 Error("SetParaboloidDimensions", "Dimensions of %s invalid: check (rlo>=0) (rhi>=0) (rlo!=rhi) dz>0",GetName());
476 return;
477 }
478 fRlo = rlo;
479 fRhi = rhi;
480 fDz = dz;
481 Double_t dd = 1./(fRhi*fRhi - fRlo*fRlo);
482 fA = 2.*fDz*dd;
483 fB = - fDz * (fRlo*fRlo + fRhi*fRhi)*dd;
484}
485
486////////////////////////////////////////////////////////////////////////////////
487/// Set paraboloid dimensions starting from an array.
488
490{
491 Double_t rlo = param[0];
492 Double_t rhi = param[1];
493 Double_t dz = param[2];
494 SetParaboloidDimensions(rlo, rhi, dz);
495}
496
497////////////////////////////////////////////////////////////////////////////////
498/// Create paraboloid mesh points.
499/// ~~~ {.cpp}
500/// Npoints = n*(n+1) + 2
501/// ifirst = 0
502/// ipoint(i,j) = 1+i*n+j; i=[0,n] j=[0,n-1]
503/// ilast = 1+n*(n+1)
504/// Nsegments = n*(2*n+3)
505/// lower: (0, j+1); j=[0,n-1]
506/// circle(i): (n*i+1+j, n*i+1+(j+1)%n); i=[0,n] j=[0,n-1]
507/// generator(i): (n*i+1+j, n*(i+1)+1+j); i,j=[0,n-1]
508/// upper: (n*n+1+j, (n+1)*n+1) j=[0,n-1]
509/// Npolygons = n*(n+2)
510/// lower: (n+j, (j+1)%n, j) j=[0,n-1]
511/// lateral(i): ((2*i+1)*n+j, 2*(i+1)*n+j, (2*i+3)*n+j, 2*(i+1)*n+(j+1)%n)
512/// i,j = [0,n-1]
513/// upper: ((2n+1)*n+j, 2*n*(n+1)+(j+1)%n, 2*n*(n+1)+j) j=[0,n-1]
514/// ~~~
515
517{
518 if (!points) return;
519 Double_t ttmin, ttmax;
520 ttmin = TMath::ATan2(-fDz, fRlo);
521 ttmax = TMath::ATan2(fDz, fRhi);
523 Double_t dtt = (ttmax-ttmin)/n;
524 Double_t dphi = 360./n;
525 Double_t tt;
526 Double_t r, z, delta;
527 Double_t phi, sph, cph;
528 Int_t indx = 0;
529 // center of the lower endcap:
530 points[indx++] = 0; // x
531 points[indx++] = 0; // y
532 points[indx++] = -fDz;
533 for (Int_t i=0; i<n+1; i++) { // nz planes = n+1
534 if (i==0) {
535 r = fRlo;
536 z = -fDz;
537 } else if (i==n) {
538 r = fRhi;
539 z = fDz;
540 } else {
541 tt = TMath::Tan(ttmin + i*dtt);
542 delta = tt*tt - 4*fA*fB; // should be always positive (a*b<0)
543 r = 0.5*(tt+TMath::Sqrt(delta))/fA;
544 z = r*tt;
545 }
546 for (Int_t j=0; j<n; j++) {
547 phi = j*dphi*TMath::DegToRad();
548 sph=TMath::Sin(phi);
549 cph=TMath::Cos(phi);
550 points[indx++] = r*cph;
551 points[indx++] = r*sph;
552 points[indx++] = z;
553 }
554 }
555 // center of the upper endcap
556 points[indx++] = 0; // x
557 points[indx++] = 0; // y
558 points[indx++] = fDz;
559}
560
561////////////////////////////////////////////////////////////////////////////////
562/// Returns numbers of vertices, segments and polygons composing the shape mesh.
563
564void TGeoParaboloid::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
565{
567 nvert = n*(n+1)+2;
568 nsegs = n*(2*n+3);
569 npols = n*(n+2);
570}
571
572////////////////////////////////////////////////////////////////////////////////
573/// Returns number of vertices on the paraboloid mesh.
574
576{
578 return (n*(n+1)+2);
579}
580
581////////////////////////////////////////////////////////////////////////////////
582/// Save a primitive as a C++ statement(s) on output stream "out".
583
584void TGeoParaboloid::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
585{
587 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
588 out << " rlo = " << fRlo << ";" << std::endl;
589 out << " rhi = " << fRhi << ";" << std::endl;
590 out << " dz = " << fDZ << ";" << std::endl;
591 out << " TGeoShape *" << GetPointerName() << " = new TGeoParaboloid(\"" << GetName() << "\", rlo,rhi,dz);" << std::endl;
593}
594
595////////////////////////////////////////////////////////////////////////////////
596/// Create paraboloid mesh points.
597
599{
600 if (!points) return;
601 Double_t ttmin, ttmax;
602 ttmin = TMath::ATan2(-fDz, fRlo);
603 ttmax = TMath::ATan2(fDz, fRhi);
605 Double_t dtt = (ttmax-ttmin)/n;
606 Double_t dphi = 360./n;
607 Double_t tt;
608 Double_t r, z, delta;
609 Double_t phi, sph, cph;
610 Int_t indx = 0;
611 // center of the lower endcap:
612 points[indx++] = 0; // x
613 points[indx++] = 0; // y
614 points[indx++] = -fDz;
615 for (Int_t i=0; i<n+1; i++) { // nz planes = n+1
616 if (i==0) {
617 r = fRlo;
618 z = -fDz;
619 } else if (i==n) {
620 r = fRhi;
621 z = fDz;
622 } else {
623 tt = TMath::Tan(ttmin + i*dtt);
624 delta = tt*tt - 4*fA*fB; // should be always positive (a*b<0)
625 r = 0.5*(tt+TMath::Sqrt(delta))/fA;
626 z = r*tt;
627 }
628 for (Int_t j=0; j<n; j++) {
629 phi = j*dphi*TMath::DegToRad();
630 sph=TMath::Sin(phi);
631 cph=TMath::Cos(phi);
632 points[indx++] = r*cph;
633 points[indx++] = r*sph;
634 points[indx++] = z;
635 }
636 }
637 // center of the upper endcap
638 points[indx++] = 0; // x
639 points[indx++] = 0; // y
640 points[indx++] = fDz;
641}
642
643////////////////////////////////////////////////////////////////////////////////
645{
646}
647
648////////////////////////////////////////////////////////////////////////////////
649/// Fills a static 3D buffer and returns a reference.
650
651const TBuffer3D & TGeoParaboloid::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
652{
653 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
654 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
655
656 if (reqSections & TBuffer3D::kRawSizes) {
658 Int_t nbPnts = n*(n+1)+2;
659 Int_t nbSegs = n*(2*n+3);
660 Int_t nbPols = n*(n+2);
661 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 2*n*5 + n*n*6)) {
663 }
664 }
665 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
666 SetPoints(buffer.fPnts);
667 if (!buffer.fLocalFrame) {
668 TransformPoints(buffer.fPnts, buffer.NbPnts());
669 }
670 SetSegsAndPols(buffer);
672 }
673
674 return buffer;
675}
676
677////////////////////////////////////////////////////////////////////////////////
678/// Check the inside status for each of the points in the array.
679/// Input: Array of point coordinates + vector size
680/// Output: Array of Booleans for the inside of each point
681
682void TGeoParaboloid::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
683{
684 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
685}
686
687////////////////////////////////////////////////////////////////////////////////
688/// Compute the normal for an array o points so that norm.dot.dir is positive
689/// Input: Arrays of point coordinates and directions + vector size
690/// Output: Array of normal directions
691
692void TGeoParaboloid::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
693{
694 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
695}
696
697////////////////////////////////////////////////////////////////////////////////
698/// Compute distance from array of input points having directions specified by dirs. Store output in dists
699
700void TGeoParaboloid::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
701{
702 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
703}
704
705////////////////////////////////////////////////////////////////////////////////
706/// Compute distance from array of input points having directions specified by dirs. Store output in dists
707
708void TGeoParaboloid::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
709{
710 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
711}
712
713////////////////////////////////////////////////////////////////////////////////
714/// Compute safe distance from each of the points in the input array.
715/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
716/// Output: Safety values
717
718void TGeoParaboloid::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
719{
720 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
721}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:114
UInt_t NbPnts() const
Definition TBuffer3D.h:80
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Int_t * fSegs
Definition TBuffer3D.h:113
Bool_t fLocalFrame
Definition TBuffer3D.h:90
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:112
Box class.
Definition TGeoBBox.h:18
Double_t fDX
Definition TGeoBBox.h:21
virtual void InspectShape() const
Prints shape parameters.
Definition TGeoBBox.cxx:819
Double_t fDY
Definition TGeoBBox.h:22
Double_t fDZ
Definition TGeoBBox.h:23
virtual void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const
Fills the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections...
Int_t GetNsegments() const
Get number of segments approximating circles.
Geometrical transformation package.
Definition TGeoMatrix.h:41
A paraboloid is defined by the revolution surface generated by a parabola and is bounded by two plane...
virtual void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const
in case shape has some negative parameters, these has to be computed in order to fit the mother
virtual Bool_t Contains(const Double_t *point) const
test if point is inside the elliptical tube
Double_t DistToParaboloid(const Double_t *point, const Double_t *dir, Bool_t in) const
Compute distance from a point to the parabola given by: z = a*rsq + b; rsq = x*x+y*y
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
Computes the closest distance from given point to this shape.
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide the paraboloid along one axis.
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const
compute distance from outside point to surface of the paraboloid and safe distance
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
virtual void SetSegsAndPols(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons.
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute closest distance from point px,py to each vertex
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
virtual void SetPoints(Double_t *points) const
Create paraboloid mesh points.
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const
compute distance from inside point to surface of the paraboloid
void SetParaboloidDimensions(Double_t rlo, Double_t rhi, Double_t dz)
Set paraboloid dimensions.
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
virtual TBuffer3D * MakeBuffer3D() const
Creates a TBuffer3D describing this shape.
virtual Int_t GetNmeshVertices() const
Returns number of vertices on the paraboloid mesh.
virtual void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
Check the inside status for each of the points in the array.
virtual void SetDimensions(Double_t *param)
Set paraboloid dimensions starting from an array.
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual void ComputeBBox()
compute bounding box of the tube
virtual void Sizeof3D() const
virtual ~TGeoParaboloid()
destructor
TGeoParaboloid()
Dummy constructor.
virtual void InspectShape() const
print shape parameters
virtual void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
Compute safe distance from each of the points in the input array.
Base abstract class for all shapes.
Definition TGeoShape.h:26
static Double_t Big()
Definition TGeoShape.h:88
Int_t GetBasicColor() const
Get the basic color (0-7).
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
const char * GetPointerName() const
Provide a pointer name containing uid.
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
virtual const char * GetName() const
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:65
@ kGeoParaboloid
Definition TGeoShape.h:62
@ kGeoRunTimeShape
Definition TGeoShape.h:41
static Double_t Tolerance()
Definition TGeoShape.h:91
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:774
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
const Int_t n
Definition legend1.C:16
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
T1 Sign(T1 a, T2 b)
Returns a value with the magnitude of a and the sign of b.
Definition TMathBase.h:175
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:644
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:79
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:660
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:592
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:586
Double_t Tan(Double_t)
Returns the tangent of an angle of x radians.
Definition TMath.h:598
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
auto * tt
Definition textangle.C:16
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345