ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Polyhedra.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4Polyhedra.cc
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 // Implementation of G4Polyhedra, a CSG polyhedra,
27 // as an inherited class of G4VCSGfaceted.
28 //
29 // Utility classes:
30 // * G4EnclosingCylinder: decided a quick check of geometry would be a
31 // good idea (for CPU speed). If the quick check fails, the regular
32 // full-blown G4VCSGfaceted version is invoked.
33 // * G4ReduciblePolygon: Really meant as a check of input parameters,
34 // this utility class also "converts" the GEANT3-like PGON/PCON
35 // arguments into the newer ones.
36 // Both these classes are implemented outside this file because they are
37 // shared with G4Polycone.
38 //
39 // Author: David C. Williams (davidw@scipp.ucsc.edu)
40 // --------------------------------------------------------------------
41 
42 #include "G4Polyhedra.hh"
43 
44 #if !defined(G4GEOM_USE_UPOLYHEDRA)
45 
46 #include "G4PolyhedraSide.hh"
47 #include "G4PolyPhiFace.hh"
48 
49 #include "G4GeomTools.hh"
50 #include "G4VoxelLimits.hh"
51 #include "G4AffineTransform.hh"
52 #include "G4BoundingEnvelope.hh"
53 
54 #include "Randomize.hh"
55 
56 #include "G4EnclosingCylinder.hh"
57 #include "G4ReduciblePolygon.hh"
58 #include "G4VPVParameterisation.hh"
59 
60 #include <sstream>
61 
62 using namespace CLHEP;
63 
64 // Constructor (GEANT3 style parameters)
65 //
66 // GEANT3 PGON radii are specified in the distance to the norm of each face.
67 //
69  G4double phiStart,
70  G4double thePhiTotal,
71  G4int theNumSide,
72  G4int numZPlanes,
73  const G4double zPlane[],
74  const G4double rInner[],
75  const G4double rOuter[] )
76  : G4VCSGfaceted( name )
77 {
78  if (theNumSide <= 0)
79  {
80  std::ostringstream message;
81  message << "Solid must have at least one side - " << GetName() << G4endl
82  << " No sides specified !";
83  G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
84  FatalErrorInArgument, message);
85  }
86 
87  //
88  // Calculate conversion factor from G3 radius to G4 radius
89  //
90  G4double phiTotal = thePhiTotal;
91  if ( (phiTotal <=0) || (phiTotal >= twopi*(1-DBL_EPSILON)) )
92  { phiTotal = twopi; }
93  G4double convertRad = std::cos(0.5*phiTotal/theNumSide);
94 
95  //
96  // Some historical stuff
97  //
99 
100  original_parameters->numSide = theNumSide;
101  original_parameters->Start_angle = phiStart;
103  original_parameters->Num_z_planes = numZPlanes;
104  original_parameters->Z_values = new G4double[numZPlanes];
105  original_parameters->Rmin = new G4double[numZPlanes];
106  original_parameters->Rmax = new G4double[numZPlanes];
107 
108  for (G4int i=0; i<numZPlanes; ++i)
109  {
110  if (( i < numZPlanes-1) && ( zPlane[i] == zPlane[i+1] ))
111  {
112  if( (rInner[i] > rOuter[i+1])
113  ||(rInner[i+1] > rOuter[i]) )
114  {
115  DumpInfo();
116  std::ostringstream message;
117  message << "Cannot create a Polyhedra with no contiguous segments."
118  << G4endl
119  << " Segments are not contiguous !" << G4endl
120  << " rMin[" << i << "] = " << rInner[i]
121  << " -- rMax[" << i+1 << "] = " << rOuter[i+1] << G4endl
122  << " rMin[" << i+1 << "] = " << rInner[i+1]
123  << " -- rMax[" << i << "] = " << rOuter[i];
124  G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
125  FatalErrorInArgument, message);
126  }
127  }
128  original_parameters->Z_values[i] = zPlane[i];
129  original_parameters->Rmin[i] = rInner[i]/convertRad;
130  original_parameters->Rmax[i] = rOuter[i]/convertRad;
131  }
132 
133 
134  //
135  // Build RZ polygon using special PCON/PGON GEANT3 constructor
136  //
137  G4ReduciblePolygon* rz =
138  new G4ReduciblePolygon( rInner, rOuter, zPlane, numZPlanes );
139  rz->ScaleA( 1/convertRad );
140 
141  //
142  // Do the real work
143  //
144  Create( phiStart, phiTotal, theNumSide, rz );
145 
146  delete rz;
147 }
148 
149 // Constructor (generic parameters)
150 //
152  G4double phiStart,
153  G4double phiTotal,
154  G4int theNumSide,
155  G4int numRZ,
156  const G4double r[],
157  const G4double z[] )
158  : G4VCSGfaceted( name ), genericPgon(true)
159 {
160  if (theNumSide <= 0)
161  {
162  std::ostringstream message;
163  message << "Solid must have at least one side - " << GetName() << G4endl
164  << " No sides specified !";
165  G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
166  FatalErrorInArgument, message);
167  }
168 
169  G4ReduciblePolygon* rz = new G4ReduciblePolygon( r, z, numRZ );
170 
171  Create( phiStart, phiTotal, theNumSide, rz );
172 
173  // Set original_parameters struct for consistency
174  //
176 
177  delete rz;
178 }
179 
180 // Create
181 //
182 // Generic create routine, called by each constructor
183 // after conversion of arguments
184 //
186  G4double phiTotal,
187  G4int theNumSide,
188  G4ReduciblePolygon* rz )
189 {
190  //
191  // Perform checks of rz values
192  //
193  if (rz->Amin() < 0.0)
194  {
195  std::ostringstream message;
196  message << "Illegal input parameters - " << GetName() << G4endl
197  << " All R values must be >= 0 !";
198  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
199  FatalErrorInArgument, message);
200  }
201 
202  G4double rzArea = rz->Area();
203  if (rzArea < -kCarTolerance)
204  {
205  rz->ReverseOrder();
206  }
207  else if (rzArea < kCarTolerance)
208  {
209  std::ostringstream message;
210  message << "Illegal input parameters - " << GetName() << G4endl
211  << " R/Z cross section is zero or near zero: " << rzArea;
212  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
213  FatalErrorInArgument, message);
214  }
215 
217  || (!rz->RemoveRedundantVertices( kCarTolerance )) )
218  {
219  std::ostringstream message;
220  message << "Illegal input parameters - " << GetName() << G4endl
221  << " Too few unique R/Z values !";
222  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
223  FatalErrorInArgument, message);
224  }
225 
226  if (rz->CrossesItself( 1/kInfinity ))
227  {
228  std::ostringstream message;
229  message << "Illegal input parameters - " << GetName() << G4endl
230  << " R/Z segments cross !";
231  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
232  FatalErrorInArgument, message);
233  }
234 
235  numCorner = rz->NumVertices();
236 
237 
238  startPhi = phiStart;
239  while( startPhi < 0 ) // Loop checking, 13.08.2015, G.Cosmo
240  startPhi += twopi;
241  //
242  // Phi opening? Account for some possible roundoff, and interpret
243  // nonsense value as representing no phi opening
244  //
245  if ( (phiTotal <= 0) || (phiTotal > twopi*(1-DBL_EPSILON)) )
246  {
247  phiIsOpen = false;
248  endPhi = phiStart+twopi;
249  }
250  else
251  {
252  phiIsOpen = true;
253 
254  //
255  // Convert phi into our convention
256  //
257  endPhi = phiStart+phiTotal;
258  while( endPhi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
259  endPhi += twopi;
260  }
261 
262  //
263  // Save number sides
264  //
265  numSide = theNumSide;
266 
267  //
268  // Allocate corner array.
269  //
271 
272  //
273  // Copy corners
274  //
275  G4ReduciblePolygonIterator iterRZ(rz);
276 
277  G4PolyhedraSideRZ *next = corners;
278  iterRZ.Begin();
279  do // Loop checking, 13.08.2015, G.Cosmo
280  {
281  next->r = iterRZ.GetA();
282  next->z = iterRZ.GetB();
283  } while( ++next, iterRZ.Next() );
284 
285  //
286  // Allocate face pointer array
287  //
289  faces = new G4VCSGface*[numFace];
290 
291  //
292  // Construct side faces
293  //
294  // To do so properly, we need to keep track of four successive RZ
295  // corners.
296  //
297  // But! Don't construct a face if both points are at zero radius!
298  //
299  G4PolyhedraSideRZ* corner = corners,
300  * prev = corners + numCorner-1,
301  * nextNext;
302  G4VCSGface** face = faces;
303  do // Loop checking, 13.08.2015, G.Cosmo
304  {
305  next = corner+1;
306  if (next >= corners+numCorner) next = corners;
307  nextNext = next+1;
308  if (nextNext >= corners+numCorner) nextNext = corners;
309 
310  if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
311 /*
312  // We must decide here if we can dare declare one of our faces
313  // as having a "valid" normal (i.e. allBehind = true). This
314  // is never possible if the face faces "inward" in r *unless*
315  // we have only one side
316  //
317  G4bool allBehind;
318  if ((corner->z > next->z) && (numSide > 1))
319  {
320  allBehind = false;
321  }
322  else
323  {
324  //
325  // Otherwise, it is only true if the line passing
326  // through the two points of the segment do not
327  // split the r/z cross section
328  //
329  allBehind = !rz->BisectedBy( corner->r, corner->z,
330  next->r, next->z, kCarTolerance );
331  }
332 */
333  *face++ = new G4PolyhedraSide( prev, corner, next, nextNext,
335  } while( prev=corner, corner=next, corner > corners );
336 
337  if (phiIsOpen)
338  {
339  //
340  // Construct phi open edges
341  //
342  *face++ = new G4PolyPhiFace( rz, startPhi, phiTotal/numSide, endPhi );
343  *face++ = new G4PolyPhiFace( rz, endPhi, phiTotal/numSide, startPhi );
344  }
345 
346  //
347  // We might have dropped a face or two: recalculate numFace
348  //
349  numFace = face-faces;
350 
351  //
352  // Make enclosingCylinder
353  //
355  new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
356 }
357 
358 // Fake default constructor - sets only member data and allocates memory
359 // for usage restricted to object persistency.
360 //
362  : G4VCSGfaceted(a), startPhi(0.), endPhi(0.)
363 {
364 }
365 
366 // Destructor
367 //
369 {
370  delete [] corners;
371  delete original_parameters;
372  delete enclosingCylinder;
373 }
374 
375 // Copy constructor
376 //
378  : G4VCSGfaceted( source )
379 {
380  CopyStuff( source );
381 }
382 
383 // Assignment operator
384 //
386 {
387  if (this == &source) return *this;
388 
389  G4VCSGfaceted::operator=( source );
390 
391  delete [] corners;
392  delete original_parameters;
393  delete enclosingCylinder;
394 
395  CopyStuff( source );
396 
397  return *this;
398 }
399 
400 // CopyStuff
401 //
402 void G4Polyhedra::CopyStuff( const G4Polyhedra& source )
403 {
404  //
405  // Simple stuff
406  //
407  numSide = source.numSide;
408  startPhi = source.startPhi;
409  endPhi = source.endPhi;
410  phiIsOpen = source.phiIsOpen;
411  numCorner = source.numCorner;
412  genericPgon= source.genericPgon;
413 
414  //
415  // The corner array
416  //
418 
419  G4PolyhedraSideRZ* corn = corners,
420  * sourceCorn = source.corners;
421  do // Loop checking, 13.08.2015, G.Cosmo
422  {
423  *corn = *sourceCorn;
424  } while( ++sourceCorn, ++corn < corners+numCorner );
425 
426  //
427  // Original parameters
428  //
429  if (source.original_parameters)
430  {
433  }
434 
435  //
436  // Enclosing cylinder
437  //
439 
440  fRebuildPolyhedron = false;
441  fpPolyhedron = nullptr;
442 }
443 
444 // Reset
445 //
446 // Recalculates and reshapes the solid, given pre-assigned scaled
447 // original_parameters.
448 //
450 {
451  if (genericPgon)
452  {
453  std::ostringstream message;
454  message << "Solid " << GetName() << " built using generic construct."
455  << G4endl << "Not applicable to the generic construct !";
456  G4Exception("G4Polyhedra::Reset()", "GeomSolids1001",
457  JustWarning, message, "Parameters NOT resetted.");
458  return true;
459  }
460 
461  //
462  // Clear old setup
463  //
465  delete [] corners;
466  delete enclosingCylinder;
467 
468  //
469  // Rebuild polyhedra
470  //
471  G4ReduciblePolygon* rz =
479  delete rz;
480 
481  return false;
482 }
483 
484 // Inside
485 //
486 // This is an override of G4VCSGfaceted::Inside, created in order
487 // to speed things up by first checking with G4EnclosingCylinder.
488 //
490 {
491  //
492  // Quick test
493  //
494  if (enclosingCylinder->MustBeOutside(p)) return kOutside;
495 
496  //
497  // Long answer
498  //
499  return G4VCSGfaceted::Inside(p);
500 }
501 
502 // DistanceToIn
503 //
504 // This is an override of G4VCSGfaceted::Inside, created in order
505 // to speed things up by first checking with G4EnclosingCylinder.
506 //
508  const G4ThreeVector& v ) const
509 {
510  //
511  // Quick test
512  //
513  if (enclosingCylinder->ShouldMiss(p,v))
514  return kInfinity;
515 
516  //
517  // Long answer
518  //
519  return G4VCSGfaceted::DistanceToIn( p, v );
520 }
521 
522 // DistanceToIn
523 //
525 {
526  return G4VCSGfaceted::DistanceToIn(p);
527 }
528 
529 // Get bounding box
530 //
532  G4ThreeVector& pMax) const
533 {
534  G4double rmin = kInfinity, rmax = -kInfinity;
535  G4double zmin = kInfinity, zmax = -kInfinity;
536  for (G4int i=0; i<GetNumRZCorner(); ++i)
537  {
538  G4PolyhedraSideRZ corner = GetCorner(i);
539  if (corner.r < rmin) rmin = corner.r;
540  if (corner.r > rmax) rmax = corner.r;
541  if (corner.z < zmin) zmin = corner.z;
542  if (corner.z > zmax) zmax = corner.z;
543  }
544 
545  G4double sphi = GetStartPhi();
546  G4double ephi = GetEndPhi();
547  G4double dphi = IsOpen() ? ephi-sphi : twopi;
548  G4int ksteps = GetNumSide();
549  G4double astep = dphi/ksteps;
550  G4double sinStep = std::sin(astep);
551  G4double cosStep = std::cos(astep);
552 
553  G4double sinCur = GetSinStartPhi();
554  G4double cosCur = GetCosStartPhi();
555  if (!IsOpen()) rmin = 0.;
556  G4double xmin = rmin*cosCur, xmax = xmin;
557  G4double ymin = rmin*sinCur, ymax = ymin;
558  for (G4int k=0; k<ksteps+1; ++k)
559  {
560  G4double x = rmax*cosCur;
561  if (x < xmin) xmin = x;
562  if (x > xmax) xmax = x;
563  G4double y = rmax*sinCur;
564  if (y < ymin) ymin = y;
565  if (y > ymax) ymax = y;
566  if (rmin > 0)
567  {
568  G4double xx = rmin*cosCur;
569  if (xx < xmin) xmin = xx;
570  if (xx > xmax) xmax = xx;
571  G4double yy = rmin*sinCur;
572  if (yy < ymin) ymin = yy;
573  if (yy > ymax) ymax = yy;
574  }
575  G4double sinTmp = sinCur;
576  sinCur = sinCur*cosStep + cosCur*sinStep;
577  cosCur = cosCur*cosStep - sinTmp*sinStep;
578  }
579  pMin.set(xmin,ymin,zmin);
580  pMax.set(xmax,ymax,zmax);
581 
582  // Check correctness of the bounding box
583  //
584  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
585  {
586  std::ostringstream message;
587  message << "Bad bounding box (min >= max) for solid: "
588  << GetName() << " !"
589  << "\npMin = " << pMin
590  << "\npMax = " << pMax;
591  G4Exception("G4Polyhedra::BoundingLimits()", "GeomMgt0001",
592  JustWarning, message);
593  DumpInfo();
594  }
595 }
596 
597 // Calculate extent under transform and specified limit
598 //
600  const G4VoxelLimits& pVoxelLimit,
601  const G4AffineTransform& pTransform,
602  G4double& pMin, G4double& pMax) const
603 {
604  G4ThreeVector bmin, bmax;
605  G4bool exist;
606 
607  // Check bounding box (bbox)
608  //
609  BoundingLimits(bmin,bmax);
610  G4BoundingEnvelope bbox(bmin,bmax);
611 #ifdef G4BBOX_EXTENT
612  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
613 #endif
614  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
615  {
616  return exist = (pMin < pMax) ? true : false;
617  }
618 
619  // To find the extent, RZ contour of the polycone is subdivided
620  // in triangles. The extent is calculated as cumulative extent of
621  // all sub-polycones formed by rotation of triangles around Z
622  //
623  G4TwoVectorList contourRZ;
624  G4TwoVectorList triangles;
625  std::vector<G4int> iout;
626  G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
627  G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
628 
629  // get RZ contour, ensure anticlockwise order of corners
630  for (G4int i=0; i<GetNumRZCorner(); ++i)
631  {
632  G4PolyhedraSideRZ corner = GetCorner(i);
633  contourRZ.push_back(G4TwoVector(corner.r,corner.z));
634  }
636  G4double area = G4GeomTools::PolygonArea(contourRZ);
637  if (area < 0.) std::reverse(contourRZ.begin(),contourRZ.end());
638 
639  // triangulate RZ countour
640  if (!G4GeomTools::TriangulatePolygon(contourRZ,triangles))
641  {
642  std::ostringstream message;
643  message << "Triangulation of RZ contour has failed for solid: "
644  << GetName() << " !"
645  << "\nExtent has been calculated using boundary box";
646  G4Exception("G4Polyhedra::CalculateExtent()",
647  "GeomMgt1002",JustWarning,message);
648  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
649  }
650 
651  // set trigonometric values
652  G4double sphi = GetStartPhi();
653  G4double ephi = GetEndPhi();
654  G4double dphi = IsOpen() ? ephi-sphi : twopi;
655  G4int ksteps = GetNumSide();
656  G4double astep = dphi/ksteps;
657  G4double sinStep = std::sin(astep);
658  G4double cosStep = std::cos(astep);
659  G4double sinStart = GetSinStartPhi();
660  G4double cosStart = GetCosStartPhi();
661 
662  // allocate vector lists
663  std::vector<const G4ThreeVectorList *> polygons;
664  polygons.resize(ksteps+1);
665  for (G4int k=0; k<ksteps+1; ++k)
666  {
667  polygons[k] = new G4ThreeVectorList(3);
668  }
669 
670  // main loop along triangles
671  pMin = kInfinity;
672  pMax = -kInfinity;
673  G4int ntria = triangles.size()/3;
674  for (G4int i=0; i<ntria; ++i)
675  {
676  G4double sinCur = sinStart;
677  G4double cosCur = cosStart;
678  G4int i3 = i*3;
679  for (G4int k=0; k<ksteps+1; ++k) // rotate triangle
680  {
681  G4ThreeVectorList* ptr = const_cast<G4ThreeVectorList*>(polygons[k]);
682  G4ThreeVectorList::iterator iter = ptr->begin();
683  iter->set(triangles[i3+0].x()*cosCur,
684  triangles[i3+0].x()*sinCur,
685  triangles[i3+0].y());
686  iter++;
687  iter->set(triangles[i3+1].x()*cosCur,
688  triangles[i3+1].x()*sinCur,
689  triangles[i3+1].y());
690  iter++;
691  iter->set(triangles[i3+2].x()*cosCur,
692  triangles[i3+2].x()*sinCur,
693  triangles[i3+2].y());
694 
695  G4double sinTmp = sinCur;
696  sinCur = sinCur*cosStep + cosCur*sinStep;
697  cosCur = cosCur*cosStep - sinTmp*sinStep;
698  }
699 
700  // set sub-envelope and adjust extent
702  G4BoundingEnvelope benv(polygons);
703  if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
704  if (emin < pMin) pMin = emin;
705  if (emax > pMax) pMax = emax;
706  if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
707  }
708  // free memory
709  for (G4int k=0; k<ksteps+1; ++k) { delete polygons[k]; polygons[k]=0;}
710  return (pMin < pMax);
711 }
712 
713 // ComputeDimensions
714 //
716  const G4int n,
717  const G4VPhysicalVolume* pRep )
718 {
719  p->ComputeDimensions(*this,n,pRep);
720 }
721 
722 // GetEntityType
723 //
725 {
726  return G4String("G4Polyhedra");
727 }
728 
729 // Make a clone of the object
730 //
732 {
733  return new G4Polyhedra(*this);
734 }
735 
736 // Stream object contents to an output stream
737 //
738 std::ostream& G4Polyhedra::StreamInfo( std::ostream& os ) const
739 {
740  G4int oldprc = os.precision(16);
741  os << "-----------------------------------------------------------\n"
742  << " *** Dump for solid - " << GetName() << " ***\n"
743  << " ===================================================\n"
744  << " Solid type: G4Polyhedra\n"
745  << " Parameters: \n"
746  << " starting phi angle : " << startPhi/degree << " degrees \n"
747  << " ending phi angle : " << endPhi/degree << " degrees \n"
748  << " number of sides : " << numSide << " \n";
749  G4int i=0;
750  if (!genericPgon)
751  {
753  os << " number of Z planes: " << numPlanes << "\n"
754  << " Z values: \n";
755  for (i=0; i<numPlanes; ++i)
756  {
757  os << " Z plane " << i << ": "
758  << original_parameters->Z_values[i] << "\n";
759  }
760  os << " Tangent distances to inner surface (Rmin): \n";
761  for (i=0; i<numPlanes; ++i)
762  {
763  os << " Z plane " << i << ": "
764  << original_parameters->Rmin[i] << "\n";
765  }
766  os << " Tangent distances to outer surface (Rmax): \n";
767  for (i=0; i<numPlanes; ++i)
768  {
769  os << " Z plane " << i << ": "
770  << original_parameters->Rmax[i] << "\n";
771  }
772  }
773  os << " number of RZ points: " << numCorner << "\n"
774  << " RZ values (corners): \n";
775  for (i=0; i<numCorner; ++i)
776  {
777  os << " "
778  << corners[i].r << ", " << corners[i].z << "\n";
779  }
780  os << "-----------------------------------------------------------\n";
781  os.precision(oldprc);
782 
783  return os;
784 }
785 
786 // GetPointOnPlane
787 //
788 // Auxiliary method for get point on surface
789 //
792  G4ThreeVector p2, G4ThreeVector p3) const
793 {
794  G4double lambda1, lambda2, chose,aOne,aTwo;
795  G4ThreeVector t, u, v, w, Area, normal;
796  aOne = 1.;
797  aTwo = 1.;
798 
799  t = p1 - p0;
800  u = p2 - p1;
801  v = p3 - p2;
802  w = p0 - p3;
803 
804  chose = G4RandFlat::shoot(0.,aOne+aTwo);
805  if( (chose>=0.) && (chose < aOne) )
806  {
807  lambda1 = G4RandFlat::shoot(0.,1.);
808  lambda2 = G4RandFlat::shoot(0.,lambda1);
809  return (p2+lambda1*v+lambda2*w);
810  }
811 
812  lambda1 = G4RandFlat::shoot(0.,1.);
813  lambda2 = G4RandFlat::shoot(0.,lambda1);
814  return (p0+lambda1*t+lambda2*u);
815 }
816 
817 // GetPointOnTriangle
818 //
819 // Auxiliary method for get point on surface
820 //
822  G4ThreeVector p2,
823  G4ThreeVector p3) const
824 {
825  G4double lambda1,lambda2;
826  G4ThreeVector v=p3-p1, w=p1-p2;
827 
828  lambda1 = G4RandFlat::shoot(0.,1.);
829  lambda2 = G4RandFlat::shoot(0.,lambda1);
830 
831  return (p2 + lambda1*w + lambda2*v);
832 }
833 
834 // GetPointOnSurface
835 //
837 {
838  if( !genericPgon ) // Polyhedra by faces
839  {
840  G4int j, numPlanes = original_parameters->Num_z_planes, Flag=0;
841  G4double chose, totArea=0., Achose1, Achose2,
842  rad1, rad2, sinphi1, sinphi2, cosphi1, cosphi2;
843  G4double a, b, l2, rang, totalPhi, ksi,
844  area, aTop=0., aBottom=0., zVal=0.;
845 
846  G4ThreeVector p0, p1, p2, p3;
847  std::vector<G4double> aVector1;
848  std::vector<G4double> aVector2;
849  std::vector<G4double> aVector3;
850 
851  totalPhi= (phiIsOpen) ? (endPhi-startPhi) : twopi;
852  ksi = totalPhi/numSide;
853  G4double cosksi = std::cos(ksi/2.);
854 
855  // Below we generate the areas relevant to our solid
856  //
857  for(j=0; j<numPlanes-1; ++j)
858  {
859  a = original_parameters->Rmax[j+1];
860  b = original_parameters->Rmax[j];
862  -original_parameters->Z_values[j+1]) + sqr(b-a);
863  area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
864  aVector1.push_back(area);
865  }
866 
867  for(j=0; j<numPlanes-1; ++j)
868  {
869  a = original_parameters->Rmin[j+1];//*cosksi;
870  b = original_parameters->Rmin[j];//*cosksi;
872  -original_parameters->Z_values[j+1]) + sqr(b-a);
873  area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
874  aVector2.push_back(area);
875  }
876 
877  for(j=0; j<numPlanes-1; ++j)
878  {
879  if(phiIsOpen == true)
880  {
881  aVector3.push_back(0.5*(original_parameters->Rmax[j]
884  -original_parameters->Rmin[j+1])
885  *std::fabs(original_parameters->Z_values[j+1]
887  }
888  else { aVector3.push_back(0.); }
889  }
890 
891  for(j=0; j<numPlanes-1; ++j)
892  {
893  totArea += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
894  }
895 
896  // Must include top and bottom areas
897  //
898  if(original_parameters->Rmax[numPlanes-1] != 0.)
899  {
900  a = original_parameters->Rmax[numPlanes-1];
901  b = original_parameters->Rmin[numPlanes-1];
902  l2 = sqr(a-b);
903  aTop = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
904  }
905 
906  if(original_parameters->Rmax[0] != 0.)
907  {
908  a = original_parameters->Rmax[0];
909  b = original_parameters->Rmin[0];
910  l2 = sqr(a-b);
911  aBottom = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
912  }
913 
914  Achose1 = 0.;
915  Achose2 = numSide*(aVector1[0]+aVector2[0])+2.*aVector3[0];
916 
917  chose = G4RandFlat::shoot(0.,totArea+aTop+aBottom);
918  if( (chose >= 0.) && (chose < aTop + aBottom) )
919  {
920  chose = G4RandFlat::shoot(startPhi,startPhi+totalPhi);
921  rang = std::floor((chose-startPhi)/ksi-0.01);
922  if(rang<0) { rang=0; }
923  rang = std::fabs(rang);
924  sinphi1 = std::sin(startPhi+rang*ksi);
925  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
926  cosphi1 = std::cos(startPhi+rang*ksi);
927  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
928  chose = G4RandFlat::shoot(0., aTop + aBottom);
929  if(chose>=0. && chose<aTop)
930  {
931  rad1 = original_parameters->Rmin[numPlanes-1];
932  rad2 = original_parameters->Rmax[numPlanes-1];
933  zVal = original_parameters->Z_values[numPlanes-1];
934  }
935  else
936  {
937  rad1 = original_parameters->Rmin[0];
938  rad2 = original_parameters->Rmax[0];
939  zVal = original_parameters->Z_values[0];
940  }
941  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
942  p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
943  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
944  p3 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
945  return GetPointOnPlane(p0,p1,p2,p3);
946  }
947  else
948  {
949  for (j=0; j<numPlanes-1; ++j)
950  {
951  if( ((chose >= Achose1) && (chose < Achose2)) || (j == numPlanes-2) )
952  {
953  Flag = j; break;
954  }
955  Achose1 += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
956  Achose2 = Achose1 + numSide*(aVector1[j+1]+aVector2[j+1])
957  + 2.*aVector3[j+1];
958  }
959  }
960 
961  // At this point we have chosen a subsection
962  // between to adjacent plane cuts...
963 
964  j = Flag;
965 
966  totArea = numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
967  chose = G4RandFlat::shoot(0.,totArea);
968 
969  if( (chose>=0.) && (chose<numSide*aVector1[j]) )
970  {
971  chose = G4RandFlat::shoot(startPhi,startPhi+totalPhi);
972  rang = std::floor((chose-startPhi)/ksi-0.01);
973  if(rang<0) { rang=0; }
974  rang = std::fabs(rang);
975  rad1 = original_parameters->Rmax[j];
976  rad2 = original_parameters->Rmax[j+1];
977  sinphi1 = std::sin(startPhi+rang*ksi);
978  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
979  cosphi1 = std::cos(startPhi+rang*ksi);
980  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
981  zVal = original_parameters->Z_values[j];
982 
983  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
984  p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
985 
986  zVal = original_parameters->Z_values[j+1];
987 
988  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
989  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
990  return GetPointOnPlane(p0,p1,p2,p3);
991  }
992  else if ( (chose >= numSide*aVector1[j])
993  && (chose <= numSide*(aVector1[j]+aVector2[j])) )
994  {
995  chose = G4RandFlat::shoot(startPhi,startPhi+totalPhi);
996  rang = std::floor((chose-startPhi)/ksi-0.01);
997  if(rang<0) { rang=0; }
998  rang = std::fabs(rang);
999  rad1 = original_parameters->Rmin[j];
1000  rad2 = original_parameters->Rmin[j+1];
1001  sinphi1 = std::sin(startPhi+rang*ksi);
1002  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
1003  cosphi1 = std::cos(startPhi+rang*ksi);
1004  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
1005  zVal = original_parameters->Z_values[j];
1006 
1007  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
1008  p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
1009 
1010  zVal = original_parameters->Z_values[j+1];
1011 
1012  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
1013  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
1014  return GetPointOnPlane(p0,p1,p2,p3);
1015  }
1016 
1017  chose = G4RandFlat::shoot(0.,2.2);
1018  if( (chose>=0.) && (chose < 1.) )
1019  {
1020  rang = startPhi;
1021  }
1022  else
1023  {
1024  rang = endPhi;
1025  }
1026 
1027  cosphi1 = std::cos(rang); rad1 = original_parameters->Rmin[j];
1028  sinphi1 = std::sin(rang); rad2 = original_parameters->Rmax[j];
1029 
1030  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
1032  p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
1034 
1035  rad1 = original_parameters->Rmax[j+1];
1036  rad2 = original_parameters->Rmin[j+1];
1037 
1038  p2 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
1040  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
1042  return GetPointOnPlane(p0,p1,p2,p3);
1043  }
1044  else // Generic polyhedra
1045  {
1046  return GetPointOnSurfaceGeneric();
1047  }
1048 }
1049 
1050 // CreatePolyhedron
1051 //
1053 {
1054  if (!genericPgon)
1055  {
1063  }
1064  else
1065  {
1066  // The following code prepares for:
1067  // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
1068  // const double xyz[][3],
1069  // const int faces_vec[][4])
1070  // Here is an extract from the header file HepPolyhedron.h:
1088  G4int nNodes;
1089  G4int nFaces;
1090  typedef G4double double3[3];
1091  double3* xyz;
1092  typedef G4int int4[4];
1093  int4* faces_vec;
1094  if (phiIsOpen)
1095  {
1096  // Triangulate open ends. Simple ear-chopping algorithm...
1097  // I'm not sure how robust this algorithm is (J.Allison).
1098  //
1099  std::vector<G4bool> chopped(numCorner, false);
1100  std::vector<G4int*> triQuads;
1101  G4int remaining = numCorner;
1102  G4int iStarter = 0;
1103  while (remaining >= 3) // Loop checking, 13.08.2015, G.Cosmo
1104  {
1105  // Find unchopped corners...
1106  //
1107  G4int A = -1, B = -1, C = -1;
1108  G4int iStepper = iStarter;
1109  do // Loop checking, 13.08.2015, G.Cosmo
1110  {
1111  if (A < 0) { A = iStepper; }
1112  else if (B < 0) { B = iStepper; }
1113  else if (C < 0) { C = iStepper; }
1114  do // Loop checking, 13.08.2015, G.Cosmo
1115  {
1116  if (++iStepper >= numCorner) iStepper = 0;
1117  }
1118  while (chopped[iStepper]);
1119  }
1120  while (C < 0 && iStepper != iStarter);
1121 
1122  // Check triangle at B is pointing outward (an "ear").
1123  // Sign of z cross product determines...
1124 
1125  G4double BAr = corners[A].r - corners[B].r;
1126  G4double BAz = corners[A].z - corners[B].z;
1127  G4double BCr = corners[C].r - corners[B].r;
1128  G4double BCz = corners[C].z - corners[B].z;
1129  if (BAr * BCz - BAz * BCr < kCarTolerance)
1130  {
1131  G4int* tq = new G4int[3];
1132  tq[0] = A + 1;
1133  tq[1] = B + 1;
1134  tq[2] = C + 1;
1135  triQuads.push_back(tq);
1136  chopped[B] = true;
1137  --remaining;
1138  }
1139  else
1140  {
1141  do // Loop checking, 13.08.2015, G.Cosmo
1142  {
1143  if (++iStarter >= numCorner) { iStarter = 0; }
1144  }
1145  while (chopped[iStarter]);
1146  }
1147  }
1148 
1149  // Transfer to faces...
1150 
1151  nNodes = (numSide + 1) * numCorner;
1152  nFaces = numSide * numCorner + 2 * triQuads.size();
1153  faces_vec = new int4[nFaces];
1154  G4int iface = 0;
1155  G4int addition = numCorner * numSide;
1156  G4int d = numCorner - 1;
1157  for (G4int iEnd = 0; iEnd < 2; ++iEnd)
1158  {
1159  for (size_t i = 0; i < triQuads.size(); ++i)
1160  {
1161  // Negative for soft/auxiliary/normally invisible edges...
1162  //
1163  G4int a, b, c;
1164  if (iEnd == 0)
1165  {
1166  a = triQuads[i][0];
1167  b = triQuads[i][1];
1168  c = triQuads[i][2];
1169  }
1170  else
1171  {
1172  a = triQuads[i][0] + addition;
1173  b = triQuads[i][2] + addition;
1174  c = triQuads[i][1] + addition;
1175  }
1176  G4int ab = std::abs(b - a);
1177  G4int bc = std::abs(c - b);
1178  G4int ca = std::abs(a - c);
1179  faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
1180  faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
1181  faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
1182  faces_vec[iface][3] = 0;
1183  ++iface;
1184  }
1185  }
1186 
1187  // Continue with sides...
1188 
1189  xyz = new double3[nNodes];
1190  const G4double dPhi = (endPhi - startPhi) / numSide;
1191  G4double phi = startPhi;
1192  G4int ixyz = 0;
1193  for (G4int iSide = 0; iSide < numSide; ++iSide)
1194  {
1195  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1196  {
1197  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1198  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1199  xyz[ixyz][2] = corners[iCorner].z;
1200  if (iCorner < numCorner - 1)
1201  {
1202  faces_vec[iface][0] = ixyz + 1;
1203  faces_vec[iface][1] = ixyz + numCorner + 1;
1204  faces_vec[iface][2] = ixyz + numCorner + 2;
1205  faces_vec[iface][3] = ixyz + 2;
1206  }
1207  else
1208  {
1209  faces_vec[iface][0] = ixyz + 1;
1210  faces_vec[iface][1] = ixyz + numCorner + 1;
1211  faces_vec[iface][2] = ixyz + 2;
1212  faces_vec[iface][3] = ixyz - numCorner + 2;
1213  }
1214  ++iface;
1215  ++ixyz;
1216  }
1217  phi += dPhi;
1218  }
1219 
1220  // Last corners...
1221 
1222  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1223  {
1224  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1225  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1226  xyz[ixyz][2] = corners[iCorner].z;
1227  ++ixyz;
1228  }
1229  }
1230  else // !phiIsOpen - i.e., a complete 360 degrees.
1231  {
1232  nNodes = numSide * numCorner;
1233  nFaces = numSide * numCorner;;
1234  xyz = new double3[nNodes];
1235  faces_vec = new int4[nFaces];
1236  // const G4double dPhi = (endPhi - startPhi) / numSide;
1237  const G4double dPhi = twopi / numSide;
1238  // !phiIsOpen endPhi-startPhi = 360 degrees.
1239  G4double phi = startPhi;
1240  G4int ixyz = 0, iface = 0;
1241  for (G4int iSide = 0; iSide < numSide; ++iSide)
1242  {
1243  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1244  {
1245  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1246  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1247  xyz[ixyz][2] = corners[iCorner].z;
1248  if (iSide < numSide - 1)
1249  {
1250  if (iCorner < numCorner - 1)
1251  {
1252  faces_vec[iface][0] = ixyz + 1;
1253  faces_vec[iface][1] = ixyz + numCorner + 1;
1254  faces_vec[iface][2] = ixyz + numCorner + 2;
1255  faces_vec[iface][3] = ixyz + 2;
1256  }
1257  else
1258  {
1259  faces_vec[iface][0] = ixyz + 1;
1260  faces_vec[iface][1] = ixyz + numCorner + 1;
1261  faces_vec[iface][2] = ixyz + 2;
1262  faces_vec[iface][3] = ixyz - numCorner + 2;
1263  }
1264  }
1265  else // Last side joins ends...
1266  {
1267  if (iCorner < numCorner - 1)
1268  {
1269  faces_vec[iface][0] = ixyz + 1;
1270  faces_vec[iface][1] = ixyz + numCorner - nFaces + 1;
1271  faces_vec[iface][2] = ixyz + numCorner - nFaces + 2;
1272  faces_vec[iface][3] = ixyz + 2;
1273  }
1274  else
1275  {
1276  faces_vec[iface][0] = ixyz + 1;
1277  faces_vec[iface][1] = ixyz - nFaces + numCorner + 1;
1278  faces_vec[iface][2] = ixyz - nFaces + 2;
1279  faces_vec[iface][3] = ixyz - numCorner + 2;
1280  }
1281  }
1282  ++ixyz;
1283  ++iface;
1284  }
1285  phi += dPhi;
1286  }
1287  }
1288  G4Polyhedron* polyhedron = new G4Polyhedron;
1289  G4int problem = polyhedron->createPolyhedron(nNodes,nFaces,xyz,faces_vec);
1290  delete [] faces_vec;
1291  delete [] xyz;
1292  if (problem)
1293  {
1294  std::ostringstream message;
1295  message << "Problem creating G4Polyhedron for: " << GetName();
1296  G4Exception("G4Polyhedra::CreatePolyhedron()", "GeomSolids1002",
1297  JustWarning, message);
1298  delete polyhedron;
1299  return nullptr;
1300  }
1301  else
1302  {
1303  return polyhedron;
1304  }
1305  }
1306 }
1307 
1308 // SetOriginalParameters
1309 //
1311 {
1312  G4int numPlanes = numCorner;
1313  G4bool isConvertible = true;
1314  G4double Zmax=rz->Bmax();
1315  rz->StartWithZMin();
1316 
1317  // Prepare vectors for storage
1318  //
1319  std::vector<G4double> Z;
1320  std::vector<G4double> Rmin;
1321  std::vector<G4double> Rmax;
1322 
1323  G4int countPlanes=1;
1324  G4int icurr=0;
1325  G4int icurl=0;
1326 
1327  // first plane Z=Z[0]
1328  //
1329  Z.push_back(corners[0].z);
1330  G4double Zprev=Z[0];
1331  if (Zprev == corners[1].z)
1332  {
1333  Rmin.push_back(corners[0].r);
1334  Rmax.push_back (corners[1].r);icurr=1;
1335  }
1336  else if (Zprev == corners[numPlanes-1].z)
1337  {
1338  Rmin.push_back(corners[numPlanes-1].r);
1339  Rmax.push_back (corners[0].r);
1340  icurl=numPlanes-1;
1341  }
1342  else
1343  {
1344  Rmin.push_back(corners[0].r);
1345  Rmax.push_back (corners[0].r);
1346  }
1347 
1348  // next planes until last
1349  //
1350  G4int inextr=0, inextl=0;
1351  for (G4int i=0; i < numPlanes-2; ++i)
1352  {
1353  inextr=1+icurr;
1354  inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1355 
1356  if((corners[inextr].z >= Zmax) & (corners[inextl].z >= Zmax)) { break; }
1357 
1358  G4double Zleft = corners[inextl].z;
1359  G4double Zright = corners[inextr].z;
1360  if(Zright>Zleft)
1361  {
1362  Z.push_back(Zleft);
1363  countPlanes++;
1364  G4double difZr=corners[inextr].z - corners[icurr].z;
1365  G4double difZl=corners[inextl].z - corners[icurl].z;
1366 
1367  if(std::fabs(difZl) < kCarTolerance)
1368  {
1369  if(std::fabs(difZr) < kCarTolerance)
1370  {
1371  Rmin.push_back(corners[inextl].r);
1372  Rmax.push_back(corners[icurr].r);
1373  }
1374  else
1375  {
1376  Rmin.push_back(corners[inextl].r);
1377  Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1378  *(corners[inextr].r - corners[icurr].r));
1379  }
1380  }
1381  else if (difZl >= kCarTolerance)
1382  {
1383  if(std::fabs(difZr) < kCarTolerance)
1384  {
1385  Rmin.push_back(corners[icurl].r);
1386  Rmax.push_back(corners[icurr].r);
1387  }
1388  else
1389  {
1390  Rmin.push_back(corners[icurl].r);
1391  Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1392  *(corners[inextr].r - corners[icurr].r));
1393  }
1394  }
1395  else
1396  {
1397  isConvertible=false; break;
1398  }
1399  icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1400  }
1401  else if(std::fabs(Zright-Zleft)<kCarTolerance) // Zright=Zleft
1402  {
1403  Z.push_back(Zleft);
1404  ++countPlanes;
1405  ++icurr;
1406 
1407  icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1408 
1409  Rmin.push_back(corners[inextl].r);
1410  Rmax.push_back (corners[inextr].r);
1411  }
1412  else // Zright<Zleft
1413  {
1414  Z.push_back(Zright);
1415  ++countPlanes;
1416 
1417  G4double difZr=corners[inextr].z - corners[icurr].z;
1418  G4double difZl=corners[inextl].z - corners[icurl].z;
1419  if(std::fabs(difZr) < kCarTolerance)
1420  {
1421  if(std::fabs(difZl) < kCarTolerance)
1422  {
1423  Rmax.push_back(corners[inextr].r);
1424  Rmin.push_back(corners[icurr].r);
1425  }
1426  else
1427  {
1428  Rmin.push_back(corners[icurl].r + (Zright-corners[icurl].z)/difZl
1429  * (corners[inextl].r - corners[icurl].r));
1430  Rmax.push_back(corners[inextr].r);
1431  }
1432  ++icurr;
1433  } // plate
1434  else if (difZr >= kCarTolerance)
1435  {
1436  if(std::fabs(difZl) < kCarTolerance)
1437  {
1438  Rmax.push_back(corners[inextr].r);
1439  Rmin.push_back (corners[icurr].r);
1440  }
1441  else
1442  {
1443  Rmax.push_back(corners[inextr].r);
1444  Rmin.push_back (corners[icurl].r+(Zright-corners[icurl].z)/difZl
1445  * (corners[inextl].r - corners[icurl].r));
1446  }
1447  ++icurr;
1448  }
1449  else
1450  {
1451  isConvertible=false; break;
1452  }
1453  }
1454  } // end for loop
1455 
1456  // last plane Z=Zmax
1457  //
1458  Z.push_back(Zmax);
1459  ++countPlanes;
1460  inextr=1+icurr;
1461  inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1462 
1463  Rmax.push_back(corners[inextr].r);
1464  Rmin.push_back(corners[inextl].r);
1465 
1466  // Set original parameters Rmin,Rmax,Z
1467  //
1468  if(isConvertible)
1469  {
1472  original_parameters->Z_values = new G4double[countPlanes];
1473  original_parameters->Rmin = new G4double[countPlanes];
1474  original_parameters->Rmax = new G4double[countPlanes];
1475 
1476  for(G4int j=0; j < countPlanes; ++j)
1477  {
1478  original_parameters->Z_values[j] = Z[j];
1479  original_parameters->Rmax[j] = Rmax[j];
1480  original_parameters->Rmin[j] = Rmin[j];
1481  }
1484  original_parameters->Num_z_planes = countPlanes;
1485 
1486  }
1487  else // Set parameters(r,z) with Rmin==0 as convention
1488  {
1489 #ifdef G4SPECSDEBUG
1490  std::ostringstream message;
1491  message << "Polyhedra " << GetName() << G4endl
1492  << "cannot be converted to Polyhedra with (Rmin,Rmaz,Z) parameters!";
1493  G4Exception("G4Polyhedra::SetOriginalParameters()",
1494  "GeomSolids0002", JustWarning, message);
1495 #endif
1498  original_parameters->Z_values = new G4double[numPlanes];
1499  original_parameters->Rmin = new G4double[numPlanes];
1500  original_parameters->Rmax = new G4double[numPlanes];
1501 
1502  for(G4int j=0; j < numPlanes; ++j)
1503  {
1505  original_parameters->Rmax[j] = corners[j].r;
1506  original_parameters->Rmin[j] = 0.0;
1507  }
1510  original_parameters->Num_z_planes = numPlanes;
1511  }
1512 }
1513 
1514 #endif