00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044 #ifndef __G4TWISTTRAPPARALLELSIDE__
00045 #define __G4TWISTTRAPPARALLELSIDE__
00046
00047 #include "G4VTwistSurface.hh"
00048
00049 #include <vector>
00050
00051 class G4TwistTrapParallelSide : public G4VTwistSurface
00052 {
00053 public:
00054
00055 G4TwistTrapParallelSide(const G4String &name,
00056 G4double PhiTwist,
00057 G4double pDz,
00058 G4double pTheta,
00059 G4double pPhi,
00060 G4double pDy1,
00061 G4double pDx1,
00062 G4double pDx2,
00063 G4double pDy2,
00064 G4double pDx3,
00065 G4double pDx4,
00066 G4double pAlph,
00067 G4double AngleSide
00068 );
00069
00070 virtual ~G4TwistTrapParallelSide();
00071
00072 virtual G4ThreeVector GetNormal(const G4ThreeVector &xx,
00073 G4bool isGlobal = false) ;
00074
00075 virtual G4int DistanceToSurface(const G4ThreeVector &gp,
00076 const G4ThreeVector &gv,
00077 G4ThreeVector gxx[],
00078 G4double distance[],
00079 G4int areacode[],
00080 G4bool isvalid[],
00081 EValidate validate = kValidateWithTol);
00082
00083 virtual G4int DistanceToSurface(const G4ThreeVector &gp,
00084 G4ThreeVector gxx[],
00085 G4double distance[],
00086 G4int areacode[]);
00087
00088 public:
00089
00090 G4TwistTrapParallelSide(__void__&);
00091
00092
00093
00094
00095 private:
00096
00097 virtual G4int GetAreaCode(const G4ThreeVector &xx,
00098 G4bool withTol = true);
00099 virtual void SetCorners();
00100 virtual void SetBoundaries();
00101
00102 void GetPhiUAtX(G4ThreeVector p, G4double &phi, G4double &u);
00103 G4ThreeVector ProjectPoint(const G4ThreeVector &p,
00104 G4bool isglobal = false);
00105
00106 virtual G4ThreeVector SurfacePoint(G4double phi, G4double u,
00107 G4bool isGlobal = false);
00108 virtual G4double GetBoundaryMin(G4double phi);
00109 virtual G4double GetBoundaryMax(G4double phi);
00110 virtual G4double GetSurfaceArea();
00111 virtual void GetFacets( G4int m, G4int n, G4double xyz[][3],
00112 G4int faces[][4], G4int iside );
00113
00114 inline G4ThreeVector NormAng(G4double phi, G4double u);
00115 inline G4double GetValueB(G4double phi) ;
00116 inline G4double Xcoef(G4double u);
00117
00118
00119 private:
00120
00121 G4double fTheta;
00122 G4double fPhi ;
00123
00124 G4double fDy1;
00125 G4double fDx1;
00126 G4double fDx2;
00127
00128 G4double fDy2;
00129 G4double fDx3;
00130 G4double fDx4;
00131
00132 G4double fDz;
00133
00134 G4double fAlph;
00135 G4double fTAlph;
00136
00137 G4double fPhiTwist;
00138
00139 G4double fAngleSide;
00140
00141 G4double fdeltaX;
00142 G4double fdeltaY;
00143
00144 G4double fDx4plus2;
00145 G4double fDx4minus2;
00146 G4double fDx3plus1;
00147 G4double fDx3minus1;
00148 G4double fDy2plus1;
00149 G4double fDy2minus1;
00150 G4double fa1md1;
00151 G4double fa2md2;
00152 };
00153
00154
00155
00156
00157
00158
00159 inline
00160 G4double G4TwistTrapParallelSide::GetValueB(G4double phi)
00161 {
00162 return ( fDy2plus1 + fDy2minus1 * ( 2 * phi ) / fPhiTwist ) ;
00163 }
00164
00165 inline
00166 G4double G4TwistTrapParallelSide::Xcoef(G4double phi)
00167 {
00168 return GetValueB(phi)/2. ;
00169 }
00170
00171 inline G4ThreeVector
00172 G4TwistTrapParallelSide::
00173 SurfacePoint( G4double phi, G4double u, G4bool isGlobal )
00174 {
00175
00176
00177 G4ThreeVector SurfPoint ( u*std::cos(phi) - Xcoef(phi)*std::sin(phi)
00178 + fdeltaX*phi/fPhiTwist,
00179 u*std::sin(phi) + Xcoef(phi)*std::cos(phi)
00180 + fdeltaY*phi/fPhiTwist,
00181 2*fDz*phi/fPhiTwist );
00182 if (isGlobal) { return (fRot * SurfPoint + fTrans); }
00183 return SurfPoint;
00184 }
00185
00186
00187 inline
00188 G4double G4TwistTrapParallelSide::GetBoundaryMin(G4double phi)
00189 {
00190 return -(fPhiTwist*(fDx2 + fDx4 - fDy2plus1*fTAlph)
00191 + 2*fDx4minus2*phi - 2*fDy2minus1*fTAlph*phi)/(2.*fPhiTwist) ;
00192 }
00193
00194 inline
00195 G4double G4TwistTrapParallelSide::GetBoundaryMax(G4double phi)
00196 {
00197 return (fDx2 + fDx4 + fDy2plus1*fTAlph)/ 2.
00198 + ((fDx4minus2 + fDy2minus1*fTAlph)*phi)/fPhiTwist ;
00199 }
00200
00201 inline
00202 G4double G4TwistTrapParallelSide::GetSurfaceArea()
00203 {
00204 return 2*fDx4plus2*fDz ;
00205 }
00206
00207 inline
00208 G4ThreeVector G4TwistTrapParallelSide::NormAng( G4double phi, G4double u )
00209 {
00210
00211
00212
00213 G4ThreeVector nvec(-2*fDz*std::sin(phi) ,
00214 2*fDz*std::cos(phi) ,
00215 -(fDy2minus1 + fPhiTwist*u + fdeltaY*std::cos(phi)
00216 -fdeltaX*std::sin(phi))) ;
00217 return nvec.unit();
00218 }
00219
00220 #endif