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 // G4Sphere 27 // 28 // Class description: 29 // 30 // A G4Sphere is, in the general case, a section of a spherical shell, 31 // between specified phi and theta angles 32 // 33 // The phi and theta segments are described by a starting angle, 34 // and the +ve delta angle for the shape. 35 // If the delta angle is >=2*pi, or >=pi the shape is treated as 36 // continuous in phi or theta respectively. 37 // 38 // Theta must lie between 0-pi (incl). 39 // 40 // Member Data: 41 // 42 // fRmin inner radius 43 // fRmax outer radius 44 // 45 // fSPhi starting angle of the segment in radians 46 // fDPhi delta angle of the segment in radians 47 // 48 // fSTheta starting angle of the segment in radians 49 // fDTheta delta angle of the segment in radians 50 // 51 // 52 // Note: 53 // Internally fSPhi & fDPhi are adjusted so that fDPhi<=2PI, 54 // and fDPhi+fSPhi<=2PI. This enables simpler comparisons to be 55 // made with (say) Phi of a point. 56 57 // 28.3.94 P.Kent: old C++ code converted to tolerant geometry 58 // 17.9.96 V.Grichine: final modifications to commit 59 // -------------------------------------------------------------------- 60 #ifndef G4SPHERE_HH 61 #define G4SPHERE_HH 62 63 #include "G4GeomTypes.hh" 64 65 #if defined(G4GEOM_USE_USOLIDS) 66 #define G4GEOM_USE_USPHERE 1 67 #endif 68 69 #if defined(G4GEOM_USE_USPHERE) 70 #define G4USphere G4Sphere 71 #include "G4USphere.hh" 72 #else 73 74 #include <CLHEP/Units/PhysicalConstants.h> 75 #include "G4CSGSolid.hh" 76 #include "G4Polyhedron.hh" 77 78 class G4VisExtent; 79 80 class G4Sphere : public G4CSGSolid 81 { 82 public: 83 84 G4Sphere(const G4String& pName, 85 G4double pRmin, G4double pRmax, 86 G4double pSPhi, G4double pDPhi, 87 G4double pSTheta, G4double pDTheta); 88 // 89 // Constructs a sphere or sphere shell section 90 // with the given name and dimensions 91 92 ~G4Sphere() override; 93 // 94 // Destructor 95 96 // Accessors 97 98 inline G4double GetInnerRadius () const; 99 inline G4double GetOuterRadius () const; 100 inline G4double GetStartPhiAngle () const; 101 inline G4double GetDeltaPhiAngle () const; 102 inline G4double GetStartThetaAngle() const; 103 inline G4double GetDeltaThetaAngle() const; 104 inline G4double GetSinStartPhi () const; 105 inline G4double GetCosStartPhi () const; 106 inline G4double GetSinEndPhi () const; 107 inline G4double GetCosEndPhi () const; 108 inline G4double GetSinStartTheta () const; 109 inline G4double GetCosStartTheta () const; 110 inline G4double GetSinEndTheta () const; 111 inline G4double GetCosEndTheta () const; 112 113 // Modifiers 114 115 inline void SetInnerRadius (G4double newRMin); 116 inline void SetOuterRadius (G4double newRmax); 117 inline void SetStartPhiAngle (G4double newSphi, G4bool trig = true); 118 inline void SetDeltaPhiAngle (G4double newDphi); 119 inline void SetStartThetaAngle(G4double newSTheta); 120 inline void SetDeltaThetaAngle(G4double newDTheta); 121 122 // Methods for solid 123 124 G4double GetCubicVolume() override; 125 G4double GetSurfaceArea() override; 126 127 void ComputeDimensions( G4VPVParameterisation* p, 128 const G4int n, 129 const G4VPhysicalVolume* pRep) override; 130 131 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const override; 132 133 G4bool CalculateExtent(const EAxis pAxis, 134 const G4VoxelLimits& pVoxelLimit, 135 const G4AffineTransform& pTransform, 136 G4double& pmin, G4double& pmax) const override; 137 138 EInside Inside(const G4ThreeVector& p) const override; 139 140 G4ThreeVector SurfaceNormal( const G4ThreeVector& p) const override; 141 142 G4double DistanceToIn(const G4ThreeVector& p, 143 const G4ThreeVector& v) const override; 144 145 G4double DistanceToIn(const G4ThreeVector& p) const override; 146 147 G4double DistanceToOut(const G4ThreeVector& p, 148 const G4ThreeVector& v, 149 const G4bool calcNorm = false, 150 G4bool* validNorm = nullptr, 151 G4ThreeVector* n = nullptr) const override; 152 153 G4double DistanceToOut(const G4ThreeVector& p) const override; 154 155 G4GeometryType GetEntityType() const override; 156 157 G4ThreeVector GetPointOnSurface() const override; 158 159 G4VSolid* Clone() const override; 160 161 std::ostream& StreamInfo(std::ostream& os) const override; 162 163 // Visualisation functions 164 165 G4VisExtent GetExtent () const override; 166 void DescribeYourselfTo(G4VGraphicsScene& scene) const override; 167 G4Polyhedron* CreatePolyhedron() const override; 168 169 G4Sphere(__void__&); 170 // 171 // Fake default constructor for usage restricted to direct object 172 // persistency for clients requiring preallocation of memory for 173 // persistifiable objects. 174 175 G4Sphere(const G4Sphere& rhs); 176 G4Sphere& operator=(const G4Sphere& rhs); 177 // Copy constructor and assignment operator. 178 179 180 private: 181 182 inline void Initialize(); 183 // 184 // Reset relevant values to zero 185 186 inline void CheckThetaAngles(G4double sTheta, G4double dTheta); 187 inline void CheckSPhiAngle(G4double sPhi); 188 inline void CheckDPhiAngle(G4double dPhi); 189 inline void CheckPhiAngles(G4double sPhi, G4double dPhi); 190 // 191 // Reset relevant flags and angle values 192 193 inline void InitializePhiTrigonometry(); 194 inline void InitializeThetaTrigonometry(); 195 // 196 // Recompute relevant trigonometric values and cache them 197 198 G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector& p) const; 199 // 200 // Algorithm for SurfaceNormal() following the original 201 // specification for points not on the surface 202 203 private: 204 205 // Used by distanceToOut 206 // 207 enum ESide {kNull,kRMin,kRMax,kSPhi,kEPhi,kSTheta,kETheta}; 208 209 // used by normal 210 // 211 enum ENorm {kNRMin,kNRMax,kNSPhi,kNEPhi,kNSTheta,kNETheta}; 212 213 G4double fRminTolerance, fRmaxTolerance, kAngTolerance, 214 kRadTolerance, fEpsilon = 2.e-11; 215 // 216 // Radial and angular tolerances 217 218 G4double fRmin, fRmax, fSPhi, fDPhi, fSTheta, fDTheta; 219 // 220 // Radial and angular dimensions 221 222 G4double sinCPhi, cosCPhi, cosHDPhi, cosHDPhiOT, cosHDPhiIT, 223 sinSPhi, cosSPhi, sinEPhi, cosEPhi, hDPhi, cPhi, ePhi; 224 // 225 // Cached trigonometric values for Phi angle 226 227 G4double sinSTheta, cosSTheta, sinETheta, cosETheta, 228 tanSTheta, tanSTheta2, tanETheta, tanETheta2, eTheta; 229 // 230 // Cached trigonometric values for Theta angle 231 232 G4bool fFullPhiSphere=false, fFullThetaSphere=false, fFullSphere=true; 233 // 234 // Flags for identification of section, shell or full sphere 235 236 G4double halfCarTolerance, halfAngTolerance; 237 // 238 // Cached half tolerance values 239 }; 240 241 #include "G4Sphere.icc" 242 243 #endif 244 245 #endif 246