File indexing completed on 2025-01-18 09:58:44
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033 #ifndef G4NeutronElectronElXsc_h
0034 #define G4NeutronElectronElXsc_h
0035
0036
0037 #include "globals.hh"
0038 #include "G4VCrossSectionDataSet.hh"
0039 #include "G4DynamicParticle.hh"
0040
0041
0042 class G4PhysicsLogVector;
0043 class G4PhysicsTable;
0044
0045 class G4NeutronElectronElXsc : public G4VCrossSectionDataSet
0046 {
0047 public:
0048
0049 G4NeutronElectronElXsc();
0050 ~G4NeutronElectronElXsc();
0051
0052 void Initialise();
0053
0054 virtual
0055 G4bool IsElementApplicable(const G4DynamicParticle*, G4int Z, const G4Material*);
0056
0057
0058 virtual
0059 G4double GetElementCrossSection(const G4DynamicParticle*,
0060 G4int Z, const G4Material*);
0061
0062 G4double GetRosenbluthXsc(const G4DynamicParticle*,
0063 G4int Z, const G4Material*);
0064
0065 G4double XscIntegrand(G4double x);
0066
0067 G4double GetElementNonRelXsc(const G4DynamicParticle*,
0068 G4int Z, const G4Material*);
0069
0070 G4double CalculateAm( G4double momentum);
0071
0072 inline G4double GetAm(){return fAm;};
0073
0074 void SetCutEnergy(G4double ec){fCutEnergy=ec;};
0075 G4double GetCutEnergy(){return fCutEnergy;};
0076
0077 void SetBiasingFactor(G4double bf){fBiasingFactor=bf;};
0078
0079 protected:
0080 G4double fM, fM2, fMv2, fme, fme2, fee, fee2;
0081 G4double fCofXsc;
0082 G4double fAm;
0083 G4int fEnergyBin;
0084 G4double fMinEnergy, fMaxEnergy, fCutEnergy;
0085 G4double fBiasingFactor;
0086
0087 G4PhysicsLogVector* fEnergyXscVector;
0088 static const G4double fXscArray[200];
0089 };
0090
0091
0092
0093
0094
0095
0096
0097 inline G4double G4NeutronElectronElXsc::CalculateAm( G4double momentum)
0098 {
0099 G4double k = momentum/CLHEP::hbarc;
0100 G4double ch = 1.13;
0101 G4double zn = 1.77*k*CLHEP::Bohr_radius;
0102 G4double zn2 = zn*zn;
0103 fAm = ch/zn2;
0104
0105 return fAm;
0106 }
0107
0108
0109
0110
0111
0112 inline G4double G4NeutronElectronElXsc::
0113 GetElementNonRelXsc(const G4DynamicParticle* aPart, G4int ZZ,
0114 const G4Material*)
0115 {
0116 G4double result(0.), te(0.), momentum(0.);
0117
0118 te = aPart->GetKineticEnergy()*fme/fM;
0119 momentum = std::sqrt( te*(te + 2.*fme) );
0120 fAm = CalculateAm(momentum);
0121
0122 result = 1. + std::log(1. +1./fAm);
0123 result *= fCofXsc;
0124 result *= ZZ;
0125
0126 return result;
0127 }
0128
0129 #endif