Geant4  10.02.p03
G4CrossSectionDataStore Class Reference

#include <G4CrossSectionDataStore.hh>

Collaboration diagram for G4CrossSectionDataStore:

Public Member Functions

 G4CrossSectionDataStore ()
 
 ~G4CrossSectionDataStore ()
 
G4double GetCrossSection (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 *)
 
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 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)
 

Private Member Functions

G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *, const G4Element *, const G4Material *aMaterial, G4int index)
 
G4CrossSectionDataStoreoperator= (const G4CrossSectionDataStore &right)
 
 G4CrossSectionDataStore (const G4CrossSectionDataStore &)
 
G4String HtmlFileName (const G4String &in) const
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Material *, G4bool requiresSlowPath)
 

Private Attributes

G4NistManagernist
 
std::vector< G4VCrossSectionDataSet * > dataSetList
 
std::vector< G4doublexsecelm
 
std::vector< G4doublexseciso
 
const G4MaterialcurrentMaterial
 
const G4ParticleDefinitionmatParticle
 
G4double matKinEnergy
 
G4double matCrossSection
 
const G4MaterialelmMaterial
 
const G4ElementcurrentElement
 
const G4ParticleDefinitionelmParticle
 
G4double elmKinEnergy
 
G4double elmCrossSection
 
G4int nDataSetList
 
G4int verboseLevel
 
G4FastPathHadronicCrossSection::controlFlag fastPathFlags
 
G4FastPathHadronicCrossSection::fastPathParameters fastPathParams
 
G4FastPathHadronicCrossSection::getCrossSectionCount counters
 
G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Cache fastPathCache
 
G4FastPathHadronicCrossSection::timing timing
 
G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Requests requests
 

Friends

struct G4FastPathHadronicCrossSection::fastPathEntry
 

Detailed Description

Definition at line 63 of file G4CrossSectionDataStore.hh.

Constructor & Destructor Documentation

◆ G4CrossSectionDataStore() [1/2]

G4CrossSectionDataStore::G4CrossSectionDataStore ( )

Definition at line 63 of file G4CrossSectionDataStore.cc.

63  :
66 {
69  currentElement = 0; //ALB 14-Aug-2012 Coverity fix.
72 }
G4FastPathHadronicCrossSection::controlFlag fastPathFlags
const G4ParticleDefinition * matParticle
G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Cache fastPathCache
G4FastPathHadronicCrossSection::getCrossSectionCount counters
static G4NistManager * Instance()
const G4ParticleDefinition * elmParticle
G4FastPathHadronicCrossSection::fastPathParameters fastPathParams
Here is the call graph for this function:

◆ ~G4CrossSectionDataStore()

G4CrossSectionDataStore::~G4CrossSectionDataStore ( )

Definition at line 76 of file G4CrossSectionDataStore.cc.

77 {}

◆ G4CrossSectionDataStore() [2/2]

G4CrossSectionDataStore::G4CrossSectionDataStore ( const G4CrossSectionDataStore )
private

Member Function Documentation

◆ ActivateFastPath()

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

Definition at line 534 of file G4CrossSectionDataStore.cc.

535 {
536  assert(pdef!=nullptr&&mat!=nullptr);
538  if ( requests.insert( { key , min_cutoff } ).second ) {
539  std::ostringstream msg;
540  msg<<"Attempting to request FastPath for couple: "<<pdef->GetParticleName()<<","<<mat->GetName();
541  msg<<" but combination already exists";
542  G4HadronicException(__FILE__,__LINE__,msg.str());
543  }
544 }
std::pair< const G4ParticleDefinition *, const G4Material * > G4CrossSectionDataStore_Key
const G4String & GetParticleName() const
G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Requests requests
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AddDataSet()

void G4CrossSectionDataStore::AddDataSet ( G4VCrossSectionDataSet p)
inline

Definition at line 160 of file G4CrossSectionDataStore.hh.

161 {
162  dataSetList.push_back(p);
163  ++nDataSetList;
164 }
std::vector< G4VCrossSectionDataSet * > dataSetList
Here is the caller graph for this function:

◆ BuildPhysicsTable()

void G4CrossSectionDataStore::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)

Definition at line 499 of file G4CrossSectionDataStore.cc.

500 {
501  if (nDataSetList == 0)
502  {
503  throw G4HadronicException(__FILE__, __LINE__,
504  "G4CrossSectionDataStore: no data sets registered");
505  return;
506  }
507  for (G4int i=0; i<nDataSetList; ++i) {
508  dataSetList[i]->BuildPhysicsTable(aParticleType);
509  }
510  //A.Dotti: if fast-path has been requested we can now create the surrogate
511  // model for fast path.
514  using my_value_type=G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Requests::value_type;
515  //Loop on all requests, if particle matches create the corresponding fsat-path
516  std::for_each( requests.begin() , requests.end() ,
517  [&aParticleType,this](const my_value_type& req) {
518  if ( aParticleType == *req.part_mat.first ) {
520  new G4FastPathHadronicCrossSection::cycleCountEntry(aParticleType.GetParticleName(),req.part_mat.second);
521  entry->fastPath =
522  new G4FastPathHadronicCrossSection::fastPathEntry(&aParticleType,req.part_mat.second,req.min_cutoff);
523  entry->fastPath->Initialize(this);
524  fastPathCache[req.part_mat] = entry;
525  }
526  }
527  );
529  }
530 }
G4FastPathHadronicCrossSection::controlFlag fastPathFlags
void Initialize(G4CrossSectionDataStore *)
G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Cache fastPathCache
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
std::vector< G4VCrossSectionDataSet * > dataSetList
fastPathEntry * fastPath
G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Requests requests
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DumpFastPath()

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

Definition at line 253 of file G4CrossSectionDataStore.cc.

254 {
256  if ( entry != nullptr ) {
257  if ( entry->fastPath != nullptr ) {
258  os<<*entry->fastPath;
259  } else {
260  os<<"#Cache entry for {"<<(pd!=nullptr?pd->GetParticleName():"UNDEFINED")<<",";
261  os<<(mat!=nullptr?mat->GetName():"UNDEFINED")<<"} found, but no fast path defined";
262  }
263  } else {
264  os<<"#Cache entry for {"<<(pd!=nullptr?pd->GetParticleName():"UNDEFINED")<<",";
265  os<<(mat!=nullptr?mat->GetName():"UNDEFINED")<<"} not found.";
266  }
267 }
G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Cache fastPathCache
const G4String & GetParticleName() const
fastPathEntry * fastPath
const G4String & GetName() const
Definition: G4Material.hh:178
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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(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;
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";
603  }
604 }
int G4int
Definition: G4Types.hh:78
std::vector< G4VCrossSectionDataSet * > dataSetList
static const double GeV
Definition: G4SIunits.hh:214
void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs) const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 }
static const G4double e2
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
std::vector< G4VCrossSectionDataSet * > dataSetList
static const G4double e1
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetCrossSection() [1/4]

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

Definition at line 155 of file G4CrossSectionDataStore.hh.

155  {
156  //By default tries to use the fast-path mechanism
157  return GetCrossSection( particle , material , false);
158 }
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
Here is the caller graph for this function:

◆ GetCrossSection() [2/4]

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

Definition at line 272 of file G4CrossSectionDataStore.cc.

275 {
276  if(mat == elmMaterial && elm == currentElement &&
277  part->GetDefinition() == elmParticle &&
278  part->GetKineticEnergy() == elmKinEnergy)
279  { return elmCrossSection; }
280 
281  elmMaterial = mat;
282  currentElement = elm;
283  elmParticle = part->GetDefinition();
284  elmKinEnergy = part->GetKineticEnergy();
285  elmCrossSection = 0.0;
286 
287  G4int i = nDataSetList-1;
288  G4int Z = G4lrint(elm->GetZ());
289  if (elm->GetNaturalAbundanceFlag() &&
290  dataSetList[i]->IsElementApplicable(part, Z, mat)) {
291 
292  // element wise cross section
293  elmCrossSection = dataSetList[i]->GetElementCrossSection(part, Z, mat);
294 
295  //G4cout << "Element wise " << elmParticle->GetParticleName()
296  // << " xsec(barn)= " << elmCrossSection/barn
297  // << " E(MeV)= " << elmKinEnergy/MeV
298  // << " Z= " << Z << " AbundFlag= " << elm->GetNaturalAbandancesFlag()
299  // <<G4endl;
300 
301  } else {
302  // isotope wise cross section
303  G4int nIso = elm->GetNumberOfIsotopes();
304  G4Isotope* iso = 0;
305 
306  // user-defined isotope abundances
307  G4IsotopeVector* isoVector = elm->GetIsotopeVector();
308  G4double* abundVector = elm->GetRelativeAbundanceVector();
309 
310  for (G4int j = 0; j<nIso; ++j) {
311  if(abundVector[j] > 0.0) {
312  iso = (*isoVector)[j];
313  elmCrossSection += abundVector[j]*
314  GetIsoCrossSection(part, Z, iso->GetN(), iso, elm, mat, i);
315  //G4cout << "Isotope wise " << elmParticle->GetParticleName()
316  // << " xsec(barn)= " << elmCrossSection/barn
317  // << " E(MeV)= " << elmKinEnergy/MeV
318  // << " Z= " << Z << " A= " << iso->GetN() << " j= " << j << G4endl;
319  }
320  }
321  }
322  //G4cout << " E(MeV)= " << elmKinEnergy/MeV
323  // << "xsec(barn)= " << elmCrossSection/barn <<G4endl;
324  return elmCrossSection;
325 }
G4bool GetNaturalAbundanceFlag() const
Definition: G4Element.hh:263
std::vector< G4Isotope * > G4IsotopeVector
int G4int
Definition: G4Types.hh:78
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
Float_t mat
G4double GetKineticEnergy() const
const G4ParticleDefinition * elmParticle
Float_t Z
std::vector< G4VCrossSectionDataSet * > dataSetList
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *, const G4Element *, const G4Material *aMaterial, G4int index)
int G4lrint(double ad)
Definition: templates.hh:163
G4int GetN() const
Definition: G4Isotope.hh:94
G4ParticleDefinition * GetDefinition() const
double G4double
Definition: G4Types.hh:76
G4double GetZ() const
Definition: G4Element.hh:131
Here is the call graph for this function:

◆ GetCrossSection() [3/4]

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

Definition at line 371 of file G4CrossSectionDataStore.cc.

376 {
377  for (G4int i = nDataSetList-1; i >= 0; --i) {
378  if (dataSetList[i]->IsIsoApplicable(part, Z, A, elm, mat) ) {
379  return dataSetList[i]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
380  }
381  }
382  G4cout << "G4CrossSectionDataStore::GetCrossSection ERROR: "
383  << " no isotope cross section found"
384  << G4endl;
385  G4cout << " for " << part->GetDefinition()->GetParticleName()
386  << " off Element " << elm->GetName()
387  << " in " << mat->GetName()
388  << " Z= " << Z << " A= " << A
389  << " E(MeV)= " << part->GetKineticEnergy()/MeV << G4endl;
390  throw G4HadronicException(__FILE__, __LINE__,
391  " no applicable data set found for the isotope");
392  return 0.0;
393 }
static const double MeV
Definition: G4SIunits.hh:211
int G4int
Definition: G4Types.hh:78
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
Float_t Z
std::vector< G4VCrossSectionDataSet * > dataSetList
const G4String & GetName() const
Definition: G4Element.hh:127
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
Definition: G4Material.hh:178
Here is the call graph for this function:

◆ GetCrossSection() [4/4]

G4double G4CrossSectionDataStore::GetCrossSection ( const G4DynamicParticle part,
const G4Material mat,
G4bool  requiresSlowPath 
)
private

Definition at line 82 of file G4CrossSectionDataStore.cc.

84 {
85  //The fast-path algorithm
86  // requiresSlowPath == true => Use slow path independently of other conditions
87  //A. Dotti: modifications to this algorithm following the studies of P. Ruth and R. Fowler
88  // on speeding up the cross-sections calculations. Their algorithm is called in the
89  // following "fast-path" while the normal approach to calcualte cross-section is
90  // referred to as "slow-path".
91  //Some details on the fast-path algorithm:
92  //The idea is to use a cached approximation of the material cross-section.
93  //Starting points:
94  //1- We need the material cross-section for navigation purposes: e.g. to calculate the PIL.
95  //2- If this interaction occurs at the end of a step we need to select the G4Element on which
96  // nucleus the interaction is actually happening. This is done calculating the element cross-section
97  // throwing a random number and selecting the appropriate nucleus (see SampleZandA function)
98  //3- To calculate the material cross section needed for 1- we use the G4Element cross sections.
99  // Since the material cross-section is simply the weighted sum of the element cross-sections.
100  //4- The slow path algorithm here accomplishes two use cases (not very good design IMHO): it
101  // calculates the material cross-section and it updates the xsecelem array that contains the
102  // relative cross-sections used by SampleZandA to select the element on which the interaction
103  // occurs.
104  //The idea of the fast-path algorithm is to replace the request for 1- with a faster calculation
105  //of the material cross-section from an approximation contained in cache.
106  //There two complications to take into account:
107  //A- If I use a fast parametrization for material cross-section, I still need to do the full
108  // calculations if an interaction occurs, because I need to calculate the element cross-sections.
109  // Since the function that updates xsecelem is the same (this one) I need to be sure that
110  // if I call this method for SampleAandZ the xsecelem is updated.
111  //B- It exists the possibility to be even fast the the fast-path algorithm: this happens when
112  // to select the element of the interaction via SampleZandI I call again this method exactly
113  // with the same conditions as when I called this method to calculate the material cross-section.
114  // In such a case xsecelem is updated and we do not need to do much more. This happens when
115  // for example a neutron undergoes an interaction at the end of the step.
116  //Dealing with B- complicates a bit the algorithm.
117  //In summary:
118  // If no fast-path algo is available (or user does not want that), go with the old plain algorithm
119  // If a fast-path algo is avilable, use it whenever possible.
120  //
121  // In general we expect user to selectively decide for which processes, materials and particle combinations
122  // we want to use the fast-path. If this is activated we also expect that the cross-section fast-path
123  // cache is created during the run initialization via calls to this method.
124  //
125  //fastPathFlags contains control flags for the fast-path algorithm:
126  // .prevCalcUsedFastPath == true => Previous call to GetCrossSection used the fast-path
127  // it is used in the decision to assess if xsecelem is
128  // correctly set-up
129  // .useFastPathIfAvailable == true => User requested the use of fast-path algorithm
130  // .initializationPhase == true => If true we are in Geant4 Init phase before the event-loop
131 
132  //Check user-request, does he want fast-path? if not
133  // OR
134  // we want fast-path and we are in initialization phase?
137  //Traditional algorithm is requested
138  requiresSlowPath=true;
139  }
140 
141  //Logging for performance calculations and counter, active only in FPDEBUG mode
143  //Measure number of cycles
145 
146  //This is the cache entry of the fast-path cross-section parametrization
148  //Did user request fast-path in first place and are we not in the initialization phase
150  //Important: if it is in initialization phase we should NOT use fast path: we are going to build it
151  //G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Key searchkey = {part->GetParticleDefinition(),mat};
152  entry = fastPathCache[{part->GetParticleDefinition(),mat}];
153  }
154 
155  //Super-fast-path: are we calling again this method for exactly the same conditions
156  //of the triplet {particle,material,energy}?
157  if(mat == currentMaterial && part->GetDefinition() == matParticle
158  && part->GetKineticEnergy() == matKinEnergy)
159  {
161  //If there is no user-request for the fast-path in first place?
162  //It means we built the xsecelem for sure, let's return immediately
164  return matCrossSection;
165  } else {
166  //Check that the last time we called this method we used the slow
167  //path: we need the data-member xsecelem to be setup correctly for the current
168  //interaction. This is ensured only if: we will do the slow path right now or we
169  //did it exactly for the same conditions of the last call.
170  if ( !fastPathFlags.prevCalcUsedFastPath && ! requiresSlowPath ) {
173  //Good everything is setup correctly, exit!
174  return matCrossSection;
175  } else {
176  //We need to follow the slow-path because
177  //xsecelem is not calculated correctly
178  requiresSlowPath = true;
179  }
180  }
181  }
182 
183  //Ok, now check if we have cached for this {particle,material,energy} the cross-section
184  //in this case let's return immediately, if we are not forced to take the slow path
185  //(e.g. as before if the xsecelem is not up-to-date we need to take the slow-path).
186  //Note that this is not equivalent to the previous ultra-fast check: we now have a map here
187  //So we can have for example a different particle.
188  if ( entry != nullptr && entry->energy == part->GetKineticEnergy() ) {
190  if ( !requiresSlowPath ) {
191  return entry->crossSection;
192  }
193  }
194 
196  matParticle = part->GetDefinition();
197  matKinEnergy = part->GetKineticEnergy();
198  matCrossSection = 0;
199 
200  //Now check if the cache entry has a fast-path cross-section calculation available
201  G4FastPathHadronicCrossSection::fastPathEntry* fast_entry = nullptr;
202  if ( entry != nullptr && ! requiresSlowPath ) {
203  fast_entry = entry->fastPath;
204  assert(fast_entry!=nullptr && !fastPathFlags.initializationPhase);
205  }
206 
207  //Each fast-path cross-section has a minimum value of validity, if energy is below
208  //that skip fast-path algorithm
209  if ( fast_entry != nullptr && part->GetKineticEnergy() < fast_entry->min_cutoff )
210  {
211  assert(requiresSlowPath==false);
212  requiresSlowPath = true;
213  }
214 
215  //Ready to use the fast-path calculation
216  if ( !requiresSlowPath && fast_entry != nullptr ) {
217  counters.FastPath();
218  //Retrieve cross-section from fast-path cache
219  matCrossSection = fast_entry->GetCrossSection(part->GetKineticEnergy());
221  } else {
222  counters.SlowPath();
223  //Remember that we are now doing the full calculation: xsecelem will
224  //be made valid
226 
227  G4int nElements = mat->GetNumberOfElements();
228  const G4double* nAtomsPerVolume = mat->GetVecNbOfAtomsPerVolume();
229 
230  if(G4int(xsecelm.size()) < nElements) { xsecelm.resize(nElements); }
231 
232  for(G4int i=0; i<nElements; ++i) {
233  matCrossSection += nAtomsPerVolume[i] *
234  GetCrossSection(part, (*mat->GetElementVector())[i], mat);
236  }
237  }
238  //Stop measurement of cpu cycles
240 
241  if ( entry != nullptr ) {
242  entry->energy = part->GetKineticEnergy();
243  entry->crossSection = matCrossSection;
244  }
245  //Some logging of timing
247  return matCrossSection;
248 }
G4FastPathHadronicCrossSection::controlFlag fastPathFlags
G4double GetCrossSection(G4double ene) const
static void logTiming(cycleCountEntry *, fastPathEntry *, timing &)
const G4ParticleDefinition * matParticle
G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Cache fastPathCache
G4double energy
const G4double min_cutoff
int G4int
Definition: G4Types.hh:78
G4FastPathHadronicCrossSection::getCrossSectionCount counters
Float_t mat
G4double GetKineticEnergy() const
static void logInvocationOneLine(cycleCountEntry *)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
fastPathEntry * fastPath
std::vector< G4double > xsecelm
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4double crossSection
G4ParticleDefinition * GetDefinition() const
static void logInvocationTriedOneLine(cycleCountEntry *)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
const G4ParticleDefinition * GetParticleDefinition() const
double G4double
Definition: G4Types.hh:76
G4FastPathHadronicCrossSection::timing timing
Here is the call graph for this function:

