Geant4  10.02.p02
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 95805 2016-02-25 11:11:49Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class source file
32 //
33 // G4ExtrudedSolid.cc
34 //
35 // Author: Ivana Hrivnacova, IPN Orsay
36 // --------------------------------------------------------------------
37 
38 #include "G4ExtrudedSolid.hh"
39 
40 #if !defined(G4GEOM_USE_UEXTRUDEDSOLID)
41 
42 #include <set>
43 #include <algorithm>
44 #include <cmath>
45 #include <iomanip>
46 
47 #include "G4GeometryTolerance.hh"
48 #include "G4PhysicalConstants.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "G4VFacet.hh"
51 #include "G4TriangularFacet.hh"
52 #include "G4QuadrangularFacet.hh"
53 
54 //_____________________________________________________________________________
55 
57  std::vector<G4TwoVector> polygon,
58  std::vector<ZSection> zsections)
59  : G4TessellatedSolid(pName),
60  fNv(polygon.size()),
61  fNz(zsections.size()),
62  fPolygon(),
63  fZSections(),
64  fTriangles(),
65  fIsConvex(false),
66  fGeometryType("G4ExtrudedSolid")
67 
68 {
69  // General constructor
70 
71  // First check input parameters
72 
73  if ( fNv < 3 )
74  {
75  std::ostringstream message;
76  message << "Number of polygon vertices < 3 - " << pName;
77  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
78  FatalErrorInArgument, message);
79  }
80 
81  if ( fNz < 2 )
82  {
83  std::ostringstream message;
84  message << "Number of z-sides < 2 - " << pName;
85  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
86  FatalErrorInArgument, message);
87  }
88 
89  for ( G4int i=0; i<fNz-1; ++i )
90  {
91  if ( zsections[i].fZ > zsections[i+1].fZ )
92  {
93  std::ostringstream message;
94  message << "Z-sections have to be ordered by z value (z0 < z1 < z2...) - "
95  << pName;
96  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
97  FatalErrorInArgument, message);
98  }
99  if ( std::fabs( zsections[i+1].fZ - zsections[i].fZ ) < kCarToleranceHalf )
100  {
101  std::ostringstream message;
102  message << "Z-sections with the same z position are not supported - "
103  << pName;
104  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0001",
105  FatalException, message);
106  }
107  }
108 
109  // Check if polygon vertices are defined clockwise
110  // (the area is positive if polygon vertices are defined anti-clockwise)
111  //
112  G4double area = 0.;
113  for ( G4int i=0; i<fNv; ++i ) {
114  G4int j = i+1;
115  if ( j == fNv ) j = 0;
116  area += 0.5 * ( polygon[i].x()*polygon[j].y() - polygon[j].x()*polygon[i].y());
117  }
118 
119  // Copy polygon
120  //
121  if ( area < 0. ) {
122  // Polygon vertices are defined clockwise, we just copy the polygon
123  for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[i]); }
124  }
125  else {
126  // Polygon vertices are defined anti-clockwise, we revert them
127  //G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
128  // JustWarning,
129  // "Polygon vertices defined anti-clockwise, reverting polygon");
130  for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[fNv-i-1]); }
131  }
132 
133 
134  // Copy z-sections
135  //
136  for ( G4int i=0; i<fNz; ++i ) { fZSections.push_back(zsections[i]); }
137 
138 
139  G4bool result = MakeFacets();
140  if (!result)
141  {
142  std::ostringstream message;
143  message << "Making facets failed - " << pName;
144  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0003",
145  FatalException, message);
146  }
147  fIsConvex = IsConvex();
148 
149 
151 }
152 
153 //_____________________________________________________________________________
154 
156  std::vector<G4TwoVector> polygon,
157  G4double dz,
158  G4TwoVector off1, G4double scale1,
159  G4TwoVector off2, G4double scale2 )
160  : G4TessellatedSolid(pName),
161  fNv(polygon.size()),
162  fNz(2),
163  fPolygon(),
164  fZSections(),
165  fTriangles(),
166  fIsConvex(false),
167  fGeometryType("G4ExtrudedSolid")
168 
169 {
170  // Special constructor for solid with 2 z-sections
171 
172  // First check input parameters
173  //
174  if ( fNv < 3 )
175  {
176  std::ostringstream message;
177  message << "Number of polygon vertices < 3 - " << pName;
178  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
179  FatalErrorInArgument, message);
180  }
181 
182  // Check if polygon vertices are defined clockwise
183  // (the area is positive if polygon vertices are defined anti-clockwise)
184 
185  G4double area = 0.;
186  for ( G4int i=0; i<fNv; ++i )
187  {
188  G4int j = i+1;
189  if ( j == fNv ) { j = 0; }
190  area += 0.5 * ( polygon[i].x()*polygon[j].y()
191  - polygon[j].x()*polygon[i].y());
192  }
193 
194  // Copy polygon
195  //
196  if ( area < 0. )
197  {
198  // Polygon vertices are defined clockwise, we just copy the polygon
199  for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[i]); }
200  }
201  else
202  {
203  // Polygon vertices are defined anti-clockwise, we revert them
204  //G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
205  // JustWarning,
206  // "Polygon vertices defined anti-clockwise, reverting polygon");
207  for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[fNv-i-1]); }
208  }
209 
210  // Copy z-sections
211  //
212  fZSections.push_back(ZSection(-dz, off1, scale1));
213  fZSections.push_back(ZSection( dz, off2, scale2));
214 
215  G4bool result = MakeFacets();
216  if (!result)
217  {
218  std::ostringstream message;
219  message << "Making facets failed - " << pName;
220  G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0003",
221  FatalException, message);
222  }
223  fIsConvex = IsConvex();
224 
226 }
227 
228 //_____________________________________________________________________________
229 
231  : G4TessellatedSolid(a), fNv(0), fNz(0), fPolygon(), fZSections(),
232  fTriangles(), fIsConvex(false), fGeometryType("G4ExtrudedSolid")
233 {
234  // Fake default constructor - sets only member data and allocates memory
235  // for usage restricted to object persistency.
236 }
237 
238 //_____________________________________________________________________________
239 
241  : G4TessellatedSolid(rhs), fNv(rhs.fNv), fNz(rhs.fNz),
242  fPolygon(rhs.fPolygon), fZSections(rhs.fZSections),
243  fTriangles(rhs.fTriangles), fIsConvex(rhs.fIsConvex),
244  fGeometryType(rhs.fGeometryType), fKScales(rhs.fKScales),
245  fScale0s(rhs.fScale0s), fKOffsets(rhs.fKOffsets), fOffset0s(rhs.fOffset0s)
246 {
247 }
248 
249 
250 //_____________________________________________________________________________
251 
253 {
254  // Check assignment to self
255  //
256  if (this == &rhs) { return *this; }
257 
258  // Copy base class data
259  //
261 
262  // Copy data
263  //
264  fNv = rhs.fNv; fNz = rhs.fNz;
265  fPolygon = rhs.fPolygon; fZSections = rhs.fZSections;
268  fScale0s = rhs.fScale0s; fKOffsets = rhs.fKOffsets;
269  fOffset0s = rhs.fOffset0s;
270 
271  return *this;
272 }
273 
274 //_____________________________________________________________________________
275 
277 {
278  // Destructor
279 }
280 
281 //_____________________________________________________________________________
282 
284 {
285  // Compute parameters for point projections p(z)
286  // to the polygon scale & offset:
287  // scale(z) = k*z + scale0
288  // offset(z) = l*z + offset0
289  // p(z) = scale(z)*p0 + offset(z)
290  // p0 = (p(z) - offset(z))/scale(z);
291  //
292 
293  for ( G4int iz=0; iz<fNz-1; ++iz)
294  {
295  G4double z1 = fZSections[iz].fZ;
296  G4double z2 = fZSections[iz+1].fZ;
297  G4double scale1 = fZSections[iz].fScale;
298  G4double scale2 = fZSections[iz+1].fScale;
299  G4TwoVector off1 = fZSections[iz].fOffset;
300  G4TwoVector off2 = fZSections[iz+1].fOffset;
301 
302  G4double kscale = (scale2 - scale1)/(z2 - z1);
303  G4double scale0 = scale2 - kscale*(z2 - z1)/2.0;
304  G4TwoVector koff = (off2 - off1)/(z2 - z1);
305  G4TwoVector off0 = off2 - koff*(z2 - z1)/2.0;
306 
307  fKScales.push_back(kscale);
308  fScale0s.push_back(scale0);
309  fKOffsets.push_back(koff);
310  fOffset0s.push_back(off0);
311  }
312 }
313 
314 
315 //_____________________________________________________________________________
316 
318 {
319  // Shift and scale vertices
320 
321  return G4ThreeVector( fPolygon[ind].x() * fZSections[iz].fScale
322  + fZSections[iz].fOffset.x(),
323  fPolygon[ind].y() * fZSections[iz].fScale
324  + fZSections[iz].fOffset.y(), fZSections[iz].fZ);
325 }
326 
327 //_____________________________________________________________________________
328 
329 
331 {
332  // Project point in the polygon scale
333  // scale(z) = k*z + scale0
334  // offset(z) = l*z + offset0
335  // p(z) = scale(z)*p0 + offset(z)
336  // p0 = (p(z) - offset(z))/scale(z);
337 
338  // Select projection (z-segment of the solid) according to p.z()
339  //
340  G4int iz = 0;
341  while ( point.z() > fZSections[iz+1].fZ && iz < fNz-2 ) { ++iz; }
342  // Loop checking, 13.08.2015, G.Cosmo
343 
344  G4double z0 = ( fZSections[iz+1].fZ + fZSections[iz].fZ )/2.0;
345  G4TwoVector p2(point.x(), point.y());
346  G4double pscale = fKScales[iz]*(point.z()-z0) + fScale0s[iz];
347  G4TwoVector poffset = fKOffsets[iz]*(point.z()-z0) + fOffset0s[iz];
348 
349  // G4cout << point << " projected to "
350  // << iz << "-th z-segment polygon as "
351  // << (p2 - poffset)/pscale << G4endl;
352 
353  // pscale is always >0 as it is an interpolation between two
354  // positive scale values
355  //
356  return (p2 - poffset)/pscale;
357 }
358 
359 //_____________________________________________________________________________
360 
362  G4TwoVector l1, G4TwoVector l2) const
363 {
364  // Return true if p is on the line through l1, l2
365 
366  if ( l1.x() == l2.x() )
367  {
368  return std::fabs(p.x() - l1.x()) < kCarToleranceHalf;
369  }
370  G4double slope= ((l2.y() - l1.y())/(l2.x() - l1.x()));
371  G4double predy= l1.y() + slope *(p.x() - l1.x());
372  G4double dy= p.y() - predy;
373 
374  // Calculate perpendicular distance
375  //
376  // G4double perpD= std::fabs(dy) / std::sqrt( 1 + slope * slope );
377  // G4bool simpleComp= (perpD<kCarToleranceHalf);
378 
379  // Check perpendicular distance vs tolerance 'directly'
380  //
381  G4bool squareComp= (dy*dy < (1+slope*slope)
383 
384  // return simpleComp;
385  return squareComp;
386 }
387 
388 //_____________________________________________________________________________
389 
391  G4TwoVector l1, G4TwoVector l2) const
392 {
393  // Return true if p is on the line through l1, l2 and lies between
394  // l1 and l2
395 
396  if ( p.x() < std::min(l1.x(), l2.x()) - kCarToleranceHalf ||
397  p.x() > std::max(l1.x(), l2.x()) + kCarToleranceHalf ||
398  p.y() < std::min(l1.y(), l2.y()) - kCarToleranceHalf ||
399  p.y() > std::max(l1.y(), l2.y()) + kCarToleranceHalf )
400  {
401  return false;
402  }
403 
404  return IsSameLine(p, l1, l2);
405 }
406 
407 //_____________________________________________________________________________
408 
410  G4TwoVector l1, G4TwoVector l2) const
411 {
412  // Return true if p1 and p2 are on the same side of the line through l1, l2
413 
414  return ( (p1.x() - l1.x()) * (l2.y() - l1.y())
415  - (l2.x() - l1.x()) * (p1.y() - l1.y()) )
416  * ( (p2.x() - l1.x()) * (l2.y() - l1.y())
417  - (l2.x() - l1.x()) * (p2.y() - l1.y()) ) > 0;
418 }
419 
420 //_____________________________________________________________________________
421 
423  G4TwoVector c, G4TwoVector p) const
424 {
425  // Return true if p is inside of triangle abc or on its edges,
426  // else returns false
427 
428  // Check extent first
429  //
430  if ( ( p.x() < a.x() && p.x() < b.x() && p.x() < c.x() ) ||
431  ( p.x() > a.x() && p.x() > b.x() && p.x() > c.x() ) ||
432  ( p.y() < a.y() && p.y() < b.y() && p.y() < c.y() ) ||
433  ( p.y() > a.y() && p.y() > b.y() && p.y() > c.y() ) ) return false;
434 
435  G4bool inside
436  = IsSameSide(p, a, b, c)
437  && IsSameSide(p, b, a, c)
438  && IsSameSide(p, c, a, b);
439 
440  G4bool onEdge
441  = IsSameLineSegment(p, a, b)
442  || IsSameLineSegment(p, b, c)
443  || IsSameLineSegment(p, c, a);
444 
445  return inside || onEdge;
446 }
447 
448 //_____________________________________________________________________________
449 
450 G4double
452 {
453  // Return the angle of the vertex in po
454 
455  G4TwoVector t1 = pa - po;
456  G4TwoVector t2 = pb - po;
457 
458  G4double result = (std::atan2(t1.y(), t1.x()) - std::atan2(t2.y(), t2.x()));
459 
460  if ( result < 0 ) result += 2*pi;
461 
462  return result;
463 }
464 
465 //_____________________________________________________________________________
466 
467 G4VFacet*
469 {
470  // Create a triangular facet from the polygon points given by indices
471  // forming the down side ( the normal goes in -z)
472 
473  std::vector<G4ThreeVector> vertices;
474  vertices.push_back(GetVertex(0, ind1));
475  vertices.push_back(GetVertex(0, ind2));
476  vertices.push_back(GetVertex(0, ind3));
477 
478  // first vertex most left
479  //
480  G4ThreeVector cross
481  = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
482 
483  if ( cross.z() > 0.0 )
484  {
485  // vertices ardered clock wise has to be reordered
486 
487  // G4cout << "G4ExtrudedSolid::MakeDownFacet: reordering vertices "
488  // << ind1 << ", " << ind2 << ", " << ind3 << G4endl;
489 
490  G4ThreeVector tmp = vertices[1];
491  vertices[1] = vertices[2];
492  vertices[2] = tmp;
493  }
494 
495  return new G4TriangularFacet(vertices[0], vertices[1],
496  vertices[2], ABSOLUTE);
497 }
498 
499 //_____________________________________________________________________________
500 
501 G4VFacet*
503 {
504  // Creates a triangular facet from the polygon points given by indices
505  // forming the upper side ( z>0 )
506 
507  std::vector<G4ThreeVector> vertices;
508  vertices.push_back(GetVertex(fNz-1, ind1));
509  vertices.push_back(GetVertex(fNz-1, ind2));
510  vertices.push_back(GetVertex(fNz-1, ind3));
511 
512  // first vertex most left
513  //
514  G4ThreeVector cross
515  = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
516 
517  if ( cross.z() < 0.0 )
518  {
519  // vertices ordered clock wise has to be reordered
520 
521  // G4cout << "G4ExtrudedSolid::MakeUpFacet: reordering vertices "
522  // << ind1 << ", " << ind2 << ", " << ind3 << G4endl;
523 
524  G4ThreeVector tmp = vertices[1];
525  vertices[1] = vertices[2];
526  vertices[2] = tmp;
527  }
528 
529  return new G4TriangularFacet(vertices[0], vertices[1],
530  vertices[2], ABSOLUTE);
531 }
532 
533 //_____________________________________________________________________________
534 
536 {
537  // Decompose polygonal sides in triangular facets
538 
539  typedef std::pair < G4TwoVector, G4int > Vertex;
540 
541  static const G4double kAngTolerance =
543 
544  // Fill one more vector
545  //
546  std::vector< Vertex > verticesToBeDone;
547  for ( G4int i=0; i<fNv; ++i )
548  {
549  verticesToBeDone.push_back(Vertex(fPolygon[i], i));
550  }
551  std::vector< Vertex > ears;
552 
553  std::vector< Vertex >::iterator c1 = verticesToBeDone.begin();
554  std::vector< Vertex >::iterator c2 = c1+1;
555  std::vector< Vertex >::iterator c3 = c1+2;
556  while ( verticesToBeDone.size()>2 ) // Loop checking, 13.08.2015, G.Cosmo
557  {
558 
559  // G4cout << "Looking at triangle : "
560  // << c1->second << " " << c2->second
561  // << " " << c3->second << G4endl;
562  //G4cout << "Looking at triangle : "
563  // << c1->first << " " << c2->first
564  // << " " << c3->first << G4endl;
565 
566  // skip concave vertices
567  //
568  G4double angle = GetAngle(c2->first, c3->first, c1->first);
569 
570  //G4cout << "angle " << angle << G4endl;
571 
572  G4int counter = 0;
573  while ( angle >= (pi-kAngTolerance) ) // Loop checking, 13.08.2015, G.Cosmo
574  {
575  // G4cout << "Skipping concave vertex " << c2->second << G4endl;
576 
577  // try next three consecutive vertices
578  //
579  c1 = c2;
580  c2 = c3;
581  ++c3;
582  if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
583 
584  //G4cout << "Looking at triangle : "
585  // << c1->first << " " << c2->first
586  // << " " << c3->first << G4endl;
587 
588  angle = GetAngle(c2->first, c3->first, c1->first);
589  //G4cout << "angle " << angle << G4endl;
590 
591  counter++;
592 
593  if ( counter > fNv) {
594  G4Exception("G4ExtrudedSolid::AddGeneralPolygonFacets",
595  "GeomSolids0003", FatalException,
596  "Triangularisation has failed.");
597  break;
598  }
599  }
600 
601  G4bool good = true;
602  std::vector< Vertex >::iterator it;
603  for ( it=verticesToBeDone.begin(); it != verticesToBeDone.end(); ++it )
604  {
605  // skip vertices of tested triangle
606  //
607  if ( it == c1 || it == c2 || it == c3 ) { continue; }
608 
609  if ( IsPointInside(c1->first, c2->first, c3->first, it->first) )
610  {
611  // G4cout << "Point " << it->second << " is inside" << G4endl;
612  good = false;
613 
614  // try next three consecutive vertices
615  //
616  c1 = c2;
617  c2 = c3;
618  ++c3;
619  if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
620  break;
621  }
622  // else
623  // { G4cout << "Point " << it->second << " is outside" << G4endl; }
624  }
625  if ( good )
626  {
627  // all points are outside triangle, we can make a facet
628 
629  // G4cout << "Found triangle : "
630  // << c1->second << " " << c2->second
631  // << " " << c3->second << G4endl;
632 
633  G4bool result;
634  result = AddFacet( MakeDownFacet(c1->second, c2->second, c3->second) );
635  if ( ! result ) { return false; }
636 
637  result = AddFacet( MakeUpFacet(c1->second, c2->second, c3->second) );
638  if ( ! result ) { return false; }
639 
640  std::vector<G4int> triangle(3);
641  triangle[0] = c1->second;
642  triangle[1] = c2->second;
643  triangle[2] = c3->second;
644  fTriangles.push_back(triangle);
645 
646  // remove the ear point from verticesToBeDone
647  //
648  verticesToBeDone.erase(c2);
649  c1 = verticesToBeDone.begin();
650  c2 = c1+1;
651  c3 = c1+2;
652  }
653  }
654  return true;
655 }
656 
657 //_____________________________________________________________________________
658 
660 {
661  // Define facets
662 
663  G4bool good;
664 
665  // Decomposition of polygonal sides in the facets
666  //
667  if ( fNv == 3 )
668  {
669  good = AddFacet( new G4TriangularFacet( GetVertex(0, 0), GetVertex(0, 1),
670  GetVertex(0, 2), ABSOLUTE) );
671  if ( ! good ) { return false; }
672 
673  good = AddFacet( new G4TriangularFacet( GetVertex(fNz-1, 2), GetVertex(fNz-1, 1),
674  GetVertex(fNz-1, 0), ABSOLUTE) );
675  if ( ! good ) { return false; }
676 
677  std::vector<G4int> triangle(3);
678  triangle[0] = 0;
679  triangle[1] = 1;
680  triangle[2] = 2;
681  fTriangles.push_back(triangle);
682  }
683 
684  else if ( fNv == 4 )
685  {
686  good = AddFacet( new G4QuadrangularFacet( GetVertex(0, 0),GetVertex(0, 1),
687  GetVertex(0, 2),GetVertex(0, 3),
688  ABSOLUTE) );
689  if ( ! good ) { return false; }
690 
691  good = AddFacet( new G4QuadrangularFacet( GetVertex(fNz-1, 3), GetVertex(fNz-1, 2),
692  GetVertex(fNz-1, 1), GetVertex(fNz-1, 0),
693  ABSOLUTE) );
694  if ( ! good ) { return false; }
695 
696  std::vector<G4int> triangle1(3);
697  triangle1[0] = 0;
698  triangle1[1] = 1;
699  triangle1[2] = 2;
700  fTriangles.push_back(triangle1);
701 
702  std::vector<G4int> triangle2(3);
703  triangle2[0] = 0;
704  triangle2[1] = 2;
705  triangle2[2] = 3;
706  fTriangles.push_back(triangle2);
707  }
708  else
709  {
710  good = AddGeneralPolygonFacets();
711  if ( ! good ) { return false; }
712  }
713 
714  // The quadrangular sides
715  //
716  for ( G4int iz = 0; iz < fNz-1; ++iz )
717  {
718  for ( G4int i = 0; i < fNv; ++i )
719  {
720  G4int j = (i+1) % fNv;
721  good = AddFacet( new G4QuadrangularFacet
722  ( GetVertex(iz, j), GetVertex(iz, i),
723  GetVertex(iz+1, i), GetVertex(iz+1, j), ABSOLUTE) );
724  if ( ! good ) { return false; }
725  }
726  }
727 
728  SetSolidClosed(true);
729 
730  return good;
731 }
732 
733 //_____________________________________________________________________________
734 
736 {
737  // Get polygon convexity (polygon is convex if all vertex angles are < pi )
738 
739  for ( G4int i=0; i< fNv; ++i )
740  {
741  G4int j = ( i + 1 ) % fNv;
742  G4int k = ( i + 2 ) % fNv;
743  G4TwoVector v1 = fPolygon[i]-fPolygon[j];
744  G4TwoVector v2 = fPolygon[k]-fPolygon[j];
745  G4double dphi = v2.phi() - v1.phi();
746  if ( dphi < 0. ) { dphi += 2.*pi; }
747 
748  if ( dphi >= pi ) { return false; }
749  }
750 
751  return true;
752 }
753 
754 //_____________________________________________________________________________
755 
757 {
758  // Return entity type
759 
760  return fGeometryType;
761 }
762 
763 //_____________________________________________________________________________
764 
766 {
767  return new G4ExtrudedSolid(*this);
768 }
769 
770 //_____________________________________________________________________________
771 
773 {
774  // Override the base class function as it fails in case of concave polygon.
775  // Project the point in the original polygon scale and check if it is inside
776  // for each triangle.
777 
778  // Check first if outside extent
779  //
780  if ( p.x() < GetMinXExtent() - kCarToleranceHalf ||
781  p.x() > GetMaxXExtent() + kCarToleranceHalf ||
782  p.y() < GetMinYExtent() - kCarToleranceHalf ||
783  p.y() > GetMaxYExtent() + kCarToleranceHalf ||
784  p.z() < GetMinZExtent() - kCarToleranceHalf ||
785  p.z() > GetMaxZExtent() + kCarToleranceHalf )
786  {
787  // G4cout << "G4ExtrudedSolid::Outside extent: " << p << G4endl;
788  return kOutside;
789  }
790 
791  // Project point p(z) to the polygon scale p0
792  //
793  G4TwoVector pscaled = ProjectPoint(p);
794 
795  // Check if on surface of polygon
796  //
797  for ( G4int i=0; i<fNv; ++i )
798  {
799  G4int j = (i+1) % fNv;
800  if ( IsSameLineSegment(pscaled, fPolygon[i], fPolygon[j]) )
801  {
802  // G4cout << "G4ExtrudedSolid::Inside return Surface (on polygon) "
803  // << G4endl;
804 
805  return kSurface;
806  }
807  }
808 
809  // Now check if inside triangles
810  //
811  std::vector< std::vector<G4int> >::const_iterator it = fTriangles.begin();
812  G4bool inside = false;
813  do // Loop checking, 13.08.2015, G.Cosmo
814  {
815  if ( IsPointInside(fPolygon[(*it)[0]], fPolygon[(*it)[1]],
816  fPolygon[(*it)[2]], pscaled) ) { inside = true; }
817  ++it;
818  } while ( (inside == false) && (it != fTriangles.end()) );
819 
820  if ( inside )
821  {
822  // Check if on surface of z sides
823  //
824  if ( std::fabs( p.z() - fZSections[0].fZ ) < kCarToleranceHalf ||
825  std::fabs( p.z() - fZSections[fNz-1].fZ ) < kCarToleranceHalf )
826  {
827  // G4cout << "G4ExtrudedSolid::Inside return Surface (on z side)"
828  // << G4endl;
829 
830  return kSurface;
831  }
832 
833  // G4cout << "G4ExtrudedSolid::Inside return Inside" << G4endl;
834 
835  return kInside;
836  }
837 
838  // G4cout << "G4ExtrudedSolid::Inside return Outside " << G4endl;
839 
840  return kOutside;
841 }
842 
843 //_____________________________________________________________________________
844 
846  const G4ThreeVector &v,
847  const G4bool calcNorm,
848  G4bool *validNorm,
849  G4ThreeVector *n) const
850 {
851  // Override the base class function to redefine validNorm
852  // (the solid can be concave)
853 
854  G4double distOut =
855  G4TessellatedSolid::DistanceToOut(p, v, calcNorm, validNorm, n);
856  if (validNorm) { *validNorm = fIsConvex; }
857 
858  return distOut;
859 }
860 
861 
862 //_____________________________________________________________________________
863 
865 {
866  // Override the overloaded base class function
867 
869 }
870 
871 //_____________________________________________________________________________
872 
873 std::ostream& G4ExtrudedSolid::StreamInfo(std::ostream &os) const
874 {
875  G4int oldprc = os.precision(16);
876  os << "-----------------------------------------------------------\n"
877  << " *** Dump for solid - " << GetName() << " ***\n"
878  << " ===================================================\n"
879  << " Solid geometry type: " << fGeometryType << G4endl;
880 
881  if ( fIsConvex)
882  { os << " Convex polygon; list of vertices:" << G4endl; }
883  else
884  { os << " Concave polygon; list of vertices:" << G4endl; }
885 
886  for ( G4int i=0; i<fNv; ++i )
887  {
888  os << std::setw(5) << "#" << i
889  << " vx = " << fPolygon[i].x()/mm << " mm"
890  << " vy = " << fPolygon[i].y()/mm << " mm" << G4endl;
891  }
892 
893  os << " Sections:" << G4endl;
894  for ( G4int iz=0; iz<fNz; ++iz )
895  {
896  os << " z = " << fZSections[iz].fZ/mm << " mm "
897  << " x0= " << fZSections[iz].fOffset.x()/mm << " mm "
898  << " y0= " << fZSections[iz].fOffset.y()/mm << " mm "
899  << " scale= " << fZSections[iz].fScale << G4endl;
900  }
901 
902 /*
903  // Triangles (for debugging)
904  os << G4endl;
905  os << " Triangles:" << G4endl;
906  os << " Triangle # vertex1 vertex2 vertex3" << G4endl;
907 
908  G4int counter = 0;
909  std::vector< std::vector<G4int> >::const_iterator it;
910  for ( it = fTriangles.begin(); it != fTriangles.end(); it++ ) {
911  std::vector<G4int> triangle = *it;
912  os << std::setw(10) << counter++
913  << std::setw(10) << triangle[0] << std::setw(10) << triangle[1] << std::setw(10) << triangle[2]
914  << G4endl;
915  }
916 */
917  os.precision(oldprc);
918 
919  return os;
920 }
921 
922 #endif
G4String GetName() const
static c2_factory< G4double > c2
G4bool IsPointInside(G4TwoVector a, G4TwoVector b, G4TwoVector c, G4TwoVector p) const
G4GeometryType GetEntityType() const
void SetSolidClosed(const G4bool t)
G4bool IsConvex() const
CLHEP::Hep3Vector G4ThreeVector
G4VSolid * Clone() const
std::vector< G4TwoVector > fOffset0s
std::vector< G4double > fScale0s
G4ExtrudedSolid(const G4String &pName, std::vector< G4TwoVector > polygon, std::vector< ZSection > zsections)
G4double GetMaxXExtent() const
std::vector< ZSection > fZSections
static G4double angle[DIM]
G4double GetMinXExtent() const
G4double a
Definition: TRTMaterials.hh:39
virtual ~G4ExtrudedSolid()
EInside Inside(const G4ThreeVector &p) const
G4double GetMaxZExtent() const
int G4int
Definition: G4Types.hh:78
G4TwoVector ProjectPoint(const G4ThreeVector &point) const
G4bool AddGeneralPolygonFacets()
G4TwoVector GetVertex(G4int index) const
G4bool IsSameLine(G4TwoVector p, G4TwoVector l1, G4TwoVector l2) const
G4double GetAngle(G4TwoVector p0, G4TwoVector pa, G4TwoVector pb) const
G4double GetMinZExtent() const
void ComputeProjectionParameters()
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
bool G4bool
Definition: G4Types.hh:79
G4double iz
Definition: TRTMaterials.hh:39
static const G4double c3
G4bool AddFacet(G4VFacet *aFacet)
std::vector< G4double > fKScales
const G4int n
static const G4double c1
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 double pi
Definition: G4SIunits.hh:74
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
G4bool IsSameLineSegment(G4TwoVector p, G4TwoVector l1, G4TwoVector l2) const
const G4double x[NPOINTSGL]
std::vector< G4TwoVector > fKOffsets
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4bool IsSameSide(G4TwoVector p1, G4TwoVector p2, G4TwoVector l1, G4TwoVector l2) const
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4VFacet * MakeUpFacet(G4int ind1, G4int ind2, G4int ind3) const
G4GeometryType fGeometryType
static const double mm
Definition: G4SIunits.hh:114
G4double GetMaxYExtent() const
G4double GetAngularTolerance() const
G4ExtrudedSolid & operator=(const G4ExtrudedSolid &rhs)
static G4GeometryTolerance * GetInstance()
G4VFacet * MakeDownFacet(G4int ind1, G4int ind2, G4int ind3) const
G4double GetMinYExtent() const
std::vector< std::vector< G4int > > fTriangles
G4TessellatedSolid & operator=(const G4TessellatedSolid &right)
std::vector< G4TwoVector > fPolygon