|
||||
File indexing completed on 2025-01-18 09:59:30
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 // author: Vladimir Grichine 0027 // 0028 // 14.10.12 V. Grichine, fFormFactor was added as the class member 0029 // 25.05.2011 first implementation 0030 // 0031 // X ray Rayleigh scattering model based on simplified form-factors 0032 // and angular distribution 0033 // 0034 0035 #ifndef G4XrayRayleighModel_h 0036 #define G4XrayRayleighModel_h 1 0037 0038 #include "G4VEmModel.hh" 0039 #include "G4ParticleChangeForGamma.hh" 0040 #include "G4Gamma.hh" 0041 0042 class G4XrayRayleighModel : public G4VEmModel 0043 { 0044 0045 public: 0046 0047 explicit G4XrayRayleighModel(const G4ParticleDefinition* p = nullptr, 0048 const G4String& nam = "XrayRayleigh"); 0049 0050 ~G4XrayRayleighModel() override; 0051 0052 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override; 0053 0054 G4double ComputeCrossSectionPerAtom( 0055 const G4ParticleDefinition*, 0056 G4double kinEnergy, 0057 G4double Z, 0058 G4double A=0, 0059 G4double cut=0, 0060 G4double emax=DBL_MAX) override; 0061 0062 void SampleSecondaries(std::vector<G4DynamicParticle*>*, 0063 const G4MaterialCutsCouple*, 0064 const G4DynamicParticle*, 0065 G4double tmin, 0066 G4double maxEnergy) override; 0067 0068 G4XrayRayleighModel & operator=(const G4XrayRayleighModel &right) = delete; 0069 G4XrayRayleighModel(const G4XrayRayleighModel&) = delete; 0070 0071 private: 0072 0073 G4ParticleChangeForGamma* fParticleChange; 0074 0075 G4double lowEnergyLimit; 0076 G4double highEnergyLimit; 0077 G4double fFormFactor; 0078 0079 G4int verboseLevel; 0080 G4bool isInitialised; 0081 0082 static const G4double fCofA; 0083 static const G4double fCofR; 0084 0085 }; 0086 0087 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 0088 0089 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |