Geant4 11.1.1
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=nullptr) final
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
 
const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final
 
void BuildPhysicsTable (const G4ParticleDefinition &) final
 
G4double IsoCrossSection (G4double ekin, G4int Z, G4int A)
 
G4double ElementCrossSection (G4double ekin, G4int Z)
 
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 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
 

Static Public Member Functions

static const char * Default_Name ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 
G4String name
 

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.

69 gamma(G4Gamma::Gamma())
70{
71 // verboseLevel = 0;
72 if(verboseLevel > 0){
73 G4cout << "G4GammaNuclearXS::G4GammaNuclearXS Initialise for Z < "
74 << MAXZGAMMAXS << G4endl;
75 }
77 if(ggXsection == nullptr) ggXsection = new G4PhotoNuclearCrossSection();
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 G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetForAllAtomsAndEnergies(G4bool val)

◆ ~G4GammaNuclearXS()

G4GammaNuclearXS::~G4GammaNuclearXS ( )
final

Definition at line 81 of file G4GammaNuclearXS.cc.

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

◆ 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");
305 FindDirectoryPath();
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 std::size_t nIso = temp.size();
325 for ( auto & elm : *table ) {
326 std::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:59
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:403
const G4String & GetParticleName() const

◆ 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";}

Referenced by G4ElectroVDNuclearModel::G4ElectroVDNuclearModel().

◆ 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

◆ 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}
double G4double
Definition: G4Types.hh:83
G4double GetKineticEnergy() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)

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}
const G4double A[17]
G4PhysicsVector * GetComponentDataByIndex(G4int Z, G4int idx)
G4double GetMaxEnergy() const
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)

Referenced by IsoCrossSection().

◆ 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

Referenced by SelectIsotope().

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

237{
238 std::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 G4int 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<(G4int)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 std::size_t nn = temp.size();
262 if(nn < nIso) { temp.resize(nIso, 0.); }
263
264 for (j=0; j<(G4int)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<(G4int)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
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetZasInt() const
Definition: G4Element.hh:132
G4double IsoCrossSection(G4double ekin, G4int Z, G4int A)
G4int GetN() const
Definition: G4Isotope.hh:93

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