Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GenericTrap.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// G4GenericTrap implementation
27//
28// Authors:
29// Tatiana Nikitina, CERN; Ivana Hrivnacova, IPN Orsay
30// Adapted from Root Arb8 implementation by Andrei Gheata, CERN
31// --------------------------------------------------------------------
32
33#include "G4GenericTrap.hh"
34
35#if !defined(G4GEOM_USE_UGENERICTRAP)
36
37#include <iomanip>
38
40#include "G4SystemOfUnits.hh"
41#include "G4TriangularFacet.hh"
43#include "G4VoxelLimits.hh"
44#include "G4AffineTransform.hh"
45#include "G4BoundingEnvelope.hh"
46#include "Randomize.hh"
47
48#include "G4VGraphicsScene.hh"
49#include "G4Polyhedron.hh"
50#include "G4VisExtent.hh"
51
52#include "G4AutoLock.hh"
53
54namespace
55{
56 G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
57}
58
59const G4int G4GenericTrap::fgkNofVertices = 8;
60const G4double G4GenericTrap::fgkTolerance = 1E-3;
61
62// --------------------------------------------------------------------
63
65 const std::vector<G4TwoVector>& vertices )
66 : G4VSolid(name), fDz(halfZ), fVertices(),
67 fMinBBoxVector(G4ThreeVector(0,0,0)), fMaxBBoxVector(G4ThreeVector(0,0,0))
68{
69 // General constructor
70 const G4double min_length=5*1.e-6;
71 G4double length = 0.;
72 G4int k=0;
73 G4String errorDescription = "InvalidSetup in \" ";
74 errorDescription += name;
75 errorDescription += "\"";
76
77 halfCarTolerance = kCarTolerance*0.5;
78
79 // Check vertices size
80
81 if ( G4int(vertices.size()) != fgkNofVertices )
82 {
83 G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids0002",
84 FatalErrorInArgument, "Number of vertices != 8");
85 }
86
87 // Check dZ
88 //
89 if (halfZ < kCarTolerance)
90 {
91 G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids0002",
92 FatalErrorInArgument, "dZ is too small or negative");
93 }
94
95 // Check Ordering and Copy vertices
96 //
97 if(CheckOrder(vertices))
98 {
99 for (G4int i=0; i<fgkNofVertices; ++i) {fVertices.push_back(vertices[i]);}
100 }
101 else
102 {
103 for (auto i=0; i <4; ++i) {fVertices.push_back(vertices[3-i]);}
104 for (auto i=0; i <4; ++i) {fVertices.push_back(vertices[7-i]);}
105 }
106
107 // Check length of segments and Adjust
108 //
109 for (auto j=0; j<2; ++j)
110 {
111 for (auto i=1; i<4; ++i)
112 {
113 k = j*4+i;
114 length = (fVertices[k]-fVertices[k-1]).mag();
115 if ( ( length < min_length) && ( length > kCarTolerance ) )
116 {
117 std::ostringstream message;
118 message << "Length segment is too small." << G4endl
119 << "Distance between " << fVertices[k-1] << " and "
120 << fVertices[k] << " is only " << length << " mm !";
121 G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids1001",
122 JustWarning, message, "Vertices will be collapsed.");
123 fVertices[k]=fVertices[k-1];
124 }
125 }
126 }
127
128 // Compute Twist
129 //
130 for(G4double & i : fTwist) { i=0.; }
131 fIsTwisted = ComputeIsTwisted();
132
133 // Compute Bounding Box
134 //
135 ComputeBBox();
136
137 // If not twisted - create tessellated solid
138 // (an alternative implementation for testing)
139 //
140#ifdef G4TESS_TEST
141 if ( !fIsTwisted ) { fTessellatedSolid = CreateTessellatedSolid(); }
142#endif
143}
144
145// --------------------------------------------------------------------
146
148 : G4VSolid(a), halfCarTolerance(0.), fDz(0.), fVertices(),
149 fMinBBoxVector(G4ThreeVector(0,0,0)), fMaxBBoxVector(G4ThreeVector(0,0,0))
150{
151 // Fake default constructor - sets only member data and allocates memory
152 // for usage restricted to object persistency.
153}
154
155// --------------------------------------------------------------------
156
158{
159 delete fTessellatedSolid;
160}
161
162// --------------------------------------------------------------------
163
165 : G4VSolid(rhs),
166 halfCarTolerance(rhs.halfCarTolerance),
167 fDz(rhs.fDz), fVertices(rhs.fVertices),
168 fIsTwisted(rhs.fIsTwisted),
169 fMinBBoxVector(rhs.fMinBBoxVector), fMaxBBoxVector(rhs.fMaxBBoxVector),
170 fVisSubdivisions(rhs.fVisSubdivisions),
171 fSurfaceArea(rhs.fSurfaceArea), fCubicVolume(rhs.fCubicVolume)
172{
173 for (auto i=0; i<4; ++i) { fTwist[i] = rhs.fTwist[i]; }
174#ifdef G4TESS_TEST
175 if (rhs.fTessellatedSolid && !fIsTwisted )
176 { fTessellatedSolid = CreateTessellatedSolid(); }
177#endif
178}
179
180// --------------------------------------------------------------------
181
183{
184 // Check assignment to self
185 //
186 if (this == &rhs) { return *this; }
187
188 // Copy base class data
189 //
191
192 // Copy data
193 //
194 halfCarTolerance = rhs.halfCarTolerance;
195 fDz = rhs.fDz; fVertices = rhs.fVertices;
196 fIsTwisted = rhs.fIsTwisted; fTessellatedSolid = nullptr;
197 fMinBBoxVector = rhs.fMinBBoxVector; fMaxBBoxVector = rhs.fMaxBBoxVector;
198 fVisSubdivisions = rhs.fVisSubdivisions;
199 fSurfaceArea = rhs.fSurfaceArea; fCubicVolume = rhs.fCubicVolume;
200
201 for (auto i=0; i<4; ++i) { fTwist[i] = rhs.fTwist[i]; }
202#ifdef G4TESS_TEST
203 if (rhs.fTessellatedSolid && !fIsTwisted )
204 { delete fTessellatedSolid; fTessellatedSolid = CreateTessellatedSolid(); }
205#endif
206 fRebuildPolyhedron = false;
207 delete fpPolyhedron; fpPolyhedron = nullptr;
208
209 return *this;
210}
211
212// --------------------------------------------------------------------
213
215G4GenericTrap::InsidePolygone(const G4ThreeVector& p,
216 const std::vector<G4TwoVector>& poly) const
217{
218 EInside in = kInside;
219 G4double cross, len2;
220 G4int count=0;
221
222 for (G4int i=0; i<4; ++i)
223 {
224 G4int j = (i+1) % 4;
225
226 cross = (p.x()-poly[i].x())*(poly[j].y()-poly[i].y())-
227 (p.y()-poly[i].y())*(poly[j].x()-poly[i].x());
228
229 len2=(poly[i]-poly[j]).mag2();
230 if (len2 > kCarTolerance)
231 {
232 if(cross*cross<=len2*halfCarTolerance*halfCarTolerance) // Surface check
233 {
234 G4double test;
235
236 // Check if p lies between the two extremes of the segment
237 //
238 G4int iMax;
239 G4int iMin;
240
241 if (poly[j].x() > poly[i].x())
242 {
243 iMax = j;
244 iMin = i;
245 }
246 else {
247 iMax = i;
248 iMin = j;
249 }
250 if ( p.x() > poly[iMax].x()+halfCarTolerance
251 || p.x() < poly[iMin].x()-halfCarTolerance )
252 {
253 return kOutside;
254 }
255
256 if (poly[j].y() > poly[i].y())
257 {
258 iMax = j;
259 iMin = i;
260 }
261 else
262 {
263 iMax = i;
264 iMin = j;
265 }
266 if ( p.y() > poly[iMax].y()+halfCarTolerance
267 || p.y() < poly[iMin].y()-halfCarTolerance )
268 {
269 return kOutside;
270 }
271
272 if ( poly[iMax].x() != poly[iMin].x() )
273 {
274 test = (p.x()-poly[iMin].x())/(poly[iMax].x()-poly[iMin].x())
275 * (poly[iMax].y()-poly[iMin].y())+poly[iMin].y();
276 }
277 else
278 {
279 test = p.y();
280 }
281
282 // Check if point is Inside Segment
283 //
284 if( (test>=(poly[iMin].y()-halfCarTolerance))
285 && (test<=(poly[iMax].y()+halfCarTolerance)) )
286 {
287 return kSurface;
288 }
289 else
290 {
291 return kOutside;
292 }
293 }
294 else if (cross<0.) { return kOutside; }
295 }
296 else
297 {
298 ++count;
299 }
300 }
301
302 // All collapsed vertices, Tet like
303 //
304 if(count==4)
305 {
306 if ( (std::fabs(p.x()-poly[0].x())
307 +std::fabs(p.y()-poly[0].y())) > halfCarTolerance )
308 {
309 in = kOutside;
310 }
311 }
312 return in;
313}
314
315// --------------------------------------------------------------------
316
318{
319 // Test if point is inside this shape
320
321#ifdef G4TESS_TEST
322 if ( fTessellatedSolid )
323 {
324 return fTessellatedSolid->Inside(p);
325 }
326#endif
327
328 EInside innew=kOutside;
329 std::vector<G4TwoVector> xy;
330
331 if (std::fabs(p.z()) <= fDz+halfCarTolerance) // First check Z range
332 {
333 // Compute intersection between Z plane containing point and the shape
334 //
335 G4double cf = 0.5*(fDz-p.z())/fDz;
336 for (auto i=0; i<4; ++i)
337 {
338 xy.push_back(fVertices[i+4]+cf*( fVertices[i]-fVertices[i+4]));
339 }
340
341 innew=InsidePolygone(p,xy);
342
343 if( (innew==kInside) || (innew==kSurface) )
344 {
345 if(std::fabs(p.z()) > fDz-halfCarTolerance) { innew=kSurface; }
346 }
347 }
348 return innew;
349}
350
351// --------------------------------------------------------------------
352
354{
355 // Calculate side nearest to p, and return normal
356 // If two sides are equidistant, sum of the Normal is returned
357
358#ifdef G4TESS_TEST
359 if ( fTessellatedSolid )
360 {
361 return fTessellatedSolid->SurfaceNormal(p);
362 }
363#endif
364
365 G4ThreeVector lnorm, sumnorm(0.,0.,0.), apprnorm(0.,0.,1.),
366 p0, p1, p2, r1, r2, r3, r4;
367 G4int noSurfaces = 0;
368 G4double distxy,distz;
369 G4bool zPlusSide=false;
370
371 distz = fDz-std::fabs(p.z());
372 if (distz < halfCarTolerance)
373 {
374 if(p.z()>0)
375 {
376 zPlusSide=true;
377 sumnorm=G4ThreeVector(0,0,1);
378 }
379 else
380 {
381 sumnorm=G4ThreeVector(0,0,-1);
382 }
383 ++noSurfaces;
384 }
385
386 // Check lateral planes
387 //
388 std:: vector<G4TwoVector> vertices;
389 G4double cf = 0.5*(fDz-p.z())/fDz;
390 for (auto i=0; i<4; ++i)
391 {
392 vertices.push_back(fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]));
393 }
394
395 // Compute distance for lateral planes
396 //
397 for (G4int q=0; q<4; ++q)
398 {
399 p0=G4ThreeVector(vertices[q].x(),vertices[q].y(),p.z());
400 if(zPlusSide)
401 {
402 p1=G4ThreeVector(fVertices[q].x(),fVertices[q].y(),-fDz);
403 }
404 else
405 {
406 p1=G4ThreeVector(fVertices[q+4].x(),fVertices[q+4].y(),fDz);
407 }
408 p2=G4ThreeVector(vertices[(q+1)%4].x(),vertices[(q+1)%4].y(),p.z());
409
410 // Collapsed vertices
411 //
412 if ( (p2-p0).mag2() < kCarTolerance )
413 {
414 if ( std::fabs(p.z()+fDz) > kCarTolerance )
415 {
416 p2=G4ThreeVector(fVertices[(q+1)%4].x(),fVertices[(q+1)%4].y(),-fDz);
417 }
418 else
419 {
420 p2=G4ThreeVector(fVertices[(q+1)%4+4].x(),fVertices[(q+1)%4+4].y(),fDz);
421 }
422 }
423 lnorm = (p1-p0).cross(p2-p0);
424 lnorm = lnorm.unit();
425 if(zPlusSide) { lnorm=-lnorm; }
426
427 // Adjust Normal for Twisted Surface
428 //
429 if ( (fIsTwisted) && (GetTwistAngle(q)!=0) )
430 {
431 G4double normP=(p2-p0).mag();
432 if(normP != 0.0)
433 {
434 G4double proj=(p-p0).dot(p2-p0)/normP;
435 if(proj<0) { proj=0; }
436 if(proj>normP) { proj=normP; }
437 G4int j=(q+1)%4;
438 r1=G4ThreeVector(fVertices[q+4].x(),fVertices[q+4].y(),fDz);
439 r2=G4ThreeVector(fVertices[j+4].x(),fVertices[j+4].y(),fDz);
440 r3=G4ThreeVector(fVertices[q].x(),fVertices[q].y(),-fDz);
441 r4=G4ThreeVector(fVertices[j].x(),fVertices[j].y(),-fDz);
442 r1=r1+proj*(r2-r1)/normP;
443 r3=r3+proj*(r4-r3)/normP;
444 r2=r1-r3;
445 r4=r2.cross(p2-p0); r4=r4.unit();
446 lnorm=r4;
447 }
448 } // End if fIsTwisted
449
450 distxy=std::fabs((p0-p).dot(lnorm));
451 if ( distxy<halfCarTolerance )
452 {
453 ++noSurfaces;
454
455 // Negative sign for Normal is taken for Outside Normal
456 //
457 sumnorm=sumnorm+lnorm;
458 }
459
460 // For ApproxSurfaceNormal
461 //
462 if (distxy<distz)
463 {
464 distz=distxy;
465 apprnorm=lnorm;
466 }
467 } // End for loop
468
469 // Calculate final Normal, add Normal in the Corners and Touching Sides
470 //
471 if ( noSurfaces == 0 )
472 {
473#ifdef G4SPECSDEBUG
474 G4Exception("G4GenericTrap::SurfaceNormal(p)", "GeomSolids1002",
475 JustWarning, "Point p is not on surface !?" );
476#endif
477 sumnorm=apprnorm;
478 // Add Approximative Surface Normal Calculation?
479 }
480 else if ( noSurfaces == 1 ) { ; }
481 else { sumnorm = sumnorm.unit(); }
482
483 return sumnorm ;
484}
485
486// --------------------------------------------------------------------
487
488G4ThreeVector G4GenericTrap::NormalToPlane( const G4ThreeVector& p,
489 const G4int ipl ) const
490{
491 // Return normal to given lateral plane ipl
492
493#ifdef G4TESS_TEST
494 if ( fTessellatedSolid )
495 {
496 return fTessellatedSolid->SurfaceNormal(p);
497 }
498#endif
499
500 G4ThreeVector lnorm, norm(0.,0.,0.), p0,p1,p2;
501
502 G4double distz = fDz-p.z();
503 G4int i=ipl; // current plane index
504
505 G4TwoVector u,v;
506 G4ThreeVector r1,r2,r3,r4;
507 G4double cf = 0.5*(fDz-p.z())/fDz;
508 G4int j=(i+1)%4;
509
510 u=fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]);
511 v=fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]);
512
513 // Compute cross product
514 //
515 p0=G4ThreeVector(u.x(),u.y(),p.z());
516
517 if (std::fabs(distz)<halfCarTolerance)
518 {
519 p1=G4ThreeVector(fVertices[i].x(),fVertices[i].y(),-fDz);
520 distz=-1;
521 }
522 else
523 {
524 p1=G4ThreeVector(fVertices[i+4].x(),fVertices[i+4].y(),fDz);
525 }
526 p2=G4ThreeVector(v.x(),v.y(),p.z());
527
528 // Collapsed vertices
529 //
530 if ( (p2-p0).mag2() < kCarTolerance )
531 {
532 if ( std::fabs(p.z()+fDz) > halfCarTolerance )
533 {
534 p2=G4ThreeVector(fVertices[j].x(),fVertices[j].y(),-fDz);
535 }
536 else
537 {
538 p2=G4ThreeVector(fVertices[j+4].x(),fVertices[j+4].y(),fDz);
539 }
540 }
541 lnorm=-(p1-p0).cross(p2-p0);
542 if (distz>-halfCarTolerance) { lnorm=-lnorm.unit(); }
543 else { lnorm=lnorm.unit(); }
544
545 // Adjust Normal for Twisted Surface
546 //
547 if( (fIsTwisted) && (GetTwistAngle(ipl)!=0) )
548 {
549 G4double normP=(p2-p0).mag();
550 if(normP != 0.0)
551 {
552 G4double proj=(p-p0).dot(p2-p0)/normP;
553 if (proj<0) { proj=0; }
554 if (proj>normP) { proj=normP; }
555
556 r1=G4ThreeVector(fVertices[i+4].x(),fVertices[i+4].y(),fDz);
557 r2=G4ThreeVector(fVertices[j+4].x(),fVertices[j+4].y(),fDz);
558 r3=G4ThreeVector(fVertices[i].x(),fVertices[i].y(),-fDz);
559 r4=G4ThreeVector(fVertices[j].x(),fVertices[j].y(),-fDz);
560 r1=r1+proj*(r2-r1)/normP;
561 r3=r3+proj*(r4-r3)/normP;
562 r2=r1-r3;
563 r4=r2.cross(p2-p0);r4=r4.unit();
564 lnorm=r4;
565 }
566 } // End if fIsTwisted
567
568 return lnorm;
569}
570
571// --------------------------------------------------------------------
572
573G4double G4GenericTrap::DistToPlane(const G4ThreeVector& p,
574 const G4ThreeVector& v,
575 const G4int ipl) const
576{
577 // Computes distance to plane ipl :
578 // ipl=0 : points 0,4,1,5
579 // ipl=1 : points 1,5,2,6
580 // ipl=2 : points 2,6,3,7
581 // ipl=3 : points 3,7,0,4
582
583 G4double xa,xb,xc,xd,ya,yb,yc,yd;
584
585 G4int j = (ipl+1)%4;
586
587 xa=fVertices[ipl].x();
588 ya=fVertices[ipl].y();
589 xb=fVertices[ipl+4].x();
590 yb=fVertices[ipl+4].y();
591 xc=fVertices[j].x();
592 yc=fVertices[j].y();
593 xd=fVertices[4+j].x();
594 yd=fVertices[4+j].y();
595
596 G4double dz2 =0.5/fDz;
597 G4double tx1 =dz2*(xb-xa);
598 G4double ty1 =dz2*(yb-ya);
599 G4double tx2 =dz2*(xd-xc);
600 G4double ty2 =dz2*(yd-yc);
601 G4double dzp =fDz+p.z();
602 G4double xs1 =xa+tx1*dzp;
603 G4double ys1 =ya+ty1*dzp;
604 G4double xs2 =xc+tx2*dzp;
605 G4double ys2 =yc+ty2*dzp;
606 G4double dxs =xs2-xs1;
607 G4double dys =ys2-ys1;
608 G4double dtx =tx2-tx1;
609 G4double dty =ty2-ty1;
610
611 G4double a = (dtx*v.y()-dty*v.x()+(tx1*ty2-tx2*ty1)*v.z())*v.z();
612 G4double b = dxs*v.y()-dys*v.x()+(dtx*p.y()-dty*p.x()+ty2*xs1-ty1*xs2
613 + tx1*ys2-tx2*ys1)*v.z();
614 G4double c=dxs*p.y()-dys*p.x()+xs1*ys2-xs2*ys1;
615 G4double q=kInfinity;
616 G4double x1,x2,y1,y2,xp,yp,zi;
617
618 if (std::fabs(a)<kCarTolerance)
619 {
620 if (std::fabs(b)<kCarTolerance) { return kInfinity; }
621 q=-c/b;
622
623 // Check if Point is on the Surface
624
625 if (q>-halfCarTolerance)
626 {
627 if (q<halfCarTolerance)
628 {
629 if (NormalToPlane(p,ipl).dot(v)<=0)
630 { if(Inside(p) != kOutside) { return 0.; } }
631 else
632 { return kInfinity; }
633 }
634
635 // Check the Intersection
636 //
637 zi=p.z()+q*v.z();
638 if (std::fabs(zi)<fDz)
639 {
640 x1=xs1+tx1*v.z()*q;
641 x2=xs2+tx2*v.z()*q;
642 xp=p.x()+q*v.x();
643 y1=ys1+ty1*v.z()*q;
644 y2=ys2+ty2*v.z()*q;
645 yp=p.y()+q*v.y();
646 zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
647 if (zi<=halfCarTolerance) { return q; }
648 }
649 }
650 return kInfinity;
651 }
652 G4double d=b*b-4*a*c;
653 if (d>=0)
654 {
655 if (a>0) { q=0.5*(-b-std::sqrt(d))/a; }
656 else { q=0.5*(-b+std::sqrt(d))/a; }
657
658 // Check if Point is on the Surface
659 //
660 if (q>-halfCarTolerance)
661 {
662 if(q<halfCarTolerance)
663 {
664 if (NormalToPlane(p,ipl).dot(v)<=0)
665 {
666 if(Inside(p)!= kOutside) { return 0.; }
667 }
668 else // Check second root; return kInfinity
669 {
670 if (a>0) { q=0.5*(-b+std::sqrt(d))/a; }
671 else { q=0.5*(-b-std::sqrt(d))/a; }
672 if (q<=halfCarTolerance) { return kInfinity; }
673 }
674 }
675 // Check the Intersection
676 //
677 zi=p.z()+q*v.z();
678 if (std::fabs(zi)<fDz)
679 {
680 x1=xs1+tx1*v.z()*q;
681 x2=xs2+tx2*v.z()*q;
682 xp=p.x()+q*v.x();
683 y1=ys1+ty1*v.z()*q;
684 y2=ys2+ty2*v.z()*q;
685 yp=p.y()+q*v.y();
686 zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
687 if (zi<=halfCarTolerance) { return q; }
688 }
689 }
690 if (a>0) { q=0.5*(-b+std::sqrt(d))/a; }
691 else { q=0.5*(-b-std::sqrt(d))/a; }
692
693 // Check if Point is on the Surface
694 //
695 if (q>-halfCarTolerance)
696 {
697 if(q<halfCarTolerance)
698 {
699 if (NormalToPlane(p,ipl).dot(v)<=0)
700 {
701 if(Inside(p) != kOutside) { return 0.; }
702 }
703 else // Check second root; return kInfinity.
704 {
705 if (a>0) { q=0.5*(-b-std::sqrt(d))/a; }
706 else { q=0.5*(-b+std::sqrt(d))/a; }
707 if (q<=halfCarTolerance) { return kInfinity; }
708 }
709 }
710 // Check the Intersection
711 //
712 zi=p.z()+q*v.z();
713 if (std::fabs(zi)<fDz)
714 {
715 x1=xs1+tx1*v.z()*q;
716 x2=xs2+tx2*v.z()*q;
717 xp=p.x()+q*v.x();
718 y1=ys1+ty1*v.z()*q;
719 y2=ys2+ty2*v.z()*q;
720 yp=p.y()+q*v.y();
721 zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
722 if (zi<=halfCarTolerance) { return q; }
723 }
724 }
725 }
726 return kInfinity;
727}
728
729// --------------------------------------------------------------------
730
732 const G4ThreeVector& v) const
733{
734#ifdef G4TESS_TEST
735 if ( fTessellatedSolid )
736 {
737 return fTessellatedSolid->DistanceToIn(p, v);
738 }
739#endif
740
741 G4double dist[5];
743
744 // Check lateral faces
745 //
746 G4int i;
747 for (i=0; i<4; ++i)
748 {
749 dist[i]=DistToPlane(p, v, i);
750 }
751
752 // Check Z planes
753 //
754 dist[4]=kInfinity;
755 if (std::fabs(p.z())>fDz-halfCarTolerance)
756 {
757 if (v.z() != 0.0)
758 {
759 G4ThreeVector pt;
760 if (p.z()>0)
761 {
762 dist[4] = (fDz-p.z())/v.z();
763 }
764 else
765 {
766 dist[4] = (-fDz-p.z())/v.z();
767 }
768 if (dist[4]<-halfCarTolerance)
769 {
770 dist[4]=kInfinity;
771 }
772 else
773 {
774 if(dist[4]<halfCarTolerance)
775 {
776 if(p.z()>0) { n=G4ThreeVector(0,0,1); }
777 else { n=G4ThreeVector(0,0,-1); }
778 if (n.dot(v)<0) { dist[4]=0.; }
779 else { dist[4]=kInfinity; }
780 }
781 pt=p+dist[4]*v;
782 if (Inside(pt)==kOutside) { dist[4]=kInfinity; }
783 }
784 }
785 }
786 G4double distmin = dist[0];
787 for (i=1; i<5 ; ++i)
788 {
789 if (dist[i] < distmin) { distmin = dist[i]; }
790 }
791
792 if (distmin<halfCarTolerance) { distmin=0.; }
793
794 return distmin;
795}
796
797// --------------------------------------------------------------------
798
800{
801 // Computes the closest distance from given point to this shape
802
803#ifdef G4TESS_TEST
804 if ( fTessellatedSolid )
805 {
806 return fTessellatedSolid->DistanceToIn(p);
807 }
808#endif
809
810 G4double safz = std::fabs(p.z())-fDz;
811 if(safz<0) { safz=0; }
812
813 G4int iseg;
814 G4double safe = safz;
815 G4double safxy = safz;
816
817 for (iseg=0; iseg<4; ++iseg)
818 {
819 safxy = SafetyToFace(p,iseg);
820 if (safxy>safe) { safe=safxy; }
821 }
822
823 return safe;
824}
825
826// --------------------------------------------------------------------
827
829G4GenericTrap::SafetyToFace(const G4ThreeVector& p, const G4int iseg) const
830{
831 // Estimate distance to lateral plane defined by segment iseg in range [0,3]
832 // Might be negative: plane seen only from inside
833
834 G4ThreeVector p1,norm;
835 G4double safe;
836
837 p1=G4ThreeVector(fVertices[iseg].x(),fVertices[iseg].y(),-fDz);
838
839 norm=NormalToPlane(p,iseg);
840 safe = (p-p1).dot(norm); // Can be negative
841
842 return safe;
843}
844
845// --------------------------------------------------------------------
846
848G4GenericTrap::DistToTriangle(const G4ThreeVector& p,
849 const G4ThreeVector& v, const G4int ipl) const
850{
851 G4double xa=fVertices[ipl].x();
852 G4double ya=fVertices[ipl].y();
853 G4double xb=fVertices[ipl+4].x();
854 G4double yb=fVertices[ipl+4].y();
855 G4int j=(ipl+1)%4;
856 G4double xc=fVertices[j].x();
857 G4double yc=fVertices[j].y();
858 G4double zab=2*fDz;
859 G4double zac=0;
860
861 if ( (std::fabs(xa-xc)+std::fabs(ya-yc)) < halfCarTolerance )
862 {
863 xc=fVertices[j+4].x();
864 yc=fVertices[j+4].y();
865 zac=2*fDz;
866 zab=2*fDz;
867
868 //Line case
869 //
870 if ( (std::fabs(xb-xc)+std::fabs(yb-yc)) < halfCarTolerance )
871 {
872 return kInfinity;
873 }
874 }
875 G4double a=(yb-ya)*zac-(yc-ya)*zab;
876 G4double b=(xc-xa)*zab-(xb-xa)*zac;
877 G4double c=(xb-xa)*(yc-ya)-(xc-xa)*(yb-ya);
878 G4double d=-xa*a-ya*b+fDz*c;
879 G4double t=a*v.x()+b*v.y()+c*v.z();
880
881 if (t!=0)
882 {
883 t=-(a*p.x()+b*p.y()+c*p.z()+d)/t;
884 }
885 if ( (t<halfCarTolerance) && (t>-halfCarTolerance) )
886 {
887 if (NormalToPlane(p,ipl).dot(v)<kCarTolerance)
888 {
889 t=kInfinity;
890 }
891 else
892 {
893 t=0;
894 }
895 }
896 if (Inside(p+v*t) != kSurface) { t=kInfinity; }
897
898 return t;
899}
900
901// --------------------------------------------------------------------
902
904 const G4ThreeVector& v,
905 const G4bool calcNorm,
906 G4bool* validNorm,
907 G4ThreeVector* n) const
908{
909#ifdef G4TESS_TEST
910 if ( fTessellatedSolid )
911 {
912 return fTessellatedSolid->DistanceToOut(p, v, calcNorm, validNorm, n);
913 }
914#endif
915
916 G4double distmin;
917 G4bool lateral_cross = false;
918 ESide side = kUndef;
919
920 if (calcNorm) { *validNorm = true; } // All normals are valid
921
922 if (v.z() < 0)
923 {
924 distmin=(-fDz-p.z())/v.z();
925 if (calcNorm) { side=kMZ; *n=G4ThreeVector(0,0,-1); }
926 }
927 else
928 {
929 if (v.z() > 0)
930 {
931 distmin = (fDz-p.z())/v.z();
932 if (calcNorm) { side=kPZ; *n=G4ThreeVector(0,0,1); }
933 }
934 else { distmin = kInfinity; }
935 }
936
937 G4double dz2 =0.5/fDz;
938 G4double xa,xb,xc,xd;
939 G4double ya,yb,yc,yd;
940
941 for (G4int ipl=0; ipl<4; ++ipl)
942 {
943 G4int j = (ipl+1)%4;
944 xa=fVertices[ipl].x();
945 ya=fVertices[ipl].y();
946 xb=fVertices[ipl+4].x();
947 yb=fVertices[ipl+4].y();
948 xc=fVertices[j].x();
949 yc=fVertices[j].y();
950 xd=fVertices[4+j].x();
951 yd=fVertices[4+j].y();
952
953 if ( ((std::fabs(xb-xd)+std::fabs(yb-yd))<halfCarTolerance)
954 || ((std::fabs(xa-xc)+std::fabs(ya-yc))<halfCarTolerance) )
955 {
956 G4double q=DistToTriangle(p,v,ipl) ;
957 if ( (q>=0) && (q<distmin) )
958 {
959 distmin=q;
960 lateral_cross=true;
961 side=ESide(ipl+1);
962 }
963 continue;
964 }
965 G4double tx1 =dz2*(xb-xa);
966 G4double ty1 =dz2*(yb-ya);
967 G4double tx2 =dz2*(xd-xc);
968 G4double ty2 =dz2*(yd-yc);
969 G4double dzp =fDz+p.z();
970 G4double xs1 =xa+tx1*dzp;
971 G4double ys1 =ya+ty1*dzp;
972 G4double xs2 =xc+tx2*dzp;
973 G4double ys2 =yc+ty2*dzp;
974 G4double dxs =xs2-xs1;
975 G4double dys =ys2-ys1;
976 G4double dtx =tx2-tx1;
977 G4double dty =ty2-ty1;
978 G4double a = (dtx*v.y()-dty*v.x()+(tx1*ty2-tx2*ty1)*v.z())*v.z();
979 G4double b = dxs*v.y()-dys*v.x()+(dtx*p.y()-dty*p.x()+ty2*xs1-ty1*xs2
980 + tx1*ys2-tx2*ys1)*v.z();
981 G4double c=dxs*p.y()-dys*p.x()+xs1*ys2-xs2*ys1;
982 G4double q=kInfinity;
983
984 if (std::fabs(a) < kCarTolerance)
985 {
986 if (std::fabs(b) < kCarTolerance) { continue; }
987 q=-c/b;
988
989 // Check for Point on the Surface
990 //
991 if ((q > -halfCarTolerance) && (q < distmin))
992 {
993 if (q < halfCarTolerance)
994 {
995 if (NormalToPlane(p,ipl).dot(v)<0.) { continue; }
996 }
997 distmin =q;
998 lateral_cross=true;
999 side=ESide(ipl+1);
1000 }
1001 continue;
1002 }
1003 G4double d=b*b-4*a*c;
1004 if (d >= 0.)
1005 {
1006 if (a > 0) { q=0.5*(-b-std::sqrt(d))/a; }
1007 else { q=0.5*(-b+std::sqrt(d))/a; }
1008
1009 // Check for Point on the Surface
1010 //
1011 if (q > -halfCarTolerance )
1012 {
1013 if (q < distmin)
1014 {
1015 if(q < halfCarTolerance)
1016 {
1017 if (NormalToPlane(p,ipl).dot(v)<0.) // Check second root
1018 {
1019 if (a > 0) { q=0.5*(-b+std::sqrt(d))/a; }
1020 else { q=0.5*(-b-std::sqrt(d))/a; }
1021 if (( q > halfCarTolerance) && (q < distmin))
1022 {
1023 distmin=q;
1024 lateral_cross = true;
1025 side=ESide(ipl+1);
1026 }
1027 continue;
1028 }
1029 }
1030 distmin = q;
1031 lateral_cross = true;
1032 side=ESide(ipl+1);
1033 }
1034 }
1035 else
1036 {
1037 if (a > 0) { q=0.5*(-b+std::sqrt(d))/a; }
1038 else { q=0.5*(-b-std::sqrt(d))/a; }
1039
1040 // Check for Point on the Surface
1041 //
1042 if ((q > -halfCarTolerance) && (q < distmin))
1043 {
1044 if (q < halfCarTolerance)
1045 {
1046 if (NormalToPlane(p,ipl).dot(v)<0.) // Check second root
1047 {
1048 if (a > 0) { q=0.5*(-b-std::sqrt(d))/a; }
1049 else { q=0.5*(-b+std::sqrt(d))/a; }
1050 if ( ( q > halfCarTolerance) && (q < distmin) )
1051 {
1052 distmin=q;
1053 lateral_cross = true;
1054 side=ESide(ipl+1);
1055 }
1056 continue;
1057 }
1058 }
1059 distmin =q;
1060 lateral_cross = true;
1061 side=ESide(ipl+1);
1062 }
1063 }
1064 }
1065 }
1066 if (!lateral_cross) // Make sure that track crosses the top or bottom
1067 {
1068 if (distmin >= kInfinity) { distmin=kCarTolerance; }
1069 G4ThreeVector pt=p+distmin*v;
1070
1071 // Check if propagated point is in the polygon
1072 //
1073 G4int i=0;
1074 if (v.z()>0.) { i=4; }
1075 std::vector<G4TwoVector> xy;
1076 for ( G4int j=0; j<4; j++) { xy.push_back(fVertices[i+j]); }
1077
1078 // Check Inside
1079 //
1080 if (InsidePolygone(pt,xy)==kOutside)
1081 {
1082 if(calcNorm)
1083 {
1084 if (v.z()>0) {side= kPZ; *n = G4ThreeVector(0,0,1);}
1085 else { side=kMZ; *n = G4ThreeVector(0,0,-1);}
1086 }
1087 return 0.;
1088 }
1089 else
1090 {
1091 if(v.z()>0) {side=kPZ;}
1092 else {side=kMZ;}
1093 }
1094 }
1095
1096 if (calcNorm)
1097 {
1098 G4ThreeVector pt=p+v*distmin;
1099 switch (side)
1100 {
1101 case kXY0:
1102 *n=NormalToPlane(pt,0);
1103 break;
1104 case kXY1:
1105 *n=NormalToPlane(pt,1);
1106 break;
1107 case kXY2:
1108 *n=NormalToPlane(pt,2);
1109 break;
1110 case kXY3:
1111 *n=NormalToPlane(pt,3);
1112 break;
1113 case kPZ:
1114 *n=G4ThreeVector(0,0,1);
1115 break;
1116 case kMZ:
1117 *n=G4ThreeVector(0,0,-1);
1118 break;
1119 default:
1120 DumpInfo();
1121 std::ostringstream message;
1122 G4long oldprc = message.precision(16);
1123 message << "Undefined side for valid surface normal to solid." << G4endl
1124 << "Position:" << G4endl
1125 << " p.x() = " << p.x()/mm << " mm" << G4endl
1126 << " p.y() = " << p.y()/mm << " mm" << G4endl
1127 << " p.z() = " << p.z()/mm << " mm" << G4endl
1128 << "Direction:" << G4endl
1129 << " v.x() = " << v.x() << G4endl
1130 << " v.y() = " << v.y() << G4endl
1131 << " v.z() = " << v.z() << G4endl
1132 << "Proposed distance :" << G4endl
1133 << " distmin = " << distmin/mm << " mm";
1134 message.precision(oldprc);
1135 G4Exception("G4GenericTrap::DistanceToOut(p,v,..)",
1136 "GeomSolids1002", JustWarning, message);
1137 break;
1138 }
1139 }
1140
1141 if (distmin<halfCarTolerance) { distmin=0.; }
1142
1143 return distmin;
1144}
1145
1146// --------------------------------------------------------------------
1147
1149{
1150
1151#ifdef G4TESS_TEST
1152 if ( fTessellatedSolid )
1153 {
1154 return fTessellatedSolid->DistanceToOut(p);
1155 }
1156#endif
1157
1158 G4double safz = fDz-std::fabs(p.z());
1159 if (safz<0) { safz = 0; }
1160
1161 G4double safe = safz;
1162 G4double safxy = safz;
1163
1164 for (G4int iseg=0; iseg<4; ++iseg)
1165 {
1166 safxy = std::fabs(SafetyToFace(p,iseg));
1167 if (safxy < safe) { safe = safxy; }
1168 }
1169
1170 return safe;
1171}
1172
1173// --------------------------------------------------------------------
1174
1176 G4ThreeVector& pMax) const
1177{
1178 pMin = GetMinimumBBox();
1179 pMax = GetMaximumBBox();
1180
1181 // Check correctness of the bounding box
1182 //
1183 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
1184 {
1185 std::ostringstream message;
1186 message << "Bad bounding box (min >= max) for solid: "
1187 << GetName() << " !"
1188 << "\npMin = " << pMin
1189 << "\npMax = " << pMax;
1190 G4Exception("G4GenericTrap::BoundingLimits()", "GeomMgt0001",
1191 JustWarning, message);
1192 DumpInfo();
1193 }
1194}
1195
1196// --------------------------------------------------------------------
1197
1198G4bool
1200 const G4VoxelLimits& pVoxelLimit,
1201 const G4AffineTransform& pTransform,
1202 G4double& pMin, G4double& pMax) const
1203{
1204 G4ThreeVector bmin, bmax;
1205 G4bool exist;
1206
1207 // Check bounding box (bbox)
1208 //
1209 BoundingLimits(bmin,bmax);
1210 G4BoundingEnvelope bbox(bmin,bmax);
1211#ifdef G4BBOX_EXTENT
1212 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
1213#endif
1214 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
1215 {
1216 return exist = pMin < pMax;
1217 }
1218
1219 // Set bounding envelope (benv) and calculate extent
1220 //
1221 // To build the bounding envelope with plane faces each side face of
1222 // the trapezoid is subdivided in triangles. Subdivision is done by
1223 // duplication of vertices in the bases in a way that the envelope be
1224 // a convex polyhedron (some faces of the envelope can be degenerate)
1225 //
1226 G4double dz = GetZHalfLength();
1227 G4ThreeVectorList baseA(8), baseB(8);
1228 for (G4int i=0; i<4; ++i)
1229 {
1230 G4TwoVector va = GetVertex(i);
1231 G4TwoVector vb = GetVertex(i+4);
1232 baseA[2*i].set(va.x(),va.y(),-dz);
1233 baseB[2*i].set(vb.x(),vb.y(), dz);
1234 }
1235 for (G4int i=0; i<4; ++i)
1236 {
1237 G4int k1=2*i, k2=(2*i+2)%8;
1238 G4double ax = (baseA[k2].x()-baseA[k1].x());
1239 G4double ay = (baseA[k2].y()-baseA[k1].y());
1240 G4double bx = (baseB[k2].x()-baseB[k1].x());
1241 G4double by = (baseB[k2].y()-baseB[k1].y());
1242 G4double znorm = ax*by - ay*bx;
1243 baseA[k1+1] = (znorm < 0.0) ? baseA[k2] : baseA[k1];
1244 baseB[k1+1] = (znorm < 0.0) ? baseB[k1] : baseB[k2];
1245 }
1246
1247 std::vector<const G4ThreeVectorList *> polygons(2);
1248 polygons[0] = &baseA;
1249 polygons[1] = &baseB;
1250
1251 G4BoundingEnvelope benv(bmin,bmax,polygons);
1252 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
1253 return exist;
1254}
1255
1256// --------------------------------------------------------------------
1257
1259{
1260 return {"G4GenericTrap"};
1261}
1262
1263// --------------------------------------------------------------------
1264
1266{
1267 return new G4GenericTrap(*this);
1268}
1269
1270// --------------------------------------------------------------------
1271
1272std::ostream& G4GenericTrap::StreamInfo(std::ostream& os) const
1273{
1274 G4long oldprc = os.precision(16);
1275 os << "-----------------------------------------------------------\n"
1276 << " *** Dump for solid - " << GetName() << " *** \n"
1277 << " =================================================== \n"
1278 << " Solid geometry type: " << GetEntityType() << G4endl
1279 << " half length Z: " << fDz/mm << " mm \n"
1280 << " list of vertices:\n";
1281
1282 for ( G4int i=0; i<fgkNofVertices; ++i )
1283 {
1284 os << std::setw(5) << "#" << i
1285 << " vx = " << fVertices[i].x()/mm << " mm"
1286 << " vy = " << fVertices[i].y()/mm << " mm" << G4endl;
1287 }
1288 os.precision(oldprc);
1289
1290 return os;
1291}
1292
1293// --------------------------------------------------------------------
1294
1296{
1297
1298#ifdef G4TESS_TEST
1299 if ( fTessellatedSolid )
1300 {
1301 return fTessellatedSolid->GetPointOnSurface();
1302 }
1303#endif
1304
1305 G4ThreeVector point;
1306 G4TwoVector u,v,w;
1307 G4double rand,area,chose,cf,lambda0,lambda1,alfa,beta,zp;
1308 G4int ipl,j;
1309
1310 std::vector<G4ThreeVector> vertices;
1311 for (auto i=0; i<4; ++i)
1312 {
1313 vertices.emplace_back(fVertices[i].x(),fVertices[i].y(),-fDz);
1314 }
1315 for (auto i=4; i<8; ++i)
1316 {
1317 vertices.emplace_back(fVertices[i].x(),fVertices[i].y(),fDz);
1318 }
1319
1320 // Surface Area of Planes
1321 //
1322 G4TwoVector A = fVertices[3] - fVertices[1];
1323 G4TwoVector B = fVertices[2] - fVertices[0];
1324 G4TwoVector C = fVertices[7] - fVertices[5];
1325 G4TwoVector D = fVertices[6] - fVertices[4];
1326 G4double Surface0 = 0.5*(A.x()*B.y() - A.y()*B.x()); //-fDz plane
1327 G4double Surface1 = GetLateralFaceArea(0);
1328 G4double Surface2 = GetLateralFaceArea(1);
1329 G4double Surface3 = GetLateralFaceArea(2);
1330 G4double Surface4 = GetLateralFaceArea(3);
1331 G4double Surface5 = 0.5*(C.x()*D.y() - C.y()*D.x()); // fDz plane
1332
1333 rand = G4UniformRand();
1334 area = Surface0+Surface1+Surface2+Surface3+Surface4+Surface5;
1335 chose = rand*area;
1336
1337 if ( ( chose < Surface0)
1338 || ( chose > (Surface0+Surface1+Surface2+Surface3+Surface4)) )
1339 { // fDz or -fDz Plane
1340 ipl = G4int(G4UniformRand()*4);
1341 j = (ipl+1)%4;
1342 if(chose < Surface0)
1343 {
1344 zp = -fDz;
1345 u = fVertices[ipl]; v = fVertices[j];
1346 w = fVertices[(ipl+3)%4];
1347 }
1348 else
1349 {
1350 zp = fDz;
1351 u = fVertices[ipl+4]; v = fVertices[j+4];
1352 w = fVertices[(ipl+3)%4+4];
1353 }
1354 alfa = G4UniformRand();
1355 beta = G4UniformRand();
1356 lambda1=alfa*beta;
1357 lambda0=alfa-lambda1;
1358 v = v-u;
1359 w = w-u;
1360 v = u+lambda0*v+lambda1*w;
1361 }
1362 else // Lateral Plane Twisted or Not
1363 {
1364 if (chose < Surface0+Surface1) { ipl=0; }
1365 else if (chose < Surface0+Surface1+Surface2) { ipl=1; }
1366 else if (chose < Surface0+Surface1+Surface2+Surface3) { ipl=2; }
1367 else { ipl=3; }
1368 j = (ipl+1)%4;
1369 zp = -fDz+G4UniformRand()*2*fDz;
1370 cf = 0.5*(fDz-zp)/fDz;
1371 u = fVertices[ipl+4]+cf*( fVertices[ipl]-fVertices[ipl+4]);
1372 v = fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]);
1373 v = u+(v-u)*G4UniformRand();
1374 }
1375 point=G4ThreeVector(v.x(),v.y(),zp);
1376
1377 return point;
1378}
1379
1380// --------------------------------------------------------------------
1381
1383{
1384 if (fCubicVolume == 0.0)
1385 {
1386 // diagonals
1387 G4TwoVector A = fVertices[3] - fVertices[1];
1388 G4TwoVector B = fVertices[2] - fVertices[0];
1389 G4TwoVector C = fVertices[7] - fVertices[5];
1390 G4TwoVector D = fVertices[6] - fVertices[4];
1391
1392 // kross products
1393 G4double AB = A.x()*B.y() - A.y()*B.x();
1394 G4double CD = C.x()*D.y() - C.y()*D.x();
1395 G4double AD = A.x()*D.y() - A.y()*D.x();
1396 G4double CB = C.x()*B.y() - C.y()*B.x();
1397
1398 fCubicVolume = ((AB + CD)/3. + (AD + CB)/6.)*fDz;
1399 }
1400 return fCubicVolume;
1401}
1402
1403// --------------------------------------------------------------------
1404
1405G4double G4GenericTrap::GetLateralFaceArea(G4int iface) const
1406{
1407 constexpr G4int NSTEP = 250;
1408 constexpr G4double dt = 1./NSTEP;
1409
1410 G4int i1 = iface, i2 = (iface + 1)%4;
1411 G4int i3 = i1 + 4, i4 = i2 + 4;
1412
1413 G4double x21 = fVertices[i2].x() - fVertices[i1].x();
1414 G4double y21 = fVertices[i2].y() - fVertices[i1].y();
1415 G4double x31 = fVertices[i3].x() - fVertices[i1].x();
1416 G4double y31 = fVertices[i3].y() - fVertices[i1].y();
1417 G4double x42 = fVertices[i4].x() - fVertices[i2].x();
1418 G4double y42 = fVertices[i4].y() - fVertices[i2].y();
1419 G4double x43 = fVertices[i4].x() - fVertices[i3].x();
1420 G4double y43 = fVertices[i4].y() - fVertices[i3].y();
1421
1422 G4double A = x21*y43 - y21*x43;
1423 G4double lmax = std::max(std::max(std::abs(x21),std::abs(y21)),
1424 std::max(std::abs(x43),std::abs(y43)));
1425 G4double eps = lmax*kCarTolerance;
1426 if (std::abs(A) < eps) // plane face
1427 {
1428 G4ThreeVector p1(fVertices[i1].x(), fVertices[i1].y(),-fDz);
1429 G4ThreeVector p2(fVertices[i2].x(), fVertices[i2].y(),-fDz);
1430 G4ThreeVector p3(fVertices[i3].x(), fVertices[i3].y(), fDz);
1431 G4ThreeVector p4(fVertices[i4].x(), fVertices[i4].y(), fDz);
1432 return ((p4 - p1).cross(p3 - p2)).mag()*0.5;
1433 }
1434
1435 // twisted face
1436 G4double B0 = x21*y31 - y21*x31;
1437 G4double B1 = x42*y31 - y42*x31;
1438 G4double HH = 4*fDz*fDz;
1439 G4double invAA = 1./(A*A);
1440 G4double sqrtAA = 2.*std::abs(A);
1441 G4double invSqrtAA = 1./sqrtAA;
1442
1443 G4double area = 0.;
1444 for (G4int i = 0; i < NSTEP; ++i)
1445 {
1446 G4double t = (i + 0.5)*dt;
1447 G4double I = y21 + (y43 - y21)*t;
1448 G4double J = x21 + (x43 - x21)*t;
1449 G4double IIJJ = HH*(I*I + J*J);
1450 G4double B = B1*t + B0;
1451
1452 G4double aa = A*A;
1453 G4double bb = 2.*A*B;
1454 G4double cc = IIJJ + B*B;
1455
1456 G4double R1 = std::sqrt(aa + bb + cc);
1457 G4double R0 = std::sqrt(cc);
1458 G4double log1 = std::log(std::abs(sqrtAA*R1 + 2.*aa + bb));
1459 G4double log0 = std::log(std::abs(sqrtAA*R0 + bb));
1460
1461 area += 0.5*R1 + 0.25*bb*invAA*(R1 - R0) + IIJJ*invSqrtAA*(log1 - log0);
1462 }
1463 return area*dt;
1464}
1465
1466// --------------------------------------------------------------------
1467
1469{
1470 if (fSurfaceArea == 0.0)
1471 {
1472 G4TwoVector A = fVertices[3] - fVertices[1];
1473 G4TwoVector B = fVertices[2] - fVertices[0];
1474 G4TwoVector C = fVertices[7] - fVertices[5];
1475 G4TwoVector D = fVertices[6] - fVertices[4];
1476 G4double S_bot = 0.5*(A.x()*B.y() - A.y()*B.x());
1477 G4double S_top = 0.5*(C.x()*D.y() - C.y()*D.x());
1478 fSurfaceArea = S_bot + S_top +
1479 GetLateralFaceArea(0) +
1480 GetLateralFaceArea(1) +
1481 GetLateralFaceArea(2) +
1482 GetLateralFaceArea(3);
1483 }
1484 return fSurfaceArea;
1485}
1486
1487// --------------------------------------------------------------------
1488
1489G4bool G4GenericTrap::ComputeIsTwisted()
1490{
1491 // Computes tangents of twist angles (angles between projections on XY plane
1492 // of corresponding -dz +dz edges).
1493
1494 G4bool twisted = false;
1495 G4double dx1, dy1, dx2, dy2;
1496 G4int nv = fgkNofVertices/2;
1497
1498 for ( G4int i=0; i<4; ++i )
1499 {
1500 dx1 = fVertices[(i+1)%nv].x()-fVertices[i].x();
1501 dy1 = fVertices[(i+1)%nv].y()-fVertices[i].y();
1502 if ( (dx1 == 0) && (dy1 == 0) ) { continue; }
1503
1504 dx2 = fVertices[nv+(i+1)%nv].x()-fVertices[nv+i].x();
1505 dy2 = fVertices[nv+(i+1)%nv].y()-fVertices[nv+i].y();
1506
1507 if ( dx2 == 0 && dy2 == 0 ) { continue; }
1508 G4double twist_angle = std::fabs(dy1*dx2 - dx1*dy2);
1509 if ( twist_angle < fgkTolerance ) { continue; }
1510 twisted = true;
1511 SetTwistAngle(i,twist_angle);
1512
1513 // Check on big angles, potentially navigation problem
1514
1515 twist_angle = std::acos( (dx1*dx2 + dy1*dy2)
1516 / (std::sqrt(dx1*dx1+dy1*dy1)
1517 * std::sqrt(dx2*dx2+dy2*dy2)) );
1518
1519 if ( std::fabs(twist_angle) > 0.5*pi+kCarTolerance )
1520 {
1521 std::ostringstream message;
1522 message << "Twisted Angle is bigger than 90 degrees - " << GetName()
1523 << G4endl
1524 << " Potential problem of malformed Solid !" << G4endl
1525 << " TwistANGLE = " << twist_angle
1526 << "*rad for lateral plane N= " << i;
1527 G4Exception("G4GenericTrap::ComputeIsTwisted()", "GeomSolids1002",
1528 JustWarning, message);
1529 }
1530 }
1531
1532 return twisted;
1533}
1534
1535// --------------------------------------------------------------------
1536
1537G4bool G4GenericTrap::CheckOrder(const std::vector<G4TwoVector>& vertices) const
1538{
1539 // Test if the vertices are in a clockwise order, if not reorder them.
1540 // Also test if they're well defined without crossing opposite segments
1541
1542 G4bool clockwise_order=true;
1543 G4double sum1 = 0.;
1544 G4double sum2 = 0.;
1545 G4int j;
1546
1547 for (G4int i=0; i<4; ++i)
1548 {
1549 j = (i+1)%4;
1550 sum1 += vertices[i].x()*vertices[j].y() - vertices[j].x()*vertices[i].y();
1551 sum2 += vertices[i+4].x()*vertices[j+4].y()
1552 - vertices[j+4].x()*vertices[i+4].y();
1553 }
1554 if (sum1*sum2 < -fgkTolerance)
1555 {
1556 std::ostringstream message;
1557 message << "Lower/upper faces defined with opposite clockwise - "
1558 << GetName();
1559 G4Exception("G4GenericTrap::CheckOrder()", "GeomSolids0002",
1560 FatalException, message);
1561 }
1562
1563 if ((sum1 > 0.)||(sum2 > 0.))
1564 {
1565 std::ostringstream message;
1566 message << "Vertices must be defined in clockwise XY planes - "
1567 << GetName();
1568 G4Exception("G4GenericTrap::CheckOrder()", "GeomSolids1001",
1569 JustWarning,message, "Re-ordering...");
1570 clockwise_order = false;
1571 }
1572
1573 // Check for illegal crossings
1574 //
1575 G4bool illegal_cross = false;
1576 illegal_cross = IsSegCrossingZ(vertices[0],vertices[4],
1577 vertices[1],vertices[5]);
1578
1579 if (!illegal_cross)
1580 {
1581 illegal_cross = IsSegCrossingZ(vertices[2],vertices[6],
1582 vertices[3],vertices[7]);
1583 }
1584 // +/- dZ planes
1585 if (!illegal_cross)
1586 {
1587 illegal_cross = IsSegCrossing(vertices[0],vertices[1],
1588 vertices[2],vertices[3]);
1589 }
1590 if (!illegal_cross)
1591 {
1592 illegal_cross = IsSegCrossing(vertices[0],vertices[3],
1593 vertices[1],vertices[2]);
1594 }
1595 if (!illegal_cross)
1596 {
1597 illegal_cross = IsSegCrossing(vertices[4],vertices[5],
1598 vertices[6],vertices[7]);
1599 }
1600 if (!illegal_cross)
1601 {
1602 illegal_cross = IsSegCrossing(vertices[4],vertices[7],
1603 vertices[5],vertices[6]);
1604 }
1605
1606 if (illegal_cross)
1607 {
1608 std::ostringstream message;
1609 message << "Malformed polygone with opposite sides - " << GetName();
1610 G4Exception("G4GenericTrap::CheckOrderAndSetup()",
1611 "GeomSolids0002", FatalException, message);
1612 }
1613 return clockwise_order;
1614}
1615
1616// --------------------------------------------------------------------
1617
1618void G4GenericTrap::ReorderVertices(std::vector<G4ThreeVector>& vertices) const
1619{
1620 // Reorder the vector of vertices
1621
1622 std::vector<G4ThreeVector> oldVertices(vertices);
1623
1624 for ( std::size_t i=0; i<oldVertices.size(); ++i )
1625 {
1626 vertices[i] = oldVertices[oldVertices.size()-1-i];
1627 }
1628}
1629
1630// --------------------------------------------------------------------
1631
1632G4bool
1633G4GenericTrap::IsSegCrossing(const G4TwoVector& a, const G4TwoVector& b,
1634 const G4TwoVector& c, const G4TwoVector& d) const
1635{
1636 // Check if segments [A,B] and [C,D] are crossing
1637
1638 G4bool stand1 = false;
1639 G4bool stand2 = false;
1640 G4double dx1,dx2,xm=0.,ym=0.,a1=0.,a2=0.,b1=0.,b2=0.;
1641 dx1=(b-a).x();
1642 dx2=(d-c).x();
1643
1644 if( std::fabs(dx1) < fgkTolerance ) { stand1 = true; }
1645 if( std::fabs(dx2) < fgkTolerance ) { stand2 = true; }
1646 if (!stand1)
1647 {
1648 a1 = (b.x()*a.y()-a.x()*b.y())/dx1;
1649 b1 = (b-a).y()/dx1;
1650 }
1651 if (!stand2)
1652 {
1653 a2 = (d.x()*c.y()-c.x()*d.y())/dx2;
1654 b2 = (d-c).y()/dx2;
1655 }
1656 if (stand1 && stand2)
1657 {
1658 // Segments parallel and vertical
1659 //
1660 if (std::fabs(a.x()-c.x())<fgkTolerance)
1661 {
1662 // Check if segments are overlapping
1663 //
1664 return ((c.y()-a.y())*(c.y()-b.y())<-fgkTolerance)
1665 || ((d.y()-a.y())*(d.y()-b.y())<-fgkTolerance)
1666 || ((a.y()-c.y())*(a.y()-d.y())<-fgkTolerance)
1667 || ((b.y()-c.y())*(b.y()-d.y())<-fgkTolerance);
1668 }
1669 // Different x values
1670 //
1671 return false;
1672 }
1673
1674 if (stand1) // First segment vertical
1675 {
1676 xm = a.x();
1677 ym = a2+b2*xm;
1678 }
1679 else
1680 {
1681 if (stand2) // Second segment vertical
1682 {
1683 xm = c.x();
1684 ym = a1+b1*xm;
1685 }
1686 else // Normal crossing
1687 {
1688 if (std::fabs(b1-b2) < fgkTolerance)
1689 {
1690 // Parallel segments, are they aligned
1691 //
1692 if (std::fabs(c.y()-(a1+b1*c.x())) > fgkTolerance) { return false; }
1693
1694 // Aligned segments, are they overlapping
1695 //
1696 if ( ((c.x()-a.x())*(c.x()-b.x())<-fgkTolerance)
1697 || ((d.x()-a.x())*(d.x()-b.x())<-fgkTolerance)
1698 || ((a.x()-c.x())*(a.x()-d.x())<-fgkTolerance)
1699 || ((b.x()-c.x())*(b.x()-d.x())<-fgkTolerance) ) { return true; }
1700
1701 return false;
1702 }
1703 xm = (a1-a2)/(b2-b1);
1704 ym = (a1*b2-a2*b1)/(b2-b1);
1705 }
1706 }
1707
1708 // Check if crossing point is both between A,B and C,D
1709 //
1710 G4double check = (xm-a.x())*(xm-b.x())+(ym-a.y())*(ym-b.y());
1711 if (check > -fgkTolerance) { return false; }
1712 check = (xm-c.x())*(xm-d.x())+(ym-c.y())*(ym-d.y());
1713 return check <= -fgkTolerance;
1714}
1715
1716// --------------------------------------------------------------------
1717
1718G4bool
1719G4GenericTrap::IsSegCrossingZ(const G4TwoVector& a, const G4TwoVector& b,
1720 const G4TwoVector& c, const G4TwoVector& d) const
1721{
1722 // Check if segments [A,B] and [C,D] are crossing when
1723 // A and C are on -dZ and B and D are on +dZ
1724
1725 // Calculate the Intersection point between two lines in 3D
1726 //
1727 G4ThreeVector temp1,temp2;
1728 G4ThreeVector v1,v2,p1,p2,p3,p4,dv;
1729 G4double q,det;
1730 p1=G4ThreeVector(a.x(),a.y(),-fDz);
1731 p2=G4ThreeVector(c.x(),c.y(),-fDz);
1732 p3=G4ThreeVector(b.x(),b.y(),fDz);
1733 p4=G4ThreeVector(d.x(),d.y(),fDz);
1734 v1=p3-p1;
1735 v2=p4-p2;
1736 dv=p2-p1;
1737
1738 // In case of Collapsed Vertices No crossing
1739 //
1740 if( (std::fabs(dv.x()) < kCarTolerance )&&
1741 (std::fabs(dv.y()) < kCarTolerance ) ) { return false; }
1742
1743 if( (std::fabs((p4-p3).x()) < kCarTolerance )&&
1744 (std::fabs((p4-p3).y()) < kCarTolerance ) ) { return false; }
1745
1746 // First estimate if Intersection is possible( if det is 0)
1747 //
1748 det = dv.x()*v1.y()*v2.z()+dv.y()*v1.z()*v2.x()
1749 - dv.x()*v1.z()*v2.y()-dv.y()*v1.x()*v2.z();
1750
1751 if (std::fabs(det)<kCarTolerance) //Intersection
1752 {
1753 temp1 = v1.cross(v2);
1754 temp2 = (p2-p1).cross(v2);
1755 if (temp1.dot(temp2) < 0) { return false; } // intersection negative
1756 q = temp1.mag();
1757
1758 if ( q < kCarTolerance ) { return false; } // parallel lines
1759 q = ((dv).cross(v2)).mag()/q;
1760
1761 if(q < 1.-kCarTolerance) { return true; }
1762 }
1763 return false;
1764}
1765
1766// --------------------------------------------------------------------
1767
1768G4VFacet*
1769G4GenericTrap::MakeDownFacet(const std::vector<G4ThreeVector>& fromVertices,
1770 G4int ind1, G4int ind2, G4int ind3) const
1771{
1772 // Create a triangular facet from the polygon points given by indices
1773 // forming the down side ( the normal goes in -z)
1774 // Do not create facet if 2 vertices are the same
1775
1776 if ( (fromVertices[ind1] == fromVertices[ind2]) ||
1777 (fromVertices[ind2] == fromVertices[ind3]) ||
1778 (fromVertices[ind1] == fromVertices[ind3]) ) { return nullptr; }
1779
1780 std::vector<G4ThreeVector> vertices;
1781 vertices.push_back(fromVertices[ind1]);
1782 vertices.push_back(fromVertices[ind2]);
1783 vertices.push_back(fromVertices[ind3]);
1784
1785 // first vertex most left
1786 //
1787 G4ThreeVector cross=(vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
1788
1789 if ( cross.z() > 0.0 )
1790 {
1791 // Should not happen, as vertices should have been reordered at this stage
1792
1793 std::ostringstream message;
1794 message << "Vertices in wrong order - " << GetName();
1795 G4Exception("G4GenericTrap::MakeDownFacet", "GeomSolids0002",
1796 FatalException, message);
1797 }
1798
1799 return new G4TriangularFacet(vertices[0], vertices[1], vertices[2], ABSOLUTE);
1800}
1801
1802// --------------------------------------------------------------------
1803
1804G4VFacet*
1805G4GenericTrap::MakeUpFacet(const std::vector<G4ThreeVector>& fromVertices,
1806 G4int ind1, G4int ind2, G4int ind3) const
1807{
1808 // Create a triangular facet from the polygon points given by indices
1809 // forming the upper side ( z>0 )
1810
1811 // Do not create facet if 2 vertices are the same
1812 //
1813 if ( (fromVertices[ind1] == fromVertices[ind2]) ||
1814 (fromVertices[ind2] == fromVertices[ind3]) ||
1815 (fromVertices[ind1] == fromVertices[ind3]) ) { return nullptr; }
1816
1817 std::vector<G4ThreeVector> vertices;
1818 vertices.push_back(fromVertices[ind1]);
1819 vertices.push_back(fromVertices[ind2]);
1820 vertices.push_back(fromVertices[ind3]);
1821
1822 // First vertex most left
1823 //
1824 G4ThreeVector cross=(vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
1825
1826 if ( cross.z() < 0.0 )
1827 {
1828 // Should not happen, as vertices should have been reordered at this stage
1829
1830 std::ostringstream message;
1831 message << "Vertices in wrong order - " << GetName();
1832 G4Exception("G4GenericTrap::MakeUpFacet", "GeomSolids0002",
1833 FatalException, message);
1834 }
1835
1836 return new G4TriangularFacet(vertices[0], vertices[1], vertices[2], ABSOLUTE);
1837}
1838
1839// --------------------------------------------------------------------
1840
1841G4VFacet*
1842G4GenericTrap::MakeSideFacet(const G4ThreeVector& downVertex0,
1843 const G4ThreeVector& downVertex1,
1844 const G4ThreeVector& upVertex1,
1845 const G4ThreeVector& upVertex0) const
1846{
1847 // Creates a triangular facet from the polygon points given by indices
1848 // forming the upper side ( z>0 )
1849
1850 if ( (downVertex0 == downVertex1) && (upVertex0 == upVertex1) )
1851 {
1852 return nullptr;
1853 }
1854
1855 if ( downVertex0 == downVertex1 )
1856 {
1857 return new G4TriangularFacet(downVertex0, upVertex1, upVertex0, ABSOLUTE);
1858 }
1859
1860 if ( upVertex0 == upVertex1 )
1861 {
1862 return new G4TriangularFacet(downVertex0, downVertex1, upVertex0, ABSOLUTE);
1863 }
1864
1865 return new G4QuadrangularFacet(downVertex0, downVertex1,
1866 upVertex1, upVertex0, ABSOLUTE);
1867}
1868
1869// --------------------------------------------------------------------
1870
1871G4TessellatedSolid* G4GenericTrap::CreateTessellatedSolid() const
1872{
1873 // 3D vertices
1874 //
1875 G4int nv = fgkNofVertices/2;
1876 std::vector<G4ThreeVector> downVertices;
1877 for ( G4int i=0; i<nv; ++i )
1878 {
1879 downVertices.emplace_back(fVertices[i].x(), fVertices[i].y(), -fDz);
1880 }
1881
1882 std::vector<G4ThreeVector> upVertices;
1883 for ( G4int i=nv; i<2*nv; ++i )
1884 {
1885 upVertices.emplace_back(fVertices[i].x(), fVertices[i].y(), fDz);
1886 }
1887
1888 // Reorder vertices if they are not ordered anti-clock wise
1889 //
1890 G4ThreeVector cross
1891 = (downVertices[1]-downVertices[0]).cross(downVertices[2]-downVertices[1]);
1892 G4ThreeVector cross1
1893 = (upVertices[1]-upVertices[0]).cross(upVertices[2]-upVertices[1]);
1894 if ( (cross.z() > 0.0) || (cross1.z() > 0.0) )
1895 {
1896 ReorderVertices(downVertices);
1897 ReorderVertices(upVertices);
1898 }
1899
1900 auto tessellatedSolid = new G4TessellatedSolid(GetName());
1901
1902 G4VFacet* facet = nullptr;
1903 facet = MakeDownFacet(downVertices, 0, 1, 2);
1904 if (facet != nullptr) { tessellatedSolid->AddFacet( facet ); }
1905 facet = MakeDownFacet(downVertices, 0, 2, 3);
1906 if (facet != nullptr) { tessellatedSolid->AddFacet( facet ); }
1907 facet = MakeUpFacet(upVertices, 0, 2, 1);
1908 if (facet != nullptr) { tessellatedSolid->AddFacet( facet ); }
1909 facet = MakeUpFacet(upVertices, 0, 3, 2);
1910 if (facet != nullptr) { tessellatedSolid->AddFacet( facet ); }
1911
1912 // The quadrangular sides
1913 //
1914 for ( G4int i = 0; i < nv; ++i )
1915 {
1916 G4int j = (i+1) % nv;
1917 facet = MakeSideFacet(downVertices[j], downVertices[i],
1918 upVertices[i], upVertices[j]);
1919
1920 if ( facet != nullptr ) { tessellatedSolid->AddFacet( facet ); }
1921 }
1922
1923 tessellatedSolid->SetSolidClosed(true);
1924
1925 return tessellatedSolid;
1926}
1927
1928// --------------------------------------------------------------------
1929
1930void G4GenericTrap::ComputeBBox()
1931{
1932 // Computes bounding box for a shape.
1933
1934 G4double minX, maxX, minY, maxY;
1935 minX = maxX = fVertices[0].x();
1936 minY = maxY = fVertices[0].y();
1937
1938 for (G4int i=1; i< fgkNofVertices; i++)
1939 {
1940 if (minX>fVertices[i].x()) { minX=fVertices[i].x(); }
1941 if (maxX<fVertices[i].x()) { maxX=fVertices[i].x(); }
1942 if (minY>fVertices[i].y()) { minY=fVertices[i].y(); }
1943 if (maxY<fVertices[i].y()) { maxY=fVertices[i].y(); }
1944 }
1945 fMinBBoxVector = G4ThreeVector(minX,minY,-fDz);
1946 fMaxBBoxVector = G4ThreeVector(maxX,maxY, fDz);
1947}
1948
1949// --------------------------------------------------------------------
1950
1952{
1953
1954#ifdef G4TESS_TEST
1955 if ( fTessellatedSolid )
1956 {
1957 return fTessellatedSolid->GetPolyhedron();
1958 }
1959#endif
1960
1961 if ( (fpPolyhedron == nullptr)
1965 {
1966 G4AutoLock l(&polyhedronMutex);
1967 delete fpPolyhedron;
1969 fRebuildPolyhedron = false;
1970 l.unlock();
1971 }
1972 return fpPolyhedron;
1973}
1974
1975// --------------------------------------------------------------------
1976
1978{
1979
1980#ifdef G4TESS_TEST
1981 if ( fTessellatedSolid )
1982 {
1983 return fTessellatedSolid->DescribeYourselfTo(scene);
1984 }
1985#endif
1986
1987 scene.AddSolid(*this);
1988}
1989
1990// --------------------------------------------------------------------
1991
1993{
1994 // Computes bounding vectors for the shape
1995
1996#ifdef G4TESS_TEST
1997 if ( fTessellatedSolid != nullptr )
1998 {
1999 return fTessellatedSolid->GetExtent();
2000 }
2001#endif
2002
2003 G4ThreeVector minVec = GetMinimumBBox();
2004 G4ThreeVector maxVec = GetMaximumBBox();
2005 return { minVec.x(), maxVec.x(),
2006 minVec.y(), maxVec.y(),
2007 minVec.z(), maxVec.z() };
2008}
2009
2010// --------------------------------------------------------------------
2011
2013{
2014
2015#ifdef G4TESS_TEST
2016 if ( fTessellatedSolid != nullptr )
2017 {
2018 return fTessellatedSolid->CreatePolyhedron();
2019 }
2020#endif
2021
2022 // Approximation of Twisted Side
2023 // Construct extra Points, if Twisted Side
2024 //
2025 G4Polyhedron* polyhedron;
2026 G4int nVertices, nFacets;
2027
2028 G4int subdivisions = 0;
2029 if (fIsTwisted)
2030 {
2031 if (GetVisSubdivisions() != 0)
2032 {
2033 subdivisions = GetVisSubdivisions();
2034 }
2035 else
2036 {
2037 // Estimation of Number of Subdivisions for smooth visualisation
2038 //
2039 G4double maxTwist = 0.;
2040 for(G4int i = 0; i < 4; ++i)
2041 {
2042 if (GetTwistAngle(i) > maxTwist) { maxTwist = GetTwistAngle(i); }
2043 }
2044
2045 // Computes bounding vectors for the shape
2046 //
2047 G4double Dx, Dy;
2048 G4ThreeVector minVec = GetMinimumBBox();
2049 G4ThreeVector maxVec = GetMaximumBBox();
2050 Dx = 0.5*(maxVec.x() - minVec.y());
2051 Dy = 0.5*(maxVec.y() - minVec.y());
2052 if (Dy > Dx) { Dx = Dy; }
2053
2054 subdivisions = 8*G4int(maxTwist/(Dx*Dx*Dx)*fDz);
2055 if (subdivisions < 4) { subdivisions = 4; }
2056 if (subdivisions > 30) { subdivisions = 30; }
2057 }
2058 }
2059 G4int sub4 = 4*subdivisions;
2060 nVertices = 8 + subdivisions*4;
2061 nFacets = 6 + subdivisions*4;
2062 G4double cf = 1./(subdivisions + 1);
2063 polyhedron = new G4Polyhedron(nVertices, nFacets);
2064
2065 // Set vertices
2066 //
2067 G4int icur = 0;
2068 for (G4int i = 0; i < 4; ++i)
2069 {
2070 G4ThreeVector v(fVertices[i].x(),fVertices[i].y(),-fDz);
2071 polyhedron->SetVertex(++icur, v);
2072 }
2073 for (G4int i = 0; i < subdivisions; ++i)
2074 {
2075 for (G4int j = 0; j < 4; ++j)
2076 {
2077 G4TwoVector u = fVertices[j]+cf*(i+1)*(fVertices[j+4]-fVertices[j]);
2078 G4ThreeVector v(u.x(),u.y(),-fDz+cf*2*fDz*(i+1));
2079 polyhedron->SetVertex(++icur, v);
2080 }
2081 }
2082 for (G4int i = 4; i < 8; ++i)
2083 {
2084 G4ThreeVector v(fVertices[i].x(),fVertices[i].y(),fDz);
2085 polyhedron->SetVertex(++icur, v);
2086 }
2087
2088 // Set facets
2089 //
2090 icur = 0;
2091 polyhedron->SetFacet(++icur, 1, 4, 3, 2); // Z-plane
2092 for (G4int i = 0; i < subdivisions + 1; ++i)
2093 {
2094 G4int is = i*4;
2095 polyhedron->SetFacet(++icur, 5+is, 8+is, 4+is, 1+is);
2096 polyhedron->SetFacet(++icur, 8+is, 7+is, 3+is, 4+is);
2097 polyhedron->SetFacet(++icur, 7+is, 6+is, 2+is, 3+is);
2098 polyhedron->SetFacet(++icur, 6+is, 5+is, 1+is, 2+is);
2099 }
2100 polyhedron->SetFacet(++icur, 5+sub4, 6+sub4, 7+sub4, 8+sub4); // Z-plane
2101
2102 polyhedron->SetReferences();
2103 polyhedron->InvertFacets();
2104
2105 return polyhedron;
2106}
2107
2108// --------------------------------------------------------------------
2109
2110#endif
std::vector< G4ThreeVector > G4ThreeVectorList
G4double C(G4double temp)
G4double B(G4double temperature)
G4double D(G4double temp)
@ JustWarning
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
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
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
#define G4UniformRand()
Definition Randomize.hh:52
double x() const
double y() const
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() const
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
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
G4GenericTrap & operator=(const G4GenericTrap &rhs)
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
G4double GetZHalfLength() const
std::ostream & StreamInfo(std::ostream &os) const override
G4TwoVector GetVertex(G4int index) const
G4GeometryType GetEntityType() const override
G4Polyhedron * GetPolyhedron() const override
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
G4bool fRebuildPolyhedron
G4double GetTwistAngle(G4int index) const
G4double GetCubicVolume() override
G4VisExtent GetExtent() const override
G4VSolid * Clone() const override
G4int GetVisSubdivisions() const
~G4GenericTrap() override
EInside Inside(const G4ThreeVector &p) const override
G4Polyhedron * CreatePolyhedron() const override
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
G4GenericTrap(const G4String &name, G4double halfZ, const std::vector< G4TwoVector > &vertices)
G4double GetSurfaceArea() override
G4ThreeVector GetPointOnSurface() const override
G4Polyhedron * fpPolyhedron
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
EInside Inside(const G4ThreeVector &p) const override
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
G4Polyhedron * GetPolyhedron() const override
G4double DistanceToOut(const G4ThreeVector &p) const override
G4ThreeVector GetPointOnSurface() const override
G4VisExtent GetExtent() const override
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
G4Polyhedron * CreatePolyhedron() const override
virtual void AddSolid(const G4Box &)=0
G4String GetName() const
void DumpInfo() const
G4double kCarTolerance
Definition G4VSolid.hh:299
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:107
static G4int GetNumberOfRotationSteps()
void SetVertex(G4int index, const G4Point3D &v)
void SetFacet(G4int index, G4int iv1, G4int iv2, G4int iv3, G4int iv4=0)
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