Geant4  10.01.p03
G4Paraboloid.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: G4Paraboloid.hh 83572 2014-09-01 15:23:27Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 // G4Paraboloid
34 //
35 // Class description:
36 //
37 // A G4Paraboloid represents a solid with parabolic profile with possible
38 // cuts along the Z axis.
39 //
40 // Member Data:
41 //
42 // dz z half lenght
43 // r1 radius at -dz
44 // r2 radius at dz
45 // r2 > r1
46 //
47 // Equation for the solid:
48 // rho^2 <= k1 * z + k2;
49 // -dz <= z <= dz
50 // r1^2 = k1 * (-dz) + k2
51 // r2^2 = k1 * ( dz) + k2
52 
53 // History:
54 // -------
55 // 10.07.2007 L.Lindroos (CERN) - First implementation
56 // --------------------------------------------------------------------
57 #ifndef G4Paraboloid_HH
58 #define G4Paraboloid_HH
59 
60 #include <CLHEP/Units/PhysicalConstants.h>
61 
62 #include "G4VSolid.hh"
63 #include "G4Polyhedron.hh"
64 
65 class G4Paraboloid : public G4VSolid
66 {
67  public: // with description
68 
69  G4Paraboloid(const G4String& pName,
70  G4double pDz,
71  G4double pR1,
72  G4double pR2);
73 
74  virtual ~G4Paraboloid();
75 
76  // Access functions
77 
78  inline G4double GetZHalfLength() const;
79  inline G4double GetRadiusMinusZ() const;
80  inline G4double GetRadiusPlusZ() const;
81 
82  inline G4double GetCubicVolume();
83  inline G4double GetSurfaceArea();
84  inline G4double CalculateSurfaceArea() const;
85 
86  // Modifiers functions
87 
88  inline void SetZHalfLength(G4double dz);
89  inline void SetRadiusMinusZ(G4double R1);
90  inline void SetRadiusPlusZ(G4double R2);
91 
92  // Solid standard methods
93 
94  //void ComputeDimensions( G4VPVParamerisation p,
95  // const G4Int n,
96  // const G4VPhysicalVolume* pRep );
97  G4bool CalculateExtent(const EAxis pAxis,
98  const G4VoxelLimits& pVoxelLimit,
99  const G4AffineTransform& pTransform,
100  G4double& pmin, G4double& pmax) const;
101  EInside Inside(const G4ThreeVector& p) const;
102  G4ThreeVector SurfaceNormal( const G4ThreeVector& p) const;
104  const G4ThreeVector& v) const;
105  G4double DistanceToIn(const G4ThreeVector& p) const;
107  const G4ThreeVector& v,
108  const G4bool calcNorm=G4bool(false),
109  G4bool *validNorm=0,
110  G4ThreeVector *n=0) const;
111  G4double DistanceToOut(const G4ThreeVector& p) const;
112 
114 
115  G4VSolid* Clone() const;
116 
117  std::ostream& StreamInfo(std::ostream& os) const;
118 
120 
121  // Visualisation functions
122 
123  void DescribeYourselfTo(G4VGraphicsScene& scene) const;
125  G4Polyhedron* GetPolyhedron () const;
126 
127  public: // without description
128 
129  G4Paraboloid(__void__&);
130  // Fake default constructor for usage restricted to direct object
131  // persistency for clients requiring preallocation of memory for
132  // persistifiable objects.
133 
134  G4Paraboloid(const G4Paraboloid& rhs);
135  G4Paraboloid& operator=(const G4Paraboloid& rhs);
136  // Copy constructor and assignment operator.
137 
138  protected: // without description
139 
142 
143  private:
144 
146  CreateRotatedVertices(const G4AffineTransform& pTransform,
147  G4int& noPolygonVertices) const;
148  // Making this mutable to allow GetPointOnSurface to have access to
149  // area function.
152 
155  // Defined to make some calculations easier to follow
156 };
157 
158 #include "G4Paraboloid.icc"
159 
160 #endif
CLHEP::Hep3Vector G4ThreeVector
G4Polyhedron * fpPolyhedron
G4Paraboloid & operator=(const G4Paraboloid &rhs)
G4double GetSurfaceArea()
G4double GetCubicVolume()
G4GeometryType GetEntityType() const
G4double GetRadiusMinusZ() const
int G4int
Definition: G4Types.hh:78
G4double fCubicVolume
G4double CalculateSurfaceArea() const
G4double fSurfaceArea
void SetRadiusPlusZ(G4double R2)
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void SetRadiusMinusZ(G4double R1)
G4double GetZHalfLength() const
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
const G4int n
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
G4VSolid * Clone() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
EInside Inside(const G4ThreeVector &p) const
std::ostream & StreamInfo(std::ostream &os) const
G4double GetRadiusPlusZ() const
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
virtual ~G4Paraboloid()
void SetZHalfLength(G4double dz)
G4ThreeVector GetPointOnSurface() const
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform, G4int &noPolygonVertices) const
double G4double
Definition: G4Types.hh:76
G4Polyhedron * GetPolyhedron() const
G4Paraboloid(const G4String &pName, G4double pDz, G4double pR1, G4double pR2)
Definition: G4Paraboloid.cc:63
G4Polyhedron * CreatePolyhedron() const
G4bool fRebuildPolyhedron