Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Hype.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// Implementation of G4Hype
27//
28// Authors:
29// Ernesto Lamanna & Francesco Safai Tehrani
30// Rome, INFN & University of Rome "La Sapienza", 9 June 1998.
31// --------------------------------------------------------------------
32
33#include "G4Hype.hh"
34
35#if !(defined(G4GEOM_USE_UHYPE) && defined(G4GEOM_USE_SYS_USOLIDS))
36
37#include "G4GeomTools.hh"
38#include "G4VoxelLimits.hh"
39#include "G4AffineTransform.hh"
40#include "G4BoundingEnvelope.hh"
41#include "G4ClippablePolygon.hh"
42#include "G4QuickRand.hh"
43
45
46#include "G4VGraphicsScene.hh"
47#include "G4VisExtent.hh"
48
49#include "G4AutoLock.hh"
50
51namespace
52{
53 G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
54 G4Mutex hypeMutex = G4MUTEX_INITIALIZER;
55}
56
57using namespace CLHEP;
58
59// Constructor - check parameters, and fills protected data members
60//
62 G4double newInnerRadius,
63 G4double newOuterRadius,
64 G4double newInnerStereo,
65 G4double newOuterStereo,
66 G4double newHalfLenZ)
67 : G4VSolid(pName)
68{
69 fHalfTol = 0.5*kCarTolerance;
70
71 // Check z-len
72 //
73 if (newHalfLenZ <= 0)
74 {
75 std::ostringstream message;
76 message << "Invalid Z half-length in solid: " << GetName() << " !"
77 << "\nZ half-length: " << newHalfLenZ/mm << " mm";
78 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
79 FatalErrorInArgument, message);
80 }
81 halfLenZ = newHalfLenZ;
82
83 // Check radii
84 //
85 if (newInnerRadius < 0 || newOuterRadius < 0 || newInnerRadius >= newOuterRadius)
86 {
87 std::ostringstream message;
88 message << "Invalid radii in solid: " << GetName() << " !"
89 << "\nInner radius: " << newInnerRadius/mm << " mm"
90 << "\nOuter radius: " << newOuterRadius/mm << " mm";
91 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
92 FatalErrorInArgument, message);
93 }
94
95 innerRadius = newInnerRadius;
96 outerRadius = newOuterRadius;
97
98 innerRadius2 = innerRadius*innerRadius;
99 outerRadius2 = outerRadius*outerRadius;
100
101 SetInnerStereo(newInnerStereo);
102 SetOuterStereo(newOuterStereo);
103
104 // Check end radii
105 //
106 if (endInnerRadius > endOuterRadius)
107 {
108 std::ostringstream message;
109 message << "Inconsistent stereo angles in solid: " << GetName() << " !"
110 << "\nEnd inner radius: " << endInnerRadius/mm << " mm"
111 << "\nEnd outer radius: " << endOuterRadius/mm << " mm";
112 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
113 FatalErrorInArgument, message);
114 }
115}
116
117// Fake default constructor - sets only member data and allocates memory
118// for usage restricted to object persistency.
119//
120G4Hype::G4Hype( __void__& a )
121 : G4VSolid(a), innerRadius(0.), outerRadius(0.), halfLenZ(0.), innerStereo(0.),
122 outerStereo(0.), tanInnerStereo(0.), tanOuterStereo(0.), tanInnerStereo2(0.),
123 tanOuterStereo2(0.), innerRadius2(0.), outerRadius2(0.), endInnerRadius2(0.),
124 endOuterRadius2(0.), endInnerRadius(0.), endOuterRadius(0.), fHalfTol(0.)
125{
126}
127
128// Destructor
129//
131{
132 delete fpPolyhedron; fpPolyhedron = nullptr;
133}
134
135// Copy constructor
136//
138 : G4VSolid(rhs), innerRadius(rhs.innerRadius),
139 outerRadius(rhs.outerRadius), halfLenZ(rhs.halfLenZ),
140 innerStereo(rhs.innerStereo), outerStereo(rhs.outerStereo),
141 tanInnerStereo(rhs.tanInnerStereo), tanOuterStereo(rhs.tanOuterStereo),
142 tanInnerStereo2(rhs.tanInnerStereo2), tanOuterStereo2(rhs.tanOuterStereo2),
143 innerRadius2(rhs.innerRadius2), outerRadius2(rhs.outerRadius2),
144 endInnerRadius2(rhs.endInnerRadius2), endOuterRadius2(rhs.endOuterRadius2),
145 endInnerRadius(rhs.endInnerRadius), endOuterRadius(rhs.endOuterRadius),
146 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
147 fInnerSurfaceArea(rhs.fInnerSurfaceArea), fOuterSurfaceArea(rhs.fOuterSurfaceArea),
148 fHalfTol(rhs.fHalfTol)
149{
150}
151
152// Assignment operator
153//
155{
156 // Check assignment to self
157 //
158 if (this == &rhs) { return *this; }
159
160 // Copy base class data
161 //
163
164 // Copy data
165 //
166 innerRadius = rhs.innerRadius; outerRadius = rhs.outerRadius;
167 halfLenZ = rhs.halfLenZ;
168 innerStereo = rhs.innerStereo; outerStereo = rhs.outerStereo;
169 tanInnerStereo = rhs.tanInnerStereo; tanOuterStereo = rhs.tanOuterStereo;
170 tanInnerStereo2 = rhs.tanInnerStereo2; tanOuterStereo2 = rhs.tanOuterStereo2;
171 innerRadius2 = rhs.innerRadius2; outerRadius2 = rhs.outerRadius2;
172 endInnerRadius2 = rhs.endInnerRadius2; endOuterRadius2 = rhs.endOuterRadius2;
173 endInnerRadius = rhs.endInnerRadius; endOuterRadius = rhs.endOuterRadius;
174 fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
175 fInnerSurfaceArea = rhs.fInnerSurfaceArea; fOuterSurfaceArea = rhs.fOuterSurfaceArea;
176 fHalfTol = rhs.fHalfTol;
177 fRebuildPolyhedron = false;
178 delete fpPolyhedron; fpPolyhedron = nullptr;
179
180 return *this;
181}
182
183// Set inner radius at z = 0
184//
186{
187 innerRadius = newIRad;
188 innerRadius2 = newIRad*newIRad;
189 endInnerRadius2 = HypeInnerRadius2(halfLenZ);
190 endInnerRadius = std::sqrt(endInnerRadius2);
191 fCubicVolume = 0.;
192 fSurfaceArea = 0.;
193 fInnerSurfaceArea =
194 G4GeomTools::HyperboloidSurfaceArea(twopi, innerRadius, tanInnerStereo, -halfLenZ, halfLenZ);
195 fRebuildPolyhedron = true;
196}
197
198// Set outer radius at z = 0
199//
201{
202 outerRadius = newORad;
203 outerRadius2 = newORad*newORad;
204 endOuterRadius2 = HypeOuterRadius2(halfLenZ);
205 endOuterRadius = std::sqrt(endOuterRadius2);
206 fCubicVolume = 0.;
207 fSurfaceArea = 0.;
208 fOuterSurfaceArea =
209 G4GeomTools::HyperboloidSurfaceArea(twopi, outerRadius, tanOuterStereo, -halfLenZ, halfLenZ);
210 fRebuildPolyhedron = true;
211}
212
213// Set half z-length
214//
216{
217 halfLenZ = newHLZ;
218 endInnerRadius2 = HypeInnerRadius2(halfLenZ);
219 endInnerRadius = std::sqrt(endInnerRadius2);
220 endOuterRadius2 = HypeOuterRadius2(halfLenZ);
221 endOuterRadius = std::sqrt(endOuterRadius2);
222 fCubicVolume = 0.;
223 fSurfaceArea = 0.;
224 fInnerSurfaceArea =
225 G4GeomTools::HyperboloidSurfaceArea(twopi, innerRadius, tanInnerStereo, -halfLenZ, halfLenZ);
226 fOuterSurfaceArea =
227 G4GeomTools::HyperboloidSurfaceArea(twopi, outerRadius, tanOuterStereo, -halfLenZ, halfLenZ);
228 fRebuildPolyhedron = true;
229}
230
231// Set inner stereo angle
232//
234{
235 innerStereo = std::abs(newISte);
236 tanInnerStereo = std::tan(innerStereo);
237 tanInnerStereo2 = tanInnerStereo*tanInnerStereo;
238 endInnerRadius2 = HypeInnerRadius2(halfLenZ);
239 endInnerRadius = std::sqrt(endInnerRadius2);
240 fCubicVolume = 0.;
241 fSurfaceArea = 0.;
242 fInnerSurfaceArea =
243 G4GeomTools::HyperboloidSurfaceArea(twopi, innerRadius, tanInnerStereo, -halfLenZ, halfLenZ);
244 fRebuildPolyhedron = true;
245}
246
247// Set outer stereo angle
248//
250{
251 outerStereo = std::abs(newOSte);
252 tanOuterStereo = std::tan(outerStereo);
253 tanOuterStereo2 = tanOuterStereo*tanOuterStereo;
254 endOuterRadius2 = HypeOuterRadius2(halfLenZ);
255 endOuterRadius = std::sqrt(endOuterRadius2);
256 fCubicVolume = 0.;
257 fSurfaceArea = 0.;
258 fOuterSurfaceArea =
259 G4GeomTools::HyperboloidSurfaceArea(twopi, outerRadius, tanOuterStereo, -halfLenZ, halfLenZ);
260 fRebuildPolyhedron = true;
261}
262
263// Dispatch to parameterisation for replication mechanism dimension
264// computation & modification.
265//
267 const G4int n,
268 const G4VPhysicalVolume* pRep)
269{
270 p->ComputeDimensions(*this,n,pRep);
271}
272
273// Get bounding box
274//
276{
277 pMin.set(-endOuterRadius,-endOuterRadius,-halfLenZ);
278 pMax.set( endOuterRadius, endOuterRadius, halfLenZ);
279
280 // Check correctness of the bounding box
281 //
282 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
283 {
284 std::ostringstream message;
285 message << "Bad bounding box (min >= max) for solid: " << GetName() << " !"
286 << "\npMin = " << pMin
287 << "\npMax = " << pMax;
288 G4Exception("G4Hype::BoundingLimits()", "GeomMgt0001",
289 JustWarning, message);
290 DumpInfo();
291 }
292}
293
294// Calculate extent under transform and specified limit
295//
297 const G4VoxelLimits& pVoxelLimit,
298 const G4AffineTransform& pTransform,
299 G4double& pMin, G4double& pMax) const
300{
301 G4ThreeVector bmin, bmax;
302
303 // Get bounding box
304 BoundingLimits(bmin,bmax);
305
306 // Find extent
307 G4BoundingEnvelope bbox(bmin,bmax);
308 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
309}
310
311// Decides whether point is inside, outside or on the surface
312//
314{
315 //
316 // Check z extents: are we outside?
317 //
318 const G4double absZ(std::fabs(p.z()));
319 if (absZ > halfLenZ + fHalfTol) { return kOutside; }
320
321 //
322 // Check outer radius
323 //
324 const G4double oRad2(HypeOuterRadius2(absZ));
325 const G4double xR2( p.x()*p.x()+p.y()*p.y() );
326
327 if (xR2 > oRad2 + kCarTolerance*endOuterRadius) { return kOutside; }
328
329 if (xR2 > oRad2 - kCarTolerance*endOuterRadius) { return kSurface; }
330
331 if (InnerSurfaceExists())
332 {
333 //
334 // Check inner radius
335 //
336 const G4double iRad2(HypeInnerRadius2(absZ));
337
338 if (xR2 < iRad2 - kCarTolerance*endInnerRadius) { return kOutside; }
339
340 if (xR2 < iRad2 + kCarTolerance*endInnerRadius) { return kSurface; }
341 }
342
343 //
344 // We are inside in radius, now check endplate surface
345 //
346 if (absZ > halfLenZ - fHalfTol) { return kSurface; }
347
348 return kInside;
349}
350
351// Returns the normal unit vector to the Hyperbolical Surface at a point
352// p on (or nearly on) the surface
353//
355{
356 //
357 // Which of the three or four surfaces are we closest to?
358 //
359 const G4double absZ(std::fabs(p.z()));
360 const G4double distZ(absZ - halfLenZ);
361 const G4double dist2Z(distZ*distZ);
362
363 const G4double xR2( p.x()*p.x()+p.y()*p.y() );
364 const G4double dist2Outer( std::fabs(xR2 - HypeOuterRadius2(absZ)) );
365
366 if (InnerSurfaceExists())
367 {
368 //
369 // Has inner surface: is this closest?
370 //
371 const G4double dist2Inner( std::fabs(xR2 - HypeInnerRadius2(absZ)) );
372 if (dist2Inner < dist2Z && dist2Inner < dist2Outer)
373 {
374 return G4ThreeVector( -p.x(), -p.y(), p.z()*tanInnerStereo2 ).unit();
375 }
376 }
377
378 //
379 // Do the "endcaps" win?
380 //
381 if (dist2Z < dist2Outer)
382 {
383 return { 0.0, 0.0, (G4double)(p.z() < 0 ? -1.0 : 1.0) };
384 }
385
386 //
387 // Outer surface wins
388 //
389 return G4ThreeVector( p.x(), p.y(), -p.z()*tanOuterStereo2 ).unit();
390}
391
392// Calculates distance to shape from outside, along normalised vector
393// - return kInfinity if no intersection,
394// or intersection distance <= tolerance
395//
396// Calculating the intersection of a line with the surfaces
397// is fairly straight forward. The difficult problem is dealing
398// with the intersections of the surfaces in a consistent manner,
399// and this accounts for the complicated logic.
400//
402 const G4ThreeVector& v ) const
403{
404 //
405 // Quick test. Beware! This assumes v is a unit vector!
406 //
407 if (std::fabs(p.x()*v.y() - p.y()*v.x()) > endOuterRadius+kCarTolerance)
408 {
409 return kInfinity;
410 }
411
412 //
413 // Take advantage of z symmetry, and reflect throught the
414 // z=0 plane so that pz is always positive
415 //
416 G4double pz(p.z()), vz(v.z());
417 if (pz < 0)
418 {
419 pz = -pz;
420 vz = -vz;
421 }
422
423 //
424 // We must be very careful if we don't want to
425 // create subtle leaks at the edges where the
426 // hyperbolic surfaces connect to the endplate.
427 // The only reliable way to do so is to make sure
428 // that the decision as to when a track passes
429 // over the edge of one surface is exactly the
430 // same decision as to when a track passes into the
431 // other surface. By "exact", we don't mean algebraicly
432 // exact, but we mean the same machine instructions
433 // should be used.
434 //
435 G4bool couldMissOuter(true),
436 couldMissInner(true),
437 cantMissInnerCylinder(false);
438
439 //
440 // Check endplate intersection
441 //
442 G4double sigz = pz-halfLenZ;
443
444 if (sigz > -fHalfTol) // equivalent to: if (pz > halfLenZ - fHalfTol)
445 {
446 //
447 // We start in front of the endplate (within roundoff)
448 // Correct direction to intersect endplate?
449 //
450 if (vz >= 0)
451 {
452 //
453 // Nope. As long as we are far enough away, we
454 // can't intersect anything
455 //
456 if (sigz > 0) { return kInfinity; }
457
458 //
459 // Otherwise, we may still hit a hyperbolic surface
460 // if the point is on the hyperbolic surface (within tolerance)
461 //
462 G4double pr2 = p.x()*p.x() + p.y()*p.y();
463 if (pr2 > endOuterRadius2 + kCarTolerance*endOuterRadius)
464 {
465 return kInfinity;
466 }
467
468 if (InnerSurfaceExists())
469 {
470 if (pr2 < endInnerRadius2 - kCarTolerance*endInnerRadius)
471 {
472 return kInfinity;
473 }
474 if ( (pr2 < endOuterRadius2 - kCarTolerance*endOuterRadius)
475 && (pr2 > endInnerRadius2 + kCarTolerance*endInnerRadius) )
476 {
477 return kInfinity;
478 }
479 }
480 else
481 {
482 if (pr2 < endOuterRadius2 - kCarTolerance*endOuterRadius)
483 {
484 return kInfinity;
485 }
486 }
487 }
488 else
489 {
490 //
491 // Where do we intersect at z = halfLenZ?
492 //
493 G4double q = -sigz/vz;
494 G4double xi = p.x() + q*v.x(),
495 yi = p.y() + q*v.y();
496
497 //
498 // Is this on the endplate? If so, return s, unless
499 // we are on the tolerant surface, in which case return 0
500 //
501 G4double pr2 = xi*xi + yi*yi;
502 if (pr2 <= endOuterRadius2)
503 {
504 if (InnerSurfaceExists())
505 {
506 if (pr2 >= endInnerRadius2) { return (sigz < fHalfTol) ? 0 : q; }
507 //
508 // This test is sufficient to ensure that the
509 // trajectory cannot miss the inner hyperbolic surface
510 // for z > 0, if the normal is correct.
511 //
512 G4double dot1 = (xi*v.x() + yi*v.y())*endInnerRadius/std::sqrt(pr2);
513 couldMissInner = (dot1 - halfLenZ*tanInnerStereo2*vz <= 0);
514
515 if (pr2 > endInnerRadius2*(1 - 2*DBL_EPSILON) )
516 {
517 //
518 // There is a potential leak if the inner
519 // surface is a cylinder
520 //
521 if ( (innerStereo < DBL_MIN)
522 && ((std::fabs(v.x()) > DBL_MIN) || (std::fabs(v.y()) > DBL_MIN)))
523 {
524 cantMissInnerCylinder = true;
525 }
526 }
527 }
528 else
529 {
530 return (sigz < fHalfTol) ? 0 : q;
531 }
532 }
533 else
534 {
535 G4double dotR( xi*v.x() + yi*v.y() );
536 if (dotR >= 0)
537 {
538 //
539 // Otherwise, if we are traveling outwards, we know
540 // we must miss the hyperbolic surfaces also, so
541 // we need not bother checking
542 //
543 return kInfinity;
544 }
545
546 //
547 // This test is sufficient to ensure that the
548 // trajectory cannot miss the outer hyperbolic surface
549 // for z > 0, if the normal is correct.
550 //
551 G4double dot1 = dotR*endOuterRadius/std::sqrt(pr2);
552 couldMissOuter = (dot1 - halfLenZ*tanOuterStereo2*vz>= 0);
553 }
554 }
555 }
556
557 //
558 // Check intersection with outer hyperbolic surface, save
559 // distance to valid intersection into "best".
560 //
561 G4double best = kInfinity;
562
563 G4double q[2];
564 G4int n = IntersectHype( p, v, outerRadius2, tanOuterStereo2, q );
565
566 if (n > 0)
567 {
568 //
569 // Potential intersection: is p on this surface?
570 //
571 if (pz < halfLenZ+fHalfTol)
572 {
573 G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeOuterRadius2(pz);
574 if (std::fabs(dr2) < kCarTolerance*endOuterRadius)
575 {
576 //
577 // Sure, but make sure we're traveling inwards at
578 // this point
579 //
580 if (p.x()*v.x() + p.y()*v.y() - pz*tanOuterStereo2*vz < 0) { return 0; }
581 }
582 }
583
584 //
585 // We are now certain that p is not on the tolerant surface.
586 // Accept only position distance q
587 //
588 for( G4int i=0; i<n; ++i )
589 {
590 if (q[i] >= 0)
591 {
592 //
593 // Check to make sure this intersection point is
594 // on the surface, but only do so if we haven't
595 // checked the endplate intersection already
596 //
597 G4double zi = pz + q[i]*vz;
598
599 if (zi < -halfLenZ) { continue; }
600 if (zi > +halfLenZ && couldMissOuter) { continue; }
601
602 //
603 // Check normal
604 //
605 G4double xi = p.x() + q[i]*v.x(),
606 yi = p.y() + q[i]*v.y();
607
608 if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz > 0) { continue; }
609
610 best = q[i];
611 break;
612 }
613 }
614 }
615
616 if (!InnerSurfaceExists()) { return best; }
617
618 //
619 // Check intersection with inner hyperbolic surface
620 //
621 n = IntersectHype( p, v, innerRadius2, tanInnerStereo2, q );
622 if (n == 0)
623 {
624 if (cantMissInnerCylinder)
625 {
626 return (sigz < fHalfTol) ? 0 : -sigz/vz;
627 }
628
629 return best;
630 }
631
632 //
633 // P on this surface?
634 //
635 if (pz < halfLenZ+fHalfTol)
636 {
637 G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeInnerRadius2(pz);
638 if (std::fabs(dr2) < kCarTolerance*endInnerRadius)
639 {
640 //
641 // Sure, but make sure we're traveling outwards at
642 // this point
643 //
644 if (p.x()*v.x() + p.y()*v.y() - pz*tanInnerStereo2*vz > 0) { return 0; }
645 }
646 }
647
648 //
649 // No, so only positive q is valid. Search for a valid intersection
650 // that is closer than the outer intersection (if it exists)
651 //
652 for( G4int i=0; i<n; ++i )
653 {
654 if (q[i] > best) { break; }
655 if (q[i] >= 0)
656 {
657 //
658 // Check to make sure this intersection point is
659 // on the surface, but only do so if we haven't
660 // checked the endplate intersection already
661 //
662 G4double zi = pz + q[i]*vz;
663
664 if (zi < -halfLenZ) { continue; }
665 if (zi > +halfLenZ && couldMissInner) { continue; }
666
667 //
668 // Check normal
669 //
670 G4double xi = p.x() + q[i]*v.x(),
671 yi = p.y() + q[i]*v.y();
672
673 if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz < 0) { continue; }
674
675 best = q[i];
676 break;
677 }
678 }
679
680 //
681 // Done
682 //
683 return best;
684}
685
686// Calculates distance to shape from outside, along perpendicular direction
687// (if one exists). May be an underestimate.
688//
689// There are five (r,z) regions:
690// 1. a point that is beyond the endcap but within the
691// endcap radii
692// 2. a point with r > outer endcap radius and with
693// a z position that is beyond the cone formed by the
694// normal of the outer hyperbolic surface at the
695// edge at which it meets the endcap.
696// 3. a point that is outside the outer surface and not in (1 or 2)
697// 4. a point that is inside the inner surface and not in (5)
698// 5. a point with radius < inner endcap radius and
699// with a z position beyond the cone formed by the
700// normal of the inner hyperbolic surface at the
701// edge at which it meets the endcap.
702// (regions 4 and 5 only exist if there is an inner surface)
703//
705{
706 G4double absZ(std::fabs(p.z()));
707
708 //
709 // Check region
710 //
711 G4double r2 = p.x()*p.x() + p.y()*p.y();
712 G4double r = std::sqrt(r2);
713
714 G4double sigz = absZ - halfLenZ;
715
716 if (r < endOuterRadius)
717 {
718 if (sigz > -fHalfTol)
719 {
720 if (InnerSurfaceExists())
721 {
722 if (r > endInnerRadius)
723 {
724 return sigz < fHalfTol ? 0 : sigz; // Region 1
725 }
726
727 G4double dr = endInnerRadius - r;
728 if (sigz > dr*tanInnerStereo2)
729 {
730 //
731 // In region 5
732 //
733 G4double answer = std::sqrt( dr*dr + sigz*sigz );
734 return answer < fHalfTol ? 0 : answer;
735 }
736 }
737 else
738 {
739 //
740 // In region 1 (no inner surface)
741 //
742 return sigz < fHalfTol ? 0 : sigz;
743 }
744 }
745 }
746 else
747 {
748 G4double dr = r - endOuterRadius;
749 if (sigz > -dr*tanOuterStereo2)
750 {
751 //
752 // In region 2
753 //
754 G4double answer = std::sqrt( dr*dr + sigz*sigz );
755 return answer < fHalfTol ? 0 : answer;
756 }
757 }
758
759 if (InnerSurfaceExists())
760 {
761 if (r2 < HypeInnerRadius2(absZ)+kCarTolerance*endInnerRadius)
762 {
763 //
764 // In region 4
765 //
766 G4double answer = ApproxDistInside( r,absZ,innerRadius,tanInnerStereo2 );
767 return answer < fHalfTol ? 0 : answer;
768 }
769 }
770
771 //
772 // We are left by elimination with region 3
773 //
774 G4double answer = ApproxDistOutside( r, absZ, outerRadius, tanOuterStereo );
775 return answer < fHalfTol ? 0 : answer;
776}
777
778// Calculates distance to surface of shape from 'inside', allowing for tolerance
779//
780// The situation here is much simplier than DistanceToIn(p,v). For
781// example, there is no need to even check whether an intersection
782// point is inside the boundary of a surface, as long as all surfaces
783// are checked and the smallest distance is used.
784//
786 const G4bool calcNorm,
787 G4bool* validNorm, G4ThreeVector* norm ) const
788{
789 static const G4ThreeVector normEnd1(0.0,0.0,+1.0);
790 static const G4ThreeVector normEnd2(0.0,0.0,-1.0);
791
792 //
793 // Keep track of closest surface
794 //
795 G4double sBest; // distance to
796 const G4ThreeVector* nBest; // normal vector
797 G4bool vBest; // whether "valid"
798
799 //
800 // Check endplate, taking advantage of symmetry.
801 // Note that the endcap is the only surface which
802 // has a "valid" normal, i.e. is a surface of which
803 // the entire solid is behind.
804 //
805 G4double pz(p.z()), vz(v.z());
806 if (vz < 0)
807 {
808 pz = -pz;
809 vz = -vz;
810 nBest = &normEnd2;
811 }
812 else
813 {
814 nBest = &normEnd1;
815 }
816
817 //
818 // Possible intercept. Are we on the surface?
819 //
820 if (pz > halfLenZ-fHalfTol)
821 {
822 if (calcNorm) { *norm = *nBest; *validNorm = true; }
823 return 0;
824 }
825
826 //
827 // Nope. Get distance. Beware of zero vz.
828 //
829 sBest = (vz > DBL_MIN) ? (halfLenZ - pz)/vz : kInfinity;
830 vBest = true;
831
832 //
833 // Check outer surface
834 //
835 G4double r2 = p.x()*p.x() + p.y()*p.y();
836
837 G4double q[2];
838 G4int n = IntersectHype( p, v, outerRadius2, tanOuterStereo2, q );
839
840 G4ThreeVector norm1, norm2;
841
842 if (n > 0)
843 {
844 //
845 // We hit somewhere. Are we on the surface?
846 //
847 G4double dr2 = r2 - HypeOuterRadius2(pz);
848 if (std::fabs(dr2) < endOuterRadius*kCarTolerance)
849 {
850 G4ThreeVector normHere( p.x(), p.y(), -p.z()*tanOuterStereo2 );
851 //
852 // Sure. But are we going the right way?
853 //
854 if (normHere.dot(v) > 0)
855 {
856 if (calcNorm) { *norm = normHere.unit(); *validNorm = false; }
857 return 0;
858 }
859 }
860
861 //
862 // Nope. Check closest positive intercept.
863 //
864 for( G4int i=0; i<n; ++i )
865 {
866 if (q[i] > sBest) { break; }
867 if (q[i] > 0)
868 {
869 //
870 // Make sure normal is correct (that this
871 // solution is an outgoing solution)
872 //
873 G4ThreeVector pk(p+q[i]*v);
874 norm1 = G4ThreeVector( pk.x(), pk.y(), -pk.z()*tanOuterStereo2 );
875 if (norm1.dot(v) > 0)
876 {
877 sBest = q[i];
878 nBest = &norm1;
879 vBest = false;
880 break;
881 }
882 }
883 }
884 }
885
886 if (InnerSurfaceExists())
887 {
888 //
889 // Check inner surface
890 //
891 n = IntersectHype( p, v, innerRadius2, tanInnerStereo2, q );
892 if (n > 0)
893 {
894 //
895 // On surface?
896 //
897 G4double dr2 = r2 - HypeInnerRadius2(pz);
898 if (std::fabs(dr2) < endInnerRadius*kCarTolerance)
899 {
900 G4ThreeVector normHere( -p.x(), -p.y(), p.z()*tanInnerStereo2 );
901 if (normHere.dot(v) > 0)
902 {
903 if (calcNorm)
904 {
905 *norm = normHere.unit();
906 *validNorm = false;
907 }
908 return 0;
909 }
910 }
911
912 //
913 // Check closest positive
914 //
915 for( G4int i=0; i<n; ++i )
916 {
917 if (q[i] > sBest) { break; }
918 if (q[i] > 0)
919 {
920 G4ThreeVector pk(p+q[i]*v);
921 norm2 = G4ThreeVector( -pk.x(), -pk.y(), pk.z()*tanInnerStereo2 );
922 if (norm2.dot(v) > 0)
923 {
924 sBest = q[i];
925 nBest = &norm2;
926 vBest = false;
927 break;
928 }
929 }
930 }
931 }
932 }
933
934 //
935 // Done!
936 //
937 if (calcNorm)
938 {
939 *validNorm = vBest;
940
941 if (nBest == &norm1 || nBest == &norm2)
942 {
943 *norm = nBest->unit();
944 }
945 else
946 {
947 *norm = *nBest;
948 }
949 }
950
951 return sBest;
952}
953
954// Calculates distance (<=actual) to closest surface of shape from inside
955//
956// May be an underestimate
957//
959{
960 //
961 // Try each surface and remember the closest
962 //
963 G4double absZ(std::fabs(p.z()));
964 G4double r(p.perp());
965
966 G4double sBest = halfLenZ - absZ;
967
968 G4double tryOuter = ApproxDistInside( r, absZ, outerRadius, tanOuterStereo2 );
969 if (tryOuter < sBest) { sBest = tryOuter; }
970
971 if (InnerSurfaceExists())
972 {
973 G4double tryInner = ApproxDistOutside( r,absZ,innerRadius,tanInnerStereo );
974 if (tryInner < sBest) { sBest = tryInner; }
975 }
976
977 return sBest < 0.5*kCarTolerance ? 0 : sBest;
978}
979
980// IntersectHype
981//
982// Decides if and where a line intersects with a hyperbolic
983// surface (of infinite extent)
984//
985// Arguments:
986// p - (in) Point on trajectory
987// v - (in) Vector along trajectory
988// r2 - (in) Square of radius at z = 0
989// tan2phi - (in) std::tan(phi)**2
990// q - (out) Up to two points of intersection, where the
991// intersection point is p + q*v, and if there are
992// two intersections, q[0] < q[1]. May be negative.
993// Returns:
994// The number of intersections. If 0, the trajectory misses.
995//
996//
997// Equation of a line:
998//
999// x = x0 + q*tx y = y0 + q*ty z = z0 + q*tz
1000//
1001// Equation of a hyperbolic surface:
1002//
1003// x**2 + y**2 = r**2 + (z*tanPhi)**2
1004//
1005// Solution is quadratic:
1006//
1007// a*q**2 + b*q + c = 0
1008//
1009// where:
1010//
1011// a = tx**2 + ty**2 - (tz*tanPhi)**2
1012//
1013// b = 2*( x0*tx + y0*ty - z0*tz*tanPhi**2 )
1014//
1015// c = x0**2 + y0**2 - r**2 - (z0*tanPhi)**2
1016//
1017//
1018G4int G4Hype::IntersectHype( const G4ThreeVector &p, const G4ThreeVector &v,
1019 G4double r2, G4double tan2Phi, G4double ss[2] ) const
1020{
1021 G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
1022 G4double tx = v.x(), ty = v.y(), tz = v.z();
1023
1024 G4double a = tx*tx + ty*ty - tz*tz*tan2Phi;
1025 G4double b = 2*( x0*tx + y0*ty - z0*tz*tan2Phi );
1026 G4double c = x0*x0 + y0*y0 - r2 - z0*z0*tan2Phi;
1027
1028 if (std::fabs(a) < DBL_MIN)
1029 {
1030 //
1031 // The trajectory is parallel to the asympotic limit of
1032 // the surface: single solution
1033 //
1034 if (std::fabs(b) < DBL_MIN) { return 0; }
1035 // Unless we travel through exact center
1036
1037 ss[0] = c/b;
1038 return 1;
1039 }
1040
1041 G4double radical = b*b - 4*a*c;
1042
1043 if (radical < -DBL_MIN) { return 0; } // No solution
1044
1045 if (radical < DBL_MIN)
1046 {
1047 //
1048 // Grazes surface
1049 //
1050 ss[0] = -b/a/2.0;
1051 return 1;
1052 }
1053
1054 radical = std::sqrt(radical);
1055
1056 G4double q = -0.5*( b + (b < 0 ? -radical : +radical) );
1057 G4double sa = q/a;
1058 G4double sb = c/q;
1059 if (sa < sb)
1060 {
1061 ss[0] = sa; ss[1] = sb;
1062 }
1063 else
1064 {
1065 ss[0] = sb; ss[1] = sa;
1066 }
1067 return 2;
1068}
1069
1070// ApproxDistOutside
1071//
1072// Finds the approximate distance of a point outside
1073// (greater radius) of a hyperbolic surface. The distance
1074// must be an underestimate. It will also be nice (although
1075// not necesary) that the estimate is always finite no
1076// matter how close the point is.
1077//
1078// Our hyperbola approaches the asymptotic limit at z = +/- infinity
1079// to the lines r = z*tanPhi. We call these lines the
1080// asymptotic limit line.
1081//
1082// We need the distance of the 2d point p(r,z) to the
1083// hyperbola r**2 = r0**2 + (z*tanPhi)**2. Find two
1084// points that bracket the true normal and use the
1085// distance to the line that connects these two points.
1086// The first such point is z=p.z. The second point is
1087// the z position on the asymptotic limit line that
1088// contains the normal on the line through the point p.
1089//
1090G4double G4Hype::ApproxDistOutside( G4double pr, G4double pz,
1091 G4double r0, G4double tanPhi ) const
1092{
1093 if (tanPhi < DBL_MIN) { return pr-r0; }
1094
1095 G4double tan2Phi = tanPhi*tanPhi;
1096
1097 //
1098 // First point
1099 //
1100 G4double z1 = pz;
1101 G4double r1 = std::sqrt( r0*r0 + z1*z1*tan2Phi );
1102
1103 //
1104 // Second point
1105 //
1106 G4double z2 = (pr*tanPhi + pz)/(1 + tan2Phi);
1107 G4double r2 = std::sqrt( r0*r0 + z2*z2*tan2Phi );
1108
1109 //
1110 // Line between them
1111 //
1112 G4double dr = r2-r1;
1113 G4double dz = z2-z1;
1114
1115 G4double len = std::sqrt(dr*dr + dz*dz);
1116 if (len < DBL_MIN)
1117 {
1118 //
1119 // The two points are the same?? I guess we
1120 // must have really bracketed the normal
1121 //
1122 dr = pr-r1;
1123 dz = pz-z1;
1124 return std::sqrt( dr*dr + dz*dz );
1125 }
1126
1127 //
1128 // Distance
1129 //
1130 return std::fabs((pr-r1)*dz - (pz-z1)*dr)/len;
1131}
1132
1133// ApproxDistInside
1134//
1135// Finds the approximate distance of a point inside
1136// of a hyperbolic surface. The distance
1137// must be an underestimate. It will also be nice (although
1138// not necesary) that the estimate is always finite no
1139// matter how close the point is.
1140//
1141// This estimate uses the distance to a line tangent to
1142// the hyperbolic function. The point of tangent is chosen
1143// by the z position point
1144//
1145// Assumes pr and pz are positive
1146//
1147G4double G4Hype::ApproxDistInside( G4double pr, G4double pz,
1148 G4double r0, G4double tan2Phi ) const
1149{
1150 if (tan2Phi < DBL_MIN) { return r0 - pr; }
1151
1152 //
1153 // Corresponding position and normal on hyperbolic
1154 //
1155 G4double rh = std::sqrt( r0*r0 + pz*pz*tan2Phi );
1156
1157 G4double dr = -rh;
1158 G4double dz = pz*tan2Phi;
1159 G4double len = std::sqrt(dr*dr + dz*dz);
1160
1161 //
1162 // Answer
1163 //
1164 return std::fabs((pr-rh)*dr)/len;
1165}
1166
1167// GetEntityType
1168//
1170{
1171 return {"G4Hype"};
1172}
1173
1174// Clone
1175//
1177{
1178 return new G4Hype(*this);
1179}
1180
1181// Streams object contents to an output stream
1182//
1183std::ostream& G4Hype::StreamInfo(std::ostream& os) const
1184{
1185 G4long oldprc = os.precision(16);
1186 os << "-----------------------------------------------------------\n"
1187 << " *** Dump for solid - " << GetName() << " ***\n"
1188 << " ===================================================\n"
1189 << " Solid type: G4Hype\n"
1190 << " Parameters: \n"
1191 << " half length Z: " << halfLenZ/mm << " mm \n"
1192 << " inner radius : " << innerRadius/mm << " mm \n"
1193 << " outer radius : " << outerRadius/mm << " mm \n"
1194 << " inner stereo angle : " << innerStereo/degree << " degrees \n"
1195 << " outer stereo angle : " << outerStereo/degree << " degrees \n"
1196 << "-----------------------------------------------------------\n";
1197 os.precision(oldprc);
1198
1199 return os;
1200}
1201
1202// Pick random point on surface
1203//
1205{
1206 G4double sbases = twopi*(endOuterRadius2 - endInnerRadius2);
1207 G4double stotal = fInnerSurfaceArea + fOuterSurfaceArea + sbases;
1208 G4double select = stotal*G4QuickRand();
1209
1210 G4double h = halfLenZ;
1211 G4double phi = twopi*G4QuickRand();
1212 G4double cosphi = std::cos(phi);
1213 G4double sinphi = std::sin(phi);
1214 if (select < sbases)
1215 {
1216 // circular bases
1217 G4double u = G4QuickRand();
1218 G4double rho = std::sqrt(u*endOuterRadius2 + (1. - u)*endInnerRadius2);
1219 G4double z = (select < 0.5*sbases) ? -h : h;
1220 return { rho*cosphi, rho*sinphi, z };
1221 }
1222
1223 // lateral surfaces (rejection sampling)
1224 G4double hh = h*h;
1225 G4double rr = (select < stotal - fInnerSurfaceArea) ? outerRadius2 : innerRadius2;
1226 G4double endrr = (select < stotal - fInnerSurfaceArea) ? endOuterRadius2 : endInnerRadius2;
1227 G4double delrr = endrr - rr;
1228 // max surface element: std::sqrt(aa*hh + (RR - aa)*(RR - aa + hh)*1.0)
1229 G4double ss_max = rr*hh + delrr*(delrr + hh);
1230 G4double u = 0.;
1231 for (auto i = 0; i < 10000; ++i)
1232 {
1233 u = 2.*G4QuickRand() - 1.;
1234 // surface element: std::sqrt(aa*hh + (RR - aa)*(RR - aa + hh)*uu)
1235 G4double ss = rr*hh + delrr*(delrr + hh)*u*u;
1236 if (ss_max*sqr(G4QuickRand()) <= ss) { break; }
1237 }
1238 G4double z = h*u;
1239 G4double rho = std::sqrt(rr + delrr*u*u);
1240 return { rho*cosphi, rho*sinphi, z };
1241}
1242
1243// GetCubicVolume
1244//
1246{
1247 if (fCubicVolume == 0)
1248 {
1249 G4AutoLock l(&hypeMutex);
1250 fCubicVolume =
1251 twopi*halfLenZ*(2.*(outerRadius2 - innerRadius2) + endOuterRadius2 - endInnerRadius2)/3.;
1252 l.unlock();
1253 }
1254 return fCubicVolume;
1255}
1256
1257// GetSurfaceArea
1258//
1260{
1261 if (fSurfaceArea == 0)
1262 {
1263 G4AutoLock l(&hypeMutex);
1264 fSurfaceArea = fInnerSurfaceArea + fOuterSurfaceArea + twopi*(endOuterRadius2 - endInnerRadius2);
1265 l.unlock();
1266 }
1267 return fSurfaceArea;
1268}
1269
1270// DescribeYourselfTo
1271//
1273{
1274 scene.AddSolid (*this);
1275}
1276
1277// GetExtent
1278//
1280{
1281 // Define the sides of the box into which the G4Tubs instance would fit.
1282 //
1283 return { -endOuterRadius, endOuterRadius,
1284 -endOuterRadius, endOuterRadius,
1285 -halfLenZ, halfLenZ };
1286}
1287
1288// CreatePolyhedron
1289//
1291{
1292 return new G4PolyhedronHype(innerRadius, outerRadius,
1293 tanInnerStereo2, tanOuterStereo2, halfLenZ);
1294}
1295
1296// GetPolyhedron
1297//
1299{
1300 if (fpPolyhedron == nullptr ||
1301 fRebuildPolyhedron ||
1302 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
1303 fpPolyhedron->GetNumberOfRotationSteps())
1304 {
1305 G4AutoLock l(&polyhedronMutex);
1306 delete fpPolyhedron;
1307 fpPolyhedron = CreatePolyhedron();
1308 fRebuildPolyhedron = false;
1309 l.unlock();
1310 }
1311 return fpPolyhedron;
1312}
1313
1314#endif // !defined(G4GEOM_USE_UHYPE) || !defined(G4GEOM_USE_SYS_USOLIDS)
G4TemplateAutoLock< G4Mutex > G4AutoLock
@ JustWarning
@ 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
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
double z() const
Hep3Vector unit() const
double x() const
double y() const
double dot(const Hep3Vector &) const
void set(double x, double y, double z)
double perp() const
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 G4double HyperboloidSurfaceArea(G4double dphi, G4double r0, G4double tanstereo, G4double zmin, G4double zmax)
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
Definition G4Hype.cc:266
G4Hype(const G4String &pName, G4double newInnerRadius, G4double newOuterRadius, G4double newInnerStereo, G4double newOuterStereo, G4double newHalfLenZ)
Definition G4Hype.cc:61
EInside Inside(const G4ThreeVector &p) const override
Definition G4Hype.cc:313
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Hype.cc:275
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
Definition G4Hype.cc:1272
G4GeometryType GetEntityType() const override
Definition G4Hype.cc:1169
G4Hype & operator=(const G4Hype &rhs)
Definition G4Hype.cc:154
G4VisExtent GetExtent() const override
Definition G4Hype.cc:1279
G4VSolid * Clone() const override
Definition G4Hype.cc:1176
void SetOuterStereo(G4double newOSte)
Definition G4Hype.cc:249
G4Polyhedron * CreatePolyhedron() const override
Definition G4Hype.cc:1290
~G4Hype() override
Definition G4Hype.cc:130
G4ThreeVector GetPointOnSurface() const override
Definition G4Hype.cc:1204
G4Polyhedron * GetPolyhedron() const override
Definition G4Hype.cc:1298
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
Definition G4Hype.cc:296
void SetOuterRadius(G4double newORad)
Definition G4Hype.cc:200
G4double GetSurfaceArea() override
Definition G4Hype.cc:1259
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
Definition G4Hype.cc:785
std::ostream & StreamInfo(std::ostream &os) const override
Definition G4Hype.cc:1183
void SetZHalfLength(G4double newHLZ)
Definition G4Hype.cc:215
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
Definition G4Hype.cc:354
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
Definition G4Hype.cc:401
void SetInnerStereo(G4double newISte)
Definition G4Hype.cc:233
void SetInnerRadius(G4double newIRad)
Definition G4Hype.cc:185
G4double GetCubicVolume() override
Definition G4Hype.cc:1245
virtual void AddSolid(const G4Box &)=0
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
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
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
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
T sqr(const T &x)
Definition templates.hh:128
#define DBL_MIN
Definition templates.hh:54
#define DBL_EPSILON
Definition templates.hh:66