Geant4-11
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 (hoshina@hepburn.s.chiba-u.ac.jp), created.
41// 13-Nov-2003 - O.Link (Oliver.Link@cern.ch), 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: // with description
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 virtual ~G4TwistedTubs();
92
94 const G4int /* n */ ,
95 const G4VPhysicalVolume* /* prep */ );
96
98
99 G4bool CalculateExtent(const EAxis pAxis,
100 const G4VoxelLimits& pVoxelLimit,
101 const G4AffineTransform& pTransform,
102 G4double& pMin,
103 G4double& pMax ) const;
104
106 const G4ThreeVector& v ) const;
107
108 G4double DistanceToIn (const G4ThreeVector& p ) const;
109
111 const G4ThreeVector& v,
112 const G4bool calcnorm = false,
113 G4bool* validnorm = nullptr,
114 G4ThreeVector* n = nullptr ) const;
115
116 G4double DistanceToOut(const G4ThreeVector& p) const;
117
118 EInside Inside (const G4ThreeVector& p) const;
119
121
122 void DescribeYourselfTo (G4VGraphicsScene& scene) const;
124 G4Polyhedron* GetPolyhedron () const;
125
126 std::ostream &StreamInfo(std::ostream& os) const;
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] ?
154 { return (fEndOuterRadius[0] > fEndOuterRadius[1] ?
156
157 G4VisExtent GetExtent () const;
159 G4VSolid* Clone() const;
160
162 // Returns an estimation of the geometrical cubic volume of the
163 // solid. Caches the computed value once computed the first time.
165 // Returns an estimation of the geometrical surface area of the
166 // solid. Caches the computed value once computed the first time.
167
169
170 public: // without description
171
172 G4TwistedTubs(__void__&);
173 // Fake default constructor for usage restricted to direct object
174 // persistency for clients requiring preallocation of memory for
175 // persistifiable objects.
176
177 G4TwistedTubs(const G4TwistedTubs& rhs);
179 // Copy constructor and assignment operator.
180
181#ifdef G4TWISTDEBUG
182 G4VTwistSurface* GetOuterHype() const { return fOuterHype; }
183#endif
184
185 private:
186
187 inline void SetFields(G4double phitwist, G4double innerrad,
188 G4double outerrad,
189 G4double negativeEndz, G4double positiveEndz);
190
191 void CreateSurfaces();
192
193 private:
194
195 G4double fPhiTwist; // Twist angle from -fZHalfLength to fZHalfLength
196 G4double fInnerRadius; // Inner-hype radius at z=0
197 G4double fOuterRadius; // Outer-hype radius at z=0
198 G4double fEndZ[2]; // z at endcaps, [0] = -ve z, [1] = +ve z
199 G4double fDPhi; // Phi-width of a segment fDPhi > 0
200 G4double fZHalfLength; // Half length along z-axis
201
202 G4double fInnerStereo; // Inner-hype stereo angle
203 G4double fOuterStereo; // Outer-hype stereo angle
204 G4double fTanInnerStereo; // std::tan(innerStereoAngle)
205 G4double fTanOuterStereo; // std::tan(outerStereoAngle)
206 G4double fKappa; // std::tan(fPhiTwist/2)/fZHalfLen;
207 G4double fEndInnerRadius[2]; // Inner-hype radii endcaps [0] -ve z, [1] +ve z
208 G4double fEndOuterRadius[2]; // Outer-hype radii endcaps [0] -ve z, [1] +ve z
209 G4double fEndPhi[2]; // Phi endcaps, [0] = -ve z, [1] = +ve z
210
211 G4double fInnerRadius2; // fInnerRadius * fInnerRadius
212 G4double fOuterRadius2; // fOuterRadius * fOuterRadius
213 G4double fTanInnerStereo2; // fInnerRadius * fInnerRadius
214 G4double fTanOuterStereo2; // fInnerRadius * fInnerRadius
215 G4double fEndZ2[2]; // fEndZ * fEndZ
216
217 G4VTwistSurface* fLowerEndcap; // Surface of -ve z
218 G4VTwistSurface* fUpperEndcap; // Surface of +ve z
219 G4VTwistSurface* fLatterTwisted; // Surface of -ve phi
220 G4VTwistSurface* fFormerTwisted; // Surface of +ve phi
221 G4VTwistSurface* fInnerHype; // Surface of -ve r
222 G4VTwistSurface* fOuterHype; // Surface of +ve r
223
224 G4double fCubicVolume = 0.0; // Cached value for cubic volume
225 G4double fSurfaceArea = 0.0; // Cached value for surface area
226
227 mutable G4bool fRebuildPolyhedron = false;
228 mutable G4Polyhedron* fpPolyhedron = nullptr; // polyhedron for vis
229
230 class LastState // last Inside result
231 {
232 public:
234 {
237 }
239 LastState(const LastState& r) : p(r.p), inside(r.inside){}
241 {
242 if (this == &r) { return *this; }
243 p = r.p; inside = r.inside;
244 return *this;
245 }
246 public:
249 };
250
251 class LastVector // last SurfaceNormal result
252 {
253 public:
255 {
258 surface = new G4VTwistSurface*[1];
259 }
261 {
262 delete [] surface;
263 }
264 LastVector(const LastVector& r) : p(r.p), vec(r.vec)
265 {
266 surface = new G4VTwistSurface*[1];
267 surface[0] = r.surface[0];
268 }
270 {
271 if (&r == this) { return *this; }
272 p = r.p; vec = r.vec;
273 delete [] surface; surface = new G4VTwistSurface*[1];
274 surface[0] = r.surface[0];
275 return *this;
276 }
277 public:
281 };
282
283 class LastValue // last G4double value
284 {
285 public:
287 {
289 value = DBL_MAX;
290 }
292 LastValue(const LastValue& r) : p(r.p), value(r.value){}
294 {
295 if (this == &r) { return *this; }
296 p = r.p; value = r.value;
297 return *this;
298 }
299 public:
302 };
303
304 class LastValueWithDoubleVector // last G4double value
305 {
306 public:
308 {
311 value = DBL_MAX;
312 }
315 : p(r.p), vec(r.vec), value(r.value){}
317 {
318 if (this == &r) { return *this; }
319 p = r.p; vec = r.vec; value = r.value;
320 return *this;
321 }
322 public:
326 };
327
334
335 };
336
337//=====================================================================
338
339//---------------------
340// inline functions
341//---------------------
342
343inline
345 G4double outerrad, G4double negativeEndz,
346 G4double positiveEndz)
347{
348 fCubicVolume = 0.;
349 fPhiTwist = phitwist;
350 fEndZ[0] = negativeEndz;
351 fEndZ[1] = positiveEndz;
352 fEndZ2[0] = fEndZ[0] * fEndZ[0];
353 fEndZ2[1] = fEndZ[1] * fEndZ[1];
354 fInnerRadius = innerrad;
355 fOuterRadius = outerrad;
358
359 if (std::fabs(fEndZ[0]) >= std::fabs(fEndZ[1]))
360 {
361 fZHalfLength = std::fabs(fEndZ[0]);
362 }
363 else
364 {
365 fZHalfLength = std::fabs(fEndZ[1]);
366 }
367
368 G4double parity = (fPhiTwist > 0 ? 1 : -1);
369 G4double tanHalfTwist = std::tan(0.5 * fPhiTwist);
370 G4double innerNumerator = std::fabs(fInnerRadius * tanHalfTwist) * parity;
371 G4double outerNumerator = std::fabs(fOuterRadius * tanHalfTwist) * parity;
372
373 fTanInnerStereo = innerNumerator / fZHalfLength;
374 fTanOuterStereo = outerNumerator / fZHalfLength;
377 fInnerStereo = std::atan2(innerNumerator, fZHalfLength);
378 fOuterStereo = std::atan2(outerNumerator, fZHalfLength);
383
384 fKappa = tanHalfTwist / fZHalfLength;
385 fEndPhi[0] = std::atan2(fEndZ[0] * tanHalfTwist, fZHalfLength);
386 fEndPhi[1] = std::atan2(fEndZ[1] * tanHalfTwist, fZHalfLength);
387
388#ifdef G4TWISTDEBUG
389 G4cout << "/********* G4TwistedTubs::SetFields() Field Parameters ***************** " << G4endl;
390 G4cout << "/* fPhiTwist : " << fPhiTwist << G4endl;
391 G4cout << "/* fEndZ(0, 1) : " << fEndZ[0] << " , " << fEndZ[1] << G4endl;
392 G4cout << "/* fEndPhi(0, 1) : " << fEndPhi[0] << " , " << fEndPhi[1] << G4endl;
393 G4cout << "/* fInnerRadius, fOuterRadius : " << fInnerRadius << " , " << fOuterRadius << G4endl;
394 G4cout << "/* fEndInnerRadius(0, 1) : " << fEndInnerRadius[0] << " , "
395 << fEndInnerRadius[1] << G4endl;
396 G4cout << "/* fEndOuterRadius(0, 1) : " << fEndOuterRadius[0] << " , "
397 << fEndOuterRadius[1] << G4endl;
398 G4cout << "/* fInnerStereo, fOuterStereo : " << fInnerStereo << " , " << fOuterStereo << G4endl;
399 G4cout << "/* tanHalfTwist, fKappa : " << tanHalfTwist << " , " << fKappa << G4endl;
400 G4cout << "/*********************************************************************** " << G4endl;
401#endif
402}
403
404#endif
static const G4double pMax
static const G4double pMin
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void set(double x, double y, double z)
LastState(const LastState &r)
LastState & operator=(const LastState &r)
LastValueWithDoubleVector & operator=(const LastValueWithDoubleVector &r)
LastValueWithDoubleVector(const LastValueWithDoubleVector &r)
LastValue & operator=(const LastValue &r)
LastValue(const LastValue &r)
G4VTwistSurface ** surface
LastVector(const LastVector &r)
LastVector & operator=(const LastVector &r)
G4double fPhiTwist
G4double fOuterRadius
G4double fEndZ2[2]
G4double GetOuterRadius() const
G4double fTanInnerStereo2
G4double fOuterStereo
LastValue fLastDistanceToOut
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4double fSurfaceArea
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
G4ThreeVector GetPointOnSurface() const
G4double fTanInnerStereo
G4double GetZHalfLength() const
LastValue fLastDistanceToIn
G4double fEndZ[2]
G4double GetPhiTwist() const
G4double GetSurfaceArea()
G4VTwistSurface * fLatterTwisted
G4double fTanOuterStereo
void DescribeYourselfTo(G4VGraphicsScene &scene) const
G4double fCubicVolume
G4Polyhedron * CreatePolyhedron() const
void SetFields(G4double phitwist, G4double innerrad, G4double outerrad, G4double negativeEndz, G4double positiveEndz)
G4double GetInnerStereo() const
G4double GetTanInnerStereo() const
void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *)
G4VTwistSurface * fFormerTwisted
G4VSolid * Clone() const
G4TwistedTubs & operator=(const G4TwistedTubs &rhs)
G4double GetEndInnerRadius() const
G4double GetEndOuterRadius() const
G4VTwistSurface * fUpperEndcap
G4double GetOuterStereo() const
void CreateSurfaces()
G4double fEndInnerRadius[2]
G4double fInnerStereo
G4double GetCubicVolume()
G4Polyhedron * GetPolyhedron() const
G4double fTanOuterStereo2
G4VTwistSurface * fLowerEndcap
G4double fInnerRadius2
G4GeometryType GetEntityType() const
virtual ~G4TwistedTubs()
LastValueWithDoubleVector fLastDistanceToOutWithV
G4double GetEndPhi(G4int i) const
G4TwistedTubs(const G4String &pname, G4double twistedangle, G4double endinnerrad, G4double endouterrad, G4double halfzlen, G4double dphi)
LastVector fLastNormal
LastValueWithDoubleVector fLastDistanceToInWithV
G4double GetTanOuterStereo() const
G4VTwistSurface * fInnerHype
G4double fInnerRadius
EInside Inside(const G4ThreeVector &p) const
G4double fZHalfLength
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=false, G4bool *validnorm=nullptr, G4ThreeVector *n=nullptr) const
G4VisExtent GetExtent() const
G4double GetTanInnerStereo2() const
G4double GetTanOuterStereo2() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4bool fRebuildPolyhedron
G4double GetEndZ(G4int i) const
G4Polyhedron * fpPolyhedron
G4double GetKappa() const
G4double GetEndInnerRadius(G4int i) const
G4double fOuterRadius2
LastState fLastInside
G4double GetInnerRadius() const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
G4double fEndPhi[2]
G4VTwistSurface * fOuterHype
G4double fEndOuterRadius[2]
std::ostream & StreamInfo(std::ostream &os) const
G4double GetDPhi() const
G4double GetEndOuterRadius(G4int i) const
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
@ kOutside
Definition: geomdefs.hh:68
static const G4double kInfinity
Definition: geomdefs.hh:41
string pname
Definition: eplot.py:33
static int parity(int x)
#define DBL_MAX
Definition: templates.hh:62