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