Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TwistTrapParallelSide.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// G4TwistTrapParallelSide
27//
28// Class description:
29//
30// Class describing a twisted boundary surface for a trapezoid.
31
32// Author: Oliver Link (CERN), 27.10.2004 - Created
33// --------------------------------------------------------------------
34#ifndef G4TWISTTRAPPARALLELSIDE_HH
35#define G4TWISTTRAPPARALLELSIDE_HH
36
37#include "G4VTwistSurface.hh"
38
39#include <vector>
40
41/**
42 * @brief G4TwistTrapParallelSide describes a twisted boundary surface for
43 * a trapezoid.
44 */
45
47{
48 public:
49
50 /**
51 * Constructs a trapezoid twisted boundary surface, given its parameters.
52 * @param[in] name The surface name.
53 * @param[in] PhiTwist The twist angle.
54 * @param[in] pDz Half z length.
55 * @param[in] pTheta Direction between end planes - polar angle.
56 * @param[in] pPhi Direction between end planes - azimuthal angle.
57 * @param[in] pDy1 Half y length at -pDz.
58 * @param[in] pDx1 Half x length at -pDz,-pDy.
59 * @param[in] pDx2 Half x length at -pDz,+pDy.
60 * @param[in] pDy2 Half y length at +pDz.
61 * @param[in] pDx3 Half x length at +pDz,-pDy.
62 * @param[in] pDx4 Half x length at +pDz,+pDy.
63 * @param[in] pAlph Tilt angle at +pDz.
64 * @param[in] AngleSide Parity.
65 */
67 G4double PhiTwist, // twist angle
68 G4double pDz, // half z lenght
69 G4double pTheta, // direction between end planes
70 G4double pPhi, // by polar and azimutal angles
71 G4double pDy1, // half y length at -pDz
72 G4double pDx1, // half x length at -pDz,-pDy
73 G4double pDx2, // half x length at -pDz,+pDy
74 G4double pDy2, // half y length at +pDz
75 G4double pDx3, // half x length at +pDz,-pDy
76 G4double pDx4, // half x length at +pDz,+pDy
77 G4double pAlph, // tilt angle at +pDz
78 G4double AngleSide // parity
79 );
80
81 /**
82 * Default destructor.
83 */
84 ~G4TwistTrapParallelSide() override = default;
85
86 /**
87 * Returns a normal vector at a surface (or very close to the surface)
88 * point at 'p'.
89 * @param[in] p The point where computing the normal.
90 * @param[in] isGlobal If true, it returns the normal in global coordinates.
91 * @returns The normal vector.
92 */
94 G4bool isGlobal = false) override ;
95
96 /**
97 * Returns the distance to surface, given point 'gp' and direction 'gv'.
98 * @param[in] gp The point from where computing the distance.
99 * @param[in] gv The direction along which computing the distance.
100 * @param[out] gxx Vector of global points based on number of solutions.
101 * @param[out] distance The distance vector based on number of solutions.
102 * @param[out] areacode The location vector based on number of solutions.
103 * @param[out] isvalid Validity vector based on number of solutions.
104 * @param[in] validate Adopted validation criteria.
105 * @returns The number of solutions.
106 */
108 const G4ThreeVector& gv,
109 G4ThreeVector gxx[],
110 G4double distance[],
111 G4int areacode[],
112 G4bool isvalid[],
113 EValidate validate = kValidateWithTol) override;
114
115 /**
116 * Returns the safety distance to surface, given point 'gp'.
117 * @param[in] gp The point from where computing the safety distance.
118 * @param[out] gxx Vector of global points based on number of solutions.
119 * @param[out] distance The distance vector based on number of solutions.
120 * @param[out] areacode The location vector based on number of solutions.
121 * @returns The number of solutions.
122 */
124 G4ThreeVector gxx[],
125 G4double distance[],
126 G4int areacode[]) override;
127
128 /**
129 * Fake default constructor for usage restricted to direct object
130 * persistency for clients requiring preallocation of memory for
131 * persistifiable objects.
132 */
133 G4TwistTrapParallelSide(__void__&);
134
135 private:
136
137 /**
138 * Returns the area code for point 'xx' using or not surface tolerance.
139 */
140 G4int GetAreaCode(const G4ThreeVector& xx,
141 G4bool withTol = true) override;
142
143 /**
144 * Setters.
145 */
146 void SetCorners() override;
147 void SetBoundaries() override;
148
149 /**
150 * Finds the closest point on surface for a given point 'p', returning
151 * 'phi' and 'u'.
152 */
153 void GetPhiUAtX(const G4ThreeVector& p, G4double& phi, G4double& u);
154
155 /**
156 * Returns projection on surface of a given point 'p'.
157 */
158 G4ThreeVector ProjectPoint(const G4ThreeVector& p,
159 G4bool isglobal = false);
160
161 /**
162 * Returns point on surface given 'phi' and 'u'.
163 */
164 inline G4ThreeVector SurfacePoint(G4double phi, G4double u,
165 G4bool isGlobal = false) override;
166
167 /**
168 * Internal accessors.
169 */
170 inline G4double GetBoundaryMin(G4double phi) override;
171 inline G4double GetBoundaryMax(G4double phi) override;
172 inline G4double GetSurfaceArea() override;
173 void GetFacets( G4int m, G4int n, G4double xyz[][3],
174 G4int faces[][4], G4int iside ) override;
175
176 inline G4ThreeVector NormAng(G4double phi, G4double u);
177 inline G4double GetValueB(G4double phi) ;
178 inline G4double Xcoef(G4double u);
179 // To calculate the w(u) function
180
181 private:
182
183 G4double fTheta;
184 G4double fPhi ;
185
186 G4double fDy1;
187 G4double fDx1;
188 G4double fDx2;
189
190 G4double fDy2;
191 G4double fDx3;
192 G4double fDx4;
193
194 G4double fDz; // Half-length along the z axis
195
196 G4double fAlph;
197 G4double fTAlph; // std::tan(fAlph)
198
199 G4double fPhiTwist; // twist angle ( dphi in surface equation)
200
201 G4double fAngleSide;
202
203 G4double fdeltaX;
204 G4double fdeltaY;
205
206 G4double fDx4plus2; // fDx4 + fDx2 == a2/2 + a1/2
207 G4double fDx4minus2; // fDx4 - fDx2 -
208 G4double fDx3plus1; // fDx3 + fDx1 == d2/2 + d1/2
209 G4double fDx3minus1; // fDx3 - fDx1 -
210 G4double fDy2plus1; // fDy2 + fDy1 == b2/2 + b1/2
211 G4double fDy2minus1; // fDy2 - fDy1 -
212 G4double fa1md1; // 2 fDx2 - 2 fDx1 == a1 - d1
213 G4double fa2md2; // 2 fDx4 - 2 fDx3
214};
215
216//========================================================
217// inline functions
218//========================================================
219
220#include "G4TwistTrapParallelSide.icc"
221
222#endif
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
~G4TwistTrapParallelSide() override=default
G4int DistanceToSurface(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector gxx[], G4double distance[], G4int areacode[], G4bool isvalid[], EValidate validate=kValidateWithTol) override
G4TwistTrapParallelSide(const G4String &name, G4double PhiTwist, G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlph, G4double AngleSide)
G4ThreeVector GetNormal(const G4ThreeVector &p, G4bool isGlobal=false) override
G4VTwistSurface(const G4String &name)
virtual G4double GetSurfaceArea()=0