Geant4-11
Public Member Functions | Static Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4GammaNuclearXS Class Referencefinal

#include <G4GammaNuclearXS.hh>

Inheritance diagram for G4GammaNuclearXS:
G4VCrossSectionDataSet

Public Member Functions

void BuildPhysicsTable (const G4ParticleDefinition &) final
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
void CrossSectionDescription (std::ostream &) const final
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
G4double ElementCrossSection (G4double ekin, G4int Z)
 
bool ForAllAtomsAndEnergies () const
 
 G4GammaNuclearXS ()
 
 G4GammaNuclearXS (const G4GammaNuclearXS &)=delete
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
 
G4double GetMaxKinEnergy () const
 
G4double GetMinKinEnergy () const
 
const G4StringGetName () const
 
virtual G4int GetVerboseLevel () const
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) final
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *mat) final
 
G4double IsoCrossSection (G4double ekin, G4int Z, G4int A)
 
G4GammaNuclearXSoperator= (const G4GammaNuclearXS &right)=delete
 
const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final
 
void SetForAllAtomsAndEnergies (G4bool val)
 
void SetMaxKinEnergy (G4double value)
 
void SetMinKinEnergy (G4double value)
 
void SetName (const G4String &nam)
 
virtual void SetVerboseLevel (G4int value)
 
 ~G4GammaNuclearXS () final
 

Static Public Member Functions

static const char * Default_Name ()
 

Protected Attributes

G4String name
 
G4int verboseLevel
 

Private Member Functions

const G4StringFindDirectoryPath ()
 
G4PhysicsVectorGetPhysicsVector (G4int Z)
 
void Initialise (G4int Z)
 
void InitialiseOnFly (G4int Z)
 
G4PhysicsVectorRetrieveVector (std::ostringstream &in, G4bool isElement, G4int Z)
 

Private Attributes

const G4int freeVectorException [11] = {4, 6, 7, 8, 27, 39, 45, 65, 67, 69, 73}
 
const G4ParticleDefinitiongamma
 
G4VCrossSectionDataSetggXsection = nullptr
 
G4bool isForAllAtomsAndEnergies
 
G4bool isMaster = false
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
G4CrossSectionDataSetRegistryregistry
 
std::vector< G4doubletemp
 

Static Private Attributes

static G4double coeff [3][3]
 
static G4ElementDatadata = nullptr
 
static G4String gDataDirectory = ""
 
static const G4int MAXZGAMMAXS = 95
 
static const G4int rTransitionBound = 150.*CLHEP::MeV
 
static G4double xs150 [MAXZGAMMAXS] = {0.0}
 

Detailed Description

Definition at line 63 of file G4GammaNuclearXS.hh.

Constructor & Destructor Documentation

◆ G4GammaNuclearXS() [1/2]

G4GammaNuclearXS::G4GammaNuclearXS ( )
explicit

Definition at line 67 of file G4GammaNuclearXS.cc.

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

References G4cout, G4endl, G4CrossSectionDataSetRegistry::GetCrossSectionDataSet(), ggXsection, G4CrossSectionDataSetRegistry::Instance(), MAXZGAMMAXS, G4VCrossSectionDataSet::SetForAllAtomsAndEnergies(), and G4VCrossSectionDataSet::verboseLevel.

◆ ~G4GammaNuclearXS()

G4GammaNuclearXS::~G4GammaNuclearXS ( )
final

Definition at line 81 of file G4GammaNuclearXS.cc.

82{
83 if(isMaster) {
84 delete data;
85 data = nullptr;
86 }
87}
static G4ElementData * data

References data, and isMaster.

◆ 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 282 of file G4GammaNuclearXS.cc.

