Logo ROOT   6.08/07
Reference Guide
TGeoEltu.cxx
Go to the documentation of this file.
1 // @(#)root/geom:$Id$
2 // Author: Mihaela Gheata 05/06/02
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 TGeoEltu
13 \ingroup Geometry_classes
14 
15 Elliptical tube class. An elliptical tube has 3 parameters
16  - A - semi-axis of the ellipse along x
17  - B - semi-axis of the ellipse along y
18  - dz - half length in z
19 
20 Begin_Macro(source)
21 {
22  TCanvas *c = new TCanvas("c", "c",0,0,600,600);
23  new TGeoManager("eltu", "poza6");
24  TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
25  TGeoMedium *med = new TGeoMedium("MED",1,mat);
26  TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
27  gGeoManager->SetTopVolume(top);
28  TGeoVolume *vol = gGeoManager->MakeEltu("ELTU",med, 30,10,40);
29  top->AddNode(vol,1);
30  gGeoManager->CloseGeometry();
31  gGeoManager->SetNsegments(50);
32  top->Draw();
33  TView *view = gPad->GetView();
34  view->ShowAxis();
35 }
36 End_Macro
37 */
38 
39 
40 #include "Riostream.h"
41 
42 #include "TGeoManager.h"
43 #include "TGeoVolume.h"
44 #include "TGeoEltu.h"
45 #include "TBuffer3D.h"
46 #include "TBuffer3DTypes.h"
47 #include "TMath.h"
48 
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 /// Dummy constructor
53 
55 {
56  SetShapeBit(TGeoShape::kGeoEltu);
57 }
58 
59 ////////////////////////////////////////////////////////////////////////////////
60 /// Default constructor specifying X and Y semiaxis length
61 
63  :TGeoTube()
64 {
66  SetEltuDimensions(a, b, dz);
67  ComputeBBox();
68 }
69 
70 ////////////////////////////////////////////////////////////////////////////////
71 /// Default constructor specifying X and Y semiaxis length
72 
74  :TGeoTube(name,0.,b,dz)
75 {
76  SetName(name);
78  SetEltuDimensions(a, b, dz);
79  ComputeBBox();
80 }
81 
82 ////////////////////////////////////////////////////////////////////////////////
83 /// Default constructor specifying minimum and maximum radius
84 /// param[0] = A
85 /// param[1] = B
86 /// param[2] = dz
87 
89 {
91  SetDimensions(param);
92  ComputeBBox();
93 }
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 /// destructor
97 
99 {
100 }
101 
102 ////////////////////////////////////////////////////////////////////////////////
103 /// Computes capacity of the shape in [length^3]
104 
106 {
107  Double_t capacity = 2.*TMath::Pi()*fDz*fRmin*fRmax;
108  return capacity;
109 }
110 
111 ////////////////////////////////////////////////////////////////////////////////
112 /// compute bounding box of the tube
113 
115 {
116  fDX = fRmin;
117  fDY = fRmax;
118  fDZ = fDz;
119 }
120 
121 ////////////////////////////////////////////////////////////////////////////////
122 /// Compute normal to closest surface from POINT.
123 
124 void TGeoEltu::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
125 {
126  Double_t a = fRmin;
127  Double_t b = fRmax;
128  Double_t safr = TMath::Abs(TMath::Sqrt(point[0]*point[0]/(a*a)+point[1]*point[1]/(b*b))-1.);
129  safr *= TMath::Min(a,b);
130  Double_t safz = TMath::Abs(fDz-TMath::Abs(point[2]));
131  if (safz<safr) {
132  norm[0] = norm[1] = 0;
133  norm[2] = TMath::Sign(1.,dir[2]);
134  return;
135  }
136  norm[2] = 0.;
137  norm[0] = point[0]*b*b;
138  norm[1] = point[1]*a*a;
139  TMath::Normalize(norm);
140 }
141 
142 ////////////////////////////////////////////////////////////////////////////////
143 /// test if point is inside the elliptical tube
144 
145 Bool_t TGeoEltu::Contains(const Double_t *point) const
146 {
147  if (TMath::Abs(point[2]) > fDz) return kFALSE;
148  Double_t r2 = (point[0]*point[0])/(fRmin*fRmin)+(point[1]*point[1])/(fRmax*fRmax);
149  if (r2>1.) return kFALSE;
150  return kTRUE;
151 }
152 
153 ////////////////////////////////////////////////////////////////////////////////
154 /// compute closest distance from point px,py to each vertex
155 
157 {
159  const Int_t numPoints=4*n;
160  return ShapeDistancetoPrimitive(numPoints, px, py);
161 }
162 
163 ////////////////////////////////////////////////////////////////////////////////
164 /// compute distance from inside point to surface of the tube
165 
166 Double_t TGeoEltu::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
167 {
168  Double_t a2=fRmin*fRmin;
169  Double_t b2=fRmax*fRmax;
170  Double_t safz1=fDz-point[2];
171  Double_t safz2=fDz+point[2];
172 
173  if (iact<3 && safe) {
174  Double_t x0=TMath::Abs(point[0]);
175  Double_t y0=TMath::Abs(point[1]);
176  Double_t x1=x0;
177  Double_t y1=TMath::Sqrt((fRmin-x0)*(fRmin+x0))*fRmax/fRmin;
178  Double_t y2=y0;
179  Double_t x2=TMath::Sqrt((fRmax-y0)*(fRmax+y0))*fRmin/fRmax;
180  Double_t d1=(x1-x0)*(x1-x0)+(y1-y0)*(y1-y0);
181  Double_t d2=(x2-x0)*(x2-x0)+(y2-y0)*(y2-y0);
182  Double_t x3,y3;
183 
184  Double_t safr=TGeoShape::Big();
185  Double_t safz = TMath::Min(safz1,safz2);
186  for (Int_t i=0; i<8; i++) {
187  if (fRmax<fRmin) {
188  x3=0.5*(x1+x2);
189  y3=TMath::Sqrt((fRmin-x3)*(fRmin+x3))*fRmax/fRmin;;
190  } else {
191  y3=0.5*(y1+y2);
192  x3=TMath::Sqrt((fRmax-y3)*(fRmax+y3))*fRmin/fRmax;
193  }
194  if (d1<d2) {
195  x2=x3;
196  y2=y3;
197  d2=(x2-x0)*(x2-x0)+(y2-y0)*(y2-y0);
198  } else {
199  x1=x3;
200  y1=y3;
201  d1=(x1-x0)*(x1-x0)+(y1-y0)*(y1-y0);
202  }
203  }
204  safr=TMath::Sqrt(d1)-1.0E-3;
205  *safe = TMath::Min(safz, safr);
206  if (iact==0) return TGeoShape::Big();
207  if ((iact==1) && (*safe>step)) return TGeoShape::Big();
208  }
209  // compute distance to surface
210  // Do Z
211  Double_t snxt = TGeoShape::Big();
212  if (dir[2]>0) {
213  snxt=safz1/dir[2];
214  } else {
215  if (dir[2]<0) snxt=-safz2/dir[2];
216  }
217  Double_t sz = snxt;
218  Double_t xz=point[0]+dir[0]*sz;
219  Double_t yz=point[1]+dir[1]*sz;
220  if ((xz*xz/a2+yz*yz/b2)<=1) return snxt;
221  // do elliptical surface
222  Double_t tolerance = TGeoShape::Tolerance();
223  Double_t u=dir[0]*dir[0]*b2+dir[1]*dir[1]*a2;
224  Double_t v=point[0]*dir[0]*b2+point[1]*dir[1]*a2;
225  Double_t w=point[0]*point[0]*b2+point[1]*point[1]*a2-a2*b2;
226  Double_t d=v*v-u*w;
227  if (d<0 || TGeoShape::IsSameWithinTolerance(u,0)) return tolerance;
228  Double_t sd=TMath::Sqrt(d);
229  snxt = (-v+sd)/u;
230 
231  if (snxt<0) return tolerance;
232  return snxt;
233 }
234 
235 ////////////////////////////////////////////////////////////////////////////////
236 /// compute distance from outside point to surface of the tube and safe distance
237 
238 Double_t TGeoEltu::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
239 {
240  Double_t safz=TMath::Abs(point[2])-fDz;
241  Double_t a2=fRmin*fRmin;
242  Double_t b2=fRmax*fRmax;
243  if (iact<3 && safe) {
244  Double_t x0=TMath::Abs(point[0]);
245  Double_t y0=TMath::Abs(point[1]);
246  *safe=0.;
247  if ((x0*x0/a2+y0*y0/b2)>=1) {
248  Double_t phi1=0;
249  Double_t phi2=0.5*TMath::Pi();
250  Double_t phi3;
251  Double_t x3=0.,y3=0.,d;
252  for (Int_t i=0; i<10; i++) {
253  phi3=(phi1+phi2)*0.5;
254  x3=fRmin*TMath::Cos(phi3);
255  y3=fRmax*TMath::Sin(phi3);
256  d=y3*a2*(x0-x3)-x3*b2*(y0-y3);
257  if (d<0) phi1=phi3;
258  else phi2=phi3;
259  }
260  *safe=TMath::Sqrt((x0-x3)*(x0-x3)+(y0-y3)*(y0-y3));
261  }
262  if (safz>0) {
263  *safe=TMath::Sqrt((*safe)*(*safe)+safz*safz);
264  }
265  if (iact==0) return TGeoShape::Big();
266  if ((iact==1) && (step<*safe)) return TGeoShape::Big();
267  }
268  // compute vector distance
269  Double_t zi, tau;
270  Double_t epsil = 10.*TGeoShape::Tolerance();
271  if (safz > -epsil) {
272  // point beyond the z limit (up or down)
273  // Check if direction is outgoing
274  if (point[2]*dir[2]>0) return TGeoShape::Big();
275  // Check if direction is perpendicular to Z axis
276  if (TGeoShape::IsSameWithinTolerance(dir[2],0)) return TGeoShape::Big();
277  // select +z or -z depending on the side of the point
278  zi = (point[2] > 0) ? fDz : -fDz;
279  // Distance to zi plane position
280  tau = (zi-point[2])/dir[2];
281  // Extrapolated coordinates at the z position of the end plane.
282  Double_t xz=point[0]+dir[0]*tau;
283  Double_t yz=point[1]+dir[1]*tau;
284  if ((xz*xz/a2+yz*yz/b2)<1) return tau;
285  }
286 
287 // Check if the bounding box is crossed within the requested distance
288  Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
289  if (sdist>=step) return TGeoShape::Big();
290  Double_t u=dir[0]*dir[0]*b2+dir[1]*dir[1]*a2; // positive
292  Double_t v=point[0]*dir[0]*b2+point[1]*dir[1]*a2;
293  Double_t w=point[0]*point[0]*b2+point[1]*point[1]*a2-a2*b2;
294  Double_t d=v*v-u*w;
295  if (d<0) return TGeoShape::Big();
296  Double_t dsq=TMath::Sqrt(d);
297  // Biggest solution - if negative, or very close to boundary
298  // no crossing (just exiting, no re-entering possible)
299  tau = (-v+dsq)/u;
300  if (tau < epsil) return TGeoShape::Big();
301  // only entering crossing must be considered (smallest)
302  tau = (-v-dsq)/u;
303  zi=point[2]+tau*dir[2];
304  // If the crossing point is not in the Z range, there is no crossing
305  if ((TMath::Abs(zi)-fDz)>0) return TGeoShape::Big();
306  // crossing is backwards (point inside the ellipse) in Z range
307  if (tau < 0) return 0.;
308  // Point is outside and crossing the elliptical tube in Z range
309  return tau;
310 }
311 
312 ////////////////////////////////////////////////////////////////////////////////
313 /// Divide the shape along one axis.
314 
315 TGeoVolume *TGeoEltu::Divide(TGeoVolume * /*voldiv*/, const char * /*divname*/, Int_t /*iaxis*/, Int_t /*ndiv*/,
316  Double_t /*start*/, Double_t /*step*/)
317 {
318  Error("Divide", "Elliptical tubes divisions not implemented");
319  return 0;
320 }
321 
322 ////////////////////////////////////////////////////////////////////////////////
323 ///--- Fill vector param[4] with the bounding cylinder parameters. The order
324 /// is the following : Rmin, Rmax, Phi1, Phi2
325 
327 {
328  param[0] = 0.; // Rmin
329  param[1] = TMath::Max(fRmin, fRmax); // Rmax
330  param[1] *= param[1];
331  param[2] = 0.; // Phi1
332  param[3] = 360.; // Phi2
333 }
334 
335 ////////////////////////////////////////////////////////////////////////////////
336 /// in case shape has some negative parameters, these has to be computed
337 /// in order to fit the mother
338 
340 {
341  if (!TestShapeBit(kGeoRunTimeShape)) return 0;
342  if (!mother->TestShapeBit(kGeoEltu)) {
343  Error("GetMakeRuntimeShape", "invalid mother");
344  return 0;
345  }
346  Double_t a, b, dz;
347  a = fRmin;
348  b = fRmax;
349  dz = fDz;
350  if (fDz<0) dz=((TGeoEltu*)mother)->GetDz();
351  if (fRmin<0)
352  a = ((TGeoEltu*)mother)->GetA();
353  if (fRmax<0)
354  a = ((TGeoEltu*)mother)->GetB();
355 
356  return (new TGeoEltu(a, b, dz));
357 }
358 
359 ////////////////////////////////////////////////////////////////////////////////
360 /// print shape parameters
361 
363 {
364  printf("*** Shape %s: TGeoEltu ***\n", GetName());
365  printf(" A = %11.5f\n", fRmin);
366  printf(" B = %11.5f\n", fRmax);
367  printf(" dz = %11.5f\n", fDz);
368  printf(" Bounding box:\n");
370 }
371 
372 ////////////////////////////////////////////////////////////////////////////////
373 /// computes the closest distance from given point to this shape, according
374 /// to option. The matching point on the shape is stored in spoint.
375 
376 Double_t TGeoEltu::Safety(const Double_t *point, Bool_t /*in*/) const
377 {
378  Double_t x0 = TMath::Abs(point[0]);
379  Double_t y0 = TMath::Abs(point[1]);
380  Double_t x1, y1, dx, dy;
381  Double_t safr, safz;
382  safr = safz = TGeoShape::Big();
383  Double_t onepls = 1.+TGeoShape::Tolerance();
384  Double_t onemin = 1.-TGeoShape::Tolerance();
385  Double_t sqdist = x0*x0/(fRmin*fRmin)+y0*y0/(fRmax*fRmax);
386  Bool_t in = kTRUE;
387  if (sqdist>onepls) in = kFALSE;
388  else if (sqdist<onemin) in = kTRUE;
389  else return 0.;
390 
391  if (in) {
392  x1 = fRmin*TMath::Sqrt(1.-(y0*y0)/(fRmax*fRmax));
393  y1 = fRmax*TMath::Sqrt(1.-(x0*x0)/(fRmin*fRmin));
394  dx = x1-x0;
395  dy = y1-y0;
396  if (TMath::Abs(dx)<TGeoShape::Tolerance()) return 0;
397  safr = dx*dy/TMath::Sqrt(dx*dx+dy*dy);
398  safz = fDz - TMath::Abs(point[2]);
399  return TMath::Min(safr,safz);
400  }
401 
402  if (TMath::Abs(x0)<TGeoShape::Tolerance()) {
403  safr = y0 - fRmax;
404  } else {
405  if (TMath::Abs(y0)<TGeoShape::Tolerance()) {
406  safr = x0 - fRmin;
407  } else {
408  Double_t f = fRmin*fRmax/TMath::Sqrt(x0*x0*fRmax*fRmax+y0*y0*fRmin*fRmin);
409  x1 = f*x0;
410  y1 = f*y0;
411  dx = x0-x1;
412  dy = y0-y1;
413  Double_t ast = fRmin*y1/fRmax;
414  Double_t bct = fRmax*x1/fRmin;
415  Double_t d = TMath::Sqrt(bct*bct+ast*ast);
416  safr = (dx*bct+dy*ast)/d;
417  }
418  }
419  safz = TMath::Abs(point[2])-fDz;
420  return TMath::Max(safr, safz);
421 }
422 
423 ////////////////////////////////////////////////////////////////////////////////
424 /// Save a primitive as a C++ statement(s) on output stream "out".
425 
426 void TGeoEltu::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
427 {
428  if (TObject::TestBit(kGeoSavePrimitive)) return;
429  out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
430  out << " a = " << fRmin << ";" << std::endl;
431  out << " b = " << fRmax << ";" << std::endl;
432  out << " dz = " << fDz << ";" << std::endl;
433  out << " TGeoShape *" << GetPointerName() << " = new TGeoEltu(\"" << GetName() << "\",a,b,dz);" << std::endl;
435 }
436 
437 ////////////////////////////////////////////////////////////////////////////////
438 /// Set dimensions of the elliptical tube.
439 
441 {
442  if ((a<=0) || (b<0) || (dz<0)) {
444  }
445  fRmin=a;
446  fRmax=b;
447  fDz=dz;
448 }
449 
450 ////////////////////////////////////////////////////////////////////////////////
451 /// Set shape dimensions starting from an array.
452 
454 {
455  Double_t a = param[0];
456  Double_t b = param[1];
457  Double_t dz = param[2];
458  SetEltuDimensions(a, b, dz);
459 }
460 
461 ////////////////////////////////////////////////////////////////////////////////
462 /// Create elliptical tube mesh points
463 
465 {
466  Double_t dz;
467  Int_t j, n;
468 
469  n = gGeoManager->GetNsegments();
470  Double_t dphi = 360./n;
471  Double_t phi = 0;
472  Double_t cph,sph;
473  dz = fDz;
474 
475  Int_t indx = 0;
476  Double_t r2,r;
477  Double_t a2=fRmin*fRmin;
478  Double_t b2=fRmax*fRmax;
479 
480  if (points) {
481  for (j = 0; j < n; j++) {
482  points[indx+6*n] = points[indx] = 0;
483  indx++;
484  points[indx+6*n] = points[indx] = 0;
485  indx++;
486  points[indx+6*n] = dz;
487  points[indx] =-dz;
488  indx++;
489  }
490  for (j = 0; j < n; j++) {
491  phi = j*dphi*TMath::DegToRad();
492  sph=TMath::Sin(phi);
493  cph=TMath::Cos(phi);
494  r2=(a2*b2)/(b2+(a2-b2)*sph*sph);
495  r=TMath::Sqrt(r2);
496  points[indx+6*n] = points[indx] = r*cph;
497  indx++;
498  points[indx+6*n] = points[indx] = r*sph;
499  indx++;
500  points[indx+6*n]= dz;
501  points[indx] =-dz;
502  indx++;
503  }
504  }
505 }
506 
507 ////////////////////////////////////////////////////////////////////////////////
508 /// Returns numbers of vertices, segments and polygons composing the shape mesh.
509 
510 void TGeoEltu::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
511 {
512  TGeoTube::GetMeshNumbers(nvert,nsegs,npols);
513 }
514 
515 ////////////////////////////////////////////////////////////////////////////////
516 /// Returns the number of vertices on the mesh.
517 
519 {
521 }
522 
523 ////////////////////////////////////////////////////////////////////////////////
524 /// Create elliptical tube mesh points
525 
527 {
528  Double_t dz;
529  Int_t j, n;
530 
531  n = gGeoManager->GetNsegments();
532  Double_t dphi = 360./n;
533  Double_t phi = 0;
534  Double_t cph,sph;
535  dz = fDz;
536 
537  Int_t indx = 0;
538  Double_t r2,r;
539  Double_t a2=fRmin*fRmin;
540  Double_t b2=fRmax*fRmax;
541 
542  if (points) {
543  for (j = 0; j < n; j++) {
544  points[indx+6*n] = points[indx] = 0;
545  indx++;
546  points[indx+6*n] = points[indx] = 0;
547  indx++;
548  points[indx+6*n] = dz;
549  points[indx] =-dz;
550  indx++;
551  }
552  for (j = 0; j < n; j++) {
553  phi = j*dphi*TMath::DegToRad();
554  sph=TMath::Sin(phi);
555  cph=TMath::Cos(phi);
556  r2=(a2*b2)/(b2+(a2-b2)*sph*sph);
557  r=TMath::Sqrt(r2);
558  points[indx+6*n] = points[indx] = r*cph;
559  indx++;
560  points[indx+6*n] = points[indx] = r*sph;
561  indx++;
562  points[indx+6*n]= dz;
563  points[indx] =-dz;
564  indx++;
565  }
566  }
567 }
568 
569 ////////////////////////////////////////////////////////////////////////////////
570 /// Fills a static 3D buffer and returns a reference.
571 
572 const TBuffer3D & TGeoEltu::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
573 {
574  static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
575  TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
576 
577  if (reqSections & TBuffer3D::kRawSizes) {
579  Int_t nbPnts = 4*n;
580  Int_t nbSegs = 8*n;
581  Int_t nbPols = 4*n;
582  if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
583  buffer.SetSectionsValid(TBuffer3D::kRawSizes);
584  }
585  }
586  if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
587  SetPoints(buffer.fPnts);
588  if (!buffer.fLocalFrame) {
589  TransformPoints(buffer.fPnts, buffer.NbPnts());
590  }
591  SetSegsAndPols(buffer);
592  buffer.SetSectionsValid(TBuffer3D::kRaw);
593  }
594 
595  return buffer;
596 }
597 
598 ////////////////////////////////////////////////////////////////////////////////
599 /// Check the inside status for each of the points in the array.
600 /// Input: Array of point coordinates + vector size
601 /// Output: Array of Booleans for the inside of each point
602 
603 void TGeoEltu::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
604 {
605  for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
606 }
607 
608 ////////////////////////////////////////////////////////////////////////////////
609 /// Compute the normal for an array o points so that norm.dot.dir is positive
610 /// Input: Arrays of point coordinates and directions + vector size
611 /// Output: Array of normal directions
612 
613 void TGeoEltu::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
614 {
615  for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
616 }
617 
618 ////////////////////////////////////////////////////////////////////////////////
619 /// Compute distance from array of input points having directions specified by dirs. Store output in dists
620 
621 void TGeoEltu::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
622 {
623  for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
624 }
625 
626 ////////////////////////////////////////////////////////////////////////////////
627 /// Compute distance from array of input points having directions specified by dirs. Store output in dists
628 
629 void TGeoEltu::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
630 {
631  for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
632 }
633 
634 ////////////////////////////////////////////////////////////////////////////////
635 /// Compute safe distance from each of the points in the input array.
636 /// Input: Array of point coordinates, array of statuses for these points, size of the arrays
637 /// Output: Safety values
638 
639 void TGeoEltu::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
640 {
641  for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
642 }
virtual void SetDimensions(Double_t *param)
Set shape dimensions starting from an array.
Definition: TGeoEltu.cxx:453
virtual void InspectShape() const
print shape parameters
Definition: TGeoEltu.cxx:362
Cylindrical tube class.
Definition: TGeoTube.h:19
Float_t Normalize(Float_t v[3])
Normalize a vector v in place.
Definition: TMath.cxx:498
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: TGeoEltu.cxx:621
virtual void SetPoints(Double_t *points) const
Create elliptical tube mesh points.
Definition: TGeoEltu.cxx:464
T1 Sign(T1 a, T2 b)
Definition: TMathBase.h:155
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
Definition: TGeoEltu.cxx:124
float Float_t
Definition: RtypesCore.h:53
virtual void SetSegsAndPols(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons.
Definition: TGeoTube.cxx:672
const char Option_t
Definition: RtypesCore.h:62
Geometrical transformation package.
Definition: TGeoMatrix.h:40
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: TGeoEltu.cxx:376
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:157
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:131
Double_t DegToRad()
Definition: TMath.h:50
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition: TGeoVolume.h:61
virtual ~TGeoEltu()
destructor
Definition: TGeoEltu.cxx:98
UInt_t NbPnts() const
Definition: TBuffer3D.h:82
Double_t fOrigin[3]
Definition: TGeoBBox.h:26
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
virtual void GetBoundingCylinder(Double_t *param) const
— Fill vector param[4] with the bounding cylinder parameters.
Definition: TGeoEltu.cxx:326
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const
in case shape has some negative parameters, these has to be computed in order to fit the mother ...
Definition: TGeoEltu.cxx:339
virtual void InspectShape() const
Prints shape parameters.
Definition: TGeoBBox.cxx:793
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: TGeoEltu.cxx:603
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
TGeoEltu()
Dummy constructor.
Definition: TGeoEltu.cxx:54
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
Definition: TGeoShape.cxx:328
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:739
static Double_t Tolerance()
Definition: TGeoShape.h:93
static const double x2[5]
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:188
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
Definition: TGeoTube.cxx:1102
Double_t fDZ
Definition: TGeoBBox.h:25
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
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: TGeoEltu.cxx:510
Double_t * fPnts
Definition: TBuffer3D.h:114
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
compute distance from outside point to surface of the tube and safe distance
Definition: TGeoEltu.cxx:238
Double_t fDz
Definition: TGeoTube.h:25
void SetSectionsValid(UInt_t mask)
Definition: TBuffer3D.h:67
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition: TGeoEltu.cxx:426
const char * GetPointerName() const
Provide a pointer name containing uid.
Definition: TGeoShape.cxx:701
Bool_t fLocalFrame
Definition: TBuffer3D.h:92
point * points
Definition: X3DBuffer.c:20
virtual const char * GetName() const
Get the shape name.
Definition: TGeoShape.cxx:250
Int_t GetNsegments() const
Get number of segments approximating circles.
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
Definition: TGeoShape.cxx:261
Bool_t SectionsValid(UInt_t mask) const
Definition: TBuffer3D.h:69
Base abstract class for all shapes.
Definition: TGeoShape.h:27
TRandom2 r(17)
virtual Bool_t Contains(const Double_t *point) const
test if point is inside the elliptical tube
Definition: TGeoEltu.cxx:145
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: TGeoTube.cxx:1113
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide the shape along one axis.
Definition: TGeoEltu.cxx:315
SVector< double, 2 > v
Definition: Dict.h:5
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: TGeoEltu.cxx:629
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition: TGeoEltu.cxx:105
void SetEltuDimensions(Double_t a, Double_t b, Double_t dz)
Set dimensions of the elliptical tube.
Definition: TGeoEltu.cxx:440
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
Definition: TGeoShape.cxx:554
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:925
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
Generic 3D primitive description class.
Definition: TBuffer3D.h:19
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: TGeoEltu.cxx:613
Double_t Cos(Double_t)
Definition: TMath.h:424
Double_t Pi()
Definition: TMath.h:44
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:279
double f(double x)
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:554
double Double_t
Definition: RtypesCore.h:55
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
Bool_t TestShapeBit(UInt_t f) const
Definition: TGeoShape.h:165
static Double_t Big()
Definition: TGeoShape.h:90
Double_t fDY
Definition: TGeoBBox.h:24
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
Definition: TGeoShape.cxx:526
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
Double_t fRmin
Definition: TGeoTube.h:23
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute closest distance from point px,py to each vertex
Definition: TGeoEltu.cxx:156
Elliptical tube class.
Definition: TGeoEltu.h:19
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
Double_t Sin(Double_t)
Definition: TMath.h:421
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
compute distance from inside point to surface of the tube
Definition: TGeoEltu.cxx:166
Double_t fDX
Definition: TGeoBBox.h:23
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: TGeoEltu.cxx:639
virtual Int_t GetNmeshVertices() const
Returns the number of vertices on the mesh.
Definition: TGeoEltu.cxx:518
virtual void ComputeBBox()
compute bounding box of the tube
Definition: TGeoEltu.cxx:114
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
const Bool_t kTRUE
Definition: Rtypes.h:91
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
Definition: TGeoEltu.cxx:572
const Int_t n
Definition: legend1.C:16
unsigned int r2[N_CITIES]
Definition: simanTSP.cxx:322
char name[80]
Definition: TGX11.cxx:109
Double_t fRmax
Definition: TGeoTube.h:24
static const double x3[11]