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