283{
284 if(verboseLevel > 0){
285 G4cout << "G4GammaNuclearXS::BuildPhysicsTable for "
286 << p.GetParticleName() << G4endl;
287 }
288 if(p.GetParticleName() != "gamma") {
290 ed << p.GetParticleName() << " is a wrong particle type -"
291 << " only gamma is allowed";
292 G4Exception("G4GammaNuclearXS::BuildPhysicsTable(..)","had012",
293 FatalException, ed, "");
294 return;
295 }
296
297 if(nullptr == data) {
298#ifdef G4MULTITHREADED
299 G4MUTEXLOCK(&gNuclearXSMutex);
300 if(nullptr == data) {
301#endif
302 isMaster = true;
303 data = new G4ElementData();
304 data->SetName("PhotoNuclear");
306#ifdef G4MULTITHREADED
307 }
308 G4MUTEXUNLOCK(&gNuclearXSMutex);
309#endif
310 }
311
312
313 // it is possible re-initialisation for the second run
314 // Upload data for elements used in geometry
316 if(isMaster) {
317 for ( auto & elm : *table ) {
318 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZGAMMAXS-1) );
319 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
320 }
321 }
322
323 // prepare isotope selection
324 size_t nIso = temp.size();
325 for ( auto & elm : *table ) {
326 size_t n = elm->GetNumberOfIsotopes();
327 if(n > nIso) { nIso = n; }
328 }
329 temp.resize(nIso, 0.0);
330}
std::vector< G4Element * > G4ElementTable
@ 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
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
G4PhysicsVector * GetElementData(G4int Z)
void SetName(const G4String &nam)
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
const G4String & FindDirectoryPath()
void Initialise(G4int Z)
std::vector< G4double > temp
const G4String & GetParticleName() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References data, FatalException, FindDirectoryPath(), G4cout, G4endl, G4Exception(), G4MUTEXLOCK, G4MUTEXUNLOCK, G4ElementData::GetElementData(), G4Element::GetElementTable(), G4ParticleDefinition::GetParticleName(), Initialise(), isMaster, G4INCL::Math::max(), MAXZGAMMAXS, G4INCL::Math::min(), CLHEP::detail::n, G4ElementData::SetName(), temp, G4VCrossSectionDataSet::verboseLevel, and Z.

◆ ComputeCrossSection()

G4double G4VCrossSectionDataSet::ComputeCrossSection ( const G4DynamicParticle part,
const G4Element elm,
const G4Material mat = nullptr 
)
inherited

Definition at line 81 of file G4VCrossSectionDataSet.cc.

84{
85 G4int Z = elm->GetZasInt();
86
87 if (IsElementApplicable(part, Z, mat)) {
88 return GetElementCrossSection(part, Z, mat);
89 }
90
91 // isotope-wise cross section making sum over available
92 // isotope cross sections, which may be incomplete, so
93 // the result is corrected
94 size_t nIso = elm->GetNumberOfIsotopes();
95 G4double fact = 0.0;
96 G4double xsec = 0.0;
97
98 // user-defined isotope abundances
99 const G4IsotopeVector* isoVector = elm->GetIsotopeVector();
100 const G4double* abundVector = elm->GetRelativeAbundanceVector();
101
102 for (size_t j=0; j<nIso; ++j) {
103 const G4Isotope* iso = (*isoVector)[j];
104 G4int A = iso->GetN();
105 if(abundVector[j] > 0.0 && IsIsoApplicable(part, Z, A, elm, mat)) {
106 fact += abundVector[j];
107 xsec += abundVector[j]*GetIsoCrossSection(part, Z, A, iso, elm, mat);
108 }
109 }
110 return (fact > 0.0) ? xsec/fact : 0.0;
111}
std::vector< G4Isotope * > G4IsotopeVector
double G4double
Definition: G4Types.hh:83
const G4double A[17]
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetZasInt() const
Definition: G4Element.hh:132
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:163
G4int GetN() const
Definition: G4Isotope.hh:93
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 G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)

