Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4TwistTubsFlatSide.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: G4TwistTubsFlatSide.hh 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 //
34 // G4TwistTubsFlatSide
35 //
36 // Class description:
37 //
38 // Class describing a flat boundary surface for a cylinder.
39 
40 // Author:
41 // 01-Aug-2002 - Kotoyo Hoshina (hoshina@hepburn.s.chiba-u.ac.jp)
42 //
43 // History:
44 // 13-Nov-2003 - O.Link (Oliver.Link@cern.ch), Integration in Geant4
45 // from original version in Jupiter-2.5.02 application.
46 // --------------------------------------------------------------------
47 #ifndef __G4TWISTTUBSFLATSIDE__
48 #define __G4TWISTTUBSFLATSIDE__
49 
50 #include "G4VTwistSurface.hh"
51 
53 {
54  public: // with description
55 
57  const G4RotationMatrix &rot,
58  const G4ThreeVector &tlate,
59  const G4ThreeVector &n,
60  const EAxis axis1 = kRho, // RHO axis !
61  const EAxis axis2 = kPhi, // PHI axis !
62  G4double axis0min = -kInfinity,
63  G4double axis1min = -kInfinity,
64  G4double axis0max = kInfinity,
65  G4double axis1max = kInfinity );
66 
67  G4TwistTubsFlatSide( const G4String &name,
68  G4double EndInnerRadius[2],
69  G4double EndOuterRadius[2],
70  G4double DPhi,
71  G4double EndPhi[2],
72  G4double EndZ[2],
73  G4int handedness ) ;
74 
75  virtual ~G4TwistTubsFlatSide();
76  virtual G4ThreeVector GetNormal(const G4ThreeVector & /* xx */ ,
77  G4bool isGlobal = false);
78  virtual G4int DistanceToSurface(const G4ThreeVector &gp,
79  const G4ThreeVector &gv,
80  G4ThreeVector gxx[],
81  G4double distance[],
82  G4int areacode[],
83  G4bool isvalid[],
84  EValidate validate = kValidateWithTol);
85 
86  virtual G4int DistanceToSurface(const G4ThreeVector &gp,
87  G4ThreeVector gxx[],
88  G4double distance[],
89  G4int areacode[]);
90 
92  G4bool isGlobal = false ) ;
93  virtual G4double GetBoundaryMin(G4double phi) ;
94  virtual G4double GetBoundaryMax(G4double phi) ;
95  virtual G4double GetSurfaceArea() { return fSurfaceArea ; }
96  virtual void GetFacets( G4int m, G4int n, G4double xyz[][3],
97  G4int faces[][4], G4int iside ) ;
98 
100 
101  public: // without description
102 
103  G4TwistTubsFlatSide(__void__&);
104  // Fake default constructor for usage restricted to direct object
105  // persistency for clients requiring preallocation of memory for
106  // persistifiable objects.
107 
108  protected: // with description
109 
110  virtual G4int GetAreaCode(const G4ThreeVector &xx,
111  G4bool withTol = true) ;
112 
113  private:
114 
115  virtual void SetCorners();
116  virtual void SetBoundaries();
117 
118 };
119 
121 SurfacePoint(G4double phi , G4double rho , G4bool isGlobal )
122 {
123  G4ThreeVector SurfPoint (rho*std::cos(phi) , rho*std::sin(phi) , 0);
124 
125  if (isGlobal) { return (fRot * SurfPoint + fTrans); }
126  return SurfPoint;
127 }
128 
129 inline
131 {
133  return std::atan2( dphimin.y(), dphimin.x() );
134 }
135 
136 inline
138 {
139  G4ThreeVector dphimax = GetCorner(sC0Max1Max);
140  return std::atan2( dphimax.y(), dphimax.x() );
141 }
142 
143 #endif
Definition: geomdefs.hh:54
const XML_Char * name
Definition: expat.h:151
virtual G4ThreeVector SurfacePoint(G4double, G4double, G4bool isGlobal=false)
Double_t xx
static const G4double kInfinity
Definition: geomdefs.hh:42
virtual G4int DistanceToSurface(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector gxx[], G4double distance[], G4int areacode[], G4bool isvalid[], EValidate validate=kValidateWithTol)
G4ThreeVector fTrans
G4RotationMatrix fRot
static constexpr double m
Definition: G4SIunits.hh:129
G4TwistTubsFlatSide(const G4String &name, const G4RotationMatrix &rot, const G4ThreeVector &tlate, const G4ThreeVector &n, const EAxis axis1=kRho, const EAxis axis2=kPhi, G4double axis0min=-kInfinity, G4double axis1min=-kInfinity, G4double axis0max=kInfinity, G4double axis1max=kInfinity)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static const G4int sC0Max1Max
virtual void GetFacets(G4int m, G4int n, G4double xyz[][3], G4int faces[][4], G4int iside)
virtual G4double GetBoundaryMax(G4double phi)
G4ThreeVector GetCorner(G4int areacode) const
int G4int
Definition: G4Types.hh:78
EAxis
Definition: geomdefs.hh:54
virtual void SetBoundaries()
virtual G4double GetBoundaryMin(G4double phi)
double x() const
virtual G4int GetAreaCode(const G4ThreeVector &xx, G4bool withTol=true)
Char_t n[5]
Definition: geomdefs.hh:54
double y() const
virtual G4ThreeVector GetNormal(const G4ThreeVector &, G4bool isGlobal=false)
static const G4int sC0Max1Min
virtual G4double GetSurfaceArea()