Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UPolyhedra.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 //
27 // $Id:$
28 //
29 // Implementation of G4UPolycone wrapper class
30 // --------------------------------------------------------------------
31 
32 #include "G4Polyhedra.hh"
33 #include "G4UPolyhedra.hh"
34 
35 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
36 
37 #include "G4GeomTools.hh"
38 #include "G4GeometryTolerance.hh"
39 #include "G4AffineTransform.hh"
40 #include "G4VPVParameterisation.hh"
41 #include "G4BoundingEnvelope.hh"
42 
43 using namespace CLHEP;
44 
46 //
47 // Constructor (GEANT3 style parameters)
48 //
49 // GEANT3 PGON radii are specified in the distance to the norm of each face.
50 //
51 G4UPolyhedra::G4UPolyhedra(const G4String& name,
52  G4double phiStart,
53  G4double phiTotal,
54  G4int numSide,
55  G4int numZPlanes,
56  const G4double zPlane[],
57  const G4double rInner[],
58  const G4double rOuter[] )
59  : Base_t(name, phiStart, phiTotal, numSide,
60  numZPlanes, zPlane, rInner, rOuter)
61 {
62  fGenericPgon = false;
63  SetOriginalParameters();
64  wrStart = phiStart;
65  while (wrStart < 0)
66  {
67  wrStart += twopi;
68  }
69  wrDelta = phiTotal;
70  if (wrDelta <= 0 || wrDelta >= twopi*(1-DBL_EPSILON))
71  {
72  wrDelta = twopi;
73  }
74  wrNumSide = numSide;
75  G4double convertRad = 1./std::cos(0.5*wrDelta/wrNumSide);
76  rzcorners.resize(0);
77  for (G4int i=0; i<numZPlanes; ++i)
78  {
79  G4double z = zPlane[i];
80  G4double r = rOuter[i]*convertRad;
81  rzcorners.push_back(G4TwoVector(r,z));
82  }
83  for (G4int i=numZPlanes-1; i>=0; --i)
84  {
85  G4double z = zPlane[i];
86  G4double r = rInner[i]*convertRad;
87  rzcorners.push_back(G4TwoVector(r,z));
88  }
89  std::vector<G4int> iout;
91 }
92 
93 
95 //
96 // Constructor (generic parameters)
97 //
98 G4UPolyhedra::G4UPolyhedra(const G4String& name,
99  G4double phiStart,
100  G4double phiTotal,
101  G4int numSide,
102  G4int numRZ,
103  const G4double r[],
104  const G4double z[] )
105  : Base_t(name, phiStart, phiTotal, numSide, numRZ, r, z)
106 {
107  fGenericPgon = true;
108  SetOriginalParameters();
109  wrStart = phiStart;
110  while (wrStart < 0)
111  {
112  wrStart += twopi;
113  }
114  wrDelta = phiTotal;
115  if (wrDelta <= 0 || wrDelta >= twopi*(1-DBL_EPSILON))
116  {
117  wrDelta = twopi;
118  }
119  wrNumSide = numSide;
120  G4double convertRad = 1./std::cos(0.5*wrDelta/wrNumSide);
121  rzcorners.resize(0);
122  for (G4int i=0; i<numRZ; ++i)
123  {
124  rzcorners.push_back(G4TwoVector(r[i]*convertRad,z[i]));
125  }
126  std::vector<G4int> iout;
128 }
129 
130 
132 //
133 // Fake default constructor - sets only member data and allocates memory
134 // for usage restricted to object persistency.
135 //
136 G4UPolyhedra::G4UPolyhedra( __void__& a )
137  : Base_t(a)
138 {
139 }
140 
141 
143 //
144 // Destructor
145 //
146 G4UPolyhedra::~G4UPolyhedra()
147 {
148 }
149 
150 
152 //
153 // Copy constructor
154 //
155 G4UPolyhedra::G4UPolyhedra( const G4UPolyhedra &source )
156  : Base_t( source )
157 {
158  fGenericPgon = source.fGenericPgon;
159  fOriginalParameters = source.fOriginalParameters;
160  wrStart = source.wrStart;
161  wrDelta = source.wrDelta;
162  wrNumSide = source.wrNumSide;
163  rzcorners = source.rzcorners;
164 }
165 
166 
168 //
169 // Assignment operator
170 //
171 G4UPolyhedra& G4UPolyhedra::operator=( const G4UPolyhedra &source )
172 {
173  if (this == &source) return *this;
174 
175  Base_t::operator=( source );
176  fGenericPgon = source.fGenericPgon;
177  fOriginalParameters = source.fOriginalParameters;
178  wrStart = source.wrStart;
179  wrDelta = source.wrDelta;
180  wrNumSide = source.wrNumSide;
181  rzcorners = source.rzcorners;
182 
183  return *this;
184 }
185 
186 
188 //
189 // Accessors & modifiers
190 //
191 G4int G4UPolyhedra::GetNumSide() const
192 {
193  return wrNumSide;
194 }
195 G4double G4UPolyhedra::GetStartPhi() const
196 {
197  return wrStart;
198 }
199 G4double G4UPolyhedra::GetEndPhi() const
200 {
201  return (wrStart + wrDelta);
202 }
203 G4double G4UPolyhedra::GetSinStartPhi() const
204 {
205  G4double phi = GetStartPhi();
206  return std::sin(phi);
207 }
208 G4double G4UPolyhedra::GetCosStartPhi() const
209 {
210  G4double phi = GetStartPhi();
211  return std::cos(phi);
212 }
213 G4double G4UPolyhedra::GetSinEndPhi() const
214 {
215  G4double phi = GetEndPhi();
216  return std::sin(phi);
217 }
218 G4double G4UPolyhedra::GetCosEndPhi() const
219 {
220  G4double phi = GetEndPhi();
221  return std::cos(phi);
222 }
223 G4bool G4UPolyhedra::IsOpen() const
224 {
225  return (wrDelta < twopi);
226 }
227 G4bool G4UPolyhedra::IsGeneric() const
228 {
229  return fGenericPgon;
230 }
231 G4int G4UPolyhedra::GetNumRZCorner() const
232 {
233  return rzcorners.size();
234 }
235 G4PolyhedraSideRZ G4UPolyhedra::GetCorner(G4int index) const
236 {
237  G4TwoVector rz = rzcorners.at(index);
238  G4PolyhedraSideRZ psiderz = { rz.x(), rz.y() };
239 
240  return psiderz;
241 }
242 G4PolyhedraHistorical* G4UPolyhedra::GetOriginalParameters() const
243 {
244  return new G4PolyhedraHistorical(fOriginalParameters);
245 }
246 void G4UPolyhedra::SetOriginalParameters()
247 {
248  G4double startPhi = GetPhiStart();
249  G4double deltaPhi = GetPhiDelta();
250  G4int numPlanes = GetZSegmentCount() + 1;
251  G4int numSides = GetSideCount();
252 
253  fOriginalParameters.Start_angle = startPhi;
254  fOriginalParameters.Opening_angle = deltaPhi;
255  fOriginalParameters.Num_z_planes = numPlanes;
256  fOriginalParameters.numSide = numSides;
257 
258  delete [] fOriginalParameters.Z_values;
259  delete [] fOriginalParameters.Rmin;
260  delete [] fOriginalParameters.Rmax;
261  fOriginalParameters.Z_values = new G4double[numPlanes];
262  fOriginalParameters.Rmin = new G4double[numPlanes];
263  fOriginalParameters.Rmax = new G4double[numPlanes];
264 
265  G4double convertRad = std::cos(0.5*deltaPhi/numSides);
266  for (G4int i=0; i<numPlanes; ++i)
267  {
268  fOriginalParameters.Z_values[i] = GetZPlanes()[i];
269  fOriginalParameters.Rmax[i] = GetRMax()[i]/convertRad;
270  fOriginalParameters.Rmin[i] = GetRMin()[i]/convertRad;
271  }
272 }
273 void G4UPolyhedra::SetOriginalParameters(G4PolyhedraHistorical* pars)
274 {
275  fOriginalParameters = *pars;
276  fRebuildPolyhedron = true;
277  Reset();
278 }
279 
281 {
282  if (fGenericPgon)
283  {
284  std::ostringstream message;
285  message << "Solid " << GetName() << " built using generic construct."
286  << G4endl << "Not applicable to the generic construct !";
287  G4Exception("G4UPolyhedra::Reset()", "GeomSolids1001",
288  JustWarning, message, "Parameters NOT resetted.");
289  return true; // error code set
290  }
291 
292  //
293  // Rebuild polyhedra based on original parameters
294  //
295  wrStart = fOriginalParameters.Start_angle;
296  while (wrStart < 0)
297  {
298  wrStart += twopi;
299  }
300  wrDelta = fOriginalParameters.Opening_angle;
301  if (wrDelta <= 0 || wrDelta >= twopi*(1-DBL_EPSILON))
302  {
303  wrDelta = twopi;
304  }
305  wrNumSide = fOriginalParameters.numSide;
306  rzcorners.resize(0);
307  for (G4int i=0; i<fOriginalParameters.Num_z_planes; ++i)
308  {
309  G4double z = fOriginalParameters.Z_values[i];
310  G4double r = fOriginalParameters.Rmax[i];
311  rzcorners.push_back(G4TwoVector(r,z));
312  }
313  for (G4int i=fOriginalParameters.Num_z_planes-1; i>=0; --i)
314  {
315  G4double z = fOriginalParameters.Z_values[i];
316  G4double r = fOriginalParameters.Rmin[i];
317  rzcorners.push_back(G4TwoVector(r,z));
318  }
319  std::vector<G4int> iout;
321 
322  return false; // error code unset
323 }
324 
325 
327 //
328 // Dispatch to parameterisation for replication mechanism dimension
329 // computation & modification.
330 //
331 void G4UPolyhedra::ComputeDimensions(G4VPVParameterisation* p,
332  const G4int n,
333  const G4VPhysicalVolume* pRep)
334 {
335  p->ComputeDimensions(*(G4Polyhedra*)this,n,pRep);
336 }
337 
338 
340 //
341 // Make a clone of the object
342 
343 G4VSolid* G4UPolyhedra::Clone() const
344 {
345  return new G4UPolyhedra(*this);
346 }
347 
348 
350 //
351 // Get bounding box
352 
353 void G4UPolyhedra::BoundingLimits(G4ThreeVector& pMin,
354  G4ThreeVector& pMax) const
355 {
356  static G4bool checkBBox = true;
357  static G4bool checkPhi = true;
358 
359  G4double rmin = kInfinity, rmax = -kInfinity;
360  G4double zmin = kInfinity, zmax = -kInfinity;
361  for (G4int i=0; i<GetNumRZCorner(); ++i)
362  {
363  G4PolyhedraSideRZ corner = GetCorner(i);
364  if (corner.r < rmin) rmin = corner.r;
365  if (corner.r > rmax) rmax = corner.r;
366  if (corner.z < zmin) zmin = corner.z;
367  if (corner.z > zmax) zmax = corner.z;
368  }
369 
370  G4double sphi = GetStartPhi();
371  G4double ephi = GetEndPhi();
372  G4double dphi = IsOpen() ? ephi-sphi : twopi;
373  G4int ksteps = GetNumSide();
374  G4double astep = dphi/ksteps;
375  G4double sinStep = std::sin(astep);
376  G4double cosStep = std::cos(astep);
377 
378  G4double sinCur = GetSinStartPhi();
379  G4double cosCur = GetCosStartPhi();
380  if (!IsOpen()) rmin = 0;
381  G4double xmin = rmin*cosCur, xmax = xmin;
382  G4double ymin = rmin*sinCur, ymax = ymin;
383  for (G4int k=0; k<ksteps+1; ++k)
384  {
385  G4double x = rmax*cosCur;
386  if (x < xmin) xmin = x;
387  if (x > xmax) xmax = x;
388  G4double y = rmax*sinCur;
389  if (y < ymin) ymin = y;
390  if (y > ymax) ymax = y;
391  if (rmin > 0)
392  {
393  G4double xx = rmin*cosCur;
394  if (xx < xmin) xmin = xx;
395  if (xx > xmax) xmax = xx;
396  G4double yy = rmin*sinCur;
397  if (yy < ymin) ymin = yy;
398  if (yy > ymax) ymax = yy;
399  }
400  G4double sinTmp = sinCur;
401  sinCur = sinCur*cosStep + cosCur*sinStep;
402  cosCur = cosCur*cosStep - sinTmp*sinStep;
403  }
404  pMin.set(xmin,ymin,zmin);
405  pMax.set(xmax,ymax,zmax);
406 
407  // Check correctness of the bounding box
408  //
409  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
410  {
411  std::ostringstream message;
412  message << "Bad bounding box (min >= max) for solid: "
413  << GetName() << " !"
414  << "\npMin = " << pMin
415  << "\npMax = " << pMax;
416  G4Exception("G4UPolyhedra::BoundingLimits()", "GeomMgt0001",
417  JustWarning, message);
418  StreamInfo(G4cout);
419  }
420 
421  // Check consistency of bounding boxes
422  //
423  if (checkBBox)
424  {
425  U3Vector vmin, vmax;
426  Extent(vmin,vmax);
427  if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
428  std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
429  std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
430  std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
431  std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
432  std::abs(pMax.z()-vmax.z()) > kCarTolerance)
433  {
434  std::ostringstream message;
435  message << "Inconsistency in bounding boxes for solid: "
436  << GetName() << " !"
437  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
438  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
439  G4Exception("G4UPolyhedra::BoundingLimits()", "GeomMgt0001",
440  JustWarning, message);
441  checkBBox = false;
442  }
443  }
444 
445  // Check consistency of angles
446  //
447  if (checkPhi)
448  {
449  if (GetStartPhi() != GetPhiStart() ||
450  GetEndPhi() != GetPhiEnd() ||
451  GetNumSide() != GetSideCount() ||
452  IsOpen() != (Base_t::GetPhiDelta() < twopi))
453  {
454  std::ostringstream message;
455  message << "Inconsistency in Phi angles or # of sides for solid: "
456  << GetName() << " !"
457  << "\nPhi start : wrapper = " << GetStartPhi()
458  << " solid = " << GetPhiStart()
459  << "\nPhi end : wrapper = " << GetEndPhi()
460  << " solid = " << GetPhiEnd()
461  << "\nPhi # sides: wrapper = " << GetNumSide()
462  << " solid = " << GetSideCount()
463  << "\nPhi is open: wrapper = " << (IsOpen() ? "true" : "false")
464  << " solid = "
465  << ((Base_t::GetPhiDelta() < twopi) ? "true" : "false");
466  G4Exception("G4UPolyhedra::BoundingLimits()", "GeomMgt0001",
467  JustWarning, message);
468  checkPhi = false;
469  }
470  }
471 }
472 
474 //
475 // Calculate extent under transform and specified limit
476 
477 G4bool
478 G4UPolyhedra::CalculateExtent(const EAxis pAxis,
479  const G4VoxelLimits& pVoxelLimit,
480  const G4AffineTransform& pTransform,
481  G4double& pMin, G4double& pMax) const
482 {
483  G4ThreeVector bmin, bmax;
484  G4bool exist;
485 
486  // Check bounding box (bbox)
487  //
488  BoundingLimits(bmin,bmax);
489  G4BoundingEnvelope bbox(bmin,bmax);
490 #ifdef G4BBOX_EXTENT
491  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
492 #endif
493  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
494  {
495  return exist = (pMin < pMax) ? true : false;
496  }
497 
498  // To find the extent, RZ contour of the polycone is subdivided
499  // in triangles. The extent is calculated as cumulative extent of
500  // all sub-polycones formed by rotation of triangles around Z
501  //
502  G4TwoVectorList contourRZ;
503  G4TwoVectorList triangles;
504  std::vector<G4int> iout;
505  G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
506  G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
507 
508  // get RZ contour, ensure anticlockwise order of corners
509  for (G4int i=0; i<GetNumRZCorner(); ++i)
510  {
511  G4PolyhedraSideRZ corner = GetCorner(i);
512  contourRZ.push_back(G4TwoVector(corner.r,corner.z));
513  }
515  G4double area = G4GeomTools::PolygonArea(contourRZ);
516  if (area < 0.) std::reverse(contourRZ.begin(),contourRZ.end());
517 
518  // triangulate RZ countour
519  if (!G4GeomTools::TriangulatePolygon(contourRZ,triangles))
520  {
521  std::ostringstream message;
522  message << "Triangulation of RZ contour has failed for solid: "
523  << GetName() << " !"
524  << "\nExtent has been calculated using boundary box";
525  G4Exception("G4UPolyhedra::CalculateExtent()",
526  "GeomMgt1002",JustWarning,message);
527  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
528  }
529 
530  // set trigonometric values
531  G4double sphi = GetStartPhi();
532  G4double ephi = GetEndPhi();
533  G4double dphi = IsOpen() ? ephi-sphi : twopi;
534  G4int ksteps = GetNumSide();
535  G4double astep = dphi/ksteps;
536  G4double sinStep = std::sin(astep);
537  G4double cosStep = std::cos(astep);
538  G4double sinStart = GetSinStartPhi();
539  G4double cosStart = GetCosStartPhi();
540 
541  // allocate vector lists
542  std::vector<const G4ThreeVectorList *> polygons;
543  polygons.resize(ksteps+1);
544  for (G4int k=0; k<ksteps+1; ++k) {
545  polygons[k] = new G4ThreeVectorList(3);
546  }
547 
548  // main loop along triangles
549  pMin = kInfinity;
550  pMax = -kInfinity;
551  G4int ntria = triangles.size()/3;
552  for (G4int i=0; i<ntria; ++i)
553  {
554  G4double sinCur = sinStart;
555  G4double cosCur = cosStart;
556  G4int i3 = i*3;
557  for (G4int k=0; k<ksteps+1; ++k) // rotate triangle
558  {
559  G4ThreeVectorList* ptr = const_cast<G4ThreeVectorList*>(polygons[k]);
560  G4ThreeVectorList::iterator iter = ptr->begin();
561  iter->set(triangles[i3+0].x()*cosCur,
562  triangles[i3+0].x()*sinCur,
563  triangles[i3+0].y());
564  iter++;
565  iter->set(triangles[i3+1].x()*cosCur,
566  triangles[i3+1].x()*sinCur,
567  triangles[i3+1].y());
568  iter++;
569  iter->set(triangles[i3+2].x()*cosCur,
570  triangles[i3+2].x()*sinCur,
571  triangles[i3+2].y());
572 
573  G4double sinTmp = sinCur;
574  sinCur = sinCur*cosStep + cosCur*sinStep;
575  cosCur = cosCur*cosStep - sinTmp*sinStep;
576  }
577 
578  // set sub-envelope and adjust extent
579  G4double emin,emax;
580  G4BoundingEnvelope benv(polygons);
581  if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
582  if (emin < pMin) pMin = emin;
583  if (emax > pMax) pMax = emax;
584  if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
585  }
586  // free memory
587  for (G4int k=0; k<ksteps+1; ++k) { delete polygons[k]; polygons[k]=0;}
588  return (pMin < pMax);
589 }
590 
591 
593 //
594 // CreatePolyhedron
595 //
596 G4Polyhedron* G4UPolyhedra::CreatePolyhedron() const
597 {
598  if (!IsGeneric())
599  {
600  return new G4PolyhedronPgon( fOriginalParameters.Start_angle,
601  fOriginalParameters.Opening_angle,
602  fOriginalParameters.numSide,
603  fOriginalParameters.Num_z_planes,
604  fOriginalParameters.Z_values,
605  fOriginalParameters.Rmin,
606  fOriginalParameters.Rmax);
607  }
608  else
609  {
610  // The following code prepares for:
611  // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
612  // const double xyz[][3],
613  // const int faces_vec[][4])
614  // Here is an extract from the header file HepPolyhedron.h:
632  G4int nNodes;
633  G4int nFaces;
634  typedef G4double double3[3];
635  double3* xyz;
636  typedef G4int int4[4];
637  int4* faces_vec;
638  if (IsOpen())
639  {
640  // Triangulate open ends. Simple ear-chopping algorithm...
641  // I'm not sure how robust this algorithm is (J.Allison).
642  //
643  std::vector<G4bool> chopped(GetNumRZCorner(), false);
644  std::vector<G4int*> triQuads;
645  G4int remaining = GetNumRZCorner();
646  G4int iStarter = 0;
647  while (remaining >= 3) // Loop checking, 13.08.2015, G.Cosmo
648  {
649  // Find unchopped corners...
650  //
651  G4int A = -1, B = -1, C = -1;
652  G4int iStepper = iStarter;
653  do // Loop checking, 13.08.2015, G.Cosmo
654  {
655  if (A < 0) { A = iStepper; }
656  else if (B < 0) { B = iStepper; }
657  else if (C < 0) { C = iStepper; }
658  do // Loop checking, 13.08.2015, G.Cosmo
659  {
660  if (++iStepper >= GetNumRZCorner()) iStepper = 0;
661  }
662  while (chopped[iStepper]);
663  }
664  while (C < 0 && iStepper != iStarter);
665 
666  // Check triangle at B is pointing outward (an "ear").
667  // Sign of z cross product determines...
668 
669  G4double BAr = GetCorner(A).r - GetCorner(B).r;
670  G4double BAz = GetCorner(A).z - GetCorner(B).z;
671  G4double BCr = GetCorner(C).r - GetCorner(B).r;
672  G4double BCz = GetCorner(C).z - GetCorner(B).z;
673  if (BAr * BCz - BAz * BCr < kCarTolerance)
674  {
675  G4int* tq = new G4int[3];
676  tq[0] = A + 1;
677  tq[1] = B + 1;
678  tq[2] = C + 1;
679  triQuads.push_back(tq);
680  chopped[B] = true;
681  --remaining;
682  }
683  else
684  {
685  do // Loop checking, 13.08.2015, G.Cosmo
686  {
687  if (++iStarter >= GetNumRZCorner()) { iStarter = 0; }
688  }
689  while (chopped[iStarter]);
690  }
691  }
692 
693  // Transfer to faces...
694  G4int numSide=GetNumSide();
695  nNodes = (numSide + 1) * GetNumRZCorner();
696  nFaces = numSide * GetNumRZCorner() + 2 * triQuads.size();
697  faces_vec = new int4[nFaces];
698  G4int iface = 0;
699  G4int addition = GetNumRZCorner() * numSide;
700  G4int d = GetNumRZCorner() - 1;
701  for (G4int iEnd = 0; iEnd < 2; ++iEnd)
702  {
703  for (size_t i = 0; i < triQuads.size(); ++i)
704  {
705  // Negative for soft/auxiliary/normally invisible edges...
706  //
707  G4int a, b, c;
708  if (iEnd == 0)
709  {
710  a = triQuads[i][0];
711  b = triQuads[i][1];
712  c = triQuads[i][2];
713  }
714  else
715  {
716  a = triQuads[i][0] + addition;
717  b = triQuads[i][2] + addition;
718  c = triQuads[i][1] + addition;
719  }
720  G4int ab = std::abs(b - a);
721  G4int bc = std::abs(c - b);
722  G4int ca = std::abs(a - c);
723  faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
724  faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
725  faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
726  faces_vec[iface][3] = 0;
727  ++iface;
728  }
729  }
730 
731  // Continue with sides...
732 
733  xyz = new double3[nNodes];
734  const G4double dPhi = (GetEndPhi() - GetStartPhi()) / numSide;
735  G4double phi = GetStartPhi();
736  G4int ixyz = 0;
737  for (G4int iSide = 0; iSide < numSide; ++iSide)
738  {
739  for (G4int iCorner = 0; iCorner < GetNumRZCorner(); ++iCorner)
740  {
741  xyz[ixyz][0] = GetCorner(iCorner).r * std::cos(phi);
742  xyz[ixyz][1] = GetCorner(iCorner).r * std::sin(phi);
743  xyz[ixyz][2] = GetCorner(iCorner).z;
744  if (iCorner < GetNumRZCorner() - 1)
745  {
746  faces_vec[iface][0] = ixyz + 1;
747  faces_vec[iface][1] = ixyz + GetNumRZCorner() + 1;
748  faces_vec[iface][2] = ixyz + GetNumRZCorner() + 2;
749  faces_vec[iface][3] = ixyz + 2;
750  }
751  else
752  {
753  faces_vec[iface][0] = ixyz + 1;
754  faces_vec[iface][1] = ixyz + GetNumRZCorner() + 1;
755  faces_vec[iface][2] = ixyz + 2;
756  faces_vec[iface][3] = ixyz - GetNumRZCorner() + 2;
757  }
758  ++iface;
759  ++ixyz;
760  }
761  phi += dPhi;
762  }
763 
764  // Last GetCorner...
765 
766  for (G4int iCorner = 0; iCorner < GetNumRZCorner(); ++iCorner)
767  {
768  xyz[ixyz][0] = GetCorner(iCorner).r * std::cos(phi);
769  xyz[ixyz][1] = GetCorner(iCorner).r * std::sin(phi);
770  xyz[ixyz][2] = GetCorner(iCorner).z;
771  ++ixyz;
772  }
773  }
774  else // !phiIsOpen - i.e., a complete 360 degrees.
775  {
776  nNodes = GetNumSide() * GetNumRZCorner();
777  nFaces = GetNumSide() * GetNumRZCorner();;
778  xyz = new double3[nNodes];
779  faces_vec = new int4[nFaces];
780  // const G4double dPhi = (endPhi - startPhi) / numSide;
781  const G4double dPhi = twopi / GetNumSide();
782  // !phiIsOpen endPhi-startPhi = 360 degrees.
783  G4double phi = GetStartPhi();
784  G4int ixyz = 0, iface = 0;
785  for (G4int iSide = 0; iSide < GetNumSide(); ++iSide)
786  {
787  for (G4int iCorner = 0; iCorner < GetNumRZCorner(); ++iCorner)
788  {
789  xyz[ixyz][0] = GetCorner(iCorner).r * std::cos(phi);
790  xyz[ixyz][1] = GetCorner(iCorner).r * std::sin(phi);
791  xyz[ixyz][2] = GetCorner(iCorner).z;
792  if (iSide < GetNumSide() - 1)
793  {
794  if (iCorner < GetNumRZCorner() - 1)
795  {
796  faces_vec[iface][0] = ixyz + 1;
797  faces_vec[iface][1] = ixyz + GetNumRZCorner() + 1;
798  faces_vec[iface][2] = ixyz + GetNumRZCorner() + 2;
799  faces_vec[iface][3] = ixyz + 2;
800  }
801  else
802  {
803  faces_vec[iface][0] = ixyz + 1;
804  faces_vec[iface][1] = ixyz + GetNumRZCorner() + 1;
805  faces_vec[iface][2] = ixyz + 2;
806  faces_vec[iface][3] = ixyz - GetNumRZCorner() + 2;
807  }
808  }
809  else // Last side joins ends...
810  {
811  if (iCorner < GetNumRZCorner() - 1)
812  {
813  faces_vec[iface][0] = ixyz + 1;
814  faces_vec[iface][1] = ixyz + GetNumRZCorner() - nFaces + 1;
815  faces_vec[iface][2] = ixyz + GetNumRZCorner() - nFaces + 2;
816  faces_vec[iface][3] = ixyz + 2;
817  }
818  else
819  {
820  faces_vec[iface][0] = ixyz + 1;
821  faces_vec[iface][1] = ixyz - nFaces + GetNumRZCorner() + 1;
822  faces_vec[iface][2] = ixyz - nFaces + 2;
823  faces_vec[iface][3] = ixyz - GetNumRZCorner() + 2;
824  }
825  }
826  ++ixyz;
827  ++iface;
828  }
829  phi += dPhi;
830  }
831  }
832  G4Polyhedron* polyhedron = new G4Polyhedron;
833  G4int prob = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec);
834  delete [] faces_vec;
835  delete [] xyz;
836  if (prob)
837  {
838  std::ostringstream message;
839  message << "Problem creating G4Polyhedron for: " << GetName();
840  G4Exception("G4Polyhedra::CreatePolyhedron()", "GeomSolids1002",
841  JustWarning, message);
842  delete polyhedron;
843  return 0;
844  }
845  else
846  {
847  return polyhedron;
848  }
849  }
850 }
851 
852 #endif // G4GEOM_USE_USOLIDS
Float_t x
Definition: compare.C:6
void set(double x, double y, double z)
const XML_Char * name
Definition: expat.h:151
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
Definition: G4GeomTools.cc:199
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
Double_t xx
G4double GetMaxExtent(const EAxis pAxis) const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
static const G4double kInfinity
Definition: geomdefs.hh:42
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
#define G4endl
Definition: G4ios.hh:61
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
Double_t z
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
double z() const
static const G4double emax
const G4double kCarTolerance
static const G4double ab
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
G4double GetMinExtent(const EAxis pAxis) const
double x() const
ntupleExperimental Reset()
static constexpr double twopi
Definition: G4SIunits.hh:76
double A(double temperature)
#define DBL_EPSILON
Definition: templates.hh:87
Float_t d
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4int
Definition: G4Types.hh:78
std::vector< G4TwoVector > G4TwoVectorList
Definition: G4GeomTools.hh:50
EAxis
Definition: geomdefs.hh:54
double C(double temp)
double y() const
G4GLOB_DLL std::ostream G4cout
double x() const
static void RemoveRedundantVertices(G4TwoVectorList &polygon, std::vector< G4int > &iout, G4double tolerance=0)
Definition: G4GeomTools.cc:311
Char_t n[5]
double y() const
double B(double temperature)
static G4double PolygonArea(const G4TwoVectorList &polygon)
Definition: G4GeomTools.cc:82