Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4CutTubs.hh
이 파일의 문서화 페이지로 가기
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: G4CutTubs.hh 105075 2017-07-11 14:22:53Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 //
34 // G4CutTubs is a tube with possible cuts in +-Z.
35 // Implementation adapted from G4Tubs (subclass of G4Tubs) and
36 // from TGEo Ctube implementation (by A.Gheata, CERN)
37 //
38 // G4CutTubs(pName,pRMin,pRMax,pDZ,pSPhi,pEPhi,pLowNorm,pHighNorm)
39 // pName,pRMin,pRMax,pDZ,pSPhi,pEPhi are the same as for G4Tubs,
40 // pLowNorm=Outside Normal at -Z
41 // pHighNorm=Outsie Normal at +Z.
42 
43 // Author: Tatiana Nikitina, CERN
44 // --------------------------------------------------------------------
45 
46 #ifndef G4CUTTUBS_HH
47 #define G4CUTTUBS_HH
48 
49 #if defined(G4GEOM_USE_USOLIDS)
50 #define G4GEOM_USE_UCTUBS 1
51 #endif
52 
53 #if defined(G4GEOM_USE_UCTUBS)
54  #define G4UCutTubs G4CutTubs
55  #include "G4UCutTubs.hh"
56 #else
57 
58 #include "G4CSGSolid.hh"
59 #include "G4Polyhedron.hh"
60 
61 class G4CutTubs : public G4CSGSolid
62 {
63  public: // with description
64 
65  G4CutTubs( const G4String& pName,
66  G4double pRMin,
67  G4double pRMax,
68  G4double pDz,
69  G4double pSPhi,
70  G4double pDPhi,
71  G4ThreeVector pLowNorm,
72  G4ThreeVector pHighNorm );
73  //
74  // Constructs a tubs with the given name and dimensions
75 
76  ~G4CutTubs();
77  //
78  // Destructor
79 
80  // Accessors
81 
82  inline G4double GetInnerRadius () const;
83  inline G4double GetOuterRadius () const;
84  inline G4double GetZHalfLength () const;
85  inline G4double GetStartPhiAngle () const;
86  inline G4double GetDeltaPhiAngle () const;
87  inline G4double GetSinStartPhi () const;
88  inline G4double GetCosStartPhi () const;
89  inline G4double GetSinEndPhi () const;
90  inline G4double GetCosEndPhi () const;
91  inline G4ThreeVector GetLowNorm () const;
92  inline G4ThreeVector GetHighNorm () const;
93 
94  // Modifiers
95 
96  inline void SetInnerRadius (G4double newRMin);
97  inline void SetOuterRadius (G4double newRMax);
98  inline void SetZHalfLength (G4double newDz);
99  inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true);
100  inline void SetDeltaPhiAngle (G4double newDPhi);
101 
102  // Methods for solid
103 
104  inline G4double GetCubicVolume();
105  inline G4double GetSurfaceArea();
106 
107  void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
108 
109  G4bool CalculateExtent( const EAxis pAxis,
110  const G4VoxelLimits& pVoxelLimit,
111  const G4AffineTransform& pTransform,
112  G4double& pmin, G4double& pmax ) const;
113 
114  EInside Inside( const G4ThreeVector& p ) const;
115 
116  G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
117 
118  G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
119  G4double DistanceToIn(const G4ThreeVector& p) const;
121  const G4bool calcNorm=G4bool(false),
122  G4bool *validNorm=0, G4ThreeVector *n=0) const;
123  G4double DistanceToOut(const G4ThreeVector& p) const;
124 
126 
128 
129  G4VSolid* Clone() const;
130 
131  std::ostream& StreamInfo( std::ostream& os ) const;
132 
133  // Visualisation functions
134 
135  void DescribeYourselfTo ( G4VGraphicsScene& scene ) const;
136  G4Polyhedron* CreatePolyhedron () const;
137 
138  public: // without description
139 
140  G4CutTubs(__void__&);
141  //
142  // Fake default constructor for usage restricted to direct object
143  // persistency for clients requiring preallocation of memory for
144  // persistifiable objects.
145 
146  G4CutTubs(const G4CutTubs& rhs);
147  G4CutTubs& operator=(const G4CutTubs& rhs);
148  // Copy constructor and assignment operator.
149 
150  // Older names for access functions
151 
152  inline G4double GetRMin() const;
153  inline G4double GetRMax() const;
154  inline G4double GetDz () const;
155  inline G4double GetSPhi() const;
156  inline G4double GetDPhi() const;
157 
158  protected:
159 
160  inline void Initialize();
161  //
162  // Reset relevant values to zero
163 
164  inline void CheckSPhiAngle(G4double sPhi);
165  inline void CheckDPhiAngle(G4double dPhi);
166  inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
167  //
168  // Reset relevant flags and angle values
169 
170  inline void InitializeTrigonometry();
171  //
172  // Recompute relevant trigonometric values and cache them
173 
175  //
176  // Algorithm for SurfaceNormal() following the original
177  // specification for points not on the surface
178 
179  G4bool IsCrossingCutPlanes() const;
180  // Check if the cutted planes are crossing.
181  // If 'true' , solid is ill defined
182 
183  G4double GetCutZ(const G4ThreeVector& p) const;
184  // Get Z value of the point on Cutted Plane
185 
186  void GetMaxMinZ(G4double& zmin,G4double& zmax)const;
187  // Get Max and Min values of Z on Cutted Plane,
188  // Used for Calculate BoundingLimits()
189 
190  private:
191 
193  //
194  // Radial and angular tolerances
195 
197  //
198  // Radial and angular dimensions
199 
202  //
203  // Cached trigonometric values
204 
206  //
207  // Flag for identification of section or full tube
208 
210  //
211  // Cached half tolerance values
212 
214  //
215  // Normals of Cut at -/+ Dz
216 };
217 
218 #include "G4CutTubs.icc"
219 
220 #endif
221 
222 #endif
G4double GetSinEndPhi() const
G4bool IsCrossingCutPlanes() const
Definition: G4CutTubs.cc:1973
void SetInnerRadius(G4double newRMin)
G4double kRadTolerance
Definition: G4CutTubs.hh:192
G4Polyhedron * CreatePolyhedron() const
Definition: G4CutTubs.cc:1913
G4double fRMax
Definition: G4CutTubs.hh:196
G4double cosHDPhiIT
Definition: G4CutTubs.hh:200
G4ThreeVector fLowNorm
Definition: G4CutTubs.hh:213
G4double GetDPhi() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4CutTubs.cc:1908
G4double GetOuterRadius() const
G4double GetRMin() const
const char * p
Definition: xmltok.h:285
G4double sinEPhi
Definition: G4CutTubs.hh:200
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4CutTubs.cc:798
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
G4double fSPhi
Definition: G4CutTubs.hh:196
G4CutTubs & operator=(const G4CutTubs &rhs)
Definition: G4CutTubs.cc:209
G4double GetCosStartPhi() const
G4double GetRMax() const
void SetOuterRadius(G4double newRMax)
G4double GetSinStartPhi() const
G4double fRMin
Definition: G4CutTubs.hh:196
G4bool fPhiFullCutTube
Definition: G4CutTubs.hh:205
G4ThreeVector GetLowNorm() const
G4double halfCarTolerance
Definition: G4CutTubs.hh:209
G4double GetCutZ(const G4ThreeVector &p) const
Definition: G4CutTubs.cc:1996
G4double cosSPhi
Definition: G4CutTubs.hh:200
G4double GetSPhi() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4CutTubs.cc:241
G4double cosCPhi
Definition: G4CutTubs.hh:200
G4ThreeVector GetPointOnSurface() const
Definition: G4CutTubs.cc:1834
G4double GetInnerRadius() const
void InitializeTrigonometry()
void Initialize()
void CheckDPhiAngle(G4double dPhi)
void CheckSPhiAngle(G4double sPhi)
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4CutTubs.cc:1809
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4CutTubs.cc:551
G4double halfAngTolerance
Definition: G4CutTubs.hh:209
G4ThreeVector fHighNorm
Definition: G4CutTubs.hh:213
G4CutTubs(const G4String &pName, G4double pRMin, G4double pRMax, G4double pDz, G4double pSPhi, G4double pDPhi, G4ThreeVector pLowNorm, G4ThreeVector pHighNorm)
Definition: G4CutTubs.cc:66
EInside Inside(const G4ThreeVector &p) const
Definition: G4CutTubs.cc:469
G4double fDz
Definition: G4CutTubs.hh:196
EInside
Definition: geomdefs.hh:58
G4GeometryType GetEntityType() const
Definition: G4CutTubs.cc:1791
G4double halfRadTolerance
Definition: G4CutTubs.hh:209
EAxis
Definition: geomdefs.hh:54
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4CutTubs.cc:358
G4double sinCPhi
Definition: G4CutTubs.hh:200
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4CutTubs.cc:654
G4double GetCosEndPhi() const
G4double sinSPhi
Definition: G4CutTubs.hh:200
G4double GetCubicVolume()
G4double fDPhi
Definition: G4CutTubs.hh:196
G4ThreeVector GetHighNorm() const
void SetDeltaPhiAngle(G4double newDPhi)
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4CutTubs.cc:1311
Char_t n[5]
G4double GetDeltaPhiAngle() const
G4double kAngTolerance
Definition: G4CutTubs.hh:192
G4double cosEPhi
Definition: G4CutTubs.hh:200
G4double GetStartPhiAngle() const
G4VSolid * Clone() const
Definition: G4CutTubs.cc:1800
G4double cosHDPhiOT
Definition: G4CutTubs.hh:200
void SetZHalfLength(G4double newDz)
void GetMaxMinZ(G4double &zmin, G4double &zmax) const
Definition: G4CutTubs.cc:2020
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4double GetSurfaceArea()
G4double GetZHalfLength() const
G4double GetDz() const