◆ GetFastPathControlFlags()

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

Definition at line 137 of file G4CrossSectionDataStore.hh.

137 { return fastPathFlags; }
G4FastPathHadronicCrossSection::controlFlag fastPathFlags
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetFastPathParameters()

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

Definition at line 135 of file G4CrossSectionDataStore.hh.

135 { return fastPathParams; }
G4FastPathHadronicCrossSection::fastPathParameters fastPathParams
Here is the caller graph for this function:

◆ GetIsoCrossSection()

G4double G4CrossSectionDataStore::GetIsoCrossSection ( const G4DynamicParticle part,
G4int  Z,
G4int  A,
const G4Isotope iso,
const G4Element elm,
const G4Material aMaterial,
G4int  index 
)
private

Definition at line 330 of file G4CrossSectionDataStore.cc.

336 {
337  // this methods is called after the check that dataSetList[idx]
338  // depend on isotopes, so for this DataSet only isotopes are checked
339 
340  // isotope-wise cross section does exist
341  if(dataSetList[idx]->IsIsoApplicable(part, Z, A, elm, mat) ) {
342  return dataSetList[idx]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
343 
344  } else {
345  // seach for other dataSet
346  for (G4int j = nDataSetList-1; j >= 0; --j) {
347  if (dataSetList[j]->IsElementApplicable(part, Z, mat)) {
348  return dataSetList[j]->GetElementCrossSection(part, Z, mat);
349  } else if (dataSetList[j]->IsIsoApplicable(part, Z, A, elm, mat)) {
350  return dataSetList[j]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
351  }
352  }
353  }
354  G4cout << "G4CrossSectionDataStore::GetCrossSection ERROR: "
355  << " no isotope cross section found"
356  << G4endl;
357  G4cout << " for " << part->GetDefinition()->GetParticleName()
358  << " off Element " << elm->GetName()
359  << " in " << mat->GetName()
360  << " Z= " << Z << " A= " << A
361  << " E(MeV)= " << part->GetKineticEnergy()/MeV << G4endl;
362  throw G4HadronicException(__FILE__, __LINE__,
363  " no applicable data set found for the isotope");
364  return 0.0;
365  //return dataSetList[idx]->ComputeCrossSection(part, elm, mat);
366 }
static const double MeV
Definition: G4SIunits.hh:211
int G4int
Definition: G4Types.hh:78
Float_t mat
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
Float_t Z
std::vector< G4VCrossSectionDataSet * > dataSetList
const G4String & GetName() const
Definition: G4Element.hh:127
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ HtmlFileName()

