|
||||
File indexing completed on 2025-01-18 09:58:38
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 // | G4LowEPComptonModel-- Geant4 Monash University | 0029 // | low energy Compton scattering model. | 0030 // | J. M. C. Brown, Monash University, Australia | 0031 // | ## Unpolarised photons only ## | 0032 // | | 0033 // | | 0034 // ********************************************************************* 0035 // | | 0036 // | The following is a Geant4 class to simulate the process of | 0037 // | bound electron Compton scattering. General code structure is | 0038 // | based on G4LowEnergyCompton.cc and G4LivermoreComptonModel.cc. | 0039 // | Algorithms for photon energy, and ejected Compton electron | 0040 // | direction taken from: | 0041 // | | 0042 // | J. M. C. Brown, M. R. Dimmock, J. E. Gillam and D. M. Paganin, | 0043 // | "A low energy bound atomic electron Compton scattering model | 0044 // | for Geant4", NIMB, Vol. 338, 77-88, 2014. | 0045 // | | 0046 // | The author acknowledges the work of the Geant4 collaboration | 0047 // | in developing the following algorithms that have been employed | 0048 // | or adapeted for the present software: | 0049 // | | 0050 // | # sampling of photon scattering angle, | 0051 // | # target element selection in composite materials, | 0052 // | # target shell selection in element, | 0053 // | # and sampling of bound electron momentum from Compton profiles. | 0054 // | | 0055 // ********************************************************************* 0056 // | | 0057 // | History: | 0058 // | -------- | 0059 // | | 0060 // | Nov. 2011 JMCB - First version | 0061 // | Feb. 2012 JMCB - Migration to Geant4 9.5 | 0062 // | Sep. 2012 JMCB - Final fixes for Geant4 9.6 | 0063 // | Feb. 2013 JMCB - Geant4 9.6 FPE fix for bug 1426 | 0064 // | Jan. 2015 JMCB - Migration to MT (Based on Livermore | 0065 // | implementation) | 0066 // | Feb. 2016 JMCB - Geant4 10.2 FPE fix for bug 1676 | 0067 // | | 0068 // ********************************************************************* 0069 0070 #ifndef G4LowEPComptonModel_h 0071 #define G4LowEPComptonModel_h 1 0072 0073 #include "G4VEmModel.hh" 0074 #include "G4PhysicsFreeVector.hh" 0075 0076 class G4ParticleChangeForGamma; 0077 class G4VAtomDeexcitation; 0078 class G4ShellData; 0079 class G4DopplerProfile; 0080 0081 class G4LowEPComptonModel : public G4VEmModel 0082 { 0083 0084 public: 0085 0086 explicit G4LowEPComptonModel(const G4ParticleDefinition* p = nullptr, 0087 const G4String& nam = "LowEPComptonModel"); 0088 virtual ~G4LowEPComptonModel(); 0089 0090 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override; 0091 void InitialiseLocal(const G4ParticleDefinition*, 0092 G4VEmModel* masterModel) override; 0093 0094 void InitialiseForElement(const G4ParticleDefinition*, G4int Z) override; 0095 0096 G4double ComputeCrossSectionPerAtom( const G4ParticleDefinition*, 0097 G4double kinEnergy, 0098 G4double Z, 0099 G4double A=0, 0100 G4double cut=0, 0101 G4double emax=DBL_MAX ) override; 0102 0103 void SampleSecondaries(std::vector<G4DynamicParticle*>*, 0104 const G4MaterialCutsCouple*, 0105 const G4DynamicParticle*, 0106 G4double tmin, 0107 G4double maxEnergy) override; 0108 0109 G4LowEPComptonModel & operator=(const G4LowEPComptonModel &right) = delete; 0110 G4LowEPComptonModel(const G4LowEPComptonModel&) = delete; 0111 0112 private: 0113 void ReadData(size_t Z, const char* path = 0); 0114 G4double ComputeScatteringFunction(G4double x, G4int Z); 0115 0116 G4ParticleChangeForGamma* fParticleChange; 0117 G4VAtomDeexcitation* fAtomDeexcitation; 0118 static G4ShellData* shellData; 0119 static G4DopplerProfile* profileData; 0120 0121 static const G4int maxZ = 99; 0122 static G4PhysicsFreeVector* data[100]; 0123 static const G4double ScatFuncFitParam[101][9]; 0124 0125 G4int verboseLevel; 0126 G4bool isInitialised; 0127 }; 0128 0129 //**************************************************************************** 0130 0131 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |