Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoHype.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Mihaela Gheata 20/11/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#include <iostream>
13
14#include "TGeoManager.h"
15#include "TGeoVolume.h"
16#include "TVirtualGeoPainter.h"
17#include "TGeoHype.h"
18#include "TBuffer3D.h"
19#include "TBuffer3DTypes.h"
20#include "TMath.h"
21
22/** \class TGeoHype
23\ingroup Shapes_classes
24
25A hyperboloid is represented as a solid limited by two planes
26perpendicular to the Z axis (top and bottom planes) and two hyperbolic
27surfaces of revolution about Z axis (inner and outer surfaces). The
28class describing hyperboloids is TGeoHype has 5 input parameters:
29
30~~~ {.cpp}
31TGeoHype(Double_t rin,Double_t stin,Double_t rout,
32Double_t stout,Double_t dz);
33~~~
34
35Begin_Macro
36{
37 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
38 new TGeoManager("hype", "hyperboloid");
39 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
40 TGeoMedium *med = new TGeoMedium("MED",1,mat);
41 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
42 gGeoManager->SetTopVolume(top);
43 TGeoVolume *vol = gGeoManager->MakeHype("HYPE",med,10, 45 ,20,45,40);
44 TGeoHype *hype = (TGeoHype*)vol->GetShape();
45 top->AddNode(vol,1);
46 gGeoManager->CloseGeometry();
47 gGeoManager->SetNsegments(80);
48 top->Draw();
49 TView *view = gPad->GetView();
50 view->ShowAxis();
51}
52End_Macro
53
54The hyperbolic surface equation is taken in the form:
55
56~~~{.cpp}
57r^2 - z^2 * tan(st)^2 = rmin^2
58~~~
59
60- `r,z:` cylindrical coordinates for a point on the surface
61- `st:` stereo angle between the hyperbola asymptotic lines and Z axis
62- `rmin:` minimum distance between hyperbola and Z axis (at `z=0`)
63
64The input parameters for a hyperboloid represent:
65
66- `rin, stin:` minimum radius and stereo angle in degrees for the inner surface
67- `rout, stout:` minimum radius and stereo angle in degrees for the outer surface
68- `dz:` half length in Z (bounding planes positions at `+/-dz`)
69
70The following conditions are mandatory in order to avoid intersections
71between the inner and outer hyperbolic surfaces in the range `+/-dz`:
72
73- `rin < rout`
74- `rout > 0`
75- `rin^2 + dz^2 * tan(stin)^2 > rout^2 + dz^2 * tan(stout)^2`
76
77Particular cases:
78
79- `rin=0, stin0:` the inner surface is conical
80- `stin=0 / stout=0:` cylindrical surface(s)
81
82*/
83
85
86////////////////////////////////////////////////////////////////////////////////
87/// Default constructor
88
90{
92 fStIn = 0.;
93 fStOut = 0.;
94 fTin = 0.;
95 fTinsq = 0.;
96 fTout = 0.;
97 fToutsq = 0.;
98}
99
100////////////////////////////////////////////////////////////////////////////////
101/// Constructor specifying hyperboloid parameters.
102
103TGeoHype::TGeoHype(Double_t rin, Double_t stin, Double_t rout, Double_t stout, Double_t dz) : TGeoTube(rin, rout, dz)
104{
106 SetHypeDimensions(rin, stin, rout, stout, dz);
107 // dz<0 can be used to force dz of hyperboloid fit the container volume
108 if (fDz < 0)
110 ComputeBBox();
111}
112////////////////////////////////////////////////////////////////////////////////
113/// Constructor specifying parameters and name.
114
115TGeoHype::TGeoHype(const char *name, Double_t rin, Double_t stin, Double_t rout, Double_t stout, Double_t dz)
116 : TGeoTube(name, rin, rout, dz)
117{
119 SetHypeDimensions(rin, stin, rout, stout, dz);
120 // dz<0 can be used to force dz of hyperboloid fit the container volume
121 if (fDz < 0)
123 ComputeBBox();
124}
125
126////////////////////////////////////////////////////////////////////////////////
127/// Default constructor specifying a list of parameters
128/// - param[0] = dz
129/// - param[1] = rin
130/// - param[2] = stin
131/// - param[3] = rout
132/// - param[4] = stout
133
134TGeoHype::TGeoHype(Double_t *param) : TGeoTube(param[1], param[3], param[0])
135{
137 SetDimensions(param);
138 // dz<0 can be used to force dz of hyperboloid fit the container volume
139 if (fDz < 0)
141 ComputeBBox();
142}
143
144////////////////////////////////////////////////////////////////////////////////
145/// destructor
146
148
149////////////////////////////////////////////////////////////////////////////////
150/// Computes capacity of the shape in [length^3]
151
153{
154 Double_t capacity = 2. * TMath::Pi() * fDz * (fRmax * fRmax - fRmin * fRmin) +
155 (2. * TMath::Pi() / 3.) * fDz * fDz * fDz * (fToutsq - fTinsq);
156 return capacity;
157}
158
159////////////////////////////////////////////////////////////////////////////////
160/// Compute bounding box of the hyperboloid
161
163{
164 if (fRmin < 0.) {
165 Warning("ComputeBBox", "Shape %s has invalid rmin=%g ! SET TO 0.", GetName(), fRmin);
166 fRmin = 0.;
167 }
168 if ((fRmin > fRmax) || (fRmin * fRmin + fTinsq * fDz * fDz > fRmax * fRmax + fToutsq * fDz * fDz)) {
170 Error("ComputeBBox", "Shape %s hyperbolic surfaces are malformed: rin=%g, stin=%g, rout=%g, stout=%g", GetName(),
172 return;
173 }
174
176 fDZ = fDz;
177}
178
179////////////////////////////////////////////////////////////////////////////////
180/// Compute normal to closest surface from POINT.
181
182void TGeoHype::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
183{
184 Double_t saf[3];
185 Double_t rsq = point[0] * point[0] + point[1] * point[1];
186 Double_t r = TMath::Sqrt(rsq);
187 Double_t rin = (HasInner()) ? (TMath::Sqrt(RadiusHypeSq(point[2], kTRUE))) : 0.;
188 Double_t rout = TMath::Sqrt(RadiusHypeSq(point[2], kFALSE));
189 saf[0] = TMath::Abs(fDz - TMath::Abs(point[2]));
190 saf[1] = (HasInner()) ? TMath::Abs(rin - r) : TGeoShape::Big();
191 saf[2] = TMath::Abs(rout - r);
192 Int_t i = TMath::LocMin(3, saf);
193 if (i == 0 || r < 1.E-10) {
194 norm[0] = norm[1] = 0.;
195 norm[2] = TMath::Sign(1., dir[2]);
196 return;
197 }
198 Double_t t = (i == 1) ? fTinsq : fToutsq;
199 ;
200 t *= -point[2] / r;
201 Double_t ct = TMath::Sqrt(1. / (1. + t * t));
202 Double_t st = t * ct;
203 Double_t phi = TMath::ATan2(point[1], point[0]);
204 Double_t cphi = TMath::Cos(phi);
205 Double_t sphi = TMath::Sin(phi);
206
207 norm[0] = ct * cphi;
208 norm[1] = ct * sphi;
209 norm[2] = st;
210 if (norm[0] * dir[0] + norm[1] * dir[1] + norm[2] * dir[2] < 0) {
211 norm[0] = -norm[0];
212 norm[1] = -norm[1];
213 norm[2] = -norm[2];
214 }
215}
216
217////////////////////////////////////////////////////////////////////////////////
218/// test if point is inside this tube
219
221{
222 if (TMath::Abs(point[2]) > fDz)
223 return kFALSE;
224 Double_t r2 = point[0] * point[0] + point[1] * point[1];
225 Double_t routsq = RadiusHypeSq(point[2], kFALSE);
226 if (r2 > routsq)
227 return kFALSE;
228 if (!HasInner())
229 return kTRUE;
230 Double_t rinsq = RadiusHypeSq(point[2], kTRUE);
231 if (r2 < rinsq)
232 return kFALSE;
233 return kTRUE;
234}
235
236////////////////////////////////////////////////////////////////////////////////
237/// compute closest distance from point px,py to each corner
238
240{
241 Int_t numPoints = GetNmeshVertices();
242 return ShapeDistancetoPrimitive(numPoints, px, py);
243}
244
245////////////////////////////////////////////////////////////////////////////////
246/// Compute distance from inside point to surface of the hyperboloid.
247
249TGeoHype::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
250{
251 if (iact < 3 && safe) {
252 *safe = Safety(point, kTRUE);
253 if (iact == 0)
254 return TGeoShape::Big();
255 if ((iact == 1) && (*safe > step))
256 return TGeoShape::Big();
257 }
258 // compute distance to surface
259 // Do Z
261 if (dir[2] > 0) {
262 sz = (fDz - point[2]) / dir[2];
263 if (sz <= 0.)
264 return 0.;
265 } else {
266 if (dir[2] < 0) {
267 sz = -(fDz + point[2]) / dir[2];
268 if (sz <= 0.)
269 return 0.;
270 }
271 }
272
273 // Do R
274 Double_t srin = TGeoShape::Big();
275 Double_t srout = TGeoShape::Big();
276 Double_t sr;
277 // inner and outer surfaces
278 Double_t s[2];
279 Int_t npos;
280 npos = DistToHype(point, dir, s, kTRUE, kTRUE);
281 if (npos)
282 srin = s[0];
283 npos = DistToHype(point, dir, s, kFALSE, kTRUE);
284 if (npos)
285 srout = s[0];
286 sr = TMath::Min(srin, srout);
287 return TMath::Min(sz, sr);
288}
289
290////////////////////////////////////////////////////////////////////////////////
291/// compute distance from outside point to surface of the hyperboloid.
292
294TGeoHype::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
295{
296 if (iact < 3 && safe) {
297 *safe = Safety(point, kFALSE);
298 if (iact == 0)
299 return TGeoShape::Big();
300 if ((iact == 1) && (step <= *safe))
301 return TGeoShape::Big();
302 }
303 // Check if the bounding box is crossed within the requested distance
304 Double_t sdist = TGeoBBox::DistFromOutside(point, dir, fDX, fDY, fDZ, fOrigin, step);
305 if (sdist >= step)
306 return TGeoShape::Big();
307 // find distance to shape
308 // Do Z
309 Double_t xi, yi, zi;
310 if (TMath::Abs(point[2]) >= fDz) {
311 // We might find Z plane crossing
312 if ((point[2] * dir[2]) < 0) {
313 // Compute distance to Z (always positive)
314 Double_t sz = (TMath::Abs(point[2]) - fDz) / TMath::Abs(dir[2]);
315 // Extrapolate
316 xi = point[0] + sz * dir[0];
317 yi = point[1] + sz * dir[1];
318 Double_t r2 = xi * xi + yi * yi;
319 Double_t rmin2 = RadiusHypeSq(fDz, kTRUE);
320 if (r2 >= rmin2) {
322 if (r2 <= rmax2)
323 return sz;
324 }
325 }
326 }
327 // We do not cross Z planes.
328 Double_t sin = TGeoShape::Big();
329 Double_t sout = TGeoShape::Big();
330 Double_t s[2];
331 Int_t npos;
332 npos = DistToHype(point, dir, s, kTRUE, kFALSE);
333 if (npos) {
334 zi = point[2] + s[0] * dir[2];
335 if (TMath::Abs(zi) <= fDz)
336 sin = s[0];
337 else if (npos == 2) {
338 zi = point[2] + s[1] * dir[2];
339 if (TMath::Abs(zi) <= fDz)
340 sin = s[1];
341 }
342 }
343 npos = DistToHype(point, dir, s, kFALSE, kFALSE);
344 if (npos) {
345 zi = point[2] + s[0] * dir[2];
346 if (TMath::Abs(zi) <= fDz)
347 sout = s[0];
348 else if (npos == 2) {
349 zi = point[2] + s[1] * dir[2];
350 if (TMath::Abs(zi) <= fDz)
351 sout = s[1];
352 }
353 }
354 return TMath::Min(sin, sout);
355}
356
357////////////////////////////////////////////////////////////////////////////////
358/// Compute distance from an arbitrary point to inner/outer surface of hyperboloid.
359/// Returns number of positive solutions. S[2] contains the solutions.
360
361Int_t TGeoHype::DistToHype(const Double_t *point, const Double_t *dir, Double_t *s, Bool_t inner, Bool_t in) const
362{
363 Double_t r0, t0, snext;
364 if (inner) {
365 if (!HasInner())
366 return 0;
367 r0 = fRmin;
368 t0 = fTinsq;
369 } else {
370 r0 = fRmax;
371 t0 = fToutsq;
372 }
373 Double_t a = dir[0] * dir[0] + dir[1] * dir[1] - t0 * dir[2] * dir[2];
374 Double_t b = t0 * point[2] * dir[2] - point[0] * dir[0] - point[1] * dir[1];
375 Double_t c = point[0] * point[0] + point[1] * point[1] - t0 * point[2] * point[2] - r0 * r0;
376
379 return 0;
380 snext = 0.5 * c / b;
381 if (snext < 0.)
382 return 0;
383 s[0] = snext;
384 return 1;
385 }
386
387 Double_t delta = b * b - a * c;
388 Double_t ainv = 1. / a;
389 Int_t npos = 0;
390 if (delta < 0.)
391 return 0;
392 delta = TMath::Sqrt(delta);
393 Double_t sone = TMath::Sign(1., ainv);
394 Int_t i = -1;
395 while (i < 2) {
396 snext = (b + i * sone * delta) * ainv;
397 i += 2;
398 if (snext < 0)
399 continue;
400 if (snext < 1.E-8) {
401 Double_t r = TMath::Sqrt(point[0] * point[0] + point[1] * point[1]);
402 Double_t t = (inner) ? fTinsq : fToutsq;
403 t *= -point[2] / r;
404 Double_t phi = TMath::ATan2(point[1], point[0]);
405 Double_t ndotd = TMath::Cos(phi) * dir[0] + TMath::Sin(phi) * dir[1] + t * dir[2];
406 if (inner)
407 ndotd *= -1;
408 if (in)
409 ndotd *= -1;
410 if (ndotd < 0)
411 s[npos++] = snext;
412 } else
413 s[npos++] = snext;
414 }
415 return npos;
416}
417
418////////////////////////////////////////////////////////////////////////////////
419/// Cannot divide hyperboloids.
420
421TGeoVolume *TGeoHype::Divide(TGeoVolume * /*voldiv*/, const char *divname, Int_t /*iaxis*/, Int_t /*ndiv*/,
422 Double_t /*start*/, Double_t /*step*/)
423{
424 Error("Divide", "Hyperboloids cannot be divided. Division volume %s not created", divname);
425 return nullptr;
426}
427
428////////////////////////////////////////////////////////////////////////////////
429/// Get range of shape for a given axis.
430
432{
433 xlo = 0;
434 xhi = 0;
435 Double_t dx = 0;
436 switch (iaxis) {
437 case 1: // R
438 xlo = fRmin;
440 dx = xhi - xlo;
441 return dx;
442 case 2: // Phi
443 xlo = 0;
444 xhi = 360;
445 dx = 360;
446 return dx;
447 case 3: // Z
448 xlo = -fDz;
449 xhi = fDz;
450 dx = xhi - xlo;
451 return dx;
452 }
453 return dx;
454}
455
456////////////////////////////////////////////////////////////////////////////////
457/// Fill vector param[4] with the bounding cylinder parameters. The order
458/// is the following : Rmin, Rmax, Phi1, Phi2, dZ
459
461{
462 param[0] = fRmin; // Rmin
463 param[0] *= param[0];
464 param[1] = TMath::Sqrt(RadiusHypeSq(fDz, kFALSE)); // Rmax
465 param[1] *= param[1];
466 param[2] = 0.; // Phi1
467 param[3] = 360.; // Phi2
468}
469
470////////////////////////////////////////////////////////////////////////////////
471/// in case shape has some negative parameters, these has to be computed
472/// in order to fit the mother
473
475{
477 return nullptr;
478 Double_t dz = fDz;
479 Double_t zmin, zmax;
480 if (fDz < 0) {
481 mother->GetAxisRange(3, zmin, zmax);
482 if (zmax < 0)
483 return nullptr;
484 dz = zmax;
485 } else {
486 Error("GetMakeRuntimeShape", "Shape %s does not have negative Z range", GetName());
487 return nullptr;
488 }
489 TGeoShape *hype = new TGeoHype(GetName(), dz, fRmax, fStOut, fRmin, fStIn);
490 return hype;
491}
492
493////////////////////////////////////////////////////////////////////////////////
494/// print shape parameters
495
497{
498 printf("*** Shape %s: TGeoHype ***\n", GetName());
499 printf(" Rin = %11.5f\n", fRmin);
500 printf(" sin = %11.5f\n", fStIn);
501 printf(" Rout = %11.5f\n", fRmax);
502 printf(" sout = %11.5f\n", fStOut);
503 printf(" dz = %11.5f\n", fDz);
504
505 printf(" Bounding box:\n");
507}
508
509////////////////////////////////////////////////////////////////////////////////
510/// Creates a TBuffer3D describing *this* shape.
511/// Coordinates are in local reference frame.
512
514{
516 Bool_t hasRmin = HasInner();
517 Int_t nbPnts = (hasRmin) ? (2 * n * n) : (n * n + 2);
518 Int_t nbSegs = (hasRmin) ? (4 * n * n) : (n * (2 * n + 1));
519 Int_t nbPols = (hasRmin) ? (2 * n * n) : (n * (n + 1));
520
521 TBuffer3D *buff =
522 new TBuffer3D(TBuffer3DTypes::kGeneric, nbPnts, 3 * nbPnts, nbSegs, 3 * nbSegs, nbPols, 6 * nbPols);
523 if (buff) {
524 SetPoints(buff->fPnts);
525 SetSegsAndPols(*buff);
526 }
527
528 return buff;
529}
530
531////////////////////////////////////////////////////////////////////////////////
532/// Fill TBuffer3D structure for segments and polygons.
533
535{
537 Int_t i, j, n;
539 Bool_t hasRmin = HasInner();
540 Int_t irin = 0;
541 Int_t irout = (hasRmin) ? (n * n) : 2;
542 // Fill segments
543 // Case hasRmin:
544 // Inner circles: [isin = 0], n (per circle) * n ( circles)
545 // iseg = isin+n*i+j , i = 0, n-1 , j = 0, n-1
546 // seg(i=1,n; j=1,n) = [irin+n*i+j] and [irin+n*i+(j+1)%n]
547 // Inner generators: [isgenin = isin+n*n], n (per circle) *(n-1) (slices)
548 // iseg = isgenin + i*n + j, i=0,n-2, j=0,n-1
549 // seg(i,j) = [irin+n*i+j] and [irin+n*(i+1)+j]
550 // Outer circles: [isout = isgenin+n*(n-1)], n (per circle) * n ( circles)
551 // iseg = isout + i*n + j , iz = 0, n-1 , j = 0, n-1
552 // seg(i=1,n; j=1,n) = [irout+n*i+j] and [irout+n*i+(j+1)%n]
553 // Outer generators: [isgenout = isout+n*n], n (per circle) *(n-1) (slices)
554 // iseg = isgenout + i*n + j, i=0,n-2, j=0,n-1
555 // seg(i,j) = [irout+n*i+j] and [irout+n*(i+1)+j]
556 // Lower cap : [islow = isgenout + n*(n-1)], n radial segments
557 // iseg = islow + j, j=0,n-1
558 // seg(j) = [irin + j] and [irout+j]
559 // Upper cap: [ishi = islow + n], nradial segments
560 // iseg = ishi + j, j=0,n-1
561 // seg[j] = [irin + n*(n-1) + j] and [irout+n*(n-1) + j]
562 //
563 // Case !hasRmin:
564 // Outer circles: [isout=0], same outer circles (n*n)
565 // Outer generators: isgenout = isout + n*n
566 // Lower cap: [islow = isgenout+n*(n-1)], n seg.
567 // iseg = islow + j, j=0,n-1
568 // seg[j] = [irin] and [irout+j]
569 // Upper cap: [ishi = islow +n]
570 // iseg = ishi + j, j=0,n-1
571 // seg[j] = [irin+1] and [irout+n*(n-1) + j]
572
573 Int_t isin = 0;
574 Int_t isgenin = (hasRmin) ? (isin + n * n) : 0;
575 Int_t isout = (hasRmin) ? (isgenin + n * (n - 1)) : 0;
576 Int_t isgenout = isout + n * n;
577 Int_t islo = isgenout + n * (n - 1);
578 Int_t ishi = islo + n;
579
580 Int_t npt = 0;
581 // Fill inner circle segments (n*n)
582 if (hasRmin) {
583 for (i = 0; i < n; i++) {
584 for (j = 0; j < n; j++) {
585 npt = 3 * (isin + n * i + j);
586 buff.fSegs[npt] = c;
587 buff.fSegs[npt + 1] = irin + n * i + j;
588 buff.fSegs[npt + 2] = irin + n * i + ((j + 1) % n);
589 }
590 }
591 // Fill inner generators (n*(n-1))
592 for (i = 0; i < n - 1; i++) {
593 for (j = 0; j < n; j++) {
594 npt = 3 * (isgenin + n * i + j);
595 buff.fSegs[npt] = c;
596 buff.fSegs[npt + 1] = irin + n * i + j;
597 buff.fSegs[npt + 2] = irin + n * (i + 1) + j;
598 }
599 }
600 }
601 // Fill outer circle segments (n*n)
602 for (i = 0; i < n; i++) {
603 for (j = 0; j < n; j++) {
604 npt = 3 * (isout + n * i + j);
605 buff.fSegs[npt] = c;
606 buff.fSegs[npt + 1] = irout + n * i + j;
607 buff.fSegs[npt + 2] = irout + n * i + ((j + 1) % n);
608 }
609 }
610 // Fill outer generators (n*(n-1))
611 for (i = 0; i < n - 1; i++) {
612 for (j = 0; j < n; j++) {
613 npt = 3 * (isgenout + n * i + j);
614 buff.fSegs[npt] = c;
615 buff.fSegs[npt + 1] = irout + n * i + j;
616 buff.fSegs[npt + 2] = irout + n * (i + 1) + j;
617 }
618 }
619 // Fill lower cap (n)
620 for (j = 0; j < n; j++) {
621 npt = 3 * (islo + j);
622 buff.fSegs[npt] = c;
623 buff.fSegs[npt + 1] = irin;
624 if (hasRmin)
625 buff.fSegs[npt + 1] += j;
626 buff.fSegs[npt + 2] = irout + j;
627 }
628 // Fill upper cap (n)
629 for (j = 0; j < n; j++) {
630 npt = 3 * (ishi + j);
631 buff.fSegs[npt] = c;
632 buff.fSegs[npt + 1] = irin + 1;
633 if (hasRmin)
634 buff.fSegs[npt + 1] += n * (n - 1) + j - 1;
635 buff.fSegs[npt + 2] = irout + n * (n - 1) + j;
636 }
637
638 // Fill polygons
639 // Inner polygons: [ipin = 0] (n-1) slices * n (edges)
640 // ipoly = ipin + n*i + j; i=0,n-2 j=0,n-1
641 // poly[i,j] = [isin+n*i+j] [isgenin+i*n+(j+1)%n] [isin+n*(i+1)+j] [isgenin+i*n+j]
642 // Outer polygons: [ipout = ipin+n*(n-1)] also (n-1)*n
643 // ipoly = ipout + n*i + j; i=0,n-2 j=0,n-1
644 // poly[i,j] = [isout+n*i+j] [isgenout+i*n+j] [isout+n*(i+1)+j] [isgenout+i*n+(j+1)%n]
645 // Lower cap: [iplow = ipout+n*(n-1): n polygons
646 // ipoly = iplow + j; j=0,n-1
647 // poly[i=0,j] = [isin+j] [islow+j] [isout+j] [islow+(j+1)%n]
648 // Upper cap: [ipup = iplow+n] : n polygons
649 // ipoly = ipup + j; j=0,n-1
650 // poly[i=n-1, j] = [isin+n*(n-1)+j] [ishi+(j+1)%n] [isout+n*(n-1)+j] [ishi+j]
651 //
652 // Case !hasRmin:
653 // ipin = 0 no inner polygons
654 // ipout = 0 same outer polygons
655 // Lower cap: iplow = ipout+n*(n-1): n polygons with 3 segments
656 // poly[i=0,j] = [isout+j] [islow+(j+1)%n] [islow+j]
657 // Upper cap: ipup = iplow+n;
658 // poly[i=n-1,j] = [isout+n*(n-1)+j] [ishi+j] [ishi+(j+1)%n]
659
660 Int_t ipin = 0;
661 Int_t ipout = (hasRmin) ? (ipin + n * (n - 1)) : 0;
662 Int_t iplo = ipout + n * (n - 1);
663 Int_t ipup = iplo + n;
664 // Inner polygons n*(n-1)
665 if (hasRmin) {
666 for (i = 0; i < n - 1; i++) {
667 for (j = 0; j < n; j++) {
668 npt = 6 * (ipin + n * i + j);
669 buff.fPols[npt] = c;
670 buff.fPols[npt + 1] = 4;
671 buff.fPols[npt + 2] = isin + n * i + j;
672 buff.fPols[npt + 3] = isgenin + i * n + ((j + 1) % n);
673 buff.fPols[npt + 4] = isin + n * (i + 1) + j;
674 buff.fPols[npt + 5] = isgenin + i * n + j;
675 }
676 }
677 }
678 // Outer polygons n*(n-1)
679 for (i = 0; i < n - 1; i++) {
680 for (j = 0; j < n; j++) {
681 npt = 6 * (ipout + n * i + j);
682 buff.fPols[npt] = c;
683 buff.fPols[npt + 1] = 4;
684 buff.fPols[npt + 2] = isout + n * i + j;
685 buff.fPols[npt + 3] = isgenout + i * n + j;
686 buff.fPols[npt + 4] = isout + n * (i + 1) + j;
687 buff.fPols[npt + 5] = isgenout + i * n + ((j + 1) % n);
688 }
689 }
690 // End caps
691 if (hasRmin) {
692 for (j = 0; j < n; j++) {
693 npt = 6 * (iplo + j);
694 buff.fPols[npt] = c + 1;
695 buff.fPols[npt + 1] = 4;
696 buff.fPols[npt + 2] = isin + j;
697 buff.fPols[npt + 3] = islo + j;
698 buff.fPols[npt + 4] = isout + j;
699 buff.fPols[npt + 5] = islo + ((j + 1) % n);
700 }
701 for (j = 0; j < n; j++) {
702 npt = 6 * (ipup + j);
703 buff.fPols[npt] = c + 2;
704 buff.fPols[npt + 1] = 4;
705 buff.fPols[npt + 2] = isin + n * (n - 1) + j;
706 buff.fPols[npt + 3] = ishi + ((j + 1) % n);
707 buff.fPols[npt + 4] = isout + n * (n - 1) + j;
708 buff.fPols[npt + 5] = ishi + j;
709 }
710 } else {
711 for (j = 0; j < n; j++) {
712 npt = 6 * iplo + 5 * j;
713 buff.fPols[npt] = c + 1;
714 buff.fPols[npt + 1] = 3;
715 buff.fPols[npt + 2] = isout + j;
716 buff.fPols[npt + 3] = islo + ((j + 1) % n);
717 buff.fPols[npt + 4] = islo + j;
718 }
719 for (j = 0; j < n; j++) {
720 npt = 6 * iplo + 5 * (n + j);
721 buff.fPols[npt] = c + 2;
722 buff.fPols[npt + 1] = 3;
723 buff.fPols[npt + 2] = isout + n * (n - 1) + j;
724 buff.fPols[npt + 3] = ishi + j;
725 buff.fPols[npt + 4] = ishi + ((j + 1) % n);
726 }
727 }
728}
729
730////////////////////////////////////////////////////////////////////////////////
731/// Compute r^2 = x^2 + y^2 at a given z coordinate, for either inner or outer hyperbolas.
732
734{
735 Double_t r0, tsq;
736 if (inner) {
737 r0 = fRmin;
738 tsq = fTinsq;
739 } else {
740 r0 = fRmax;
741 tsq = fToutsq;
742 }
743 return (r0 * r0 + tsq * z * z);
744}
745
746////////////////////////////////////////////////////////////////////////////////
747/// Compute z^2 at a given r^2, for either inner or outer hyperbolas.
748
750{
751 Double_t r0, tsq;
752 if (inner) {
753 r0 = fRmin;
754 tsq = fTinsq;
755 } else {
756 r0 = fRmax;
757 tsq = fToutsq;
758 }
759 if (TMath::Abs(tsq) < TGeoShape::Tolerance())
760 return TGeoShape::Big();
761 return ((r * r - r0 * r0) / tsq);
762}
763
764////////////////////////////////////////////////////////////////////////////////
765/// computes the closest distance from given point to this shape, according
766/// to option. The matching point on the shape is stored in spoint.
767
769{
770 Double_t safe, safrmin, safrmax;
771 if (in) {
772 safe = fDz - TMath::Abs(point[2]);
773 safrmin = SafetyToHype(point, kTRUE, in);
774 if (safrmin < safe)
775 safe = safrmin;
776 safrmax = SafetyToHype(point, kFALSE, in);
777 if (safrmax < safe)
778 safe = safrmax;
779 } else {
780 safe = -fDz + TMath::Abs(point[2]);
781 safrmin = SafetyToHype(point, kTRUE, in);
782 if (safrmin > safe)
783 safe = safrmin;
784 safrmax = SafetyToHype(point, kFALSE, in);
785 if (safrmax > safe)
786 safe = safrmax;
787 }
788 return safe;
789}
790
791////////////////////////////////////////////////////////////////////////////////
792/// Compute an underestimate of the closest distance from a point to inner or
793/// outer infinite hyperbolas.
794
796{
797 Double_t r, rsq, rhsq, rh, dr, tsq, saf;
798 if (inner && !HasInner())
799 return (in) ? TGeoShape::Big() : -TGeoShape::Big();
800 rsq = point[0] * point[0] + point[1] * point[1];
801 r = TMath::Sqrt(rsq);
802 rhsq = RadiusHypeSq(point[2], inner);
803 rh = TMath::Sqrt(rhsq);
804 dr = r - rh;
805 if (inner) {
806 if (!in && dr > 0)
807 return -TGeoShape::Big();
809 return TMath::Abs(dr);
811 return TMath::Abs(dr / TMath::Sqrt(1. + fTinsq));
812 tsq = fTinsq;
813 } else {
814 if (!in && dr < 0)
815 return -TGeoShape::Big();
817 return TMath::Abs(dr);
818 tsq = fToutsq;
819 }
821 return 0.;
822 // 1. dr<0 => approximate safety with distance to tangent to hyperbola in z = |point[2]|
823 Double_t m;
824 if (dr < 0) {
825 m = rh / (tsq * TMath::Abs(point[2]));
826 saf = -m * dr / TMath::Sqrt(1. + m * m);
827 return saf;
828 }
829 // 2. dr>0 => approximate safety with distance from point to segment P1(r(z0),z0) and P2(r0, z(r0))
830 m = (TMath::Sqrt(ZHypeSq(r, inner)) - TMath::Abs(point[2])) / dr;
831 saf = m * dr / TMath::Sqrt(1. + m * m);
832 return saf;
833}
834
835////////////////////////////////////////////////////////////////////////////////
836/// Save a primitive as a C++ statement(s) on output stream "out".
837
838void TGeoHype::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
839{
841 return;
842 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
843 out << " rin = " << fRmin << ";" << std::endl;
844 out << " stin = " << fStIn << ";" << std::endl;
845 out << " rout = " << fRmax << ";" << std::endl;
846 out << " stout = " << fStOut << ";" << std::endl;
847 out << " dz = " << fDz << ";" << std::endl;
848 out << " TGeoShape *" << GetPointerName() << " = new TGeoHype(\"" << GetName() << "\",rin,stin,rout,stout,dz);"
849 << std::endl;
851}
852
853////////////////////////////////////////////////////////////////////////////////
854/// Set dimensions of the hyperboloid.
855
857{
858 fRmin = rin;
859 fRmax = rout;
860 fDz = dz;
861 fStIn = stin;
862 fStOut = stout;
864 fTinsq = fTin * fTin;
866 fToutsq = fTout * fTout;
867 if ((fRmin == 0) && (fStIn == 0))
869 else
871}
872
873////////////////////////////////////////////////////////////////////////////////
874/// Set dimensions of the hyperboloid starting from an array.
875/// - param[0] = dz
876/// - param[1] = rin
877/// - param[2] = stin
878/// - param[3] = rout
879/// - param[4] = stout
880
882{
883 Double_t dz = param[0];
884 Double_t rin = param[1];
885 Double_t stin = param[2];
886 Double_t rout = param[3];
887 Double_t stout = param[4];
888 SetHypeDimensions(rin, stin, rout, stout, dz);
889}
890
891////////////////////////////////////////////////////////////////////////////////
892/// create tube mesh points
893
895{
896 Double_t z, dz, r;
897 Int_t i, j, n;
898 if (!points)
899 return;
901 Double_t dphi = 360. / n;
902 Double_t phi = 0;
903 dz = 2. * fDz / (n - 1);
904
905 Int_t indx = 0;
906
907 if (HasInner()) {
908 // Inner surface points
909 for (i = 0; i < n; i++) {
910 z = -fDz + i * dz;
912 for (j = 0; j < n; j++) {
913 phi = j * dphi * TMath::DegToRad();
914 points[indx++] = r * TMath::Cos(phi);
915 points[indx++] = r * TMath::Sin(phi);
916 points[indx++] = z;
917 }
918 }
919 } else {
920 points[indx++] = 0.;
921 points[indx++] = 0.;
922 points[indx++] = -fDz;
923 points[indx++] = 0.;
924 points[indx++] = 0.;
925 points[indx++] = fDz;
926 }
927 // Outer surface points
928 for (i = 0; i < n; i++) {
929 z = -fDz + i * dz;
931 for (j = 0; j < n; j++) {
932 phi = j * dphi * TMath::DegToRad();
933 points[indx++] = r * TMath::Cos(phi);
934 points[indx++] = r * TMath::Sin(phi);
935 points[indx++] = z;
936 }
937 }
938}
939
940////////////////////////////////////////////////////////////////////////////////
941/// create tube mesh points
942
944{
945 Double_t z, dz, r;
946 Int_t i, j, n;
947 if (!points)
948 return;
950 Double_t dphi = 360. / n;
951 Double_t phi = 0;
952 dz = 2. * fDz / (n - 1);
953
954 Int_t indx = 0;
955
956 if (HasInner()) {
957 // Inner surface points
958 for (i = 0; i < n; i++) {
959 z = -fDz + i * dz;
961 for (j = 0; j < n; j++) {
962 phi = j * dphi * TMath::DegToRad();
963 points[indx++] = r * TMath::Cos(phi);
964 points[indx++] = r * TMath::Sin(phi);
965 points[indx++] = z;
966 }
967 }
968 } else {
969 points[indx++] = 0.;
970 points[indx++] = 0.;
971 points[indx++] = -fDz;
972 points[indx++] = 0.;
973 points[indx++] = 0.;
974 points[indx++] = fDz;
975 }
976 // Outer surface points
977 for (i = 0; i < n; i++) {
978 z = -fDz + i * dz;
980 for (j = 0; j < n; j++) {
981 phi = j * dphi * TMath::DegToRad();
982 points[indx++] = r * TMath::Cos(phi);
983 points[indx++] = r * TMath::Sin(phi);
984 points[indx++] = z;
985 }
986 }
987}
988
989////////////////////////////////////////////////////////////////////////////////
990/// Returns numbers of vertices, segments and polygons composing the shape mesh.
991
992void TGeoHype::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
993{
995 Bool_t hasRmin = HasInner();
996 nvert = (hasRmin) ? (2 * n * n) : (n * n + 2);
997 nsegs = (hasRmin) ? (4 * n * n) : (n * (2 * n + 1));
998 npols = (hasRmin) ? (2 * n * n) : (n * (n + 1));
999}
1000
1001////////////////////////////////////////////////////////////////////////////////
1002/// Return number of vertices of the mesh representation
1003
1005{
1007 Int_t numPoints = (HasRmin()) ? (2 * n * n) : (n * n + 2);
1008 return numPoints;
1009}
1010
1011////////////////////////////////////////////////////////////////////////////////
1012/// fill size of this 3-D object
1013
1014void TGeoHype::Sizeof3D() const {}
1015
1016////////////////////////////////////////////////////////////////////////////////
1017/// Fills a static 3D buffer and returns a reference.
1018
1019const TBuffer3D &TGeoHype::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
1020{
1021 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
1022
1023 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
1024
1025 if (reqSections & TBuffer3D::kRawSizes) {
1027 Bool_t hasRmin = HasInner();
1028 Int_t nbPnts = (hasRmin) ? (2 * n * n) : (n * n + 2);
1029 Int_t nbSegs = (hasRmin) ? (4 * n * n) : (n * (2 * n + 1));
1030 Int_t nbPols = (hasRmin) ? (2 * n * n) : (n * (n + 1));
1031 if (buffer.SetRawSizes(nbPnts, 3 * nbPnts, nbSegs, 3 * nbSegs, nbPols, 6 * nbPols)) {
1033 }
1034 }
1035 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
1036 SetPoints(buffer.fPnts);
1037 if (!buffer.fLocalFrame) {
1038 TransformPoints(buffer.fPnts, buffer.NbPnts());
1039 }
1040
1041 SetSegsAndPols(buffer);
1043 }
1044
1045 return buffer;
1046}
1047
1048////////////////////////////////////////////////////////////////////////////////
1049/// Check the inside status for each of the points in the array.
1050/// Input: Array of point coordinates + vector size
1051/// Output: Array of Booleans for the inside of each point
1052
1053void TGeoHype::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
1054{
1055 for (Int_t i = 0; i < vecsize; i++)
1056 inside[i] = Contains(&points[3 * i]);
1057}
1058
1059////////////////////////////////////////////////////////////////////////////////
1060/// Compute the normal for an array o points so that norm.dot.dir is positive
1061/// Input: Arrays of point coordinates and directions + vector size
1062/// Output: Array of normal directions
1063
1064void TGeoHype::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
1065{
1066 for (Int_t i = 0; i < vecsize; i++)
1067 ComputeNormal(&points[3 * i], &dirs[3 * i], &norms[3 * i]);
1068}
1069
1070////////////////////////////////////////////////////////////////////////////////
1071/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1072
1073void TGeoHype::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
1074 Double_t *step) const
1075{
1076 for (Int_t i = 0; i < vecsize; i++)
1077 dists[i] = DistFromInside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1078}
1079
1080////////////////////////////////////////////////////////////////////////////////
1081/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1082
1083void TGeoHype::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
1084 Double_t *step) const
1085{
1086 for (Int_t i = 0; i < vecsize; i++)
1087 dists[i] = DistFromOutside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1088}
1089
1090////////////////////////////////////////////////////////////////////////////////
1091/// Compute safe distance from each of the points in the input array.
1092/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1093/// Output: Safety values
1094
1095void TGeoHype::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
1096{
1097 for (Int_t i = 0; i < vecsize; i++)
1098 safe[i] = Safety(&points[3 * i], inside[i]);
1099}
#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
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
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:229
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
#define isin(address, start, length)
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:115
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:114
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
A hyperboloid is represented as a solid limited by two planes perpendicular to the Z axis (top and bo...
Definition TGeoHype.h:17
Int_t GetNmeshVertices() const override
Return number of vertices of the mesh representation.
const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const override
Fills a static 3D buffer and returns a reference.
Double_t SafetyToHype(const Double_t *point, Bool_t inner, Bool_t in) const
Compute an underestimate of the closest distance from a point to inner or outer infinite hyperbolas.
Definition TGeoHype.cxx:795
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 TGeoHype.cxx:474
Double_t ZHypeSq(Double_t r, Bool_t inner) const
Compute z^2 at a given r^2, for either inner or outer hyperbolas.
Definition TGeoHype.cxx:749
~TGeoHype() override
destructor
Definition TGeoHype.cxx:147
void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) override
Compute normal to closest surface from POINT.
Definition TGeoHype.cxx:182
void SetSegsAndPols(TBuffer3D &buff) const override
Fill TBuffer3D structure for segments and polygons.
Definition TGeoHype.cxx:534
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...
Double_t RadiusHypeSq(Double_t z, Bool_t inner) const
Compute r^2 = x^2 + y^2 at a given z coordinate, for either inner or outer hyperbolas.
Definition TGeoHype.cxx:733
TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step) override
Cannot divide hyperboloids.
Definition TGeoHype.cxx:421
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 TGeoHype.cxx:768
void ComputeBBox() override
Compute bounding box of the hyperboloid.
Definition TGeoHype.cxx:162
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save a primitive as a C++ statement(s) on output stream "out".
Definition TGeoHype.cxx:838
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
compute closest distance from point px,py to each corner
Definition TGeoHype.cxx:239
Double_t Capacity() const override
Computes capacity of the shape in [length^3].
Definition TGeoHype.cxx:152
void Sizeof3D() const override
fill size of this 3-D object
Int_t DistToHype(const Double_t *point, const Double_t *dir, Double_t *s, Bool_t inner, Bool_t in) const
Compute distance from an arbitrary point to inner/outer surface of hyperboloid.
Definition TGeoHype.cxx:361
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...
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.
void SetDimensions(Double_t *param) override
Set dimensions of the hyperboloid starting from an array.
Definition TGeoHype.cxx:881
void SetPoints(Double_t *points) const override
create tube mesh points
Definition TGeoHype.cxx:894
TGeoHype()
Default constructor.
Definition TGeoHype.cxx:89
Double_t fStIn
Definition TGeoHype.h:23
Double_t fToutsq
Definition TGeoHype.h:31
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 hyperboloid.
Definition TGeoHype.cxx:249
void GetBoundingCylinder(Double_t *param) const override
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoHype.cxx:460
void SetHypeDimensions(Double_t rin, Double_t stin, Double_t rout, Double_t stout, Double_t dz)
Set dimensions of the hyperboloid.
Definition TGeoHype.cxx:856
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 hyperboloid.
Definition TGeoHype.cxx:294
void InspectShape() const override
print shape parameters
Definition TGeoHype.cxx:496
Double_t fStOut
Definition TGeoHype.h:24
Bool_t Contains(const Double_t *point) const override
test if point is inside this tube
Definition TGeoHype.cxx:220
Double_t fTinsq
Definition TGeoHype.h:30
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 TGeoHype.cxx:992
Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const override
Get range of shape for a given axis.
Definition TGeoHype.cxx:431
TBuffer3D * MakeBuffer3D() const override
Creates a TBuffer3D describing this shape.
Definition TGeoHype.cxx:513
Double_t fTout
Definition TGeoHype.h:29
Double_t fTin
Definition TGeoHype.h:28
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...
Bool_t HasInner() const
Definition TGeoHype.h:74
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.
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
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const =0
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.
const char * GetName() const override
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:64
@ kGeoInvalidShape
Definition TGeoShape.h:41
@ 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
Double_t fRmin
Definition TGeoTube.h:20
Double_t fDz
Definition TGeoTube.h:22
Double_t fRmax
Definition TGeoTube.h:21
Bool_t HasRmin() const
Definition TGeoTube.h:75
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
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
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Definition TMath.h:982
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:646
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
Double_t Tan(Double_t)
Returns the tangent of an angle of x radians.
Definition TMath.h:600
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
TMarker m
Definition textangle.C:8