Geant4  10.00.p02
G4Cons.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: G4Cons.cc 81636 2014-06-04 09:06:08Z gcosmo $
28 // GEANT4 tag $Name: $
29 //
30 //
31 // class G4Cons
32 //
33 // Implementation for G4Cons class
34 //
35 // History:
36 //
37 // 05.04.12 M.Kelsey: GetPointOnSurface() throw flat in sqrt(r)
38 // 12.10.09 T.Nikitina: Added to DistanceToIn(p,v) check on the direction in
39 // case of point on surface
40 // 03.05.05 V.Grichine: SurfaceNormal(p) according to J. Apostolakis proposal
41 // 13.09.96 V.Grichine: Review and final modifications
42 // ~1994 P.Kent: Created, as main part of the geometry prototype
43 // --------------------------------------------------------------------
44 
45 #include "G4Cons.hh"
46 
47 #if !defined(G4GEOM_USE_UCONS)
48 
49 #include "G4VoxelLimits.hh"
50 #include "G4AffineTransform.hh"
51 #include "G4GeometryTolerance.hh"
52 
53 #include "G4VPVParameterisation.hh"
54 
55 #include "meshdefs.hh"
56 
57 #include "Randomize.hh"
58 
59 #include "G4VGraphicsScene.hh"
60 
61 using namespace CLHEP;
62 
64 //
65 // Private enum: Not for external use - used by distanceToOut
66 
68 
69 // used by normal
70 
72 
74 //
75 // constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
76 // - note if pDPhi>2PI then reset to 2PI
77 
78 G4Cons::G4Cons( const G4String& pName,
79  G4double pRmin1, G4double pRmax1,
80  G4double pRmin2, G4double pRmax2,
81  G4double pDz,
82  G4double pSPhi, G4double pDPhi)
83  : G4CSGSolid(pName), fRmin1(pRmin1), fRmin2(pRmin2),
84  fRmax1(pRmax1), fRmax2(pRmax2), fDz(pDz), fSPhi(0.), fDPhi(0.)
85 {
88 
92 
93  // Check z-len
94  //
95  if ( pDz < 0 )
96  {
97  std::ostringstream message;
98  message << "Invalid Z half-length for Solid: " << GetName() << G4endl
99  << " hZ = " << pDz;
100  G4Exception("G4Cons::G4Cons()", "GeomSolids0002",
101  FatalException, message);
102  }
103 
104  // Check radii
105  //
106  if (((pRmin1>=pRmax1) || (pRmin2>=pRmax2) || (pRmin1<0)) && (pRmin2<0))
107  {
108  std::ostringstream message;
109  message << "Invalid values of radii for Solid: " << GetName() << G4endl
110  << " pRmin1 = " << pRmin1 << ", pRmin2 = " << pRmin2
111  << ", pRmax1 = " << pRmax1 << ", pRmax2 = " << pRmax2;
112  G4Exception("G4Cons::G4Cons()", "GeomSolids0002",
113  FatalException, message) ;
114  }
115  if( (pRmin1 == 0.0) && (pRmin2 > 0.0) ) { fRmin1 = 1e3*kRadTolerance ; }
116  if( (pRmin2 == 0.0) && (pRmin1 > 0.0) ) { fRmin2 = 1e3*kRadTolerance ; }
117 
118  // Check angles
119  //
120  CheckPhiAngles(pSPhi, pDPhi);
121 }
122 
124 //
125 // Fake default constructor - sets only member data and allocates memory
126 // for usage restricted to object persistency.
127 //
128 G4Cons::G4Cons( __void__& a )
129  : G4CSGSolid(a), kRadTolerance(0.), kAngTolerance(0.),
130  fRmin1(0.), fRmin2(0.), fRmax1(0.), fRmax2(0.), fDz(0.),
131  fSPhi(0.), fDPhi(0.), sinCPhi(0.), cosCPhi(0.), cosHDPhiOT(0.),
132  cosHDPhiIT(0.), sinSPhi(0.), cosSPhi(0.), sinEPhi(0.), cosEPhi(0.),
133  fPhiFullCone(false), halfCarTolerance(0.), halfRadTolerance(0.),
134  halfAngTolerance(0.)
135 
136 {
137 }
138 
140 //
141 // Destructor
142 
144 {
145 }
146 
148 //
149 // Copy constructor
150 
152  : G4CSGSolid(rhs), kRadTolerance(rhs.kRadTolerance),
153  kAngTolerance(rhs.kAngTolerance), fRmin1(rhs.fRmin1), fRmin2(rhs.fRmin2),
154  fRmax1(rhs.fRmax1), fRmax2(rhs.fRmax2), fDz(rhs.fDz), fSPhi(rhs.fSPhi),
155  fDPhi(rhs.fDPhi), sinCPhi(rhs.sinCPhi), cosCPhi(rhs.cosCPhi),
156  cosHDPhiOT(rhs.cosHDPhiOT), cosHDPhiIT(rhs.cosHDPhiIT),
157  sinSPhi(rhs.sinSPhi), cosSPhi(rhs.cosSPhi), sinEPhi(rhs.sinEPhi),
158  cosEPhi(rhs.cosEPhi), fPhiFullCone(rhs.fPhiFullCone),
159  halfCarTolerance(rhs.halfCarTolerance),
160  halfRadTolerance(rhs.halfRadTolerance),
161  halfAngTolerance(rhs.halfAngTolerance)
162 {
164 }
165 
167 //
168 // Assignment operator
169 
171 {
172  // Check assignment to self
173  //
174  if (this == &rhs) { return *this; }
175 
176  // Copy base class data
177  //
179 
180  // Copy data
181  //
184  fRmin1 = rhs.fRmin1; fRmin2 = rhs.fRmin2;
185  fRmax1 = rhs.fRmax1; fRmax2 = rhs.fRmax2;
186  fDz = rhs.fDz; fSPhi = rhs.fSPhi; fDPhi = rhs.fDPhi;
187  sinCPhi = rhs.sinCPhi; cosCPhi = rhs.cosCPhi;
189  sinSPhi = rhs.sinSPhi; cosSPhi = rhs.cosSPhi;
190  sinEPhi = rhs.sinEPhi; cosEPhi = rhs.cosEPhi;
196 
197  return *this;
198 }
199 
201 //
202 // Return whether point inside/outside/on surface
203 
205 {
206  G4double r2, rl, rh, pPhi, tolRMin, tolRMax; // rh2, rl2 ;
207  EInside in;
208 
209  if (std::fabs(p.z()) > fDz + halfCarTolerance ) { return in = kOutside; }
210  else if(std::fabs(p.z()) >= fDz - halfCarTolerance ) { in = kSurface; }
211  else { in = kInside; }
212 
213  r2 = p.x()*p.x() + p.y()*p.y() ;
214  rl = 0.5*(fRmin2*(p.z() + fDz) + fRmin1*(fDz - p.z()))/fDz ;
215  rh = 0.5*(fRmax2*(p.z()+fDz)+fRmax1*(fDz-p.z()))/fDz;
216 
217  // rh2 = rh*rh;
218 
219  tolRMin = rl - halfRadTolerance;
220  if ( tolRMin < 0 ) { tolRMin = 0; }
221  tolRMax = rh + halfRadTolerance;
222 
223  if ( (r2<tolRMin*tolRMin) || (r2>tolRMax*tolRMax) ) { return in = kOutside; }
224 
225  if (rl) { tolRMin = rl + halfRadTolerance; }
226  else { tolRMin = 0.0; }
227  tolRMax = rh - halfRadTolerance;
228 
229  if (in == kInside) // else it's kSurface already
230  {
231  if ( (r2 < tolRMin*tolRMin) || (r2 >= tolRMax*tolRMax) ) { in = kSurface; }
232  }
233  if ( !fPhiFullCone && ((p.x() != 0.0) || (p.y() != 0.0)) )
234  {
235  pPhi = std::atan2(p.y(),p.x()) ;
236 
237  if ( pPhi < fSPhi - halfAngTolerance ) { pPhi += twopi; }
238  else if ( pPhi > fSPhi + fDPhi + halfAngTolerance ) { pPhi -= twopi; }
239 
240  if ( (pPhi < fSPhi - halfAngTolerance) ||
241  (pPhi > fSPhi + fDPhi + halfAngTolerance) ) { return in = kOutside; }
242 
243  else if (in == kInside) // else it's kSurface anyway already
244  {
245  if ( (pPhi < fSPhi + halfAngTolerance) ||
246  (pPhi > fSPhi + fDPhi - halfAngTolerance) ) { in = kSurface; }
247  }
248  }
249  else if ( !fPhiFullCone ) { in = kSurface; }
250 
251  return in ;
252 }
253 
255 //
256 // Dispatch to parameterisation for replication mechanism dimension
257 // computation & modification.
258 
260  const G4int n,
261  const G4VPhysicalVolume* pRep )
262 {
263  p->ComputeDimensions(*this,n,pRep) ;
264 }
265 
266 
268 //
269 // Calculate extent under transform and specified limit
270 
272  const G4VoxelLimits& pVoxelLimit,
273  const G4AffineTransform& pTransform,
274  G4double& pMin,
275  G4double& pMax ) const
276 {
277  if ( !pTransform.IsRotated() && (fDPhi == twopi)
278  && (fRmin1 == 0) && (fRmin2 == 0) )
279  {
280  // Special case handling for unrotated solid cones
281  // Compute z/x/y mins and maxs for bounding box respecting limits,
282  // with early returns if outside limits. Then switch() on pAxis,
283  // and compute exact x and y limit for x/y case
284 
285  G4double xoffset, xMin, xMax ;
286  G4double yoffset, yMin, yMax ;
287  G4double zoffset, zMin, zMax ;
288 
289  G4double diff1, diff2, delta, maxDiff, newMin, newMax, RMax ;
290  G4double xoff1, xoff2, yoff1, yoff2 ;
291 
292  zoffset = pTransform.NetTranslation().z();
293  zMin = zoffset - fDz ;
294  zMax = zoffset + fDz ;
295 
296  if (pVoxelLimit.IsZLimited())
297  {
298  if( (zMin > pVoxelLimit.GetMaxZExtent() + kCarTolerance) ||
299  (zMax < pVoxelLimit.GetMinZExtent() - kCarTolerance) )
300  {
301  return false ;
302  }
303  else
304  {
305  if ( zMin < pVoxelLimit.GetMinZExtent() )
306  {
307  zMin = pVoxelLimit.GetMinZExtent() ;
308  }
309  if ( zMax > pVoxelLimit.GetMaxZExtent() )
310  {
311  zMax = pVoxelLimit.GetMaxZExtent() ;
312  }
313  }
314  }
315  xoffset = pTransform.NetTranslation().x() ;
316  RMax = (fRmax2 >= fRmax1) ? zMax : zMin ;
317  xMax = xoffset + (fRmax1 + fRmax2)*0.5 +
318  (RMax - zoffset)*(fRmax2 - fRmax1)/(2*fDz) ;
319  xMin = 2*xoffset-xMax ;
320 
321  if (pVoxelLimit.IsXLimited())
322  {
323  if ( (xMin > pVoxelLimit.GetMaxXExtent() + kCarTolerance) ||
324  (xMax < pVoxelLimit.GetMinXExtent() - kCarTolerance) )
325  {
326  return false ;
327  }
328  else
329  {
330  if ( xMin < pVoxelLimit.GetMinXExtent() )
331  {
332  xMin = pVoxelLimit.GetMinXExtent() ;
333  }
334  if ( xMax > pVoxelLimit.GetMaxXExtent() )
335  {
336  xMax=pVoxelLimit.GetMaxXExtent() ;
337  }
338  }
339  }
340  yoffset = pTransform.NetTranslation().y() ;
341  yMax = yoffset + (fRmax1 + fRmax2)*0.5 +
342  (RMax - zoffset)*(fRmax2 - fRmax1)/(2*fDz) ;
343  yMin = 2*yoffset-yMax ;
344  RMax = yMax - yoffset ; // = max radius due to Zmax/Zmin cuttings
345 
346  if (pVoxelLimit.IsYLimited())
347  {
348  if ( (yMin > pVoxelLimit.GetMaxYExtent() + kCarTolerance) ||
349  (yMax < pVoxelLimit.GetMinYExtent() - kCarTolerance) )
350  {
351  return false ;
352  }
353  else
354  {
355  if ( yMin < pVoxelLimit.GetMinYExtent() )
356  {
357  yMin = pVoxelLimit.GetMinYExtent() ;
358  }
359  if ( yMax > pVoxelLimit.GetMaxYExtent() )
360  {
361  yMax = pVoxelLimit.GetMaxYExtent() ;
362  }
363  }
364  }
365  switch (pAxis) // Known to cut cones
366  {
367  case kXAxis:
368  yoff1 = yoffset - yMin ;
369  yoff2 = yMax - yoffset ;
370 
371  if ((yoff1 >= 0) && (yoff2 >= 0)) // Y limits cross max/min x
372  { // => no change
373  pMin = xMin ;
374  pMax = xMax ;
375  }
376  else
377  {
378  // Y limits don't cross max/min x => compute max delta x,
379  // hence new mins/maxs
380  delta=RMax*RMax-yoff1*yoff1;
381  diff1=(delta>0.) ? std::sqrt(delta) : 0.;
382  delta=RMax*RMax-yoff2*yoff2;
383  diff2=(delta>0.) ? std::sqrt(delta) : 0.;
384  maxDiff = (diff1>diff2) ? diff1:diff2 ;
385  newMin = xoffset - maxDiff ;
386  newMax = xoffset + maxDiff ;
387  pMin = ( newMin < xMin ) ? xMin : newMin ;
388  pMax = ( newMax > xMax) ? xMax : newMax ;
389  }
390  break ;
391 
392  case kYAxis:
393  xoff1 = xoffset - xMin ;
394  xoff2 = xMax - xoffset ;
395 
396  if ((xoff1 >= 0) && (xoff2 >= 0) ) // X limits cross max/min y
397  { // => no change
398  pMin = yMin ;
399  pMax = yMax ;
400  }
401  else
402  {
403  // X limits don't cross max/min y => compute max delta y,
404  // hence new mins/maxs
405  delta=RMax*RMax-xoff1*xoff1;
406  diff1=(delta>0.) ? std::sqrt(delta) : 0.;
407  delta=RMax*RMax-xoff2*xoff2;
408  diff2=(delta>0.) ? std::sqrt(delta) : 0.;
409  maxDiff = (diff1 > diff2) ? diff1:diff2 ;
410  newMin = yoffset - maxDiff ;
411  newMax = yoffset + maxDiff ;
412  pMin = (newMin < yMin) ? yMin : newMin ;
413  pMax = (newMax > yMax) ? yMax : newMax ;
414  }
415  break ;
416 
417  case kZAxis:
418  pMin = zMin ;
419  pMax = zMax ;
420  break ;
421 
422  default:
423  break ;
424  }
425  pMin -= kCarTolerance ;
426  pMax += kCarTolerance ;
427 
428  return true ;
429  }
430  else // Calculate rotated vertex coordinates
431  {
432  G4int i, noEntries, noBetweenSections4 ;
433  G4bool existsAfterClip = false ;
434  G4ThreeVectorList* vertices = CreateRotatedVertices(pTransform) ;
435 
436  pMin = +kInfinity ;
437  pMax = -kInfinity ;
438 
439  noEntries = vertices->size() ;
440  noBetweenSections4 = noEntries-4 ;
441 
442  for ( i = 0 ; i < noEntries ; i += 4 )
443  {
444  ClipCrossSection(vertices, i, pVoxelLimit, pAxis, pMin, pMax) ;
445  }
446  for ( i = 0 ; i < noBetweenSections4 ; i += 4 )
447  {
448  ClipBetweenSections(vertices, i, pVoxelLimit, pAxis, pMin, pMax) ;
449  }
450  if ( (pMin != kInfinity) || (pMax != -kInfinity) )
451  {
452  existsAfterClip = true ;
453 
454  // Add 2*tolerance to avoid precision troubles
455 
456  pMin -= kCarTolerance ;
457  pMax += kCarTolerance ;
458  }
459  else
460  {
461  // Check for case where completely enveloping clipping volume
462  // If point inside then we are confident that the solid completely
463  // envelopes the clipping volume. Hence set min/max extents according
464  // to clipping volume extents along the specified axis.
465 
466  G4ThreeVector clipCentre(
467  (pVoxelLimit.GetMinXExtent() + pVoxelLimit.GetMaxXExtent())*0.5,
468  (pVoxelLimit.GetMinYExtent() + pVoxelLimit.GetMaxYExtent())*0.5,
469  (pVoxelLimit.GetMinZExtent() + pVoxelLimit.GetMaxZExtent())*0.5 ) ;
470 
471  if (Inside(pTransform.Inverse().TransformPoint(clipCentre)) != kOutside)
472  {
473  existsAfterClip = true ;
474  pMin = pVoxelLimit.GetMinExtent(pAxis) ;
475  pMax = pVoxelLimit.GetMaxExtent(pAxis) ;
476  }
477  }
478  delete vertices ;
479  return existsAfterClip ;
480  }
481 }
482 
484 //
485 // Return unit normal of surface closest to p
486 // - note if point on z axis, ignore phi divided sides
487 // - unsafe if point close to z axis a rmin=0 - no explicit checks
488 
490 {
491  G4int noSurfaces = 0;
492  G4double rho, pPhi;
493  G4double distZ, distRMin, distRMax;
494  G4double distSPhi = kInfinity, distEPhi = kInfinity;
495  G4double tanRMin, secRMin, pRMin, widRMin;
496  G4double tanRMax, secRMax, pRMax, widRMax;
497 
498  G4ThreeVector norm, sumnorm(0.,0.,0.), nZ = G4ThreeVector(0.,0.,1.);
499  G4ThreeVector nR, nr(0.,0.,0.), nPs, nPe;
500 
501  distZ = std::fabs(std::fabs(p.z()) - fDz);
502  rho = std::sqrt(p.x()*p.x() + p.y()*p.y());
503 
504  tanRMin = (fRmin2 - fRmin1)*0.5/fDz;
505  secRMin = std::sqrt(1 + tanRMin*tanRMin);
506  pRMin = rho - p.z()*tanRMin;
507  widRMin = fRmin2 - fDz*tanRMin;
508  distRMin = std::fabs(pRMin - widRMin)/secRMin;
509 
510  tanRMax = (fRmax2 - fRmax1)*0.5/fDz;
511  secRMax = std::sqrt(1+tanRMax*tanRMax);
512  pRMax = rho - p.z()*tanRMax;
513  widRMax = fRmax2 - fDz*tanRMax;
514  distRMax = std::fabs(pRMax - widRMax)/secRMax;
515 
516  if (!fPhiFullCone) // Protected against (0,0,z)
517  {
518  if ( rho )
519  {
520  pPhi = std::atan2(p.y(),p.x());
521 
522  if (pPhi < fSPhi-halfCarTolerance) { pPhi += twopi; }
523  else if (pPhi > fSPhi+fDPhi+halfCarTolerance) { pPhi -= twopi; }
524 
525  distSPhi = std::fabs( pPhi - fSPhi );
526  distEPhi = std::fabs( pPhi - fSPhi - fDPhi );
527  }
528  else if( !(fRmin1) || !(fRmin2) )
529  {
530  distSPhi = 0.;
531  distEPhi = 0.;
532  }
533  nPs = G4ThreeVector(std::sin(fSPhi), -std::cos(fSPhi), 0);
534  nPe = G4ThreeVector(-std::sin(fSPhi+fDPhi), std::cos(fSPhi+fDPhi), 0);
535  }
536  if ( rho > halfCarTolerance )
537  {
538  nR = G4ThreeVector(p.x()/rho/secRMax, p.y()/rho/secRMax, -tanRMax/secRMax);
539  if (fRmin1 || fRmin2)
540  {
541  nr = G4ThreeVector(-p.x()/rho/secRMin,-p.y()/rho/secRMin,tanRMin/secRMin);
542  }
543  }
544 
545  if( distRMax <= halfCarTolerance )
546  {
547  noSurfaces ++;
548  sumnorm += nR;
549  }
550  if( (fRmin1 || fRmin2) && (distRMin <= halfCarTolerance) )
551  {
552  noSurfaces ++;
553  sumnorm += nr;
554  }
555  if( !fPhiFullCone )
556  {
557  if (distSPhi <= halfAngTolerance)
558  {
559  noSurfaces ++;
560  sumnorm += nPs;
561  }
562  if (distEPhi <= halfAngTolerance)
563  {
564  noSurfaces ++;
565  sumnorm += nPe;
566  }
567  }
568  if (distZ <= halfCarTolerance)
569  {
570  noSurfaces ++;
571  if ( p.z() >= 0.) { sumnorm += nZ; }
572  else { sumnorm -= nZ; }
573  }
574  if ( noSurfaces == 0 )
575  {
576 #ifdef G4CSGDEBUG
577  G4Exception("G4Cons::SurfaceNormal(p)", "GeomSolids1002",
578  JustWarning, "Point p is not on surface !?" );
579 #endif
580  norm = ApproxSurfaceNormal(p);
581  }
582  else if ( noSurfaces == 1 ) { norm = sumnorm; }
583  else { norm = sumnorm.unit(); }
584 
585  return norm ;
586 }
587 
589 //
590 // Algorithm for SurfaceNormal() following the original specification
591 // for points not on the surface
592 
594 {
595  ENorm side ;
596  G4ThreeVector norm ;
597  G4double rho, phi ;
598  G4double distZ, distRMin, distRMax, distSPhi, distEPhi, distMin ;
599  G4double tanRMin, secRMin, pRMin, widRMin ;
600  G4double tanRMax, secRMax, pRMax, widRMax ;
601 
602  distZ = std::fabs(std::fabs(p.z()) - fDz) ;
603  rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
604 
605  tanRMin = (fRmin2 - fRmin1)*0.5/fDz ;
606  secRMin = std::sqrt(1 + tanRMin*tanRMin) ;
607  pRMin = rho - p.z()*tanRMin ;
608  widRMin = fRmin2 - fDz*tanRMin ;
609  distRMin = std::fabs(pRMin - widRMin)/secRMin ;
610 
611  tanRMax = (fRmax2 - fRmax1)*0.5/fDz ;
612  secRMax = std::sqrt(1+tanRMax*tanRMax) ;
613  pRMax = rho - p.z()*tanRMax ;
614  widRMax = fRmax2 - fDz*tanRMax ;
615  distRMax = std::fabs(pRMax - widRMax)/secRMax ;
616 
617  if (distRMin < distRMax) // First minimum
618  {
619  if (distZ < distRMin)
620  {
621  distMin = distZ ;
622  side = kNZ ;
623  }
624  else
625  {
626  distMin = distRMin ;
627  side = kNRMin ;
628  }
629  }
630  else
631  {
632  if (distZ < distRMax)
633  {
634  distMin = distZ ;
635  side = kNZ ;
636  }
637  else
638  {
639  distMin = distRMax ;
640  side = kNRMax ;
641  }
642  }
643  if ( !fPhiFullCone && rho ) // Protected against (0,0,z)
644  {
645  phi = std::atan2(p.y(),p.x()) ;
646 
647  if (phi < 0) { phi += twopi; }
648 
649  if (fSPhi < 0) { distSPhi = std::fabs(phi - (fSPhi + twopi))*rho; }
650  else { distSPhi = std::fabs(phi - fSPhi)*rho; }
651 
652  distEPhi = std::fabs(phi - fSPhi - fDPhi)*rho ;
653 
654  // Find new minimum
655 
656  if (distSPhi < distEPhi)
657  {
658  if (distSPhi < distMin) { side = kNSPhi; }
659  }
660  else
661  {
662  if (distEPhi < distMin) { side = kNEPhi; }
663  }
664  }
665  switch (side)
666  {
667  case kNRMin: // Inner radius
668  rho *= secRMin ;
669  norm = G4ThreeVector(-p.x()/rho, -p.y()/rho, tanRMin/secRMin) ;
670  break ;
671  case kNRMax: // Outer radius
672  rho *= secRMax ;
673  norm = G4ThreeVector(p.x()/rho, p.y()/rho, -tanRMax/secRMax) ;
674  break ;
675  case kNZ: // +/- dz
676  if (p.z() > 0) { norm = G4ThreeVector(0,0,1); }
677  else { norm = G4ThreeVector(0,0,-1); }
678  break ;
679  case kNSPhi:
680  norm = G4ThreeVector(std::sin(fSPhi), -std::cos(fSPhi), 0) ;
681  break ;
682  case kNEPhi:
683  norm=G4ThreeVector(-std::sin(fSPhi+fDPhi), std::cos(fSPhi+fDPhi), 0) ;
684  break ;
685  default: // Should never reach this case...
686  DumpInfo();
687  G4Exception("G4Cons::ApproxSurfaceNormal()",
688  "GeomSolids1002", JustWarning,
689  "Undefined side for valid surface normal to solid.");
690  break ;
691  }
692  return norm ;
693 }
694 
696 //
697 // Calculate distance to shape from outside, along normalised vector
698 // - return kInfinity if no intersection, or intersection distance <= tolerance
699 //
700 // - Compute the intersection with the z planes
701 // - if at valid r, phi, return
702 //
703 // -> If point is outside cone, compute intersection with rmax1*0.5
704 // - if at valid phi,z return
705 // - if inside outer cone, handle case when on tolerant outer cone
706 // boundary and heading inwards(->0 to in)
707 //
708 // -> Compute intersection with inner cone, taking largest +ve root
709 // - if valid (in z,phi), save intersction
710 //
711 // -> If phi segmented, compute intersections with phi half planes
712 // - return smallest of valid phi intersections and
713 // inner radius intersection
714 //
715 // NOTE:
716 // - `if valid' implies tolerant checking of intersection points
717 // - z, phi intersection from Tubs
718 
720  const G4ThreeVector& v ) const
721 {
722  G4double snxt = kInfinity ; // snxt = default return value
723  const G4double dRmax = 50*(fRmax1+fRmax2);// 100*(Rmax1+Rmax2)/2.
724 
725  G4double tanRMax,secRMax,rMaxAv,rMaxOAv ; // Data for cones
726  G4double tanRMin,secRMin,rMinAv,rMinOAv ;
727  G4double rout,rin ;
728 
729  G4double tolORMin,tolORMin2,tolIRMin,tolIRMin2 ; // `generous' radii squared
730  G4double tolORMax2,tolIRMax,tolIRMax2 ;
731  G4double tolODz,tolIDz ;
732 
733  G4double Dist,sd,xi,yi,zi,ri=0.,risec,rhoi2,cosPsi ; // Intersection point vars
734 
735  G4double t1,t2,t3,b,c,d ; // Quadratic solver variables
736  G4double nt1,nt2,nt3 ;
737  G4double Comp ;
738 
739  G4ThreeVector Normal;
740 
741  // Cone Precalcs
742 
743  tanRMin = (fRmin2 - fRmin1)*0.5/fDz ;
744  secRMin = std::sqrt(1.0 + tanRMin*tanRMin) ;
745  rMinAv = (fRmin1 + fRmin2)*0.5 ;
746 
747  if (rMinAv > halfRadTolerance)
748  {
749  rMinOAv = rMinAv - halfRadTolerance ;
750  }
751  else
752  {
753  rMinOAv = 0.0 ;
754  }
755  tanRMax = (fRmax2 - fRmax1)*0.5/fDz ;
756  secRMax = std::sqrt(1.0 + tanRMax*tanRMax) ;
757  rMaxAv = (fRmax1 + fRmax2)*0.5 ;
758  rMaxOAv = rMaxAv + halfRadTolerance ;
759 
760  // Intersection with z-surfaces
761 
762  tolIDz = fDz - halfCarTolerance ;
763  tolODz = fDz + halfCarTolerance ;
764 
765  if (std::fabs(p.z()) >= tolIDz)
766  {
767  if ( p.z()*v.z() < 0 ) // at +Z going in -Z or visa versa
768  {
769  sd = (std::fabs(p.z()) - fDz)/std::fabs(v.z()) ; // Z intersect distance
770 
771  if( sd < 0.0 ) { sd = 0.0; } // negative dist -> zero
772 
773  xi = p.x() + sd*v.x() ; // Intersection coords
774  yi = p.y() + sd*v.y() ;
775  rhoi2 = xi*xi + yi*yi ;
776 
777  // Check validity of intersection
778  // Calculate (outer) tolerant radi^2 at intersecion
779 
780  if (v.z() > 0)
781  {
782  tolORMin = fRmin1 - halfRadTolerance*secRMin ;
783  tolIRMin = fRmin1 + halfRadTolerance*secRMin ;
784  tolIRMax = fRmax1 - halfRadTolerance*secRMin ;
785  tolORMax2 = (fRmax1 + halfRadTolerance*secRMax)*
786  (fRmax1 + halfRadTolerance*secRMax) ;
787  }
788  else
789  {
790  tolORMin = fRmin2 - halfRadTolerance*secRMin ;
791  tolIRMin = fRmin2 + halfRadTolerance*secRMin ;
792  tolIRMax = fRmax2 - halfRadTolerance*secRMin ;
793  tolORMax2 = (fRmax2 + halfRadTolerance*secRMax)*
794  (fRmax2 + halfRadTolerance*secRMax) ;
795  }
796  if ( tolORMin > 0 )
797  {
798  tolORMin2 = tolORMin*tolORMin ;
799  tolIRMin2 = tolIRMin*tolIRMin ;
800  }
801  else
802  {
803  tolORMin2 = 0.0 ;
804  tolIRMin2 = 0.0 ;
805  }
806  if ( tolIRMax > 0 ) { tolIRMax2 = tolIRMax*tolIRMax; }
807  else { tolIRMax2 = 0.0; }
808 
809  if ( (tolIRMin2 <= rhoi2) && (rhoi2 <= tolIRMax2) )
810  {
811  if ( !fPhiFullCone && rhoi2 )
812  {
813  // Psi = angle made with central (average) phi of shape
814 
815  cosPsi = (xi*cosCPhi + yi*sinCPhi)/std::sqrt(rhoi2) ;
816 
817  if (cosPsi >= cosHDPhiIT) { return sd; }
818  }
819  else
820  {
821  return sd;
822  }
823  }
824  }
825  else // On/outside extent, and heading away -> cannot intersect
826  {
827  return snxt ;
828  }
829  }
830 
831 // ----> Can not intersect z surfaces
832 
833 
834 // Intersection with outer cone (possible return) and
835 // inner cone (must also check phi)
836 //
837 // Intersection point (xi,yi,zi) on line x=p.x+t*v.x etc.
838 //
839 // Intersects with x^2+y^2=(a*z+b)^2
840 //
841 // where a=tanRMax or tanRMin
842 // b=rMaxAv or rMinAv
843 //
844 // (vx^2+vy^2-(a*vz)^2)t^2+2t(pxvx+pyvy-a*vz(a*pz+b))+px^2+py^2-(a*pz+b)^2=0 ;
845 // t1 t2 t3
846 //
847 // \--------u-------/ \-----------v----------/ \---------w--------/
848 //
849 
850  t1 = 1.0 - v.z()*v.z() ;
851  t2 = p.x()*v.x() + p.y()*v.y() ;
852  t3 = p.x()*p.x() + p.y()*p.y() ;
853  rin = tanRMin*p.z() + rMinAv ;
854  rout = tanRMax*p.z() + rMaxAv ;
855 
856  // Outer Cone Intersection
857  // Must be outside/on outer cone for valid intersection
858 
859  nt1 = t1 - (tanRMax*v.z())*(tanRMax*v.z()) ;
860  nt2 = t2 - tanRMax*v.z()*rout ;
861  nt3 = t3 - rout*rout ;
862 
863  if (std::fabs(nt1) > kRadTolerance) // Equation quadratic => 2 roots
864  {
865  b = nt2/nt1;
866  c = nt3/nt1;
867  d = b*b-c ;
868  if ( (nt3 > rout*rout*kRadTolerance*kRadTolerance*secRMax*secRMax)
869  || (rout < 0) )
870  {
871  // If outside real cone (should be rho-rout>kRadTolerance*0.5
872  // NOT rho^2 etc) saves a std::sqrt() at expense of accuracy
873 
874  if (d >= 0)
875  {
876 
877  if ((rout < 0) && (nt3 <= 0))
878  {
879  // Inside `shadow cone' with -ve radius
880  // -> 2nd root could be on real cone
881 
882  if (b>0) { sd = c/(-b-std::sqrt(d)); }
883  else { sd = -b + std::sqrt(d); }
884  }
885  else
886  {
887  if ((b <= 0) && (c >= 0)) // both >=0, try smaller root
888  {
889  sd=c/(-b+std::sqrt(d));
890  }
891  else
892  {
893  if ( c <= 0 ) // second >=0
894  {
895  sd = -b + std::sqrt(d) ;
896  if((sd<0) & (sd>-halfRadTolerance)) sd=0;
897  }
898  else // both negative, travel away
899  {
900  return kInfinity ;
901  }
902  }
903  }
904  if ( sd > 0 ) // If 'forwards'. Check z intersection
905  {
906  if ( sd>dRmax ) // Avoid rounding errors due to precision issues on
907  { // 64 bits systems. Split long distances and recompute
908  G4double fTerm = sd-std::fmod(sd,dRmax);
909  sd = fTerm + DistanceToIn(p+fTerm*v,v);
910  }
911  zi = p.z() + sd*v.z() ;
912 
913  if (std::fabs(zi) <= tolODz)
914  {
915  // Z ok. Check phi intersection if reqd
916 
917  if ( fPhiFullCone ) { return sd; }
918  else
919  {
920  xi = p.x() + sd*v.x() ;
921  yi = p.y() + sd*v.y() ;
922  ri = rMaxAv + zi*tanRMax ;
923  cosPsi = (xi*cosCPhi + yi*sinCPhi)/ri ;
924 
925  if ( cosPsi >= cosHDPhiIT ) { return sd; }
926  }
927  }
928  } // end if (sd>0)
929  }
930  }
931  else
932  {
933  // Inside outer cone
934  // check not inside, and heading through G4Cons (-> 0 to in)
935 
936  if ( ( t3 > (rin + halfRadTolerance*secRMin)*
937  (rin + halfRadTolerance*secRMin) )
938  && (nt2 < 0) && (d >= 0) && (std::fabs(p.z()) <= tolIDz) )
939  {
940  // Inside cones, delta r -ve, inside z extent
941  // Point is on the Surface => check Direction using Normal.dot(v)
942 
943  xi = p.x() ;
944  yi = p.y() ;
945  risec = std::sqrt(xi*xi + yi*yi)*secRMax ;
946  Normal = G4ThreeVector(xi/risec,yi/risec,-tanRMax/secRMax) ;
947  if ( !fPhiFullCone )
948  {
949  cosPsi = (p.x()*cosCPhi + p.y()*sinCPhi)/std::sqrt(t3) ;
950  if ( cosPsi >= cosHDPhiIT )
951  {
952  if ( Normal.dot(v) <= 0 ) { return 0.0; }
953  }
954  }
955  else
956  {
957  if ( Normal.dot(v) <= 0 ) { return 0.0; }
958  }
959  }
960  }
961  }
962  else // Single root case
963  {
964  if ( std::fabs(nt2) > kRadTolerance )
965  {
966  sd = -0.5*nt3/nt2 ;
967 
968  if ( sd < 0 ) { return kInfinity; } // travel away
969  else // sd >= 0, If 'forwards'. Check z intersection
970  {
971  zi = p.z() + sd*v.z() ;
972 
973  if ((std::fabs(zi) <= tolODz) && (nt2 < 0))
974  {
975  // Z ok. Check phi intersection if reqd
976 
977  if ( fPhiFullCone ) { return sd; }
978  else
979  {
980  xi = p.x() + sd*v.x() ;
981  yi = p.y() + sd*v.y() ;
982  ri = rMaxAv + zi*tanRMax ;
983  cosPsi = (xi*cosCPhi + yi*sinCPhi)/ri ;
984 
985  if (cosPsi >= cosHDPhiIT) { return sd; }
986  }
987  }
988  }
989  }
990  else // travel || cone surface from its origin
991  {
992  sd = kInfinity ;
993  }
994  }
995 
996  // Inner Cone Intersection
997  // o Space is divided into 3 areas:
998  // 1) Radius greater than real inner cone & imaginary cone & outside
999  // tolerance
1000  // 2) Radius less than inner or imaginary cone & outside tolarance
1001  // 3) Within tolerance of real or imaginary cones
1002  // - Extra checks needed for 3's intersections
1003  // => lots of duplicated code
1004 
1005  if (rMinAv)
1006  {
1007  nt1 = t1 - (tanRMin*v.z())*(tanRMin*v.z()) ;
1008  nt2 = t2 - tanRMin*v.z()*rin ;
1009  nt3 = t3 - rin*rin ;
1010 
1011  if ( nt1 )
1012  {
1013  if ( nt3 > rin*kRadTolerance*secRMin )
1014  {
1015  // At radius greater than real & imaginary cones
1016  // -> 2nd root, with zi check
1017 
1018  b = nt2/nt1 ;
1019  c = nt3/nt1 ;
1020  d = b*b-c ;
1021  if (d >= 0) // > 0
1022  {
1023  if(b>0){sd = c/( -b-std::sqrt(d));}
1024  else {sd = -b + std::sqrt(d) ;}
1025 
1026  if ( sd >= 0 ) // > 0
1027  {
1028  if ( sd>dRmax ) // Avoid rounding errors due to precision issues on
1029  { // 64 bits systems. Split long distance and recompute
1030  G4double fTerm = sd-std::fmod(sd,dRmax);
1031  sd = fTerm + DistanceToIn(p+fTerm*v,v);
1032  }
1033  zi = p.z() + sd*v.z() ;
1034 
1035  if ( std::fabs(zi) <= tolODz )
1036  {
1037  if ( !fPhiFullCone )
1038  {
1039  xi = p.x() + sd*v.x() ;
1040  yi = p.y() + sd*v.y() ;
1041  ri = rMinAv + zi*tanRMin ;
1042  cosPsi = (xi*cosCPhi + yi*sinCPhi)/ri ;
1043 
1044  if (cosPsi >= cosHDPhiIT)
1045  {
1046  if ( sd > halfRadTolerance ) { snxt=sd; }
1047  else
1048  {
1049  // Calculate a normal vector in order to check Direction
1050 
1051  risec = std::sqrt(xi*xi + yi*yi)*secRMin ;
1052  Normal = G4ThreeVector(-xi/risec,-yi/risec,tanRMin/secRMin);
1053  if ( Normal.dot(v) <= 0 ) { snxt = sd; }
1054  }
1055  }
1056  }
1057  else
1058  {
1059  if ( sd > halfRadTolerance ) { return sd; }
1060  else
1061  {
1062  // Calculate a normal vector in order to check Direction
1063 
1064  xi = p.x() + sd*v.x() ;
1065  yi = p.y() + sd*v.y() ;
1066  risec = std::sqrt(xi*xi + yi*yi)*secRMin ;
1067  Normal = G4ThreeVector(-xi/risec,-yi/risec,tanRMin/secRMin) ;
1068  if ( Normal.dot(v) <= 0 ) { return sd; }
1069  }
1070  }
1071  }
1072  }
1073  }
1074  }
1075  else if ( nt3 < -rin*kRadTolerance*secRMin )
1076  {
1077  // Within radius of inner cone (real or imaginary)
1078  // -> Try 2nd root, with checking intersection is with real cone
1079  // -> If check fails, try 1st root, also checking intersection is
1080  // on real cone
1081 
1082  b = nt2/nt1 ;
1083  c = nt3/nt1 ;
1084  d = b*b - c ;
1085 
1086  if ( d >= 0 ) // > 0
1087  {
1088  if (b>0) { sd = c/(-b-std::sqrt(d)); }
1089  else { sd = -b + std::sqrt(d); }
1090  zi = p.z() + sd*v.z() ;
1091  ri = rMinAv + zi*tanRMin ;
1092 
1093  if ( ri > 0 )
1094  {
1095  if ( (sd >= 0) && (std::fabs(zi) <= tolODz) ) // sd > 0
1096  {
1097  if ( sd>dRmax ) // Avoid rounding errors due to precision issues
1098  { // seen on 64 bits systems. Split and recompute
1099  G4double fTerm = sd-std::fmod(sd,dRmax);
1100  sd = fTerm + DistanceToIn(p+fTerm*v,v);
1101  }
1102  if ( !fPhiFullCone )
1103  {
1104  xi = p.x() + sd*v.x() ;
1105  yi = p.y() + sd*v.y() ;
1106  cosPsi = (xi*cosCPhi + yi*sinCPhi)/ri ;
1107 
1108  if (cosPsi >= cosHDPhiOT)
1109  {
1110  if ( sd > halfRadTolerance ) { snxt=sd; }
1111  else
1112  {
1113  // Calculate a normal vector in order to check Direction
1114 
1115  risec = std::sqrt(xi*xi + yi*yi)*secRMin ;
1116  Normal = G4ThreeVector(-xi/risec,-yi/risec,tanRMin/secRMin);
1117  if ( Normal.dot(v) <= 0 ) { snxt = sd; }
1118  }
1119  }
1120  }
1121  else
1122  {
1123  if( sd > halfRadTolerance ) { return sd; }
1124  else
1125  {
1126  // Calculate a normal vector in order to check Direction
1127 
1128  xi = p.x() + sd*v.x() ;
1129  yi = p.y() + sd*v.y() ;
1130  risec = std::sqrt(xi*xi + yi*yi)*secRMin ;
1131  Normal = G4ThreeVector(-xi/risec,-yi/risec,tanRMin/secRMin) ;
1132  if ( Normal.dot(v) <= 0 ) { return sd; }
1133  }
1134  }
1135  }
1136  }
1137  else
1138  {
1139  if (b>0) { sd = -b - std::sqrt(d); }
1140  else { sd = c/(-b+std::sqrt(d)); }
1141  zi = p.z() + sd*v.z() ;
1142  ri = rMinAv + zi*tanRMin ;
1143 
1144  if ( (sd >= 0) && (ri > 0) && (std::fabs(zi) <= tolODz) ) // sd>0
1145  {
1146  if ( sd>dRmax ) // Avoid rounding errors due to precision issues
1147  { // seen on 64 bits systems. Split and recompute
1148  G4double fTerm = sd-std::fmod(sd,dRmax);
1149  sd = fTerm + DistanceToIn(p+fTerm*v,v);
1150  }
1151  if ( !fPhiFullCone )
1152  {
1153  xi = p.x() + sd*v.x() ;
1154  yi = p.y() + sd*v.y() ;
1155  cosPsi = (xi*cosCPhi + yi*sinCPhi)/ri ;
1156 
1157  if (cosPsi >= cosHDPhiIT)
1158  {
1159  if ( sd > halfRadTolerance ) { snxt=sd; }
1160  else
1161  {
1162  // Calculate a normal vector in order to check Direction
1163 
1164  risec = std::sqrt(xi*xi + yi*yi)*secRMin ;
1165  Normal = G4ThreeVector(-xi/risec,-yi/risec,tanRMin/secRMin);
1166  if ( Normal.dot(v) <= 0 ) { snxt = sd; }
1167  }
1168  }
1169  }
1170  else
1171  {
1172  if ( sd > halfRadTolerance ) { return sd; }
1173  else
1174  {
1175  // Calculate a normal vector in order to check Direction
1176 
1177  xi = p.x() + sd*v.x() ;
1178  yi = p.y() + sd*v.y() ;
1179  risec = std::sqrt(xi*xi + yi*yi)*secRMin ;
1180  Normal = G4ThreeVector(-xi/risec,-yi/risec,tanRMin/secRMin) ;
1181  if ( Normal.dot(v) <= 0 ) { return sd; }
1182  }
1183  }
1184  }
1185  }
1186  }
1187  }
1188  else
1189  {
1190  // Within kRadTol*0.5 of inner cone (real OR imaginary)
1191  // ----> Check not travelling through (=>0 to in)
1192  // ----> if not:
1193  // -2nd root with validity check
1194 
1195  if ( std::fabs(p.z()) <= tolODz )
1196  {
1197  if ( nt2 > 0 )
1198  {
1199  // Inside inner real cone, heading outwards, inside z range
1200 
1201  if ( !fPhiFullCone )
1202  {
1203  cosPsi = (p.x()*cosCPhi + p.y()*sinCPhi)/std::sqrt(t3) ;
1204 
1205  if (cosPsi >= cosHDPhiIT) { return 0.0; }
1206  }
1207  else { return 0.0; }
1208  }
1209  else
1210  {
1211  // Within z extent, but not travelling through
1212  // -> 2nd root or kInfinity if 1st root on imaginary cone
1213 
1214  b = nt2/nt1 ;
1215  c = nt3/nt1 ;
1216  d = b*b - c ;
1217 
1218  if ( d >= 0 ) // > 0
1219  {
1220  if (b>0) { sd = -b - std::sqrt(d); }
1221  else { sd = c/(-b+std::sqrt(d)); }
1222  zi = p.z() + sd*v.z() ;
1223  ri = rMinAv + zi*tanRMin ;
1224 
1225  if ( ri > 0 ) // 2nd root
1226  {
1227  if (b>0) { sd = c/(-b-std::sqrt(d)); }
1228  else { sd = -b + std::sqrt(d); }
1229 
1230  zi = p.z() + sd*v.z() ;
1231 
1232  if ( (sd >= 0) && (std::fabs(zi) <= tolODz) ) // sd>0
1233  {
1234  if ( sd>dRmax ) // Avoid rounding errors due to precision issue
1235  { // seen on 64 bits systems. Split and recompute
1236  G4double fTerm = sd-std::fmod(sd,dRmax);
1237  sd = fTerm + DistanceToIn(p+fTerm*v,v);
1238  }
1239  if ( !fPhiFullCone )
1240  {
1241  xi = p.x() + sd*v.x() ;
1242  yi = p.y() + sd*v.y() ;
1243  ri = rMinAv + zi*tanRMin ;
1244  cosPsi = (xi*cosCPhi + yi*sinCPhi)/ri ;
1245 
1246  if ( cosPsi >= cosHDPhiIT ) { snxt = sd; }
1247  }
1248  else { return sd; }
1249  }
1250  }
1251  else { return kInfinity; }
1252  }
1253  }
1254  }
1255  else // 2nd root
1256  {
1257  b = nt2/nt1 ;
1258  c = nt3/nt1 ;
1259  d = b*b - c ;
1260 
1261  if ( d > 0 )
1262  {
1263  if (b>0) { sd = c/(-b-std::sqrt(d)); }
1264  else { sd = -b + std::sqrt(d) ; }
1265  zi = p.z() + sd*v.z() ;
1266 
1267  if ( (sd >= 0) && (std::fabs(zi) <= tolODz) ) // sd>0
1268  {
1269  if ( sd>dRmax ) // Avoid rounding errors due to precision issues
1270  { // seen on 64 bits systems. Split and recompute
1271  G4double fTerm = sd-std::fmod(sd,dRmax);
1272  sd = fTerm + DistanceToIn(p+fTerm*v,v);
1273  }
1274  if ( !fPhiFullCone )
1275  {
1276  xi = p.x() + sd*v.x();
1277  yi = p.y() + sd*v.y();
1278  ri = rMinAv + zi*tanRMin ;
1279  cosPsi = (xi*cosCPhi + yi*sinCPhi)/ri;
1280 
1281  if (cosPsi >= cosHDPhiIT) { snxt = sd; }
1282  }
1283  else { return sd; }
1284  }
1285  }
1286  }
1287  }
1288  }
1289  }
1290 
1291  // Phi segment intersection
1292  //
1293  // o Tolerant of points inside phi planes by up to kCarTolerance*0.5
1294  //
1295  // o NOTE: Large duplication of code between sphi & ephi checks
1296  // -> only diffs: sphi -> ephi, Comp -> -Comp and half-plane
1297  // intersection check <=0 -> >=0
1298  // -> Should use some form of loop Construct
1299 
1300  if ( !fPhiFullCone )
1301  {
1302  // First phi surface (starting phi)
1303 
1304  Comp = v.x()*sinSPhi - v.y()*cosSPhi ;
1305 
1306  if ( Comp < 0 ) // Component in outwards normal dirn
1307  {
1308  Dist = (p.y()*cosSPhi - p.x()*sinSPhi) ;
1309 
1310  if (Dist < halfCarTolerance)
1311  {
1312  sd = Dist/Comp ;
1313 
1314  if ( sd < snxt )
1315  {
1316  if ( sd < 0 ) { sd = 0.0; }
1317 
1318  zi = p.z() + sd*v.z() ;
1319 
1320  if ( std::fabs(zi) <= tolODz )
1321  {
1322  xi = p.x() + sd*v.x() ;
1323  yi = p.y() + sd*v.y() ;
1324  rhoi2 = xi*xi + yi*yi ;
1325  tolORMin2 = (rMinOAv + zi*tanRMin)*(rMinOAv + zi*tanRMin) ;
1326  tolORMax2 = (rMaxOAv + zi*tanRMax)*(rMaxOAv + zi*tanRMax) ;
1327 
1328  if ( (rhoi2 >= tolORMin2) && (rhoi2 <= tolORMax2) )
1329  {
1330  // z and r intersections good - check intersecting with
1331  // correct half-plane
1332 
1333  if ((yi*cosCPhi - xi*sinCPhi) <= 0 ) { snxt = sd; }
1334  }
1335  }
1336  }
1337  }
1338  }
1339 
1340  // Second phi surface (Ending phi)
1341 
1342  Comp = -(v.x()*sinEPhi - v.y()*cosEPhi) ;
1343 
1344  if ( Comp < 0 ) // Component in outwards normal dirn
1345  {
1346  Dist = -(p.y()*cosEPhi - p.x()*sinEPhi) ;
1347  if (Dist < halfCarTolerance)
1348  {
1349  sd = Dist/Comp ;
1350 
1351  if ( sd < snxt )
1352  {
1353  if ( sd < 0 ) { sd = 0.0; }
1354 
1355  zi = p.z() + sd*v.z() ;
1356 
1357  if (std::fabs(zi) <= tolODz)
1358  {
1359  xi = p.x() + sd*v.x() ;
1360  yi = p.y() + sd*v.y() ;
1361  rhoi2 = xi*xi + yi*yi ;
1362  tolORMin2 = (rMinOAv + zi*tanRMin)*(rMinOAv + zi*tanRMin) ;
1363  tolORMax2 = (rMaxOAv + zi*tanRMax)*(rMaxOAv + zi*tanRMax) ;
1364 
1365  if ( (rhoi2 >= tolORMin2) && (rhoi2 <= tolORMax2) )
1366  {
1367  // z and r intersections good - check intersecting with
1368  // correct half-plane
1369 
1370  if ( (yi*cosCPhi - xi*sinCPhi) >= 0.0 ) { snxt = sd; }
1371  }
1372  }
1373  }
1374  }
1375  }
1376  }
1377  if (snxt < halfCarTolerance) { snxt = 0.; }
1378 
1379  return snxt ;
1380 }
1381 
1383 //
1384 // Calculate distance (<= actual) to closest surface of shape from outside
1385 // - Calculate distance to z, radial planes
1386 // - Only to phi planes if outside phi extent
1387 // - Return 0 if point inside
1388 
1390 {
1391  G4double safe=0.0, rho, safeR1, safeR2, safeZ, safePhi, cosPsi ;
1392  G4double tanRMin, secRMin, pRMin ;
1393  G4double tanRMax, secRMax, pRMax ;
1394 
1395  rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
1396  safeZ = std::fabs(p.z()) - fDz ;
1397 
1398  if ( fRmin1 || fRmin2 )
1399  {
1400  tanRMin = (fRmin2 - fRmin1)*0.5/fDz ;
1401  secRMin = std::sqrt(1.0 + tanRMin*tanRMin) ;
1402  pRMin = tanRMin*p.z() + (fRmin1 + fRmin2)*0.5 ;
1403  safeR1 = (pRMin - rho)/secRMin ;
1404 
1405  tanRMax = (fRmax2 - fRmax1)*0.5/fDz ;
1406  secRMax = std::sqrt(1.0 + tanRMax*tanRMax) ;
1407  pRMax = tanRMax*p.z() + (fRmax1 + fRmax2)*0.5 ;
1408  safeR2 = (rho - pRMax)/secRMax ;
1409 
1410  if ( safeR1 > safeR2) { safe = safeR1; }
1411  else { safe = safeR2; }
1412  }
1413  else
1414  {
1415  tanRMax = (fRmax2 - fRmax1)*0.5/fDz ;
1416  secRMax = std::sqrt(1.0 + tanRMax*tanRMax) ;
1417  pRMax = tanRMax*p.z() + (fRmax1 + fRmax2)*0.5 ;
1418  safe = (rho - pRMax)/secRMax ;
1419  }
1420  if ( safeZ > safe ) { safe = safeZ; }
1421 
1422  if ( !fPhiFullCone && rho )
1423  {
1424  // Psi=angle from central phi to point
1425 
1426  cosPsi = (p.x()*cosCPhi + p.y()*sinCPhi)/rho ;
1427 
1428  if ( cosPsi < std::cos(fDPhi*0.5) ) // Point lies outside phi range
1429  {
1430  if ( (p.y()*cosCPhi - p.x()*sinCPhi) <= 0.0 )
1431  {
1432  safePhi = std::fabs(p.x()*std::sin(fSPhi)-p.y()*std::cos(fSPhi));
1433  }
1434  else
1435  {
1436  safePhi = std::fabs(p.x()*sinEPhi-p.y()*cosEPhi);
1437  }
1438  if ( safePhi > safe ) { safe = safePhi; }
1439  }
1440  }
1441  if ( safe < 0.0 ) { safe = 0.0; }
1442 
1443  return safe ;
1444 }
1445 
1447 //
1448 // Calculate distance to surface of shape from 'inside', allowing for tolerance
1449 // - Only Calc rmax intersection if no valid rmin intersection
1450 
1452  const G4ThreeVector& v,
1453  const G4bool calcNorm,
1454  G4bool *validNorm,
1455  G4ThreeVector *n) const
1456 {
1457  ESide side = kNull, sider = kNull, sidephi = kNull;
1458 
1459  G4double snxt,srd,sphi,pdist ;
1460 
1461  G4double tanRMax, secRMax, rMaxAv ; // Data for outer cone
1462  G4double tanRMin, secRMin, rMinAv ; // Data for inner cone
1463 
1464  G4double t1, t2, t3, rout, rin, nt1, nt2, nt3 ;
1465  G4double b, c, d, sr2, sr3 ;
1466 
1467  // Vars for intersection within tolerance
1468 
1469  ESide sidetol = kNull ;
1470  G4double slentol = kInfinity ;
1471 
1472  // Vars for phi intersection:
1473 
1474  G4double pDistS, compS, pDistE, compE, sphi2, xi, yi, risec, vphi ;
1475  G4double zi, ri, deltaRoi2 ;
1476 
1477  // Z plane intersection
1478 
1479  if ( v.z() > 0.0 )
1480  {
1481  pdist = fDz - p.z() ;
1482 
1483  if (pdist > halfCarTolerance)
1484  {
1485  snxt = pdist/v.z() ;
1486  side = kPZ ;
1487  }
1488  else
1489  {
1490  if (calcNorm)
1491  {
1492  *n = G4ThreeVector(0,0,1) ;
1493  *validNorm = true ;
1494  }
1495  return snxt = 0.0;
1496  }
1497  }
1498  else if ( v.z() < 0.0 )
1499  {
1500  pdist = fDz + p.z() ;
1501 
1502  if ( pdist > halfCarTolerance)
1503  {
1504  snxt = -pdist/v.z() ;
1505  side = kMZ ;
1506  }
1507  else
1508  {
1509  if ( calcNorm )
1510  {
1511  *n = G4ThreeVector(0,0,-1) ;
1512  *validNorm = true ;
1513  }
1514  return snxt = 0.0 ;
1515  }
1516  }
1517  else // Travel perpendicular to z axis
1518  {
1519  snxt = kInfinity ;
1520  side = kNull ;
1521  }
1522 
1523  // Radial Intersections
1524  //
1525  // Intersection with outer cone (possible return) and
1526  // inner cone (must also check phi)
1527  //
1528  // Intersection point (xi,yi,zi) on line x=p.x+t*v.x etc.
1529  //
1530  // Intersects with x^2+y^2=(a*z+b)^2
1531  //
1532  // where a=tanRMax or tanRMin
1533  // b=rMaxAv or rMinAv
1534  //
1535  // (vx^2+vy^2-(a*vz)^2)t^2+2t(pxvx+pyvy-a*vz(a*pz+b))+px^2+py^2-(a*pz+b)^2=0 ;
1536  // t1 t2 t3
1537  //
1538  // \--------u-------/ \-----------v----------/ \---------w--------/
1539 
1540  tanRMax = (fRmax2 - fRmax1)*0.5/fDz ;
1541  secRMax = std::sqrt(1.0 + tanRMax*tanRMax) ;
1542  rMaxAv = (fRmax1 + fRmax2)*0.5 ;
1543 
1544 
1545  t1 = 1.0 - v.z()*v.z() ; // since v normalised
1546  t2 = p.x()*v.x() + p.y()*v.y() ;
1547  t3 = p.x()*p.x() + p.y()*p.y() ;
1548  rout = tanRMax*p.z() + rMaxAv ;
1549 
1550  nt1 = t1 - (tanRMax*v.z())*(tanRMax*v.z()) ;
1551  nt2 = t2 - tanRMax*v.z()*rout ;
1552  nt3 = t3 - rout*rout ;
1553 
1554  if (v.z() > 0.0)
1555  {
1556  deltaRoi2 = snxt*snxt*t1 + 2*snxt*t2 + t3
1557  - fRmax2*(fRmax2 + kRadTolerance*secRMax);
1558  }
1559  else if ( v.z() < 0.0 )
1560  {
1561  deltaRoi2 = snxt*snxt*t1 + 2*snxt*t2 + t3
1562  - fRmax1*(fRmax1 + kRadTolerance*secRMax);
1563  }
1564  else
1565  {
1566  deltaRoi2 = 1.0;
1567  }
1568 
1569  if ( nt1 && (deltaRoi2 > 0.0) )
1570  {
1571  // Equation quadratic => 2 roots : second root must be leaving
1572 
1573  b = nt2/nt1 ;
1574  c = nt3/nt1 ;
1575  d = b*b - c ;
1576 
1577  if ( d >= 0 )
1578  {
1579  // Check if on outer cone & heading outwards
1580  // NOTE: Should use rho-rout>-kRadTolerance*0.5
1581 
1582  if (nt3 > -halfRadTolerance && nt2 >= 0 )
1583  {
1584  if (calcNorm)
1585  {
1586  risec = std::sqrt(t3)*secRMax ;
1587  *validNorm = true ;
1588  *n = G4ThreeVector(p.x()/risec,p.y()/risec,-tanRMax/secRMax);
1589  }
1590  return snxt=0 ;
1591  }
1592  else
1593  {
1594  sider = kRMax ;
1595  if (b>0) { srd = -b - std::sqrt(d); }
1596  else { srd = c/(-b+std::sqrt(d)) ; }
1597 
1598  zi = p.z() + srd*v.z() ;
1599  ri = tanRMax*zi + rMaxAv ;
1600 
1601  if ((ri >= 0) && (-halfRadTolerance <= srd) && (srd <= halfRadTolerance))
1602  {
1603  // An intersection within the tolerance
1604  // we will Store it in case it is good -
1605  //
1606  slentol = srd ;
1607  sidetol = kRMax ;
1608  }
1609  if ( (ri < 0) || (srd < halfRadTolerance) )
1610  {
1611  // Safety: if both roots -ve ensure that srd cannot `win'
1612  // distance to out
1613 
1614  if (b>0) { sr2 = c/(-b-std::sqrt(d)); }
1615  else { sr2 = -b + std::sqrt(d); }
1616  zi = p.z() + sr2*v.z() ;
1617  ri = tanRMax*zi + rMaxAv ;
1618 
1619  if ((ri >= 0) && (sr2 > halfRadTolerance))
1620  {
1621  srd = sr2;
1622  }
1623  else
1624  {
1625  srd = kInfinity ;
1626 
1627  if( (-halfRadTolerance <= sr2) && ( sr2 <= halfRadTolerance) )
1628  {
1629  // An intersection within the tolerance.
1630  // Storing it in case it is good.
1631 
1632  slentol = sr2 ;
1633  sidetol = kRMax ;
1634  }
1635  }
1636  }
1637  }
1638  }
1639  else
1640  {
1641  // No intersection with outer cone & not parallel
1642  // -> already outside, no intersection
1643 
1644  if ( calcNorm )
1645  {
1646  risec = std::sqrt(t3)*secRMax;
1647  *validNorm = true;
1648  *n = G4ThreeVector(p.x()/risec,p.y()/risec,-tanRMax/secRMax);
1649  }
1650  return snxt = 0.0 ;
1651  }
1652  }
1653  else if ( nt2 && (deltaRoi2 > 0.0) )
1654  {
1655  // Linear case (only one intersection) => point outside outer cone
1656 
1657  if ( calcNorm )
1658  {
1659  risec = std::sqrt(t3)*secRMax;
1660  *validNorm = true;
1661  *n = G4ThreeVector(p.x()/risec,p.y()/risec,-tanRMax/secRMax);
1662  }
1663  return snxt = 0.0 ;
1664  }
1665  else
1666  {
1667  // No intersection -> parallel to outer cone
1668  // => Z or inner cone intersection
1669 
1670  srd = kInfinity ;
1671  }
1672 
1673  // Check possible intersection within tolerance
1674 
1675  if ( slentol <= halfCarTolerance )
1676  {
1677  // An intersection within the tolerance was found.
1678  // We must accept it only if the momentum points outwards.
1679  //
1680  // G4ThreeVector ptTol ; // The point of the intersection
1681  // ptTol= p + slentol*v ;
1682  // ri=tanRMax*zi+rMaxAv ;
1683  //
1684  // Calculate a normal vector, as below
1685 
1686  xi = p.x() + slentol*v.x();
1687  yi = p.y() + slentol*v.y();
1688  risec = std::sqrt(xi*xi + yi*yi)*secRMax;
1689  G4ThreeVector Normal = G4ThreeVector(xi/risec,yi/risec,-tanRMax/secRMax);
1690 
1691  if ( Normal.dot(v) > 0 ) // We will leave the Cone immediatelly
1692  {
1693  if ( calcNorm )
1694  {
1695  *n = Normal.unit() ;
1696  *validNorm = true ;
1697  }
1698  return snxt = 0.0 ;
1699  }
1700  else // On the surface, but not heading out so we ignore this intersection
1701  { // (as it is within tolerance).
1702  slentol = kInfinity ;
1703  }
1704  }
1705 
1706  // Inner Cone intersection
1707 
1708  if ( fRmin1 || fRmin2 )
1709  {
1710  tanRMin = (fRmin2 - fRmin1)*0.5/fDz ;
1711  nt1 = t1 - (tanRMin*v.z())*(tanRMin*v.z()) ;
1712 
1713  if ( nt1 )
1714  {
1715  secRMin = std::sqrt(1.0 + tanRMin*tanRMin) ;
1716  rMinAv = (fRmin1 + fRmin2)*0.5 ;
1717  rin = tanRMin*p.z() + rMinAv ;
1718  nt2 = t2 - tanRMin*v.z()*rin ;
1719  nt3 = t3 - rin*rin ;
1720 
1721  // Equation quadratic => 2 roots : first root must be leaving
1722 
1723  b = nt2/nt1 ;
1724  c = nt3/nt1 ;
1725  d = b*b - c ;
1726 
1727  if ( d >= 0.0 )
1728  {
1729  // NOTE: should be rho-rin<kRadTolerance*0.5,
1730  // but using squared versions for efficiency
1731 
1732  if (nt3 < kRadTolerance*(rin + kRadTolerance*0.25))
1733  {
1734  if ( nt2 < 0.0 )
1735  {
1736  if (calcNorm) { *validNorm = false; }
1737  return snxt = 0.0;
1738  }
1739  }
1740  else
1741  {
1742  if (b>0) { sr2 = -b - std::sqrt(d); }
1743  else { sr2 = c/(-b+std::sqrt(d)); }
1744  zi = p.z() + sr2*v.z() ;
1745  ri = tanRMin*zi + rMinAv ;
1746 
1747  if( (ri>=0.0)&&(-halfRadTolerance<=sr2)&&(sr2<=halfRadTolerance) )
1748  {
1749  // An intersection within the tolerance
1750  // storing it in case it is good.
1751 
1752  slentol = sr2 ;
1753  sidetol = kRMax ;
1754  }
1755  if( (ri<0) || (sr2 < halfRadTolerance) )
1756  {
1757  if (b>0) { sr3 = c/(-b-std::sqrt(d)); }
1758  else { sr3 = -b + std::sqrt(d) ; }
1759 
1760  // Safety: if both roots -ve ensure that srd cannot `win'
1761  // distancetoout
1762 
1763  if ( sr3 > halfRadTolerance )
1764  {
1765  if( sr3 < srd )
1766  {
1767  zi = p.z() + sr3*v.z() ;
1768  ri = tanRMin*zi + rMinAv ;
1769 
1770  if ( ri >= 0.0 )
1771  {
1772  srd=sr3 ;
1773  sider=kRMin ;
1774  }
1775  }
1776  }
1777  else if ( sr3 > -halfRadTolerance )
1778  {
1779  // Intersection in tolerance. Store to check if it's good
1780 
1781  slentol = sr3 ;
1782  sidetol = kRMin ;
1783  }
1784  }
1785  else if ( (sr2 < srd) && (sr2 > halfCarTolerance) )
1786  {
1787  srd = sr2 ;
1788  sider = kRMin ;
1789  }
1790  else if (sr2 > -halfCarTolerance)
1791  {
1792  // Intersection in tolerance. Store to check if it's good
1793 
1794  slentol = sr2 ;
1795  sidetol = kRMin ;
1796  }
1797  if( slentol <= halfCarTolerance )
1798  {
1799  // An intersection within the tolerance was found.
1800  // We must accept it only if the momentum points outwards.
1801 
1802  G4ThreeVector Normal ;
1803 
1804  // Calculate a normal vector, as below
1805 
1806  xi = p.x() + slentol*v.x() ;
1807  yi = p.y() + slentol*v.y() ;
1808  if( sidetol==kRMax )
1809  {
1810  risec = std::sqrt(xi*xi + yi*yi)*secRMax ;
1811  Normal = G4ThreeVector(xi/risec,yi/risec,-tanRMax/secRMax) ;
1812  }
1813  else
1814  {
1815  risec = std::sqrt(xi*xi + yi*yi)*secRMin ;
1816  Normal = G4ThreeVector(-xi/risec,-yi/risec,tanRMin/secRMin) ;
1817  }
1818  if( Normal.dot(v) > 0 )
1819  {
1820  // We will leave the cone immediately
1821 
1822  if( calcNorm )
1823  {
1824  *n = Normal.unit() ;
1825  *validNorm = true ;
1826  }
1827  return snxt = 0.0 ;
1828  }
1829  else
1830  {
1831  // On the surface, but not heading out so we ignore this
1832  // intersection (as it is within tolerance).
1833 
1834  slentol = kInfinity ;
1835  }
1836  }
1837  }
1838  }
1839  }
1840  }
1841 
1842  // Linear case => point outside inner cone ---> outer cone intersect
1843  //
1844  // Phi Intersection
1845 
1846  if ( !fPhiFullCone )
1847  {
1848  // add angle calculation with correction
1849  // of the difference in domain of atan2 and Sphi
1850 
1851  vphi = std::atan2(v.y(),v.x()) ;
1852 
1853  if ( vphi < fSPhi - halfAngTolerance ) { vphi += twopi; }
1854  else if ( vphi > fSPhi + fDPhi + halfAngTolerance ) { vphi -= twopi; }
1855 
1856  if ( p.x() || p.y() ) // Check if on z axis (rho not needed later)
1857  {
1858  // pDist -ve when inside
1859 
1860  pDistS = p.x()*sinSPhi - p.y()*cosSPhi ;
1861  pDistE = -p.x()*sinEPhi + p.y()*cosEPhi ;
1862 
1863  // Comp -ve when in direction of outwards normal
1864 
1865  compS = -sinSPhi*v.x() + cosSPhi*v.y() ;
1866  compE = sinEPhi*v.x() - cosEPhi*v.y() ;
1867 
1868  sidephi = kNull ;
1869 
1870  if( ( (fDPhi <= pi) && ( (pDistS <= halfCarTolerance)
1871  && (pDistE <= halfCarTolerance) ) )
1872  || ( (fDPhi > pi) && !((pDistS > halfCarTolerance)
1873  && (pDistE > halfCarTolerance) ) ) )
1874  {
1875  // Inside both phi *full* planes
1876  if ( compS < 0 )
1877  {
1878  sphi = pDistS/compS ;
1879  if (sphi >= -halfCarTolerance)
1880  {
1881  xi = p.x() + sphi*v.x() ;
1882  yi = p.y() + sphi*v.y() ;
1883 
1884  // Check intersecting with correct half-plane
1885  // (if not -> no intersect)
1886  //
1887  if ( (std::fabs(xi)<=kCarTolerance)
1888  && (std::fabs(yi)<=kCarTolerance) )
1889  {
1890  sidephi= kSPhi;
1891  if ( ( fSPhi-halfAngTolerance <= vphi )
1892  && ( fSPhi+fDPhi+halfAngTolerance >=vphi ) )
1893  {
1894  sphi = kInfinity;
1895  }
1896  }
1897  else
1898  if ( (yi*cosCPhi-xi*sinCPhi)>=0 )
1899  {
1900  sphi = kInfinity ;
1901  }
1902  else
1903  {
1904  sidephi = kSPhi ;
1905  if ( pDistS > -halfCarTolerance )
1906  {
1907  sphi = 0.0 ; // Leave by sphi immediately
1908  }
1909  }
1910  }
1911  else
1912  {
1913  sphi = kInfinity ;
1914  }
1915  }
1916  else
1917  {
1918  sphi = kInfinity ;
1919  }
1920 
1921  if ( compE < 0 )
1922  {
1923  sphi2 = pDistE/compE ;
1924 
1925  // Only check further if < starting phi intersection
1926  //
1927  if ( (sphi2 > -halfCarTolerance) && (sphi2 < sphi) )
1928  {
1929  xi = p.x() + sphi2*v.x() ;
1930  yi = p.y() + sphi2*v.y() ;
1931 
1932  // Check intersecting with correct half-plane
1933 
1934  if ( (std::fabs(xi)<=kCarTolerance)
1935  && (std::fabs(yi)<=kCarTolerance) )
1936  {
1937  // Leaving via ending phi
1938 
1939  if(!( (fSPhi-halfAngTolerance <= vphi)
1940  && (fSPhi+fDPhi+halfAngTolerance >= vphi) ) )
1941  {
1942  sidephi = kEPhi ;
1943  if ( pDistE <= -halfCarTolerance ) { sphi = sphi2; }
1944  else { sphi = 0.0; }
1945  }
1946  }
1947  else // Check intersecting with correct half-plane
1948  if ( yi*cosCPhi-xi*sinCPhi >= 0 )
1949  {
1950  // Leaving via ending phi
1951 
1952  sidephi = kEPhi ;
1953  if ( pDistE <= -halfCarTolerance ) { sphi = sphi2; }
1954  else { sphi = 0.0; }
1955  }
1956  }
1957  }
1958  }
1959  else
1960  {
1961  sphi = kInfinity ;
1962  }
1963  }
1964  else
1965  {
1966  // On z axis + travel not || to z axis -> if phi of vector direction
1967  // within phi of shape, Step limited by rmax, else Step =0
1968 
1969  if ( (fSPhi-halfAngTolerance <= vphi)
1970  && (vphi <= fSPhi+fDPhi+halfAngTolerance) )
1971  {
1972  sphi = kInfinity ;
1973  }
1974  else
1975  {
1976  sidephi = kSPhi ; // arbitrary
1977  sphi = 0.0 ;
1978  }
1979  }
1980  if ( sphi < snxt ) // Order intersecttions
1981  {
1982  snxt=sphi ;
1983  side=sidephi ;
1984  }
1985  }
1986  if ( srd < snxt ) // Order intersections
1987  {
1988  snxt = srd ;
1989  side = sider ;
1990  }
1991  if (calcNorm)
1992  {
1993  switch(side)
1994  { // Note: returned vector not normalised
1995  case kRMax: // (divide by frmax for unit vector)
1996  xi = p.x() + snxt*v.x() ;
1997  yi = p.y() + snxt*v.y() ;
1998  risec = std::sqrt(xi*xi + yi*yi)*secRMax ;
1999  *n = G4ThreeVector(xi/risec,yi/risec,-tanRMax/secRMax) ;
2000  *validNorm = true ;
2001  break ;
2002  case kRMin:
2003  *validNorm = false ; // Rmin is inconvex
2004  break ;
2005  case kSPhi:
2006  if ( fDPhi <= pi )
2007  {
2008  *n = G4ThreeVector(sinSPhi, -cosSPhi, 0);
2009  *validNorm = true ;
2010  }
2011  else
2012  {
2013  *validNorm = false ;
2014  }
2015  break ;
2016  case kEPhi:
2017  if ( fDPhi <= pi )
2018  {
2019  *n = G4ThreeVector(-sinEPhi, cosEPhi, 0);
2020  *validNorm = true ;
2021  }
2022  else
2023  {
2024  *validNorm = false ;
2025  }
2026  break ;
2027  case kPZ:
2028  *n = G4ThreeVector(0,0,1) ;
2029  *validNorm = true ;
2030  break ;
2031  case kMZ:
2032  *n = G4ThreeVector(0,0,-1) ;
2033  *validNorm = true ;
2034  break ;
2035  default:
2036  G4cout << G4endl ;
2037  DumpInfo();
2038  std::ostringstream message;
2039  G4int oldprc = message.precision(16) ;
2040  message << "Undefined side for valid surface normal to solid."
2041  << G4endl
2042  << "Position:" << G4endl << G4endl
2043  << "p.x() = " << p.x()/mm << " mm" << G4endl
2044  << "p.y() = " << p.y()/mm << " mm" << G4endl
2045  << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl
2046  << "pho at z = " << std::sqrt( p.x()*p.x()+p.y()*p.y() )/mm
2047  << " mm" << G4endl << G4endl ;
2048  if( p.x() != 0. || p.y() != 0.)
2049  {
2050  message << "point phi = " << std::atan2(p.y(),p.x())/degree
2051  << " degree" << G4endl << G4endl ;
2052  }
2053  message << "Direction:" << G4endl << G4endl
2054  << "v.x() = " << v.x() << G4endl
2055  << "v.y() = " << v.y() << G4endl
2056  << "v.z() = " << v.z() << G4endl<< G4endl
2057  << "Proposed distance :" << G4endl<< G4endl
2058  << "snxt = " << snxt/mm << " mm" << G4endl ;
2059  message.precision(oldprc) ;
2060  G4Exception("G4Cons::DistanceToOut(p,v,..)","GeomSolids1002",
2061  JustWarning, message) ;
2062  break ;
2063  }
2064  }
2065  if (snxt < halfCarTolerance) { snxt = 0.; }
2066 
2067  return snxt ;
2068 }
2069 
2071 //
2072 // Calculate distance (<=actual) to closest surface of shape from inside
2073 
2075 {
2076  G4double safe=0.0, rho, safeR1, safeR2, safeZ, safePhi;
2077  G4double tanRMin, secRMin, pRMin;
2078  G4double tanRMax, secRMax, pRMax;
2079 
2080 #ifdef G4CSGDEBUG
2081  if( Inside(p) == kOutside )
2082  {
2083  G4int oldprc=G4cout.precision(16) ;
2084  G4cout << G4endl ;
2085  DumpInfo();
2086  G4cout << "Position:" << G4endl << G4endl ;
2087  G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ;
2088  G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ;
2089  G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ;
2090  G4cout << "pho at z = " << std::sqrt( p.x()*p.x()+p.y()*p.y() )/mm
2091  << " mm" << G4endl << G4endl ;
2092  if( (p.x() != 0.) || (p.x() != 0.) )
2093  {
2094  G4cout << "point phi = " << std::atan2(p.y(),p.x())/degree
2095  << " degree" << G4endl << G4endl ;
2096  }
2097  G4cout.precision(oldprc) ;
2098  G4Exception("G4Cons::DistanceToOut(p)", "GeomSolids1002",
2099  JustWarning, "Point p is outside !?" );
2100  }
2101 #endif
2102 
2103  rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
2104  safeZ = fDz - std::fabs(p.z()) ;
2105 
2106  if (fRmin1 || fRmin2)
2107  {
2108  tanRMin = (fRmin2 - fRmin1)*0.5/fDz ;
2109  secRMin = std::sqrt(1.0 + tanRMin*tanRMin) ;
2110  pRMin = tanRMin*p.z() + (fRmin1 + fRmin2)*0.5 ;
2111  safeR1 = (rho - pRMin)/secRMin ;
2112  }
2113  else
2114  {
2115  safeR1 = kInfinity ;
2116  }
2117 
2118  tanRMax = (fRmax2 - fRmax1)*0.5/fDz ;
2119  secRMax = std::sqrt(1.0 + tanRMax*tanRMax) ;
2120  pRMax = tanRMax*p.z() + (fRmax1+fRmax2)*0.5 ;
2121  safeR2 = (pRMax - rho)/secRMax ;
2122 
2123  if (safeR1 < safeR2) { safe = safeR1; }
2124  else { safe = safeR2; }
2125  if (safeZ < safe) { safe = safeZ ; }
2126 
2127  // Check if phi divided, Calc distances closest phi plane
2128 
2129  if (!fPhiFullCone)
2130  {
2131  // Above/below central phi of G4Cons?
2132 
2133  if ( (p.y()*cosCPhi - p.x()*sinCPhi) <= 0 )
2134  {
2135  safePhi = -(p.x()*sinSPhi - p.y()*cosSPhi) ;
2136  }
2137  else
2138  {
2139  safePhi = (p.x()*sinEPhi - p.y()*cosEPhi) ;
2140  }
2141  if (safePhi < safe) { safe = safePhi; }
2142  }
2143  if ( safe < 0 ) { safe = 0; }
2144 
2145  return safe ;
2146 }
2147 
2149 //
2150 // Create a List containing the transformed vertices
2151 // Ordering [0-3] -fDz cross section
2152 // [4-7] +fDz cross section such that [0] is below [4],
2153 // [1] below [5] etc.
2154 // Note:
2155 // Caller has deletion resposibility
2156 // Potential improvement: For last slice, use actual ending angle
2157 // to avoid rounding error problems.
2158 
2161 {
2162  G4ThreeVectorList* vertices ;
2163  G4ThreeVector vertex0, vertex1, vertex2, vertex3 ;
2164  G4double meshAngle, meshRMax1, meshRMax2, crossAngle;
2165  G4double cosCrossAngle, sinCrossAngle, sAngle ;
2166  G4double rMaxX1, rMaxX2, rMaxY1, rMaxY2, rMinX1, rMinX2, rMinY1, rMinY2 ;
2167  G4int crossSection, noCrossSections ;
2168 
2169  // Compute no of cross-sections necessary to mesh cone
2170 
2171  noCrossSections = G4int(fDPhi/kMeshAngleDefault) + 1 ;
2172 
2173  if (noCrossSections < kMinMeshSections)
2174  {
2175  noCrossSections = kMinMeshSections ;
2176  }
2177  else if (noCrossSections > kMaxMeshSections)
2178  {
2179  noCrossSections = kMaxMeshSections ;
2180  }
2181  meshAngle = fDPhi/(noCrossSections - 1) ;
2182 
2183  meshRMax1 = fRmax1/std::cos(meshAngle*0.5) ;
2184  meshRMax2 = fRmax2/std::cos(meshAngle*0.5) ;
2185 
2186  // If complete in phi, set start angle such that mesh will be at RMax
2187  // on the x axis. Will give better extent calculations when not rotated.
2188 
2189  if ( fPhiFullCone && (fSPhi == 0.0) )
2190  {
2191  sAngle = -meshAngle*0.5 ;
2192  }
2193  else
2194  {
2195  sAngle = fSPhi ;
2196  }
2197  vertices = new G4ThreeVectorList();
2198 
2199  if (vertices)
2200  {
2201  vertices->reserve(noCrossSections*4) ;
2202  for (crossSection = 0 ; crossSection < noCrossSections ; crossSection++)
2203  {
2204  // Compute coordinates of cross section at section crossSection
2205 
2206  crossAngle = sAngle + crossSection*meshAngle ;
2207  cosCrossAngle = std::cos(crossAngle) ;
2208  sinCrossAngle = std::sin(crossAngle) ;
2209 
2210  rMaxX1 = meshRMax1*cosCrossAngle ;
2211  rMaxY1 = meshRMax1*sinCrossAngle ;
2212  rMaxX2 = meshRMax2*cosCrossAngle ;
2213  rMaxY2 = meshRMax2*sinCrossAngle ;
2214 
2215  rMinX1 = fRmin1*cosCrossAngle ;
2216  rMinY1 = fRmin1*sinCrossAngle ;
2217  rMinX2 = fRmin2*cosCrossAngle ;
2218  rMinY2 = fRmin2*sinCrossAngle ;
2219 
2220  vertex0 = G4ThreeVector(rMinX1,rMinY1,-fDz) ;
2221  vertex1 = G4ThreeVector(rMaxX1,rMaxY1,-fDz) ;
2222  vertex2 = G4ThreeVector(rMaxX2,rMaxY2,+fDz) ;
2223  vertex3 = G4ThreeVector(rMinX2,rMinY2,+fDz) ;
2224 
2225  vertices->push_back(pTransform.TransformPoint(vertex0)) ;
2226  vertices->push_back(pTransform.TransformPoint(vertex1)) ;
2227  vertices->push_back(pTransform.TransformPoint(vertex2)) ;
2228  vertices->push_back(pTransform.TransformPoint(vertex3)) ;
2229  }
2230  }
2231  else
2232  {
2233  DumpInfo();
2234  G4Exception("G4Cons::CreateRotatedVertices()",
2235  "GeomSolids0003", FatalException,
2236  "Error in allocation of vertices. Out of memory !");
2237  }
2238 
2239  return vertices ;
2240 }
2241 
2243 //
2244 // GetEntityType
2245 
2247 {
2248  return G4String("G4Cons");
2249 }
2250 
2252 //
2253 // Make a clone of the object
2254 //
2256 {
2257  return new G4Cons(*this);
2258 }
2259 
2261 //
2262 // Stream object contents to an output stream
2263 
2264 std::ostream& G4Cons::StreamInfo(std::ostream& os) const
2265 {
2266  G4int oldprc = os.precision(16);
2267  os << "-----------------------------------------------------------\n"
2268  << " *** Dump for solid - " << GetName() << " ***\n"
2269  << " ===================================================\n"
2270  << " Solid type: G4Cons\n"
2271  << " Parameters: \n"
2272  << " inside -fDz radius: " << fRmin1/mm << " mm \n"
2273  << " outside -fDz radius: " << fRmax1/mm << " mm \n"
2274  << " inside +fDz radius: " << fRmin2/mm << " mm \n"
2275  << " outside +fDz radius: " << fRmax2/mm << " mm \n"
2276  << " half length in Z : " << fDz/mm << " mm \n"
2277  << " starting angle of segment: " << fSPhi/degree << " degrees \n"
2278  << " delta angle of segment : " << fDPhi/degree << " degrees \n"
2279  << "-----------------------------------------------------------\n";
2280  os.precision(oldprc);
2281 
2282  return os;
2283 }
2284 
2285 
2286 
2288 //
2289 // GetPointOnSurface
2290 
2292 {
2293  // declare working variables
2294  //
2295  G4double Aone, Atwo, Athree, Afour, Afive, slin, slout, phi;
2296  G4double zRand, cosu, sinu, rRand1, rRand2, chose, rone, rtwo, qone, qtwo;
2297  rone = (fRmax1-fRmax2)/(2.*fDz);
2298  rtwo = (fRmin1-fRmin2)/(2.*fDz);
2299  qone=0.; qtwo=0.;
2300  if(fRmax1!=fRmax2) { qone = fDz*(fRmax1+fRmax2)/(fRmax1-fRmax2); }
2301  if(fRmin1!=fRmin2) { qtwo = fDz*(fRmin1+fRmin2)/(fRmin1-fRmin2); }
2302  slin = std::sqrt(sqr(fRmin1-fRmin2)+sqr(2.*fDz));
2303  slout = std::sqrt(sqr(fRmax1-fRmax2)+sqr(2.*fDz));
2304  Aone = 0.5*fDPhi*(fRmax2 + fRmax1)*slout;
2305  Atwo = 0.5*fDPhi*(fRmin2 + fRmin1)*slin;
2306  Athree = 0.5*fDPhi*(fRmax1*fRmax1-fRmin1*fRmin1);
2307  Afour = 0.5*fDPhi*(fRmax2*fRmax2-fRmin2*fRmin2);
2308  Afive = fDz*(fRmax1-fRmin1+fRmax2-fRmin2);
2309 
2311  cosu = std::cos(phi); sinu = std::sin(phi);
2312  rRand1 = GetRadiusInRing(fRmin1, fRmin2);
2313  rRand2 = GetRadiusInRing(fRmax1, fRmax2);
2314 
2315  if ( (fSPhi == 0.) && fPhiFullCone ) { Afive = 0.; }
2316  chose = RandFlat::shoot(0.,Aone+Atwo+Athree+Afour+2.*Afive);
2317 
2318  if( (chose >= 0.) && (chose < Aone) )
2319  {
2320  if(fRmin1 != fRmin2)
2321  {
2322  zRand = RandFlat::shoot(-1.*fDz,fDz);
2323  return G4ThreeVector (rtwo*cosu*(qtwo-zRand),
2324  rtwo*sinu*(qtwo-zRand), zRand);
2325  }
2326  else
2327  {
2328  return G4ThreeVector(fRmin1*cosu, fRmin2*sinu,
2329  RandFlat::shoot(-1.*fDz,fDz));
2330  }
2331  }
2332  else if( (chose >= Aone) && (chose <= Aone + Atwo) )
2333  {
2334  if(fRmax1 != fRmax2)
2335  {
2336  zRand = RandFlat::shoot(-1.*fDz,fDz);
2337  return G4ThreeVector (rone*cosu*(qone-zRand),
2338  rone*sinu*(qone-zRand), zRand);
2339  }
2340  else
2341  {
2342  return G4ThreeVector(fRmax1*cosu, fRmax2*sinu,
2343  RandFlat::shoot(-1.*fDz,fDz));
2344  }
2345  }
2346  else if( (chose >= Aone + Atwo) && (chose < Aone + Atwo + Athree) )
2347  {
2348  return G4ThreeVector (rRand1*cosu, rRand1*sinu, -1*fDz);
2349  }
2350  else if( (chose >= Aone + Atwo + Athree)
2351  && (chose < Aone + Atwo + Athree + Afour) )
2352  {
2353  return G4ThreeVector (rRand2*cosu,rRand2*sinu,fDz);
2354  }
2355  else if( (chose >= Aone + Atwo + Athree + Afour)
2356  && (chose < Aone + Atwo + Athree + Afour + Afive) )
2357  {
2358  zRand = RandFlat::shoot(-1.*fDz,fDz);
2359  rRand1 = RandFlat::shoot(fRmin2-((zRand-fDz)/(2.*fDz))*(fRmin1-fRmin2),
2360  fRmax2-((zRand-fDz)/(2.*fDz))*(fRmax1-fRmax2));
2361  return G4ThreeVector (rRand1*std::cos(fSPhi),
2362  rRand1*std::sin(fSPhi), zRand);
2363  }
2364  else
2365  {
2366  zRand = RandFlat::shoot(-1.*fDz,fDz);
2367  rRand1 = RandFlat::shoot(fRmin2-((zRand-fDz)/(2.*fDz))*(fRmin1-fRmin2),
2368  fRmax2-((zRand-fDz)/(2.*fDz))*(fRmax1-fRmax2));
2369  return G4ThreeVector (rRand1*std::cos(fSPhi+fDPhi),
2370  rRand1*std::sin(fSPhi+fDPhi), zRand);
2371  }
2372 }
2373 
2375 //
2376 // Methods for visualisation
2377 
2379 {
2380  scene.AddSolid (*this);
2381 }
2382 
2384 {
2386 }
2387 
2388 #endif
G4Cons(const G4String &pName, G4double pRmin1, G4double pRmax1, G4double pRmin2, G4double pRmax2, G4double pDz, G4double pSPhi, G4double pDPhi)
Definition: G4Cons.cc:78
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)
Definition: G4Cons.cc:67
G4double halfRadTolerance
Definition: G4Cons.hh:233
G4double fDz
Definition: G4Cons.hh:220
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double GetMinYExtent() const
G4bool fPhiFullCone
Definition: G4Cons.hh:229
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Cons.cc:2160
Definition: G4Cons.cc:67
G4AffineTransform Inverse() const
G4Cons & operator=(const G4Cons &rhs)
Definition: G4Cons.cc:170
G4double sinEPhi
Definition: G4Cons.hh:224
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Cons.cc:719
G4bool IsYLimited() const
G4double GetRadiusInRing(G4double rmin, G4double rmax) const
Definition: G4CSGSolid.cc:101
EInside Inside(const G4ThreeVector &p) const
Definition: G4Cons.cc:204
const G4double pi
G4double fRmin1
Definition: G4Cons.hh:220
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Cons.cc:1451
G4bool IsRotated() const
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Cons.cc:593
G4ThreeVector NetTranslation() const
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Cons.cc:2264
G4bool IsXLimited() const
G4double a
Definition: TRTMaterials.hh:39
const G4int kMinMeshSections
Definition: meshdefs.hh:45
virtual void AddSolid(const G4Box &)=0
int G4int
Definition: G4Types.hh:78
G4GeometryType GetEntityType() const
Definition: G4Cons.cc:2246
G4double fRmax2
Definition: G4Cons.hh:220
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Cons.cc:259
Definition: G4Cons.cc:67
void DumpInfo() const
ENorm
Definition: G4Cons.cc:71
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4double cosCPhi
Definition: G4Cons.hh:224
G4double kRadTolerance
Definition: G4Cons.hh:216
G4double GetMaxXExtent() const
G4double GetMinZExtent() const
G4GLOB_DLL std::ostream G4cout
G4double fRmin2
Definition: G4Cons.hh:220
Definition: G4Cons.cc:67
G4double sinSPhi
Definition: G4Cons.hh:224
G4Polyhedron * fpPolyhedron
Definition: G4CSGSolid.hh:80
G4double fSPhi
Definition: G4Cons.hh:220
G4double GetRadialTolerance() const
bool G4bool
Definition: G4Types.hh:79
G4double halfCarTolerance
Definition: G4Cons.hh:233
Definition: G4Cons.hh:83
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4CSGSolid.cc:124
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
const G4int n
G4double cosSPhi
Definition: G4Cons.hh:224
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4double GetMinXExtent() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double GetMaxZExtent() const
G4double fDPhi
Definition: G4Cons.hh:220
~G4Cons()
Definition: G4Cons.cc:143
G4double halfAngTolerance
Definition: G4Cons.hh:233
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
G4Polyhedron * CreatePolyhedron() const
Definition: G4Cons.cc:2383
Definition: G4Cons.cc:67
G4double sinCPhi
Definition: G4Cons.hh:224
static const double degree
Definition: G4SIunits.hh:125
G4double kAngTolerance
Definition: G4Cons.hh:216
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Cons.cc:2378
G4double cosEPhi
Definition: G4Cons.hh:224
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Cons.cc:271
#define G4endl
Definition: G4ios.hh:61
G4double GetMaxYExtent() const
G4double kCarTolerance
Definition: G4VSolid.hh:305
Definition: G4Cons.cc:71
ESide
Definition: G4Cons.cc:67
T sqr(const T &x)
Definition: templates.hh:145
G4double cosHDPhiOT
Definition: G4Cons.hh:224
G4ThreeVector GetPointOnSurface() const
Definition: G4Cons.cc:2291
Definition: G4Cons.cc:67
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Cons.cc:489
G4double cosHDPhiIT
Definition: G4Cons.hh:224
double G4double
Definition: G4Types.hh:76
G4VSolid * Clone() const
Definition: G4Cons.cc:2255
G4double GetMaxExtent(const EAxis pAxis) const
static const G4double e3
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition: G4CSGSolid.cc:82
G4bool IsZLimited() const
static const double mm
Definition: G4SIunits.hh:102
Definition: G4Cons.cc:71
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
Definition: G4Cons.cc:67
const G4double kMeshAngleDefault
Definition: meshdefs.hh:42
static G4GeometryTolerance * GetInstance()
G4double fRmax1
Definition: G4Cons.hh:220
Definition: G4Cons.cc:71
const G4int kMaxMeshSections
Definition: meshdefs.hh:46
Definition: G4Cons.cc:71
Definition: G4Cons.cc:71