ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ExtrudedSolid.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4ExtrudedSolid.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 // G4ExtrudedSolid implementation
27 //
28 // Author: Ivana Hrivnacova, IPN Orsay
29 //
30 // CHANGE HISTORY
31 // --------------
32 //
33 // 31.10.2017 E.Tcherniaev: added implementation for a non-convex
34 // right prism
35 // 08.09.2017 E.Tcherniaev: added implementation for a convex
36 // right prism
37 // 21.10.2016 E.Tcherniaev: reimplemented CalculateExtent(),
38 // used G4GeomTools::PolygonArea() to calculate area,
39 // replaced IsConvex() with G4GeomTools::IsConvex()
40 // 02.03.2016 E.Tcherniaev: added CheckPolygon() to remove
41 // collinear and coincident points from polygon
42 // --------------------------------------------------------------------
43 
44 #include "G4ExtrudedSolid.hh"
45 
46 #if !defined(G4GEOM_USE_UEXTRUDEDSOLID)
47 
48 #include <set>
49 #include <algorithm>
50 #include <cmath>
51 #include <iomanip>
52 
53 #include "G4GeomTools.hh"
54 #include "G4VoxelLimits.hh"
55 #include "G4AffineTransform.hh"
56 #include "G4BoundingEnvelope.hh"
57 
58 #include "G4GeometryTolerance.hh"
59 #include "G4PhysicalConstants.hh"
60 #include "G4SystemOfUnits.hh"
61 #include "G4TriangularFacet.hh"
62 #include "G4QuadrangularFacet.hh"
63 
64 //_____________________________________________________________________________
65 
67  const std::vector<G4TwoVector>& polygon,
68  const std::vector<ZSection>& zsections)
69  : G4TessellatedSolid(pName),
70  fNv(polygon.size()),
71  fNz(zsections.size()),
72  fIsConvex(false),
73  fGeometryType("G4ExtrudedSolid"),
74  fSolidType(0)
75 {
76  // General constructor
77 
78  // First check input parameters
79 
80  if (fNv < 3)
81  {
82  std::ostringstream message;
83  message << "Number of vertices in polygon < 3 - " << pName;
84  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
85  FatalErrorInArgument, message);
86  }
87 
88  if (fNz < 2)
89  {
90  std::ostringstream message;
91  message << "Number of z-sides < 2 - " << pName;
92  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
93  FatalErrorInArgument, message);
94  }
95 
96  for ( G4int i=0; i<fNz-1; ++i )
97  {
98  if ( zsections[i].fZ > zsections[i+1].fZ )
99  {
100  std::ostringstream message;
101  message << "Z-sections have to be ordered by z value (z0 < z1 < z2...) - "
102  << pName;
103  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
104  FatalErrorInArgument, message);
105  }
106  if ( std::fabs( zsections[i+1].fZ - zsections[i].fZ ) < kCarToleranceHalf )
107  {
108  std::ostringstream message;
109  message << "Z-sections with the same z position are not supported - "
110  << pName;
111  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0001",
112  FatalException, message);
113  }
114  }
115 
116  // Copy polygon
117  //
118  fPolygon = polygon;
119 
120  // Remove collinear and coincident vertices, if any
121  //
122  std::vector<G4int> removedVertices;
124  2*kCarTolerance);
125  if (removedVertices.size() != 0)
126  {
127  G4int nremoved = removedVertices.size();
128  std::ostringstream message;
129  message << "The following "<< nremoved
130  << " vertices have been removed from polygon in " << pName
131  << "\nas collinear or coincident with other vertices: "
132  << removedVertices[0];
133  for (G4int i=1; i<nremoved; ++i) message << ", " << removedVertices[i];
134  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
135  JustWarning, message);
136  }
137 
138  fNv = fPolygon.size();
139  if (fNv < 3)
140  {
141  std::ostringstream message;
142  message << "Number of vertices in polygon after removal < 3 - " << pName;
143  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
144  FatalErrorInArgument, message);
145  }
146 
147  // Check if polygon vertices are defined clockwise
148  // (the area is positive if polygon vertices are defined anti-clockwise)
149  //
151  {
152  // Polygon vertices are defined anti-clockwise, we revert them
153  // G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
154  // JustWarning,
155  // "Polygon vertices defined anti-clockwise, reverting polygon");
156  std::reverse(fPolygon.begin(),fPolygon.end());
157  }
158 
159  // Copy z-sections
160  //
161  fZSections = zsections;
162 
163  G4bool result = MakeFacets();
164  if (!result)
165  {
166  std::ostringstream message;
167  message << "Making facets failed - " << pName;
168  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0003",
169  FatalException, message);
170  }
172 
174 
175  // Check if the solid is a right prism, if so then set lateral planes
176  //
177  if ((fNz == 2)
178  && (fZSections[0].fScale == 1) && (fZSections[1].fScale == 1)
179  && (fZSections[0].fOffset == G4TwoVector(0,0))
180  && (fZSections[1].fOffset == G4TwoVector(0,0)))
181  {
182  fSolidType = (fIsConvex) ? 1 : 2; // 1 - convex, 2 - non-convex right prism
184  }
185 }
186 
187 //_____________________________________________________________________________
188 
190  const std::vector<G4TwoVector>& polygon,
191  G4double dz,
192  const G4TwoVector& off1, G4double scale1,
193  const G4TwoVector& off2, G4double scale2 )
194  : G4TessellatedSolid(pName),
195  fNv(polygon.size()),
196  fNz(2),
197  fGeometryType("G4ExtrudedSolid")
198 {
199  // Special constructor for solid with 2 z-sections
200 
201  // First check input parameters
202  //
203  if (fNv < 3)
204  {
205  std::ostringstream message;
206  message << "Number of vertices in polygon < 3 - " << pName;
207  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
208  FatalErrorInArgument, message);
209  }
210 
211  // Copy polygon
212  //
213  fPolygon = polygon;
214 
215  // Remove collinear and coincident vertices, if any
216  //
217  std::vector<G4int> removedVertices;
219  2*kCarTolerance);
220  if (removedVertices.size() != 0)
221  {
222  G4int nremoved = removedVertices.size();
223  std::ostringstream message;
224  message << "The following "<< nremoved
225  << " vertices have been removed from polygon in " << pName
226  << "\nas collinear or coincident with other vertices: "
227  << removedVertices[0];
228  for (G4int i=1; i<nremoved; ++i) message << ", " << removedVertices[i];
229  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
230  JustWarning, message);
231  }
232 
233  fNv = fPolygon.size();
234  if (fNv < 3)
235  {
236  std::ostringstream message;
237  message << "Number of vertices in polygon after removal < 3 - " << pName;
238  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
239  FatalErrorInArgument, message);
240  }
241 
242  // Check if polygon vertices are defined clockwise
243  // (the area is positive if polygon vertices are defined anti-clockwise)
244  //
246  {
247  // Polygon vertices are defined anti-clockwise, we revert them
248  // G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
249  // JustWarning,
250  // "Polygon vertices defined anti-clockwise, reverting polygon");
251  std::reverse(fPolygon.begin(),fPolygon.end());
252  }
253 
254  // Copy z-sections
255  //
256  fZSections.push_back(ZSection(-dz, off1, scale1));
257  fZSections.push_back(ZSection( dz, off2, scale2));
258 
259  G4bool result = MakeFacets();
260  if (!result)
261  {
262  std::ostringstream message;
263  message << "Making facets failed - " << pName;
264  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0003",
265  FatalException, message);
266  }
268 
270 
271  // Check if the solid is a right prism, if so then set lateral planes
272  //
273  if ((scale1 == 1) && (scale2 == 1)
274  && (off1 == G4TwoVector(0,0)) && (off2 == G4TwoVector(0,0)))
275  {
276  fSolidType = (fIsConvex) ? 1 : 2; // 1 - convex, 2 - non-convex right prism
278  }
279 }
280 
281 //_____________________________________________________________________________
282 
284  : G4TessellatedSolid(a), fNv(0), fNz(0),
285  fGeometryType("G4ExtrudedSolid")
286 {
287  // Fake default constructor - sets only member data and allocates memory
288  // for usage restricted to object persistency.
289 }
290 
291 //_____________________________________________________________________________
292 
294  : G4TessellatedSolid(rhs), fNv(rhs.fNv), fNz(rhs.fNz),
295  fPolygon(rhs.fPolygon), fZSections(rhs.fZSections),
296  fTriangles(rhs.fTriangles), fIsConvex(rhs.fIsConvex),
297  fGeometryType(rhs.fGeometryType),
298  fSolidType(rhs.fSolidType), fPlanes(rhs.fPlanes),
299  fLines(rhs.fLines), fLengths(rhs.fLengths),
300  fKScales(rhs.fKScales), fScale0s(rhs.fScale0s),
301  fKOffsets(rhs.fKOffsets), fOffset0s(rhs.fOffset0s)
302 {
303 }
304 
305 //_____________________________________________________________________________
306 
308 {
309  // Check assignment to self
310  //
311  if (this == &rhs) { return *this; }
312 
313  // Copy base class data
314  //
316 
317  // Copy data
318  //
319  fNv = rhs.fNv; fNz = rhs.fNz;
320  fPolygon = rhs.fPolygon; fZSections = rhs.fZSections;
323  fSolidType = rhs.fSolidType; fPlanes = rhs.fPlanes;
324  fLines = rhs.fLines; fLengths = rhs.fLengths;
325  fKScales = rhs.fKScales; fScale0s = rhs.fScale0s;
326  fKOffsets = rhs.fKOffsets; fOffset0s = rhs.fOffset0s;
327 
328  return *this;
329 }
330 
331 //_____________________________________________________________________________
332 
334 {
335  // Destructor
336 }
337 
338 //_____________________________________________________________________________
339 
341 {
342  // Compute parameters for point projections p(z)
343  // to the polygon scale & offset:
344  // scale(z) = k*z + scale0
345  // offset(z) = l*z + offset0
346  // p(z) = scale(z)*p0 + offset(z)
347  // p0 = (p(z) - offset(z))/scale(z);
348  //
349 
350  for ( G4int iz=0; iz<fNz-1; ++iz)
351  {
352  G4double z1 = fZSections[iz].fZ;
353  G4double z2 = fZSections[iz+1].fZ;
354  G4double scale1 = fZSections[iz].fScale;
355  G4double scale2 = fZSections[iz+1].fScale;
356  G4TwoVector off1 = fZSections[iz].fOffset;
357  G4TwoVector off2 = fZSections[iz+1].fOffset;
358 
359  G4double kscale = (scale2 - scale1)/(z2 - z1);
360  G4double scale0 = scale2 - kscale*(z2 - z1)/2.0;
361  G4TwoVector koff = (off2 - off1)/(z2 - z1);
362  G4TwoVector off0 = off2 - koff*(z2 - z1)/2.0;
363 
364  fKScales.push_back(kscale);
365  fScale0s.push_back(scale0);
366  fKOffsets.push_back(koff);
367  fOffset0s.push_back(off0);
368  }
369 }
370 
371 //_____________________________________________________________________________
372 
374 {
375  // Compute lateral planes: a*x + b*y + c*z + d = 0
376  //
377  G4int Nv = fPolygon.size();
378  fPlanes.resize(Nv);
379  for (G4int i=0, k=Nv-1; i<Nv; k=i++)
380  {
381  G4TwoVector norm = (fPolygon[i] - fPolygon[k]).unit();
382  fPlanes[i].a = -norm.y();
383  fPlanes[i].b = norm.x();
384  fPlanes[i].c = 0;
385  fPlanes[i].d = norm.y()*fPolygon[i].x() - norm.x()*fPolygon[i].y();
386  }
387 
388  // Compute edge equations: x = k*y + m
389  // and edge lengths
390  //
391  fLines.resize(Nv);
392  fLengths.resize(Nv);
393  for (G4int i=0, k=Nv-1; i<Nv; k=i++)
394  {
395  if (fPolygon[k].y() == fPolygon[i].y())
396  {
397  fLines[i].k = 0;
398  fLines[i].m = fPolygon[i].x();
399  }
400  else
401  {
402  G4double ctg = (fPolygon[k].x()-fPolygon[i].x())/(fPolygon[k].y()-fPolygon[i].y());
403  fLines[i].k = ctg;
404  fLines[i].m = fPolygon[i].x() - ctg*fPolygon[i].y();
405  }
406  fLengths[i] = (fPolygon[i] - fPolygon[k]).mag();
407  }
408 }
409 
410 //_____________________________________________________________________________
411 
413 {
414  // Shift and scale vertices
415 
416  return G4ThreeVector( fPolygon[ind].x() * fZSections[iz].fScale
417  + fZSections[iz].fOffset.x(),
418  fPolygon[ind].y() * fZSections[iz].fScale
419  + fZSections[iz].fOffset.y(), fZSections[iz].fZ);
420 }
421 
422 //_____________________________________________________________________________
423 
425 {
426  // Project point in the polygon scale
427  // scale(z) = k*z + scale0
428  // offset(z) = l*z + offset0
429  // p(z) = scale(z)*p0 + offset(z)
430  // p0 = (p(z) - offset(z))/scale(z);
431 
432  // Select projection (z-segment of the solid) according to p.z()
433  //
434  G4int iz = 0;
435  while ( point.z() > fZSections[iz+1].fZ && iz < fNz-2 ) { ++iz; }
436  // Loop checking, 13.08.2015, G.Cosmo
437 
438  G4double z0 = ( fZSections[iz+1].fZ + fZSections[iz].fZ )/2.0;
439  G4TwoVector p2(point.x(), point.y());
440  G4double pscale = fKScales[iz]*(point.z()-z0) + fScale0s[iz];
441  G4TwoVector poffset = fKOffsets[iz]*(point.z()-z0) + fOffset0s[iz];
442 
443  // G4cout << point << " projected to "
444  // << iz << "-th z-segment polygon as "
445  // << (p2 - poffset)/pscale << G4endl;
446 
447  // pscale is always >0 as it is an interpolation between two
448  // positive scale values
449  //
450  return (p2 - poffset)/pscale;
451 }
452 
453 //_____________________________________________________________________________
454 
456  const G4TwoVector& l1,
457  const G4TwoVector& l2) const
458 {
459  // Return true if p is on the line through l1, l2
460 
461  if ( l1.x() == l2.x() )
462  {
463  return std::fabs(p.x() - l1.x()) < kCarToleranceHalf;
464  }
465  G4double slope= ((l2.y() - l1.y())/(l2.x() - l1.x()));
466  G4double predy= l1.y() + slope *(p.x() - l1.x());
467  G4double dy= p.y() - predy;
468 
469  // Calculate perpendicular distance
470  //
471  // G4double perpD= std::fabs(dy) / std::sqrt( 1 + slope * slope );
472  // G4bool simpleComp= (perpD<kCarToleranceHalf);
473 
474  // Check perpendicular distance vs tolerance 'directly'
475  //
476  G4bool squareComp = (dy*dy < (1+slope*slope)
478 
479  // return simpleComp;
480  return squareComp;
481 }
482 
483 //_____________________________________________________________________________
484 
486  const G4TwoVector& l1,
487  const G4TwoVector& l2) const
488 {
489  // Return true if p is on the line through l1, l2 and lies between
490  // l1 and l2
491 
492  if ( p.x() < std::min(l1.x(), l2.x()) - kCarToleranceHalf ||
493  p.x() > std::max(l1.x(), l2.x()) + kCarToleranceHalf ||
494  p.y() < std::min(l1.y(), l2.y()) - kCarToleranceHalf ||
495  p.y() > std::max(l1.y(), l2.y()) + kCarToleranceHalf )
496  {
497  return false;
498  }
499 
500  return IsSameLine(p, l1, l2);
501 }
502 
503 //_____________________________________________________________________________
504 
506  const G4TwoVector& p2,
507  const G4TwoVector& l1,
508  const G4TwoVector& l2) const
509 {
510  // Return true if p1 and p2 are on the same side of the line through l1, l2
511 
512  return ( (p1.x() - l1.x()) * (l2.y() - l1.y())
513  - (l2.x() - l1.x()) * (p1.y() - l1.y()) )
514  * ( (p2.x() - l1.x()) * (l2.y() - l1.y())
515  - (l2.x() - l1.x()) * (p2.y() - l1.y()) ) > 0;
516 }
517 
518 //_____________________________________________________________________________
519 
521  const G4TwoVector& b,
522  const G4TwoVector& c,
523  const G4TwoVector& p) const
524 {
525  // Return true if p is inside of triangle abc or on its edges,
526  // else returns false
527 
528  // Check extent first
529  //
530  if ( ( p.x() < a.x() && p.x() < b.x() && p.x() < c.x() ) ||
531  ( p.x() > a.x() && p.x() > b.x() && p.x() > c.x() ) ||
532  ( p.y() < a.y() && p.y() < b.y() && p.y() < c.y() ) ||
533  ( p.y() > a.y() && p.y() > b.y() && p.y() > c.y() ) ) return false;
534 
535  G4bool inside
536  = IsSameSide(p, a, b, c)
537  && IsSameSide(p, b, a, c)
538  && IsSameSide(p, c, a, b);
539 
540  G4bool onEdge
541  = IsSameLineSegment(p, a, b)
542  || IsSameLineSegment(p, b, c)
543  || IsSameLineSegment(p, c, a);
544 
545  return inside || onEdge;
546 }
547 
548 //_____________________________________________________________________________
549 
550 G4double
552  const G4TwoVector& pa,
553  const G4TwoVector& pb) const
554 {
555  // Return the angle of the vertex in po
556 
557  G4TwoVector t1 = pa - po;
558  G4TwoVector t2 = pb - po;
559 
560  G4double result = (std::atan2(t1.y(), t1.x()) - std::atan2(t2.y(), t2.x()));
561 
562  if ( result < 0 ) result += 2*pi;
563 
564  return result;
565 }
566 
567 //_____________________________________________________________________________
568 
569 G4VFacet*
571 {
572  // Create a triangular facet from the polygon points given by indices
573  // forming the down side ( the normal goes in -z)
574 
575  std::vector<G4ThreeVector> vertices;
576  vertices.push_back(GetVertex(0, ind1));
577  vertices.push_back(GetVertex(0, ind2));
578  vertices.push_back(GetVertex(0, ind3));
579 
580  // first vertex most left
581  //
583  = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
584 
585  if ( cross.z() > 0.0 )
586  {
587  // vertices ordered clock wise has to be reordered
588 
589  // G4cout << "G4ExtrudedSolid::MakeDownFacet: reordering vertices "
590  // << ind1 << ", " << ind2 << ", " << ind3 << G4endl;
591 
592  G4ThreeVector tmp = vertices[1];
593  vertices[1] = vertices[2];
594  vertices[2] = tmp;
595  }
596 
597  return new G4TriangularFacet(vertices[0], vertices[1],
598  vertices[2], ABSOLUTE);
599 }
600 
601 //_____________________________________________________________________________
602 
603 G4VFacet*
605 {
606  // Creates a triangular facet from the polygon points given by indices
607  // forming the upper side ( z>0 )
608 
609  std::vector<G4ThreeVector> vertices;
610  vertices.push_back(GetVertex(fNz-1, ind1));
611  vertices.push_back(GetVertex(fNz-1, ind2));
612  vertices.push_back(GetVertex(fNz-1, ind3));
613 
614  // first vertex most left
615  //
617  = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
618 
619  if ( cross.z() < 0.0 )
620  {
621  // vertices ordered clock wise has to be reordered
622 
623  // G4cout << "G4ExtrudedSolid::MakeUpFacet: reordering vertices "
624  // << ind1 << ", " << ind2 << ", " << ind3 << G4endl;
625 
626  G4ThreeVector tmp = vertices[1];
627  vertices[1] = vertices[2];
628  vertices[2] = tmp;
629  }
630 
631  return new G4TriangularFacet(vertices[0], vertices[1],
632  vertices[2], ABSOLUTE);
633 }
634 
635 //_____________________________________________________________________________
636 
638 {
639  // Decompose polygonal sides in triangular facets
640 
641  typedef std::pair < G4TwoVector, G4int > Vertex;
642 
643  static const G4double kAngTolerance =
645 
646  // Fill one more vector
647  //
648  std::vector< Vertex > verticesToBeDone;
649  for ( G4int i=0; i<fNv; ++i )
650  {
651  verticesToBeDone.push_back(Vertex(fPolygon[i], i));
652  }
653  std::vector< Vertex > ears;
654 
655  std::vector< Vertex >::iterator c1 = verticesToBeDone.begin();
656  std::vector< Vertex >::iterator c2 = c1+1;
657  std::vector< Vertex >::iterator c3 = c1+2;
658  while ( verticesToBeDone.size()>2 ) // Loop checking, 13.08.2015, G.Cosmo
659  {
660 
661  // G4cout << "Looking at triangle : "
662  // << c1->second << " " << c2->second
663  // << " " << c3->second << G4endl;
664  //G4cout << "Looking at triangle : "
665  // << c1->first << " " << c2->first
666  // << " " << c3->first << G4endl;
667 
668  // skip concave vertices
669  //
670  G4double angle = GetAngle(c2->first, c3->first, c1->first);
671 
672  //G4cout << "angle " << angle << G4endl;
673 
674  G4int counter = 0;
675  while ( angle >= (pi-kAngTolerance) ) // Loop checking, 13.08.2015, G.Cosmo
676  {
677  // G4cout << "Skipping concave vertex " << c2->second << G4endl;
678 
679  // try next three consecutive vertices
680  //
681  c1 = c2;
682  c2 = c3;
683  ++c3;
684  if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
685 
686  //G4cout << "Looking at triangle : "
687  // << c1->first << " " << c2->first
688  // << " " << c3->first << G4endl;
689 
690  angle = GetAngle(c2->first, c3->first, c1->first);
691  //G4cout << "angle " << angle << G4endl;
692 
693  ++counter;
694 
695  if ( counter > fNv)
696  {
697  G4Exception("G4ExtrudedSolid::AddGeneralPolygonFacets",
698  "GeomSolids0003", FatalException,
699  "Triangularisation has failed.");
700  break;
701  }
702  }
703 
704  G4bool good = true;
705  for ( auto it=verticesToBeDone.cbegin(); it!=verticesToBeDone.cend(); ++it )
706  {
707  // skip vertices of tested triangle
708  //
709  if ( it == c1 || it == c2 || it == c3 ) { continue; }
710 
711  if ( IsPointInside(c1->first, c2->first, c3->first, it->first) )
712  {
713  // G4cout << "Point " << it->second << " is inside" << G4endl;
714  good = false;
715 
716  // try next three consecutive vertices
717  //
718  c1 = c2;
719  c2 = c3;
720  ++c3;
721  if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
722  break;
723  }
724  // else
725  // { G4cout << "Point " << it->second << " is outside" << G4endl; }
726  }
727  if ( good )
728  {
729  // all points are outside triangle, we can make a facet
730 
731  // G4cout << "Found triangle : "
732  // << c1->second << " " << c2->second
733  // << " " << c3->second << G4endl;
734 
735  G4bool result;
736  result = AddFacet( MakeDownFacet(c1->second, c2->second, c3->second) );
737  if ( ! result ) { return false; }
738 
739  result = AddFacet( MakeUpFacet(c1->second, c2->second, c3->second) );
740  if ( ! result ) { return false; }
741 
742  std::vector<G4int> triangle(3);
743  triangle[0] = c1->second;
744  triangle[1] = c2->second;
745  triangle[2] = c3->second;
746  fTriangles.push_back(triangle);
747 
748  // remove the ear point from verticesToBeDone
749  //
750  verticesToBeDone.erase(c2);
751  c1 = verticesToBeDone.begin();
752  c2 = c1+1;
753  c3 = c1+2;
754  }
755  }
756  return true;
757 }
758 
759 //_____________________________________________________________________________
760 
762 {
763  // Define facets
764 
765  G4bool good;
766 
767  // Decomposition of polygonal sides in the facets
768  //
769  if ( fNv == 3 )
770  {
771  good = AddFacet( new G4TriangularFacet( GetVertex(0, 0), GetVertex(0, 1),
772  GetVertex(0, 2), ABSOLUTE) );
773  if ( ! good ) { return false; }
774 
775  good = AddFacet( new G4TriangularFacet( GetVertex(fNz-1, 2),
776  GetVertex(fNz-1, 1),
777  GetVertex(fNz-1, 0),
778  ABSOLUTE) );
779  if ( ! good ) { return false; }
780 
781  std::vector<G4int> triangle(3);
782  triangle[0] = 0;
783  triangle[1] = 1;
784  triangle[2] = 2;
785  fTriangles.push_back(triangle);
786  }
787 
788  else if ( fNv == 4 )
789  {
790  good = AddFacet( new G4QuadrangularFacet( GetVertex(0, 0),GetVertex(0, 1),
791  GetVertex(0, 2),GetVertex(0, 3),
792  ABSOLUTE) );
793  if ( ! good ) { return false; }
794 
795  good = AddFacet( new G4QuadrangularFacet( GetVertex(fNz-1, 3),
796  GetVertex(fNz-1, 2),
797  GetVertex(fNz-1, 1),
798  GetVertex(fNz-1, 0),
799  ABSOLUTE) );
800  if ( ! good ) { return false; }
801 
802  std::vector<G4int> triangle1(3);
803  triangle1[0] = 0;
804  triangle1[1] = 1;
805  triangle1[2] = 2;
806  fTriangles.push_back(triangle1);
807 
808  std::vector<G4int> triangle2(3);
809  triangle2[0] = 0;
810  triangle2[1] = 2;
811  triangle2[2] = 3;
812  fTriangles.push_back(triangle2);
813  }
814  else
815  {
816  good = AddGeneralPolygonFacets();
817  if ( ! good ) { return false; }
818  }
819 
820  // The quadrangular sides
821  //
822  for ( G4int iz = 0; iz < fNz-1; ++iz )
823  {
824  for ( G4int i = 0; i < fNv; ++i )
825  {
826  G4int j = (i+1) % fNv;
827  good = AddFacet( new G4QuadrangularFacet
828  ( GetVertex(iz, j), GetVertex(iz, i),
829  GetVertex(iz+1, i), GetVertex(iz+1, j), ABSOLUTE) );
830  if ( ! good ) { return false; }
831  }
832  }
833 
834  SetSolidClosed(true);
835 
836  return good;
837 }
838 
839 //_____________________________________________________________________________
840 
842 {
843  // Return entity type
844 
845  return fGeometryType;
846 }
847 
848 //_____________________________________________________________________________
849 
851 {
852  return new G4ExtrudedSolid(*this);
853 }
854 
855 //_____________________________________________________________________________
856 
858 {
859  switch (fSolidType)
860  {
861  case 1: // convex right prism
862  {
863  G4double dist = std::max(fZSections[0].fZ-p.z(),p.z()-fZSections[1].fZ);
864  if (dist > kCarToleranceHalf) { return kOutside; }
865 
866  G4int np = fPlanes.size();
867  for (G4int i=0; i<np; ++i)
868  {
869  G4double dd = fPlanes[i].a*p.x() + fPlanes[i].b*p.y() + fPlanes[i].d;
870  if (dd > dist) { dist = dd; }
871  }
872  if (dist > kCarToleranceHalf) { return kOutside; }
873  return (dist > -kCarToleranceHalf) ? kSurface : kInside;
874  }
875  case 2: // non-convex right prism
876  {
877  G4double distz = std::max(fZSections[0].fZ-p.z(),p.z()-fZSections[1].fZ);
878  if (distz > kCarToleranceHalf) { return kOutside; }
879 
880  G4bool in = PointInPolygon(p);
881  if (distz > -kCarToleranceHalf && in) { return kSurface; }
882 
884  if (in)
885  {
886  return (dd >= 0) ? kInside : kSurface;
887  }
888  else
889  {
890  return (dd > 0) ? kOutside : kSurface;
891  }
892  }
893  }
894 
895  // Override the base class function as it fails in case of concave polygon.
896  // Project the point in the original polygon scale and check if it is inside
897  // for each triangle.
898 
899  // Check first if outside extent
900  //
901  if ( p.x() < GetMinXExtent() - kCarToleranceHalf ||
902  p.x() > GetMaxXExtent() + kCarToleranceHalf ||
903  p.y() < GetMinYExtent() - kCarToleranceHalf ||
904  p.y() > GetMaxYExtent() + kCarToleranceHalf ||
905  p.z() < GetMinZExtent() - kCarToleranceHalf ||
906  p.z() > GetMaxZExtent() + kCarToleranceHalf )
907  {
908  // G4cout << "G4ExtrudedSolid::Outside extent: " << p << G4endl;
909  return kOutside;
910  }
911 
912  // Project point p(z) to the polygon scale p0
913  //
914  G4TwoVector pscaled = ProjectPoint(p);
915 
916  // Check if on surface of polygon
917  //
918  for ( G4int i=0; i<fNv; ++i )
919  {
920  G4int j = (i+1) % fNv;
921  if ( IsSameLineSegment(pscaled, fPolygon[i], fPolygon[j]) )
922  {
923  // G4cout << "G4ExtrudedSolid::Inside return Surface (on polygon) "
924  // << G4endl;
925 
926  return kSurface;
927  }
928  }
929 
930  // Now check if inside triangles
931  //
932  auto it = fTriangles.cbegin();
933  G4bool inside = false;
934  do // Loop checking, 13.08.2015, G.Cosmo
935  {
936  if ( IsPointInside(fPolygon[(*it)[0]], fPolygon[(*it)[1]],
937  fPolygon[(*it)[2]], pscaled) ) { inside = true; }
938  ++it;
939  } while ( (inside == false) && (it != fTriangles.cend()) );
940 
941  if ( inside )
942  {
943  // Check if on surface of z sides
944  //
945  if ( std::fabs( p.z() - fZSections[0].fZ ) < kCarToleranceHalf ||
946  std::fabs( p.z() - fZSections[fNz-1].fZ ) < kCarToleranceHalf )
947  {
948  // G4cout << "G4ExtrudedSolid::Inside return Surface (on z side)"
949  // << G4endl;
950 
951  return kSurface;
952  }
953 
954  // G4cout << "G4ExtrudedSolid::Inside return Inside" << G4endl;
955 
956  return kInside;
957  }
958 
959  // G4cout << "G4ExtrudedSolid::Inside return Outside " << G4endl;
960 
961  return kOutside;
962 }
963 
964 //_____________________________________________________________________________
965 
967 {
968  G4int nsurf = 0;
969  G4double nx = 0, ny = 0, nz = 0;
970  switch (fSolidType)
971  {
972  case 1: // convex right prism
973  {
974  if (std::abs(p.z() - fZSections[0].fZ) <= kCarToleranceHalf)
975  {
976  nz = -1; ++nsurf;
977  }
978  if (std::abs(p.z() - fZSections[1].fZ) <= kCarToleranceHalf)
979  {
980  nz = 1; ++nsurf;
981  }
982  for (G4int i=0; i<fNv; ++i)
983  {
984  G4double dd = fPlanes[i].a*p.x() + fPlanes[i].b*p.y() + fPlanes[i].d;
985  if (std::abs(dd) > kCarToleranceHalf) continue;
986  nx += fPlanes[i].a;
987  ny += fPlanes[i].b;
988  ++nsurf;
989  }
990  break;
991  }
992  case 2: // non-convex right prism
993  {
994  if (std::abs(p.z() - fZSections[0].fZ) <= kCarToleranceHalf)
995  {
996  nz = -1; ++nsurf;
997  }
998  if (std::abs(p.z() - fZSections[1].fZ) <= kCarToleranceHalf)
999  {
1000  nz = 1; ++nsurf;
1001  }
1002 
1003  G4double sqrCarToleranceHalf = kCarToleranceHalf*kCarToleranceHalf;
1004  for (G4int i=0, k=fNv-1; i<fNv; k=i++)
1005  {
1006  G4double ix = p.x() - fPolygon[i].x();
1007  G4double iy = p.y() - fPolygon[i].y();
1008  G4double u = fPlanes[i].a*iy - fPlanes[i].b*ix;
1009  if (u < 0)
1010  {
1011  if (ix*ix + iy*iy > sqrCarToleranceHalf) continue;
1012  }
1013  else if (u > fLengths[i])
1014  {
1015  G4double kx = p.x() - fPolygon[k].x();
1016  G4double ky = p.y() - fPolygon[k].y();
1017  if (kx*kx + ky*ky > sqrCarToleranceHalf) continue;
1018  }
1019  else
1020  {
1021  G4double dd = fPlanes[i].a*p.x() + fPlanes[i].b*p.y() + fPlanes[i].d;
1022  if (dd*dd > sqrCarToleranceHalf) continue;
1023  }
1024  nx += fPlanes[i].a;
1025  ny += fPlanes[i].b;
1026  ++nsurf;
1027  }
1028  break;
1029  }
1030  default:
1031  {
1033  }
1034  }
1035 
1036  // Return normal (right prism)
1037  //
1038  if (nsurf == 1)
1039  {
1040  return G4ThreeVector(nx,ny,nz);
1041  }
1042  else if (nsurf != 0) // edge or corner
1043  {
1044  return G4ThreeVector(nx,ny,nz).unit();
1045  }
1046  else
1047  {
1048  // Point is not on the surface, compute approximate normal
1049  //
1050 #ifdef G4CSGDEBUG
1051  std::ostringstream message;
1052  G4int oldprc = message.precision(16);
1053  message << "Point p is not on surface (!?) of solid: "
1054  << GetName() << G4endl;
1055  message << "Position:\n";
1056  message << " p.x() = " << p.x()/mm << " mm\n";
1057  message << " p.y() = " << p.y()/mm << " mm\n";
1058  message << " p.z() = " << p.z()/mm << " mm";
1059  G4cout.precision(oldprc) ;
1060  G4Exception("G4TesselatedSolid::SurfaceNormal(p)", "GeomSolids1002",
1061  JustWarning, message );
1062  DumpInfo();
1063 #endif
1064  return ApproxSurfaceNormal(p);
1065  }
1066 }
1067 
1068 //_____________________________________________________________________________
1069 
1071 {
1072  // This method is valid only for right prisms and
1073  // normally should not be called
1074 
1075  if (fSolidType == 1 || fSolidType == 2)
1076  {
1077  // Find distances to z-planes
1078  //
1079  G4double dz0 = fZSections[0].fZ - p.z();
1080  G4double dz1 = p.z() - fZSections[1].fZ;
1081  G4double ddz0 = dz0*dz0;
1082  G4double ddz1 = dz1*dz1;
1083 
1084  // Find nearest lateral side and distance to it
1085  //
1086  G4int iside = 0;
1087  G4double dd = DBL_MAX;
1088  for (G4int i=0, k=fNv-1; i<fNv; k=i++)
1089  {
1090  G4double ix = p.x() - fPolygon[i].x();
1091  G4double iy = p.y() - fPolygon[i].y();
1092  G4double u = fPlanes[i].a*iy - fPlanes[i].b*ix;
1093  if (u < 0)
1094  {
1095  G4double tmp = ix*ix + iy*iy;
1096  if (tmp < dd) { dd = tmp; iside = i; }
1097  }
1098  else if (u > fLengths[i])
1099  {
1100  G4double kx = p.x() - fPolygon[k].x();
1101  G4double ky = p.y() - fPolygon[k].y();
1102  G4double tmp = kx*kx + ky*ky;
1103  if (tmp < dd) { dd = tmp; iside = i; }
1104  }
1105  else
1106  {
1107  G4double tmp = fPlanes[i].a*p.x() + fPlanes[i].b*p.y() + fPlanes[i].d;
1108  tmp *= tmp;
1109  if (tmp < dd) { dd = tmp; iside = i; }
1110  }
1111  }
1112 
1113  // Find region
1114  //
1115  // 3 | 1 | 3
1116  // ----+-------+----
1117  // 2 | 0 | 2
1118  // ----+-------+----
1119  // 3 | 1 | 3
1120  //
1121  G4int iregion = 0;
1122  if (std::max(dz0,dz1) > 0) iregion = 1;
1123 
1124  G4bool in = PointInPolygon(p);
1125  if (!in) iregion += 2;
1126 
1127  // Return normal
1128  //
1129  switch (iregion)
1130  {
1131  case 0:
1132  {
1133  if (ddz0 <= ddz1 && ddz0 <= dd) return G4ThreeVector(0, 0,-1);
1134  if (ddz1 <= ddz0 && ddz1 <= dd) return G4ThreeVector(0, 0, 1);
1135  return G4ThreeVector(fPlanes[iside].a, fPlanes[iside].b, 0);
1136  }
1137  case 1:
1138  {
1139  return G4ThreeVector(0, 0, (dz0 > dz1) ? -1 : 1);
1140  }
1141  case 2:
1142  {
1143  return G4ThreeVector(fPlanes[iside].a, fPlanes[iside].b, 0);
1144  }
1145  case 3:
1146  {
1147  G4double dzmax = std::max(dz0,dz1);
1148  if (dzmax*dzmax > dd) return G4ThreeVector(0,0,(dz0 > dz1) ? -1 : 1);
1149  return G4ThreeVector(fPlanes[iside].a,fPlanes[iside].b, 0);
1150  }
1151  }
1152  }
1153  return G4ThreeVector(0,0,0);
1154 }
1155 
1156 //_____________________________________________________________________________
1157 
1159  const G4ThreeVector& v) const
1160 {
1161  G4double z0 = fZSections[0].fZ;
1162  G4double z1 = fZSections[fNz-1].fZ;
1163  if ((p.z() <= z0 + kCarToleranceHalf) && v.z() <= 0) return kInfinity;
1164  if ((p.z() >= z1 - kCarToleranceHalf) && v.z() >= 0) return kInfinity;
1165 
1166  switch (fSolidType)
1167  {
1168  case 1: // convex right prism
1169  {
1170  // Intersection with Z planes
1171  //
1172  G4double dz = (z1 - z0)*0.5;
1173  G4double pz = p.z() - dz - z0;
1174 
1175  G4double invz = (v.z() == 0) ? DBL_MAX : -1./v.z();
1176  G4double ddz = (invz < 0) ? dz : -dz;
1177  G4double tzmin = (pz + ddz)*invz;
1178  G4double tzmax = (pz - ddz)*invz;
1179 
1180  // Intersection with lateral planes
1181  //
1182  G4int np = fPlanes.size();
1183  G4double txmin = tzmin, txmax = tzmax;
1184  for (G4int i=0; i<np; ++i)
1185  {
1186  G4double cosa = fPlanes[i].a*v.x()+fPlanes[i].b*v.y();
1187  G4double dist = fPlanes[i].a*p.x()+fPlanes[i].b*p.y()+fPlanes[i].d;
1188  if (dist >= -kCarToleranceHalf)
1189  {
1190  if (cosa >= 0) { return kInfinity; }
1191  G4double tmp = -dist/cosa;
1192  if (txmin < tmp) { txmin = tmp; }
1193  }
1194  else if (cosa > 0)
1195  {
1196  G4double tmp = -dist/cosa;
1197  if (txmax > tmp) { txmax = tmp; }
1198  }
1199  }
1200 
1201  // Find distance
1202  //
1203  G4double tmin = txmin, tmax = txmax;
1204  if (tmax <= tmin + kCarToleranceHalf) // touch or no hit
1205  {
1206  return kInfinity;
1207  }
1208  return (tmin < kCarToleranceHalf) ? 0. : tmin;
1209  }
1210  case 2: // non-convex right prism
1211  {
1212  }
1213  }
1215 }
1216 
1217 //_____________________________________________________________________________
1218 
1220 {
1221  switch (fSolidType)
1222  {
1223  case 1: // convex right prism
1224  {
1225  G4double dist = std::max(fZSections[0].fZ-p.z(),p.z()-fZSections[1].fZ);
1226  G4int np = fPlanes.size();
1227  for (G4int i=0; i<np; ++i)
1228  {
1229  G4double dd = fPlanes[i].a*p.x() + fPlanes[i].b*p.y() + fPlanes[i].d;
1230  if (dd > dist) dist = dd;
1231  }
1232  return (dist > 0) ? dist : 0.;
1233  }
1234  case 2: // non-convex right prism
1235  {
1236  G4bool in = PointInPolygon(p);
1237  if (in)
1238  {
1239  G4double distz= std::max(fZSections[0].fZ-p.z(),p.z()-fZSections[1].fZ);
1240  return (distz > 0) ? distz : 0;
1241  }
1242  else
1243  {
1244  G4double distz= std::max(fZSections[0].fZ-p.z(),p.z()-fZSections[1].fZ);
1246  if (distz > 0) dd += distz*distz;
1247  return std::sqrt(dd);
1248  }
1249  }
1250  }
1251 
1252  // General case: use tessellated solid
1254 }
1255 
1256 //_____________________________________________________________________________
1257 
1259  const G4ThreeVector &v,
1260  const G4bool calcNorm,
1261  G4bool* validNorm,
1262  G4ThreeVector* n) const
1263 {
1264  G4bool getnorm = calcNorm;
1265  if (getnorm) *validNorm = true;
1266 
1267  G4double z0 = fZSections[0].fZ;
1268  G4double z1 = fZSections[fNz-1].fZ;
1269  if ((p.z() <= z0 + kCarToleranceHalf) && v.z() < 0)
1270  {
1271  if (getnorm) n->set(0,0,-1);
1272  return 0;
1273  }
1274  if ((p.z() >= z1 - kCarToleranceHalf) && v.z() > 0)
1275  {
1276  if (getnorm) n->set(0,0,1);
1277  return 0;
1278  }
1279 
1280  switch (fSolidType)
1281  {
1282  case 1: // convex right prism
1283  {
1284  // Intersection with Z planes
1285  //
1286  G4double dz = (z1 - z0)*0.5;
1287  G4double pz = p.z() - 0.5 * (z0 + z1);
1288 
1289  G4double vz = v.z();
1290  G4double tmax = (vz == 0) ? DBL_MAX : (std::copysign(dz,vz) - pz)/vz;
1291  G4int iside = (vz < 0) ? -4 : -2; // little trick: (-4+3)=-1, (-2+3)=+1
1292 
1293  // Intersection with lateral planes
1294  //
1295  G4int np = fPlanes.size();
1296  for (G4int i=0; i<np; ++i)
1297  {
1298  G4double cosa = fPlanes[i].a*v.x()+fPlanes[i].b*v.y();
1299  if (cosa > 0)
1300  {
1301  G4double dist = fPlanes[i].a*p.x()+fPlanes[i].b*p.y()+fPlanes[i].d;
1302  if (dist >= -kCarToleranceHalf)
1303  {
1304  if (getnorm) n->set(fPlanes[i].a, fPlanes[i].b, fPlanes[i].c);
1305  return 0;
1306  }
1307  G4double tmp = -dist/cosa;
1308  if (tmax > tmp) { tmax = tmp; iside = i; }
1309  }
1310  }
1311 
1312  // Set normal, if required, and return distance
1313  //
1314  if (getnorm)
1315  {
1316  if (iside < 0)
1317  { n->set(0, 0, iside + 3); } // (-4+3)=-1, (-2+3)=+1
1318  else
1319  { n->set(fPlanes[iside].a, fPlanes[iside].b, fPlanes[iside].c); }
1320  }
1321  return tmax;
1322  }
1323  case 2: // non-convex right prism
1324  {
1325  }
1326  }
1327 
1328  // Override the base class function to redefine validNorm
1329  // (the solid can be concave)
1330 
1331  G4double distOut =
1332  G4TessellatedSolid::DistanceToOut(p, v, calcNorm, validNorm, n);
1333  if (validNorm) { *validNorm = fIsConvex; }
1334 
1335  return distOut;
1336 }
1337 
1338 //_____________________________________________________________________________
1339 
1341 {
1342  switch (fSolidType)
1343  {
1344  case 1: // convex right prism
1345  {
1346  G4double dist = std::max(fZSections[0].fZ-p.z(),p.z()-fZSections[1].fZ);
1347  G4int np = fPlanes.size();
1348  for (G4int i=0; i<np; ++i)
1349  {
1350  G4double dd = fPlanes[i].a*p.x() + fPlanes[i].b*p.y() + fPlanes[i].d;
1351  if (dd > dist) dist = dd;
1352  }
1353  return (dist < 0) ? -dist : 0.;
1354  }
1355  case 2: // non-convex right prism
1356  {
1357  G4double distz = std::max(fZSections[0].fZ-p.z(),p.z()-fZSections[1].fZ);
1358  G4bool in = PointInPolygon(p);
1359  if (distz >= 0 || (!in)) return 0; // point is outside
1360  return std::min(-distz,std::sqrt(DistanceToPolygonSqr(p)));
1361  }
1362  }
1363 
1364  // General case: use tessellated solid
1366 }
1367 
1368 //_____________________________________________________________________________
1369 // Get bounding box
1370 
1372  G4ThreeVector& pMax) const
1373 {
1374  G4double xmin0 = kInfinity, xmax0 = -kInfinity;
1375  G4double ymin0 = kInfinity, ymax0 = -kInfinity;
1376 
1377  for (G4int i=0; i<GetNofVertices(); ++i)
1378  {
1379  G4double x = fPolygon[i].x();
1380  if (x < xmin0) xmin0 = x;
1381  if (x > xmax0) xmax0 = x;
1382  G4double y = fPolygon[i].y();
1383  if (y < ymin0) ymin0 = y;
1384  if (y > ymax0) ymax0 = y;
1385  }
1386 
1389 
1390  G4int nsect = GetNofZSections();
1391  for (G4int i=0; i<nsect; ++i)
1392  {
1393  ZSection zsect = GetZSection(i);
1394  G4double dx = zsect.fOffset.x();
1395  G4double dy = zsect.fOffset.y();
1396  G4double scale = zsect.fScale;
1397  xmin = std::min(xmin,xmin0*scale+dx);
1398  xmax = std::max(xmax,xmax0*scale+dx);
1399  ymin = std::min(ymin,ymin0*scale+dy);
1400  ymax = std::max(ymax,ymax0*scale+dy);
1401  }
1402 
1403  G4double zmin = GetZSection(0).fZ;
1404  G4double zmax = GetZSection(nsect-1).fZ;
1405 
1406  pMin.set(xmin,ymin,zmin);
1407  pMax.set(xmax,ymax,zmax);
1408 
1409  // Check correctness of the bounding box
1410  //
1411  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
1412  {
1413  std::ostringstream message;
1414  message << "Bad bounding box (min >= max) for solid: "
1415  << GetName() << " !"
1416  << "\npMin = " << pMin
1417  << "\npMax = " << pMax;
1418  G4Exception("G4ExtrudedSolid::BoundingLimits()",
1419  "GeomMgt0001", JustWarning, message);
1420  DumpInfo();
1421  }
1422 }
1423 
1424 //_____________________________________________________________________________
1425 // Calculate extent under transform and specified limit
1426 
1427 G4bool
1429  const G4VoxelLimits& pVoxelLimit,
1430  const G4AffineTransform& pTransform,
1431  G4double& pMin, G4double& pMax) const
1432 {
1433  G4ThreeVector bmin, bmax;
1434  G4bool exist;
1435 
1436  // Check bounding box (bbox)
1437  //
1438  BoundingLimits(bmin,bmax);
1439  G4BoundingEnvelope bbox(bmin,bmax);
1440 #ifdef G4BBOX_EXTENT
1441  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
1442 #endif
1443  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
1444  {
1445  return exist = (pMin < pMax) ? true : false;
1446  }
1447 
1448  // To find the extent, the base polygon is subdivided in triangles.
1449  // The extent is calculated as cumulative extent of the parts
1450  // formed by extrusion of the triangles
1451  //
1452  G4TwoVectorList triangles;
1453  G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
1454  G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
1455 
1456  // triangulate the base polygon
1457  if (!G4GeomTools::TriangulatePolygon(fPolygon,triangles))
1458  {
1459  std::ostringstream message;
1460  message << "Triangulation of the base polygon has failed for solid: "
1461  << GetName() << " !"
1462  << "\nExtent has been calculated using boundary box";
1463  G4Exception("G4ExtrudedSolid::CalculateExtent()",
1464  "GeomMgt1002",JustWarning,message);
1465  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
1466  }
1467 
1468  // allocate vector lists
1469  G4int nsect = GetNofZSections();
1470  std::vector<const G4ThreeVectorList *> polygons;
1471  polygons.resize(nsect);
1472  for (G4int k=0; k<nsect; ++k) { polygons[k] = new G4ThreeVectorList(3); }
1473 
1474  // main loop along triangles
1475  pMin = kInfinity;
1476  pMax = -kInfinity;
1477  G4int ntria = triangles.size()/3;
1478  for (G4int i=0; i<ntria; ++i)
1479  {
1480  G4int i3 = i*3;
1481  for (G4int k=0; k<nsect; ++k) // extrude triangle
1482  {
1483  ZSection zsect = GetZSection(k);
1484  G4double z = zsect.fZ;
1485  G4double dx = zsect.fOffset.x();
1486  G4double dy = zsect.fOffset.y();
1487  G4double scale = zsect.fScale;
1488 
1489  G4ThreeVectorList* ptr = const_cast<G4ThreeVectorList*>(polygons[k]);
1490  G4ThreeVectorList::iterator iter = ptr->begin();
1491  G4double x0 = triangles[i3+0].x()*scale+dx;
1492  G4double y0 = triangles[i3+0].y()*scale+dy;
1493  iter->set(x0,y0,z);
1494  iter++;
1495  G4double x1 = triangles[i3+1].x()*scale+dx;
1496  G4double y1 = triangles[i3+1].y()*scale+dy;
1497  iter->set(x1,y1,z);
1498  iter++;
1499  G4double x2 = triangles[i3+2].x()*scale+dx;
1500  G4double y2 = triangles[i3+2].y()*scale+dy;
1501  iter->set(x2,y2,z);
1502  }
1503 
1504  // set sub-envelope and adjust extent
1505  G4double emin,emax;
1506  G4BoundingEnvelope benv(polygons);
1507  if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
1508  if (emin < pMin) pMin = emin;
1509  if (emax > pMax) pMax = emax;
1510  if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
1511  }
1512  // free memory
1513  for (G4int k=0; k<nsect; ++k) { delete polygons[k]; polygons[k]=0;}
1514  return (pMin < pMax);
1515 }
1516 
1517 //_____________________________________________________________________________
1518 
1519 std::ostream& G4ExtrudedSolid::StreamInfo(std::ostream &os) const
1520 {
1521  G4int oldprc = os.precision(16);
1522  os << "-----------------------------------------------------------\n"
1523  << " *** Dump for solid - " << GetName() << " ***\n"
1524  << " ===================================================\n"
1525  << " Solid geometry type: " << fGeometryType << G4endl;
1526 
1527  if ( fIsConvex)
1528  { os << " Convex polygon; list of vertices:" << G4endl; }
1529  else
1530  { os << " Concave polygon; list of vertices:" << G4endl; }
1531 
1532  for ( G4int i=0; i<fNv; ++i )
1533  {
1534  os << std::setw(5) << "#" << i
1535  << " vx = " << fPolygon[i].x()/mm << " mm"
1536  << " vy = " << fPolygon[i].y()/mm << " mm" << G4endl;
1537  }
1538 
1539  os << " Sections:" << G4endl;
1540  for ( G4int iz=0; iz<fNz; ++iz )
1541  {
1542  os << " z = " << fZSections[iz].fZ/mm << " mm "
1543  << " x0= " << fZSections[iz].fOffset.x()/mm << " mm "
1544  << " y0= " << fZSections[iz].fOffset.y()/mm << " mm "
1545  << " scale= " << fZSections[iz].fScale << G4endl;
1546  }
1547 
1548 /*
1549  // Triangles (for debugging)
1550  os << G4endl;
1551  os << " Triangles:" << G4endl;
1552  os << " Triangle # vertex1 vertex2 vertex3" << G4endl;
1553 
1554  G4int counter = 0;
1555  std::vector< std::vector<G4int> >::const_iterator it;
1556  for ( it = fTriangles.begin(); it != fTriangles.end(); it++ ) {
1557  std::vector<G4int> triangle = *it;
1558  os << std::setw(10) << counter++
1559  << std::setw(10) << triangle[0] << std::setw(10) << triangle[1]
1560  << std::setw(10) << triangle[2]
1561  << G4endl;
1562  }
1563 */
1564  os.precision(oldprc);
1565 
1566  return os;
1567 }
1568 
1569 #endif