References A, G4VCrossSectionDataSet::GetElementCrossSection(), G4VCrossSectionDataSet::GetIsoCrossSection(), G4Element::GetIsotopeVector(), G4Isotope::GetN(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), G4Element::GetZasInt(), G4VCrossSectionDataSet::IsElementApplicable(), G4VCrossSectionDataSet::IsIsoApplicable(), and Z.

Referenced by G4VCrossSectionDataSet::GetCrossSection().

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 89 of file G4GammaNuclearXS.cc.

90{
91 outFile << "G4GammaNuclearXS calculates the gamma nuclear\n"
92 << "cross-section for GDR energy region on nuclei using data from the high precision\n"
93 << "IAEA photonuclear database (2019). Then liniear connection\n"
94 <<"implemented with previous CHIPS photonuclear model\n";
95}

◆ Default_Name()

static const char * G4GammaNuclearXS::Default_Name ( )
inlinestatic

Definition at line 71 of file G4GammaNuclearXS.hh.

71{return "GammaNuclearXS";}

◆ DumpPhysicsTable()

void G4VCrossSectionDataSet::DumpPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ ElementCrossSection()

G4double G4GammaNuclearXS::ElementCrossSection ( G4double  ekin,
G4int  Z 
)

Definition at line 145 of file G4GammaNuclearXS.cc.

146{
147 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), ekin);
148 return GetElementCrossSection(&theGamma, ZZ);
149}
CLHEP::Hep3Vector G4ThreeVector
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final

References gamma, and GetElementCrossSection().

◆ FindDirectoryPath()

const G4String & G4GammaNuclearXS::FindDirectoryPath ( )
private

Definition at line 332 of file G4GammaNuclearXS.cc.

333{
334 // check environment variable
335 // build the complete string identifying the file with the data set
336 if(gDataDirectory.empty()) {
337 char* path = std::getenv("G4PARTICLEXSDATA");
338 if (nullptr != path) {
339 std::ostringstream ost;
340 ost << path << "/gamma/inel";
341 gDataDirectory = ost.str();
342 } else {
343 G4Exception("G4GammaNuclearXS::Initialise(..)","had013",
345 "Environment variable G4PARTICLEXSDATA is not defined");
346 }
347 }
348 return gDataDirectory;
349}
static G4String gDataDirectory

References FatalException, G4Exception(), and gDataDirectory.

Referenced by BuildPhysicsTable(), and Initialise().

◆ ForAllAtomsAndEnergies()

bool G4VCrossSectionDataSet::ForAllAtomsAndEnergies ( ) const
inlineinherited

◆ GetCrossSection()

G4double G4VCrossSectionDataSet::GetCrossSection ( const G4DynamicParticle dp,
const G4Element elm,
const G4Material mat = nullptr 
)
inlineinherited

Definition at line 187 of file G4VCrossSectionDataSet.hh.

190{
191 return ComputeCrossSection(dp, elm, mat);
192}
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)

References G4VCrossSectionDataSet::ComputeCrossSection().

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 112 of file G4GammaNuclearXS.cc.

114{
115 const G4int Z = (ZZ >= MAXZGAMMAXS) ? MAXZGAMMAXS - 1 : ZZ;
116 auto pv = GetPhysicsVector(Z);
117
118 if(pv == nullptr) {
119 return ggXsection->GetElementCrossSection(aParticle, Z, mat);
120 }
121 const G4double emax = pv->GetMaxEnergy();
122 const G4double ekin = aParticle->GetKineticEnergy();
123 G4double xs = 0.0;
124 if(ekin <= emax) {
125 xs = pv->Value(ekin);
126 } else if(ekin >= rTransitionBound){
127 xs = ggXsection->GetElementCrossSection(aParticle, Z, mat);
128 } else {
129 const G4double rxs = xs150[Z];
130 const G4double lxs = pv->Value(emax);
131 xs = lxs + (ekin - emax)*(rxs - lxs)/(rTransitionBound-emax);
132 }
133
134#ifdef G4VERBOSE
135 if(verboseLevel > 1) {
136 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
137 << ", nElmXS(b)= " << xs/CLHEP::barn
138 << G4endl;
139 }
140#endif
141 return xs;
142}
static const G4double emax
G4double GetKineticEnergy() const
static G4double xs150[MAXZGAMMAXS]
static const G4int rTransitionBound
G4PhysicsVector * GetPhysicsVector(G4int Z)
static constexpr double barn
Definition: SystemOfUnits.h:86
static constexpr double MeV

References CLHEP::barn, emax, G4cout, G4endl, G4VCrossSectionDataSet::GetElementCrossSection(), G4DynamicParticle::GetKineticEnergy(), GetPhysicsVector(), ggXsection, MAXZGAMMAXS, CLHEP::MeV, rTransitionBound, G4VCrossSectionDataSet::verboseLevel, xs150, and Z.

Referenced by ElementCrossSection().

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 158 of file G4GammaNuclearXS.cc.

163{
164 const G4int Z = (ZZ >= MAXZGAMMAXS) ? MAXZGAMMAXS - 1 : ZZ;
165 /*
166 G4cout << "IsoCrossSection Z= " << Z << " A= " << A
167 << " Amin= " << amin[Z] << " Amax= " << amax[Z]
168 << " E(MeV)= " << ekin << G4endl;
169 */
170 auto pv = GetPhysicsVector(Z);
171 if(pv == nullptr) {
172 return ggXsection->GetIsoCrossSection(aParticle, Z, A);
173 }
174 const G4double ekin = aParticle->GetKineticEnergy();
175 const G4double emax = pv->GetMaxEnergy();
176 G4double xs = 0.0;
177
178 // compute isotope cross section if applicable
179 if(amin[Z] < amax[Z] && A >= amin[Z] && A <= amax[Z] &&
180 ekin < rTransitionBound) {
181 auto pviso = data->GetComponentDataByIndex(Z, A - amin[Z]);
182 // isotope file exists
183 if(nullptr != pviso) {
184 const G4double emaxiso = pviso->GetMaxEnergy();
185 if(ekin <= emaxiso) {
186 xs = pviso->Value(ekin);
187 } else {
189 theGamma(gamma, G4ThreeVector(0,0,1.), rTransitionBound);
190 const G4double rxs = ggXsection->GetIsoCrossSection(&theGamma, Z, A);
191 const G4double lxs = pviso->Value(emaxiso);
192 xs = lxs + (ekin - emaxiso)*(rxs - lxs)/(rTransitionBound-emaxiso);
193 }
194#ifdef G4VERBOSE
195 if(verboseLevel > 1) {
196 G4cout << "G4GammaNuclearXS::IsoXS: Z= " << Z << " A= " << A
197 << " Ekin(MeV)= " << ekin/CLHEP::MeV
198 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
199 }
200#endif
201 return xs;
202 }
203 }
204
205 // use element x-section
206 // for the hydrogen target there is no element data
207 if(ekin <= emax && Z != 1) {
208 xs = pv->Value(ekin)*A/aeff[Z];
209
210 // CHIPS for high energy and for the hydrogen target
211 } else if(ekin >= rTransitionBound || Z == 1) {
212 if(Z <= 2 && ekin > 10.*GeV) {
213 xs = coeff[Z][A - amin[Z]]*
214 ggXsection->GetElementCrossSection(aParticle, Z, 0);
215 } else {
216 xs = ggXsection->GetIsoCrossSection(aParticle, Z, A);
217 }
218
219 // transition GDR to CHIPS
220 } else {
221 const G4double rxs = xs150[Z];
222 const G4double lxs = pv->Value(emax)*A/aeff[Z];
223 xs = lxs + (ekin - emax)*(rxs - lxs)/(rTransitionBound-emax);
224 }
225#ifdef G4VERBOSE
226 if(verboseLevel > 1) {
227 G4cout << "G4GammaNuclearXS::IsoXS: Z= " << Z << " A= " << A
228 << " Ekin(MeV)= " << ekin/CLHEP::MeV
229 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
230 }
231#endif
232 return xs;
233}
static const G4int amax[]
static const G4double aeff[]
static const G4int amin[]
static constexpr double GeV
Definition: G4SIunits.hh:203
G4PhysicsVector * GetComponentDataByIndex(G4int Z, G4int idx)
static G4double coeff[3][3]
G4double GetMaxEnergy() const

