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