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