source: trunk/source/geometry/solids/specific/include/G4VTwistedFaceted.hh @ 850

Last change on this file since 850 was 850, checked in by garnier, 16 years ago

geant4.8.2 beta

File size: 10.3 KB
Line 
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,v 1.10 2006/10/20 13:45:20 gcosmo Exp $
28// GEANT4 tag $Name: HEAD $
29//
30// --------------------------------------------------------------------
31// GEANT 4 class header file
32//
33//
34// G4VTwistedFaceted
35//
36// Class description:
37//
38//  G4VTwistedFaceted is an abstract base class for twisted boxoids:
39//  G4TwistedTrd, G4TwistedTrap and G4TwistedBox
40
41// Author:
42//
43//   27-Oct-2004 - O.Link (Oliver.Link@cern.ch)
44//
45// --------------------------------------------------------------------
46
47#ifndef __G4VTWISTEDFACETED__
48#define __G4VTWISTEDFACETED__
49
50#include "G4VSolid.hh"
51#include "G4TwistTrapAlphaSide.hh"
52#include "G4TwistTrapParallelSide.hh"
53#include "G4TwistBoxSide.hh"
54#include "G4TwistTrapFlatSide.hh"
55
56class G4SolidExtentList;
57class G4ClippablePolygon;
58
59class G4VTwistedFaceted: public G4VSolid
60{
61 public:  // with description
62 
63  G4VTwistedFaceted(const G4String &pname,    // Name of instance
64                          G4double PhiTwist,  // twist angle
65                          G4double pDz,       // half z lenght
66                          G4double pTheta,  // direction between end planes
67                          G4double pPhi,    // defined by polar & azim. angles
68                          G4double pDy1,    // half y length at -pDz
69                          G4double pDx1,    // half x length at -pDz,-pDy
70                          G4double pDx2,    // half x length at -pDz,+pDy
71                          G4double pDy2,    // half y length at +pDz
72                          G4double pDx3,    // half x length at +pDz,-pDy
73                          G4double pDx4,    // half x length at +pDz,+pDy
74                          G4double pAlph    // tilt angle at +pDz
75                   );
76 
77  virtual ~G4VTwistedFaceted();
78             
79  virtual void ComputeDimensions(G4VPVParameterisation*,
80                                 const G4int,
81                                 const G4VPhysicalVolume*  );
82 
83  virtual G4bool CalculateExtent(const EAxis               pAxis,
84                                 const G4VoxelLimits      &pVoxelLimit,
85                                 const G4AffineTransform  &pTransform,
86                                       G4double           &pMin,
87                                       G4double           &pMax ) const;
88
89  virtual G4double DistanceToIn (const G4ThreeVector &p,
90                                 const G4ThreeVector &v ) const;
91
92  virtual G4double DistanceToIn (const G4ThreeVector &p ) const;
93   
94  virtual G4double DistanceToOut(const G4ThreeVector &p, 
95                                 const G4ThreeVector &v,
96                                 const G4bool         calcnorm  = false,
97                                       G4bool        *validnorm = 0, 
98                                       G4ThreeVector *n=0 ) const;
99
100  virtual G4double DistanceToOut(const G4ThreeVector &p) const;
101 
102  virtual EInside Inside (const G4ThreeVector &p) const;
103
104  virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const;
105
106  G4ThreeVector GetPointOnSurface() const;
107  G4ThreeVector GetPointInSolid(G4double z) const;
108 
109  virtual inline G4double GetCubicVolume() ;
110  virtual inline G4double GetSurfaceArea() ;
111
112  virtual void            DescribeYourselfTo (G4VGraphicsScene &scene) const;
113  virtual G4Polyhedron   *CreatePolyhedron   () const ;
114  virtual G4NURBS        *CreateNURBS        () 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
151 protected:  // with description
152
153  G4ThreeVectorList*
154    CreateRotatedVertices(const G4AffineTransform& pTransform) const;
155      // Create the List of transformed vertices in the format required
156      // for G4VSolid:: ClipCrossSection and ClipBetweenSections.
157
158 private:
159 
160  void CreateSurfaces();
161
162 private:
163 
164  G4double fTheta;   
165  G4double fPhi ;
166
167  G4double fDy1;   
168  G4double fDx1;     
169  G4double fDx2;     
170 
171  G4double fDy2;   
172  G4double fDx3;     
173  G4double fDx4;     
174
175  G4double fDz;        // Half-length along the z axis
176
177  G4double fDx ;       // maximum side in x
178  G4double fDy ;       // maximum side in y
179
180  G4double fAlph ;
181  G4double fTAlph ;    // std::tan(fAlph)
182
183  G4double fdeltaX ;
184  G4double fdeltaY ;
185   
186  G4double fPhiTwist;  // twist angle ( dphi in surface equation)
187
188  G4double fAngleSide;
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  mutable G4Polyhedron* fpPolyhedron;  // pointer to polyhedron for vis
202
203  class LastState              // last Inside result
204  {
205    public:
206      LastState()
207      {
208        p.set(kInfinity,kInfinity,kInfinity);
209        inside = kOutside;
210      }
211      ~LastState(){}
212    public:
213      G4ThreeVector p;
214        EInside       inside;
215  };
216             
217  class LastVector             // last SurfaceNormal result
218  {
219    public:
220      LastVector()
221      {
222        p.set(kInfinity,kInfinity,kInfinity);
223        vec.set(kInfinity,kInfinity,kInfinity);
224        surface = new G4VTwistSurface*[1];
225      }
226      ~LastVector()
227      {
228        delete [] surface;
229      }
230    public:
231      G4ThreeVector   p;
232      G4ThreeVector   vec;
233      G4VTwistSurface    **surface;
234  };
235
236  class LastValue              // last G4double value
237  {
238    public:
239      LastValue()
240      {
241        p.set(kInfinity,kInfinity,kInfinity);
242        value = DBL_MAX;
243      }
244      ~LastValue(){}
245    public:
246      G4ThreeVector p;
247      G4double      value;
248  };
249             
250  class LastValueWithDoubleVector   // last G4double value
251  {
252    public:
253      LastValueWithDoubleVector()
254      {
255        p.set(kInfinity,kInfinity,kInfinity);
256        vec.set(kInfinity,kInfinity,kInfinity);
257        value = DBL_MAX;
258      }
259      ~LastValueWithDoubleVector(){}
260      public:
261        G4ThreeVector p;
262        G4ThreeVector vec;
263        G4double      value;
264  };
265             
266  LastState    fLastInside;
267  LastVector   fLastNormal;
268  LastValue    fLastDistanceToIn;
269  LastValue    fLastDistanceToOut;
270  LastValueWithDoubleVector   fLastDistanceToInWithV;
271  LastValueWithDoubleVector   fLastDistanceToOutWithV;
272
273 };
274
275//=====================================================================
276
277inline
278G4double G4VTwistedFaceted::GetCubicVolume()
279{
280  if(fCubicVolume != 0.) ;
281  else   fCubicVolume = 2 * fDz
282                      * ( ( fDx1 + fDx2 ) * fDy1 + ( fDx3 + fDx4 ) * fDy2  );
283  return fCubicVolume;
284}
285
286inline
287G4double G4VTwistedFaceted::GetSurfaceArea()
288{
289  if(fSurfaceArea != 0.) ;
290  else   fSurfaceArea = G4VSolid::GetSurfaceArea();
291  return fSurfaceArea;
292}
293
294inline
295G4double G4VTwistedFaceted::GetValueA(G4double phi) const
296{
297  return ( fDx4 + fDx2  + ( fDx4 - fDx2 ) * ( 2 * phi ) / fPhiTwist  ) ;
298}
299
300inline
301G4double G4VTwistedFaceted::GetValueD(G4double phi) const
302{
303  return ( fDx3 + fDx1  + ( fDx3 - fDx1 ) * ( 2 * phi ) / fPhiTwist  ) ;
304} 
305
306inline 
307G4double G4VTwistedFaceted::GetValueB(G4double phi) const
308{
309  return ( fDy2 + fDy1  + ( fDy2 - fDy1 ) * ( 2 * phi ) / fPhiTwist ) ;
310}
311
312inline
313G4double G4VTwistedFaceted::Xcoef(G4double u, G4double phi, G4double ftg) const 
314{
315  return GetValueA(phi)/2. + (GetValueD(phi)-GetValueA(phi))/4. 
316    - u*( ( GetValueD(phi)-GetValueA(phi) ) / ( 2 * GetValueB(phi) ) - ftg );
317}
318
319#endif
Note: See TracBrowser for help on using the repository browser.