ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TGeoHype.cxx
Go to the documentation of this file.
1 // @(#)root/geom:$Id$
2 // Author: Mihaela Gheata 20/11/04
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 
13 #include "Riostream.h"
14 
15 #include "TGeoManager.h"
16 #include "TGeoVolume.h"
17 #include "TVirtualGeoPainter.h"
18 #include "TGeoHype.h"
19 #include "TVirtualPad.h"
20 #include "TBuffer3D.h"
21 #include "TBuffer3DTypes.h"
22 #include "TMath.h"
23 
24 //_____________________________________________________________________________
25 // TGeoHype - Hyperboloid class defined by 5 parameters. Bounded by:
26 // - Two z planes at z=+/-dz
27 // - Inner and outer lateral surfaces. These represent the surfaces
28 // described by the revolution of 2 hyperbolas about the Z axis:
29 // r^2 - (t*z)^2 = a^2
30 //
31 // r = distance between hyperbola and Z axis at coordinate z
32 // t = tangent of the stereo angle (angle made by hyperbola
33 // asimptotic lines and Z axis). t=0 means cylindrical surface.
34 // a = distance between hyperbola and Z axis at z=0
35 //
36 // The inner hyperbolic surface is described by:
37 // r^2 - (tin*z)^2 = rin^2
38 // - absence of the inner surface (filled hyperboloid can be forced
39 // by rin=0 and sin=0
40 // The outer hyperbolic surface is described by:
41 // r^2 - (tout*z)^2 = rout^2
42 // TGeoHype parameters: dz[cm], rin[cm], sin[deg], rout[cm], sout[deg].
43 // MANDATORY conditions:
44 // - rin < rout
45 // - rout > 0
46 // - rin^2 + (tin*dz)^2 > rout^2 + (tout*dz)^2
47 // SUPPORTED CASES:
48 // - rin = 0, tin != 0 => inner surface conical
49 // - tin=0 AND/OR tout=0 => corresponding surface(s) cyllindrical
50 // e.g. tin=0 AND tout=0 => shape becomes a tube with: rmin,rmax,dz
51 //
52 //_____________________________________________________________________________
53 
54 
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 /// Default constructor
59 
61 {
62  SetShapeBit(TGeoShape::kGeoHype);
63  fStIn = 0.;
64  fStOut = 0.;
65  fTin = 0.;
66  fTinsq = 0.;
67  fTout = 0.;
68  fToutsq = 0.;
69 }
70 
71 
72 ////////////////////////////////////////////////////////////////////////////////
73 /// Constructor specifying hyperboloid parameters.
74 
76  :TGeoTube(rin, rout, dz)
77 {
79  SetHypeDimensions(rin, stin, rout, stout, dz);
80  // dz<0 can be used to force dz of hyperboloid fit the container volume
82  ComputeBBox();
83 }
84 ////////////////////////////////////////////////////////////////////////////////
85 /// Constructor specifying parameters and name.
86 
87 TGeoHype::TGeoHype(const char *name,Double_t rin, Double_t stin, Double_t rout, Double_t stout, Double_t dz)
88  :TGeoTube(name, rin, rout, dz)
89 {
91  SetHypeDimensions(rin, stin, rout, stout, dz);
92  // dz<0 can be used to force dz of hyperboloid fit the container volume
94  ComputeBBox();
95 }
96 
97 ////////////////////////////////////////////////////////////////////////////////
98 /// Default constructor specifying a list of parameters
99 /// param[0] = dz
100 /// param[1] = rin
101 /// param[2] = stin
102 /// param[3] = rout
103 /// param[4] = stout
104 
106  :TGeoTube(param[1],param[3],param[0])
107 {
109  SetDimensions(param);
110  // dz<0 can be used to force dz of hyperboloid fit the container volume
112  ComputeBBox();
113 }
114 
115 ////////////////////////////////////////////////////////////////////////////////
116 /// destructor
117 
119 {
120 }
121 
122 ////////////////////////////////////////////////////////////////////////////////
123 /// Computes capacity of the shape in [length^3]
124 
126 {
127  Double_t capacity = 2.*TMath::Pi()*fDz*(fRmax*fRmax-fRmin*fRmin) +
128  (2.*TMath::Pi()/3.)*fDz*fDz*fDz*(fToutsq-fTinsq);
129  return capacity;
130 }
131 
132 ////////////////////////////////////////////////////////////////////////////////
133 /// Compute bounding box of the hyperboloid
134 
136 {
137  if (fRmin<0.) {
138  Warning("ComputeBBox", "Shape %s has invalid rmin=%g ! SET TO 0.", GetName(),fRmin);
139  fRmin = 0.;
140  }
143  Error("ComputeBBox", "Shape %s hyperbolic surfaces are malformed: rin=%g, stin=%g, rout=%g, stout=%g",
144  GetName(), fRmin, fStIn, fRmax, fStOut);
145  return;
146  }
147 
149  fDZ = fDz;
150 }
151 
152 ////////////////////////////////////////////////////////////////////////////////
153 /// Compute normal to closest surface from POINT.
154 
156 {
157  Double_t saf[3];
158  Double_t rsq = point[0]*point[0]+point[1]*point[1];
159  Double_t r = TMath::Sqrt(rsq);
160  Double_t rin = (HasInner())?(TMath::Sqrt(RadiusHypeSq(point[2],kTRUE))):0.;
161  Double_t rout = TMath::Sqrt(RadiusHypeSq(point[2],kFALSE));
162  saf[0] = TMath::Abs(fDz-TMath::Abs(point[2]));
163  saf[1] = (HasInner())?TMath::Abs(rin-r):TGeoShape::Big();
164  saf[2] = TMath::Abs(rout-r);
165  Int_t i = TMath::LocMin(3,saf);
166  if (i==0 || r<1.E-10) {
167  norm[0] = norm[1] = 0.;
168  norm[2] = TMath::Sign(1.,dir[2]);
169  return;
170  }
171  Double_t t = (i==1)?fTinsq:fToutsq;;
172  t *= -point[2]/r;
173  Double_t ct = TMath::Sqrt(1./(1.+t*t));
174  Double_t st = t * ct;
175  Double_t phi = TMath::ATan2(point[1], point[0]);
176  Double_t cphi = TMath::Cos(phi);
177  Double_t sphi = TMath::Sin(phi);
178 
179  norm[0] = ct*cphi;
180  norm[1] = ct*sphi;
181  norm[2] = st;
182  if (norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2]<0) {
183  norm[0] = -norm[0];
184  norm[1] = -norm[1];
185  norm[2] = -norm[2];
186  }
187 }
188 
189 ////////////////////////////////////////////////////////////////////////////////
190 /// test if point is inside this tube
191 
192 Bool_t TGeoHype::Contains(const Double_t *point) const
193 {
194  if (TMath::Abs(point[2]) > fDz) return kFALSE;
195  Double_t r2 = point[0]*point[0]+point[1]*point[1];
196  Double_t routsq = RadiusHypeSq(point[2], kFALSE);
197  if (r2>routsq) return kFALSE;
198  if (!HasInner()) return kTRUE;
199  Double_t rinsq = RadiusHypeSq(point[2], kTRUE);
200  if (r2<rinsq) return kFALSE;
201  return kTRUE;
202 }
203 
204 ////////////////////////////////////////////////////////////////////////////////
205 /// compute closest distance from point px,py to each corner
206 
208 {
209  Int_t numPoints = GetNmeshVertices();
210  return ShapeDistancetoPrimitive(numPoints, px, py);
211 }
212 
213 ////////////////////////////////////////////////////////////////////////////////
214 /// Compute distance from inside point to surface of the hyperboloid.
215 
216 Double_t TGeoHype::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
217 {
218  if (iact<3 && safe) {
219  *safe = Safety(point, kTRUE);
220  if (iact==0) return TGeoShape::Big();
221  if ((iact==1) && (*safe>step)) return TGeoShape::Big();
222  }
223  // compute distance to surface
224  // Do Z
225  Double_t sz = TGeoShape::Big();
226  if (dir[2]>0) {
227  sz = (fDz-point[2])/dir[2];
228  if (sz<=0.) return 0.;
229  } else {
230  if (dir[2]<0) {
231  sz = -(fDz+point[2])/dir[2];
232  if (sz<=0.) return 0.;
233  }
234  }
235 
236 
237  // Do R
238  Double_t srin = TGeoShape::Big();
239  Double_t srout = TGeoShape::Big();
240  Double_t sr;
241  // inner and outer surfaces
242  Double_t s[2];
243  Int_t npos;
244  npos = DistToHype(point, dir, s, kTRUE, kTRUE);
245  if (npos) srin = s[0];
246  npos = DistToHype(point, dir, s, kFALSE, kTRUE);
247  if (npos) srout = s[0];
248  sr = TMath::Min(srin, srout);
249  return TMath::Min(sz,sr);
250 }
251 
252 
253 ////////////////////////////////////////////////////////////////////////////////
254 /// compute distance from outside point to surface of the hyperboloid.
255 
256 Double_t TGeoHype::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
257 {
258  if (iact<3 && safe) {
259  *safe = Safety(point, kFALSE);
260  if (iact==0) return TGeoShape::Big();
261  if ((iact==1) && (step<=*safe)) return TGeoShape::Big();
262  }
263 // Check if the bounding box is crossed within the requested distance
264  Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
265  if (sdist>=step) return TGeoShape::Big();
266  // find distance to shape
267  // Do Z
268  Double_t xi, yi, zi;
269  Double_t sz = TGeoShape::Big();
270  if (TMath::Abs(point[2])>=fDz) {
271  // We might find Z plane crossing
272  if ((point[2]*dir[2]) < 0) {
273  // Compute distance to Z (always positive)
274  sz = (TMath::Abs(point[2])-fDz)/TMath::Abs(dir[2]);
275  // Extrapolate
276  xi = point[0]+sz*dir[0];
277  yi = point[1]+sz*dir[1];
278  Double_t r2 = xi*xi + yi*yi;
279  Double_t rmin2 = RadiusHypeSq(fDz, kTRUE);
280  if (r2 >= rmin2) {
281  Double_t rmax2 = RadiusHypeSq(fDz, kFALSE);
282  if (r2 <= rmax2) return sz;
283  }
284  }
285  }
286  // We do not cross Z planes.
288  Double_t sout = TGeoShape::Big();
289  Double_t s[2];
290  Int_t npos;
291  npos = DistToHype(point, dir, s, kTRUE, kFALSE);
292  if (npos) {
293  zi = point[2] + s[0]*dir[2];
294  if (TMath::Abs(zi) <= fDz) sin = s[0];
295  else if (npos==2) {
296  zi = point[2] + s[1]*dir[2];
297  if (TMath::Abs(zi) <= fDz) sin = s[1];
298  }
299  }
300  npos = DistToHype(point, dir, s, kFALSE, kFALSE);
301  if (npos) {
302  zi = point[2] + s[0]*dir[2];
303  if (TMath::Abs(zi) <= fDz) sout = s[0];
304  else if (npos==2) {
305  zi = point[2] + s[1]*dir[2];
306  if (TMath::Abs(zi) <= fDz) sout = s[1];
307  }
308  }
309  return TMath::Min(sin, sout);
310 }
311 
312 ////////////////////////////////////////////////////////////////////////////////
313 /// Compute distance from an arbitrary point to inner/outer surface of hyperboloid.
314 /// Returns number of positive solutions. S[2] contains the solutions.
315 
316 Int_t TGeoHype::DistToHype(const Double_t *point, const Double_t *dir, Double_t *s, Bool_t inner, Bool_t in) const
317 {
318  Double_t r0, t0, snext;
319  if (inner) {
320  if (!HasInner()) return 0;
321  r0 = fRmin;
322  t0 = fTinsq;
323  } else {
324  r0 = fRmax;
325  t0 = fToutsq;
326  }
327  Double_t a = dir[0]*dir[0] + dir[1]*dir[1] - t0*dir[2]*dir[2];
328  Double_t b = t0*point[2]*dir[2] - point[0]*dir[0] - point[1]*dir[1];
329  Double_t c = point[0]*point[0] + point[1]*point[1] - t0*point[2]*point[2] - r0*r0;
330 
331  if (TMath::Abs(a) < TGeoShape::Tolerance()) {
332  if (TMath::Abs(b) < TGeoShape::Tolerance()) return 0;
333  snext = 0.5*c/b;
334  if (snext < 0.) return 0;
335  s[0] = snext;
336  return 1;
337  }
338 
339  Double_t delta = b*b - a*c;
340  Double_t ainv = 1./a;
341  Int_t npos = 0;
342  if (delta < 0.) return 0;
343  delta = TMath::Sqrt(delta);
344  Double_t sone = TMath::Sign(1.,ainv);
345  Int_t i = -1;
346  while (i<2) {
347  snext = (b + i*sone*delta)*ainv;
348  i += 2;
349  if (snext<0) continue;
350  if (snext<1.E-8) {
351  Double_t r = TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
352  Double_t t = (inner)?fTinsq:fToutsq;
353  t *= -point[2]/r;
354  Double_t phi = TMath::ATan2(point[1], point[0]);
355  Double_t ndotd = TMath::Cos(phi)*dir[0]+TMath::Sin(phi)*dir[1]+t*dir[2];
356  if (inner) ndotd *= -1;
357  if (in) ndotd *= -1;
358  if (ndotd<0) s[npos++] = snext;
359  } else s[npos++] = snext;
360  }
361  return npos;
362 }
363 
364 ////////////////////////////////////////////////////////////////////////////////
365 /// Cannot divide hyperboloids.
366 
367 TGeoVolume *TGeoHype::Divide(TGeoVolume * /*voldiv*/, const char *divname, Int_t /*iaxis*/, Int_t /*ndiv*/,
368  Double_t /*start*/, Double_t /*step*/)
369 {
370  Error("Divide", "Hyperboloids cannot be divided. Division volume %s not created", divname);
371  return 0;
372 }
373 
374 ////////////////////////////////////////////////////////////////////////////////
375 /// Get range of shape for a given axis.
376 
378 {
379  xlo = 0;
380  xhi = 0;
381  Double_t dx = 0;
382  switch (iaxis) {
383  case 1: // R
384  xlo = fRmin;
386  dx = xhi-xlo;
387  return dx;
388  case 2: // Phi
389  xlo = 0;
390  xhi = 360;
391  dx = 360;
392  return dx;
393  case 3: // Z
394  xlo = -fDz;
395  xhi = fDz;
396  dx = xhi-xlo;
397  return dx;
398  }
399  return dx;
400 }
401 
402 ////////////////////////////////////////////////////////////////////////////////
403 ///--- Fill vector param[4] with the bounding cylinder parameters. The order
404 /// is the following : Rmin, Rmax, Phi1, Phi2, dZ
405 
407 {
408  param[0] = fRmin; // Rmin
409  param[0] *= param[0];
410  param[1] = TMath::Sqrt(RadiusHypeSq(fDz, kFALSE)); // Rmax
411  param[1] *= param[1];
412  param[2] = 0.; // Phi1
413  param[3] = 360.; // Phi1
414 }
415 
416 ////////////////////////////////////////////////////////////////////////////////
417 /// in case shape has some negative parameters, these has to be computed
418 /// in order to fit the mother
419 
421 {
422  if (!TestShapeBit(kGeoRunTimeShape)) return 0;
423  Double_t dz;
424  Double_t zmin,zmax;
425  dz = fDz;
426  if (fDz<0) {
427  mother->GetAxisRange(3,zmin,zmax);
428  if (zmax<0) return 0;
429  dz=zmax;
430  } else {
431  Error("GetMakeRuntimeShape", "Shape %s does not have negative Z range", GetName());
432  return 0;
433  }
434  TGeoShape *hype = new TGeoHype(GetName(), dz, fRmax, fStOut, fRmin, fStIn);
435  return hype;
436 }
437 
438 ////////////////////////////////////////////////////////////////////////////////
439 /// print shape parameters
440 
442 {
443  printf("*** Shape %s: TGeoHype ***\n", GetName());
444  printf(" Rin = %11.5f\n", fRmin);
445  printf(" sin = %11.5f\n", fStIn);
446  printf(" Rout = %11.5f\n", fRmax);
447  printf(" sout = %11.5f\n", fStOut);
448  printf(" dz = %11.5f\n", fDz);
449 
450  printf(" Bounding box:\n");
452 }
453 
454 ////////////////////////////////////////////////////////////////////////////////
455 /// Creates a TBuffer3D describing *this* shape.
456 /// Coordinates are in local reference frame.
457 
459 {
461  Bool_t hasRmin = HasInner();
462  Int_t nbPnts = (hasRmin)?(2*n*n):(n*n+2);
463  Int_t nbSegs = (hasRmin)?(4*n*n):(n*(2*n+1));
464  Int_t nbPols = (hasRmin)?(2*n*n):(n*(n+1));
465 
467  nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols);
468  if (buff)
469  {
470  SetPoints(buff->fPnts);
471  SetSegsAndPols(*buff);
472  }
473 
474  return buff;
475 }
476 
477 ////////////////////////////////////////////////////////////////////////////////
478 /// Fill TBuffer3D structure for segments and polygons.
479 
481 {
482  Int_t c = GetBasicColor();
483  Int_t i, j, n;
484  n = gGeoManager->GetNsegments();
485  Bool_t hasRmin = HasInner();
486  Int_t irin = 0;
487  Int_t irout = (hasRmin)?(n*n):2;
488  // Fill segments
489  // Case hasRmin:
490  // Inner circles: [isin = 0], n (per circle) * n ( circles)
491  // iseg = isin+n*i+j , i = 0, n-1 , j = 0, n-1
492  // seg(i=1,n; j=1,n) = [irin+n*i+j] and [irin+n*i+(j+1)%n]
493  // Inner generators: [isgenin = isin+n*n], n (per circle) *(n-1) (slices)
494  // iseg = isgenin + i*n + j, i=0,n-2, j=0,n-1
495  // seg(i,j) = [irin+n*i+j] and [irin+n*(i+1)+j]
496  // Outer circles: [isout = isgenin+n*(n-1)], n (per circle) * n ( circles)
497  // iseg = isout + i*n + j , iz = 0, n-1 , j = 0, n-1
498  // seg(i=1,n; j=1,n) = [irout+n*i+j] and [irout+n*i+(j+1)%n]
499  // Outer generators: [isgenout = isout+n*n], n (per circle) *(n-1) (slices)
500  // iseg = isgenout + i*n + j, i=0,n-2, j=0,n-1
501  // seg(i,j) = [irout+n*i+j] and [irout+n*(i+1)+j]
502  // Lower cap : [islow = isgenout + n*(n-1)], n radial segments
503  // iseg = islow + j, j=0,n-1
504  // seg(j) = [irin + j] and [irout+j]
505  // Upper cap: [ishi = islow + n], nradial segments
506  // iseg = ishi + j, j=0,n-1
507  // seg[j] = [irin + n*(n-1) + j] and [irout+n*(n-1) + j]
508  //
509  // Case !hasRmin:
510  // Outer circles: [isout=0], same outer circles (n*n)
511  // Outer generators: isgenout = isout + n*n
512  // Lower cap: [islow = isgenout+n*(n-1)], n seg.
513  // iseg = islow + j, j=0,n-1
514  // seg[j] = [irin] and [irout+j]
515  // Upper cap: [ishi = islow +n]
516  // iseg = ishi + j, j=0,n-1
517  // seg[j] = [irin+1] and [irout+n*(n-1) + j]
518 
519  Int_t isin = 0;
520  Int_t isgenin = (hasRmin)?(isin+n*n):0;
521  Int_t isout = (hasRmin)?(isgenin+n*(n-1)):0;
522  Int_t isgenout = isout+n*n;
523  Int_t islo = isgenout+n*(n-1);
524  Int_t ishi = islo + n;
525 
526  Int_t npt = 0;
527  // Fill inner circle segments (n*n)
528  if (hasRmin) {
529  for (i=0; i<n; i++) {
530  for (j=0; j<n; j++) {
531  npt = 3*(isin+n*i+j);
532  buff.fSegs[npt] = c;
533  buff.fSegs[npt+1] = irin+n*i+j;
534  buff.fSegs[npt+2] = irin+n*i+((j+1)%n);
535  }
536  }
537  // Fill inner generators (n*(n-1))
538  for (i=0; i<n-1; i++) {
539  for (j=0; j<n; j++) {
540  npt = 3*(isgenin+n*i+j);
541  buff.fSegs[npt] = c;
542  buff.fSegs[npt+1] = irin+n*i+j;
543  buff.fSegs[npt+2] = irin+n*(i+1)+j;
544  }
545  }
546  }
547  // Fill outer circle segments (n*n)
548  for (i=0; i<n; i++) {
549  for (j=0; j<n; j++) {
550  npt = 3*(isout + n*i+j);
551  buff.fSegs[npt] = c;
552  buff.fSegs[npt+1] = irout+n*i+j;
553  buff.fSegs[npt+2] = irout+n*i+((j+1)%n);
554  }
555  }
556  // Fill outer generators (n*(n-1))
557  for (i=0; i<n-1; i++) {
558  for (j=0; j<n; j++) {
559  npt = 3*(isgenout+n*i+j);
560  buff.fSegs[npt] = c;
561  buff.fSegs[npt+1] = irout+n*i+j;
562  buff.fSegs[npt+2] = irout+n*(i+1)+j;
563  }
564  }
565  // Fill lower cap (n)
566  for (j=0; j<n; j++) {
567  npt = 3*(islo+j);
568  buff.fSegs[npt] = c;
569  buff.fSegs[npt+1] = irin;
570  if (hasRmin) buff.fSegs[npt+1] += j;
571  buff.fSegs[npt+2] = irout + j;
572  }
573  // Fill upper cap (n)
574  for (j=0; j<n; j++) {
575  npt = 3*(ishi+j);
576  buff.fSegs[npt] = c;
577  buff.fSegs[npt+1] = irin+1;
578  if (hasRmin) buff.fSegs[npt+1] += n*(n-1)+j-1;
579  buff.fSegs[npt+2] = irout + n*(n-1)+j;
580  }
581 
582  // Fill polygons
583  // Inner polygons: [ipin = 0] (n-1) slices * n (edges)
584  // ipoly = ipin + n*i + j; i=0,n-2 j=0,n-1
585  // poly[i,j] = [isin+n*i+j] [isgenin+i*n+(j+1)%n] [isin+n*(i+1)+j] [isgenin+i*n+j]
586  // Outer polygons: [ipout = ipin+n*(n-1)] also (n-1)*n
587  // ipoly = ipout + n*i + j; i=0,n-2 j=0,n-1
588  // poly[i,j] = [isout+n*i+j] [isgenout+i*n+j] [isout+n*(i+1)+j] [isgenout+i*n+(j+1)%n]
589  // Lower cap: [iplow = ipout+n*(n-1): n polygons
590  // ipoly = iplow + j; j=0,n-1
591  // poly[i=0,j] = [isin+j] [islow+j] [isout+j] [islow+(j+1)%n]
592  // Upper cap: [ipup = iplow+n] : n polygons
593  // ipoly = ipup + j; j=0,n-1
594  // poly[i=n-1, j] = [isin+n*(n-1)+j] [ishi+(j+1)%n] [isout+n*(n-1)+j] [ishi+j]
595  //
596  // Case !hasRmin:
597  // ipin = 0 no inner polygons
598  // ipout = 0 same outer polygons
599  // Lower cap: iplow = ipout+n*(n-1): n polygons with 3 segments
600  // poly[i=0,j] = [isout+j] [islow+(j+1)%n] [islow+j]
601  // Upper cap: ipup = iplow+n;
602  // poly[i=n-1,j] = [isout+n*(n-1)+j] [ishi+j] [ishi+(j+1)%n]
603 
604  Int_t ipin = 0;
605  Int_t ipout = (hasRmin)?(ipin+n*(n-1)):0;
606  Int_t iplo = ipout+n*(n-1);
607  Int_t ipup = iplo+n;
608  // Inner polygons n*(n-1)
609  if (hasRmin) {
610  for (i=0; i<n-1; i++) {
611  for (j=0; j<n; j++) {
612  npt = 6*(ipin+n*i+j);
613  buff.fPols[npt] = c;
614  buff.fPols[npt+1] = 4;
615  buff.fPols[npt+2] = isin+n*i+j;
616  buff.fPols[npt+3] = isgenin+i*n+((j+1)%n);
617  buff.fPols[npt+4] = isin+n*(i+1)+j;
618  buff.fPols[npt+5] = isgenin+i*n+j;
619  }
620  }
621  }
622  // Outer polygons n*(n-1)
623  for (i=0; i<n-1; i++) {
624  for (j=0; j<n; j++) {
625  npt = 6*(ipout+n*i+j);
626  buff.fPols[npt] = c;
627  buff.fPols[npt+1] = 4;
628  buff.fPols[npt+2] = isout+n*i+j;
629  buff.fPols[npt+3] = isgenout+i*n+j;
630  buff.fPols[npt+4] = isout+n*(i+1)+j;
631  buff.fPols[npt+5] = isgenout+i*n+((j+1)%n);
632  }
633  }
634  // End caps
635  if (hasRmin) {
636  for (j=0; j<n; j++) {
637  npt = 6*(iplo+j);
638  buff.fPols[npt] = c+1;
639  buff.fPols[npt+1] = 4;
640  buff.fPols[npt+2] = isin+j;
641  buff.fPols[npt+3] = islo+j;
642  buff.fPols[npt+4] = isout+j;
643  buff.fPols[npt+5] = islo+((j+1)%n);
644  }
645  for (j=0; j<n; j++) {
646  npt = 6*(ipup+j);
647  buff.fPols[npt] = c+2;
648  buff.fPols[npt+1] = 4;
649  buff.fPols[npt+2] = isin+n*(n-1)+j;
650  buff.fPols[npt+3] = ishi+((j+1)%n);
651  buff.fPols[npt+4] = isout+n*(n-1)+j;
652  buff.fPols[npt+5] = ishi+j;
653  }
654  } else {
655  for (j=0; j<n; j++) {
656  npt = 6*iplo+5*j;
657  buff.fPols[npt] = c+1;
658  buff.fPols[npt+1] = 3;
659  buff.fPols[npt+2] = isout+j;
660  buff.fPols[npt+3] = islo+((j+1)%n);
661  buff.fPols[npt+4] = islo+j;
662  }
663  for (j=0; j<n; j++) {
664  npt = 6*iplo+5*(n+j);
665  buff.fPols[npt] = c+2;
666  buff.fPols[npt+1] = 3;
667  buff.fPols[npt+2] = isout+n*(n-1)+j;
668  buff.fPols[npt+3] = ishi+j;
669  buff.fPols[npt+4] = ishi+((j+1)%n);
670  }
671  }
672 }
673 
674 ////////////////////////////////////////////////////////////////////////////////
675 /// Compute r^2 = x^2 + y^2 at a given z coordinate, for either inner or outer hyperbolas.
676 
678 {
679  Double_t r0, tsq;
680  if (inner) {
681  r0 = fRmin;
682  tsq = fTinsq;
683  } else {
684  r0 = fRmax;
685  tsq = fToutsq;
686  }
687  return (r0*r0+tsq*z*z);
688 }
689 
690 ////////////////////////////////////////////////////////////////////////////////
691 /// Compute z^2 at a given r^2, for either inner or outer hyperbolas.
692 
694 {
695  Double_t r0, tsq;
696  if (inner) {
697  r0 = fRmin;
698  tsq = fTinsq;
699  } else {
700  r0 = fRmax;
701  tsq = fToutsq;
702  }
703  if (TMath::Abs(tsq) < TGeoShape::Tolerance()) return TGeoShape::Big();
704  return ((r*r-r0*r0)/tsq);
705 }
706 
707 ////////////////////////////////////////////////////////////////////////////////
708 /// computes the closest distance from given point to this shape, according
709 /// to option. The matching point on the shape is stored in spoint.
710 
711 Double_t TGeoHype::Safety(const Double_t *point, Bool_t in) const
712 {
713  Double_t safe, safrmin, safrmax;
714  if (in) {
715  safe = fDz-TMath::Abs(point[2]);
716  safrmin = SafetyToHype(point, kTRUE, in);
717  if (safrmin < safe) safe = safrmin;
718  safrmax = SafetyToHype(point, kFALSE,in);
719  if (safrmax < safe) safe = safrmax;
720  } else {
721  safe = -fDz+TMath::Abs(point[2]);
722  safrmin = SafetyToHype(point, kTRUE, in);
723  if (safrmin > safe) safe = safrmin;
724  safrmax = SafetyToHype(point, kFALSE,in);
725  if (safrmax > safe) safe = safrmax;
726  }
727  return safe;
728 }
729 
730 ////////////////////////////////////////////////////////////////////////////////
731 /// Compute an underestimate of the closest distance from a point to inner or
732 /// outer infinite hyperbolas.
733 
734 Double_t TGeoHype::SafetyToHype(const Double_t *point, Bool_t inner, Bool_t in) const
735 {
736  Double_t r, rsq, rhsq, rh, dr, tsq, saf;
737  if (inner && !HasInner()) return (in)?TGeoShape::Big():-TGeoShape::Big();
738  rsq = point[0]*point[0]+point[1]*point[1];
739  r = TMath::Sqrt(rsq);
740  rhsq = RadiusHypeSq(point[2], inner);
741  rh = TMath::Sqrt(rhsq);
742  dr = r - rh;
743  if (inner) {
744  if (!in && dr>0) return -TGeoShape::Big();
745  if (TMath::Abs(fStIn) < TGeoShape::Tolerance()) return TMath::Abs(dr);
746  if (fRmin<TGeoShape::Tolerance()) return TMath::Abs(dr/TMath::Sqrt(1.+ fTinsq));
747  tsq = fTinsq;
748  } else {
749  if (!in && dr<0) return -TGeoShape::Big();
750  if (TMath::Abs(fStOut) < TGeoShape::Tolerance()) return TMath::Abs(dr);
751  tsq = fToutsq;
752  }
753  if (TMath::Abs(dr)<TGeoShape::Tolerance()) return 0.;
754  // 1. dr<0 => approximate safety with distance to tangent to hyperbola in z = |point[2]|
755  Double_t m;
756  if (dr<0) {
757  m = rh/(tsq*TMath::Abs(point[2]));
758  saf = -m*dr/TMath::Sqrt(1.+m*m);
759  return saf;
760  }
761  // 2. dr>0 => approximate safety with distance from point to segment P1(r(z0),z0) and P2(r0, z(r0))
762  m = (TMath::Sqrt(ZHypeSq(r,inner)) - TMath::Abs(point[2]))/dr;
763  saf = m*dr/TMath::Sqrt(1.+m*m);
764  return saf;
765 }
766 
767 ////////////////////////////////////////////////////////////////////////////////
768 /// Save a primitive as a C++ statement(s) on output stream "out".
769 
770 void TGeoHype::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
771 {
772  if (TObject::TestBit(kGeoSavePrimitive)) return;
773  out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
774  out << " rin = " << fRmin << ";" << std::endl;
775  out << " stin = " << fStIn << ";" << std::endl;
776  out << " rout = " << fRmax << ";" << std::endl;
777  out << " stout = " << fStOut << ";" << std::endl;
778  out << " dz = " << fDz << ";" << std::endl;
779  out << " TGeoShape *" << GetPointerName() << " = new TGeoHype(\"" << GetName() << "\",rin,stin,rout,stout,dz);" << std::endl;
781 }
782 
783 ////////////////////////////////////////////////////////////////////////////////
784 /// Set dimensions of the hyperboloid.
785 
787 {
788  fRmin = rin;
789  fRmax = rout;
790  fDz = dz;
791  fStIn = stin;
792  fStOut = stout;
794  fTinsq = fTin*fTin;
796  fToutsq = fTout*fTout;
797  if ((fRmin==0) && (fStIn==0)) SetShapeBit(kGeoRSeg, kTRUE);
798  else SetShapeBit(kGeoRSeg, kFALSE);
799 }
800 
801 ////////////////////////////////////////////////////////////////////////////////
802 /// Set dimensions of the hyperboloid starting from an array.
803 /// param[0] = dz
804 /// param[1] = rin
805 /// param[2] = stin
806 /// param[3] = rout
807 /// param[4] = stout
808 
810 {
811  Double_t dz = param[0];
812  Double_t rin = param[1];
813  Double_t stin = param[2];
814  Double_t rout = param[3];
815  Double_t stout = param[4];
816  SetHypeDimensions(rin, stin, rout, stout, dz);
817 }
818 
819 ////////////////////////////////////////////////////////////////////////////////
820 /// create tube mesh points
821 
823 {
824  Double_t z,dz,r;
825  Int_t i,j, n;
826  if (!points) return;
827  n = gGeoManager->GetNsegments();
828  Double_t dphi = 360./n;
829  Double_t phi = 0;
830  dz = 2.*fDz/(n-1);
831 
832  Int_t indx = 0;
833 
834  if (HasInner()) {
835  // Inner surface points
836  for (i=0; i<n; i++) {
837  z = -fDz+i*dz;
838  r = TMath::Sqrt(RadiusHypeSq(z, kTRUE));
839  for (j=0; j<n; j++) {
840  phi = j*dphi*TMath::DegToRad();
841  points[indx++] = r * TMath::Cos(phi);
842  points[indx++] = r * TMath::Sin(phi);
843  points[indx++] = z;
844  }
845  }
846  } else {
847  points[indx++] = 0.;
848  points[indx++] = 0.;
849  points[indx++] = -fDz;
850  points[indx++] = 0.;
851  points[indx++] = 0.;
852  points[indx++] = fDz;
853  }
854  // Outer surface points
855  for (i=0; i<n; i++) {
856  z = -fDz + i*dz;
858  for (j=0; j<n; j++) {
859  phi = j*dphi*TMath::DegToRad();
860  points[indx++] = r * TMath::Cos(phi);
861  points[indx++] = r * TMath::Sin(phi);
862  points[indx++] = z;
863  }
864  }
865 }
866 
867 ////////////////////////////////////////////////////////////////////////////////
868 /// create tube mesh points
869 
871 {
872  Double_t z,dz,r;
873  Int_t i,j, n;
874  if (!points) return;
875  n = gGeoManager->GetNsegments();
876  Double_t dphi = 360./n;
877  Double_t phi = 0;
878  dz = 2.*fDz/(n-1);
879 
880  Int_t indx = 0;
881 
882  if (HasInner()) {
883  // Inner surface points
884  for (i=0; i<n; i++) {
885  z = -fDz+i*dz;
886  r = TMath::Sqrt(RadiusHypeSq(z, kTRUE));
887  for (j=0; j<n; j++) {
888  phi = j*dphi*TMath::DegToRad();
889  points[indx++] = r * TMath::Cos(phi);
890  points[indx++] = r * TMath::Sin(phi);
891  points[indx++] = z;
892  }
893  }
894  } else {
895  points[indx++] = 0.;
896  points[indx++] = 0.;
897  points[indx++] = -fDz;
898  points[indx++] = 0.;
899  points[indx++] = 0.;
900  points[indx++] = fDz;
901  }
902  // Outer surface points
903  for (i=0; i<n; i++) {
904  z = -fDz + i*dz;
906  for (j=0; j<n; j++) {
907  phi = j*dphi*TMath::DegToRad();
908  points[indx++] = r * TMath::Cos(phi);
909  points[indx++] = r * TMath::Sin(phi);
910  points[indx++] = z;
911  }
912  }
913 }
914 
915 ////////////////////////////////////////////////////////////////////////////////
916 /// Returns numbers of vertices, segments and polygons composing the shape mesh.
917 
918 void TGeoHype::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
919 {
921  Bool_t hasRmin = HasInner();
922  nvert = (hasRmin)?(2*n*n):(n*n+2);
923  nsegs = (hasRmin)?(4*n*n):(n*(2*n+1));
924  npols = (hasRmin)?(2*n*n):(n*(n+1));
925 }
926 
927 ////////////////////////////////////////////////////////////////////////////////
928 /// Return number of vertices of the mesh representation
929 
931 {
933  Int_t numPoints = (HasRmin())?(2*n*n):(n*n+2);
934  return numPoints;
935 }
936 
937 ////////////////////////////////////////////////////////////////////////////////
938 ////// fill size of this 3-D object
939 //// TVirtualGeoPainter *painter = gGeoManager->GetGeomPainter();
940 //// if (!painter) return;
941 //// Int_t n = gGeoManager->GetNsegments();
942 //// Int_t numPoints = n*4;
943 //// Int_t numSegs = n*8;
944 //// Int_t numPolys = n*4;
945 //// painter->AddSize3D(numPoints, numSegs, numPolys);
946 
947 void TGeoHype::Sizeof3D() const
948 {
949 }
950 
951 ////////////////////////////////////////////////////////////////////////////////
952 /// Fills a static 3D buffer and returns a reference.
953 
954 const TBuffer3D & TGeoHype::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
955 {
957 
958  TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
959 
960  if (reqSections & TBuffer3D::kRawSizes) {
962  Bool_t hasRmin = HasInner();
963  Int_t nbPnts = (hasRmin)?(2*n*n):(n*n+2);
964  Int_t nbSegs = (hasRmin)?(4*n*n):(n*(2*n+1));
965  Int_t nbPols = (hasRmin)?(2*n*n):(n*(n+1));
966  if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
967  buffer.SetSectionsValid(TBuffer3D::kRawSizes);
968  }
969  }
970  if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
971  SetPoints(buffer.fPnts);
972  if (!buffer.fLocalFrame) {
973  TransformPoints(buffer.fPnts, buffer.NbPnts());
974  }
975 
976  SetSegsAndPols(buffer);
977  buffer.SetSectionsValid(TBuffer3D::kRaw);
978  }
979 
980  return buffer;
981 }
982 
983 ////////////////////////////////////////////////////////////////////////////////
984 /// Check the inside status for each of the points in the array.
985 /// Input: Array of point coordinates + vector size
986 /// Output: Array of Booleans for the inside of each point
987 
988 void TGeoHype::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
989 {
990  for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
991 }
992 
993 ////////////////////////////////////////////////////////////////////////////////
994 /// Compute the normal for an array o points so that norm.dot.dir is positive
995 /// Input: Arrays of point coordinates and directions + vector size
996 /// Output: Array of normal directions
997 
998 void TGeoHype::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
999 {
1000  for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
1001 }
1002 
1003 ////////////////////////////////////////////////////////////////////////////////
1004 /// Compute distance from array of input points having directions specisied by dirs. Store output in dists
1005 
1006 void TGeoHype::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1007 {
1008  for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
1009 }
1010 
1011 ////////////////////////////////////////////////////////////////////////////////
1012 /// Compute distance from array of input points having directions specisied by dirs. Store output in dists
1013 
1014 void TGeoHype::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1015 {
1016  for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
1017 }
1018 
1019 ////////////////////////////////////////////////////////////////////////////////
1020 /// Compute safe distance from each of the points in the input array.
1021 /// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1022 /// Output: Safety values
1023 
1024 void TGeoHype::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
1025 {
1026  for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
1027 }
Double_t fStOut
Definition: TGeoHype.h:55
virtual ~TGeoHype()
destructor
Definition: TGeoHype.cxx:118
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition: TGeoHype.cxx:125
tuple buffer
Definition: tree.py:99
#define snext(osub1, osub2)
Definition: triangle.c:1167
Int_t GetNsegments() const
Get number of segments approximating circles.
T1 Sign(T1 a, T2 b)
Definition: TMathBase.h:155
virtual TBuffer3D * MakeBuffer3D() const
Creates a TBuffer3D describing this shape.
Definition: TGeoHype.cxx:458
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Definition: TGeoHype.cxx:377
float Float_t
Definition: RtypesCore.h:53
return c
ClassImp(TGeoHype) TGeoHype
Default constructor.
Definition: TGeoHype.cxx:55
const char Option_t
Definition: RtypesCore.h:62
Int_t GetBasicColor() const
Get the basic color (0-7).
Definition: TGeoShape.cxx:672
virtual void SetPoints(Double_t *points) const
create tube mesh points
Definition: TGeoHype.cxx:822
Double_t DegToRad()
Definition: TMath.h:50
void hype()
Definition: geodemo.C:1056
Double_t fOrigin[3]
Definition: TGeoBBox.h:36
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition: TGeoHype.cxx:770
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
Double_t SafetyToHype(const Double_t *point, Bool_t inner, Bool_t in) const
Compute an underestimate of the closest distance from a point to inner or outer infinite hyperbolas...
Definition: TGeoHype.cxx:734
Float_t py
Definition: hprod.C:33
Double_t fTout
Definition: TGeoHype.h:60
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:386
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
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: TGeoHype.cxx:711
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:732
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: TGeoHype.cxx:420
Double_t fTin
Definition: TGeoHype.h:59
static Double_t Tolerance()
Definition: TGeoShape.h:101
virtual Bool_t Contains(const Double_t *point) const
test if point is inside this tube
Definition: TGeoHype.cxx:192
virtual void ComputeBBox()
Compute bounding box of the hyperboloid.
Definition: TGeoHype.cxx:135
void SetHypeDimensions(Double_t rin, Double_t stin, Double_t rout, Double_t stout, Double_t dz)
Set dimensions of the hyperboloid.
Definition: TGeoHype.cxx:786
Double_t fDZ
Definition: TGeoBBox.h:35
Double_t * fPnts
Definition: TBuffer3D.h:114
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
compute closest distance from point px,py to each corner
Definition: TGeoHype.cxx:207
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:454
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
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 hyperboloid.
Definition: TGeoHype.cxx:256
double sin(double)
Float_t z[5]
Definition: Ifit.C:16
Double_t fDz
Definition: TGeoTube.h:34
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Cannot divide hyperboloids.
Definition: TGeoHype.cxx:367
char * out
Definition: TBase64.cxx:29
void SetSectionsValid(UInt_t mask)
Definition: TBuffer3D.h:67
Int_t * fPols
Definition: TBuffer3D.h:116
Bool_t fLocalFrame
Definition: TBuffer3D.h:92
point * points
Definition: X3DBuffer.c:20
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
Definition: TGeoShape.cxx:551
Double_t ZHypeSq(Double_t r, Bool_t inner) const
Compute z^2 at a given r^2, for either inner or outer hyperbolas.
Definition: TGeoHype.cxx:693
TThread * t[5]
Definition: threadsh1.C:13
ROOT::R::TRInterface & r
Definition: Object.C:4
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 specisied by dirs. Store output in dist...
Definition: TGeoHype.cxx:1006
virtual void SetDimensions(Double_t *param)
Set dimensions of the hyperboloid starting from an array.
Definition: TGeoHype.cxx:809
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:187
virtual void InspectShape() const
print shape parameters
Definition: TGeoHype.cxx:441
Double_t fTinsq
Definition: TGeoHype.h:61
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:173
TMarker * m
Definition: textangle.C:8
#define isin(address, start, length)
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: TGeoHype.cxx:918
Double_t E()
Definition: TMath.h:54
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:357
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
Definition: TGeoHype.cxx:954
Generic 3D primitive description class.
Definition: TBuffer3D.h:19
Double_t fStIn
Definition: TGeoHype.h:54
virtual void Sizeof3D() const
Definition: TGeoHype.cxx:947
virtual const char * GetName() const
Get the shape name.
Definition: TGeoShape.cxx:247
virtual void SetSegsAndPols(TBuffer3D &buff) const
Fill TBuffer3D structure for segments and polygons.
Definition: TGeoHype.cxx:480
Double_t Cos(Double_t)
Definition: TMath.h:424
Double_t Pi()
Definition: TMath.h:44
virtual Int_t GetNmeshVertices() const
Return number of vertices of the mesh representation.
Definition: TGeoHype.cxx:930
Float_t phi
Definition: shapesAnim.C:6
Bool_t HasInner() const
Definition: TGeoHype.h:100
virtual void InspectShape() const
Prints shape parameters.
Definition: TGeoBBox.cxx:749
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:556
double Double_t
Definition: RtypesCore.h:55
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:989
Int_t DistToHype(const Double_t *point, const Double_t *dir, Double_t *s, Bool_t inner, Bool_t in) const
Compute distance from an arbitrary point to inner/outer surface of hyperboloid.
Definition: TGeoHype.cxx:316
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
Definition: TGeoHype.cxx:155
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
void dir(char *path=0)
Definition: rootalias.C:30
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
Definition: TGeoShape.cxx:258
Bool_t HasRmin() const
Definition: TGeoTube.h:80
const char * GetPointerName() const
Provide a pointer name containing uid.
Definition: TGeoShape.cxx:698
static Double_t Big()
Definition: TGeoShape.h:98
Double_t fToutsq
Definition: TGeoHype.h:62
Double_t fDY
Definition: TGeoBBox.h:34
#define name(a, b)
Definition: linkTestLib0.cpp:5
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
Definition: TGeoShape.cxx:523
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 hyperboloid.
Definition: TGeoHype.cxx:216
Float_t px
Definition: hprod.C:33
Int_t * fSegs
Definition: TBuffer3D.h:115
Double_t fRmin
Definition: TGeoTube.h:32
Bool_t TestShapeBit(UInt_t f) const
Definition: TGeoShape.h:172
UInt_t NbPnts() const
Definition: TBuffer3D.h:82
Double_t Sin(Double_t)
Definition: TMath.h:421
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const =0
Double_t fDX
Definition: TGeoBBox.h:33
Double_t dr
Definition: parallelcoord.C:14
virtual void GetBoundingCylinder(Double_t *param) const
— Fill vector param[4] with the bounding cylinder parameters.
Definition: TGeoHype.cxx:406
tuple ct
Definition: tornado.py:53
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: TGeoHype.cxx:998
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
Long64_t LocMin(Long64_t n, const T *a)
Definition: TMath.h:695
const Bool_t kTRUE
Definition: Rtypes.h:91
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 specisied by dirs. Store output in dist...
Definition: TGeoHype.cxx:1014
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
const Int_t n
Definition: legend1.C:16
Double_t Tan(Double_t)
Definition: TMath.h:427
Double_t RadiusHypeSq(Double_t z, Bool_t inner) const
Compute r^2 = x^2 + y^2 at a given z coordinate, for either inner or outer hyperbolas.
Definition: TGeoHype.cxx:677
unsigned int r2[N_CITIES]
Definition: simanTSP.cxx:322
Double_t fRmax
Definition: TGeoTube.h:33
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: TGeoHype.cxx:1024
Bool_t SectionsValid(UInt_t mask) const
Definition: TBuffer3D.h:69
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: TGeoHype.cxx:988
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904