Geant4 11.1.1
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//
27// -------------------------------------------------------------------
28// File name: G4CrossSectionDataStore
29//
30// Modifications:
31// 23.01.2009 V.Ivanchenko move constructor and destructor to source,
32// use STL vector instead of C-array
33//
34// August 2011 Re-designed
35// by G. Folger, V. Ivantchenko, T. Koi and D.H. Wright
36
37// Class Description
38// This is the class to which cross section data sets may be registered.
39// An instance of it is contained in each hadronic process, allowing
40// the use of the AddDataSet() method to tailor the cross sections to
41// your application.
42// Class Description - End
43
44#ifndef G4CrossSectionDataStore_h
45#define G4CrossSectionDataStore_h 1
46
47#include "globals.hh"
49#include "G4DynamicParticle.hh"
50#include "G4PhysicsVector.hh"
51#include <vector>
52#include <iostream>
53
54class G4Nucleus;
56class G4Isotope;
57class G4Element;
58class G4Material;
59class G4NistManager;
60
62{
63public:
64
66
68
69 // Run time cross section per unit volume
72
73 // Cross section per element
75 const G4Element*, const G4Material*);
76
77 // Cross section per isotope
79 const G4Isotope*,
80 const G4Element*, const G4Material*);
81
82 // Sample Z and A of a target nucleus and upload into G4Nucleus
84 G4Nucleus& target);
85
86 // Initialisation before run
88
89 // Dump store to G4cout
91
92 // Dump store as html
93 void DumpHtml(const G4ParticleDefinition&, std::ofstream&) const;
95
97 void AddDataSet(G4VCrossSectionDataSet*, std::size_t);
98 inline const std::vector<G4VCrossSectionDataSet*>& GetDataSetList() const;
99
100 inline void SetVerboseLevel(G4int value);
101
102 // may be used by special processes
103 inline void SetForcedElement(const G4Element*);
104
106 (const G4CrossSectionDataStore &right) = delete;
108
109private:
110
111 G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
112 const G4Isotope*, const G4Element*,
113 const G4Material*, const G4int index);
114
115 G4String HtmlFileName(const G4String & in) const;
116
117 G4NistManager* nist;
118 const G4Material* currentMaterial = nullptr;
119 const G4ParticleDefinition* matParticle = nullptr;
120 const G4Element* forcedElement = nullptr;
121 G4double matKinEnergy = 0.0;
122 G4double matCrossSection = 0.0;
123
124 G4int nDataSetList = 0;
125 G4int verboseLevel = 1;
126
127 std::vector<G4VCrossSectionDataSet*> dataSetList;
128 std::vector<G4double> xsecelm;
129 std::vector<G4double> xseciso;
130};
131
133{
134 verboseLevel = value;
135}
136
138{
139 forcedElement = ptr;
140}
141
142inline const std::vector<G4VCrossSectionDataSet*>&
144{
145 return dataSetList;
146}
147
148inline G4double
150 const G4Material* mat)
151{
152 if(dp->GetKineticEnergy() != matKinEnergy || mat != currentMaterial ||
153 dp->GetDefinition() != matParticle) {
154 ComputeCrossSection(dp, mat);
155 }
156 return matCrossSection;
157}
158
159#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
void DumpHtml(const G4ParticleDefinition &, std::ofstream &) const
~G4CrossSectionDataStore()=default
const std::vector< G4VCrossSectionDataSet * > & GetDataSetList() const
void BuildPhysicsTable(const G4ParticleDefinition &)
void AddDataSet(G4VCrossSectionDataSet *)
G4CrossSectionDataStore(const G4CrossSectionDataStore &)=delete
void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs) const
void DumpPhysicsTable(const G4ParticleDefinition &)
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
void SetForcedElement(const G4Element *)
const G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const