|
||||
File indexing completed on 2025-01-18 09:58:48
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 // 0028 // 0029 //////////////////////////////////////////////////////////////////////// 0030 // Optical Photon Rayleigh Scattering Class Definition 0031 //////////////////////////////////////////////////////////////////////// 0032 // 0033 // File: G4OpRayleigh.hh 0034 // Description: Discrete Process -- Rayleigh scattering of optical photons 0035 // Version: 1.0 0036 // Created: 1996-05-31 0037 // Author: Juliet Armstrong 0038 // Updated: 2014-08-20 allow for more material types 0039 // 2005-07-28 add G4ProcessType to constructor 0040 // 1999-10-29 add method and class descriptors 0041 // 1997-04-09 by Peter Gumplinger 0042 // > new physics/tracking scheme 0043 // 0044 //////////////////////////////////////////////////////////////////////// 0045 0046 #ifndef G4OpRayleigh_h 0047 #define G4OpRayleigh_h 1 0048 0049 #include "G4VDiscreteProcess.hh" 0050 #include "G4OpticalPhoton.hh" 0051 #include "G4PhysicsTable.hh" 0052 0053 class G4OpRayleigh : public G4VDiscreteProcess 0054 { 0055 public: 0056 explicit G4OpRayleigh(const G4String& processName = "OpRayleigh", 0057 G4ProcessType type = fOptical); 0058 virtual ~G4OpRayleigh(); 0059 0060 virtual G4bool IsApplicable( 0061 const G4ParticleDefinition& aParticleType) override; 0062 // Returns true -> 'is applicable' only for an optical photon. 0063 0064 virtual void BuildPhysicsTable( 0065 const G4ParticleDefinition& aParticleType) override; 0066 // Build thePhysicsTable at a right time 0067 0068 virtual G4double GetMeanFreePath(const G4Track& aTrack, G4double, 0069 G4ForceCondition*) override; 0070 // Returns the mean free path for Rayleigh scattering 0071 0072 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 0073 const G4Step& aStep) override; 0074 // This is the method implementing Rayleigh scattering. 0075 0076 virtual G4PhysicsTable* GetPhysicsTable() const; 0077 // Returns the address of the physics table. 0078 0079 virtual void DumpPhysicsTable() const; 0080 // Prints the physics table. 0081 0082 virtual void PreparePhysicsTable(const G4ParticleDefinition&) override; 0083 virtual void Initialise(); 0084 0085 void SetVerboseLevel(G4int); 0086 0087 protected: 0088 G4PhysicsTable* thePhysicsTable; 0089 0090 private: 0091 G4OpRayleigh(const G4OpRayleigh& right) = delete; 0092 G4OpRayleigh& operator=(const G4OpRayleigh& right) = delete; 0093 0094 /// Calculates the mean free paths for a material as a function of 0095 /// photon energy 0096 G4PhysicsFreeVector* CalculateRayleighMeanFreePaths( 0097 const G4Material* material) const; 0098 0099 size_t idx_rslength = 0; 0100 }; 0101 0102 //////////////////// 0103 // Inline methods 0104 //////////////////// 0105 0106 inline G4bool G4OpRayleigh::IsApplicable( 0107 const G4ParticleDefinition& aParticleType) 0108 { 0109 return (&aParticleType == G4OpticalPhoton::OpticalPhoton()); 0110 } 0111 0112 inline void G4OpRayleigh::DumpPhysicsTable() const 0113 { 0114 for(size_t i = 0; i < thePhysicsTable->entries(); ++i) 0115 { 0116 ((G4PhysicsFreeVector*) (*thePhysicsTable)[i])->DumpValues(); 0117 } 0118 } 0119 0120 inline G4PhysicsTable* G4OpRayleigh::GetPhysicsTable() const 0121 { 0122 return thePhysicsTable; 0123 } 0124 0125 #endif /* G4OpRayleigh_h */
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |