Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GammaNuclearXS Class Referencefinal

#include <G4GammaNuclearXS.hh>

+ Inheritance diagram for G4GammaNuclearXS:

Public Member Functions

 G4GammaNuclearXS ()
 
 ~G4GammaNuclearXS () final
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) final
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *mat) final
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat) final
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
 
const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final
 
void BuildPhysicsTable (const G4ParticleDefinition &) final
 
void CrossSectionDescription (std::ostream &) const final
 
G4GammaNuclearXSoperator= (const G4GammaNuclearXS &right)=delete
 
 G4GammaNuclearXS (const G4GammaNuclearXS &)=delete
 
- 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 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 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 G4int GetVerboseLevel () 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
 

Static Public Member Functions

static const char * Default_Name ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 57 of file G4GammaNuclearXS.hh.

Constructor & Destructor Documentation

◆ G4GammaNuclearXS() [1/2]

G4GammaNuclearXS::G4GammaNuclearXS ( )
explicit

Definition at line 67 of file G4GammaNuclearXS.cc.

69 ggXsection(nullptr),
70 gamma(G4Gamma::Gamma()),
71 isMaster(false)
72{
73 // verboseLevel = 0;
74 if(verboseLevel > 0){
75 G4cout << "G4GammaNuclearXS::G4GammaNuclearXS Initialise for Z < "
76 << MAXZEL << G4endl;
77 }
79 if(ggXsection == nullptr) ggXsection = new G4PhotoNuclearCrossSection();
81}
const G4int MAXZEL
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
static G4CrossSectionDataSetRegistry * Instance()
static const char * Default_Name()
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetForAllAtomsAndEnergies(G4bool val)

◆ ~G4GammaNuclearXS()

G4GammaNuclearXS::~G4GammaNuclearXS ( )
final

Definition at line 83 of file G4GammaNuclearXS.cc.

84{
85 if(isMaster) {
86 for(G4int i=0; i<MAXZEL; ++i) {
87 delete data[i];
88 data[i] = nullptr;
89 }
90 }
91}
int G4int
Definition: G4Types.hh:85

◆ G4GammaNuclearXS() [2/2]

G4GammaNuclearXS::G4GammaNuclearXS ( const G4GammaNuclearXS )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4GammaNuclearXS::BuildPhysicsTable ( const G4ParticleDefinition p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 181 of file G4GammaNuclearXS.cc.

182{
183 if(verboseLevel > 0){
184 G4cout << "G4GammaNuclearXS::BuildPhysicsTable for "
185 << p.GetParticleName() << G4endl;
186 }
187 if(p.GetParticleName() != "gamma") {
189 ed << p.GetParticleName() << " is a wrong particle type -"
190 << " only gamma is allowed";
191 G4Exception("G4GammaNuclearXS::BuildPhysicsTable(..)","had012",
192 FatalException, ed, "");
193 return;
194 }
195 if(0. == coeff[0]) {
196#ifdef G4MULTITHREADED
197 G4MUTEXLOCK(&gNuclearXSMutex);
198 if(0. == coeff[0]) {
199#endif
200 coeff[0] = 1.0;
201 isMaster = true;
202 FindDirectoryPath();
203#ifdef G4MULTITHREADED
204 }
205 G4MUTEXUNLOCK(&gNuclearXSMutex);
206#endif
207 }
208
209 // it is possible re-initialisation for the second run
210 if(isMaster) {
211
212 // Access to elements
213 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
214 size_t numOfCouples = theCoupleTable->GetTableSize();
215 for(size_t j=0; j<numOfCouples; ++j) {
216 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
217 auto elmVec = mat->GetElementVector();
218 size_t numOfElem = mat->GetNumberOfElements();
219 for (size_t ie = 0; ie < numOfElem; ++ie) {
220 G4int Z = std::max(1,std::min(((*elmVec)[ie])->GetZasInt(), MAXZEL-1));
221 if(data[Z] == nullptr) { Initialise(Z); }
222 }
223 }
224 }
225}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
const G4String & GetParticleName() const
static G4ProductionCutsTable * GetProductionCutsTable()

◆ CrossSectionDescription()

void G4GammaNuclearXS::CrossSectionDescription ( std::ostream &  outFile) const
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 93 of file G4GammaNuclearXS.cc.

94{
95 outFile << "G4GammaNuclearXS calculates the gamma nuclear\n"
96 << "cross section on nuclei using data from the high precision\n"
97 << "LEND gamma database. The data are simplified and smoothed over\n"
98 << "the resonance region in order to reduce CPU time.\n"
99 << "For high energies Glauber-Gribiv cross section is used.\n";
100}

◆ Default_Name()

static const char * G4GammaNuclearXS::Default_Name ( )
inlinestatic

Definition at line 65 of file G4GammaNuclearXS.hh.

65{return "G4GammaNuclearXS";}

◆ GetElementCrossSection()

G4double G4GammaNuclearXS::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material mat 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 117 of file G4GammaNuclearXS.cc.

119{
120 G4double xs = 0.0;
121 G4double ekin = aParticle->GetKineticEnergy();
122
123 G4int Z = (ZZ >= MAXZEL) ? MAXZEL - 1 : ZZ;
124
125 auto pv = GetPhysicsVector(Z);
126 if(pv == nullptr) { return xs; }
127 // G4cout << "G4GammaNuclearXS::GetCrossSection e= " << ekin
128 // << " Z= " << Z << G4endl;
129
130 if(ekin <= pv->GetMaxEnergy()) {
131 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
132 } else {
133 xs = coeff[Z]*ggXsection->GetElementCrossSection(aParticle, Z, mat);
134 }
135
136#ifdef G4VERBOSE
137 if(verboseLevel > 1) {
138 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
139 << ", nElmXS(b)= " << xs/CLHEP::barn
140 << G4endl;
141 }
142#endif
143 return xs;
144}
double G4double
Definition: G4Types.hh:83
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)

Referenced by GetIsoCrossSection().

◆ GetIsoCrossSection()

G4double G4GammaNuclearXS::GetIsoCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
G4int  A,
const G4Isotope iso,
const G4Element elm,
const G4Material mat 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 146 of file G4GammaNuclearXS.cc.

151{
152 return GetElementCrossSection(aParticle, Z, mat) * A/aeff[Z];
153}
double A(double temperature)
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat) final

◆ IsElementApplicable()

G4bool G4GammaNuclearXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material  
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 103 of file G4GammaNuclearXS.cc.

105{
106 return true;
107}

◆ IsIsoApplicable()

G4bool G4GammaNuclearXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element ,
const G4Material mat 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 109 of file G4GammaNuclearXS.cc.

112{
113 return true;
114}

◆ operator=()

G4GammaNuclearXS & G4GammaNuclearXS::operator= ( const G4GammaNuclearXS right)
delete

◆ SelectIsotope()

const G4Isotope * G4GammaNuclearXS::SelectIsotope ( const G4Element anElement,
G4double  kinEnergy,
G4double  logE 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 155 of file G4GammaNuclearXS.cc.

157{
158 size_t nIso = anElement->GetNumberOfIsotopes();
159 const G4Isotope* iso = anElement->GetIsotope(0);
160
161 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
162 if(1 == nIso) { return iso; }
163
164 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
166 G4double sum = 0.0;
167 size_t j;
168
169 // isotope wise cross section not used
170 for (j=0; j<nIso; ++j) {
171 sum += abundVector[j];
172 if(q <= sum) {
173 iso = anElement->GetIsotope(j);
174 break;
175 }
176 }
177 return iso;
178}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158

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