File indexing completed on 2025-01-18 09:59:15
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034 #ifndef G4TWISTTRAPPARALLELSIDE_HH
0035 #define G4TWISTTRAPPARALLELSIDE_HH
0036
0037 #include "G4VTwistSurface.hh"
0038
0039 #include <vector>
0040
0041 class G4TwistTrapParallelSide : public G4VTwistSurface
0042 {
0043 public:
0044
0045 G4TwistTrapParallelSide(const G4String& name,
0046 G4double PhiTwist,
0047 G4double pDz,
0048 G4double pTheta,
0049 G4double pPhi,
0050 G4double pDy1,
0051 G4double pDx1,
0052 G4double pDx2,
0053 G4double pDy2,
0054 G4double pDx3,
0055 G4double pDx4,
0056 G4double pAlph,
0057 G4double AngleSide
0058 );
0059
0060 ~G4TwistTrapParallelSide() override;
0061
0062 G4ThreeVector GetNormal(const G4ThreeVector& xx,
0063 G4bool isGlobal = false) override ;
0064
0065 G4int DistanceToSurface(const G4ThreeVector& gp,
0066 const G4ThreeVector& gv,
0067 G4ThreeVector gxx[],
0068 G4double distance[],
0069 G4int areacode[],
0070 G4bool isvalid[],
0071 EValidate validate = kValidateWithTol) override;
0072
0073 G4int DistanceToSurface(const G4ThreeVector& gp,
0074 G4ThreeVector gxx[],
0075 G4double distance[],
0076 G4int areacode[]) override;
0077
0078 G4TwistTrapParallelSide(__void__&);
0079
0080
0081
0082
0083 private:
0084
0085 G4int GetAreaCode(const G4ThreeVector& xx,
0086 G4bool withTol = true) override;
0087 void SetCorners() override;
0088 void SetBoundaries() override;
0089
0090 void GetPhiUAtX(const G4ThreeVector& p, G4double& phi, G4double& u);
0091 G4ThreeVector ProjectPoint(const G4ThreeVector& p,
0092 G4bool isglobal = false);
0093
0094 inline G4ThreeVector SurfacePoint(G4double phi, G4double u,
0095 G4bool isGlobal = false) override;
0096 inline G4double GetBoundaryMin(G4double phi) override;
0097 inline G4double GetBoundaryMax(G4double phi) override;
0098 inline G4double GetSurfaceArea() override;
0099 void GetFacets( G4int m, G4int n, G4double xyz[][3],
0100 G4int faces[][4], G4int iside ) override;
0101
0102 inline G4ThreeVector NormAng(G4double phi, G4double u);
0103 inline G4double GetValueB(G4double phi) ;
0104 inline G4double Xcoef(G4double u);
0105
0106
0107 private:
0108
0109 G4double fTheta;
0110 G4double fPhi ;
0111
0112 G4double fDy1;
0113 G4double fDx1;
0114 G4double fDx2;
0115
0116 G4double fDy2;
0117 G4double fDx3;
0118 G4double fDx4;
0119
0120 G4double fDz;
0121
0122 G4double fAlph;
0123 G4double fTAlph;
0124
0125 G4double fPhiTwist;
0126
0127 G4double fAngleSide;
0128
0129 G4double fdeltaX;
0130 G4double fdeltaY;
0131
0132 G4double fDx4plus2;
0133 G4double fDx4minus2;
0134 G4double fDx3plus1;
0135 G4double fDx3minus1;
0136 G4double fDy2plus1;
0137 G4double fDy2minus1;
0138 G4double fa1md1;
0139 G4double fa2md2;
0140 };
0141
0142
0143
0144
0145
0146
0147 inline
0148 G4double G4TwistTrapParallelSide::GetValueB(G4double phi)
0149 {
0150 return ( fDy2plus1 + fDy2minus1 * ( 2 * phi ) / fPhiTwist ) ;
0151 }
0152
0153 inline
0154 G4double G4TwistTrapParallelSide::Xcoef(G4double phi)
0155 {
0156 return GetValueB(phi)/2. ;
0157 }
0158
0159 inline G4ThreeVector
0160 G4TwistTrapParallelSide::
0161 SurfacePoint( G4double phi, G4double u, G4bool isGlobal )
0162 {
0163
0164
0165 G4ThreeVector SurfPoint ( u*std::cos(phi) - Xcoef(phi)*std::sin(phi)
0166 + fdeltaX*phi/fPhiTwist,
0167 u*std::sin(phi) + Xcoef(phi)*std::cos(phi)
0168 + fdeltaY*phi/fPhiTwist,
0169 2*fDz*phi/fPhiTwist );
0170 if (isGlobal) { return (fRot * SurfPoint + fTrans); }
0171 return SurfPoint;
0172 }
0173
0174
0175 inline
0176 G4double G4TwistTrapParallelSide::GetBoundaryMin(G4double phi)
0177 {
0178 return -(fPhiTwist*(fDx2 + fDx4 - fDy2plus1*fTAlph)
0179 + 2*fDx4minus2*phi - 2*fDy2minus1*fTAlph*phi)/(2.*fPhiTwist) ;
0180 }
0181
0182 inline
0183 G4double G4TwistTrapParallelSide::GetBoundaryMax(G4double phi)
0184 {
0185 return (fDx2 + fDx4 + fDy2plus1*fTAlph)/ 2.
0186 + ((fDx4minus2 + fDy2minus1*fTAlph)*phi)/fPhiTwist ;
0187 }
0188
0189 inline
0190 G4double G4TwistTrapParallelSide::GetSurfaceArea()
0191 {
0192 return 2*fDx4plus2*fDz ;
0193 }
0194
0195 inline
0196 G4ThreeVector G4TwistTrapParallelSide::NormAng( G4double phi, G4double u )
0197 {
0198
0199
0200
0201 G4ThreeVector nvec(-2*fDz*std::sin(phi) ,
0202 2*fDz*std::cos(phi) ,
0203 -(fDy2minus1 + fPhiTwist*u + fdeltaY*std::cos(phi)
0204 -fdeltaX*std::sin(phi))) ;
0205 return nvec.unit();
0206 }
0207
0208 #endif