G4String G4CrossSectionDataStore::HtmlFileName ( const G4String in) const
private

Definition at line 631 of file G4CrossSectionDataStore.cc.

632 {
633  G4String str(in);
634  // replace blanks by _ C++11 version:
635 #ifdef G4USE_STD11
636  std::transform(str.begin(), str.end(), str.begin(), [](char ch) {
637  return ch == ' ' ? '_' : ch;
638  });
639 #else
640  // and now in ancient language
641  for(std::string::iterator it = str.begin(); it != str.end(); ++it) {
642  if(*it == ' ') *it = '_';
643  }
644 #endif
645  str=str + ".html";
646  return str;
647 }
Here is the caller graph for this function:

◆ operator=()

G4CrossSectionDataStore& G4CrossSectionDataStore::operator= ( const G4CrossSectionDataStore right)
private

◆ PrintCrossSectionHtml()

void G4CrossSectionDataStore::PrintCrossSectionHtml ( const G4VCrossSectionDataSet cs) const

Definition at line 608 of file G4CrossSectionDataStore.cc.

609 {
610  G4String dirName(getenv("G4PhysListDocDir"));
611  G4String physListName(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 }
G4String HtmlFileName(const G4String &in) const
virtual void CrossSectionDescription(std::ostream &) const
const G4String & GetName() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleZandA()

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

Definition at line 398 of file G4CrossSectionDataStore.cc.

401 {
403 
404  G4int nElements = mat->GetNumberOfElements();
405  const G4ElementVector* theElementVector = mat->GetElementVector();
406  G4Element* anElement = (*theElementVector)[0];
407 
408  G4double cross = GetCrossSection(part, mat , true);
409  // select element from a compound
410  if(1 < nElements) {
411  cross *= G4UniformRand();
412  for(G4int i=0; i<nElements; ++i) {
413  if(cross <= xsecelm[i]) {
414  anElement = (*theElementVector)[i];
415  break;
416  }
417  }
418  }
419 
420  G4int Z = G4lrint(anElement->GetZ());
421  G4Isotope* iso = 0;
422 
423  G4int i = nDataSetList-1;
424  if (dataSetList[i]->IsElementApplicable(part, Z, mat)) {
425 
426  //----------------------------------------------------------------
427  // element-wise cross section
428  // isotope cross section is not computed
429  //----------------------------------------------------------------
430  G4int nIso = anElement->GetNumberOfIsotopes();
431  if (0 >= nIso) {
432  G4cout << " Element " << anElement->GetName() << " Z= " << Z
433  << " has no isotopes " << G4endl;
434  throw G4HadronicException(__FILE__, __LINE__,
435  " Isotope vector is not defined");
436  return anElement;
437  }
438  // isotope abundances
439  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
440  iso = (*isoVector)[0];
441 
442  // more than 1 isotope
443  if(1 < nIso) {
444  iso = dataSetList[i]->SelectIsotope(anElement, part->GetKineticEnergy());
445  }
446 
447  } else {
448 
449  //----------------------------------------------------------------
450  // isotope-wise cross section
451  // isotope cross section is computed
452  //----------------------------------------------------------------
453  G4int nIso = anElement->GetNumberOfIsotopes();
454  cross = 0.0;
455 
456  if (0 >= nIso) {
457  G4cout << " Element " << anElement->GetName() << " Z= " << Z
458  << " has no isotopes " << G4endl;
459  throw G4HadronicException(__FILE__, __LINE__,
460  " Isotope vector is not defined");
461  return anElement;
462  }
463 
464  // user-defined isotope abundances
465  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
466  iso = (*isoVector)[0];
467 
468  // more than 1 isotope
469  if(1 < nIso) {
470  G4double* abundVector = anElement->GetRelativeAbundanceVector();
471  if(G4int(xseciso.size()) < nIso) { xseciso.resize(nIso); }
472 
473  for (G4int j = 0; j<nIso; ++j) {
474  G4double xsec = 0.0;
475  if(abundVector[j] > 0.0) {
476  iso = (*isoVector)[j];
477  xsec = abundVector[j]*
478  GetIsoCrossSection(part, Z, iso->GetN(), iso, anElement, mat, i);
479  }
480  cross += xsec;
481  xseciso[j] = cross;
482  }
483  cross *= G4UniformRand();
484  for (G4int j = 0; j<nIso; ++j) {
485  if(cross <= xseciso[j]) {
486  iso = (*isoVector)[j];
487  break;
488  }
489  }
490  }
491  }
492  target.SetIsotope(iso);
493  return anElement;
494 }
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
std::vector< G4Isotope * > G4IsotopeVector
std::vector< G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:78
G4FastPathHadronicCrossSection::getCrossSectionCount counters
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
Float_t mat
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
std::vector< G4double > xseciso
Float_t Z
std::vector< G4VCrossSectionDataSet * > dataSetList
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *, const G4Element *, const G4Material *aMaterial, G4int index)
const G4String & GetName() const
Definition: G4Element.hh:127
int G4lrint(double ad)
Definition: templates.hh:163
std::vector< G4double > xsecelm
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
#define G4endl
Definition: G4ios.hh:61
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
double G4double
Definition: G4Types.hh:76
G4double GetZ() const
Definition: G4Element.hh:131
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetVerboseLevel()

