Geant4  10.00.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 81689 2014-06-04 15:09:00Z 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 
468 }
469 
470 
471 //
472 // Reset
473 //
474 // Recalculates and reshapes the solid, given pre-assigned scaled
475 // original_parameters.
476 //
478 {
479  if (genericPgon)
480  {
481  std::ostringstream message;
482  message << "Solid " << GetName() << " built using generic construct."
483  << G4endl << "Not applicable to the generic construct !";
484  G4Exception("G4Polyhedra::Reset()", "GeomSolids1001",
485  JustWarning, message, "Parameters NOT resetted.");
486  return 1;
487  }
488 
489  //
490  // Clear old setup
491  //
493  delete [] corners;
494  delete enclosingCylinder;
495 
496  //
497  // Rebuild polyhedra
498  //
499  G4ReduciblePolygon *rz =
507  delete rz;
508 
509  return 0;
510 }
511 
512 
513 //
514 // Inside
515 //
516 // This is an override of G4VCSGfaceted::Inside, created in order
517 // to speed things up by first checking with G4EnclosingCylinder.
518 //
520 {
521  //
522  // Quick test
523  //
524  if (enclosingCylinder->MustBeOutside(p)) return kOutside;
525 
526  //
527  // Long answer
528  //
529  return G4VCSGfaceted::Inside(p);
530 }
531 
532 
533 //
534 // DistanceToIn
535 //
536 // This is an override of G4VCSGfaceted::Inside, created in order
537 // to speed things up by first checking with G4EnclosingCylinder.
538 //
540  const G4ThreeVector &v ) const
541 {
542  //
543  // Quick test
544  //
545  if (enclosingCylinder->ShouldMiss(p,v))
546  return kInfinity;
547 
548  //
549  // Long answer
550  //
551  return G4VCSGfaceted::DistanceToIn( p, v );
552 }
553 
554 
555 //
556 // DistanceToIn
557 //
559 {
560  return G4VCSGfaceted::DistanceToIn(p);
561 }
562 
563 
564 //
565 // ComputeDimensions
566 //
568  const G4int n,
569  const G4VPhysicalVolume* pRep )
570 {
571  p->ComputeDimensions(*this,n,pRep);
572 }
573 
574 
575 //
576 // GetEntityType
577 //
579 {
580  return G4String("G4Polyhedra");
581 }
582 
583 
584 //
585 // Make a clone of the object
586 //
588 {
589  return new G4Polyhedra(*this);
590 }
591 
592 
593 //
594 // Stream object contents to an output stream
595 //
596 std::ostream& G4Polyhedra::StreamInfo( std::ostream& os ) const
597 {
598  G4int oldprc = os.precision(16);
599  os << "-----------------------------------------------------------\n"
600  << " *** Dump for solid - " << GetName() << " ***\n"
601  << " ===================================================\n"
602  << " Solid type: G4Polyhedra\n"
603  << " Parameters: \n"
604  << " starting phi angle : " << startPhi/degree << " degrees \n"
605  << " ending phi angle : " << endPhi/degree << " degrees \n";
606  G4int i=0;
607  if (!genericPgon)
608  {
610  os << " number of Z planes: " << numPlanes << "\n"
611  << " Z values: \n";
612  for (i=0; i<numPlanes; i++)
613  {
614  os << " Z plane " << i << ": "
615  << original_parameters->Z_values[i] << "\n";
616  }
617  os << " Tangent distances to inner surface (Rmin): \n";
618  for (i=0; i<numPlanes; i++)
619  {
620  os << " Z plane " << i << ": "
621  << original_parameters->Rmin[i] << "\n";
622  }
623  os << " Tangent distances to outer surface (Rmax): \n";
624  for (i=0; i<numPlanes; i++)
625  {
626  os << " Z plane " << i << ": "
627  << original_parameters->Rmax[i] << "\n";
628  }
629  }
630  os << " number of RZ points: " << numCorner << "\n"
631  << " RZ values (corners): \n";
632  for (i=0; i<numCorner; i++)
633  {
634  os << " "
635  << corners[i].r << ", " << corners[i].z << "\n";
636  }
637  os << "-----------------------------------------------------------\n";
638  os.precision(oldprc);
639 
640  return os;
641 }
642 
643 
644 //
645 // GetPointOnPlane
646 //
647 // Auxiliary method for get point on surface
648 //
650  G4ThreeVector p2, G4ThreeVector p3) const
651 {
652  G4double lambda1, lambda2, chose,aOne,aTwo;
653  G4ThreeVector t, u, v, w, Area, normal;
654  aOne = 1.;
655  aTwo = 1.;
656 
657  t = p1 - p0;
658  u = p2 - p1;
659  v = p3 - p2;
660  w = p0 - p3;
661 
662  chose = RandFlat::shoot(0.,aOne+aTwo);
663  if( (chose>=0.) && (chose < aOne) )
664  {
665  lambda1 = RandFlat::shoot(0.,1.);
666  lambda2 = RandFlat::shoot(0.,lambda1);
667  return (p2+lambda1*v+lambda2*w);
668  }
669 
670  lambda1 = RandFlat::shoot(0.,1.);
671  lambda2 = RandFlat::shoot(0.,lambda1);
672  return (p0+lambda1*t+lambda2*u);
673 }
674 
675 
676 //
677 // GetPointOnTriangle
678 //
679 // Auxiliary method for get point on surface
680 //
682  G4ThreeVector p2,
683  G4ThreeVector p3) const
684 {
685  G4double lambda1,lambda2;
686  G4ThreeVector v=p3-p1, w=p1-p2;
687 
688  lambda1 = RandFlat::shoot(0.,1.);
689  lambda2 = RandFlat::shoot(0.,lambda1);
690 
691  return (p2 + lambda1*w + lambda2*v);
692 }
693 
694 
695 //
696 // GetPointOnSurface
697 //
699 {
700  if( !genericPgon ) // Polyhedra by faces
701  {
702  G4int j, numPlanes = original_parameters->Num_z_planes, Flag=0;
703  G4double chose, totArea=0., Achose1, Achose2,
704  rad1, rad2, sinphi1, sinphi2, cosphi1, cosphi2;
705  G4double a, b, l2, rang, totalPhi, ksi,
706  area, aTop=0., aBottom=0., zVal=0.;
707 
708  G4ThreeVector p0, p1, p2, p3;
709  std::vector<G4double> aVector1;
710  std::vector<G4double> aVector2;
711  std::vector<G4double> aVector3;
712 
713  totalPhi= (phiIsOpen) ? (endPhi-startPhi) : twopi;
714  ksi = totalPhi/numSide;
715  G4double cosksi = std::cos(ksi/2.);
716 
717  // Below we generate the areas relevant to our solid
718  //
719  for(j=0; j<numPlanes-1; j++)
720  {
721  a = original_parameters->Rmax[j+1];
722  b = original_parameters->Rmax[j];
724  -original_parameters->Z_values[j+1]) + sqr(b-a);
725  area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
726  aVector1.push_back(area);
727  }
728 
729  for(j=0; j<numPlanes-1; j++)
730  {
731  a = original_parameters->Rmin[j+1];//*cosksi;
732  b = original_parameters->Rmin[j];//*cosksi;
734  -original_parameters->Z_values[j+1]) + sqr(b-a);
735  area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
736  aVector2.push_back(area);
737  }
738 
739  for(j=0; j<numPlanes-1; j++)
740  {
741  if(phiIsOpen == true)
742  {
743  aVector3.push_back(0.5*(original_parameters->Rmax[j]
746  -original_parameters->Rmin[j+1])
747  *std::fabs(original_parameters->Z_values[j+1]
749  }
750  else { aVector3.push_back(0.); }
751  }
752 
753  for(j=0; j<numPlanes-1; j++)
754  {
755  totArea += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
756  }
757 
758  // Must include top and bottom areas
759  //
760  if(original_parameters->Rmax[numPlanes-1] != 0.)
761  {
762  a = original_parameters->Rmax[numPlanes-1];
763  b = original_parameters->Rmin[numPlanes-1];
764  l2 = sqr(a-b);
765  aTop = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
766  }
767 
768  if(original_parameters->Rmax[0] != 0.)
769  {
770  a = original_parameters->Rmax[0];
771  b = original_parameters->Rmin[0];
772  l2 = sqr(a-b);
773  aBottom = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
774  }
775 
776  Achose1 = 0.;
777  Achose2 = numSide*(aVector1[0]+aVector2[0])+2.*aVector3[0];
778 
779  chose = RandFlat::shoot(0.,totArea+aTop+aBottom);
780  if( (chose >= 0.) && (chose < aTop + aBottom) )
781  {
782  chose = RandFlat::shoot(startPhi,startPhi+totalPhi);
783  rang = std::floor((chose-startPhi)/ksi-0.01);
784  if(rang<0) { rang=0; }
785  rang = std::fabs(rang);
786  sinphi1 = std::sin(startPhi+rang*ksi);
787  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
788  cosphi1 = std::cos(startPhi+rang*ksi);
789  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
790  chose = RandFlat::shoot(0., aTop + aBottom);
791  if(chose>=0. && chose<aTop)
792  {
793  rad1 = original_parameters->Rmin[numPlanes-1];
794  rad2 = original_parameters->Rmax[numPlanes-1];
795  zVal = original_parameters->Z_values[numPlanes-1];
796  }
797  else
798  {
799  rad1 = original_parameters->Rmin[0];
800  rad2 = original_parameters->Rmax[0];
801  zVal = original_parameters->Z_values[0];
802  }
803  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
804  p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
805  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
806  p3 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
807  return GetPointOnPlane(p0,p1,p2,p3);
808  }
809  else
810  {
811  for (j=0; j<numPlanes-1; j++)
812  {
813  if( ((chose >= Achose1) && (chose < Achose2)) || (j == numPlanes-2) )
814  {
815  Flag = j; break;
816  }
817  Achose1 += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
818  Achose2 = Achose1 + numSide*(aVector1[j+1]+aVector2[j+1])
819  + 2.*aVector3[j+1];
820  }
821  }
822 
823  // At this point we have chosen a subsection
824  // between to adjacent plane cuts...
825 
826  j = Flag;
827 
828  totArea = numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
829  chose = RandFlat::shoot(0.,totArea);
830 
831  if( (chose>=0.) && (chose<numSide*aVector1[j]) )
832  {
833  chose = RandFlat::shoot(startPhi,startPhi+totalPhi);
834  rang = std::floor((chose-startPhi)/ksi-0.01);
835  if(rang<0) { rang=0; }
836  rang = std::fabs(rang);
837  rad1 = original_parameters->Rmax[j];
838  rad2 = original_parameters->Rmax[j+1];
839  sinphi1 = std::sin(startPhi+rang*ksi);
840  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
841  cosphi1 = std::cos(startPhi+rang*ksi);
842  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
843  zVal = original_parameters->Z_values[j];
844 
845  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
846  p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
847 
848  zVal = original_parameters->Z_values[j+1];
849 
850  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
851  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
852  return GetPointOnPlane(p0,p1,p2,p3);
853  }
854  else if ( (chose >= numSide*aVector1[j])
855  && (chose <= numSide*(aVector1[j]+aVector2[j])) )
856  {
857  chose = RandFlat::shoot(startPhi,startPhi+totalPhi);
858  rang = std::floor((chose-startPhi)/ksi-0.01);
859  if(rang<0) { rang=0; }
860  rang = std::fabs(rang);
861  rad1 = original_parameters->Rmin[j];
862  rad2 = original_parameters->Rmin[j+1];
863  sinphi1 = std::sin(startPhi+rang*ksi);
864  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
865  cosphi1 = std::cos(startPhi+rang*ksi);
866  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
867  zVal = original_parameters->Z_values[j];
868 
869  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
870  p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
871 
872  zVal = original_parameters->Z_values[j+1];
873 
874  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
875  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
876  return GetPointOnPlane(p0,p1,p2,p3);
877  }
878 
879  chose = RandFlat::shoot(0.,2.2);
880  if( (chose>=0.) && (chose < 1.) )
881  {
882  rang = startPhi;
883  }
884  else
885  {
886  rang = endPhi;
887  }
888 
889  cosphi1 = std::cos(rang); rad1 = original_parameters->Rmin[j];
890  sinphi1 = std::sin(rang); rad2 = original_parameters->Rmax[j];
891 
892  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
894  p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
896 
897  rad1 = original_parameters->Rmax[j+1];
898  rad2 = original_parameters->Rmin[j+1];
899 
900  p2 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
902  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
904  return GetPointOnPlane(p0,p1,p2,p3);
905  }
906  else // Generic polyhedra
907  {
908  return GetPointOnSurfaceGeneric();
909  }
910 }
911 
912 //
913 // CreatePolyhedron
914 //
916 {
917  if (!genericPgon)
918  {
926  }
927  else
928  {
929  // The following code prepares for:
930  // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
931  // const double xyz[][3],
932  // const int faces_vec[][4])
933  // Here is an extract from the header file HepPolyhedron.h:
951  G4int nNodes;
952  G4int nFaces;
953  typedef G4double double3[3];
954  double3* xyz;
955  typedef G4int int4[4];
956  int4* faces_vec;
957  if (phiIsOpen)
958  {
959  // Triangulate open ends. Simple ear-chopping algorithm...
960  // I'm not sure how robust this algorithm is (J.Allison).
961  //
962  std::vector<G4bool> chopped(numCorner, false);
963  std::vector<G4int*> triQuads;
964  G4int remaining = numCorner;
965  G4int iStarter = 0;
966  while (remaining >= 3)
967  {
968  // Find unchopped corners...
969  //
970  G4int A = -1, B = -1, C = -1;
971  G4int iStepper = iStarter;
972  do
973  {
974  if (A < 0) { A = iStepper; }
975  else if (B < 0) { B = iStepper; }
976  else if (C < 0) { C = iStepper; }
977  do
978  {
979  if (++iStepper >= numCorner) iStepper = 0;
980  }
981  while (chopped[iStepper]);
982  }
983  while (C < 0 && iStepper != iStarter);
984 
985  // Check triangle at B is pointing outward (an "ear").
986  // Sign of z cross product determines...
987 
988  G4double BAr = corners[A].r - corners[B].r;
989  G4double BAz = corners[A].z - corners[B].z;
990  G4double BCr = corners[C].r - corners[B].r;
991  G4double BCz = corners[C].z - corners[B].z;
992  if (BAr * BCz - BAz * BCr < kCarTolerance)
993  {
994  G4int* tq = new G4int[3];
995  tq[0] = A + 1;
996  tq[1] = B + 1;
997  tq[2] = C + 1;
998  triQuads.push_back(tq);
999  chopped[B] = true;
1000  --remaining;
1001  }
1002  else
1003  {
1004  do
1005  {
1006  if (++iStarter >= numCorner) { iStarter = 0; }
1007  }
1008  while (chopped[iStarter]);
1009  }
1010  }
1011 
1012  // Transfer to faces...
1013 
1014  nNodes = (numSide + 1) * numCorner;
1015  nFaces = numSide * numCorner + 2 * triQuads.size();
1016  faces_vec = new int4[nFaces];
1017  G4int iface = 0;
1018  G4int addition = numCorner * numSide;
1019  G4int d = numCorner - 1;
1020  for (G4int iEnd = 0; iEnd < 2; ++iEnd)
1021  {
1022  for (size_t i = 0; i < triQuads.size(); ++i)
1023  {
1024  // Negative for soft/auxiliary/normally invisible edges...
1025  //
1026  G4int a, b, c;
1027  if (iEnd == 0)
1028  {
1029  a = triQuads[i][0];
1030  b = triQuads[i][1];
1031  c = triQuads[i][2];
1032  }
1033  else
1034  {
1035  a = triQuads[i][0] + addition;
1036  b = triQuads[i][2] + addition;
1037  c = triQuads[i][1] + addition;
1038  }
1039  G4int ab = std::abs(b - a);
1040  G4int bc = std::abs(c - b);
1041  G4int ca = std::abs(a - c);
1042  faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
1043  faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
1044  faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
1045  faces_vec[iface][3] = 0;
1046  ++iface;
1047  }
1048  }
1049 
1050  // Continue with sides...
1051 
1052  xyz = new double3[nNodes];
1053  const G4double dPhi = (endPhi - startPhi) / numSide;
1054  G4double phi = startPhi;
1055  G4int ixyz = 0;
1056  for (G4int iSide = 0; iSide < numSide; ++iSide)
1057  {
1058  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1059  {
1060  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1061  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1062  xyz[ixyz][2] = corners[iCorner].z;
1063  if (iCorner < numCorner - 1)
1064  {
1065  faces_vec[iface][0] = ixyz + 1;
1066  faces_vec[iface][1] = ixyz + numCorner + 1;
1067  faces_vec[iface][2] = ixyz + numCorner + 2;
1068  faces_vec[iface][3] = ixyz + 2;
1069  }
1070  else
1071  {
1072  faces_vec[iface][0] = ixyz + 1;
1073  faces_vec[iface][1] = ixyz + numCorner + 1;
1074  faces_vec[iface][2] = ixyz + 2;
1075  faces_vec[iface][3] = ixyz - numCorner + 2;
1076  }
1077  ++iface;
1078  ++ixyz;
1079  }
1080  phi += dPhi;
1081  }
1082 
1083  // Last corners...
1084 
1085  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1086  {
1087  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1088  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1089  xyz[ixyz][2] = corners[iCorner].z;
1090  ++ixyz;
1091  }
1092  }
1093  else // !phiIsOpen - i.e., a complete 360 degrees.
1094  {
1095  nNodes = numSide * numCorner;
1096  nFaces = numSide * numCorner;;
1097  xyz = new double3[nNodes];
1098  faces_vec = new int4[nFaces];
1099  // const G4double dPhi = (endPhi - startPhi) / numSide;
1100  const G4double dPhi = twopi / numSide; // !phiIsOpen endPhi-startPhi = 360 degrees.
1101  G4double phi = startPhi;
1102  G4int ixyz = 0, iface = 0;
1103  for (G4int iSide = 0; iSide < numSide; ++iSide)
1104  {
1105  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1106  {
1107  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1108  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1109  xyz[ixyz][2] = corners[iCorner].z;
1110  if (iSide < numSide - 1)
1111  {
1112  if (iCorner < numCorner - 1)
1113  {
1114  faces_vec[iface][0] = ixyz + 1;
1115  faces_vec[iface][1] = ixyz + numCorner + 1;
1116  faces_vec[iface][2] = ixyz + numCorner + 2;
1117  faces_vec[iface][3] = ixyz + 2;
1118  }
1119  else
1120  {
1121  faces_vec[iface][0] = ixyz + 1;
1122  faces_vec[iface][1] = ixyz + numCorner + 1;
1123  faces_vec[iface][2] = ixyz + 2;
1124  faces_vec[iface][3] = ixyz - numCorner + 2;
1125  }
1126  }
1127  else // Last side joins ends...
1128  {
1129  if (iCorner < numCorner - 1)
1130  {
1131  faces_vec[iface][0] = ixyz + 1;
1132  faces_vec[iface][1] = ixyz + numCorner - nFaces + 1;
1133  faces_vec[iface][2] = ixyz + numCorner - nFaces + 2;
1134  faces_vec[iface][3] = ixyz + 2;
1135  }
1136  else
1137  {
1138  faces_vec[iface][0] = ixyz + 1;
1139  faces_vec[iface][1] = ixyz - nFaces + numCorner + 1;
1140  faces_vec[iface][2] = ixyz - nFaces + 2;
1141  faces_vec[iface][3] = ixyz - numCorner + 2;
1142  }
1143  }
1144  ++ixyz;
1145  ++iface;
1146  }
1147  phi += dPhi;
1148  }
1149  }
1150  G4Polyhedron* polyhedron = new G4Polyhedron;
1151  G4int problem = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec);
1152  delete [] faces_vec;
1153  delete [] xyz;
1154  if (problem)
1155  {
1156  std::ostringstream message;
1157  message << "Problem creating G4Polyhedron for: " << GetName();
1158  G4Exception("G4Polyhedra::CreatePolyhedron()", "GeomSolids1002",
1159  JustWarning, message);
1160  delete polyhedron;
1161  return 0;
1162  }
1163  else
1164  {
1165  return polyhedron;
1166  }
1167  }
1168 }
1169 
1170 
1172 {
1173  G4int numPlanes = (G4int)numCorner;
1174  G4bool isConvertible=true;
1175  G4double Zmax=rz->Bmax();
1176  rz->StartWithZMin();
1177 
1178  // Prepare vectors for storage
1179  //
1180  std::vector<G4double> Z;
1181  std::vector<G4double> Rmin;
1182  std::vector<G4double> Rmax;
1183 
1184  G4int countPlanes=1;
1185  G4int icurr=0;
1186  G4int icurl=0;
1187 
1188  // first plane Z=Z[0]
1189  //
1190  Z.push_back(corners[0].z);
1191  G4double Zprev=Z[0];
1192  if (Zprev == corners[1].z)
1193  {
1194  Rmin.push_back(corners[0].r);
1195  Rmax.push_back (corners[1].r);icurr=1;
1196  }
1197  else if (Zprev == corners[numPlanes-1].z)
1198  {
1199  Rmin.push_back(corners[numPlanes-1].r);
1200  Rmax.push_back (corners[0].r);
1201  icurl=numPlanes-1;
1202  }
1203  else
1204  {
1205  Rmin.push_back(corners[0].r);
1206  Rmax.push_back (corners[0].r);
1207  }
1208 
1209  // next planes until last
1210  //
1211  G4int inextr=0, inextl=0;
1212  for (G4int i=0; i < numPlanes-2; i++)
1213  {
1214  inextr=1+icurr;
1215  inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1216 
1217  if((corners[inextr].z >= Zmax) & (corners[inextl].z >= Zmax)) { break; }
1218 
1219  G4double Zleft = corners[inextl].z;
1220  G4double Zright = corners[inextr].z;
1221  if(Zright>Zleft)
1222  {
1223  Z.push_back(Zleft);
1224  countPlanes++;
1225  G4double difZr=corners[inextr].z - corners[icurr].z;
1226  G4double difZl=corners[inextl].z - corners[icurl].z;
1227 
1228  if(std::fabs(difZl) < kCarTolerance)
1229  {
1230  if(corners[inextl].r >= corners[icurl].r)
1231  {
1232  Rmin.push_back(corners[icurl].r);
1233  Rmax.push_back(Rmax[countPlanes-2]);
1234  Rmax[countPlanes-2]=corners[icurl].r;
1235  }
1236  else
1237  {
1238  Rmin.push_back(corners[inextl].r);
1239  Rmax.push_back(corners[icurl].r);
1240  }
1241  }
1242  else if (difZl >= kCarTolerance)
1243  {
1244  Rmin.push_back(corners[inextl].r);
1245  Rmax.push_back (corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1246  *(corners[inextr].r - corners[icurr].r));
1247  }
1248  else
1249  {
1250  isConvertible=false; break;
1251  }
1252  icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1253  }
1254  else if(std::fabs(Zright-Zleft)<kCarTolerance) // Zright=Zleft
1255  {
1256  Z.push_back(Zleft);
1257  countPlanes++;
1258  icurr++;
1259 
1260  icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1261 
1262  Rmin.push_back(corners[inextl].r);
1263  Rmax.push_back (corners[inextr].r);
1264  }
1265  else // Zright<Zleft
1266  {
1267  Z.push_back(Zright);
1268  countPlanes++;
1269 
1270  G4double difZr=corners[inextr].z - corners[icurr].z;
1271  G4double difZl=corners[inextl].z - corners[icurl].z;
1272  if(std::fabs(difZr) < kCarTolerance)
1273  {
1274  if(corners[inextr].r >= corners[icurr].r)
1275  {
1276  Rmin.push_back(corners[icurr].r);
1277  Rmax.push_back(corners[inextr].r);
1278  }
1279  else
1280  {
1281  Rmin.push_back(corners[inextr].r);
1282  Rmax.push_back(corners[icurr].r);
1283  Rmax[countPlanes-2]=corners[inextr].r;
1284  }
1285  icurr++;
1286  } // plate
1287  else if (difZr >= kCarTolerance)
1288  {
1289  if(std::fabs(difZl)<kCarTolerance)
1290  {
1291  Rmax.push_back(corners[inextr].r);
1292  Rmin.push_back (corners[icurr].r);
1293  }
1294  else
1295  {
1296  Rmax.push_back(corners[inextr].r);
1297  Rmin.push_back (corners[icurl].r+(Zright-corners[icurl].z)/difZl
1298  * (corners[inextl].r - corners[icurl].r));
1299  }
1300  icurr++;
1301  }
1302  else
1303  {
1304  isConvertible=false; break;
1305  }
1306  }
1307  } // end for loop
1308 
1309  // last plane Z=Zmax
1310  //
1311  Z.push_back(Zmax);
1312  countPlanes++;
1313  inextr=1+icurr;
1314  inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1315 
1316  if(corners[inextr].z==corners[inextl].z)
1317  {
1318  Rmax.push_back(corners[inextr].r);
1319  Rmin.push_back(corners[inextl].r);
1320  }
1321  else
1322  {
1323  Rmax.push_back(corners[inextr].r);
1324  Rmin.push_back(corners[inextl].r);
1325  }
1326 
1327  // Set original parameters Rmin,Rmax,Z
1328  //
1329  if(isConvertible)
1330  {
1333  original_parameters->Z_values = new G4double[countPlanes];
1334  original_parameters->Rmin = new G4double[countPlanes];
1335  original_parameters->Rmax = new G4double[countPlanes];
1336 
1337  for(G4int j=0; j < countPlanes; j++)
1338  {
1339  original_parameters->Z_values[j] = Z[j];
1340  original_parameters->Rmax[j] = Rmax[j];
1341  original_parameters->Rmin[j] = Rmin[j];
1342  }
1345  original_parameters->Num_z_planes = countPlanes;
1346 
1347  }
1348  else // Set parameters(r,z) with Rmin==0 as convention
1349  {
1350 #ifdef G4SPECSDEBUG
1351  std::ostringstream message;
1352  message << "Polyhedra " << GetName() << G4endl
1353  << "cannot be converted to Polyhedra with (Rmin,Rmaz,Z) parameters!";
1354  G4Exception("G4Polyhedra::SetOriginalParameters()",
1355  "GeomSolids0002", JustWarning, message);
1356 #endif
1359  original_parameters->Z_values = new G4double[numPlanes];
1360  original_parameters->Rmin = new G4double[numPlanes];
1361  original_parameters->Rmax = new G4double[numPlanes];
1362 
1363  for(G4int j=0; j < numPlanes; j++)
1364  {
1366  original_parameters->Rmax[j] = corners[j].r;
1367  original_parameters->Rmin[j] = 0.0;
1368  }
1371  original_parameters->Num_z_planes = numPlanes;
1372  }
1373  //return isConvertible;
1374 }
1375 
1376 #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
G4Polyhedron * fpPolyhedron
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
G4GeometryType GetEntityType() const
Definition: G4Polyhedra.cc:578
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:649
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:477
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:698
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:567
G4int NumVertices() const
G4double Bmax() const
void ScaleA(G4double scale)
bool G4bool
Definition: G4Types.hh:79
G4double startPhi
Definition: G4Polyhedra.hh:181
#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:539
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
EInside Inside(const G4ThreeVector &p) const
Definition: G4Polyhedra.cc:519
static const G4double ab
G4Polyhedron * CreatePolyhedron() const
Definition: G4Polyhedra.cc:915
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:681
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:596
G4VSolid * Clone() const
Definition: G4Polyhedra.cc:587
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
virtual G4Polyhedron * GetPolyhedron() const