Geant4 Cross Reference |
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 // G4Trap 27 // 28 // Class description: 29 // 30 // A G4Trap is a general trapezoid: The faces perpendicular to the 31 // z planes are trapezia, and their centres are not necessarily on 32 // a line parallel to the z axis. 33 // 34 // Note that of the 11 parameters described below, only 9 are really 35 // independent - a check for planarity is made in the calculation of the 36 // equation for each plane. If the planes are not parallel, a call to 37 // G4Exception is made. 38 // 39 // pDz Half-length along the z-axis 40 // pTheta Polar angle of the line joining the centres of the faces 41 // at -/+pDz 42 // pPhi Azimuthal angle of the line joing the centre of the face at 43 // -pDz to the centre of the face at +pDz 44 // pDy1 Half-length along y of the face at -pDz 45 // pDx1 Half-length along x of the side at y=-pDy1 of the face at -pDz 46 // pDx2 Half-length along x of the side at y=+pDy1 of the face at -pDz 47 // pAlp1 Angle with respect to the y axis from the centre of the side 48 // at y=-pDy1 to the centre at y=+pDy1 of the face at -pDz 49 // 50 // pDy2 Half-length along y of the face at +pDz 51 // pDx3 Half-length along x of the side at y=-pDy2 of the face at +pDz 52 // pDx4 Half-length along x of the side at y=+pDy2 of the face at +pDz 53 // pAlp2 Angle with respect to the y axis from the centre of the side 54 // at y=-pDy2 to the centre at y=+pDy2 of the face at +pDz 55 // 56 // 57 // Member Data: 58 // 59 // fDz Half-length along the z axis 60 // fTthetaCphi = std::tan(pTheta)*std::cos(pPhi) 61 // fTthetaSphi = std::tan(pTheta)*std::sin(pPhi) 62 // These combinations are suitable for creation of the trapezoid corners 63 // 64 // fDy1 Half-length along y of the face at -fDz 65 // fDx1 Half-length along x of the side at y=-fDy1 of the face at -fDz 66 // fDx2 Half-length along x of the side at y=+fDy1 of the face at -fDz 67 // fTalpha1 Tan of Angle with respect to the y axis from the centre of 68 // the side at y=-fDy1 to the centre at y=+fDy1 of the face 69 // at -fDz 70 // 71 // fDy2 Half-length along y of the face at +fDz 72 // fDx3 Half-length along x of the side at y=-fDy2 of the face at +fDz 73 // fDx4 Half-length along x of the side at y=+fDy2 of the face at +fDz 74 // fTalpha2 Tan of Angle with respect to the y axis from the centre of 75 // the side at y=-fDy2 to the centre at y=+fDy2 of the face 76 // at +fDz 77 // 78 // TrapSidePlane fPlanes[4] Plane equations of the faces not at +/-fDz 79 // NOTE: order is important !!! 80 81 // 23.3.94 P.Kent: Old C++ code converted to tolerant geometry 82 // 9.9.96 V.Grichine: Final modifications before to commit 83 // 8.12.97 J.Allison: Added "nominal" contructor and method SetAllParameters 84 // -------------------------------------------------------------------- 85 #ifndef G4TRAP_HH 86 #define G4TRAP_HH 87 88 #include "G4Types.hh" 89 90 struct TrapSidePlane 91 { 92 G4double a,b,c,d; // Normal unit vector (a,b,c) and offset (d) 93 // => Ax+By+Cz+D=0 94 }; 95 96 #include "G4GeomTypes.hh" 97 98 #if defined(G4GEOM_USE_USOLIDS) 99 #define G4GEOM_USE_UTRAP 1 100 #endif 101 102 #if defined(G4GEOM_USE_UTRAP) 103 #define G4UTrap G4Trap 104 #include "G4UTrap.hh" 105 #else 106 107 #include "G4CSGSolid.hh" 108 109 class G4Trap : public G4CSGSolid 110 { 111 112 public: 113 114 G4Trap( const G4String& pName, 115 G4double pDz, 116 G4double pTheta, G4double pPhi, 117 G4double pDy1, G4double pDx1, G4double pDx2, 118 G4double pAlp1, 119 G4double pDy2, G4double pDx3, G4double pDx4, 120 G4double pAlp2 ); 121 // 122 // The most general constructor for G4Trap which prepares plane 123 // equations and corner coordinates from parameters 124 125 G4Trap( const G4String& pName, 126 const G4ThreeVector pt[8] ) ; 127 // 128 // Prepares plane equations and parameters from corner coordinates 129 130 G4Trap( const G4String& pName, 131 G4double pZ, 132 G4double pY, 133 G4double pX, G4double pLTX ); 134 // 135 // Constructor for Right Angular Wedge from STEP (assumes pLTX<=pX) 136 137 G4Trap( const G4String& pName, 138 G4double pDx1, G4double pDx2, 139 G4double pDy1, G4double pDy2, 140 G4double pDz ); 141 // 142 // Constructor for G4Trd 143 144 G4Trap(const G4String& pName, 145 G4double pDx, G4double pDy, G4double pDz, 146 G4double pAlpha, G4double pTheta, G4double pPhi ); 147 // 148 // Constructor for G4Para 149 150 G4Trap( const G4String& pName ); 151 // 152 // Constructor for "nominal" G4Trap whose parameters are to be set 153 // by a G4VPVParamaterisation later 154 155 ~G4Trap() override ; 156 // 157 // Destructor 158 159 // Accessors 160 161 inline G4double GetZHalfLength() const; 162 inline G4double GetYHalfLength1() const; 163 inline G4double GetXHalfLength1() const; 164 inline G4double GetXHalfLength2() const; 165 inline G4double GetTanAlpha1() const; 166 inline G4double GetYHalfLength2() const; 167 inline G4double GetXHalfLength3() const; 168 inline G4double GetXHalfLength4() const; 169 inline G4double GetTanAlpha2() const; 170 // 171 // Returns coordinates of unit vector along straight 172 // line joining centers of -/+fDz planes 173 174 inline TrapSidePlane GetSidePlane( G4int n ) const; 175 inline G4ThreeVector GetSymAxis() const; 176 177 inline G4double GetPhi() const; 178 inline G4double GetTheta() const; 179 inline G4double GetAlpha1() const; 180 inline G4double GetAlpha2() const; 181 // Obtain (re)computed values of original parameters 182 183 // Modifiers 184 185 void SetAllParameters ( G4double pDz, 186 G4double pTheta, 187 G4double pPhi, 188 G4double pDy1, 189 G4double pDx1, 190 G4double pDx2, 191 G4double pAlp1, 192 G4double pDy2, 193 G4double pDx3, 194 G4double pDx4, 195 G4double pAlp2 ); 196 197 // Methods for solid 198 199 G4double GetCubicVolume() override; 200 G4double GetSurfaceArea() override; 201 202 void ComputeDimensions( G4VPVParameterisation* p, 203 const G4int n, 204 const G4VPhysicalVolume* pRep ) override; 205 206 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const override; 207 208 G4bool CalculateExtent(const EAxis pAxis, 209 const G4VoxelLimits& pVoxelLimit, 210 const G4AffineTransform& pTransform, 211 G4double& pMin, G4double& pMax) const override; 212 213 EInside Inside( const G4ThreeVector& p ) const override; 214 215 G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const override; 216 217 G4double DistanceToIn(const G4ThreeVector& p, 218 const G4ThreeVector& v) const override; 219 220 G4double DistanceToIn( const G4ThreeVector& p ) const override; 221 222 G4double DistanceToOut(const G4ThreeVector& p, const G4ThreeVector& v, 223 const G4bool calcNorm = false, 224 G4bool* validNorm = nullptr, 225 G4ThreeVector* n = nullptr) const override; 226 227 G4double DistanceToOut( const G4ThreeVector& p ) const override; 228 229 G4GeometryType GetEntityType() const override; 230 231 G4ThreeVector GetPointOnSurface() const override; 232 233 G4bool IsFaceted() const override; 234 235 G4VSolid* Clone() const override; 236 237 std::ostream& StreamInfo( std::ostream& os ) const override; 238 239 // Visualisation functions 240 241 void DescribeYourselfTo (G4VGraphicsScene& scene) const override; 242 G4Polyhedron* CreatePolyhedron () const override; 243 244 G4Trap(__void__&); 245 // Fake default constructor for usage restricted to direct object 246 // persistency for clients requiring preallocation of memory for 247 // persistifiable objects. 248 249 G4Trap(const G4Trap& rhs); 250 G4Trap& operator=(const G4Trap& rhs); 251 // Copy constructor and assignment operator. 252 253 protected: 254 255 void MakePlanes(); 256 void MakePlanes( const G4ThreeVector pt[8] ); 257 G4bool MakePlane( const G4ThreeVector& p1, 258 const G4ThreeVector& p2, 259 const G4ThreeVector& p3, 260 const G4ThreeVector& p4, 261 TrapSidePlane& plane ) ; 262 void SetCachedValues(); 263 264 private: 265 266 void CheckParameters(); 267 // Check parameters 268 269 void GetVertices(G4ThreeVector pt[8]) const; 270 // Compute coordinates of the trap vertices from planes 271 272 G4ThreeVector ApproxSurfaceNormal( const G4ThreeVector& p ) const; 273 // Algorithm for SurfaceNormal() following the original 274 // specification for points not on the surface 275 276 private: 277 278 G4double halfCarTolerance; 279 G4double fDz,fTthetaCphi,fTthetaSphi; 280 G4double fDy1,fDx1,fDx2,fTalpha1; 281 G4double fDy2,fDx3,fDx4,fTalpha2; 282 TrapSidePlane fPlanes[4]; 283 G4double fAreas[6]; 284 G4int fTrapType; 285 }; 286 287 #include "G4Trap.icc" 288 289 #endif 290 291 #endif 292