Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoTube.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 24/10/01
3// TGeoTube::Contains() and DistFromInside/In() 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 TGeoTube
14\ingroup Tubes
15
16Cylindrical tube class. A tube has 3 parameters :
17 - `rmin` : minimum radius
18 - `rmax` : maximum radius
19 - `dz` : half length
20
21~~~ {.cpp}
22TGeoTube(Double_t rmin,Double_t rmax,Double_t dz);
23~~~
24
25Begin_Macro
26{
27 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
28 new TGeoManager("tube", "poza2");
29 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
30 TGeoMedium *med = new TGeoMedium("MED",1,mat);
31 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
32 gGeoManager->SetTopVolume(top);
33 TGeoVolume *vol = gGeoManager->MakeTube("TUBE",med, 20,30,40);
34 vol->SetLineWidth(2);
35 top->AddNode(vol,1);
36 gGeoManager->CloseGeometry();
37 gGeoManager->SetNsegments(80);
38 top->Draw();
39 TView *view = gPad->GetView();
40 view->ShowAxis();
41}
42End_Macro
43*/
44
45/** \class TGeoTubeSeg
46\ingroup Tubes
47
48A tube segment is a tube having a range in phi. The tube segment class
49derives from **`TGeoTube`**, having 2 extra parameters: `phi1` and
50`phi2`.
51
52~~~ {.cpp}
53TGeoTubeSeg(Double_t rmin,Double_t rmax,Double_t dz,
54Double_t phi1,Double_t phi2);
55~~~
56
57Here `phi1` and `phi2 `are the starting and ending `phi `values in
58degrees. The `general phi convention` is that the shape ranges from
59`phi1` to `phi2` going counterclockwise. The angles can be defined with
60either negative or positive values. They are stored such that `phi1` is
61converted to `[0,360]` and `phi2 > phi1`.
62
63Begin_Macro
64{
65 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
66 new TGeoManager("tubeseg", "poza3");
67 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
68 TGeoMedium *med = new TGeoMedium("MED",1,mat);
69 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
70 gGeoManager->SetTopVolume(top);
71 TGeoVolume *vol = gGeoManager->MakeTubs("TUBESEG",med, 20,30,40,-30,270);
72 vol->SetLineWidth(2);
73 top->AddNode(vol,1);
74 gGeoManager->CloseGeometry();
75 gGeoManager->SetNsegments(80);
76 top->Draw();
77 TView *view = gPad->GetView();
78 view->ShowAxis();
79}
80End_Macro
81*/
82
83/** \class TGeoCtub
84\ingroup Tubes
85
86The cut tubes constructor has the form:
87
88~~~ {.cpp}
89TGeoCtub(Double_t rmin,Double_t rmax,Double_t dz,
90Double_t phi1,Double_t phi2,
91Double_t nxlow,Double_t nylow,Double_t nzlow, Double_t nxhi,
92Double_t nyhi,Double_t nzhi);
93~~~
94
95Begin_Macro
96{
97 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
98 new TGeoManager("ctub", "poza3");
99 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
100 TGeoMedium *med = new TGeoMedium("MED",1,mat);
101 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
102 gGeoManager->SetTopVolume(top);
103 Double_t theta = 160.*TMath::Pi()/180.;
104 Double_t phi = 30.*TMath::Pi()/180.;
105 Double_t nlow[3];
106 nlow[0] = TMath::Sin(theta)*TMath::Cos(phi);
107 nlow[1] = TMath::Sin(theta)*TMath::Sin(phi);
108 nlow[2] = TMath::Cos(theta);
109 theta = 20.*TMath::Pi()/180.;
110 phi = 60.*TMath::Pi()/180.;
111 Double_t nhi[3];
112 nhi[0] = TMath::Sin(theta)*TMath::Cos(phi);
113 nhi[1] = TMath::Sin(theta)*TMath::Sin(phi);
114 nhi[2] = TMath::Cos(theta);
115 TGeoVolume *vol = gGeoManager->MakeCtub("CTUB",med, 20,30,40,-30,250, nlow[0], nlow[1], nlow[2], nhi[0],nhi[1],nhi[2]);
116 vol->SetLineWidth(2);
117 top->AddNode(vol,1);
118 gGeoManager->CloseGeometry();
119 gGeoManager->SetNsegments(80);
120 top->Draw();
121 TView *view = gPad->GetView();
122 view->ShowAxis();
123}
124End_Macro
125
126A cut tube is a tube segment cut with two planes. The centers of the 2
127sections are positioned at `dZ`. Each cut plane is therefore defined by
128a point `(0,0,dZ)` and its normal unit vector pointing outside the
129shape:
130
131`Nlow=(Nx,Ny,Nz<0)`, `Nhigh=(Nx',Ny',Nz'>0)`.
132*/
133
134#include <iostream>
135
136#include "TGeoManager.h"
137#include "TGeoVolume.h"
138#include "TVirtualGeoPainter.h"
139#include "TGeoTube.h"
140#include "TBuffer3D.h"
141#include "TBuffer3DTypes.h"
142#include "TMath.h"
143
145
146////////////////////////////////////////////////////////////////////////////////
147/// Default constructor
148
150{
152 fRmin = 0.0;
153 fRmax = 0.0;
154 fDz = 0.0;
155}
156
157
158////////////////////////////////////////////////////////////////////////////////
159/// Default constructor specifying minimum and maximum radius
160
162 :TGeoBBox(0, 0, 0)
163{
165 SetTubeDimensions(rmin, rmax, dz);
166 if ((fDz<0) || (fRmin<0) || (fRmax<0)) {
168// if (fRmax<=fRmin) SetShapeBit(kGeoInvalidShape);
169// printf("tube : dz=%f rmin=%f rmax=%f\n", dz, rmin, rmax);
170 }
171 ComputeBBox();
172}
173////////////////////////////////////////////////////////////////////////////////
174/// Default constructor specifying minimum and maximum radius
175
177 :TGeoBBox(name, 0, 0, 0)
178{
180 SetTubeDimensions(rmin, rmax, dz);
181 if ((fDz<0) || (fRmin<0) || (fRmax<0)) {
183// if (fRmax<=fRmin) SetShapeBit(kGeoInvalidShape);
184// printf("tube : dz=%f rmin=%f rmax=%f\n", dz, rmin, rmax);
185 }
186 ComputeBBox();
187}
188
189////////////////////////////////////////////////////////////////////////////////
190/// Default constructor specifying minimum and maximum radius
191/// - param[0] = Rmin
192/// - param[1] = Rmax
193/// - param[2] = dz
194
196 :TGeoBBox(0, 0, 0)
197{
199 SetDimensions(param);
200 if ((fDz<0) || (fRmin<0) || (fRmax<0)) SetShapeBit(kGeoRunTimeShape);
201 ComputeBBox();
202}
203
204////////////////////////////////////////////////////////////////////////////////
205/// destructor
206
208{
209}
210
211////////////////////////////////////////////////////////////////////////////////
212/// Computes capacity of the shape in [length^3]
213
215{
217}
218
219////////////////////////////////////////////////////////////////////////////////
220/// Computes capacity of the shape in [length^3]
221
223{
224 Double_t capacity = 2.*TMath::Pi()*(rmax*rmax-rmin*rmin)*dz;
225 return capacity;
226}
227
228////////////////////////////////////////////////////////////////////////////////
229/// compute bounding box of the tube
230
232{
233 fDX = fDY = fRmax;
234 fDZ = fDz;
235}
236
237////////////////////////////////////////////////////////////////////////////////
238/// Compute normal to closest surface from POINT.
239
240void TGeoTube::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
241{
242 Double_t saf[3];
243 Double_t rsq = point[0]*point[0]+point[1]*point[1];
244 Double_t r = TMath::Sqrt(rsq);
245 saf[0] = TMath::Abs(fDz-TMath::Abs(point[2]));
246 saf[1] = (fRmin>1E-10)?TMath::Abs(r-fRmin):TGeoShape::Big();
247 saf[2] = TMath::Abs(fRmax-r);
248 Int_t i = TMath::LocMin(3,saf);
249 if (i==0) {
250 norm[0] = norm[1] = 0.;
251 norm[2] = TMath::Sign(1.,dir[2]);
252 return;
253 }
254 norm[2] = 0;
255 Double_t phi = TMath::ATan2(point[1], point[0]);
256 norm[0] = TMath::Cos(phi);
257 norm[1] = TMath::Sin(phi);
258 if (norm[0]*dir[0]+norm[1]*dir[1]<0) {
259 norm[0] = -norm[0];
260 norm[1] = -norm[1];
261 }
262}
263
264////////////////////////////////////////////////////////////////////////////////
265/// Compute normal to closest surface from POINT.
266
267void TGeoTube::ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm,
268 Double_t /*rmin*/, Double_t /*rmax*/, Double_t /*dz*/)
269{
270 norm[2] = 0;
271 Double_t phi = TMath::ATan2(point[1], point[0]);
272 norm[0] = TMath::Cos(phi);
273 norm[1] = TMath::Sin(phi);
274 if (norm[0]*dir[0]+norm[1]*dir[1]<0) {
275 norm[0] = -norm[0];
276 norm[1] = -norm[1];
277 }
278}
279
280////////////////////////////////////////////////////////////////////////////////
281/// test if point is inside this tube
282
284{
285 if (TMath::Abs(point[2]) > fDz) return kFALSE;
286 Double_t r2 = point[0]*point[0]+point[1]*point[1];
287 if ((r2<fRmin*fRmin) || (r2>fRmax*fRmax)) return kFALSE;
288 return kTRUE;
289}
290
291////////////////////////////////////////////////////////////////////////////////
292/// compute closest distance from point px,py to each corner
293
295{
297 Int_t numPoints = 4*n;
298 if (!HasRmin()) numPoints = 2*(n+1);
299 return ShapeDistancetoPrimitive(numPoints, px, py);
300}
301
302////////////////////////////////////////////////////////////////////////////////
303/// Compute distance from inside point to surface of the tube (static)
304/// Boundary safe algorithm.
305/// compute distance to surface
306/// Do Z
307
309{
311 if (dir[2]) {
312 sz = (TMath::Sign(dz, dir[2])-point[2])/dir[2];
313 if (sz<=0) return 0.0;
314 }
315 // Do R
316 Double_t nsq=dir[0]*dir[0]+dir[1]*dir[1];
317 if (TMath::Abs(nsq)<TGeoShape::Tolerance()) return sz;
318 Double_t rsq=point[0]*point[0]+point[1]*point[1];
319 Double_t rdotn=point[0]*dir[0]+point[1]*dir[1];
320 Double_t b,d;
321 // inner cylinder
322 if (rmin>0) {
323 // Protection in case point is actually outside the tube
324 if (rsq <= rmin*rmin+TGeoShape::Tolerance()) {
325 if (rdotn<0) return 0.0;
326 } else {
327 if (rdotn<0) {
328 DistToTube(rsq,nsq,rdotn,rmin,b,d);
329 if (d>0) {
330 Double_t sr = -b-d;
331 if (sr>0) return TMath::Min(sz,sr);
332 }
333 }
334 }
335 }
336 // outer cylinder
337 if (rsq >= rmax*rmax-TGeoShape::Tolerance()) {
338 if (rdotn>=0) return 0.0;
339 }
340 DistToTube(rsq,nsq,rdotn,rmax,b,d);
341 if (d>0) {
342 Double_t sr = -b+d;
343 if (sr>0) return TMath::Min(sz,sr);
344 }
345 return 0.;
346}
347
348////////////////////////////////////////////////////////////////////////////////
349/// Compute distance from inside point to surface of the tube
350/// Boundary safe algorithm.
351
352Double_t TGeoTube::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
353{
354 if (iact<3 && safe) {
355 *safe = Safety(point, kTRUE);
356 if (iact==0) return TGeoShape::Big();
357 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
358 }
359 // compute distance to surface
360 return DistFromInsideS(point, dir, fRmin, fRmax, fDz);
361}
362
363////////////////////////////////////////////////////////////////////////////////
364/// Static method to compute distance from outside point to a tube with given parameters
365/// Boundary safe algorithm.
366/// check Z planes
367
369{
370 Double_t xi,yi,zi;
371 Double_t rmaxsq = rmax*rmax;
372 Double_t rminsq = rmin*rmin;
373 zi = dz - TMath::Abs(point[2]);
374 Bool_t in = kFALSE;
375 Bool_t inz = (zi<0)?kFALSE:kTRUE;
376 if (!inz) {
377 if (point[2]*dir[2]>=0) return TGeoShape::Big();
378 Double_t s = -zi/TMath::Abs(dir[2]);
379 xi = point[0]+s*dir[0];
380 yi = point[1]+s*dir[1];
381 Double_t r2=xi*xi+yi*yi;
382 if ((rminsq<=r2) && (r2<=rmaxsq)) return s;
383 }
384
385 Double_t rsq = point[0]*point[0]+point[1]*point[1];
386 // check outer cyl. surface
387 Double_t nsq=dir[0]*dir[0]+dir[1]*dir[1];
388 Double_t rdotn=point[0]*dir[0]+point[1]*dir[1];
389 Double_t b,d;
390 Bool_t inrmax = kFALSE;
391 Bool_t inrmin = kFALSE;
392 if (rsq<=rmaxsq+TGeoShape::Tolerance()) inrmax = kTRUE;
393 if (rsq>=rminsq-TGeoShape::Tolerance()) inrmin = kTRUE;
394 in = inz & inrmin & inrmax;
395 // If inside, we are most likely on a boundary within machine precision.
396 if (in) {
397 Bool_t checkout = kFALSE;
398 Double_t r = TMath::Sqrt(rsq);
399 if (zi<rmax-r) {
400 if ((TGeoShape::IsSameWithinTolerance(rmin,0)) || (zi<r-rmin)) {
401 if (point[2]*dir[2]<0) return 0.0;
402 return TGeoShape::Big();
403 }
404 }
405 if ((rmaxsq-rsq) < (rsq-rminsq)) checkout = kTRUE;
406 if (checkout) {
407 if (rdotn>=0) return TGeoShape::Big();
408 return 0.0;
409 }
410 if (TGeoShape::IsSameWithinTolerance(rmin,0)) return 0.0;
411 if (rdotn>=0) return 0.0;
412 // Ray exiting rmin -> check (+) solution for inner tube
414 DistToTube(rsq, nsq, rdotn, rmin, b, d);
415 if (d>0) {
416 Double_t s=-b+d;
417 if (s>0) {
418 zi=point[2]+s*dir[2];
419 if (TMath::Abs(zi)<=dz) return s;
420 }
421 }
422 return TGeoShape::Big();
423 }
424 // Check outer cylinder (only r>rmax has to be considered)
426 if (!inrmax) {
427 DistToTube(rsq, nsq, rdotn, rmax, b, d);
428 if (d>0) {
429 Double_t s=-b-d;
430 if (s>0) {
431 zi=point[2]+s*dir[2];
432 if (TMath::Abs(zi)<=dz) return s;
433 }
434 }
435 }
436 // check inner cylinder
437 if (rmin>0) {
438 DistToTube(rsq, nsq, rdotn, rmin, b, d);
439 if (d>0) {
440 Double_t s=-b+d;
441 if (s>0) {
442 zi=point[2]+s*dir[2];
443 if (TMath::Abs(zi)<=dz) return s;
444 }
445 }
446 }
447 return TGeoShape::Big();
448}
449
450////////////////////////////////////////////////////////////////////////////////
451/// Compute distance from outside point to surface of the tube and safe distance
452/// Boundary safe algorithm.
453/// fist localize point w.r.t tube
454
455Double_t TGeoTube::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
456{
457 if (iact<3 && safe) {
458 *safe = Safety(point, kFALSE);
459 if (iact==0) return TGeoShape::Big();
460 if ((iact==1) && (step<=*safe)) return TGeoShape::Big();
461 }
462// Check if the bounding box is crossed within the requested distance
463 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
464 if (sdist>=step) return TGeoShape::Big();
465 // find distance to shape
466 return DistFromOutsideS(point, dir, fRmin, fRmax, fDz);
467}
468
469////////////////////////////////////////////////////////////////////////////////
470/// Static method computing the distance to a tube with given radius, starting from
471/// POINT along DIR director cosines. The distance is computed as :
472/// RSQ = point[0]*point[0]+point[1]*point[1]
473/// NSQ = dir[0]*dir[0]+dir[1]*dir[1] ---> should NOT be 0 !!!
474/// RDOTN = point[0]*dir[0]+point[1]*dir[1]
475/// The distance can be computed as :
476/// D = -B +/- DELTA
477/// where DELTA.GT.0 and D.GT.0
478
480{
481 Double_t t1 = 1./nsq;
482 Double_t t3=rsq-(radius*radius);
483 b = t1*rdotn;
484 Double_t c =t1*t3;
485 delta = b*b-c;
486 if (delta>0) {
487 delta=TMath::Sqrt(delta);
488 } else {
489 delta = -1;
490 }
491}
492
493////////////////////////////////////////////////////////////////////////////////
494/// Divide this tube shape belonging to volume "voldiv" into ndiv volumes
495/// called divname, from start position with the given step. Returns pointer
496/// to created division cell volume in case of Z divisions. For radial division
497/// creates all volumes with different shapes and returns pointer to volume that
498/// was divided. In case a wrong division axis is supplied, returns pointer to
499/// volume that was divided.
500
501TGeoVolume *TGeoTube::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
502 Double_t start, Double_t step)
503{
504 TGeoShape *shape; //--- shape to be created
505 TGeoVolume *vol; //--- division volume to be created
506 TGeoVolumeMulti *vmulti; //--- generic divided volume
507 TGeoPatternFinder *finder; //--- finder to be attached
508 TString opt = ""; //--- option to be attached
509 Int_t id;
510 Double_t end = start+ndiv*step;
511 switch (iaxis) {
512 case 1: //--- R division
513 finder = new TGeoPatternCylR(voldiv, ndiv, start, end);
514 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
515 voldiv->SetFinder(finder);
516 finder->SetDivIndex(voldiv->GetNdaughters());
517 for (id=0; id<ndiv; id++) {
518 shape = new TGeoTube(start+id*step, start+(id+1)*step, fDz);
519 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
520 vmulti->AddVolume(vol);
521 opt = "R";
522 voldiv->AddNodeOffset(vol, id, 0, opt.Data());
523 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
524 }
525 return vmulti;
526 case 2: //--- Phi division
527 finder = new TGeoPatternCylPhi(voldiv, ndiv, start, end);
528 voldiv->SetFinder(finder);
529 finder->SetDivIndex(voldiv->GetNdaughters());
530 shape = new TGeoTubeSeg(fRmin, fRmax, fDz, -step/2, step/2);
531 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
532 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
533 vmulti->AddVolume(vol);
534 opt = "Phi";
535 for (id=0; id<ndiv; id++) {
536 voldiv->AddNodeOffset(vol, id, start+id*step+step/2, opt.Data());
537 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
538 }
539 return vmulti;
540 case 3: //--- Z division
541 finder = new TGeoPatternZ(voldiv, ndiv, start, start+ndiv*step);
542 voldiv->SetFinder(finder);
543 finder->SetDivIndex(voldiv->GetNdaughters());
544 shape = new TGeoTube(fRmin, fRmax, step/2);
545 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
546 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
547 vmulti->AddVolume(vol);
548 opt = "Z";
549 for (id=0; id<ndiv; id++) {
550 voldiv->AddNodeOffset(vol, id, start+step/2+id*step, opt.Data());
551 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
552 }
553 return vmulti;
554 default:
555 Error("Divide", "In shape %s wrong axis type for division", GetName());
556 return 0;
557 }
558}
559
560////////////////////////////////////////////////////////////////////////////////
561/// Returns name of axis IAXIS.
562
563const char *TGeoTube::GetAxisName(Int_t iaxis) const
564{
565 switch (iaxis) {
566 case 1:
567 return "R";
568 case 2:
569 return "PHI";
570 case 3:
571 return "Z";
572 default:
573 return "UNDEFINED";
574 }
575}
576
577////////////////////////////////////////////////////////////////////////////////
578/// Get range of shape for a given axis.
579
581{
582 xlo = 0;
583 xhi = 0;
584 Double_t dx = 0;
585 switch (iaxis) {
586 case 1:
587 xlo = fRmin;
588 xhi = fRmax;
589 dx = xhi-xlo;
590 return dx;
591 case 2:
592 xlo = 0;
593 xhi = 360;
594 dx = 360;
595 return dx;
596 case 3:
597 xlo = -fDz;
598 xhi = fDz;
599 dx = xhi-xlo;
600 return dx;
601 }
602 return dx;
603}
604
605////////////////////////////////////////////////////////////////////////////////
606/// Fill vector param[4] with the bounding cylinder parameters. The order
607/// is the following : Rmin, Rmax, Phi1, Phi2, dZ
608
610{
611 param[0] = fRmin; // Rmin
612 param[0] *= param[0];
613 param[1] = fRmax; // Rmax
614 param[1] *= param[1];
615 param[2] = 0.; // Phi1
616 param[3] = 360.; // Phi2
617}
618
619////////////////////////////////////////////////////////////////////////////////
620/// in case shape has some negative parameters, these has to be computed
621/// in order to fit the mother
622
624{
625 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
626 Double_t rmin, rmax, dz;
628 rmin = fRmin;
629 rmax = fRmax;
630 dz = fDz;
631 if (fDz<0) {
632 mother->GetAxisRange(3,xmin,xmax);
633 if (xmax<0) return 0;
634 dz=xmax;
635 }
636 mother->GetAxisRange(1,xmin,xmax);
637 if (fRmin<0) {
638 if (xmin<0) return 0;
639 rmin = xmin;
640 }
641 if (fRmax<0) {
642 if (xmax<=0) return 0;
643 rmax = xmax;
644 }
645
646 return (new TGeoTube(GetName(), rmin, rmax, dz));
647}
648
649////////////////////////////////////////////////////////////////////////////////
650/// print shape parameters
651
653{
654 printf("*** Shape %s: TGeoTube ***\n", GetName());
655 printf(" Rmin = %11.5f\n", fRmin);
656 printf(" Rmax = %11.5f\n", fRmax);
657 printf(" dz = %11.5f\n", fDz);
658 printf(" Bounding box:\n");
660}
661
662////////////////////////////////////////////////////////////////////////////////
663/// Creates a TBuffer3D describing *this* shape.
664/// Coordinates are in local reference frame.
665
667{
669 Int_t nbPnts = 4*n;
670 Int_t nbSegs = 8*n;
671 Int_t nbPols = 4*n;
672 if (!HasRmin()) {
673 nbPnts = 2*(n+1);
674 nbSegs = 5*n;
675 nbPols = 3*n;
676 }
678 nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols);
679 if (buff)
680 {
681 SetPoints(buff->fPnts);
682 SetSegsAndPols(*buff);
683 }
684
685 return buff;
686}
687
688////////////////////////////////////////////////////////////////////////////////
689/// Fill TBuffer3D structure for segments and polygons.
690
692{
693 Int_t i, j,indx;
695 Int_t c = (((buffer.fColor) %8) -1) * 4;
696 if (c < 0) c = 0;
697
698 if (HasRmin()) {
699 // circle segments:
700 // lower rmin circle: i=0, (0, n-1)
701 // lower rmax circle: i=1, (n, 2n-1)
702 // upper rmin circle: i=2, (2n, 3n-1)
703 // upper rmax circle: i=1, (3n, 4n-1)
704 for (i = 0; i < 4; i++) {
705 for (j = 0; j < n; j++) {
706 indx = 3*(i*n+j);
707 buffer.fSegs[indx ] = c;
708 buffer.fSegs[indx+1] = i*n+j;
709 buffer.fSegs[indx+2] = i*n+(j+1)%n;
710 }
711 }
712 // Z-parallel segments
713 // inner: i=4, (4n, 5n-1)
714 // outer: i=5, (5n, 6n-1)
715 for (i = 4; i < 6; i++) {
716 for (j = 0; j < n; j++) {
717 indx = 3*(i*n+j);
718 buffer.fSegs[indx ] = c+1;
719 buffer.fSegs[indx+1] = (i-4)*n+j;
720 buffer.fSegs[indx+2] = (i-2)*n+j;
721 }
722 }
723 // Radial segments
724 // lower: i=6, (6n, 7n-1)
725 // upper: i=7, (7n, 8n-1)
726 for (i = 6; i < 8; i++) {
727 for (j = 0; j < n; j++) {
728 indx = 3*(i*n+j);
729 buffer.fSegs[indx ] = c;
730 buffer.fSegs[indx+1] = 2*(i-6)*n+j;
731 buffer.fSegs[indx+2] = (2*(i-6)+1)*n+j;
732 }
733 }
734 // Polygons
735 i=0;
736 // Inner lateral (0, n-1)
737 for (j = 0; j < n; j++) {
738 indx = 6*(i*n+j);
739 buffer.fPols[indx ] = c;
740 buffer.fPols[indx+1] = 4;
741 buffer.fPols[indx+2] = j;
742 buffer.fPols[indx+3] = 4*n+(j+1)%n;
743 buffer.fPols[indx+4] = 2*n+j;
744 buffer.fPols[indx+5] = 4*n+j;
745 }
746 i=1;
747 // Outer lateral (n,2n-1)
748 for (j = 0; j < n; j++) {
749 indx = 6*(i*n+j);
750 buffer.fPols[indx ] = c+1;
751 buffer.fPols[indx+1] = 4;
752 buffer.fPols[indx+2] = n+j;
753 buffer.fPols[indx+3] = 5*n+j;
754 buffer.fPols[indx+4] = 3*n+j;
755 buffer.fPols[indx+5] = 5*n+(j+1)%n;
756 }
757 i=2;
758 // lower disc (2n, 3n-1)
759 for (j = 0; j < n; j++) {
760 indx = 6*(i*n+j);
761 buffer.fPols[indx ] = c;
762 buffer.fPols[indx+1] = 4;
763 buffer.fPols[indx+2] = j;
764 buffer.fPols[indx+3] = 6*n+j;
765 buffer.fPols[indx+4] = n+j;
766 buffer.fPols[indx+5] = 6*n+(j+1)%n;
767 }
768 i=3;
769 // upper disc (3n, 4n-1)
770 for (j = 0; j < n; j++) {
771 indx = 6*(i*n+j);
772 buffer.fPols[indx ] = c;
773 buffer.fPols[indx+1] = 4;
774 buffer.fPols[indx+2] = 2*n+j;
775 buffer.fPols[indx+3] = 7*n+(j+1)%n;
776 buffer.fPols[indx+4] = 3*n+j;
777 buffer.fPols[indx+5] = 7*n+j;
778 }
779 return;
780 }
781 // Rmin=0 tubes
782 // circle segments
783 // lower rmax circle: i=0, (0, n-1)
784 // upper rmax circle: i=1, (n, 2n-1)
785 for (i = 0; i < 2; i++) {
786 for (j = 0; j < n; j++) {
787 indx = 3*(i*n+j);
788 buffer.fSegs[indx ] = c;
789 buffer.fSegs[indx+1] = 2+i*n+j;
790 buffer.fSegs[indx+2] = 2+i*n+(j+1)%n;
791 }
792 }
793 // Z-parallel segments (2n,3n-1)
794 for (j = 0; j < n; j++) {
795 indx = 3*(2*n+j);
796 buffer.fSegs[indx ] = c+1;
797 buffer.fSegs[indx+1] = 2+j;
798 buffer.fSegs[indx+2] = 2+n+j;
799 }
800 // Radial segments
801 // Lower circle: i=3, (3n,4n-1)
802 // Upper circle: i=4, (4n,5n-1)
803 for (i = 3; i < 5; i++) {
804 for (j = 0; j < n; j++) {
805 indx = 3*(i*n+j);
806 buffer.fSegs[indx ] = c;
807 buffer.fSegs[indx+1] = i-3;
808 buffer.fSegs[indx+2] = 2+(i-3)*n+j;
809 }
810 }
811 // Polygons
812 // lateral (0,n-1)
813 for (j = 0; j < n; j++) {
814 indx = 6*j;
815 buffer.fPols[indx ] = c+1;
816 buffer.fPols[indx+1] = 4;
817 buffer.fPols[indx+2] = j;
818 buffer.fPols[indx+3] = 2*n+j;
819 buffer.fPols[indx+4] = n+j;
820 buffer.fPols[indx+5] = 2*n+(j+1)%n;
821 }
822 // bottom triangles (n,2n-1)
823 for (j = 0; j < n; j++) {
824 indx = 6*n + 5*j;
825 buffer.fPols[indx ] = c;
826 buffer.fPols[indx+1] = 3;
827 buffer.fPols[indx+2] = j;
828 buffer.fPols[indx+3] = 3*n+(j+1)%n;
829 buffer.fPols[indx+4] = 3*n+j;
830 }
831 // top triangles (2n,3n-1)
832 for (j = 0; j < n; j++) {
833 indx = 6*n + 5*n + 5*j;
834 buffer.fPols[indx ] = c;
835 buffer.fPols[indx+1] = 3;
836 buffer.fPols[indx+2] = n+j;
837 buffer.fPols[indx+3] = 4*n+j;
838 buffer.fPols[indx+4] = 4*n+(j+1)%n;
839 }
840}
841
842////////////////////////////////////////////////////////////////////////////////
843/// computes the closest distance from given point to this shape, according
844/// to option. The matching point on the shape is stored in spoint.
845
847{
848#ifndef NEVER
849 Double_t r = TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
850 Double_t safe, safrmin, safrmax;
851 if (in) {
852 safe = fDz-TMath::Abs(point[2]); // positive if inside
853 if (fRmin>1E-10) {
854 safrmin = r-fRmin;
855 if (safrmin < safe) safe = safrmin;
856 }
857 safrmax = fRmax-r;
858 if (safrmax < safe) safe = safrmax;
859 } else {
860 safe = -fDz+TMath::Abs(point[2]);
861 if (fRmin>1E-10) {
862 safrmin = -r+fRmin;
863 if (safrmin > safe) safe = safrmin;
864 }
865 safrmax = -fRmax+r;
866 if (safrmax > safe) safe = safrmax;
867 }
868 return safe;
869#else
870 Double_t saf[3];
871 Double_t rsq = point[0]*point[0]+point[1]*point[1];
872 Double_t r = TMath::Sqrt(rsq);
873 saf[0] = fDz-TMath::Abs(point[2]); // positive if inside
874 saf[1] = (fRmin>1E-10)?(r-fRmin):TGeoShape::Big();
875 saf[2] = fRmax-r;
876 if (in) return saf[TMath::LocMin(3,saf)];
877 for (Int_t i=0; i<3; i++) saf[i]=-saf[i];
878 return saf[TMath::LocMax(3,saf)];
879#endif
880}
881
882////////////////////////////////////////////////////////////////////////////////
883/// computes the closest distance from given point to this shape, according
884/// to option. The matching point on the shape is stored in spoint.
885
886Double_t TGeoTube::SafetyS(const Double_t *point, Bool_t in, Double_t rmin, Double_t rmax, Double_t dz, Int_t skipz)
887{
888 Double_t saf[3];
889 Double_t rsq = point[0]*point[0]+point[1]*point[1];
890 Double_t r = TMath::Sqrt(rsq);
891 switch (skipz) {
892 case 1: // skip lower Z plane
893 saf[0] = dz - point[2];
894 break;
895 case 2: // skip upper Z plane
896 saf[0] = dz + point[2];
897 break;
898 case 3: // skip both
899 saf[0] = TGeoShape::Big();
900 break;
901 default:
902 saf[0] = dz-TMath::Abs(point[2]);
903 }
904 saf[1] = (rmin>1E-10)?(r-rmin):TGeoShape::Big();
905 saf[2] = rmax-r;
906// printf("saf0=%g saf1=%g saf2=%g in=%d skipz=%d\n", saf[0],saf[1],saf[2],in,skipz);
907 if (in) return saf[TMath::LocMin(3,saf)];
908 for (Int_t i=0; i<3; i++) saf[i]=-saf[i];
909 return saf[TMath::LocMax(3,saf)];
910}
911
912////////////////////////////////////////////////////////////////////////////////
913/// Save a primitive as a C++ statement(s) on output stream "out".
914
915void TGeoTube::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
916{
918 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
919 out << " rmin = " << fRmin << ";" << std::endl;
920 out << " rmax = " << fRmax << ";" << std::endl;
921 out << " dz = " << fDz << ";" << std::endl;
922 out << " TGeoShape *" << GetPointerName() << " = new TGeoTube(\"" << GetName() << "\",rmin,rmax,dz);" << std::endl;
924}
925
926////////////////////////////////////////////////////////////////////////////////
927/// Set tube dimensions.
928
930{
931 fRmin = rmin;
932 fRmax = rmax;
933 fDz = dz;
934 if (fRmin>0 && fRmax>0 && fRmin>=fRmax)
935 Error("SetTubeDimensions", "In shape %s wrong rmin=%g rmax=%g", GetName(), rmin,rmax);
936}
937
938////////////////////////////////////////////////////////////////////////////////
939/// Set tube dimensions starting from a list.
940
942{
943 Double_t rmin = param[0];
944 Double_t rmax = param[1];
945 Double_t dz = param[2];
946 SetTubeDimensions(rmin, rmax, dz);
947}
948
949////////////////////////////////////////////////////////////////////////////////
950/// Fills array with n random points located on the line segments of the shape mesh.
951/// The output array must be provided with a length of minimum 3*npoints. Returns
952/// true if operation is implemented.
953
955{
956 if (npoints > (npoints/2)*2) {
957 Error("GetPointsOnSegments","Npoints must be even number");
958 return kFALSE;
959 }
960 Int_t nc = 0;
961 if (HasRmin()) nc = (Int_t)TMath::Sqrt(0.5*npoints);
962 else nc = (Int_t)TMath::Sqrt(1.*npoints);
963 Double_t dphi = TMath::TwoPi()/nc;
964 Double_t phi = 0;
965 Int_t ntop = 0;
966 if (HasRmin()) ntop = npoints/2 - nc*(nc-1);
967 else ntop = npoints - nc*(nc-1);
968 Double_t dz = 2*fDz/(nc-1);
969 Double_t z = 0;
970 Int_t icrt = 0;
971 Int_t nphi = nc;
972 // loop z sections
973 for (Int_t i=0; i<nc; i++) {
974 if (i == (nc-1)) nphi = ntop;
975 z = -fDz + i*dz;
976 // loop points on circle sections
977 for (Int_t j=0; j<nphi; j++) {
978 phi = j*dphi;
979 if (HasRmin()) {
980 array[icrt++] = fRmin * TMath::Cos(phi);
981 array[icrt++] = fRmin * TMath::Sin(phi);
982 array[icrt++] = z;
983 }
984 array[icrt++] = fRmax * TMath::Cos(phi);
985 array[icrt++] = fRmax * TMath::Sin(phi);
986 array[icrt++] = z;
987 }
988 }
989 return kTRUE;
990}
991
992////////////////////////////////////////////////////////////////////////////////
993/// create tube mesh points
994
996{
997 Double_t dz;
998 Int_t j, n;
1000 Double_t dphi = 360./n;
1001 Double_t phi = 0;
1002 dz = fDz;
1003 Int_t indx = 0;
1004 if (points) {
1005 if (HasRmin()) {
1006 // 4*n points
1007 // (0,n-1) lower rmin circle
1008 // (2n, 3n-1) upper rmin circle
1009 for (j = 0; j < n; j++) {
1010 phi = j*dphi*TMath::DegToRad();
1011 points[indx+6*n] = points[indx] = fRmin * TMath::Cos(phi);
1012 indx++;
1013 points[indx+6*n] = points[indx] = fRmin * TMath::Sin(phi);
1014 indx++;
1015 points[indx+6*n] = dz;
1016 points[indx] =-dz;
1017 indx++;
1018 }
1019 // (n, 2n-1) lower rmax circle
1020 // (3n, 4n-1) upper rmax circle
1021 for (j = 0; j < n; j++) {
1022 phi = j*dphi*TMath::DegToRad();
1023 points[indx+6*n] = points[indx] = fRmax * TMath::Cos(phi);
1024 indx++;
1025 points[indx+6*n] = points[indx] = fRmax * TMath::Sin(phi);
1026 indx++;
1027 points[indx+6*n]= dz;
1028 points[indx] =-dz;
1029 indx++;
1030 }
1031 } else {
1032 // centers of lower/upper circles (0,1)
1033 points[indx++] = 0.;
1034 points[indx++] = 0.;
1035 points[indx++] = -dz;
1036 points[indx++] = 0.;
1037 points[indx++] = 0.;
1038 points[indx++] = dz;
1039 // lower rmax circle (2, 2+n-1)
1040 // upper rmax circle (2+n, 2+2n-1)
1041 for (j = 0; j < n; j++) {
1042 phi = j*dphi*TMath::DegToRad();
1043 points[indx+3*n] = points[indx] = fRmax * TMath::Cos(phi);
1044 indx++;
1045 points[indx+3*n] = points[indx] = fRmax * TMath::Sin(phi);
1046 indx++;
1047 points[indx+3*n]= dz;
1048 points[indx] =-dz;
1049 indx++;
1050 }
1051 }
1052 }
1053}
1054
1055////////////////////////////////////////////////////////////////////////////////
1056/// create tube mesh points
1057
1059{
1060 Double_t dz;
1061 Int_t j, n;
1063 Double_t dphi = 360./n;
1064 Double_t phi = 0;
1065 dz = fDz;
1066 Int_t indx = 0;
1067 if (points) {
1068 if (HasRmin()) {
1069 // 4*n points
1070 // (0,n-1) lower rmin circle
1071 // (2n, 3n-1) upper rmin circle
1072 for (j = 0; j < n; j++) {
1073 phi = j*dphi*TMath::DegToRad();
1074 points[indx+6*n] = points[indx] = fRmin * TMath::Cos(phi);
1075 indx++;
1076 points[indx+6*n] = points[indx] = fRmin * TMath::Sin(phi);
1077 indx++;
1078 points[indx+6*n] = dz;
1079 points[indx] =-dz;
1080 indx++;
1081 }
1082 // (n, 2n-1) lower rmax circle
1083 // (3n, 4n-1) upper rmax circle
1084 for (j = 0; j < n; j++) {
1085 phi = j*dphi*TMath::DegToRad();
1086 points[indx+6*n] = points[indx] = fRmax * TMath::Cos(phi);
1087 indx++;
1088 points[indx+6*n] = points[indx] = fRmax * TMath::Sin(phi);
1089 indx++;
1090 points[indx+6*n]= dz;
1091 points[indx] =-dz;
1092 indx++;
1093 }
1094 } else {
1095 // centers of lower/upper circles (0,1)
1096 points[indx++] = 0.;
1097 points[indx++] = 0.;
1098 points[indx++] = -dz;
1099 points[indx++] = 0.;
1100 points[indx++] = 0.;
1101 points[indx++] = dz;
1102 // lower rmax circle (2, 2+n-1)
1103 // upper rmax circle (2+n, 2+2n-1)
1104 for (j = 0; j < n; j++) {
1105 phi = j*dphi*TMath::DegToRad();
1106 points[indx+3*n] = points[indx] = fRmax * TMath::Cos(phi);
1107 indx++;
1108 points[indx+3*n] = points[indx] = fRmax * TMath::Sin(phi);
1109 indx++;
1110 points[indx+3*n]= dz;
1111 points[indx] =-dz;
1112 indx++;
1113 }
1114 }
1115 }
1116}
1117
1118////////////////////////////////////////////////////////////////////////////////
1119/// Return number of vertices of the mesh representation
1120
1122{
1124 Int_t numPoints = n*4;
1125 if (!HasRmin()) numPoints = 2*(n+1);
1126 return numPoints;
1127}
1128
1129////////////////////////////////////////////////////////////////////////////////
1130/// Returns numbers of vertices, segments and polygons composing the shape mesh.
1131
1132void TGeoTube::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
1133{
1135 nvert = n*4;
1136 nsegs = n*8;
1137 npols = n*4;
1138 if (!HasRmin()) {
1139 nvert = 2*(n+1);
1140 nsegs = 5*n;
1141 npols = 3*n;
1142 } else {
1143 nvert = n*4;
1144 nsegs = n*8;
1145 npols = n*4;
1146 }
1147}
1148
1149////////////////////////////////////////////////////////////////////////////////
1150/// fill size of this 3-D object
1151
1153{
1154}
1155
1156////////////////////////////////////////////////////////////////////////////////
1157/// Fills a static 3D buffer and returns a reference.
1158
1159const TBuffer3D & TGeoTube::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
1160{
1161 static TBuffer3DTube buffer;
1162 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
1163
1164 if (reqSections & TBuffer3D::kShapeSpecific) {
1165 buffer.fRadiusInner = fRmin;
1166 buffer.fRadiusOuter = fRmax;
1167 buffer.fHalfLength = fDz;
1169 }
1170 if (reqSections & TBuffer3D::kRawSizes) {
1172 Int_t nbPnts = 4*n;
1173 Int_t nbSegs = 8*n;
1174 Int_t nbPols = 4*n;
1175 if (!HasRmin()) {
1176 nbPnts = 2*(n+1);
1177 nbSegs = 5*n;
1178 nbPols = 3*n;
1179 }
1180 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
1182 }
1183 }
1184 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
1185 SetPoints(buffer.fPnts);
1186 if (!buffer.fLocalFrame) {
1187 TransformPoints(buffer.fPnts, buffer.NbPnts());
1188 }
1189 SetSegsAndPols(buffer);
1191 }
1192
1193 return buffer;
1194}
1195
1196////////////////////////////////////////////////////////////////////////////////
1197/// Check the inside status for each of the points in the array.
1198/// Input: Array of point coordinates + vector size
1199/// Output: Array of Booleans for the inside of each point
1200
1201void TGeoTube::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
1202{
1203 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
1204}
1205
1206////////////////////////////////////////////////////////////////////////////////
1207/// Compute the normal for an array o points so that norm.dot.dir is positive
1208/// Input: Arrays of point coordinates and directions + vector size
1209/// Output: Array of normal directions
1210
1211void TGeoTube::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
1212{
1213 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
1214}
1215
1216////////////////////////////////////////////////////////////////////////////////
1217/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1218
1219void TGeoTube::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1220{
1221 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
1222}
1223
1224////////////////////////////////////////////////////////////////////////////////
1225/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1226
1227void TGeoTube::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1228{
1229 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
1230}
1231
1232////////////////////////////////////////////////////////////////////////////////
1233/// Compute safe distance from each of the points in the input array.
1234/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1235/// Output: Safety values
1236
1237void TGeoTube::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
1238{
1239 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
1240}
1241
1243
1244////////////////////////////////////////////////////////////////////////////////
1245/// Default constructor
1246
1248 :TGeoTube(),
1249 fPhi1(0.), fPhi2(0.), fS1(0.), fC1(0.), fS2(0.), fC2(0.), fSm(0.), fCm(0.), fCdfi(0.)
1250{
1252}
1253
1254////////////////////////////////////////////////////////////////////////////////
1255/// Default constructor specifying minimum and maximum radius.
1256/// The segment will be from phiStart to phiEnd expressed in degree.
1257
1259 Double_t phiStart, Double_t phiEnd)
1260 :TGeoTube(rmin, rmax, dz),
1261 fPhi1(0.), fPhi2(0.), fS1(0.), fC1(0.), fS2(0.), fC2(0.), fSm(0.), fCm(0.), fCdfi(0.)
1262{
1264 SetTubsDimensions(rmin, rmax, dz, phiStart, phiEnd);
1265 ComputeBBox();
1266}
1267
1268////////////////////////////////////////////////////////////////////////////////
1269/// Default constructor specifying minimum and maximum radius
1270/// The segment will be from phiStart to phiEnd expressed in degree.
1271
1273 Double_t phiStart, Double_t phiEnd)
1274 :TGeoTube(name, rmin, rmax, dz)
1275{
1277 SetTubsDimensions(rmin, rmax, dz, phiStart, phiEnd);
1278 ComputeBBox();
1279}
1280
1281////////////////////////////////////////////////////////////////////////////////
1282/// Default constructor specifying minimum and maximum radius
1283/// - param[0] = Rmin
1284/// - param[1] = Rmax
1285/// - param[2] = dz
1286/// - param[3] = phi1
1287/// - param[4] = phi2
1288
1290 :TGeoTube(0, 0, 0)
1291{
1293 SetDimensions(param);
1294 ComputeBBox();
1295}
1296
1297////////////////////////////////////////////////////////////////////////////////
1298/// destructor
1299
1301{
1302}
1303
1304////////////////////////////////////////////////////////////////////////////////
1305/// Function called after streaming an object of this class.
1306
1308{
1310}
1311
1312////////////////////////////////////////////////////////////////////////////////
1313/// Init frequently used trigonometric values
1314
1316{
1319 fC1 = TMath::Cos(phi1);
1320 fS1 = TMath::Sin(phi1);
1321 fC2 = TMath::Cos(phi2);
1322 fS2 = TMath::Sin(phi2);
1323 Double_t fio = 0.5*(phi1+phi2);
1324 fCm = TMath::Cos(fio);
1325 fSm = TMath::Sin(fio);
1326 Double_t dfi = 0.5*(phi2-phi1);
1327 fCdfi = TMath::Cos(dfi);
1328}
1329
1330////////////////////////////////////////////////////////////////////////////////
1331/// Computes capacity of the shape in [length^3]
1332
1334{
1336}
1337
1338////////////////////////////////////////////////////////////////////////////////
1339/// Computes capacity of the shape in [length^3]
1340
1342{
1343 Double_t capacity = TMath::Abs(phiEnd-phiStart)*TMath::DegToRad()*(rmax*rmax-rmin*rmin)*dz;
1344 return capacity;
1345}
1346
1347////////////////////////////////////////////////////////////////////////////////
1348/// compute bounding box of the tube segment
1349
1351{
1352 Double_t xc[4];
1353 Double_t yc[4];
1354 xc[0] = fRmax*fC1;
1355 yc[0] = fRmax*fS1;
1356 xc[1] = fRmax*fC2;
1357 yc[1] = fRmax*fS2;
1358 xc[2] = fRmin*fC1;
1359 yc[2] = fRmin*fS1;
1360 xc[3] = fRmin*fC2;
1361 yc[3] = fRmin*fS2;
1362
1363 Double_t xmin = xc[TMath::LocMin(4, &xc[0])];
1364 Double_t xmax = xc[TMath::LocMax(4, &xc[0])];
1365 Double_t ymin = yc[TMath::LocMin(4, &yc[0])];
1366 Double_t ymax = yc[TMath::LocMax(4, &yc[0])];
1367
1368 Double_t dp = fPhi2-fPhi1;
1369 if (dp<0) dp+=360;
1370 Double_t ddp = -fPhi1;
1371 if (ddp<0) ddp+= 360;
1372 if (ddp>360) ddp-=360;
1373 if (ddp<=dp) xmax = fRmax;
1374 ddp = 90-fPhi1;
1375 if (ddp<0) ddp+= 360;
1376 if (ddp>360) ddp-=360;
1377 if (ddp<=dp) ymax = fRmax;
1378 ddp = 180-fPhi1;
1379 if (ddp<0) ddp+= 360;
1380 if (ddp>360) ddp-=360;
1381 if (ddp<=dp) xmin = -fRmax;
1382 ddp = 270-fPhi1;
1383 if (ddp<0) ddp+= 360;
1384 if (ddp>360) ddp-=360;
1385 if (ddp<=dp) ymin = -fRmax;
1386 fOrigin[0] = (xmax+xmin)/2;
1387 fOrigin[1] = (ymax+ymin)/2;
1388 fOrigin[2] = 0;
1389 fDX = (xmax-xmin)/2;
1390 fDY = (ymax-ymin)/2;
1391 fDZ = fDz;
1392}
1393
1394////////////////////////////////////////////////////////////////////////////////
1395/// Compute normal to closest surface from POINT.
1396
1397void TGeoTubeSeg::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
1398{
1399 Double_t saf[3];
1400 Double_t rsq = point[0]*point[0]+point[1]*point[1];
1401 Double_t r = TMath::Sqrt(rsq);
1402 saf[0] = TMath::Abs(fDz-TMath::Abs(point[2]));
1403 saf[1] = (fRmin>1E-10)?TMath::Abs(r-fRmin):TGeoShape::Big();
1404 saf[2] = TMath::Abs(fRmax-r);
1405 Int_t i = TMath::LocMin(3,saf);
1406 if (((fPhi2-fPhi1)<360.) && TGeoShape::IsCloseToPhi(saf[i], point,fC1,fS1,fC2,fS2)) {
1407 TGeoShape::NormalPhi(point,dir,norm,fC1,fS1,fC2,fS2);
1408 return;
1409 }
1410 if (i==0) {
1411 norm[0] = norm[1] = 0.;
1412 norm[2] = TMath::Sign(1.,dir[2]);
1413 return;
1414 };
1415 norm[2] = 0;
1416 Double_t phi = TMath::ATan2(point[1], point[0]);
1417 norm[0] = TMath::Cos(phi);
1418 norm[1] = TMath::Sin(phi);
1419 if (norm[0]*dir[0]+norm[1]*dir[1]<0) {
1420 norm[0] = -norm[0];
1421 norm[1] = -norm[1];
1422 }
1423}
1424
1425////////////////////////////////////////////////////////////////////////////////
1426/// Compute normal to closest surface from POINT.
1427
1428void TGeoTubeSeg::ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm,
1429 Double_t rmin, Double_t rmax, Double_t /*dz*/,
1431{
1432 Double_t saf[2];
1433 Double_t rsq = point[0]*point[0]+point[1]*point[1];
1434 Double_t r = TMath::Sqrt(rsq);
1435 saf[0] = (rmin>1E-10)?TMath::Abs(r-rmin):TGeoShape::Big();
1436 saf[1] = TMath::Abs(rmax-r);
1437 Int_t i = TMath::LocMin(2,saf);
1438 if (TGeoShape::IsCloseToPhi(saf[i], point,c1,s1,c2,s2)) {
1439 TGeoShape::NormalPhi(point,dir,norm,c1,s1,c2,s2);
1440 return;
1441 }
1442 norm[2] = 0;
1443 Double_t phi = TMath::ATan2(point[1], point[0]);
1444 norm[0] = TMath::Cos(phi);
1445 norm[1] = TMath::Sin(phi);
1446 if (norm[0]*dir[0]+norm[1]*dir[1]<0) {
1447 norm[0] = -norm[0];
1448 norm[1] = -norm[1];
1449 }
1450}
1451
1452////////////////////////////////////////////////////////////////////////////////
1453/// test if point is inside this tube segment
1454/// first check if point is inside the tube
1455
1457{
1458 if (!TGeoTube::Contains(point)) return kFALSE;
1459 return IsInPhiRange(point, fPhi1, fPhi2);
1460}
1461
1462////////////////////////////////////////////////////////////////////////////////
1463/// compute closest distance from point px,py to each corner
1464
1466{
1468 const Int_t numPoints = 4*n;
1469 return ShapeDistancetoPrimitive(numPoints, px, py);
1470}
1471
1472////////////////////////////////////////////////////////////////////////////////
1473/// Compute distance from inside point to surface of the tube segment (static)
1474/// Boundary safe algorithm.
1475/// Do Z
1476
1479{
1480 Double_t stube = TGeoTube::DistFromInsideS(point,dir,rmin,rmax,dz);
1481 if (stube<=0) return 0.0;
1482 Double_t rsq = point[0]*point[0]+point[1]*point[1];
1483 Double_t r = TMath::Sqrt(rsq);
1484 Double_t cpsi=point[0]*cm+point[1]*sm;
1485 if (cpsi>r*cdfi+TGeoShape::Tolerance()) {
1486 Double_t sfmin = TGeoShape::DistToPhiMin(point, dir, s1, c1, s2, c2, sm, cm);
1487 return TMath::Min(stube,sfmin);
1488 }
1489 // Point on the phi boundary or outside
1490 // which one: phi1 or phi2
1491 Double_t ddotn, xi, yi;
1492 if (TMath::Abs(point[1]-s1*r) < TMath::Abs(point[1]-s2*r)) {
1493 ddotn = s1*dir[0]-c1*dir[1];
1494 if (ddotn>=0) return 0.0;
1495 ddotn = -s2*dir[0]+c2*dir[1];
1496 if (ddotn<=0) return stube;
1497 Double_t sfmin = s2*point[0]-c2*point[1];
1498 if (sfmin<=0) return stube;
1499 sfmin /= ddotn;
1500 if (sfmin >= stube) return stube;
1501 xi = point[0]+sfmin*dir[0];
1502 yi = point[1]+sfmin*dir[1];
1503 if (yi*cm-xi*sm<0) return stube;
1504 return sfmin;
1505 }
1506 ddotn = -s2*dir[0]+c2*dir[1];
1507 if (ddotn>=0) return 0.0;
1508 ddotn = s1*dir[0]-c1*dir[1];
1509 if (ddotn<=0) return stube;
1510 Double_t sfmin = -s1*point[0]+c1*point[1];
1511 if (sfmin<=0) return stube;
1512 sfmin /= ddotn;
1513 if (sfmin >= stube) return stube;
1514 xi = point[0]+sfmin*dir[0];
1515 yi = point[1]+sfmin*dir[1];
1516 if (yi*cm-xi*sm>0) return stube;
1517 return sfmin;
1518}
1519
1520////////////////////////////////////////////////////////////////////////////////
1521/// Compute distance from inside point to surface of the tube segment
1522/// Boundary safe algorithm.
1523
1524Double_t TGeoTubeSeg::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
1525{
1526 if (iact<3 && safe) {
1527 *safe = SafetyS(point, kTRUE, fRmin, fRmax, fDz, fPhi1, fPhi2);
1528 if (iact==0) return TGeoShape::Big();
1529 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
1530 }
1531 if ((fPhi2-fPhi1)>=360.) return TGeoTube::DistFromInsideS(point,dir,fRmin,fRmax,fDz);
1532
1533 // compute distance to surface
1535}
1536
1537////////////////////////////////////////////////////////////////////////////////
1538/// Static method to compute distance to arbitrary tube segment from outside point
1539/// Boundary safe algorithm.
1540
1543 Double_t cm, Double_t sm, Double_t cdfi)
1544{
1545 Double_t r2, cpsi, s;
1546 // check Z planes
1547 Double_t xi, yi, zi;
1548 zi = dz - TMath::Abs(point[2]);
1549 Double_t rmaxsq = rmax*rmax;
1550 Double_t rminsq = rmin*rmin;
1551 Double_t snxt=TGeoShape::Big();
1552 Bool_t in = kFALSE;
1553 Bool_t inz = (zi<0)?kFALSE:kTRUE;
1554 if (!inz) {
1555 if (point[2]*dir[2]>=0) return TGeoShape::Big();
1556 s = -zi/TMath::Abs(dir[2]);
1557 xi = point[0]+s*dir[0];
1558 yi = point[1]+s*dir[1];
1559 r2=xi*xi+yi*yi;
1560 if ((rminsq<=r2) && (r2<=rmaxsq)) {
1561 cpsi=(xi*cm+yi*sm)/TMath::Sqrt(r2);
1562 if (cpsi>=cdfi) return s;
1563 }
1564 }
1565
1566 // check outer cyl. surface
1567 Double_t rsq = point[0]*point[0]+point[1]*point[1];
1568 Double_t r = TMath::Sqrt(rsq);
1569 Double_t nsq=dir[0]*dir[0]+dir[1]*dir[1];
1570 Double_t rdotn=point[0]*dir[0]+point[1]*dir[1];
1571 Double_t b,d;
1572 Bool_t inrmax = kFALSE;
1573 Bool_t inrmin = kFALSE;
1574 Bool_t inphi = kFALSE;
1575 if (rsq<=rmaxsq+TGeoShape::Tolerance()) inrmax = kTRUE;
1576 if (rsq>=rminsq-TGeoShape::Tolerance()) inrmin = kTRUE;
1577 cpsi=point[0]*cm+point[1]*sm;
1578 if (cpsi>r*cdfi-TGeoShape::Tolerance()) inphi = kTRUE;
1579 in = inz & inrmin & inrmax & inphi;
1580 // If inside, we are most likely on a boundary within machine precision.
1581 if (in) {
1582 Bool_t checkout = kFALSE;
1583 Double_t safphi = (cpsi-r*cdfi)*TMath::Sqrt(1.-cdfi*cdfi);
1584// Double_t sch, cch;
1585 // check if on Z boundaries
1586 if (zi<rmax-r) {
1587 if (TGeoShape::IsSameWithinTolerance(rmin,0) || (zi<r-rmin)) {
1588 if (zi<safphi) {
1589 if (point[2]*dir[2]<0) return 0.0;
1590 return TGeoShape::Big();
1591 }
1592 }
1593 }
1594 if ((rmaxsq-rsq) < (rsq-rminsq)) checkout = kTRUE;
1595 // check if on Rmax boundary
1596 if (checkout && (rmax-r<safphi)) {
1597 if (rdotn>=0) return TGeoShape::Big();
1598 return 0.0;
1599 }
1600 if (TMath::Abs(nsq)<TGeoShape::Tolerance()) return TGeoShape::Big();
1601 // check if on phi boundary
1602 if (TGeoShape::IsSameWithinTolerance(rmin,0) || (safphi<r-rmin)) {
1603 // We may cross again a phi of rmin boundary
1604 // check first if we are on phi1 or phi2
1605 Double_t un;
1606 if (point[0]*c1 + point[1]*s1 > point[0]*c2 + point[1]*s2) {
1607 un = dir[0]*s1-dir[1]*c1;
1608 if (un < 0) return 0.0;
1609 if (cdfi>=0) return TGeoShape::Big();
1610 un = -dir[0]*s2+dir[1]*c2;
1611 if (un<0) {
1612 s = -point[0]*s2+point[1]*c2;
1613 if (s>0) {
1614 s /= (-un);
1615 zi = point[2]+s*dir[2];
1616 if (TMath::Abs(zi)<=dz) {
1617 xi = point[0]+s*dir[0];
1618 yi = point[1]+s*dir[1];
1619 r2=xi*xi+yi*yi;
1620 if ((rminsq<=r2) && (r2<=rmaxsq)) {
1621 if ((yi*cm-xi*sm)>0) return s;
1622 }
1623 }
1624 }
1625 }
1626 } else {
1627 un = -dir[0]*s2+dir[1]*c2;
1628 if (un < 0) return 0.0;
1629 if (cdfi>=0) return TGeoShape::Big();
1630 un = dir[0]*s1-dir[1]*c1;
1631 if (un<0) {
1632 s = point[0]*s1-point[1]*c1;
1633 if (s>0) {
1634 s /= (-un);
1635 zi = point[2]+s*dir[2];
1636 if (TMath::Abs(zi)<=dz) {
1637 xi = point[0]+s*dir[0];
1638 yi = point[1]+s*dir[1];
1639 r2=xi*xi+yi*yi;
1640 if ((rminsq<=r2) && (r2<=rmaxsq)) {
1641 if ((yi*cm-xi*sm)<0) return s;
1642 }
1643 }
1644 }
1645 }
1646 }
1647 // We may also cross rmin, (+) solution
1648 if (rdotn>=0) return TGeoShape::Big();
1649 if (cdfi>=0) return TGeoShape::Big();
1650 DistToTube(rsq, nsq, rdotn, rmin, b, d);
1651 if (d>0) {
1652 s=-b+d;
1653 if (s>0) {
1654 zi=point[2]+s*dir[2];
1655 if (TMath::Abs(zi)<=dz) {
1656 xi=point[0]+s*dir[0];
1657 yi=point[1]+s*dir[1];
1658 if ((xi*cm+yi*sm) >= rmin*cdfi) return s;
1659 }
1660 }
1661 }
1662 return TGeoShape::Big();
1663 }
1664 // we are on rmin boundary: we may cross again rmin or a phi facette
1665 if (rdotn>=0) return 0.0;
1666 DistToTube(rsq, nsq, rdotn, rmin, b, d);
1667 if (d>0) {
1668 s=-b+d;
1669 if (s>0) {
1670 zi=point[2]+s*dir[2];
1671 if (TMath::Abs(zi)<=dz) {
1672 // now check phi range
1673 xi=point[0]+s*dir[0];
1674 yi=point[1]+s*dir[1];
1675 if ((xi*cm+yi*sm) >= rmin*cdfi) return s;
1676 // now we really have to check any phi crossing
1677 Double_t un=-dir[0]*s1+dir[1]*c1;
1678 if (un > 0) {
1679 s=point[0]*s1-point[1]*c1;
1680 if (s>=0) {
1681 s /= un;
1682 zi=point[2]+s*dir[2];
1683 if (TMath::Abs(zi)<=dz) {
1684 xi=point[0]+s*dir[0];
1685 yi=point[1]+s*dir[1];
1686 r2=xi*xi+yi*yi;
1687 if ((rminsq<=r2) && (r2<=rmaxsq)) {
1688 if ((yi*cm-xi*sm)<=0) {
1689 if (s<snxt) snxt=s;
1690 }
1691 }
1692 }
1693 }
1694 }
1695 un=dir[0]*s2-dir[1]*c2;
1696 if (un > 0) {
1697 s=(point[1]*c2-point[0]*s2)/un;
1698 if (s>=0 && s<snxt) {
1699 zi=point[2]+s*dir[2];
1700 if (TMath::Abs(zi)<=dz) {
1701 xi=point[0]+s*dir[0];
1702 yi=point[1]+s*dir[1];
1703 r2=xi*xi+yi*yi;
1704 if ((rminsq<=r2) && (r2<=rmaxsq)) {
1705 if ((yi*cm-xi*sm)>=0) {
1706 return s;
1707 }
1708 }
1709 }
1710 }
1711 }
1712 return snxt;
1713 }
1714 }
1715 }
1716 return TGeoShape::Big();
1717 }
1718 // only r>rmax has to be considered
1719 if (TMath::Abs(nsq)<TGeoShape::Tolerance()) return TGeoShape::Big();
1720 if (rsq>=rmax*rmax) {
1721 if (rdotn>=0) return TGeoShape::Big();
1722 TGeoTube::DistToTube(rsq, nsq, rdotn, rmax, b, d);
1723 if (d>0) {
1724 s=-b-d;
1725 if (s>0) {
1726 zi=point[2]+s*dir[2];
1727 if (TMath::Abs(zi)<=dz) {
1728 xi=point[0]+s*dir[0];
1729 yi=point[1]+s*dir[1];
1730 cpsi = xi*cm+yi*sm;
1731 if (cpsi>=rmax*cdfi) return s;
1732 }
1733 }
1734 }
1735 }
1736 // check inner cylinder
1737 if (rmin>0) {
1738 TGeoTube::DistToTube(rsq, nsq, rdotn, rmin, b, d);
1739 if (d>0) {
1740 s=-b+d;
1741 if (s>0) {
1742 zi=point[2]+s*dir[2];
1743 if (TMath::Abs(zi)<=dz) {
1744 xi=point[0]+s*dir[0];
1745 yi=point[1]+s*dir[1];
1746 cpsi = xi*cm+yi*sm;
1747 if (cpsi>=rmin*cdfi) snxt=s;
1748 }
1749 }
1750 }
1751 }
1752 // check phi planes
1753 Double_t un=-dir[0]*s1+dir[1]*c1;
1754 if (un > 0) {
1755 s=point[0]*s1-point[1]*c1;
1756 if (s>=0) {
1757 s /= un;
1758 zi=point[2]+s*dir[2];
1759 if (TMath::Abs(zi)<=dz) {
1760 xi=point[0]+s*dir[0];
1761 yi=point[1]+s*dir[1];
1762 r2=xi*xi+yi*yi;
1763 if ((rminsq<=r2) && (r2<=rmaxsq)) {
1764 if ((yi*cm-xi*sm)<=0) {
1765 if (s<snxt) snxt=s;
1766 }
1767 }
1768 }
1769 }
1770 }
1771 un=dir[0]*s2-dir[1]*c2;
1772 if (un > 0) {
1773 s=point[1]*c2-point[0]*s2;
1774 if (s>=0) {
1775 s /= un;
1776 zi=point[2]+s*dir[2];
1777 if (TMath::Abs(zi)<=dz) {
1778 xi=point[0]+s*dir[0];
1779 yi=point[1]+s*dir[1];
1780 r2=xi*xi+yi*yi;
1781 if ((rminsq<=r2) && (r2<=rmaxsq)) {
1782 if ((yi*cm-xi*sm)>=0) {
1783 if (s<snxt) snxt=s;
1784 }
1785 }
1786 }
1787 }
1788 }
1789 return snxt;
1790}
1791
1792////////////////////////////////////////////////////////////////////////////////
1793/// compute distance from outside point to surface of the tube segment
1794/// fist localize point w.r.t tube
1795
1796Double_t TGeoTubeSeg::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
1797{
1798 if (iact<3 && safe) {
1799 *safe = SafetyS(point, kFALSE, fRmin, fRmax, fDz, fPhi1, fPhi2);
1800 if (iact==0) return TGeoShape::Big();
1801 if ((iact==1) && (step<=*safe)) return TGeoShape::Big();
1802 }
1803// Check if the bounding box is crossed within the requested distance
1804 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
1805 if (sdist>=step) return TGeoShape::Big();
1806 if ((fPhi2-fPhi1)>=360.) return TGeoTube::DistFromOutsideS(point,dir,fRmin,fRmax,fDz);
1807
1808 // find distance to shape
1809 return TGeoTubeSeg::DistFromOutsideS(point, dir, fRmin, fRmax, fDz, fC1, fS1, fC2, fS2, fCm, fSm, fCdfi);
1810}
1811
1812////////////////////////////////////////////////////////////////////////////////
1813/// Divide this tube segment shape belonging to volume "voldiv" into ndiv volumes
1814/// called divname, from start position with the given step. Returns pointer
1815/// to created division cell volume in case of Z divisions. For radialdivision
1816/// creates all volumes with different shapes and returns pointer to volume that
1817/// was divided. In case a wrong division axis is supplied, returns pointer to
1818/// volume that was divided.
1819
1820TGeoVolume *TGeoTubeSeg::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
1821 Double_t start, Double_t step)
1822{
1823 TGeoShape *shape; //--- shape to be created
1824 TGeoVolume *vol; //--- division volume to be created
1825 TGeoVolumeMulti *vmulti; //--- generic divided volume
1826 TGeoPatternFinder *finder; //--- finder to be attached
1827 TString opt = ""; //--- option to be attached
1828 Double_t dphi;
1829 Int_t id;
1830 Double_t end = start+ndiv*step;
1831 switch (iaxis) {
1832 case 1: //--- R division
1833 finder = new TGeoPatternCylR(voldiv, ndiv, start, end);
1834 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
1835 voldiv->SetFinder(finder);
1836 finder->SetDivIndex(voldiv->GetNdaughters());
1837 for (id=0; id<ndiv; id++) {
1838 shape = new TGeoTubeSeg(start+id*step, start+(id+1)*step, fDz, fPhi1, fPhi2);
1839 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
1840 vmulti->AddVolume(vol);
1841 opt = "R";
1842 voldiv->AddNodeOffset(vol, id, 0, opt.Data());
1843 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
1844 }
1845 return vmulti;
1846 case 2: //--- Phi division
1847 dphi = fPhi2-fPhi1;
1848 if (dphi<0) dphi+=360.;
1849 if (step<=0) {step=dphi/ndiv; start=fPhi1; end=fPhi2;}
1850 finder = new TGeoPatternCylPhi(voldiv, ndiv, start, end);
1851 voldiv->SetFinder(finder);
1852 finder->SetDivIndex(voldiv->GetNdaughters());
1853 shape = new TGeoTubeSeg(fRmin, fRmax, fDz, -step/2, step/2);
1854 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
1855 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
1856 vmulti->AddVolume(vol);
1857 opt = "Phi";
1858 for (id=0; id<ndiv; id++) {
1859 voldiv->AddNodeOffset(vol, id, start+id*step+step/2, opt.Data());
1860 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
1861 }
1862 return vmulti;
1863 case 3: //--- Z division
1864 finder = new TGeoPatternZ(voldiv, ndiv, start, end);
1865 voldiv->SetFinder(finder);
1866 finder->SetDivIndex(voldiv->GetNdaughters());
1867 shape = new TGeoTubeSeg(fRmin, fRmax, step/2, fPhi1, fPhi2);
1868 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
1869 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
1870 vmulti->AddVolume(vol);
1871 opt = "Z";
1872 for (id=0; id<ndiv; id++) {
1873 voldiv->AddNodeOffset(vol, id, start+step/2+id*step, opt.Data());
1874 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
1875 }
1876 return vmulti;
1877 default:
1878 Error("Divide", "In shape %s wrong axis type for division", GetName());
1879 return 0;
1880 }
1881}
1882
1883////////////////////////////////////////////////////////////////////////////////
1884/// Get range of shape for a given axis.
1885
1887{
1888 xlo = 0;
1889 xhi = 0;
1890 Double_t dx = 0;
1891 switch (iaxis) {
1892 case 1:
1893 xlo = fRmin;
1894 xhi = fRmax;
1895 dx = xhi-xlo;
1896 return dx;
1897 case 2:
1898 xlo = fPhi1;
1899 xhi = fPhi2;
1900 dx = xhi-xlo;
1901 return dx;
1902 case 3:
1903 xlo = -fDz;
1904 xhi = fDz;
1905 dx = xhi-xlo;
1906 return dx;
1907 }
1908 return dx;
1909}
1910
1911////////////////////////////////////////////////////////////////////////////////
1912/// Fill vector param[4] with the bounding cylinder parameters. The order
1913/// is the following : Rmin, Rmax, Phi1, Phi2
1914
1916{
1917 param[0] = fRmin;
1918 param[0] *= param[0];
1919 param[1] = fRmax;
1920 param[1] *= param[1];
1921 param[2] = fPhi1;
1922 param[3] = fPhi2;
1923}
1924
1925////////////////////////////////////////////////////////////////////////////////
1926/// in case shape has some negative parameters, these has to be computed
1927/// in order to fit the mother
1928
1930{
1931 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
1932 if (!mother->TestShapeBit(kGeoTube)) {
1933 Error("GetMakeRuntimeShape", "Invalid mother for shape %s", GetName());
1934 return 0;
1935 }
1936 Double_t rmin, rmax, dz;
1937 rmin = fRmin;
1938 rmax = fRmax;
1939 dz = fDz;
1940 if (fDz<0) dz=((TGeoTube*)mother)->GetDz();
1941 if (fRmin<0)
1942 rmin = ((TGeoTube*)mother)->GetRmin();
1943 if ((fRmax<0) || (fRmax<=fRmin))
1944 rmax = ((TGeoTube*)mother)->GetRmax();
1945
1946 return (new TGeoTubeSeg(GetName(),rmin, rmax, dz, fPhi1, fPhi2));
1947}
1948
1949////////////////////////////////////////////////////////////////////////////////
1950/// print shape parameters
1951
1953{
1954 printf("*** Shape %s: TGeoTubeSeg ***\n", GetName());
1955 printf(" Rmin = %11.5f\n", fRmin);
1956 printf(" Rmax = %11.5f\n", fRmax);
1957 printf(" dz = %11.5f\n", fDz);
1958 printf(" phi1 = %11.5f\n", fPhi1);
1959 printf(" phi2 = %11.5f\n", fPhi2);
1960 printf(" Bounding box:\n");
1962}
1963
1964////////////////////////////////////////////////////////////////////////////////
1965/// Creates a TBuffer3D describing *this* shape.
1966/// Coordinates are in local reference frame.
1967
1969{
1971 Int_t nbPnts = 4*n;
1972 Int_t nbSegs = 2*nbPnts;
1973 Int_t nbPols = nbPnts-2;
1974
1976 nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols);
1977 if (buff)
1978 {
1979 SetPoints(buff->fPnts);
1980 SetSegsAndPols(*buff);
1981 }
1982
1983 return buff;
1984}
1985
1986////////////////////////////////////////////////////////////////////////////////
1987/// Fill TBuffer3D structure for segments and polygons.
1988
1990{
1991 Int_t i, j;
1993 Int_t c = GetBasicColor();
1994
1995 memset(buff.fSegs, 0, buff.NbSegs()*3*sizeof(Int_t));
1996 for (i = 0; i < 4; i++) {
1997 for (j = 1; j < n; j++) {
1998 buff.fSegs[(i*n+j-1)*3 ] = c;
1999 buff.fSegs[(i*n+j-1)*3+1] = i*n+j-1;
2000 buff.fSegs[(i*n+j-1)*3+2] = i*n+j;
2001 }
2002 }
2003 for (i = 4; i < 6; i++) {
2004 for (j = 0; j < n; j++) {
2005 buff.fSegs[(i*n+j)*3 ] = c+1;
2006 buff.fSegs[(i*n+j)*3+1] = (i-4)*n+j;
2007 buff.fSegs[(i*n+j)*3+2] = (i-2)*n+j;
2008 }
2009 }
2010 for (i = 6; i < 8; i++) {
2011 for (j = 0; j < n; j++) {
2012 buff.fSegs[(i*n+j)*3 ] = c;
2013 buff.fSegs[(i*n+j)*3+1] = 2*(i-6)*n+j;
2014 buff.fSegs[(i*n+j)*3+2] = (2*(i-6)+1)*n+j;
2015 }
2016 }
2017
2018 Int_t indx = 0;
2019 memset(buff.fPols, 0, buff.NbPols()*6*sizeof(Int_t));
2020 i = 0;
2021 for (j = 0; j < n-1; j++) {
2022 buff.fPols[indx++] = c;
2023 buff.fPols[indx++] = 4;
2024 buff.fPols[indx++] = (4+i)*n+j+1;
2025 buff.fPols[indx++] = (2+i)*n+j;
2026 buff.fPols[indx++] = (4+i)*n+j;
2027 buff.fPols[indx++] = i*n+j;
2028 }
2029 i = 1;
2030 for (j = 0; j < n-1; j++) {
2031 buff.fPols[indx++] = c;
2032 buff.fPols[indx++] = 4;
2033 buff.fPols[indx++] = i*n+j;
2034 buff.fPols[indx++] = (4+i)*n+j;
2035 buff.fPols[indx++] = (2+i)*n+j;
2036 buff.fPols[indx++] = (4+i)*n+j+1;
2037 }
2038 i = 2;
2039 for (j = 0; j < n-1; j++) {
2040 buff.fPols[indx++] = c+i;
2041 buff.fPols[indx++] = 4;
2042 buff.fPols[indx++] = (i-2)*2*n+j;
2043 buff.fPols[indx++] = (4+i)*n+j;
2044 buff.fPols[indx++] = ((i-2)*2+1)*n+j;
2045 buff.fPols[indx++] = (4+i)*n+j+1;
2046 }
2047 i = 3;
2048 for (j = 0; j < n-1; j++) {
2049 buff.fPols[indx++] = c+i;
2050 buff.fPols[indx++] = 4;
2051 buff.fPols[indx++] = (4+i)*n+j+1;
2052 buff.fPols[indx++] = ((i-2)*2+1)*n+j;
2053 buff.fPols[indx++] = (4+i)*n+j;
2054 buff.fPols[indx++] = (i-2)*2*n+j;
2055 }
2056 buff.fPols[indx++] = c+2;
2057 buff.fPols[indx++] = 4;
2058 buff.fPols[indx++] = 6*n;
2059 buff.fPols[indx++] = 4*n;
2060 buff.fPols[indx++] = 7*n;
2061 buff.fPols[indx++] = 5*n;
2062 buff.fPols[indx++] = c+2;
2063 buff.fPols[indx++] = 4;
2064 buff.fPols[indx++] = 6*n-1;
2065 buff.fPols[indx++] = 8*n-1;
2066 buff.fPols[indx++] = 5*n-1;
2067 buff.fPols[indx++] = 7*n-1;
2068}
2069
2070////////////////////////////////////////////////////////////////////////////////
2071/// computes the closest distance from given point InitTrigonometry();to this shape, according
2072/// to option. The matching point on the shape is stored in spoint.
2073
2075{
2076 Double_t saf[3];
2077 Double_t rsq = point[0]*point[0]+point[1]*point[1];
2078 Double_t r = TMath::Sqrt(rsq);
2079 if (in) {
2080 saf[0] = fDz-TMath::Abs(point[2]);
2081 saf[1] = r-fRmin;
2082 saf[2] = fRmax-r;
2083 Double_t safe = saf[TMath::LocMin(3,saf)];
2084 if ((fPhi2-fPhi1)>=360.) return safe;
2085 Double_t safphi = TGeoShape::SafetyPhi(point,in,fPhi1,fPhi2);
2086 return TMath::Min(safe, safphi);
2087 }
2088 // Point expected to be outside
2089 Bool_t inphi = kFALSE;
2090 Double_t cpsi=point[0]*fCm+point[1]*fSm;
2091 saf[0] = TMath::Abs(point[2])-fDz;
2092 if (cpsi>r*fCdfi-TGeoShape::Tolerance()) inphi = kTRUE;
2093 if (inphi) {
2094 saf[1] = fRmin-r;
2095 saf[2] = r-fRmax;
2096 Double_t safe = saf[TMath::LocMax(3,saf)];
2097 safe = TMath::Max(0., safe);
2098 return safe;
2099 }
2100 // Point outside the phi range
2101 // Compute projected radius of the (r,phi) position vector onto
2102 // phi1 and phi2 edges and take the maximum for choosing the side.
2103 Double_t rproj = TMath::Max(point[0]*fC1+point[1]*fS1, point[0]*fC2+point[1]*fS2);
2104 saf[1] = fRmin-rproj;
2105 saf[2] = rproj-fRmax;
2106 Double_t safe = TMath::Max(saf[1], saf[2]);
2107 if ((fPhi2-fPhi1)>=360.) return TMath::Max(safe,saf[0]);
2108 if (safe>0) {
2109 // rproj not within (rmin,rmax) - > no need to calculate safphi
2110 safe = TMath::Sqrt(rsq-rproj*rproj+safe*safe);
2111 return (saf[0]<0) ? safe : TMath::Sqrt(safe*safe+saf[0]*saf[0]);
2112 }
2113 Double_t safphi = TGeoShape::SafetyPhi(point,in,fPhi1,fPhi2);
2114 return (saf[0]<0) ? safphi : TMath::Sqrt(saf[0]*saf[0]+safphi*safphi);
2115}
2116
2117////////////////////////////////////////////////////////////////////////////////
2118/// Static method to compute the closest distance from given point to this shape.
2119
2121 Double_t phi1d, Double_t phi2d, Int_t skipz)
2122{
2123 Double_t saf[3];
2124 Double_t rsq = point[0]*point[0]+point[1]*point[1];
2125 Double_t r = TMath::Sqrt(rsq);
2126
2127 switch (skipz) {
2128 case 1: // skip lower Z plane
2129 saf[0] = dz - point[2];
2130 break;
2131 case 2: // skip upper Z plane
2132 saf[0] = dz + point[2];
2133 break;
2134 case 3: // skip both
2135 saf[0] = TGeoShape::Big();
2136 break;
2137 default:
2138 saf[0] = dz-TMath::Abs(point[2]);
2139 }
2140
2141 if (in) {
2142 saf[1] = r-rmin;
2143 saf[2] = rmax-r;
2144 Double_t safe = saf[TMath::LocMin(3,saf)];
2145 if ((phi2d-phi1d)>=360.) return safe;
2146 Double_t safphi = TGeoShape::SafetyPhi(point,in,phi1d,phi2d);
2147 return TMath::Min(safe, safphi);
2148 }
2149 // Point expected to be outside
2150 saf[0] = -saf[0];
2151 Bool_t inphi = kFALSE;
2152 Double_t phi1 = phi1d*TMath::DegToRad();
2153 Double_t phi2 = phi2d*TMath::DegToRad();
2154
2155 Double_t fio = 0.5*(phi1+phi2);
2156 Double_t cm = TMath::Cos(fio);
2157 Double_t sm = TMath::Sin(fio);
2158 Double_t cpsi=point[0]*cm+point[1]*sm;
2159 Double_t dfi = 0.5*(phi2-phi1);
2160 Double_t cdfi = TMath::Cos(dfi);
2161 if (cpsi>r*cdfi-TGeoShape::Tolerance()) inphi = kTRUE;
2162 if (inphi) {
2163 saf[1] = rmin-r;
2164 saf[2] = r-rmax;
2165 Double_t safe = saf[TMath::LocMax(3,saf)];
2166 safe = TMath::Max(0., safe);
2167 return safe;
2168 }
2169 // Point outside the phi range
2170 // Compute projected radius of the (r,phi) position vector onto
2171 // phi1 and phi2 edges and take the maximum for choosing the side.
2172 Double_t c1 = TMath::Cos(phi1);
2173 Double_t s1 = TMath::Sin(phi1);
2174 Double_t c2 = TMath::Cos(phi2);
2175 Double_t s2 = TMath::Sin(phi2);
2176
2177 Double_t rproj = TMath::Max(point[0]*c1+point[1]*s1, point[0]*c2+point[1]*s2);
2178 saf[1] = rmin-rproj;
2179 saf[2] = rproj-rmax;
2180 Double_t safe = TMath::Max(saf[1], saf[2]);
2181 if ((phi2d-phi1d)>=360.) return TMath::Max(safe,saf[0]);
2182 if (safe>0) {
2183 // rproj not within (rmin,rmax) - > no need to calculate safphi
2184 safe = TMath::Sqrt(rsq-rproj*rproj+safe*safe);
2185 return (saf[0]<0) ? safe : TMath::Sqrt(safe*safe+saf[0]*saf[0]);
2186 }
2187 Double_t safphi = TGeoShape::SafetyPhi(point,in,phi1d,phi2d);
2188 return (saf[0]<0) ? safphi : TMath::Sqrt(saf[0]*saf[0]+safphi*safphi);
2189}
2190
2191////////////////////////////////////////////////////////////////////////////////
2192/// Save a primitive as a C++ statement(s) on output stream "out".
2193
2194void TGeoTubeSeg::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
2195{
2197 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
2198 out << " rmin = " << fRmin << ";" << std::endl;
2199 out << " rmax = " << fRmax << ";" << std::endl;
2200 out << " dz = " << fDz << ";" << std::endl;
2201 out << " phi1 = " << fPhi1 << ";" << std::endl;
2202 out << " phi2 = " << fPhi2 << ";" << std::endl;
2203 out << " TGeoShape *" << GetPointerName() << " = new TGeoTubeSeg(\"" << GetName() << "\",rmin,rmax,dz,phi1,phi2);" << std::endl;
2205}
2206
2207////////////////////////////////////////////////////////////////////////////////
2208/// Set dimensions of the tube segment.
2209/// The segment will be from phiStart to phiEnd expressed in degree.
2210
2212 Double_t phiStart, Double_t phiEnd)
2213{
2214 fRmin = rmin;
2215 fRmax = rmax;
2216 fDz = dz;
2217 fPhi1 = phiStart;
2218 if (fPhi1 < 0) fPhi1 += 360.;
2219 fPhi2 = phiEnd;
2220 while (fPhi2<=fPhi1) fPhi2+=360.;
2221 if (TGeoShape::IsSameWithinTolerance(fPhi1,fPhi2)) Fatal("SetTubsDimensions", "In shape %s invalid phi1=%g, phi2=%g\n", GetName(), fPhi1, fPhi2);
2223}
2224
2225////////////////////////////////////////////////////////////////////////////////
2226/// Set dimensions of the tube segment starting from a list.
2227
2229{
2230 Double_t rmin = param[0];
2231 Double_t rmax = param[1];
2232 Double_t dz = param[2];
2233 Double_t phi1 = param[3];
2234 Double_t phi2 = param[4];
2235 SetTubsDimensions(rmin, rmax, dz, phi1, phi2);
2236}
2237
2238////////////////////////////////////////////////////////////////////////////////
2239/// Fills array with n random points located on the line segments of the shape mesh.
2240/// The output array must be provided with a length of minimum 3*npoints. Returns
2241/// true if operation is implemented.
2242
2244{
2245 if (npoints > (npoints/2)*2) {
2246 Error("GetPointsOnSegments","Npoints must be even number");
2247 return kFALSE;
2248 }
2249 Int_t nc = (Int_t)TMath::Sqrt(0.5*npoints);
2250 Double_t dphi = (fPhi2-fPhi1)*TMath::DegToRad()/(nc-1);
2251 Double_t phi = 0;
2252 Double_t phi1 = fPhi1 * TMath::DegToRad();
2253 Int_t ntop = npoints/2 - nc*(nc-1);
2254 Double_t dz = 2*fDz/(nc-1);
2255 Double_t z = 0;
2256 Int_t icrt = 0;
2257 Int_t nphi = nc;
2258 // loop z sections
2259 for (Int_t i=0; i<nc; i++) {
2260 if (i == (nc-1)) {
2261 nphi = ntop;
2262 dphi = (fPhi2-fPhi1)*TMath::DegToRad()/(nphi-1);
2263 }
2264 z = -fDz + i*dz;
2265 // loop points on circle sections
2266 for (Int_t j=0; j<nphi; j++) {
2267 phi = phi1 + j*dphi;
2268 array[icrt++] = fRmin * TMath::Cos(phi);
2269 array[icrt++] = fRmin * TMath::Sin(phi);
2270 array[icrt++] = z;
2271 array[icrt++] = fRmax * TMath::Cos(phi);
2272 array[icrt++] = fRmax * TMath::Sin(phi);
2273 array[icrt++] = z;
2274 }
2275 }
2276 return kTRUE;
2277}
2278
2279////////////////////////////////////////////////////////////////////////////////
2280/// Create tube segment mesh points.
2281
2283{
2284 Double_t dz;
2285 Int_t j, n;
2286 Double_t phi, phi1, phi2, dphi;
2287 phi1 = fPhi1;
2288 phi2 = fPhi2;
2289 if (phi2<phi1) phi2+=360.;
2290 n = gGeoManager->GetNsegments()+1;
2291
2292 dphi = (phi2-phi1)/(n-1);
2293 dz = fDz;
2294
2295 if (points) {
2296 Int_t indx = 0;
2297
2298 for (j = 0; j < n; j++) {
2299 phi = (phi1+j*dphi)*TMath::DegToRad();
2300 points[indx+6*n] = points[indx] = fRmin * TMath::Cos(phi);
2301 indx++;
2302 points[indx+6*n] = points[indx] = fRmin * TMath::Sin(phi);
2303 indx++;
2304 points[indx+6*n] = dz;
2305 points[indx] =-dz;
2306 indx++;
2307 }
2308 for (j = 0; j < n; j++) {
2309 phi = (phi1+j*dphi)*TMath::DegToRad();
2310 points[indx+6*n] = points[indx] = fRmax * TMath::Cos(phi);
2311 indx++;
2312 points[indx+6*n] = points[indx] = fRmax * TMath::Sin(phi);
2313 indx++;
2314 points[indx+6*n]= dz;
2315 points[indx] =-dz;
2316 indx++;
2317 }
2318 }
2319}
2320
2321////////////////////////////////////////////////////////////////////////////////
2322/// Create tube segment mesh points.
2323
2325{
2326 Double_t dz;
2327 Int_t j, n;
2328 Double_t phi, phi1, phi2, dphi;
2329 phi1 = fPhi1;
2330 phi2 = fPhi2;
2331 if (phi2<phi1) phi2+=360.;
2332 n = gGeoManager->GetNsegments()+1;
2333
2334 dphi = (phi2-phi1)/(n-1);
2335 dz = fDz;
2336
2337 if (points) {
2338 Int_t indx = 0;
2339
2340 for (j = 0; j < n; j++) {
2341 phi = (phi1+j*dphi)*TMath::DegToRad();
2342 points[indx+6*n] = points[indx] = fRmin * TMath::Cos(phi);
2343 indx++;
2344 points[indx+6*n] = points[indx] = fRmin * TMath::Sin(phi);
2345 indx++;
2346 points[indx+6*n] = dz;
2347 points[indx] =-dz;
2348 indx++;
2349 }
2350 for (j = 0; j < n; j++) {
2351 phi = (phi1+j*dphi)*TMath::DegToRad();
2352 points[indx+6*n] = points[indx] = fRmax * TMath::Cos(phi);
2353 indx++;
2354 points[indx+6*n] = points[indx] = fRmax * TMath::Sin(phi);
2355 indx++;
2356 points[indx+6*n]= dz;
2357 points[indx] =-dz;
2358 indx++;
2359 }
2360 }
2361}
2362
2363////////////////////////////////////////////////////////////////////////////////
2364/// Returns numbers of vertices, segments and polygons composing the shape mesh.
2365
2366void TGeoTubeSeg::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
2367{
2369 nvert = n*4;
2370 nsegs = n*8;
2371 npols = n*4 - 2;
2372}
2373
2374////////////////////////////////////////////////////////////////////////////////
2375/// Return number of vertices of the mesh representation
2376
2378{
2380 Int_t numPoints = n*4;
2381 return numPoints;
2382}
2383
2384////////////////////////////////////////////////////////////////////////////////
2385/// fill size of this 3-D object
2386
2388{
2389}
2390
2391////////////////////////////////////////////////////////////////////////////////
2392/// Fills a static 3D buffer and returns a reference.
2393
2394const TBuffer3D & TGeoTubeSeg::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
2395{
2396 static TBuffer3DTubeSeg buffer;
2397 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
2398
2399 if (reqSections & TBuffer3D::kShapeSpecific) {
2400 // These from TBuffer3DTube / TGeoTube
2401 buffer.fRadiusInner = fRmin;
2402 buffer.fRadiusOuter = fRmax;
2403 buffer.fHalfLength = fDz;
2404 buffer.fPhiMin = fPhi1;
2405 buffer.fPhiMax = fPhi2;
2407 }
2408 if (reqSections & TBuffer3D::kRawSizes) {
2410 Int_t nbPnts = 4*n;
2411 Int_t nbSegs = 2*nbPnts;
2412 Int_t nbPols = nbPnts-2;
2413 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
2415 }
2416 }
2417 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
2418 SetPoints(buffer.fPnts);
2419 if (!buffer.fLocalFrame) {
2420 TransformPoints(buffer.fPnts, buffer.NbPnts());
2421 }
2422 SetSegsAndPols(buffer);
2424 }
2425
2426 return buffer;
2427}
2428
2429////////////////////////////////////////////////////////////////////////////////
2430/// Check the inside status for each of the points in the array.
2431/// Input: Array of point coordinates + vector size
2432/// Output: Array of Booleans for the inside of each point
2433
2434void TGeoTubeSeg::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
2435{
2436 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
2437}
2438
2439////////////////////////////////////////////////////////////////////////////////
2440/// Compute the normal for an array o points so that norm.dot.dir is positive
2441/// Input: Arrays of point coordinates and directions + vector size
2442/// Output: Array of normal directions
2443
2444void TGeoTubeSeg::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
2445{
2446 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
2447}
2448
2449////////////////////////////////////////////////////////////////////////////////
2450/// Compute distance from array of input points having directions specified by dirs. Store output in dists
2451
2452void TGeoTubeSeg::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
2453{
2454 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
2455}
2456
2457////////////////////////////////////////////////////////////////////////////////
2458/// Compute distance from array of input points having directions specified by dirs. Store output in dists
2459
2460void TGeoTubeSeg::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
2461{
2462 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
2463}
2464
2465////////////////////////////////////////////////////////////////////////////////
2466/// Compute safe distance from each of the points in the input array.
2467/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
2468/// Output: Safety values
2469
2470void TGeoTubeSeg::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
2471{
2472 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
2473}
2474
2476
2478{
2479// default ctor
2480 fNlow[0] = fNlow[1] = fNhigh[0] = fNhigh[1] = 0.;
2481 fNlow[2] = -1;
2482 fNhigh[2] = 1;
2483}
2484
2485////////////////////////////////////////////////////////////////////////////////
2486/// constructor
2487
2489 Double_t lx, Double_t ly, Double_t lz, Double_t tx, Double_t ty, Double_t tz)
2490 :TGeoTubeSeg(rmin, rmax, dz, phi1, phi2)
2491{
2492 fNlow[0] = lx;
2493 fNlow[1] = ly;
2494 fNlow[2] = lz;
2495 fNhigh[0] = tx;
2496 fNhigh[1] = ty;
2497 fNhigh[2] = tz;
2499 ComputeBBox();
2500}
2501
2502////////////////////////////////////////////////////////////////////////////////
2503/// constructor
2504
2505TGeoCtub::TGeoCtub(const char *name, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2,
2506 Double_t lx, Double_t ly, Double_t lz, Double_t tx, Double_t ty, Double_t tz)
2507 :TGeoTubeSeg(name, rmin, rmax, dz, phi1, phi2)
2508{
2509 fNlow[0] = lx;
2510 fNlow[1] = ly;
2511 fNlow[2] = lz;
2512 fNhigh[0] = tx;
2513 fNhigh[1] = ty;
2514 fNhigh[2] = tz;
2516 ComputeBBox();
2517}
2518
2519////////////////////////////////////////////////////////////////////////////////
2520/// ctor with parameters
2521
2523 :TGeoTubeSeg(0,0,0,0,0)
2524{
2525 SetCtubDimensions(params[0], params[1], params[2], params[3], params[4], params[5],
2526 params[6], params[7], params[8], params[9], params[10]);
2528}
2529
2530////////////////////////////////////////////////////////////////////////////////
2531/// destructor
2532
2534{
2535}
2536
2537////////////////////////////////////////////////////////////////////////////////
2538/// Computes capacity of the shape in [length^3]
2539
2541{
2542 Double_t capacity = TGeoTubeSeg::Capacity();
2543 return capacity;
2544}
2545
2546////////////////////////////////////////////////////////////////////////////////
2547/// compute minimum bounding box of the ctub
2548
2550{
2552 if ((fNlow[2]>-(1E-10)) || (fNhigh[2]<1E-10)) {
2553 Error("ComputeBBox", "In shape %s wrong definition of cut planes", GetName());
2554 return;
2555 }
2556 Double_t xc=0, yc=0;
2557 Double_t zmin=0, zmax=0;
2558 Double_t z1;
2559 Double_t z[8];
2560 // check if nxy is in the phi range
2561 Double_t phi_low = TMath::ATan2(fNlow[1], fNlow[0]) *TMath::RadToDeg();
2563 Bool_t in_range_low = kFALSE;
2564 Bool_t in_range_hi = kFALSE;
2565
2566 Int_t i;
2567 for (i=0; i<2; i++) {
2568 if (phi_low<0) phi_low+=360.;
2569 Double_t dphi = fPhi2 -fPhi1;
2570 if (dphi < 0) dphi+=360.;
2571 Double_t ddp = phi_low-fPhi1;
2572 if (ddp<0) ddp += 360.;
2573 if (ddp <= dphi) {
2574 xc = fRmin*TMath::Cos(phi_low*TMath::DegToRad());
2575 yc = fRmin*TMath::Sin(phi_low*TMath::DegToRad());
2576 z1 = GetZcoord(xc, yc, -fDz);
2577 xc = fRmax*TMath::Cos(phi_low*TMath::DegToRad());
2578 yc = fRmax*TMath::Sin(phi_low*TMath::DegToRad());
2579 z1 = TMath::Min(z1, GetZcoord(xc, yc, -fDz));
2580 if (in_range_low)
2581 zmin = TMath::Min(zmin, z1);
2582 else
2583 zmin = z1;
2584 in_range_low = kTRUE;
2585 }
2586 phi_low += 180;
2587 if (phi_low>360) phi_low-=360.;
2588 }
2589
2590 for (i=0; i<2; i++) {
2591 if (phi_hi<0) phi_hi+=360.;
2592 Double_t dphi = fPhi2 -fPhi1;
2593 if (dphi < 0) dphi+=360.;
2594 Double_t ddp = phi_hi-fPhi1;
2595 if (ddp<0) ddp += 360.;
2596 if (ddp <= dphi) {
2597 xc = fRmin*TMath::Cos(phi_hi*TMath::DegToRad());
2598 yc = fRmin*TMath::Sin(phi_hi*TMath::DegToRad());
2599 z1 = GetZcoord(xc, yc, fDz);
2600 xc = fRmax*TMath::Cos(phi_hi*TMath::DegToRad());
2601 yc = fRmax*TMath::Sin(phi_hi*TMath::DegToRad());
2602 z1 = TMath::Max(z1, GetZcoord(xc, yc, fDz));
2603 if (in_range_hi)
2604 zmax = TMath::Max(zmax, z1);
2605 else
2606 zmax = z1;
2607 in_range_hi = kTRUE;
2608 }
2609 phi_hi += 180;
2610 if (phi_hi>360) phi_hi-=360.;
2611 }
2612
2613
2614 xc = fRmin*fC1;
2615 yc = fRmin*fS1;
2616 z[0] = GetZcoord(xc, yc, -fDz);
2617 z[4] = GetZcoord(xc, yc, fDz);
2618
2619 xc = fRmin*fC2;
2620 yc = fRmin*fS2;
2621 z[1] = GetZcoord(xc, yc, -fDz);
2622 z[5] = GetZcoord(xc, yc, fDz);
2623
2624 xc = fRmax*fC1;
2625 yc = fRmax*fS1;
2626 z[2] = GetZcoord(xc, yc, -fDz);
2627 z[6] = GetZcoord(xc, yc, fDz);
2628
2629 xc = fRmax*fC2;
2630 yc = fRmax*fS2;
2631 z[3] = GetZcoord(xc, yc, -fDz);
2632 z[7] = GetZcoord(xc, yc, fDz);
2633
2634 z1 = z[TMath::LocMin(4, &z[0])];
2635 if (in_range_low)
2636 zmin = TMath::Min(zmin, z1);
2637 else
2638 zmin = z1;
2639
2640 z1 = z[TMath::LocMax(4, &z[4])+4];
2641 if (in_range_hi)
2642 zmax = TMath::Max(zmax, z1);
2643 else
2644 zmax = z1;
2645
2646 fDZ = 0.5*(zmax-zmin);
2647 fOrigin[2] = 0.5*(zmax+zmin);
2648}
2649
2650////////////////////////////////////////////////////////////////////////////////
2651/// Compute normal to closest surface from POINT.
2652
2653void TGeoCtub::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
2654{
2655 Double_t saf[4];
2656 Bool_t isseg = kTRUE;
2657 if (TMath::Abs(fPhi2-fPhi1-360.)<1E-8) isseg=kFALSE;
2658 Double_t rsq = point[0]*point[0]+point[1]*point[1];
2659 Double_t r = TMath::Sqrt(rsq);
2660
2661 saf[0] = TMath::Abs(point[0]*fNlow[0] + point[1]*fNlow[1] + (fDz+point[2])*fNlow[2]);
2662 saf[1] = TMath::Abs(point[0]*fNhigh[0] + point[1]*fNhigh[1] - (fDz-point[2])*fNhigh[2]);
2663 saf[2] = (fRmin>1E-10)?TMath::Abs(r-fRmin):TGeoShape::Big();
2664 saf[3] = TMath::Abs(fRmax-r);
2665 Int_t i = TMath::LocMin(4,saf);
2666 if (isseg) {
2667 if (TGeoShape::IsCloseToPhi(saf[i], point,fC1,fS1,fC2,fS2)) {
2668 TGeoShape::NormalPhi(point,dir,norm,fC1,fS1,fC2,fS2);
2669 return;
2670 }
2671 }
2672 if (i==0) {
2673 memcpy(norm, fNlow, 3*sizeof(Double_t));
2674 if (norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2]<0) {
2675 norm[0] = -norm[0];
2676 norm[1] = -norm[1];
2677 norm[2] = -norm[2];
2678 }
2679 return;
2680 }
2681 if (i==1) {
2682 memcpy(norm, fNhigh, 3*sizeof(Double_t));
2683 if (norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2]<0) {
2684 norm[0] = -norm[0];
2685 norm[1] = -norm[1];
2686 norm[2] = -norm[2];
2687 }
2688 return;
2689 }
2690
2691 norm[2] = 0;
2692 Double_t phi = TMath::ATan2(point[1], point[0]);
2693 norm[0] = TMath::Cos(phi);
2694 norm[1] = TMath::Sin(phi);
2695 if (norm[0]*dir[0]+norm[1]*dir[1]<0) {
2696 norm[0] = -norm[0];
2697 norm[1] = -norm[1];
2698 }
2699}
2700
2701////////////////////////////////////////////////////////////////////////////////
2702/// check if point is contained in the cut tube
2703/// check the lower cut plane
2704
2706{
2707 Double_t zin = point[0]*fNlow[0]+point[1]*fNlow[1]+(point[2]+fDz)*fNlow[2];
2708 if (zin>0) return kFALSE;
2709 // check the higher cut plane
2710 zin = point[0]*fNhigh[0]+point[1]*fNhigh[1]+(point[2]-fDz)*fNhigh[2];
2711 if (zin>0) return kFALSE;
2712 // check radius
2713 Double_t r2 = point[0]*point[0]+point[1]*point[1];
2714 if ((r2<fRmin*fRmin) || (r2>fRmax*fRmax)) return kFALSE;
2715 // check phi
2716 Double_t phi = TMath::ATan2(point[1], point[0]) * TMath::RadToDeg();
2717 if (phi < 0 ) phi+=360.;
2718 Double_t dphi = fPhi2 -fPhi1;
2719 Double_t ddp = phi-fPhi1;
2720 if (ddp<0) ddp += 360.;
2721// if (ddp>360) ddp-=360;
2722 if (ddp > dphi) return kFALSE;
2723 return kTRUE;
2724}
2725
2726////////////////////////////////////////////////////////////////////////////////
2727/// Get range of shape for a given axis.
2728
2730{
2731 xlo = 0;
2732 xhi = 0;
2733 Double_t dx = 0;
2734 switch (iaxis) {
2735 case 1:
2736 xlo = fRmin;
2737 xhi = fRmax;
2738 dx = xhi-xlo;
2739 return dx;
2740 case 2:
2741 xlo = fPhi1;
2742 xhi = fPhi2;
2743 dx = xhi-xlo;
2744 return dx;
2745 }
2746 return dx;
2747}
2748
2749////////////////////////////////////////////////////////////////////////////////
2750/// compute real Z coordinate of a point belonging to either lower or
2751/// higher caps (z should be either +fDz or -fDz)
2752
2754{
2755 Double_t newz = 0;
2756 if (zc<0) newz = -fDz-(xc*fNlow[0]+yc*fNlow[1])/fNlow[2];
2757 else newz = fDz-(xc*fNhigh[0]+yc*fNhigh[1])/fNhigh[2];
2758 return newz;
2759}
2760
2761////////////////////////////////////////////////////////////////////////////////
2762/// compute distance from outside point to surface of the cut tube
2763
2764Double_t TGeoCtub::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
2765{
2766 if (iact<3 && safe) {
2767 *safe = Safety(point, kFALSE);
2768 if (iact==0) return TGeoShape::Big();
2769 if ((iact==1) && (step<=*safe)) return TGeoShape::Big();
2770 }
2771// Check if the bounding box is crossed within the requested distance
2772 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
2773 if (sdist>=step) return TGeoShape::Big();
2774 Double_t saf[2];
2775 saf[0] = point[0]*fNlow[0] + point[1]*fNlow[1] + (fDz+point[2])*fNlow[2];
2776 saf[1] = point[0]*fNhigh[0] + point[1]*fNhigh[1] + (point[2]-fDz)*fNhigh[2];
2777 Double_t rsq = point[0]*point[0]+point[1]*point[1];
2778 Double_t r = TMath::Sqrt(rsq);
2779 Double_t cpsi=0;
2780 Bool_t tub = kFALSE;
2781 if (TMath::Abs(fPhi2-fPhi1-360.)<1E-8) tub = kTRUE;
2782
2783 // find distance to shape
2784 Double_t r2;
2785 Double_t calf = dir[0]*fNlow[0]+dir[1]*fNlow[1]+dir[2]*fNlow[2];
2786 // check Z planes
2787 Double_t xi, yi, zi, s;
2788 if (saf[0]>0) {
2789 if (calf<0) {
2790 s = -saf[0]/calf;
2791 xi = point[0]+s*dir[0];
2792 yi = point[1]+s*dir[1];
2793 r2=xi*xi+yi*yi;
2794 if (((fRmin*fRmin)<=r2) && (r2<=(fRmax*fRmax))) {
2795 if (tub) return s;
2796 cpsi=(xi*fCm+yi*fSm)/TMath::Sqrt(r2);
2797 if (cpsi>=fCdfi) return s;
2798 }
2799 }
2800 }
2801 calf = dir[0]*fNhigh[0]+dir[1]*fNhigh[1]+dir[2]*fNhigh[2];
2802 if (saf[1]>0) {
2803 if (calf<0) {
2804 s = -saf[1]/calf;
2805 xi = point[0]+s*dir[0];
2806 yi = point[1]+s*dir[1];
2807 r2=xi*xi+yi*yi;
2808 if (((fRmin*fRmin)<=r2) && (r2<=(fRmax*fRmax))) {
2809 if (tub) return s;
2810 cpsi=(xi*fCm+yi*fSm)/TMath::Sqrt(r2);
2811 if (cpsi>=fCdfi) return s;
2812 }
2813 }
2814 }
2815
2816 // check outer cyl. surface
2817 Double_t nsq=dir[0]*dir[0]+dir[1]*dir[1];
2818 if (TMath::Abs(nsq)<1E-10) return TGeoShape::Big();
2819 Double_t rdotn=point[0]*dir[0]+point[1]*dir[1];
2820 Double_t b,d;
2821 // only r>fRmax coming inwards has to be considered
2822 if (r>fRmax && rdotn<0) {
2823 TGeoTube::DistToTube(rsq, nsq, rdotn, fRmax, b, d);
2824 if (d>0) {
2825 s=-b-d;
2826 if (s>0) {
2827 xi=point[0]+s*dir[0];
2828 yi=point[1]+s*dir[1];
2829 zi=point[2]+s*dir[2];
2830 if ((-xi*fNlow[0]-yi*fNlow[1]-(zi+fDz)*fNlow[2])>0) {
2831 if ((-xi*fNhigh[0]-yi*fNhigh[1]+(fDz-zi)*fNhigh[2])>0) {
2832 if (tub) return s;
2833 cpsi=(xi*fCm+yi*fSm)/fRmax;
2834 if (cpsi>=fCdfi) return s;
2835 }
2836 }
2837 }
2838 }
2839 }
2840 // check inner cylinder
2841 Double_t snxt=TGeoShape::Big();
2842 if (fRmin>0) {
2843 TGeoTube::DistToTube(rsq, nsq, rdotn, fRmin, b, d);
2844 if (d>0) {
2845 s=-b+d;
2846 if (s>0) {
2847 xi=point[0]+s*dir[0];
2848 yi=point[1]+s*dir[1];
2849 zi=point[2]+s*dir[2];
2850 if ((-xi*fNlow[0]-yi*fNlow[1]-(zi+fDz)*fNlow[2])>0) {
2851 if ((-xi*fNhigh[0]-yi*fNhigh[1]+(fDz-zi)*fNhigh[2])>0) {
2852 if (tub) return s;
2853 cpsi=(xi*fCm+yi*fSm)/fRmin;
2854 if (cpsi>=fCdfi) snxt=s;
2855 }
2856 }
2857 }
2858 }
2859 }
2860 // check phi planes
2861 if (tub) return snxt;
2862 Double_t un=dir[0]*fS1-dir[1]*fC1;
2863 if (un<-TGeoShape::Tolerance()) {
2864 s=(point[1]*fC1-point[0]*fS1)/un;
2865 if (s>=0) {
2866 xi=point[0]+s*dir[0];
2867 yi=point[1]+s*dir[1];
2868 zi=point[2]+s*dir[2];
2869 if ((-xi*fNlow[0]-yi*fNlow[1]-(zi+fDz)*fNlow[2])>0) {
2870 if ((-xi*fNhigh[0]-yi*fNhigh[1]+(fDz-zi)*fNhigh[2])>0) {
2871 r2=xi*xi+yi*yi;
2872 if ((fRmin*fRmin<=r2) && (r2<=fRmax*fRmax)) {
2873 if ((yi*fCm-xi*fSm)<=0) {
2874 if (s<snxt) snxt=s;
2875 }
2876 }
2877 }
2878 }
2879 }
2880 }
2881 un=dir[0]*fS2-dir[1]*fC2;
2882 if (un>TGeoShape::Tolerance()) {
2883 s=(point[1]*fC2-point[0]*fS2)/un;
2884 if (s>=0) {
2885 xi=point[0]+s*dir[0];
2886 yi=point[1]+s*dir[1];
2887 zi=point[2]+s*dir[2];
2888 if ((-xi*fNlow[0]-yi*fNlow[1]-(zi+fDz)*fNlow[2])>0) {
2889 if ((-xi*fNhigh[0]-yi*fNhigh[1]+(fDz-zi)*fNhigh[2])>0) {
2890 r2=xi*xi+yi*yi;
2891 if ((fRmin*fRmin<=r2) && (r2<=fRmax*fRmax)) {
2892 if ((yi*fCm-xi*fSm)>=0) {
2893 if (s<snxt) snxt=s;
2894 }
2895 }
2896 }
2897 }
2898 }
2899 }
2900 return snxt;
2901}
2902
2903////////////////////////////////////////////////////////////////////////////////
2904/// compute distance from inside point to surface of the cut tube
2905
2906Double_t TGeoCtub::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
2907{
2908 if (iact<3 && safe) *safe = Safety(point, kTRUE);
2909 if (iact==0) return TGeoShape::Big();
2910 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
2911 Double_t rsq = point[0]*point[0]+point[1]*point[1];
2912 Bool_t tub = kFALSE;
2913 if (TMath::Abs(fPhi2-fPhi1-360.)<1E-8) tub = kTRUE;
2914 // compute distance to surface
2915 // Do Z
2916 Double_t sz = TGeoShape::Big();
2917 Double_t saf[2];
2918 saf[0] = -point[0]*fNlow[0] - point[1]*fNlow[1] - (fDz+point[2])*fNlow[2];
2919 saf[1] = -point[0]*fNhigh[0] - point[1]*fNhigh[1] + (fDz-point[2])*fNhigh[2];
2920 Double_t calf = dir[0]*fNlow[0]+dir[1]*fNlow[1]+dir[2]*fNlow[2];
2921 if (calf>0) sz = saf[0]/calf;
2922
2923 calf = dir[0]*fNhigh[0]+dir[1]*fNhigh[1]+dir[2]*fNhigh[2];
2924 if (calf>0) {
2925 Double_t sz1 = saf[1]/calf;
2926 if (sz1<sz) sz = sz1;
2927 }
2928
2929 // Do R
2930 Double_t nsq=dir[0]*dir[0]+dir[1]*dir[1];
2931 // track parallel to Z
2932 if (TMath::Abs(nsq)<1E-10) return sz;
2933 Double_t rdotn=point[0]*dir[0]+point[1]*dir[1];
2935 Double_t b, d;
2936 Bool_t skip_outer = kFALSE;
2937 // inner cylinder
2938 if (fRmin>1E-10) {
2939 TGeoTube::DistToTube(rsq, nsq, rdotn, fRmin, b, d);
2940 if (d>0) {
2941 sr=-b-d;
2942 if (sr>0) skip_outer = kTRUE;
2943 }
2944 }
2945 // outer cylinder
2946 if (!skip_outer) {
2947 TGeoTube::DistToTube(rsq, nsq, rdotn, fRmax, b, d);
2948 if (d>0) {
2949 sr=-b+d;
2950 if (sr<0) sr=TGeoShape::Big();
2951 } else {
2952 return 0.; // already outside
2953 }
2954 }
2955 // phi planes
2956 Double_t sfmin = TGeoShape::Big();
2957 if (!tub) sfmin=TGeoShape::DistToPhiMin(point, dir, fS1, fC1, fS2, fC2, fSm, fCm);
2958 return TMath::Min(TMath::Min(sz,sr), sfmin);
2959}
2960
2961////////////////////////////////////////////////////////////////////////////////
2962/// Divide the tube along one axis.
2963
2964TGeoVolume *TGeoCtub::Divide(TGeoVolume * /*voldiv*/, const char * /*divname*/, Int_t /*iaxis*/, Int_t /*ndiv*/,
2965 Double_t /*start*/, Double_t /*step*/)
2966{
2967 Warning("Divide", "In shape %s division of a cut tube not implemented", GetName());
2968 return 0;
2969}
2970
2971////////////////////////////////////////////////////////////////////////////////
2972/// in case shape has some negative parameters, these has to be computed
2973/// in order to fit the mother
2974
2976{
2977 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
2978 if (!mother->TestShapeBit(kGeoTube)) {
2979 Error("GetMakeRuntimeShape", "Invalid mother for shape %s", GetName());
2980 return 0;
2981 }
2982 Double_t rmin, rmax, dz;
2983 rmin = fRmin;
2984 rmax = fRmax;
2985 dz = fDz;
2986 if (fDz<0) dz=((TGeoTube*)mother)->GetDz();
2987 if (fRmin<0)
2988 rmin = ((TGeoTube*)mother)->GetRmin();
2989 if ((fRmax<0) || (fRmax<=fRmin))
2990 rmax = ((TGeoTube*)mother)->GetRmax();
2991
2992 return (new TGeoCtub(rmin, rmax, dz, fPhi1, fPhi2, fNlow[0], fNlow[1], fNlow[2],
2993 fNhigh[0], fNhigh[1], fNhigh[2]));
2994}
2995
2996////////////////////////////////////////////////////////////////////////////////
2997/// print shape parameters
2998
3000{
3001 printf("*** Shape %s: TGeoCtub ***\n", GetName());
3002 printf(" lx = %11.5f\n", fNlow[0]);
3003 printf(" ly = %11.5f\n", fNlow[1]);
3004 printf(" lz = %11.5f\n", fNlow[2]);
3005 printf(" tx = %11.5f\n", fNhigh[0]);
3006 printf(" ty = %11.5f\n", fNhigh[1]);
3007 printf(" tz = %11.5f\n", fNhigh[2]);
3009}
3010
3011////////////////////////////////////////////////////////////////////////////////
3012/// computes the closest distance from given point to this shape, according
3013/// to option. The matching point on the shape is stored in spoint.
3014
3016{
3017 Double_t saf[4];
3018 Double_t rsq = point[0]*point[0]+point[1]*point[1];
3019 Double_t r = TMath::Sqrt(rsq);
3020 Bool_t isseg = kTRUE;
3021 if (TMath::Abs(fPhi2-fPhi1-360.)<1E-8) isseg=kFALSE;
3022
3023 saf[0] = -point[0]*fNlow[0] - point[1]*fNlow[1] - (fDz+point[2])*fNlow[2];
3024 saf[1] = -point[0]*fNhigh[0] - point[1]*fNhigh[1] + (fDz-point[2])*fNhigh[2];
3025 saf[2] = (fRmin<1E-10 && !isseg)?TGeoShape::Big():(r-fRmin);
3026 saf[3] = fRmax-r;
3027 Double_t safphi = TGeoShape::Big();
3028 if (isseg) safphi = TGeoShape::SafetyPhi(point, in, fPhi1, fPhi2);
3029
3030 if (in) {
3031 Double_t safe = saf[TMath::LocMin(4,saf)];
3032 return TMath::Min(safe, safphi);
3033 }
3034 for (Int_t i=0; i<4; i++) saf[i]=-saf[i];
3035 Double_t safe = saf[TMath::LocMax(4,saf)];
3036 if (isseg) return TMath::Max(safe, safphi);
3037 return safe;
3038}
3039
3040////////////////////////////////////////////////////////////////////////////////
3041/// set dimensions of a cut tube
3042
3044 Double_t lx, Double_t ly, Double_t lz, Double_t tx, Double_t ty, Double_t tz)
3045{
3046 SetTubsDimensions(rmin, rmax, dz, phi1, phi2);
3047 fNlow[0] = lx;
3048 fNlow[1] = ly;
3049 fNlow[2] = lz;
3050 fNhigh[0] = tx;
3051 fNhigh[1] = ty;
3052 fNhigh[2] = tz;
3053 ComputeBBox();
3054}
3055
3056////////////////////////////////////////////////////////////////////////////////
3057/// Save a primitive as a C++ statement(s) on output stream "out".
3058
3059void TGeoCtub::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
3060{
3062 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
3063 out << " rmin = " << fRmin << ";" << std::endl;
3064 out << " rmax = " << fRmax << ";" << std::endl;
3065 out << " dz = " << fDz << ";" << std::endl;
3066 out << " phi1 = " << fPhi1 << ";" << std::endl;
3067 out << " phi2 = " << fPhi2 << ";" << std::endl;
3068 out << " lx = " << fNlow[0] << ";" << std::endl;
3069 out << " ly = " << fNlow[1] << ";" << std::endl;
3070 out << " lz = " << fNlow[2] << ";" << std::endl;
3071 out << " tx = " << fNhigh[0] << ";" << std::endl;
3072 out << " ty = " << fNhigh[1] << ";" << std::endl;
3073 out << " tz = " << fNhigh[2] << ";" << std::endl;
3074 out << " TGeoShape *" << GetPointerName() << " = new TGeoCtub(\"" << GetName() << "\",rmin,rmax,dz,phi1,phi2,lx,ly,lz,tx,ty,tz);" << std::endl; TObject::SetBit(TGeoShape::kGeoSavePrimitive);
3075}
3076
3077////////////////////////////////////////////////////////////////////////////////
3078/// Set dimensions of the cut tube starting from a list.
3079
3081{
3082 SetCtubDimensions(param[0], param[1], param[2], param[3], param[4], param[5],
3083 param[6], param[7], param[8], param[9], param[10]);
3084 ComputeBBox();
3085}
3086
3087////////////////////////////////////////////////////////////////////////////////
3088/// Fills array with n random points located on the line segments of the shape mesh.
3089/// The output array must be provided with a length of minimum 3*npoints. Returns
3090/// true if operation is implemented.
3091
3093{
3094 return kFALSE;
3095}
3096
3097////////////////////////////////////////////////////////////////////////////////
3098/// Create mesh points for the cut tube.
3099
3101{
3102 Double_t dz;
3103 Int_t j, n;
3104 Double_t phi, phi1, phi2, dphi;
3105 phi1 = fPhi1;
3106 phi2 = fPhi2;
3107 if (phi2<phi1) phi2+=360.;
3108 n = gGeoManager->GetNsegments()+1;
3109
3110 dphi = (phi2-phi1)/(n-1);
3111 dz = fDz;
3112
3113 if (points) {
3114 Int_t indx = 0;
3115
3116 for (j = 0; j < n; j++) {
3117 phi = (phi1+j*dphi)*TMath::DegToRad();
3118 points[indx+6*n] = points[indx] = fRmin * TMath::Cos(phi);
3119 indx++;
3120 points[indx+6*n] = points[indx] = fRmin * TMath::Sin(phi);
3121 indx++;
3122 points[indx+6*n] = GetZcoord(points[indx-2], points[indx-1], dz);
3123 points[indx] = GetZcoord(points[indx-2], points[indx-1], -dz);
3124 indx++;
3125 }
3126 for (j = 0; j < n; j++) {
3127 phi = (phi1+j*dphi)*TMath::DegToRad();
3128 points[indx+6*n] = points[indx] = fRmax * TMath::Cos(phi);
3129 indx++;
3130 points[indx+6*n] = points[indx] = fRmax * TMath::Sin(phi);
3131 indx++;
3132 points[indx+6*n]= GetZcoord(points[indx-2], points[indx-1], dz);
3133 points[indx] = GetZcoord(points[indx-2], points[indx-1], -dz);
3134 indx++;
3135 }
3136 }
3137}
3138
3139////////////////////////////////////////////////////////////////////////////////
3140/// Create mesh points for the cut tube.
3141
3143{
3144 Double_t dz;
3145 Int_t j, n;
3146 Double_t phi, phi1, phi2, dphi;
3147 phi1 = fPhi1;
3148 phi2 = fPhi2;
3149 if (phi2<phi1) phi2+=360.;
3150 n = gGeoManager->GetNsegments()+1;
3151
3152 dphi = (phi2-phi1)/(n-1);
3153 dz = fDz;
3154
3155 if (points) {
3156 Int_t indx = 0;
3157
3158 for (j = 0; j < n; j++) {
3159 phi = (phi1+j*dphi)*TMath::DegToRad();
3160 points[indx+6*n] = points[indx] = fRmin * TMath::Cos(phi);
3161 indx++;
3162 points[indx+6*n] = points[indx] = fRmin * TMath::Sin(phi);
3163 indx++;
3164 points[indx+6*n] = GetZcoord(points[indx-2], points[indx-1], dz);
3165 points[indx] = GetZcoord(points[indx-2], points[indx-1], -dz);
3166 indx++;
3167 }
3168 for (j = 0; j < n; j++) {
3169 phi = (phi1+j*dphi)*TMath::DegToRad();
3170 points[indx+6*n] = points[indx] = fRmax * TMath::Cos(phi);
3171 indx++;
3172 points[indx+6*n] = points[indx] = fRmax * TMath::Sin(phi);
3173 indx++;
3174 points[indx+6*n]= GetZcoord(points[indx-2], points[indx-1], dz);
3175 points[indx] = GetZcoord(points[indx-2], points[indx-1], -dz);
3176 indx++;
3177 }
3178 }
3179}
3180
3181////////////////////////////////////////////////////////////////////////////////
3182/// Returns numbers of vertices, segments and polygons composing the shape mesh.
3183
3184void TGeoCtub::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
3185{
3186 TGeoTubeSeg::GetMeshNumbers(nvert,nsegs,npols);
3187}
3188
3189////////////////////////////////////////////////////////////////////////////////
3190/// Return number of vertices of the mesh representation
3191
3193{
3195 Int_t numPoints = n*4;
3196 return numPoints;
3197}
3198
3199////////////////////////////////////////////////////////////////////////////////
3200/// Fills a static 3D buffer and returns a reference.
3201
3202const TBuffer3D & TGeoCtub::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
3203{
3204 static TBuffer3DCutTube buffer;
3205
3206 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
3207
3208 if (reqSections & TBuffer3D::kShapeSpecific) {
3209 // These from TBuffer3DCutTube / TGeoCtub
3210 buffer.fRadiusInner = fRmin;
3211 buffer.fRadiusOuter = fRmax;
3212 buffer.fHalfLength = fDz;
3213 buffer.fPhiMin = fPhi1;
3214 buffer.fPhiMax = fPhi2;
3215
3216 for (UInt_t i = 0; i < 3; i++ ) {
3217 buffer.fLowPlaneNorm[i] = fNlow[i];
3218 buffer.fHighPlaneNorm[i] = fNhigh[i];
3219 }
3221 }
3222 if (reqSections & TBuffer3D::kRawSizes) {
3224 Int_t nbPnts = 4*n;
3225 Int_t nbSegs = 2*nbPnts;
3226 Int_t nbPols = nbPnts-2;
3227 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
3229 }
3230 }
3231 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
3232 SetPoints(buffer.fPnts);
3233 if (!buffer.fLocalFrame) {
3234 TransformPoints(buffer.fPnts, buffer.NbPnts());
3235 }
3236 SetSegsAndPols(buffer);
3238 }
3239
3240 return buffer;
3241}
3242
3243////////////////////////////////////////////////////////////////////////////////
3244/// Check the inside status for each of the points in the array.
3245/// Input: Array of point coordinates + vector size
3246/// Output: Array of Booleans for the inside of each point
3247
3248void TGeoCtub::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
3249{
3250 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
3251}
3252
3253////////////////////////////////////////////////////////////////////////////////
3254/// Compute the normal for an array o points so that norm.dot.dir is positive
3255/// Input: Arrays of point coordinates and directions + vector size
3256/// Output: Array of normal directions
3257
3258void TGeoCtub::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
3259{
3260 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
3261}
3262
3263////////////////////////////////////////////////////////////////////////////////
3264/// Compute distance from array of input points having directions specified by dirs. Store output in dists
3265
3266void TGeoCtub::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
3267{
3268 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
3269}
3270
3271////////////////////////////////////////////////////////////////////////////////
3272/// Compute distance from array of input points having directions specified by dirs. Store output in dists
3273
3274void TGeoCtub::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
3275{
3276 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
3277}
3278
3279////////////////////////////////////////////////////////////////////////////////
3280/// Compute safe distance from each of the points in the input array.
3281/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
3282/// Output: Safety values
3283
3284void TGeoCtub::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
3285{
3286 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
3287}
ROOT::R::TRInterface & r
Definition Object.C:4
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define s1(x)
Definition RSha256.hxx:91
int Int_t
Definition RtypesCore.h:45
const Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:187
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:231
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:245
XFontStruct * id
Definition TGX11.cxx:109
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
float xmin
float ymin
float xmax
float ymax
point * points
Definition X3DBuffer.c:22
Cut tube segment description class - see TBuffer3DTypes for producer classes.
Definition TBuffer3D.h:211
Double_t fLowPlaneNorm[3]
Definition TBuffer3D.h:223
Double_t fHighPlaneNorm[3]
Definition TBuffer3D.h:224
Tube segment description class - see TBuffer3DTypes for producer classes.
Definition TBuffer3D.h:184
Double_t fPhiMax
Definition TBuffer3D.h:203
Double_t fPhiMin
Definition TBuffer3D.h:202
Complete tube description class - see TBuffer3DTypes for producer classes.
Definition TBuffer3D.h:156
Double_t fRadiusInner
Definition TBuffer3D.h:174
Double_t fRadiusOuter
Definition TBuffer3D.h:175
Double_t fHalfLength
Definition TBuffer3D.h:176
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:114
UInt_t NbPols() const
Definition TBuffer3D.h:82
UInt_t NbPnts() const
Definition TBuffer3D.h:80
UInt_t NbSegs() const
Definition TBuffer3D.h:81
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
@ kShapeSpecific
Definition TBuffer3D.h:52
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Int_t * fSegs
Definition TBuffer3D.h:113
Bool_t fLocalFrame
Definition TBuffer3D.h:90
Int_t fColor
Definition TBuffer3D.h:88
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:112
Box class.
Definition TGeoBBox.h:18
Double_t fDX
Definition TGeoBBox.h:21
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
Compute distance from outside point to surface of the box.
Definition TGeoBBox.cxx:458
virtual void InspectShape() const
Prints shape parameters.
Definition TGeoBBox.cxx:819
Double_t fOrigin[3]
Definition TGeoBBox.h:24
Double_t fDY
Definition TGeoBBox.h:22
Double_t fDZ
Definition TGeoBBox.h:23
virtual void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const
Fills the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections...
The cut tubes constructor has the form:
Definition TGeoTube.h:172
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
computes the closest distance from given point to this shape, according to option.
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
compute distance from inside point to surface of the cut tube
virtual Bool_t Contains(const Double_t *point) const
check if point is contained in the cut tube check the lower cut plane
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
virtual void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual void ComputeBBox()
compute minimum bounding box of the ctub
virtual ~TGeoCtub()
destructor
virtual void SetPoints(Double_t *points) const
Create mesh points for the cut tube.
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
virtual TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const
in case shape has some negative parameters, these has to be computed in order to fit the mother
virtual void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
Compute safe distance from each of the points in the input array.
virtual void SetDimensions(Double_t *param)
Set dimensions of the cut tube starting from a list.
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
virtual void InspectShape() const
print shape parameters
Double_t GetZcoord(Double_t xc, Double_t yc, Double_t zc) const
compute real Z coordinate of a point belonging to either lower or higher caps (z should be either +fD...
virtual void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
Check the inside status for each of the points in the array.
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
compute distance from outside point to surface of the cut tube
virtual Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const
Fills array with n random points located on the line segments of the shape mesh.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
void SetCtubDimensions(Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2, Double_t lx, Double_t ly, Double_t lz, Double_t tx, Double_t ty, Double_t tz)
set dimensions of a cut tube
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide the tube along one axis.
Double_t fNlow[3]
Definition TGeoTube.h:175
Double_t fNhigh[3]
Definition TGeoTube.h:176
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
Int_t GetNsegments() const
Get number of segments approximating circles.
Geometrical transformation package.
Definition TGeoMatrix.h:41
Node containing an offset.
Definition TGeoNode.h:184
Base finder class for patterns.
void SetDivIndex(Int_t index)
Base abstract class for all shapes.
Definition TGeoShape.h:26
static Double_t Big()
Definition TGeoShape.h:88
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const =0
Int_t GetBasicColor() const
Get the basic color (0-7).
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
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.
static Bool_t IsInPhiRange(const Double_t *point, Double_t phi1, Double_t phi2)
Static method to check if a point is in the phi range (phi1, phi2) [degrees].
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
virtual const char * GetName() const
Get the shape name.
static void NormalPhi(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Static method to compute normal to phi planes.
@ kGeoSavePrimitive
Definition TGeoShape.h:65
@ kGeoTubeSeg
Definition TGeoShape.h:48
@ kGeoRunTimeShape
Definition TGeoShape.h:41
static Double_t Tolerance()
Definition TGeoShape.h:91
static Bool_t IsCloseToPhi(Double_t epsil, const Double_t *point, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
True if point is closer than epsil to one of the phi planes defined by c1,s1 or c2,...
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:163
A tube segment is a tube having a range in phi.
Definition TGeoTube.h:92
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this tube segment shape belonging to volume "voldiv" into ndiv volumes called divname,...
virtual void Sizeof3D() const
fill size of this 3-D object
virtual void AfterStreamer()
Function called after streaming an object of this class.
virtual ~TGeoTubeSeg()
destructor
TGeoTubeSeg()
Default constructor.
virtual Bool_t Contains(const Double_t *point) const
test if point is inside this tube segment first check if point is inside the tube
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
Double_t fPhi1
Definition TGeoTube.h:95
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
virtual void ComputeBBox()
compute bounding box of the tube segment
Double_t fPhi2
Definition TGeoTube.h:96
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
static Double_t DistFromInsideS(const Double_t *point, const Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz, 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) Boundary safe algorithm.
Double_t fC2
Definition TGeoTube.h:101
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
compute distance from outside point to surface of the tube segment fist localize point w....
static Double_t DistFromOutsideS(const Double_t *point, const Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz, Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi)
Static method to compute distance to arbitrary tube segment from outside point Boundary safe algorith...
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Double_t fCdfi
Definition TGeoTube.h:104
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
virtual void SetPoints(Double_t *points) const
Create tube segment mesh points.
Double_t fCm
Definition TGeoTube.h:103
static Double_t SafetyS(const Double_t *point, Bool_t in, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2, Int_t skipz=0)
Static method to compute the closest distance from given point to this shape.
virtual Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const
Fills array with n random points located on the line segments of the shape mesh.
void InitTrigonometry()
Init frequently used trigonometric values.
void SetTubsDimensions(Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2)
Set dimensions of the tube segment.
virtual TBuffer3D * MakeBuffer3D() const
Creates a TBuffer3D describing this shape.
Double_t fC1
Definition TGeoTube.h:99
Double_t fS1
Definition TGeoTube.h:98
Double_t fS2
Definition TGeoTube.h:100
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
computes the closest distance from given point InitTrigonometry();to this shape, according to option.
virtual void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
Check the inside status for each of the points in the array.
virtual void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
Compute safe distance from each of the points in the input array.
virtual void SetSegsAndPols(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons.
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
virtual TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const
in case shape has some negative parameters, these has to be computed in order to fit the mother
virtual void InspectShape() const
print shape parameters
static void ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t rmin, Double_t rmax, Double_t dz, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Compute normal to closest surface from POINT.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute closest distance from point px,py to each corner
Double_t fSm
Definition TGeoTube.h:102
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
Compute distance from inside point to surface of the tube segment Boundary safe algorithm.
virtual void SetDimensions(Double_t *param)
Set dimensions of the tube segment starting from a list.
Cylindrical tube class.
Definition TGeoTube.h:18
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
Compute distance from inside point to surface of the tube Boundary safe algorithm.
Definition TGeoTube.cxx:352
virtual void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
Check the inside status for each of the points in the array.
virtual void InspectShape() const
print shape parameters
Definition TGeoTube.cxx:652
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition TGeoTube.cxx:214
virtual const char * GetAxisName(Int_t iaxis) const
Returns name of axis IAXIS.
Definition TGeoTube.cxx:563
static void DistToTube(Double_t rsq, Double_t nsq, Double_t rdotn, Double_t radius, Double_t &b, Double_t &delta)
Static method computing the distance to a tube with given radius, starting from POINT along DIR direc...
Definition TGeoTube.cxx:479
void SetTubeDimensions(Double_t rmin, Double_t rmax, Double_t dz)
Set tube dimensions.
Definition TGeoTube.cxx:929
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoTube.cxx:609
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition TGeoTube.cxx:915
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
virtual void SetSegsAndPols(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons.
Definition TGeoTube.cxx:691
Double_t fRmin
Definition TGeoTube.h:21
virtual void SetDimensions(Double_t *param)
Set tube dimensions starting from a list.
Definition TGeoTube.cxx:941
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
Compute distance from outside point to surface of the tube and safe distance Boundary safe algorithm.
Definition TGeoTube.cxx:455
Double_t fDz
Definition TGeoTube.h:23
static Double_t DistFromOutsideS(const Double_t *point, const Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz)
Static method to compute distance from outside point to a tube with given parameters Boundary safe al...
Definition TGeoTube.cxx:368
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
computes the closest distance from given point to this shape, according to option.
Definition TGeoTube.cxx:846
virtual void ComputeBBox()
compute bounding box of the tube
Definition TGeoTube.cxx:231
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
Definition TGeoTube.cxx:240
static Double_t SafetyS(const Double_t *point, Bool_t in, Double_t rmin, Double_t rmax, Double_t dz, Int_t skipz=0)
computes the closest distance from given point to this shape, according to option.
Definition TGeoTube.cxx:886
Double_t fRmax
Definition TGeoTube.h:22
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute closest distance from point px,py to each corner
Definition TGeoTube.cxx:294
virtual Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const
Fills array with n random points located on the line segments of the shape mesh.
Definition TGeoTube.cxx:954
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Definition TGeoTube.cxx:580
virtual ~TGeoTube()
destructor
Definition TGeoTube.cxx:207
virtual void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
virtual void Sizeof3D() const
fill size of this 3-D object
static void ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t rmin, Double_t rmax, Double_t dz)
Compute normal to closest surface from POINT.
Definition TGeoTube.cxx:267
virtual void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
Compute safe distance from each of the points in the input array.
virtual Bool_t Contains(const Double_t *point) const
test if point is inside this tube
Definition TGeoTube.cxx:283
virtual TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const
in case shape has some negative parameters, these has to be computed in order to fit the mother
Definition TGeoTube.cxx:623
virtual void SetPoints(Double_t *points) const
create tube mesh points
Definition TGeoTube.cxx:995
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
Bool_t HasRmin() const
Definition TGeoTube.h:72
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this tube shape belonging to volume "voldiv" into ndiv volumes called divname,...
Definition TGeoTube.cxx:501
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
virtual TBuffer3D * MakeBuffer3D() const
Creates a TBuffer3D describing this shape.
Definition TGeoTube.cxx:666
static Double_t DistFromInsideS(const Double_t *point, const Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz)
Compute distance from inside point to surface of the tube (static) Boundary safe algorithm.
Definition TGeoTube.cxx:308
TGeoTube()
Default constructor.
Definition TGeoTube.cxx:149
Volume families.
Definition TGeoVolume.h:255
void AddVolume(TGeoVolume *vol)
Add a volume with valid shape to the list of volumes.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
void AddNodeOffset(TGeoVolume *vol, Int_t copy_no, Double_t offset=0, Option_t *option="")
Add a division node to the list of nodes.
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:174
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoVolume.h:232
Int_t GetNdaughters() const
Definition TGeoVolume.h:351
TObjArray * GetNodes()
Definition TGeoVolume.h:168
TObject * At(Int_t idx) const
Definition TObjArray.h:164
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:200
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:766
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
Basic string class.
Definition TString.h:136
const char * Data() const
Definition TString.h:369
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)
Return index of array with the minimum element.
Definition TMath.h:922
Short_t Max(Short_t a, Short_t b)
Definition TMathBase.h:208
T1 Sign(T1 a, T2 b)
Definition TMathBase.h:161
Double_t ATan2(Double_t y, Double_t x)
Definition TMath.h:629
Long64_t LocMax(Long64_t n, const T *a)
Return index of array with the maximum element.
Definition TMath.h:950
constexpr Double_t DegToRad()
Conversion from degree to radian:
Definition TMath.h:81
Double_t Sqrt(Double_t x)
Definition TMath.h:641
Short_t Min(Short_t a, Short_t b)
Definition TMathBase.h:176
Double_t Cos(Double_t)
Definition TMath.h:593
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Definition TMath.h:589
constexpr Double_t RadToDeg()
Conversion from radian to degree:
Definition TMath.h:73
Short_t Abs(Short_t d)
Definition TMathBase.h:120
constexpr Double_t TwoPi()
Definition TMath.h:44
auto * t1
Definition textangle.C:20