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