Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Trd.hh
이 파일의 문서화 페이지로 가기
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 104894 2017-06-26 13:30:00Z 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 #include "G4Polyhedron.hh"
71 
72 class G4Trd : public G4CSGSolid
73 {
74  public: // with description
75 
76  G4Trd( const G4String& pName,
77  G4double pdx1, G4double pdx2,
78  G4double pdy1, G4double pdy2,
79  G4double pdz );
80  //
81  // Constructs a trapezoid with name, and half lengths
82 
83  ~G4Trd();
84  //
85  // Destructor
86 
87  // Accessors
88 
89  inline G4double GetXHalfLength1() const;
90  inline G4double GetXHalfLength2() const;
91  inline G4double GetYHalfLength1() const;
92  inline G4double GetYHalfLength2() const;
93  inline G4double GetZHalfLength() const;
94 
95  // Modifiers
96 
97  inline void SetXHalfLength1(G4double val);
98  inline void SetXHalfLength2(G4double val);
99  inline void SetYHalfLength1(G4double val);
100  inline void SetYHalfLength2(G4double val);
101  inline void SetZHalfLength(G4double val);
102 
103  void SetAllParameters ( G4double pdx1, G4double pdx2,
104  G4double pdy1, G4double pdy2,
105  G4double pdz );
106 
107  // Methods of solid
108 
111 
113  const G4int n,
114  const G4VPhysicalVolume* pRep );
115 
116  void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
117 
118  G4bool CalculateExtent( const EAxis pAxis,
119  const G4VoxelLimits& pVoxelLimit,
120  const G4AffineTransform& pTransform,
121  G4double& pMin, G4double& pMax ) const;
122 
123  EInside Inside( const G4ThreeVector& p ) const;
124 
125  G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const;
126 
128  const G4ThreeVector& v ) const;
129 
130  G4double DistanceToIn( const G4ThreeVector& p ) const;
131 
133  const G4ThreeVector& v,
134  const G4bool calcNorm=false,
135  G4bool *validNorm=0,
136  G4ThreeVector *n=0 ) const;
137 
138  G4double DistanceToOut( const G4ThreeVector& p ) const;
139 
141 
143 
144  G4VSolid* Clone() const;
145 
146  std::ostream& StreamInfo( std::ostream& os ) const;
147 
148  // Visualisation functions
149 
150  void DescribeYourselfTo (G4VGraphicsScene& scene) const;
151  G4Polyhedron* CreatePolyhedron () const;
152 
153  public: // without description
154 
155  G4Trd(__void__&);
156  // Fake default constructor for usage restricted to direct object
157  // persistency for clients requiring preallocation of memory for
158  // persistifiable objects.
159 
160  G4Trd(const G4Trd& rhs);
161  G4Trd& operator=(const G4Trd& rhs);
162  // Copy constructor and assignment operator
163 
164  private:
165 
166  void CheckParameters();
167  // Check parameters
168 
169  void MakePlanes();
170  // Set side planes
171 
173  // Algorithm for SurfaceNormal() following the original
174  // specification for points not on the surface
175 
176  private:
177 
180  struct { G4double a,b,c,d; } fPlanes[4];
181 };
182 
183 #include "G4Trd.icc"
184 
185 #endif
186 
187 #endif
G4double fDy2
Definition: G4Trd.hh:179
void SetXHalfLength2(G4double val)
G4double a
Definition: G4Trd.hh:180
void CheckParameters()
Definition: G4Trd.cc:159
void SetYHalfLength1(G4double val)
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Trd.cc:260
void SetYHalfLength2(G4double val)
G4double halfCarTolerance
Definition: G4Trd.hh:178
G4double GetCubicVolume()
Definition: G4Trd.cc:218
const char * p
Definition: xmltok.h:285
void SetXHalfLength1(G4double val)
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Trd.cc:440
G4double fDz
Definition: G4Trd.hh:179
G4double GetSurfaceArea()
Definition: G4Trd.cc:232
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Trd.cc:575
G4double fDy1
Definition: G4Trd.hh:179
G4Trd(const G4String &pName, G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:66
G4double GetYHalfLength2() const
G4VSolid * Clone() const
Definition: G4Trd.cc:698
G4double GetZHalfLength() const
G4double GetXHalfLength1() const
EInside Inside(const G4ThreeVector &p) const
Definition: G4Trd.cc:342
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Trd.cc:707
void SetZHalfLength(G4double val)
void MakePlanes()
Definition: G4Trd.cc:181
G4double GetYHalfLength1() const
~G4Trd()
Definition: G4Trd.cc:93
G4Polyhedron * CreatePolyhedron() const
Definition: G4Trd.cc:793
G4double fDx1
Definition: G4Trd.hh:179
int G4int
Definition: G4Types.hh:78
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
Definition: G4Trd.hh:72
G4double GetXHalfLength2() const
G4GeometryType GetEntityType() const
Definition: G4Trd.cc:689
struct G4Trd::@34 fPlanes[4]
Char_t n[5]
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Trd.cc:359
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Trd.cc:249
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Trd.cc:788
G4double fDx2
Definition: G4Trd.hh:179
G4double b
Definition: G4Trd.hh:180
G4double c
Definition: G4Trd.hh:180
void SetAllParameters(G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:138
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Trd.cc:464
G4Trd & operator=(const G4Trd &rhs)
Definition: G4Trd.cc:113
G4ThreeVector GetPointOnSurface() const
Definition: G4Trd.cc:730
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
Definition: G4Trd.cc:291
G4double d
Definition: G4Trd.hh:180