Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VTwistedFaceted.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// G4VTwistedFaceted implementation
27//
28// Author: Oliver Link (CERN), 27.10.2004 - Created
29// --------------------------------------------------------------------
30
31#include "G4VTwistedFaceted.hh"
32
34#include "G4SystemOfUnits.hh"
35#include "G4VoxelLimits.hh"
36#include "G4AffineTransform.hh"
37#include "G4BoundingEnvelope.hh"
38#include "G4SolidExtentList.hh"
39#include "G4ClippablePolygon.hh"
42
43#include "G4VGraphicsScene.hh"
44#include "G4Polyhedron.hh"
45#include "G4VisExtent.hh"
46
47#include "G4QuickRand.hh"
48
49#include "G4AutoLock.hh"
50
51namespace
52{
53 G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
55}
56
57
58//=====================================================================
59//* constructors ------------------------------------------------------
60
62G4VTwistedFaceted( const G4String& pname, // Name of instance
63 G4double PhiTwist, // twist angle
64 G4double pDz, // half z length
65 G4double pTheta, // direction between end planes
66 G4double pPhi, // defined by polar and azim. angles
67 G4double pDy1, // half y length at -pDz
68 G4double pDx1, // half x length at -pDz,-pDy
69 G4double pDx2, // half x length at -pDz,+pDy
70 G4double pDy2, // half y length at +pDz
71 G4double pDx3, // half x length at +pDz,-pDy
72 G4double pDx4, // half x length at +pDz,+pDy
73 G4double pAlph ) // tilt angle
74 : G4VSolid(pname),
75 fLowerEndcap(nullptr), fUpperEndcap(nullptr), fSide0(nullptr),
76 fSide90(nullptr), fSide180(nullptr), fSide270(nullptr)
77{
78
79 G4double pDytmp ;
80 G4double fDxUp ;
81 G4double fDxDown ;
82
83 fDx1 = pDx1 ;
84 fDx2 = pDx2 ;
85 fDx3 = pDx3 ;
86 fDx4 = pDx4 ;
87 fDy1 = pDy1 ;
88 fDy2 = pDy2 ;
89 fDz = pDz ;
90
91 G4double kAngTolerance
93
94 // maximum values
95 //
96 fDxDown = ( fDx1 > fDx2 ? fDx1 : fDx2 ) ;
97 fDxUp = ( fDx3 > fDx4 ? fDx3 : fDx4 ) ;
98 fDx = ( fDxUp > fDxDown ? fDxUp : fDxDown ) ;
99 fDy = ( fDy1 > fDy2 ? fDy1 : fDy2 ) ;
100
101 // planarity check
102 //
103 if ( fDx1 != fDx2 && fDx3 != fDx4 )
104 {
105 pDytmp = fDy1 * ( fDx3 - fDx4 ) / ( fDx1 - fDx2 ) ;
106 if ( std::fabs(pDytmp - fDy2) > kCarTolerance )
107 {
108 std::ostringstream message;
109 message << "Not planar surface in untwisted Trapezoid: "
110 << GetName() << G4endl
111 << "fDy2 is " << fDy2 << " but should be "
112 << pDytmp << ".";
113 G4Exception("G4VTwistedFaceted::G4VTwistedFaceted()", "GeomSolids0002",
114 FatalErrorInArgument, message);
115 }
116 }
117
118#ifdef G4TWISTDEBUG
119 if ( fDx1 == fDx2 && fDx3 == fDx4 )
120 {
121 G4cout << "Trapezoid is a box" << G4endl ;
122 }
123
124#endif
125
126 if ( ( fDx1 == fDx2 && fDx3 != fDx4 ) || ( fDx1 != fDx2 && fDx3 == fDx4 ) )
127 {
128 std::ostringstream message;
129 message << "Not planar surface in untwisted Trapezoid: "
130 << GetName() << G4endl
131 << "One endcap is rectangular, the other is a trapezoid." << G4endl
132 << "For planarity reasons they have to be rectangles or trapezoids "
133 << "on both sides.";
134 G4Exception("G4VTwistedFaceted::G4VTwistedFaceted()", "GeomSolids0002",
135 FatalErrorInArgument, message);
136 }
137
138 // twist angle
139 //
140 fPhiTwist = PhiTwist ;
141
142 // tilt angle
143 //
144 fAlph = pAlph ;
145 fTAlph = std::tan(fAlph) ;
146
147 fTheta = pTheta ;
148 fPhi = pPhi ;
149
150 // dx in surface equation
151 //
152 fdeltaX = 2 * fDz * std::tan(fTheta) * std::cos(fPhi) ;
153
154 // dy in surface equation
155 //
156 fdeltaY = 2 * fDz * std::tan(fTheta) * std::sin(fPhi) ;
157
158 if ( ( fDx1 <= 2*kCarTolerance)
159 || ( fDx2 <= 2*kCarTolerance)
160 || ( fDx3 <= 2*kCarTolerance)
161 || ( fDx4 <= 2*kCarTolerance)
162 || ( fDy1 <= 2*kCarTolerance)
163 || ( fDy2 <= 2*kCarTolerance)
164 || ( fDz <= 2*kCarTolerance)
165 || ( std::fabs(fPhiTwist) <= 2*kAngTolerance )
166 || ( std::fabs(fPhiTwist) >= pi/2 )
167 || ( std::fabs(fAlph) >= pi/2 )
168 || fTheta >= pi/2 || fTheta < 0
169 )
170 {
171 std::ostringstream message;
172 message << "Invalid dimensions. Too small, or twist angle too big: "
173 << GetName() << G4endl
174 << "fDx 1-4 = " << fDx1/cm << ", " << fDx2/cm << ", "
175 << fDx3/cm << ", " << fDx4/cm << " cm" << G4endl
176 << "fDy 1-2 = " << fDy1/cm << ", " << fDy2/cm << ", "
177 << " cm" << G4endl
178 << "fDz = " << fDz/cm << " cm" << G4endl
179 << " twistangle " << fPhiTwist/deg << " deg" << G4endl
180 << " phi,theta = " << fPhi/deg << ", " << fTheta/deg << " deg";
181 G4Exception("G4TwistedTrap::G4VTwistedFaceted()",
182 "GeomSolids0002", FatalErrorInArgument, message);
183 }
184 CreateSurfaces();
185}
186
187
188//=====================================================================
189//* Fake default constructor ------------------------------------------
190
192 : G4VSolid(a),
193 fLowerEndcap(nullptr), fUpperEndcap(nullptr),
194 fSide0(nullptr), fSide90(nullptr), fSide180(nullptr), fSide270(nullptr)
195{
196}
197
198
199//=====================================================================
200//* destructor --------------------------------------------------------
201
203{
204 delete fLowerEndcap ;
205 delete fUpperEndcap ;
206
207 delete fSide0 ;
208 delete fSide90 ;
209 delete fSide180 ;
210 delete fSide270 ;
211 delete fpPolyhedron; fpPolyhedron = nullptr;
212}
213
214
215//=====================================================================
216//* Copy constructor --------------------------------------------------
217
219 : G4VSolid(rhs),
221 fTheta(rhs.fTheta), fPhi(rhs.fPhi),
222 fDy1(rhs.fDy1), fDx1(rhs.fDx1), fDx2(rhs.fDx2), fDy2(rhs.fDy2),
223 fDx3(rhs.fDx3), fDx4(rhs.fDx4), fDz(rhs.fDz), fDx(rhs.fDx), fDy(rhs.fDy),
224 fAlph(rhs.fAlph), fTAlph(rhs.fTAlph), fdeltaX(rhs.fdeltaX),
225 fdeltaY(rhs.fdeltaY), fPhiTwist(rhs.fPhiTwist), fLowerEndcap(nullptr),
226 fUpperEndcap(nullptr), fSide0(nullptr), fSide90(nullptr), fSide180(nullptr), fSide270(nullptr)
227{
228 CreateSurfaces();
229}
230
231
232//=====================================================================
233//* Assignment operator -----------------------------------------------
234
236{
237 // Check assignment to self
238 //
239 if (this == &rhs) { return *this; }
240
241 // Copy base class data
242 //
244
245 // Copy data
246 //
247 fTheta = rhs.fTheta; fPhi = rhs.fPhi;
248 fDy1= rhs.fDy1; fDx1= rhs.fDx1; fDx2= rhs.fDx2; fDy2= rhs.fDy2;
249 fDx3= rhs.fDx3; fDx4= rhs.fDx4; fDz= rhs.fDz; fDx= rhs.fDx; fDy= rhs.fDy;
250 fAlph= rhs.fAlph; fTAlph= rhs.fTAlph; fdeltaX= rhs.fdeltaX;
251 fdeltaY= rhs.fdeltaY; fPhiTwist= rhs.fPhiTwist; fLowerEndcap= nullptr;
252 fUpperEndcap= nullptr; fSide0= nullptr; fSide90= nullptr; fSide180= nullptr; fSide270= nullptr;
254 fRebuildPolyhedron = false;
255 delete fpPolyhedron; fpPolyhedron = nullptr;
256
257 CreateSurfaces();
258
259 return *this;
260}
261
262
263//=====================================================================
264//* ComputeDimensions -------------------------------------------------
265
267 const G4int ,
268 const G4VPhysicalVolume* )
269{
270 G4Exception("G4VTwistedFaceted::ComputeDimensions()",
271 "GeomSolids0001", FatalException,
272 "G4VTwistedFaceted does not support Parameterisation.");
273}
274
275
276//=====================================================================
277//* Extent ------------------------------------------------------------
278
280 G4ThreeVector& pMax) const
281{
282 G4double cosPhi = std::cos(fPhi);
283 G4double sinPhi = std::sin(fPhi);
284 G4double tanTheta = std::tan(fTheta);
285 G4double tanAlpha = fTAlph;
286
287 G4double xmid1 = fDy1*tanAlpha;
288 G4double x1 = std::abs(xmid1 + fDx1);
289 G4double x2 = std::abs(xmid1 - fDx1);
290 G4double x3 = std::abs(xmid1 + fDx2);
291 G4double x4 = std::abs(xmid1 - fDx2);
292 G4double xmax1 = std::max(std::max(std::max(x1, x2), x3), x4);
293 G4double rmax1 = std::sqrt(xmax1*xmax1 + fDy1*fDy1);
294
295 G4double xmid2 = fDy2*tanAlpha;
296 G4double x5 = std::abs(xmid2 + fDx3);
297 G4double x6 = std::abs(xmid2 - fDx3);
298 G4double x7 = std::abs(xmid2 + fDx4);
299 G4double x8 = std::abs(xmid2 - fDx4);
300 G4double xmax2 = std::max(std::max(std::max(x5, x6), x7), x8);
301 G4double rmax2 = std::sqrt(xmax2*xmax2 + fDy2*fDy2);
302
303 G4double x0 = fDz*tanTheta*cosPhi;
304 G4double y0 = fDz*tanTheta*sinPhi;
305 G4double xmin = std::min(-x0 - rmax1, x0 - rmax2);
306 G4double ymin = std::min(-y0 - rmax1, y0 - rmax2);
307 G4double xmax = std::max(-x0 + rmax1, x0 + rmax2);
308 G4double ymax = std::max(-y0 + rmax1, y0 + rmax2);
309 pMin.set(xmin, ymin,-fDz);
310 pMax.set(xmax, ymax, fDz);
311}
312
313
314//=====================================================================
315//* CalculateExtent ---------------------------------------------------
316
317G4bool
319 const G4VoxelLimits& pVoxelLimit,
320 const G4AffineTransform& pTransform,
321 G4double& pMin,
322 G4double& pMax ) const
323{
324 G4ThreeVector bmin, bmax;
325
326 // Get bounding box
327 BoundingLimits(bmin,bmax);
328
329 // Find extent
330 G4BoundingEnvelope bbox(bmin,bmax);
331 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
332}
333
334
335//=====================================================================
336//* Inside ------------------------------------------------------------
337
339{
340
341 EInside tmpin = kOutside ;
342
343 G4double phi = p.z()/(2*fDz) * fPhiTwist ; // rotate the point to z=0
344 G4double cphi = std::cos(-phi) ;
345 G4double sphi = std::sin(-phi) ;
346
347 G4double px = p.x() + fdeltaX * ( -phi/fPhiTwist) ; // shift
348 G4double py = p.y() + fdeltaY * ( -phi/fPhiTwist) ;
349 G4double pz = p.z() ;
350
351 G4double posx = px * cphi - py * sphi ; // rotation
352 G4double posy = px * sphi + py * cphi ;
353 G4double posz = pz ;
354
355 G4double xMin = Xcoef(posy,phi,fTAlph) - 2*Xcoef(posy,phi,0.) ;
356 G4double xMax = Xcoef(posy,phi,fTAlph) ;
357
358 G4double yMax = GetValueB(phi)/2. ; // b(phi)/2 is limit
359 G4double yMin = -yMax ;
360
361#ifdef G4TWISTDEBUG
362
363 G4cout << "inside called: p = " << p << G4endl ;
364 G4cout << "fDx1 = " << fDx1 << G4endl ;
365 G4cout << "fDx2 = " << fDx2 << G4endl ;
366 G4cout << "fDx3 = " << fDx3 << G4endl ;
367 G4cout << "fDx4 = " << fDx4 << G4endl ;
368
369 G4cout << "fDy1 = " << fDy1 << G4endl ;
370 G4cout << "fDy2 = " << fDy2 << G4endl ;
371
372 G4cout << "fDz = " << fDz << G4endl ;
373
374 G4cout << "Tilt angle alpha = " << fAlph << G4endl ;
375 G4cout << "phi,theta = " << fPhi << " , " << fTheta << G4endl ;
376
377 G4cout << "Twist angle = " << fPhiTwist << G4endl ;
378
379 G4cout << "posx = " << posx << G4endl ;
380 G4cout << "posy = " << posy << G4endl ;
381 G4cout << "xMin = " << xMin << G4endl ;
382 G4cout << "xMax = " << xMax << G4endl ;
383 G4cout << "yMin = " << yMin << G4endl ;
384 G4cout << "yMax = " << yMax << G4endl ;
385
386#endif
387
388
389 if ( posx <= xMax - kCarTolerance*0.5
390 && posx >= xMin + kCarTolerance*0.5 )
391 {
392 if ( posy <= yMax - kCarTolerance*0.5
393 && posy >= yMin + kCarTolerance*0.5 )
394 {
395 if (std::fabs(posz) <= fDz - kCarTolerance*0.5 )
396 {
397 tmpin = kInside ;
398 }
399 else if (std::fabs(posz) <= fDz + kCarTolerance*0.5 )
400 {
401 tmpin = kSurface ;
402 }
403 }
404 else if ( posy <= yMax + kCarTolerance*0.5
405 && posy >= yMin - kCarTolerance*0.5 )
406 {
407 if (std::fabs(posz) <= fDz + kCarTolerance*0.5 ) { tmpin = kSurface ; }
408 }
409 }
410 else if ( posx <= xMax + kCarTolerance*0.5
411 && posx >= xMin - kCarTolerance*0.5 )
412 {
413 if ( posy <= yMax + kCarTolerance*0.5
414 && posy >= yMin - kCarTolerance*0.5 )
415 {
416 if (std::fabs(posz) <= fDz + kCarTolerance*0.5) { tmpin = kSurface ; }
417 }
418 }
419
420#ifdef G4TWISTDEBUG
421 G4cout << "inside = " << tmpin << G4endl ;
422#endif
423
424 return tmpin;
425
426}
427
428
429//=====================================================================
430//* SurfaceNormal -----------------------------------------------------
431
433{
434 //
435 // return the normal unit vector to the Hyperbolical Surface at a point
436 // p on (or nearly on) the surface
437 //
438 // Which of the three or four surfaces are we closest to?
439 //
440
441
442 G4double distance = kInfinity;
443
444 G4VTwistSurface* surfaces[6];
445
446 surfaces[0] = fSide0 ;
447 surfaces[1] = fSide90 ;
448 surfaces[2] = fSide180 ;
449 surfaces[3] = fSide270 ;
450 surfaces[4] = fLowerEndcap;
451 surfaces[5] = fUpperEndcap;
452
453 G4ThreeVector xx;
454 G4ThreeVector bestxx;
455 G4int i;
456 G4int besti = -1;
457 for (i=0; i< 6; i++)
458 {
459 G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
460 if (tmpdistance < distance)
461 {
462 distance = tmpdistance;
463 bestxx = xx;
464 besti = i;
465 }
466 }
467
468 return surfaces[besti]->GetNormal(bestxx, true);
469}
470
471
472//=====================================================================
473//* DistanceToIn (p, v) -----------------------------------------------
474
476 const G4ThreeVector& v ) const
477{
478
479 // DistanceToIn (p, v):
480 // Calculate distance to surface of shape from `outside'
481 // along with the v, allowing for tolerance.
482 // The function returns kInfinity if no intersection or
483 // just grazing within tolerance.
484
485 //
486 // Calculate DistanceToIn(p,v)
487 //
488
489 EInside currentside = Inside(p);
490
491 if (currentside == kInside)
492 {
493 }
494 else if (currentside == kSurface)
495 {
496 // particle is just on a boundary.
497 // if the particle is entering to the volume, return 0
498 //
499 G4ThreeVector normal = SurfaceNormal(p);
500 if (normal*v < 0)
501 {
502 return 0;
503 }
504 }
505
506 // now, we can take smallest positive distance.
507
508 // Initialize
509 //
510 G4double distance = kInfinity;
511
512 // Find intersections and choose nearest one
513 //
514 G4VTwistSurface *surfaces[6];
515
516 surfaces[0] = fSide0;
517 surfaces[1] = fSide90 ;
518 surfaces[2] = fSide180 ;
519 surfaces[3] = fSide270 ;
520 surfaces[4] = fLowerEndcap;
521 surfaces[5] = fUpperEndcap;
522
523 G4ThreeVector xx;
524 G4ThreeVector bestxx;
525 for (const auto & surface : surfaces)
526 {
527#ifdef G4TWISTDEBUG
528 G4cout << G4endl << "surface " << &surface - &*surfaces << ": " << G4endl << G4endl ;
529#endif
530 G4double tmpdistance = surface->DistanceToIn(p, v, xx);
531#ifdef G4TWISTDEBUG
532 G4cout << "Solid DistanceToIn : distance = " << tmpdistance << G4endl ;
533 G4cout << "intersection point = " << xx << G4endl ;
534#endif
535 if (tmpdistance < distance)
536 {
537 distance = tmpdistance;
538 bestxx = xx;
539 }
540 }
541
542#ifdef G4TWISTDEBUG
543 G4cout << "best distance = " << distance << G4endl ;
544#endif
545
546 // timer.Stop();
547 return distance;
548}
549
550
551//=====================================================================
552//* DistanceToIn (p) --------------------------------------------------
553
555{
556 // DistanceToIn(p):
557 // Calculate distance to surface of shape from `outside',
558 // allowing for tolerance
559 //
560
561 //
562 // Calculate DistanceToIn(p)
563 //
564
565 EInside currentside = Inside(p);
566
567 switch (currentside)
568 {
569 case (kInside) :
570 {
571 }
572
573 case (kSurface) :
574 {
575 return 0;
576 }
577
578 case (kOutside) :
579 {
580 // Initialize
581 //
582 G4double distance = kInfinity;
583
584 // Find intersections and choose nearest one
585 //
586 G4VTwistSurface* surfaces[6];
587
588 surfaces[0] = fSide0;
589 surfaces[1] = fSide90 ;
590 surfaces[2] = fSide180 ;
591 surfaces[3] = fSide270 ;
592 surfaces[4] = fLowerEndcap;
593 surfaces[5] = fUpperEndcap;
594
595 G4ThreeVector xx;
596 G4ThreeVector bestxx;
597 for (const auto & surface : surfaces)
598 {
599 G4double tmpdistance = surface->DistanceTo(p, xx);
600 if (tmpdistance < distance)
601 {
602 distance = tmpdistance;
603 bestxx = xx;
604 }
605 }
606 return distance;
607 }
608
609 default:
610 {
611 G4Exception("G4VTwistedFaceted::DistanceToIn(p)", "GeomSolids0003",
612 FatalException, "Unknown point location!");
613 }
614 } // switch end
615
616 return 0.;
617}
618
619
620//=====================================================================
621//* DistanceToOut (p, v) ----------------------------------------------
622
625 const G4ThreeVector& v,
626 const G4bool calcNorm,
627 G4bool* validNorm,
628 G4ThreeVector* norm ) const
629{
630 // DistanceToOut (p, v):
631 // Calculate distance to surface of shape from `inside'
632 // along with the v, allowing for tolerance.
633 // The function returns kInfinity if no intersection or
634 // just grazing within tolerance.
635
636 //
637 // Calculate DistanceToOut(p,v)
638 //
639
640 EInside currentside = Inside(p);
641
642 if (currentside == kOutside)
643 {
644 }
645 else if (currentside == kSurface)
646 {
647 // particle is just on a boundary.
648 // if the particle is exiting from the volume, return 0
649 //
650 G4ThreeVector normal = SurfaceNormal(p);
651 if (normal*v > 0)
652 {
653 if (calcNorm)
654 {
655 *norm = normal;
656 *validNorm = true;
657 }
658 // timer.Stop();
659 return 0;
660 }
661 }
662
663 // now, we can take smallest positive distance.
664
665 // Initialize
666 G4double distance = kInfinity;
667
668 // find intersections and choose nearest one.
669 G4VTwistSurface *surfaces[6];
670
671 surfaces[0] = fSide0;
672 surfaces[1] = fSide90 ;
673 surfaces[2] = fSide180 ;
674 surfaces[3] = fSide270 ;
675 surfaces[4] = fLowerEndcap;
676 surfaces[5] = fUpperEndcap;
677
678 G4int besti = -1;
679 G4ThreeVector xx;
680 G4ThreeVector bestxx;
681 for (auto i=0; i<6 ; ++i)
682 {
683 G4double tmpdistance = surfaces[i]->DistanceToOut(p, v, xx);
684 if (tmpdistance < distance)
685 {
686 distance = tmpdistance;
687 bestxx = xx;
688 besti = i;
689 }
690 }
691
692 if (calcNorm)
693 {
694 if (besti != -1)
695 {
696 *norm = (surfaces[besti]->GetNormal(p, true));
697 *validNorm = surfaces[besti]->IsValidNorm();
698 }
699 }
700
701 return distance;
702}
703
704
705//=====================================================================
706//* DistanceToOut (p) -------------------------------------------------
707
709{
710 // DistanceToOut(p):
711 // Calculate distance to surface of shape from `inside',
712 // allowing for tolerance
713
714 //
715 // Calculate DistanceToOut(p)
716 //
717
718 EInside currentside = Inside(p);
719 G4double retval = kInfinity;
720
721 switch (currentside)
722 {
723 case (kOutside) :
724 {
725#ifdef G4SPECSDEBUG
726 G4int oldprc = G4cout.precision(16) ;
727 G4cout << G4endl ;
728 DumpInfo();
729 G4cout << "Position:" << G4endl << G4endl ;
730 G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ;
731 G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ;
732 G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ;
733 G4cout.precision(oldprc) ;
734 G4Exception("G4VTwistedFaceted::DistanceToOut(p)", "GeomSolids1002",
735 JustWarning, "Point p is outside !?" );
736#endif
737 break;
738 }
739 case (kSurface) :
740 {
741 retval = 0;
742 break;
743 }
744
745 case (kInside) :
746 {
747 // Initialize
748 //
749 G4double distance = kInfinity;
750
751 // find intersections and choose nearest one
752 //
753 G4VTwistSurface* surfaces[6];
754
755 surfaces[0] = fSide0;
756 surfaces[1] = fSide90 ;
757 surfaces[2] = fSide180 ;
758 surfaces[3] = fSide270 ;
759 surfaces[4] = fLowerEndcap;
760 surfaces[5] = fUpperEndcap;
761
762 G4ThreeVector xx;
763 G4ThreeVector bestxx;
764 for (const auto & surface : surfaces)
765 {
766 G4double tmpdistance = surface->DistanceTo(p, xx);
767 if (tmpdistance < distance)
768 {
769 distance = tmpdistance;
770 bestxx = xx;
771 }
772 }
773 retval = distance;
774 break;
775 }
776
777 default :
778 {
779 G4Exception("G4VTwistedFaceted::DistanceToOut(p)", "GeomSolids0003",
780 FatalException, "Unknown point location!");
781 break;
782 }
783 } // switch end
784
785 return retval;
786}
787
788
789//=====================================================================
790//* StreamInfo --------------------------------------------------------
791
792std::ostream& G4VTwistedFaceted::StreamInfo(std::ostream& os) const
793{
794 //
795 // Stream object contents to an output stream
796 //
797 G4long oldprc = os.precision(16);
798 os << "-----------------------------------------------------------\n"
799 << " *** Dump for solid - " << GetName() << " ***\n"
800 << " ===================================================\n"
801 << " Solid type: G4VTwistedFaceted\n"
802 << " Parameters: \n"
803 << " polar angle theta = " << fTheta/degree << " deg" << G4endl
804 << " azimuthal angle phi = " << fPhi/degree << " deg" << G4endl
805 << " tilt angle alpha = " << fAlph/degree << " deg" << G4endl
806 << " TWIST angle = " << fPhiTwist/degree << " deg" << G4endl
807 << " Half length along y (lower endcap) = " << fDy1/cm << " cm"
808 << G4endl
809 << " Half length along x (lower endcap, bottom) = " << fDx1/cm << " cm"
810 << G4endl
811 << " Half length along x (lower endcap, top) = " << fDx2/cm << " cm"
812 << G4endl
813 << " Half length along y (upper endcap) = " << fDy2/cm << " cm"
814 << G4endl
815 << " Half length along x (upper endcap, bottom) = " << fDx3/cm << " cm"
816 << G4endl
817 << " Half length along x (upper endcap, top) = " << fDx4/cm << " cm"
818 << G4endl
819 << "-----------------------------------------------------------\n";
820 os.precision(oldprc);
821
822 return os;
823}
824
825
826//=====================================================================
827//* DiscribeYourselfTo ------------------------------------------------
828
830{
831 scene.AddSolid (*this);
832}
833
834
835//=====================================================================
836//* GetExtent ---------------------------------------------------------
837
839{
840 G4double maxRad = std::sqrt( fDx*fDx + fDy*fDy);
841
842 return { -maxRad, maxRad ,
843 -maxRad, maxRad ,
844 -fDz, fDz };
845}
846
847
848//=====================================================================
849//* CreateSurfaces ----------------------------------------------------
850
851void G4VTwistedFaceted::CreateSurfaces()
852{
853
854 // create 6 surfaces of TwistedTub.
855
856 if ( fDx1 == fDx2 && fDx3 == fDx4 ) // special case : Box
857 {
858 fSide0 = new G4TwistBoxSide("0deg", fPhiTwist, fDz, fTheta, fPhi,
859 fDy1, fDx1, fDx1, fDy2, fDx3, fDx3, fAlph, 0.*deg);
860 fSide180 = new G4TwistBoxSide("180deg", fPhiTwist, fDz, fTheta, fPhi+pi,
861 fDy1, fDx1, fDx1, fDy2, fDx3, fDx3, fAlph, 180.*deg);
862 }
863 else // default general case
864 {
865 fSide0 = new G4TwistTrapAlphaSide("0deg" ,fPhiTwist, fDz, fTheta,
866 fPhi, fDy1, fDx1, fDx2, fDy2, fDx3, fDx4, fAlph, 0.*deg);
867 fSide180 = new G4TwistTrapAlphaSide("180deg", fPhiTwist, fDz, fTheta,
868 fPhi+pi, fDy1, fDx2, fDx1, fDy2, fDx4, fDx3, fAlph, 180.*deg);
869 }
870
871 // create parallel sides
872 //
873 fSide90 = new G4TwistTrapParallelSide("90deg", fPhiTwist, fDz, fTheta,
874 fPhi, fDy1, fDx1, fDx2, fDy2, fDx3, fDx4, fAlph, 0.*deg);
875 fSide270 = new G4TwistTrapParallelSide("270deg", fPhiTwist, fDz, fTheta,
876 fPhi+pi, fDy1, fDx2, fDx1, fDy2, fDx4, fDx3, fAlph, 180.*deg);
877
878 // create endcaps
879 //
880 fUpperEndcap = new G4TwistTrapFlatSide("UpperCap",fPhiTwist, fDx3, fDx4, fDy2,
881 fDz, fAlph, fPhi, fTheta, 1 );
882 fLowerEndcap = new G4TwistTrapFlatSide("LowerCap",fPhiTwist, fDx1, fDx2, fDy1,
883 fDz, fAlph, fPhi, fTheta, -1 );
884
885 // Set neighbour surfaces
886
887 fSide0->SetNeighbours( fSide270 , fLowerEndcap , fSide90 , fUpperEndcap );
888 fSide90->SetNeighbours( fSide0 , fLowerEndcap , fSide180 , fUpperEndcap );
889 fSide180->SetNeighbours(fSide90 , fLowerEndcap , fSide270 , fUpperEndcap );
890 fSide270->SetNeighbours(fSide180 , fLowerEndcap , fSide0 , fUpperEndcap );
891 fUpperEndcap->SetNeighbours( fSide180, fSide270 , fSide0 , fSide90 );
892 fLowerEndcap->SetNeighbours( fSide180, fSide270 , fSide0 , fSide90 );
893
894}
895
896//=====================================================================
897//* GetLateralFaceArea ------------------------------------------------
898
900G4VTwistedFaceted::GetLateralFaceArea(const G4TwoVector& p1,
901 const G4TwoVector& p2,
902 const G4TwoVector& p3,
903 const G4TwoVector& p4) const
904{
905 constexpr G4int NSTEP = 100;
906 constexpr G4double dt = 1./NSTEP;
907
908 G4double h = 2.*fDz;
909 G4double hh = h*h;
910 G4double hTanTheta = h*std::tan(fTheta);
911 G4double x1 = p1.x();
912 G4double y1 = p1.y();
913 G4double x21 = p2.x() - p1.x();
914 G4double y21 = p2.y() - p1.y();
915 G4double x31 = p3.x() - p1.x();
916 G4double y31 = p3.y() - p1.y();
917 G4double x42 = p4.x() - p2.x();
918 G4double y42 = p4.y() - p2.y();
919 G4double x43 = p4.x() - p3.x();
920 G4double y43 = p4.y() - p3.y();
921
922 // check if face is plane (just in case)
923 G4double lmax = std::max(std::max(std::abs(x21),std::abs(y21)),
924 std::max(std::abs(x43),std::abs(y43)));
925 G4double eps = lmax*kCarTolerance;
926 if (std::abs(fPhiTwist) < kCarTolerance &&
927 std::abs(x21*y43 - y21*x43) < eps)
928 {
929 G4double x0 = hTanTheta*std::cos(fPhi);
930 G4double y0 = hTanTheta*std::sin(fPhi);
931 G4ThreeVector A(p4.x() - p1.x() + x0, p4.y() - p1.y() + y0, h);
932 G4ThreeVector B(p3.x() - p2.x() + x0, p3.y() - p2.y() + y0, h);
933 return (A.cross(B)).mag()*0.5;
934 }
935
936 // twisted face
937 G4double area = 0.;
938 for (G4int i = 0; i < NSTEP; ++i)
939 {
940 G4double t = (i + 0.5)*dt;
941 G4double I = x21 + (x42 - x31)*t;
942 G4double J = y21 + (y42 - y31)*t;
943 G4double II = I*I;
944 G4double JJ = J*J;
945 G4double IIJJ = hh*(I*I + J*J);
946
947 G4double ang = fPhi + fPhiTwist*(0.5 - t);
948 G4double A = fPhiTwist*(II + JJ) + x21*y43 - x43*y21;
949 G4double B = fPhiTwist*(I*(x1 + x31*t) + J*(y1 + y31*t)) +
950 hTanTheta*(I*std::sin(ang) - J*std::cos(ang)) +
951 (I*y31 - J*x31);
952
953 G4double invAA = 1./(A*A);
954 G4double sqrtAA = 2.*std::abs(A);
955 G4double invSqrtAA = 1./sqrtAA;
956
957 G4double aa = A*A;
958 G4double bb = 2.*A*B;
959 G4double cc = IIJJ + B*B;
960
961 G4double R1 = std::sqrt(aa + bb + cc);
962 G4double R0 = std::sqrt(cc);
963 G4double log1 = std::log(std::abs(sqrtAA*R1 + 2.*aa + bb));
964 G4double log0 = std::log(std::abs(sqrtAA*R0 + bb));
965
966 area += 0.5*R1 + 0.25*bb*invAA*(R1 - R0) + IIJJ*invSqrtAA*(log1 - log0);
967 }
968 return area*dt;
969}
970
971//=====================================================================
972//* GetCubicVolume ----------------------------------------------------
973
975{
976 if(fCubicVolume == 0)
977 {
978 G4AutoLock l(&vtwMutex);
979 fCubicVolume = ((fDx1 + fDx2 + fDx3 + fDx4)*(fDy1 + fDy2) +
980 (fDx4 + fDx3 - fDx2 - fDx1)*(fDy2 - fDy1)/3)*fDz;
981 l.unlock();
982 }
983 return fCubicVolume;
984}
985
986//=====================================================================
987//* GetSurfaceArea ----------------------------------------------------
988
990{
991 if (fSurfaceArea == 0)
992 {
993 G4AutoLock l(&vtwMutex);
994 G4TwoVector vv[8];
995 vv[0] = G4TwoVector(-fDx1 - fDy1*fTAlph,-fDy1);
996 vv[1] = G4TwoVector( fDx1 - fDy1*fTAlph,-fDy1);
997 vv[2] = G4TwoVector(-fDx2 + fDy1*fTAlph, fDy1);
998 vv[3] = G4TwoVector( fDx2 + fDy1*fTAlph, fDy1);
999 vv[4] = G4TwoVector(-fDx3 - fDy2*fTAlph,-fDy2);
1000 vv[5] = G4TwoVector( fDx3 - fDy2*fTAlph,-fDy2);
1001 vv[6] = G4TwoVector(-fDx4 + fDy2*fTAlph, fDy2);
1002 vv[7] = G4TwoVector( fDx4 + fDy2*fTAlph, fDy2);
1003 fSurfaceArea = 2.*(fDy1*(fDx1 + fDx2) + fDy2*(fDx3 + fDx4)) +
1004 GetLateralFaceArea(vv[0], vv[1], vv[4], vv[5]) +
1005 GetLateralFaceArea(vv[1], vv[3], vv[5], vv[7]) +
1006 GetLateralFaceArea(vv[3], vv[2], vv[7], vv[6]) +
1007 GetLateralFaceArea(vv[2], vv[0], vv[6], vv[4]);
1008 l.unlock();
1009 }
1010 return fSurfaceArea;
1011}
1012
1013//=====================================================================
1014//* GetEntityType -----------------------------------------------------
1015
1017{
1018 return {"G4VTwistedFaceted"};
1019}
1020
1021
1022//=====================================================================
1023//* GetPolyhedron -----------------------------------------------------
1024
1026{
1027 if (fpPolyhedron == nullptr ||
1029 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
1030 fpPolyhedron->GetNumberOfRotationSteps())
1031 {
1032 G4AutoLock l(&polyhedronMutex);
1033 delete fpPolyhedron;
1035 fRebuildPolyhedron = false;
1036 l.unlock();
1037 }
1038
1039 return fpPolyhedron;
1040}
1041
1042
1043//=====================================================================
1044//* GetPointOnSurface -------------------------------------------------
1045
1047{
1048
1049 G4double phi = fPhiTwist*(G4QuickRand() - 0.5);
1050 G4double u , umin, umax ; // variable for twisted surfaces
1051 G4double y ; // variable for flat surface (top and bottom)
1052
1053 // Compute the areas. Attention: Only correct for trapezoids
1054 // where the twisting is done along the z-axis. In the general case
1055 // the computed surface area is more difficult. However this simplification
1056 // does not affect the tracking through the solid.
1057
1058 G4double a1 = fSide0->GetSurfaceArea();
1059 G4double a2 = fSide90->GetSurfaceArea();
1060 G4double a3 = fSide180->GetSurfaceArea() ;
1061 G4double a4 = fSide270->GetSurfaceArea() ;
1062 G4double a5 = fLowerEndcap->GetSurfaceArea() ;
1063 G4double a6 = fUpperEndcap->GetSurfaceArea() ;
1064
1065#ifdef G4TWISTDEBUG
1066 G4cout << "Surface 0 deg = " << a1 << G4endl ;
1067 G4cout << "Surface 90 deg = " << a2 << G4endl ;
1068 G4cout << "Surface 180 deg = " << a3 << G4endl ;
1069 G4cout << "Surface 270 deg = " << a4 << G4endl ;
1070 G4cout << "Surface Lower = " << a5 << G4endl ;
1071 G4cout << "Surface Upper = " << a6 << G4endl ;
1072#endif
1073
1074 G4double chose = (a1 + a2 + a3 + a4 + a5 + a6)*G4QuickRand() ;
1075
1076 if(chose < a1)
1077 {
1078 umin = fSide0->GetBoundaryMin(phi) ;
1079 umax = fSide0->GetBoundaryMax(phi) ;
1080 u = umin + (umax - umin)*G4QuickRand();
1081 return fSide0->SurfacePoint(phi, u, true) ; // point on 0 deg surface
1082 }
1083 if( (chose >= a1) && (chose < a1 + a2 ) )
1084 {
1085 umin = fSide90->GetBoundaryMin(phi) ;
1086 umax = fSide90->GetBoundaryMax(phi) ;
1087 u = umin + (umax - umin)*G4QuickRand();
1088 return fSide90->SurfacePoint(phi, u, true); // point on 90 deg surface
1089 }
1090 if( (chose >= a1 + a2 ) && (chose < a1 + a2 + a3 ) )
1091 {
1092 umin = fSide180->GetBoundaryMin(phi) ;
1093 umax = fSide180->GetBoundaryMax(phi) ;
1094 u = umin + (umax - umin)*G4QuickRand();
1095 return fSide180->SurfacePoint(phi, u, true); // point on 180 deg surface
1096 }
1097 if( (chose >= a1 + a2 + a3 ) && (chose < a1 + a2 + a3 + a4 ) )
1098 {
1099 umin = fSide270->GetBoundaryMin(phi) ;
1100 umax = fSide270->GetBoundaryMax(phi) ;
1101 u = umin + (umax - umin)*G4QuickRand();
1102 return fSide270->SurfacePoint(phi, u, true); // point on 270 deg surface
1103 }
1104 if( (chose >= a1 + a2 + a3 + a4 ) && (chose < a1 + a2 + a3 + a4 + a5 ) )
1105 {
1106 y = fDy1*(2.*G4QuickRand() - 1.);
1107 umin = fLowerEndcap->GetBoundaryMin(y) ;
1108 umax = fLowerEndcap->GetBoundaryMax(y) ;
1109 u = umin + (umax - umin)*G4QuickRand();
1110 return fLowerEndcap->SurfacePoint(u,y,true); // point on lower endcap
1111 }
1112 y = fDy2*(2.*G4QuickRand() - 1.);
1113 umin = fUpperEndcap->GetBoundaryMin(y) ;
1114 umax = fUpperEndcap->GetBoundaryMax(y) ;
1115 u = umin + (umax - umin)*G4QuickRand();
1116 return fUpperEndcap->SurfacePoint(u,y,true); // point on upper endcap
1117}
1118
1119
1120//=====================================================================
1121//* CreatePolyhedron --------------------------------------------------
1122
1124{
1125 // number of meshes
1126 const G4int k =
1128 std::abs(fPhiTwist) / twopi) + 2;
1129 const G4int n = k;
1130
1131 const G4int nnodes = 4*(k-1)*(n-2) + 2*k*k ;
1132 const G4int nfaces = 4*(k-1)*(n-1) + 2*(k-1)*(k-1) ;
1133
1134 auto ph = new G4Polyhedron;
1135 typedef G4double G4double3[3];
1136 typedef G4int G4int4[4];
1137 auto xyz = new G4double3[nnodes]; // number of nodes
1138 auto faces = new G4int4[nfaces] ; // number of faces
1139
1140 fLowerEndcap->GetFacets(k,k,xyz,faces,0) ;
1141 fUpperEndcap->GetFacets(k,k,xyz,faces,1) ;
1142 fSide270->GetFacets(k,n,xyz,faces,2) ;
1143 fSide0->GetFacets(k,n,xyz,faces,3) ;
1144 fSide90->GetFacets(k,n,xyz,faces,4) ;
1145 fSide180->GetFacets(k,n,xyz,faces,5) ;
1146
1147 ph->createPolyhedron(nnodes,nfaces,xyz,faces);
1148
1149 return ph;
1150}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4double B(G4double temperature)
@ JustWarning
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4QuickRand(uint32_t seed=0)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
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
G4String G4GeometryType
Definition G4VSolid.hh:70
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double x() const
double y() const
double z() 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 CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
static G4GeometryTolerance * GetInstance()
G4double GetAngularTolerance() const
G4TwistBoxSide describes a twisted boundary surface for a trapezoid.
virtual void AddSolid(const G4Box &)=0
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4String GetName() const
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
void DumpInfo() const
G4double kCarTolerance
Definition G4VSolid.hh:418
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108
G4VTwistSurface is a base class for boundary surface of a G4VSolid.
virtual G4ThreeVector GetNormal(const G4ThreeVector &p, G4bool isGlobal)=0
G4bool IsValidNorm() const
virtual G4double DistanceTo(const G4ThreeVector &gp, G4ThreeVector &gxx)
virtual G4double DistanceToOut(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector &gxxbest)
G4Polyhedron * GetPolyhedron() const override
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
G4ThreeVector GetPointOnSurface() const override
G4VTwistedFaceted & operator=(const G4VTwistedFaceted &rhs)
G4Polyhedron * fpPolyhedron
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4double GetCubicVolume() override
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
EInside Inside(const G4ThreeVector &p) const override
void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *) override
G4Polyhedron * CreatePolyhedron() const override
std::ostream & StreamInfo(std::ostream &os) const override
G4GeometryType GetEntityType() const override
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=false, G4bool *validnorm=nullptr, G4ThreeVector *n=nullptr) const override
G4VTwistedFaceted(const G4String &pname, G4double PhiTwist, G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlph)
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
G4VisExtent GetExtent() const override
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
G4double GetValueB(G4double phi) const
G4double Xcoef(G4double u, G4double phi, G4double ftg) const
G4double GetSurfaceArea() override
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
static G4int GetNumberOfRotationSteps()
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