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