void G4CrossSectionDataStore::SetVerboseLevel ( G4int  value)
inline

Definition at line 166 of file G4CrossSectionDataStore.hh.

167 {
168  verboseLevel = value;
169 }

Friends And Related Function Documentation

◆ G4FastPathHadronicCrossSection::fastPathEntry

Definition at line 141 of file G4CrossSectionDataStore.hh.

Member Data Documentation

◆ counters

G4FastPathHadronicCrossSection::getCrossSectionCount G4CrossSectionDataStore::counters
private

Definition at line 148 of file G4CrossSectionDataStore.hh.

◆ currentElement

const G4Element* G4CrossSectionDataStore::currentElement
private

Definition at line 125 of file G4CrossSectionDataStore.hh.

◆ currentMaterial

const G4Material* G4CrossSectionDataStore::currentMaterial
private

Definition at line 119 of file G4CrossSectionDataStore.hh.

◆ dataSetList

std::vector<G4VCrossSectionDataSet*> G4CrossSectionDataStore::dataSetList
private

Definition at line 115 of file G4CrossSectionDataStore.hh.

◆ elmCrossSection

G4double G4CrossSectionDataStore::elmCrossSection
private

Definition at line 128 of file G4CrossSectionDataStore.hh.

◆ elmKinEnergy

