Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Trd.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: G4Trd.hh 76263 2013-11-08 11:41:52Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 //
34 // G4Trd
35 //
36 // Class description:
37 //
38 // A G4Trd is a trapezoid with the x and y dimensions varying along z
39 // functions:
40 //
41 // Member Data:
42 //
43 // fDx1 Half-length along x at the surface positioned at -dz
44 // fDx2 Half-length along x at the surface positioned at +dz
45 // fDy1 Half-length along y at the surface positioned at -dz
46 // fDy2 Half-length along y at the surface positioned at +dz
47 // fDz Half-length along z axis
48 
49 // History:
50 // 12.01.95 P.Kent: Old prototype code converted to thick geometry
51 // 17.02.95 P.Kent: Exiting normal return
52 // 19.08.96 P.Kent, V.Grichine: Fs in accordance with G4Box
53 // 21.04.97 J.Apostolakis: Added Set Methods
54 // 19.11.99 V.Grichine: kUndefined was added to Eside enum
55 // --------------------------------------------------------------------
56 
57 #ifndef G4TRD_HH
58 #define G4TRD_HH
59 
60 #if defined(G4GEOM_USE_USOLIDS)
61 #define G4GEOM_USE_UTRD 1
62 #endif
63 
64 #if defined(G4GEOM_USE_UTRD)
65  #define G4UTrd G4Trd
66  #include "G4UTrd.hh"
67 #else
68 
69 #include "G4CSGSolid.hh"
70 
71 class G4Trd : public G4CSGSolid
72 {
73  public: // with description
74 
75  G4Trd( const G4String& pName,
76  G4double pdx1, G4double pdx2,
77  G4double pdy1, G4double pdy2,
78  G4double pdz );
79  //
80  // Constructs a trapezoid with name, and half lengths
81 
82  ~G4Trd();
83  //
84  // Destructor
85 
86  // Accessors
87 
88  inline G4double GetXHalfLength1() const;
89  inline G4double GetXHalfLength2() const;
90  inline G4double GetYHalfLength1() const;
91  inline G4double GetYHalfLength2() const;
92  inline G4double GetZHalfLength() const;
93 
94  // Modifiers
95 
96  inline void SetXHalfLength1(G4double val);
97  inline void SetXHalfLength2(G4double val);
98  inline void SetYHalfLength1(G4double val);
99  inline void SetYHalfLength2(G4double val);
100  inline void SetZHalfLength(G4double val);
101 
102  // Methods of solid
103 
104  inline G4double GetCubicVolume();
105  inline G4double GetSurfaceArea();
106 
108  const G4int n,
109  const G4VPhysicalVolume* pRep );
110 
111  G4bool CalculateExtent( const EAxis pAxis,
112  const G4VoxelLimits& pVoxelLimit,
113  const G4AffineTransform& pTransform,
114  G4double& pMin, G4double& pMax ) const;
115 
116  EInside Inside( const G4ThreeVector& p ) const;
117 
118  G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
119 
121  const G4ThreeVector& v ) const;
122 
123  G4double DistanceToIn( const G4ThreeVector& p ) const;
124 
126  const G4ThreeVector& v,
127  const G4bool calcNorm=false,
128  G4bool *validNorm=0,
129  G4ThreeVector *n=0 ) const;
130 
131  G4double DistanceToOut( const G4ThreeVector& p ) const;
132 
133  void CheckAndSetAllParameters ( G4double pdx1, G4double pdx2,
134  G4double pdy1, G4double pdy2,
135  G4double pdz );
136 
137  void SetAllParameters ( G4double pdx1, G4double pdx2,
138  G4double pdy1, G4double pdy2,
139  G4double pdz );
140 
142 
144 
145  G4VSolid* Clone() const;
146 
147  std::ostream& StreamInfo( std::ostream& os ) const;
148 
149  // Visualisation functions
150 
151  void DescribeYourselfTo (G4VGraphicsScene& scene) const;
152  G4Polyhedron* CreatePolyhedron () const;
153 
154  public: // without description
155 
157  // Codes for faces (kPX=plus x face,kMY= minus y face etc)
158 
159  G4Trd(__void__&);
160  // Fake default constructor for usage restricted to direct object
161  // persistency for clients requiring preallocation of memory for
162  // persistifiable objects.
163 
164  G4Trd(const G4Trd& rhs);
165  G4Trd& operator=(const G4Trd& rhs);
166  // Copy constructor and assignment operator.
167 
169  CreateRotatedVertices( const G4AffineTransform& pTransform ) const;
170  //
171  // Creates the List of transformed vertices in the format required
172  // for G4CSGSolid:: ClipCrossSection and ClipBetweenSections
173 
175  // Algorithm for SurfaceNormal() following the original
176  // specification for points not on the surface
177 
178  private:
179 
180  G4double fDx1,fDx2,fDy1,fDy2,fDz;
181 };
182 
183 #include "G4Trd.icc"
184 
185 #endif
186 
187 #endif
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Trd.cc:1466
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:169
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Trd.cc:1319
~G4Trd()
Definition: G4Trd.cc:128
ESide
Definition: G4Trd.hh:156
G4double GetYHalfLength1() const
void SetXHalfLength2(G4double val)
const char * p
Definition: xmltok.h:285
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Trd.cc:1376
G4GeometryType GetEntityType() const
Definition: G4Trd.cc:1358
Definition: G4Trd.hh:71
void SetZHalfLength(G4double val)
G4double GetZHalfLength() const
int G4int
Definition: G4Types.hh:78
G4VSolid * Clone() const
Definition: G4Trd.cc:1367
G4ThreeVector GetPointOnSurface() const
Definition: G4Trd.cc:1403
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Trd.cc:503
G4double GetXHalfLength2() const
void SetYHalfLength1(G4double val)
G4Polyhedron * CreatePolyhedron() const
Definition: G4Trd.cc:1471
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Trd.cc:596
G4double GetCubicVolume()
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
EInside Inside(const G4ThreeVector &p) const
Definition: G4Trd.cc:371
G4double GetYHalfLength2() const
const G4int n
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Trd.cc:869
void CheckAndSetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:73
G4Trd & operator=(const G4Trd &rhs)
Definition: G4Trd.cc:146
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Trd.cc:181
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Trd.cc:434
G4double GetXHalfLength1() const
void SetXHalfLength1(G4double val)
G4Trd(const G4String &pName, G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:60
void SetYHalfLength2(G4double val)
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
Definition: G4Trd.cc:193
double G4double
Definition: G4Types.hh:76
G4double GetSurfaceArea()