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