Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4CrossSectionDataStore.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4CrossSectionDataStore.cc 98827 2016-08-12 12:37:39Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4CrossSectionDataStore
34 //
35 // Modifications:
36 // 23.01.2009 V.Ivanchenko add destruction of data sets
37 // 29.04.2010 G.Folger modifictaions for integer A & Z
38 // 14.03.2011 V.Ivanchenko fixed DumpPhysicsTable
39 // 15.08.2011 G.Folger, V.Ivanchenko, T.Koi, D.Wright redesign the class
40 // 07.03.2013 M.Maire cosmetic in DumpPhysicsTable
41 //
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
44 
46 #include "G4SystemOfUnits.hh"
47 #include "G4UnitsTable.hh"
48 #include "G4HadronicException.hh"
49 #include "G4HadTmpUtil.hh"
50 #include "Randomize.hh"
51 #include "G4Nucleus.hh"
52 
53 #include "G4DynamicParticle.hh"
54 #include "G4Isotope.hh"
55 #include "G4Element.hh"
56 #include "G4Material.hh"
57 #include "G4NistManager.hh"
58 #include <algorithm>
59 
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
62 
64  nDataSetList(0), verboseLevel(0),fastPathFlags(),fastPathParams(),
65  counters(),fastPathCache()
66 {
67  nist = G4NistManager::Instance();
68  currentMaterial = elmMaterial = 0;
69  currentElement = 0; //ALB 14-Aug-2012 Coverity fix.
70  matParticle = elmParticle = 0;
71  matKinEnergy = elmKinEnergy = matCrossSection = elmCrossSection = 0.0;
72 }
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
75 
77 {}
78 
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
83  const G4Material* mat , G4bool requiresSlowPath)
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?
135  if ( !fastPathFlags.useFastPathIfAvailable
136  || (fastPathFlags.useFastPathIfAvailable&&fastPathFlags.initializationPhase) ) {
137  //Traditional algorithm is requested
138  requiresSlowPath=true;
139  }
140 
141  //Logging for performance calculations and counter, active only in FPDEBUG mode
142  counters.MethodCalled();
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
149  if ( fastPathFlags.useFastPathIfAvailable && !fastPathFlags.initializationPhase ) {
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
163  if ( !fastPathFlags.useFastPathIfAvailable ) {
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 ) {
171  counters.HitOneLine();
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 
195  currentMaterial = mat;
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());
220  fastPathFlags.prevCalcUsedFastPath=true;
221  } else {
222  counters.SlowPath();
223  //Remember that we are now doing the full calculation: xsecelem will
224  //be made valid
225  fastPathFlags.prevCalcUsedFastPath=false;
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);
235  xsecelm[i] = matCrossSection;
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
246  G4FastPathHadronicCrossSection::logTiming(entry,fast_entry,timing);
247  return matCrossSection;
248 }
249 
250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
251 
252 void
254 {
255  const G4FastPathHadronicCrossSection::cycleCountEntry* entry = fastPathCache[{pd,mat}];
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 }
268 
269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
270 
271 G4double
273  const G4Element* elm,
274  const G4Material* mat)
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 }
326 
327 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
328 
329 G4double
330 G4CrossSectionDataStore::GetIsoCrossSection(const G4DynamicParticle* part,
331  G4int Z, G4int A,
332  const G4Isotope* iso,
333  const G4Element* elm,
334  const G4Material* mat,
335  G4int idx)
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 }
365 
366 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
367 
368 G4double
370  G4int Z, G4int A,
371  const G4Isotope* iso,
372  const G4Element* elm,
373  const G4Material* mat)
374 {
375  for (G4int i = nDataSetList-1; i >= 0; --i) {
376  if (dataSetList[i]->IsIsoApplicable(part, Z, A, elm, mat) ) {
377  return dataSetList[i]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
378  }
379  }
380  G4cout << "G4CrossSectionDataStore::GetCrossSection ERROR: "
381  << " no isotope cross section found"
382  << G4endl;
383  G4cout << " for " << part->GetDefinition()->GetParticleName()
384  << " off Element " << elm->GetName()
385  << " in " << mat->GetName()
386  << " Z= " << Z << " A= " << A
387  << " E(MeV)= " << part->GetKineticEnergy()/MeV << G4endl;
388  throw G4HadronicException(__FILE__, __LINE__,
389  " no applicable data set found for the isotope");
390 }
391 
392 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
393 
394 G4Element*
396  const G4Material* mat,
397  G4Nucleus& target)
398 {
399  counters.SampleZandA();
400 
401  G4int nElements = mat->GetNumberOfElements();
402  const G4ElementVector* theElementVector = mat->GetElementVector();
403  G4Element* anElement = (*theElementVector)[0];
404 
405  G4double cross = GetCrossSection(part, mat , true);
406  // select element from a compound
407  if(1 < nElements) {
408  cross *= G4UniformRand();
409  for(G4int i=0; i<nElements; ++i) {
410  if(cross <= xsecelm[i]) {
411  anElement = (*theElementVector)[i];
412  break;
413  }
414  }
415  }
416 
417  G4int Z = G4lrint(anElement->GetZ());
418  G4Isotope* iso = 0;
419 
420  G4int i = nDataSetList-1;
421  if (dataSetList[i]->IsElementApplicable(part, Z, mat)) {
422 
423  //----------------------------------------------------------------
424  // element-wise cross section
425  // isotope cross section is not computed
426  //----------------------------------------------------------------
427  G4int nIso = anElement->GetNumberOfIsotopes();
428  if (0 >= nIso) {
429  G4cout << " Element " << anElement->GetName() << " Z= " << Z
430  << " has no isotopes " << G4endl;
431  throw G4HadronicException(__FILE__, __LINE__,
432  " Isotope vector is not defined");
433  return anElement;
434  }
435  // isotope abundances
436  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
437  iso = (*isoVector)[0];
438 
439  // more than 1 isotope
440  if(1 < nIso) {
441  iso = dataSetList[i]->SelectIsotope(anElement, part->GetKineticEnergy());
442  }
443 
444  } else {
445 
446  //----------------------------------------------------------------
447  // isotope-wise cross section
448  // isotope cross section is computed
449  //----------------------------------------------------------------
450  G4int nIso = anElement->GetNumberOfIsotopes();
451  cross = 0.0;
452 
453  if (0 >= nIso) {
454  G4cout << " Element " << anElement->GetName() << " Z= " << Z
455  << " has no isotopes " << G4endl;
456  throw G4HadronicException(__FILE__, __LINE__,
457  " Isotope vector is not defined");
458  return anElement;
459  }
460 
461  // user-defined isotope abundances
462  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
463  iso = (*isoVector)[0];
464 
465  // more than 1 isotope
466  if(1 < nIso) {
467  G4double* abundVector = anElement->GetRelativeAbundanceVector();
468  if(G4int(xseciso.size()) < nIso) { xseciso.resize(nIso); }
469 
470  for (G4int j = 0; j<nIso; ++j) {
471  G4double xsec = 0.0;
472  if(abundVector[j] > 0.0) {
473  iso = (*isoVector)[j];
474  xsec = abundVector[j]*
475  GetIsoCrossSection(part, Z, iso->GetN(), iso, anElement, mat, i);
476  }
477  cross += xsec;
478  xseciso[j] = cross;
479  }
480  cross *= G4UniformRand();
481  for (G4int j = 0; j<nIso; ++j) {
482  if(cross <= xseciso[j]) {
483  iso = (*isoVector)[j];
484  break;
485  }
486  }
487  }
488  }
489  target.SetIsotope(iso);
490  return anElement;
491 }
492 
493 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
494 
495 void
497 {
498  if (nDataSetList == 0)
499  {
500  throw G4HadronicException(__FILE__, __LINE__,
501  "G4CrossSectionDataStore: no data sets registered");
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 ) {
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 }
528 
529 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
530 
532 {
533  assert(pdef!=nullptr&&mat!=nullptr);
535  if ( requests.insert( { key , min_cutoff } ).second ) {
536  std::ostringstream msg;
537  msg<<"Attempting to request FastPath for couple: "<<pdef->GetParticleName()<<","<<mat->GetName();
538  msg<<" but combination already exists";
539  throw G4HadronicException(__FILE__,__LINE__,msg.str());
540  }
541 }
542 
543 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
544 
545 void
547 {
548  // Print out all cross section data sets used and the energies at
549  // which they apply
550 
551  if (nDataSetList == 0) {
552  G4cout << "WARNING - G4CrossSectionDataStore::DumpPhysicsTable: "
553  << " no data sets registered" << G4endl;
554  return;
555  }
556 
557  for (G4int i = nDataSetList-1; i >= 0; --i) {
558  G4double e1 = dataSetList[i]->GetMinKinEnergy();
559  G4double e2 = dataSetList[i]->GetMaxKinEnergy();
560  G4cout
561  << " Cr_sctns: " << std::setw(25) << dataSetList[i]->GetName() << ": "
562  << G4BestUnit(e1, "Energy")
563  << " ---> "
564  << G4BestUnit(e2, "Energy") << "\n";
565  if (dataSetList[i]->GetName() == "G4CrossSectionPairGG") {
566  dataSetList[i]->DumpPhysicsTable(aParticleType);
567  }
568  }
569 }
570 
571 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
572 #include <typeinfo>
574  std::ofstream& outFile) const
575 {
576  // Write cross section data set info to html physics list
577  // documentation page
578 
579  G4double ehi = 0;
580  G4double elo = 0;
581  G4String physListName(getenv("G4PhysListName"));
582  for (G4int i = nDataSetList-1; i > 0; i--) {
583  elo = dataSetList[i]->GetMinKinEnergy()/GeV;
584  ehi = dataSetList[i]->GetMaxKinEnergy()/GeV;
585  outFile << " <li><b><a href=\"" << physListName << "_"
586  << dataSetList[i]->GetName() << ".html\"> "
587  << dataSetList[i]->GetName() << "</a> from "
588  << elo << " GeV to " << ehi << " GeV </b></li>\n";
589  //G4cerr << i << ": XS for " << pD.GetParticleName() << " : " << dataSetList[i]->GetName()
590  // << " typeid : " << typeid(dataSetList[i]).name()<< G4endl;
591  PrintCrossSectionHtml(dataSetList[i]);
592  }
593 
594  G4double defaultHi = dataSetList[0]->GetMaxKinEnergy()/GeV;
595  if (ehi < defaultHi) {
596  outFile << " <li><b><a href=\"" << dataSetList[0]->GetName() << ".html\"> "
597  << dataSetList[0]->GetName() << "</a> from "
598  << ehi << " GeV to " << defaultHi << " GeV </b></li>\n";
599  PrintCrossSectionHtml(dataSetList[0]);
600  }
601 }
602 
603 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
604 
606 {
607  G4String dirName(getenv("G4PhysListDocDir"));
608  G4String physListName(getenv("G4PhysListName"));
609 
610  G4String pathName = dirName + "/" + physListName + "_" + HtmlFileName(cs->GetName());
611  std::ofstream outCS;
612  outCS.open(pathName);
613  outCS << "<html>\n";
614  outCS << "<head>\n";
615  outCS << "<title>Description of " << cs->GetName()
616  << "</title>\n";
617  outCS << "</head>\n";
618  outCS << "<body>\n";
619 
620  cs->CrossSectionDescription(outCS);
621 
622  outCS << "</body>\n";
623  outCS << "</html>\n";
624 
625 }
626 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
627 
628 G4String G4CrossSectionDataStore::HtmlFileName(const G4String & in) const
629 {
630  G4String str(in);
631  // replace blanks by _ C++11 version:
632 #ifdef G4USE_STD11
633  std::transform(str.begin(), str.end(), str.begin(), [](char ch) {
634  return ch == ' ' ? '_' : ch;
635  });
636 #else
637  // and now in ancient language
638  for(std::string::iterator it = str.begin(); it != str.end(); ++it) {
639  if(*it == ' ') *it = '_';
640  }
641 #endif
642  str=str + ".html";
643  return str;
644 }
645 
646 
647 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
G4bool GetNaturalAbundanceFlag() const
Definition: G4Element.hh:262
void Initialize(G4CrossSectionDataStore *)
std::vector< G4Isotope * > G4IsotopeVector
const XML_Char * target
Definition: expat.h:268
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
virtual void CrossSectionDescription(std::ostream &) const
static void logTiming(cycleCountEntry *, fastPathEntry *, timing &)
const G4String & GetName() const
Definition: G4Material.hh:178
G4double GetZ() const
Definition: G4Element.hh:131
G4ParticleDefinition * GetDefinition() const
static constexpr double second
Definition: G4SIunits.hh:157
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
const G4String & GetName() const
G4double energy
G4double GetCrossSection(G4double ene) const
const G4double min_cutoff
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
const G4String & GetParticleName() const
void ActivateFastPath(const G4ParticleDefinition *, const G4Material *, G4double)
std::pair< const G4ParticleDefinition *, const G4Material * > G4CrossSectionDataStore_Key
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
static void logInvocationOneLine(cycleCountEntry *)
bool G4bool
Definition: G4Types.hh:79
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
const G4ParticleDefinition * GetParticleDefinition() const
void DumpFastPath(const G4ParticleDefinition *, const G4Material *, std::ostream &os)
fastPathEntry * fastPath
void BuildPhysicsTable(const G4ParticleDefinition &)
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:163
void DumpPhysicsTable(const G4ParticleDefinition &)
void DumpHtml(const G4ParticleDefinition &, std::ofstream &) const
int G4lrint(double ad)
Definition: templates.hh:163
G4double crossSection
static constexpr double GeV
Definition: G4SIunits.hh:217
void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs) const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
static void logInvocationTriedOneLine(cycleCountEntry *)
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127