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 // G4Cons 27 // 28 // Class description: 29 // 30 // A G4Cons is, in the general case, a Phi segment of a cone, with 31 // half-length fDz, inner and outer radii specified at -fDz and +fDz. 32 // The Phi segment is described by a starting fSPhi angle, and the 33 // +fDPhi delta angle for the shape. 34 // If the delta angle is >=2*pi, the shape is treated as continuous 35 // in Phi 36 // 37 // Member Data: 38 // 39 // fRmin1 inside radius at -fDz 40 // fRmin2 inside radius at +fDz 41 // fRmax1 outside radius at -fDz 42 // fRmax2 outside radius at +fDz 43 // fDz half length in z 44 // 45 // fSPhi starting angle of the segment in radians 46 // fDPhi delta angle of the segment in radians 47 // 48 // fPhiFullCone Boolean variable used for indicate the Phi Section 49 // 50 // Note: 51 // Internally fSPhi & fDPhi are adjusted so that fDPhi<=2PI, 52 // and fDPhi+fSPhi<=2PI. This enables simpler comparisons to be 53 // made with (say) Phi of a point. 54 55 // 19.3.94 P.Kent: Old C++ code converted to tolerant geometry 56 // 13.9.96 V.Grichine: Final modifications to commit 57 // -------------------------------------------------------------------- 58 #ifndef G4CONS_HH 59 #define G4CONS_HH 60 61 #include "G4GeomTypes.hh" 62 63 #if defined(G4GEOM_USE_USOLIDS) 64 #define G4GEOM_USE_UCONS 1 65 #endif 66 67 #if defined(G4GEOM_USE_UCONS) 68 #define G4UCons G4Cons 69 #include "G4UCons.hh" 70 #else 71 72 #include <CLHEP/Units/PhysicalConstants.h> 73 74 #include "G4CSGSolid.hh" 75 #include "G4Polyhedron.hh" 76 77 class G4Cons : public G4CSGSolid 78 { 79 public: 80 81 G4Cons(const G4String& pName, 82 G4double pRmin1, G4double pRmax1, 83 G4double pRmin2, G4double pRmax2, 84 G4double pDz, 85 G4double pSPhi, G4double pDPhi); 86 // 87 // Constructs a cone with the given name and dimensions 88 89 ~G4Cons() override ; 90 // 91 // Destructor 92 93 // Accessors 94 95 inline G4double GetInnerRadiusMinusZ() const; 96 inline G4double GetOuterRadiusMinusZ() const; 97 inline G4double GetInnerRadiusPlusZ() const; 98 inline G4double GetOuterRadiusPlusZ() const; 99 inline G4double GetZHalfLength() const; 100 inline G4double GetStartPhiAngle() const; 101 inline G4double GetDeltaPhiAngle() const; 102 inline G4double GetSinStartPhi() const; 103 inline G4double GetCosStartPhi() const; 104 inline G4double GetSinEndPhi() const; 105 inline G4double GetCosEndPhi() const; 106 107 // Modifiers 108 109 inline void SetInnerRadiusMinusZ (G4double Rmin1 ); 110 inline void SetOuterRadiusMinusZ (G4double Rmax1 ); 111 inline void SetInnerRadiusPlusZ (G4double Rmin2 ); 112 inline void SetOuterRadiusPlusZ (G4double Rmax2 ); 113 inline void SetZHalfLength (G4double newDz ); 114 inline void SetStartPhiAngle (G4double newSPhi, G4bool trig=true); 115 inline void SetDeltaPhiAngle (G4double newDPhi); 116 117 // Other methods for solid 118 119 inline G4double GetCubicVolume() override; 120 inline G4double GetSurfaceArea() override; 121 122 void ComputeDimensions( G4VPVParameterisation* p, 123 const G4int n, 124 const G4VPhysicalVolume* pRep) override; 125 126 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const override; 127 128 G4bool CalculateExtent(const EAxis pAxis, 129 const G4VoxelLimits& pVoxelLimit, 130 const G4AffineTransform& pTransform, 131 G4double& pMin, G4double& pMax) const override; 132 133 EInside Inside( const G4ThreeVector& p ) const override; 134 135 G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const override; 136 137 G4double DistanceToIn (const G4ThreeVector& p, 138 const G4ThreeVector& v) const override; 139 G4double DistanceToIn (const G4ThreeVector& p) const override; 140 G4double DistanceToOut(const G4ThreeVector& p, 141 const G4ThreeVector& v, 142 const G4bool calcNorm = false, 143 G4bool* validNorm = nullptr, 144 G4ThreeVector* n = nullptr) const override; 145 G4double DistanceToOut(const G4ThreeVector& p) const override; 146 147 G4GeometryType GetEntityType() const override; 148 149 G4ThreeVector GetPointOnSurface() const override; 150 151 G4VSolid* Clone() const override; 152 153 std::ostream& StreamInfo(std::ostream& os) const override; 154 155 // Visualisation functions 156 157 void DescribeYourselfTo( G4VGraphicsScene& scene ) const override; 158 G4Polyhedron* CreatePolyhedron() const override; 159 160 G4Cons(__void__&); 161 // 162 // Fake default constructor for usage restricted to direct object 163 // persistency for clients requiring preallocation of memory for 164 // persistifiable objects. 165 166 G4Cons(const G4Cons& rhs); 167 G4Cons& operator=(const G4Cons& rhs); 168 // Copy constructor and assignment operator. 169 170 private: 171 172 inline void Initialize(); 173 // 174 // Reset relevant values to zero 175 176 inline void CheckSPhiAngle(G4double sPhi); 177 inline void CheckDPhiAngle(G4double dPhi); 178 inline void CheckPhiAngles(G4double sPhi, G4double dPhi); 179 // 180 // Reset relevant flags and angle values 181 182 inline void InitializeTrigonometry(); 183 // 184 // Recompute relevant trigonometric values and cache them 185 186 G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector& p) const; 187 // 188 // Algorithm for SurfaceNormal() following the original 189 // specification for points not on the surface 190 191 private: 192 193 // Used by distanceToOut 194 // 195 enum ESide {kNull,kRMin,kRMax,kSPhi,kEPhi,kPZ,kMZ}; 196 197 // used by normal 198 // 199 enum ENorm {kNRMin,kNRMax,kNSPhi,kNEPhi,kNZ}; 200 201 G4double kRadTolerance, kAngTolerance; 202 // 203 // Radial and angular tolerances 204 205 G4double fRmin1, fRmin2, fRmax1, fRmax2, fDz, fSPhi, fDPhi; 206 // 207 // Radial and angular dimensions 208 209 G4double sinCPhi, cosCPhi, cosHDPhi, cosHDPhiOT, cosHDPhiIT, 210 sinSPhi, cosSPhi, sinEPhi, cosEPhi; 211 // 212 // Cached trigonometric values 213 214 G4bool fPhiFullCone = false; 215 // 216 // Flag for identification of section or full cone 217 218 G4double halfCarTolerance, halfRadTolerance, halfAngTolerance; 219 // 220 // Cached half tolerance values 221 }; 222 223 #include "G4Cons.icc" 224 225 #endif 226 227 #endif 228