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