Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronElectronElXsc Class Reference

#include <G4NeutronElectronElXsc.hh>

+ Inheritance diagram for G4NeutronElectronElXsc:

Public Member Functions

 G4NeutronElectronElXsc ()
 
 ~G4NeutronElectronElXsc ()
 
void Initialise ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *)
 
G4double GetRosenbluthXsc (const G4DynamicParticle *, G4int Z, const G4Material *)
 
G4double XscIntegrand (G4double x)
 
G4double GetElementNonRelXsc (const G4DynamicParticle *, G4int Z, const G4Material *)
 
G4double CalculateAm (G4double momentum)
 
G4double GetAm ()
 
void SetCutEnergy (G4double ec)
 
G4double GetCutEnergy ()
 
void SetBiasingFactor (G4double bf)
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
bool ForAllAtomsAndEnergies () const
 
void SetForAllAtomsAndEnergies (G4bool val)
 
const G4StringGetName () const
 
void SetName (const G4String &nam)
 
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete
 

Protected Attributes

G4double fM
 
G4double fM2
 
G4double fMv2
 
G4double fme
 
G4double fme2
 
G4double fee
 
G4double fee2
 
G4double fCofXsc
 
G4double fAm
 
G4int fEnergyBin
 
G4double fMinEnergy
 
G4double fMaxEnergy
 
G4double fCutEnergy
 
G4double fBiasingFactor
 
G4PhysicsLogVectorfEnergyXscVector
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 
G4String name
 

Static Protected Attributes

static const G4double fXscArray [200]
 

Detailed Description

Definition at line 45 of file G4NeutronElectronElXsc.hh.

Constructor & Destructor Documentation

◆ G4NeutronElectronElXsc()

G4NeutronElectronElXsc::G4NeutronElectronElXsc ( )

Definition at line 50 of file G4NeutronElectronElXsc.cc.

51 : G4VCrossSectionDataSet("NuElectronCcXsc")
52{
53 // neutron magneton squared
54
55 fM = neutron_mass_c2; // neutron mass
56 fM2 = fM*fM;
57 fme = electron_mass_c2;
58 fme2 = fme*fme;
59 fMv2 = 0.7056*GeV*GeV;
60 fee = fme;
61 fee2 = fee*fee;
62 fAm = 0.001;
63
64 fCofXsc = pi*fine_structure_const*fine_structure_const*hbarc*hbarc;
65 fCofXsc *= 3.6481; // neutron Fm^2(0)
66 fCofXsc /= fM*fM;
67
68 // G4cout<<"fCofXsc = "<<fCofXsc*GeV/cm2<<" cm2/GeV"<<G4endl;
69
70 // G4cout<<"hbarc = "<<hbarc/MeV/fermi<<" MeV*fermi"<<G4endl;
71
72 fCutEnergy = 0.; // default value
73
74 fEnergyBin = 200;
75 fMinEnergy = 1.*MeV;
76 fMaxEnergy = 10000.*GeV;
77
79
80 for( G4int iTkin = 0; iTkin < fEnergyBin; iTkin++) fEnergyXscVector->PutValue(iTkin, fXscArray[iTkin]*microbarn);
81
82 fBiasingFactor = 1.;
83
84 // Initialise();
85}
int G4int
Definition: G4Types.hh:85
static const G4double fXscArray[200]
G4PhysicsLogVector * fEnergyXscVector
void PutValue(const std::size_t index, const G4double value)
const G4double pi

◆ ~G4NeutronElectronElXsc()

G4NeutronElectronElXsc::~G4NeutronElectronElXsc ( )

Definition at line 87 of file G4NeutronElectronElXsc.cc.

88{
89 if( fEnergyXscVector )
90 {
91 delete fEnergyXscVector;
93 }
94}

Member Function Documentation

◆ CalculateAm()

G4double G4NeutronElectronElXsc::CalculateAm ( G4double  momentum)
inline

Definition at line 97 of file G4NeutronElectronElXsc.hh.

