Geant4  10.02.p02
G4GenericPolycone.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: G4GenericPolycone.cc 72131 2013-07-11 07:12:24Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class source file
32 //
33 //
34 // G4GenericPolycone.cc
35 //
36 // Implementation of a CSG polycone
37 //
38 // --------------------------------------------------------------------
39 
40 #include "G4GenericPolycone.hh"
41 
42 #if !defined(G4GEOM_USE_UGENERICPOLYCONE)
43 
44 #include "G4PolyconeSide.hh"
45 #include "G4PolyPhiFace.hh"
46 
47 #include "Randomize.hh"
48 
49 #include "G4Polyhedron.hh"
50 #include "G4EnclosingCylinder.hh"
51 #include "G4ReduciblePolygon.hh"
52 #include "G4VPVParameterisation.hh"
53 
54 using namespace CLHEP;
55 
56 //
57 // Constructor (generic parameters)
58 //
60  G4double phiStart,
61  G4double phiTotal,
62  G4int numRZ,
63  const G4double r[],
64  const G4double z[] )
65  : G4VCSGfaceted( name )
66 {
67 
68  G4ReduciblePolygon *rz = new G4ReduciblePolygon( r, z, numRZ );
69 
70  Create( phiStart, phiTotal, rz );
71 
72  // Set original_parameters struct for consistency
73  //
74  //SetOriginalParameters(rz);
75 
76  delete rz;
77 }
78 
79 //
80 // Create
81 //
82 // Generic create routine, called by each constructor after
83 // conversion of arguments
84 //
86  G4double phiTotal,
87  G4ReduciblePolygon *rz )
88 {
89  //
90  // Perform checks of rz values
91  //
92  if (rz->Amin() < 0.0)
93  {
94  std::ostringstream message;
95  message << "Illegal input parameters - " << GetName() << G4endl
96  << " All R values must be >= 0 !";
97  G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
98  FatalErrorInArgument, message);
99  }
100 
101  G4double rzArea = rz->Area();
102  if (rzArea < -kCarTolerance)
103  rz->ReverseOrder();
104 
105  else if (rzArea < -kCarTolerance)
106  {
107  std::ostringstream message;
108  message << "Illegal input parameters - " << GetName() << G4endl
109  << " R/Z cross section is zero or near zero: " << rzArea;
110  G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
111  FatalErrorInArgument, message);
112  }
113 
115  || (!rz->RemoveRedundantVertices( kCarTolerance )) )
116  {
117  std::ostringstream message;
118  message << "Illegal input parameters - " << GetName() << G4endl
119  << " Too few unique R/Z values !";
120  G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
121  FatalErrorInArgument, message);
122  }
123 
124  if (rz->CrossesItself(1/kInfinity))
125  {
126  std::ostringstream message;
127  message << "Illegal input parameters - " << GetName() << G4endl
128  << " R/Z segments cross !";
129  G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
130  FatalErrorInArgument, message);
131  }
132 
133  numCorner = rz->NumVertices();
134 
135  //
136  // Phi opening? Account for some possible roundoff, and interpret
137  // nonsense value as representing no phi opening
138  //
139  if (phiTotal <= 0 || phiTotal > twopi-1E-10)
140  {
141  phiIsOpen = false;
142  startPhi = 0;
143  endPhi = twopi;
144  }
145  else
146  {
147  phiIsOpen = true;
148 
149  //
150  // Convert phi into our convention
151  //
152  startPhi = phiStart;
153  while( startPhi < 0 ) // Loop checking, 13.08.2015, G.Cosmo
154  startPhi += twopi;
155 
156  endPhi = phiStart+phiTotal;
157  while( endPhi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
158  endPhi += twopi;
159  }
160 
161  //
162  // Allocate corner array.
163  //
165 
166  //
167  // Copy corners
168  //
169  G4ReduciblePolygonIterator iterRZ(rz);
170 
171  G4PolyconeSideRZ *next = corners;
172  iterRZ.Begin();
173  do // Loop checking, 13.08.2015, G.Cosmo
174  {
175  next->r = iterRZ.GetA();
176  next->z = iterRZ.GetB();
177  } while( ++next, iterRZ.Next() );
178 
179  //
180  // Allocate face pointer array
181  //
183  faces = new G4VCSGface*[numFace];
184 
185  //
186  // Construct conical faces
187  //
188  // But! Don't construct a face if both points are at zero radius!
189  //
190  G4PolyconeSideRZ *corner = corners,
191  *prev = corners + numCorner-1,
192  *nextNext;
193  G4VCSGface **face = faces;
194  do // Loop checking, 13.08.2015, G.Cosmo
195  {
196  next = corner+1;
197  if (next >= corners+numCorner) next = corners;
198  nextNext = next+1;
199  if (nextNext >= corners+numCorner) nextNext = corners;
200 
201  if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
202 
203  //
204  // We must decide here if we can dare declare one of our faces
205  // as having a "valid" normal (i.e. allBehind = true). This
206  // is never possible if the face faces "inward" in r.
207  //
208  G4bool allBehind;
209  if (corner->z > next->z)
210  {
211  allBehind = false;
212  }
213  else
214  {
215  //
216  // Otherwise, it is only true if the line passing
217  // through the two points of the segment do not
218  // split the r/z cross section
219  //
220  allBehind = !rz->BisectedBy( corner->r, corner->z,
221  next->r, next->z, kCarTolerance );
222  }
223 
224  *face++ = new G4PolyconeSide( prev, corner, next, nextNext,
225  startPhi, endPhi-startPhi, phiIsOpen, allBehind );
226  } while( prev=corner, corner=next, corner > corners );
227 
228  if (phiIsOpen)
229  {
230  //
231  // Construct phi open edges
232  //
233  *face++ = new G4PolyPhiFace( rz, startPhi, 0, endPhi );
234  *face++ = new G4PolyPhiFace( rz, endPhi, 0, startPhi );
235  }
236 
237  //
238  // We might have dropped a face or two: recalculate numFace
239  //
240  numFace = face-faces;
241 
242  //
243  // Make enclosingCylinder
244  //
246  new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
247 }
248 
249 
250 //
251 // Fake default constructor - sets only member data and allocates memory
252 // for usage restricted to object persistency.
253 //
255  : G4VCSGfaceted(a), startPhi(0.), endPhi(0.), phiIsOpen(false),
256  numCorner(0), corners(0), enclosingCylinder(0)
257 {
258 }
259 
260 
261 //
262 // Destructor
263 //
265 {
266  delete [] corners;
267  delete enclosingCylinder;
268 }
269 
270 
271 //
272 // Copy constructor
273 //
275  : G4VCSGfaceted( source )
276 {
277  CopyStuff( source );
278 }
279 
280 
281 //
282 // Assignment operator
283 //
285 {
286  if (this == &source) return *this;
287 
288  G4VCSGfaceted::operator=( source );
289 
290  delete [] corners;
291  // if (original_parameters) delete original_parameters;
292 
293  delete enclosingCylinder;
294 
295  CopyStuff( source );
296 
297  return *this;
298 }
299 
300 
301 //
302 // CopyStuff
303 //
305 {
306  //
307  // Simple stuff
308  //
309  startPhi = source.startPhi;
310  endPhi = source.endPhi;
311  phiIsOpen = source.phiIsOpen;
312  numCorner = source.numCorner;
313 
314  //
315  // The corner array
316  //
318 
319  G4PolyconeSideRZ *corn = corners,
320  *sourceCorn = source.corners;
321  do // Loop checking, 13.08.2015, G.Cosmo
322  {
323  *corn = *sourceCorn;
324  } while( ++sourceCorn, ++corn < corners+numCorner );
325 
326  //
327  // Enclosing cylinder
328  //
330 
331  fRebuildPolyhedron = false;
332  fpPolyhedron = 0;
333 }
334 
335 
336 //
337 // Reset
338 //
340 {
341 
342  std::ostringstream message;
343  message << "Solid " << GetName() << " built using generic construct."
344  << G4endl << "Not applicable to the generic construct !";
345  G4Exception("G4GenericPolycone::Reset()", "GeomSolids1001",
346  JustWarning, message, "Parameters NOT resetted.");
347  return 1;
348 
349 }
350 
351 
352 //
353 // Inside
354 //
355 // This is an override of G4VCSGfaceted::Inside, created in order
356 // to speed things up by first checking with G4EnclosingCylinder.
357 //
359 {
360  //
361  // Quick test
362  //
363  if (enclosingCylinder->MustBeOutside(p)) return kOutside;
364 
365  //
366  // Long answer
367  //
368  return G4VCSGfaceted::Inside(p);
369 }
370 
371 
372 //
373 // DistanceToIn
374 //
375 // This is an override of G4VCSGfaceted::Inside, created in order
376 // to speed things up by first checking with G4EnclosingCylinder.
377 //
379  const G4ThreeVector &v ) const
380 {
381  //
382  // Quick test
383  //
384  if (enclosingCylinder->ShouldMiss(p,v))
385  return kInfinity;
386 
387  //
388  // Long answer
389  //
390  return G4VCSGfaceted::DistanceToIn( p, v );
391 }
392 
393 
394 //
395 // DistanceToIn
396 //
398 {
399  return G4VCSGfaceted::DistanceToIn(p);
400 }
401 
402 
403 //
404 // ComputeDimensions
405 //
406 /*void G4GenericPolycone::ComputeDimensions( G4VPVParameterisation* p,
407  const G4int n,
408  const G4VPhysicalVolume* pRep )
409 {
410  p->ComputeDimensions(*this,n,pRep);
411 }
412 */
413 //
414 // GetEntityType
415 //
417 {
418  return G4String("G4GenericPolycone");
419 }
420 
421 
422 //
423 // Make a clone of the object
424 //
426 {
427  return new G4GenericPolycone(*this);
428 }
429 
430 //
431 // Stream object contents to an output stream
432 //
433 std::ostream& G4GenericPolycone::StreamInfo( std::ostream& os ) const
434 {
435  G4int oldprc = os.precision(16);
436  os << "-----------------------------------------------------------\n"
437  << " *** Dump for solid - " << GetName() << " ***\n"
438  << " ===================================================\n"
439  << " Solid type: G4GenericPolycone\n"
440  << " Parameters: \n"
441  << " starting phi angle : " << startPhi/degree << " degrees \n"
442  << " ending phi angle : " << endPhi/degree << " degrees \n";
443  G4int i=0;
444 
445  os << " number of RZ points: " << numCorner << "\n"
446  << " RZ values (corners): \n";
447  for (i=0; i<numCorner; i++)
448  {
449  os << " "
450  << corners[i].r << ", " << corners[i].z << "\n";
451  }
452  os << "-----------------------------------------------------------\n";
453  os.precision(oldprc);
454 
455  return os;
456 }
457 
458 
459 
460 //
461 // GetPointOnSurface
462 //
464 {
465  return GetPointOnSurfaceGeneric();
466 
467 }
468 
469 //
470 // CreatePolyhedron
471 //
473 {
474  // The following code prepares for:
475  // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
476  // const double xyz[][3],
477  // const int faces_vec[][4])
478  // Here is an extract from the header file HepPolyhedron.h:
496  const G4int numSide =
497  G4int(G4Polyhedron::GetNumberOfRotationSteps()
498  * (endPhi - startPhi) / twopi) + 1;
499  G4int nNodes;
500  G4int nFaces;
501  typedef G4double double3[3];
502  double3* xyz;
503  typedef G4int int4[4];
504  int4* faces_vec;
505  if (phiIsOpen)
506  {
507  // Triangulate open ends. Simple ear-chopping algorithm...
508  // I'm not sure how robust this algorithm is (J.Allison).
509  //
510  std::vector<G4bool> chopped(numCorner, false);
511  std::vector<G4int*> triQuads;
512  G4int remaining = numCorner;
513  G4int iStarter = 0;
514  while (remaining >= 3) // Loop checking, 13.08.2015, G.Cosmo
515  {
516  // Find unchopped corners...
517  //
518  G4int A = -1, B = -1, C = -1;
519  G4int iStepper = iStarter;
520  do // Loop checking, 13.08.2015, G.Cosmo
521  {
522  if (A < 0) { A = iStepper; }
523  else if (B < 0) { B = iStepper; }
524  else if (C < 0) { C = iStepper; }
525  do // Loop checking, 13.08.2015, G.Cosmo
526  {
527  if (++iStepper >= numCorner) { iStepper = 0; }
528  }
529  while (chopped[iStepper]);
530  }
531  while (C < 0 && iStepper != iStarter);
532 
533  // Check triangle at B is pointing outward (an "ear").
534  // Sign of z cross product determines...
535  //
536  G4double BAr = corners[A].r - corners[B].r;
537  G4double BAz = corners[A].z - corners[B].z;
538  G4double BCr = corners[C].r - corners[B].r;
539  G4double BCz = corners[C].z - corners[B].z;
540  if (BAr * BCz - BAz * BCr < kCarTolerance)
541  {
542  G4int* tq = new G4int[3];
543  tq[0] = A + 1;
544  tq[1] = B + 1;
545  tq[2] = C + 1;
546  triQuads.push_back(tq);
547  chopped[B] = true;
548  --remaining;
549  }
550  else
551  {
552  do // Loop checking, 13.08.2015, G.Cosmo
553  {
554  if (++iStarter >= numCorner) { iStarter = 0; }
555  }
556  while (chopped[iStarter]);
557  }
558  }
559  // Transfer to faces...
560  //
561  nNodes = (numSide + 1) * numCorner;
562  nFaces = numSide * numCorner + 2 * triQuads.size();
563  faces_vec = new int4[nFaces];
564  G4int iface = 0;
565  G4int addition = numCorner * numSide;
566  G4int d = numCorner - 1;
567  for (G4int iEnd = 0; iEnd < 2; ++iEnd)
568  {
569  for (size_t i = 0; i < triQuads.size(); ++i)
570  {
571  // Negative for soft/auxiliary/normally invisible edges...
572  //
573  G4int a, b, c;
574  if (iEnd == 0)
575  {
576  a = triQuads[i][0];
577  b = triQuads[i][1];
578  c = triQuads[i][2];
579  }
580  else
581  {
582  a = triQuads[i][0] + addition;
583  b = triQuads[i][2] + addition;
584  c = triQuads[i][1] + addition;
585  }
586  G4int ab = std::abs(b - a);
587  G4int bc = std::abs(c - b);
588  G4int ca = std::abs(a - c);
589  faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
590  faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
591  faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
592  faces_vec[iface][3] = 0;
593  ++iface;
594  }
595  }
596 
597  // Continue with sides...
598 
599  xyz = new double3[nNodes];
600  const G4double dPhi = (endPhi - startPhi) / numSide;
601  G4double phi = startPhi;
602  G4int ixyz = 0;
603  for (G4int iSide = 0; iSide < numSide; ++iSide)
604  {
605  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
606  {
607  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
608  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
609  xyz[ixyz][2] = corners[iCorner].z;
610  if (iSide == 0) // startPhi
611  {
612  if (iCorner < numCorner - 1)
613  {
614  faces_vec[iface][0] = ixyz + 1;
615  faces_vec[iface][1] = -(ixyz + numCorner + 1);
616  faces_vec[iface][2] = ixyz + numCorner + 2;
617  faces_vec[iface][3] = ixyz + 2;
618  }
619  else
620  {
621  faces_vec[iface][0] = ixyz + 1;
622  faces_vec[iface][1] = -(ixyz + numCorner + 1);
623  faces_vec[iface][2] = ixyz + 2;
624  faces_vec[iface][3] = ixyz - numCorner + 2;
625  }
626  }
627  else if (iSide == numSide - 1) // endPhi
628  {
629  if (iCorner < numCorner - 1)
630  {
631  faces_vec[iface][0] = ixyz + 1;
632  faces_vec[iface][1] = ixyz + numCorner + 1;
633  faces_vec[iface][2] = ixyz + numCorner + 2;
634  faces_vec[iface][3] = -(ixyz + 2);
635  }
636  else
637  {
638  faces_vec[iface][0] = ixyz + 1;
639  faces_vec[iface][1] = ixyz + numCorner + 1;
640  faces_vec[iface][2] = ixyz + 2;
641  faces_vec[iface][3] = -(ixyz - numCorner + 2);
642  }
643  }
644  else
645  {
646  if (iCorner < numCorner - 1)
647  {
648  faces_vec[iface][0] = ixyz + 1;
649  faces_vec[iface][1] = -(ixyz + numCorner + 1);
650  faces_vec[iface][2] = ixyz + numCorner + 2;
651  faces_vec[iface][3] = -(ixyz + 2);
652  }
653  else
654  {
655  faces_vec[iface][0] = ixyz + 1;
656  faces_vec[iface][1] = -(ixyz + numCorner + 1);
657  faces_vec[iface][2] = ixyz + 2;
658  faces_vec[iface][3] = -(ixyz - numCorner + 2);
659  }
660  }
661  ++iface;
662  ++ixyz;
663  }
664  phi += dPhi;
665  }
666 
667  // Last corners...
668 
669  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
670  {
671  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
672  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
673  xyz[ixyz][2] = corners[iCorner].z;
674  ++ixyz;
675  }
676  }
677  else // !phiIsOpen - i.e., a complete 360 degrees.
678  {
679  nNodes = numSide * numCorner;
680  nFaces = numSide * numCorner;;
681  xyz = new double3[nNodes];
682  faces_vec = new int4[nFaces];
683  const G4double dPhi = (endPhi - startPhi) / numSide;
684  G4double phi = startPhi;
685  G4int ixyz = 0, iface = 0;
686  for (G4int iSide = 0; iSide < numSide; ++iSide)
687  {
688  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
689  {
690  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
691  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
692  xyz[ixyz][2] = corners[iCorner].z;
693 
694  if (iSide < numSide - 1)
695  {
696  if (iCorner < numCorner - 1)
697  {
698  faces_vec[iface][0] = ixyz + 1;
699  faces_vec[iface][1] = -(ixyz + numCorner + 1);
700  faces_vec[iface][2] = ixyz + numCorner + 2;
701  faces_vec[iface][3] = -(ixyz + 2);
702  }
703  else
704  {
705  faces_vec[iface][0] = ixyz + 1;
706  faces_vec[iface][1] = -(ixyz + numCorner + 1);
707  faces_vec[iface][2] = ixyz + 2;
708  faces_vec[iface][3] = -(ixyz - numCorner + 2);
709  }
710  }
711  else // Last side joins ends...
712  {
713  if (iCorner < numCorner - 1)
714  {
715  faces_vec[iface][0] = ixyz + 1;
716  faces_vec[iface][1] = -(ixyz + numCorner - nFaces + 1);
717  faces_vec[iface][2] = ixyz + numCorner - nFaces + 2;
718  faces_vec[iface][3] = -(ixyz + 2);
719  }
720  else
721  {
722  faces_vec[iface][0] = ixyz + 1;
723  faces_vec[iface][1] = -(ixyz - nFaces + numCorner + 1);
724  faces_vec[iface][2] = ixyz - nFaces + 2;
725  faces_vec[iface][3] = -(ixyz - numCorner + 2);
726  }
727  }
728  ++ixyz;
729  ++iface;
730  }
731  phi += dPhi;
732  }
733  }
734  G4Polyhedron* polyhedron = new G4Polyhedron;
735  G4int problem = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec);
736  delete [] faces_vec;
737  delete [] xyz;
738  if (problem)
739  {
740  std::ostringstream message;
741  message << "Problem creating G4Polyhedron for: " << GetName();
742  G4Exception("G4GenericPolycone::CreatePolyhedron()", "GeomSolids1002",
743  JustWarning, message);
744  delete polyhedron;
745  return 0;
746  }
747  else
748  {
749  return polyhedron;
750  }
751 }
752 
753 #endif
G4String GetName() const
G4bool CrossesItself(G4double tolerance)
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)
G4bool fRebuildPolyhedron
G4Polyhedron * fpPolyhedron
static const G4double kInfinity
Definition: geomdefs.hh:42
G4PolyconeSideRZ * corners
CLHEP::Hep3Vector G4ThreeVector
G4double Amin() const
G4double z
Definition: TRTMaterials.hh:39
G4String name
Definition: TRTMaterials.hh:40
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4double a
Definition: TRTMaterials.hh:39
double C(double temp)
G4bool MustBeOutside(const G4ThreeVector &p) const
double B(double temperature)
int G4int
Definition: G4Types.hh:78
G4GenericPolycone(const G4String &name, G4double phiStart, G4double phiTotal, G4int numRZ, const G4double r[], const G4double z[])
G4GeometryType GetEntityType() const
G4Polyhedron * CreatePolyhedron() const
G4bool RemoveDuplicateVertices(G4double tolerance)
G4ThreeVector GetPointOnSurface() const
G4EnclosingCylinder * enclosingCylinder
void Create(G4double phiStart, G4double phiTotal, G4ReduciblePolygon *rz)
G4bool RemoveRedundantVertices(G4double tolerance)
G4ThreeVector GetPointOnSurfaceGeneric() const
double A(double temperature)
G4int NumVertices() const
void CopyStuff(const G4GenericPolycone &source)
bool G4bool
Definition: G4Types.hh:79
std::ostream & StreamInfo(std::ostream &os) const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
static const double twopi
Definition: G4SIunits.hh:75
G4VCSGface ** faces
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4GenericPolycone & operator=(const G4GenericPolycone &source)
static const G4double ab
EInside Inside(const G4ThreeVector &p) const
EInside
Definition: geomdefs.hh:58
static const double degree
Definition: G4SIunits.hh:143
G4bool BisectedBy(G4double a1, G4double b1, G4double a2, G4double b2, G4double tolerance)
#define G4endl
Definition: G4ios.hh:61
G4double kCarTolerance
Definition: G4VSolid.hh:304
double G4double
Definition: G4Types.hh:76
G4VSolid * Clone() const
virtual EInside Inside(const G4ThreeVector &p) const
G4bool ShouldMiss(const G4ThreeVector &p, const G4ThreeVector &v) const