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 // 27 // 28 // 29 //////////////////////////////////////////////////////////////////////// 30 // Optical Photon Rayleigh Scattering Class Definition 31 //////////////////////////////////////////////////////////////////////// 32 // 33 // File: G4OpRayleigh.hh 34 // Description: Discrete Process -- Rayleigh scattering of optical photons 35 // Version: 1.0 36 // Created: 1996-05-31 37 // Author: Juliet Armstrong 38 // Updated: 2014-08-20 allow for more material types 39 // 2005-07-28 add G4ProcessType to constructor 40 // 1999-10-29 add method and class descriptors 41 // 1997-04-09 by Peter Gumplinger 42 // > new physics/tracking scheme 43 // 44 //////////////////////////////////////////////////////////////////////// 45 46 #ifndef G4OpRayleigh_h 47 #define G4OpRayleigh_h 1 48 49 #include "G4VDiscreteProcess.hh" 50 #include "G4OpticalPhoton.hh" 51 #include "G4PhysicsTable.hh" 52 53 class G4OpRayleigh : public G4VDiscreteProcess 54 { 55 public: 56 explicit G4OpRayleigh(const G4String& processName = "OpRayleigh", 57 G4ProcessType type = fOptical); 58 virtual ~G4OpRayleigh(); 59 60 virtual G4bool IsApplicable( 61 const G4ParticleDefinition& aParticleType) override; 62 // Returns true -> 'is applicable' only for an optical photon. 63 64 virtual void BuildPhysicsTable( 65 const G4ParticleDefinition& aParticleType) override; 66 // Build thePhysicsTable at a right time 67 68 virtual G4double GetMeanFreePath(const G4Track& aTrack, G4double, 69 G4ForceCondition*) override; 70 // Returns the mean free path for Rayleigh scattering 71 72 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 73 const G4Step& aStep) override; 74 // This is the method implementing Rayleigh scattering. 75 76 virtual G4PhysicsTable* GetPhysicsTable() const; 77 // Returns the address of the physics table. 78 79 virtual void DumpPhysicsTable() const; 80 // Prints the physics table. 81 82 virtual void PreparePhysicsTable(const G4ParticleDefinition&) override; 83 virtual void Initialise(); 84 85 void SetVerboseLevel(G4int); 86 87 protected: 88 G4PhysicsTable* thePhysicsTable; 89 90 private: 91 G4OpRayleigh(const G4OpRayleigh& right) = delete; 92 G4OpRayleigh& operator=(const G4OpRayleigh& right) = delete; 93 94 /// Calculates the mean free paths for a material as a function of 95 /// photon energy 96 G4PhysicsFreeVector* CalculateRayleighMeanFreePaths( 97 const G4Material* material) const; 98 99 size_t idx_rslength = 0; 100 }; 101 102 //////////////////// 103 // Inline methods 104 //////////////////// 105 106 inline G4bool G4OpRayleigh::IsApplicable( 107 const G4ParticleDefinition& aParticleType) 108 { 109 return (&aParticleType == G4OpticalPhoton::OpticalPhoton()); 110 } 111 112 inline void G4OpRayleigh::DumpPhysicsTable() const 113 { 114 for(size_t i = 0; i < thePhysicsTable->entries(); ++i) 115 { 116 ((G4PhysicsFreeVector*) (*thePhysicsTable)[i])->DumpValues(); 117 } 118 } 119 120 inline G4PhysicsTable* G4OpRayleigh::GetPhysicsTable() const 121 { 122 return thePhysicsTable; 123 } 124 125 #endif /* G4OpRayleigh_h */ 126