Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// GEANT4 tag $Name: not supported by cvs2svn $
28//
29// -------------------------------------------------------------------
30// File name: G4CrossSectionDataStore
31//
32// Modifications:
33// 23.01.2009 V.Ivanchenko move constructor and destructor to source,
34// use STL vector instead of C-array
35//
36// August 2011 Re-designed
37// by G. Folger, V. Ivantchenko, T. Koi and D.H. Wright
38
39// Class Description
40// This is the class to which cross section data sets may be registered.
41// An instance of it is contained in each hadronic process, allowing
42// the use of the AddDataSet() method to tailor the cross sections to
43// your application.
44// Class Description - End
45
46#ifndef G4CrossSectionDataStore_h
47#define G4CrossSectionDataStore_h 1
48
49#include "globals.hh"
51#include <vector>
52
53class G4Nucleus;
56class G4Isotope;
57class G4Element;
58class G4Material;
59class G4NistManager;
60
62{
63public:
64
66
68
69 // Cross section per unit volume is computed (inverse mean free path)
71
72 // Cross section per element is computed
74 const G4Element*, const G4Material*);
75
76 // Cross section per isotope is computed
78 const G4Isotope*,
79 const G4Element*, const G4Material*);
80
81 // Sample Z and A of a target nucleus and upload into G4Nucleus
83 G4Nucleus& target);
84
85 // Initialisation before run
87
88 // Dump store to G4cout
90
91 // Dump store as html
92 void DumpHtml(const G4ParticleDefinition&, std::ofstream&);
93
95
96 inline void SetVerboseLevel(G4int value);
97
98private:
99
100 G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
101 const G4Isotope*,
102 const G4Element*, const G4Material* aMaterial,
103 G4int index);
104
105 G4CrossSectionDataStore & operator=(const G4CrossSectionDataStore &right);
107
108 G4NistManager* nist;
109
110 std::vector<G4VCrossSectionDataSet*> dataSetList;
111 std::vector<G4double> xsecelm;
112 std::vector<G4double> xseciso;
113
114 const G4Material* currentMaterial;
115 const G4ParticleDefinition* matParticle;
116 G4double matKinEnergy;
117 G4double matCrossSection;
118
119 const G4Material* elmMaterial;
120 const G4Element* currentElement;
121 const G4ParticleDefinition* elmParticle;
122 G4double elmKinEnergy;
123 G4double elmCrossSection;
124
125 G4int nDataSetList;
126 G4int verboseLevel;
127};
128
130{
131 dataSetList.push_back(p);
132 ++nDataSetList;
133}
134
136{
137 verboseLevel = value;
138}
139
140#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void BuildPhysicsTable(const G4ParticleDefinition &)
void AddDataSet(G4VCrossSectionDataSet *)
void DumpPhysicsTable(const G4ParticleDefinition &)
G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
void DumpHtml(const G4ParticleDefinition &, std::ofstream &)