Geant4  10.01.p03
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 79491 2014-03-05 15:24:29Z 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 
109  // Modifiers
110 
111  inline void SetInnerRadiusMinusZ (G4double Rmin1 );
112  inline void SetOuterRadiusMinusZ (G4double Rmax1 );
113  inline void SetInnerRadiusPlusZ (G4double Rmin2 );
114  inline void SetOuterRadiusPlusZ (G4double Rmax2 );
115  inline void SetZHalfLength (G4double newDz );
116  inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true);
117  inline void SetDeltaPhiAngle (G4double newDPhi);
118 
119  // Other methods for solid
120 
121  inline G4double GetCubicVolume();
122  inline G4double GetSurfaceArea();
123 
125  const G4int n,
126  const G4VPhysicalVolume* pRep );
127 
128  G4bool CalculateExtent( const EAxis pAxis,
129  const G4VoxelLimits& pVoxelLimit,
130  const G4AffineTransform& pTransform,
131  G4double& pmin, G4double& pmax ) const;
132 
133  EInside Inside( const G4ThreeVector& p ) const;
134 
135  G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
136 
138  const G4ThreeVector& v) const;
139  G4double DistanceToIn (const G4ThreeVector& p) const;
141  const G4ThreeVector& v,
142  const G4bool calcNorm=G4bool(false),
143  G4bool *validNorm=0,
144  G4ThreeVector *n=0) const;
145  G4double DistanceToOut(const G4ThreeVector& p) const;
146 
148 
150 
151  G4VSolid* Clone() const;
152 
153  std::ostream& StreamInfo(std::ostream& os) const;
154 
155  // Visualisation functions
156 
157  void DescribeYourselfTo( G4VGraphicsScene& scene ) const;
159 
160  public: // without description
161 
162  G4Cons(__void__&);
163  //
164  // Fake default constructor for usage restricted to direct object
165  // persistency for clients requiring preallocation of memory for
166  // persistifiable objects.
167 
168  G4Cons(const G4Cons& rhs);
169  G4Cons& operator=(const G4Cons& rhs);
170  // Copy constructor and assignment operator.
171 
172  // Old access functions
173 
174  inline G4double GetRmin1() const;
175  inline G4double GetRmax1() const;
176  inline G4double GetRmin2() const;
177  inline G4double GetRmax2() const;
178  inline G4double GetDz() const;
179  inline G4double GetSPhi() const;
180  inline G4double GetDPhi() const;
181 
182  private:
183 
185  CreateRotatedVertices(const G4AffineTransform& pTransform) const;
186 
187  inline void Initialize();
188  //
189  // Reset relevant values to zero
190 
191  inline void CheckSPhiAngle(G4double sPhi);
192  inline void CheckDPhiAngle(G4double dPhi);
193  inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
194  //
195  // Reset relevant flags and angle values
196 
197  inline void InitializeTrigonometry();
198  //
199  // Recompute relevant trigonometric values and cache them
200 
202  //
203  // Algorithm for SurfaceNormal() following the original
204  // specification for points not on the surface
205 
206  private:
207 
208  // Used by distanceToOut
209  //
211 
212  // used by normal
213  //
215 
217  //
218  // Radial and angular tolerances
219 
221  //
222  // Radial and angular dimensions
223 
226  //
227  // Cached trigonometric values
228 
230  //
231  // Flag for identification of section or full cone
232 
234  //
235  // Cached half tolerance values
236 };
237 
238 #include "G4Cons.icc"
239 
240 #endif
241 
242 #endif
G4Cons(const G4String &pName, G4double pRmin1, G4double pRmax1, G4double pRmin2, G4double pRmax2, G4double pDz, G4double pSPhi, G4double pDPhi)
Definition: G4Cons.cc:82
void SetZHalfLength(G4double newDz)
G4double halfRadTolerance
Definition: G4Cons.hh:233
G4double fDz
Definition: G4Cons.hh:220
void Initialize()
void SetInnerRadiusMinusZ(G4double Rmin1)
G4bool fPhiFullCone
Definition: G4Cons.hh:229
CLHEP::Hep3Vector G4ThreeVector
G4double GetDz() const
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Cons.cc:2161
G4Cons & operator=(const G4Cons &rhs)
Definition: G4Cons.cc:172
G4double sinEPhi
Definition: G4Cons.hh:224
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Cons.cc:720
EInside Inside(const G4ThreeVector &p) const
Definition: G4Cons.cc:205
G4double fRmin1
Definition: G4Cons.hh:220
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:1452
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Cons.cc:594
void SetOuterRadiusPlusZ(G4double Rmax2)
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Cons.cc:2265
G4double GetOuterRadiusMinusZ() const
int G4int
Definition: G4Types.hh:78
G4GeometryType GetEntityType() const
Definition: G4Cons.cc:2247
G4double fRmax2
Definition: G4Cons.hh:220
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Cons.cc:260
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4double cosCPhi
Definition: G4Cons.hh:224
void CheckSPhiAngle(G4double sPhi)
G4double kRadTolerance
Definition: G4Cons.hh:216
G4double GetSPhi() const
G4double fRmin2
Definition: G4Cons.hh:220
G4double sinSPhi
Definition: G4Cons.hh:224
G4double fSPhi
Definition: G4Cons.hh:220
void SetDeltaPhiAngle(G4double newDPhi)
bool G4bool
Definition: G4Types.hh:79
G4double halfCarTolerance
Definition: G4Cons.hh:233
Definition: G4Cons.hh:83
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
G4double GetStartPhiAngle() const
void InitializeTrigonometry()
const G4int n
G4double GetInnerRadiusPlusZ() const
G4double cosSPhi
Definition: G4Cons.hh:224
G4double GetCubicVolume()
G4double fDPhi
Definition: G4Cons.hh:220
~G4Cons()
Definition: G4Cons.cc:146
G4double halfAngTolerance
Definition: G4Cons.hh:233
EInside
Definition: geomdefs.hh:58
G4double GetRmin1() const
G4double GetSurfaceArea()
EAxis
Definition: geomdefs.hh:54
G4Polyhedron * CreatePolyhedron() const
Definition: G4Cons.cc:2384
G4double sinCPhi
Definition: G4Cons.hh:224
G4double kAngTolerance
Definition: G4Cons.hh:216
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Cons.cc:2379
G4double cosEPhi
Definition: G4Cons.hh:224
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Cons.cc:272
G4double GetDPhi() const
void SetOuterRadiusMinusZ(G4double Rmax1)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4double cosHDPhiOT
Definition: G4Cons.hh:224
G4ThreeVector GetPointOnSurface() const
Definition: G4Cons.cc:2292
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Cons.cc:490
G4double cosHDPhiIT
Definition: G4Cons.hh:224
double G4double
Definition: G4Types.hh:76
G4double GetInnerRadiusMinusZ() const
G4VSolid * Clone() const
Definition: G4Cons.cc:2256
void SetInnerRadiusPlusZ(G4double Rmin2)
G4double GetOuterRadiusPlusZ() const
G4double GetRmax2() const
G4double GetZHalfLength() const
void CheckDPhiAngle(G4double dPhi)
G4double fRmax1
Definition: G4Cons.hh:220
G4double GetDeltaPhiAngle() const
G4double GetRmax1() const