G4double G4CrossSectionDataStore::elmKinEnergy
private

Definition at line 127 of file G4CrossSectionDataStore.hh.

◆ elmMaterial

const G4Material* G4CrossSectionDataStore::elmMaterial
private

Definition at line 124 of file G4CrossSectionDataStore.hh.

◆ elmParticle

const G4ParticleDefinition* G4CrossSectionDataStore::elmParticle
private

Definition at line 126 of file G4CrossSectionDataStore.hh.

◆ fastPathCache

G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Cache G4CrossSectionDataStore::fastPathCache
private

Definition at line 150 of file G4CrossSectionDataStore.hh.

◆ fastPathFlags

G4FastPathHadronicCrossSection::controlFlag G4CrossSectionDataStore::fastPathFlags
private

Definition at line 145 of file G4CrossSectionDataStore.hh.

◆ fastPathParams

G4FastPathHadronicCrossSection::fastPathParameters G4CrossSectionDataStore::fastPathParams
private

Definition at line 146 of file G4CrossSectionDataStore.hh.

◆ matCrossSection

G4double G4CrossSectionDataStore::matCrossSection
private

Definition at line 122 of file G4CrossSectionDataStore.hh.

◆ matKinEnergy

G4double G4CrossSectionDataStore::matKinEnergy
private

Definition at line 121 of file G4CrossSectionDataStore.hh.

◆ matParticle

const G4ParticleDefinition* G4CrossSectionDataStore::matParticle
private

Definition at line 120 of file G4CrossSectionDataStore.hh.

◆ nDataSetList

G4int G4CrossSectionDataStore::nDataSetList
private

Definition at line 130 of file G4CrossSectionDataStore.hh.

◆ nist

G4NistManager* G4CrossSectionDataStore::nist
private

Definition at line 113 of file G4CrossSectionDataStore.hh.

◆ requests

G4FastPathHadronicCrossSection::G4CrossSectionDataStore_Requests G4CrossSectionDataStore::requests
private

Definition at line 152 of file G4CrossSectionDataStore.hh.

◆ timing

G4FastPathHadronicCrossSection::timing G4CrossSectionDataStore::timing
private

Definition at line 151 of file G4CrossSectionDataStore.hh.

◆ verboseLevel

G4int G4CrossSectionDataStore::verboseLevel
private

Definition at line 131 of file G4CrossSectionDataStore.hh.

◆ xsecelm

std::vector<G4double> G4CrossSectionDataStore::xsecelm
private

Definition at line 116 of file G4CrossSectionDataStore.hh.

◆ xseciso

std::vector<G4double> G4CrossSectionDataStore::xseciso
private

Definition at line 117 of file G4CrossSectionDataStore.hh.


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