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