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