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