Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoEltu.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Mihaela Gheata 05/06/02
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 TGeoEltu
13\ingroup Tubes
14
15An elliptical tube is defined by the two semi-axes `A` and `B`. It ranges
16from `-dZ` to `+dZ` as all other tubes:
17
18~~~ {.cpp}
19TGeoEltu(Double_t a,Double_t b,Double_t dz);
20~~~
21
22Begin_Macro
23{
24 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
25 new TGeoManager("eltu", "poza6");
26 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
27 TGeoMedium *med = new TGeoMedium("MED",1,mat);
28 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
29 gGeoManager->SetTopVolume(top);
30 TGeoVolume *vol = gGeoManager->MakeEltu("ELTU",med, 30,10,40);
31 top->AddNode(vol,1);
32 gGeoManager->CloseGeometry();
33 gGeoManager->SetNsegments(50);
34 top->Draw();
35 TView *view = gPad->GetView();
36 view->ShowAxis();
37}
38End_Macro
39*/
40
41
42#include <iostream>
43
44#include "TGeoManager.h"
45#include "TGeoVolume.h"
46#include "TGeoEltu.h"
47#include "TBuffer3D.h"
48#include "TBuffer3DTypes.h"
49#include "TMath.h"
50
52
53////////////////////////////////////////////////////////////////////////////////
54/// Dummy constructor
55
57{
59}
60
61////////////////////////////////////////////////////////////////////////////////
62/// Default constructor specifying X and Y semiaxis length
63
65 :TGeoTube(0, 0, 0)
66{
68 SetEltuDimensions(a, b, dz);
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Default constructor specifying X and Y semiaxis length
74
76 :TGeoTube(name,0.,b,dz)
77{
80 SetEltuDimensions(a, b, dz);
82}
83
84////////////////////////////////////////////////////////////////////////////////
85/// Default constructor specifying minimum and maximum radius
86/// param[0] = A
87/// param[1] = B
88/// param[2] = dz
89
91{
93 SetDimensions(param);
95}
96
97////////////////////////////////////////////////////////////////////////////////
98/// destructor
99
101{
102}
103
104////////////////////////////////////////////////////////////////////////////////
105/// Computes capacity of the shape in [length^3]
106
108{
109 Double_t capacity = 2.*TMath::Pi()*fDz*fRmin*fRmax;
110 return capacity;
111}
112
113////////////////////////////////////////////////////////////////////////////////
114/// compute bounding box of the tube
115
117{
118 fDX = fRmin;
119 fDY = fRmax;
120 fDZ = fDz;
121}
122
123////////////////////////////////////////////////////////////////////////////////
124/// Compute normal to closest surface from POINT.
125
126void TGeoEltu::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
127{
128 Double_t a = fRmin;
129 Double_t b = fRmax;
130 Double_t safr = TMath::Abs(TMath::Sqrt(point[0]*point[0]/(a*a)+point[1]*point[1]/(b*b))-1.);
131 safr *= TMath::Min(a,b);
132 Double_t safz = TMath::Abs(fDz-TMath::Abs(point[2]));
133 if (safz<safr) {
134 norm[0] = norm[1] = 0;
135 norm[2] = TMath::Sign(1.,dir[2]);
136 return;
137 }
138 norm[2] = 0.;
139 norm[0] = point[0]*b*b;
140 norm[1] = point[1]*a*a;
141 TMath::Normalize(norm);
142}
143
144////////////////////////////////////////////////////////////////////////////////
145/// test if point is inside the elliptical tube
146
148{
149 if (TMath::Abs(point[2]) > fDz) return kFALSE;
150 Double_t r2 = (point[0]*point[0])/(fRmin*fRmin)+(point[1]*point[1])/(fRmax*fRmax);
151 if (r2>1.) return kFALSE;
152 return kTRUE;
153}
154
155////////////////////////////////////////////////////////////////////////////////
156/// compute closest distance from point px,py to each vertex
157
159{
161 const Int_t numPoints=4*n;
162 return ShapeDistancetoPrimitive(numPoints, px, py);
163}
164
165////////////////////////////////////////////////////////////////////////////////
166/// compute distance from inside point to surface of the tube
167
168Double_t TGeoEltu::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
169{
172 Double_t safz1=fDz-point[2];
173 Double_t safz2=fDz+point[2];
174
175 if (iact<3 && safe) {
176 Double_t x0=TMath::Abs(point[0]);
177 Double_t y0=TMath::Abs(point[1]);
178 Double_t x1=x0;
180 Double_t y2=y0;
182 Double_t d1=(x1-x0)*(x1-x0)+(y1-y0)*(y1-y0);
183 Double_t d2=(x2-x0)*(x2-x0)+(y2-y0)*(y2-y0);
184 Double_t x3,y3;
185
186 Double_t safz = TMath::Min(safz1,safz2);
187 for (Int_t i=0; i<8; i++) {
188 if (fRmax<fRmin) {
189 x3=0.5*(x1+x2);
191 } else {
192 y3=0.5*(y1+y2);
193 x3=TMath::Sqrt((fRmax-y3)*(fRmax+y3))*fRmin/fRmax;
194 }
195 if (d1<d2) {
196 x2=x3;
197 y2=y3;
198 d2=(x2-x0)*(x2-x0)+(y2-y0)*(y2-y0);
199 } else {
200 x1=x3;
201 y1=y3;
202 d1=(x1-x0)*(x1-x0)+(y1-y0)*(y1-y0);
203 }
204 }
205 Double_t safr = TMath::Sqrt(d1)-1.0E-3;
206 *safe = TMath::Min(safz, safr);
207 if (iact==0) return TGeoShape::Big();
208 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
209 }
210 // compute distance to surface
211 // Do Z
212 Double_t snxt = TGeoShape::Big();
213 if (dir[2]>0) {
214 snxt=safz1/dir[2];
215 } else {
216 if (dir[2]<0) snxt=-safz2/dir[2];
217 }
218 Double_t sz = snxt;
219 Double_t xz=point[0]+dir[0]*sz;
220 Double_t yz=point[1]+dir[1]*sz;
221 if ((xz*xz/a2+yz*yz/b2)<=1) return snxt;
222 // do elliptical surface
223 Double_t tolerance = TGeoShape::Tolerance();
224 Double_t u=dir[0]*dir[0]*b2+dir[1]*dir[1]*a2;
225 Double_t v=point[0]*dir[0]*b2+point[1]*dir[1]*a2;
226 Double_t w=point[0]*point[0]*b2+point[1]*point[1]*a2-a2*b2;
227 Double_t d=v*v-u*w;
228 if (d<0 || TGeoShape::IsSameWithinTolerance(u,0)) return tolerance;
230 snxt = (-v+sd)/u;
231
232 if (snxt<0) return tolerance;
233 return snxt;
234}
235
236////////////////////////////////////////////////////////////////////////////////
237/// compute distance from outside point to surface of the tube and safe distance
238
239Double_t TGeoEltu::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
240{
241 Double_t safz=TMath::Abs(point[2])-fDz;
244 if (iact<3 && safe) {
245 Double_t x0=TMath::Abs(point[0]);
246 Double_t y0=TMath::Abs(point[1]);
247 *safe=0.;
248 if ((x0*x0/a2+y0*y0/b2)>=1) {
249 Double_t phi1=0;
250 Double_t phi2=0.5*TMath::Pi();
251 Double_t phi3;
252 Double_t x3=0.,y3=0.,d;
253 for (Int_t i=0; i<10; i++) {
254 phi3=(phi1+phi2)*0.5;
255 x3=fRmin*TMath::Cos(phi3);
256 y3=fRmax*TMath::Sin(phi3);
257 d=y3*a2*(x0-x3)-x3*b2*(y0-y3);
258 if (d<0) phi1=phi3;
259 else phi2=phi3;
260 }
261 *safe=TMath::Sqrt((x0-x3)*(x0-x3)+(y0-y3)*(y0-y3));
262 }
263 if (safz>0) {
264 *safe=TMath::Sqrt((*safe)*(*safe)+safz*safz);
265 }
266 if (iact==0) return TGeoShape::Big();
267 if ((iact==1) && (step<*safe)) return TGeoShape::Big();
268 }
269 // compute vector distance
270 Double_t zi, tau;
271 Double_t epsil = 10.*TGeoShape::Tolerance();
272 if (safz > -epsil) {
273 // point beyond the z limit (up or down)
274 // Check if direction is outgoing
275 if (point[2]*dir[2]>0) return TGeoShape::Big();
276 // Check if direction is perpendicular to Z axis
277 if (TGeoShape::IsSameWithinTolerance(dir[2],0)) return TGeoShape::Big();
278 // select +z or -z depending on the side of the point
279 zi = (point[2] > 0) ? fDz : -fDz;
280 // Distance to zi plane position
281 tau = (zi-point[2])/dir[2];
282 // Extrapolated coordinates at the z position of the end plane.
283 Double_t xz=point[0]+dir[0]*tau;
284 Double_t yz=point[1]+dir[1]*tau;
285 if ((xz*xz/a2+yz*yz/b2)<1) return tau;
286 }
287
288// Check if the bounding box is crossed within the requested distance
289 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
290 if (sdist>=step) return TGeoShape::Big();
291 Double_t u=dir[0]*dir[0]*b2+dir[1]*dir[1]*a2; // positive
293 Double_t v=point[0]*dir[0]*b2+point[1]*dir[1]*a2;
294 Double_t w=point[0]*point[0]*b2+point[1]*point[1]*a2-a2*b2;
295 Double_t d=v*v-u*w;
296 if (d<0) return TGeoShape::Big();
298 // Biggest solution - if negative, or very close to boundary
299 // no crossing (just exiting, no re-entering possible)
300 tau = (-v+dsq)/u;
301 if (tau < epsil) return TGeoShape::Big();
302 // only entering crossing must be considered (smallest)
303 tau = (-v-dsq)/u;
304 zi=point[2]+tau*dir[2];
305 // If the crossing point is not in the Z range, there is no crossing
306 if ((TMath::Abs(zi)-fDz)>0) return TGeoShape::Big();
307 // crossing is backwards (point inside the ellipse) in Z range
308 if (tau < 0) return 0.;
309 // Point is outside and crossing the elliptical tube in Z range
310 return tau;
311}
312
313////////////////////////////////////////////////////////////////////////////////
314/// Divide the shape along one axis.
315
316TGeoVolume *TGeoEltu::Divide(TGeoVolume * /*voldiv*/, const char * /*divname*/, Int_t /*iaxis*/, Int_t /*ndiv*/,
317 Double_t /*start*/, Double_t /*step*/)
318{
319 Error("Divide", "Elliptical tubes divisions not implemented");
320 return 0;
321}
322
323////////////////////////////////////////////////////////////////////////////////
324/// Fill vector param[4] with the bounding cylinder parameters. The order
325/// is the following : Rmin, Rmax, Phi1, Phi2
326
328{
329 param[0] = 0.; // Rmin
330 param[1] = TMath::Max(fRmin, fRmax); // Rmax
331 param[1] *= param[1];
332 param[2] = 0.; // Phi1
333 param[3] = 360.; // Phi2
334}
335
336////////////////////////////////////////////////////////////////////////////////
337/// in case shape has some negative parameters, these has to be computed
338/// in order to fit the mother
339
341{
342 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
343 if (!mother->TestShapeBit(kGeoEltu)) {
344 Error("GetMakeRuntimeShape", "invalid mother");
345 return 0;
346 }
347 Double_t a, b, dz;
348 a = fRmin;
349 b = fRmax;
350 dz = fDz;
351 if (fDz<0) dz=((TGeoEltu*)mother)->GetDz();
352 if (fRmin<0)
353 a = ((TGeoEltu*)mother)->GetA();
354 if (fRmax<0)
355 a = ((TGeoEltu*)mother)->GetB();
356
357 return (new TGeoEltu(a, b, dz));
358}
359
360////////////////////////////////////////////////////////////////////////////////
361/// print shape parameters
362
364{
365 printf("*** Shape %s: TGeoEltu ***\n", GetName());
366 printf(" A = %11.5f\n", fRmin);
367 printf(" B = %11.5f\n", fRmax);
368 printf(" dz = %11.5f\n", fDz);
369 printf(" Bounding box:\n");
371}
372
373////////////////////////////////////////////////////////////////////////////////
374/// computes the closest distance from given point to this shape, according
375/// to option. The matching point on the shape is stored in spoint.
376
377Double_t TGeoEltu::Safety(const Double_t *point, Bool_t /*in*/) const
378{
379 Double_t x0 = TMath::Abs(point[0]);
380 Double_t y0 = TMath::Abs(point[1]);
381 Double_t x1, y1, dx, dy;
382 Double_t safr, safz;
383 Double_t onepls = 1.+TGeoShape::Tolerance();
384 Double_t onemin = 1.-TGeoShape::Tolerance();
385 Double_t sqdist = x0*x0/(fRmin*fRmin)+y0*y0/(fRmax*fRmax);
386 Bool_t in = kTRUE;
387 if (sqdist>onepls) in = kFALSE;
388 else if (sqdist<onemin) in = kTRUE;
389 else return 0.;
390
391 if (in) {
392 x1 = fRmin*TMath::Sqrt(1.-(y0*y0)/(fRmax*fRmax));
393 y1 = fRmax*TMath::Sqrt(1.-(x0*x0)/(fRmin*fRmin));
394 dx = x1-x0;
395 dy = y1-y0;
396 if (TMath::Abs(dx)<TGeoShape::Tolerance()) return 0;
397 safr = dx*dy/TMath::Sqrt(dx*dx+dy*dy);
398 safz = fDz - TMath::Abs(point[2]);
399 return TMath::Min(safr,safz);
400 }
401
403 safr = y0 - fRmax;
404 } else {
406 safr = x0 - fRmin;
407 } else {
409 x1 = f*x0;
410 y1 = f*y0;
411 dx = x0-x1;
412 dy = y0-y1;
413 Double_t ast = fRmin*y1/fRmax;
414 Double_t bct = fRmax*x1/fRmin;
415 Double_t d = TMath::Sqrt(bct*bct+ast*ast);
416 safr = (dx*bct+dy*ast)/d;
417 }
418 }
419 safz = TMath::Abs(point[2])-fDz;
420 return TMath::Max(safr, safz);
421}
422
423////////////////////////////////////////////////////////////////////////////////
424/// Save a primitive as a C++ statement(s) on output stream "out".
425
426void TGeoEltu::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
427{
429 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
430 out << " a = " << fRmin << ";" << std::endl;
431 out << " b = " << fRmax << ";" << std::endl;
432 out << " dz = " << fDz << ";" << std::endl;
433 out << " TGeoShape *" << GetPointerName() << " = new TGeoEltu(\"" << GetName() << "\",a,b,dz);" << std::endl;
435}
436
437////////////////////////////////////////////////////////////////////////////////
438/// Set dimensions of the elliptical tube.
439
441{
442 if ((a<=0) || (b<0) || (dz<0)) {
444 }
445 fRmin=a;
446 fRmax=b;
447 fDz=dz;
448}
449
450////////////////////////////////////////////////////////////////////////////////
451/// Set shape dimensions starting from an array.
452
454{
455 Double_t a = param[0];
456 Double_t b = param[1];
457 Double_t dz = param[2];
458 SetEltuDimensions(a, b, dz);
459}
460
461////////////////////////////////////////////////////////////////////////////////
462/// Create elliptical tube mesh points
463
465{
466 Double_t dz;
467 Int_t j, n;
468
470 Double_t dphi = 360./n;
471 Double_t phi = 0;
472 Double_t cph,sph;
473 dz = fDz;
474
475 Int_t indx = 0;
476 Double_t r2,r;
479
480 if (points) {
481 for (j = 0; j < n; j++) {
482 points[indx+6*n] = points[indx] = 0;
483 indx++;
484 points[indx+6*n] = points[indx] = 0;
485 indx++;
486 points[indx+6*n] = dz;
487 points[indx] =-dz;
488 indx++;
489 }
490 for (j = 0; j < n; j++) {
491 phi = j*dphi*TMath::DegToRad();
492 sph=TMath::Sin(phi);
493 cph=TMath::Cos(phi);
494 r2=(a2*b2)/(b2+(a2-b2)*sph*sph);
495 r=TMath::Sqrt(r2);
496 points[indx+6*n] = points[indx] = r*cph;
497 indx++;
498 points[indx+6*n] = points[indx] = r*sph;
499 indx++;
500 points[indx+6*n]= dz;
501 points[indx] =-dz;
502 indx++;
503 }
504 }
505}
506
507////////////////////////////////////////////////////////////////////////////////
508/// Returns numbers of vertices, segments and polygons composing the shape mesh.
509
510void TGeoEltu::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
511{
512 TGeoTube::GetMeshNumbers(nvert,nsegs,npols);
513}
514
515////////////////////////////////////////////////////////////////////////////////
516/// Returns the number of vertices on the mesh.
517
519{
521}
522
523////////////////////////////////////////////////////////////////////////////////
524/// Create elliptical tube mesh points
525
527{
528 Double_t dz;
529 Int_t j, n;
530
532 Double_t dphi = 360./n;
533 Double_t phi = 0;
534 Double_t cph,sph;
535 dz = fDz;
536
537 Int_t indx = 0;
538 Double_t r2,r;
541
542 if (points) {
543 for (j = 0; j < n; j++) {
544 points[indx+6*n] = points[indx] = 0;
545 indx++;
546 points[indx+6*n] = points[indx] = 0;
547 indx++;
548 points[indx+6*n] = dz;
549 points[indx] =-dz;
550 indx++;
551 }
552 for (j = 0; j < n; j++) {
553 phi = j*dphi*TMath::DegToRad();
554 sph=TMath::Sin(phi);
555 cph=TMath::Cos(phi);
556 r2=(a2*b2)/(b2+(a2-b2)*sph*sph);
557 r=TMath::Sqrt(r2);
558 points[indx+6*n] = points[indx] = r*cph;
559 indx++;
560 points[indx+6*n] = points[indx] = r*sph;
561 indx++;
562 points[indx+6*n]= dz;
563 points[indx] =-dz;
564 indx++;
565 }
566 }
567}
568
569////////////////////////////////////////////////////////////////////////////////
570/// Fills a static 3D buffer and returns a reference.
571
572const TBuffer3D & TGeoEltu::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
573{
574 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
575 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
576
577 if (reqSections & TBuffer3D::kRawSizes) {
579 Int_t nbPnts = 4*n;
580 Int_t nbSegs = 8*n;
581 Int_t nbPols = 4*n;
582 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
584 }
585 }
586 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
587 SetPoints(buffer.fPnts);
588 if (!buffer.fLocalFrame) {
589 TransformPoints(buffer.fPnts, buffer.NbPnts());
590 }
591 SetSegsAndPols(buffer);
593 }
594
595 return buffer;
596}
597
598////////////////////////////////////////////////////////////////////////////////
599/// Check the inside status for each of the points in the array.
600/// Input: Array of point coordinates + vector size
601/// Output: Array of Booleans for the inside of each point
602
603void TGeoEltu::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
604{
605 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
606}
607
608////////////////////////////////////////////////////////////////////////////////
609/// Compute the normal for an array o points so that norm.dot.dir is positive
610/// Input: Arrays of point coordinates and directions + vector size
611/// Output: Array of normal directions
612
613void TGeoEltu::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
614{
615 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
616}
617
618////////////////////////////////////////////////////////////////////////////////
619/// Compute distance from array of input points having directions specified by dirs. Store output in dists
620
621void TGeoEltu::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
622{
623 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
624}
625
626////////////////////////////////////////////////////////////////////////////////
627/// Compute distance from array of input points having directions specified by dirs. Store output in dists
628
629void TGeoEltu::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
630{
631 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
632}
633
634////////////////////////////////////////////////////////////////////////////////
635/// Compute safe distance from each of the points in the input array.
636/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
637/// Output: Safety values
638
639void TGeoEltu::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
640{
641 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
642}
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define a(i)
Definition RSha256.hxx:99
static const double x3[11]
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
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:197
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 x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
Generic 3D primitive description class.
Definition TBuffer3D.h:18
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
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
Double_t fDX
Definition TGeoBBox.h:21
virtual void InspectShape() const
Prints shape parameters.
Definition TGeoBBox.cxx:819
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 box.
Definition TGeoBBox.cxx:458
Double_t fOrigin[3]
Definition TGeoBBox.h:24
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...
An elliptical tube is defined by the two semi-axes A and B.
Definition TGeoEltu.h:18
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
Definition TGeoEltu.cxx:340
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition TGeoEltu.cxx:426
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
Definition TGeoEltu.cxx:510
virtual Bool_t Contains(const Double_t *point) const
test if point is inside the elliptical tube
Definition TGeoEltu.cxx:147
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
Definition TGeoEltu.cxx:126
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoEltu.cxx:327
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
Definition TGeoEltu.cxx:572
virtual ~TGeoEltu()
destructor
Definition TGeoEltu.cxx:100
virtual void SetDimensions(Double_t *param)
Set shape dimensions starting from an array.
Definition TGeoEltu.cxx:453
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.
Definition TGeoEltu.cxx:603
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...
Definition TGeoEltu.cxx:629
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...
Definition TGeoEltu.cxx:621
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
computes the closest distance from given point to this shape, according to option.
Definition TGeoEltu.cxx:377
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 tube
Definition TGeoEltu.cxx:168
TGeoEltu()
Dummy constructor.
Definition TGeoEltu.cxx:56
virtual void InspectShape() const
print shape parameters
Definition TGeoEltu.cxx:363
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute closest distance from point px,py to each vertex
Definition TGeoEltu.cxx:158
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.
Definition TGeoEltu.cxx:639
virtual void ComputeBBox()
compute bounding box of the tube
Definition TGeoEltu.cxx:116
virtual Int_t GetNmeshVertices() const
Returns the number of vertices on the mesh.
Definition TGeoEltu.cxx:518
void SetEltuDimensions(Double_t a, Double_t b, Double_t dz)
Set dimensions of the elliptical tube.
Definition TGeoEltu.cxx:440
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition TGeoEltu.cxx:107
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide the shape along one axis.
Definition TGeoEltu.cxx:316
virtual void SetPoints(Double_t *points) const
Create elliptical tube mesh points.
Definition TGeoEltu.cxx:464
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 tube and safe distance
Definition TGeoEltu.cxx:239
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...
Definition TGeoEltu.cxx:613
Int_t GetNsegments() const
Get number of segments approximating circles.
Geometrical transformation package.
Definition TGeoMatrix.h:41
Base abstract class for all shapes.
Definition TGeoShape.h:26
static Double_t Big()
Definition TGeoShape.h:88
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.
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
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
@ kGeoRunTimeShape
Definition TGeoShape.h:41
static Double_t Tolerance()
Definition TGeoShape.h:91
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:163
Cylindrical tube class.
Definition TGeoTube.h:18
virtual void SetSegsAndPols(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons.
Definition TGeoTube.cxx:691
Double_t fRmin
Definition TGeoTube.h:21
Double_t fDz
Definition TGeoTube.h:23
Double_t fRmax
Definition TGeoTube.h:22
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
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
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
Float_t Normalize(Float_t v[3])
Normalize a vector v in place.
Definition TMath.cxx:518
T1 Sign(T1 a, T2 b)
Returns a value with the magnitude of a and the sign of b.
Definition TMathBase.h:175
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
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123