/home/docs/checkouts/readthedocs.org/user_builds/ratpac/checkouts/latest/src/geo/include/RAT/GLG4TorusStack.hh Source File

Ratpac-two: /home/docs/checkouts/readthedocs.org/user_builds/ratpac/checkouts/latest/src/geo/include/RAT/GLG4TorusStack.hh Source File
Ratpac-two
GLG4TorusStack.hh
1 // This file is part of the GenericLAND software library.
2 // $Id: GLG4TorusStack.hh,v 1.1 2005/08/30 19:55:22 volsung Exp $
3 //
4 // This code partly derived from the intellectual property of
5 // the RD44 GEANT4 collaboration.
6 //
7 // By copying, distributing or modifying the Program (or any work
8 // based on the Program) you indicate your acceptance of this statement,
9 // and all its terms, whatever they may be.
10 //
11 //
12 // class GLG4TorusStack
13 //
14 // A stack of torii sliced normal to the Z axis, with curved
15 // sides parallel to the z-axis. Each torus has a certain swept
16 // radius about which it is centered, and a certain cross-sectional
17 // radius. A negative cross-sectional radius means to use the "inner"
18 // surface of the torus as the boundary of the toroid stack.
19 // A zero cross-sectional radius may be used to indicate a cylinder.
20 // The complete stack is filled and complete in phi angle.
21 // The torii are actually specified in terms of the z and "rho" coordinates
22 // of the edges of the torii and the z coordinate of each segment's plane of
23 // symmetry; the swept radii ("a") and torii radii ("b")
24 // are calculated accordingly.
25 //
26 //
27 // Member functions:
28 //
29 // As inherited from G4CSGSolid+
30 //
31 // GLG4TorusStack(const G4String &pName )
32 //
33 // Construct an empty TorusStack with the given name.
34 //
35 //
36 // void SetAllParameters
37 // (G4int n_, // number of Z-segments
38 // const G4double z_edge_[ ], // n+1 edges of Z-segments
39 // const double rho_edge_[ ], // n+1 dist. from Z-axis
40 // const G4double z_o_[ ] ) // z-origins (n total)
41 //
42 // Define number of torii in stack and the dimensions of each.
43 //
44 // ****************************************************************/
45 //
46 // Protected:
47 //
48 // G4ThreeVectorList*
49 // CreateRotatedVertices(const G4AffineTransform& pTransform) const
50 //
51 // Create the List of transformed vertices in the format required
52 // for G4VSolid:: ClipCrossSection and ClipBetweenSections.
53 //
54 //
55 // 1999.11.22 G.Horton-Smith First version of GLG4TorusStack
56 
57 #ifndef GLG4TorusStack_HH
58 #define GLG4TorusStack_HH
59 
60 #include "G4CSGSolid.hh"
61 
62 class GLG4TorusStack : public G4CSGSolid {
63  public:
64  GLG4TorusStack(const G4String &pName);
65 
66  virtual ~GLG4TorusStack();
67 
68  void SetAllParameters(G4int n, // number of Z-segments
69  const G4double z_edge[], // n+1 edges of Z-segments
70  const G4double rho_edge[], // n+1 dist. from Z-axis
71  const G4double z_o[], // tor z-origins (n total)
72  GLG4TorusStack *inner = 0); // inner TorusStack (opt.)
73 
74  void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep);
75 
76  G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform,
77  G4double &pmin, G4double &pmax) const;
78 
79  G4int GetN() const { return n; }
80  G4double GetZEdge(int i) const { return z_edge[i]; }
81  G4double GetRhoEdge(int i) const { return rho_edge[i]; }
82  G4double GetZo(int i) const { return z_o[i]; }
83  G4double GetRo(int i) const { return r_o[i]; }
84  G4double GetRadius(int i) const { return radius[i]; }
85  GLG4TorusStack *GetInner() const { return inner; }
86 
87  EInside Inside(const G4ThreeVector &p) const;
88 
89  G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const;
90 
91  G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const;
92  G4double DistanceToIn(const G4ThreeVector &p) const;
93  G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm = G4bool(false),
94  G4bool *validNorm = 0, G4ThreeVector *n = 0) const;
95  G4double DistanceToOut(const G4ThreeVector &p) const;
96 
97  // Naming method (pseudo-RTTI : run-time type identification)
98  virtual G4GeometryType GetEntityType() const { return G4String("GLG4TorusStack"); }
99 
100  // Visualisation functions
101  void DescribeYourselfTo(G4VGraphicsScene &scene) const;
102  G4VisExtent GetExtent() const;
103  G4Polyhedron *CreatePolyhedron() const;
104 
105  // other generally useful functions -- public so others can use them!
106 
107  // find first torus intersection -- s == distance from p along v
108  static G4int FindFirstTorusRoot(G4double a, // swept radius
109  G4double b, // cross-section radius
110  const G4ThreeVector &p, // start point of ray
111  const G4ThreeVector &v, // direction of ray
112  G4double smin, // lower bracket on root
113  G4double smax, // upper bracket on root
114  G4bool fEntering, // true if looking for out->in
115  G4double &sout); // distance to root, if found
116 
117  // find first root of arbitrary function in bracketed interval
118  class RootFinder {
119  public:
120  virtual void f_and_Df(G4double x, G4double &f, G4double &Df) = 0;
121  //
122  // f_and_Df must be overridden by user to be function which sets
123  // f and Df to the value of the function and the derivative w.r.t.
124  // x of the function, respectively, for a given x
125  //
126  int FindRoot(G4double smin, // lower bound on root
127  G4double smax, // uppper bound on root
128  G4double tol, // tolerance on root
129  G4bool fFindFallingRoot, // true == "want dF<0 at root"
130  // false == "want dF>0 at root"
131  G4double &sout); // returned root
132  //
133  // FindRoot returns 1 if root found, 0 if no root found.
134  //
135  };
136 
137  // return integer i such that z[i] <= z_lu < z[i+1]
138  // or -1 if z_lu < z[0] or z_lu >= z[n-1]
139  static G4int FindInOrderedList(double z_lu, const double *z, int n);
140 
141  // a helpful function that doesn't hurt anything to call:
142  G4int FindNearestSegment(G4double pr, G4double pz) const;
143 
144  protected:
145  G4ThreeVectorList *CreateRotatedVertices(const G4AffineTransform &pTransform, G4int &noPolygonVertices) const;
146 
147  EInside Inside1(const G4ThreeVector &p) const;
148 
149  void CheckABRho();
150 
151  int n; // number of Z-segments
152  double *z_edge; // n+1 edges of Z-segments
153  double *rho_edge; // n+1 2-d distance from Z-axis at each edge
154  double *z_o; // z-origins, one for each toroid segment (n total)
155  double *r_o; // swept radii, one for each toroid segment (n total)
156  double *radius; // torus radii, one for each toroid segment (n total)
157  double max_rho; // maxium distance of surface from Z axis
158  double myRadTolerance; // because Geant4.1.0 default is too small
159  static double surfaceTolerance; // in GEANT4.9.0 it is not a global const
160  static double radTolerance; // in GEANT4.9.0 it is not a global const
161  static double angTolerance; // in GEANT4.9.0 it is not a global const
162  GLG4TorusStack *inner; // because G4SubtractionSolid is bad
163 };
164 
165 #endif
Definition: GLG4TorusStack.hh:118
Definition: GLG4TorusStack.hh:62