File indexing completed on 2025-01-18 09:59:27
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
0035 #ifndef G4VTWISTEDFACETED_HH
0036 #define G4VTWISTEDFACETED_HH
0037
0038 #include "G4VSolid.hh"
0039 #include "G4TwoVector.hh"
0040 #include "G4TwistTrapAlphaSide.hh"
0041 #include "G4TwistTrapParallelSide.hh"
0042 #include "G4TwistBoxSide.hh"
0043 #include "G4TwistTrapFlatSide.hh"
0044
0045 class G4SolidExtentList;
0046 class G4ClippablePolygon;
0047
0048 class G4VTwistedFaceted: public G4VSolid
0049 {
0050 public:
0051
0052 G4VTwistedFaceted(const G4String& pname,
0053 G4double PhiTwist,
0054 G4double pDz,
0055 G4double pTheta,
0056 G4double pPhi,
0057 G4double pDy1,
0058 G4double pDx1,
0059 G4double pDx2,
0060 G4double pDy2,
0061 G4double pDx3,
0062 G4double pDx4,
0063 G4double pAlph
0064 );
0065
0066 virtual ~G4VTwistedFaceted();
0067
0068 void ComputeDimensions( G4VPVParameterisation*,
0069 const G4int,
0070 const G4VPhysicalVolume* ) override;
0071
0072 void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override;
0073
0074 G4bool CalculateExtent(const EAxis pAxis,
0075 const G4VoxelLimits& pVoxelLimit,
0076 const G4AffineTransform& pTransform,
0077 G4double& pMin,
0078 G4double& pMax ) const override;
0079
0080 G4double DistanceToIn (const G4ThreeVector& p,
0081 const G4ThreeVector& v ) const override;
0082
0083 G4double DistanceToIn (const G4ThreeVector& p ) const override;
0084
0085 G4double DistanceToOut(const G4ThreeVector& p,
0086 const G4ThreeVector& v,
0087 const G4bool calcnorm = false,
0088 G4bool* validnorm = nullptr,
0089 G4ThreeVector* n = nullptr ) const override;
0090
0091 G4double DistanceToOut(const G4ThreeVector& p) const override;
0092
0093 EInside Inside (const G4ThreeVector& p) const override;
0094
0095 G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const override;
0096
0097 G4ThreeVector GetPointOnSurface() const override;
0098 G4ThreeVector GetPointInSolid(G4double z) const;
0099
0100 G4double GetCubicVolume() override;
0101 G4double GetSurfaceArea() override;
0102
0103 void DescribeYourselfTo (G4VGraphicsScene& scene) const override;
0104 G4Polyhedron* CreatePolyhedron () const override;
0105 G4Polyhedron* GetPolyhedron () const override;
0106
0107 std::ostream &StreamInfo(std::ostream& os) const override;
0108
0109
0110
0111 inline G4double GetTwistAngle() const { return fPhiTwist; }
0112
0113 inline G4double GetDx1 () const { return fDx1 ; }
0114 inline G4double GetDx2 () const { return fDx2 ; }
0115 inline G4double GetDx3 () const { return fDx3 ; }
0116 inline G4double GetDx4 () const { return fDx4 ; }
0117 inline G4double GetDy1 () const { return fDy1 ; }
0118 inline G4double GetDy2 () const { return fDy2 ; }
0119 inline G4double GetDz () const { return fDz ; }
0120 inline G4double GetPhi () const { return fPhi ; }
0121 inline G4double GetTheta () const { return fTheta ; }
0122 inline G4double GetAlpha () const { return fAlph ; }
0123
0124 inline G4double Xcoef(G4double u, G4double phi, G4double ftg) const;
0125
0126
0127 inline G4double GetValueA(G4double phi) const;
0128 inline G4double GetValueB(G4double phi) const;
0129 inline G4double GetValueD(G4double phi) const;
0130
0131 G4VisExtent GetExtent () const override;
0132 G4GeometryType GetEntityType() const override;
0133
0134 G4VTwistedFaceted(__void__&);
0135
0136
0137
0138
0139 G4VTwistedFaceted(const G4VTwistedFaceted& rhs);
0140 G4VTwistedFaceted& operator=(const G4VTwistedFaceted& rhs);
0141
0142
0143 protected:
0144
0145 mutable G4bool fRebuildPolyhedron = false;
0146 mutable G4Polyhedron* fpPolyhedron = nullptr;
0147
0148 G4double fCubicVolume = 0.0;
0149 G4double fSurfaceArea = 0.0;
0150
0151 private:
0152
0153 double GetLateralFaceArea(const G4TwoVector& p1,
0154 const G4TwoVector& p2,
0155 const G4TwoVector& p3,
0156 const G4TwoVector& p4) const;
0157 void CreateSurfaces();
0158
0159 private:
0160
0161 G4double fTheta;
0162 G4double fPhi ;
0163
0164 G4double fDy1;
0165 G4double fDx1;
0166 G4double fDx2;
0167
0168 G4double fDy2;
0169 G4double fDx3;
0170 G4double fDx4;
0171
0172 G4double fDz;
0173
0174 G4double fDx ;
0175 G4double fDy ;
0176
0177 G4double fAlph ;
0178 G4double fTAlph ;
0179
0180 G4double fdeltaX ;
0181 G4double fdeltaY ;
0182
0183 G4double fPhiTwist;
0184
0185 G4VTwistSurface* fLowerEndcap ;
0186 G4VTwistSurface* fUpperEndcap ;
0187
0188 G4VTwistSurface* fSide0 ;
0189 G4VTwistSurface* fSide90 ;
0190 G4VTwistSurface* fSide180 ;
0191 G4VTwistSurface* fSide270 ;
0192
0193 };
0194
0195
0196
0197 inline
0198 G4double G4VTwistedFaceted::GetValueA(G4double phi) const
0199 {
0200 return ( fDx4 + fDx2 + ( fDx4 - fDx2 ) * ( 2 * phi ) / fPhiTwist ) ;
0201 }
0202
0203 inline
0204 G4double G4VTwistedFaceted::GetValueD(G4double phi) const
0205 {
0206 return ( fDx3 + fDx1 + ( fDx3 - fDx1 ) * ( 2 * phi ) / fPhiTwist ) ;
0207 }
0208
0209 inline
0210 G4double G4VTwistedFaceted::GetValueB(G4double phi) const
0211 {
0212 return ( fDy2 + fDy1 + ( fDy2 - fDy1 ) * ( 2 * phi ) / fPhiTwist ) ;
0213 }
0214
0215 inline
0216 G4double G4VTwistedFaceted::Xcoef(G4double u, G4double phi, G4double ftg) const
0217 {
0218 return GetValueA(phi)/2. + (GetValueD(phi)-GetValueA(phi))/4.
0219 - u*( ( GetValueD(phi)-GetValueA(phi) ) / ( 2 * GetValueB(phi) ) - ftg );
0220 }
0221
0222 #endif