References A, aeff, amax, amin, CLHEP::barn, coeff, data, emax, G4cout, G4endl, gamma, G4ElementData::GetComponentDataByIndex(), G4VCrossSectionDataSet::GetElementCrossSection(), G4VCrossSectionDataSet::GetIsoCrossSection(), G4DynamicParticle::GetKineticEnergy(), G4PhysicsVector::GetMaxEnergy(), GetPhysicsVector(), GeV, ggXsection, MAXZGAMMAXS, CLHEP::MeV, rTransitionBound, G4VCrossSectionDataSet::verboseLevel, xs150, and Z.

Referenced by IsoCrossSection().

◆ GetMaxKinEnergy()

G4double G4VCrossSectionDataSet::GetMaxKinEnergy ( ) const
inlineinherited

◆ GetMinKinEnergy()

G4double G4VCrossSectionDataSet::GetMinKinEnergy ( ) const
inlineinherited

◆ GetName()

const G4String & G4VCrossSectionDataSet::GetName ( ) const
inlineinherited

◆ GetPhysicsVector()

G4PhysicsVector * G4GammaNuclearXS::GetPhysicsVector ( G4int  Z)
inlineprivate

Definition at line 140 of file G4GammaNuclearXS.hh.

141{
143 if(pv == nullptr) {
145 pv = data->GetElementData(Z);
146 }
147 return pv;
148}
void InitialiseOnFly(G4int Z)

References data, G4ElementData::GetElementData(), InitialiseOnFly(), and Z.

Referenced by GetElementCrossSection(), and GetIsoCrossSection().

◆ GetVerboseLevel()

G4int G4VCrossSectionDataSet::GetVerboseLevel ( ) const
inlinevirtualinherited

◆ Initialise()

void G4GammaNuclearXS::Initialise ( G4int  Z)
private

Definition at line 364 of file G4GammaNuclearXS.cc.

