Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VTwistedFaceted.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// G4VTwistedFaceted
27//
28// Class description:
29//
30// G4VTwistedFaceted is an abstract base class for twisted boxoids:
31// G4TwistedTrd, G4TwistedTrap and G4TwistedBox
32
33// Author: 27-Oct-2004 - O.Link (Oliver.Link@cern.ch)
34// --------------------------------------------------------------------
35#ifndef G4VTWISTEDFACETED_HH
36#define G4VTWISTEDFACETED_HH
37
38#include "G4VSolid.hh"
39#include "G4TwoVector.hh"
42#include "G4TwistBoxSide.hh"
44
47
49{
50 public: // with description
51
52 G4VTwistedFaceted(const G4String& pname, // Name of instance
53 G4double PhiTwist, // twist angle
54 G4double pDz, // half z lenght
55 G4double pTheta, // direction between end planes
56 G4double pPhi, // defined by polar & azim. angles
57 G4double pDy1, // half y length at -pDz
58 G4double pDx1, // half x length at -pDz,-pDy
59 G4double pDx2, // half x length at -pDz,+pDy
60 G4double pDy2, // half y length at +pDz
61 G4double pDx3, // half x length at +pDz,-pDy
62 G4double pDx4, // half x length at +pDz,+pDy
63 G4double pAlph // tilt angle at +pDz
64 );
65
66 virtual ~G4VTwistedFaceted();
67
69 const G4int,
70 const G4VPhysicalVolume* );
71
72 virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const;
73
74 virtual G4bool CalculateExtent(const EAxis pAxis,
75 const G4VoxelLimits& pVoxelLimit,
76 const G4AffineTransform& pTransform,
77 G4double& pMin,
78 G4double& pMax ) const;
79
80 virtual G4double DistanceToIn (const G4ThreeVector& p,
81 const G4ThreeVector& v ) const;
82
83 virtual G4double DistanceToIn (const G4ThreeVector& p ) const;
84
85 virtual G4double DistanceToOut(const G4ThreeVector& p,
86 const G4ThreeVector& v,
87 const G4bool calcnorm = false,
88 G4bool* validnorm = nullptr,
89 G4ThreeVector* n = nullptr ) const;
90
91 virtual G4double DistanceToOut(const G4ThreeVector& p) const;
92
93 virtual EInside Inside (const G4ThreeVector& p) const;
94
95 virtual G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const;
96
99
100 virtual G4double GetCubicVolume();
101 virtual G4double GetSurfaceArea();
102
103 virtual void DescribeYourselfTo (G4VGraphicsScene& scene) const;
104 virtual G4Polyhedron* CreatePolyhedron () const;
105 virtual G4Polyhedron* GetPolyhedron () const;
106
107 virtual std::ostream &StreamInfo(std::ostream& os) const;
108
109 // accessors
110
111 inline G4double GetTwistAngle() const { return fPhiTwist; }
112
113 inline G4double GetDx1 () const { return fDx1 ; }
114 inline G4double GetDx2 () const { return fDx2 ; }
115 inline G4double GetDx3 () const { return fDx3 ; }
116 inline G4double GetDx4 () const { return fDx4 ; }
117 inline G4double GetDy1 () const { return fDy1 ; }
118 inline G4double GetDy2 () const { return fDy2 ; }
119 inline G4double GetDz () const { return fDz ; }
120 inline G4double GetPhi () const { return fPhi ; }
121 inline G4double GetTheta () const { return fTheta ; }
122 inline G4double GetAlpha () const { return fAlph ; }
123
124 inline G4double Xcoef(G4double u, G4double phi, G4double ftg) const;
125 // For calculating the w(u) function
126
127 inline G4double GetValueA(G4double phi) const;
128 inline G4double GetValueB(G4double phi) const;
129 inline G4double GetValueD(G4double phi) const;
130
131 virtual G4VisExtent GetExtent () const;
132 virtual G4GeometryType GetEntityType() const;
133
134 public: // without description
135
136 G4VTwistedFaceted(__void__&);
137 // Fake default constructor for usage restricted to direct object
138 // persistency for clients requiring preallocation of memory for
139 // persistifiable objects.
140
143 // Copy constructor and assignment operator.
144
145 protected: // with description
146
147 mutable G4bool fRebuildPolyhedron = false;
148 mutable G4Polyhedron* fpPolyhedron = nullptr; // polyhedron for vis
149
150 private:
151
152 double GetLateralFaceArea(const G4TwoVector& p1,
153 const G4TwoVector& p2,
154 const G4TwoVector& p3,
155 const G4TwoVector& p4) const;
156 void CreateSurfaces();
157
158 private:
159
160 G4double fTheta;
161 G4double fPhi ;
162
163 G4double fDy1;
164 G4double fDx1;
165 G4double fDx2;
166
167 G4double fDy2;
168 G4double fDx3;
169 G4double fDx4;
170
171 G4double fDz; // Half-length along the z axis
172
173 G4double fDx ; // maximum side in x
174 G4double fDy ; // maximum side in y
175
176 G4double fAlph ;
177 G4double fTAlph ; // std::tan(fAlph)
178
179 G4double fdeltaX ;
180 G4double fdeltaY ;
181
182 G4double fPhiTwist; // twist angle ( dphi in surface equation)
183
184 G4VTwistSurface* fLowerEndcap ; // surface of -ve z
185 G4VTwistSurface* fUpperEndcap ; // surface of +ve z
186
187 G4VTwistSurface* fSide0 ; // Twisted Side at phi = 0 deg
188 G4VTwistSurface* fSide90 ; // Twisted Side at phi = 90 deg
189 G4VTwistSurface* fSide180 ; // Twisted Side at phi = 180 deg
190 G4VTwistSurface* fSide270 ; // Twisted Side at phi = 270 deg
191
192 protected:
193
194 G4double fCubicVolume = 0.0; // volume of the solid
195 G4double fSurfaceArea = 0.0; // area of the solid
196
197 private:
198
199 class LastState // last Inside result
200 {
201 public:
202 LastState()
203 {
204 p.set(kInfinity,kInfinity,kInfinity); inside = kOutside;
205 }
206 ~LastState(){}
207 LastState(const LastState& r) : p(r.p), inside(r.inside){}
208 LastState& operator=(const LastState& r)
209 {
210 if (this == &r) { return *this; }
211 p = r.p; inside = r.inside;
212 return *this;
213 }
214 public:
216 EInside inside;
217 };
218
219 class LastVector // last SurfaceNormal result
220 {
221 public:
222 LastVector()
223 {
224 p.set(kInfinity,kInfinity,kInfinity);
225 vec.set(kInfinity,kInfinity,kInfinity);
226 surface = new G4VTwistSurface*[1];
227 }
228 ~LastVector()
229 {
230 delete [] surface;
231 }
232 LastVector(const LastVector& r) : p(r.p), vec(r.vec)
233 {
234 surface = new G4VTwistSurface*[1];
235 surface[0] = r.surface[0];
236 }
237 LastVector& operator=(const LastVector& r)
238 {
239 if (&r == this) { return *this; }
240 p = r.p; vec = r.vec;
241 delete [] surface; surface = new G4VTwistSurface*[1];
242 surface[0] = r.surface[0];
243 return *this;
244 }
245 public:
247 G4ThreeVector vec;
248 G4VTwistSurface **surface;
249 };
250
251 class LastValue // last G4double value
252 {
253 public:
254 LastValue()
255 {
256 p.set(kInfinity,kInfinity,kInfinity);
257 value = DBL_MAX;
258 }
259 ~LastValue(){}
260 LastValue(const LastValue& r) : p(r.p), value(r.value){}
261 LastValue& operator=(const LastValue& r)
262 {
263 if (this == &r) { return *this; }
264 p = r.p; value = r.value;
265 return *this;
266 }
267 public:
269 G4double value;
270 };
271
272 class LastValueWithDoubleVector // last G4double value
273 {
274 public:
275 LastValueWithDoubleVector()
276 {
277 p.set(kInfinity,kInfinity,kInfinity);
278 vec.set(kInfinity,kInfinity,kInfinity);
279 value = DBL_MAX;
280 }
281 ~LastValueWithDoubleVector(){}
282 LastValueWithDoubleVector(const LastValueWithDoubleVector& r)
283 : p(r.p), vec(r.vec), value(r.value){}
284 LastValueWithDoubleVector& operator=(const LastValueWithDoubleVector& r)
285 {
286 if (this == &r) { return *this; }
287 p = r.p; vec = r.vec; value = r.value;
288 return *this;
289 }
290 public:
292 G4ThreeVector vec;
293 G4double value;
294 };
295
296 LastState fLastInside;
297 LastVector fLastNormal;
298 LastValue fLastDistanceToIn;
299 LastValue fLastDistanceToOut;
300 LastValueWithDoubleVector fLastDistanceToInWithV;
301 LastValueWithDoubleVector fLastDistanceToOutWithV;
302 };
303
304//=====================================================================
305
306inline
308{
309 return ( fDx4 + fDx2 + ( fDx4 - fDx2 ) * ( 2 * phi ) / fPhiTwist ) ;
310}
311
312inline
314{
315 return ( fDx3 + fDx1 + ( fDx3 - fDx1 ) * ( 2 * phi ) / fPhiTwist ) ;
316}
317
318inline
320{
321 return ( fDy2 + fDy1 + ( fDy2 - fDy1 ) * ( 2 * phi ) / fPhiTwist ) ;
322}
323
324inline
326{
327 return GetValueA(phi)/2. + (GetValueD(phi)-GetValueA(phi))/4.
328 - u*( ( GetValueD(phi)-GetValueA(phi) ) / ( 2 * GetValueB(phi) ) - ftg );
329}
330
331#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void set(double x, double y, double z)
G4double GetValueD(G4double phi) const
G4double GetDy1() const
virtual G4double GetCubicVolume()
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
G4double GetValueA(G4double phi) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4double GetDy2() const
G4double GetTheta() const
G4VTwistedFaceted & operator=(const G4VTwistedFaceted &rhs)
G4Polyhedron * fpPolyhedron
G4double GetPhi() const
virtual G4Polyhedron * GetPolyhedron() const
G4ThreeVector GetPointOnSurface() const
virtual G4GeometryType GetEntityType() const
G4double GetDx3() const
G4double GetTwistAngle() const
G4ThreeVector GetPointInSolid(G4double z) const
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
G4double GetDx4() const
G4double GetAlpha() const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=false, G4bool *validnorm=nullptr, G4ThreeVector *n=nullptr) const
virtual EInside Inside(const G4ThreeVector &p) const
virtual std::ostream & StreamInfo(std::ostream &os) const
virtual void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *)
G4double GetDz() const
virtual G4VisExtent GetExtent() const
virtual G4Polyhedron * CreatePolyhedron() const
G4double GetDx1() const
G4double GetDx2() const
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const
G4double GetValueB(G4double phi) const
G4double Xcoef(G4double u, G4double phi, G4double ftg) const
virtual G4double GetSurfaceArea()
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
@ kOutside
Definition: geomdefs.hh:68
yystype & operator=(const yystype &right)
Definition: G4UItokenNum.hh:78
#define DBL_MAX
Definition: templates.hh:62