Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CutTubs.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// G4CutTubs
27//
28// Class description:
29//
30// G4CutTubs is a tube with possible cuts in +-Z.
31// Implementation adapted from G4Tubs (subclass of G4Tubs) and
32// from TGEo Ctube implementation (by A.Gheata, CERN)
33//
34// G4CutTubs(pName,pRMin,pRMax,pDZ,pSPhi,pEPhi,pLowNorm,pHighNorm)
35// pName,pRMin,pRMax,pDZ,pSPhi,pEPhi are the same as for G4Tubs,
36// pLowNorm=Outside Normal at -Z
37// pHighNorm=Outsie Normal at +Z.
38
39// Author: Tatiana Nikitina, CERN
40// --------------------------------------------------------------------
41
42#ifndef G4CUTTUBS_HH
43#define G4CUTTUBS_HH
44
45#include "G4GeomTypes.hh"
46
47#if defined(G4GEOM_USE_USOLIDS)
48#define G4GEOM_USE_UCTUBS 1
49#endif
50
51#if defined(G4GEOM_USE_UCTUBS)
52 #define G4UCutTubs G4CutTubs
53 #include "G4UCutTubs.hh"
54#else
55
56#include "G4CSGSolid.hh"
57#include "G4Polyhedron.hh"
58
59class G4CutTubs : public G4CSGSolid
60{
61 public: // with description
62
63 G4CutTubs( const G4String& pName,
64 G4double pRMin,
65 G4double pRMax,
66 G4double pDz,
67 G4double pSPhi,
68 G4double pDPhi,
69 G4ThreeVector pLowNorm,
70 G4ThreeVector pHighNorm );
71 //
72 // Constructs a tubs with the given name and dimensions
73
74 ~G4CutTubs();
75 //
76 // Destructor
77
78 // Accessors
79
80 inline G4double GetInnerRadius () const;
81 inline G4double GetOuterRadius () const;
82 inline G4double GetZHalfLength () const;
83 inline G4double GetStartPhiAngle () const;
84 inline G4double GetDeltaPhiAngle () const;
85 inline G4double GetSinStartPhi () const;
86 inline G4double GetCosStartPhi () const;
87 inline G4double GetSinEndPhi () const;
88 inline G4double GetCosEndPhi () const;
89 inline G4ThreeVector GetLowNorm () const;
90 inline G4ThreeVector GetHighNorm () const;
91
92 // Modifiers
93
94 inline void SetInnerRadius (G4double newRMin);
95 inline void SetOuterRadius (G4double newRMax);
96 inline void SetZHalfLength (G4double newDz);
97 inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true);
98 inline void SetDeltaPhiAngle (G4double newDPhi);
99
100 // Methods for solid
101
104
105 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
106
107 G4bool CalculateExtent( const EAxis pAxis,
108 const G4VoxelLimits& pVoxelLimit,
109 const G4AffineTransform& pTransform,
110 G4double& pmin, G4double& pmax ) const;
111
112 EInside Inside( const G4ThreeVector& p ) const;
113
114 G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
115
116 G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
117 G4double DistanceToIn(const G4ThreeVector& p) const;
119 const G4bool calcNorm = false,
120 G4bool* validNorm = nullptr,
121 G4ThreeVector* n = nullptr) const;
122 G4double DistanceToOut(const G4ThreeVector& p) const;
123
125
127
128 G4VSolid* Clone() const;
129
130 std::ostream& StreamInfo( std::ostream& os ) const;
131
132 // Visualisation functions
133
134 void DescribeYourselfTo ( G4VGraphicsScene& scene ) const;
136
137 public: // without description
138
139 G4CutTubs(__void__&);
140 //
141 // Fake default constructor for usage restricted to direct object
142 // persistency for clients requiring preallocation of memory for
143 // persistifiable objects.
144
145 G4CutTubs(const G4CutTubs& rhs);
146 G4CutTubs& operator=(const G4CutTubs& rhs);
147 // Copy constructor and assignment operator.
148
149 // Older names for access functions
150
151 inline G4double GetRMin() const;
152 inline G4double GetRMax() const;
153 inline G4double GetDz () const;
154 inline G4double GetSPhi() const;
155 inline G4double GetDPhi() const;
156
157 protected:
158
159 inline void Initialize();
160 //
161 // Reset relevant values to zero
162
163 inline void CheckSPhiAngle(G4double sPhi);
164 inline void CheckDPhiAngle(G4double dPhi);
165 inline void CheckPhiAngles(G4double sPhi, G4double dPhi);
166 //
167 // Reset relevant flags and angle values
168
170 //
171 // Recompute relevant trigonometric values and cache them
172
174 //
175 // Algorithm for SurfaceNormal() following the original
176 // specification for points not on the surface
177
179 // Check if the cutted planes are crossing.
180 // If 'true' , solid is ill defined
181
182 G4double GetCutZ(const G4ThreeVector& p) const;
183 // Get Z value of the point on Cutted Plane
184
185 private:
186
187 G4double kRadTolerance, kAngTolerance;
188 //
189 // Radial and angular tolerances
190
191 G4double fRMin, fRMax, fDz, fSPhi, fDPhi;
192 mutable G4double fZMin, fZMax;
193 //
194 // Radial and angular dimensions
195
196 G4double sinCPhi, cosCPhi, cosHDPhi, cosHDPhiOT, cosHDPhiIT,
197 sinSPhi, cosSPhi, sinEPhi, cosEPhi;
198 //
199 // Cached trigonometric values
200
201 G4bool fPhiFullCutTube = false;
202 //
203 // Flag for identification of section or full tube
204
205 G4double halfCarTolerance, halfRadTolerance, halfAngTolerance;
206 //
207 // Cached half tolerance values
208
209 G4ThreeVector fLowNorm, fHighNorm;
210 //
211 // Normals of Cut at -/+ Dz
212};
213
214#include "G4CutTubs.icc"
215
216#endif
217
218#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4ThreeVector GetHighNorm() const
G4double GetCubicVolume()
Definition: G4CutTubs.cc:245
G4double GetStartPhiAngle() const
void Initialize()
G4double GetZHalfLength() const
G4double GetInnerRadius() const
EInside Inside(const G4ThreeVector &p) const
Definition: G4CutTubs.cc:578
void SetInnerRadius(G4double newRMin)
G4ThreeVector GetLowNorm() const
void CheckDPhiAngle(G4double dPhi)
void SetStartPhiAngle(G4double newSPhi, G4bool trig=true)
void SetDeltaPhiAngle(G4double newDPhi)
G4double GetCutZ(const G4ThreeVector &p) const
Definition: G4CutTubs.cc:2163
G4Polyhedron * CreatePolyhedron() const
Definition: G4CutTubs.cc:2066
G4double GetDeltaPhiAngle() const
G4double GetOuterRadius() const
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4CutTubs.cc:763
G4bool IsCrossingCutPlanes() const
Definition: G4CutTubs.cc:2128
G4double GetSinStartPhi() const
G4GeometryType GetEntityType() const
Definition: G4CutTubs.cc:1900
void CheckPhiAngles(G4double sPhi, G4double dPhi)
void SetZHalfLength(G4double newDz)
G4CutTubs & operator=(const G4CutTubs &rhs)
Definition: G4CutTubs.cc:212
G4double GetRMin() const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4CutTubs.cc:350
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
Definition: G4CutTubs.cc:1420
G4VSolid * Clone() const
Definition: G4CutTubs.cc:1909
G4double GetSurfaceArea()
Definition: G4CutTubs.cc:296
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4CutTubs.cc:467
G4ThreeVector GetPointOnSurface() const
Definition: G4CutTubs.cc:1943
void InitializeTrigonometry()
G4double GetSinEndPhi() const
G4double GetSPhi() const
G4double GetDz() const
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4CutTubs.cc:1918
void CheckSPhiAngle(G4double sPhi)
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4CutTubs.cc:660
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4CutTubs.cc:907
G4double GetDPhi() const
G4double GetCosStartPhi() const
G4double GetCosEndPhi() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4CutTubs.cc:2061
G4double GetRMax() const
void SetOuterRadius(G4double newRMax)
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67