365{
366 if(nullptr != data->GetElementData(Z)) { return; }
367
368 // upload data from file
369 std::ostringstream ost;
370 ost << FindDirectoryPath() << Z ;
371 G4PhysicsVector* v = RetrieveVector(ost, true, Z);
372
374 /*
375 G4cout << "G4GammaNuclearXS::Initialise for Z= " << Z
376 << " A= " << Amean << " Amin= " << amin[Z]
377 << " Amax= " << amax[Z] << G4endl;
378 */
379 // upload isotope data
381 xs150[Z] = ggXsection->GetElementCrossSection(&theGamma, Z, 0);
382 if(amax[Z] > amin[Z]) {
383 size_t nmax = (size_t)(amax[Z]-amin[Z]+1);
385 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
386 std::ostringstream ost1;
387 ost1 << gDataDirectory << Z << "_" << A;
388 G4PhysicsVector* v1 = RetrieveVector(ost1, false, Z);
389 data->AddComponent(Z, A, v1);
390 if(Z<=2){
391 theGamma.SetKineticEnergy(10.*GeV);
392 G4double sig1 = ggXsection->GetIsoCrossSection(&theGamma, Z, A);
393 G4double sig2 = ggXsection->GetElementCrossSection(&theGamma, Z, 0);
394 if(sig2 > 0.) coeff[Z][A-amin[Z]]=(sig1/sig2);
395 else coeff[Z][A-amin[Z]]=1.;
396 }
397 }
398 }
399}
void InitialiseForComponent(G4int Z, G4int nComponents=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
G4PhysicsVector * RetrieveVector(std::ostringstream &in, G4bool isElement, G4int Z)

References A, G4ElementData::AddComponent(), amax, amin, coeff, data, FindDirectoryPath(), gamma, gDataDirectory, G4VCrossSectionDataSet::GetElementCrossSection(), G4ElementData::GetElementData(), G4VCrossSectionDataSet::GetIsoCrossSection(), GeV, ggXsection, G4ElementData::InitialiseForComponent(), G4ElementData::InitialiseForElement(), RetrieveVector(), rTransitionBound, G4DynamicParticle::SetKineticEnergy(), xs150, and Z.

Referenced by BuildPhysicsTable(), and InitialiseOnFly().

◆ InitialiseOnFly()

void G4GammaNuclearXS::InitialiseOnFly ( G4int  Z)
private

Definition at line 351 of file G4GammaNuclearXS.cc.

352{
353#ifdef G4MULTITHREADED
354 G4MUTEXLOCK(&gNuclearXSMutex);
355 if(nullptr == data->GetElementData(Z)) {
356#endif
357 Initialise(Z);
358#ifdef G4MULTITHREADED
359 }
360 G4MUTEXUNLOCK(&gNuclearXSMutex);
361#endif
362}

References data, G4MUTEXLOCK, G4MUTEXUNLOCK, G4ElementData::GetElementData(), Initialise(), and Z.

Referenced by GetPhysicsVector().

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 98 of file G4GammaNuclearXS.cc.

100{
101 return true;
102}

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 104 of file G4GammaNuclearXS.cc.

107{
108 return true;
109}

◆ IsoCrossSection()

G4double G4GammaNuclearXS::IsoCrossSection ( G4double  ekin,
G4int  Z,
G4int  A 
)

Definition at line 152 of file G4GammaNuclearXS.cc.

153{
154 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), ekin);
155 return GetIsoCrossSection(&theGamma, ZZ, A);
156}
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final

References A, gamma, and GetIsoCrossSection().

Referenced by SelectIsotope().

◆ operator=()

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

◆ RetrieveVector()

G4PhysicsVector * G4GammaNuclearXS::RetrieveVector ( std::ostringstream &  in,
G4bool  isElement,
G4int  Z 
)
private

Definition at line 402 of file G4GammaNuclearXS.cc.

403{
404 G4PhysicsVector* v = nullptr;
405
406 std::ifstream filein(ost.str().c_str());
407 if (!filein.is_open()) {
408 if(isElement) {
410 ed << "Data file <" << ost.str().c_str()
411 << "> is not opened!";
412 G4Exception("G4GammaNuclearXS::RetrieveVector(..)","had014",
413 FatalException, ed, "Check G4PARTICLEXSDATA");
414 }
415 } else {
416 if(verboseLevel > 1) {
417 G4cout << "File " << ost.str()
418 << " is opened by G4GammaNuclearXS" << G4endl;
419 }
420 // retrieve data from DB
421 if(std::find(std::begin(freeVectorException), std::end(freeVectorException), Z ) == std::end(freeVectorException) && isElement) {
422 v = new G4PhysicsLinearVector();
423 } else {
424 v = new G4PhysicsVector();
425 }
426 if(!v->Retrieve(filein, true)) {
428 ed << "Data file <" << ost.str().c_str()
429 << "> is not retrieved!";
430 G4Exception("G4GammaNuclearXS::RetrieveVector(..)","had015",
431 FatalException, ed, "Check G4PARTICLEXSDATA");
432 }
433 }
434 return v;
435}
const G4int freeVectorException[11]
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)

References FatalException, freeVectorException, G4cout, G4endl, G4Exception(), G4PhysicsVector::Retrieve(), G4VCrossSectionDataSet::verboseLevel, and Z.

Referenced by Initialise().

◆ SelectIsotope()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 235 of file G4GammaNuclearXS.cc.

237{
238 size_t nIso = anElement->GetNumberOfIsotopes();
239 const G4Isotope* iso = anElement->GetIsotope(0);
240
241 if(1 == nIso) { return iso; }
242
243 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
245 G4double sum = 0.0;
246 size_t j;
247 G4int Z = anElement->GetZasInt();
248
249 // condition to use only isotope abundance
250 if(amax[Z] == amin[Z] || kinEnergy > rTransitionBound || Z >= MAXZGAMMAXS ) {
251 for (j=0; j<nIso; ++j) {
252 sum += abundVector[j];
253 if(q <= sum) {
254 iso = anElement->GetIsotope(j);
255 break;
256 }
257 }
258 return iso;
259 }
260 // use isotope cross sections
261 size_t nn = temp.size();
262 if(nn < nIso) { temp.resize(nIso, 0.); }
263
264 for (j=0; j<nIso; ++j) {
265 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
266 // << " abund= " << abundVector[j] << G4endl;
267 sum += abundVector[j]*
268 IsoCrossSection(kinEnergy, Z, anElement->GetIsotope(j)->GetN());
269 temp[j] = sum;
270 }
271 sum *= q;
272 for (j = 0; j<nIso; ++j) {
273 if(temp[j] >= sum) {
274 iso = anElement->GetIsotope(j);
275 break;
276 }
277 }
278 return iso;
279}
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
G4double IsoCrossSection(G4double ekin, G4int Z, G4int A)

References amax, amin, G4UniformRand, G4Element::GetIsotope(), G4Isotope::GetN(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), G4Element::GetZasInt(), IsoCrossSection(), MAXZGAMMAXS, G4InuclParticleNames::nn, rTransitionBound, temp, and Z.

◆ SetForAllAtomsAndEnergies()

void G4VCrossSectionDataSet::SetForAllAtomsAndEnergies ( G4bool  val)
inlineinherited

◆ SetMaxKinEnergy()

void G4VCrossSectionDataSet::SetMaxKinEnergy ( G4double  value)
inlineinherited

◆ SetMinKinEnergy()

void G4VCrossSectionDataSet::SetMinKinEnergy ( G4double  value)
inlineinherited

◆ SetName()

void G4VCrossSectionDataSet::SetName ( const G4String nam)
inlineinherited

Definition at line 240 of file G4VCrossSectionDataSet.hh.

241{
242 name = nam;
243}

References G4VCrossSectionDataSet::name.

Referenced by G4ParticleHPInelasticData::G4ParticleHPInelasticData().

◆ SetVerboseLevel()

void G4VCrossSectionDataSet::SetVerboseLevel ( G4int  value)
inlinevirtualinherited

Field Documentation

◆ coeff

G4double G4GammaNuclearXS::coeff
staticprivate

Definition at line 130 of file G4GammaNuclearXS.hh.

Referenced by GetIsoCrossSection(), and Initialise().

◆ data

G4ElementData * G4GammaNuclearXS::data = nullptr
staticprivate

◆ freeVectorException

const G4int G4GammaNuclearXS::freeVectorException[11] = {4, 6, 7, 8, 27, 39, 45, 65, 67, 69, 73}
private

Definition at line 126 of file G4GammaNuclearXS.hh.

Referenced by RetrieveVector().

◆ gamma

const G4ParticleDefinition* G4GammaNuclearXS::gamma
private

◆ gDataDirectory

G4String G4GammaNuclearXS::gDataDirectory = ""
staticprivate

Definition at line 132 of file G4GammaNuclearXS.hh.

Referenced by FindDirectoryPath(), and Initialise().

◆ ggXsection

G4VCrossSectionDataSet* G4GammaNuclearXS::ggXsection = nullptr
private

◆ isForAllAtomsAndEnergies

G4bool G4VCrossSectionDataSet::isForAllAtomsAndEnergies
privateinherited

◆ isMaster

G4bool G4GammaNuclearXS::isMaster = false
private

Definition at line 119 of file G4GammaNuclearXS.hh.

Referenced by BuildPhysicsTable(), and ~G4GammaNuclearXS().

◆ maxKinEnergy

G4double G4VCrossSectionDataSet::maxKinEnergy
privateinherited

◆ MAXZGAMMAXS

const G4int G4GammaNuclearXS::MAXZGAMMAXS = 95
staticprivate

◆ minKinEnergy

G4double G4VCrossSectionDataSet::minKinEnergy
privateinherited

◆ name

G4String G4VCrossSectionDataSet::name
protectedinherited

◆ registry

G4CrossSectionDataSetRegistry* G4VCrossSectionDataSet::registry
privateinherited

◆ rTransitionBound

const G4int G4GammaNuclearXS::rTransitionBound = 150.*CLHEP::MeV
staticprivate

◆ temp

std::vector<G4double> G4GammaNuclearXS::temp
private

Definition at line 116 of file G4GammaNuclearXS.hh.

Referenced by BuildPhysicsTable(), and SelectIsotope().

◆ verboseLevel

G4int G4VCrossSectionDataSet::verboseLevel
protectedinherited

Definition at line 168 of file G4VCrossSectionDataSet.hh.

Referenced by G4BGGNucleonElasticXS::BuildPhysicsTable(), G4BGGPionElasticXS::BuildPhysicsTable(), G4BGGPionInelasticXS::BuildPhysicsTable(), BuildPhysicsTable(), G4NeutronCaptureXS::BuildPhysicsTable(), G4NeutronElasticXS::BuildPhysicsTable(), G4NeutronInelasticXS::BuildPhysicsTable(), G4ParticleInelasticXS::BuildPhysicsTable(), G4BGGNucleonInelasticXS::BuildPhysicsTable(), G4LENDCrossSection::create_used_target_map(), G4BGGNucleonElasticXS::G4BGGNucleonElasticXS(), G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS(), G4BGGPionElasticXS::G4BGGPionElasticXS(), G4BGGPionInelasticXS::G4BGGPionInelasticXS(), G4GammaNuclearXS(), G4NeutronCaptureXS::G4NeutronCaptureXS(), G4NeutronElasticXS::G4NeutronElasticXS(), G4NeutronInelasticXS::G4NeutronInelasticXS(), G4ParticleInelasticXS::G4ParticleInelasticXS(), G4NeutronCaptureXS::GetElementCrossSection(), G4NeutronElasticXS::GetElementCrossSection(), G4NeutronInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetElementCrossSection(), G4BGGPionElasticXS::GetElementCrossSection(), G4BGGPionInelasticXS::GetElementCrossSection(), G4BGGNucleonInelasticXS::GetElementCrossSection(), GetElementCrossSection(), G4ParticleInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetIsoCrossSection(), G4BGGPionElasticXS::GetIsoCrossSection(), G4BGGPionInelasticXS::GetIsoCrossSection(), GetIsoCrossSection(), G4BGGNucleonInelasticXS::GetIsoCrossSection(), G4VCrossSectionDataSet::GetVerboseLevel(), G4NeutronElasticXS::Initialise(), G4ParticleInelasticXS::IsoCrossSection(), G4NeutronCaptureXS::IsoCrossSection(), G4NeutronInelasticXS::IsoCrossSection(), RetrieveVector(), G4NeutronCaptureXS::RetrieveVector(), G4NeutronInelasticXS::RetrieveVector(), G4ParticleInelasticXS::RetrieveVector(), and G4VCrossSectionDataSet::SetVerboseLevel().

◆ xs150

G4double G4GammaNuclearXS::xs150 = {0.0}
staticprivate

Definition at line 131 of file G4GammaNuclearXS.hh.

Referenced by GetElementCrossSection(), GetIsoCrossSection(), and Initialise().


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