Geant4  10.03
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 100820 2016-11-02 15:18:48Z 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 #ifndef G4Cons_HH
67 #define G4Cons_HH
68 
69 #if defined(G4GEOM_USE_USOLIDS)
70 #define G4GEOM_USE_UCONS 1
71 #endif
72 
73 #if defined(G4GEOM_USE_UCONS)
74  #define G4UCons G4Cons
75  #include "G4UCons.hh"
76 #else
77 
78 #include <CLHEP/Units/PhysicalConstants.h>
79 
80 #include "G4CSGSolid.hh"
81 #include "G4Polyhedron.hh"
82 
83 class G4Cons : public G4CSGSolid
84 {
85  public: // with description
86 
87  G4Cons(const G4String& pName,
88  G4double pRmin1, G4double pRmax1,
89  G4double pRmin2, G4double pRmax2,
90  G4double pDz,
91  G4double pSPhi, G4double pDPhi);
92  //
93  // Constructs a cone with the given name and dimensions
94 
95  ~G4Cons() ;
96  //
97  // Destructor
98 
99  // Accessors
100 
101  inline G4double GetInnerRadiusMinusZ() const;
102  inline G4double GetOuterRadiusMinusZ() const;
103  inline G4double GetInnerRadiusPlusZ() const;
104  inline G4double GetOuterRadiusPlusZ() const;
105  inline G4double GetZHalfLength() const;
106  inline G4double GetStartPhiAngle() const;
107  inline G4double GetDeltaPhiAngle() const;
108  inline G4double GetSinStartPhi() const;
109  inline G4double GetCosStartPhi() const;
110  inline G4double GetSinEndPhi() const;
111  inline G4double GetCosEndPhi() const;
112 
113  // Modifiers
114 
115  inline void SetInnerRadiusMinusZ (G4double Rmin1 );
116  inline void SetOuterRadiusMinusZ (G4double Rmax1 );
117  inline void SetInnerRadiusPlusZ (G4double Rmin2 );
118  inline void SetOuterRadiusPlusZ (G4double Rmax2 );
119  inline void SetZHalfLength (G4double newDz );
120  inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true);
121  inline void SetDeltaPhiAngle (G4double newDPhi);
122 
123  // Other methods for solid
124 
125  inline G4double GetCubicVolume();
126  inline G4double GetSurfaceArea();
127 
129  const G4int n,
130  const G4VPhysicalVolume* pRep );
131 
132  void Extent(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
133 
134  G4bool CalculateExtent( const EAxis pAxis,
135  const G4VoxelLimits& pVoxelLimit,
136  const G4AffineTransform& pTransform,
137  G4double& pMin, G4double& pMax ) const;
138 
139  EInside Inside( const G4ThreeVector& p ) const;
140 
141  G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
142 
144  const G4ThreeVector& v) const;
145  G4double DistanceToIn (const G4ThreeVector& p) const;
147  const G4ThreeVector& v,
148  const G4bool calcNorm=G4bool(false),
149  G4bool *validNorm=0,
150  G4ThreeVector *n=0) const;
151  G4double DistanceToOut(const G4ThreeVector& p) const;
152 
154 
156 
157  G4VSolid* Clone() const;
158 
159  std::ostream& StreamInfo(std::ostream& os) const;
160 
161  // Visualisation functions
162 
163  void DescribeYourselfTo( G4VGraphicsScene& scene ) const;
165 
166  public: // without description
167 
168  G4Cons(__void__&);
169  //
170  // Fake default constructor for usage restricted to direct object
171  // persistency for clients requiring preallocation of memory for
172  // persistifiable objects.
173 
174  G4Cons(const G4Cons& rhs);
175  G4Cons& operator=(const G4Cons& rhs);
176  // Copy constructor and assignment operator.
177 
178  // Old access functions
179 
180  inline G4double GetRmin1() const;
181  inline G4double GetRmax1() const;
182  inline G4double GetRmin2() const;
183  inline G4double GetRmax2() const;
184  inline G4double GetDz() const;
185  inline G4double GetSPhi() const;
186  inline G4double GetDPhi() const;
187 
188  private:
189 
190  inline void Initialize();
191  //
192  // Reset relevant values to zero
193 
194  inline void CheckSPhiAngle(G4double sPhi);
195  inline void CheckDPhiAngle(G4double dPhi);
196  inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
197  //
198  // Reset relevant flags and angle values
199 
200  inline void InitializeTrigonometry();
201  //
202  // Recompute relevant trigonometric values and cache them
203 
205  //
206  // Algorithm for SurfaceNormal() following the original
207  // specification for points not on the surface
208 
209  private:
210 
211  // Used by distanceToOut
212  //
214 
215  // used by normal
216  //
218 
220  //
221  // Radial and angular tolerances
222 
224  //
225  // Radial and angular dimensions
226 
229  //
230  // Cached trigonometric values
231 
233  //
234  // Flag for identification of section or full cone
235 
237  //
238  // Cached half tolerance values
239 };
240 
241 #include "G4Cons.icc"
242 
243 #endif
244 
245 #endif
G4Cons(const G4String &pName, G4double pRmin1, G4double pRmax1, G4double pRmin2, G4double pRmax2, G4double pDz, G4double pSPhi, G4double pDPhi)
Definition: G4Cons.cc:87
void SetZHalfLength(G4double newDz)
G4double halfRadTolerance
Definition: G4Cons.hh:236
G4double fDz
Definition: G4Cons.hh:223
G4double GetSinEndPhi() const
void Initialize()
void SetInnerRadiusMinusZ(G4double Rmin1)
G4bool fPhiFullCone
Definition: G4Cons.hh:232
CLHEP::Hep3Vector G4ThreeVector
G4double GetDz() const
G4Cons & operator=(const G4Cons &rhs)
Definition: G4Cons.cc:177
G4double sinEPhi
Definition: G4Cons.hh:227
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Cons.cc:663
EInside Inside(const G4ThreeVector &p) const
Definition: G4Cons.cc:210
G4double fRmin1
Definition: G4Cons.hh:223
G4double GetRmin2() const
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Cons.cc:1395
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Cons.cc:537
void SetOuterRadiusPlusZ(G4double Rmax2)
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Cons.cc:2114
G4double GetOuterRadiusMinusZ() const
int G4int
Definition: G4Types.hh:78
G4GeometryType GetEntityType() const
Definition: G4Cons.cc:2096
G4double fRmax2
Definition: G4Cons.hh:223
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Cons.cc:265
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4double cosCPhi
Definition: G4Cons.hh:227
void CheckSPhiAngle(G4double sPhi)
G4double GetSinStartPhi() const
G4double kRadTolerance
Definition: G4Cons.hh:219
G4double GetSPhi() const
G4double fRmin2
Definition: G4Cons.hh:223
G4double sinSPhi
Definition: G4Cons.hh:227
G4double fSPhi
Definition: G4Cons.hh:223
void SetDeltaPhiAngle(G4double newDPhi)
bool G4bool
Definition: G4Types.hh:79
G4double halfCarTolerance
Definition: G4Cons.hh:236
Definition: G4Cons.hh:83
G4double GetStartPhiAngle() const
void InitializeTrigonometry()
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
Definition: G4Cons.cc:318
const G4int n
G4double GetInnerRadiusPlusZ() const
G4double cosSPhi
Definition: G4Cons.hh:227
G4double GetCubicVolume()
G4double GetCosEndPhi() const
void Extent(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Cons.cc:276
G4double fDPhi
Definition: G4Cons.hh:223
~G4Cons()
Definition: G4Cons.cc:151
G4double GetCosStartPhi() const
G4double halfAngTolerance
Definition: G4Cons.hh:236
EInside
Definition: geomdefs.hh:58
G4double GetRmin1() const
G4double GetSurfaceArea()
EAxis
Definition: geomdefs.hh:54
G4Polyhedron * CreatePolyhedron() const
Definition: G4Cons.cc:2233
G4double sinCPhi
Definition: G4Cons.hh:227
G4double kAngTolerance
Definition: G4Cons.hh:219
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Cons.cc:2228
G4double cosEPhi
Definition: G4Cons.hh:227
G4double GetDPhi() const
void SetOuterRadiusMinusZ(G4double Rmax1)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4double cosHDPhiOT
Definition: G4Cons.hh:227
G4ThreeVector GetPointOnSurface() const
Definition: G4Cons.cc:2141
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Cons.cc:433
G4double cosHDPhiIT
Definition: G4Cons.hh:227
double G4double
Definition: G4Types.hh:76
G4double GetInnerRadiusMinusZ() const
G4VSolid * Clone() const
Definition: G4Cons.cc:2105
void SetInnerRadiusPlusZ(G4double Rmin2)
G4double GetOuterRadiusPlusZ() const
G4double GetRmax2() const
G4double GetZHalfLength() const
void CheckDPhiAngle(G4double dPhi)
G4double fRmax1
Definition: G4Cons.hh:223
G4double GetDeltaPhiAngle() const
G4double GetRmax1() const