Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoPcon.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 24/10/01
3// TGeoPcon::Contains() 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 TGeoPcon
14\ingroup Shapes_classes
15
16A polycone is represented by a sequence of tubes/cones, glued together
17at defined Z planes. The polycone might have a phi segmentation, which
18globally applies to all the pieces. It has to be defined in two steps:
19
201. First call the TGeoPcon constructor to define a polycone:
21
22 ~~~{.cpp}
23 TGeoPcon(Double_t phi1,Double_t dphi,Int_t nz
24 ~~~
25
26 - `phi1:` starting phi angle in degrees
27 - `dphi:` total phi range
28 - `nz:` number of Z planes defining polycone sections (minimum 2)
29
302. Define one by one all sections [0, nz-1]
31
32~~~{.cpp}
33void TGeoPcon::DefineSection(Int_t i,Double_t z,
34Double_t rmin, Double_t rmax);
35~~~
36
37 - `i:` section index [0, nz-1]
38 - `z:` z coordinate of the section
39 - `rmin:` minimum radius corresponding too this section
40 - `rmax:` maximum radius.
41
42The first section (`i=0`) has to be positioned always the lowest Z
43coordinate. It defines the radii of the first cone/tube segment at its
44lower Z. The next section defines the end-cap of the first segment, but
45it can represent also the beginning of the next one. Any discontinuity
46in the radius has to be represented by a section defined at the same Z
47coordinate as the previous one. The Z coordinates of all sections must
48be sorted in increasing order. Any radius or Z coordinate of a given
49plane have corresponding getters:
50
51~~~{.cpp}
52Double_t TGeoPcon::GetRmin(Int_t i);
53Double_t TGeoPcon::GetRmax(Int_t i);
54Double_t TGeoPcon::GetZ(Int_t i);
55~~~
56
57Note that the last section should be defined last, since it triggers the
58computation of the bounding box of the polycone.
59
60Begin_Macro
61{
62 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
63 new TGeoManager("pcon", "poza10");
64 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
65 TGeoMedium *med = new TGeoMedium("MED",1,mat);
66 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
67 gGeoManager->SetTopVolume(top);
68 TGeoVolume *vol = gGeoManager->MakePcon("PCON",med, -30.0,300,4);
69 TGeoPcon *pcon = (TGeoPcon*)(vol->GetShape());
70 pcon->DefineSection(0,0,15,20);
71 pcon->DefineSection(1,20,15,20);
72 pcon->DefineSection(2,20,15,25);
73 pcon->DefineSection(3,50,15,20);
74 vol->SetLineWidth(2);
75 top->AddNode(vol,1);
76 gGeoManager->CloseGeometry();
77 gGeoManager->SetNsegments(30);
78 top->Draw();
79 TView *view = gPad->GetView();
80 view->ShowAxis();
81}
82End_Macro
83*/
84
85#include "TGeoPcon.h"
86
87#include <iostream>
88
89#include "TBuffer.h"
90#include "TGeoManager.h"
91#include "TGeoVolume.h"
92#include "TVirtualGeoPainter.h"
93#include "TGeoTube.h"
94#include "TGeoCone.h"
95#include "TBuffer3D.h"
96#include "TBuffer3DTypes.h"
97#include "TMath.h"
98
100
101////////////////////////////////////////////////////////////////////////////////
102/// dummy ctor
103
105 :TGeoBBox(),
106 fNz(0),
107 fPhi1(0.),
108 fDphi(0.),
109 fRmin(nullptr),
110 fRmax(nullptr),
111 fZ(nullptr),
112 fFullPhi(kFALSE),
113 fC1(0.),
114 fS1(0.),
115 fC2(0.),
116 fS2(0.),
117 fCm(0.),
118 fSm(0.),
119 fCdphi(0.)
120{
122}
123
124////////////////////////////////////////////////////////////////////////////////
125/// Default constructor
126
128 :TGeoBBox(0, 0, 0),
129 fNz(nz),
130 fPhi1(phi),
131 fDphi(dphi),
132 fRmin(nullptr),
133 fRmax(nullptr),
134 fZ(nullptr),
135 fFullPhi(kFALSE),
136 fC1(0.),
137 fS1(0.),
138 fC2(0.),
139 fS2(0.),
140 fCm(0.),
141 fSm(0.),
142 fCdphi(0.)
143{
145 while (fPhi1<0) fPhi1+=360.;
146 fRmin = new Double_t [nz];
147 fRmax = new Double_t [nz];
148 fZ = new Double_t [nz];
149 memset(fRmin, 0, nz*sizeof(Double_t));
150 memset(fRmax, 0, nz*sizeof(Double_t));
151 memset(fZ, 0, nz*sizeof(Double_t));
153 Double_t phi1 = fPhi1;
154 Double_t phi2 = phi1+fDphi;
155 Double_t phim = 0.5*(phi1+phi2);
163}
164
165////////////////////////////////////////////////////////////////////////////////
166/// Default constructor
167
168TGeoPcon::TGeoPcon(const char *name, Double_t phi, Double_t dphi, Int_t nz)
169 :TGeoBBox(name, 0, 0, 0),
170 fNz(nz),
171 fPhi1(phi),
172 fDphi(dphi),
173 fRmin(nullptr),
174 fRmax(nullptr),
175 fZ(nullptr),
176 fFullPhi(kFALSE),
177 fC1(0.),
178 fS1(0.),
179 fC2(0.),
180 fS2(0.),
181 fCm(0.),
182 fSm(0.),
183 fCdphi(0.)
184{
186 while (fPhi1<0) fPhi1+=360.;
187 fRmin = new Double_t [nz];
188 fRmax = new Double_t [nz];
189 fZ = new Double_t [nz];
190 memset(fRmin, 0, nz*sizeof(Double_t));
191 memset(fRmax, 0, nz*sizeof(Double_t));
192 memset(fZ, 0, nz*sizeof(Double_t));
194 Double_t phi1 = fPhi1;
195 Double_t phi2 = phi1+fDphi;
196 Double_t phim = 0.5*(phi1+phi2);
204}
205
206////////////////////////////////////////////////////////////////////////////////
207/// Default constructor in GEANT3 style
208/// - param[0] = phi1
209/// - param[1] = dphi
210/// - param[2] = nz
211/// - param[3] = z1
212/// - param[4] = Rmin1
213/// - param[5] = Rmax1
214/// ...
215
217 :TGeoBBox(0, 0, 0),
218 fNz(0),
219 fPhi1(0.),
220 fDphi(0.),
221 fRmin(nullptr),
222 fRmax(nullptr),
223 fZ(nullptr),
224 fFullPhi(kFALSE),
225 fC1(0.),
226 fS1(0.),
227 fC2(0.),
228 fS2(0.),
229 fCm(0.),
230 fSm(0.),
231 fCdphi(0.)
232{
234 SetDimensions(param);
235 ComputeBBox();
236}
237
238////////////////////////////////////////////////////////////////////////////////
239/// destructor
240
242{
243 if (fRmin) { delete[] fRmin; fRmin = nullptr; }
244 if (fRmax) { delete[] fRmax; fRmax = nullptr; }
245 if (fZ) { delete[] fZ; fZ = nullptr; }
246}
247
248////////////////////////////////////////////////////////////////////////////////
249/// Computes capacity of the shape in [length^3]
250
252{
253 Int_t ipl;
254 Double_t rmin1, rmax1, rmin2, rmax2, phi1, phi2, dz;
255 Double_t capacity = 0.;
256 phi1 = fPhi1;
257 phi2 = fPhi1 + fDphi;
258 for (ipl=0; ipl<fNz-1; ipl++) {
259 dz = 0.5*(fZ[ipl+1]-fZ[ipl]);
260 if (dz < TGeoShape::Tolerance()) continue;
261 rmin1 = fRmin[ipl];
262 rmax1 = fRmax[ipl];
263 rmin2 = fRmin[ipl+1];
264 rmax2 = fRmax[ipl+1];
265 capacity += TGeoConeSeg::Capacity(dz,rmin1,rmax1,rmin2,rmax2,phi1,phi2);
266 }
267 return capacity;
268}
269
270////////////////////////////////////////////////////////////////////////////////
271/// compute bounding box of the pcon
272/// Check if the sections are in increasing Z order
273
275{
276 for (Int_t isec=0; isec<fNz-1; isec++) {
277 if (TMath::Abs(fZ[isec]-fZ[isec+1]) < TGeoShape::Tolerance()) {
278 fZ[isec+1]=fZ[isec];
279 if (IsSameWithinTolerance(fRmin[isec], fRmin[isec+1]) &&
280 IsSameWithinTolerance(fRmax[isec], fRmax[isec+1])) {
281 InspectShape();
282 Error("ComputeBBox", "Duplicated section %d/%d for shape %s", isec, isec+1, GetName());
283 }
284 }
285 if (fZ[isec]>fZ[isec+1]) {
286 InspectShape();
287 Fatal("ComputeBBox", "Wrong section order");
288 }
289 }
290 // Check if the last sections are valid
291 if (TMath::Abs(fZ[1]-fZ[0]) < TGeoShape::Tolerance() ||
293 InspectShape();
294 Fatal("ComputeBBox","Shape %s at index %d: Not allowed first two or last two sections at same Z",
296 }
297 Double_t zmin = TMath::Min(fZ[0], fZ[fNz-1]);
298 Double_t zmax = TMath::Max(fZ[0], fZ[fNz-1]);
299 // find largest rmax an smallest rmin
300 Double_t rmin, rmax;
301 rmin = fRmin[TMath::LocMin(fNz, fRmin)];
302 rmax = fRmax[TMath::LocMax(fNz, fRmax)];
303
304 Double_t xc[4];
305 Double_t yc[4];
306 xc[0] = rmax*fC1;
307 yc[0] = rmax*fS1;
308 xc[1] = rmax*fC2;
309 yc[1] = rmax*fS2;
310 xc[2] = rmin*fC1;
311 yc[2] = rmin*fS1;
312 xc[3] = rmin*fC2;
313 yc[3] = rmin*fS2;
314
315 Double_t xmin = xc[TMath::LocMin(4, &xc[0])];
316 Double_t xmax = xc[TMath::LocMax(4, &xc[0])];
317 Double_t ymin = yc[TMath::LocMin(4, &yc[0])];
318 Double_t ymax = yc[TMath::LocMax(4, &yc[0])];
319
320 Double_t ddp = -fPhi1;
321 if (ddp<0) ddp+= 360;
322 if (ddp<=fDphi) xmax = rmax;
323 ddp = 90-fPhi1;
324 if (ddp<0) ddp+= 360;
325 if (ddp<=fDphi) ymax = rmax;
326 ddp = 180-fPhi1;
327 if (ddp<0) ddp+= 360;
328 if (ddp<=fDphi) xmin = -rmax;
329 ddp = 270-fPhi1;
330 if (ddp<0) ddp+= 360;
331 if (ddp<=fDphi) ymin = -rmax;
332 fOrigin[0] = (xmax+xmin)/2;
333 fOrigin[1] = (ymax+ymin)/2;
334 fOrigin[2] = (zmax+zmin)/2;
335 fDX = (xmax-xmin)/2;
336 fDY = (ymax-ymin)/2;
337 fDZ = (zmax-zmin)/2;
339}
340
341////////////////////////////////////////////////////////////////////////////////
342/// Compute normal to closest surface from POINT.
343
344void TGeoPcon::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
345{
346 memset(norm,0,3*sizeof(Double_t));
347 Double_t r;
348 Double_t ptnew[3];
349 Double_t dz, rmin1, rmax1, rmin2, rmax2;
350 Bool_t is_tube;
351 Int_t ipl = TMath::BinarySearch(fNz, fZ, point[2]);
352 if (ipl==(fNz-1) || ipl<0) {
353 // point outside Z range
354 norm[2] = TMath::Sign(1., dir[2]);
355 return;
356 }
357 Int_t iplclose = ipl;
358 if ((fZ[ipl+1]-point[2])<(point[2]-fZ[ipl])) iplclose++;
359 dz = TMath::Abs(fZ[iplclose]-point[2]);
360 if (dz<1E-5) {
361 if (iplclose==0 || iplclose==(fNz-1)) {
362 norm[2] = TMath::Sign(1., dir[2]);
363 return;
364 }
365 if (iplclose==ipl && TGeoShape::IsSameWithinTolerance(fZ[ipl],fZ[ipl-1])) {
366 r = TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
367 if (r<TMath::Max(fRmin[ipl],fRmin[ipl-1]) || r>TMath::Min(fRmax[ipl],fRmax[ipl-1])) {
368 norm[2] = TMath::Sign(1., dir[2]);
369 return;
370 }
371 } else {
372 if (TGeoShape::IsSameWithinTolerance(fZ[iplclose],fZ[iplclose+1])) {
373 r = TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
374 if (r<TMath::Max(fRmin[iplclose],fRmin[iplclose+1]) || r>TMath::Min(fRmax[iplclose],fRmax[iplclose+1])) {
375 norm[2] = TMath::Sign(1., dir[2]);
376 return;
377 }
378 }
379 }
380 } //-> Z done
381 memcpy(ptnew, point, 3*sizeof(Double_t));
382 dz = 0.5*(fZ[ipl+1]-fZ[ipl]);
384 norm[2] = TMath::Sign(1., dir[2]);
385 return;
386 }
387 ptnew[2] -= 0.5*(fZ[ipl]+fZ[ipl+1]);
388 rmin1 = fRmin[ipl];
389 rmax1 = fRmax[ipl];
390 rmin2 = fRmin[ipl+1];
391 rmax2 = fRmax[ipl+1];
392 is_tube = (TGeoShape::IsSameWithinTolerance(rmin1,rmin2) && TGeoShape::IsSameWithinTolerance(rmax1,rmax2))?kTRUE:kFALSE;
393 if (!fFullPhi) {
394 if (is_tube) TGeoTubeSeg::ComputeNormalS(ptnew,dir,norm,rmin1,rmax1,dz,fC1,fS1,fC2,fS2);
395 else TGeoConeSeg::ComputeNormalS(ptnew,dir,norm,dz,rmin1,rmax1,rmin2,rmax2,fC1,fS1,fC2,fS2);
396 } else {
397 if (is_tube) TGeoTube::ComputeNormalS(ptnew,dir,norm,rmin1,rmax1,dz);
398 else TGeoCone::ComputeNormalS(ptnew,dir,norm,dz,rmin1,rmax1,rmin2,rmax2);
399 }
400}
401
402////////////////////////////////////////////////////////////////////////////////
403/// test if point is inside this shape
404/// check total z range
405
407{
408 if ((point[2]<fZ[0]) || (point[2]>fZ[fNz-1])) return kFALSE;
409 // check R squared
410 Double_t r2 = point[0]*point[0]+point[1]*point[1];
411
412 Int_t izl = 0;
413 Int_t izh = fNz-1;
414 Int_t izt = (fNz-1)/2;
415 while ((izh-izl)>1) {
416 if (point[2] > fZ[izt]) izl = izt;
417 else izh = izt;
418 izt = (izl+izh)>>1;
419 }
420 // the point is in the section bounded by izl and izh Z planes
421
422 // compute Rmin and Rmax and test the value of R squared
423 Double_t rmin, rmax;
425 rmin = TMath::Min(fRmin[izl], fRmin[izh]);
426 rmax = TMath::Max(fRmax[izl], fRmax[izh]);
427 } else {
428 Double_t dz = fZ[izh] - fZ[izl];
429 Double_t dz1 = point[2] - fZ[izl];
430 rmin = (fRmin[izl]*(dz-dz1)+fRmin[izh]*dz1)/dz;
431 rmax = (fRmax[izl]*(dz-dz1)+fRmax[izh]*dz1)/dz;
432 }
433 if ((r2<rmin*rmin) || (r2>rmax*rmax)) return kFALSE;
434 // now check phi
436 if (r2<1E-10) return kTRUE;
437 Double_t phi = TMath::ATan2(point[1], point[0]) * TMath::RadToDeg();
438 if (phi < 0) phi+=360.0;
439 Double_t ddp = phi-fPhi1;
440 if (ddp<0) ddp+=360.;
441 if (ddp<=fDphi) return kTRUE;
442 return kFALSE;
443}
444
445////////////////////////////////////////////////////////////////////////////////
446/// compute closest distance from point px,py to each corner
447
449{
451 const Int_t numPoints = 2*n*fNz;
452 return ShapeDistancetoPrimitive(numPoints, px, py);
453}
454
455////////////////////////////////////////////////////////////////////////////////
456/// compute distance from inside point to surface of the polycone
457
458Double_t TGeoPcon::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
459{
460 if (iact<3 && safe) {
461 *safe = Safety(point, kTRUE);
462 if (iact==0) return TGeoShape::Big();
463 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
464 }
465 Double_t snxt = TGeoShape::Big();
466 Double_t sstep = 1E-6;
467 Double_t point_new[3];
468 // determine which z segment contains the point
469 Int_t ipl = TMath::BinarySearch(fNz, fZ, point[2]+TMath::Sign(1.E-10,dir[2]));
470 if (ipl<0) ipl=0;
471 if (ipl==(fNz-1)) ipl--;
472 Double_t dz = 0.5*(fZ[ipl+1]-fZ[ipl]);
473 Bool_t special_case = kFALSE;
474 if (dz<1e-9) {
475 // radius changing segment, make sure track is not in the XY plane
476 if (TGeoShape::IsSameWithinTolerance(dir[2], 0)) {
477 special_case = kTRUE;
478 } else {
479 //check if a close point is still contained
480 point_new[0] = point[0]+sstep*dir[0];
481 point_new[1] = point[1]+sstep*dir[1];
482 point_new[2] = point[2]+sstep*dir[2];
483 if (!Contains(point_new)) return 0.;
484 return (DistFromInside(point_new,dir,iact,step,safe)+sstep);
485 }
486 }
487 // determine if the current segment is a tube or a cone
488 Bool_t intub = kTRUE;
489 if (!TGeoShape::IsSameWithinTolerance(fRmin[ipl],fRmin[ipl+1])) intub=kFALSE;
490 else if (!TGeoShape::IsSameWithinTolerance(fRmax[ipl],fRmax[ipl+1])) intub=kFALSE;
491 // determine phi segmentation
492 memcpy(point_new, point, 2*sizeof(Double_t));
493 // new point in reference system of the current segment
494 point_new[2] = point[2]-0.5*(fZ[ipl]+fZ[ipl+1]);
495
496 if (special_case) {
497 if (!fFullPhi)
498 snxt = TGeoTubeSeg::DistFromInsideS(point_new, dir,
499 TMath::Min(fRmin[ipl],fRmin[ipl+1]), TMath::Max(fRmax[ipl],fRmax[ipl+1]),
500 dz, fC1,fS1,fC2,fS2,fCm,fSm,fCdphi);
501 else
502 snxt = TGeoTube::DistFromInsideS(point_new, dir,
503 TMath::Min(fRmin[ipl],fRmin[ipl+1]), TMath::Max(fRmax[ipl],fRmax[ipl+1]),dz);
504 return snxt;
505 }
506 if (intub) {
507 if (!fFullPhi)
508 snxt = TGeoTubeSeg::DistFromInsideS(point_new, dir, fRmin[ipl], fRmax[ipl],dz, fC1,fS1,fC2,fS2,fCm,fSm,fCdphi);
509 else
510 snxt = TGeoTube::DistFromInsideS(point_new, dir, fRmin[ipl], fRmax[ipl],dz);
511 } else {
512 if (!fFullPhi)
513 snxt = TGeoConeSeg::DistFromInsideS(point_new,dir,dz,fRmin[ipl],fRmax[ipl],fRmin[ipl+1],fRmax[ipl+1],fC1,fS1,fC2,fS2,fCm,fSm,fCdphi);
514 else
515 snxt=TGeoCone::DistFromInsideS(point_new,dir,dz,fRmin[ipl],fRmax[ipl],fRmin[ipl+1], fRmax[ipl+1]);
516 }
517
518 for (Int_t i=0; i<3; i++) point_new[i]=point[i]+(snxt+1E-6)*dir[i];
519 if (!Contains(&point_new[0])) return snxt;
520
521 snxt += DistFromInside(&point_new[0], dir, 3) + 1E-6;
522 return snxt;
523}
524
525////////////////////////////////////////////////////////////////////////////////
526/// compute distance to a pcon Z slice. Segment iz must be valid
527
528Double_t TGeoPcon::DistToSegZ(const Double_t *point, const Double_t *dir, Int_t &iz) const
529{
530 Double_t zmin=fZ[iz];
531 Double_t zmax=fZ[iz+1];
532 if (TGeoShape::IsSameWithinTolerance(zmin,zmax)) {
533 if (TGeoShape::IsSameWithinTolerance(dir[2],0)) return TGeoShape::Big();
534 Int_t istep=(dir[2]>0)?1:-1;
535 iz+=istep;
536 if (iz<0 || iz>(fNz-2)) return TGeoShape::Big();
537 return DistToSegZ(point,dir,iz);
538 }
539 Double_t dz=0.5*(zmax-zmin);
540 Double_t local[3];
541 memcpy(&local[0], point, 3*sizeof(Double_t));
542 local[2]=point[2]-0.5*(zmin+zmax);
543 Double_t snxt;
544 Double_t rmin1=fRmin[iz];
545 Double_t rmax1=fRmax[iz];
546 Double_t rmin2=fRmin[iz+1];
547 Double_t rmax2=fRmax[iz+1];
548
550 if (fFullPhi) snxt=TGeoTube::DistFromOutsideS(local, dir, rmin1, rmax1, dz);
551 else snxt=TGeoTubeSeg::DistFromOutsideS(local,dir,rmin1,rmax1,dz,fC1,fS1,fC2,fS2,fCm,fSm,fCdphi);
552 } else {
553 if (fFullPhi) snxt=TGeoCone::DistFromOutsideS(local,dir,dz,rmin1, rmax1,rmin2,rmax2);
554 else snxt=TGeoConeSeg::DistFromOutsideS(local,dir,dz,rmin1,rmax1,rmin2,rmax2,fC1,fS1,fC2,fS2,fCm,fSm,fCdphi);
555 }
556 if (snxt<1E20) return snxt;
557 // check next segment
558 if (TGeoShape::IsSameWithinTolerance(dir[2],0)) return TGeoShape::Big();
559 Int_t istep=(dir[2]>0)?1:-1;
560 iz+=istep;
561 if (iz<0 || iz>(fNz-2)) return TGeoShape::Big();
562 return DistToSegZ(point,dir,iz);
563}
564
565////////////////////////////////////////////////////////////////////////////////
566/// compute distance from outside point to surface of the tube
567
568Double_t TGeoPcon::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
569{
570 if ((iact<3) && safe) {
571 *safe = Safety(point, kFALSE);
572 if ((iact==1) && (*safe>step)) return TGeoShape::Big();
573 if (iact==0) return TGeoShape::Big();
574 }
575 // check if ray intersect outscribed cylinder
576 if ((point[2]<fZ[0]) && (dir[2]<=0)) return TGeoShape::Big();
577 if ((point[2]>fZ[fNz-1]) && (dir[2]>=0)) return TGeoShape::Big();
578// Check if the bounding box is crossed within the requested distance
579 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
580 if (sdist>=step) return TGeoShape::Big();
581
582 Double_t r2 = point[0]*point[0]+point[1]*point[1];
583 Double_t radmax=0;
584 radmax=fRmax[TMath::LocMax(fNz, fRmax)];
585 if (r2>(radmax*radmax)) {
586 Double_t rpr=-point[0]*dir[0]-point[1]*dir[1];
587 Double_t nxy=dir[0]*dir[0]+dir[1]*dir[1];
588 if (rpr<TMath::Sqrt((r2-radmax*radmax)*nxy)) return TGeoShape::Big();
589 }
590
591 // find in which Z segment we are
592 Int_t ipl = TMath::BinarySearch(fNz, fZ, point[2]);
593 Int_t ifirst = ipl;
594 if (ifirst<0) {
595 ifirst = 0;
596 } else if (ifirst>=(fNz-1)) {
597 ifirst = fNz-2;
598 }
599 // find if point is in the phi gap
600 Double_t phi=0;
601 if (!fFullPhi) {
602 phi = TMath::ATan2(point[1], point[0]);
603 if (phi<0) phi+=2.*TMath::Pi();
604 }
605
606 // compute distance to boundary
607 return DistToSegZ(point,dir,ifirst);
608}
609
610////////////////////////////////////////////////////////////////////////////////
611/// Defines z position of a section plane, rmin and rmax at this z. Sections
612/// should be defined in increasing or decreasing Z order and the last section
613/// HAS to be snum = fNz-1
614
616{
617 if ((snum<0) || (snum>=fNz)) return;
618 fZ[snum] = z;
619 fRmin[snum] = rmin;
620 fRmax[snum] = rmax;
621 if (rmin>rmax)
622 Warning("DefineSection", "Shape %s: invalid rmin=%g rmax=%g", GetName(), rmin, rmax);
623 if (snum==(fNz-1)) {
624 // Reorder sections in increasing Z order
625 if (fZ[0] > fZ[snum]) {
626 Int_t iz = 0;
627 Int_t izi = fNz-1;
628 Double_t temp;
629 while (iz<izi) {
630 temp = fZ[iz];
631 fZ[iz] = fZ[izi];
632 fZ[izi] = temp;
633 temp = fRmin[iz];
634 fRmin[iz] = fRmin[izi];
635 fRmin[izi] = temp;
636 temp = fRmax[iz];
637 fRmax[iz] = fRmax[izi];
638 fRmax[izi] = temp;
639 iz++;
640 izi--;
641 }
642 }
643 ComputeBBox();
644 }
645}
646
647////////////////////////////////////////////////////////////////////////////////
648/// Returns number of segments on each mesh circle segment.
649
651{
652 return gGeoManager->GetNsegments();
653}
654
655////////////////////////////////////////////////////////////////////////////////
656/// Divide this polycone shape belonging to volume "voldiv" into ndiv volumes
657/// called divname, from start position with the given step. Returns pointer
658/// to created division cell volume in case of Z divisions. Z divisions can be
659/// performed if the divided range is in between two consecutive Z planes.
660/// In case a wrong division axis is supplied, returns pointer to
661/// volume that was divided.
662
663TGeoVolume *TGeoPcon::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
664 Double_t start, Double_t step)
665{
666 TGeoShape *shape; //--- shape to be created
667 TGeoVolume *vol; //--- division volume to be created
668 TGeoVolumeMulti *vmulti; //--- generic divided volume
669 TGeoPatternFinder *finder; //--- finder to be attached
670 TString opt = ""; //--- option to be attached
671 Double_t zmin = start;
672 Double_t zmax = start+ndiv*step;
673 Int_t isect = -1;
674 Int_t is, id, ipl;
675 switch (iaxis) {
676 case 1: //--- R division
677 Error("Divide", "Shape %s: cannot divide a pcon on radius", GetName());
678 return 0;
679 case 2: //--- Phi division
680 finder = new TGeoPatternCylPhi(voldiv, ndiv, start, start+ndiv*step);
681 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
682 voldiv->SetFinder(finder);
683 finder->SetDivIndex(voldiv->GetNdaughters());
684 shape = new TGeoPcon(-step/2, step, fNz);
685 for (is=0; is<fNz; is++)
686 ((TGeoPcon*)shape)->DefineSection(is, fZ[is], fRmin[is], fRmax[is]);
687 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
688 vmulti->AddVolume(vol);
689 opt = "Phi";
690 for (id=0; id<ndiv; id++) {
691 voldiv->AddNodeOffset(vol, id, start+id*step+step/2, opt.Data());
692 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
693 }
694 return vmulti;
695 case 3: //--- Z division
696 // find start plane
697 for (ipl=0; ipl<fNz-1; ipl++) {
698 if (start<fZ[ipl]) continue;
699 else {
700 if ((start+ndiv*step)>fZ[ipl+1]) continue;
701 }
702 isect = ipl;
703 zmin = fZ[isect];
704 zmax= fZ[isect+1];
705 break;
706 }
707 if (isect<0) {
708 Error("Divide", "Shape %s: cannot divide pcon on Z if divided region is not between 2 planes", GetName());
709 return 0;
710 }
711 finder = new TGeoPatternZ(voldiv, ndiv, start, start+ndiv*step);
712 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
713 voldiv->SetFinder(finder);
714 finder->SetDivIndex(voldiv->GetNdaughters());
715 opt = "Z";
716 for (id=0; id<ndiv; id++) {
717 Double_t z1 = start+id*step;
718 Double_t z2 = start+(id+1)*step;
719 Double_t rmin1 = (fRmin[isect]*(zmax-z1)-fRmin[isect+1]*(zmin-z1))/(zmax-zmin);
720 Double_t rmax1 = (fRmax[isect]*(zmax-z1)-fRmax[isect+1]*(zmin-z1))/(zmax-zmin);
721 Double_t rmin2 = (fRmin[isect]*(zmax-z2)-fRmin[isect+1]*(zmin-z2))/(zmax-zmin);
722 Double_t rmax2 = (fRmax[isect]*(zmax-z2)-fRmax[isect+1]*(zmin-z2))/(zmax-zmin);
724 Bool_t is_seg = (fDphi<360)?kTRUE:kFALSE;
725 if (is_seg) {
726 if (is_tube) shape=new TGeoTubeSeg(fRmin[isect],fRmax[isect],step/2, fPhi1, fPhi1+fDphi);
727 else shape=new TGeoConeSeg(step/2, rmin1, rmax1, rmin2, rmax2, fPhi1, fPhi1+fDphi);
728 } else {
729 if (is_tube) shape=new TGeoTube(fRmin[isect],fRmax[isect],step/2);
730 else shape = new TGeoCone(step/2,rmin1,rmax1,rmin2,rmax2);
731 }
732 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
733 vmulti->AddVolume(vol);
734 voldiv->AddNodeOffset(vol, id, start+id*step+step/2, opt.Data());
735 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
736 }
737 return vmulti;
738 default:
739 Error("Divide", "Shape %s: Wrong axis %d for division",GetName(), iaxis);
740 return 0;
741 }
742}
743
744////////////////////////////////////////////////////////////////////////////////
745/// Returns name of axis IAXIS.
746
747const char *TGeoPcon::GetAxisName(Int_t iaxis) const
748{
749 switch (iaxis) {
750 case 1:
751 return "R";
752 case 2:
753 return "PHI";
754 case 3:
755 return "Z";
756 default:
757 return "UNDEFINED";
758 }
759}
760
761////////////////////////////////////////////////////////////////////////////////
762/// Get range of shape for a given axis.
763
765{
766 xlo = 0;
767 xhi = 0;
768 Double_t dx = 0;
769 switch (iaxis) {
770 case 2:
771 xlo = fPhi1;
772 xhi = fPhi1 + fDphi;
773 dx = fDphi;
774 return dx;
775 case 3:
776 xlo = fZ[0];
777 xhi = fZ[fNz-1];
778 dx = xhi-xlo;
779 return dx;
780 }
781 return dx;
782}
783
784////////////////////////////////////////////////////////////////////////////////
785/// Fill vector param[4] with the bounding cylinder parameters. The order
786/// is the following : Rmin, Rmax, Phi1, Phi2
787
789{
790 param[0] = fRmin[0]; // Rmin
791 param[1] = fRmax[0]; // Rmax
792 for (Int_t i=1; i<fNz; i++) {
793 if (fRmin[i] < param[0]) param[0] = fRmin[i];
794 if (fRmax[i] > param[1]) param[1] = fRmax[i];
795 }
796 param[0] *= param[0];
797 param[1] *= param[1];
799 param[2] = 0.;
800 param[3] = 360.;
801 return;
802 }
803 param[2] = (fPhi1<0)?(fPhi1+360.):fPhi1; // Phi1
804 param[3] = param[2]+fDphi; // Phi2
805}
806
807////////////////////////////////////////////////////////////////////////////////
808/// Returns Rmin for Z segment IPL.
809
811{
812 if (ipl<0 || ipl>(fNz-1)) {
813 Error("GetRmin","ipl=%i out of range (0,%i) in shape %s",ipl,fNz-1,GetName());
814 return 0.;
815 }
816 return fRmin[ipl];
817}
818
819////////////////////////////////////////////////////////////////////////////////
820/// Returns Rmax for Z segment IPL.
821
823{
824 if (ipl<0 || ipl>(fNz-1)) {
825 Error("GetRmax","ipl=%i out of range (0,%i) in shape %s",ipl,fNz-1,GetName());
826 return 0.;
827 }
828 return fRmax[ipl];
829}
830
831////////////////////////////////////////////////////////////////////////////////
832/// Returns Z for segment IPL.
833
835{
836 if (ipl<0 || ipl>(fNz-1)) {
837 Error("GetZ","ipl=%i out of range (0,%i) in shape %s",ipl,fNz-1,GetName());
838 return 0.;
839 }
840 return fZ[ipl];
841}
842
843////////////////////////////////////////////////////////////////////////////////
844/// print shape parameters
845
847{
848 printf("*** Shape %s: TGeoPcon ***\n", GetName());
849 printf(" Nz = %i\n", fNz);
850 printf(" phi1 = %11.5f\n", fPhi1);
851 printf(" dphi = %11.5f\n", fDphi);
852 for (Int_t ipl=0; ipl<fNz; ipl++)
853 printf(" plane %i: z=%11.5f Rmin=%11.5f Rmax=%11.5f\n", ipl, fZ[ipl], fRmin[ipl], fRmax[ipl]);
854 printf(" Bounding box:\n");
856}
857
858////////////////////////////////////////////////////////////////////////////////
859/// Creates a TBuffer3D describing *this* shape.
860/// Coordinates are in local reference frame.
861
863{
864 Int_t nbPnts, nbSegs, nbPols;
865 GetMeshNumbers(nbPnts, nbSegs, nbPols);
866 if (nbPnts <= 0) return nullptr;
867
869 nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols);
870 if (buff)
871 {
872 SetPoints(buff->fPnts);
873 SetSegsAndPols(*buff);
874 }
875
876 return buff;
877}
878
879////////////////////////////////////////////////////////////////////////////////
880/// Fill TBuffer3D structure for segments and polygons.
881
883{
884 if (!HasInsideSurface()) {
886 return;
887 }
888
889 Int_t i, j;
890 const Int_t n = gGeoManager->GetNsegments()+1;
891 Int_t nz = GetNz();
892 if (nz < 2) return;
893 Int_t nbPnts = nz*2*n;
894 if (nbPnts <= 0) return;
895 Double_t dphi = GetDphi();
896
897 Bool_t specialCase = TGeoShape::IsSameWithinTolerance(dphi,360);
899
900 Int_t indx = 0, indx2, k;
901
902 //inside & outside circles, number of segments: 2*nz*(n-1)
903 // special case number of segments: 2*nz*n
904 for (i = 0; i < nz*2; i++) {
905 indx2 = i*n;
906 for (j = 1; j < n; j++) {
907 buff.fSegs[indx++] = c;
908 buff.fSegs[indx++] = indx2+j-1;
909 buff.fSegs[indx++] = indx2+j;
910 }
911 if (specialCase) {
912 buff.fSegs[indx++] = c;
913 buff.fSegs[indx++] = indx2+j-1;
914 buff.fSegs[indx++] = indx2;
915 }
916 }
917
918 //bottom & top lines, number of segments: 2*n
919 for (i = 0; i < 2; i++) {
920 indx2 = i*(nz-1)*2*n;
921 for (j = 0; j < n; j++) {
922 buff.fSegs[indx++] = c;
923 buff.fSegs[indx++] = indx2+j;
924 buff.fSegs[indx++] = indx2+n+j;
925 }
926 }
927
928 //inside & outside cylinders, number of segments: 2*(nz-1)*n
929 for (i = 0; i < (nz-1); i++) {
930 //inside cylinder
931 indx2 = i*n*2;
932 for (j = 0; j < n; j++) {
933 buff.fSegs[indx++] = c+2;
934 buff.fSegs[indx++] = indx2+j;
935 buff.fSegs[indx++] = indx2+n*2+j;
936 }
937 //outside cylinder
938 indx2 = i*n*2+n;
939 for (j = 0; j < n; j++) {
940 buff.fSegs[indx++] = c+3;
941 buff.fSegs[indx++] = indx2+j;
942 buff.fSegs[indx++] = indx2+n*2+j;
943 }
944 }
945
946 //left & right sections, number of segments: 2*(nz-2)
947 // special case number of segments: 0
948 if (!specialCase) {
949 for (i = 1; i < (nz-1); i++) {
950 for (j = 0; j < 2; j++) {
951 buff.fSegs[indx++] = c;
952 buff.fSegs[indx++] = 2*i * n + j*(n-1);
953 buff.fSegs[indx++] = (2*i+1) * n + j*(n-1);
954 }
955 }
956 }
957
958 Int_t m = n - 1 + (specialCase ? 1 : 0);
959 indx = 0;
960
961 //bottom & top, number of polygons: 2*(n-1)
962 // special case number of polygons: 2*n
963 for (j = 0; j < n-1; j++) {
964 buff.fPols[indx++] = c+3;
965 buff.fPols[indx++] = 4;
966 buff.fPols[indx++] = 2*nz*m+j;
967 buff.fPols[indx++] = m+j;
968 buff.fPols[indx++] = 2*nz*m+j+1;
969 buff.fPols[indx++] = j;
970 }
971 for (j = 0; j < n-1; j++) {
972 buff.fPols[indx++] = c+3;
973 buff.fPols[indx++] = 4;
974 buff.fPols[indx++] = 2*nz*m+n+j;
975 buff.fPols[indx++] = (nz*2-2)*m+j;
976 buff.fPols[indx++] = 2*nz*m+n+j+1;
977 buff.fPols[indx++] = (nz*2-2)*m+m+j;
978 }
979 if (specialCase) {
980 buff.fPols[indx++] = c+3;
981 buff.fPols[indx++] = 4;
982 buff.fPols[indx++] = 2*nz*m+j;
983 buff.fPols[indx++] = m+j;
984 buff.fPols[indx++] = 2*nz*m;
985 buff.fPols[indx++] = j;
986
987 buff.fPols[indx++] = c+3;
988 buff.fPols[indx++] = 4;
989 buff.fPols[indx++] = 2*nz*m+n+j;
990 buff.fPols[indx++] = (nz*2-2)*m+m+j;
991 buff.fPols[indx++] = 2*nz*m+n;
992 buff.fPols[indx++] = (nz*2-2)*m+j;
993 }
994
995 //inside & outside, number of polygons: (nz-1)*2*(n-1)
996 for (k = 0; k < (nz-1); k++) {
997 for (j = 0; j < n-1; j++) {
998 buff.fPols[indx++] = c;
999 buff.fPols[indx++] = 4;
1000 buff.fPols[indx++] = 2*k*m+j;
1001 buff.fPols[indx++] = nz*2*m+(2*k+2)*n+j+1;
1002 buff.fPols[indx++] = (2*k+2)*m+j;
1003 buff.fPols[indx++] = nz*2*m+(2*k+2)*n+j;
1004 }
1005 for (j = 0; j < n-1; j++) {
1006 buff.fPols[indx++] = c+1;
1007 buff.fPols[indx++] = 4;
1008 buff.fPols[indx++] = (2*k+1)*m+j;
1009 buff.fPols[indx++] = nz*2*m+(2*k+3)*n+j;
1010 buff.fPols[indx++] = (2*k+3)*m+j;
1011 buff.fPols[indx++] = nz*2*m+(2*k+3)*n+j+1;
1012 }
1013 if (specialCase) {
1014 buff.fPols[indx++] = c;
1015 buff.fPols[indx++] = 4;
1016 buff.fPols[indx++] = 2*k*m+j;
1017 buff.fPols[indx++] = nz*2*m+(2*k+2)*n;
1018 buff.fPols[indx++] = (2*k+2)*m+j;
1019 buff.fPols[indx++] = nz*2*m+(2*k+2)*n+j;
1020
1021 buff.fPols[indx++] = c+1;
1022 buff.fPols[indx++] = 4;
1023 buff.fPols[indx++] = (2*k+1)*m+j;
1024 buff.fPols[indx++] = nz*2*m+(2*k+3)*n+j;
1025 buff.fPols[indx++] = (2*k+3)*m+j;
1026 buff.fPols[indx++] = nz*2*m+(2*k+3)*n;
1027 }
1028 }
1029
1030 //left & right sections, number of polygons: 2*(nz-1)
1031 // special case number of polygons: 0
1032 if (!specialCase) {
1033 indx2 = nz*2*(n-1);
1034 for (k = 0; k < (nz-1); k++) {
1035 buff.fPols[indx++] = c+2;
1036 buff.fPols[indx++] = 4;
1037 buff.fPols[indx++] = k==0 ? indx2 : indx2+2*nz*n+2*(k-1);
1038 buff.fPols[indx++] = indx2+2*(k+1)*n;
1039 buff.fPols[indx++] = indx2+2*nz*n+2*k;
1040 buff.fPols[indx++] = indx2+(2*k+3)*n;
1041
1042 buff.fPols[indx++] = c+2;
1043 buff.fPols[indx++] = 4;
1044 buff.fPols[indx++] = k==0 ? indx2+n-1 : indx2+2*nz*n+2*(k-1)+1;
1045 buff.fPols[indx++] = indx2+(2*k+3)*n+n-1;
1046 buff.fPols[indx++] = indx2+2*nz*n+2*k+1;
1047 buff.fPols[indx++] = indx2+2*(k+1)*n+n-1;
1048 }
1049 buff.fPols[indx-8] = indx2+n;
1050 buff.fPols[indx-2] = indx2+2*n-1;
1051 }
1052}
1053
1054
1055////////////////////////////////////////////////////////////////////////////////
1056/// Fill TBuffer3D structure for segments and polygons, when no inner surface exists
1057
1059{
1060 const Int_t n = gGeoManager->GetNsegments() + 1;
1061 const Int_t nz = GetNz();
1062 const Int_t nbPnts = nz * n + 2;
1063
1064 if ((nz < 2) || (nbPnts <= 0) || (n < 2)) return;
1065
1066 Int_t c = GetBasicColor();
1067
1068 Int_t indx = 0, indx1 = 0, indx2 = 0, i, j;
1069
1070 // outside circles, number of segments: nz*n
1071 for (i = 0; i < nz; i++) {
1072 indx2 = i * n;
1073 for (j = 1; j < n; j++) {
1074 buff.fSegs[indx++] = c;
1075 buff.fSegs[indx++] = indx2 + j - 1;
1076 buff.fSegs[indx++] = indx2 + j % (n-1);
1077 }
1078 }
1079
1080 indx2 = 0;
1081 // bottom lines
1082 for (j = 0; j < n; j++) {
1083 buff.fSegs[indx++] = c;
1084 buff.fSegs[indx++] = indx2 + j % (n-1);
1085 buff.fSegs[indx++] = nbPnts - 2;
1086 }
1087
1088 indx2 = (nz-1) * n;
1089 // top lines
1090 for (j = 0; j < n; j++) {
1091 buff.fSegs[indx++] = c;
1092 buff.fSegs[indx++] = indx2 + j % (n-1);
1093 buff.fSegs[indx++] = nbPnts - 1;
1094 }
1095
1096 // outside cylinders, number of segments: (nz-1)*n
1097 for (i = 0; i < (nz - 1); i++) {
1098 // outside cylinder
1099 indx2 = i * n;
1100 for (j = 0; j < n; j++) {
1101 buff.fSegs[indx++] = c;
1102 buff.fSegs[indx++] = indx2 + j % (n-1);
1103 buff.fSegs[indx++] = indx2 + n + j % (n-1);
1104 }
1105 }
1106
1107 indx = 0;
1108
1109 // bottom cap
1110 indx1 = 0; // start of first z layer
1111 indx2 = nz*(n-1);
1112 for (j = 0; j < n - 1; j++) {
1113 buff.fPols[indx++] = c;
1114 buff.fPols[indx++] = 3;
1115 buff.fPols[indx++] = indx1 + j;
1116 buff.fPols[indx++] = indx2 + j + 1;
1117 buff.fPols[indx++] = indx2 + j;
1118 }
1119
1120 // top cap
1121 indx1 = (nz-1)*(n-1); // start last z layer
1122 indx2 = nz*(n-1) + n;
1123 for (j = 0; j < n - 1; j++) {
1124 buff.fPols[indx++] = c;
1125 buff.fPols[indx++] = 3;
1126 buff.fPols[indx++] = indx1 + j; // last z layer
1127 buff.fPols[indx++] = indx2 + j;
1128 buff.fPols[indx++] = indx2 + j + 1;
1129 }
1130
1131 // outside, number of polygons: (nz-1)*(n-1)
1132 for (Int_t k = 0; k < (nz - 1); k++) {
1133 indx1 = k*(n-1);
1134 indx2 = nz*(n-1) + n*2 + k*n;
1135 for (j = 0; j < n-1; j++) {
1136 buff.fPols[indx++] = c;
1137 buff.fPols[indx++] = 4;
1138 buff.fPols[indx++] = indx1 + j;
1139 buff.fPols[indx++] = indx2 + j;
1140 buff.fPols[indx++] = indx1 + j + (n-1);
1141 buff.fPols[indx++] = indx2 + j + 1;
1142 }
1143 }
1144}
1145
1146////////////////////////////////////////////////////////////////////////////////
1147/// Compute safety from POINT to segment between planes ipl, ipl+1 within safmin.
1148
1150{
1151 if (ipl<0 || ipl>fNz-2) return (safmin+1.); // error in input plane
1152// Get info about segment.
1153 Double_t dz = 0.5*(fZ[ipl+1]-fZ[ipl]);
1154 if (dz<1E-9) return 1E9; // radius-changing segment
1155 Double_t ptnew[3];
1156 memcpy(ptnew, point, 3*sizeof(Double_t));
1157 ptnew[2] -= 0.5*(fZ[ipl]+fZ[ipl+1]);
1158 Double_t safe = TMath::Abs(ptnew[2])-dz;
1159 if (safe>safmin) return TGeoShape::Big(); // means: stop checking further segments
1160 Double_t rmin1 = fRmin[ipl];
1161 Double_t rmax1 = fRmax[ipl];
1162 Double_t rmin2 = fRmin[ipl+1];
1163 Double_t rmax2 = fRmax[ipl+1];
1165 if (!fFullPhi) {
1166 if (is_tube) safe = TGeoTubeSeg::SafetyS(ptnew,in,rmin1,rmax1, dz,fPhi1,fPhi1+fDphi,0);
1167 else safe = TGeoConeSeg::SafetyS(ptnew,in,dz,rmin1,rmax1,rmin2,rmax2,fPhi1,fPhi1+fDphi,0);
1168 } else {
1169 if (is_tube) safe = TGeoTube::SafetyS(ptnew,in,rmin1,rmax1,dz,0);
1170 else safe = TGeoCone::SafetyS(ptnew,in,dz,rmin1,rmax1,rmin2,rmax2,0);
1171 }
1172 if (safe<0) safe=0;
1173 return safe;
1174}
1175
1176////////////////////////////////////////////////////////////////////////////////
1177/// computes the closest distance from given point to this shape, according
1178/// to option. The matching point on the shape is stored in spoint.
1179/// localize the Z segment
1180
1182{
1183 Double_t safmin, saftmp;
1184 Double_t dz;
1185 Int_t ipl, iplane;
1186
1187 if (in) {
1188 //---> point is inside pcon
1189 ipl = TMath::BinarySearch(fNz, fZ, point[2]);
1190 if (ipl==(fNz-1)) return 0; // point on last Z boundary
1191 if (ipl<0) return 0; // point on first Z boundary
1192 if (ipl>0 && TGeoShape::IsSameWithinTolerance(fZ[ipl-1],fZ[ipl]) && TGeoShape::IsSameWithinTolerance(point[2],fZ[ipl-1])) ipl--;
1193 dz = 0.5*(fZ[ipl+1]-fZ[ipl]);
1194 if (dz<1E-8) {
1195 // Point on a segment-changing plane
1196 safmin = TMath::Min(point[2]-fZ[ipl-1],fZ[ipl+2]-point[2]);
1197 saftmp = TGeoShape::Big();
1198 if (fDphi<360) saftmp = TGeoShape::SafetyPhi(point,in,fPhi1,fPhi1+fDphi);
1199 if (saftmp<safmin) safmin = saftmp;
1200 Double_t radius = TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
1201 if (fRmin[ipl]>0) safmin = TMath::Min(safmin, radius-fRmin[ipl]);
1202 if (fRmin[ipl+1]>0) safmin = TMath::Min(safmin, radius-fRmin[ipl+1]);
1203 safmin = TMath::Min(safmin, fRmax[ipl]-radius);
1204 safmin = TMath::Min(safmin, fRmax[ipl+1]-radius);
1205 if (safmin<0) safmin = 0;
1206 return safmin;
1207 }
1208 // Check safety for current segment
1209 safmin = SafetyToSegment(point, ipl);
1210 if (safmin>1E10) {
1211 // something went wrong - point is not inside current segment
1212 return 0.;
1213 }
1214 if (safmin<1E-6) return TMath::Abs(safmin); // point on radius-changing plane
1215 // check increasing iplanes
1216/*
1217 iplane = ipl+1;
1218 saftmp = 0.;
1219 while ((iplane<fNz-1) && saftmp<1E10) {
1220 saftmp = TMath::Abs(SafetyToSegment(point,iplane,kFALSE,safmin));
1221 if (saftmp<safmin) safmin=saftmp;
1222 iplane++;
1223 }
1224 // now decreasing nplanes
1225 iplane = ipl-1;
1226 saftmp = 0.;
1227 while ((iplane>=0) && saftmp<1E10) {
1228 saftmp = TMath::Abs(SafetyToSegment(point,iplane,kFALSE,safmin));
1229 if (saftmp<safmin) safmin=saftmp;
1230 iplane--;
1231 }
1232*/
1233 return safmin;
1234 }
1235 //---> point is outside pcon
1236 ipl = TMath::BinarySearch(fNz, fZ, point[2]);
1237 if (ipl<0) ipl=0;
1238 else if (ipl==fNz-1) ipl=fNz-2;
1239 dz = 0.5*(fZ[ipl+1]-fZ[ipl]);
1240 if (dz<1E-8 && (ipl+2<fNz)) {
1241 ipl++;
1242 dz = 0.5*(fZ[ipl+1]-fZ[ipl]);
1243 }
1244 // Check safety for current segment
1245 safmin = SafetyToSegment(point, ipl, kFALSE);
1246 if (safmin<1E-6) return TMath::Abs(safmin); // point on radius-changing plane
1247 saftmp = 0.;
1248 // check increasing iplanes
1249 iplane = ipl+1;
1250 saftmp = 0.;
1251 while ((iplane<fNz-1) && saftmp<1E10) {
1252 saftmp = TMath::Abs(SafetyToSegment(point,iplane,kFALSE,safmin));
1253 if (saftmp<safmin) safmin=saftmp;
1254 iplane++;
1255 }
1256 // now decreasing nplanes
1257 iplane = ipl-1;
1258 saftmp = 0.;
1259 while ((iplane>=0) && saftmp<1E10) {
1260 saftmp = TMath::Abs(SafetyToSegment(point,iplane,kFALSE,safmin));
1261 if (saftmp<safmin) safmin=saftmp;
1262 iplane--;
1263 }
1264 return safmin;
1265}
1266
1267////////////////////////////////////////////////////////////////////////////////
1268/// Save a primitive as a C++ statement(s) on output stream "out".
1269
1270void TGeoPcon::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
1271{
1273 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
1274 out << " phi1 = " << fPhi1 << ";" << std::endl;
1275 out << " dphi = " << fDphi << ";" << std::endl;
1276 out << " nz = " << fNz << ";" << std::endl;
1277 out << " auto " << GetPointerName() << " = new TGeoPcon(\"" << GetName() << "\", phi1, dphi, nz);" << std::endl;
1278 for (Int_t i=0; i<fNz; i++) {
1279 out << " z = " << fZ[i] << ";" << std::endl;
1280 out << " rmin = " << fRmin[i] << ";" << std::endl;
1281 out << " rmax = " << fRmax[i] << ";" << std::endl;
1282 out << " " << GetPointerName() << "->DefineSection(" << i << ", z,rmin,rmax);" << std::endl;
1283 }
1285}
1286
1287////////////////////////////////////////////////////////////////////////////////
1288/// Set polycone dimensions starting from an array.
1289
1291{
1292 fPhi1 = param[0];
1293 while (fPhi1<0) fPhi1 += 360.;
1294 fDphi = param[1];
1295 fNz = (Int_t)param[2];
1296 if (fNz<2) {
1297 Error("SetDimensions","Pcon %s: Number of Z sections must be > 2", GetName());
1298 return;
1299 }
1300 if (fRmin) delete [] fRmin;
1301 if (fRmax) delete [] fRmax;
1302 if (fZ) delete [] fZ;
1303 fRmin = new Double_t [fNz];
1304 fRmax = new Double_t [fNz];
1305 fZ = new Double_t [fNz];
1306 memset(fRmin, 0, fNz*sizeof(Double_t));
1307 memset(fRmax, 0, fNz*sizeof(Double_t));
1308 memset(fZ, 0, fNz*sizeof(Double_t));
1310 Double_t phi1 = fPhi1;
1311 Double_t phi2 = phi1+fDphi;
1312 Double_t phim = 0.5*(phi1+phi2);
1313 fC1 = TMath::Cos(phi1*TMath::DegToRad());
1314 fS1 = TMath::Sin(phi1*TMath::DegToRad());
1315 fC2 = TMath::Cos(phi2*TMath::DegToRad());
1316 fS2 = TMath::Sin(phi2*TMath::DegToRad());
1317 fCm = TMath::Cos(phim*TMath::DegToRad());
1318 fSm = TMath::Sin(phim*TMath::DegToRad());
1320
1321 for (Int_t i=0; i<fNz; i++)
1322 DefineSection(i, param[3+3*i], param[4+3*i], param[5+3*i]);
1323}
1324
1325////////////////////////////////////////////////////////////////////////////////
1326/// create polycone mesh points
1327
1329{
1330 Double_t phi, dphi;
1332 dphi = fDphi/(n-1);
1333 Int_t i, j;
1334 Int_t indx = 0;
1335
1336 Bool_t hasInside = HasInsideSurface();
1337
1338 if (points) {
1339 for (i = 0; i < fNz; i++) {
1340 if (hasInside)
1341 for (j = 0; j < n; j++) {
1342 phi = (fPhi1+j*dphi)*TMath::DegToRad();
1343 points[indx++] = fRmin[i] * TMath::Cos(phi);
1344 points[indx++] = fRmin[i] * TMath::Sin(phi);
1345 points[indx++] = fZ[i];
1346 }
1347 for (j = 0; j < n; j++) {
1348 phi = (fPhi1+j*dphi)*TMath::DegToRad();
1349 points[indx++] = fRmax[i] * TMath::Cos(phi);
1350 points[indx++] = fRmax[i] * TMath::Sin(phi);
1351 points[indx++] = fZ[i];
1352 }
1353 }
1354 if (!hasInside) {
1355 points[indx++] = 0;
1356 points[indx++] = 0;
1357 points[indx++] = fZ[0];
1358
1359 points[indx++] = 0;
1360 points[indx++] = 0;
1361 points[indx++] = fZ[GetNz()-1];
1362 }
1363 }
1364}
1365
1366////////////////////////////////////////////////////////////////////////////////
1367/// create polycone mesh points
1368
1370{
1371 Double_t phi, dphi;
1373 dphi = fDphi/(n-1);
1374 Int_t i, j;
1375 Int_t indx = 0;
1376
1377 Bool_t hasInside = HasInsideSurface();
1378
1379 if (points) {
1380 for (i = 0; i < fNz; i++) {
1381 if (hasInside)
1382 for (j = 0; j < n; j++) {
1383 phi = (fPhi1+j*dphi)*TMath::DegToRad();
1384 points[indx++] = fRmin[i] * TMath::Cos(phi);
1385 points[indx++] = fRmin[i] * TMath::Sin(phi);
1386 points[indx++] = fZ[i];
1387 }
1388 for (j = 0; j < n; j++) {
1389 phi = (fPhi1+j*dphi)*TMath::DegToRad();
1390 points[indx++] = fRmax[i] * TMath::Cos(phi);
1391 points[indx++] = fRmax[i] * TMath::Sin(phi);
1392 points[indx++] = fZ[i];
1393 }
1394 }
1395 if (!hasInside) {
1396 points[indx++] = 0;
1397 points[indx++] = 0;
1398 points[indx++] = fZ[0];
1399
1400 points[indx++] = 0;
1401 points[indx++] = 0;
1402 points[indx++] = fZ[GetNz()-1];
1403 }
1404 }
1405}
1406////////////////////////////////////////////////////////////////////////////////
1407/// Return number of vertices of the mesh representation
1408
1410{
1411 Int_t nvert, nsegs, npols;
1412 GetMeshNumbers(nvert, nsegs, npols);
1413 return nvert;
1414}
1415
1416////////////////////////////////////////////////////////////////////////////////
1417/// fill size of this 3-D object
1418
1420{
1421}
1422
1423////////////////////////////////////////////////////////////////////////////////
1424/// Returns true when pgon has internal surface
1425/// It will be only disabled when all Rmin values are 0
1426
1428{
1429 // only when full 360 is used, internal part can be excluded
1430 Bool_t specialCase = TGeoShape::IsSameWithinTolerance(GetDphi(), 360);
1431 if (!specialCase) return kTRUE;
1432
1433 for (Int_t i = 0; i < GetNz(); i++)
1434 if (fRmin[i] > 0.)
1435 return kTRUE;
1436
1437 return kFALSE;
1438}
1439
1440
1441////////////////////////////////////////////////////////////////////////////////
1442/// Returns numbers of vertices, segments and polygons composing the shape mesh.
1443
1444void TGeoPcon::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
1445{
1446 nvert = nsegs = npols = 0;
1447
1449 Int_t nz = GetNz();
1450 if (nz < 2) return;
1451
1452 if (HasInsideSurface()) {
1453 Bool_t specialCase = TGeoShape::IsSameWithinTolerance(GetDphi(),360);
1454 nvert = nz*2*n;
1455 nsegs = 4*(nz*n-1+(specialCase ? 1 : 0));
1456 npols = 2*(nz*n-1+(specialCase ? 1 : 0));
1457 } else {
1458 nvert = nz * n + 2;
1459 nsegs = nz * (n - 1) + n * 2 + (nz - 1) * n;
1460 npols = 2 * (n - 1) + (nz - 1) * (n - 1);
1461 }
1462}
1463
1464////////////////////////////////////////////////////////////////////////////////
1465/// Fills a static 3D buffer and returns a reference.
1466
1467const TBuffer3D & TGeoPcon::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
1468{
1469 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
1470
1471 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
1472
1473 if (reqSections & TBuffer3D::kRawSizes) {
1474 Int_t nbPnts, nbSegs, nbPols;
1475 GetMeshNumbers(nbPnts, nbSegs, nbPols);
1476 if (nbPnts > 0) {
1477 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
1479 }
1480 }
1481 }
1482 // TODO: Push down to TGeoShape?? Would have to do raw sizes set first..
1483 // can rest of TGeoShape be deferred until after this?
1484 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
1485 SetPoints(buffer.fPnts);
1486 if (!buffer.fLocalFrame) {
1487 TransformPoints(buffer.fPnts, buffer.NbPnts());
1488 }
1489
1490 SetSegsAndPols(buffer);
1492 }
1493
1494 return buffer;
1495}
1496
1497////////////////////////////////////////////////////////////////////////////////
1498/// Stream an object of class TGeoPcon.
1499
1501{
1502 if (R__b.IsReading()) {
1503 R__b.ReadClassBuffer(TGeoPcon::Class(),this);
1505 Double_t phi1 = fPhi1;
1506 Double_t phi2 = phi1+fDphi;
1507 Double_t phim = 0.5*(phi1+phi2);
1508 fC1 = TMath::Cos(phi1*TMath::DegToRad());
1509 fS1 = TMath::Sin(phi1*TMath::DegToRad());
1510 fC2 = TMath::Cos(phi2*TMath::DegToRad());
1511 fS2 = TMath::Sin(phi2*TMath::DegToRad());
1512 fCm = TMath::Cos(phim*TMath::DegToRad());
1513 fSm = TMath::Sin(phim*TMath::DegToRad());
1515 } else {
1516 R__b.WriteClassBuffer(TGeoPcon::Class(),this);
1517 }
1518}
1519
1520////////////////////////////////////////////////////////////////////////////////
1521/// Check the inside status for each of the points in the array.
1522/// Input: Array of point coordinates + vector size
1523/// Output: Array of Booleans for the inside of each point
1524
1525void TGeoPcon::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
1526{
1527 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
1528}
1529
1530////////////////////////////////////////////////////////////////////////////////
1531/// Compute the normal for an array o points so that norm.dot.dir is positive
1532/// Input: Arrays of point coordinates and directions + vector size
1533/// Output: Array of normal directions
1534
1535void TGeoPcon::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
1536{
1537 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
1538}
1539
1540////////////////////////////////////////////////////////////////////////////////
1541/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1542
1543void TGeoPcon::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1544{
1545 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
1546}
1547
1548////////////////////////////////////////////////////////////////////////////////
1549/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1550
1551void TGeoPcon::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1552{
1553 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
1554}
1555
1556////////////////////////////////////////////////////////////////////////////////
1557/// Compute safe distance from each of the points in the input array.
1558/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1559/// Output: Safety values
1560
1561void TGeoPcon::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
1562{
1563 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
1564}
#define c(i)
Definition RSha256.hxx:101
#define e(i)
Definition RSha256.hxx:103
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize id
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
float xmin
float ymin
float xmax
float ymax
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:114
UInt_t NbPnts() const
Definition TBuffer3D.h:80
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Int_t * fSegs
Definition TBuffer3D.h:113
Bool_t fLocalFrame
Definition TBuffer3D.h:90
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:112
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Box class.
Definition TGeoBBox.h:18
Double_t fDX
Definition TGeoBBox.h:21
virtual void InspectShape() const
Prints shape parameters.
Definition TGeoBBox.cxx:819
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const
Compute distance from outside point to surface of the box.
Definition TGeoBBox.cxx:458
Double_t fOrigin[3]
Definition TGeoBBox.h:24
Double_t fDY
Definition TGeoBBox.h:22
Double_t fDZ
Definition TGeoBBox.h:23
virtual void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const
Fills the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections...
A cone segment is a cone having a range in phi.
Definition TGeoCone.h:102
static void ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Compute normal to closest surface from POINT.
static Double_t DistFromOutsideS(const Double_t *point, const Double_t *dir, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi)
compute distance from outside point to surface of arbitrary tube
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
static Double_t SafetyS(const Double_t *point, Bool_t in, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2, Int_t skipz=0)
Static method to compute the closest distance from given point to this shape.
static Double_t DistFromInsideS(const Double_t *point, const Double_t *dir, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi)
compute distance from inside point to surface of the tube segment
The cones are defined by 5 parameters:
Definition TGeoCone.h:18
static void ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
Compute normal to closest surface from POINT.
Definition TGeoCone.cxx:240
static Double_t DistFromInsideS(const Double_t *point, const Double_t *dir, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
Compute distance from inside point to surface of the cone (static) Boundary safe algorithm.
Definition TGeoCone.cxx:291
static Double_t DistFromOutsideS(const Double_t *point, const Double_t *dir, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
Compute distance from outside point to surface of the tube Boundary safe algorithm.
Definition TGeoCone.cxx:377
static Double_t SafetyS(const Double_t *point, Bool_t in, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Int_t skipz=0)
computes the closest distance from given point to this shape, according to option.
Definition TGeoCone.cxx:901
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
TObjArray * GetListOfShapes() const
Int_t GetNsegments() const
Get number of segments approximating circles.
Node containing an offset.
Definition TGeoNode.h:184
Base finder class for patterns.
void SetDivIndex(Int_t index)
A polycone is represented by a sequence of tubes/cones, glued together at defined Z planes.
Definition TGeoPcon.h:18
Double_t fSm
Cosine of (phi1+phi2)/2.
Definition TGeoPcon.h:33
virtual const char * GetAxisName(Int_t iaxis) const
Returns name of axis IAXIS.
Definition TGeoPcon.cxx:747
Double_t * GetRmax() const
Definition TGeoPcon.h:81
virtual void SetDimensions(Double_t *param)
Set polycone dimensions starting from an array.
Double_t GetDphi() const
Definition TGeoPcon.h:76
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const
compute distance from outside point to surface of the tube
Definition TGeoPcon.cxx:568
Double_t SafetyToSegment(const Double_t *point, Int_t ipl, Bool_t in=kTRUE, Double_t safmin=TGeoShape::Big()) const
Compute safety from POINT to segment between planes ipl, ipl+1 within safmin.
Double_t * fRmax
Definition TGeoPcon.h:25
virtual void Sizeof3D() const
fill size of this 3-D object
virtual void SetPoints(Double_t *points) const
create polycone mesh points
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const
compute distance from inside point to surface of the polycone
Definition TGeoPcon.cxx:458
virtual void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
Check the inside status for each of the points in the array.
virtual void ComputeBBox()
compute bounding box of the pcon Check if the sections are in increasing Z order
Definition TGeoPcon.cxx:274
Double_t * fRmin
Definition TGeoPcon.h:24
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
Definition TGeoPcon.cxx:344
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute closest distance from point px,py to each corner
Definition TGeoPcon.cxx:448
Int_t fNz
Definition TGeoPcon.h:21
Double_t * fZ
Definition TGeoPcon.h:26
virtual void DefineSection(Int_t snum, Double_t z, Double_t rmin, Double_t rmax)
Defines z position of a section plane, rmin and rmax at this z.
Definition TGeoPcon.cxx:615
Double_t fC1
Full phi range flag.
Definition TGeoPcon.h:28
Double_t fCdphi
Sine of (phi1+phi2)/2.
Definition TGeoPcon.h:34
Bool_t fFullPhi
Definition TGeoPcon.h:27
Double_t fS1
Cosine of phi1.
Definition TGeoPcon.h:29
virtual void InspectShape() const
print shape parameters
Definition TGeoPcon.cxx:846
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this polycone shape belonging to volume "voldiv" into ndiv volumes called divname,...
Definition TGeoPcon.cxx:663
virtual void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
Compute safe distance from each of the points in the input array.
virtual void SetSegsAndPols(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons.
Definition TGeoPcon.cxx:882
virtual Bool_t Contains(const Double_t *point) const
test if point is inside this shape check total z range
Definition TGeoPcon.cxx:406
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
Double_t DistToSegZ(const Double_t *point, const Double_t *dir, Int_t &iz) const
compute distance to a pcon Z slice. Segment iz must be valid
Definition TGeoPcon.cxx:528
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition TGeoPcon.cxx:251
Double_t fCm
Sine of phi1+dphi.
Definition TGeoPcon.h:32
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
Double_t * GetZ() const
Definition TGeoPcon.h:83
virtual void Streamer(TBuffer &)
Stream an object of class TGeoPcon.
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
static TClass * Class()
virtual TBuffer3D * MakeBuffer3D() const
Creates a TBuffer3D describing this shape.
Definition TGeoPcon.cxx:862
virtual void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
Bool_t HasInsideSurface() const
Returns true when pgon has internal surface It will be only disabled when all Rmin values are 0.
virtual Int_t GetNsegments() const
Returns number of segments on each mesh circle segment.
Definition TGeoPcon.cxx:650
Double_t fPhi1
Definition TGeoPcon.h:22
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Definition TGeoPcon.cxx:764
Double_t fS2
Cosine of phi1+dphi.
Definition TGeoPcon.h:31
virtual ~TGeoPcon()
destructor
Definition TGeoPcon.cxx:241
void SetSegsAndPolsNoInside(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons, when no inner surface exists.
Double_t fDphi
Definition TGeoPcon.h:23
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
Double_t fC2
Sine of phi1.
Definition TGeoPcon.h:30
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Int_t GetNz() const
Definition TGeoPcon.h:77
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoPcon.cxx:788
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
computes the closest distance from given point to this shape, according to option.
Double_t * GetRmin() const
Definition TGeoPcon.h:79
TGeoPcon()
dummy ctor
Definition TGeoPcon.cxx:104
Base abstract class for all shapes.
Definition TGeoShape.h:26
static Double_t Big()
Definition TGeoShape.h:88
Int_t GetBasicColor() const
Get the basic color (0-7).
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
static Double_t SafetyPhi(const Double_t *point, Bool_t in, Double_t phi1, Double_t phi2)
Static method to compute safety w.r.t a phi corner defined by cosines/sines of the angles phi1,...
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
const char * GetPointerName() const
Provide a pointer name containing uid.
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
virtual const char * GetName() const
Get the shape name.
@ kGeoClosedShape
Definition TGeoShape.h:60
@ kGeoSavePrimitive
Definition TGeoShape.h:65
static Double_t Tolerance()
Definition TGeoShape.h:91
A tube segment is a tube having a range in phi.
Definition TGeoTube.h:92
static Double_t DistFromInsideS(const Double_t *point, const Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz, Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi)
Compute distance from inside point to surface of the tube segment (static) Boundary safe algorithm.
static Double_t DistFromOutsideS(const Double_t *point, const Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz, Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi)
Static method to compute distance to arbitrary tube segment from outside point Boundary safe algorith...
static Double_t SafetyS(const Double_t *point, Bool_t in, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2, Int_t skipz=0)
Static method to compute the closest distance from given point to this shape.
static void ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t rmin, Double_t rmax, Double_t dz, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Compute normal to closest surface from POINT.
Cylindrical tube class.
Definition TGeoTube.h:18
static Double_t DistFromOutsideS(const Double_t *point, const Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz)
Static method to compute distance from outside point to a tube with given parameters Boundary safe al...
Definition TGeoTube.cxx:368
static Double_t SafetyS(const Double_t *point, Bool_t in, Double_t rmin, Double_t rmax, Double_t dz, Int_t skipz=0)
computes the closest distance from given point to this shape, according to option.
Definition TGeoTube.cxx:886
static void ComputeNormalS(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t rmin, Double_t rmax, Double_t dz)
Compute normal to closest surface from POINT.
Definition TGeoTube.cxx:267
static Double_t DistFromInsideS(const Double_t *point, const Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz)
Compute distance from inside point to surface of the tube (static) Boundary safe algorithm.
Definition TGeoTube.cxx:308
Volume families.
Definition TGeoVolume.h:255
void AddVolume(TGeoVolume *vol)
Add a volume with valid shape to the list of volumes.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
void AddNodeOffset(TGeoVolume *vol, Int_t copy_no, Double_t offset=0, Option_t *option="")
Add a division node to the list of nodes.
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:174
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoVolume.h:232
Int_t GetNdaughters() const
Definition TGeoVolume.h:351
TObjArray * GetNodes()
Definition TGeoVolume.h:168
Int_t IndexOf(const TObject *obj) const override
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:956
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:774
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:998
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:380
const Int_t n
Definition legend1.C:16
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Definition TMath.h:980
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
T1 Sign(T1 a, T2 b)
Returns a value with the magnitude of a and the sign of b.
Definition TMathBase.h:175
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:644
Long64_t LocMax(Long64_t n, const T *a)
Returns index of array with the maximum element.
Definition TMath.h:1010
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:79
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:660
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:592
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:586
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:347
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
TMarker m
Definition textangle.C:8