Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
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 //
27 // $Id: G4TwistedTubs.hh 104316 2017-05-24 13:04:23Z 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  void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
108 
109  G4bool CalculateExtent(const EAxis pAxis,
110  const G4VoxelLimits &pVoxelLimit,
111  const G4AffineTransform &pTransform,
112  G4double &pMin,
113  G4double &pMax ) const;
114 
116  const G4ThreeVector &v ) const;
117 
118  G4double DistanceToIn (const G4ThreeVector &p ) const;
119 
121  const G4ThreeVector &v,
122  const G4bool calcnorm=G4bool(false),
123  G4bool *validnorm=0,
124  G4ThreeVector *n=0 ) const;
125 
126  G4double DistanceToOut(const G4ThreeVector &p) const;
127 
128  EInside Inside (const G4ThreeVector &p) const;
129 
130  G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const;
131 
132  void DescribeYourselfTo (G4VGraphicsScene &scene) const;
133  G4Polyhedron *CreatePolyhedron () const;
134  G4Polyhedron *GetPolyhedron () const;
135 
136  std::ostream &StreamInfo(std::ostream& os) const;
137 
138  // accessors
139 
140  inline G4double GetDPhi () const { return fDPhi ; }
141  inline G4double GetPhiTwist () const { return fPhiTwist ; }
142  inline G4double GetInnerRadius () const { return fInnerRadius; }
143  inline G4double GetOuterRadius () const { return fOuterRadius; }
144  inline G4double GetInnerStereo () const { return fInnerStereo; }
145  inline G4double GetOuterStereo () const { return fOuterStereo; }
146  inline G4double GetZHalfLength () const { return fZHalfLength; }
147  inline G4double GetKappa () const { return fKappa ; }
148 
149  inline G4double GetTanInnerStereo () const { return fTanInnerStereo ; }
150  inline G4double GetTanInnerStereo2() const { return fTanInnerStereo2 ; }
151  inline G4double GetTanOuterStereo () const { return fTanOuterStereo ; }
152  inline G4double GetTanOuterStereo2() const { return fTanOuterStereo2 ; }
153 
154  inline G4double GetEndZ (G4int i) const { return fEndZ[i] ; }
155  inline G4double GetEndPhi (G4int i) const { return fEndPhi[i]; }
157  { return fEndInnerRadius[i]; }
159  { return fEndOuterRadius[i]; }
160  inline G4double GetEndInnerRadius () const
161  { return (fEndInnerRadius[0] > fEndInnerRadius[1] ?
162  fEndInnerRadius[0] : fEndInnerRadius[1]); }
163  inline G4double GetEndOuterRadius () const
164  { return (fEndOuterRadius[0] > fEndOuterRadius[1] ?
165  fEndOuterRadius[0] : fEndOuterRadius[1]); }
166 
167  G4VisExtent GetExtent () const;
169  G4VSolid* Clone() const;
170 
172  // Returns an estimation of the geometrical cubic volume of the
173  // solid. Caches the computed value once computed the first time.
175  // Returns an estimation of the geometrical surface area of the
176  // solid. Caches the computed value once computed the first time.
177 
179 
180  public: // without description
181 
182  G4TwistedTubs(__void__&);
183  // Fake default constructor for usage restricted to direct object
184  // persistency for clients requiring preallocation of memory for
185  // persistifiable objects.
186 
187  G4TwistedTubs(const G4TwistedTubs& rhs);
188  G4TwistedTubs& operator=(const G4TwistedTubs& rhs);
189  // Copy constructor and assignment operator.
190 
191 #ifdef G4TWISTDEBUG
192  G4VTwistSurface * GetOuterHype() const { return fOuterHype; }
193 #endif
194 
195  private:
196 
197  inline void SetFields(G4double phitwist, G4double innerrad,
198  G4double outerrad,
199  G4double negativeEndz, G4double positiveEndz);
200 
201  void CreateSurfaces();
202 
203  private:
204 
205  G4double fPhiTwist; // Twist angle from -fZHalfLength to fZHalfLength
206  G4double fInnerRadius; // Inner-hype radius at z=0
207  G4double fOuterRadius; // Outer-hype radius at z=0
208  G4double fEndZ[2]; // z at endcaps, [0] = -ve z, [1] = +ve z
209  G4double fDPhi; // Phi-width of a segment fDPhi > 0
210  G4double fZHalfLength; // Half length along z-axis
211 
212  G4double fInnerStereo; // Inner-hype stereo angle
213  G4double fOuterStereo; // Outer-hype stereo angle
214  G4double fTanInnerStereo; // std::tan(innerStereoAngle)
215  G4double fTanOuterStereo; // std::tan(outerStereoAngle)
216  G4double fKappa; // std::tan(fPhiTwist/2)/fZHalfLen;
217  G4double fEndInnerRadius[2]; // Inner-hype radii endcaps [0] -ve z, [1] +ve z
218  G4double fEndOuterRadius[2]; // Outer-hype radii endcaps [0] -ve z, [1] +ve z
219  G4double fEndPhi[2]; // Phi endcaps, [0] = -ve z, [1] = +ve z
220 
221  G4double fInnerRadius2; // fInnerRadius * fInnerRadius
222  G4double fOuterRadius2; // fOuterRadius * fOuterRadius
223  G4double fTanInnerStereo2; // fInnerRadius * fInnerRadius
224  G4double fTanOuterStereo2; // fInnerRadius * fInnerRadius
225  G4double fEndZ2[2]; // fEndZ * fEndZ
226 
227  G4VTwistSurface *fLowerEndcap; // Surface of -ve z
228  G4VTwistSurface *fUpperEndcap; // Surface of +ve z
229  G4VTwistSurface *fLatterTwisted; // Surface of -ve phi
230  G4VTwistSurface *fFormerTwisted; // Surface of +ve phi
231  G4VTwistSurface *fInnerHype; // Surface of -ve r
232  G4VTwistSurface *fOuterHype; // Surface of +ve r
233 
234  G4double fCubicVolume; // Cached value for cubic volume
235  G4double fSurfaceArea; // Cached value for surface area
236 
238  mutable G4Polyhedron* fpPolyhedron; // pointer to polyhedron for vis
239 
240  class LastState // last Inside result
241  {
242  public:
244  {
246  inside = kOutside;
247  }
249  LastState(const LastState& r) : p(r.p), inside(r.inside){}
251  {
252  if (this == &r) { return *this; }
253  p = r.p; inside = r.inside;
254  return *this;
255  }
256  public:
259  };
260 
261  class LastVector // last SurfaceNormal result
262  {
263  public:
265  {
268  surface = new G4VTwistSurface*[1];
269  }
271  {
272  delete [] surface;
273  }
274  LastVector(const LastVector& r) : p(r.p), vec(r.vec)
275  {
276  surface = new G4VTwistSurface*[1];
277  surface[0] = r.surface[0];
278  }
280  {
281  if (&r == this) { return *this; }
282  p = r.p; vec = r.vec;
283  delete [] surface; surface = new G4VTwistSurface*[1];
284  surface[0] = r.surface[0];
285  return *this;
286  }
287  public:
291  };
292 
293  class LastValue // last G4double value
294  {
295  public:
297  {
299  value = DBL_MAX;
300  }
302  LastValue(const LastValue& r) : p(r.p), value(r.value){}
304  {
305  if (this == &r) { return *this; }
306  p = r.p; value = r.value;
307  return *this;
308  }
309  public:
312  };
313 
314  class LastValueWithDoubleVector // last G4double value
315  {
316  public:
318  {
321  value = DBL_MAX;
322  }
325  : p(r.p), vec(r.vec), value(r.value){}
327  {
328  if (this == &r) { return *this; }
329  p = r.p; vec = r.vec; value = r.value;
330  return *this;
331  }
332  public:
336  };
337 
344 
345  };
346 
347 //=====================================================================
348 
349 //---------------------
350 // inline functions
351 //---------------------
352 
353 inline
354 void G4TwistedTubs::SetFields(G4double phitwist, G4double innerrad,
355  G4double outerrad, G4double negativeEndz,
356  G4double positiveEndz)
357 {
358  fCubicVolume = 0.;
359  fPhiTwist = phitwist;
360  fEndZ[0] = negativeEndz;
361  fEndZ[1] = positiveEndz;
362  fEndZ2[0] = fEndZ[0] * fEndZ[0];
363  fEndZ2[1] = fEndZ[1] * fEndZ[1];
364  fInnerRadius = innerrad;
365  fOuterRadius = outerrad;
368 
369  if (std::fabs(fEndZ[0]) >= std::fabs(fEndZ[1])) {
370  fZHalfLength = std::fabs(fEndZ[0]);
371  } else {
372  fZHalfLength = std::fabs(fEndZ[1]);
373  }
374 
375  G4double parity = (fPhiTwist > 0 ? 1 : -1);
376  G4double tanHalfTwist = std::tan(0.5 * fPhiTwist);
377  G4double innerNumerator = std::fabs(fInnerRadius * tanHalfTwist) * parity;
378  G4double outerNumerator = std::fabs(fOuterRadius * tanHalfTwist) * parity;
379 
380  fTanInnerStereo = innerNumerator / fZHalfLength;
381  fTanOuterStereo = outerNumerator / fZHalfLength;
384  fInnerStereo = std::atan2(innerNumerator, fZHalfLength);
385  fOuterStereo = std::atan2(outerNumerator, fZHalfLength);
386  fEndInnerRadius[0] = std::sqrt(fInnerRadius2 + fEndZ2[0] * fTanInnerStereo2);
387  fEndInnerRadius[1] = std::sqrt(fInnerRadius2 + fEndZ2[1] * fTanInnerStereo2);
388  fEndOuterRadius[0] = std::sqrt(fOuterRadius2 + fEndZ2[0] * fTanOuterStereo2);
389  fEndOuterRadius[1] = std::sqrt(fOuterRadius2 + fEndZ2[1] * fTanOuterStereo2);
390 
391  fKappa = tanHalfTwist / fZHalfLength;
392  fEndPhi[0] = std::atan2(fEndZ[0] * tanHalfTwist, fZHalfLength);
393  fEndPhi[1] = std::atan2(fEndZ[1] * tanHalfTwist, fZHalfLength);
394 
395 #ifdef G4TWISTDEBUG
396  G4cout << "/********* G4TwistedTubs::SetFields() Field Parameters ***************** " << G4endl;
397  G4cout << "/* fPhiTwist : " << fPhiTwist << G4endl;
398  G4cout << "/* fEndZ(0, 1) : " << fEndZ[0] << " , " << fEndZ[1] << G4endl;
399  G4cout << "/* fEndPhi(0, 1) : " << fEndPhi[0] << " , " << fEndPhi[1] << G4endl;
400  G4cout << "/* fInnerRadius, fOuterRadius : " << fInnerRadius << " , " << fOuterRadius << G4endl;
401  G4cout << "/* fEndInnerRadius(0, 1) : " << fEndInnerRadius[0] << " , "
402  << fEndInnerRadius[1] << G4endl;
403  G4cout << "/* fEndOuterRadius(0, 1) : " << fEndOuterRadius[0] << " , "
404  << fEndOuterRadius[1] << G4endl;
405  G4cout << "/* fInnerStereo, fOuterStereo : " << fInnerStereo << " , " << fOuterStereo << G4endl;
406  G4cout << "/* tanHalfTwist, fKappa : " << tanHalfTwist << " , " << fKappa << G4endl;
407  G4cout << "/*********************************************************************** " << G4endl;
408 #endif
409 }
410 
411 #endif
G4double fTanOuterStereo
void set(double x, double y, double z)
G4double GetZHalfLength() const
G4double fEndPhi[2]
void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *)
G4Polyhedron * fpPolyhedron
G4double GetEndOuterRadius() const
LastValue fLastDistanceToOut
G4double fPhiTwist
LastVector & operator=(const LastVector &r)
G4Polyhedron * GetPolyhedron() const
G4double fOuterRadius
G4double GetTanOuterStereo2() const
G4double GetSurfaceArea()
static const G4double kInfinity
Definition: geomdefs.hh:42
G4VTwistSurface * fInnerHype
void DescribeYourselfTo(G4VGraphicsScene &scene) const
G4double fInnerStereo
G4double fTanOuterStereo2
LastState(const LastState &r)
G4double GetDPhi() const
G4double GetKappa() const
#define G4endl
Definition: G4ios.hh:61
const char * p
Definition: xmltok.h:285
G4double GetOuterRadius() const
G4bool fRebuildPolyhedron
G4TwistedTubs(const G4String &pname, G4double twistedangle, G4double endinnerrad, G4double endouterrad, G4double halfzlen, G4double dphi)
G4VTwistSurface * fLatterTwisted
G4double fTanInnerStereo2
G4double fTanInnerStereo
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=G4bool(false), G4bool *validnorm=0, G4ThreeVector *n=0) const
G4VTwistSurface * fLowerEndcap
G4GeometryType GetEntityType() const
G4double GetInnerStereo() const
LastValue(const LastValue &r)
G4double fCubicVolume
G4VisExtent GetExtent() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4VTwistSurface * fFormerTwisted
LastState fLastInside
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double fOuterStereo
LastValueWithDoubleVector & operator=(const LastValueWithDoubleVector &r)
G4double GetEndZ(G4int i) const
const XML_Char int const XML_Char * value
Definition: expat.h:331
G4VSolid * Clone() const
LastVector fLastNormal
G4double GetCubicVolume()
EInside Inside(const G4ThreeVector &p) const
G4VTwistSurface * fUpperEndcap
G4TwistedTubs & operator=(const G4TwistedTubs &rhs)
LastValue fLastDistanceToIn
LastValue & operator=(const LastValue &r)
G4double fEndOuterRadius[2]
G4Polyhedron * CreatePolyhedron() const
G4double fInnerRadius2
G4double GetPhiTwist() const
G4double GetEndInnerRadius() const
G4VTwistSurface ** surface
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
std::ostream & StreamInfo(std::ostream &os) const
int G4int
Definition: G4Types.hh:78
EInside
Definition: geomdefs.hh:58
G4double GetEndOuterRadius(G4int i) const
LastValueWithDoubleVector fLastDistanceToInWithV
G4VTwistSurface * fOuterHype
EAxis
Definition: geomdefs.hh:54
G4ThreeVector GetPointOnSurface() const
LastState & operator=(const LastState &r)
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
LastVector(const LastVector &r)
virtual ~G4TwistedTubs()
void SetFields(G4double phitwist, G4double innerrad, G4double outerrad, G4double negativeEndz, G4double positiveEndz)
LastValueWithDoubleVector(const LastValueWithDoubleVector &r)
G4double fEndInnerRadius[2]
G4GLOB_DLL std::ostream G4cout
G4double GetEndInnerRadius(G4int i) const
Char_t n[5]
G4double fInnerRadius
G4double fZHalfLength
static int parity(int x)
G4double fSurfaceArea
G4double GetEndPhi(G4int i) const
#define DBL_MAX
Definition: templates.hh:83
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
G4double GetTanInnerStereo2() const
G4double fEndZ[2]
G4double GetInnerRadius() const
G4double fEndZ2[2]
void CreateSurfaces()
G4double GetTanOuterStereo() const
LastValueWithDoubleVector fLastDistanceToOutWithV
G4double GetOuterStereo() const
G4double fOuterRadius2
G4double GetTanInnerStereo() const