Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronicProcessStore.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//
29// GEANT4 Class header file
30//
31//
32// File name: G4HadronicProcessStore
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 09.05.2008
37//
38// Modifications:
39//
40//
41// Class Description:
42//
43
44// -------------------------------------------------------------------
45//
46
47#ifndef G4HadronicProcessStore_h
48#define G4HadronicProcessStore_h 1
49
50
51#include "globals.hh"
52#include "G4DynamicParticle.hh"
53#include "G4ThreeVector.hh"
54#include "G4HadronicProcess.hh"
59#include <map>
60#include <vector>
61#include <iostream>
62
63class G4Element;
66
67class G4HadronicProcessStore
68{
69
70friend class G4ThreadLocalSingleton<G4HadronicProcessStore>;
71
72public:
73
74 static G4HadronicProcessStore* Instance();
75
77
79 const G4ParticleDefinition* particle,
80 G4double kineticEnergy,
81 const G4VProcess* process,
82 const G4Element* element,
83 const G4Material* material=nullptr);
84
86 const G4ParticleDefinition* particle,
87 G4double kineticEnergy,
88 const G4VProcess* process,
89 const G4Material* material);
90
92 const G4ParticleDefinition *aParticle,
93 G4double kineticEnergy,
94 const G4Material *material);
95
97 const G4ParticleDefinition *aParticle,
98 G4double kineticEnergy,
99 const G4Element *anElement,
100 const G4Material* mat=nullptr);
101
103 const G4ParticleDefinition *aParticle,
104 G4double kineticEnergy,
105 G4int Z, G4int A);
106
108 const G4ParticleDefinition *aParticle,
109 G4double kineticEnergy,
110 const G4Material *material);
111
113 const G4ParticleDefinition *aParticle,
114 G4double kineticEnergy,
115 const G4Element *anElement, const G4Material* mat=0);
116
118 const G4ParticleDefinition *aParticle,
119 G4double kineticEnergy,
120 G4int Z, G4int A);
121
123 const G4ParticleDefinition *aParticle,
124 G4double kineticEnergy,
125 const G4Material *material);
126
128 const G4ParticleDefinition *aParticle,
129 G4double kineticEnergy,
130 const G4Element *anElement,
131 const G4Material* mat=nullptr);
132
134 const G4ParticleDefinition *aParticle,
135 G4double kineticEnergy,
136 G4int Z, G4int A);
137
139 const G4ParticleDefinition *aParticle,
140 G4double kineticEnergy,
141 const G4Material *material);
142
144 const G4ParticleDefinition *aParticle,
145 G4double kineticEnergy,
146 const G4Element *anElement,
147 const G4Material* mat=nullptr);
148
150 const G4ParticleDefinition *aParticle,
151 G4double kineticEnergy,
152 G4int Z, G4int A);
153
155 const G4ParticleDefinition *aParticle,
156 G4double kineticEnergy,
157 const G4Material *material);
158
160 const G4ParticleDefinition *aParticle,
161 G4double kineticEnergy,
162 const G4Element *anElement,
163 const G4Material* mat=nullptr);
164
166 const G4ParticleDefinition *aParticle,
167 G4double kineticEnergy,
168 G4int Z, G4int A);
169
170 // register/deregister processes following G4HadronicProcess interface
172
174 const G4ParticleDefinition*);
175
178
180
181 // register/deregister processes following only G4VProcess interface
183
185 const G4ParticleDefinition*);
186
188
189 void SetBuildXSTable(G4bool val);
190
191 G4bool GetBuildXSTable() const;
192
193 void PrintInfo(const G4ParticleDefinition*);
194
195 void Dump(G4int level);
196 void DumpHtml();
197 void PrintHtml(const G4ParticleDefinition*, std::ofstream&);
198 void PrintModelHtml(const G4HadronicInteraction * model) const;
199
200 void SetVerbose(G4int val);
202 // these methods are obsolete and will be removed
203
205 G4HadronicProcessType subType);
206
207 // Energy-momentum non-conservation limits and reporting
208 void SetEpReportLevel(G4int level);
209
210 void SetProcessAbsLevel(G4double absoluteLevel);
211
212 void SetProcessRelLevel(G4double relativeLevel);
213
214private:
215
216 // constructor
217 G4HadronicProcessStore();
218
219 // print process info
220 void Print(G4int idxProcess, G4int idxParticle);
221
222 G4String HtmlFileName(const G4String &) const;
223
224 static G4ThreadLocal G4HadronicProcessStore* instance;
225
226 typedef const G4ParticleDefinition* PD;
227 typedef G4HadronicProcess* HP;
228 typedef G4HadronicInteraction* HI;
229
230 // hadronic processes following G4HadronicProcess interface
231 std::vector<G4HadronicProcess*> process;
232 std::vector<G4HadronicInteraction*> model;
233 std::vector<G4String> modelName;
234 std::vector<PD> particle;
235 std::vector<G4int> wasPrinted;
236
237 std::multimap<PD,HP> p_map;
238 std::multimap<HP,HI> m_map;
239
240 // hadronic processes following only G4VProcess interface
241 std::vector<G4VProcess*> extraProcess;
242 std::multimap<PD,G4VProcess*> ep_map;
243
245
246 // counters and options
247 G4int n_proc = 0;
248 G4int n_model = 0;
249 G4int n_part = 0;
250 G4int n_extra = 0;
251
252 G4bool buildTableStart = true;
253 G4bool buildXSTable = false;
254
255 // cache
256 HP currentProcess = nullptr;
257 PD currentParticle = nullptr;
258 PD theGenericIon;
259
260 G4DynamicParticle localDP;
261
262 G4HadronicEPTestMessenger* theEPTestMessenger;
263};
264
265#endif
G4HadronicProcessType
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
void DeRegister(G4HadronicProcess *)
G4double GetCaptureCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=nullptr)
G4double GetCaptureCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
G4double GetCrossSectionPerVolume(const G4ParticleDefinition *particle, G4double kineticEnergy, const G4VProcess *process, const G4Material *material)
G4double GetCaptureCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
G4HadronicProcess * FindProcess(const G4ParticleDefinition *, G4HadronicProcessType subType)
void RegisterParticle(G4HadronicProcess *, const G4ParticleDefinition *)
G4double GetChargeExchangeCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
void PrintHtml(const G4ParticleDefinition *, std::ofstream &)
void SetProcessAbsLevel(G4double absoluteLevel)
G4double GetChargeExchangeCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
G4double GetFissionCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
G4double GetInelasticCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=nullptr)
G4double GetInelasticCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
G4double GetInelasticCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
void SetProcessRelLevel(G4double relativeLevel)
void DeRegisterExtraProcess(G4VProcess *)
G4double GetFissionCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=nullptr)
G4double GetChargeExchangeCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=nullptr)
void RegisterExtraProcess(G4VProcess *)
G4double GetElasticCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
static G4HadronicProcessStore * Instance()
G4double GetElasticCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
void RegisterInteraction(G4HadronicProcess *, G4HadronicInteraction *)
G4double GetCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double kineticEnergy, const G4VProcess *process, const G4Element *element, const G4Material *material=nullptr)
void Register(G4HadronicProcess *)
void PrintModelHtml(const G4HadronicInteraction *model) const
G4double GetElasticCrossSectionPerIsotope(const G4ParticleDefinition *aParticle, G4double kineticEnergy, G4int Z, G4int A)
void PrintInfo(const G4ParticleDefinition *)
G4double GetFissionCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
#define G4ThreadLocal
Definition tls.hh:77