Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TwistedTubs.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// G4TwistedTubs
27//
28// Class description:
29//
30// G4TwistedTubs is a sector of a twisted hollow cylinder.
31// A twisted cylinder which is placed along with z-axis and is
32// separated into phi-segments should become a hyperboloid, and
33// its each segmented piece should be tilted with a stereo angle.
34// G4TwistedTubs is a G4VSolid.
35//
36// Details of the implementation: "Development of a Geant4 solid
37// for stereo mini-jet cells in a cylindrical drift chamber",
38// Computer Physics Communications 153 (2003) pp.373–391
39
40// 01-Aug-2002 - Kotoyo Hoshina ([email protected]), created.
41// 13-Nov-2003 - O.Link ([email protected]), Integration in Geant4
42// from original version in Jupiter-2.5.02 application.
43// --------------------------------------------------------------------
44#ifndef G4TWISTEDTUBS_HH
45#define G4TWISTEDTUBS_HH
46
47#include "G4VSolid.hh"
49#include "G4TwistTubsSide.hh"
51
54
55class G4TwistedTubs : public G4VSolid
56{
57 public:
58
59 G4TwistedTubs(const G4String& pname, // Name of instance
60 G4double twistedangle, // Twisted angle
61 G4double endinnerrad, // Inner radius at endcap
62 G4double endouterrad, // Outer radius at endcap
63 G4double halfzlen, // half z length
64 G4double dphi); // Phi angle of a segment
65
66 G4TwistedTubs(const G4String& pname, // Name of instance
67 G4double twistedangle, // Stereo angle
68 G4double endinnerrad, // Inner radius at endcap
69 G4double endouterrad, // Outer radius at endcap
70 G4double halfzlen, // half z length
71 G4int nseg, // Number of segments in totalPhi
72 G4double totphi); // Total angle of all segments
73
74 G4TwistedTubs(const G4String& pname, // Name of instance
75 G4double twistedangle, // Twisted angle
76 G4double innerrad, // Inner radius at z=0
77 G4double outerrad, // Outer radius at z=0
78 G4double negativeEndz, // -ve z endplate
79 G4double positiveEndz, // +ve z endplate
80 G4double dphi); // Phi angle of a segment
81
82 G4TwistedTubs(const G4String& pname, // Name of instance
83 G4double twistedangle, // Stereo angle
84 G4double innerrad, // Inner radius at z=0
85 G4double outerrad, // Outer radius at z=0
86 G4double negativeEndz, // -ve z endplate
87 G4double positiveEndz, // +ve z endplate
88 G4int nseg, // Number of segments in totalPhi
89 G4double totphi); // Total angle of all segments
90
91 ~G4TwistedTubs() override;
92
94 const G4int /* n */ ,
95 const G4VPhysicalVolume* /* prep */ ) override;
96
97 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const override;
98
99 G4bool CalculateExtent(const EAxis pAxis,
100 const G4VoxelLimits& pVoxelLimit,
101 const G4AffineTransform& pTransform,
102 G4double& pMin,
103 G4double& pMax ) const override;
104
106 const G4ThreeVector& v ) const override;
107
108 G4double DistanceToIn (const G4ThreeVector& p ) const override;
109
111 const G4ThreeVector& v,
112 const G4bool calcnorm = false,
113 G4bool* validnorm = nullptr,
114 G4ThreeVector* n = nullptr ) const override;
115
116 G4double DistanceToOut(const G4ThreeVector& p) const override;
117
118 EInside Inside (const G4ThreeVector& p) const override;
119
120 G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const override;
121
122 void DescribeYourselfTo (G4VGraphicsScene& scene) const override;
123 G4Polyhedron* CreatePolyhedron () const override;
124 G4Polyhedron* GetPolyhedron () const override;
125
126 std::ostream &StreamInfo(std::ostream& os) const override;
127
128 // accessors
129
130 inline G4double GetDPhi () const { return fDPhi ; }
131 inline G4double GetPhiTwist () const { return fPhiTwist ; }
132 inline G4double GetInnerRadius () const { return fInnerRadius; }
133 inline G4double GetOuterRadius () const { return fOuterRadius; }
134 inline G4double GetInnerStereo () const { return fInnerStereo; }
135 inline G4double GetOuterStereo () const { return fOuterStereo; }
136 inline G4double GetZHalfLength () const { return fZHalfLength; }
137 inline G4double GetKappa () const { return fKappa ; }
138
139 inline G4double GetTanInnerStereo () const { return fTanInnerStereo ; }
140 inline G4double GetTanInnerStereo2() const { return fTanInnerStereo2 ; }
141 inline G4double GetTanOuterStereo () const { return fTanOuterStereo ; }
142 inline G4double GetTanOuterStereo2() const { return fTanOuterStereo2 ; }
143
144 inline G4double GetEndZ (G4int i) const { return fEndZ[i] ; }
145 inline G4double GetEndPhi (G4int i) const { return fEndPhi[i]; }
147 { return fEndInnerRadius[i]; }
149 { return fEndOuterRadius[i]; }
151 { return (fEndInnerRadius[0] > fEndInnerRadius[1] ?
152 fEndInnerRadius[0] : fEndInnerRadius[1]); }
154 { return (fEndOuterRadius[0] > fEndOuterRadius[1] ?
155 fEndOuterRadius[0] : fEndOuterRadius[1]); }
156
157 G4VisExtent GetExtent () const override;
158 G4GeometryType GetEntityType() const override;
159 G4VSolid* Clone() const override;
160
161 G4double GetCubicVolume() override;
162 // Returns an estimation of the geometrical cubic volume of the
163 // solid. Caches the computed value once computed the first time.
164 G4double GetSurfaceArea() override;
165 // Returns the geometrical surface area of the solid.
166 // Caches the computed value once computed the first time.
167
168 G4ThreeVector GetPointOnSurface() const override ;
169
170 G4TwistedTubs(__void__&);
171 // Fake default constructor for usage restricted to direct object
172 // persistency for clients requiring preallocation of memory for
173 // persistifiable objects.
174
175 G4TwistedTubs(const G4TwistedTubs& rhs);
177 // Copy constructor and assignment operator.
178
179#ifdef G4TWISTDEBUG
180 G4VTwistSurface* GetOuterHype() const { return fOuterHype; }
181#endif
182
183 private:
184
185 inline void SetFields(G4double phitwist, G4double innerrad,
186 G4double outerrad,
187 G4double negativeEndz, G4double positiveEndz);
188 void CreateSurfaces();
189 G4double GetLateralArea(G4double a, G4double r, G4double z) const;
190 G4double GetPhiCutArea(G4double a, G4double r, G4double z) const;
191
192 private:
193
194 G4double fPhiTwist; // Twist angle from -fZHalfLength to fZHalfLength
195 G4double fInnerRadius; // Inner-hype radius at z=0
196 G4double fOuterRadius; // Outer-hype radius at z=0
197 G4double fEndZ[2]; // z at endcaps, [0] = -ve z, [1] = +ve z
198 G4double fDPhi; // Phi-width of a segment fDPhi > 0
199 G4double fZHalfLength; // Half length along z-axis
200
201 G4double fInnerStereo; // Inner-hype stereo angle
202 G4double fOuterStereo; // Outer-hype stereo angle
203 G4double fTanInnerStereo; // std::tan(innerStereoAngle)
204 G4double fTanOuterStereo; // std::tan(outerStereoAngle)
205 G4double fKappa; // std::tan(fPhiTwist/2)/fZHalfLen;
206 G4double fEndInnerRadius[2]; // Inner-hype radii endcaps [0] -ve z, [1] +ve z
207 G4double fEndOuterRadius[2]; // Outer-hype radii endcaps [0] -ve z, [1] +ve z
208 G4double fEndPhi[2]; // Phi endcaps, [0] = -ve z, [1] = +ve z
209
210 G4double fInnerRadius2; // fInnerRadius * fInnerRadius
211 G4double fOuterRadius2; // fOuterRadius * fOuterRadius
212 G4double fTanInnerStereo2; // fInnerRadius * fInnerRadius
213 G4double fTanOuterStereo2; // fInnerRadius * fInnerRadius
214 G4double fEndZ2[2]; // fEndZ * fEndZ
215
216 G4VTwistSurface* fLowerEndcap; // Surface of -ve z
217 G4VTwistSurface* fUpperEndcap; // Surface of +ve z
218 G4VTwistSurface* fLatterTwisted; // Surface of -ve phi
219 G4VTwistSurface* fFormerTwisted; // Surface of +ve phi
220 G4VTwistSurface* fInnerHype; // Surface of -ve r
221 G4VTwistSurface* fOuterHype; // Surface of +ve r
222
223 G4double fCubicVolume = 0.0; // Cached value for cubic volume
224 G4double fSurfaceArea = 0.0; // Cached value for surface area
225
226 mutable G4bool fRebuildPolyhedron = false;
227 mutable G4Polyhedron* fpPolyhedron = nullptr; // polyhedron for vis
228
229 class LastState // last Inside result
230 {
231 public:
232 LastState()
233 {
234 p.set(kInfinity,kInfinity,kInfinity);
235 inside = kOutside;
236 }
237 ~LastState()= default;
238 LastState(const LastState& r) = default;
239 LastState& operator=(const LastState& r)
240 {
241 if (this == &r) { return *this; }
242 p = r.p; inside = r.inside;
243 return *this;
244 }
245 public:
247 EInside inside;
248 };
249
250 class LastVector // last SurfaceNormal result
251 {
252 public:
253 LastVector()
254 {
255 p.set(kInfinity,kInfinity,kInfinity);
256 vec.set(kInfinity,kInfinity,kInfinity);
257 surface = new G4VTwistSurface*[1];
258 }
259 ~LastVector()
260 {
261 delete [] surface;
262 }
263 LastVector(const LastVector& r) : p(r.p), vec(r.vec)
264 {
265 surface = new G4VTwistSurface*[1];
266 surface[0] = r.surface[0];
267 }
268 LastVector& operator=(const LastVector& r)
269 {
270 if (&r == this) { return *this; }
271 p = r.p; vec = r.vec;
272 delete [] surface; surface = new G4VTwistSurface*[1];
273 surface[0] = r.surface[0];
274 return *this;
275 }
276 public:
278 G4ThreeVector vec;
279 G4VTwistSurface **surface;
280 };
281
282 class LastValue // last G4double value
283 {
284 public:
285 LastValue()
286 {
287 p.set(kInfinity,kInfinity,kInfinity);
288 value = DBL_MAX;
289 }
290 ~LastValue()= default;
291 LastValue(const LastValue& r) = default;
292 LastValue& operator=(const LastValue& r)
293 {
294 if (this == &r) { return *this; }
295 p = r.p; value = r.value;
296 return *this;
297 }
298 public:
300 G4double value;
301 };
302
303 class LastValueWithDoubleVector // last G4double value
304 {
305 public:
306 LastValueWithDoubleVector()
307 {
308 p.set(kInfinity,kInfinity,kInfinity);
309 vec.set(kInfinity,kInfinity,kInfinity);
310 value = DBL_MAX;
311 }
312 ~LastValueWithDoubleVector()= default;
313 LastValueWithDoubleVector(const LastValueWithDoubleVector& r) = default;
314 LastValueWithDoubleVector& operator=(const LastValueWithDoubleVector& r)
315 {
316 if (this == &r) { return *this; }
317 p = r.p; vec = r.vec; value = r.value;
318 return *this;
319 }
320 public:
322 G4ThreeVector vec;
323 G4double value;
324 };
325
326 LastState fLastInside;
327 LastVector fLastNormal;
328 LastValue fLastDistanceToIn;
329 LastValue fLastDistanceToOut;
330 LastValueWithDoubleVector fLastDistanceToInWithV;
331 LastValueWithDoubleVector fLastDistanceToOutWithV;
332};
333
334//=====================================================================
335
336//---------------------
337// inline functions
338//---------------------
339
340inline
341void G4TwistedTubs::SetFields(G4double phitwist, G4double innerrad,
342 G4double outerrad, G4double negativeEndz,
343 G4double positiveEndz)
344{
345 fCubicVolume = 0.;
346 fPhiTwist = phitwist;
347 fEndZ[0] = negativeEndz;
348 fEndZ[1] = positiveEndz;
349 fEndZ2[0] = fEndZ[0] * fEndZ[0];
350 fEndZ2[1] = fEndZ[1] * fEndZ[1];
351 fInnerRadius = innerrad;
352 fOuterRadius = outerrad;
353 fInnerRadius2 = fInnerRadius * fInnerRadius;
354 fOuterRadius2 = fOuterRadius * fOuterRadius;
355
356 if (std::fabs(fEndZ[0]) >= std::fabs(fEndZ[1]))
357 {
358 fZHalfLength = std::fabs(fEndZ[0]);
359 }
360 else
361 {
362 fZHalfLength = std::fabs(fEndZ[1]);
363 }
364
365 G4double parity = (fPhiTwist > 0 ? 1 : -1);
366 G4double tanHalfTwist = std::tan(0.5 * fPhiTwist);
367 G4double innerNumerator = std::fabs(fInnerRadius * tanHalfTwist) * parity;
368 G4double outerNumerator = std::fabs(fOuterRadius * tanHalfTwist) * parity;
369
370 fTanInnerStereo = innerNumerator / fZHalfLength;
371 fTanOuterStereo = outerNumerator / fZHalfLength;
372 fTanInnerStereo2 = fTanInnerStereo * fTanInnerStereo;
373 fTanOuterStereo2 = fTanOuterStereo * fTanOuterStereo;
374 fInnerStereo = std::atan2(innerNumerator, fZHalfLength);
375 fOuterStereo = std::atan2(outerNumerator, fZHalfLength);
376 fEndInnerRadius[0] = std::sqrt(fInnerRadius2 + fEndZ2[0] * fTanInnerStereo2);
377 fEndInnerRadius[1] = std::sqrt(fInnerRadius2 + fEndZ2[1] * fTanInnerStereo2);
378 fEndOuterRadius[0] = std::sqrt(fOuterRadius2 + fEndZ2[0] * fTanOuterStereo2);
379 fEndOuterRadius[1] = std::sqrt(fOuterRadius2 + fEndZ2[1] * fTanOuterStereo2);
380
381 fKappa = tanHalfTwist / fZHalfLength;
382 fEndPhi[0] = std::atan2(fEndZ[0] * tanHalfTwist, fZHalfLength);
383 fEndPhi[1] = std::atan2(fEndZ[1] * tanHalfTwist, fZHalfLength);
384
385#ifdef G4TWISTDEBUG
386 G4cout << "/********* G4TwistedTubs::SetFields() Field Parameters ***************** " << G4endl;
387 G4cout << "/* fPhiTwist : " << fPhiTwist << G4endl;
388 G4cout << "/* fEndZ(0, 1) : " << fEndZ[0] << " , " << fEndZ[1] << G4endl;
389 G4cout << "/* fEndPhi(0, 1) : " << fEndPhi[0] << " , " << fEndPhi[1] << G4endl;
390 G4cout << "/* fInnerRadius, fOuterRadius : " << fInnerRadius << " , " << fOuterRadius << G4endl;
391 G4cout << "/* fEndInnerRadius(0, 1) : " << fEndInnerRadius[0] << " , "
392 << fEndInnerRadius[1] << G4endl;
393 G4cout << "/* fEndOuterRadius(0, 1) : " << fEndOuterRadius[0] << " , "
394 << fEndOuterRadius[1] << G4endl;
395 G4cout << "/* fInnerStereo, fOuterStereo : " << fInnerStereo << " , " << fOuterStereo << G4endl;
396 G4cout << "/* tanHalfTwist, fKappa : " << tanHalfTwist << " , " << fKappa << G4endl;
397 G4cout << "/*********************************************************************** " << G4endl;
398#endif
399}
400
401#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void set(double x, double y, double z)
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=false, G4bool *validnorm=nullptr, G4ThreeVector *n=nullptr) const override
G4double GetOuterRadius() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
G4double GetZHalfLength() const
G4double GetPhiTwist() const
G4Polyhedron * CreatePolyhedron() const override
G4double GetInnerStereo() const
G4double GetTanInnerStereo() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
G4TwistedTubs & operator=(const G4TwistedTubs &rhs)
G4double GetEndInnerRadius() const
G4double GetEndOuterRadius() const
G4double GetOuterStereo() const
G4VSolid * Clone() const override
G4Polyhedron * GetPolyhedron() const override
~G4TwistedTubs() override
G4GeometryType GetEntityType() const override
G4double GetEndPhi(G4int i) const
G4TwistedTubs(const G4String &pname, G4double twistedangle, G4double endinnerrad, G4double endouterrad, G4double halfzlen, G4double dphi)
G4double GetTanOuterStereo() const
std::ostream & StreamInfo(std::ostream &os) const override
G4double GetSurfaceArea() override
G4double GetCubicVolume() override
void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *) override
G4double GetTanInnerStereo2() const
G4double GetTanOuterStereo2() const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4VisExtent GetExtent() const override
G4ThreeVector GetPointOnSurface() const override
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
G4double GetEndZ(G4int i) const
G4double GetKappa() const
G4double GetEndInnerRadius(G4int i) const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
G4double GetInnerRadius() const
EInside Inside(const G4ThreeVector &p) const override
G4double GetDPhi() const
G4double GetEndOuterRadius(G4int i) const
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67
@ kOutside
Definition geomdefs.hh:68
#define DBL_MAX
Definition templates.hh:62