Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
27//
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4HadronicProcessStore
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 09.05.2008
39//
40// Modifications:
41//
42//
43// Class Description:
44//
45
46// -------------------------------------------------------------------
47//
48
49#ifndef G4HadronicProcessStore_h
50#define G4HadronicProcessStore_h 1
51
52
53#include "globals.hh"
54#include "G4DynamicParticle.hh"
55#include "G4ThreeVector.hh"
56#include "G4HadronicProcess.hh"
60#include <map>
61#include <vector>
62#include <iostream>
63
64class G4Element;
66
67
69{
70
71public:
72
74
76
77 void Clean();
79 const G4ParticleDefinition* particle,
80 G4double kineticEnergy,
81 const G4VProcess* process,
82 const G4Element* element);
83
85 const G4ParticleDefinition* particle,
86 G4double kineticEnergy,
87 const G4VProcess* process,
88 const G4Material* material);
89
91 const G4ParticleDefinition *aParticle,
92 G4double kineticEnergy,
93 const G4Material *material);
94
96 const G4ParticleDefinition *aParticle,
97 G4double kineticEnergy,
98 const G4Element *anElement, const G4Material* mat=0);
99
101 const G4ParticleDefinition *aParticle,
102 G4double kineticEnergy,
103 G4int Z, G4int A);
104
106 const G4ParticleDefinition *aParticle,
107 G4double kineticEnergy,
108 const G4Material *material);
109
111 const G4ParticleDefinition *aParticle,
112 G4double kineticEnergy,
113 const G4Element *anElement, const G4Material* mat=0);
114
116 const G4ParticleDefinition *aParticle,
117 G4double kineticEnergy,
118 G4int Z, G4int A);
119
121 const G4ParticleDefinition *aParticle,
122 G4double kineticEnergy,
123 const G4Material *material);
124
126 const G4ParticleDefinition *aParticle,
127 G4double kineticEnergy,
128 const G4Element *anElement, const G4Material* mat=0);
129
131 const G4ParticleDefinition *aParticle,
132 G4double kineticEnergy,
133 G4int Z, G4int A);
134
136 const G4ParticleDefinition *aParticle,
137 G4double kineticEnergy,
138 const G4Material *material);
139
141 const G4ParticleDefinition *aParticle,
142 G4double kineticEnergy,
143 const G4Element *anElement, const G4Material* mat=0);
144
146 const G4ParticleDefinition *aParticle,
147 G4double kineticEnergy,
148 G4int Z, G4int A);
149
151 const G4ParticleDefinition *aParticle,
152 G4double kineticEnergy,
153 const G4Material *material);
154
156 const G4ParticleDefinition *aParticle,
157 G4double kineticEnergy,
158 const G4Element *anElement, const G4Material* mat=0);
159
161 const G4ParticleDefinition *aParticle,
162 G4double kineticEnergy,
163 G4int Z, G4int A);
164
165 // register/deregister processes following G4HadronicProcess interface
167
169 const G4ParticleDefinition*);
170
173
175
176 // register/deregister processes following only G4VProcess interface
178
180 const G4ParticleDefinition*);
181
183
184 void PrintInfo(const G4ParticleDefinition*);
185
186 void Dump(G4int level);
187 void DumpHtml();
188 void PrintHtml(const G4ParticleDefinition*, std::ofstream&);
189 void PrintModelHtml(const G4HadronicInteraction * model) const;
190
191 void SetVerbose(G4int val);
192
194
196 G4HadronicProcessType subType);
197
198 // Energy-momentum non-conservation limits and reporting
199 void SetEpReportLevel(G4int level);
200
201 void SetProcessAbsLevel(G4double absoluteLevel);
202
203 void SetProcessRelLevel(G4double relativeLevel);
204
205private:
206
207 // constructor
209
210 // print process info
211 void Print(G4int idxProcess, G4int idxParticle);
212
213 static G4HadronicProcessStore* theInstance;
214
215 typedef const G4ParticleDefinition* PD;
216 typedef G4HadronicProcess* HP;
217 typedef G4HadronicInteraction* HI;
218
219 // hadronic processes following G4HadronicProcess interface
220 std::vector<G4HadronicProcess*> process;
221 std::vector<G4HadronicInteraction*> model;
222 std::vector<G4String> modelName;
223 std::vector<PD> particle;
224 std::vector<G4int> wasPrinted;
225
226 std::multimap<PD,HP> p_map;
227 std::multimap<HP,HI> m_map;
228
229 // hadronic processes following only G4VProcess interface
230 std::vector<G4VProcess*> extraProcess;
231 std::multimap<PD,G4VProcess*> ep_map;
232
233 // counters and options
234 G4int n_proc;
235 G4int n_model;
236 G4int n_part;
237 G4int n_extra;
238
239 G4int verbose;
240 G4bool buildTableStart;
241
242 // cash
243 HP currentProcess;
244 PD currentParticle;
245
246 G4DynamicParticle localDP;
247
248 G4HadronicEPTestMessenger* theEPTestMessenger;
249};
250
251
252#endif
253
G4HadronicProcessType
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void DeRegister(G4HadronicProcess *)
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 GetFissionCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
G4double GetCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double kineticEnergy, const G4VProcess *process, const G4Element *element)
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 *)
void RegisterExtraProcess(G4VProcess *)
G4double GetElasticCrossSectionPerVolume(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Material *material)
G4double GetCaptureCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
void SetEpReportLevel(G4int level)
static G4HadronicProcessStore * Instance()
G4double GetElasticCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
void RegisterInteraction(G4HadronicProcess *, G4HadronicInteraction *)
G4double GetChargeExchangeCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
G4double GetInelasticCrossSectionPerAtom(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4Element *anElement, const G4Material *mat=0)
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)