Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ExtrudedSolid.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4ExtrudedSolid.cc 101118 2016-11-07 09:10:59Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class source file
32 //
33 // G4ExtrudedSolid.cc
34 //
35 // Author: Ivana Hrivnacova, IPN Orsay
36 //
37 // CHANGE HISTORY
38 // --------------
39 //
40 // 21.10.2016 E.Tcherniaev: added Extent() and CalculateExtent(),
41 // used G4GeomTools::PolygonArea() to calculate area,
42 // replaced IsConvex() with G4GeomTools::IsConvex()
43 // 02.03.2016 E.Tcherniaev: added CheckPolygon() to remove
44 // collinear and coincident points from polygon
45 // --------------------------------------------------------------------
46 
47 #include "G4ExtrudedSolid.hh"
48 
49 #if !defined(G4GEOM_USE_UEXTRUDEDSOLID)
50 
51 #include <set>
52 #include <algorithm>
53 #include <cmath>
54 #include <iomanip>
55 
56 #include "G4GeomTools.hh"
57 #include "G4VoxelLimits.hh"
58 #include "G4AffineTransform.hh"
59 #include "G4BoundingEnvelope.hh"
60 
61 #include "G4GeometryTolerance.hh"
62 #include "G4PhysicalConstants.hh"
63 #include "G4SystemOfUnits.hh"
64 #include "G4VFacet.hh"
65 #include "G4TriangularFacet.hh"
66 #include "G4QuadrangularFacet.hh"
67 
68 //_____________________________________________________________________________
69 
71  const std::vector<G4TwoVector>& polygon,
72  const std::vector<ZSection>& zsections)
73  : G4TessellatedSolid(pName),
74  fNv(polygon.size()),
75  fNz(zsections.size()),
76  fPolygon(),
77  fZSections(),
78  fTriangles(),
79  fIsConvex(false),
80  fGeometryType("G4ExtrudedSolid")
81 
82 {
83  // General constructor
84 
85  // First check input parameters
86 
87  if (fNv < 3)
88  {
89  std::ostringstream message;
90  message << "Number of vertices in polygon < 3 - " << pName;
91  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
92  FatalErrorInArgument, message);
93  }
94 
95  if (fNz < 2)
96  {
97  std::ostringstream message;
98  message << "Number of z-sides < 2 - " << pName;
99  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
100  FatalErrorInArgument, message);
101  }
102 
103  for ( G4int i=0; i<fNz-1; ++i )
104  {
105  if ( zsections[i].fZ > zsections[i+1].fZ )
106  {
107  std::ostringstream message;
108  message << "Z-sections have to be ordered by z value (z0 < z1 < z2...) - "
109  << pName;
110  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
111  FatalErrorInArgument, message);
112  }
113  if ( std::fabs( zsections[i+1].fZ - zsections[i].fZ ) < kCarToleranceHalf )
114  {
115  std::ostringstream message;
116  message << "Z-sections with the same z position are not supported - "
117  << pName;
118  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0001",
119  FatalException, message);
120  }
121  }
122 
123  // Copy polygon
124  //
125  fPolygon = polygon;
126 
127  // Remove collinear and coincident vertices, if any
128  //
129  std::vector<G4int> removedVertices;
130  G4GeomTools::RemoveRedundantVertices(fPolygon,removedVertices,
131  2*kCarTolerance);
132  if (removedVertices.size() != 0)
133  {
134  G4int nremoved = removedVertices.size();
135  std::ostringstream message;
136  message << "The following "<< nremoved
137  << " vertices have been removed from polygon in " << pName
138  << "\nas collinear or coincident with other vertices: "
139  << removedVertices[0];
140  for (G4int i=1; i<nremoved; ++i) message << ", " << removedVertices[i];
141  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
142  JustWarning, message);
143  }
144 
145  fNv = fPolygon.size();
146  if (fNv < 3)
147  {
148  std::ostringstream message;
149  message << "Number of vertices in polygon after removal < 3 - " << pName;
150  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
151  FatalErrorInArgument, message);
152  }
153 
154  // Check if polygon vertices are defined clockwise
155  // (the area is positive if polygon vertices are defined anti-clockwise)
156  //
157  if (G4GeomTools::PolygonArea(fPolygon) > 0.)
158  {
159  // Polygon vertices are defined anti-clockwise, we revert them
160  // G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
161  // JustWarning,
162  // "Polygon vertices defined anti-clockwise, reverting polygon");
163  std::reverse(fPolygon.begin(),fPolygon.end());
164  }
165 
166  // Copy z-sections
167  //
168  fZSections = zsections;
169 
170  G4bool result = MakeFacets();
171  if (!result)
172  {
173  std::ostringstream message;
174  message << "Making facets failed - " << pName;
175  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0003",
176  FatalException, message);
177  }
178  fIsConvex = G4GeomTools::IsConvex(fPolygon);
179 
180  ComputeProjectionParameters();
181 }
182 
183 //_____________________________________________________________________________
184 
186  const std::vector<G4TwoVector>& polygon,
187  G4double dz,
188  const G4TwoVector& off1, G4double scale1,
189  const G4TwoVector& off2, G4double scale2 )
190  : G4TessellatedSolid(pName),
191  fNv(polygon.size()),
192  fNz(2),
193  fPolygon(),
194  fZSections(),
195  fTriangles(),
196  fIsConvex(false),
197  fGeometryType("G4ExtrudedSolid")
198 
199 {
200  // Special constructor for solid with 2 z-sections
201 
202  // First check input parameters
203  //
204  if (fNv < 3)
205  {
206  std::ostringstream message;
207  message << "Number of vertices in polygon < 3 - " << pName;
208  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
209  FatalErrorInArgument, message);
210  }
211 
212  // Copy polygon
213  //
214  fPolygon = polygon;
215 
216  // Remove collinear and coincident vertices, if any
217  //
218  std::vector<G4int> removedVertices;
219  G4GeomTools::RemoveRedundantVertices(fPolygon,removedVertices,
220  2*kCarTolerance);
221  if (removedVertices.size() != 0)
222  {
223  G4int nremoved = removedVertices.size();
224  std::ostringstream message;
225  message << "The following "<< nremoved
226  << " vertices have been removed from polygon in " << pName
227  << "\nas collinear or coincident with other vertices: "
228  << removedVertices[0];
229  for (G4int i=1; i<nremoved; ++i) message << ", " << removedVertices[i];
230  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
231  JustWarning, message);
232  }
233 
234  fNv = fPolygon.size();
235  if (fNv < 3)
236  {
237  std::ostringstream message;
238  message << "Number of vertices in polygon after removal < 3 - " << pName;
239  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
240  FatalErrorInArgument, message);
241  }
242 
243  // Check if polygon vertices are defined clockwise
244  // (the area is positive if polygon vertices are defined anti-clockwise)
245  //
246  if (G4GeomTools::PolygonArea(fPolygon) > 0.)
247  {
248  // Polygon vertices are defined anti-clockwise, we revert them
249  // G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
250  // JustWarning,
251  // "Polygon vertices defined anti-clockwise, reverting polygon");
252  std::reverse(fPolygon.begin(),fPolygon.end());
253  }
254 
255  // Copy z-sections
256  //
257  fZSections.push_back(ZSection(-dz, off1, scale1));
258  fZSections.push_back(ZSection( dz, off2, scale2));
259 
260  G4bool result = MakeFacets();
261  if (!result)
262  {
263  std::ostringstream message;
264  message << "Making facets failed - " << pName;
265  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0003",
266  FatalException, message);
267  }
268  fIsConvex = G4GeomTools::IsConvex(fPolygon);
269 
270  ComputeProjectionParameters();
271 }
272 
273 //_____________________________________________________________________________
274 
276  : G4TessellatedSolid(a), fNv(0), fNz(0), fPolygon(), fZSections(),
277  fTriangles(), fIsConvex(false), fGeometryType("G4ExtrudedSolid")
278 {
279  // Fake default constructor - sets only member data and allocates memory
280  // for usage restricted to object persistency.
281 }
282 
283 //_____________________________________________________________________________
284 
286  : G4TessellatedSolid(rhs), fNv(rhs.fNv), fNz(rhs.fNz),
287  fPolygon(rhs.fPolygon), fZSections(rhs.fZSections),
288  fTriangles(rhs.fTriangles), fIsConvex(rhs.fIsConvex),
289  fGeometryType(rhs.fGeometryType), fKScales(rhs.fKScales),
290  fScale0s(rhs.fScale0s), fKOffsets(rhs.fKOffsets), fOffset0s(rhs.fOffset0s)
291 {
292 }
293 
294 
295 //_____________________________________________________________________________
296 
298 {
299  // Check assignment to self
300  //
301  if (this == &rhs) { return *this; }
302 
303  // Copy base class data
304  //
306 
307  // Copy data
308  //
309  fNv = rhs.fNv; fNz = rhs.fNz;
310  fPolygon = rhs.fPolygon; fZSections = rhs.fZSections;
311  fTriangles = rhs.fTriangles; fIsConvex = rhs.fIsConvex;
312  fGeometryType = rhs.fGeometryType; fKScales = rhs.fKScales;
313  fScale0s = rhs.fScale0s; fKOffsets = rhs.fKOffsets;
314  fOffset0s = rhs.fOffset0s;
315 
316  return *this;
317 }
318 
319 //_____________________________________________________________________________
320 
322 {
323  // Destructor
324 }
325 
326 //_____________________________________________________________________________
327 
328 void G4ExtrudedSolid::ComputeProjectionParameters()
329 {
330  // Compute parameters for point projections p(z)
331  // to the polygon scale & offset:
332  // scale(z) = k*z + scale0
333  // offset(z) = l*z + offset0
334  // p(z) = scale(z)*p0 + offset(z)
335  // p0 = (p(z) - offset(z))/scale(z);
336  //
337 
338  for ( G4int iz=0; iz<fNz-1; ++iz)
339  {
340  G4double z1 = fZSections[iz].fZ;
341  G4double z2 = fZSections[iz+1].fZ;
342  G4double scale1 = fZSections[iz].fScale;
343  G4double scale2 = fZSections[iz+1].fScale;
344  G4TwoVector off1 = fZSections[iz].fOffset;
345  G4TwoVector off2 = fZSections[iz+1].fOffset;
346 
347  G4double kscale = (scale2 - scale1)/(z2 - z1);
348  G4double scale0 = scale2 - kscale*(z2 - z1)/2.0;
349  G4TwoVector koff = (off2 - off1)/(z2 - z1);
350  G4TwoVector off0 = off2 - koff*(z2 - z1)/2.0;
351 
352  fKScales.push_back(kscale);
353  fScale0s.push_back(scale0);
354  fKOffsets.push_back(koff);
355  fOffset0s.push_back(off0);
356  }
357 }
358 
359 
360 //_____________________________________________________________________________
361 
363 {
364  // Shift and scale vertices
365 
366  return G4ThreeVector( fPolygon[ind].x() * fZSections[iz].fScale
367  + fZSections[iz].fOffset.x(),
368  fPolygon[ind].y() * fZSections[iz].fScale
369  + fZSections[iz].fOffset.y(), fZSections[iz].fZ);
370 }
371 
372 //_____________________________________________________________________________
373 
374 
375 G4TwoVector G4ExtrudedSolid::ProjectPoint(const G4ThreeVector& point) const
376 {
377  // Project point in the polygon scale
378  // scale(z) = k*z + scale0
379  // offset(z) = l*z + offset0
380  // p(z) = scale(z)*p0 + offset(z)
381  // p0 = (p(z) - offset(z))/scale(z);
382 
383  // Select projection (z-segment of the solid) according to p.z()
384  //
385  G4int iz = 0;
386  while ( point.z() > fZSections[iz+1].fZ && iz < fNz-2 ) { ++iz; }
387  // Loop checking, 13.08.2015, G.Cosmo
388 
389  G4double z0 = ( fZSections[iz+1].fZ + fZSections[iz].fZ )/2.0;
390  G4TwoVector p2(point.x(), point.y());
391  G4double pscale = fKScales[iz]*(point.z()-z0) + fScale0s[iz];
392  G4TwoVector poffset = fKOffsets[iz]*(point.z()-z0) + fOffset0s[iz];
393 
394  // G4cout << point << " projected to "
395  // << iz << "-th z-segment polygon as "
396  // << (p2 - poffset)/pscale << G4endl;
397 
398  // pscale is always >0 as it is an interpolation between two
399  // positive scale values
400  //
401  return (p2 - poffset)/pscale;
402 }
403 
404 //_____________________________________________________________________________
405 
406 G4bool G4ExtrudedSolid::IsSameLine(const G4TwoVector& p,
407  const G4TwoVector& l1,
408  const G4TwoVector& l2) const
409 {
410  // Return true if p is on the line through l1, l2
411 
412  if ( l1.x() == l2.x() )
413  {
414  return std::fabs(p.x() - l1.x()) < kCarToleranceHalf;
415  }
416  G4double slope= ((l2.y() - l1.y())/(l2.x() - l1.x()));
417  G4double predy= l1.y() + slope *(p.x() - l1.x());
418  G4double dy= p.y() - predy;
419 
420  // Calculate perpendicular distance
421  //
422  // G4double perpD= std::fabs(dy) / std::sqrt( 1 + slope * slope );
423  // G4bool simpleComp= (perpD<kCarToleranceHalf);
424 
425  // Check perpendicular distance vs tolerance 'directly'
426  //
427  G4bool squareComp = (dy*dy < (1+slope*slope)
429 
430  // return simpleComp;
431  return squareComp;
432 }
433 
434 //_____________________________________________________________________________
435 
436 G4bool G4ExtrudedSolid::IsSameLineSegment(const G4TwoVector& p,
437  const G4TwoVector& l1,
438  const G4TwoVector& l2) const
439 {
440  // Return true if p is on the line through l1, l2 and lies between
441  // l1 and l2
442 
443  if ( p.x() < std::min(l1.x(), l2.x()) - kCarToleranceHalf ||
444  p.x() > std::max(l1.x(), l2.x()) + kCarToleranceHalf ||
445  p.y() < std::min(l1.y(), l2.y()) - kCarToleranceHalf ||
446  p.y() > std::max(l1.y(), l2.y()) + kCarToleranceHalf )
447  {
448  return false;
449  }
450 
451  return IsSameLine(p, l1, l2);
452 }
453 
454 //_____________________________________________________________________________
455 
456 G4bool G4ExtrudedSolid::IsSameSide(const G4TwoVector& p1,
457  const G4TwoVector& p2,
458  const G4TwoVector& l1,
459  const G4TwoVector& l2) const
460 {
461  // Return true if p1 and p2 are on the same side of the line through l1, l2
462 
463  return ( (p1.x() - l1.x()) * (l2.y() - l1.y())
464  - (l2.x() - l1.x()) * (p1.y() - l1.y()) )
465  * ( (p2.x() - l1.x()) * (l2.y() - l1.y())
466  - (l2.x() - l1.x()) * (p2.y() - l1.y()) ) > 0;
467 }
468 
469 //_____________________________________________________________________________
470 
471 G4bool G4ExtrudedSolid::IsPointInside(const G4TwoVector& a,
472  const G4TwoVector& b,
473  const G4TwoVector& c,
474  const G4TwoVector& p) const
475 {
476  // Return true if p is inside of triangle abc or on its edges,
477  // else returns false
478 
479  // Check extent first
480  //
481  if ( ( p.x() < a.x() && p.x() < b.x() && p.x() < c.x() ) ||
482  ( p.x() > a.x() && p.x() > b.x() && p.x() > c.x() ) ||
483  ( p.y() < a.y() && p.y() < b.y() && p.y() < c.y() ) ||
484  ( p.y() > a.y() && p.y() > b.y() && p.y() > c.y() ) ) return false;
485 
486  G4bool inside
487  = IsSameSide(p, a, b, c)
488  && IsSameSide(p, b, a, c)
489  && IsSameSide(p, c, a, b);
490 
491  G4bool onEdge
492  = IsSameLineSegment(p, a, b)
493  || IsSameLineSegment(p, b, c)
494  || IsSameLineSegment(p, c, a);
495 
496  return inside || onEdge;
497 }
498 
499 //_____________________________________________________________________________
500 
501 G4double
502 G4ExtrudedSolid::GetAngle(const G4TwoVector& po,
503  const G4TwoVector& pa,
504  const G4TwoVector& pb) const
505 {
506  // Return the angle of the vertex in po
507 
508  G4TwoVector t1 = pa - po;
509  G4TwoVector t2 = pb - po;
510 
511  G4double result = (std::atan2(t1.y(), t1.x()) - std::atan2(t2.y(), t2.x()));
512 
513  if ( result < 0 ) result += 2*pi;
514 
515  return result;
516 }
517 
518 //_____________________________________________________________________________
519 
520 G4VFacet*
521 G4ExtrudedSolid::MakeDownFacet(G4int ind1, G4int ind2, G4int ind3) const
522 {
523  // Create a triangular facet from the polygon points given by indices
524  // forming the down side ( the normal goes in -z)
525 
526  std::vector<G4ThreeVector> vertices;
527  vertices.push_back(GetVertex(0, ind1));
528  vertices.push_back(GetVertex(0, ind2));
529  vertices.push_back(GetVertex(0, ind3));
530 
531  // first vertex most left
532  //
533  G4ThreeVector cross
534  = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
535 
536  if ( cross.z() > 0.0 )
537  {
538  // vertices ordered clock wise has to be reordered
539 
540  // G4cout << "G4ExtrudedSolid::MakeDownFacet: reordering vertices "
541  // << ind1 << ", " << ind2 << ", " << ind3 << G4endl;
542 
543  G4ThreeVector tmp = vertices[1];
544  vertices[1] = vertices[2];
545  vertices[2] = tmp;
546  }
547 
548  return new G4TriangularFacet(vertices[0], vertices[1],
549  vertices[2], ABSOLUTE);
550 }
551 
552 //_____________________________________________________________________________
553 
554 G4VFacet*
555 G4ExtrudedSolid::MakeUpFacet(G4int ind1, G4int ind2, G4int ind3) const
556 {
557  // Creates a triangular facet from the polygon points given by indices
558  // forming the upper side ( z>0 )
559 
560  std::vector<G4ThreeVector> vertices;
561  vertices.push_back(GetVertex(fNz-1, ind1));
562  vertices.push_back(GetVertex(fNz-1, ind2));
563  vertices.push_back(GetVertex(fNz-1, ind3));
564 
565  // first vertex most left
566  //
567  G4ThreeVector cross
568  = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
569 
570  if ( cross.z() < 0.0 )
571  {
572  // vertices ordered clock wise has to be reordered
573 
574  // G4cout << "G4ExtrudedSolid::MakeUpFacet: reordering vertices "
575  // << ind1 << ", " << ind2 << ", " << ind3 << G4endl;
576 
577  G4ThreeVector tmp = vertices[1];
578  vertices[1] = vertices[2];
579  vertices[2] = tmp;
580  }
581 
582  return new G4TriangularFacet(vertices[0], vertices[1],
583  vertices[2], ABSOLUTE);
584 }
585 
586 //_____________________________________________________________________________
587 
588 G4bool G4ExtrudedSolid::AddGeneralPolygonFacets()
589 {
590  // Decompose polygonal sides in triangular facets
591 
592  typedef std::pair < G4TwoVector, G4int > Vertex;
593 
594  static const G4double kAngTolerance =
596 
597  // Fill one more vector
598  //
599  std::vector< Vertex > verticesToBeDone;
600  for ( G4int i=0; i<fNv; ++i )
601  {
602  verticesToBeDone.push_back(Vertex(fPolygon[i], i));
603  }
604  std::vector< Vertex > ears;
605 
606  std::vector< Vertex >::iterator c1 = verticesToBeDone.begin();
607  std::vector< Vertex >::iterator c2 = c1+1;
608  std::vector< Vertex >::iterator c3 = c1+2;
609  while ( verticesToBeDone.size()>2 ) // Loop checking, 13.08.2015, G.Cosmo
610  {
611 
612  // G4cout << "Looking at triangle : "
613  // << c1->second << " " << c2->second
614  // << " " << c3->second << G4endl;
615  //G4cout << "Looking at triangle : "
616  // << c1->first << " " << c2->first
617  // << " " << c3->first << G4endl;
618 
619  // skip concave vertices
620  //
621  G4double angle = GetAngle(c2->first, c3->first, c1->first);
622 
623  //G4cout << "angle " << angle << G4endl;
624 
625  G4int counter = 0;
626  while ( angle >= (pi-kAngTolerance) ) // Loop checking, 13.08.2015, G.Cosmo
627  {
628  // G4cout << "Skipping concave vertex " << c2->second << G4endl;
629 
630  // try next three consecutive vertices
631  //
632  c1 = c2;
633  c2 = c3;
634  ++c3;
635  if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
636 
637  //G4cout << "Looking at triangle : "
638  // << c1->first << " " << c2->first
639  // << " " << c3->first << G4endl;
640 
641  angle = GetAngle(c2->first, c3->first, c1->first);
642  //G4cout << "angle " << angle << G4endl;
643 
644  counter++;
645 
646  if ( counter > fNv) {
647  G4Exception("G4ExtrudedSolid::AddGeneralPolygonFacets",
648  "GeomSolids0003", FatalException,
649  "Triangularisation has failed.");
650  break;
651  }
652  }
653 
654  G4bool good = true;
655  std::vector< Vertex >::iterator it;
656  for ( it=verticesToBeDone.begin(); it != verticesToBeDone.end(); ++it )
657  {
658  // skip vertices of tested triangle
659  //
660  if ( it == c1 || it == c2 || it == c3 ) { continue; }
661 
662  if ( IsPointInside(c1->first, c2->first, c3->first, it->first) )
663  {
664  // G4cout << "Point " << it->second << " is inside" << G4endl;
665  good = false;
666 
667  // try next three consecutive vertices
668  //
669  c1 = c2;
670  c2 = c3;
671  ++c3;
672  if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
673  break;
674  }
675  // else
676  // { G4cout << "Point " << it->second << " is outside" << G4endl; }
677  }
678  if ( good )
679  {
680  // all points are outside triangle, we can make a facet
681 
682  // G4cout << "Found triangle : "
683  // << c1->second << " " << c2->second
684  // << " " << c3->second << G4endl;
685 
686  G4bool result;
687  result = AddFacet( MakeDownFacet(c1->second, c2->second, c3->second) );
688  if ( ! result ) { return false; }
689 
690  result = AddFacet( MakeUpFacet(c1->second, c2->second, c3->second) );
691  if ( ! result ) { return false; }
692 
693  std::vector<G4int> triangle(3);
694  triangle[0] = c1->second;
695  triangle[1] = c2->second;
696  triangle[2] = c3->second;
697  fTriangles.push_back(triangle);
698 
699  // remove the ear point from verticesToBeDone
700  //
701  verticesToBeDone.erase(c2);
702  c1 = verticesToBeDone.begin();
703  c2 = c1+1;
704  c3 = c1+2;
705  }
706  }
707  return true;
708 }
709 
710 //_____________________________________________________________________________
711 
712 G4bool G4ExtrudedSolid::MakeFacets()
713 {
714  // Define facets
715 
716  G4bool good;
717 
718  // Decomposition of polygonal sides in the facets
719  //
720  if ( fNv == 3 )
721  {
722  good = AddFacet( new G4TriangularFacet( GetVertex(0, 0), GetVertex(0, 1),
723  GetVertex(0, 2), ABSOLUTE) );
724  if ( ! good ) { return false; }
725 
726  good = AddFacet( new G4TriangularFacet( GetVertex(fNz-1, 2),
727  GetVertex(fNz-1, 1),
728  GetVertex(fNz-1, 0),
729  ABSOLUTE) );
730  if ( ! good ) { return false; }
731 
732  std::vector<G4int> triangle(3);
733  triangle[0] = 0;
734  triangle[1] = 1;
735  triangle[2] = 2;
736  fTriangles.push_back(triangle);
737  }
738 
739  else if ( fNv == 4 )
740  {
741  good = AddFacet( new G4QuadrangularFacet( GetVertex(0, 0),GetVertex(0, 1),
742  GetVertex(0, 2),GetVertex(0, 3),
743  ABSOLUTE) );
744  if ( ! good ) { return false; }
745 
746  good = AddFacet( new G4QuadrangularFacet( GetVertex(fNz-1, 3),
747  GetVertex(fNz-1, 2),
748  GetVertex(fNz-1, 1),
749  GetVertex(fNz-1, 0),
750  ABSOLUTE) );
751  if ( ! good ) { return false; }
752 
753  std::vector<G4int> triangle1(3);
754  triangle1[0] = 0;
755  triangle1[1] = 1;
756  triangle1[2] = 2;
757  fTriangles.push_back(triangle1);
758 
759  std::vector<G4int> triangle2(3);
760  triangle2[0] = 0;
761  triangle2[1] = 2;
762  triangle2[2] = 3;
763  fTriangles.push_back(triangle2);
764  }
765  else
766  {
767  good = AddGeneralPolygonFacets();
768  if ( ! good ) { return false; }
769  }
770 
771  // The quadrangular sides
772  //
773  for ( G4int iz = 0; iz < fNz-1; ++iz )
774  {
775  for ( G4int i = 0; i < fNv; ++i )
776  {
777  G4int j = (i+1) % fNv;
778  good = AddFacet( new G4QuadrangularFacet
779  ( GetVertex(iz, j), GetVertex(iz, i),
780  GetVertex(iz+1, i), GetVertex(iz+1, j), ABSOLUTE) );
781  if ( ! good ) { return false; }
782  }
783  }
784 
785  SetSolidClosed(true);
786 
787  return good;
788 }
789 
790 //_____________________________________________________________________________
791 
793 {
794  // Return entity type
795 
796  return fGeometryType;
797 }
798 
799 //_____________________________________________________________________________
800 
802 {
803  return new G4ExtrudedSolid(*this);
804 }
805 
806 //_____________________________________________________________________________
807 
809 {
810  // Override the base class function as it fails in case of concave polygon.
811  // Project the point in the original polygon scale and check if it is inside
812  // for each triangle.
813 
814  // Check first if outside extent
815  //
816  if ( p.x() < GetMinXExtent() - kCarToleranceHalf ||
817  p.x() > GetMaxXExtent() + kCarToleranceHalf ||
818  p.y() < GetMinYExtent() - kCarToleranceHalf ||
819  p.y() > GetMaxYExtent() + kCarToleranceHalf ||
820  p.z() < GetMinZExtent() - kCarToleranceHalf ||
821  p.z() > GetMaxZExtent() + kCarToleranceHalf )
822  {
823  // G4cout << "G4ExtrudedSolid::Outside extent: " << p << G4endl;
824  return kOutside;
825  }
826 
827  // Project point p(z) to the polygon scale p0
828  //
829  G4TwoVector pscaled = ProjectPoint(p);
830 
831  // Check if on surface of polygon
832  //
833  for ( G4int i=0; i<fNv; ++i )
834  {
835  G4int j = (i+1) % fNv;
836  if ( IsSameLineSegment(pscaled, fPolygon[i], fPolygon[j]) )
837  {
838  // G4cout << "G4ExtrudedSolid::Inside return Surface (on polygon) "
839  // << G4endl;
840 
841  return kSurface;
842  }
843  }
844 
845  // Now check if inside triangles
846  //
847  std::vector< std::vector<G4int> >::const_iterator it = fTriangles.begin();
848  G4bool inside = false;
849  do // Loop checking, 13.08.2015, G.Cosmo
850  {
851  if ( IsPointInside(fPolygon[(*it)[0]], fPolygon[(*it)[1]],
852  fPolygon[(*it)[2]], pscaled) ) { inside = true; }
853  ++it;
854  } while ( (inside == false) && (it != fTriangles.end()) );
855 
856  if ( inside )
857  {
858  // Check if on surface of z sides
859  //
860  if ( std::fabs( p.z() - fZSections[0].fZ ) < kCarToleranceHalf ||
861  std::fabs( p.z() - fZSections[fNz-1].fZ ) < kCarToleranceHalf )
862  {
863  // G4cout << "G4ExtrudedSolid::Inside return Surface (on z side)"
864  // << G4endl;
865 
866  return kSurface;
867  }
868 
869  // G4cout << "G4ExtrudedSolid::Inside return Inside" << G4endl;
870 
871  return kInside;
872  }
873 
874  // G4cout << "G4ExtrudedSolid::Inside return Outside " << G4endl;
875 
876  return kOutside;
877 }
878 
879 //_____________________________________________________________________________
880 
882  const G4ThreeVector &v,
883  const G4bool calcNorm,
884  G4bool *validNorm,
885  G4ThreeVector *n) const
886 {
887  // Override the base class function to redefine validNorm
888  // (the solid can be concave)
889 
890  G4double distOut =
891  G4TessellatedSolid::DistanceToOut(p, v, calcNorm, validNorm, n);
892  if (validNorm) { *validNorm = fIsConvex; }
893 
894  return distOut;
895 }
896 
897 
898 //_____________________________________________________________________________
899 
901 {
902  // Override the overloaded base class function
903 
905 }
906 
908 //
909 // Get bounding box
910 
912 {
913  G4double xmin0 = kInfinity, xmax0 = -kInfinity;
914  G4double ymin0 = kInfinity, ymax0 = -kInfinity;
915 
916  for (G4int i=0; i<GetNofVertices(); ++i)
917  {
918  G4double x = fPolygon[i].x();
919  if (x < xmin0) xmin0 = x;
920  if (x > xmax0) xmax0 = x;
921  G4double y = fPolygon[i].y();
922  if (y < ymin0) ymin0 = y;
923  if (y > ymax0) ymax0 = y;
924  }
925 
926  G4double xmin = kInfinity, xmax = -kInfinity;
927  G4double ymin = kInfinity, ymax = -kInfinity;
928 
929  G4int nsect = GetNofZSections();
930  for (G4int i=0; i<nsect; ++i)
931  {
932  ZSection zsect = GetZSection(i);
933  G4double dx = zsect.fOffset.x();
934  G4double dy = zsect.fOffset.y();
935  G4double scale = zsect.fScale;
936  xmin = std::min(xmin,xmin0*scale+dx);
937  xmax = std::max(xmax,xmax0*scale+dx);
938  ymin = std::min(ymin,ymin0*scale+dy);
939  ymax = std::max(ymax,ymax0*scale+dy);
940  }
941 
942  G4double zmin = GetZSection(0).fZ;
943  G4double zmax = GetZSection(nsect-1).fZ;
944 
945  pMin.set(xmin,ymin,zmin);
946  pMax.set(xmax,ymax,zmax);
947 
948  // Check correctness of the bounding box
949  //
950  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
951  {
952  std::ostringstream message;
953  message << "Bad bounding box (min >= max) for solid: "
954  << GetName() << " !"
955  << "\npMin = " << pMin
956  << "\npMax = " << pMax;
957  G4Exception("G4ExtrudedSolid::Extent()",
958  "GeomMgt0001", JustWarning, message);
959  DumpInfo();
960  }
961 }
962 
964 //
965 // Calculate extent under transform and specified limit
966 
967 G4bool
969  const G4VoxelLimits& pVoxelLimit,
970  const G4AffineTransform& pTransform,
971  G4double& pMin, G4double& pMax) const
972 {
973  G4ThreeVector bmin, bmax;
974  G4bool exist;
975 
976  // Check bounding box (bbox)
977  //
978  Extent(bmin,bmax);
979  G4BoundingEnvelope bbox(bmin,bmax);
980 #ifdef G4BBOX_EXTENT
981  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
982 #endif
983  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
984  {
985  return exist = (pMin < pMax) ? true : false;
986  }
987 
988  // To find the extent, the base polygon is subdivided in triangles.
989  // The extent is calculated as cumulative extent of the parts
990  // formed by extrusion of the triangles
991  //
992  G4TwoVectorList triangles;
993  G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
994  G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
995 
996  // triangulate the base polygon
997  if (!G4GeomTools::TriangulatePolygon(fPolygon,triangles))
998  {
999  std::ostringstream message;
1000  message << "Triangulation of the base polygon has failed for solid: "
1001  << GetName() << " !"
1002  << "\nExtent has been calculated using boundary box";
1003  G4Exception("G4ExtrudedSolid::CalculateExtent()",
1004  "GeomMgt1002",JustWarning,message);
1005  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
1006  }
1007 
1008  // allocate vector lists
1009  G4int nsect = GetNofZSections();
1010  std::vector<const G4ThreeVectorList *> polygons;
1011  polygons.resize(nsect);
1012  for (G4int k=0; k<nsect; ++k) { polygons[k] = new G4ThreeVectorList(3); }
1013 
1014  // main loop along triangles
1015  pMin = kInfinity;
1016  pMax = -kInfinity;
1017  G4int ntria = triangles.size()/3;
1018  for (G4int i=0; i<ntria; ++i)
1019  {
1020  G4int i3 = i*3;
1021  for (G4int k=0; k<nsect; ++k) // extrude triangle
1022  {
1023  ZSection zsect = GetZSection(k);
1024  G4double z = zsect.fZ;
1025  G4double dx = zsect.fOffset.x();
1026  G4double dy = zsect.fOffset.y();
1027  G4double scale = zsect.fScale;
1028 
1029  G4ThreeVectorList* ptr = const_cast<G4ThreeVectorList*>(polygons[k]);
1030  G4ThreeVectorList::iterator iter = ptr->begin();
1031  G4double x0 = triangles[i3+0].x()*scale+dx;
1032  G4double y0 = triangles[i3+0].y()*scale+dy;
1033  iter->set(x0,y0,z);
1034  iter++;
1035  G4double x1 = triangles[i3+1].x()*scale+dx;
1036  G4double y1 = triangles[i3+1].y()*scale+dy;
1037  iter->set(x1,y1,z);
1038  iter++;
1039  G4double x2 = triangles[i3+2].x()*scale+dx;
1040  G4double y2 = triangles[i3+2].y()*scale+dy;
1041  iter->set(x2,y2,z);
1042  }
1043 
1044  // set sub-envelope and adjust extent
1045  G4double emin,emax;
1046  G4BoundingEnvelope benv(polygons);
1047  if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
1048  if (emin < pMin) pMin = emin;
1049  if (emax > pMax) pMax = emax;
1050  if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
1051  }
1052  // free memory
1053  for (G4int k=0; k<nsect; ++k) { delete polygons[k]; polygons[k]=0;}
1054  return (pMin < pMax);
1055 }
1056 
1057 //_____________________________________________________________________________
1058 
1059 std::ostream& G4ExtrudedSolid::StreamInfo(std::ostream &os) const
1060 {
1061  G4int oldprc = os.precision(16);
1062  os << "-----------------------------------------------------------\n"
1063  << " *** Dump for solid - " << GetName() << " ***\n"
1064  << " ===================================================\n"
1065  << " Solid geometry type: " << fGeometryType << G4endl;
1066 
1067  if ( fIsConvex)
1068  { os << " Convex polygon; list of vertices:" << G4endl; }
1069  else
1070  { os << " Concave polygon; list of vertices:" << G4endl; }
1071 
1072  for ( G4int i=0; i<fNv; ++i )
1073  {
1074  os << std::setw(5) << "#" << i
1075  << " vx = " << fPolygon[i].x()/mm << " mm"
1076  << " vy = " << fPolygon[i].y()/mm << " mm" << G4endl;
1077  }
1078 
1079  os << " Sections:" << G4endl;
1080  for ( G4int iz=0; iz<fNz; ++iz )
1081  {
1082  os << " z = " << fZSections[iz].fZ/mm << " mm "
1083  << " x0= " << fZSections[iz].fOffset.x()/mm << " mm "
1084  << " y0= " << fZSections[iz].fOffset.y()/mm << " mm "
1085  << " scale= " << fZSections[iz].fScale << G4endl;
1086  }
1087 
1088 /*
1089  // Triangles (for debugging)
1090  os << G4endl;
1091  os << " Triangles:" << G4endl;
1092  os << " Triangle # vertex1 vertex2 vertex3" << G4endl;
1093 
1094  G4int counter = 0;
1095  std::vector< std::vector<G4int> >::const_iterator it;
1096  for ( it = fTriangles.begin(); it != fTriangles.end(); it++ ) {
1097  std::vector<G4int> triangle = *it;
1098  os << std::setw(10) << counter++
1099  << std::setw(10) << triangle[0] << std::setw(10) << triangle[1]
1100  << std::setw(10) << triangle[2]
1101  << G4endl;
1102  }
1103 */
1104  os.precision(oldprc);
1105 
1106  return os;
1107 }
1108 
1109 #endif
G4double G4ParticleHPJENDLHEData::G4double result
void set(double x, double y, double z)
G4String GetName() const
double y() const
G4GeometryType GetEntityType() const
void SetSolidClosed(const G4bool t)
void Extent(G4ThreeVector &pMin, G4ThreeVector &pMax) const
double x() const
static constexpr double mm
Definition: G4SIunits.hh:115
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
G4VSolid * Clone() const
double x() const
std::vector< G4TwoVector > G4TwoVectorList
Definition: G4GeomTools.hh:50
const char * p
Definition: xmltok.h:285
G4double GetMaxXExtent() const
static G4double angle[DIM]
G4ExtrudedSolid(const G4String &pName, const std::vector< G4TwoVector > &polygon, const std::vector< ZSection > &zsections)
static G4double PolygonArea(const G4TwoVectorList &polygon)
Definition: G4GeomTools.cc:82
G4double GetMinXExtent() const
virtual ~G4ExtrudedSolid()
EInside Inside(const G4ThreeVector &p) const
G4double GetMaxZExtent() const
int G4int
Definition: G4Types.hh:78
double z() const
G4TwoVector GetVertex(G4int index) const
void DumpInfo() const
G4double GetMinZExtent() const
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
virtual G4double DistanceToOut(const G4ThreeVector &p) const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
bool G4bool
Definition: G4Types.hh:79
G4int GetNofVertices() const
G4bool AddFacet(G4VFacet *aFacet)
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
Definition: G4GeomTools.cc:178
std::vector< G4ThreeVector > G4ThreeVectorList
G4int GetNofZSections() const
std::ostream & StreamInfo(std::ostream &os) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double emax
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
double y() const
ZSection GetZSection(G4int index) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4bool IsConvex(const G4TwoVectorList &polygon)
Definition: G4GeomTools.cc:150
#define G4endl
Definition: G4ios.hh:61
G4double kCarTolerance
Definition: G4VSolid.hh:307
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
G4double GetMaxExtent(const EAxis pAxis) const
G4double GetMaxYExtent() const
G4double GetAngularTolerance() const
G4ExtrudedSolid & operator=(const G4ExtrudedSolid &rhs)
G4double GetMinExtent(const EAxis pAxis) const
static G4GeometryTolerance * GetInstance()
G4double GetMinYExtent() const
G4TessellatedSolid & operator=(const G4TessellatedSolid &right)
static void RemoveRedundantVertices(G4TwoVectorList &polygon, std::vector< G4int > &iout, G4double tolerance=0)
Definition: G4GeomTools.cc:293