ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4TwistedTubs.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4TwistedTubs.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 // G4TwistedTubs
27 //
28 // Class description:
29 //
30 // G4TwistedTubs is a sort of twisted 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 // It can have inner & outer surfaces as well as G4TwistedTubs,
36 // but cannot has different stereo angles between the inner surface
37 // and outer surface.
38 
39 // 01-Aug-2002 - Kotoyo Hoshina (hoshina@hepburn.s.chiba-u.ac.jp), created.
40 // 13-Nov-2003 - O.Link (Oliver.Link@cern.ch), Integration in Geant4
41 // from original version in Jupiter-2.5.02 application.
42 // --------------------------------------------------------------------
43 #ifndef G4TWISTEDTUBS_HH
44 #define G4TWISTEDTUBS_HH
45 
46 #include "G4VSolid.hh"
47 #include "G4TwistTubsFlatSide.hh"
48 #include "G4TwistTubsSide.hh"
49 #include "G4TwistTubsHypeSide.hh"
50 
51 class G4SolidExtentList;
52 class G4ClippablePolygon;
53 
54 class G4TwistedTubs : public G4VSolid
55 {
56  public: // with description
57 
58  G4TwistedTubs(const G4String& pname, // Name of instance
59  G4double twistedangle, // Twisted angle
60  G4double endinnerrad, // Inner radius at endcap
61  G4double endouterrad, // Outer radius at endcap
62  G4double halfzlen, // half z length
63  G4double dphi); // Phi angle of a segment
64 
65  G4TwistedTubs(const G4String& pname, // Name of instance
66  G4double twistedangle, // Stereo angle
67  G4double endinnerrad, // Inner radius at endcap
68  G4double endouterrad, // Outer radius at endcap
69  G4double halfzlen, // half z length
70  G4int nseg, // Number of segments in totalPhi
71  G4double totphi); // Total angle of all segments
72 
73  G4TwistedTubs(const G4String& pname, // Name of instance
74  G4double twistedangle, // Twisted angle
75  G4double innerrad, // Inner radius at z=0
76  G4double outerrad, // Outer radius at z=0
77  G4double negativeEndz, // -ve z endplate
78  G4double positiveEndz, // +ve z endplate
79  G4double dphi); // Phi angle of a segment
80 
81  G4TwistedTubs(const G4String& pname, // Name of instance
82  G4double twistedangle, // Stereo angle
83  G4double innerrad, // Inner radius at z=0
84  G4double outerrad, // Outer radius at z=0
85  G4double negativeEndz, // -ve z endplate
86  G4double positiveEndz, // +ve z endplate
87  G4int nseg, // Number of segments in totalPhi
88  G4double totphi); // Total angle of all segments
89 
90  virtual ~G4TwistedTubs();
91 
93  const G4int /* n */ ,
94  const G4VPhysicalVolume* /* prep */ );
95 
97 
98  G4bool CalculateExtent(const EAxis pAxis,
99  const G4VoxelLimits& pVoxelLimit,
100  const G4AffineTransform& pTransform,
101  G4double& pMin,
102  G4double& pMax ) const;
103 
105  const G4ThreeVector& v ) const;
106 
107  G4double DistanceToIn (const G4ThreeVector& p ) const;
108 
110  const G4ThreeVector& v,
111  const G4bool calcnorm = false,
112  G4bool* validnorm = nullptr,
113  G4ThreeVector* n = nullptr ) const;
114 
115  G4double DistanceToOut(const G4ThreeVector& p) const;
116 
117  EInside Inside (const G4ThreeVector& p) const;
118 
119  G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const;
120 
121  void DescribeYourselfTo (G4VGraphicsScene& scene) const;
122  G4Polyhedron* CreatePolyhedron () const;
123  G4Polyhedron* GetPolyhedron () const;
124 
125  std::ostream &StreamInfo(std::ostream& os) const;
126 
127  // accessors
128 
129  inline G4double GetDPhi () const { return fDPhi ; }
130  inline G4double GetPhiTwist () const { return fPhiTwist ; }
131  inline G4double GetInnerRadius () const { return fInnerRadius; }
132  inline G4double GetOuterRadius () const { return fOuterRadius; }
133  inline G4double GetInnerStereo () const { return fInnerStereo; }
134  inline G4double GetOuterStereo () const { return fOuterStereo; }
135  inline G4double GetZHalfLength () const { return fZHalfLength; }
136  inline G4double GetKappa () const { return fKappa ; }
137 
138  inline G4double GetTanInnerStereo () const { return fTanInnerStereo ; }
139  inline G4double GetTanInnerStereo2() const { return fTanInnerStereo2 ; }
140  inline G4double GetTanOuterStereo () const { return fTanOuterStereo ; }
141  inline G4double GetTanOuterStereo2() const { return fTanOuterStereo2 ; }
142 
143  inline G4double GetEndZ (G4int i) const { return fEndZ[i] ; }
144  inline G4double GetEndPhi (G4int i) const { return fEndPhi[i]; }
146  { return fEndInnerRadius[i]; }
148  { return fEndOuterRadius[i]; }
149  inline G4double GetEndInnerRadius () const
150  { return (fEndInnerRadius[0] > fEndInnerRadius[1] ?
151  fEndInnerRadius[0] : fEndInnerRadius[1]); }
152  inline G4double GetEndOuterRadius () const
153  { return (fEndOuterRadius[0] > fEndOuterRadius[1] ?
154  fEndOuterRadius[0] : fEndOuterRadius[1]); }
155 
156  G4VisExtent GetExtent () const;
158  G4VSolid* Clone() const;
159 
161  // Returns an estimation of the geometrical cubic volume of the
162  // solid. Caches the computed value once computed the first time.
164  // Returns an estimation of the geometrical surface area of the
165  // solid. Caches the computed value once computed the first time.
166 
168 
169  public: // without description
170 
171  G4TwistedTubs(__void__&);
172  // Fake default constructor for usage restricted to direct object
173  // persistency for clients requiring preallocation of memory for
174  // persistifiable objects.
175 
176  G4TwistedTubs(const G4TwistedTubs& rhs);
177  G4TwistedTubs& operator=(const G4TwistedTubs& rhs);
178  // Copy constructor and assignment operator.
179 
180 #ifdef G4TWISTDEBUG
181  G4VTwistSurface* GetOuterHype() const { return fOuterHype; }
182 #endif
183 
184  private:
185 
186  inline void SetFields(G4double phitwist, G4double innerrad,
187  G4double outerrad,
188  G4double negativeEndz, G4double positiveEndz);
189 
190  void CreateSurfaces();
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:
233  {
235  inside = kOutside;
236  }
238  LastState(const LastState& r) : p(r.p), inside(r.inside){}
240  {
241  if (this == &r) { return *this; }
242  p = r.p; inside = r.inside;
243  return *this;
244  }
245  public:
248  };
249 
250  class LastVector // last SurfaceNormal result
251  {
252  public:
254  {
257  surface = new G4VTwistSurface*[1];
258  }
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  }
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:
280  };
281 
282  class LastValue // last G4double value
283  {
284  public:
286  {
288  value = DBL_MAX;
289  }
291  LastValue(const LastValue& r) : p(r.p), value(r.value){}
293  {
294  if (this == &r) { return *this; }
295  p = r.p; value = r.value;
296  return *this;
297  }
298  public:
301  };
302 
303  class LastValueWithDoubleVector // last G4double value
304  {
305  public:
307  {
310  value = DBL_MAX;
311  }
314  : p(r.p), vec(r.vec), value(r.value){}
316  {
317  if (this == &r) { return *this; }
318  p = r.p; vec = r.vec; value = r.value;
319  return *this;
320  }
321  public:
325  };
326 
333 
334  };
335 
336 //=====================================================================
337 
338 //---------------------
339 // inline functions
340 //---------------------
341 
342 inline
344  G4double outerrad, G4double negativeEndz,
345  G4double positiveEndz)
346 {
347  fCubicVolume = 0.;
348  fPhiTwist = phitwist;
349  fEndZ[0] = negativeEndz;
350  fEndZ[1] = positiveEndz;
351  fEndZ2[0] = fEndZ[0] * fEndZ[0];
352  fEndZ2[1] = fEndZ[1] * fEndZ[1];
357 
358  if (std::fabs(fEndZ[0]) >= std::fabs(fEndZ[1]))
359  {
360  fZHalfLength = std::fabs(fEndZ[0]);
361  }
362  else
363  {
364  fZHalfLength = std::fabs(fEndZ[1]);
365  }
366 
367  G4double parity = (fPhiTwist > 0 ? 1 : -1);
368  G4double tanHalfTwist = std::tan(0.5 * fPhiTwist);
369  G4double innerNumerator = std::fabs(fInnerRadius * tanHalfTwist) * parity;
370  G4double outerNumerator = std::fabs(fOuterRadius * tanHalfTwist) * parity;
371 
372  fTanInnerStereo = innerNumerator / fZHalfLength;
373  fTanOuterStereo = outerNumerator / fZHalfLength;
376  fInnerStereo = std::atan2(innerNumerator, fZHalfLength);
377  fOuterStereo = std::atan2(outerNumerator, fZHalfLength);
378  fEndInnerRadius[0] = std::sqrt(fInnerRadius2 + fEndZ2[0] * fTanInnerStereo2);
379  fEndInnerRadius[1] = std::sqrt(fInnerRadius2 + fEndZ2[1] * fTanInnerStereo2);
380  fEndOuterRadius[0] = std::sqrt(fOuterRadius2 + fEndZ2[0] * fTanOuterStereo2);
381  fEndOuterRadius[1] = std::sqrt(fOuterRadius2 + fEndZ2[1] * fTanOuterStereo2);
382 
383  fKappa = tanHalfTwist / fZHalfLength;
384  fEndPhi[0] = std::atan2(fEndZ[0] * tanHalfTwist, fZHalfLength);
385  fEndPhi[1] = std::atan2(fEndZ[1] * tanHalfTwist, fZHalfLength);
386 
387 #ifdef G4TWISTDEBUG
388  G4cout << "/********* G4TwistedTubs::SetFields() Field Parameters ***************** " << G4endl;
389  G4cout << "/* fPhiTwist : " << fPhiTwist << G4endl;
390  G4cout << "/* fEndZ(0, 1) : " << fEndZ[0] << " , " << fEndZ[1] << G4endl;
391  G4cout << "/* fEndPhi(0, 1) : " << fEndPhi[0] << " , " << fEndPhi[1] << G4endl;
392  G4cout << "/* fInnerRadius, fOuterRadius : " << fInnerRadius << " , " << fOuterRadius << G4endl;
393  G4cout << "/* fEndInnerRadius(0, 1) : " << fEndInnerRadius[0] << " , "
394  << fEndInnerRadius[1] << G4endl;
395  G4cout << "/* fEndOuterRadius(0, 1) : " << fEndOuterRadius[0] << " , "
396  << fEndOuterRadius[1] << G4endl;
397  G4cout << "/* fInnerStereo, fOuterStereo : " << fInnerStereo << " , " << fOuterStereo << G4endl;
398  G4cout << "/* tanHalfTwist, fKappa : " << tanHalfTwist << " , " << fKappa << G4endl;
399  G4cout << "/*********************************************************************** " << G4endl;
400 #endif
401 }
402 
403 #endif