Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoCone.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 31/01/02
3// TGeoCone::Contains() and DistFromInside() implemented by Mihaela Gheata
4
5/*************************************************************************
6 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
7 * All rights reserved. *
8 * *
9 * For the licensing terms see $ROOTSYS/LICENSE. *
10 * For the list of contributors see $ROOTSYS/README/CREDITS. *
11 *************************************************************************/
12
13/** \class TGeoCone
14\ingroup Cones
15
16The cones are defined by 5 parameters:
17
18~~~{.cpp}
19TGeoCone(Double_t dz,Double_t rmin1,Double_t rmax1,
20Double_t rmin2,Double_t rmax2);
21~~~
22
23 - `rmin1:` internal radius at Z is `-dz`
24 - `rmax1:` external radius at Z is `-dz`
25 - `rmin2:` internal radius at Z is `+dz`
26 - `rmax2:` external radius at Z is `+dz`
27 - `dz:` half length in Z (a cone ranges from `-dz` to +`dz`)
28
29A cone has Z-axis as its symmetry axis.
30
31Begin_Macro
32{
33 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
34 new TGeoManager("cone", "poza4");
35 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
36 TGeoMedium *med = new TGeoMedium("MED",1,mat);
37 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
38 gGeoManager->SetTopVolume(top);
39 TGeoVolume *vol = gGeoManager->MakeCone("CONE",med, 40,10,20,35,45);
40 vol->SetLineWidth(2);
41 top->AddNode(vol,1);
42 gGeoManager->CloseGeometry();
43 gGeoManager->SetNsegments(30);
44 top->Draw();
45 TView *view = gPad->GetView();
46 view->ShowAxis();
47}
48End_Macro
49*/
50
51
52/** \class TGeoConeSeg
53\ingroup Cones
54
55A cone segment is a cone having a range in `phi`. The cone segment class
56derives from **`TGeoCone`**, having two extra parameters: `phi1` and
57`phi2`.
58
59~~~{.cpp}
60TGeoConeSeg(Double_t dz,Double_t rmin1,Double_t rmax1,
61Double_t rmin2,Double_t rmax2,Double_t phi1,Double_t phi2);
62~~~
63
64Parameters `phi1` and `phi2` have the same meaning and convention as for
65tube segments.
66
67Begin_Macro
68{
69 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
70 new TGeoManager("coneseg", "poza5");
71 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
72 TGeoMedium *med = new TGeoMedium("MED",1,mat);
73 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
74 gGeoManager->SetTopVolume(top);
75 TGeoVolume *vol = gGeoManager->MakeCons("CONESEG",med, 40,30,40,10,20,-30,250);
76 top->AddNode(vol,1);
77 gGeoManager->CloseGeometry();
78 gGeoManager->SetNsegments(30);
79 top->Draw();
80 TView *view = gPad->GetView();
81 view->ShowAxis();
82}
83End_Macro
84*/
85
86#include <iostream>
87
88#include "TGeoManager.h"
89#include "TGeoVolume.h"
90#include "TVirtualGeoPainter.h"
91#include "TGeoCone.h"
92#include "TBuffer3D.h"
93#include "TBuffer3DTypes.h"
94#include "TMath.h"
95
97
98////////////////////////////////////////////////////////////////////////////////
99/// Default constructor
100
102{
104 fDz = 0.0;
105 fRmin1 = 0.0;
106 fRmax1 = 0.0;
107 fRmin2 = 0.0;
108 fRmax2 = 0.0;
109}
110
111////////////////////////////////////////////////////////////////////////////////
112/// Default constructor specifying minimum and maximum radius
113
115 Double_t rmin2, Double_t rmax2)
116 :TGeoBBox(0, 0, 0)
117{
119 SetConeDimensions(dz, rmin1, rmax1, rmin2, rmax2);
120 if ((dz<0) || (rmin1<0) || (rmax1<0) || (rmin2<0) || (rmax2<0)) {
122 }
123 else ComputeBBox();
124}
125
126////////////////////////////////////////////////////////////////////////////////
127/// Default constructor specifying minimum and maximum radius
128
129TGeoCone::TGeoCone(const char *name, Double_t dz, Double_t rmin1, Double_t rmax1,
130 Double_t rmin2, Double_t rmax2)
131 :TGeoBBox(name, 0, 0, 0)
132{
134 SetConeDimensions(dz, rmin1, rmax1, rmin2, rmax2);
135 if ((dz<0) || (rmin1<0) || (rmax1<0) || (rmin2<0) || (rmax2<0)) {
137 }
138 else ComputeBBox();
139}
140
141////////////////////////////////////////////////////////////////////////////////
142/// Default constructor specifying minimum and maximum radius
143/// - param[0] = dz
144/// - param[1] = Rmin1
145/// - param[2] = Rmax1
146/// - param[3] = Rmin2
147/// - param[4] = Rmax2
148
150 :TGeoBBox(0, 0, 0)
151{
153 SetDimensions(param);
154 if ((fDz<0) || (fRmin1<0) || (fRmax1<0) || (fRmin2<0) || (fRmax2<0))
156 else ComputeBBox();
157}
158
159////////////////////////////////////////////////////////////////////////////////
160/// Computes capacity of the shape in [length^3]
161
163{
165}
166
167////////////////////////////////////////////////////////////////////////////////
168/// Computes capacity of the shape in [length^3]
169
171{
172 Double_t capacity = (2.*dz*TMath::Pi()/3.)*(rmax1*rmax1+rmax2*rmax2+rmax1*rmax2-
173 rmin1*rmin1-rmin2*rmin2-rmin1*rmin2);
174 return capacity;
175}
176
177////////////////////////////////////////////////////////////////////////////////
178/// destructor
179
181{
182}
183
184////////////////////////////////////////////////////////////////////////////////
185/// compute bounding box of the sphere
186
188{
189 TGeoBBox *box = (TGeoBBox*)this;
191 memset(fOrigin, 0, 3*sizeof(Double_t));
192}
193
194////////////////////////////////////////////////////////////////////////////////
195/// Compute normal to closest surface from POINT.
196
197void TGeoCone::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
198{
199 Double_t safr,safe,phi;
200 memset(norm,0,3*sizeof(Double_t));
201 phi = TMath::ATan2(point[1],point[0]);
202 Double_t cphi = TMath::Cos(phi);
203 Double_t sphi = TMath::Sin(phi);
204 Double_t ro1 = 0.5*(fRmin1+fRmin2);
205 Double_t tg1 = 0.5*(fRmin2-fRmin1)/fDz;
206 Double_t cr1 = 1./TMath::Sqrt(1.+tg1*tg1);
207 Double_t ro2 = 0.5*(fRmax1+fRmax2);
208 Double_t tg2 = 0.5*(fRmax2-fRmax1)/fDz;
209 Double_t cr2 = 1./TMath::Sqrt(1.+tg2*tg2);
210
211 Double_t r=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
212 Double_t rin = tg1*point[2]+ro1;
213 Double_t rout = tg2*point[2]+ro2;
214 safe = TMath::Abs(fDz-TMath::Abs(point[2]));
215 norm[2] = 1;
216
217 safr = (ro1>0)?(TMath::Abs((r-rin)*cr1)):TGeoShape::Big();
218 if (safr<safe) {
219 safe = safr;
220 norm[0] = cr1*cphi;
221 norm[1] = cr1*sphi;
222 norm[2] = -tg1*cr1;
223 }
224 safr = TMath::Abs((rout-r)*cr2);
225 if (safr<safe) {
226 norm[0] = cr2*cphi;
227 norm[1] = cr2*sphi;
228 norm[2] = -tg2*cr2;
229 }
230 if (norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2]<0) {
231 norm[0] = -norm[0];
232 norm[1] = -norm[1];
233 norm[2] = -norm[2];
234 }
235}
236
237////////////////////////////////////////////////////////////////////////////////
238/// Compute normal to closest surface from POINT.
239
240void TGeoCone::ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm,
241 Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
242{
243 Double_t safe,phi;
244 memset(norm,0,3*sizeof(Double_t));
245 phi = TMath::ATan2(point[1],point[0]);
246 Double_t cphi = TMath::Cos(phi);
247 Double_t sphi = TMath::Sin(phi);
248 Double_t ro1 = 0.5*(rmin1+rmin2);
249 Double_t tg1 = 0.5*(rmin2-rmin1)/dz;
250 Double_t cr1 = 1./TMath::Sqrt(1.+tg1*tg1);
251 Double_t ro2 = 0.5*(rmax1+rmax2);
252 Double_t tg2 = 0.5*(rmax2-rmax1)/dz;
253 Double_t cr2 = 1./TMath::Sqrt(1.+tg2*tg2);
254
255 Double_t r=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
256 Double_t rin = tg1*point[2]+ro1;
257 Double_t rout = tg2*point[2]+ro2;
258 safe = (ro1>0)?(TMath::Abs((r-rin)*cr1)):TGeoShape::Big();
259 norm[0] = cr1*cphi;
260 norm[1] = cr1*sphi;
261 norm[2] = -tg1*cr1;
262 if (TMath::Abs((rout-r)*cr2)<safe) {
263 norm[0] = cr2*cphi;
264 norm[1] = cr2*sphi;
265 norm[2] = -tg2*cr2;
266 }
267 if (norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2]<0) {
268 norm[0] = -norm[0];
269 norm[1] = -norm[1];
270 norm[2] = -norm[2];
271 }
272}
273
274////////////////////////////////////////////////////////////////////////////////
275/// test if point is inside this cone
276
278{
279 if (TMath::Abs(point[2]) > fDz) return kFALSE;
280 Double_t r2 = point[0]*point[0]+point[1]*point[1];
281 Double_t rl = 0.5*(fRmin2*(point[2]+fDz)+fRmin1*(fDz-point[2]))/fDz;
282 Double_t rh = 0.5*(fRmax2*(point[2]+fDz)+fRmax1*(fDz-point[2]))/fDz;
283 if ((r2<rl*rl) || (r2>rh*rh)) return kFALSE;
284 return kTRUE;
285}
286
287////////////////////////////////////////////////////////////////////////////////
288/// Compute distance from inside point to surface of the cone (static)
289/// Boundary safe algorithm.
290
292 Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
293{
294 if (dz<=0) return TGeoShape::Big();
295 // compute distance to surface
296 // Do Z
298 if (dir[2]) {
299 sz = (TMath::Sign(dz, dir[2])-point[2])/dir[2];
300 if (sz<=0) return 0.0;
301 }
302 Double_t rsq=point[0]*point[0]+point[1]*point[1];
303 Double_t zinv = 1./dz;
304 Double_t rin = 0.5*(rmin1+rmin2+(rmin2-rmin1)*point[2]*zinv);
305 // Do Rmin
306 Double_t b,delta,zi;
307 if (rin>0) {
308 // Protection in case point is actually outside the cone
309 if (rsq < rin*(rin+TGeoShape::Tolerance())) {
310 Double_t ddotn = point[0]*dir[0]+point[1]*dir[1]+0.5*(rmin1-rmin2)*dir[2]*zinv*TMath::Sqrt(rsq);
311 if (ddotn<=0) return 0.0;
312 } else {
313 TGeoCone::DistToCone(point, dir, dz, rmin1, rmin2, b, delta);
314 if (delta>0) {
315 Double_t sr = -b-delta;
316 if (sr>0) {
317 zi = point[2]+sr*dir[2];
318 if (TMath::Abs(zi)<=dz) return TMath::Min(sz,sr);
319 }
320 sr = -b+delta;
321 if (sr>0) {
322 zi = point[2]+sr*dir[2];
323 if (TMath::Abs(zi)<=dz) return TMath::Min(sz,sr);
324 }
325 }
326 }
327 }
328 // Do Rmax
329 Double_t rout = 0.5*(rmax1+rmax2+(rmax2-rmax1)*point[2]*zinv);
330 if (rsq > rout*(rout-TGeoShape::Tolerance())) {
331 Double_t ddotn = point[0]*dir[0]+point[1]*dir[1]+0.5*(rmax1-rmax2)*dir[2]*zinv*TMath::Sqrt(rsq);
332 if (ddotn>=0) return 0.0;
333 TGeoCone::DistToCone(point, dir, dz, rmax1, rmax2, b, delta);
334 if (delta<0) return 0.0;
335 Double_t sr = -b+delta;
336 if (sr<0) return sz;
337 if (TMath::Abs(-b-delta)>sr) return sz;
338 zi = point[2]+sr*dir[2];
339 if (TMath::Abs(zi)<=dz) return TMath::Min(sz,sr);
340 return sz;
341 }
342 TGeoCone::DistToCone(point, dir, dz, rmax1, rmax2, b, delta);
343 if (delta>0) {
344 Double_t sr = -b-delta;
345 if (sr>0) {
346 zi = point[2]+sr*dir[2];
347 if (TMath::Abs(zi)<=dz) return TMath::Min(sz,sr);
348 }
349 sr = -b+delta;
350 if (sr>TGeoShape::Tolerance()) {
351 zi = point[2]+sr*dir[2];
352 if (TMath::Abs(zi)<=dz) return TMath::Min(sz,sr);
353 }
354 }
355 return sz;
356}
357
358////////////////////////////////////////////////////////////////////////////////
359/// Compute distance from inside point to surface of the cone
360/// Boundary safe algorithm.
361
362Double_t TGeoCone::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
363{
364 if (iact<3 && safe) {
365 *safe = Safety(point, kTRUE);
366 if (iact==0) return TGeoShape::Big();
367 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
368 }
369 // compute distance to surface
370 return TGeoCone::DistFromInsideS(point, dir, fDz, fRmin1, fRmax1, fRmin2, fRmax2);
371}
372
373////////////////////////////////////////////////////////////////////////////////
374/// Compute distance from outside point to surface of the tube
375/// Boundary safe algorithm.
376
378 Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
379{
380 // compute distance to Z planes
381 if (dz<=0) return TGeoShape::Big();
382 Double_t snxt;
383 Double_t xp, yp, zp;
384 Bool_t inz = kTRUE;
385
386 if (point[2]<=-dz) {
387 if (dir[2]<=0) return TGeoShape::Big();
388 snxt = (-dz-point[2])/dir[2];
389 xp = point[0]+snxt*dir[0];
390 yp = point[1]+snxt*dir[1];
391 Double_t r2 = xp*xp+yp*yp;
392 if ((r2>=rmin1*rmin1) && (r2<=rmax1*rmax1)) return snxt;
393 inz = kFALSE;
394 } else {
395 if (point[2]>=dz) {
396 if (dir[2]>=0) return TGeoShape::Big();
397 snxt = (dz-point[2])/dir[2];
398 xp = point[0]+snxt*dir[0];
399 yp = point[1]+snxt*dir[1];
400 Double_t r2 = xp*xp+yp*yp;
401 if ((r2>=rmin2*rmin2) && (r2<=rmax2*rmax2)) return snxt;
402 inz = kFALSE;
403 }
404 }
405
406 Double_t rsq = point[0]*point[0]+point[1]*point[1];
407 Double_t dzinv = 1./dz;
408 Double_t ro1=0.5*(rmin1+rmin2);
409 Bool_t hasrmin = (ro1>0)?kTRUE:kFALSE;
410 Double_t tg1 = 0.;
411 Double_t rin = 0.;
412 Bool_t inrmin = kTRUE; // r>=rmin
413 if (hasrmin) {
414 tg1=0.5*(rmin2-rmin1)*dzinv;
415 rin=ro1+tg1*point[2];
416 if (rin>0 && rsq<rin*(rin-TGeoShape::Tolerance())) inrmin=kFALSE;
417 }
418 Double_t ro2=0.5*(rmax1+rmax2);
419 Double_t tg2=0.5*(rmax2-rmax1)*dzinv;
420 Double_t rout=tg2*point[2]+ro2;
421 Bool_t inrmax = kFALSE;
422 if (rout>0 && rsq<rout*(rout+TGeoShape::Tolerance())) inrmax=kTRUE;
423 Bool_t in = inz & inrmin & inrmax;
424 Double_t b,delta;
425 // If inside cone, we are most likely on a boundary within machine precision.
426 if (in) {
428 Double_t safz = dz-TMath::Abs(point[2]); // positive
429 Double_t safrmin = (hasrmin)?(r-rin):TGeoShape::Big();
430 Double_t safrmax = rout - r;
431 if (safz<=safrmin && safz<=safrmax) {
432 // on Z boundary
433 if (point[2]*dir[2]<0) return 0.0;
434 return TGeoShape::Big();
435 }
436 if (safrmax<safrmin) {
437 // on rmax boundary
438 Double_t ddotn = point[0]*dir[0]+point[1]*dir[1]-tg2*dir[2]*r;
439 if (ddotn<=0) return 0.0;
440 return TGeoShape::Big();
441 }
442 // on rmin boundary
443 Double_t ddotn = point[0]*dir[0]+point[1]*dir[1]-tg1*dir[2]*r;
444 if (ddotn>=0) return 0.0;
445 // we can cross (+) solution of rmin
446 TGeoCone::DistToCone(point, dir, dz, rmin1, rmin2, b, delta);
447
448 if (delta<0) return 0.0;
449 snxt = -b+delta;
450 if (snxt<0) return TGeoShape::Big();
451 if (TMath::Abs(-b-delta)>snxt) return TGeoShape::Big();
452 zp = point[2]+snxt*dir[2];
453 if (TMath::Abs(zp)<=dz) return snxt;
454 return TGeoShape::Big();
455 }
456
457 // compute distance to inner cone
458 snxt = TGeoShape::Big();
459 if (!inrmin) {
460 // ray can cross inner cone (but not only!)
461 TGeoCone::DistToCone(point, dir, dz, rmin1, rmin2, b, delta);
462 if (delta<0) return TGeoShape::Big();
463 snxt = -b+delta;
464 if (snxt>0) {
465 zp = point[2]+snxt*dir[2];
466 if (TMath::Abs(zp)<=dz) return snxt;
467 }
468 snxt = -b-delta;
469 if (snxt>0) {
470 zp = point[2]+snxt*dir[2];
471 if (TMath::Abs(zp)<=dz) return snxt;
472 }
473 snxt = TGeoShape::Big();
474 } else {
475 if (hasrmin) {
476 TGeoCone::DistToCone(point, dir, dz, rmin1, rmin2, b, delta);
477 if (delta>0) {
478 Double_t din = -b+delta;
479 if (din>0) {
480 zp = point[2]+din*dir[2];
481 if (TMath::Abs(zp)<=dz) snxt = din;
482 }
483 }
484 }
485 }
486
487 if (inrmax) return snxt;
488 // We can cross outer cone, both solutions possible
489 // compute distance to outer cone
490 TGeoCone::DistToCone(point, dir, dz, rmax1, rmax2, b, delta);
491 if (delta<0) return snxt;
492 Double_t dout = -b-delta;
493 if (dout>0 && dout<snxt) {
494 zp = point[2]+dout*dir[2];
495 if (TMath::Abs(zp)<=dz) return dout;
496 }
497 dout = -b+delta;
498 if (dout<=0 || dout>snxt) return snxt;
499 zp = point[2]+dout*dir[2];
500 if (TMath::Abs(zp)<=dz) return dout;
501 return snxt;
502}
503
504////////////////////////////////////////////////////////////////////////////////
505/// compute distance from outside point to surface of the tube
506
507Double_t TGeoCone::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
508{
509 // compute safe radius
510 if (iact<3 && safe) {
511 *safe = Safety(point, kFALSE);
512 if (iact==0) return TGeoShape::Big();
513 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
514 }
515 // Check if the bounding box is crossed within the requested distance
516 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
517 if (sdist>=step) return TGeoShape::Big();
518 // compute distance to Z planes
520}
521
522////////////////////////////////////////////////////////////////////////////////
523/// Static method to compute distance to a conical surface with :
524/// - r1, z1 - radius and Z position of lower base
525/// - r2, z2 - radius and Z position of upper base
526
527void TGeoCone::DistToCone(const Double_t *point, const Double_t *dir, Double_t dz, Double_t r1, Double_t r2,
528 Double_t &b, Double_t &delta)
529{
530 b = 0;
531 delta = -1.;
532 if (dz<0) return;
533 Double_t ro0 = 0.5*(r1+r2);
534 Double_t tz = 0.5*(r2-r1)/dz;
535 Double_t rsq = point[0]*point[0] + point[1]*point[1];
536 Double_t rc = ro0 + point[2]*tz;
537
538 Double_t a = dir[0]*dir[0] + dir[1]*dir[1] - tz*tz*dir[2]*dir[2];
539 b = point[0]*dir[0] + point[1]*dir[1] - tz*rc*dir[2];
540 Double_t c = rsq - rc*rc;
541
543 if (TMath::Abs(b)<TGeoShape::Tolerance()) return;
544 b = 0.5*c/b;
545 delta = 0.;
546 return;
547 }
548 a = 1./a;
549 b *= a;
550 c *= a;
551 delta = b*b - c;
552 if (delta>0) {
553 delta = TMath::Sqrt(delta);
554 } else {
555 delta = -1.;
556 }
557}
558
559////////////////////////////////////////////////////////////////////////////////
560/// compute closest distance from point px,py to each corner
561
563{
565 const Int_t numPoints = 4*n;
566 return ShapeDistancetoPrimitive(numPoints, px, py);
567}
568
569////////////////////////////////////////////////////////////////////////////////
570/// Divide this cone shape belonging to volume "voldiv" into ndiv volumes
571/// called divname, from start position with the given step. Returns pointer
572/// to created division cell volume in case of Z divisions. For Z division
573/// creates all volumes with different shapes and returns pointer to volume that
574/// was divided. In case a wrong division axis is supplied, returns pointer to
575/// volume that was divided.
576
577TGeoVolume *TGeoCone::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
578 Double_t start, Double_t step)
579{
580 TGeoShape *shape; //--- shape to be created
581 TGeoVolume *vol; //--- division volume to be created
582 TGeoVolumeMulti *vmulti; //--- generic divided volume
583 TGeoPatternFinder *finder; //--- finder to be attached
584 TString opt = ""; //--- option to be attached
585 Int_t id;
586 Double_t end = start+ndiv*step;
587 switch (iaxis) {
588 case 1: //--- R division
589 Error("Divide","division of a cone on R not implemented");
590 return 0;
591 case 2: // --- Phi division
592 finder = new TGeoPatternCylPhi(voldiv, ndiv, start, end);
593 voldiv->SetFinder(finder);
594 finder->SetDivIndex(voldiv->GetNdaughters());
595 shape = new TGeoConeSeg(fDz, fRmin1, fRmax1, fRmin2, fRmax2, -step/2, step/2);
596 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
597 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
598 vmulti->AddVolume(vol);
599 opt = "Phi";
600 for (id=0; id<ndiv; id++) {
601 voldiv->AddNodeOffset(vol, id, start+id*step+step/2, opt.Data());
602 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
603 }
604 return vmulti;
605 case 3: //--- Z division
606 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
607 finder = new TGeoPatternZ(voldiv, ndiv, start, end);
608 voldiv->SetFinder(finder);
609 finder->SetDivIndex(voldiv->GetNdaughters());
610 for (id=0; id<ndiv; id++) {
611 Double_t z1 = start+id*step;
612 Double_t z2 = start+(id+1)*step;
613 Double_t rmin1n = 0.5*(fRmin1*(fDz-z1)+fRmin2*(fDz+z1))/fDz;
614 Double_t rmax1n = 0.5*(fRmax1*(fDz-z1)+fRmax2*(fDz+z1))/fDz;
615 Double_t rmin2n = 0.5*(fRmin1*(fDz-z2)+fRmin2*(fDz+z2))/fDz;
616 Double_t rmax2n = 0.5*(fRmax1*(fDz-z2)+fRmax2*(fDz+z2))/fDz;
617 shape = new TGeoCone(0.5*step,rmin1n, rmax1n, rmin2n, rmax2n);
618 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
619 vmulti->AddVolume(vol);
620 opt = "Z";
621 voldiv->AddNodeOffset(vol, id, start+id*step+step/2, opt.Data());
622 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
623 }
624 return vmulti;
625 default:
626 Error("Divide", "Wrong axis type for division");
627 return 0;
628 }
629}
630
631////////////////////////////////////////////////////////////////////////////////
632/// Returns name of axis IAXIS.
633
634const char *TGeoCone::GetAxisName(Int_t iaxis) const
635{
636 switch (iaxis) {
637 case 1:
638 return "R";
639 case 2:
640 return "PHI";
641 case 3:
642 return "Z";
643 default:
644 return "undefined";
645 }
646}
647
648////////////////////////////////////////////////////////////////////////////////
649/// Get range of shape for a given axis.
650
652{
653 xlo = 0;
654 xhi = 0;
655 Double_t dx = 0;
656 switch (iaxis) {
657 case 2:
658 xlo = 0.;
659 xhi = 360.;
660 return 360.;
661 case 3:
662 xlo = -fDz;
663 xhi = fDz;
664 dx = xhi-xlo;
665 return dx;
666 }
667 return dx;
668}
669
670////////////////////////////////////////////////////////////////////////////////
671/// Fill vector param[4] with the bounding cylinder parameters. The order
672/// is the following : Rmin, Rmax, Phi1, Phi2, dZ
673
675{
676 param[0] = TMath::Min(fRmin1, fRmin2); // Rmin
677 param[0] *= param[0];
678 param[1] = TMath::Max(fRmax1, fRmax2); // Rmax
679 param[1] *= param[1];
680 param[2] = 0.; // Phi1
681 param[3] = 360.; // Phi2
682}
683
684////////////////////////////////////////////////////////////////////////////////
685/// in case shape has some negative parameters, these has to be computed
686/// in order to fit the mother
687
689{
690 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
691 if (!mother->TestShapeBit(kGeoCone)) {
692 Error("GetMakeRuntimeShape", "invalid mother");
693 return 0;
694 }
695 Double_t rmin1, rmax1, rmin2, rmax2, dz;
696 rmin1 = fRmin1;
697 rmax1 = fRmax1;
698 rmin2 = fRmin2;
699 rmax2 = fRmax2;
700 dz = fDz;
701 if (fDz<0) dz=((TGeoCone*)mother)->GetDz();
702 if (fRmin1<0)
703 rmin1 = ((TGeoCone*)mother)->GetRmin1();
704 if (fRmax1<0)
705 rmax1 = ((TGeoCone*)mother)->GetRmax1();
706 if (fRmin2<0)
707 rmin2 = ((TGeoCone*)mother)->GetRmin2();
708 if (fRmax2<0)
709 rmax2 = ((TGeoCone*)mother)->GetRmax2();
710
711 return (new TGeoCone(GetName(), dz, rmin1, rmax1, rmin2, rmax2));
712}
713
714////////////////////////////////////////////////////////////////////////////////
715/// Fills array with n random points located on the line segments of the shape mesh.
716/// The output array must be provided with a length of minimum 3*npoints. Returns
717/// true if operation is implemented.
718
720{
721 if (npoints > (npoints/2)*2) {
722 Error("GetPointsOnSegments","Npoints must be even number");
723 return kFALSE;
724 }
725 Bool_t hasrmin = (fRmin1>0 || fRmin2>0)?kTRUE:kFALSE;
726 Int_t nc = 0;
727 if (hasrmin) nc = (Int_t)TMath::Sqrt(0.5*npoints);
728 else nc = (Int_t)TMath::Sqrt(1.*npoints);
729 Double_t dphi = TMath::TwoPi()/nc;
730 Double_t phi = 0;
731 Int_t ntop = 0;
732 if (hasrmin) ntop = npoints/2 - nc*(nc-1);
733 else ntop = npoints - nc*(nc-1);
734 Double_t dz = 2*fDz/(nc-1);
735 Double_t z = 0;
736 Int_t icrt = 0;
737 Int_t nphi = nc;
738 Double_t rmin = 0.;
739 Double_t rmax = 0.;
740 // loop z sections
741 for (Int_t i=0; i<nc; i++) {
742 if (i == (nc-1)) nphi = ntop;
743 z = -fDz + i*dz;
744 if (hasrmin) rmin = 0.5*(fRmin1+fRmin2) + 0.5*(fRmin2-fRmin1)*z/fDz;
745 rmax = 0.5*(fRmax1+fRmax2) + 0.5*(fRmax2-fRmax1)*z/fDz;
746 // loop points on circle sections
747 for (Int_t j=0; j<nphi; j++) {
748 phi = j*dphi;
749 if (hasrmin) {
750 array[icrt++] = rmin * TMath::Cos(phi);
751 array[icrt++] = rmin * TMath::Sin(phi);
752 array[icrt++] = z;
753 }
754 array[icrt++] = rmax * TMath::Cos(phi);
755 array[icrt++] = rmax * TMath::Sin(phi);
756 array[icrt++] = z;
757 }
758 }
759 return kTRUE;
760}
761
762
763////////////////////////////////////////////////////////////////////////////////
764/// print shape parameters
765
767{
768 printf("*** Shape %s TGeoCone ***\n", GetName());
769 printf(" dz =: %11.5f\n", fDz);
770 printf(" Rmin1 = %11.5f\n", fRmin1);
771 printf(" Rmax1 = %11.5f\n", fRmax1);
772 printf(" Rmin2 = %11.5f\n", fRmin2);
773 printf(" Rmax2 = %11.5f\n", fRmax2);
774 printf(" Bounding box:\n");
776}
777
778////////////////////////////////////////////////////////////////////////////////
779/// Creates a TBuffer3D describing *this* shape.
780/// Coordinates are in local reference frame.
781
783{
785 Int_t nbPnts = 4*n;
786 Int_t nbSegs = 8*n;
787 Int_t nbPols = 4*n;
789 nbPnts, 3*nbPnts,
790 nbSegs, 3*nbSegs,
791 nbPols, 6*nbPols);
792
793 if (buff)
794 {
795 SetPoints(buff->fPnts);
796 SetSegsAndPols(*buff);
797 }
798 return buff;
799}
800
801////////////////////////////////////////////////////////////////////////////////
802/// Fill TBuffer3D structure for segments and polygons.
803
805{
806 Int_t i,j;
809
810 for (i = 0; i < 4; i++) {
811 for (j = 0; j < n; j++) {
812 buffer.fSegs[(i*n+j)*3 ] = c;
813 buffer.fSegs[(i*n+j)*3+1] = i*n+j;
814 buffer.fSegs[(i*n+j)*3+2] = i*n+j+1;
815 }
816 buffer.fSegs[(i*n+j-1)*3+2] = i*n;
817 }
818 for (i = 4; i < 6; i++) {
819 for (j = 0; j < n; j++) {
820 buffer.fSegs[(i*n+j)*3 ] = c+1;
821 buffer.fSegs[(i*n+j)*3+1] = (i-4)*n+j;
822 buffer.fSegs[(i*n+j)*3+2] = (i-2)*n+j;
823 }
824 }
825 for (i = 6; i < 8; i++) {
826 for (j = 0; j < n; j++) {
827 buffer.fSegs[(i*n+j)*3 ] = c;
828 buffer.fSegs[(i*n+j)*3+1] = 2*(i-6)*n+j;
829 buffer.fSegs[(i*n+j)*3+2] = (2*(i-6)+1)*n+j;
830 }
831 }
832
833 Int_t indx = 0;
834 i=0;
835 for (j = 0; j < n; j++) {
836 indx = 6*(i*n+j);
837 buffer.fPols[indx ] = c;
838 buffer.fPols[indx+1] = 4;
839 buffer.fPols[indx+5] = i*n+j;
840 buffer.fPols[indx+4] = (4+i)*n+j;
841 buffer.fPols[indx+3] = (2+i)*n+j;
842 buffer.fPols[indx+2] = (4+i)*n+j+1;
843 }
844 buffer.fPols[indx+2] = (4+i)*n;
845 i=1;
846 for (j = 0; j < n; j++) {
847 indx = 6*(i*n+j);
848 buffer.fPols[indx ] = c;
849 buffer.fPols[indx+1] = 4;
850 buffer.fPols[indx+2] = i*n+j;
851 buffer.fPols[indx+3] = (4+i)*n+j;
852 buffer.fPols[indx+4] = (2+i)*n+j;
853 buffer.fPols[indx+5] = (4+i)*n+j+1;
854 }
855 buffer.fPols[indx+5] = (4+i)*n;
856 i=2;
857 for (j = 0; j < n; j++) {
858 indx = 6*(i*n+j);
859 buffer.fPols[indx ] = c+i;
860 buffer.fPols[indx+1] = 4;
861 buffer.fPols[indx+2] = (i-2)*2*n+j;
862 buffer.fPols[indx+3] = (4+i)*n+j;
863 buffer.fPols[indx+4] = ((i-2)*2+1)*n+j;
864 buffer.fPols[indx+5] = (4+i)*n+j+1;
865 }
866 buffer.fPols[indx+5] = (4+i)*n;
867 i=3;
868 for (j = 0; j < n; j++) {
869 indx = 6*(i*n+j);
870 buffer.fPols[indx ] = c+i;
871 buffer.fPols[indx+1] = 4;
872 buffer.fPols[indx+5] = (i-2)*2*n+j;
873 buffer.fPols[indx+4] = (4+i)*n+j;
874 buffer.fPols[indx+3] = ((i-2)*2+1)*n+j;
875 buffer.fPols[indx+2] = (4+i)*n+j+1;
876 }
877 buffer.fPols[indx+2] = (4+i)*n;
878}
879
880////////////////////////////////////////////////////////////////////////////////
881/// computes the closest distance from given point to this shape, according
882/// to option. The matching point on the shape is stored in spoint.
883
885{
886 Double_t saf[4];
887 Double_t r=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
888 saf[0] = TGeoShape::SafetySeg(r,point[2], fRmin1, -fDz, fRmax1, -fDz, !in);
889 saf[1] = TGeoShape::SafetySeg(r,point[2], fRmax2, fDz, fRmin2, fDz, !in);
890 saf[2] = TGeoShape::SafetySeg(r,point[2], fRmin2, fDz, fRmin1, -fDz, !in);
891 saf[3] = TGeoShape::SafetySeg(r,point[2], fRmax1, -fDz, fRmax2, fDz, !in);
892 Double_t safety = saf[TMath::LocMin(4,saf)];
893 if (safety>1.E20) safety = 0.;
894 return safety;
895}
896
897////////////////////////////////////////////////////////////////////////////////
898/// computes the closest distance from given point to this shape, according
899/// to option. The matching point on the shape is stored in spoint.
900
902 Double_t rmin2, Double_t rmax2, Int_t skipz)
903{
904 Double_t saf[4];
905 Double_t r=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
906// Double_t rin = tg1*point[2]+ro1;
907// Double_t rout = tg2*point[2]+ro2;
908 switch (skipz) {
909 case 1: // skip lower Z plane
910 saf[0] = TGeoShape::Big();
911 saf[1] = TGeoShape::SafetySeg(r,point[2], rmax2, dz, rmin2, dz, !in);
912 break;
913 case 2: // skip upper Z plane
914 saf[0] = TGeoShape::SafetySeg(r,point[2], rmin1, -dz, rmax1, -dz, !in);
915 saf[1] = TGeoShape::Big();
916 break;
917 case 3: // skip both
918 saf[0] = saf[1] = TGeoShape::Big();
919 break;
920 default:
921 saf[0] = TGeoShape::SafetySeg(r,point[2], rmin1, -dz, rmax1, -dz, !in);
922 saf[1] = TGeoShape::SafetySeg(r,point[2], rmax2, dz, rmin2, dz, !in);
923 }
924 // Safety to inner part
925 if (rmin1>0 || rmin2>0)
926 saf[2] = TGeoShape::SafetySeg(r,point[2], rmin2, dz, rmin1, -dz, !in);
927 else
928 saf[2] = TGeoShape::Big();
929 saf[3] = TGeoShape::SafetySeg(r,point[2], rmax1, -dz, rmax2, dz, !in);
930 return saf[TMath::LocMin(4,saf)];
931}
932
933////////////////////////////////////////////////////////////////////////////////
934/// Save a primitive as a C++ statement(s) on output stream "out".
935
936void TGeoCone::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
937{
939 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
940 out << " dz = " << fDz << ";" << std::endl;
941 out << " rmin1 = " << fRmin1 << ";" << std::endl;
942 out << " rmax1 = " << fRmax1 << ";" << std::endl;
943 out << " rmin2 = " << fRmin2 << ";" << std::endl;
944 out << " rmax2 = " << fRmax2 << ";" << std::endl;
945 out << " TGeoShape *" << GetPointerName() << " = new TGeoCone(\"" << GetName() << "\", dz,rmin1,rmax1,rmin2,rmax2);" << std::endl;
947}
948
949////////////////////////////////////////////////////////////////////////////////
950/// Set cone dimensions.
951
953 Double_t rmin2, Double_t rmax2)
954{
955 if (rmin1>=0) {
956 if (rmax1>0) {
957 if (rmin1<=rmax1) {
958 // normal rmin/rmax
959 fRmin1 = rmin1;
960 fRmax1 = rmax1;
961 } else {
962 fRmin1 = rmax1;
963 fRmax1 = rmin1;
964 Warning("SetConeDimensions", "rmin1>rmax1 Switch rmin1<->rmax1");
966 }
967 } else {
968 // run-time
969 fRmin1 = rmin1;
970 fRmax1 = rmax1;
971 }
972 } else {
973 // run-time
974 fRmin1 = rmin1;
975 fRmax1 = rmax1;
976 }
977 if (rmin2>=0) {
978 if (rmax2>0) {
979 if (rmin2<=rmax2) {
980 // normal rmin/rmax
981 fRmin2 = rmin2;
982 fRmax2 = rmax2;
983 } else {
984 fRmin2 = rmax2;
985 fRmax2 = rmin2;
986 Warning("SetConeDimensions", "rmin2>rmax2 Switch rmin2<->rmax2");
988 }
989 } else {
990 // run-time
991 fRmin2 = rmin2;
992 fRmax2 = rmax2;
993 }
994 } else {
995 // run-time
996 fRmin2 = rmin2;
997 fRmax2 = rmax2;
998 }
999
1000 fDz = dz;
1001}
1002
1003////////////////////////////////////////////////////////////////////////////////
1004/// Set cone dimensions from an array.
1005
1007{
1008 Double_t dz = param[0];
1009 Double_t rmin1 = param[1];
1010 Double_t rmax1 = param[2];
1011 Double_t rmin2 = param[3];
1012 Double_t rmax2 = param[4];
1013 SetConeDimensions(dz, rmin1, rmax1, rmin2, rmax2);
1014}
1015
1016////////////////////////////////////////////////////////////////////////////////
1017/// Create cone mesh points.
1018
1020{
1021 Double_t dz, phi, dphi;
1022 Int_t j, n;
1023
1025 dphi = 360./n;
1026 dz = fDz;
1027 Int_t indx = 0;
1028
1029 if (points) {
1030 for (j = 0; j < n; j++) {
1031 phi = j*dphi*TMath::DegToRad();
1032 points[indx++] = fRmin1 * TMath::Cos(phi);
1033 points[indx++] = fRmin1 * TMath::Sin(phi);
1034 points[indx++] = -dz;
1035 }
1036
1037 for (j = 0; j < n; j++) {
1038 phi = j*dphi*TMath::DegToRad();
1039 points[indx++] = fRmax1 * TMath::Cos(phi);
1040 points[indx++] = fRmax1 * TMath::Sin(phi);
1041 points[indx++] = -dz;
1042 }
1043
1044 for (j = 0; j < n; j++) {
1045 phi = j*dphi*TMath::DegToRad();
1046 points[indx++] = fRmin2 * TMath::Cos(phi);
1047 points[indx++] = fRmin2 * TMath::Sin(phi);
1048 points[indx++] = dz;
1049 }
1050
1051 for (j = 0; j < n; j++) {
1052 phi = j*dphi*TMath::DegToRad();
1053 points[indx++] = fRmax2 * TMath::Cos(phi);
1054 points[indx++] = fRmax2 * TMath::Sin(phi);
1055 points[indx++] = dz;
1056 }
1057 }
1058}
1059
1060////////////////////////////////////////////////////////////////////////////////
1061/// Create cone mesh points.
1062
1064{
1065 Double_t dz, phi, dphi;
1066 Int_t j, n;
1067
1069 dphi = 360./n;
1070 dz = fDz;
1071 Int_t indx = 0;
1072
1073 if (points) {
1074 for (j = 0; j < n; j++) {
1075 phi = j*dphi*TMath::DegToRad();
1076 points[indx++] = fRmin1 * TMath::Cos(phi);
1077 points[indx++] = fRmin1 * TMath::Sin(phi);
1078 points[indx++] = -dz;
1079 }
1080
1081 for (j = 0; j < n; j++) {
1082 phi = j*dphi*TMath::DegToRad();
1083 points[indx++] = fRmax1 * TMath::Cos(phi);
1084 points[indx++] = fRmax1 * TMath::Sin(phi);
1085 points[indx++] = -dz;
1086 }
1087
1088 for (j = 0; j < n; j++) {
1089 phi = j*dphi*TMath::DegToRad();
1090 points[indx++] = fRmin2 * TMath::Cos(phi);
1091 points[indx++] = fRmin2 * TMath::Sin(phi);
1092 points[indx++] = dz;
1093 }
1094
1095 for (j = 0; j < n; j++) {
1096 phi = j*dphi*TMath::DegToRad();
1097 points[indx++] = fRmax2 * TMath::Cos(phi);
1098 points[indx++] = fRmax2 * TMath::Sin(phi);
1099 points[indx++] = dz;
1100 }
1101 }
1102}
1103
1104////////////////////////////////////////////////////////////////////////////////
1105/// Returns numbers of vertices, segments and polygons composing the shape mesh.
1106
1107void TGeoCone::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
1108{
1110 nvert = n*4;
1111 nsegs = n*8;
1112 npols = n*4;
1113}
1114
1115////////////////////////////////////////////////////////////////////////////////
1116/// Return number of vertices of the mesh representation
1117
1119{
1121 Int_t numPoints = n*4;
1122 return numPoints;
1123}
1124
1125////////////////////////////////////////////////////////////////////////////////
1126/// Fill size of this 3-D object
1127
1129{
1130}
1131
1132////////////////////////////////////////////////////////////////////////////////
1133/// Fills a static 3D buffer and returns a reference.
1134
1135const TBuffer3D & TGeoCone::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
1136{
1137 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
1138
1139 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
1140
1141 if (reqSections & TBuffer3D::kRawSizes) {
1143 Int_t nbPnts = 4*n;
1144 Int_t nbSegs = 8*n;
1145 Int_t nbPols = 4*n;
1146 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
1148 }
1149 }
1150
1151 // TODO: Can we push this as common down to TGeoShape?
1152 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
1153 SetPoints(buffer.fPnts);
1154 if (!buffer.fLocalFrame) {
1155 TransformPoints(buffer.fPnts, buffer.NbPnts());
1156 }
1157
1158 SetSegsAndPols(buffer);
1160 }
1161
1162 return buffer;
1163}
1164
1165////////////////////////////////////////////////////////////////////////////////
1166/// Check the inside status for each of the points in the array.
1167/// Input: Array of point coordinates + vector size
1168/// Output: Array of Booleans for the inside of each point
1169
1170void TGeoCone::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
1171{
1172 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
1173}
1174
1175////////////////////////////////////////////////////////////////////////////////
1176/// Compute the normal for an array o points so that norm.dot.dir is positive
1177/// Input: Arrays of point coordinates and directions + vector size
1178/// Output: Array of normal directions
1179
1180void TGeoCone::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
1181{
1182 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
1183}
1184
1185////////////////////////////////////////////////////////////////////////////////
1186/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1187
1188void TGeoCone::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1189{
1190 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
1191}
1192
1193////////////////////////////////////////////////////////////////////////////////
1194/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1195
1196void TGeoCone::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1197{
1198 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
1199}
1200
1201////////////////////////////////////////////////////////////////////////////////
1202/// Compute safe distance from each of the points in the input array.
1203/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1204/// Output: Safety values
1205
1206void TGeoCone::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
1207{
1208 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
1209}
1210
1212
1213////////////////////////////////////////////////////////////////////////////////
1214/// Default constructor
1215
1217 :TGeoCone(),
1218 fPhi1(0.), fPhi2(0.), fS1(0.), fC1(0.), fS2(0.), fC2(0.), fSm(0.), fCm(0.), fCdfi(0.)
1219{
1221 fPhi1 = fPhi2 = 0.0;
1222}
1223
1224////////////////////////////////////////////////////////////////////////////////
1225/// Default constructor specifying minimum and maximum radius
1226
1228 Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2)
1229 :TGeoCone(dz, rmin1, rmax1, rmin2, rmax2),
1230 fPhi1(0.), fPhi2(0.), fS1(0.), fC1(0.), fS2(0.), fC2(0.), fSm(0.), fCm(0.), fCdfi(0.)
1231
1232{
1234 SetConsDimensions(dz, rmin1, rmax1, rmin2, rmax2, phi1, phi2);
1235 ComputeBBox();
1236}
1237
1238////////////////////////////////////////////////////////////////////////////////
1239/// Default constructor specifying minimum and maximum radius
1240
1242 Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2)
1243 :TGeoCone(name, dz, rmin1, rmax1, rmin2, rmax2),
1244 fPhi1(0.), fPhi2(0.), fS1(0.), fC1(0.), fS2(0.), fC2(0.), fSm(0.), fCm(0.), fCdfi(0.)
1245{
1247 SetConsDimensions(dz, rmin1, rmax1, rmin2, rmax2, phi1, phi2);
1248 ComputeBBox();
1249}
1250
1251////////////////////////////////////////////////////////////////////////////////
1252/// Default constructor specifying minimum and maximum radius
1253/// - param[0] = dz
1254/// - param[1] = Rmin1
1255/// - param[2] = Rmax1
1256/// - param[3] = Rmin2
1257/// - param[4] = Rmax2
1258/// - param[5] = phi1
1259/// - param[6] = phi2
1260
1262 :TGeoCone(0,0,0,0,0),
1263 fPhi1(0.), fPhi2(0.), fS1(0.), fC1(0.), fS2(0.), fC2(0.), fSm(0.), fCm(0.), fCdfi(0.)
1264{
1266 SetDimensions(param);
1267 ComputeBBox();
1268}
1269
1270////////////////////////////////////////////////////////////////////////////////
1271/// destructor
1272
1274{
1275}
1276
1277////////////////////////////////////////////////////////////////////////////////
1278/// Function called after streaming an object of this class.
1279
1281{
1283}
1284
1285////////////////////////////////////////////////////////////////////////////////
1286/// Init frequently used trigonometric values
1287
1289{
1292 fC1 = TMath::Cos(phi1);
1293 fS1 = TMath::Sin(phi1);
1294 fC2 = TMath::Cos(phi2);
1295 fS2 = TMath::Sin(phi2);
1296 Double_t fio = 0.5*(phi1+phi2);
1297 fCm = TMath::Cos(fio);
1298 fSm = TMath::Sin(fio);
1299 Double_t dfi = 0.5*(phi2-phi1);
1300 fCdfi = TMath::Cos(dfi);
1301}
1302
1303////////////////////////////////////////////////////////////////////////////////
1304/// Computes capacity of the shape in [length^3]
1305
1307{
1309}
1310
1311////////////////////////////////////////////////////////////////////////////////
1312/// Computes capacity of the shape in [length^3]
1313
1315{
1316 Double_t capacity = (TMath::Abs(phi2-phi1)*TMath::DegToRad()*dz/3.)*
1317 (rmax1*rmax1+rmax2*rmax2+rmax1*rmax2-
1318 rmin1*rmin1-rmin2*rmin2-rmin1*rmin2);
1319 return capacity;
1320}
1321
1322////////////////////////////////////////////////////////////////////////////////
1323/// compute bounding box of the tube segment
1324
1326{
1327 Double_t rmin, rmax;
1328 rmin = TMath::Min(fRmin1, fRmin2);
1329 rmax = TMath::Max(fRmax1, fRmax2);
1330
1331 Double_t xc[4];
1332 Double_t yc[4];
1333 xc[0] = rmax*fC1;
1334 yc[0] = rmax*fS1;
1335 xc[1] = rmax*fC2;
1336 yc[1] = rmax*fS2;
1337 xc[2] = rmin*fC1;
1338 yc[2] = rmin*fS1;
1339 xc[3] = rmin*fC2;
1340 yc[3] = rmin*fS2;
1341
1342 Double_t xmin = xc[TMath::LocMin(4, &xc[0])];
1343 Double_t xmax = xc[TMath::LocMax(4, &xc[0])];
1344 Double_t ymin = yc[TMath::LocMin(4, &yc[0])];
1345 Double_t ymax = yc[TMath::LocMax(4, &yc[0])];
1346
1347 Double_t dp = fPhi2-fPhi1;
1348 Double_t ddp = -fPhi1;
1349 if (ddp<0) ddp+= 360;
1350 if (ddp<=dp) xmax = rmax;
1351 ddp = 90-fPhi1;
1352 if (ddp<0) ddp+= 360;
1353 if (ddp<=dp) ymax = rmax;
1354 ddp = 180-fPhi1;
1355 if (ddp<0) ddp+= 360;
1356 if (ddp<=dp) xmin = -rmax;
1357 ddp = 270-fPhi1;
1358 if (ddp<0) ddp+= 360;
1359 if (ddp<=dp) ymin = -rmax;
1360 fOrigin[0] = (xmax+xmin)/2;
1361 fOrigin[1] = (ymax+ymin)/2;
1362 fOrigin[2] = 0;
1363 fDX = (xmax-xmin)/2;
1364 fDY = (ymax-ymin)/2;
1365 fDZ = fDz;
1366}
1367
1368////////////////////////////////////////////////////////////////////////////////
1369/// Compute normal to closest surface from POINT.
1370
1371void TGeoConeSeg::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
1372{
1373 Double_t saf[3];
1374 Double_t ro1 = 0.5*(fRmin1+fRmin2);
1375 Double_t tg1 = 0.5*(fRmin2-fRmin1)/fDz;
1376 Double_t cr1 = 1./TMath::Sqrt(1.+tg1*tg1);
1377 Double_t ro2 = 0.5*(fRmax1+fRmax2);
1378 Double_t tg2 = 0.5*(fRmax2-fRmax1)/fDz;
1379 Double_t cr2 = 1./TMath::Sqrt(1.+tg2*tg2);
1380
1381 Double_t r=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
1382 Double_t rin = tg1*point[2]+ro1;
1383 Double_t rout = tg2*point[2]+ro2;
1384 saf[0] = TMath::Abs(fDz-TMath::Abs(point[2]));
1385 saf[1] = (ro1>0)?(TMath::Abs((r-rin)*cr1)):TGeoShape::Big();
1386 saf[2] = TMath::Abs((rout-r)*cr2);
1387 Int_t i = TMath::LocMin(3,saf);
1388 if (((fPhi2-fPhi1)<360.) && TGeoShape::IsCloseToPhi(saf[i], point,fC1,fS1,fC2,fS2)) {
1389 TGeoShape::NormalPhi(point,dir,norm,fC1,fS1,fC2,fS2);
1390 return;
1391 }
1392 if (i==0) {
1393 norm[0] = norm[1] = 0.;
1394 norm[2] = TMath::Sign(1.,dir[2]);
1395 return;
1396 }
1397
1398 Double_t phi = TMath::ATan2(point[1],point[0]);
1399 Double_t cphi = TMath::Cos(phi);
1400 Double_t sphi = TMath::Sin(phi);
1401
1402 if (i==1) {
1403 norm[0] = cr1*cphi;
1404 norm[1] = cr1*sphi;
1405 norm[2] = -tg1*cr1;
1406 } else {
1407 norm[0] = cr2*cphi;
1408 norm[1] = cr2*sphi;
1409 norm[2] = -tg2*cr2;
1410 }
1411
1412 if (norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2]<0) {
1413 norm[0] = -norm[0];
1414 norm[1] = -norm[1];
1415 norm[2] = -norm[2];
1416 }
1417}
1418
1419////////////////////////////////////////////////////////////////////////////////
1420/// Compute normal to closest surface from POINT.
1421
1422void TGeoConeSeg::ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm,
1423 Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2,
1425{
1426 Double_t saf[2];
1427 Double_t ro1 = 0.5*(rmin1+rmin2);
1428 Double_t tg1 = 0.5*(rmin2-rmin1)/dz;
1429 Double_t cr1 = 1./TMath::Sqrt(1.+tg1*tg1);
1430 Double_t ro2 = 0.5*(rmax1+rmax2);
1431 Double_t tg2 = 0.5*(rmax2-rmax1)/dz;
1432 Double_t cr2 = 1./TMath::Sqrt(1.+tg2*tg2);
1433
1434 Double_t r=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
1435 Double_t rin = tg1*point[2]+ro1;
1436 Double_t rout = tg2*point[2]+ro2;
1437 saf[0] = (ro1>0)?(TMath::Abs((r-rin)*cr1)):TGeoShape::Big();
1438 saf[1] = TMath::Abs((rout-r)*cr2);
1439 Int_t i = TMath::LocMin(2,saf);
1440 if (TGeoShape::IsCloseToPhi(saf[i], point,c1,s1,c2,s2)) {
1441 TGeoShape::NormalPhi(point,dir,norm,c1,s1,c2,s2);
1442 return;
1443 }
1444
1445 Double_t phi = TMath::ATan2(point[1],point[0]);
1446 Double_t cphi = TMath::Cos(phi);
1447 Double_t sphi = TMath::Sin(phi);
1448
1449 if (i==0) {
1450 norm[0] = cr1*cphi;
1451 norm[1] = cr1*sphi;
1452 norm[2] = -tg1*cr1;
1453 } else {
1454 norm[0] = cr2*cphi;
1455 norm[1] = cr2*sphi;
1456 norm[2] = -tg2*cr2;
1457 }
1458
1459 if (norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2]<0) {
1460 norm[0] = -norm[0];
1461 norm[1] = -norm[1];
1462 norm[2] = -norm[2];
1463 }
1464}
1465
1466////////////////////////////////////////////////////////////////////////////////
1467/// test if point is inside this sphere
1468
1470{
1471 if (!TGeoCone::Contains(point)) return kFALSE;
1472 Double_t dphi = fPhi2 - fPhi1;
1473 if (dphi >= 360.) return kTRUE;
1474 Double_t phi = TMath::ATan2(point[1], point[0]) * TMath::RadToDeg();
1475 if (phi < 0 ) phi+=360.;
1476 Double_t ddp = phi-fPhi1;
1477 if (ddp < 0) ddp+=360.;
1478// if (ddp > 360) ddp-=360;
1479 if (ddp > dphi) return kFALSE;
1480 return kTRUE;
1481}
1482
1483////////////////////////////////////////////////////////////////////////////////
1484/// Static method to compute distance to a conical surface with :
1485/// - r1, z1 - radius and Z position of lower base
1486/// - r2, z2 - radius and Z position of upper base
1487/// - phi1, phi2 - phi limits
1488
1490{
1491 Double_t dz = z2-z1;
1492 if (dz<=0) {
1493 return TGeoShape::Big();
1494 }
1495
1496 Double_t dphi = phi2 - phi1;
1497 Bool_t hasphi = kTRUE;
1498 if (dphi >= 360.) hasphi=kFALSE;
1499 if (dphi < 0) dphi+=360.;
1500// printf("phi1=%f phi2=%f dphi=%f\n", phi1, phi2, dphi);
1501
1502 Double_t ro0 = 0.5*(r1+r2);
1503 Double_t fz = (r2-r1)/dz;
1504 Double_t r0sq = point[0]*point[0] + point[1]*point[1];
1505 Double_t rc = ro0 + fz*(point[2]-0.5*(z1+z2));
1506
1507 Double_t a = dir[0]*dir[0] + dir[1]*dir[1] - fz*fz*dir[2]*dir[2];
1508 Double_t b = point[0]*dir[0] + point[1]*dir[1] - fz*rc*dir[2];
1509 Double_t c = r0sq - rc*rc;
1510
1511 if (a==0) return TGeoShape::Big();
1512 a = 1./a;
1513 b *= a;
1514 c *= a;
1515 Double_t delta = b*b - c;
1516 if (delta<0) return TGeoShape::Big();
1517 delta = TMath::Sqrt(delta);
1518
1519 Double_t snxt = -b-delta;
1520 Double_t ptnew[3];
1521 Double_t ddp, phi;
1522 if (snxt>0) {
1523 // check Z range
1524 ptnew[2] = point[2] + snxt*dir[2];
1525 if (((ptnew[2]-z1)*(ptnew[2]-z2)) < 0) {
1526 // check phi range
1527 if (!hasphi) return snxt;
1528 ptnew[0] = point[0] + snxt*dir[0];
1529 ptnew[1] = point[1] + snxt*dir[1];
1530 phi = TMath::ATan2(ptnew[1], ptnew[0]) * TMath::RadToDeg();
1531 if (phi < 0 ) phi+=360.;
1532 ddp = phi-phi1;
1533 if (ddp < 0) ddp+=360.;
1534 // printf("snxt1=%f phi=%f ddp=%f\n", snxt, phi, ddp);
1535 if (ddp<=dphi) return snxt;
1536 }
1537 }
1538 snxt = -b+delta;
1539 if (snxt>0) {
1540 // check Z range
1541 ptnew[2] = point[2] + snxt*dir[2];
1542 if (((ptnew[2]-z1)*(ptnew[2]-z2)) < 0) {
1543 // check phi range
1544 if (!hasphi) return snxt;
1545 ptnew[0] = point[0] + snxt*dir[0];
1546 ptnew[1] = point[1] + snxt*dir[1];
1547 phi = TMath::ATan2(ptnew[1], ptnew[0]) * TMath::RadToDeg();
1548 if (phi < 0 ) phi+=360.;
1549 ddp = phi-phi1;
1550 if (ddp < 0) ddp+=360.;
1551 // printf("snxt2=%f phi=%f ddp=%f\n", snxt, phi, ddp);
1552 if (ddp<=dphi) return snxt;
1553 }
1554 }
1555 return TGeoShape::Big();
1556}
1557
1558////////////////////////////////////////////////////////////////////////////////
1559/// compute distance from inside point to surface of the tube segment
1560
1562 Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2,
1564{
1565 if (dz<=0) return TGeoShape::Big();
1566 // Do Z
1567 Double_t scone = TGeoCone::DistFromInsideS(point,dir,dz,rmin1,rmax1,rmin2,rmax2);
1568 if (scone<=0) return 0.0;
1569 Double_t sfmin = TGeoShape::Big();
1570 Double_t rsq = point[0]*point[0]+point[1]*point[1];
1571 Double_t r = TMath::Sqrt(rsq);
1572 Double_t cpsi=point[0]*cm+point[1]*sm;
1573 if (cpsi>r*cdfi+TGeoShape::Tolerance()) {
1574 sfmin = TGeoShape::DistToPhiMin(point, dir, s1, c1, s2, c2, sm, cm);
1575 return TMath::Min(scone,sfmin);
1576 }
1577 // Point on the phi boundary or outside
1578 // which one: phi1 or phi2
1579 Double_t ddotn, xi, yi;
1580 if (TMath::Abs(point[1]-s1*r) < TMath::Abs(point[1]-s2*r)) {
1581 ddotn = s1*dir[0]-c1*dir[1];
1582 if (ddotn>=0) return 0.0;
1583 ddotn = -s2*dir[0]+c2*dir[1];
1584 if (ddotn<=0) return scone;
1585 sfmin = s2*point[0]-c2*point[1];
1586 if (sfmin<=0) return scone;
1587 sfmin /= ddotn;
1588 if (sfmin >= scone) return scone;
1589 xi = point[0]+sfmin*dir[0];
1590 yi = point[1]+sfmin*dir[1];
1591 if (yi*cm-xi*sm<0) return scone;
1592 return sfmin;
1593 }
1594 ddotn = -s2*dir[0]+c2*dir[1];
1595 if (ddotn>=0) return 0.0;
1596 ddotn = s1*dir[0]-c1*dir[1];
1597 if (ddotn<=0) return scone;
1598 sfmin = -s1*point[0]+c1*point[1];
1599 if (sfmin<=0) return scone;
1600 sfmin /= ddotn;
1601 if (sfmin >= scone) return scone;
1602 xi = point[0]+sfmin*dir[0];
1603 yi = point[1]+sfmin*dir[1];
1604 if (yi*cm-xi*sm>0) return scone;
1605 return sfmin;
1606}
1607
1608////////////////////////////////////////////////////////////////////////////////
1609/// compute distance from inside point to surface of the tube segment
1610
1611Double_t TGeoConeSeg::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
1612{
1613 if (iact<3 && safe) {
1615 if (iact==0) return TGeoShape::Big();
1616 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
1617 }
1618 if ((fPhi2-fPhi1)>=360.) return TGeoCone::DistFromInsideS(point,dir,fDz,fRmin1,fRmax1,fRmin2,fRmax2);
1619
1620 // compute distance to surface
1622}
1623
1624////////////////////////////////////////////////////////////////////////////////
1625/// compute distance from outside point to surface of arbitrary tube
1626
1628 Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2,
1630{
1631 if (dz<=0) return TGeoShape::Big();
1632 Double_t r2, cpsi;
1633 // check Z planes
1634 Double_t xi, yi, zi;
1635 Double_t b,delta;
1636 zi = dz - TMath::Abs(point[2]);
1637 Double_t rin,rout;
1639 Double_t snxt=TGeoShape::Big();
1640 Bool_t in = kFALSE;
1641 Bool_t inz = (zi<0)?kFALSE:kTRUE;
1642 if (!inz) {
1643 if (point[2]*dir[2]>=0) return TGeoShape::Big();
1644 s = -zi/TMath::Abs(dir[2]);
1645 xi = point[0]+s*dir[0];
1646 yi = point[1]+s*dir[1];
1647 r2=xi*xi+yi*yi;
1648 if (dir[2]>0) {
1649 rin = rmin1;
1650 rout = rmax1;
1651 } else {
1652 rin = rmin2;
1653 rout = rmax2;
1654 }
1655 if ((rin*rin<=r2) && (r2<=rout*rout)) {
1656 cpsi=xi*cm+yi*sm;
1657 if (cpsi>=(cdfi*TMath::Sqrt(r2))) return s;
1658 }
1659 }
1660 Double_t zinv = 1./dz;
1661 Double_t rsq = point[0]*point[0]+point[1]*point[1];
1662 Double_t r = TMath::Sqrt(rsq);
1663 Double_t ro1=0.5*(rmin1+rmin2);
1664 Bool_t hasrmin = (ro1>0)?kTRUE:kFALSE;
1665 Double_t tg1 = 0.0;
1666 Bool_t inrmin = kFALSE;
1667 rin = 0.0;
1668 if (hasrmin) {
1669 tg1=0.5*(rmin2-rmin1)*zinv;
1670 rin = ro1+tg1*point[2];
1671 if (rsq > rin*(rin-TGeoShape::Tolerance())) inrmin=kTRUE;
1672 } else {
1673 inrmin = kTRUE;
1674 }
1675 Double_t ro2=0.5*(rmax1+rmax2);
1676 Double_t tg2=0.5*(rmax2-rmax1)*zinv;
1677 rout = ro2+tg2*point[2];
1678 Bool_t inrmax = kFALSE;
1679 if (r < rout+TGeoShape::Tolerance()) inrmax = kTRUE;
1680 Bool_t inphi = kFALSE;
1681 cpsi=point[0]*cm+point[1]*sm;
1682 if (cpsi>r*cdfi-TGeoShape::Tolerance()) inphi = kTRUE;
1683 in = inz & inrmin & inrmax & inphi;
1684 // If inside, we are most likely on a boundary within machine precision.
1685 if (in) {
1686 Double_t safphi = (cpsi-r*cdfi)*TMath::Sqrt(1.-cdfi*cdfi);
1687 Double_t safrmin = (hasrmin)?TMath::Abs(r-rin):(TGeoShape::Big());
1688 Double_t safrmax = TMath::Abs(r-rout);
1689 // check if on Z boundaries
1690 if (zi<safrmax && zi<safrmin && zi<safphi) {
1691 if (point[2]*dir[2]<0) return 0.0;
1692 return TGeoShape::Big();
1693 }
1694 // check if on Rmax boundary
1695 if (safrmax<safrmin && safrmax<safphi) {
1696 Double_t ddotn = point[0]*dir[0]+point[1]*dir[1]-tg2*dir[2]*r;
1697 if (ddotn<=0) return 0.0;
1698 return TGeoShape::Big();
1699 }
1700 // check if on phi boundary
1701 if (safphi<safrmin) {
1702 // We may cross again a phi of rmin boundary
1703 // check first if we are on phi1 or phi2
1704 Double_t un;
1705 if (point[0]*c1 + point[1]*s1 > point[0]*c2 + point[1]*s2) {
1706 un = dir[0]*s1-dir[1]*c1;
1707 if (un < 0) return 0.0;
1708 if (cdfi>=0) return TGeoShape::Big();
1709 un = -dir[0]*s2+dir[1]*c2;
1710 if (un<0) {
1711 s = -point[0]*s2+point[1]*c2;
1712 if (s>0) {
1713 s /= (-un);
1714 zi = point[2]+s*dir[2];
1715 if (TMath::Abs(zi)<=dz) {
1716 xi = point[0]+s*dir[0];
1717 yi = point[1]+s*dir[1];
1718 if ((yi*cm-xi*sm)>0) {
1719 r2=xi*xi+yi*yi;
1720 rin = ro1+tg1*zi;
1721 rout = ro2+tg2*zi;
1722 if ((rin*rin<=r2) && (rout*rout>=r2)) return s;
1723 }
1724 }
1725 }
1726 }
1727 } else {
1728 un = -dir[0]*s2+dir[1]*c2;
1729 if (un < 0) return 0.0;
1730 if (cdfi>=0) return TGeoShape::Big();
1731 un = dir[0]*s1-dir[1]*c1;
1732 if (un<0) {
1733 s = point[0]*s1-point[1]*c1;
1734 if (s>0) {
1735 s /= (-un);
1736 zi = point[2]+s*dir[2];
1737 if (TMath::Abs(zi)<=dz) {
1738 xi = point[0]+s*dir[0];
1739 yi = point[1]+s*dir[1];
1740 if ((yi*cm-xi*sm)<0) {
1741 r2=xi*xi+yi*yi;
1742 rin = ro1+tg1*zi;
1743 rout = ro2+tg2*zi;
1744 if ((rin*rin<=r2) && (rout*rout>=r2)) return s;
1745 }
1746 }
1747 }
1748 }
1749 }
1750 // We may also cross rmin, second solution coming from outside
1751 Double_t ddotn = point[0]*dir[0]+point[1]*dir[1]-tg1*dir[2]*r;
1752 if (ddotn>=0) return TGeoShape::Big();
1753 if (cdfi>=0) return TGeoShape::Big();
1754 TGeoCone::DistToCone(point, dir, dz, rmin1, rmin2, b, delta);
1755 if (delta<0) return TGeoShape::Big();
1756 snxt = -b-delta;
1757 if (snxt<0) return TGeoShape::Big();
1758 snxt = -b+delta;
1759 zi = point[2]+snxt*dir[2];
1760 if (TMath::Abs(zi)>dz) return TGeoShape::Big();
1761 xi = point[0]+snxt*dir[0];
1762 yi = point[1]+snxt*dir[1];
1763 r2=xi*xi+yi*yi;
1764 cpsi=xi*cm+yi*sm;
1765 if (cpsi>=(cdfi*TMath::Sqrt(r2))) return snxt;
1766 return TGeoShape::Big();
1767 }
1768 // We are on rmin boundary: we may cross again rmin or a phi facette
1769 Double_t ddotn = point[0]*dir[0]+point[1]*dir[1]-tg1*dir[2]*r;
1770 if (ddotn>=0) return 0.0;
1771 TGeoCone::DistToCone(point, dir, dz, rmin1, rmin2, b, delta);
1772 if (delta<0) return 0.0;
1773 snxt = -b+delta;
1774 if (snxt<0) return TGeoShape::Big();
1775 if (TMath::Abs(-b-delta)>snxt) return TGeoShape::Big();
1776 zi = point[2]+snxt*dir[2];
1777 if (TMath::Abs(zi)>dz) return TGeoShape::Big();
1778 // OK, we cross rmin at snxt - check if within phi range
1779 xi = point[0]+snxt*dir[0];
1780 yi = point[1]+snxt*dir[1];
1781 r2=xi*xi+yi*yi;
1782 cpsi=xi*cm+yi*sm;
1783 if (cpsi>=(cdfi*TMath::Sqrt(r2))) return snxt;
1784 // we cross rmin in the phi gap - we may cross a phi facette
1785 if (cdfi>=0) return TGeoShape::Big();
1786 Double_t un=-dir[0]*s1+dir[1]*c1;
1787 if (un > 0) {
1788 s=point[0]*s1-point[1]*c1;
1789 if (s>=0) {
1790 s /= un;
1791 zi=point[2]+s*dir[2];
1792 if (TMath::Abs(zi)<=dz) {
1793 xi=point[0]+s*dir[0];
1794 yi=point[1]+s*dir[1];
1795 if ((yi*cm-xi*sm)<=0) {
1796 r2=xi*xi+yi*yi;
1797 rin = ro1+tg1*zi;
1798 rout = ro2+tg2*zi;
1799 if ((rin*rin<=r2) && (rout*rout>=r2)) return s;
1800 }
1801 }
1802 }
1803 }
1804 un=dir[0]*s2-dir[1]*c2;
1805 if (un > 0) {
1806 s=(point[1]*c2-point[0]*s2)/un;
1807 if (s>=0) {
1808 zi=point[2]+s*dir[2];
1809 if (TMath::Abs(zi)<=dz) {
1810 xi=point[0]+s*dir[0];
1811 yi=point[1]+s*dir[1];
1812 if ((yi*cm-xi*sm)>=0) {
1813 r2=xi*xi+yi*yi;
1814 rin = ro1+tg1*zi;
1815 rout = ro2+tg2*zi;
1816 if ((rin*rin<=r2) && (rout*rout>=r2)) return s;
1817 }
1818 }
1819 }
1820 }
1821 return TGeoShape::Big();
1822 }
1823
1824 // The point is really outside
1825 Double_t sr1 = TGeoShape::Big();
1826 if (!inrmax) {
1827 // check crossing with outer cone
1828 TGeoCone::DistToCone(point, dir, dz, rmax1, rmax2, b, delta);
1829 if (delta>=0) {
1830 s = -b-delta;
1831 if (s>0) {
1832 zi=point[2]+s*dir[2];
1833 if (TMath::Abs(zi)<=dz) {
1834 xi=point[0]+s*dir[0];
1835 yi=point[1]+s*dir[1];
1836 r2=xi*xi+yi*yi;
1837 cpsi=xi*cm+yi*sm;
1838 if (cpsi>=(cdfi*TMath::Sqrt(r2))) return s; // rmax crossing
1839 }
1840 }
1841 s = -b+delta;
1842 if (s>0) {
1843 zi=point[2]+s*dir[2];
1844 if (TMath::Abs(zi)<=dz) {
1845 xi=point[0]+s*dir[0];
1846 yi=point[1]+s*dir[1];
1847 r2=xi*xi+yi*yi;
1848 cpsi=xi*cm+yi*sm;
1849 if (cpsi>=(cdfi*TMath::Sqrt(r2))) sr1=s;
1850 }
1851 }
1852 }
1853 }
1854 // check crossing with inner cone
1855 Double_t sr2 = TGeoShape::Big();
1856 TGeoCone::DistToCone(point, dir, dz, rmin1, rmin2, b, delta);
1857 if (delta>=0) {
1858 s = -b-delta;
1859 if (s>0) {
1860 zi=point[2]+s*dir[2];
1861 if (TMath::Abs(zi)<=dz) {
1862 xi=point[0]+s*dir[0];
1863 yi=point[1]+s*dir[1];
1864 r2=xi*xi+yi*yi;
1865 cpsi=xi*cm+yi*sm;
1866 if (cpsi>=(cdfi*TMath::Sqrt(r2))) sr2=s;
1867 }
1868 }
1869 if (sr2>1E10) {
1870 s = -b+delta;
1871 if (s>0) {
1872 zi=point[2]+s*dir[2];
1873 if (TMath::Abs(zi)<=dz) {
1874 xi=point[0]+s*dir[0];
1875 yi=point[1]+s*dir[1];
1876 r2=xi*xi+yi*yi;
1877 cpsi=xi*cm+yi*sm;
1878 if (cpsi>=(cdfi*TMath::Sqrt(r2))) sr2=s;
1879 }
1880 }
1881 }
1882 }
1883 snxt = TMath::Min(sr1,sr2);
1884 // Check phi crossing
1885 s = TGeoShape::DistToPhiMin(point,dir,s1,c1,s2,c2,sm,cm,kFALSE);
1886 if (s>snxt) return snxt;
1887 zi=point[2]+s*dir[2];
1888 if (TMath::Abs(zi)>dz) return snxt;
1889 xi=point[0]+s*dir[0];
1890 yi=point[1]+s*dir[1];
1891 r2=xi*xi+yi*yi;
1892 rout = ro2+tg2*zi;
1893 if (r2>rout*rout) return snxt;
1894 rin = ro1+tg1*zi;
1895 if (r2>=rin*rin) return s; // phi crossing
1896 return snxt;
1897}
1898
1899////////////////////////////////////////////////////////////////////////////////
1900/// compute distance from outside point to surface of the tube
1901
1902Double_t TGeoConeSeg::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
1903{
1904 // compute safe radius
1905 if (iact<3 && safe) {
1906 *safe = Safety(point, kFALSE);
1907 if (iact==0) return TGeoShape::Big();
1908 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
1909 }
1910 // Check if the bounding box is crossed within the requested distance
1911 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
1912 if (sdist>=step) return TGeoShape::Big();
1913 if ((fPhi2-fPhi1)>=360.) return TGeoCone::DistFromOutsideS(point,dir,fDz,fRmin1,fRmax1,fRmin2,fRmax2);
1915}
1916
1917////////////////////////////////////////////////////////////////////////////////
1918/// compute closest distance from point px,py to each corner
1919
1921{
1923 const Int_t numPoints = 4*n;
1924 return ShapeDistancetoPrimitive(numPoints, px, py);
1925}
1926
1927////////////////////////////////////////////////////////////////////////////////
1928/// Divide this cone segment shape belonging to volume "voldiv" into ndiv volumes
1929/// called divname, from start position with the given step. Returns pointer
1930/// to created division cell volume in case of Z divisions. For Z division
1931/// creates all volumes with different shapes and returns pointer to volume that
1932/// was divided. In case a wrong division axis is supplied, returns pointer to
1933/// volume that was divided.
1934
1935TGeoVolume *TGeoConeSeg::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
1936 Double_t start, Double_t step)
1937{
1938 TGeoShape *shape; //--- shape to be created
1939 TGeoVolume *vol; //--- division volume to be created
1940 TGeoVolumeMulti *vmulti; //--- generic divided volume
1941 TGeoPatternFinder *finder; //--- finder to be attached
1942 TString opt = ""; //--- option to be attached
1943 Double_t dphi;
1944 Int_t id;
1945 Double_t end = start+ndiv*step;
1946 switch (iaxis) {
1947 case 1: //--- R division
1948 Error("Divide","division of a cone segment on R not implemented");
1949 return 0;
1950 case 2: //--- Phi division
1951 dphi = fPhi2-fPhi1;
1952 if (dphi<0) dphi+=360.;
1953 finder = new TGeoPatternCylPhi(voldiv, ndiv, start, end);
1954 voldiv->SetFinder(finder);
1955 finder->SetDivIndex(voldiv->GetNdaughters());
1956 shape = new TGeoConeSeg(fDz, fRmin1, fRmax1, fRmin2, fRmax2, -step/2, step/2);
1957 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
1958 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
1959 vmulti->AddVolume(vol);
1960 opt = "Phi";
1961 for (id=0; id<ndiv; id++) {
1962 voldiv->AddNodeOffset(vol, id, start+id*step+step/2, opt.Data());
1963 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
1964 }
1965 return vmulti;
1966 case 3: //--- Z division
1967 finder = new TGeoPatternZ(voldiv, ndiv, start, end);
1968 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
1969 voldiv->SetFinder(finder);
1970 finder->SetDivIndex(voldiv->GetNdaughters());
1971 for (id=0; id<ndiv; id++) {
1972 Double_t z1 = start+id*step;
1973 Double_t z2 = start+(id+1)*step;
1974 Double_t rmin1n = 0.5*(fRmin1*(fDz-z1)+fRmin2*(fDz+z1))/fDz;
1975 Double_t rmax1n = 0.5*(fRmax1*(fDz-z1)+fRmax2*(fDz+z1))/fDz;
1976 Double_t rmin2n = 0.5*(fRmin1*(fDz-z2)+fRmin2*(fDz+z2))/fDz;
1977 Double_t rmax2n = 0.5*(fRmax1*(fDz-z2)+fRmax2*(fDz+z2))/fDz;
1978 shape = new TGeoConeSeg(step/2, rmin1n, rmax1n, rmin2n, rmax2n, fPhi1, fPhi2);
1979 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
1980 vmulti->AddVolume(vol);
1981 opt = "Z";
1982 voldiv->AddNodeOffset(vol, id, start+id*step+step/2, opt.Data());
1983 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
1984 }
1985 return vmulti;
1986 default:
1987 Error("Divide", "Wrong axis type for division");
1988 return 0;
1989 }
1990}
1991
1992////////////////////////////////////////////////////////////////////////////////
1993/// Get range of shape for a given axis.
1994
1996{
1997 xlo = 0;
1998 xhi = 0;
1999 Double_t dx = 0;
2000 switch (iaxis) {
2001 case 2:
2002 xlo = fPhi1;
2003 xhi = fPhi2;
2004 dx = xhi-xlo;
2005 return dx;
2006 case 3:
2007 xlo = -fDz;
2008 xhi = fDz;
2009 dx = xhi-xlo;
2010 return dx;
2011 }
2012 return dx;
2013}
2014
2015////////////////////////////////////////////////////////////////////////////////
2016/// Fill vector param[4] with the bounding cylinder parameters. The order
2017/// is the following : Rmin, Rmax, Phi1, Phi2
2018
2020{
2021 param[0] = TMath::Min(fRmin1, fRmin2); // Rmin
2022 param[0] *= param[0];
2023 param[1] = TMath::Max(fRmax1, fRmax2); // Rmax
2024 param[1] *= param[1];
2025 param[2] = (fPhi1<0)?(fPhi1+360.):fPhi1; // Phi1
2026 param[3] = fPhi2; // Phi2
2027 while (param[3]<param[2]) param[3]+=360.;
2028}
2029
2030////////////////////////////////////////////////////////////////////////////////
2031/// in case shape has some negative parameters, these has to be computed
2032/// in order to fit the mother
2033
2035{
2036 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
2037 if (!mother->TestShapeBit(kGeoConeSeg)) {
2038 Error("GetMakeRuntimeShape", "invalid mother");
2039 return 0;
2040 }
2041 Double_t rmin1, rmax1, rmin2, rmax2, dz;
2042 rmin1 = fRmin1;
2043 rmax1 = fRmax1;
2044 rmin2 = fRmin2;
2045 rmax2 = fRmax2;
2046 dz = fDz;
2047 if (fDz<0) dz=((TGeoCone*)mother)->GetDz();
2048 if (fRmin1<0)
2049 rmin1 = ((TGeoCone*)mother)->GetRmin1();
2050 if ((fRmax1<0) || (fRmax1<fRmin1))
2051 rmax1 = ((TGeoCone*)mother)->GetRmax1();
2052 if (fRmin2<0)
2053 rmin2 = ((TGeoCone*)mother)->GetRmin2();
2054 if ((fRmax2<0) || (fRmax2<fRmin2))
2055 rmax2 = ((TGeoCone*)mother)->GetRmax2();
2056
2057 return (new TGeoConeSeg(GetName(), dz, rmin1, rmax1, rmin2, rmax2, fPhi1, fPhi2));
2058}
2059
2060////////////////////////////////////////////////////////////////////////////////
2061/// print shape parameters
2062
2064{
2065 printf("*** Shape %s: TGeoConeSeg ***\n", GetName());
2066 printf(" dz = %11.5f\n", fDz);
2067 printf(" Rmin1 = %11.5f\n", fRmin1);
2068 printf(" Rmax1 = %11.5f\n", fRmax1);
2069 printf(" Rmin2 = %11.5f\n", fRmin2);
2070 printf(" Rmax2 = %11.5f\n", fRmax2);
2071 printf(" phi1 = %11.5f\n", fPhi1);
2072 printf(" phi2 = %11.5f\n", fPhi2);
2073 printf(" Bounding box:\n");
2075}
2076
2077 ///////////////////////////////////////////////////////////////////////////////
2078 /// Creates a TBuffer3D describing *this* shape.
2079 /// Coordinates are in local reference frame.
2080
2082{
2084 Int_t nbPnts = 4*n;
2085 Int_t nbSegs = 2*nbPnts;
2086 Int_t nbPols = nbPnts-2;
2088 nbPnts, 3*nbPnts,
2089 nbSegs, 3*nbSegs,
2090 nbPols, 6*nbPols);
2091
2092 if (buff)
2093 {
2094 SetPoints(buff->fPnts);
2095 SetSegsAndPols(*buff);
2096 }
2097
2098 return buff;
2099}
2100
2101////////////////////////////////////////////////////////////////////////////////
2102/// Fill TBuffer3D structure for segments and polygons.
2103
2105{
2106 Int_t i, j;
2108 Int_t c = GetBasicColor();
2109
2110 memset(buffer.fSegs, 0, buffer.NbSegs()*3*sizeof(Int_t));
2111 for (i = 0; i < 4; i++) {
2112 for (j = 1; j < n; j++) {
2113 buffer.fSegs[(i*n+j-1)*3 ] = c;
2114 buffer.fSegs[(i*n+j-1)*3+1] = i*n+j-1;
2115 buffer.fSegs[(i*n+j-1)*3+2] = i*n+j;
2116 }
2117 }
2118 for (i = 4; i < 6; i++) {
2119 for (j = 0; j < n; j++) {
2120 buffer.fSegs[(i*n+j)*3 ] = c+1;
2121 buffer.fSegs[(i*n+j)*3+1] = (i-4)*n+j;
2122 buffer.fSegs[(i*n+j)*3+2] = (i-2)*n+j;
2123 }
2124 }
2125 for (i = 6; i < 8; i++) {
2126 for (j = 0; j < n; j++) {
2127 buffer.fSegs[(i*n+j)*3 ] = c;
2128 buffer.fSegs[(i*n+j)*3+1] = 2*(i-6)*n+j;
2129 buffer.fSegs[(i*n+j)*3+2] = (2*(i-6)+1)*n+j;
2130 }
2131 }
2132
2133 Int_t indx = 0;
2134 memset(buffer.fPols, 0, buffer.NbPols()*6*sizeof(Int_t));
2135 i = 0;
2136 for (j = 0; j < n-1; j++) {
2137 buffer.fPols[indx++] = c;
2138 buffer.fPols[indx++] = 4;
2139 buffer.fPols[indx++] = (4+i)*n+j+1;
2140 buffer.fPols[indx++] = (2+i)*n+j;
2141 buffer.fPols[indx++] = (4+i)*n+j;
2142 buffer.fPols[indx++] = i*n+j;
2143 }
2144 i = 1;
2145 for (j = 0; j < n-1; j++) {
2146 buffer.fPols[indx++] = c;
2147 buffer.fPols[indx++] = 4;
2148 buffer.fPols[indx++] = i*n+j;
2149 buffer.fPols[indx++] = (4+i)*n+j;
2150 buffer.fPols[indx++] = (2+i)*n+j;
2151 buffer.fPols[indx++] = (4+i)*n+j+1;
2152 }
2153 i = 2;
2154 for (j = 0; j < n-1; j++) {
2155 buffer.fPols[indx++] = c+i;
2156 buffer.fPols[indx++] = 4;
2157 buffer.fPols[indx++] = (i-2)*2*n+j;
2158 buffer.fPols[indx++] = (4+i)*n+j;
2159 buffer.fPols[indx++] = ((i-2)*2+1)*n+j;
2160 buffer.fPols[indx++] = (4+i)*n+j+1;
2161 }
2162 i = 3;
2163 for (j = 0; j < n-1; j++) {
2164 buffer.fPols[indx++] = c+i;
2165 buffer.fPols[indx++] = 4;
2166 buffer.fPols[indx++] = (4+i)*n+j+1;
2167 buffer.fPols[indx++] = ((i-2)*2+1)*n+j;
2168 buffer.fPols[indx++] = (4+i)*n+j;
2169 buffer.fPols[indx++] = (i-2)*2*n+j;
2170 }
2171 buffer.fPols[indx++] = c+2;
2172 buffer.fPols[indx++] = 4;
2173 buffer.fPols[indx++] = 6*n;
2174 buffer.fPols[indx++] = 4*n;
2175 buffer.fPols[indx++] = 7*n;
2176 buffer.fPols[indx++] = 5*n;
2177 buffer.fPols[indx++] = c+2;
2178 buffer.fPols[indx++] = 4;
2179 buffer.fPols[indx++] = 6*n-1;
2180 buffer.fPols[indx++] = 8*n-1;
2181 buffer.fPols[indx++] = 5*n-1;
2182 buffer.fPols[indx++] = 7*n-1;
2183}
2184
2185////////////////////////////////////////////////////////////////////////////////
2186/// computes the closest distance from given point to this shape, according
2187/// to option. The matching point on the shape is stored in spoint.
2188
2190{
2191 Double_t safe = TGeoCone::Safety(point,in);
2192 if ((fPhi2-fPhi1)>=360.) return safe;
2193 Double_t safphi = TGeoShape::SafetyPhi(point, in, fPhi1, fPhi2);
2194 if (in) return TMath::Min(safe, safphi);
2195 if (safe>1.E10) return safphi;
2196 return TMath::Max(safe, safphi);
2197}
2198
2199////////////////////////////////////////////////////////////////////////////////
2200/// Static method to compute the closest distance from given point to this shape.
2201
2203 Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2, Int_t skipz)
2204{
2205 Double_t safe = TGeoCone::SafetyS(point,in,dz,rmin1,rmax1,rmin2,rmax2,skipz);
2206 if ((phi2-phi1)>=360.) return safe;
2207 Double_t safphi = TGeoShape::SafetyPhi(point,in,phi1,phi2);
2208 if (in) return TMath::Min(safe, safphi);
2209 if (safe>1.E10) return safphi;
2210 return TMath::Max(safe, safphi);
2211}
2212
2213////////////////////////////////////////////////////////////////////////////////
2214/// Save a primitive as a C++ statement(s) on output stream "out".
2215
2216void TGeoConeSeg::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
2217{
2219 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
2220 out << " dz = " << fDz << ";" << std::endl;
2221 out << " rmin1 = " << fRmin1 << ";" << std::endl;
2222 out << " rmax1 = " << fRmax1 << ";" << std::endl;
2223 out << " rmin2 = " << fRmin2 << ";" << std::endl;
2224 out << " rmax2 = " << fRmax2 << ";" << std::endl;
2225 out << " phi1 = " << fPhi1 << ";" << std::endl;
2226 out << " phi2 = " << fPhi2 << ";" << std::endl;
2227 out << " TGeoShape *" << GetPointerName() << " = new TGeoConeSeg(\"" << GetName() << "\", dz,rmin1,rmax1,rmin2,rmax2,phi1,phi2);" << std::endl;
2229}
2230
2231////////////////////////////////////////////////////////////////////////////////
2232/// Set dimensions of the cone segment.
2233
2235 Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2)
2236{
2237 fDz = dz;
2238 fRmin1 = rmin1;
2239 fRmax1 = rmax1;
2240 fRmin2 = rmin2;
2241 fRmax2 = rmax2;
2242 fPhi1 = phi1;
2243 while (fPhi1<0) fPhi1+=360.;
2244 fPhi2 = phi2;
2245 while (fPhi2<=fPhi1) fPhi2+=360.;
2246 if (TGeoShape::IsSameWithinTolerance(fPhi1,fPhi2)) Fatal("SetConsDimensions", "In shape %s invalid phi1=%g, phi2=%g\n", GetName(), fPhi1, fPhi2);
2248}
2249
2250////////////////////////////////////////////////////////////////////////////////
2251/// Set dimensions of the cone segment from an array.
2252
2254{
2255 Double_t dz = param[0];
2256 Double_t rmin1 = param[1];
2257 Double_t rmax1 = param[2];
2258 Double_t rmin2 = param[3];
2259 Double_t rmax2 = param[4];
2260 Double_t phi1 = param[5];
2261 Double_t phi2 = param[6];
2262 SetConsDimensions(dz, rmin1, rmax1,rmin2, rmax2, phi1, phi2);
2263}
2264
2265////////////////////////////////////////////////////////////////////////////////
2266/// Create cone segment mesh points.
2267
2269{
2270 Int_t j, n;
2271 Float_t dphi,phi,phi1, phi2,dz;
2272
2273 n = gGeoManager->GetNsegments()+1;
2274 dz = fDz;
2275 phi1 = fPhi1;
2276 phi2 = fPhi2;
2277
2278 dphi = (phi2-phi1)/(n-1);
2279
2280 Int_t indx = 0;
2281
2282 if (points) {
2283 for (j = 0; j < n; j++) {
2284 phi = (fPhi1+j*dphi)*TMath::DegToRad();
2285 points[indx++] = fRmin1 * TMath::Cos(phi);
2286 points[indx++] = fRmin1 * TMath::Sin(phi);
2287 points[indx++] = -dz;
2288 }
2289 for (j = 0; j < n; j++) {
2290 phi = (fPhi1+j*dphi)*TMath::DegToRad();
2291 points[indx++] = fRmax1 * TMath::Cos(phi);
2292 points[indx++] = fRmax1 * TMath::Sin(phi);
2293 points[indx++] = -dz;
2294 }
2295 for (j = 0; j < n; j++) {
2296 phi = (fPhi1+j*dphi)*TMath::DegToRad();
2297 points[indx++] = fRmin2 * TMath::Cos(phi);
2298 points[indx++] = fRmin2 * TMath::Sin(phi);
2299 points[indx++] = dz;
2300 }
2301 for (j = 0; j < n; j++) {
2302 phi = (fPhi1+j*dphi)*TMath::DegToRad();
2303 points[indx++] = fRmax2 * TMath::Cos(phi);
2304 points[indx++] = fRmax2 * TMath::Sin(phi);
2305 points[indx++] = dz;
2306 }
2307 }
2308}
2309
2310////////////////////////////////////////////////////////////////////////////////
2311/// Create cone segment mesh points.
2312
2314{
2315 Int_t j, n;
2316 Float_t dphi,phi,phi1, phi2,dz;
2317
2318 n = gGeoManager->GetNsegments()+1;
2319 dz = fDz;
2320 phi1 = fPhi1;
2321 phi2 = fPhi2;
2322
2323 dphi = (phi2-phi1)/(n-1);
2324
2325 Int_t indx = 0;
2326
2327 if (points) {
2328 for (j = 0; j < n; j++) {
2329 phi = (fPhi1+j*dphi)*TMath::DegToRad();
2330 points[indx++] = fRmin1 * TMath::Cos(phi);
2331 points[indx++] = fRmin1 * TMath::Sin(phi);
2332 points[indx++] = -dz;
2333 }
2334 for (j = 0; j < n; j++) {
2335 phi = (fPhi1+j*dphi)*TMath::DegToRad();
2336 points[indx++] = fRmax1 * TMath::Cos(phi);
2337 points[indx++] = fRmax1 * TMath::Sin(phi);
2338 points[indx++] = -dz;
2339 }
2340 for (j = 0; j < n; j++) {
2341 phi = (fPhi1+j*dphi)*TMath::DegToRad();
2342 points[indx++] = fRmin2 * TMath::Cos(phi);
2343 points[indx++] = fRmin2 * TMath::Sin(phi);
2344 points[indx++] = dz;
2345 }
2346 for (j = 0; j < n; j++) {
2347 phi = (fPhi1+j*dphi)*TMath::DegToRad();
2348 points[indx++] = fRmax2 * TMath::Cos(phi);
2349 points[indx++] = fRmax2 * TMath::Sin(phi);
2350 points[indx++] = dz;
2351 }
2352 }
2353}
2354
2355////////////////////////////////////////////////////////////////////////////////
2356/// Returns numbers of vertices, segments and polygons composing the shape mesh.
2357
2358void TGeoConeSeg::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
2359{
2361 nvert = n*4;
2362 nsegs = n*8;
2363 npols = n*4-2;
2364}
2365
2366////////////////////////////////////////////////////////////////////////////////
2367/// Return number of vertices of the mesh representation
2368
2370{
2372 Int_t numPoints = n*4;
2373 return numPoints;
2374}
2375
2376////////////////////////////////////////////////////////////////////////////////
2377/// Fill size of this 3-D object
2378
2380{
2381}
2382
2383////////////////////////////////////////////////////////////////////////////////
2384/// Fills a static 3D buffer and returns a reference.
2385
2386const TBuffer3D & TGeoConeSeg::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
2387{
2388 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
2389
2390 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
2391
2392 if (reqSections & TBuffer3D::kRawSizes) {
2394 Int_t nbPnts = 4*n;
2395 Int_t nbSegs = 2*nbPnts;
2396 Int_t nbPols = nbPnts-2;
2397 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
2399 }
2400 }
2401 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
2402 SetPoints(buffer.fPnts);
2403 if (!buffer.fLocalFrame) {
2404 TransformPoints(buffer.fPnts, buffer.NbPnts());
2405 }
2406
2407 SetSegsAndPols(buffer);
2409 }
2410
2411 return buffer;
2412}
2413
2414////////////////////////////////////////////////////////////////////////////////
2415/// Fills array with n random points located on the line segments of the shape mesh.
2416/// The output array must be provided with a length of minimum 3*npoints. Returns
2417/// true if operation is implemented.
2418
2420{
2421 if (npoints > (npoints/2)*2) {
2422 Error("GetPointsOnSegments","Npoints must be even number");
2423 return kFALSE;
2424 }
2425 Int_t nc = (Int_t)TMath::Sqrt(0.5*npoints);
2426 Double_t dphi = (fPhi2-fPhi1)*TMath::DegToRad()/(nc-1);
2427 Double_t phi = 0;
2428 Double_t phi1 = fPhi1 * TMath::DegToRad();
2429 Int_t ntop = npoints/2 - nc*(nc-1);
2430 Double_t dz = 2*fDz/(nc-1);
2431 Double_t z = 0;
2432 Double_t rmin = 0.;
2433 Double_t rmax = 0.;
2434 Int_t icrt = 0;
2435 Int_t nphi = nc;
2436 // loop z sections
2437 for (Int_t i=0; i<nc; i++) {
2438 if (i == (nc-1)) {
2439 nphi = ntop;
2440 dphi = (fPhi2-fPhi1)*TMath::DegToRad()/(nphi-1);
2441 }
2442 z = -fDz + i*dz;
2443 rmin = 0.5*(fRmin1+fRmin2) + 0.5*(fRmin2-fRmin1)*z/fDz;
2444 rmax = 0.5*(fRmax1+fRmax2) + 0.5*(fRmax2-fRmax1)*z/fDz;
2445 // loop points on circle sections
2446 for (Int_t j=0; j<nphi; j++) {
2447 phi = phi1 + j*dphi;
2448 array[icrt++] = rmin * TMath::Cos(phi);
2449 array[icrt++] = rmin * TMath::Sin(phi);
2450 array[icrt++] = z;
2451 array[icrt++] = rmax * TMath::Cos(phi);
2452 array[icrt++] = rmax * TMath::Sin(phi);
2453 array[icrt++] = z;
2454 }
2455 }
2456 return kTRUE;
2457}
2458
2459////////////////////////////////////////////////////////////////////////////////
2460/// Check the inside status for each of the points in the array.
2461/// Input: Array of point coordinates + vector size
2462/// Output: Array of Booleans for the inside of each point
2463
2464void TGeoConeSeg::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
2465{
2466 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
2467}
2468
2469////////////////////////////////////////////////////////////////////////////////
2470/// Compute the normal for an array o points so that norm.dot.dir is positive
2471/// Input: Arrays of point coordinates and directions + vector size
2472/// Output: Array of normal directions
2473
2474void TGeoConeSeg::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
2475{
2476 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
2477}
2478
2479////////////////////////////////////////////////////////////////////////////////
2480/// Compute distance from array of input points having directions specified by dirs. Store output in dists
2481
2482void TGeoConeSeg::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
2483{
2484 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
2485}
2486
2487////////////////////////////////////////////////////////////////////////////////
2488/// Compute distance from array of input points having directions specified by dirs. Store output in dists
2489
2490void TGeoConeSeg::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
2491{
2492 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
2493}
2494
2495////////////////////////////////////////////////////////////////////////////////
2496/// Compute safe distance from each of the points in the input array.
2497/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
2498/// Output: Safety values
2499
2500void TGeoConeSeg::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
2501{
2502 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
2503}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define s1(x)
Definition RSha256.hxx:91
int Int_t
Definition RtypesCore.h:45
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
void Fatal(const char *location, const char *msgfmt,...)
Use this function in case of a fatal error. It will abort the program.
Definition TError.cxx:256
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 id
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
float xmin
float ymin
float xmax
float ymax
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:114
UInt_t NbPols() const
Definition TBuffer3D.h:82
UInt_t NbPnts() const
Definition TBuffer3D.h:80
UInt_t NbSegs() const
Definition TBuffer3D.h:81
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Int_t * fSegs
Definition TBuffer3D.h:113
Bool_t fLocalFrame
Definition TBuffer3D.h:90
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:112
Box class.
Definition TGeoBBox.h:18
Double_t fDX
Definition TGeoBBox.h:21
void SetBoxDimensions(Double_t dx, Double_t dy, Double_t dz, Double_t *origin=nullptr)
Set parameters of the box.
Definition TGeoBBox.cxx:926
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...
A cone segment is a cone having a range in phi.
Definition TGeoCone.h:102
static void ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Compute normal to closest surface from POINT.
Double_t fCdfi
Definition TGeoCone.h:114
virtual Bool_t Contains(const Double_t *point) const
test if point is inside this sphere
virtual TBuffer3D * MakeBuffer3D() const
Creates a TBuffer3D describing this shape.
virtual ~TGeoConeSeg()
destructor
Double_t fCm
Definition TGeoCone.h:113
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
static Double_t DistFromOutsideS(const Double_t *point, const Double_t *dir, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi)
compute distance from outside point to surface of arbitrary tube
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
TGeoConeSeg()
Default constructor.
virtual Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const
Fills array with n random points located on the line segments of the shape mesh.
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 segment
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
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
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
virtual void Sizeof3D() const
Fill size of this 3-D object.
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...
Double_t fC1
Definition TGeoCone.h:109
void InitTrigonometry()
Init frequently used trigonometric values.
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
static Double_t SafetyS(const Double_t *point, Bool_t in, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2, Int_t skipz=0)
Static method to compute the closest distance from given point to this shape.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute closest distance from point px,py to each corner
virtual void SetDimensions(Double_t *param)
Set dimensions of the cone segment from an array.
Double_t fPhi1
Definition TGeoCone.h:105
virtual void InspectShape() const
print shape parameters
void SetConsDimensions(Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2)
Set dimensions of the cone segment.
Double_t fS2
Definition TGeoCone.h:110
virtual void AfterStreamer()
Function called after streaming an object of this class.
Double_t fS1
Definition TGeoCone.h:108
virtual void SetPoints(Double_t *points) const
Create cone segment mesh points.
static Double_t DistFromInsideS(const Double_t *point, const Double_t *dir, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi)
compute distance from inside point to surface of the tube segment
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this cone segment shape belonging to volume "voldiv" into ndiv volumes called divname,...
static Double_t DistToCons(const Double_t *point, const Double_t *dir, Double_t r1, Double_t z1, Double_t r2, Double_t z2, Double_t phi1, Double_t phi2)
Static method to compute distance to a conical surface with :
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Double_t fSm
Definition TGeoCone.h:112
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.
virtual void ComputeBBox()
compute bounding box of the tube segment
virtual void SetSegsAndPols(TBuffer3D &buffer) const
Fill TBuffer3D structure for segments and polygons.
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.
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
Check the inside status for each of the points in the array.
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
Double_t fPhi2
Definition TGeoCone.h:106
Double_t fC2
Definition TGeoCone.h:111
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
The cones are defined by 5 parameters:
Definition TGeoCone.h:18
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
virtual const char * GetAxisName(Int_t iaxis) const
Returns name of axis IAXIS.
Definition TGeoCone.cxx:634
virtual void Sizeof3D() const
Fill size of this 3-D object.
virtual void SetDimensions(Double_t *param)
Set cone dimensions from an array.
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoCone.cxx:674
static void ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
Compute normal to closest surface from POINT.
Definition TGeoCone.cxx:240
void SetConeDimensions(Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
Set cone dimensions.
Definition TGeoCone.cxx:952
static Double_t DistFromInsideS(const Double_t *point, const Double_t *dir, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
Compute distance from inside point to surface of the cone (static) Boundary safe algorithm.
Definition TGeoCone.cxx:291
virtual ~TGeoCone()
destructor
Definition TGeoCone.cxx:180
Double_t fRmax1
Definition TGeoCone.h:23
Double_t fRmax2
Definition TGeoCone.h:25
virtual void ComputeBBox()
compute bounding box of the sphere
Definition TGeoCone.cxx:187
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition TGeoCone.cxx:162
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.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute closest distance from point px,py to each corner
Definition TGeoCone.cxx:562
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual Bool_t Contains(const Double_t *point) const
test if point is inside this cone
Definition TGeoCone.cxx:277
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Definition TGeoCone.cxx:651
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
Definition TGeoCone.cxx:197
static Double_t DistFromOutsideS(const Double_t *point, const Double_t *dir, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
Compute distance from outside point to surface of the tube Boundary safe algorithm.
Definition TGeoCone.cxx:377
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...
Double_t fDz
Definition TGeoCone.h:21
TGeoCone()
Default constructor.
Definition TGeoCone.cxx:101
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 TGeoCone.cxx:884
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
Definition TGeoCone.cxx:507
static void DistToCone(const Double_t *point, const Double_t *dir, Double_t dz, Double_t r1, Double_t r2, Double_t &b, Double_t &delta)
Static method to compute distance to a conical surface with :
Definition TGeoCone.cxx:527
static Double_t SafetyS(const Double_t *point, Bool_t in, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Int_t skipz=0)
computes the closest distance from given point to this shape, according to option.
Definition TGeoCone.cxx:901
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition TGeoCone.cxx:936
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this cone shape belonging to volume "voldiv" into ndiv volumes called divname,...
Definition TGeoCone.cxx:577
virtual TBuffer3D * MakeBuffer3D() const
Creates a TBuffer3D describing this shape.
Definition TGeoCone.cxx:782
virtual Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const
Fills array with n random points located on the line segments of the shape mesh.
Definition TGeoCone.cxx:719
Double_t fRmin1
Definition TGeoCone.h:22
Double_t fRmin2
Definition TGeoCone.h:24
virtual void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
Check the inside status for each of the points in the array.
virtual void SetSegsAndPols(TBuffer3D &buffer) const
Fill TBuffer3D structure for segments and polygons.
Definition TGeoCone.cxx:804
virtual void SetPoints(Double_t *points) const
Create cone mesh points.
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 TGeoCone.cxx:688
virtual void InspectShape() const
print shape parameters
Definition TGeoCone.cxx:766
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 cone Boundary safe algorithm.
Definition TGeoCone.cxx:362
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
Int_t GetNsegments() const
Get number of segments approximating circles.
Geometrical transformation package.
Definition TGeoMatrix.h:41
Node containing an offset.
Definition TGeoNode.h:184
Base finder class for patterns.
void SetDivIndex(Int_t index)
Base abstract class for all shapes.
Definition TGeoShape.h:26
static Double_t Big()
Definition TGeoShape.h:88
Int_t GetBasicColor() const
Get the basic color (0-7).
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
static Double_t DistToPhiMin(const Double_t *point, const Double_t *dir, Double_t s1, Double_t c1, Double_t s2, Double_t c2, Double_t sm, Double_t cm, Bool_t in=kTRUE)
compute distance from point (inside phi) to both phi planes. Return minimum.
static Double_t SafetyPhi(const Double_t *point, Bool_t in, Double_t phi1, Double_t phi2)
Static method to compute safety w.r.t a phi corner defined by cosines/sines of the angles phi1,...
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.
static void NormalPhi(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Static method to compute normal to phi planes.
static Double_t SafetySeg(Double_t r, Double_t z, Double_t r1, Double_t z1, Double_t r2, Double_t z2, Bool_t outer)
Compute distance from point of coordinates (r,z) to segment (r1,z1):(r2,z2)
@ kGeoSavePrimitive
Definition TGeoShape.h:65
@ kGeoRunTimeShape
Definition TGeoShape.h:41
@ kGeoConeSeg
Definition TGeoShape.h:50
static Double_t Tolerance()
Definition TGeoShape.h:91
static Bool_t IsCloseToPhi(Double_t epsil, const Double_t *point, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
True if point is closer than epsil to one of the phi planes defined by c1,s1 or c2,...
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:163
Volume families.
Definition TGeoVolume.h:255
void AddVolume(TGeoVolume *vol)
Add a volume with valid shape to the list of volumes.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
void AddNodeOffset(TGeoVolume *vol, Int_t copy_no, Double_t offset=0, Option_t *option="")
Add a division node to the list of nodes.
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:174
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoVolume.h:232
Int_t GetNdaughters() const
Definition TGeoVolume.h:351
TObjArray * GetNodes()
Definition TGeoVolume.h:168
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
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
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:956
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:774
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:380
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
return c1
Definition legend1.C:41
const Int_t n
Definition legend1.C:16
return c2
Definition legend2.C:14
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Definition TMath.h:980
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
T1 Sign(T1 a, T2 b)
Returns a value with the magnitude of a and the sign of b.
Definition TMathBase.h:175
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:644
Long64_t LocMax(Long64_t n, const T *a)
Returns index of array with the maximum element.
Definition TMath.h:1010
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
constexpr Double_t RadToDeg()
Conversion from radian to degree: .
Definition TMath.h:72
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
constexpr Double_t TwoPi()
Definition TMath.h:44