Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Cons.hh
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: G4Cons.hh 69788 2013-05-15 12:06:57Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 // G4Cons
34 //
35 // Class description:
36 //
37 // A G4Cons is, in the general case, a Phi segment of a cone, with
38 // half-length fDz, inner and outer radii specified at -fDz and +fDz.
39 // The Phi segment is described by a starting fSPhi angle, and the
40 // +fDPhi delta angle for the shape.
41 // If the delta angle is >=2*pi, the shape is treated as continuous
42 // in Phi
43 //
44 // Member Data:
45 //
46 // fRmin1 inside radius at -fDz
47 // fRmin2 inside radius at +fDz
48 // fRmax1 outside radius at -fDz
49 // fRmax2 outside radius at +fDz
50 // fDz half length in z
51 //
52 // fSPhi starting angle of the segment in radians
53 // fDPhi delta angle of the segment in radians
54 //
55 // fPhiFullCone Boolean variable used for indicate the Phi Section
56 //
57 // Note:
58 // Internally fSPhi & fDPhi are adjusted so that fDPhi<=2PI,
59 // and fDPhi+fSPhi<=2PI. This enables simpler comparisons to be
60 // made with (say) Phi of a point.
61 
62 // History:
63 // 19.3.94 P.Kent: Old C++ code converted to tolerant geometry
64 // 13.9.96 V.Grichine: Final modifications to commit
65 // --------------------------------------------------------------------
66 
67 #ifndef G4Cons_HH
68 #define G4Cons_HH
69 
71 
72 #include "G4CSGSolid.hh"
73 
74 class G4Cons : public G4CSGSolid
75 {
76  public: // with description
77 
78  G4Cons(const G4String& pName,
79  G4double pRmin1, G4double pRmax1,
80  G4double pRmin2, G4double pRmax2,
81  G4double pDz,
82  G4double pSPhi, G4double pDPhi);
83  //
84  // Constructs a cone with the given name and dimensions
85 
86  ~G4Cons() ;
87  //
88  // Destructor
89 
90  // Accessors
91 
92  inline G4double GetInnerRadiusMinusZ() const;
93  inline G4double GetOuterRadiusMinusZ() const;
94  inline G4double GetInnerRadiusPlusZ() const;
95  inline G4double GetOuterRadiusPlusZ() const;
96  inline G4double GetZHalfLength() const;
97  inline G4double GetStartPhiAngle() const;
98  inline G4double GetDeltaPhiAngle() const;
99 
100  // Modifiers
101 
102  inline void SetInnerRadiusMinusZ (G4double Rmin1 );
103  inline void SetOuterRadiusMinusZ (G4double Rmax1 );
104  inline void SetInnerRadiusPlusZ (G4double Rmin2 );
105  inline void SetOuterRadiusPlusZ (G4double Rmax2 );
106  inline void SetZHalfLength (G4double newDz );
107  inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true);
108  inline void SetDeltaPhiAngle (G4double newDPhi);
109 
110  // Other methods for solid
111 
112  inline G4double GetCubicVolume();
113  inline G4double GetSurfaceArea();
114 
116  const G4int n,
117  const G4VPhysicalVolume* pRep );
118 
119  G4bool CalculateExtent( const EAxis pAxis,
120  const G4VoxelLimits& pVoxelLimit,
121  const G4AffineTransform& pTransform,
122  G4double& pmin, G4double& pmax ) const;
123 
124  EInside Inside( const G4ThreeVector& p ) const;
125 
126  G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
127 
129  const G4ThreeVector& v) const;
130  G4double DistanceToIn (const G4ThreeVector& p) const;
132  const G4ThreeVector& v,
133  const G4bool calcNorm=G4bool(false),
134  G4bool *validNorm=0,
135  G4ThreeVector *n=0) const;
136  G4double DistanceToOut(const G4ThreeVector& p) const;
137 
139 
141 
142  G4VSolid* Clone() const;
143 
144  std::ostream& StreamInfo(std::ostream& os) const;
145 
146  // Visualisation functions
147 
148  void DescribeYourselfTo( G4VGraphicsScene& scene ) const;
150  G4NURBS* CreateNURBS() const;
151 
152  public: // without description
153 
154  G4Cons(__void__&);
155  //
156  // Fake default constructor for usage restricted to direct object
157  // persistency for clients requiring preallocation of memory for
158  // persistifiable objects.
159 
160  G4Cons(const G4Cons& rhs);
161  G4Cons& operator=(const G4Cons& rhs);
162  // Copy constructor and assignment operator.
163 
164  // Old access functions
165 
166  inline G4double GetRmin1() const;
167  inline G4double GetRmax1() const;
168  inline G4double GetRmin2() const;
169  inline G4double GetRmax2() const;
170  inline G4double GetDz() const;
171  inline G4double GetSPhi() const;
172  inline G4double GetDPhi() const;
173 
174  private:
175 
177  CreateRotatedVertices(const G4AffineTransform& pTransform) const;
178 
179  inline void Initialize();
180  //
181  // Reset relevant values to zero
182 
183  inline void CheckSPhiAngle(G4double sPhi);
184  inline void CheckDPhiAngle(G4double dPhi);
185  inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
186  //
187  // Reset relevant flags and angle values
188 
189  inline void InitializeTrigonometry();
190  //
191  // Recompute relevant trigonometric values and cache them
192 
193  G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector& p) const;
194  //
195  // Algorithm for SurfaceNormal() following the original
196  // specification for points not on the surface
197 
198  private:
199 
200  // Used by distanceToOut
201  //
202  enum ESide {kNull,kRMin,kRMax,kSPhi,kEPhi,kPZ,kMZ};
203 
204  // used by normal
205  //
206  enum ENorm {kNRMin,kNRMax,kNSPhi,kNEPhi,kNZ};
207 
208  G4double kRadTolerance, kAngTolerance;
209  //
210  // Radial and angular tolerances
211 
212  G4double fRmin1, fRmin2, fRmax1, fRmax2, fDz, fSPhi, fDPhi;
213  //
214  // Radial and angular dimensions
215 
216  G4double sinCPhi, cosCPhi, cosHDPhiOT, cosHDPhiIT,
217  sinSPhi, cosSPhi, sinEPhi, cosEPhi;
218  //
219  // Cached trigonometric values
220 
221  G4bool fPhiFullCone;
222  //
223  // Flag for identification of section or full cone
224 };
225 
226 #include "G4Cons.icc"
227 
228 #endif