Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoTorus.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 28/07/03
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TGeoTorus
13\ingroup Geometry_classes
14
15Torus segment class. A torus has 5 parameters :
16 - R - axial radius
17 - Rmin - inner radius
18 - Rmax - outer radius
19 - Phi1 - starting phi
20 - Dphi - phi extent
21
22
23Begin_Macro(source)
24{
25 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
26 new TGeoManager("torus", "poza2");
27 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
28 TGeoMedium *med = new TGeoMedium("MED",1,mat);
29 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
30 gGeoManager->SetTopVolume(top);
31 TGeoVolume *vol = gGeoManager->MakeTorus("TORUS",med, 40,20,25,0,270);
32 top->AddNode(vol,1);
33 gGeoManager->CloseGeometry();
34 gGeoManager->SetNsegments(30);
35 top->Draw();
36 TView *view = gPad->GetView();
37 view->ShowAxis();
38}
39End_Macro
40*/
41
42#include <iostream>
43
44#include "TGeoManager.h"
45#include "TGeoVolume.h"
46#include "TGeoTube.h"
47#include "TVirtualGeoPainter.h"
48#include "TGeoTorus.h"
49#include "TBuffer3D.h"
50#include "TBuffer3DTypes.h"
51#include "TMath.h"
52
54
55////////////////////////////////////////////////////////////////////////////////
56/// Default constructor
57
59{
61 fR = 0.0;
62 fRmin = 0.0;
63 fRmax = 0.0;
64 fPhi1 = 0.0;
65 fDphi = 0.0;
66}
67
68////////////////////////////////////////////////////////////////////////////////
69/// Constructor without name.
70
72 :TGeoBBox(0, 0, 0)
73{
75 SetTorusDimensions(r, rmin, rmax, phi1, dphi);
76 if ((fRmin<0) || (fRmax<0))
79}
80
81////////////////////////////////////////////////////////////////////////////////
82/// Constructor with name.
83
85 :TGeoBBox(name, 0, 0, 0)
86{
88 SetTorusDimensions(r, rmin, rmax, phi1, dphi);
89 if ((fRmin<0) || (fRmax<0))
92}
93
94////////////////////////////////////////////////////////////////////////////////
95/// Constructor based on an array of parameters.
96/// - param[0] = R
97/// - param[1] = Rmin
98/// - param[2] = Rmax
99/// - param[3] = Phi1
100/// - param[4] = Dphi
101
103 :TGeoBBox(0, 0, 0)
104{
106 SetDimensions(param);
108 ComputeBBox();
109}
110
111////////////////////////////////////////////////////////////////////////////////
112/// Computes capacity of the shape in [length^3]
113
115{
116 Double_t capacity = (fDphi/180.)*TMath::Pi()*TMath::Pi()*fR*(fRmax*fRmax-fRmin*fRmin);
117 return capacity;
118}
119
120////////////////////////////////////////////////////////////////////////////////
121/// Compute bounding box of the torus.
122
124{
125 fDZ = fRmax;
127 fDX = fDY = fR+fRmax;
128 return;
129 }
130 Double_t xc[4];
131 Double_t yc[4];
140
141 Double_t xmin = xc[TMath::LocMin(4, &xc[0])];
142 Double_t xmax = xc[TMath::LocMax(4, &xc[0])];
143 Double_t ymin = yc[TMath::LocMin(4, &yc[0])];
144 Double_t ymax = yc[TMath::LocMax(4, &yc[0])];
145 Double_t ddp = -fPhi1;
146 if (ddp<0) ddp+= 360;
147 if (ddp<=fDphi) xmax = fR+fRmax;
148 ddp = 90-fPhi1;
149 if (ddp<0) ddp+= 360;
150 if (ddp>360) ddp-=360;
151 if (ddp<=fDphi) ymax = fR+fRmax;
152 ddp = 180-fPhi1;
153 if (ddp<0) ddp+= 360;
154 if (ddp>360) ddp-=360;
155 if (ddp<=fDphi) xmin = -(fR+fRmax);
156 ddp = 270-fPhi1;
157 if (ddp<0) ddp+= 360;
158 if (ddp>360) ddp-=360;
159 if (ddp<=fDphi) ymin = -(fR+fRmax);
160 fOrigin[0] = (xmax+xmin)/2;
161 fOrigin[1] = (ymax+ymin)/2;
162 fOrigin[2] = 0;
163 fDX = (xmax-xmin)/2;
164 fDY = (ymax-ymin)/2;
165}
166
167////////////////////////////////////////////////////////////////////////////////
168/// Compute normal to closest surface from POINT.
169
170void TGeoTorus::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
171{
172 Double_t phi = TMath::ATan2(point[1],point[0]);
173 if (fDphi<360) {
176 Double_t c1 = TMath::Cos(phi1);
177 Double_t s1 = TMath::Sin(phi1);
178 Double_t c2 = TMath::Cos(phi2);
179 Double_t s2 = TMath::Sin(phi2);
180
181 Double_t daxis = Daxis(point,dir,0);
182 if ((fRmax-daxis)>1E-5) {
183 if (TGeoShape::IsSameWithinTolerance(fRmin,0) || (daxis-fRmin)>1E-5) {
184 TGeoShape::NormalPhi(point,dir,norm,c1,s1,c2,s2);
185 return;
186 }
187 }
188 }
189 Double_t r0[3];
190 r0[0] = fR*TMath::Cos(phi);
191 r0[1] = fR*TMath::Sin(phi);
192 r0[2] = 0;
193 Double_t normsq = 0;
194 for (Int_t i=0; i<3; i++) {
195 norm[i] = point[i] - r0[i];
196 normsq += norm[i]*norm[i];
197 }
198
199 normsq = TMath::Sqrt(normsq);
200 norm[0] /= normsq;
201 norm[1] /= normsq;
202 norm[2] /= normsq;
203 if (dir[0]*norm[0]+dir[1]*norm[1]+dir[2]*norm[2] < 0) {
204 norm[0] = -norm[0];
205 norm[1] = -norm[1];
206 norm[2] = -norm[2];
207 }
208}
209
210////////////////////////////////////////////////////////////////////////////////
211/// Test if point is inside the torus.
212/// check phi range
213
215{
217 Double_t phi = TMath::ATan2(point[1], point[0]) * TMath::RadToDeg();
218 if (phi < 0) phi+=360.0;
219 Double_t ddp = phi-fPhi1;
220 if (ddp<0) ddp+=360.;
221 if (ddp>fDphi) return kFALSE;
222 }
223 //check radius
224 Double_t rxy = TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
225 Double_t radsq = (rxy-fR)*(rxy-fR) + point[2]*point[2];
226 if (radsq<fRmin*fRmin) return kFALSE;
227 if (radsq>fRmax*fRmax) return kFALSE;
228 return kTRUE;
229}
230
231////////////////////////////////////////////////////////////////////////////////
232/// Compute closest distance from point px,py to each vertex.
233
235{
237 Int_t numPoints = n*(n-1);
238 if (fRmin>0) numPoints *= 2;
239 else if (fDphi<360) numPoints += 2;
240 return ShapeDistancetoPrimitive(numPoints, px, py);
241}
242
243////////////////////////////////////////////////////////////////////////////////
244/// Computes distance to axis of the torus from point pt + t*dir;
245
247{
248 Double_t p[3];
249 for (Int_t i=0; i<3; i++) p[i] = pt[i]+t*dir[i];
250 Double_t rxy = TMath::Sqrt(p[0]*p[0]+p[1]*p[1]);
251 return TMath::Sqrt((rxy-fR)*(rxy-fR)+p[2]*p[2]);
252}
253
254////////////////////////////////////////////////////////////////////////////////
255/// Computes derivative w.r.t. t of the distance to axis of the torus from point pt + t*dir;
256
258{
259 Double_t p[3];
260 for (Int_t i=0; i<3; i++) p[i] = pt[i]+t*dir[i];
261 Double_t rxy = TMath::Sqrt(p[0]*p[0]+p[1]*p[1]);
262 if (rxy<1E-4) return ((p[2]*dir[2]-fR*TMath::Sqrt(dir[0]*dir[0]+dir[1]*dir[1]))/TMath::Sqrt(fR*fR+p[2]*p[2]));
263 Double_t d = TMath::Sqrt((rxy-fR)*(rxy-fR)+p[2]*p[2]);
264 if (TGeoShape::IsSameWithinTolerance(d,0)) return 0.;
265 Double_t dd = (p[0]*dir[0]+p[1]*dir[1]+p[2]*dir[2] - (p[0]*dir[0]+p[1]*dir[1])*fR/rxy)/d;
266 return dd;
267}
268
269////////////////////////////////////////////////////////////////////////////////
270/// Second derivative of distance to torus axis w.r.t t.
271
273{
274 Double_t p[3];
275 for (Int_t i=0; i<3; i++) p[i] = pt[i]+t*dir[i];
276 Double_t rxy = TMath::Sqrt(p[0]*p[0]+p[1]*p[1]);
277 if (rxy<1E-6) return 0;
278 Double_t daxis = TMath::Sqrt((rxy-fR)*(rxy-fR)+p[2]*p[2]);
279 if (TGeoShape::IsSameWithinTolerance(daxis,0)) return 0;
280 Double_t ddaxis = (p[0]*dir[0]+p[1]*dir[1]+p[2]*dir[2] - (p[0]*dir[0]+p[1]*dir[1])*fR/rxy)/daxis;
281 Double_t dddaxis = 1 - ddaxis*ddaxis - (1-dir[2]*dir[2])*fR/rxy +
282 fR*(p[0]*dir[0]+p[1]*dir[1])*(p[0]*dir[0]+p[1]*dir[1])/(rxy*rxy*rxy);
283 dddaxis /= daxis;
284 return dddaxis;
285}
286
287////////////////////////////////////////////////////////////////////////////////
288/// Compute distance from inside point to surface of the torus.
289
290Double_t TGeoTorus::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
291{
292 if (iact<3 && safe) {
293 *safe = Safety(point, kTRUE);
294 if (iact==0) return TGeoShape::Big();
295 if ((iact==1) && (step<=*safe)) return TGeoShape::Big();
296 }
298 Bool_t hasphi = (fDphi<360)?kTRUE:kFALSE;
299 Bool_t hasrmin = (fRmin>0)?kTRUE:kFALSE;
300 Double_t dout = ToBoundary(point,dir,fRmax,kTRUE);
301// Double_t dax = Daxis(point,dir,dout);
302 Double_t din = (hasrmin)?ToBoundary(point,dir,fRmin,kTRUE):TGeoShape::Big();
303 snext = TMath::Min(dout,din);
304 if (snext>1E10) return TGeoShape::Tolerance();
305 Double_t dphi = TGeoShape::Big();
306 if (hasphi) {
307 // Torus segment case.
308 Double_t c1,s1,c2,s2,cm,sm,cdfi;
311 c1=TMath::Cos(phi1);
312 s1=TMath::Sin(phi1);
313 c2=TMath::Cos(phi2);
314 s2=TMath::Sin(phi2);
315 Double_t fio=0.5*(phi1+phi2);
316 cm=TMath::Cos(fio);
317 sm=TMath::Sin(fio);
318 cdfi = TMath::Cos(0.5*(phi2-phi1));
319 dphi = TGeoTubeSeg::DistFromInsideS(point,dir,fR-fRmax,fR+fRmax, fRmax, c1,s1,c2,s2,cm,sm,cdfi);
320 Double_t daxis = Daxis(point,dir,dphi);
321 if (daxis>=fRmin+1.E-8 && daxis<=fRmax-1.E-8) snext=TMath::Min(snext,dphi);
322 }
323 return snext;
324}
325
326////////////////////////////////////////////////////////////////////////////////
327/// Compute distance from outside point to surface of the torus.
328
329Double_t TGeoTorus::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
330{
331 if (iact<3 && safe) {
332 *safe = Safety(point, kFALSE);
333 if (iact==0) return TGeoShape::Big();
334 if ((iact==1) && (step<=*safe)) return TGeoShape::Big();
335 }
336// Check if the bounding box is crossed within the requested distance
337 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
338 if (sdist>=step) return TGeoShape::Big();
339 Double_t daxis;
340 Bool_t hasphi = (fDphi<360)?kTRUE:kFALSE;
341// Bool_t hasrmin = (fRmin>0)?kTRUE:kFALSE;
342 Double_t c1=0,s1=0,c2=0,s2=0,cm=0,sm=0,cdfi=0;
343 Bool_t inphi = kFALSE;
344 Double_t phi, ddp, phi1,phi2,fio;
345 Double_t rxy2,dd;
347 Double_t pt[3];
348 Int_t i;
349
350 if (hasphi) {
351 // Torus segment case.
352 phi=TMath::ATan2(point[1], point[0])*TMath::RadToDeg();;
353 if (phi<0) phi+=360;
354 ddp = phi-fPhi1;
355 if (ddp<0) ddp+=360;;
356 if (ddp<=fDphi) inphi=kTRUE;
357 phi1=fPhi1*TMath::DegToRad();
358 phi2=(fPhi1+fDphi)*TMath::DegToRad();
359 c1=TMath::Cos(phi1);
360 s1=TMath::Sin(phi1);
361 c2=TMath::Cos(phi2);
362 s2=TMath::Sin(phi2);
363 fio=0.5*(phi1+phi2);
364 cm=TMath::Cos(fio);
365 sm=TMath::Sin(fio);
366 cdfi=TMath::Cos(0.5*(phi2-phi1));
367 }
368 // Check if we are inside or outside the bounding ring.
369 Bool_t inbring = kFALSE;
370 if (TMath::Abs(point[2]) <= fRmax) {
371 rxy2 = point[0]*point[0]+point[1]*point[1];
372 if ((rxy2>=(fR-fRmax)*(fR-fRmax)) && (rxy2<=(fR+fRmax)*(fR+fRmax))) {
373 if (!hasphi || inphi) inbring=kTRUE;
374 }
375 }
376
377 // If outside the ring, compute distance to it.
378 Double_t dring = TGeoShape::Big();
379 Double_t eps = 1.E-8;
380 snext = 0;
381 daxis = -1;
382 memcpy(pt,point,3*sizeof(Double_t));
383 if (!inbring) {
384 if (hasphi) dring = TGeoTubeSeg::DistFromOutsideS(point,dir,TMath::Max(0.,fR-fRmax-eps),fR+fRmax+eps, fRmax+eps, c1,s1,c2,s2,cm,sm,cdfi);
385 else dring = TGeoTube::DistFromOutsideS(point,dir,TMath::Max(0.,fR-fRmax-eps),fR+fRmax+eps, fRmax+eps);
386 // If not crossing it, return BIG.
387 if (dring>1E10) return TGeoShape::Big();
388 snext = dring;
389 // Check if the crossing is due to phi.
390 daxis = Daxis(point,dir,snext);
391 if (daxis>=fRmin && daxis<fRmax) return snext;
392 // Not a phi crossing -> propagate until we cross the ring.
393 for (i=0; i<3; i++) pt[i] = point[i]+snext*dir[i];
394 }
395 // Point pt is inside the bounding ring, no phi crossing so far.
396 // Check if we are in the hole.
397 if (daxis<0) daxis = Daxis(pt,dir,0);
398 if (daxis<fRmin+1.E-8) {
399 // We are in the hole. Check if we came from outside.
400 if (snext>0) {
401 // we can cross either the inner torus or exit the other hole.
402 snext += 0.1*eps;
403 for (i=0; i<3; i++) pt[i] += 0.1*eps*dir[i];
404 }
405 // We are in the hole from the beginning.
406 // find first crossing with inner torus
407 dd = ToBoundary(pt,dir, fRmin,kFALSE);
408 // find exit distance from inner bounding ring
409 if (hasphi) dring = TGeoTubeSeg::DistFromInsideS(pt,dir,fR-fRmin,fR+fRmin, fRmin, c1,s1,c2,s2,cm,sm,cdfi);
410 else dring = TGeoTube::DistFromInsideS(pt,dir,fR-fRmin,fR+fRmin, fRmin);
411 if (dd<dring) return (snext+dd);
412 // we were exiting a hole inside phi hole
413 snext += dring+ eps;
414 for (i=0; i<3; i++) pt[i] = point[i] + snext*dir[i];
415 snext += DistFromOutside(pt,dir,3);
416 return snext;
417 }
418 // We are inside the outer ring, having daxis>fRmax
419 // Compute distance to exit the bounding ring (again)
420 if (snext>0) {
421 // we can cross either the inner torus or exit the other hole.
422 snext += 0.1*eps;
423 for (i=0; i<3; i++) pt[i] += 0.1*eps*dir[i];
424 }
425 // Check intersection with outer torus
426 dd = ToBoundary(pt, dir, fRmax, kFALSE);
427 if (hasphi) dring = TGeoTubeSeg::DistFromInsideS(pt,dir,TMath::Max(0.,fR-fRmax-eps),fR+fRmax+eps, fRmax+eps, c1,s1,c2,s2,cm,sm,cdfi);
428 else dring = TGeoTube::DistFromInsideS(pt,dir,TMath::Max(0.,fR-fRmax-eps),fR+fRmax+eps, fRmax+eps);
429 if (dd<dring) {
430 snext += dd;
431 return snext;
432 }
433 // We are exiting the bounding ring before crossing the torus -> propagate
434 snext += dring+eps;
435 for (i=0; i<3; i++) pt[i] = point[i] + snext*dir[i];
436 snext += DistFromOutside(pt,dir,3);
437 return snext;
438}
439
440////////////////////////////////////////////////////////////////////////////////
441/// Divide this torus shape belonging to volume "voldiv" into ndiv volumes
442/// called divname, from start position with the given step.
443
444TGeoVolume *TGeoTorus::Divide(TGeoVolume * /*voldiv*/, const char * /*divname*/, Int_t /*iaxis*/, Int_t /*ndiv*/,
445 Double_t /*start*/, Double_t /*step*/)
446{
447 return 0;
448}
449
450////////////////////////////////////////////////////////////////////////////////
451/// Returns name of axis IAXIS.
452
453const char *TGeoTorus::GetAxisName(Int_t iaxis) const
454{
455 switch (iaxis) {
456 case 1:
457 return "R";
458 case 2:
459 return "PHI";
460 case 3:
461 return "Z";
462 default:
463 return "UNDEFINED";
464 }
465}
466
467////////////////////////////////////////////////////////////////////////////////
468/// Get range of shape for a given axis.
469
471{
472 xlo = 0;
473 xhi = 0;
474 Double_t dx = 0;
475 switch (iaxis) {
476 case 1:
477 xlo = fRmin;
478 xhi = fRmax;
479 dx = xhi-xlo;
480 return dx;
481 case 2:
482 xlo = fPhi1;
483 xhi = fPhi1+fDphi;
484 dx = fDphi;
485 return dx;
486 case 3:
487 dx = 0;
488 return dx;
489 }
490 return dx;
491}
492
493////////////////////////////////////////////////////////////////////////////////
494/// Fill vector param[4] with the bounding cylinder parameters. The order
495/// is the following : Rmin, Rmax, Phi1, Phi2, dZ
496
498{
499 param[0] = (fR-fRmax); // Rmin
500 param[1] = (fR+fRmax); // Rmax
501 param[2] = fPhi1; // Phi1
502 param[3] = fPhi1+fDphi; // Phi2
503}
504
505////////////////////////////////////////////////////////////////////////////////
506/// Create a shape fitting the mother.
507
509{
510 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
511 Error("GetMakeRuntimeShape", "parametrized toruses not supported");
512 return 0;
513}
514
515////////////////////////////////////////////////////////////////////////////////
516/// print shape parameters
517
519{
520 printf("*** Shape %s: TGeoTorus ***\n", GetName());
521 printf(" R = %11.5f\n", fR);
522 printf(" Rmin = %11.5f\n", fRmin);
523 printf(" Rmax = %11.5f\n", fRmax);
524 printf(" Phi1 = %11.5f\n", fPhi1);
525 printf(" Dphi = %11.5f\n", fDphi);
526 printf(" Bounding box:\n");
528}
529
530////////////////////////////////////////////////////////////////////////////////
531/// Creates a TBuffer3D describing *this* shape.
532/// Coordinates are in local reference frame.
533
535{
537 Int_t nbPnts = n*(n-1);
538 Bool_t hasrmin = (GetRmin()>0)?kTRUE:kFALSE;
539 Bool_t hasphi = (GetDphi()<360)?kTRUE:kFALSE;
540 if (hasrmin) nbPnts *= 2;
541 else if (hasphi) nbPnts += 2;
542
543 Int_t nbSegs = (2*n-1)*(n-1);
544 Int_t nbPols = (n-1)*(n-1);
545 if (hasrmin) {
546 nbSegs += (2*n-1)*(n-1);
547 nbPols += (n-1)*(n-1);
548 }
549 if (hasphi) {
550 nbSegs += 2*(n-1);
551 nbPols += 2*(n-1);
552 }
553
555 nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols);
556 if (buff)
557 {
558 SetPoints(buff->fPnts);
559 SetSegsAndPols(*buff);
560 }
561
562 return buff;
563}
564
565////////////////////////////////////////////////////////////////////////////////
566/// Fill TBuffer3D structure for segments and polygons.
567
569{
570 Int_t i, j;
572 Int_t nbPnts = n*(n-1);
573 Int_t indx, indp, startcap=0;
574 Bool_t hasrmin = (GetRmin()>0)?kTRUE:kFALSE;
575 Bool_t hasphi = (GetDphi()<360)?kTRUE:kFALSE;
576 if (hasrmin) nbPnts *= 2;
577 else if (hasphi) nbPnts += 2;
579
580 indp = n*(n-1); // start index for points on inner surface
581 memset(buff.fSegs, 0, buff.NbSegs()*3*sizeof(Int_t));
582
583 // outer surface phi circles = n*(n-1) -> [0, n*(n-1) -1]
584 // connect point j with point j+1 on same row
585 indx = 0;
586 for (i = 0; i < n; i++) { // rows [0,n-1]
587 for (j = 0; j < n-1; j++) { // points on a row [0, n-2]
588 buff.fSegs[indx+(i*(n-1)+j)*3] = c;
589 buff.fSegs[indx+(i*(n-1)+j)*3+1] = i*(n-1)+j; // j on row i
590 buff.fSegs[indx+(i*(n-1)+j)*3+2] = i*(n-1)+((j+1)%(n-1)); // j+1 on row i
591 }
592 }
593 indx += 3*n*(n-1);
594 // outer surface generators = (n-1)*(n-1) -> [n*(n-1), (2*n-1)*(n-1) -1]
595 // connect point j on row i with point j on row i+1
596 for (i = 0; i < n-1; i++) { // rows [0, n-2]
597 for (j = 0; j < n-1; j++) { // points on a row [0, n-2]
598 buff.fSegs[indx+(i*(n-1)+j)*3] = c;
599 buff.fSegs[indx+(i*(n-1)+j)*3+1] = i*(n-1)+j; // j on row i
600 buff.fSegs[indx+(i*(n-1)+j)*3+2] = (i+1)*(n-1)+j; // j on row i+1
601 }
602 }
603 indx += 3*(n-1)*(n-1);
604 startcap = (2*n-1)*(n-1);
605
606 if (hasrmin) {
607 // inner surface phi circles = n*(n-1) -> [(2*n-1)*(n-1), (3*n-1)*(n-1) -1]
608 // connect point j with point j+1 on same row
609 for (i = 0; i < n; i++) { // rows [0, n-1]
610 for (j = 0; j < n-1; j++) { // points on a row [0, n-2]
611 buff.fSegs[indx+(i*(n-1)+j)*3] = c; // lighter color
612 buff.fSegs[indx+(i*(n-1)+j)*3+1] = indp + i*(n-1)+j; // j on row i
613 buff.fSegs[indx+(i*(n-1)+j)*3+2] = indp + i*(n-1)+((j+1)%(n-1)); // j+1 on row i
614 }
615 }
616 indx += 3*n*(n-1);
617 // inner surface generators = (n-1)*n -> [(3*n-1)*(n-1), (4*n-2)*(n-1) -1]
618 // connect point j on row i with point j on row i+1
619 for (i = 0; i < n-1; i++) { // rows [0, n-2]
620 for (j = 0; j < n-1; j++) { // points on a row [0, n-2]
621 buff.fSegs[indx+(i*(n-1)+j)*3] = c; // lighter color
622 buff.fSegs[indx+(i*(n-1)+j)*3+1] = indp + i*(n-1)+j; // j on row i
623 buff.fSegs[indx+(i*(n-1)+j)*3+2] = indp + (i+1)*(n-1)+j; // j on row i+1
624 }
625 }
626 indx += 3*(n-1)*(n-1);
627 startcap = (4*n-2)*(n-1);
628 }
629
630 if (hasphi) {
631 if (hasrmin) {
632 // endcaps = 2*(n-1) -> [(4*n-2)*(n-1), 4*n*(n-1)-1]
633 i = 0;
634 for (j = 0; j < n-1; j++) {
635 buff.fSegs[indx+j*3] = c+1;
636 buff.fSegs[indx+j*3+1] = (n-1)*i+j; // outer j on row 0
637 buff.fSegs[indx+j*3+2] = indp+(n-1)*i+j; // inner j on row 0
638 }
639 indx += 3*(n-1);
640 i = n-1;
641 for (j = 0; j < n-1; j++) {
642 buff.fSegs[indx+j*3] = c+1;
643 buff.fSegs[indx+j*3+1] = (n-1)*i+j; // outer j on row n-1
644 buff.fSegs[indx+j*3+2] = indp+(n-1)*i+j; // inner j on row n-1
645 }
646 indx += 3*(n-1);
647 } else {
648 i = 0;
649 for (j = 0; j < n-1; j++) {
650 buff.fSegs[indx+j*3] = c+1;
651 buff.fSegs[indx+j*3+1] = (n-1)*i+j; // outer j on row 0
652 buff.fSegs[indx+j*3+2] = n*(n-1); // center of first endcap
653 }
654 indx += 3*(n-1);
655 i = n-1;
656 for (j = 0; j < n-1; j++) {
657 buff.fSegs[indx+j*3] = c+1;
658 buff.fSegs[indx+j*3+1] = (n-1)*i+j; // outer j on row n-1
659 buff.fSegs[indx+j*3+2] = n*(n-1)+1; // center of second endcap
660 }
661 indx += 3*(n-1);
662 }
663 }
664
665 indx = 0;
666 memset(buff.fPols, 0, buff.NbPols()*6*sizeof(Int_t));
667
668 // outer surface = (n-1)*(n-1) -> [0, (n-1)*(n-1)-1]
669 // normal pointing out
670 for (i=0; i<n-1; i++) {
671 for (j=0; j<n-1; j++) {
672 buff.fPols[indx++] = c;
673 buff.fPols[indx++] = 4;
674 buff.fPols[indx++] = n*(n-1)+(n-1)*i+((j+1)%(n-1)); // generator j+1 on outer row i
675 buff.fPols[indx++] = (n-1)*(i+1)+j; // seg j on outer row i+1
676 buff.fPols[indx++] = n*(n-1)+(n-1)*i+j; // generator j on outer row i
677 buff.fPols[indx++] = (n-1)*i+j; // seg j on outer row i
678 }
679 }
680 if (hasrmin) {
681 indp = (2*n-1)*(n-1); // start index of inner segments
682 // inner surface = (n-1)*(n-1) -> [(n-1)*(n-1), 2*(n-1)*(n-1)-1]
683 // normal pointing out
684 for (i=0; i<n-1; i++) {
685 for (j=0; j<n-1; j++) {
686 buff.fPols[indx++] = c;
687 buff.fPols[indx++] = 4;
688 buff.fPols[indx++] = indp+n*(n-1)+(n-1)*i+j; // generator j on inner row i
689 buff.fPols[indx++] = indp+(n-1)*(i+1)+j; // seg j on inner row i+1
690 buff.fPols[indx++] = indp+n*(n-1)+(n-1)*i+((j+1)%(n-1)); // generator j+1 on inner r>
691 buff.fPols[indx++] = indp+(n-1)*i+j; // seg j on inner row i
692 }
693 }
694 }
695 if (hasphi) {
696 // endcaps = 2*(n-1) -> [2*(n-1)*(n-1), 2*n*(n-1)-1]
697 i=0; // row 0
698 Int_t np = (hasrmin)?4:3;
699 for (j=0; j<n-1; j++) {
700 buff.fPols[indx++] = c+1;
701 buff.fPols[indx++] = np;
702 buff.fPols[indx++] = j; // seg j on outer row 0 a
703 buff.fPols[indx++] = startcap+j; // endcap j on row 0 d
704 if(hasrmin)
705 buff.fPols[indx++] = indp+j; // seg j on inner row 0 c
706 buff.fPols[indx++] = startcap+((j+1)%(n-1)); // endcap j+1 on row 0 b
707 }
708
709 i=n-1; // row n-1
710 for (j=0; j<n-1; j++) {
711 buff.fPols[indx++] = c+1;
712 buff.fPols[indx++] = np;
713 buff.fPols[indx++] = (n-1)*i+j; // seg j on outer row n-1 a
714 buff.fPols[indx++] = startcap+(n-1)+((j+1)%(n-1)); // endcap j+1 on row n-1 d
715 if (hasrmin)
716 buff.fPols[indx++] = indp+(n-1)*i+j; // seg j on inner row n-1 c
717 buff.fPols[indx++] = startcap+(n-1)+j; // endcap j on row n-1 b
718 }
719 }
720}
721
722////////////////////////////////////////////////////////////////////////////////
723/// computes the closest distance from given point to this shape, according
724/// to option. The matching point on the shape is stored in spoint.
725
727{
728 Double_t saf[2];
729 Int_t i;
730 Double_t rxy = TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
731 Double_t rad = TMath::Sqrt((rxy-fR)*(rxy-fR) + point[2]*point[2]);
732 saf[0] = rad-fRmin;
733 saf[1] = fRmax-rad;
735 if (in) return TMath::Min(saf[0],saf[1]);
736 for (i=0; i<2; i++) saf[i]=-saf[i];
737 return TMath::Max(saf[0], saf[1]);
738 }
739
740 Double_t safphi = TGeoShape::SafetyPhi(point,in,fPhi1, fPhi1+fDphi);
741 Double_t safe = TGeoShape::Big();
742 if (in) {
743 safe = TMath::Min(saf[0], saf[1]);
744 return TMath::Min(safe, safphi);
745 }
746 for (i=0; i<2; i++) saf[i]=-saf[i];
747 safe = TMath::Max(saf[0], saf[1]);
748 return TMath::Max(safe, safphi);
749}
750
751////////////////////////////////////////////////////////////////////////////////
752/// Save a primitive as a C++ statement(s) on output stream "out".
753
754void TGeoTorus::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
755{
757 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
758 out << " r = " << fR << ";" << std::endl;
759 out << " rmin = " << fRmin << ";" << std::endl;
760 out << " rmax = " << fRmax << ";" << std::endl;
761 out << " phi1 = " << fPhi1 << ";" << std::endl;
762 out << " dphi = " << fDphi << ";" << std::endl;
763 out << " TGeoShape *" << GetPointerName() << " = new TGeoTorus(\"" << GetName() << "\",r,rmin,rmax,phi1,dphi);" << std::endl;
765}
766
767////////////////////////////////////////////////////////////////////////////////
768/// Set torus dimensions.
769
771 Double_t phi1, Double_t dphi)
772{
773 fR = r;
774 fRmin = rmin;
775 fRmax = rmax;
776 fPhi1 = phi1;
777 if (fPhi1<0) fPhi1+=360.;
778 fDphi = dphi;
779}
780
781////////////////////////////////////////////////////////////////////////////////
782/// Set torus dimensions starting from a list.
783
785{
786 SetTorusDimensions(param[0], param[1], param[2], param[3], param[4]);
787}
788
789////////////////////////////////////////////////////////////////////////////////
790/// Create torus mesh points
791
793{
794 if (!points) return;
796 Double_t phin, phout;
797 Double_t dpin = 360./(n-1);
798 Double_t dpout = fDphi/(n-1);
799 Double_t co,so,ci,si;
801 Int_t i,j;
802 Int_t indx = 0;
803 // loop outer mesh -> n*(n-1) points [0, 3*n*(n-1)-1]
804 for (i=0; i<n; i++) {
805 phout = (fPhi1+i*dpout)*TMath::DegToRad();
806 co = TMath::Cos(phout);
807 so = TMath::Sin(phout);
808 for (j=0; j<n-1; j++) {
809 phin = j*dpin*TMath::DegToRad();
810 ci = TMath::Cos(phin);
811 si = TMath::Sin(phin);
812 points[indx++] = (fR+fRmax*ci)*co;
813 points[indx++] = (fR+fRmax*ci)*so;
814 points[indx++] = fRmax*si;
815 }
816 }
817
818 if (havermin) {
819 // loop inner mesh -> n*(n-1) points [3*n*(n-1), 6*n*(n-1)]
820 for (i=0; i<n; i++) {
821 phout = (fPhi1+i*dpout)*TMath::DegToRad();
822 co = TMath::Cos(phout);
823 so = TMath::Sin(phout);
824 for (j=0; j<n-1; j++) {
825 phin = j*dpin*TMath::DegToRad();
826 ci = TMath::Cos(phin);
827 si = TMath::Sin(phin);
828 points[indx++] = (fR+fRmin*ci)*co;
829 points[indx++] = (fR+fRmin*ci)*so;
830 points[indx++] = fRmin*si;
831 }
832 }
833 } else {
834 if (fDphi<360.) {
835 // just add extra 2 points on the centers of the 2 phi cuts [3*n*n, 3*n*n+1]
838 points[indx++] = 0;
841 points[indx++] = 0;
842 }
843 }
844}
845
846////////////////////////////////////////////////////////////////////////////////
847/// Create torus mesh points
848
850{
851 if (!points) return;
853 Double_t phin, phout;
854 Double_t dpin = 360./(n-1);
855 Double_t dpout = fDphi/(n-1);
856 Double_t co,so,ci,si;
858 Int_t i,j;
859 Int_t indx = 0;
860 // loop outer mesh -> n*(n-1) points [0, 3*n*(n-1)-1]
861 // plane i = 0, n-1 point j = 0, n-1 ipoint = n*i + j
862 for (i=0; i<n; i++) {
863 phout = (fPhi1+i*dpout)*TMath::DegToRad();
864 co = TMath::Cos(phout);
865 so = TMath::Sin(phout);
866 for (j=0; j<n-1; j++) {
867 phin = j*dpin*TMath::DegToRad();
868 ci = TMath::Cos(phin);
869 si = TMath::Sin(phin);
870 points[indx++] = (fR+fRmax*ci)*co;
871 points[indx++] = (fR+fRmax*ci)*so;
872 points[indx++] = fRmax*si;
873 }
874 }
875
876 if (havermin) {
877 // loop inner mesh -> n*(n-1) points [3*n*(n-1), 6*n*(n-1)]
878 // plane i = 0, n-1 point j = 0, n-1 ipoint = n*n + n*i + j
879 for (i=0; i<n; i++) {
880 phout = (fPhi1+i*dpout)*TMath::DegToRad();
881 co = TMath::Cos(phout);
882 so = TMath::Sin(phout);
883 for (j=0; j<n-1; j++) {
884 phin = j*dpin*TMath::DegToRad();
885 ci = TMath::Cos(phin);
886 si = TMath::Sin(phin);
887 points[indx++] = (fR+fRmin*ci)*co;
888 points[indx++] = (fR+fRmin*ci)*so;
889 points[indx++] = fRmin*si;
890 }
891 }
892 } else {
893 if (fDphi<360.) {
894 // just add extra 2 points on the centers of the 2 phi cuts [n*n, n*n+1]
895 // ip1 = n*(n-1) + 0;
896 // ip2 = n*(n-1) + 1
899 points[indx++] = 0;
902 points[indx++] = 0;
903 }
904 }
905}
906
907////////////////////////////////////////////////////////////////////////////////
908/// Return number of vertices of the mesh representation
909
911{
913 Int_t numPoints = n*(n-1);
914 if (fRmin>TGeoShape::Tolerance()) numPoints *= 2;
915 else if (fDphi<360.) numPoints += 2;
916 return numPoints;
917}
918
919////////////////////////////////////////////////////////////////////////////////
920/// fill size of this 3-D object
921
923{
924}
925
926////////////////////////////////////////////////////////////////////////////////
927/// Find real solutions of the cubic equation : x^3 + a*x^2 + b*x + c = 0
928/// Input: a,b,c
929/// Output: x[3] real solutions
930/// Returns number of real solutions (1 or 3)
931
933{
934 const Double_t ott = 1./3.;
935 const Double_t sq3 = TMath::Sqrt(3.);
936 Int_t ireal = 1;
937 Double_t p = b-a*a*ott;
938 Double_t q = c-a*b*ott+2.*a*a*a*ott*ott*ott;
939 Double_t delta = 4*p*p*p+27*q*q;
940// Double_t y1r, y1i, y2r, y2i;
941 Double_t t,u;
942 if (delta>=0) {
943 delta = TMath::Sqrt(delta);
944 t = (-3*q*sq3+delta)/(6*sq3);
945 u = (3*q*sq3+delta)/(6*sq3);
946 x[0] = TMath::Sign(1.,t)*TMath::Power(TMath::Abs(t),ott)-
947 TMath::Sign(1.,u)*TMath::Power(TMath::Abs(u),ott)-a*ott;
948 } else {
949 delta = TMath::Sqrt(-delta);
950 t = -0.5*q;
951 u = delta/(6*sq3);
952 x[0] = 2.*TMath::Power(t*t+u*u,0.5*ott) * TMath::Cos(ott*TMath::ATan2(u,t));
953 x[0] -= a*ott;
954 }
955
956 t = x[0]*x[0]+a*x[0]+b;
957 u = a+x[0];
958 delta = u*u-4.*t;
959 if (delta>=0) {
960 ireal = 3;
961 delta = TMath::Sqrt(delta);
962 x[1] = 0.5*(-u-delta);
963 x[2] = 0.5*(-u+delta);
964 }
965 return ireal;
966}
967
968////////////////////////////////////////////////////////////////////////////////
969/// Find real solutions of the quartic equation : x^4 + a*x^3 + b*x^2 + c*x + d = 0
970/// Input: a,b,c,d
971/// Output: x[4] - real solutions
972/// Returns number of real solutions (0 to 3)
973
975{
976 Double_t e = b-3.*a*a/8.;
977 Double_t f = c+a*a*a/8.-0.5*a*b;
978 Double_t g = d-3.*a*a*a*a/256. + a*a*b/16. - a*c/4.;
979 Double_t xx[4];
980 Int_t ind[4];
981 Double_t delta;
982 Double_t h=0;
983 Int_t ireal = 0;
984 Int_t i;
986 delta = e*e-4.*g;
987 if (delta<0) return 0;
988 delta = TMath::Sqrt(delta);
989 h = 0.5*(-e-delta);
990 if (h>=0) {
991 h = TMath::Sqrt(h);
992 x[ireal++] = -h-0.25*a;
993 x[ireal++] = h-0.25*a;
994 }
995 h = 0.5*(-e+delta);
996 if (h>=0) {
997 h = TMath::Sqrt(h);
998 x[ireal++] = -h-0.25*a;
999 x[ireal++] = h-0.25*a;
1000 }
1001 if (ireal>0) {
1002 TMath::Sort(ireal, x, ind,kFALSE);
1003 for (i=0; i<ireal; i++) xx[i] = x[ind[i]];
1004 memcpy(x,xx,ireal*sizeof(Double_t));
1005 }
1006 return ireal;
1007 }
1008
1010 x[ireal++] = -0.25*a;
1011 ind[0] = SolveCubic(0,e,f,xx);
1012 for (i=0; i<ind[0]; i++) x[ireal++] = xx[i]-0.25*a;
1013 if (ireal>0) {
1014 TMath::Sort(ireal, x, ind,kFALSE);
1015 for (i=0; i<ireal; i++) xx[i] = x[ind[i]];
1016 memcpy(x,xx,ireal*sizeof(Double_t));
1017 }
1018 return ireal;
1019 }
1020
1021
1022 ireal = SolveCubic(2.*e, e*e-4.*g, -f*f, xx);
1023 if (ireal==1) {
1024 if (xx[0]<=0) return 0;
1025 h = TMath::Sqrt(xx[0]);
1026 } else {
1027 // 3 real solutions of the cubic
1028 for (i=0; i<3; i++) {
1029 h = xx[i];
1030 if (h>=0) break;
1031 }
1032 if (h<=0) return 0;
1033 h = TMath::Sqrt(h);
1034 }
1035 Double_t j = 0.5*(e+h*h-f/h);
1036 ireal = 0;
1037 delta = h*h-4.*j;
1038 if (delta>=0) {
1039 delta = TMath::Sqrt(delta);
1040 x[ireal++] = 0.5*(-h-delta)-0.25*a;
1041 x[ireal++] = 0.5*(-h+delta)-0.25*a;
1042 }
1043 delta = h*h-4.*g/j;
1044 if (delta>=0) {
1045 delta = TMath::Sqrt(delta);
1046 x[ireal++] = 0.5*(h-delta)-0.25*a;
1047 x[ireal++] = 0.5*(h+delta)-0.25*a;
1048 }
1049 if (ireal>0) {
1050 TMath::Sort(ireal, x, ind,kFALSE);
1051 for (i=0; i<ireal; i++) xx[i] = x[ind[i]];
1052 memcpy(x,xx,ireal*sizeof(Double_t));
1053 }
1054 return ireal;
1055}
1056
1057////////////////////////////////////////////////////////////////////////////////
1058/// Returns distance to the surface or the torus (fR,r) from a point, along
1059/// a direction. Point is close enough to the boundary so that the distance
1060/// to the torus is decreasing while moving along the given direction.
1061
1063{
1064 // Compute coefficients of the quartic
1067 Double_t r0sq = pt[0]*pt[0]+pt[1]*pt[1]+pt[2]*pt[2];
1068 Double_t rdotn = pt[0]*dir[0]+pt[1]*dir[1]+pt[2]*dir[2];
1069 Double_t rsumsq = fR*fR+r*r;
1070 Double_t a = 4.*rdotn;
1071 Double_t b = 2.*(r0sq+2.*rdotn*rdotn-rsumsq+2.*fR*fR*dir[2]*dir[2]);
1072 Double_t c = 4.*(r0sq*rdotn-rsumsq*rdotn+2.*fR*fR*pt[2]*dir[2]);
1073 Double_t d = r0sq*r0sq-2.*r0sq*rsumsq+4.*fR*fR*pt[2]*pt[2]+(fR*fR-r*r)*(fR*fR-r*r);
1074
1075 Double_t x[4],y[4];
1076 Int_t nsol = 0;
1077
1078 if (TMath::Abs(dir[2])<1E-3 && TMath::Abs(pt[2])<0.1*r) {
1079 Double_t r0 = fR - TMath::Sqrt((r-pt[2])*(r+pt[2]));
1080 Double_t b0 = (pt[0]*dir[0]+pt[1]*dir[1])/(dir[0]*dir[0]+dir[1]*dir[1]);
1081 Double_t c0 = (pt[0]*pt[0] + (pt[1]-r0)*(pt[1]+r0))/(dir[0]*dir[0]+dir[1]*dir[1]);
1082 Double_t delta = b0*b0-c0;
1083 if (delta>0) {
1084 y[nsol] = -b0-TMath::Sqrt(delta);
1085 if (y[nsol]>-tol) nsol++;
1086 y[nsol] = -b0+TMath::Sqrt(delta);
1087 if (y[nsol]>-tol) nsol++;
1088 }
1089 r0 = fR + TMath::Sqrt((r-pt[2])*(r+pt[2]));
1090 c0 = (pt[0]*pt[0] + (pt[1]-r0)*(pt[1]+r0))/(dir[0]*dir[0]+dir[1]*dir[1]);
1091 delta = b0*b0-c0;
1092 if (delta>0) {
1093 y[nsol] = -b0-TMath::Sqrt(delta);
1094 if (y[nsol]>-tol) nsol++;
1095 y[nsol] = -b0+TMath::Sqrt(delta);
1096 if (y[nsol]>-tol) nsol++;
1097 }
1098 if (nsol) {
1099 // Sort solutions
1100 Int_t ind[4];
1101 TMath::Sort(nsol, y, ind,kFALSE);
1102 for (Int_t j=0; j<nsol; j++) x[j] = y[ind[j]];
1103 }
1104 } else {
1105 nsol = SolveQuartic(a,b,c,d,x);
1106 }
1107 if (!nsol) return TGeoShape::Big();
1108 // look for first positive solution
1109 Double_t phi, ndotd;
1110 Double_t r0[3], norm[3];
1112 for (Int_t i=0; i<nsol; i++) {
1113 if (x[i]<-10) continue;
1114 phi = TMath::ATan2(pt[1]+x[i]*dir[1],pt[0]+x[i]*dir[0]);
1115 r0[0] = fR*TMath::Cos(phi);
1116 r0[1] = fR*TMath::Sin(phi);
1117 r0[2] = 0;
1118 for (Int_t ipt=0; ipt<3; ipt++) norm[ipt] = pt[ipt]+x[i]*dir[ipt] - r0[ipt];
1119 ndotd = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
1120 if (inner^in) {
1121 if (ndotd<0) continue;
1122 } else {
1123 if (ndotd>0) continue;
1124 }
1125 s = x[i];
1126 Double_t eps = TGeoShape::Big();
1127 Double_t delta = s*s*s*s + a*s*s*s + b*s*s + c*s + d;
1128 Double_t eps0 = -delta/(4.*s*s*s + 3.*a*s*s + 2.*b*s + c);
1129 while (TMath::Abs(eps)>TGeoShape::Tolerance()) {
1130 if (TMath::Abs(eps0)>100) break;
1131 s += eps0;
1132 if (TMath::Abs(s+eps0)<TGeoShape::Tolerance()) break;
1133 delta = s*s*s*s + a*s*s*s + b*s*s + c*s + d;
1134 eps = -delta/(4.*s*s*s + 3.*a*s*s + 2.*b*s + c);
1135 if (TMath::Abs(eps)>TMath::Abs(eps0)) break;
1136 eps0 = eps;
1137 }
1138 if (s<-TGeoShape::Tolerance()) continue;
1139 return TMath::Max(0.,s);
1140 }
1141 return TGeoShape::Big();
1142}
1143
1144////////////////////////////////////////////////////////////////////////////////
1145/// Returns numbers of vertices, segments and polygons composing the shape mesh.
1146
1147void TGeoTorus::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
1148{
1150 nvert = n*(n-1);
1151 Bool_t hasrmin = (GetRmin()>0)?kTRUE:kFALSE;
1152 Bool_t hasphi = (GetDphi()<360)?kTRUE:kFALSE;
1153 if (hasrmin) nvert *= 2;
1154 else if (hasphi) nvert += 2;
1155 nsegs = (2*n-1)*(n-1);
1156 npols = (n-1)*(n-1);
1157 if (hasrmin) {
1158 nsegs += (2*n-1)*(n-1);
1159 npols += (n-1)*(n-1);
1160 }
1161 if (hasphi) {
1162 nsegs += 2*(n-1);
1163 npols += 2*(n-1);
1164 }
1165}
1166
1167////////////////////////////////////////////////////////////////////////////////
1168/// Fills a static 3D buffer and returns a reference.
1169
1170const TBuffer3D & TGeoTorus::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
1171{
1172 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
1173
1174 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
1175
1176 if (reqSections & TBuffer3D::kRawSizes) {
1178 Int_t nbPnts = n*(n-1);
1179 Bool_t hasrmin = (GetRmin()>0)?kTRUE:kFALSE;
1180 Bool_t hasphi = (GetDphi()<360)?kTRUE:kFALSE;
1181 if (hasrmin) nbPnts *= 2;
1182 else if (hasphi) nbPnts += 2;
1183
1184 Int_t nbSegs = (2*n-1)*(n-1);
1185 Int_t nbPols = (n-1)*(n-1);
1186 if (hasrmin) {
1187 nbSegs += (2*n-1)*(n-1);
1188 nbPols += (n-1)*(n-1);
1189 }
1190 if (hasphi) {
1191 nbSegs += 2*(n-1);
1192 nbPols += 2*(n-1);
1193 }
1194
1195 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
1197 }
1198 }
1199 // TODO: Push down to TGeoShape?? But would have to do raw sizes set first..
1200 // can rest of TGeoShape be deferred until after
1201 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
1202 SetPoints(buffer.fPnts);
1203 if (!buffer.fLocalFrame) {
1204 TransformPoints(buffer.fPnts, buffer.NbPnts());
1205 }
1206
1207 SetSegsAndPols(buffer);
1209 }
1210
1211 return buffer;
1212}
1213
1214////////////////////////////////////////////////////////////////////////////////
1215/// Check the inside status for each of the points in the array.
1216/// Input: Array of point coordinates + vector size
1217/// Output: Array of Booleans for the inside of each point
1218
1219void TGeoTorus::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
1220{
1221 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
1222}
1223
1224////////////////////////////////////////////////////////////////////////////////
1225/// Compute the normal for an array o points so that norm.dot.dir is positive
1226/// Input: Arrays of point coordinates and directions + vector size
1227/// Output: Array of normal directions
1228
1229void TGeoTorus::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
1230{
1231 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
1232}
1233
1234////////////////////////////////////////////////////////////////////////////////
1235/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1236
1237void TGeoTorus::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1238{
1239 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
1240}
1241
1242////////////////////////////////////////////////////////////////////////////////
1243/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1244
1245void TGeoTorus::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1246{
1247 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
1248}
1249
1250////////////////////////////////////////////////////////////////////////////////
1251/// Compute safe distance from each of the points in the input array.
1252/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1253/// Output: Safety values
1254
1255void TGeoTorus::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
1256{
1257 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
1258}
ROOT::R::TRInterface & r
Definition Object.C:4
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define g(i)
Definition RSha256.hxx:105
#define a(i)
Definition RSha256.hxx:99
#define s1(x)
Definition RSha256.hxx:91
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
const Bool_t kFALSE
Definition RtypesCore.h:92
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:91
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
float xmin
float * q
float ymin
float xmax
float ymax
point * points
Definition X3DBuffer.c:22
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:114
UInt_t NbPols() const
Definition TBuffer3D.h:82
UInt_t NbPnts() const
Definition TBuffer3D.h:80
UInt_t NbSegs() const
Definition TBuffer3D.h:81
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
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
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:429
virtual void InspectShape() const
Prints shape parameters.
Definition TGeoBBox.cxx:792
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...
Int_t GetNsegments() const
Get number of segments approximating circles.
Geometrical transformation package.
Definition TGeoMatrix.h:41
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.
static void NormalPhi(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Static method to compute normal to phi planes.
@ kGeoSavePrimitive
Definition TGeoShape.h:65
@ kGeoRunTimeShape
Definition TGeoShape.h:41
static Double_t Tolerance()
Definition TGeoShape.h:91
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:163
Torus segment class.
Definition TGeoTorus.h:18
Int_t SolveQuartic(Double_t a, Double_t b, Double_t c, Double_t d, Double_t *x) const
Find real solutions of the quartic equation : x^4 + a*x^3 + b*x^2 + c*x + d = 0 Input: a,...
virtual Bool_t Contains(const Double_t *point) const
Test if point is inside the torus.
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
virtual TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const
Create a shape fitting the mother.
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 torus.
Double_t ToBoundary(const Double_t *pt, const Double_t *dir, Double_t r, Bool_t in) const
Returns distance to the surface or the torus (fR,r) from a point, along a direction.
virtual void InspectShape() const
print shape parameters
Double_t DDDaxis(const Double_t *pt, const Double_t *dir, Double_t t) const
Second derivative of distance to torus axis w.r.t t.
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
virtual 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 torus.
Double_t Daxis(const Double_t *pt, const Double_t *dir, Double_t t) const
Computes distance to axis of the torus from point pt + t*dir;.
virtual void SetDimensions(Double_t *param)
Set torus dimensions starting from a list.
void SetTorusDimensions(Double_t r, Double_t rmin, Double_t rmax, Double_t phi1, Double_t dphi)
Set torus dimensions.
virtual void SetSegsAndPols(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons.
Double_t GetRmin() const
Definition TGeoTorus.h:70
Double_t fR
Definition TGeoTorus.h:21
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
virtual void SetPoints(Double_t *points) const
Create torus mesh points.
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
Double_t fRmax
Definition TGeoTorus.h:23
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Double_t fDphi
Definition TGeoTorus.h:25
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute closest distance from point px,py to each vertex.
Double_t DDaxis(const Double_t *pt, const Double_t *dir, Double_t t) const
Computes derivative w.r.t. t of the distance to axis of the torus from point pt + t*dir;.
virtual void ComputeBBox()
Compute bounding box of the torus.
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this torus shape belonging to volume "voldiv" into ndiv volumes called divname,...
virtual void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
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 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 fPhi1
Definition TGeoTorus.h:24
TGeoTorus()
Default constructor.
Definition TGeoTorus.cxx:58
virtual TBuffer3D * MakeBuffer3D() const
Creates a TBuffer3D describing this shape.
Int_t SolveCubic(Double_t a, Double_t b, Double_t c, Double_t *x) const
Find real solutions of the cubic equation : x^3 + a*x^2 + b*x + c = 0 Input: a,b,c Output: x[3] real ...
virtual const char * GetAxisName(Int_t iaxis) const
Returns name of axis IAXIS.
Double_t GetDphi() const
Definition TGeoTorus.h:73
Double_t fRmin
Definition TGeoTorus.h:22
virtual void Sizeof3D() const
fill size of this 3-D object
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
virtual void 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.
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 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:347
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:286
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
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:130
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:696
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:893
TPaveText * pt
return c1
Definition legend1.C:41
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
return c2
Definition legend2.C:14
Long64_t LocMin(Long64_t n, const T *a)
Return index of array with the minimum element.
Definition TMath.h:972
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:679
Long64_t LocMax(Long64_t n, const T *a)
Return index of array with the maximum element.
Definition TMath.h:1000
constexpr Double_t DegToRad()
Conversion from degree to radian:
Definition TMath.h:81
Double_t Sqrt(Double_t x)
Definition TMath.h:691
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition TMath.h:735
Short_t Min(Short_t a, Short_t b)
Definition TMathBase.h:180
Double_t Cos(Double_t)
Definition TMath.h:643
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Definition TMath.h:639
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Definition TMathBase.h:362
constexpr Double_t RadToDeg()
Conversion from radian to degree:
Definition TMath.h:73
Short_t Abs(Short_t d)
Definition TMathBase.h:120
#define snext(osub1, osub2)
Definition triangle.c:1167