#include <G4ElNucleusSFcs.hh>
|
| G4ElNucleusSFcs () |
|
virtual | ~G4ElNucleusSFcs () |
|
virtual void | CrossSectionDescription (std::ostream &) const |
|
virtual G4bool | IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) |
|
virtual G4double | GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) |
|
G4double | ThresholdEnergy () |
|
G4double | GetRatio (G4int Z, G4int A) |
|
| G4VCrossSectionDataSet (const G4String &nam="") |
|
virtual | ~G4VCrossSectionDataSet () |
|
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 | 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 G4Isotope * | SelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) |
|
virtual void | BuildPhysicsTable (const G4ParticleDefinition &) |
|
virtual void | DumpPhysicsTable (const G4ParticleDefinition &) |
|
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 G4String & | GetName () const |
|
void | SetName (const G4String &nam) |
|
G4VCrossSectionDataSet & | operator= (const G4VCrossSectionDataSet &right)=delete |
|
| G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete |
|
Definition at line 51 of file G4ElNucleusSFcs.hh.
◆ G4ElNucleusSFcs()
G4ElNucleusSFcs::G4ElNucleusSFcs |
( |
| ) |
|
Definition at line 42 of file G4ElNucleusSFcs.cc.
43{
45}
static const char * Default_Name()
G4VCrossSectionDataSet(const G4String &nam="")
◆ ~G4ElNucleusSFcs()
G4ElNucleusSFcs::~G4ElNucleusSFcs |
( |
| ) |
|
|
virtual |
Definition at line 47 of file G4ElNucleusSFcs.cc.
48{
49 if( fCHIPScs != nullptr ) delete fCHIPScs;
50}
◆ CrossSectionDescription()
void G4ElNucleusSFcs::CrossSectionDescription |
( |
std::ostream & | outFile | ) |
const |
|
virtual |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 53 of file G4ElNucleusSFcs.cc.
54{
55 outFile << "G4ElNucleusSFcs provides the inelastic\n"
56 << "cross section for e- and e+ interactions with nuclei according to the structure function approach."
57 << "all energies.\n";
58}
◆ Default_Name()
static const char * G4ElNucleusSFcs::Default_Name |
( |
| ) |
|
|
inlinestatic |
◆ GetIsoCrossSection()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 70 of file G4ElNucleusSFcs.cc.
72{
74
76
78
79 if( ratio > 0. ) xsc /= ratio;
80
81 return xsc;
82}
G4double GetRatio(G4int Z, G4int A)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat)
◆ GetRatio()
Definition at line 88 of file G4ElNucleusSFcs.cc.
89{
91
92 if ( Z == 1 &&
A == 1 )
return 1.51;
93 else if ( Z == 1 &&
A == 2 )
return 0.33;
94 else if ( Z == 1 &&
A == 3 )
return 0.27;
95 else if ( Z == 2 &&
A == 4 )
return 1.81;
96 else if ( Z == 6 &&
A == 12 )
return 2.26;
97 else if ( Z == 7 &&
A == 14 )
return 2.47;
98 else if ( Z == 8 &&
A == 16 )
return 2.61;
99 else if ( Z == 13 &&
A == 27 )
return 2.57;
100 else if ( Z == 14 &&
A == 28 )
return 2.49;
101 else if ( Z == 18 &&
A == 40 )
return 2.72;
102 else if ( Z == 22 &&
A == 48 )
return 2.71;
103 else if ( Z == 26 &&
A == 56 )
return 2.79;
104 else if ( Z == 29 &&
A == 64 )
return 2.78;
105 else if ( Z == 32 &&
A == 73 )
return 2.87;
106 else if ( Z == 42 &&
A == 96 )
return 3.02;
107 else if ( Z == 46 &&
A == 106 )
return 3.02;
108 else if ( Z == 47 &&
A == 108 )
return 2.99;
109 else if ( Z == 48 &&
A == 112 )
return 3.00;
110 else if ( Z == 74 &&
A == 184 )
return 3.44;
111 else if ( Z == 79 &&
A == 200 )
return 3.49;
112 else if ( Z == 82 &&
A == 207 )
return 3.48;
113 else if ( Z == 92 &&
A == 238 )
return 3.88;
114 else
115 {
116 G4int it(0), iMax(19);
118
119 for ( it = 0; it < iMax; ++it ) if ( zz <= fZZ[it] ) break;
120
121 if ( it == 0 ) return fRR[0];
122 else if( it == iMax ) return fRR[iMax-1];
123 else
124 {
129
130 if( x1 >= x2 ) return fRR[it];
131 else
132 {
134 ratio = y1 + ( zz - x1 )*angle;
135 }
136 }
137 }
138 return ratio;
139}
Referenced by GetIsoCrossSection().
◆ IsElementApplicable()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 60 of file G4ElNucleusSFcs.cc.
61{
62
65
66 return ( Z > 0 && Z < 120 && eKin > thTkin );
67}
G4double GetKineticEnergy() const
G4double ThresholdEnergy()
◆ ThresholdEnergy()
G4double G4ElNucleusSFcs::ThresholdEnergy |
( |
| ) |
|
The documentation for this class was generated from the following files: