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

#include <G4CrossSectionDataStore.hh>

Public Member Functions

 G4CrossSectionDataStore ()
 
 ~G4CrossSectionDataStore ()
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Material *)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Material *)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *)
 
G4double GetCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *, const G4Element *, const G4Material *)
 
const G4ElementSampleZandA (const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DumpPhysicsTable (const G4ParticleDefinition &)
 
void DumpHtml (const G4ParticleDefinition &, std::ofstream &) const
 
void PrintCrossSectionHtml (const G4VCrossSectionDataSet *cs) const
 
void AddDataSet (G4VCrossSectionDataSet *)
 
void AddDataSet (G4VCrossSectionDataSet *, size_t)
 
void SetVerboseLevel (G4int value)
 
const G4FastPathHadronicCrossSection::fastPathParametersGetFastPathParameters () const
 
const G4FastPathHadronicCrossSection::controlFlagGetFastPathControlFlags () const
 
void DumpFastPath (const G4ParticleDefinition *, const G4Material *, std::ostream &os)
 
void ActivateFastPath (const G4ParticleDefinition *, const G4Material *, G4double)
 

Friends

struct G4FastPathHadronicCrossSection::fastPathEntry
 

Detailed Description

Definition at line 62 of file G4CrossSectionDataStore.hh.

Constructor & Destructor Documentation

◆ G4CrossSectionDataStore()

G4CrossSectionDataStore::G4CrossSectionDataStore ( )

Definition at line 60 of file G4CrossSectionDataStore.cc.

60 :
61 nDataSetList(0), verboseLevel(0),fastPathFlags(),fastPathParams(),
62 counters(),fastPathCache()
63{
65 currentMaterial = elmMaterial = nullptr;
66 currentElement = nullptr; //ALB 14-Aug-2012 Coverity fix.
67 matParticle = elmParticle = nullptr;
68 matKinEnergy = elmKinEnergy = matCrossSection = elmCrossSection = 0.0;
69}
static G4NistManager * Instance()

◆ ~G4CrossSectionDataStore()

G4CrossSectionDataStore::~G4CrossSectionDataStore ( )

Definition at line 73 of file G4CrossSectionDataStore.cc.

74{}

Member Function Documentation

◆ ActivateFastPath()

void G4CrossSectionDataStore::ActivateFastPath ( const G4ParticleDefinition pdef,
const G4Material mat,
G4double  min_cutoff 
)

Definition at line 531 of file G4CrossSectionDataStore.cc.

533{
534 assert(pdef!=nullptr&&mat!=nullptr);
536 if ( requests.insert( { key , min_cutoff } ).second ) {
538 ed << "Attempting to request FastPath for couple: <"
539 << pdef->GetParticleName() << ", " <<mat->GetName()
540 << "> but combination already exists" << G4endl;
541 G4Exception("G4CrossSectionDataStore::ActivateFastPath", "had001",
542 FatalException, ed);
543 }
544}
@ 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 G4endl
Definition: G4ios.hh:57
const G4String & GetName() const
Definition: G4Material.hh:175
const G4String & GetParticleName() const
std::pair< const G4ParticleDefinition *, const G4Material * > G4CrossSectionDataStore_Key

◆ AddDataSet() [1/2]

void G4CrossSectionDataStore::AddDataSet ( G4VCrossSectionDataSet p)

Definition at line 644 of file G4CrossSectionDataStore.cc.

645{
646 if(p->ForAllAtomsAndEnergies()) {
647 dataSetList.clear();
648 nDataSetList = 0;
649 }
650 dataSetList.push_back(p);
651 ++nDataSetList;
652}

Referenced by G4HadronicProcess::AddDataSet(), G4ElectronNuclearProcess::G4ElectronNuclearProcess(), G4PhotoNuclearProcess::G4PhotoNuclearProcess(), G4PositronNuclearProcess::G4PositronNuclearProcess(), and G4HadronPhysicsShielding::Neutron().

◆ AddDataSet() [2/2]

void G4CrossSectionDataStore::AddDataSet ( G4VCrossSectionDataSet p,
size_t  i 
)

Definition at line 658 of file G4CrossSectionDataStore.cc.

659{
660 if(p->ForAllAtomsAndEnergies()) {
661 dataSetList.clear();
662 dataSetList.push_back(p);
663 nDataSetList = 1;
664 } else {
665 if ( i > dataSetList.size() ) i = dataSetList.size();
666 std::vector< G4VCrossSectionDataSet* >::iterator it = dataSetList.end() - i;
667 dataSetList.insert(it , p);
668 ++nDataSetList;
669 }
670}

◆ BuildPhysicsTable()

void G4CrossSectionDataStore::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)

Definition at line 494 of file G4CrossSectionDataStore.cc.

495{
496 if (nDataSetList == 0) {
498 ed << "No cross section is registered for "
499 << aParticleType.GetParticleName() << G4endl;
500 G4Exception("G4CrossSectionDataStore::BuildPhysicsTable", "had001",
501 FatalException, ed);
502 return;
503 }
504 for (G4int i=0; i<nDataSetList; ++i) {
505 dataSetList[i]->BuildPhysicsTable(aParticleType);
506 }
507 //A.Dotti: if fast-path has been requested we can now create the surrogate
508 // model for fast path.
509 if ( fastPathFlags.useFastPathIfAvailable ) {
510 fastPathFlags.initializationPhase = true;
511 using my_value_type=G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Requests::value_type;
512 //Loop on all requests, if particle matches create the corresponding fsat-path
513 std::for_each( requests.begin() , requests.end() ,
514 [&aParticleType,this](const my_value_type& req) {
515 if ( aParticleType == *req.part_mat.first ) {
516 G4FastPathHadronicCrossSection::cycleCountEntry* entry =
517 new G4FastPathHadronicCrossSection::cycleCountEntry(aParticleType.GetParticleName(),req.part_mat.second);
518 entry->fastPath =
519 new G4FastPathHadronicCrossSection::fastPathEntry(&aParticleType,req.part_mat.second,req.min_cutoff);
520 entry->fastPath->Initialize(this);
521 fastPathCache[req.part_mat] = entry;
522 }
523 }
524 );
525 fastPathFlags.initializationPhase = false;
526 }
527}
int G4int
Definition: G4Types.hh:85

Referenced by G4HadronicProcess::BuildPhysicsTable().

◆ ComputeCrossSection()

G4double G4CrossSectionDataStore::ComputeCrossSection ( const G4DynamicParticle part,
const G4Material mat 
)

Definition at line 268 of file G4CrossSectionDataStore.cc.

270{
271 if(mat == currentMaterial && part->GetDefinition() == matParticle
272 && part->GetKineticEnergy() == matKinEnergy) {
273 return matCrossSection;
274 }
275 currentMaterial = mat;
276 matParticle = part->GetDefinition();
277 matKinEnergy = part->GetKineticEnergy();
278 matCrossSection = 0.0;
279
280 size_t nElements = mat->GetNumberOfElements();
281 const G4double* nAtomsPerVolume = mat->GetVecNbOfAtomsPerVolume();
282
283 if(xsecelm.size() < nElements) { xsecelm.resize(nElements); }
284
285 for(size_t i=0; i<nElements; ++i) {
286 matCrossSection += nAtomsPerVolume[i] *
287 GetCrossSection(part, mat->GetElement(i), mat);
288 xsecelm[i] = matCrossSection;
289 }
290 return matCrossSection;
291}
double G4double
Definition: G4Types.hh:83
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204

Referenced by G4GammaGeneralProcess::BuildPhysicsTable(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4HadronicProcess::GetMeanFreePath(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4HadronXSDataTable::Initialise(), G4HadronicProcess::PostStepDoIt(), and G4GammaGeneralProcess::SampleHadSecondaries().

◆ DumpFastPath()

void G4CrossSectionDataStore::DumpFastPath ( const G4ParticleDefinition pd,
const G4Material mat,
std::ostream &  os 
)

Definition at line 249 of file G4CrossSectionDataStore.cc.

250{
251 const G4FastPathHadronicCrossSection::cycleCountEntry* entry = fastPathCache[{pd,mat}];
252 if ( entry != nullptr ) {
253 if ( entry->fastPath != nullptr ) {
254 os<<*entry->fastPath;
255 } else {
256 os<<"#Cache entry for {"<<(pd!=nullptr?pd->GetParticleName():"UNDEFINED")<<",";
257 os<<(mat!=nullptr?mat->GetName():"UNDEFINED")<<"} found, but no fast path defined";
258 }
259 } else {
260 os<<"#Cache entry for {"<<(pd!=nullptr?pd->GetParticleName():"UNDEFINED")<<",";
261 os<<(mat!=nullptr?mat->GetName():"UNDEFINED")<<"} not found.";
262 }
263}

◆ DumpHtml()

void G4CrossSectionDataStore::DumpHtml ( const G4ParticleDefinition ,
std::ofstream &  outFile 
) const

Definition at line 576 of file G4CrossSectionDataStore.cc.

578{
579 // Write cross section data set info to html physics list
580 // documentation page
581
582 G4double ehi = 0;
583 G4double elo = 0;
584 G4String physListName(std::getenv("G4PhysListName"));
585 for (G4int i = nDataSetList-1; i > 0; i--) {
586 elo = dataSetList[i]->GetMinKinEnergy()/GeV;
587 ehi = dataSetList[i]->GetMaxKinEnergy()/GeV;
588 outFile << " <li><b><a href=\"" << physListName << "_"
589 << dataSetList[i]->GetName() << ".html\"> "
590 << dataSetList[i]->GetName() << "</a> from "
591 << elo << " GeV to " << ehi << " GeV </b></li>\n";
592 //G4cerr << i << ": XS for " << pD.GetParticleName() << " : " << dataSetList[i]->GetName()
593 // << " typeid : " << typeid(dataSetList[i]).name()<< G4endl;
594 PrintCrossSectionHtml(dataSetList[i]);
595 }
596
597 G4double defaultHi = dataSetList[0]->GetMaxKinEnergy()/GeV;
598 if (ehi < defaultHi) {
599 outFile << " <li><b><a href=\"" << dataSetList[0]->GetName() << ".html\"> "
600 << dataSetList[0]->GetName() << "</a> from "
601 << ehi << " GeV to " << defaultHi << " GeV </b></li>\n";
602 PrintCrossSectionHtml(dataSetList[0]);
603 }
604}
void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs) const

Referenced by G4HadronicProcessStore::PrintHtml().

◆ DumpPhysicsTable()

void G4CrossSectionDataStore::DumpPhysicsTable ( const G4ParticleDefinition aParticleType)

Definition at line 549 of file G4CrossSectionDataStore.cc.

550{
551 // Print out all cross section data sets used and the energies at
552 // which they apply
553
554 if (nDataSetList == 0) {
555 G4cout << "WARNING - G4CrossSectionDataStore::DumpPhysicsTable: "
556 << " no data sets registered" << G4endl;
557 return;
558 }
559
560 for (G4int i = nDataSetList-1; i >= 0; --i) {
561 G4double e1 = dataSetList[i]->GetMinKinEnergy();
562 G4double e2 = dataSetList[i]->GetMaxKinEnergy();
563 G4cout
564 << " Cr_sctns: " << std::setw(25) << dataSetList[i]->GetName() << ": "
565 << G4BestUnit(e1, "Energy")
566 << " ---> "
567 << G4BestUnit(e2, "Energy") << "\n";
568 if (dataSetList[i]->GetName() == "G4CrossSectionPairGG") {
569 dataSetList[i]->DumpPhysicsTable(aParticleType);
570 }
571 }
572}
#define G4BestUnit(a, b)
G4GLOB_DLL std::ostream G4cout

Referenced by G4ChargeExchangeProcess::DumpPhysicsTable(), and G4HadronicProcess::DumpPhysicsTable().

◆ GetCrossSection() [1/3]

G4double G4CrossSectionDataStore::GetCrossSection ( const G4DynamicParticle part,
const G4Element elm,
const G4Material mat 
)

Definition at line 296 of file G4CrossSectionDataStore.cc.

299{
300 if(mat == elmMaterial && elm == currentElement &&
301 part->GetDefinition() == elmParticle &&
302 part->GetKineticEnergy() == elmKinEnergy)
303 { return elmCrossSection; }
304
305 elmMaterial = mat;
306 currentElement = elm;
307 elmParticle = part->GetDefinition();
308 elmKinEnergy = part->GetKineticEnergy();
309 elmCrossSection = 0.0;
310
311 G4int i = nDataSetList-1;
312 G4int Z = elm->GetZasInt();
313 if (elm->GetNaturalAbundanceFlag() &&
314 dataSetList[i]->IsElementApplicable(part, Z, mat)) {
315
316 // element wise cross section
317 elmCrossSection = dataSetList[i]->GetElementCrossSection(part, Z, mat);
318
319 //G4cout << "Element wise " << elmParticle->GetParticleName()
320 // << " xsec(barn)= " << elmCrossSection/barn
321 // << " E(MeV)= " << elmKinEnergy/MeV
322 // << " Z= " << Z << " AbundFlag= " << elm->GetNaturalAbandancesFlag()
323 // <<G4endl;
324
325 } else {
326 // isotope wise cross section
327 size_t nIso = elm->GetNumberOfIsotopes();
328
329 // user-defined isotope abundances
330 const G4double* abundVector = elm->GetRelativeAbundanceVector();
331
332 for (size_t j=0; j<nIso; ++j) {
333 if(abundVector[j] > 0.0) {
334 const G4Isotope* iso = elm->GetIsotope(j);
335 elmCrossSection += abundVector[j]*
336 GetIsoCrossSection(part, Z, iso->GetN(), iso, elm, mat, i);
337 //G4cout << "Isotope wise " << elmParticle->GetParticleName()
338 // << " xsec(barn)= " << elmCrossSection/barn
339 // << " E(MeV)= " << elmKinEnergy/MeV
340 // << " Z= " << Z << " A= " << iso->GetN() << " j= " << j << G4endl;
341 }
342 }
343 }
344 //G4cout << " E(MeV)= " << elmKinEnergy/MeV
345 // << "xsec(barn)= " << elmCrossSection/barn <<G4endl;
346 return elmCrossSection;
347}
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4bool GetNaturalAbundanceFlag() const
Definition: G4Element.hh:261
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4int GetZasInt() const
Definition: G4Element.hh:131
G4int GetN() const
Definition: G4Isotope.hh:93

◆ GetCrossSection() [2/3]

G4double G4CrossSectionDataStore::GetCrossSection ( const G4DynamicParticle particle,
const G4Material material 
)
inline

◆ GetCrossSection() [3/3]

G4double G4CrossSectionDataStore::GetCrossSection ( const G4DynamicParticle part,
G4int  Z,
G4int  A,
const G4Isotope iso,
const G4Element elm,
const G4Material mat 
)

Definition at line 390 of file G4CrossSectionDataStore.cc.

395{
396 for (G4int i = nDataSetList-1; i >= 0; --i) {
397 if (dataSetList[i]->IsIsoApplicable(part, Z, A, elm, mat) ) {
398 return dataSetList[i]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
399 }
400 }
402 ed << "No isotope cross section found for "
403 << part->GetDefinition()->GetParticleName()
404 << " off Element " << elm->GetName()
405 << " in " << mat->GetName() << " Z= " << Z << " A= " << A
406 << " E(MeV)= " << part->GetKineticEnergy()/MeV << G4endl;
407 G4Exception("G4CrossSectionDataStore::GetCrossSection", "had001",
408 FatalException, ed);
409 return 0.0;
410}
double A(double temperature)
const G4String & GetName() const
Definition: G4Element.hh:126

◆ GetFastPathControlFlags()

const G4FastPathHadronicCrossSection::controlFlag & G4CrossSectionDataStore::GetFastPathControlFlags ( ) const
inline

Definition at line 138 of file G4CrossSectionDataStore.hh.

138{ return fastPathFlags; }

Referenced by G4FastPathHadronicCrossSection::fastPathEntry::Initialize().

◆ GetFastPathParameters()

const G4FastPathHadronicCrossSection::fastPathParameters & G4CrossSectionDataStore::GetFastPathParameters ( ) const
inline

Definition at line 136 of file G4CrossSectionDataStore.hh.

136{ return fastPathParams; }

Referenced by G4FastPathHadronicCrossSection::fastPathEntry::Initialize().

◆ PrintCrossSectionHtml()

void G4CrossSectionDataStore::PrintCrossSectionHtml ( const G4VCrossSectionDataSet cs) const

Definition at line 608 of file G4CrossSectionDataStore.cc.

609{
610 G4String dirName(std::getenv("G4PhysListDocDir"));
611 G4String physListName(std::getenv("G4PhysListName"));
612
613 G4String pathName = dirName + "/" + physListName + "_" + HtmlFileName(cs->GetName());
614 std::ofstream outCS;
615 outCS.open(pathName);
616 outCS << "<html>\n";
617 outCS << "<head>\n";
618 outCS << "<title>Description of " << cs->GetName()
619 << "</title>\n";
620 outCS << "</head>\n";
621 outCS << "<body>\n";
622
623 cs->CrossSectionDescription(outCS);
624
625 outCS << "</body>\n";
626 outCS << "</html>\n";
627
628}
const G4String & GetName() const
virtual void CrossSectionDescription(std::ostream &) const

Referenced by DumpHtml().

◆ SampleZandA()

const G4Element * G4CrossSectionDataStore::SampleZandA ( const G4DynamicParticle part,
const G4Material mat,
G4Nucleus target 
)

Definition at line 415 of file G4CrossSectionDataStore.cc.

418{
419 size_t nElements = mat->GetNumberOfElements();
420 const G4Element* anElement = mat->GetElement(0);
421
422 // select element from a compound
423 if(1 < nElements) {
424 G4double cross = matCrossSection*G4UniformRand();
425 for(size_t i=0; i<nElements; ++i) {
426 if(cross <= xsecelm[i]) {
427 anElement = mat->GetElement(i);
428 break;
429 }
430 }
431 }
432
433 G4int Z = anElement->GetZasInt();
434 const G4Isotope* iso = nullptr;
435
436 G4int i = nDataSetList-1;
437 if (dataSetList[i]->IsElementApplicable(part, Z, mat)) {
438
439 //----------------------------------------------------------------
440 // element-wise cross section
441 // isotope cross section is not computed
442 //----------------------------------------------------------------
443 size_t nIso = anElement->GetNumberOfIsotopes();
444 iso = anElement->GetIsotope(0);
445
446 // more than 1 isotope
447 if(1 < nIso) {
448 iso = dataSetList[i]->SelectIsotope(anElement,
449 part->GetKineticEnergy(),
450 part->GetLogKineticEnergy());
451 }
452 } else {
453
454 //----------------------------------------------------------------
455 // isotope-wise cross section
456 // isotope cross section is computed
457 //----------------------------------------------------------------
458 size_t nIso = anElement->GetNumberOfIsotopes();
459 iso = anElement->GetIsotope(0);
460
461 // more than 1 isotope
462 if(1 < nIso) {
463 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
464 if(xseciso.size() < nIso) { xseciso.resize(nIso); }
465
466 G4double cross = 0.0;
467 size_t j;
468 for (j = 0; j<nIso; ++j) {
469 G4double xsec = 0.0;
470 if(abundVector[j] > 0.0) {
471 iso = anElement->GetIsotope(j);
472 xsec = abundVector[j]*
473 GetIsoCrossSection(part, Z, iso->GetN(), iso, anElement, mat, i);
474 }
475 cross += xsec;
476 xseciso[j] = cross;
477 }
478 cross *= G4UniformRand();
479 for (j = 0; j<nIso; ++j) {
480 if(cross <= xseciso[j]) {
481 iso = anElement->GetIsotope(j);
482 break;
483 }
484 }
485 }
486 }
487 target.SetIsotope(iso);
488 return anElement;
489}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetLogKineticEnergy() const
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122

Referenced by G4NeutrinoElectronProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), and G4HadronElasticProcess::PostStepDoIt().

◆ SetVerboseLevel()

void G4CrossSectionDataStore::SetVerboseLevel ( G4int  value)
inline

Definition at line 161 of file G4CrossSectionDataStore.hh.

162{
163 verboseLevel = value;
164}

Friends And Related Function Documentation

◆ G4FastPathHadronicCrossSection::fastPathEntry

Definition at line 142 of file G4CrossSectionDataStore.hh.


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