98{
99 G4double k = momentum/CLHEP::hbarc;
100 G4double ch = 1.13;
101 G4double zn = 1.77*k*CLHEP::Bohr_radius;
102 G4double zn2 = zn*zn;
103 fAm = ch/zn2;
104
105 return fAm;
106}
double G4double
Definition: G4Types.hh:83

Referenced by GetElementNonRelXsc(), and GetRosenbluthXsc().

◆ GetAm()

G4double G4NeutronElectronElXsc::GetAm ( )
inline

Definition at line 72 of file G4NeutronElectronElXsc.hh.

72{return fAm;};

◆ GetCutEnergy()

G4double G4NeutronElectronElXsc::GetCutEnergy ( )
inline

Definition at line 75 of file G4NeutronElectronElXsc.hh.

75{return fCutEnergy;};

◆ GetElementCrossSection()

G4double G4NeutronElectronElXsc::GetElementCrossSection ( const G4DynamicParticle aPart,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 141 of file G4NeutronElectronElXsc.cc.

144{
145 G4double result = 0., Tkin;
146
147 Tkin = aPart->GetKineticEnergy();
148
149 result = fEnergyXscVector->Value(Tkin);
150
151 result *= ZZ; // incoherent sum over all element electrons
152
153 result *= fBiasingFactor;
154
155 return result;
156}
G4double GetKineticEnergy() const
G4double Value(const G4double energy, std::size_t &lastidx) const

◆ GetElementNonRelXsc()

G4double G4NeutronElectronElXsc::GetElementNonRelXsc ( const G4DynamicParticle aPart,
G4int  Z,
const G4Material  
)
inline

Definition at line 112 of file G4NeutronElectronElXsc.hh.

115{
116 G4double result(0.), te(0.), momentum(0.);
117
118 te = aPart->GetKineticEnergy()*fme/fM;
119 momentum = std::sqrt( te*(te + 2.*fme) );
120 fAm = CalculateAm(momentum);
121
122 result = 1. + std::log(1. +1./fAm);
123 result *= fCofXsc; //*energy;
124 result *= ZZ; // incoherent sum over all element electrons
125
126 return result;
127}
G4double CalculateAm(G4double momentum)

◆ GetRosenbluthXsc()

G4double G4NeutronElectronElXsc::GetRosenbluthXsc ( const G4DynamicParticle aPart,
G4int  Z,
const G4Material  
)

Definition at line 162 of file G4NeutronElectronElXsc.cc.

165{
166 G4double result = 0., momentum;
167
168 fee = aPart->GetTotalEnergy()*fme/fM;
169 fee2 = fee*fee;
170 momentum = sqrt( fee2 - fme2 );
171 fAm = CalculateAm(momentum);
172
173 // G4Integrator<G4NeutronElectronElXsc, G4double(G4NeutronElectronElXsc::*)(G4double)> integral;
174
175 // result = integral.Legendre96( this, &G4NeutronElectronElXsc::XscIntegrand, 0., 1. );
176
177 result *= fCofXsc;
178
179 result *= ZZ; // incoherent sum over all element electrons
180
181 return result;
182}
G4double GetTotalEnergy() const

Referenced by Initialise().

◆ Initialise()

void G4NeutronElectronElXsc::Initialise ( )

Definition at line 116 of file G4NeutronElectronElXsc.cc.

117{
118 G4int iTkin;
119 G4double Tkin, rosxsc, xsc, delta, err=1.e-5;
120 const G4ThreeVector mDir = G4ThreeVector(0.,0.,1.);
123
125
126 for( iTkin = 0; iTkin < fEnergyBin; iTkin++)
127 {
128 Tkin = fEnergyXscVector->GetLowEdgeEnergy(iTkin);
129 dP = G4DynamicParticle(pD, mDir, Tkin);
130 rosxsc = GetRosenbluthXsc(&dP, 1, mat);
131 fEnergyXscVector->PutValue(iTkin, rosxsc); // xscV.PutValue(evt, rosxsc); //
132 xsc= fEnergyXscVector->Value(Tkin); // xsc= xscV.GetLowEdgeEnergy(evt); //
133 delta = 0.5*std::abs( (rosxsc-xsc) )/(rosxsc+xsc);
134 if(delta > err) G4cout<<Tkin/GeV<<" GeV, rosxsc = "<<rosxsc/microbarn<<"umb, v-xsc = "<<xsc<<" umb"<<G4endl;
135 }
136 return;
137}
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetRosenbluthXsc(const G4DynamicParticle *, G4int Z, const G4Material *)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
G4double GetLowEdgeEnergy(const std::size_t index) const

◆ IsElementApplicable()

G4bool G4NeutronElectronElXsc::IsElementApplicable ( const G4DynamicParticle aPart,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 101 of file G4NeutronElectronElXsc.cc.

102{
103 G4bool result = false;
104 G4String pName = aPart->GetDefinition()->GetParticleName();
105 G4double Tkin = aPart->GetKineticEnergy();
106
107 if( pName == "neutron" &&
108 Tkin >= fMinEnergy &&
109 Tkin <= fMaxEnergy ) result = true;
110
111 return result;
112}
bool G4bool
Definition: G4Types.hh:86
G4ParticleDefinition * GetDefinition() const
const G4String & GetParticleName() const

◆ SetBiasingFactor()

void G4NeutronElectronElXsc::SetBiasingFactor ( G4double  bf)
inline

Definition at line 77 of file G4NeutronElectronElXsc.hh.

77{fBiasingFactor=bf;};

◆ SetCutEnergy()

void G4NeutronElectronElXsc::SetCutEnergy ( G4double  ec)
inline

Definition at line 74 of file G4NeutronElectronElXsc.hh.

74{fCutEnergy=ec;};

◆ XscIntegrand()

G4double G4NeutronElectronElXsc::XscIntegrand ( G4double  x)

Definition at line 190 of file G4NeutronElectronElXsc.cc.

191{
192 G4double result = 1., q2, znq2, znf, znf2, znf4;
193
194 znq2 = 1. + 2.*fee*x/fM;
195
196 q2 = 4.*fee2*x/znq2;
197
198 znf = 1 + q2/fMv2;
199 znf2 = znf*znf;
200 znf4 = znf2*znf2;
201
202 result /= ( x + fAm )*znq2*znq2*znf4;
203
204 result *= ( 1 - x )/( 1 + q2/4./fM2 ) + 2.*x;
205
206 return result;
207}

Member Data Documentation

◆ fAm

G4double G4NeutronElectronElXsc::fAm
protected

◆ fBiasingFactor

G4double G4NeutronElectronElXsc::fBiasingFactor
protected

◆ fCofXsc

G4double G4NeutronElectronElXsc::fCofXsc
protected

◆ fCutEnergy

G4double G4NeutronElectronElXsc::fCutEnergy
protected

Definition at line 84 of file G4NeutronElectronElXsc.hh.

Referenced by G4NeutronElectronElXsc(), GetCutEnergy(), and SetCutEnergy().

◆ fee

G4double G4NeutronElectronElXsc::fee
protected

◆ fee2

G4double G4NeutronElectronElXsc::fee2
protected

◆ fEnergyBin

G4int G4NeutronElectronElXsc::fEnergyBin
protected

Definition at line 83 of file G4NeutronElectronElXsc.hh.

Referenced by G4NeutronElectronElXsc(), and Initialise().

◆ fEnergyXscVector

G4PhysicsLogVector* G4NeutronElectronElXsc::fEnergyXscVector
protected

◆ fM

G4double G4NeutronElectronElXsc::fM
protected

◆ fM2

G4double G4NeutronElectronElXsc::fM2
protected

Definition at line 80 of file G4NeutronElectronElXsc.hh.

Referenced by G4NeutronElectronElXsc(), and XscIntegrand().

◆ fMaxEnergy

G4double G4NeutronElectronElXsc::fMaxEnergy
protected

Definition at line 84 of file G4NeutronElectronElXsc.hh.

Referenced by G4NeutronElectronElXsc(), and IsElementApplicable().

◆ fme

G4double G4NeutronElectronElXsc::fme
protected

◆ fme2

G4double G4NeutronElectronElXsc::fme2
protected

Definition at line 80 of file G4NeutronElectronElXsc.hh.

Referenced by G4NeutronElectronElXsc(), and GetRosenbluthXsc().

◆ fMinEnergy

G4double G4NeutronElectronElXsc::fMinEnergy
protected

Definition at line 84 of file G4NeutronElectronElXsc.hh.

Referenced by G4NeutronElectronElXsc(), and IsElementApplicable().

◆ fMv2

G4double G4NeutronElectronElXsc::fMv2
protected

Definition at line 80 of file G4NeutronElectronElXsc.hh.

Referenced by G4NeutronElectronElXsc(), and XscIntegrand().

◆ fXscArray

const G4double G4NeutronElectronElXsc::fXscArray
staticprotected
Initial value:
= {
1.52681, 1.54903, 1.57123, 1.59341, 1.61556, 1.63769, 1.6598, 1.68189, 1.70396,
1.72601, 1.74805, 1.77007, 1.79208, 1.81407, 1.83605, 1.85801, 1.87997, 1.90192,
1.92385, 1.94578, 1.96771, 1.98962, 2.01154, 2.03345, 2.05535, 2.07725, 2.09915,
2.12105, 2.14295, 2.16485, 2.18675, 2.20865, 2.23055, 2.25244, 2.27433, 2.29621,
2.31807, 2.33992, 2.36173, 2.38351, 2.40524, 2.42691, 2.4485, 2.47, 2.49138,
2.51262, 2.53369, 2.55457, 2.57524, 2.59565, 2.61577, 2.63559, 2.65505, 2.67414,
2.69281, 2.71104, 2.72881, 2.74607, 2.76282, 2.77903, 2.79467, 2.80974, 2.82422,
2.83811, 2.85139, 2.86408, 2.87616, 2.88764, 2.89854, 2.90885, 2.91859, 2.92777,
2.93641, 2.94453, 2.95213, 2.95924, 2.96588, 2.97207, 2.97782, 2.98316, 2.98811,
2.99268, 2.9969, 3.00078, 3.00435, 3.00761, 3.01059, 3.01331, 3.01578, 3.01801,
3.02003, 3.02185, 3.02347, 3.02491, 3.02619, 3.02732, 3.0283, 3.02915, 3.02988,
3.03049, 3.03099, 3.03139, 3.03169, 3.03191, 3.03203, 3.03208, 3.03205, 3.03195,
3.03177, 3.03152, 3.0312, 3.03081, 3.03034, 3.0298, 3.02919, 3.02849, 3.02771,
3.02684, 3.02588, 3.02482, 3.02365, 3.02237, 3.02097, 3.01943, 3.01775, 3.0159,
3.01389, 3.01169, 3.00929, 3.00666, 3.00379, 3.00065, 2.99722, 2.99347, 2.98936,
2.98487, 2.97996, 2.97459, 2.9687, 2.96226, 2.9552, 2.94748, 2.93903, 2.92977,
2.91965, 2.90858, 2.89649, 2.88329, 2.86889, 2.85321, 2.83615, 2.81764, 2.7976,
2.77594, 2.7526, 2.72754, 2.70071, 2.67209, 2.64171, 2.60957, 2.57575, 2.54031,
2.50336, 2.46504, 2.42548, 2.38484, 2.34328, 2.30099, 2.2581, 2.21478, 2.17115,
2.12735, 2.08345, 2.03954, 1.99569, 1.95191, 1.90825, 1.86471, 1.82129, 1.77799,
1.7348, 1.69171, 1.64869, 1.60575, 1.56286, 1.52, 1.47718, 1.43437, 1.39157,
1.34877, 1.30596, 1.26314, 1.22031, 1.17746, 1.13459, 1.0917, 1.04879, 1.00585,
0.962892, 0.919908 }

Definition at line 88 of file G4NeutronElectronElXsc.hh.

Referenced by G4NeutronElectronElXsc().


The documentation for this class was generated from the following files: