Geant4 9.6.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//
27// $Id$
28// $Original: G4Hype.cc,v 1.0 1998/06/09 16:57:50 safai Exp $
29//
30//
31// --------------------------------------------------------------------
32// GEANT 4 class source file
33//
34//
35// G4Hype.cc
36//
37// --------------------------------------------------------------------
38//
39// Authors:
40// Ernesto Lamanna ([email protected]) &
41// Francesco Safai Tehrani ([email protected])
42// Rome, INFN & University of Rome "La Sapienza", 9 June 1998.
43//
44// --------------------------------------------------------------------
45
46#include "G4Hype.hh"
47
48#include "G4VoxelLimits.hh"
49#include "G4AffineTransform.hh"
50#include "G4SolidExtentList.hh"
51#include "G4ClippablePolygon.hh"
52
54
55#include "meshdefs.hh"
56
57#include <cmath>
58
59#include "Randomize.hh"
60
61#include "G4VGraphicsScene.hh"
62#include "G4Polyhedron.hh"
63#include "G4VisExtent.hh"
64#include "G4NURBS.hh"
65#include "G4NURBStube.hh"
66#include "G4NURBScylinder.hh"
67#include "G4NURBStubesector.hh"
68
69using namespace CLHEP;
70
71// Constructor - check parameters, and fills protected data members
73 G4double newInnerRadius,
74 G4double newOuterRadius,
75 G4double newInnerStereo,
76 G4double newOuterStereo,
77 G4double newHalfLenZ)
78 : G4VSolid(pName), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
79{
80 // Check z-len
81 //
82 if (newHalfLenZ<=0)
83 {
84 std::ostringstream message;
85 message << "Invalid Z half-length - " << GetName() << G4endl
86 << " Invalid Z half-length: "
87 << newHalfLenZ/mm << " mm";
88 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
89 FatalErrorInArgument, message);
90 }
91 halfLenZ=newHalfLenZ;
92
93 // Check radii
94 //
95 if (newInnerRadius<0 || newOuterRadius<0)
96 {
97 std::ostringstream message;
98 message << "Invalid radii - " << GetName() << G4endl
99 << " Invalid radii ! Inner radius: "
100 << newInnerRadius/mm << " mm" << G4endl
101 << " Outer radius: "
102 << newOuterRadius/mm << " mm";
103 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
104 FatalErrorInArgument, message);
105 }
106 if (newInnerRadius >= newOuterRadius)
107 {
108 std::ostringstream message;
109 message << "Outer > inner radius - " << GetName() << G4endl
110 << " Invalid radii ! Inner radius: "
111 << newInnerRadius/mm << " mm" << G4endl
112 << " Outer radius: "
113 << newOuterRadius/mm << " mm";
114 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
115 FatalErrorInArgument, message);
116 }
117
118 innerRadius=newInnerRadius;
119 outerRadius=newOuterRadius;
120
123
124 SetInnerStereo( newInnerStereo );
125 SetOuterStereo( newOuterStereo );
126}
127
128
129//
130// Fake default constructor - sets only member data and allocates memory
131// for usage restricted to object persistency.
132//
133G4Hype::G4Hype( __void__& a )
134 : G4VSolid(a), innerRadius(0.), outerRadius(0.), halfLenZ(0.), innerStereo(0.),
135 outerStereo(0.), tanInnerStereo(0.), tanOuterStereo(0.), tanInnerStereo2(0.),
136 tanOuterStereo2(0.), innerRadius2(0.), outerRadius2(0.), endInnerRadius2(0.),
137 endOuterRadius2(0.), endInnerRadius(0.), endOuterRadius(0.),
138 fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
139{
140}
141
142
143//
144// Destructor
145//
147{
148 delete fpPolyhedron;
149}
150
151
152//
153// Copy constructor
154//
156 : G4VSolid(rhs), innerRadius(rhs.innerRadius),
157 outerRadius(rhs.outerRadius), halfLenZ(rhs.halfLenZ),
158 innerStereo(rhs.innerStereo), outerStereo(rhs.outerStereo),
159 tanInnerStereo(rhs.tanInnerStereo), tanOuterStereo(rhs.tanOuterStereo),
160 tanInnerStereo2(rhs.tanInnerStereo2), tanOuterStereo2(rhs.tanOuterStereo2),
161 innerRadius2(rhs.innerRadius2), outerRadius2(rhs.outerRadius2),
162 endInnerRadius2(rhs.endInnerRadius2), endOuterRadius2(rhs.endOuterRadius2),
163 endInnerRadius(rhs.endInnerRadius), endOuterRadius(rhs.endOuterRadius),
164 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
165 fpPolyhedron(0)
166{
167}
168
169
170//
171// Assignment operator
172//
174{
175 // Check assignment to self
176 //
177 if (this == &rhs) { return *this; }
178
179 // Copy base class data
180 //
182
183 // Copy data
184 //
186 halfLenZ = rhs.halfLenZ;
193 fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
194 fpPolyhedron = 0;
195
196 return *this;
197}
198
199
200//
201// Dispatch to parameterisation for replication mechanism dimension
202// computation & modification.
203//
205 const G4int n,
206 const G4VPhysicalVolume* pRep)
207{
208 p->ComputeDimensions(*this,n,pRep);
209}
210
211
212//
213// CalculateExtent
214//
216 const G4VoxelLimits &voxelLimit,
217 const G4AffineTransform &transform,
218 G4double &min, G4double &max ) const
219{
220 G4SolidExtentList extentList( axis, voxelLimit );
221
222 //
223 // Choose phi size of our segment(s) based on constants as
224 // defined in meshdefs.hh
225 //
226 G4int numPhi = kMaxMeshSections;
227 G4double sigPhi = twopi/numPhi;
228 G4double rFudge = 1.0/std::cos(0.5*sigPhi);
229
230 //
231 // We work around in phi building polygons along the way.
232 // As a reasonable compromise between accuracy and
233 // complexity (=cpu time), the following facets are chosen:
234 //
235 // 1. If outerRadius/endOuterRadius > 0.95, approximate
236 // the outer surface as a cylinder, and use one
237 // rectangular polygon (0-1) to build its mesh.
238 //
239 // Otherwise, use two trapazoidal polygons that
240 // meet at z = 0 (0-4-1)
241 //
242 // 2. If there is no inner surface, then use one
243 // polygon for each entire endcap. (0) and (1)
244 //
245 // Otherwise, use a trapazoidal polygon for each
246 // phi segment of each endcap. (0-2) and (1-3)
247 //
248 // 3. For the inner surface, if innerRadius/endInnerRadius > 0.95,
249 // approximate the inner surface as a cylinder of
250 // radius innerRadius and use one rectangular polygon
251 // to build each phi segment of its mesh. (2-3)
252 //
253 // Otherwise, use one rectangular polygon centered
254 // at z = 0 (5-6) and two connecting trapazoidal polygons
255 // for each phi segment (2-5) and (3-6).
256 //
257
258 G4bool splitOuter = (outerRadius/endOuterRadius < 0.95);
259 G4bool splitInner = 0;
260 if (InnerSurfaceExists())
261 {
262 splitInner = (innerRadius/endInnerRadius < 0.95);
263 }
264
265 //
266 // Vertex assignments (v and w arrays)
267 // [0] and [1] are mandatory
268 // the rest are optional
269 //
270 // + -
271 // [0]------[4]------[1] <--- outer radius
272 // | |
273 // | |
274 // [2]---[5]---[6]---[3] <--- inner radius
275 //
276
277
278 G4ClippablePolygon endPoly1, endPoly2;
279
280 G4double phi = 0,
281 cosPhi = std::cos(phi),
282 sinPhi = std::sin(phi);
283 G4ThreeVector v0( rFudge*endOuterRadius*cosPhi,
284 rFudge*endOuterRadius*sinPhi,
285 +halfLenZ ),
286 v1( rFudge*endOuterRadius*cosPhi,
287 rFudge*endOuterRadius*sinPhi,
288 -halfLenZ ),
289 v2, v3, v4, v5, v6,
290 w0, w1, w2, w3, w4, w5, w6;
291 transform.ApplyPointTransform( v0 );
292 transform.ApplyPointTransform( v1 );
293
294 G4double zInnerSplit=0.;
295 if (InnerSurfaceExists())
296 {
297 if (splitInner)
298 {
299 v2 = transform.TransformPoint(
301 endInnerRadius*sinPhi, +halfLenZ ) );
302 v3 = transform.TransformPoint(
304 endInnerRadius*sinPhi, -halfLenZ ) );
305 //
306 // Find intersection of line normal to inner
307 // surface at z = halfLenZ and line r=innerRadius
308 //
311
312 zInnerSplit = halfLenZ + (innerRadius - endInnerRadius)*zn/rn;
313
314 //
315 // Build associated vertices
316 //
317 v5 = transform.TransformPoint(
319 innerRadius*sinPhi, +zInnerSplit ) );
320 v6 = transform.TransformPoint(
322 innerRadius*sinPhi, -zInnerSplit ) );
323 }
324 else
325 {
326 v2 = transform.TransformPoint(
328 innerRadius*sinPhi, +halfLenZ ) );
329 v3 = transform.TransformPoint(
331 innerRadius*sinPhi, -halfLenZ ) );
332 }
333 }
334
335 if (splitOuter)
336 {
337 v4 = transform.TransformPoint(
338 G4ThreeVector( rFudge*outerRadius*cosPhi,
339 rFudge*outerRadius*sinPhi, 0 ) );
340 }
341
342 //
343 // Loop over phi segments
344 //
345 do
346 {
347 phi += sigPhi;
348 if (numPhi == 1) phi = 0; // Try to avoid roundoff
349 cosPhi = std::cos(phi),
350 sinPhi = std::sin(phi);
351
352 G4double r(rFudge*endOuterRadius);
353 w0 = G4ThreeVector( r*cosPhi, r*sinPhi, +halfLenZ );
354 w1 = G4ThreeVector( r*cosPhi, r*sinPhi, -halfLenZ );
355 transform.ApplyPointTransform( w0 );
356 transform.ApplyPointTransform( w1 );
357
358 //
359 // Outer hyperbolic surface
360 //
361 if (splitOuter)
362 {
363 r = rFudge*outerRadius;
364 w4 = G4ThreeVector( r*cosPhi, r*sinPhi, 0 );
365 transform.ApplyPointTransform( w4 );
366
367 AddPolyToExtent( v0, v4, w4, w0, voxelLimit, axis, extentList );
368 AddPolyToExtent( v4, v1, w1, w4, voxelLimit, axis, extentList );
369 }
370 else
371 {
372 AddPolyToExtent( v0, v1, w1, w0, voxelLimit, axis, extentList );
373 }
374
375 if (InnerSurfaceExists())
376 {
377 //
378 // Inner hyperbolic surface
379 //
380 if (splitInner)
381 {
382 w2 = G4ThreeVector( endInnerRadius*cosPhi,
383 endInnerRadius*sinPhi, +halfLenZ );
384 w3 = G4ThreeVector( endInnerRadius*cosPhi,
385 endInnerRadius*sinPhi, -halfLenZ );
386 transform.ApplyPointTransform( w2 );
387 transform.ApplyPointTransform( w3 );
388
389 w5 = G4ThreeVector( innerRadius*cosPhi,
390 innerRadius*sinPhi, +zInnerSplit );
391 w6 = G4ThreeVector( innerRadius*cosPhi,
392 innerRadius*sinPhi, -zInnerSplit );
393 transform.ApplyPointTransform( w5 );
394 transform.ApplyPointTransform( w6 );
395 AddPolyToExtent( v3, v6, w6, w3, voxelLimit, axis, extentList );
396 AddPolyToExtent( v6, v5, w5, w6, voxelLimit, axis, extentList );
397 AddPolyToExtent( v5, v2, w2, w5, voxelLimit, axis, extentList );
398 }
399 else
400 {
401 w2 = G4ThreeVector( innerRadius*cosPhi,
402 innerRadius*sinPhi, +halfLenZ );
403 w3 = G4ThreeVector( innerRadius*cosPhi,
404 innerRadius*sinPhi, -halfLenZ );
405 transform.ApplyPointTransform( w2 );
406 transform.ApplyPointTransform( w3 );
407
408 AddPolyToExtent( v3, v2, w2, w3, voxelLimit, axis, extentList );
409 }
410
411 //
412 // Endplate segments
413 //
414 AddPolyToExtent( v1, v3, w3, w1, voxelLimit, axis, extentList );
415 AddPolyToExtent( v2, v0, w0, w2, voxelLimit, axis, extentList );
416 }
417 else
418 {
419 //
420 // Continue building endplate polygons
421 //
422 endPoly1.AddVertexInOrder( v0 );
423 endPoly2.AddVertexInOrder( v1 );
424 }
425
426 //
427 // Next phi segments
428 //
429 v0 = w0;
430 v1 = w1;
431 if (InnerSurfaceExists())
432 {
433 v2 = w2;
434 v3 = w3;
435 if (splitInner)
436 {
437 v5 = w5;
438 v6 = w6;
439 }
440 }
441 if (splitOuter) v4 = w4;
442
443 } while( --numPhi > 0 );
444
445
446 //
447 // Don't forget about the endplate polygons, if
448 // we use them
449 //
450 if (!InnerSurfaceExists())
451 {
452 if (endPoly1.PartialClip( voxelLimit, axis ))
453 {
454 static const G4ThreeVector normal(0,0,+1);
455 endPoly1.SetNormal( transform.TransformAxis(normal) );
456 extentList.AddSurface( endPoly1 );
457 }
458
459 if (endPoly2.PartialClip( voxelLimit, axis ))
460 {
461 static const G4ThreeVector normal(0,0,-1);
462 endPoly2.SetNormal( transform.TransformAxis(normal) );
463 extentList.AddSurface( endPoly2 );
464 }
465 }
466
467 //
468 // Return min/max value
469 //
470 return extentList.GetExtent( min, max );
471}
472
473
474//
475// AddPolyToExtent (static)
476//
477// Utility function for CalculateExtent
478//
480 const G4ThreeVector &v1,
481 const G4ThreeVector &w1,
482 const G4ThreeVector &w0,
483 const G4VoxelLimits &voxelLimit,
484 const EAxis axis,
485 G4SolidExtentList &extentList )
486{
487 G4ClippablePolygon phiPoly;
488
489 phiPoly.AddVertexInOrder( v0 );
490 phiPoly.AddVertexInOrder( v1 );
491 phiPoly.AddVertexInOrder( w1 );
492 phiPoly.AddVertexInOrder( w0 );
493
494 if (phiPoly.PartialClip( voxelLimit, axis ))
495 {
496 phiPoly.SetNormal( (v1-v0).cross(w0-v0).unit() );
497 extentList.AddSurface( phiPoly );
498 }
499}
500
501
502//
503// Decides whether point is inside,outside or on the surface
504//
506{
507 static const G4double halfTol = 0.5*kCarTolerance;
508
509 //
510 // Check z extents: are we outside?
511 //
512 const G4double absZ(std::fabs(p.z()));
513 if (absZ > halfLenZ + halfTol) return kOutside;
514
515 //
516 // Check outer radius
517 //
518 const G4double oRad2(HypeOuterRadius2(absZ));
519 const G4double xR2( p.x()*p.x()+p.y()*p.y() );
520
521 if (xR2 > oRad2 + kCarTolerance*endOuterRadius) return kOutside;
522
523 if (xR2 > oRad2 - kCarTolerance*endOuterRadius) return kSurface;
524
525 if (InnerSurfaceExists())
526 {
527 //
528 // Check inner radius
529 //
530 const G4double iRad2(HypeInnerRadius2(absZ));
531
532 if (xR2 < iRad2 - kCarTolerance*endInnerRadius) return kOutside;
533
534 if (xR2 < iRad2 + kCarTolerance*endInnerRadius) return kSurface;
535 }
536
537 //
538 // We are inside in radius, now check endplate surface
539 //
540 if (absZ > halfLenZ - halfTol) return kSurface;
541
542 return kInside;
543}
544
545
546
547//
548// return the normal unit vector to the Hyperbolical Surface at a point
549// p on (or nearly on) the surface
550//
552{
553 //
554 // Which of the three or four surfaces are we closest to?
555 //
556 const G4double absZ(std::fabs(p.z()));
557 const G4double distZ(absZ - halfLenZ);
558 const G4double dist2Z(distZ*distZ);
559
560 const G4double xR2( p.x()*p.x()+p.y()*p.y() );
561 const G4double dist2Outer( std::fabs(xR2 - HypeOuterRadius2(absZ)) );
562
563 if (InnerSurfaceExists())
564 {
565 //
566 // Has inner surface: is this closest?
567 //
568 const G4double dist2Inner( std::fabs(xR2 - HypeInnerRadius2(absZ)) );
569 if (dist2Inner < dist2Z && dist2Inner < dist2Outer)
570 return G4ThreeVector( -p.x(), -p.y(), p.z()*tanInnerStereo2 ).unit();
571 }
572
573 //
574 // Do the "endcaps" win?
575 //
576 if (dist2Z < dist2Outer)
577 return G4ThreeVector( 0.0, 0.0, p.z() < 0 ? -1.0 : 1.0 );
578
579
580 //
581 // Outer surface wins
582 //
583 return G4ThreeVector( p.x(), p.y(), -p.z()*tanOuterStereo2 ).unit();
584}
585
586
587//
588// Calculate distance to shape from outside, along normalised vector
589// - return kInfinity if no intersection,
590// or intersection distance <= tolerance
591//
592// Calculating the intersection of a line with the surfaces
593// is fairly straight forward. The difficult problem is dealing
594// with the intersections of the surfaces in a consistent manner,
595// and this accounts for the complicated logic.
596//
598 const G4ThreeVector& v ) const
599{
600 static const G4double halfTol = 0.5*kCarTolerance;
601
602 //
603 // Quick test. Beware! This assumes v is a unit vector!
604 //
605 if (std::fabs(p.x()*v.y() - p.y()*v.x()) > endOuterRadius+kCarTolerance)
606 return kInfinity;
607
608 //
609 // Take advantage of z symmetry, and reflect throught the
610 // z=0 plane so that pz is always positive
611 //
612 G4double pz(p.z()), vz(v.z());
613 if (pz < 0)
614 {
615 pz = -pz;
616 vz = -vz;
617 }
618
619 //
620 // We must be very careful if we don't want to
621 // create subtle leaks at the edges where the
622 // hyperbolic surfaces connect to the endplate.
623 // The only reliable way to do so is to make sure
624 // that the decision as to when a track passes
625 // over the edge of one surface is exactly the
626 // same decision as to when a track passes into the
627 // other surface. By "exact", we don't mean algebraicly
628 // exact, but we mean the same machine instructions
629 // should be used.
630 //
631 G4bool couldMissOuter(true),
632 couldMissInner(true),
633 cantMissInnerCylinder(false);
634
635 //
636 // Check endplate intersection
637 //
638 G4double sigz = pz-halfLenZ;
639
640 if (sigz > -halfTol) // equivalent to: if (pz > halfLenZ - halfTol)
641 {
642 //
643 // We start in front of the endplate (within roundoff)
644 // Correct direction to intersect endplate?
645 //
646 if (vz >= 0)
647 {
648 //
649 // Nope. As long as we are far enough away, we
650 // can't intersect anything
651 //
652 if (sigz > 0) return kInfinity;
653
654 //
655 // Otherwise, we may still hit a hyperbolic surface
656 // if the point is on the hyperbolic surface (within tolerance)
657 //
658 G4double pr2 = p.x()*p.x() + p.y()*p.y();
660 return kInfinity;
661
662 if (InnerSurfaceExists())
663 {
665 return kInfinity;
668 return kInfinity;
669 }
670 else
671 {
673 return kInfinity;
674 }
675 }
676 else
677 {
678 //
679 // Where do we intersect at z = halfLenZ?
680 //
681 G4double q = -sigz/vz;
682 G4double xi = p.x() + q*v.x(),
683 yi = p.y() + q*v.y();
684
685 //
686 // Is this on the endplate? If so, return s, unless
687 // we are on the tolerant surface, in which case return 0
688 //
689 G4double pr2 = xi*xi + yi*yi;
690 if (pr2 <= endOuterRadius2)
691 {
692 if (InnerSurfaceExists())
693 {
694 if (pr2 >= endInnerRadius2) return (sigz < halfTol) ? 0 : q;
695 //
696 // This test is sufficient to ensure that the
697 // trajectory cannot miss the inner hyperbolic surface
698 // for z > 0, if the normal is correct.
699 //
700 G4double dot1 = (xi*v.x() + yi*v.y())*endInnerRadius/std::sqrt(pr2);
701 couldMissInner = (dot1 - halfLenZ*tanInnerStereo2*vz <= 0);
702
703 if (pr2 > endInnerRadius2*(1 - 2*DBL_EPSILON) )
704 {
705 //
706 // There is a potential leak if the inner
707 // surface is a cylinder
708 //
709 if ( (innerStereo < DBL_MIN)
710 && ((std::fabs(v.x()) > DBL_MIN) || (std::fabs(v.y()) > DBL_MIN)) )
711 cantMissInnerCylinder = true;
712 }
713 }
714 else
715 {
716 return (sigz < halfTol) ? 0 : q;
717 }
718 }
719 else
720 {
721 G4double dotR( xi*v.x() + yi*v.y() );
722 if (dotR >= 0)
723 {
724 //
725 // Otherwise, if we are traveling outwards, we know
726 // we must miss the hyperbolic surfaces also, so
727 // we need not bother checking
728 //
729 return kInfinity;
730 }
731 else
732 {
733 //
734 // This test is sufficient to ensure that the
735 // trajectory cannot miss the outer hyperbolic surface
736 // for z > 0, if the normal is correct.
737 //
738 G4double dot1 = dotR*endOuterRadius/std::sqrt(pr2);
739 couldMissOuter = (dot1 - halfLenZ*tanOuterStereo2*vz>= 0);
740 }
741 }
742 }
743 }
744
745 //
746 // Check intersection with outer hyperbolic surface, save
747 // distance to valid intersection into "best".
748 //
749 G4double best = kInfinity;
750
751 G4double q[2];
753
754 if (n > 0)
755 {
756 //
757 // Potential intersection: is p on this surface?
758 //
759 if (pz < halfLenZ+halfTol)
760 {
761 G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeOuterRadius2(pz);
762 if (std::fabs(dr2) < kCarTolerance*endOuterRadius)
763 {
764 //
765 // Sure, but make sure we're traveling inwards at
766 // this point
767 //
768 if (p.x()*v.x() + p.y()*v.y() - pz*tanOuterStereo2*vz < 0)
769 return 0;
770 }
771 }
772
773 //
774 // We are now certain that p is not on the tolerant surface.
775 // Accept only position distance q
776 //
777 G4int i;
778 for( i=0; i<n; i++ )
779 {
780 if (q[i] >= 0)
781 {
782 //
783 // Check to make sure this intersection point is
784 // on the surface, but only do so if we haven't
785 // checked the endplate intersection already
786 //
787 G4double zi = pz + q[i]*vz;
788
789 if (zi < -halfLenZ) continue;
790 if (zi > +halfLenZ && couldMissOuter) continue;
791
792 //
793 // Check normal
794 //
795 G4double xi = p.x() + q[i]*v.x(),
796 yi = p.y() + q[i]*v.y();
797
798 if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz > 0) continue;
799
800 best = q[i];
801 break;
802 }
803 }
804 }
805
806 if (!InnerSurfaceExists()) return best;
807
808 //
809 // Check intersection with inner hyperbolic surface
810 //
812 if (n == 0)
813 {
814 if (cantMissInnerCylinder) return (sigz < halfTol) ? 0 : -sigz/vz;
815
816 return best;
817 }
818
819 //
820 // P on this surface?
821 //
822 if (pz < halfLenZ+halfTol)
823 {
824 G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeInnerRadius2(pz);
825 if (std::fabs(dr2) < kCarTolerance*endInnerRadius)
826 {
827 //
828 // Sure, but make sure we're traveling outwards at
829 // this point
830 //
831 if (p.x()*v.x() + p.y()*v.y() - pz*tanInnerStereo2*vz > 0) return 0;
832 }
833 }
834
835 //
836 // No, so only positive q is valid. Search for a valid intersection
837 // that is closer than the outer intersection (if it exists)
838 //
839 G4int i;
840 for( i=0; i<n; i++ )
841 {
842 if (q[i] > best) break;
843 if (q[i] >= 0)
844 {
845 //
846 // Check to make sure this intersection point is
847 // on the surface, but only do so if we haven't
848 // checked the endplate intersection already
849 //
850 G4double zi = pz + q[i]*vz;
851
852 if (zi < -halfLenZ) continue;
853 if (zi > +halfLenZ && couldMissInner) continue;
854
855 //
856 // Check normal
857 //
858 G4double xi = p.x() + q[i]*v.x(),
859 yi = p.y() + q[i]*v.y();
860
861 if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz < 0) continue;
862
863 best = q[i];
864 break;
865 }
866 }
867
868 //
869 // Done
870 //
871 return best;
872}
873
874
875//
876// Calculate distance to shape from outside, along perpendicular direction
877// (if one exists). May be an underestimate.
878//
879// There are five (r,z) regions:
880// 1. a point that is beyond the endcap but within the
881// endcap radii
882// 2. a point with r > outer endcap radius and with
883// a z position that is beyond the cone formed by the
884// normal of the outer hyperbolic surface at the
885// edge at which it meets the endcap.
886// 3. a point that is outside the outer surface and not in (1 or 2)
887// 4. a point that is inside the inner surface and not in (5)
888// 5. a point with radius < inner endcap radius and
889// with a z position beyond the cone formed by the
890// normal of the inner hyperbolic surface at the
891// edge at which it meets the endcap.
892// (regions 4 and 5 only exist if there is an inner surface)
893//
895{
896 static const G4double halfTol(0.5*kCarTolerance);
897
898 G4double absZ(std::fabs(p.z()));
899
900 //
901 // Check region
902 //
903 G4double r2 = p.x()*p.x() + p.y()*p.y();
904 G4double r = std::sqrt(r2);
905
906 G4double sigz = absZ - halfLenZ;
907
908 if (r < endOuterRadius)
909 {
910 if (sigz > -halfTol)
911 {
912 if (InnerSurfaceExists())
913 {
914 if (r > endInnerRadius)
915 return sigz < halfTol ? 0 : sigz; // Region 1
916
917 G4double dr = endInnerRadius - r;
918 if (sigz > dr*tanInnerStereo2)
919 {
920 //
921 // In region 5
922 //
923 G4double answer = std::sqrt( dr*dr + sigz*sigz );
924 return answer < halfTol ? 0 : answer;
925 }
926 }
927 else
928 {
929 //
930 // In region 1 (no inner surface)
931 //
932 return sigz < halfTol ? 0 : sigz;
933 }
934 }
935 }
936 else
937 {
938 G4double dr = r - endOuterRadius;
939 if (sigz > -dr*tanOuterStereo2)
940 {
941 //
942 // In region 2
943 //
944 G4double answer = std::sqrt( dr*dr + sigz*sigz );
945 return answer < halfTol ? 0 : answer;
946 }
947 }
948
949 if (InnerSurfaceExists())
950 {
952 {
953 //
954 // In region 4
955 //
957 return answer < halfTol ? 0 : answer;
958 }
959 }
960
961 //
962 // We are left by elimination with region 3
963 //
965 return answer < halfTol ? 0 : answer;
966}
967
968
969//
970// Calculate distance to surface of shape from `inside', allowing for tolerance
971//
972// The situation here is much simplier than DistanceToIn(p,v). For
973// example, there is no need to even check whether an intersection
974// point is inside the boundary of a surface, as long as all surfaces
975// are checked and the smallest distance is used.
976//
978 const G4bool calcNorm,
979 G4bool *validNorm, G4ThreeVector *norm ) const
980{
981 static const G4double halfTol = 0.5*kCarTolerance;
982
983
984 static const G4ThreeVector normEnd1(0.0,0.0,+1.0);
985 static const G4ThreeVector normEnd2(0.0,0.0,-1.0);
986
987 //
988 // Keep track of closest surface
989 //
990 G4double sBest; // distance to
991 const G4ThreeVector *nBest; // normal vector
992 G4bool vBest; // whether "valid"
993
994 //
995 // Check endplate, taking advantage of symmetry.
996 // Note that the endcap is the only surface which
997 // has a "valid" normal, i.e. is a surface of which
998 // the entire solid is behind.
999 //
1000 G4double pz(p.z()), vz(v.z());
1001 if (vz < 0)
1002 {
1003 pz = -pz;
1004 vz = -vz;
1005 nBest = &normEnd2;
1006 }
1007 else
1008 nBest = &normEnd1;
1009
1010 //
1011 // Possible intercept. Are we on the surface?
1012 //
1013 if (pz > halfLenZ-halfTol)
1014 {
1015 if (calcNorm) { *norm = *nBest; *validNorm = true; }
1016 return 0;
1017 }
1018
1019 //
1020 // Nope. Get distance. Beware of zero vz.
1021 //
1022 sBest = (vz > DBL_MIN) ? (halfLenZ - pz)/vz : kInfinity;
1023 vBest = true;
1024
1025 //
1026 // Check outer surface
1027 //
1028 G4double r2 = p.x()*p.x() + p.y()*p.y();
1029
1030 G4double q[2];
1032
1033 G4ThreeVector norm1, norm2;
1034
1035 if (n > 0)
1036 {
1037 //
1038 // We hit somewhere. Are we on the surface?
1039 //
1040 G4double dr2 = r2 - HypeOuterRadius2(pz);
1041 if (std::fabs(dr2) < endOuterRadius*kCarTolerance)
1042 {
1043 G4ThreeVector normHere( p.x(), p.y(), -p.z()*tanOuterStereo2 );
1044 //
1045 // Sure. But are we going the right way?
1046 //
1047 if (normHere.dot(v) > 0)
1048 {
1049 if (calcNorm) { *norm = normHere.unit(); *validNorm = false; }
1050 return 0;
1051 }
1052 }
1053
1054 //
1055 // Nope. Check closest positive intercept.
1056 //
1057 G4int i;
1058 for( i=0; i<n; i++ )
1059 {
1060 if (q[i] > sBest) break;
1061 if (q[i] > 0)
1062 {
1063 //
1064 // Make sure normal is correct (that this
1065 // solution is an outgoing solution)
1066 //
1067 G4ThreeVector pk(p+q[i]*v);
1068 norm1 = G4ThreeVector( pk.x(), pk.y(), -pk.z()*tanOuterStereo2 );
1069 if (norm1.dot(v) > 0)
1070 {
1071 sBest = q[i];
1072 nBest = &norm1;
1073 vBest = false;
1074 break;
1075 }
1076 }
1077 }
1078 }
1079
1080 if (InnerSurfaceExists())
1081 {
1082 //
1083 // Check inner surface
1084 //
1086 if (n > 0)
1087 {
1088 //
1089 // On surface?
1090 //
1091 G4double dr2 = r2 - HypeInnerRadius2(pz);
1092 if (std::fabs(dr2) < endInnerRadius*kCarTolerance)
1093 {
1094 G4ThreeVector normHere( -p.x(), -p.y(), p.z()*tanInnerStereo2 );
1095 if (normHere.dot(v) > 0)
1096 {
1097 if (calcNorm)
1098 {
1099 *norm = normHere.unit();
1100 *validNorm = false;
1101 }
1102 return 0;
1103 }
1104 }
1105
1106 //
1107 // Check closest positive
1108 //
1109 G4int i;
1110 for( i=0; i<n; i++ )
1111 {
1112 if (q[i] > sBest) break;
1113 if (q[i] > 0)
1114 {
1115 G4ThreeVector pk(p+q[i]*v);
1116 norm2 = G4ThreeVector( -pk.x(), -pk.y(), pk.z()*tanInnerStereo2 );
1117 if (norm2.dot(v) > 0)
1118 {
1119 sBest = q[i];
1120 nBest = &norm2;
1121 vBest = false;
1122 break;
1123 }
1124 }
1125 }
1126 }
1127 }
1128
1129 //
1130 // Done!
1131 //
1132 if (calcNorm)
1133 {
1134 *validNorm = vBest;
1135
1136 if (nBest == &norm1 || nBest == &norm2)
1137 *norm = nBest->unit();
1138 else
1139 *norm = *nBest;
1140 }
1141
1142 return sBest;
1143}
1144
1145
1146//
1147// Calculate distance (<=actual) to closest surface of shape from inside
1148//
1149// May be an underestimate
1150//
1152{
1153 //
1154 // Try each surface and remember the closest
1155 //
1156 G4double absZ(std::fabs(p.z()));
1157 G4double r(p.perp());
1158
1159 G4double sBest = halfLenZ - absZ;
1160
1161 G4double tryOuter = ApproxDistInside( r, absZ, outerRadius, tanOuterStereo2 );
1162 if (tryOuter < sBest)
1163 sBest = tryOuter;
1164
1165 if (InnerSurfaceExists())
1166 {
1168 if (tryInner < sBest) sBest = tryInner;
1169 }
1170
1171 return sBest < 0.5*kCarTolerance ? 0 : sBest;
1172}
1173
1174
1175//
1176// IntersectHype (static)
1177//
1178// Decide if and where a line intersects with a hyperbolic
1179// surface (of infinite extent)
1180//
1181// Arguments:
1182// p - (in) Point on trajectory
1183// v - (in) Vector along trajectory
1184// r2 - (in) Square of radius at z = 0
1185// tan2phi - (in) std::tan(phi)**2
1186// q - (out) Up to two points of intersection, where the
1187// intersection point is p + q*v, and if there are
1188// two intersections, q[0] < q[1]. May be negative.
1189// Returns:
1190// The number of intersections. If 0, the trajectory misses.
1191//
1192//
1193// Equation of a line:
1194//
1195// x = x0 + q*tx y = y0 + q*ty z = z0 + q*tz
1196//
1197// Equation of a hyperbolic surface:
1198//
1199// x**2 + y**2 = r**2 + (z*tanPhi)**2
1200//
1201// Solution is quadratic:
1202//
1203// a*q**2 + b*q + c = 0
1204//
1205// where:
1206//
1207// a = tx**2 + ty**2 - (tz*tanPhi)**2
1208//
1209// b = 2*( x0*tx + y0*ty - z0*tz*tanPhi**2 )
1210//
1211// c = x0**2 + y0**2 - r**2 - (z0*tanPhi)**2
1212//
1213//
1215 G4double r2, G4double tan2Phi, G4double ss[2] )
1216{
1217 G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
1218 G4double tx = v.x(), ty = v.y(), tz = v.z();
1219
1220 G4double a = tx*tx + ty*ty - tz*tz*tan2Phi;
1221 G4double b = 2*( x0*tx + y0*ty - z0*tz*tan2Phi );
1222 G4double c = x0*x0 + y0*y0 - r2 - z0*z0*tan2Phi;
1223
1224 if (std::fabs(a) < DBL_MIN)
1225 {
1226 //
1227 // The trajectory is parallel to the asympotic limit of
1228 // the surface: single solution
1229 //
1230 if (std::fabs(b) < DBL_MIN) return 0; // Unless we travel through exact center
1231
1232 ss[0] = c/b;
1233 return 1;
1234 }
1235
1236
1237 G4double radical = b*b - 4*a*c;
1238
1239 if (radical < -DBL_MIN) return 0; // No solution
1240
1241 if (radical < DBL_MIN)
1242 {
1243 //
1244 // Grazes surface
1245 //
1246 ss[0] = -b/a/2.0;
1247 return 1;
1248 }
1249
1250 radical = std::sqrt(radical);
1251
1252 G4double q = -0.5*( b + (b < 0 ? -radical : +radical) );
1253 G4double sa = q/a;
1254 G4double sb = c/q;
1255 if (sa < sb) { ss[0] = sa; ss[1] = sb; } else { ss[0] = sb; ss[1] = sa; }
1256 return 2;
1257}
1258
1259
1260//
1261// ApproxDistOutside (static)
1262//
1263// Find the approximate distance of a point outside
1264// (greater radius) of a hyperbolic surface. The distance
1265// must be an underestimate. It will also be nice (although
1266// not necesary) that the estimate is always finite no
1267// matter how close the point is.
1268//
1269// Our hyperbola approaches the asymptotic limit at z = +/- infinity
1270// to the lines r = z*tanPhi. We call these lines the
1271// asymptotic limit line.
1272//
1273// We need the distance of the 2d point p(r,z) to the
1274// hyperbola r**2 = r0**2 + (z*tanPhi)**2. Find two
1275// points that bracket the true normal and use the
1276// distance to the line that connects these two points.
1277// The first such point is z=p.z. The second point is
1278// the z position on the asymptotic limit line that
1279// contains the normal on the line through the point p.
1280//
1282 G4double r0, G4double tanPhi )
1283{
1284 if (tanPhi < DBL_MIN) return pr-r0;
1285
1286 G4double tan2Phi = tanPhi*tanPhi;
1287
1288 //
1289 // First point
1290 //
1291 G4double z1 = pz;
1292 G4double r1 = std::sqrt( r0*r0 + z1*z1*tan2Phi );
1293
1294 //
1295 // Second point
1296 //
1297 G4double z2 = (pr*tanPhi + pz)/(1 + tan2Phi);
1298 G4double r2 = std::sqrt( r0*r0 + z2*z2*tan2Phi );
1299
1300 //
1301 // Line between them
1302 //
1303 G4double dr = r2-r1;
1304 G4double dz = z2-z1;
1305
1306 G4double len = std::sqrt(dr*dr + dz*dz);
1307 if (len < DBL_MIN)
1308 {
1309 //
1310 // The two points are the same?? I guess we
1311 // must have really bracketed the normal
1312 //
1313 dr = pr-r1;
1314 dz = pz-z1;
1315 return std::sqrt( dr*dr + dz*dz );
1316 }
1317
1318 //
1319 // Distance
1320 //
1321 return std::fabs((pr-r1)*dz - (pz-z1)*dr)/len;
1322}
1323
1324//
1325// ApproxDistInside (static)
1326//
1327// Find the approximate distance of a point inside
1328// of a hyperbolic surface. The distance
1329// must be an underestimate. It will also be nice (although
1330// not necesary) that the estimate is always finite no
1331// matter how close the point is.
1332//
1333// This estimate uses the distance to a line tangent to
1334// the hyperbolic function. The point of tangent is chosen
1335// by the z position point
1336//
1337// Assumes pr and pz are positive
1338//
1340 G4double r0, G4double tan2Phi )
1341{
1342 if (tan2Phi < DBL_MIN) return r0 - pr;
1343
1344 //
1345 // Corresponding position and normal on hyperbolic
1346 //
1347 G4double rh = std::sqrt( r0*r0 + pz*pz*tan2Phi );
1348
1349 G4double dr = -rh;
1350 G4double dz = pz*tan2Phi;
1351 G4double len = std::sqrt(dr*dr + dz*dz);
1352
1353 //
1354 // Answer
1355 //
1356 return std::fabs((pr-rh)*dr)/len;
1357}
1358
1359
1360//
1361// GetEntityType
1362//
1364{
1365 return G4String("G4Hype");
1366}
1367
1368
1369//
1370// Clone
1371//
1373{
1374 return new G4Hype(*this);
1375}
1376
1377
1378//
1379// GetCubicVolume
1380//
1382{
1383 if(fCubicVolume != 0.) {;}
1384 else { fCubicVolume = G4VSolid::GetCubicVolume(); }
1385 return fCubicVolume;
1386}
1387
1388
1389//
1390// GetSurfaceArea
1391//
1393{
1394 if(fSurfaceArea != 0.) {;}
1395 else { fSurfaceArea = G4VSolid::GetSurfaceArea(); }
1396 return fSurfaceArea;
1397}
1398
1399
1400//
1401// Stream object contents to an output stream
1402//
1403std::ostream& G4Hype::StreamInfo(std::ostream& os) const
1404{
1405 G4int oldprc = os.precision(16);
1406 os << "-----------------------------------------------------------\n"
1407 << " *** Dump for solid - " << GetName() << " ***\n"
1408 << " ===================================================\n"
1409 << " Solid type: G4Hype\n"
1410 << " Parameters: \n"
1411 << " half length Z: " << halfLenZ/mm << " mm \n"
1412 << " inner radius : " << innerRadius/mm << " mm \n"
1413 << " outer radius : " << outerRadius/mm << " mm \n"
1414 << " inner stereo angle : " << innerStereo/degree << " degrees \n"
1415 << " outer stereo angle : " << outerStereo/degree << " degrees \n"
1416 << "-----------------------------------------------------------\n";
1417 os.precision(oldprc);
1418
1419 return os;
1420}
1421
1422
1423
1424//
1425// GetPointOnSurface
1426//
1428{
1429 G4double xRand, yRand, zRand, r2 , aOne, aTwo, aThree, chose, sinhu;
1430 G4double phi, cosphi, sinphi, rBar2Out, rBar2In, alpha, t, rOut, rIn2, rOut2;
1431
1432 // we use the formula of the area of a surface of revolution to compute
1433 // the areas, using the equation of the hyperbola:
1434 // x^2 + y^2 = (z*tanphi)^2 + r^2
1435
1436 rBar2Out = outerRadius2;
1437 alpha = 2.*pi*rBar2Out*std::cos(outerStereo)/tanOuterStereo;
1439 t = std::log(t+std::sqrt(sqr(t)+1));
1440 aOne = std::fabs(2.*alpha*(std::sinh(2.*t)/4.+t/2.));
1441
1442
1443 rBar2In = innerRadius2;
1444 alpha = 2.*pi*rBar2In*std::cos(innerStereo)/tanInnerStereo;
1446 t = std::log(t+std::sqrt(sqr(t)+1));
1447 aTwo = std::fabs(2.*alpha*(std::sinh(2.*t)/4.+t/2.));
1448
1451
1452 if(outerStereo == 0.) {aOne = std::fabs(2.*pi*outerRadius*2.*halfLenZ);}
1453 if(innerStereo == 0.) {aTwo = std::fabs(2.*pi*innerRadius*2.*halfLenZ);}
1454
1455 phi = RandFlat::shoot(0.,2.*pi);
1456 cosphi = std::cos(phi);
1457 sinphi = std::sin(phi);
1460
1461 chose = RandFlat::shoot(0.,aOne+aTwo+2.*aThree);
1462 if(chose>=0. && chose < aOne)
1463 {
1464 if(outerStereo != 0.)
1465 {
1466 zRand = outerRadius*sinhu/tanOuterStereo;
1467 xRand = std::sqrt(sqr(sinhu)+1)*outerRadius*cosphi;
1468 yRand = std::sqrt(sqr(sinhu)+1)*outerRadius*sinphi;
1469 return G4ThreeVector (xRand, yRand, zRand);
1470 }
1471 else
1472 {
1473 return G4ThreeVector(outerRadius*cosphi,outerRadius*sinphi,
1475 }
1476 }
1477 else if(chose>=aOne && chose<aOne+aTwo)
1478 {
1479 if(innerStereo != 0.)
1480 {
1483 zRand = innerRadius*sinhu/tanInnerStereo;
1484 xRand = std::sqrt(sqr(sinhu)+1)*innerRadius*cosphi;
1485 yRand = std::sqrt(sqr(sinhu)+1)*innerRadius*sinphi;
1486 return G4ThreeVector (xRand, yRand, zRand);
1487 }
1488 else
1489 {
1490 return G4ThreeVector(innerRadius*cosphi,innerRadius*sinphi,
1492 }
1493 }
1494 else if(chose>=aOne+aTwo && chose<aOne+aTwo+aThree)
1495 {
1498 rOut = std::sqrt(rOut2) ;
1499
1500 do {
1501 xRand = RandFlat::shoot(-rOut,rOut) ;
1502 yRand = RandFlat::shoot(-rOut,rOut) ;
1503 r2 = xRand*xRand + yRand*yRand ;
1504 } while ( ! ( r2 >= rIn2 && r2 <= rOut2 ) ) ;
1505
1506 zRand = halfLenZ;
1507 return G4ThreeVector (xRand, yRand, zRand);
1508 }
1509 else
1510 {
1513 rOut = std::sqrt(rOut2) ;
1514
1515 do {
1516 xRand = RandFlat::shoot(-rOut,rOut) ;
1517 yRand = RandFlat::shoot(-rOut,rOut) ;
1518 r2 = xRand*xRand + yRand*yRand ;
1519 } while ( ! ( r2 >= rIn2 && r2 <= rOut2 ) ) ;
1520
1521 zRand = -1.*halfLenZ;
1522 return G4ThreeVector (xRand, yRand, zRand);
1523 }
1524}
1525
1526
1527//
1528// DescribeYourselfTo
1529//
1531{
1532 scene.AddSolid (*this);
1533}
1534
1535
1536//
1537// GetExtent
1538//
1540{
1541 // Define the sides of the box into which the G4Tubs instance would fit.
1542 //
1545 -halfLenZ, halfLenZ );
1546}
1547
1548
1549//
1550// CreatePolyhedron
1551//
1553{
1556}
1557
1558
1559//
1560// GetPolyhedron
1561//
1563{
1564 if (!fpPolyhedron ||
1566 fpPolyhedron->GetNumberOfRotationSteps())
1567 {
1568 delete fpPolyhedron;
1569 fpPolyhedron = CreatePolyhedron();
1570 }
1571 return fpPolyhedron;
1572}
1573
1574
1575//
1576// CreateNURBS
1577//
1579{
1580 // Tube for now!!!
1581 //
1583}
1584
1585
1586//
1587// asinh
1588//
1589G4double G4Hype::asinh(G4double arg)
1590{
1591 return std::log(arg+std::sqrt(sqr(arg)+1));
1592}
@ FatalErrorInArgument
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
double z() const
Hep3Vector unit() const
double x() const
double y() const
double dot(const Hep3Vector &) const
double perp() const
static double shoot()
Definition: RandFlat.cc:59
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
void ApplyPointTransform(G4ThreeVector &vec) const
virtual G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
virtual void AddVertexInOrder(const G4ThreeVector vertex)
void SetNormal(const G4ThreeVector &newNormal)
Definition: G4Hype.hh:67
G4ThreeVector GetPointOnSurface() const
Definition: G4Hype.cc:1427
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Hype.cc:977
G4double endInnerRadius
Definition: G4Hype.hh:184
G4double tanOuterStereo
Definition: G4Hype.hh:177
G4Hype(const G4String &pName, G4double newInnerRadius, G4double newOuterRadius, G4double newInnerStereo, G4double newOuterStereo, G4double newHalfLenZ)
Definition: G4Hype.cc:72
G4Polyhedron * CreatePolyhedron() const
Definition: G4Hype.cc:1552
G4VisExtent GetExtent() const
Definition: G4Hype.cc:1539
static void AddPolyToExtent(const G4ThreeVector &v0, const G4ThreeVector &v1, const G4ThreeVector &w1, const G4ThreeVector &w0, const G4VoxelLimits &voxelLimit, const EAxis axis, G4SolidExtentList &extentList)
Definition: G4Hype.cc:479
G4Polyhedron * GetPolyhedron() const
Definition: G4Hype.cc:1562
G4double tanOuterStereo2
Definition: G4Hype.hh:179
static G4int IntersectHype(const G4ThreeVector &p, const G4ThreeVector &v, G4double r2, G4double tan2Phi, G4double s[2])
Definition: G4Hype.cc:1214
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Hype.cc:551
G4NURBS * CreateNURBS() const
Definition: G4Hype.cc:1578
G4double tanInnerStereo2
Definition: G4Hype.hh:178
G4double halfLenZ
Definition: G4Hype.hh:170
G4Hype & operator=(const G4Hype &rhs)
Definition: G4Hype.cc:173
G4double endInnerRadius2
Definition: G4Hype.hh:182
EInside Inside(const G4ThreeVector &p) const
Definition: G4Hype.cc:505
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Hype.cc:204
G4VSolid * Clone() const
Definition: G4Hype.cc:1372
G4double innerRadius2
Definition: G4Hype.hh:180
static G4double ApproxDistOutside(G4double pr, G4double pz, G4double r0, G4double tanPhi)
Definition: G4Hype.cc:1281
void SetOuterStereo(G4double newOSte)
G4double endOuterRadius
Definition: G4Hype.hh:185
G4double outerStereo
Definition: G4Hype.hh:172
G4double outerRadius
Definition: G4Hype.hh:169
G4double GetCubicVolume()
Definition: G4Hype.cc:1381
G4double tanInnerStereo
Definition: G4Hype.hh:176
G4double endOuterRadius2
Definition: G4Hype.hh:183
G4bool InnerSurfaceExists() const
G4double outerRadius2
Definition: G4Hype.hh:181
G4double HypeOuterRadius2(G4double zVal) const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Hype.cc:597
G4double innerRadius
Definition: G4Hype.hh:168
virtual ~G4Hype()
Definition: G4Hype.cc:146
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Hype.cc:1403
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Hype.cc:215
G4double innerStereo
Definition: G4Hype.hh:171
static G4double ApproxDistInside(G4double pr, G4double pz, G4double r0, G4double tan2Phi)
Definition: G4Hype.cc:1339
G4GeometryType GetEntityType() const
Definition: G4Hype.cc:1363
G4double HypeInnerRadius2(G4double zVal) const
G4double GetSurfaceArea()
Definition: G4Hype.cc:1392
void SetInnerStereo(G4double newISte)
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Hype.cc:1530
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
void AddSurface(const G4ClippablePolygon &surface)
G4bool GetExtent(G4double &min, G4double &max) const
virtual void AddSolid(const G4Box &)=0
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4String GetName() const
G4double kCarTolerance
Definition: G4VSolid.hh:307
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:248
static G4int GetNumberOfRotationSteps()
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:58
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4int kMaxMeshSections
Definition: meshdefs.hh:46
Definition: DoubConv.h:17
T sqr(const T &x)
Definition: templates.hh:145
#define DBL_MIN
Definition: templates.hh:75
#define DBL_EPSILON
Definition: templates.hh:87