Geant4  10.00.p02
G4Tubs.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4Tubs.cc 81636 2014-06-04 09:06:08Z gcosmo $
28 //
29 //
30 // class G4Tubs
31 //
32 // History:
33 //
34 // 05.04.12 M.Kelsey: Use sqrt(r) in GetPointOnSurface() for uniform points
35 // 02.08.07 T.Nikitina: bug fixed in DistanceToOut(p,v,..) for negative value under sqrt
36 // for the case: p on the surface and v is tangent to the surface
37 // 11.05.07 T.Nikitina: bug fixed in DistanceToOut(p,v,..) for phi < 2pi
38 // 03.05.05 V.Grichine: SurfaceNormal(p) according to J. Apostolakis proposal
39 // 16.03.05 V.Grichine: SurfaceNormal(p) with edges/corners for boolean
40 // 20.07.01 V.Grichine: bug fixed in Inside(p)
41 // 20.02.01 V.Grichine: bug fixed in Inside(p) and CalculateExtent was
42 // simplified base on G4Box::CalculateExtent
43 // 07.12.00 V.Grichine: phi-section algorithm was changed in Inside(p)
44 // 28.11.00 V.Grichine: bug fixed in Inside(p)
45 // 31.10.00 V.Grichine: assign srd, sphi in Distance ToOut(p,v,...)
46 // 08.08.00 V.Grichine: more stable roots of 2-equation in DistanceToOut(p,v,..)
47 // 02.08.00 V.Grichine: point is outside check in Distance ToOut(p)
48 // 17.05.00 V.Grichine: bugs (#76,#91) fixed in Distance ToOut(p,v,...)
49 // 31.03.00 V.Grichine: bug fixed in Inside(p)
50 // 19.11.99 V.Grichine: side = kNull in DistanceToOut(p,v,...)
51 // 13.10.99 V.Grichine: bugs fixed in DistanceToIn(p,v)
52 // 28.05.99 V.Grichine: bugs fixed in DistanceToOut(p,v,...)
53 // 25.05.99 V.Grichine: bugs fixed in DistanceToIn(p,v)
54 // 23.03.99 V.Grichine: bug fixed in DistanceToIn(p,v)
55 // 09.10.98 V.Grichine: modifications in DistanceToOut(p,v,...)
56 // 18.06.98 V.Grichine: n-normalisation in DistanceToOut(p,v)
57 //
58 // 1994-95 P.Kent: implementation
59 //
61 
62 #include "G4Tubs.hh"
63 
64 #if !defined(G4GEOM_USE_UTUBS)
65 
66 #include "G4VoxelLimits.hh"
67 #include "G4AffineTransform.hh"
68 #include "G4GeometryTolerance.hh"
69 
70 #include "G4VPVParameterisation.hh"
71 
72 #include "Randomize.hh"
73 
74 #include "meshdefs.hh"
75 
76 #include "G4VGraphicsScene.hh"
77 
78 using namespace CLHEP;
79 
81 //
82 // Constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
83 // - note if pdphi>2PI then reset to 2PI
84 
85 G4Tubs::G4Tubs( const G4String &pName,
86  G4double pRMin, G4double pRMax,
87  G4double pDz,
88  G4double pSPhi, G4double pDPhi )
89  : G4CSGSolid(pName), fRMin(pRMin), fRMax(pRMax), fDz(pDz), fSPhi(0), fDPhi(0)
90 {
91 
94 
98 
99  if (pDz<=0) // Check z-len
100  {
101  std::ostringstream message;
102  message << "Negative Z half-length (" << pDz << ") in solid: " << GetName();
103  G4Exception("G4Tubs::G4Tubs()", "GeomSolids0002", FatalException, message);
104  }
105  if ( (pRMin >= pRMax) || (pRMin < 0) ) // Check radii
106  {
107  std::ostringstream message;
108  message << "Invalid values for radii in solid: " << GetName()
109  << G4endl
110  << " pRMin = " << pRMin << ", pRMax = " << pRMax;
111  G4Exception("G4Tubs::G4Tubs()", "GeomSolids0002", FatalException, message);
112  }
113 
114  // Check angles
115  //
116  CheckPhiAngles(pSPhi, pDPhi);
117 }
118 
120 //
121 // Fake default constructor - sets only member data and allocates memory
122 // for usage restricted to object persistency.
123 //
124 G4Tubs::G4Tubs( __void__& a )
125  : G4CSGSolid(a), kRadTolerance(0.), kAngTolerance(0.),
126  fRMin(0.), fRMax(0.), fDz(0.), fSPhi(0.), fDPhi(0.),
127  sinCPhi(0.), cosCPhi(0.), cosHDPhiOT(0.), cosHDPhiIT(0.),
128  sinSPhi(0.), cosSPhi(0.), sinEPhi(0.), cosEPhi(0.),
129  fPhiFullTube(false), halfCarTolerance(0.), halfRadTolerance(0.),
130  halfAngTolerance(0.)
131 
132 {
133 }
134 
136 //
137 // Destructor
138 
140 {
141 }
142 
144 //
145 // Copy constructor
146 
148  : G4CSGSolid(rhs),
149  kRadTolerance(rhs.kRadTolerance), kAngTolerance(rhs.kAngTolerance),
150  fRMin(rhs.fRMin), fRMax(rhs.fRMax), fDz(rhs.fDz),
151  fSPhi(rhs.fSPhi), fDPhi(rhs.fDPhi),
152  sinCPhi(rhs.sinCPhi), cosCPhi(rhs.sinCPhi),
153  cosHDPhiOT(rhs.cosHDPhiOT), cosHDPhiIT(rhs.cosHDPhiOT),
154  sinSPhi(rhs.sinSPhi), cosSPhi(rhs.cosSPhi),
155  sinEPhi(rhs.sinEPhi), cosEPhi(rhs.cosEPhi), fPhiFullTube(rhs.fPhiFullTube),
156  halfCarTolerance(rhs.halfCarTolerance),
157  halfRadTolerance(rhs.halfRadTolerance),
158  halfAngTolerance(rhs.halfAngTolerance)
159 {
161 }
162 
164 //
165 // Assignment operator
166 
168 {
169  // Check assignment to self
170  //
171  if (this == &rhs) { return *this; }
172 
173  // Copy base class data
174  //
176 
177  // Copy data
178  //
180  fRMin = rhs.fRMin; fRMax = rhs.fRMax; fDz = rhs.fDz;
181  fSPhi = rhs.fSPhi; fDPhi = rhs.fDPhi;
182  sinCPhi = rhs.sinCPhi; cosCPhi = rhs.sinCPhi;
184  sinSPhi = rhs.sinSPhi; cosSPhi = rhs.cosSPhi;
185  sinEPhi = rhs.sinEPhi; cosEPhi = rhs.cosEPhi;
191 
192  return *this;
193 }
194 
196 //
197 // Dispatch to parameterisation for replication mechanism dimension
198 // computation & modification.
199 
201  const G4int n,
202  const G4VPhysicalVolume* pRep )
203 {
204  p->ComputeDimensions(*this,n,pRep) ;
205 }
206 
208 //
209 // Calculate extent under transform and specified limit
210 
212  const G4VoxelLimits& pVoxelLimit,
213  const G4AffineTransform& pTransform,
214  G4double& pMin,
215  G4double& pMax ) const
216 {
217 
218  if ( (!pTransform.IsRotated()) && (fDPhi == twopi) && (fRMin == 0) )
219  {
220  // Special case handling for unrotated solid tubes
221  // Compute x/y/z mins and maxs fro bounding box respecting limits,
222  // with early returns if outside limits. Then switch() on pAxis,
223  // and compute exact x and y limit for x/y case
224 
225  G4double xoffset, xMin, xMax;
226  G4double yoffset, yMin, yMax;
227  G4double zoffset, zMin, zMax;
228 
229  G4double diff1, diff2, maxDiff, newMin, newMax;
230  G4double xoff1, xoff2, yoff1, yoff2, delta;
231 
232  xoffset = pTransform.NetTranslation().x();
233  xMin = xoffset - fRMax;
234  xMax = xoffset + fRMax;
235 
236  if (pVoxelLimit.IsXLimited())
237  {
238  if ( (xMin > pVoxelLimit.GetMaxXExtent())
239  || (xMax < pVoxelLimit.GetMinXExtent()) )
240  {
241  return false;
242  }
243  else
244  {
245  if (xMin < pVoxelLimit.GetMinXExtent())
246  {
247  xMin = pVoxelLimit.GetMinXExtent();
248  }
249  if (xMax > pVoxelLimit.GetMaxXExtent())
250  {
251  xMax = pVoxelLimit.GetMaxXExtent();
252  }
253  }
254  }
255  yoffset = pTransform.NetTranslation().y();
256  yMin = yoffset - fRMax;
257  yMax = yoffset + fRMax;
258 
259  if ( pVoxelLimit.IsYLimited() )
260  {
261  if ( (yMin > pVoxelLimit.GetMaxYExtent())
262  || (yMax < pVoxelLimit.GetMinYExtent()) )
263  {
264  return false;
265  }
266  else
267  {
268  if (yMin < pVoxelLimit.GetMinYExtent())
269  {
270  yMin = pVoxelLimit.GetMinYExtent();
271  }
272  if (yMax > pVoxelLimit.GetMaxYExtent())
273  {
274  yMax=pVoxelLimit.GetMaxYExtent();
275  }
276  }
277  }
278  zoffset = pTransform.NetTranslation().z();
279  zMin = zoffset - fDz;
280  zMax = zoffset + fDz;
281 
282  if ( pVoxelLimit.IsZLimited() )
283  {
284  if ( (zMin > pVoxelLimit.GetMaxZExtent())
285  || (zMax < pVoxelLimit.GetMinZExtent()) )
286  {
287  return false;
288  }
289  else
290  {
291  if (zMin < pVoxelLimit.GetMinZExtent())
292  {
293  zMin = pVoxelLimit.GetMinZExtent();
294  }
295  if (zMax > pVoxelLimit.GetMaxZExtent())
296  {
297  zMax = pVoxelLimit.GetMaxZExtent();
298  }
299  }
300  }
301  switch ( pAxis ) // Known to cut cylinder
302  {
303  case kXAxis :
304  {
305  yoff1 = yoffset - yMin;
306  yoff2 = yMax - yoffset;
307 
308  if ( (yoff1 >= 0) && (yoff2 >= 0) ) // Y limits cross max/min x
309  { // => no change
310  pMin = xMin;
311  pMax = xMax;
312  }
313  else
314  {
315  // Y limits don't cross max/min x => compute max delta x,
316  // hence new mins/maxs
317 
318  delta = fRMax*fRMax - yoff1*yoff1;
319  diff1 = (delta>0.) ? std::sqrt(delta) : 0.;
320  delta = fRMax*fRMax - yoff2*yoff2;
321  diff2 = (delta>0.) ? std::sqrt(delta) : 0.;
322  maxDiff = (diff1 > diff2) ? diff1:diff2;
323  newMin = xoffset - maxDiff;
324  newMax = xoffset + maxDiff;
325  pMin = (newMin < xMin) ? xMin : newMin;
326  pMax = (newMax > xMax) ? xMax : newMax;
327  }
328  break;
329  }
330  case kYAxis :
331  {
332  xoff1 = xoffset - xMin;
333  xoff2 = xMax - xoffset;
334 
335  if ( (xoff1 >= 0) && (xoff2 >= 0) ) // X limits cross max/min y
336  { // => no change
337  pMin = yMin;
338  pMax = yMax;
339  }
340  else
341  {
342  // X limits don't cross max/min y => compute max delta y,
343  // hence new mins/maxs
344 
345  delta = fRMax*fRMax - xoff1*xoff1;
346  diff1 = (delta>0.) ? std::sqrt(delta) : 0.;
347  delta = fRMax*fRMax - xoff2*xoff2;
348  diff2 = (delta>0.) ? std::sqrt(delta) : 0.;
349  maxDiff = (diff1 > diff2) ? diff1 : diff2;
350  newMin = yoffset - maxDiff;
351  newMax = yoffset + maxDiff;
352  pMin = (newMin < yMin) ? yMin : newMin;
353  pMax = (newMax > yMax) ? yMax : newMax;
354  }
355  break;
356  }
357  case kZAxis:
358  {
359  pMin = zMin;
360  pMax = zMax;
361  break;
362  }
363  default:
364  break;
365  }
366  pMin -= kCarTolerance;
367  pMax += kCarTolerance;
368  return true;
369  }
370  else // Calculate rotated vertex coordinates
371  {
372  G4int i, noEntries, noBetweenSections4;
373  G4bool existsAfterClip = false;
374  G4ThreeVectorList* vertices = CreateRotatedVertices(pTransform);
375 
376  pMin = kInfinity;
377  pMax = -kInfinity;
378 
379  noEntries = vertices->size();
380  noBetweenSections4 = noEntries - 4;
381 
382  for ( i = 0 ; i < noEntries ; i += 4 )
383  {
384  ClipCrossSection(vertices, i, pVoxelLimit, pAxis, pMin, pMax);
385  }
386  for ( i = 0 ; i < noBetweenSections4 ; i += 4 )
387  {
388  ClipBetweenSections(vertices, i, pVoxelLimit, pAxis, pMin, pMax);
389  }
390  if ( (pMin != kInfinity) || (pMax != -kInfinity) )
391  {
392  existsAfterClip = true;
393  pMin -= kCarTolerance; // Add 2*tolerance to avoid precision troubles
394  pMax += kCarTolerance;
395  }
396  else
397  {
398  // Check for case where completely enveloping clipping volume
399  // If point inside then we are confident that the solid completely
400  // envelopes the clipping volume. Hence set min/max extents according
401  // to clipping volume extents along the specified axis.
402 
403  G4ThreeVector clipCentre(
404  (pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
405  (pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
406  (pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5 );
407 
408  if ( Inside(pTransform.Inverse().TransformPoint(clipCentre)) != kOutside )
409  {
410  existsAfterClip = true;
411  pMin = pVoxelLimit.GetMinExtent(pAxis);
412  pMax = pVoxelLimit.GetMaxExtent(pAxis);
413  }
414  }
415  delete vertices;
416  return existsAfterClip;
417  }
418 }
419 
420 
422 //
423 // Return whether point inside/outside/on surface
424 
426 {
427  G4double r2,pPhi,tolRMin,tolRMax;
428  EInside in = kOutside ;
429 
430  if (std::fabs(p.z()) <= fDz - halfCarTolerance)
431  {
432  r2 = p.x()*p.x() + p.y()*p.y() ;
433 
434  if (fRMin) { tolRMin = fRMin + halfRadTolerance ; }
435  else { tolRMin = 0 ; }
436 
437  tolRMax = fRMax - halfRadTolerance ;
438 
439  if ((r2 >= tolRMin*tolRMin) && (r2 <= tolRMax*tolRMax))
440  {
441  if ( fPhiFullTube )
442  {
443  in = kInside ;
444  }
445  else
446  {
447  // Try inner tolerant phi boundaries (=>inside)
448  // if not inside, try outer tolerant phi boundaries
449 
450  if ( (tolRMin==0) && (std::fabs(p.x())<=halfCarTolerance)
451  && (std::fabs(p.y())<=halfCarTolerance) )
452  {
453  in=kSurface;
454  }
455  else
456  {
457  pPhi = std::atan2(p.y(),p.x()) ;
458  if ( pPhi < -halfAngTolerance ) { pPhi += twopi; } // 0<=pPhi<2pi
459 
460  if ( fSPhi >= 0 )
461  {
462  if ( (std::fabs(pPhi) < halfAngTolerance)
463  && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
464  {
465  pPhi += twopi ; // 0 <= pPhi < 2pi
466  }
467  if ( (pPhi >= fSPhi + halfAngTolerance)
468  && (pPhi <= fSPhi + fDPhi - halfAngTolerance) )
469  {
470  in = kInside ;
471  }
472  else if ( (pPhi >= fSPhi - halfAngTolerance)
473  && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
474  {
475  in = kSurface ;
476  }
477  }
478  else // fSPhi < 0
479  {
480  if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
481  && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;} //kOutside
482  else if ( (pPhi <= fSPhi + twopi + halfAngTolerance)
483  && (pPhi >= fSPhi + fDPhi - halfAngTolerance) )
484  {
485  in = kSurface ;
486  }
487  else
488  {
489  in = kInside ;
490  }
491  }
492  }
493  }
494  }
495  else // Try generous boundaries
496  {
497  tolRMin = fRMin - halfRadTolerance ;
498  tolRMax = fRMax + halfRadTolerance ;
499 
500  if ( tolRMin < 0 ) { tolRMin = 0; }
501 
502  if ( (r2 >= tolRMin*tolRMin) && (r2 <= tolRMax*tolRMax) )
503  {
504  if (fPhiFullTube || (r2 <=halfRadTolerance*halfRadTolerance) )
505  { // Continuous in phi or on z-axis
506  in = kSurface ;
507  }
508  else // Try outer tolerant phi boundaries only
509  {
510  pPhi = std::atan2(p.y(),p.x()) ;
511 
512  if ( pPhi < -halfAngTolerance) { pPhi += twopi; } // 0<=pPhi<2pi
513  if ( fSPhi >= 0 )
514  {
515  if ( (std::fabs(pPhi) < halfAngTolerance)
516  && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
517  {
518  pPhi += twopi ; // 0 <= pPhi < 2pi
519  }
520  if ( (pPhi >= fSPhi - halfAngTolerance)
521  && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
522  {
523  in = kSurface ;
524  }
525  }
526  else // fSPhi < 0
527  {
528  if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
529  && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;} // kOutside
530  else
531  {
532  in = kSurface ;
533  }
534  }
535  }
536  }
537  }
538  }
539  else if (std::fabs(p.z()) <= fDz + halfCarTolerance)
540  { // Check within tolerant r limits
541  r2 = p.x()*p.x() + p.y()*p.y() ;
542  tolRMin = fRMin - halfRadTolerance ;
543  tolRMax = fRMax + halfRadTolerance ;
544 
545  if ( tolRMin < 0 ) { tolRMin = 0; }
546 
547  if ( (r2 >= tolRMin*tolRMin) && (r2 <= tolRMax*tolRMax) )
548  {
549  if (fPhiFullTube || (r2 <=halfRadTolerance*halfRadTolerance))
550  { // Continuous in phi or on z-axis
551  in = kSurface ;
552  }
553  else // Try outer tolerant phi boundaries
554  {
555  pPhi = std::atan2(p.y(),p.x()) ;
556 
557  if ( pPhi < -halfAngTolerance ) { pPhi += twopi; } // 0<=pPhi<2pi
558  if ( fSPhi >= 0 )
559  {
560  if ( (std::fabs(pPhi) < halfAngTolerance)
561  && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
562  {
563  pPhi += twopi ; // 0 <= pPhi < 2pi
564  }
565  if ( (pPhi >= fSPhi - halfAngTolerance)
566  && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
567  {
568  in = kSurface;
569  }
570  }
571  else // fSPhi < 0
572  {
573  if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
574  && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;}
575  else
576  {
577  in = kSurface ;
578  }
579  }
580  }
581  }
582  }
583  return in;
584 }
585 
587 //
588 // Return unit normal of surface closest to p
589 // - note if point on z axis, ignore phi divided sides
590 // - unsafe if point close to z axis a rmin=0 - no explicit checks
591 
593 {
594  G4int noSurfaces = 0;
595  G4double rho, pPhi;
596  G4double distZ, distRMin, distRMax;
597  G4double distSPhi = kInfinity, distEPhi = kInfinity;
598 
599  G4ThreeVector norm, sumnorm(0.,0.,0.);
600  G4ThreeVector nZ = G4ThreeVector(0, 0, 1.0);
601  G4ThreeVector nR, nPs, nPe;
602 
603  rho = std::sqrt(p.x()*p.x() + p.y()*p.y());
604 
605  distRMin = std::fabs(rho - fRMin);
606  distRMax = std::fabs(rho - fRMax);
607  distZ = std::fabs(std::fabs(p.z()) - fDz);
608 
609  if (!fPhiFullTube) // Protected against (0,0,z)
610  {
611  if ( rho > halfCarTolerance )
612  {
613  pPhi = std::atan2(p.y(),p.x());
614 
615  if(pPhi < fSPhi- halfCarTolerance) { pPhi += twopi; }
616  else if(pPhi > fSPhi+fDPhi+ halfCarTolerance) { pPhi -= twopi; }
617 
618  distSPhi = std::fabs(pPhi - fSPhi);
619  distEPhi = std::fabs(pPhi - fSPhi - fDPhi);
620  }
621  else if( !fRMin )
622  {
623  distSPhi = 0.;
624  distEPhi = 0.;
625  }
626  nPs = G4ThreeVector(std::sin(fSPhi),-std::cos(fSPhi),0);
627  nPe = G4ThreeVector(-std::sin(fSPhi+fDPhi),std::cos(fSPhi+fDPhi),0);
628  }
629  if ( rho > halfCarTolerance ) { nR = G4ThreeVector(p.x()/rho,p.y()/rho,0); }
630 
631  if( distRMax <= halfCarTolerance )
632  {
633  noSurfaces ++;
634  sumnorm += nR;
635  }
636  if( fRMin && (distRMin <= halfCarTolerance) )
637  {
638  noSurfaces ++;
639  sumnorm -= nR;
640  }
641  if( fDPhi < twopi )
642  {
643  if (distSPhi <= halfAngTolerance)
644  {
645  noSurfaces ++;
646  sumnorm += nPs;
647  }
648  if (distEPhi <= halfAngTolerance)
649  {
650  noSurfaces ++;
651  sumnorm += nPe;
652  }
653  }
654  if (distZ <= halfCarTolerance)
655  {
656  noSurfaces ++;
657  if ( p.z() >= 0.) { sumnorm += nZ; }
658  else { sumnorm -= nZ; }
659  }
660  if ( noSurfaces == 0 )
661  {
662 #ifdef G4CSGDEBUG
663  G4Exception("G4Tubs::SurfaceNormal(p)", "GeomSolids1002",
664  JustWarning, "Point p is not on surface !?" );
665  G4int oldprc = G4cout.precision(20);
666  G4cout<< "G4Tubs::SN ( "<<p.x()<<", "<<p.y()<<", "<<p.z()<<" ); "
667  << G4endl << G4endl;
668  G4cout.precision(oldprc) ;
669 #endif
670  norm = ApproxSurfaceNormal(p);
671  }
672  else if ( noSurfaces == 1 ) { norm = sumnorm; }
673  else { norm = sumnorm.unit(); }
674 
675  return norm;
676 }
677 
679 //
680 // Algorithm for SurfaceNormal() following the original specification
681 // for points not on the surface
682 
684 {
685  ENorm side ;
686  G4ThreeVector norm ;
687  G4double rho, phi ;
688  G4double distZ, distRMin, distRMax, distSPhi, distEPhi, distMin ;
689 
690  rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
691 
692  distRMin = std::fabs(rho - fRMin) ;
693  distRMax = std::fabs(rho - fRMax) ;
694  distZ = std::fabs(std::fabs(p.z()) - fDz) ;
695 
696  if (distRMin < distRMax) // First minimum
697  {
698  if ( distZ < distRMin )
699  {
700  distMin = distZ ;
701  side = kNZ ;
702  }
703  else
704  {
705  distMin = distRMin ;
706  side = kNRMin ;
707  }
708  }
709  else
710  {
711  if ( distZ < distRMax )
712  {
713  distMin = distZ ;
714  side = kNZ ;
715  }
716  else
717  {
718  distMin = distRMax ;
719  side = kNRMax ;
720  }
721  }
722  if (!fPhiFullTube && rho ) // Protected against (0,0,z)
723  {
724  phi = std::atan2(p.y(),p.x()) ;
725 
726  if ( phi < 0 ) { phi += twopi; }
727 
728  if ( fSPhi < 0 )
729  {
730  distSPhi = std::fabs(phi - (fSPhi + twopi))*rho ;
731  }
732  else
733  {
734  distSPhi = std::fabs(phi - fSPhi)*rho ;
735  }
736  distEPhi = std::fabs(phi - fSPhi - fDPhi)*rho ;
737 
738  if (distSPhi < distEPhi) // Find new minimum
739  {
740  if ( distSPhi < distMin )
741  {
742  side = kNSPhi ;
743  }
744  }
745  else
746  {
747  if ( distEPhi < distMin )
748  {
749  side = kNEPhi ;
750  }
751  }
752  }
753  switch ( side )
754  {
755  case kNRMin : // Inner radius
756  {
757  norm = G4ThreeVector(-p.x()/rho, -p.y()/rho, 0) ;
758  break ;
759  }
760  case kNRMax : // Outer radius
761  {
762  norm = G4ThreeVector(p.x()/rho, p.y()/rho, 0) ;
763  break ;
764  }
765  case kNZ : // + or - dz
766  {
767  if ( p.z() > 0 ) { norm = G4ThreeVector(0,0,1) ; }
768  else { norm = G4ThreeVector(0,0,-1); }
769  break ;
770  }
771  case kNSPhi:
772  {
773  norm = G4ThreeVector(std::sin(fSPhi), -std::cos(fSPhi), 0) ;
774  break ;
775  }
776  case kNEPhi:
777  {
778  norm = G4ThreeVector(-std::sin(fSPhi+fDPhi), std::cos(fSPhi+fDPhi), 0) ;
779  break;
780  }
781  default: // Should never reach this case ...
782  {
783  DumpInfo();
784  G4Exception("G4Tubs::ApproxSurfaceNormal()",
785  "GeomSolids1002", JustWarning,
786  "Undefined side for valid surface normal to solid.");
787  break ;
788  }
789  }
790  return norm;
791 }
792 
794 //
795 //
796 // Calculate distance to shape from outside, along normalised vector
797 // - return kInfinity if no intersection, or intersection distance <= tolerance
798 //
799 // - Compute the intersection with the z planes
800 // - if at valid r, phi, return
801 //
802 // -> If point is outer outer radius, compute intersection with rmax
803 // - if at valid phi,z return
804 //
805 // -> Compute intersection with inner radius, taking largest +ve root
806 // - if valid (in z,phi), save intersction
807 //
808 // -> If phi segmented, compute intersections with phi half planes
809 // - return smallest of valid phi intersections and
810 // inner radius intersection
811 //
812 // NOTE:
813 // - 'if valid' implies tolerant checking of intersection points
814 
816  const G4ThreeVector& v ) const
817 {
818  G4double snxt = kInfinity ; // snxt = default return value
819  G4double tolORMin2, tolIRMax2 ; // 'generous' radii squared
820  G4double tolORMax2, tolIRMin2, tolODz, tolIDz ;
821  const G4double dRmax = 100.*fRMax;
822 
823  // Intersection point variables
824  //
825  G4double Dist, sd, xi, yi, zi, rho2, inum, iden, cosPsi, Comp ;
826  G4double t1, t2, t3, b, c, d ; // Quadratic solver variables
827 
828  // Calculate tolerant rmin and rmax
829 
830  if (fRMin > kRadTolerance)
831  {
832  tolORMin2 = (fRMin - halfRadTolerance)*(fRMin - halfRadTolerance) ;
833  tolIRMin2 = (fRMin + halfRadTolerance)*(fRMin + halfRadTolerance) ;
834  }
835  else
836  {
837  tolORMin2 = 0.0 ;
838  tolIRMin2 = 0.0 ;
839  }
840  tolORMax2 = (fRMax + halfRadTolerance)*(fRMax + halfRadTolerance) ;
841  tolIRMax2 = (fRMax - halfRadTolerance)*(fRMax - halfRadTolerance) ;
842 
843  // Intersection with Z surfaces
844 
845  tolIDz = fDz - halfCarTolerance ;
846  tolODz = fDz + halfCarTolerance ;
847 
848  if (std::fabs(p.z()) >= tolIDz)
849  {
850  if ( p.z()*v.z() < 0 ) // at +Z going in -Z or visa versa
851  {
852  sd = (std::fabs(p.z()) - fDz)/std::fabs(v.z()) ; // Z intersect distance
853 
854  if(sd < 0.0) { sd = 0.0; }
855 
856  xi = p.x() + sd*v.x() ; // Intersection coords
857  yi = p.y() + sd*v.y() ;
858  rho2 = xi*xi + yi*yi ;
859 
860  // Check validity of intersection
861 
862  if ((tolIRMin2 <= rho2) && (rho2 <= tolIRMax2))
863  {
864  if (!fPhiFullTube && rho2)
865  {
866  // Psi = angle made with central (average) phi of shape
867  //
868  inum = xi*cosCPhi + yi*sinCPhi ;
869  iden = std::sqrt(rho2) ;
870  cosPsi = inum/iden ;
871  if (cosPsi >= cosHDPhiIT) { return sd ; }
872  }
873  else
874  {
875  return sd ;
876  }
877  }
878  }
879  else
880  {
881  if ( snxt<halfCarTolerance ) { snxt=0; }
882  return snxt ; // On/outside extent, and heading away
883  // -> cannot intersect
884  }
885  }
886 
887  // -> Can not intersect z surfaces
888  //
889  // Intersection with rmax (possible return) and rmin (must also check phi)
890  //
891  // Intersection point (xi,yi,zi) on line x=p.x+t*v.x etc.
892  //
893  // Intersects with x^2+y^2=R^2
894  //
895  // Hence (v.x^2+v.y^2)t^2+ 2t(p.x*v.x+p.y*v.y)+p.x^2+p.y^2-R^2=0
896  // t1 t2 t3
897 
898  t1 = 1.0 - v.z()*v.z() ;
899  t2 = p.x()*v.x() + p.y()*v.y() ;
900  t3 = p.x()*p.x() + p.y()*p.y() ;
901 
902  if ( t1 > 0 ) // Check not || to z axis
903  {
904  b = t2/t1 ;
905  c = t3 - fRMax*fRMax ;
906  if ((t3 >= tolORMax2) && (t2<0)) // This also handles the tangent case
907  {
908  // Try outer cylinder intersection
909  // c=(t3-fRMax*fRMax)/t1;
910 
911  c /= t1 ;
912  d = b*b - c ;
913 
914  if (d >= 0) // If real root
915  {
916  sd = c/(-b+std::sqrt(d));
917  if (sd >= 0) // If 'forwards'
918  {
919  if ( sd>dRmax ) // Avoid rounding errors due to precision issues on
920  { // 64 bits systems. Split long distances and recompute
921  G4double fTerm = sd-std::fmod(sd,dRmax);
922  sd = fTerm + DistanceToIn(p+fTerm*v,v);
923  }
924  // Check z intersection
925  //
926  zi = p.z() + sd*v.z() ;
927  if (std::fabs(zi)<=tolODz)
928  {
929  // Z ok. Check phi intersection if reqd
930  //
931  if (fPhiFullTube)
932  {
933  return sd ;
934  }
935  else
936  {
937  xi = p.x() + sd*v.x() ;
938  yi = p.y() + sd*v.y() ;
939  cosPsi = (xi*cosCPhi + yi*sinCPhi)/fRMax ;
940  if (cosPsi >= cosHDPhiIT) { return sd ; }
941  }
942  } // end if std::fabs(zi)
943  } // end if (sd>=0)
944  } // end if (d>=0)
945  } // end if (r>=fRMax)
946  else
947  {
948  // Inside outer radius :
949  // check not inside, and heading through tubs (-> 0 to in)
950 
951  if ((t3 > tolIRMin2) && (t2 < 0) && (std::fabs(p.z()) <= tolIDz))
952  {
953  // Inside both radii, delta r -ve, inside z extent
954 
955  if (!fPhiFullTube)
956  {
957  inum = p.x()*cosCPhi + p.y()*sinCPhi ;
958  iden = std::sqrt(t3) ;
959  cosPsi = inum/iden ;
960  if (cosPsi >= cosHDPhiIT)
961  {
962  // In the old version, the small negative tangent for the point
963  // on surface was not taken in account, and returning 0.0 ...
964  // New version: check the tangent for the point on surface and
965  // if no intersection, return kInfinity, if intersection instead
966  // return sd.
967  //
968  c = t3-fRMax*fRMax;
969  if ( c<=0.0 )
970  {
971  return 0.0;
972  }
973  else
974  {
975  c = c/t1 ;
976  d = b*b-c;
977  if ( d>=0.0 )
978  {
979  snxt = c/(-b+std::sqrt(d)); // using safe solution
980  // for quadratic equation
981  if ( snxt < halfCarTolerance ) { snxt=0; }
982  return snxt ;
983  }
984  else
985  {
986  return kInfinity;
987  }
988  }
989  }
990  }
991  else
992  {
993  // In the old version, the small negative tangent for the point
994  // on surface was not taken in account, and returning 0.0 ...
995  // New version: check the tangent for the point on surface and
996  // if no intersection, return kInfinity, if intersection instead
997  // return sd.
998  //
999  c = t3 - fRMax*fRMax;
1000  if ( c<=0.0 )
1001  {
1002  return 0.0;
1003  }
1004  else
1005  {
1006  c = c/t1 ;
1007  d = b*b-c;
1008  if ( d>=0.0 )
1009  {
1010  snxt= c/(-b+std::sqrt(d)); // using safe solution
1011  // for quadratic equation
1012  if ( snxt < halfCarTolerance ) { snxt=0; }
1013  return snxt ;
1014  }
1015  else
1016  {
1017  return kInfinity;
1018  }
1019  }
1020  } // end if (!fPhiFullTube)
1021  } // end if (t3>tolIRMin2)
1022  } // end if (Inside Outer Radius)
1023  if ( fRMin ) // Try inner cylinder intersection
1024  {
1025  c = (t3 - fRMin*fRMin)/t1 ;
1026  d = b*b - c ;
1027  if ( d >= 0.0 ) // If real root
1028  {
1029  // Always want 2nd root - we are outside and know rmax Hit was bad
1030  // - If on surface of rmin also need farthest root
1031 
1032  sd =( b > 0. )? c/(-b - std::sqrt(d)) : (-b + std::sqrt(d));
1033  if (sd >= -halfCarTolerance) // check forwards
1034  {
1035  // Check z intersection
1036  //
1037  if(sd < 0.0) { sd = 0.0; }
1038  if ( sd>dRmax ) // Avoid rounding errors due to precision issues seen
1039  { // 64 bits systems. Split long distances and recompute
1040  G4double fTerm = sd-std::fmod(sd,dRmax);
1041  sd = fTerm + DistanceToIn(p+fTerm*v,v);
1042  }
1043  zi = p.z() + sd*v.z() ;
1044  if (std::fabs(zi) <= tolODz)
1045  {
1046  // Z ok. Check phi
1047  //
1048  if ( fPhiFullTube )
1049  {
1050  return sd ;
1051  }
1052  else
1053  {
1054  xi = p.x() + sd*v.x() ;
1055  yi = p.y() + sd*v.y() ;
1056  cosPsi = (xi*cosCPhi + yi*sinCPhi)/fRMin ;
1057  if (cosPsi >= cosHDPhiIT)
1058  {
1059  // Good inner radius isect
1060  // - but earlier phi isect still possible
1061 
1062  snxt = sd ;
1063  }
1064  }
1065  } // end if std::fabs(zi)
1066  } // end if (sd>=0)
1067  } // end if (d>=0)
1068  } // end if (fRMin)
1069  }
1070 
1071  // Phi segment intersection
1072  //
1073  // o Tolerant of points inside phi planes by up to kCarTolerance*0.5
1074  //
1075  // o NOTE: Large duplication of code between sphi & ephi checks
1076  // -> only diffs: sphi -> ephi, Comp -> -Comp and half-plane
1077  // intersection check <=0 -> >=0
1078  // -> use some form of loop Construct ?
1079  //
1080  if ( !fPhiFullTube )
1081  {
1082  // First phi surface (Starting phi)
1083  //
1084  Comp = v.x()*sinSPhi - v.y()*cosSPhi ;
1085 
1086  if ( Comp < 0 ) // Component in outwards normal dirn
1087  {
1088  Dist = (p.y()*cosSPhi - p.x()*sinSPhi) ;
1089 
1090  if ( Dist < halfCarTolerance )
1091  {
1092  sd = Dist/Comp ;
1093 
1094  if (sd < snxt)
1095  {
1096  if ( sd < 0 ) { sd = 0.0; }
1097  zi = p.z() + sd*v.z() ;
1098  if ( std::fabs(zi) <= tolODz )
1099  {
1100  xi = p.x() + sd*v.x() ;
1101  yi = p.y() + sd*v.y() ;
1102  rho2 = xi*xi + yi*yi ;
1103 
1104  if ( ( (rho2 >= tolIRMin2) && (rho2 <= tolIRMax2) )
1105  || ( (rho2 > tolORMin2) && (rho2 < tolIRMin2)
1106  && ( v.y()*cosSPhi - v.x()*sinSPhi > 0 )
1107  && ( v.x()*cosSPhi + v.y()*sinSPhi >= 0 ) )
1108  || ( (rho2 > tolIRMax2) && (rho2 < tolORMax2)
1109  && (v.y()*cosSPhi - v.x()*sinSPhi > 0)
1110  && (v.x()*cosSPhi + v.y()*sinSPhi < 0) ) )
1111  {
1112  // z and r intersections good
1113  // - check intersecting with correct half-plane
1114  //
1115  if ((yi*cosCPhi-xi*sinCPhi) <= halfCarTolerance) { snxt = sd; }
1116  }
1117  }
1118  }
1119  }
1120  }
1121 
1122  // Second phi surface (Ending phi)
1123 
1124  Comp = -(v.x()*sinEPhi - v.y()*cosEPhi) ;
1125 
1126  if (Comp < 0 ) // Component in outwards normal dirn
1127  {
1128  Dist = -(p.y()*cosEPhi - p.x()*sinEPhi) ;
1129 
1130  if ( Dist < halfCarTolerance )
1131  {
1132  sd = Dist/Comp ;
1133 
1134  if (sd < snxt)
1135  {
1136  if ( sd < 0 ) { sd = 0; }
1137  zi = p.z() + sd*v.z() ;
1138  if ( std::fabs(zi) <= tolODz )
1139  {
1140  xi = p.x() + sd*v.x() ;
1141  yi = p.y() + sd*v.y() ;
1142  rho2 = xi*xi + yi*yi ;
1143  if ( ( (rho2 >= tolIRMin2) && (rho2 <= tolIRMax2) )
1144  || ( (rho2 > tolORMin2) && (rho2 < tolIRMin2)
1145  && (v.x()*sinEPhi - v.y()*cosEPhi > 0)
1146  && (v.x()*cosEPhi + v.y()*sinEPhi >= 0) )
1147  || ( (rho2 > tolIRMax2) && (rho2 < tolORMax2)
1148  && (v.x()*sinEPhi - v.y()*cosEPhi > 0)
1149  && (v.x()*cosEPhi + v.y()*sinEPhi < 0) ) )
1150  {
1151  // z and r intersections good
1152  // - check intersecting with correct half-plane
1153  //
1154  if ( (yi*cosCPhi-xi*sinCPhi) >= 0 ) { snxt = sd; }
1155  } //?? >=-halfCarTolerance
1156  }
1157  }
1158  }
1159  } // Comp < 0
1160  } // !fPhiFullTube
1161  if ( snxt<halfCarTolerance ) { snxt=0; }
1162  return snxt ;
1163 }
1164 
1166 //
1167 // Calculate distance to shape from outside, along normalised vector
1168 // - return kInfinity if no intersection, or intersection distance <= tolerance
1169 //
1170 // - Compute the intersection with the z planes
1171 // - if at valid r, phi, return
1172 //
1173 // -> If point is outer outer radius, compute intersection with rmax
1174 // - if at valid phi,z return
1175 //
1176 // -> Compute intersection with inner radius, taking largest +ve root
1177 // - if valid (in z,phi), save intersction
1178 //
1179 // -> If phi segmented, compute intersections with phi half planes
1180 // - return smallest of valid phi intersections and
1181 // inner radius intersection
1182 //
1183 // NOTE:
1184 // - Precalculations for phi trigonometry are Done `just in time'
1185 // - `if valid' implies tolerant checking of intersection points
1186 // Calculate distance (<= actual) to closest surface of shape from outside
1187 // - Calculate distance to z, radial planes
1188 // - Only to phi planes if outside phi extent
1189 // - Return 0 if point inside
1190 
1192 {
1193  G4double safe=0.0, rho, safe1, safe2, safe3 ;
1194  G4double safePhi, cosPsi ;
1195 
1196  rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
1197  safe1 = fRMin - rho ;
1198  safe2 = rho - fRMax ;
1199  safe3 = std::fabs(p.z()) - fDz ;
1200 
1201  if ( safe1 > safe2 ) { safe = safe1; }
1202  else { safe = safe2; }
1203  if ( safe3 > safe ) { safe = safe3; }
1204 
1205  if ( (!fPhiFullTube) && (rho) )
1206  {
1207  // Psi=angle from central phi to point
1208  //
1209  cosPsi = (p.x()*cosCPhi + p.y()*sinCPhi)/rho ;
1210 
1211  if ( cosPsi < std::cos(fDPhi*0.5) )
1212  {
1213  // Point lies outside phi range
1214 
1215  if ( (p.y()*cosCPhi - p.x()*sinCPhi) <= 0 )
1216  {
1217  safePhi = std::fabs(p.x()*sinSPhi - p.y()*cosSPhi) ;
1218  }
1219  else
1220  {
1221  safePhi = std::fabs(p.x()*sinEPhi - p.y()*cosEPhi) ;
1222  }
1223  if ( safePhi > safe ) { safe = safePhi; }
1224  }
1225  }
1226  if ( safe < 0 ) { safe = 0; }
1227  return safe ;
1228 }
1229 
1231 //
1232 // Calculate distance to surface of shape from `inside', allowing for tolerance
1233 // - Only Calc rmax intersection if no valid rmin intersection
1234 
1236  const G4ThreeVector& v,
1237  const G4bool calcNorm,
1238  G4bool *validNorm,
1239  G4ThreeVector *n ) const
1240 {
1241  ESide side=kNull , sider=kNull, sidephi=kNull ;
1242  G4double snxt, srd=kInfinity, sphi=kInfinity, pdist ;
1243  G4double deltaR, t1, t2, t3, b, c, d2, roMin2 ;
1244 
1245  // Vars for phi intersection:
1246 
1247  G4double pDistS, compS, pDistE, compE, sphi2, xi, yi, vphi, roi2 ;
1248 
1249  // Z plane intersection
1250 
1251  if (v.z() > 0 )
1252  {
1253  pdist = fDz - p.z() ;
1254  if ( pdist > halfCarTolerance )
1255  {
1256  snxt = pdist/v.z() ;
1257  side = kPZ ;
1258  }
1259  else
1260  {
1261  if (calcNorm)
1262  {
1263  *n = G4ThreeVector(0,0,1) ;
1264  *validNorm = true ;
1265  }
1266  return snxt = 0 ;
1267  }
1268  }
1269  else if ( v.z() < 0 )
1270  {
1271  pdist = fDz + p.z() ;
1272 
1273  if ( pdist > halfCarTolerance )
1274  {
1275  snxt = -pdist/v.z() ;
1276  side = kMZ ;
1277  }
1278  else
1279  {
1280  if (calcNorm)
1281  {
1282  *n = G4ThreeVector(0,0,-1) ;
1283  *validNorm = true ;
1284  }
1285  return snxt = 0.0 ;
1286  }
1287  }
1288  else
1289  {
1290  snxt = kInfinity ; // Travel perpendicular to z axis
1291  side = kNull;
1292  }
1293 
1294  // Radial Intersections
1295  //
1296  // Find intersection with cylinders at rmax/rmin
1297  // Intersection point (xi,yi,zi) on line x=p.x+t*v.x etc.
1298  //
1299  // Intersects with x^2+y^2=R^2
1300  //
1301  // Hence (v.x^2+v.y^2)t^2+ 2t(p.x*v.x+p.y*v.y)+p.x^2+p.y^2-R^2=0
1302  //
1303  // t1 t2 t3
1304 
1305  t1 = 1.0 - v.z()*v.z() ; // since v normalised
1306  t2 = p.x()*v.x() + p.y()*v.y() ;
1307  t3 = p.x()*p.x() + p.y()*p.y() ;
1308 
1309  if ( snxt > 10*(fDz+fRMax) ) { roi2 = 2*fRMax*fRMax; }
1310  else { roi2 = snxt*snxt*t1 + 2*snxt*t2 + t3; } // radius^2 on +-fDz
1311 
1312  if ( t1 > 0 ) // Check not parallel
1313  {
1314  // Calculate srd, r exit distance
1315 
1316  if ( (t2 >= 0.0) && (roi2 > fRMax*(fRMax + kRadTolerance)) )
1317  {
1318  // Delta r not negative => leaving via rmax
1319 
1320  deltaR = t3 - fRMax*fRMax ;
1321 
1322  // NOTE: Should use rho-fRMax<-kRadTolerance*0.5
1323  // - avoid sqrt for efficiency
1324 
1325  if ( deltaR < -kRadTolerance*fRMax )
1326  {
1327  b = t2/t1 ;
1328  c = deltaR/t1 ;
1329  d2 = b*b-c;
1330  if( d2 >= 0 ) { srd = c/( -b - std::sqrt(d2)); }
1331  else { srd = 0.; }
1332  sider = kRMax ;
1333  }
1334  else
1335  {
1336  // On tolerant boundary & heading outwards (or perpendicular to)
1337  // outer radial surface -> leaving immediately
1338 
1339  if ( calcNorm )
1340  {
1341  *n = G4ThreeVector(p.x()/fRMax,p.y()/fRMax,0) ;
1342  *validNorm = true ;
1343  }
1344  return snxt = 0 ; // Leaving by rmax immediately
1345  }
1346  }
1347  else if ( t2 < 0. ) // i.e. t2 < 0; Possible rmin intersection
1348  {
1349  roMin2 = t3 - t2*t2/t1 ; // min ro2 of the plane of movement
1350 
1351  if ( fRMin && (roMin2 < fRMin*(fRMin - kRadTolerance)) )
1352  {
1353  deltaR = t3 - fRMin*fRMin ;
1354  b = t2/t1 ;
1355  c = deltaR/t1 ;
1356  d2 = b*b - c ;
1357 
1358  if ( d2 >= 0 ) // Leaving via rmin
1359  {
1360  // NOTE: SHould use rho-rmin>kRadTolerance*0.5
1361  // - avoid sqrt for efficiency
1362 
1363  if (deltaR > kRadTolerance*fRMin)
1364  {
1365  srd = c/(-b+std::sqrt(d2));
1366  sider = kRMin ;
1367  }
1368  else
1369  {
1370  if ( calcNorm ) { *validNorm = false; } // Concave side
1371  return snxt = 0.0;
1372  }
1373  }
1374  else // No rmin intersect -> must be rmax intersect
1375  {
1376  deltaR = t3 - fRMax*fRMax ;
1377  c = deltaR/t1 ;
1378  d2 = b*b-c;
1379  if( d2 >=0. )
1380  {
1381  srd = -b + std::sqrt(d2) ;
1382  sider = kRMax ;
1383  }
1384  else // Case: On the border+t2<kRadTolerance
1385  // (v is perpendicular to the surface)
1386  {
1387  if (calcNorm)
1388  {
1389  *n = G4ThreeVector(p.x()/fRMax,p.y()/fRMax,0) ;
1390  *validNorm = true ;
1391  }
1392  return snxt = 0.0;
1393  }
1394  }
1395  }
1396  else if ( roi2 > fRMax*(fRMax + kRadTolerance) )
1397  // No rmin intersect -> must be rmax intersect
1398  {
1399  deltaR = t3 - fRMax*fRMax ;
1400  b = t2/t1 ;
1401  c = deltaR/t1;
1402  d2 = b*b-c;
1403  if( d2 >= 0 )
1404  {
1405  srd = -b + std::sqrt(d2) ;
1406  sider = kRMax ;
1407  }
1408  else // Case: On the border+t2<kRadTolerance
1409  // (v is perpendicular to the surface)
1410  {
1411  if (calcNorm)
1412  {
1413  *n = G4ThreeVector(p.x()/fRMax,p.y()/fRMax,0) ;
1414  *validNorm = true ;
1415  }
1416  return snxt = 0.0;
1417  }
1418  }
1419  }
1420 
1421  // Phi Intersection
1422 
1423  if ( !fPhiFullTube )
1424  {
1425  // add angle calculation with correction
1426  // of the difference in domain of atan2 and Sphi
1427  //
1428  vphi = std::atan2(v.y(),v.x()) ;
1429 
1430  if ( vphi < fSPhi - halfAngTolerance ) { vphi += twopi; }
1431  else if ( vphi > fSPhi + fDPhi + halfAngTolerance ) { vphi -= twopi; }
1432 
1433 
1434  if ( p.x() || p.y() ) // Check if on z axis (rho not needed later)
1435  {
1436  // pDist -ve when inside
1437 
1438  pDistS = p.x()*sinSPhi - p.y()*cosSPhi ;
1439  pDistE = -p.x()*sinEPhi + p.y()*cosEPhi ;
1440 
1441  // Comp -ve when in direction of outwards normal
1442 
1443  compS = -sinSPhi*v.x() + cosSPhi*v.y() ;
1444  compE = sinEPhi*v.x() - cosEPhi*v.y() ;
1445 
1446  sidephi = kNull;
1447 
1448  if( ( (fDPhi <= pi) && ( (pDistS <= halfCarTolerance)
1449  && (pDistE <= halfCarTolerance) ) )
1450  || ( (fDPhi > pi) && !((pDistS > halfCarTolerance)
1451  && (pDistE > halfCarTolerance) ) ) )
1452  {
1453  // Inside both phi *full* planes
1454 
1455  if ( compS < 0 )
1456  {
1457  sphi = pDistS/compS ;
1458 
1459  if (sphi >= -halfCarTolerance)
1460  {
1461  xi = p.x() + sphi*v.x() ;
1462  yi = p.y() + sphi*v.y() ;
1463 
1464  // Check intersecting with correct half-plane
1465  // (if not -> no intersect)
1466  //
1467  if( (std::fabs(xi)<=kCarTolerance)&&(std::fabs(yi)<=kCarTolerance) )
1468  {
1469  sidephi = kSPhi;
1470  if (((fSPhi-halfAngTolerance)<=vphi)
1471  &&((fSPhi+fDPhi+halfAngTolerance)>=vphi))
1472  {
1473  sphi = kInfinity;
1474  }
1475  }
1476  else if ( yi*cosCPhi-xi*sinCPhi >=0 )
1477  {
1478  sphi = kInfinity ;
1479  }
1480  else
1481  {
1482  sidephi = kSPhi ;
1483  if ( pDistS > -halfCarTolerance )
1484  {
1485  sphi = 0.0 ; // Leave by sphi immediately
1486  }
1487  }
1488  }
1489  else
1490  {
1491  sphi = kInfinity ;
1492  }
1493  }
1494  else
1495  {
1496  sphi = kInfinity ;
1497  }
1498 
1499  if ( compE < 0 )
1500  {
1501  sphi2 = pDistE/compE ;
1502 
1503  // Only check further if < starting phi intersection
1504  //
1505  if ( (sphi2 > -halfCarTolerance) && (sphi2 < sphi) )
1506  {
1507  xi = p.x() + sphi2*v.x() ;
1508  yi = p.y() + sphi2*v.y() ;
1509 
1510  if ((std::fabs(xi)<=kCarTolerance)&&(std::fabs(yi)<=kCarTolerance))
1511  {
1512  // Leaving via ending phi
1513  //
1514  if( !((fSPhi-halfAngTolerance <= vphi)
1515  &&(fSPhi+fDPhi+halfAngTolerance >= vphi)) )
1516  {
1517  sidephi = kEPhi ;
1518  if ( pDistE <= -halfCarTolerance ) { sphi = sphi2 ; }
1519  else { sphi = 0.0 ; }
1520  }
1521  }
1522  else // Check intersecting with correct half-plane
1523 
1524  if ( (yi*cosCPhi-xi*sinCPhi) >= 0)
1525  {
1526  // Leaving via ending phi
1527  //
1528  sidephi = kEPhi ;
1529  if ( pDistE <= -halfCarTolerance ) { sphi = sphi2 ; }
1530  else { sphi = 0.0 ; }
1531  }
1532  }
1533  }
1534  }
1535  else
1536  {
1537  sphi = kInfinity ;
1538  }
1539  }
1540  else
1541  {
1542  // On z axis + travel not || to z axis -> if phi of vector direction
1543  // within phi of shape, Step limited by rmax, else Step =0
1544 
1545  if ( (fSPhi - halfAngTolerance <= vphi)
1546  && (vphi <= fSPhi + fDPhi + halfAngTolerance ) )
1547  {
1548  sphi = kInfinity ;
1549  }
1550  else
1551  {
1552  sidephi = kSPhi ; // arbitrary
1553  sphi = 0.0 ;
1554  }
1555  }
1556  if (sphi < snxt) // Order intersecttions
1557  {
1558  snxt = sphi ;
1559  side = sidephi ;
1560  }
1561  }
1562  if (srd < snxt) // Order intersections
1563  {
1564  snxt = srd ;
1565  side = sider ;
1566  }
1567  }
1568  if (calcNorm)
1569  {
1570  switch(side)
1571  {
1572  case kRMax:
1573  // Note: returned vector not normalised
1574  // (divide by fRMax for unit vector)
1575  //
1576  xi = p.x() + snxt*v.x() ;
1577  yi = p.y() + snxt*v.y() ;
1578  *n = G4ThreeVector(xi/fRMax,yi/fRMax,0) ;
1579  *validNorm = true ;
1580  break ;
1581 
1582  case kRMin:
1583  *validNorm = false ; // Rmin is inconvex
1584  break ;
1585 
1586  case kSPhi:
1587  if ( fDPhi <= pi )
1588  {
1589  *n = G4ThreeVector(sinSPhi,-cosSPhi,0) ;
1590  *validNorm = true ;
1591  }
1592  else
1593  {
1594  *validNorm = false ;
1595  }
1596  break ;
1597 
1598  case kEPhi:
1599  if (fDPhi <= pi)
1600  {
1601  *n = G4ThreeVector(-sinEPhi,cosEPhi,0) ;
1602  *validNorm = true ;
1603  }
1604  else
1605  {
1606  *validNorm = false ;
1607  }
1608  break ;
1609 
1610  case kPZ:
1611  *n = G4ThreeVector(0,0,1) ;
1612  *validNorm = true ;
1613  break ;
1614 
1615  case kMZ:
1616  *n = G4ThreeVector(0,0,-1) ;
1617  *validNorm = true ;
1618  break ;
1619 
1620  default:
1621  G4cout << G4endl ;
1622  DumpInfo();
1623  std::ostringstream message;
1624  G4int oldprc = message.precision(16);
1625  message << "Undefined side for valid surface normal to solid."
1626  << G4endl
1627  << "Position:" << G4endl << G4endl
1628  << "p.x() = " << p.x()/mm << " mm" << G4endl
1629  << "p.y() = " << p.y()/mm << " mm" << G4endl
1630  << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl
1631  << "Direction:" << G4endl << G4endl
1632  << "v.x() = " << v.x() << G4endl
1633  << "v.y() = " << v.y() << G4endl
1634  << "v.z() = " << v.z() << G4endl << G4endl
1635  << "Proposed distance :" << G4endl << G4endl
1636  << "snxt = " << snxt/mm << " mm" << G4endl ;
1637  message.precision(oldprc) ;
1638  G4Exception("G4Tubs::DistanceToOut(p,v,..)", "GeomSolids1002",
1639  JustWarning, message);
1640  break ;
1641  }
1642  }
1643  if ( snxt<halfCarTolerance ) { snxt=0 ; }
1644 
1645  return snxt ;
1646 }
1647 
1649 //
1650 // Calculate distance (<=actual) to closest surface of shape from inside
1651 
1653 {
1654  G4double safe=0.0, rho, safeR1, safeR2, safeZ, safePhi ;
1655  rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
1656 
1657 #ifdef G4CSGDEBUG
1658  if( Inside(p) == kOutside )
1659  {
1660  G4int oldprc = G4cout.precision(16) ;
1661  G4cout << G4endl ;
1662  DumpInfo();
1663  G4cout << "Position:" << G4endl << G4endl ;
1664  G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ;
1665  G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ;
1666  G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ;
1667  G4cout.precision(oldprc) ;
1668  G4Exception("G4Tubs::DistanceToOut(p)", "GeomSolids1002",
1669  JustWarning, "Point p is outside !?");
1670  }
1671 #endif
1672 
1673  if ( fRMin )
1674  {
1675  safeR1 = rho - fRMin ;
1676  safeR2 = fRMax - rho ;
1677 
1678  if ( safeR1 < safeR2 ) { safe = safeR1 ; }
1679  else { safe = safeR2 ; }
1680  }
1681  else
1682  {
1683  safe = fRMax - rho ;
1684  }
1685  safeZ = fDz - std::fabs(p.z()) ;
1686 
1687  if ( safeZ < safe ) { safe = safeZ ; }
1688 
1689  // Check if phi divided, Calc distances closest phi plane
1690  //
1691  if ( !fPhiFullTube )
1692  {
1693  if ( p.y()*cosCPhi-p.x()*sinCPhi <= 0 )
1694  {
1695  safePhi = -(p.x()*sinSPhi - p.y()*cosSPhi) ;
1696  }
1697  else
1698  {
1699  safePhi = (p.x()*sinEPhi - p.y()*cosEPhi) ;
1700  }
1701  if (safePhi < safe) { safe = safePhi ; }
1702  }
1703  if ( safe < 0 ) { safe = 0 ; }
1704 
1705  return safe ;
1706 }
1707 
1709 //
1710 // Create a List containing the transformed vertices
1711 // Ordering [0-3] -fDz cross section
1712 // [4-7] +fDz cross section such that [0] is below [4],
1713 // [1] below [5] etc.
1714 // Note:
1715 // Caller has deletion resposibility
1716 // Potential improvement: For last slice, use actual ending angle
1717 // to avoid rounding error problems.
1718 
1721 {
1722  G4ThreeVectorList* vertices ;
1723  G4ThreeVector vertex0, vertex1, vertex2, vertex3 ;
1724  G4double meshAngle, meshRMax, crossAngle,
1725  cosCrossAngle, sinCrossAngle, sAngle;
1726  G4double rMaxX, rMaxY, rMinX, rMinY, meshRMin ;
1727  G4int crossSection, noCrossSections;
1728 
1729  // Compute no of cross-sections necessary to mesh tube
1730  //
1731  noCrossSections = G4int(fDPhi/kMeshAngleDefault) + 1 ;
1732 
1733  if ( noCrossSections < kMinMeshSections )
1734  {
1735  noCrossSections = kMinMeshSections ;
1736  }
1737  else if (noCrossSections>kMaxMeshSections)
1738  {
1739  noCrossSections = kMaxMeshSections ;
1740  }
1741  // noCrossSections = 4 ;
1742 
1743  meshAngle = fDPhi/(noCrossSections - 1) ;
1744  // meshAngle = fDPhi/(noCrossSections) ;
1745 
1746  meshRMax = (fRMax+100*kCarTolerance)/std::cos(meshAngle*0.5) ;
1747  meshRMin = fRMin - 100*kCarTolerance ;
1748 
1749  // If complete in phi, set start angle such that mesh will be at fRMax
1750  // on the x axis. Will give better extent calculations when not rotated.
1751 
1752  if (fPhiFullTube && (fSPhi == 0) ) { sAngle = -meshAngle*0.5 ; }
1753  else { sAngle = fSPhi ; }
1754 
1755  vertices = new G4ThreeVectorList();
1756 
1757  if ( vertices )
1758  {
1759  vertices->reserve(noCrossSections*4);
1760  for (crossSection = 0 ; crossSection < noCrossSections ; crossSection++ )
1761  {
1762  // Compute coordinates of cross section at section crossSection
1763 
1764  crossAngle = sAngle + crossSection*meshAngle ;
1765  cosCrossAngle = std::cos(crossAngle) ;
1766  sinCrossAngle = std::sin(crossAngle) ;
1767 
1768  rMaxX = meshRMax*cosCrossAngle ;
1769  rMaxY = meshRMax*sinCrossAngle ;
1770 
1771  if(meshRMin <= 0.0)
1772  {
1773  rMinX = 0.0 ;
1774  rMinY = 0.0 ;
1775  }
1776  else
1777  {
1778  rMinX = meshRMin*cosCrossAngle ;
1779  rMinY = meshRMin*sinCrossAngle ;
1780  }
1781  vertex0 = G4ThreeVector(rMinX,rMinY,-fDz) ;
1782  vertex1 = G4ThreeVector(rMaxX,rMaxY,-fDz) ;
1783  vertex2 = G4ThreeVector(rMaxX,rMaxY,+fDz) ;
1784  vertex3 = G4ThreeVector(rMinX,rMinY,+fDz) ;
1785 
1786  vertices->push_back(pTransform.TransformPoint(vertex0)) ;
1787  vertices->push_back(pTransform.TransformPoint(vertex1)) ;
1788  vertices->push_back(pTransform.TransformPoint(vertex2)) ;
1789  vertices->push_back(pTransform.TransformPoint(vertex3)) ;
1790  }
1791  }
1792  else
1793  {
1794  DumpInfo();
1795  G4Exception("G4Tubs::CreateRotatedVertices()",
1796  "GeomSolids0003", FatalException,
1797  "Error in allocation of vertices. Out of memory !");
1798  }
1799  return vertices ;
1800 }
1801 
1803 //
1804 // Stream object contents to an output stream
1805 
1807 {
1808  return G4String("G4Tubs");
1809 }
1810 
1812 //
1813 // Make a clone of the object
1814 //
1816 {
1817  return new G4Tubs(*this);
1818 }
1819 
1821 //
1822 // Stream object contents to an output stream
1823 
1824 std::ostream& G4Tubs::StreamInfo( std::ostream& os ) const
1825 {
1826  G4int oldprc = os.precision(16);
1827  os << "-----------------------------------------------------------\n"
1828  << " *** Dump for solid - " << GetName() << " ***\n"
1829  << " ===================================================\n"
1830  << " Solid type: G4Tubs\n"
1831  << " Parameters: \n"
1832  << " inner radius : " << fRMin/mm << " mm \n"
1833  << " outer radius : " << fRMax/mm << " mm \n"
1834  << " half length Z: " << fDz/mm << " mm \n"
1835  << " starting phi : " << fSPhi/degree << " degrees \n"
1836  << " delta phi : " << fDPhi/degree << " degrees \n"
1837  << "-----------------------------------------------------------\n";
1838  os.precision(oldprc);
1839 
1840  return os;
1841 }
1842 
1844 //
1845 // GetPointOnSurface
1846 
1848 {
1849  G4double xRand, yRand, zRand, phi, cosphi, sinphi, chose,
1850  aOne, aTwo, aThr, aFou;
1851  G4double rRand;
1852 
1853  aOne = 2.*fDz*fDPhi*fRMax;
1854  aTwo = 2.*fDz*fDPhi*fRMin;
1855  aThr = 0.5*fDPhi*(fRMax*fRMax-fRMin*fRMin);
1856  aFou = 2.*fDz*(fRMax-fRMin);
1857 
1858  phi = RandFlat::shoot(fSPhi, fSPhi+fDPhi);
1859  cosphi = std::cos(phi);
1860  sinphi = std::sin(phi);
1861 
1862  rRand = GetRadiusInRing(fRMin,fRMax);
1863 
1864  if( (fSPhi == 0) && (fDPhi == twopi) ) { aFou = 0; }
1865 
1866  chose = RandFlat::shoot(0.,aOne+aTwo+2.*aThr+2.*aFou);
1867 
1868  if( (chose >=0) && (chose < aOne) )
1869  {
1870  xRand = fRMax*cosphi;
1871  yRand = fRMax*sinphi;
1872  zRand = RandFlat::shoot(-1.*fDz,fDz);
1873  return G4ThreeVector (xRand, yRand, zRand);
1874  }
1875  else if( (chose >= aOne) && (chose < aOne + aTwo) )
1876  {
1877  xRand = fRMin*cosphi;
1878  yRand = fRMin*sinphi;
1879  zRand = RandFlat::shoot(-1.*fDz,fDz);
1880  return G4ThreeVector (xRand, yRand, zRand);
1881  }
1882  else if( (chose >= aOne + aTwo) && (chose < aOne + aTwo + aThr) )
1883  {
1884  xRand = rRand*cosphi;
1885  yRand = rRand*sinphi;
1886  zRand = fDz;
1887  return G4ThreeVector (xRand, yRand, zRand);
1888  }
1889  else if( (chose >= aOne + aTwo + aThr) && (chose < aOne + aTwo + 2.*aThr) )
1890  {
1891  xRand = rRand*cosphi;
1892  yRand = rRand*sinphi;
1893  zRand = -1.*fDz;
1894  return G4ThreeVector (xRand, yRand, zRand);
1895  }
1896  else if( (chose >= aOne + aTwo + 2.*aThr)
1897  && (chose < aOne + aTwo + 2.*aThr + aFou) )
1898  {
1899  xRand = rRand*std::cos(fSPhi);
1900  yRand = rRand*std::sin(fSPhi);
1901  zRand = RandFlat::shoot(-1.*fDz,fDz);
1902  return G4ThreeVector (xRand, yRand, zRand);
1903  }
1904  else
1905  {
1906  xRand = rRand*std::cos(fSPhi+fDPhi);
1907  yRand = rRand*std::sin(fSPhi+fDPhi);
1908  zRand = RandFlat::shoot(-1.*fDz,fDz);
1909  return G4ThreeVector (xRand, yRand, zRand);
1910  }
1911 }
1912 
1914 //
1915 // Methods for visualisation
1916 
1918 {
1919  scene.AddSolid (*this) ;
1920 }
1921 
1923 {
1924  return new G4PolyhedronTubs (fRMin, fRMax, fDz, fSPhi, fDPhi) ;
1925 }
1926 #endif
G4String GetName() const
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:347
ThreeVector shoot(const G4int Ap, const G4int Af)
G4bool fPhiFullTube
Definition: G4Tubs.hh:226
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double GetMinYExtent() const
CLHEP::Hep3Vector G4ThreeVector
G4double fDPhi
Definition: G4Tubs.hh:217
G4double halfCarTolerance
Definition: G4Tubs.hh:230
G4AffineTransform Inverse() const
G4bool IsYLimited() const
G4double GetRadiusInRing(G4double rmin, G4double rmax) const
Definition: G4CSGSolid.cc:101
G4Tubs(const G4String &pName, G4double pRMin, G4double pRMax, G4double pDz, G4double pSPhi, G4double pDPhi)
Definition: G4Tubs.cc:85
const G4double pi
Definition: G4Tubs.hh:85
G4double cosEPhi
Definition: G4Tubs.hh:221
G4bool IsRotated() const
G4ThreeVector NetTranslation() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tubs.cc:815
G4bool IsXLimited() const
G4double a
Definition: TRTMaterials.hh:39
const G4int kMinMeshSections
Definition: meshdefs.hh:45
G4ThreeVector GetPointOnSurface() const
Definition: G4Tubs.cc:1847
virtual void AddSolid(const G4Box &)=0
virtual G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tubs.cc:683
int G4int
Definition: G4Types.hh:78
G4double cosHDPhiOT
Definition: G4Tubs.hh:221
void DumpInfo() const
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Tubs.cc:1824
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Tubs.cc:1917
G4double GetMaxXExtent() const
G4double sinSPhi
Definition: G4Tubs.hh:221
G4double GetMinZExtent() const
G4GeometryType GetEntityType() const
Definition: G4Tubs.cc:1806
G4GLOB_DLL std::ostream G4cout
G4Polyhedron * fpPolyhedron
Definition: G4CSGSolid.hh:80
G4double fRMin
Definition: G4Tubs.hh:217
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tubs.cc:1922
G4double GetRadialTolerance() const
bool G4bool
Definition: G4Types.hh:79
G4VSolid * Clone() const
Definition: G4Tubs.cc:1815
G4double cosCPhi
Definition: G4Tubs.hh:221
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4CSGSolid.cc:124
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
G4double fSPhi
Definition: G4Tubs.hh:217
const G4int n
virtual ~G4Tubs()
Definition: G4Tubs.cc:139
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Tubs.cc:211
G4double sinEPhi
Definition: G4Tubs.hh:221
G4double halfAngTolerance
Definition: G4Tubs.hh:230
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Tubs.cc:1235
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4double GetMinXExtent() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double GetMaxZExtent() const
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tubs.cc:425
G4double halfRadTolerance
Definition: G4Tubs.hh:230
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
static const double degree
Definition: G4SIunits.hh:125
G4double fRMax
Definition: G4Tubs.hh:217
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tubs.cc:592
#define G4endl
Definition: G4ios.hh:61
G4double GetMaxYExtent() const
G4double kCarTolerance
Definition: G4VSolid.hh:305
G4double kRadTolerance
Definition: G4Tubs.hh:213
G4double fDz
Definition: G4Tubs.hh:217
G4double cosSPhi
Definition: G4Tubs.hh:221
double G4double
Definition: G4Types.hh:76
G4double GetMaxExtent(const EAxis pAxis) const
static const G4double d2
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition: G4CSGSolid.cc:82
G4bool IsZLimited() const
G4double cosHDPhiIT
Definition: G4Tubs.hh:221
static const double mm
Definition: G4SIunits.hh:102
G4double GetAngularTolerance() const
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:378
G4double sinCPhi
Definition: G4Tubs.hh:221
G4double GetMinExtent(const EAxis pAxis) const
const G4double kMeshAngleDefault
Definition: meshdefs.hh:42
static G4GeometryTolerance * GetInstance()
G4double kAngTolerance
Definition: G4Tubs.hh:213
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Tubs.cc:200
G4Tubs & operator=(const G4Tubs &rhs)
Definition: G4Tubs.cc:167
const G4int kMaxMeshSections
Definition: meshdefs.hh:46
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Tubs.cc:1720