Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 //
27 // $Id: G4VTwistedFaceted.hh 99781 2016-10-05 10:18:54Z gcosmo $
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 (Oliver.Link@cern.ch)
43 //
44 // --------------------------------------------------------------------
45 
46 #ifndef __G4VTWISTEDFACETED__
47 #define __G4VTWISTEDFACETED__
48 
49 #include "G4VSolid.hh"
50 #include "G4TwistTrapAlphaSide.hh"
52 #include "G4TwistBoxSide.hh"
53 #include "G4TwistTrapFlatSide.hh"
54 
55 class G4SolidExtentList;
56 class G4ClippablePolygon;
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 void Extent(G4ThreeVector &pMin, G4ThreeVector &pMax) const;
83 
84  virtual G4bool CalculateExtent(const EAxis pAxis,
85  const G4VoxelLimits &pVoxelLimit,
86  const G4AffineTransform &pTransform,
87  G4double &pMin,
88  G4double &pMax ) const;
89 
90  virtual G4double DistanceToIn (const G4ThreeVector &p,
91  const G4ThreeVector &v ) const;
92 
93  virtual G4double DistanceToIn (const G4ThreeVector &p ) const;
94 
95  virtual G4double DistanceToOut(const G4ThreeVector &p,
96  const G4ThreeVector &v,
97  const G4bool calcnorm = false,
98  G4bool *validnorm = 0,
99  G4ThreeVector *n=0 ) const;
100 
101  virtual G4double DistanceToOut(const G4ThreeVector &p) const;
102 
103  virtual EInside Inside (const G4ThreeVector &p) const;
104 
105  virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const;
106 
109 
110  virtual inline G4double GetCubicVolume() ;
111  virtual inline G4double GetSurfaceArea() ;
112 
113  virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const;
114  virtual G4Polyhedron *CreatePolyhedron () const ;
115  virtual G4Polyhedron *GetPolyhedron () const;
116 
117  virtual std::ostream &StreamInfo(std::ostream& os) const;
118 
119  // accessors
120 
121  inline G4double GetTwistAngle () const { return fPhiTwist; }
122 
123  inline G4double GetDx1 () const { return fDx1 ; }
124  inline G4double GetDx2 () const { return fDx2 ; }
125  inline G4double GetDx3 () const { return fDx3 ; }
126  inline G4double GetDx4 () const { return fDx4 ; }
127  inline G4double GetDy1 () const { return fDy1 ; }
128  inline G4double GetDy2 () const { return fDy2 ; }
129  inline G4double GetDz () const { return fDz ; }
130  inline G4double GetPhi () const { return fPhi ; }
131  inline G4double GetTheta () const { return fTheta ; }
132  inline G4double GetAlpha () const { return fAlph ; }
133 
134  inline G4double Xcoef(G4double u,G4double phi, G4double ftg) const ;
135  // For calculating the w(u) function
136 
137  inline G4double GetValueA(G4double phi) const;
138  inline G4double GetValueB(G4double phi) const;
139  inline G4double GetValueD(G4double phi) const;
140 
141  virtual G4VisExtent GetExtent () const;
142  virtual G4GeometryType GetEntityType() const;
143 
144  public: // without description
145 
146  G4VTwistedFaceted(__void__&);
147  // Fake default constructor for usage restricted to direct object
148  // persistency for clients requiring preallocation of memory for
149  // persistifiable objects.
150 
153  // Copy constructor and assignment operator.
154 
155  protected: // with description
156 
158  mutable G4Polyhedron* fpPolyhedron; // pointer to polyhedron for vis
159 
160  private:
161 
162  void CreateSurfaces();
163 
164  private:
165 
166  G4double fTheta;
167  G4double fPhi ;
168 
169  G4double fDy1;
170  G4double fDx1;
171  G4double fDx2;
172 
173  G4double fDy2;
174  G4double fDx3;
175  G4double fDx4;
176 
177  G4double fDz; // Half-length along the z axis
178 
179  G4double fDx ; // maximum side in x
180  G4double fDy ; // maximum side in y
181 
182  G4double fAlph ;
183  G4double fTAlph ; // std::tan(fAlph)
184 
185  G4double fdeltaX ;
186  G4double fdeltaY ;
187 
188  G4double fPhiTwist; // twist angle ( dphi in surface equation)
189 
190  G4VTwistSurface *fLowerEndcap ; // surface of -ve z
191  G4VTwistSurface *fUpperEndcap ; // surface of +ve z
192 
193  G4VTwistSurface *fSide0 ; // Twisted Side at phi = 0 deg
194  G4VTwistSurface *fSide90 ; // Twisted Side at phi = 90 deg
195  G4VTwistSurface *fSide180 ; // Twisted Side at phi = 180 deg
196  G4VTwistSurface *fSide270 ; // Twisted Side at phi = 270 deg
197 
198  G4double fCubicVolume ; // volume of the solid
199  G4double fSurfaceArea ; // area of the solid
200 
201  class LastState // last Inside result
202  {
203  public:
204  LastState()
205  {
207  }
208  ~LastState(){}
209  LastState(const LastState& r) : p(r.p), inside(r.inside){}
210  LastState& operator=(const LastState& r)
211  {
212  if (this == &r) { return *this; }
213  p = r.p; inside = r.inside;
214  return *this;
215  }
216  public:
218  EInside inside;
219  };
220 
221  class LastVector // last SurfaceNormal result
222  {
223  public:
224  LastVector()
225  {
227  vec.set(kInfinity,kInfinity,kInfinity);
228  surface = new G4VTwistSurface*[1];
229  }
230  ~LastVector()
231  {
232  delete [] surface;
233  }
234  LastVector(const LastVector& r) : p(r.p), vec(r.vec)
235  {
236  surface = new G4VTwistSurface*[1];
237  surface[0] = r.surface[0];
238  }
239  LastVector& operator=(const LastVector& r)
240  {
241  if (&r == this) { return *this; }
242  p = r.p; vec = r.vec;
243  delete [] surface; surface = new G4VTwistSurface*[1];
244  surface[0] = r.surface[0];
245  return *this;
246  }
247  public:
249  G4ThreeVector vec;
250  G4VTwistSurface **surface;
251  };
252 
253  class LastValue // last G4double value
254  {
255  public:
256  LastValue()
257  {
259  value = DBL_MAX;
260  }
261  ~LastValue(){}
262  LastValue(const LastValue& r) : p(r.p), value(r.value){}
263  LastValue& operator=(const LastValue& r)
264  {
265  if (this == &r) { return *this; }
266  p = r.p; value = r.value;
267  return *this;
268  }
269  public:
271  G4double value;
272  };
273 
274  class LastValueWithDoubleVector // last G4double value
275  {
276  public:
277  LastValueWithDoubleVector()
278  {
280  vec.set(kInfinity,kInfinity,kInfinity);
281  value = DBL_MAX;
282  }
283  ~LastValueWithDoubleVector(){}
284  LastValueWithDoubleVector(const LastValueWithDoubleVector& r)
285  : p(r.p), vec(r.vec), value(r.value){}
286  LastValueWithDoubleVector& operator=(const LastValueWithDoubleVector& r)
287  {
288  if (this == &r) { return *this; }
289  p = r.p; vec = r.vec; value = r.value;
290  return *this;
291  }
292  public:
294  G4ThreeVector vec;
295  G4double value;
296  };
297 
298  LastState fLastInside;
299  LastVector fLastNormal;
300  LastValue fLastDistanceToIn;
301  LastValue fLastDistanceToOut;
302  LastValueWithDoubleVector fLastDistanceToInWithV;
303  LastValueWithDoubleVector fLastDistanceToOutWithV;
304 
305  };
306 
307 //=====================================================================
308 
309 inline
311 {
312  if(fCubicVolume != 0.) ;
313  else fCubicVolume = 2 * fDz
314  * ( ( fDx1 + fDx2 ) * fDy1 + ( fDx3 + fDx4 ) * fDy2 );
315  return fCubicVolume;
316 }
317 
318 inline
320 {
321  if(fSurfaceArea != 0.) ;
322  else fSurfaceArea = G4VSolid::GetSurfaceArea();
323  return fSurfaceArea;
324 }
325 
326 inline
328 {
329  return ( fDx4 + fDx2 + ( fDx4 - fDx2 ) * ( 2 * phi ) / fPhiTwist ) ;
330 }
331 
332 inline
334 {
335  return ( fDx3 + fDx1 + ( fDx3 - fDx1 ) * ( 2 * phi ) / fPhiTwist ) ;
336 }
337 
338 inline
340 {
341  return ( fDy2 + fDy1 + ( fDy2 - fDy1 ) * ( 2 * phi ) / fPhiTwist ) ;
342 }
343 
344 inline
346 {
347  return GetValueA(phi)/2. + (GetValueD(phi)-GetValueA(phi))/4.
348  - u*( ( GetValueD(phi)-GetValueA(phi) ) / ( 2 * GetValueB(phi) ) - ftg );
349 }
350 
351 #endif
void set(double x, double y, double z)
virtual void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *)
virtual G4double GetSurfaceArea()
virtual G4double GetCubicVolume()
static const G4double kInfinity
Definition: geomdefs.hh:42
G4VTwistedFaceted(const G4String &pname, G4double PhiTwist, G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlph)
virtual G4GeometryType GetEntityType() const
const char * p
Definition: xmltok.h:285
G4ThreeVector GetPointOnSurface() const
G4double GetDy2() const
G4ThreeVector GetPointInSolid(G4double z) const
virtual G4Polyhedron * CreatePolyhedron() const
G4double GetDy1() const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4double GetDx3() const
virtual void Extent(G4ThreeVector &pMin, G4ThreeVector &pMax) const
G4VTwistedFaceted & operator=(const G4VTwistedFaceted &rhs)
G4double GetTwistAngle() const
int G4int
Definition: G4Types.hh:78
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4double Xcoef(G4double u, G4double phi, G4double ftg) const
G4double GetValueD(G4double phi) const
G4double GetTheta() const
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
const XML_Char int const XML_Char * value
Definition: expat.h:331
bool G4bool
Definition: G4Types.hh:79
G4double GetValueB(G4double phi) const
const G4int n
tuple v
Definition: test.py:18
string pname
Definition: eplot.py:33
G4double GetDx2() const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=false, G4bool *validnorm=0, G4ThreeVector *n=0) const
G4double GetPhi() const
virtual G4Polyhedron * GetPolyhedron() const
EInside
Definition: geomdefs.hh:58
virtual G4VisExtent GetExtent() const
virtual EInside Inside(const G4ThreeVector &p) const
EAxis
Definition: geomdefs.hh:54
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const
tuple z
Definition: test.py:28
G4double GetAlpha() const
G4double GetDx1() const
double G4double
Definition: G4Types.hh:76
virtual std::ostream & StreamInfo(std::ostream &os) const
G4double GetDx4() const
#define DBL_MAX
Definition: templates.hh:83
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:251
G4Polyhedron * fpPolyhedron
G4double GetValueA(G4double phi) const
G4double GetDz() const