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