Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CrossSectionDataStore.hh
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: $
27 //
28 // -------------------------------------------------------------------
29 // File name: G4CrossSectionDataStore
30 //
31 // Modifications:
32 // 23.01.2009 V.Ivanchenko move constructor and destructor to source,
33 // use STL vector instead of C-array
34 //
35 // August 2011 Re-designed
36 // by G. Folger, V. Ivantchenko, T. Koi and D.H. Wright
37 
38 // Class Description
39 // This is the class to which cross section data sets may be registered.
40 // An instance of it is contained in each hadronic process, allowing
41 // the use of the AddDataSet() method to tailor the cross sections to
42 // your application.
43 // Class Description - End
44 
45 #ifndef G4CrossSectionDataStore_h
46 #define G4CrossSectionDataStore_h 1
47 
48 #include "globals.hh"
51 #include "G4DynamicParticle.hh"
52 #include "G4PhysicsVector.hh"
53 #include <vector>
54 #include <iostream>
55 
56 class G4Nucleus;
58 class G4Isotope;
59 class G4Element;
60 class G4Material;
61 class G4NistManager;
62 
64 {
65 public:
66 
68 
70 
71  // Cross section per unit volume is computed (inverse mean free path)
73 
74  // Cross section per element is computed
76  const G4Element*, const G4Material*);
77 
78  // Cross section per isotope is computed
80  const G4Isotope*,
81  const G4Element*, const G4Material*);
82 
83  // Sample Z and A of a target nucleus and upload into G4Nucleus
85  G4Nucleus& target);
86 
87  // Initialisation before run
89 
90  // Dump store to G4cout
92 
93  // Dump store as html
94  void DumpHtml(const G4ParticleDefinition&, std::ofstream&) const;
95  void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs) const;
96 
97  inline void AddDataSet(G4VCrossSectionDataSet*);
98 
99  inline void SetVerboseLevel(G4int value);
100 
101 private:
102 
103  G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
104  const G4Isotope*,
105  const G4Element*, const G4Material* aMaterial,
106  G4int index);
107 
108  G4CrossSectionDataStore & operator=(const G4CrossSectionDataStore &right);
110 
111  G4String HtmlFileName(const G4String & in) const;
112 
113  G4NistManager* nist;
114 
115  std::vector<G4VCrossSectionDataSet*> dataSetList;
116  std::vector<G4double> xsecelm;
117  std::vector<G4double> xseciso;
118 
119  const G4Material* currentMaterial;
120  const G4ParticleDefinition* matParticle;
121  G4double matKinEnergy;
122  G4double matCrossSection;
123 
124  const G4Material* elmMaterial;
125  const G4Element* currentElement;
126  const G4ParticleDefinition* elmParticle;
127  G4double elmKinEnergy;
128  G4double elmCrossSection;
129 
130  G4int nDataSetList;
131  G4int verboseLevel;
132  //Fast path: caching
133 public:
135  GetFastPathParameters() const { return fastPathParams; }
137  GetFastPathControlFlags() const { return fastPathFlags; }
138  void DumpFastPath( const G4ParticleDefinition* , const G4Material* , std::ostream& os);
140 private:
142  //The following method is called by the public one GetCrossSection(const G4DynamicParticle*, const G4Material*)
143  //The third parameter is used to force the calculation of cross-sections skipping the fast-path mechanism
144  G4double GetCrossSection(const G4DynamicParticle*, const G4Material*, G4bool requiresSlowPath);
147  //Counters
149  //TODO: share this among threads
153 };
154 
156  //By default tries to use the fast-path mechanism
157  return GetCrossSection( particle , material , false);
158 }
159 
161 {
162  dataSetList.push_back(p);
163  ++nDataSetList;
164 }
165 
167 {
168  verboseLevel = value;
169 }
170 
171 #endif
const G4FastPathHadronicCrossSection::controlFlag & GetFastPathControlFlags() const
void AddDataSet(G4VCrossSectionDataSet *)
const XML_Char * target
Definition: expat.h:268
const char * p
Definition: xmltok.h:285
const G4ParticleDefinition *const particle
int G4int
Definition: G4Types.hh:78
const G4Material *const material
const G4FastPathHadronicCrossSection::fastPathParameters & GetFastPathParameters() const
std::set< fastPathRequestConfig_t, fastPathRequestConfig_Less > G4CrossSectionDataStore_Requests
void ActivateFastPath(const G4ParticleDefinition *, const G4Material *, G4double)
std::unordered_map< G4CrossSectionDataStore_Key, cycleCountEntry *, G4CrossSectionDataStore_Key_Hash, G4CrossSectionDataStore_Key_EqualTo > G4CrossSectionDataStore_Cache
double A(double temperature)
const XML_Char int const XML_Char * value
Definition: expat.h:331
bool G4bool
Definition: G4Types.hh:79
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
void DumpFastPath(const G4ParticleDefinition *, const G4Material *, std::ostream &os)
void BuildPhysicsTable(const G4ParticleDefinition &)
void DumpPhysicsTable(const G4ParticleDefinition &)
void DumpHtml(const G4ParticleDefinition &, std::ofstream &) const
void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs) const
G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
double G4double
Definition: G4Types.hh:76