Geant4  10.01.p02
G4Polyhedra.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: G4Polyhedra.cc 88948 2015-03-16 16:26:50Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class source file
32 //
33 //
34 // G4Polyhedra.cc
35 //
36 // Implementation of a CSG polyhedra, as an inherited class of G4VCSGfaceted.
37 //
38 // To be done:
39 // * Cracks: there are probably small cracks in the seams between the
40 // phi face (G4PolyPhiFace) and sides (G4PolyhedraSide) that are not
41 // entirely leakproof. Also, I am not sure all vertices are leak proof.
42 // * Many optimizations are possible, but not implemented.
43 // * Visualization needs to be updated outside of this routine.
44 //
45 // Utility classes:
46 // * G4EnclosingCylinder: I decided a quick check of geometry would be a
47 // good idea (for CPU speed). If the quick check fails, the regular
48 // full-blown G4VCSGfaceted version is invoked.
49 // * G4ReduciblePolygon: Really meant as a check of input parameters,
50 // this utility class also "converts" the GEANT3-like PGON/PCON
51 // arguments into the newer ones.
52 // Both these classes are implemented outside this file because they are
53 // shared with G4Polycone.
54 //
55 // --------------------------------------------------------------------
56 
57 #include "G4Polyhedra.hh"
58 
59 #if !defined(G4GEOM_USE_UPOLYHEDRA)
60 
61 #include "G4PolyhedraSide.hh"
62 #include "G4PolyPhiFace.hh"
63 
64 #include "Randomize.hh"
65 
66 #include "G4EnclosingCylinder.hh"
67 #include "G4ReduciblePolygon.hh"
68 #include "G4VPVParameterisation.hh"
69 
70 #include <sstream>
71 
72 using namespace CLHEP;
73 
74 //
75 // Constructor (GEANT3 style parameters)
76 //
77 // GEANT3 PGON radii are specified in the distance to the norm of each face.
78 //
80  G4double phiStart,
81  G4double thePhiTotal,
82  G4int theNumSide,
83  G4int numZPlanes,
84  const G4double zPlane[],
85  const G4double rInner[],
86  const G4double rOuter[] )
87  : G4VCSGfaceted( name ), genericPgon(false)
88 {
89  if (theNumSide <= 0)
90  {
91  std::ostringstream message;
92  message << "Solid must have at least one side - " << GetName() << G4endl
93  << " No sides specified !";
94  G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
95  FatalErrorInArgument, message);
96  }
97 
98  //
99  // Calculate conversion factor from G3 radius to G4 radius
100  //
101  G4double phiTotal = thePhiTotal;
102  if ( (phiTotal <=0) || (phiTotal >= twopi*(1-DBL_EPSILON)) )
103  { phiTotal = twopi; }
104  G4double convertRad = std::cos(0.5*phiTotal/theNumSide);
105 
106  //
107  // Some historical stuff
108  //
110 
111  original_parameters->numSide = theNumSide;
112  original_parameters->Start_angle = phiStart;
114  original_parameters->Num_z_planes = numZPlanes;
115  original_parameters->Z_values = new G4double[numZPlanes];
116  original_parameters->Rmin = new G4double[numZPlanes];
117  original_parameters->Rmax = new G4double[numZPlanes];
118 
119  G4int i;
120  for (i=0; i<numZPlanes; i++)
121  {
122  if (( i < numZPlanes-1) && ( zPlane[i] == zPlane[i+1] ))
123  {
124  if( (rInner[i] > rOuter[i+1])
125  ||(rInner[i+1] > rOuter[i]) )
126  {
127  DumpInfo();
128  std::ostringstream message;
129  message << "Cannot create a Polyhedra with no contiguous segments."
130  << G4endl
131  << " Segments are not contiguous !" << G4endl
132  << " rMin[" << i << "] = " << rInner[i]
133  << " -- rMax[" << i+1 << "] = " << rOuter[i+1] << G4endl
134  << " rMin[" << i+1 << "] = " << rInner[i+1]
135  << " -- rMax[" << i << "] = " << rOuter[i];
136  G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
137  FatalErrorInArgument, message);
138  }
139  }
140  original_parameters->Z_values[i] = zPlane[i];
141  original_parameters->Rmin[i] = rInner[i]/convertRad;
142  original_parameters->Rmax[i] = rOuter[i]/convertRad;
143  }
144 
145 
146  //
147  // Build RZ polygon using special PCON/PGON GEANT3 constructor
148  //
149  G4ReduciblePolygon *rz =
150  new G4ReduciblePolygon( rInner, rOuter, zPlane, numZPlanes );
151  rz->ScaleA( 1/convertRad );
152 
153  //
154  // Do the real work
155  //
156  Create( phiStart, phiTotal, theNumSide, rz );
157 
158  delete rz;
159 }
160 
161 
162 //
163 // Constructor (generic parameters)
164 //
166  G4double phiStart,
167  G4double phiTotal,
168  G4int theNumSide,
169  G4int numRZ,
170  const G4double r[],
171  const G4double z[] )
172  : G4VCSGfaceted( name ), genericPgon(true)
173 {
174  if (theNumSide <= 0)
175  {
176  std::ostringstream message;
177  message << "Solid must have at least one side - " << GetName() << G4endl
178  << " No sides specified !";
179  G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
180  FatalErrorInArgument, message);
181  }
182 
183  G4ReduciblePolygon *rz = new G4ReduciblePolygon( r, z, numRZ );
184 
185  Create( phiStart, phiTotal, theNumSide, rz );
186 
187  // Set original_parameters struct for consistency
188  //
190 
191  delete rz;
192 }
193 
194 
195 //
196 // Create
197 //
198 // Generic create routine, called by each constructor
199 // after conversion of arguments
200 //
202  G4double phiTotal,
203  G4int theNumSide,
204  G4ReduciblePolygon *rz )
205 {
206  //
207  // Perform checks of rz values
208  //
209  if (rz->Amin() < 0.0)
210  {
211  std::ostringstream message;
212  message << "Illegal input parameters - " << GetName() << G4endl
213  << " All R values must be >= 0 !";
214  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
215  FatalErrorInArgument, message);
216  }
217 
218  G4double rzArea = rz->Area();
219  if (rzArea < -kCarTolerance)
220  rz->ReverseOrder();
221 
222  else if (rzArea < -kCarTolerance)
223  {
224  std::ostringstream message;
225  message << "Illegal input parameters - " << GetName() << G4endl
226  << " R/Z cross section is zero or near zero: " << rzArea;
227  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
228  FatalErrorInArgument, message);
229  }
230 
232  || (!rz->RemoveRedundantVertices( kCarTolerance )) )
233  {
234  std::ostringstream message;
235  message << "Illegal input parameters - " << GetName() << G4endl
236  << " Too few unique R/Z values !";
237  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
238  FatalErrorInArgument, message);
239  }
240 
241  if (rz->CrossesItself( 1/kInfinity ))
242  {
243  std::ostringstream message;
244  message << "Illegal input parameters - " << GetName() << G4endl
245  << " R/Z segments cross !";
246  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
247  FatalErrorInArgument, message);
248  }
249 
250  numCorner = rz->NumVertices();
251 
252 
253  startPhi = phiStart;
254  while( startPhi < 0 ) startPhi += twopi;
255  //
256  // Phi opening? Account for some possible roundoff, and interpret
257  // nonsense value as representing no phi opening
258  //
259  if ( (phiTotal <= 0) || (phiTotal > twopi*(1-DBL_EPSILON)) )
260  {
261  phiIsOpen = false;
262  endPhi = phiStart+twopi;
263  }
264  else
265  {
266  phiIsOpen = true;
267 
268  //
269  // Convert phi into our convention
270  //
271  endPhi = phiStart+phiTotal;
272  while( endPhi < startPhi ) endPhi += twopi;
273  }
274 
275  //
276  // Save number sides
277  //
278  numSide = theNumSide;
279 
280  //
281  // Allocate corner array.
282  //
284 
285  //
286  // Copy corners
287  //
288  G4ReduciblePolygonIterator iterRZ(rz);
289 
290  G4PolyhedraSideRZ *next = corners;
291  iterRZ.Begin();
292  do
293  {
294  next->r = iterRZ.GetA();
295  next->z = iterRZ.GetB();
296  } while( ++next, iterRZ.Next() );
297 
298  //
299  // Allocate face pointer array
300  //
302  faces = new G4VCSGface*[numFace];
303 
304  //
305  // Construct side faces
306  //
307  // To do so properly, we need to keep track of four successive RZ
308  // corners.
309  //
310  // But! Don't construct a face if both points are at zero radius!
311  //
312  G4PolyhedraSideRZ *corner = corners,
313  *prev = corners + numCorner-1,
314  *nextNext;
315  G4VCSGface **face = faces;
316  do
317  {
318  next = corner+1;
319  if (next >= corners+numCorner) next = corners;
320  nextNext = next+1;
321  if (nextNext >= corners+numCorner) nextNext = corners;
322 
323  if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
324 /*
325  // We must decide here if we can dare declare one of our faces
326  // as having a "valid" normal (i.e. allBehind = true). This
327  // is never possible if the face faces "inward" in r *unless*
328  // we have only one side
329  //
330  G4bool allBehind;
331  if ((corner->z > next->z) && (numSide > 1))
332  {
333  allBehind = false;
334  }
335  else
336  {
337  //
338  // Otherwise, it is only true if the line passing
339  // through the two points of the segment do not
340  // split the r/z cross section
341  //
342  allBehind = !rz->BisectedBy( corner->r, corner->z,
343  next->r, next->z, kCarTolerance );
344  }
345 */
346  *face++ = new G4PolyhedraSide( prev, corner, next, nextNext,
348  } while( prev=corner, corner=next, corner > corners );
349 
350  if (phiIsOpen)
351  {
352  //
353  // Construct phi open edges
354  //
355  *face++ = new G4PolyPhiFace( rz, startPhi, phiTotal/numSide, endPhi );
356  *face++ = new G4PolyPhiFace( rz, endPhi, phiTotal/numSide, startPhi );
357  }
358 
359  //
360  // We might have dropped a face or two: recalculate numFace
361  //
362  numFace = face-faces;
363 
364  //
365  // Make enclosingCylinder
366  //
368  new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
369 }
370 
371 
372 //
373 // Fake default constructor - sets only member data and allocates memory
374 // for usage restricted to object persistency.
375 //
377  : G4VCSGfaceted(a), numSide(0), startPhi(0.), endPhi(0.),
378  phiIsOpen(false), genericPgon(false), numCorner(0), corners(0),
379  original_parameters(0), enclosingCylinder(0)
380 {
381 }
382 
383 
384 //
385 // Destructor
386 //
388 {
389  delete [] corners;
391 
392  delete enclosingCylinder;
393 }
394 
395 
396 //
397 // Copy constructor
398 //
400  : G4VCSGfaceted( source )
401 {
402  CopyStuff( source );
403 }
404 
405 
406 //
407 // Assignment operator
408 //
410 {
411  if (this == &source) return *this;
412 
413  G4VCSGfaceted::operator=( source );
414 
415  delete [] corners;
417 
418  delete enclosingCylinder;
419 
420  CopyStuff( source );
421 
422  return *this;
423 }
424 
425 
426 //
427 // CopyStuff
428 //
429 void G4Polyhedra::CopyStuff( const G4Polyhedra &source )
430 {
431  //
432  // Simple stuff
433  //
434  numSide = source.numSide;
435  startPhi = source.startPhi;
436  endPhi = source.endPhi;
437  phiIsOpen = source.phiIsOpen;
438  numCorner = source.numCorner;
439  genericPgon= source.genericPgon;
440 
441  //
442  // The corner array
443  //
445 
446  G4PolyhedraSideRZ *corn = corners,
447  *sourceCorn = source.corners;
448  do
449  {
450  *corn = *sourceCorn;
451  } while( ++sourceCorn, ++corn < corners+numCorner );
452 
453  //
454  // Original parameters
455  //
456  if (source.original_parameters)
457  {
460  }
461 
462  //
463  // Enclosing cylinder
464  //
466 
467  fRebuildPolyhedron = false;
468  fpPolyhedron = 0;
469 }
470 
471 
472 //
473 // Reset
474 //
475 // Recalculates and reshapes the solid, given pre-assigned scaled
476 // original_parameters.
477 //
479 {
480  if (genericPgon)
481  {
482  std::ostringstream message;
483  message << "Solid " << GetName() << " built using generic construct."
484  << G4endl << "Not applicable to the generic construct !";
485  G4Exception("G4Polyhedra::Reset()", "GeomSolids1001",
486  JustWarning, message, "Parameters NOT resetted.");
487  return 1;
488  }
489 
490  //
491  // Clear old setup
492  //
494  delete [] corners;
495  delete enclosingCylinder;
496 
497  //
498  // Rebuild polyhedra
499  //
500  G4ReduciblePolygon *rz =
508  delete rz;
509 
510  return 0;
511 }
512 
513 
514 //
515 // Inside
516 //
517 // This is an override of G4VCSGfaceted::Inside, created in order
518 // to speed things up by first checking with G4EnclosingCylinder.
519 //
521 {
522  //
523  // Quick test
524  //
525  if (enclosingCylinder->MustBeOutside(p)) return kOutside;
526 
527  //
528  // Long answer
529  //
530  return G4VCSGfaceted::Inside(p);
531 }
532 
533 
534 //
535 // DistanceToIn
536 //
537 // This is an override of G4VCSGfaceted::Inside, created in order
538 // to speed things up by first checking with G4EnclosingCylinder.
539 //
541  const G4ThreeVector &v ) const
542 {
543  //
544  // Quick test
545  //
546  if (enclosingCylinder->ShouldMiss(p,v))
547  return kInfinity;
548 
549  //
550  // Long answer
551  //
552  return G4VCSGfaceted::DistanceToIn( p, v );
553 }
554 
555 
556 //
557 // DistanceToIn
558 //
560 {
561  return G4VCSGfaceted::DistanceToIn(p);
562 }
563 
564 
565 //
566 // ComputeDimensions
567 //
569  const G4int n,
570  const G4VPhysicalVolume* pRep )
571 {
572  p->ComputeDimensions(*this,n,pRep);
573 }
574 
575 
576 //
577 // GetEntityType
578 //
580 {
581  return G4String("G4Polyhedra");
582 }
583 
584 
585 //
586 // Make a clone of the object
587 //
589 {
590  return new G4Polyhedra(*this);
591 }
592 
593 
594 //
595 // Stream object contents to an output stream
596 //
597 std::ostream& G4Polyhedra::StreamInfo( std::ostream& os ) const
598 {
599  G4int oldprc = os.precision(16);
600  os << "-----------------------------------------------------------\n"
601  << " *** Dump for solid - " << GetName() << " ***\n"
602  << " ===================================================\n"
603  << " Solid type: G4Polyhedra\n"
604  << " Parameters: \n"
605  << " starting phi angle : " << startPhi/degree << " degrees \n"
606  << " ending phi angle : " << endPhi/degree << " degrees \n"
607  << " number of sides : " << numSide << " \n";
608  G4int i=0;
609  if (!genericPgon)
610  {
612  os << " number of Z planes: " << numPlanes << "\n"
613  << " Z values: \n";
614  for (i=0; i<numPlanes; i++)
615  {
616  os << " Z plane " << i << ": "
617  << original_parameters->Z_values[i] << "\n";
618  }
619  os << " Tangent distances to inner surface (Rmin): \n";
620  for (i=0; i<numPlanes; i++)
621  {
622  os << " Z plane " << i << ": "
623  << original_parameters->Rmin[i] << "\n";
624  }
625  os << " Tangent distances to outer surface (Rmax): \n";
626  for (i=0; i<numPlanes; i++)
627  {
628  os << " Z plane " << i << ": "
629  << original_parameters->Rmax[i] << "\n";
630  }
631  }
632  os << " number of RZ points: " << numCorner << "\n"
633  << " RZ values (corners): \n";
634  for (i=0; i<numCorner; i++)
635  {
636  os << " "
637  << corners[i].r << ", " << corners[i].z << "\n";
638  }
639  os << "-----------------------------------------------------------\n";
640  os.precision(oldprc);
641 
642  return os;
643 }
644 
645 
646 //
647 // GetPointOnPlane
648 //
649 // Auxiliary method for get point on surface
650 //
652  G4ThreeVector p2, G4ThreeVector p3) const
653 {
654  G4double lambda1, lambda2, chose,aOne,aTwo;
655  G4ThreeVector t, u, v, w, Area, normal;
656  aOne = 1.;
657  aTwo = 1.;
658 
659  t = p1 - p0;
660  u = p2 - p1;
661  v = p3 - p2;
662  w = p0 - p3;
663 
664  chose = RandFlat::shoot(0.,aOne+aTwo);
665  if( (chose>=0.) && (chose < aOne) )
666  {
667  lambda1 = RandFlat::shoot(0.,1.);
668  lambda2 = RandFlat::shoot(0.,lambda1);
669  return (p2+lambda1*v+lambda2*w);
670  }
671 
672  lambda1 = RandFlat::shoot(0.,1.);
673  lambda2 = RandFlat::shoot(0.,lambda1);
674  return (p0+lambda1*t+lambda2*u);
675 }
676 
677 
678 //
679 // GetPointOnTriangle
680 //
681 // Auxiliary method for get point on surface
682 //
685  G4ThreeVector p3) const
686 {
687  G4double lambda1,lambda2;
688  G4ThreeVector v=p3-p1, w=p1-p2;
689 
690  lambda1 = RandFlat::shoot(0.,1.);
691  lambda2 = RandFlat::shoot(0.,lambda1);
692 
693  return (p2 + lambda1*w + lambda2*v);
694 }
695 
696 
697 //
698 // GetPointOnSurface
699 //
701 {
702  if( !genericPgon ) // Polyhedra by faces
703  {
704  G4int j, numPlanes = original_parameters->Num_z_planes, Flag=0;
705  G4double chose, totArea=0., Achose1, Achose2,
706  rad1, rad2, sinphi1, sinphi2, cosphi1, cosphi2;
707  G4double a, b, l2, rang, totalPhi, ksi,
708  area, aTop=0., aBottom=0., zVal=0.;
709 
710  G4ThreeVector p0, p1, p2, p3;
711  std::vector<G4double> aVector1;
712  std::vector<G4double> aVector2;
713  std::vector<G4double> aVector3;
714 
715  totalPhi= (phiIsOpen) ? (endPhi-startPhi) : twopi;
716  ksi = totalPhi/numSide;
717  G4double cosksi = std::cos(ksi/2.);
718 
719  // Below we generate the areas relevant to our solid
720  //
721  for(j=0; j<numPlanes-1; j++)
722  {
723  a = original_parameters->Rmax[j+1];
724  b = original_parameters->Rmax[j];
726  -original_parameters->Z_values[j+1]) + sqr(b-a);
727  area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
728  aVector1.push_back(area);
729  }
730 
731  for(j=0; j<numPlanes-1; j++)
732  {
733  a = original_parameters->Rmin[j+1];//*cosksi;
734  b = original_parameters->Rmin[j];//*cosksi;
736  -original_parameters->Z_values[j+1]) + sqr(b-a);
737  area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
738  aVector2.push_back(area);
739  }
740 
741  for(j=0; j<numPlanes-1; j++)
742  {
743  if(phiIsOpen == true)
744  {
745  aVector3.push_back(0.5*(original_parameters->Rmax[j]
748  -original_parameters->Rmin[j+1])
749  *std::fabs(original_parameters->Z_values[j+1]
751  }
752  else { aVector3.push_back(0.); }
753  }
754 
755  for(j=0; j<numPlanes-1; j++)
756  {
757  totArea += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
758  }
759 
760  // Must include top and bottom areas
761  //
762  if(original_parameters->Rmax[numPlanes-1] != 0.)
763  {
764  a = original_parameters->Rmax[numPlanes-1];
765  b = original_parameters->Rmin[numPlanes-1];
766  l2 = sqr(a-b);
767  aTop = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
768  }
769 
770  if(original_parameters->Rmax[0] != 0.)
771  {
772  a = original_parameters->Rmax[0];
773  b = original_parameters->Rmin[0];
774  l2 = sqr(a-b);
775  aBottom = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
776  }
777 
778  Achose1 = 0.;
779  Achose2 = numSide*(aVector1[0]+aVector2[0])+2.*aVector3[0];
780 
781  chose = RandFlat::shoot(0.,totArea+aTop+aBottom);
782  if( (chose >= 0.) && (chose < aTop + aBottom) )
783  {
784  chose = RandFlat::shoot(startPhi,startPhi+totalPhi);
785  rang = std::floor((chose-startPhi)/ksi-0.01);
786  if(rang<0) { rang=0; }
787  rang = std::fabs(rang);
788  sinphi1 = std::sin(startPhi+rang*ksi);
789  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
790  cosphi1 = std::cos(startPhi+rang*ksi);
791  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
792  chose = RandFlat::shoot(0., aTop + aBottom);
793  if(chose>=0. && chose<aTop)
794  {
795  rad1 = original_parameters->Rmin[numPlanes-1];
796  rad2 = original_parameters->Rmax[numPlanes-1];
797  zVal = original_parameters->Z_values[numPlanes-1];
798  }
799  else
800  {
801  rad1 = original_parameters->Rmin[0];
802  rad2 = original_parameters->Rmax[0];
803  zVal = original_parameters->Z_values[0];
804  }
805  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
806  p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
807  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
808  p3 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
809  return GetPointOnPlane(p0,p1,p2,p3);
810  }
811  else
812  {
813  for (j=0; j<numPlanes-1; j++)
814  {
815  if( ((chose >= Achose1) && (chose < Achose2)) || (j == numPlanes-2) )
816  {
817  Flag = j; break;
818  }
819  Achose1 += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
820  Achose2 = Achose1 + numSide*(aVector1[j+1]+aVector2[j+1])
821  + 2.*aVector3[j+1];
822  }
823  }
824 
825  // At this point we have chosen a subsection
826  // between to adjacent plane cuts...
827 
828  j = Flag;
829 
830  totArea = numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
831  chose = RandFlat::shoot(0.,totArea);
832 
833  if( (chose>=0.) && (chose<numSide*aVector1[j]) )
834  {
835  chose = RandFlat::shoot(startPhi,startPhi+totalPhi);
836  rang = std::floor((chose-startPhi)/ksi-0.01);
837  if(rang<0) { rang=0; }
838  rang = std::fabs(rang);
839  rad1 = original_parameters->Rmax[j];
840  rad2 = original_parameters->Rmax[j+1];
841  sinphi1 = std::sin(startPhi+rang*ksi);
842  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
843  cosphi1 = std::cos(startPhi+rang*ksi);
844  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
845  zVal = original_parameters->Z_values[j];
846 
847  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
848  p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
849 
850  zVal = original_parameters->Z_values[j+1];
851 
852  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
853  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
854  return GetPointOnPlane(p0,p1,p2,p3);
855  }
856  else if ( (chose >= numSide*aVector1[j])
857  && (chose <= numSide*(aVector1[j]+aVector2[j])) )
858  {
859  chose = RandFlat::shoot(startPhi,startPhi+totalPhi);
860  rang = std::floor((chose-startPhi)/ksi-0.01);
861  if(rang<0) { rang=0; }
862  rang = std::fabs(rang);
863  rad1 = original_parameters->Rmin[j];
864  rad2 = original_parameters->Rmin[j+1];
865  sinphi1 = std::sin(startPhi+rang*ksi);
866  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
867  cosphi1 = std::cos(startPhi+rang*ksi);
868  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
869  zVal = original_parameters->Z_values[j];
870 
871  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
872  p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
873 
874  zVal = original_parameters->Z_values[j+1];
875 
876  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
877  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
878  return GetPointOnPlane(p0,p1,p2,p3);
879  }
880 
881  chose = RandFlat::shoot(0.,2.2);
882  if( (chose>=0.) && (chose < 1.) )
883  {
884  rang = startPhi;
885  }
886  else
887  {
888  rang = endPhi;
889  }
890 
891  cosphi1 = std::cos(rang); rad1 = original_parameters->Rmin[j];
892  sinphi1 = std::sin(rang); rad2 = original_parameters->Rmax[j];
893 
894  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
896  p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
898 
899  rad1 = original_parameters->Rmax[j+1];
900  rad2 = original_parameters->Rmin[j+1];
901 
902  p2 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
904  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
906  return GetPointOnPlane(p0,p1,p2,p3);
907  }
908  else // Generic polyhedra
909  {
910  return GetPointOnSurfaceGeneric();
911  }
912 }
913 
914 //
915 // CreatePolyhedron
916 //
918 {
919  if (!genericPgon)
920  {
928  }
929  else
930  {
931  // The following code prepares for:
932  // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
933  // const double xyz[][3],
934  // const int faces_vec[][4])
935  // Here is an extract from the header file HepPolyhedron.h:
953  G4int nNodes;
954  G4int nFaces;
955  typedef G4double double3[3];
956  double3* xyz;
957  typedef G4int int4[4];
958  int4* faces_vec;
959  if (phiIsOpen)
960  {
961  // Triangulate open ends. Simple ear-chopping algorithm...
962  // I'm not sure how robust this algorithm is (J.Allison).
963  //
964  std::vector<G4bool> chopped(numCorner, false);
965  std::vector<G4int*> triQuads;
966  G4int remaining = numCorner;
967  G4int iStarter = 0;
968  while (remaining >= 3)
969  {
970  // Find unchopped corners...
971  //
972  G4int A = -1, B = -1, C = -1;
973  G4int iStepper = iStarter;
974  do
975  {
976  if (A < 0) { A = iStepper; }
977  else if (B < 0) { B = iStepper; }
978  else if (C < 0) { C = iStepper; }
979  do
980  {
981  if (++iStepper >= numCorner) iStepper = 0;
982  }
983  while (chopped[iStepper]);
984  }
985  while (C < 0 && iStepper != iStarter);
986 
987  // Check triangle at B is pointing outward (an "ear").
988  // Sign of z cross product determines...
989 
990  G4double BAr = corners[A].r - corners[B].r;
991  G4double BAz = corners[A].z - corners[B].z;
992  G4double BCr = corners[C].r - corners[B].r;
993  G4double BCz = corners[C].z - corners[B].z;
994  if (BAr * BCz - BAz * BCr < kCarTolerance)
995  {
996  G4int* tq = new G4int[3];
997  tq[0] = A + 1;
998  tq[1] = B + 1;
999  tq[2] = C + 1;
1000  triQuads.push_back(tq);
1001  chopped[B] = true;
1002  --remaining;
1003  }
1004  else
1005  {
1006  do
1007  {
1008  if (++iStarter >= numCorner) { iStarter = 0; }
1009  }
1010  while (chopped[iStarter]);
1011  }
1012  }
1013 
1014  // Transfer to faces...
1015 
1016  nNodes = (numSide + 1) * numCorner;
1017  nFaces = numSide * numCorner + 2 * triQuads.size();
1018  faces_vec = new int4[nFaces];
1019  G4int iface = 0;
1020  G4int addition = numCorner * numSide;
1021  G4int d = numCorner - 1;
1022  for (G4int iEnd = 0; iEnd < 2; ++iEnd)
1023  {
1024  for (size_t i = 0; i < triQuads.size(); ++i)
1025  {
1026  // Negative for soft/auxiliary/normally invisible edges...
1027  //
1028  G4int a, b, c;
1029  if (iEnd == 0)
1030  {
1031  a = triQuads[i][0];
1032  b = triQuads[i][1];
1033  c = triQuads[i][2];
1034  }
1035  else
1036  {
1037  a = triQuads[i][0] + addition;
1038  b = triQuads[i][2] + addition;
1039  c = triQuads[i][1] + addition;
1040  }
1041  G4int ab = std::abs(b - a);
1042  G4int bc = std::abs(c - b);
1043  G4int ca = std::abs(a - c);
1044  faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
1045  faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
1046  faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
1047  faces_vec[iface][3] = 0;
1048  ++iface;
1049  }
1050  }
1051 
1052  // Continue with sides...
1053 
1054  xyz = new double3[nNodes];
1055  const G4double dPhi = (endPhi - startPhi) / numSide;
1056  G4double phi = startPhi;
1057  G4int ixyz = 0;
1058  for (G4int iSide = 0; iSide < numSide; ++iSide)
1059  {
1060  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1061  {
1062  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1063  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1064  xyz[ixyz][2] = corners[iCorner].z;
1065  if (iCorner < numCorner - 1)
1066  {
1067  faces_vec[iface][0] = ixyz + 1;
1068  faces_vec[iface][1] = ixyz + numCorner + 1;
1069  faces_vec[iface][2] = ixyz + numCorner + 2;
1070  faces_vec[iface][3] = ixyz + 2;
1071  }
1072  else
1073  {
1074  faces_vec[iface][0] = ixyz + 1;
1075  faces_vec[iface][1] = ixyz + numCorner + 1;
1076  faces_vec[iface][2] = ixyz + 2;
1077  faces_vec[iface][3] = ixyz - numCorner + 2;
1078  }
1079  ++iface;
1080  ++ixyz;
1081  }
1082  phi += dPhi;
1083  }
1084 
1085  // Last corners...
1086 
1087  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1088  {
1089  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1090  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1091  xyz[ixyz][2] = corners[iCorner].z;
1092  ++ixyz;
1093  }
1094  }
1095  else // !phiIsOpen - i.e., a complete 360 degrees.
1096  {
1097  nNodes = numSide * numCorner;
1098  nFaces = numSide * numCorner;;
1099  xyz = new double3[nNodes];
1100  faces_vec = new int4[nFaces];
1101  // const G4double dPhi = (endPhi - startPhi) / numSide;
1102  const G4double dPhi = twopi / numSide; // !phiIsOpen endPhi-startPhi = 360 degrees.
1103  G4double phi = startPhi;
1104  G4int ixyz = 0, iface = 0;
1105  for (G4int iSide = 0; iSide < numSide; ++iSide)
1106  {
1107  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1108  {
1109  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1110  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1111  xyz[ixyz][2] = corners[iCorner].z;
1112  if (iSide < numSide - 1)
1113  {
1114  if (iCorner < numCorner - 1)
1115  {
1116  faces_vec[iface][0] = ixyz + 1;
1117  faces_vec[iface][1] = ixyz + numCorner + 1;
1118  faces_vec[iface][2] = ixyz + numCorner + 2;
1119  faces_vec[iface][3] = ixyz + 2;
1120  }
1121  else
1122  {
1123  faces_vec[iface][0] = ixyz + 1;
1124  faces_vec[iface][1] = ixyz + numCorner + 1;
1125  faces_vec[iface][2] = ixyz + 2;
1126  faces_vec[iface][3] = ixyz - numCorner + 2;
1127  }
1128  }
1129  else // Last side joins ends...
1130  {
1131  if (iCorner < numCorner - 1)
1132  {
1133  faces_vec[iface][0] = ixyz + 1;
1134  faces_vec[iface][1] = ixyz + numCorner - nFaces + 1;
1135  faces_vec[iface][2] = ixyz + numCorner - nFaces + 2;
1136  faces_vec[iface][3] = ixyz + 2;
1137  }
1138  else
1139  {
1140  faces_vec[iface][0] = ixyz + 1;
1141  faces_vec[iface][1] = ixyz - nFaces + numCorner + 1;
1142  faces_vec[iface][2] = ixyz - nFaces + 2;
1143  faces_vec[iface][3] = ixyz - numCorner + 2;
1144  }
1145  }
1146  ++ixyz;
1147  ++iface;
1148  }
1149  phi += dPhi;
1150  }
1151  }
1152  G4Polyhedron* polyhedron = new G4Polyhedron;
1153  G4int problem = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec);
1154  delete [] faces_vec;
1155  delete [] xyz;
1156  if (problem)
1157  {
1158  std::ostringstream message;
1159  message << "Problem creating G4Polyhedron for: " << GetName();
1160  G4Exception("G4Polyhedra::CreatePolyhedron()", "GeomSolids1002",
1161  JustWarning, message);
1162  delete polyhedron;
1163  return 0;
1164  }
1165  else
1166  {
1167  return polyhedron;
1168  }
1169  }
1170 }
1171 
1172 
1174 {
1175  G4int numPlanes = (G4int)numCorner;
1176  G4bool isConvertible=true;
1177  G4double Zmax=rz->Bmax();
1178  rz->StartWithZMin();
1179 
1180  // Prepare vectors for storage
1181  //
1182  std::vector<G4double> Z;
1183  std::vector<G4double> Rmin;
1184  std::vector<G4double> Rmax;
1185 
1186  G4int countPlanes=1;
1187  G4int icurr=0;
1188  G4int icurl=0;
1189 
1190  // first plane Z=Z[0]
1191  //
1192  Z.push_back(corners[0].z);
1193  G4double Zprev=Z[0];
1194  if (Zprev == corners[1].z)
1195  {
1196  Rmin.push_back(corners[0].r);
1197  Rmax.push_back (corners[1].r);icurr=1;
1198  }
1199  else if (Zprev == corners[numPlanes-1].z)
1200  {
1201  Rmin.push_back(corners[numPlanes-1].r);
1202  Rmax.push_back (corners[0].r);
1203  icurl=numPlanes-1;
1204  }
1205  else
1206  {
1207  Rmin.push_back(corners[0].r);
1208  Rmax.push_back (corners[0].r);
1209  }
1210 
1211  // next planes until last
1212  //
1213  G4int inextr=0, inextl=0;
1214  for (G4int i=0; i < numPlanes-2; i++)
1215  {
1216  inextr=1+icurr;
1217  inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1218 
1219  if((corners[inextr].z >= Zmax) & (corners[inextl].z >= Zmax)) { break; }
1220 
1221  G4double Zleft = corners[inextl].z;
1222  G4double Zright = corners[inextr].z;
1223  if(Zright>Zleft)
1224  {
1225  Z.push_back(Zleft);
1226  countPlanes++;
1227  G4double difZr=corners[inextr].z - corners[icurr].z;
1228  G4double difZl=corners[inextl].z - corners[icurl].z;
1229 
1230  if(std::fabs(difZl) < kCarTolerance)
1231  {
1232  if(std::fabs(difZr) < kCarTolerance)
1233  {
1234  Rmin.push_back(corners[inextl].r);
1235  Rmax.push_back(corners[icurr].r);
1236  }
1237  else
1238  {
1239  Rmin.push_back(corners[inextl].r);
1240  Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1241  *(corners[inextr].r - corners[icurr].r));
1242  }
1243  }
1244  else if (difZl >= kCarTolerance)
1245  {
1246  if(std::fabs(difZr) < kCarTolerance)
1247  {
1248  Rmin.push_back(corners[icurl].r);
1249  Rmax.push_back(corners[icurr].r);
1250  }
1251  else
1252  {
1253  Rmin.push_back(corners[icurl].r);
1254  Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1255  *(corners[inextr].r - corners[icurr].r));
1256  }
1257  }
1258  else
1259  {
1260  isConvertible=false; break;
1261  }
1262  icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1263  }
1264  else if(std::fabs(Zright-Zleft)<kCarTolerance) // Zright=Zleft
1265  {
1266  Z.push_back(Zleft);
1267  countPlanes++;
1268  icurr++;
1269 
1270  icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1271 
1272  Rmin.push_back(corners[inextl].r);
1273  Rmax.push_back (corners[inextr].r);
1274  }
1275  else // Zright<Zleft
1276  {
1277  Z.push_back(Zright);
1278  countPlanes++;
1279 
1280  G4double difZr=corners[inextr].z - corners[icurr].z;
1281  G4double difZl=corners[inextl].z - corners[icurl].z;
1282  if(std::fabs(difZr) < kCarTolerance)
1283  {
1284  if(std::fabs(difZl) < kCarTolerance)
1285  {
1286  Rmax.push_back(corners[inextr].r);
1287  Rmin.push_back(corners[icurr].r);
1288  }
1289  else
1290  {
1291  Rmin.push_back(corners[icurl].r + (Zright-corners[icurl].z)/difZl
1292  * (corners[inextl].r - corners[icurl].r));
1293  Rmax.push_back(corners[inextr].r);
1294  }
1295  icurr++;
1296  } // plate
1297  else if (difZr >= kCarTolerance)
1298  {
1299  if(std::fabs(difZl) < kCarTolerance)
1300  {
1301  Rmax.push_back(corners[inextr].r);
1302  Rmin.push_back (corners[icurr].r);
1303  }
1304  else
1305  {
1306  Rmax.push_back(corners[inextr].r);
1307  Rmin.push_back (corners[icurl].r+(Zright-corners[icurl].z)/difZl
1308  * (corners[inextl].r - corners[icurl].r));
1309  }
1310  icurr++;
1311  }
1312  else
1313  {
1314  isConvertible=false; break;
1315  }
1316  }
1317  } // end for loop
1318 
1319  // last plane Z=Zmax
1320  //
1321  Z.push_back(Zmax);
1322  countPlanes++;
1323  inextr=1+icurr;
1324  inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1325 
1326  if(corners[inextr].z==corners[inextl].z)
1327  {
1328  Rmax.push_back(corners[inextr].r);
1329  Rmin.push_back(corners[inextl].r);
1330  }
1331  else
1332  {
1333  Rmax.push_back(corners[inextr].r);
1334  Rmin.push_back(corners[inextl].r);
1335  }
1336 
1337  // Set original parameters Rmin,Rmax,Z
1338  //
1339  if(isConvertible)
1340  {
1343  original_parameters->Z_values = new G4double[countPlanes];
1344  original_parameters->Rmin = new G4double[countPlanes];
1345  original_parameters->Rmax = new G4double[countPlanes];
1346 
1347  for(G4int j=0; j < countPlanes; j++)
1348  {
1349  original_parameters->Z_values[j] = Z[j];
1350  original_parameters->Rmax[j] = Rmax[j];
1351  original_parameters->Rmin[j] = Rmin[j];
1352  }
1355  original_parameters->Num_z_planes = countPlanes;
1356 
1357  }
1358  else // Set parameters(r,z) with Rmin==0 as convention
1359  {
1360 #ifdef G4SPECSDEBUG
1361  std::ostringstream message;
1362  message << "Polyhedra " << GetName() << G4endl
1363  << "cannot be converted to Polyhedra with (Rmin,Rmaz,Z) parameters!";
1364  G4Exception("G4Polyhedra::SetOriginalParameters()",
1365  "GeomSolids0002", JustWarning, message);
1366 #endif
1369  original_parameters->Z_values = new G4double[numPlanes];
1370  original_parameters->Rmin = new G4double[numPlanes];
1371  original_parameters->Rmax = new G4double[numPlanes];
1372 
1373  for(G4int j=0; j < numPlanes; j++)
1374  {
1376  original_parameters->Rmax[j] = corners[j].r;
1377  original_parameters->Rmin[j] = 0.0;
1378  }
1381  original_parameters->Num_z_planes = numPlanes;
1382  }
1383  //return isConvertible;
1384 }
1385 
1386 #endif
void CopyStuff(const G4Polyhedra &source)
Definition: G4Polyhedra.cc:429
G4String GetName() const
G4bool CrossesItself(G4double tolerance)
ThreeVector shoot(const G4int Ap, const G4int Af)
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)
virtual ~G4Polyhedra()
Definition: G4Polyhedra.cc:387
G4bool fRebuildPolyhedron
G4Polyhedron * fpPolyhedron
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
G4GeometryType GetEntityType() const
Definition: G4Polyhedra.cc:579
G4double Amin() const
G4double z
Definition: TRTMaterials.hh:39
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4String name
Definition: TRTMaterials.hh:40
G4ThreeVector GetPointOnPlane(G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3) const
Definition: G4Polyhedra.cc:651
G4bool genericPgon
Definition: G4Polyhedra.hh:184
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
G4double a
Definition: TRTMaterials.hh:39
G4bool MustBeOutside(const G4ThreeVector &p) const
G4bool Reset()
Definition: G4Polyhedra.cc:478
int G4int
Definition: G4Types.hh:78
void Create(G4double phiStart, G4double phiTotal, G4int numSide, G4ReduciblePolygon *rz)
Definition: G4Polyhedra.cc:201
G4ThreeVector GetPointOnSurface() const
Definition: G4Polyhedra.cc:700
void DumpInfo() const
G4double endPhi
Definition: G4Polyhedra.hh:182
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4bool RemoveDuplicateVertices(G4double tolerance)
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:189
G4bool RemoveRedundantVertices(G4double tolerance)
G4ThreeVector GetPointOnSurfaceGeneric() const
G4int numCorner
Definition: G4Polyhedra.hh:185
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Polyhedra.cc:568
G4int NumVertices() const
G4double Bmax() const
void ScaleA(G4double scale)
bool G4bool
Definition: G4Types.hh:79
G4double startPhi
Definition: G4Polyhedra.hh:181
const G4double p2
const G4double p1
#define DBL_EPSILON
Definition: templates.hh:87
const G4int n
static const G4double A[nN]
G4VCSGface ** faces
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Polyhedra.cc:540
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
EInside Inside(const G4ThreeVector &p) const
Definition: G4Polyhedra.cc:520
static const G4double ab
const G4double p0
G4Polyhedron * CreatePolyhedron() const
Definition: G4Polyhedra.cc:917
EInside
Definition: geomdefs.hh:58
static const double degree
Definition: G4SIunits.hh:125
G4Polyhedra & operator=(const G4Polyhedra &source)
Definition: G4Polyhedra.cc:409
#define G4endl
Definition: G4ios.hh:61
G4double kCarTolerance
Definition: G4VSolid.hh:305
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
void SetOriginalParameters(G4PolyhedraHistorical *pars)
G4ThreeVector GetPointOnTriangle(G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2) const
Definition: G4Polyhedra.cc:683
virtual EInside Inside(const G4ThreeVector &p) const
G4bool phiIsOpen
Definition: G4Polyhedra.hh:183
G4bool ShouldMiss(const G4ThreeVector &p, const G4ThreeVector &v) const
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Polyhedra.cc:597
G4VSolid * Clone() const
Definition: G4Polyhedra.cc:588
G4Polyhedra(const G4String &name, G4double phiStart, G4double phiTotal, G4int numSide, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])
Definition: G4Polyhedra.cc:79