source: trunk/source/processes/hadronic/management/include/G4HadronicProcessStore.hh@ 1344

Last change on this file since 1344 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 6.7 KB
RevLine 
[966]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//
[1340]26// $Id: G4HadronicProcessStore.hh,v 1.7 2010/07/05 14:50:15 vnivanch Exp $
27// GEANT4 tag $Name: hadr-man-V09-03-04 $
[966]28//
29//
30// -------------------------------------------------------------------
31//
32// GEANT4 Class header file
33//
34//
35// File name: G4HadronicProcessStore
36//
37// Author: Vladimir Ivanchenko
38//
39// Creation date: 09.05.2008
40//
41// Modifications:
42//
43//
44// Class Description:
45//
46
47// -------------------------------------------------------------------
48//
49
50#ifndef G4HadronicProcessStore_h
51#define G4HadronicProcessStore_h 1
52
53
54#include "globals.hh"
55#include "G4DynamicParticle.hh"
[1340]56#include "G4ThreeVector.hh"
[966]57#include "G4HadronicProcess.hh"
58#include "G4HadronicInteraction.hh"
59#include "G4ParticleDefinition.hh"
60#include "G4HadronicProcessType.hh"
61#include <map>
62#include <vector>
63
64class G4Element;
[1315]65class G4HadronicEPTestMessenger;
[966]66
[1315]67
[966]68class G4HadronicProcessStore
69{
70
71public:
72
73 static G4HadronicProcessStore* Instance();
74
75 ~G4HadronicProcessStore();
76
[1055]77 void Clean();
78
[966]79 G4double GetInelasticCrossSectionPerVolume(
80 const G4ParticleDefinition *aParticle,
81 G4double kineticEnergy,
82 const G4Material *material);
83
84 G4double GetInelasticCrossSectionPerAtom(
85 const G4ParticleDefinition *aParticle,
86 G4double kineticEnergy,
87 const G4Element *anElement);
88
89 G4double GetInelasticCrossSectionPerIsotope(
90 const G4ParticleDefinition *aParticle,
91 G4double kineticEnergy,
92 G4int Z, G4int A);
93
94 G4double GetElasticCrossSectionPerVolume(
95 const G4ParticleDefinition *aParticle,
96 G4double kineticEnergy,
97 const G4Material *material);
98
99 G4double GetElasticCrossSectionPerAtom(
100 const G4ParticleDefinition *aParticle,
101 G4double kineticEnergy,
102 const G4Element *anElement);
103
104 G4double GetElasticCrossSectionPerIsotope(
105 const G4ParticleDefinition *aParticle,
106 G4double kineticEnergy,
107 G4int Z, G4int A);
108
109 G4double GetCaptureCrossSectionPerVolume(
110 const G4ParticleDefinition *aParticle,
111 G4double kineticEnergy,
112 const G4Material *material);
113
114 G4double GetCaptureCrossSectionPerAtom(
115 const G4ParticleDefinition *aParticle,
116 G4double kineticEnergy,
117 const G4Element *anElement);
118
119 G4double GetCaptureCrossSectionPerIsotope(
120 const G4ParticleDefinition *aParticle,
121 G4double kineticEnergy,
122 G4int Z, G4int A);
123
124 G4double GetFissionCrossSectionPerVolume(
125 const G4ParticleDefinition *aParticle,
126 G4double kineticEnergy,
127 const G4Material *material);
128
129 G4double GetFissionCrossSectionPerAtom(
130 const G4ParticleDefinition *aParticle,
131 G4double kineticEnergy,
132 const G4Element *anElement);
133
134 G4double GetFissionCrossSectionPerIsotope(
135 const G4ParticleDefinition *aParticle,
136 G4double kineticEnergy,
137 G4int Z, G4int A);
138
139 G4double GetChargeExchangeCrossSectionPerVolume(
140 const G4ParticleDefinition *aParticle,
141 G4double kineticEnergy,
142 const G4Material *material);
143
144 G4double GetChargeExchangeCrossSectionPerAtom(
145 const G4ParticleDefinition *aParticle,
146 G4double kineticEnergy,
147 const G4Element *anElement);
148
149 G4double GetChargeExchangeCrossSectionPerIsotope(
150 const G4ParticleDefinition *aParticle,
151 G4double kineticEnergy,
152 G4int Z, G4int A);
153
154 // register/deregister processes following G4HadronicProcess interface
155 void Register(G4HadronicProcess*);
156
157 void RegisterParticle(G4HadronicProcess*,
158 const G4ParticleDefinition*);
159
160 void RegisterInteraction(G4HadronicProcess*,
161 G4HadronicInteraction*);
162
163 void DeRegister(G4HadronicProcess*);
164
165 // register/deregister processes following only G4VProcess interface
166 void RegisterExtraProcess(G4VProcess*);
167
168 void RegisterParticleForExtraProcess(G4VProcess*,
169 const G4ParticleDefinition*);
170
171 void DeRegisterExtraProcess(G4VProcess*);
172
173 void PrintInfo(const G4ParticleDefinition*);
174
175 void Dump(G4int level);
176
177 void SetVerbose(G4int val);
178
179 G4int GetVerbose();
180
181 G4HadronicProcess* FindProcess(const G4ParticleDefinition*,
182 G4HadronicProcessType subType);
183
[1315]184 // Energy-momentum non-conservation limits and reporting
185 void SetEpReportLevel(G4int level);
186
187 void SetProcessAbsLevel(G4double absoluteLevel);
188
189 void SetProcessRelLevel(G4double relativeLevel);
190
[966]191private:
192
193 // constructor
194 G4HadronicProcessStore();
195
196 // print process info
197 void Print(G4int idxProcess, G4int idxParticle);
198
199 static G4HadronicProcessStore* theInstance;
200
201 typedef const G4ParticleDefinition* PD;
202 typedef G4HadronicProcess* HP;
203 typedef G4HadronicInteraction* HI;
204
205 // hadronic processes following G4HadronicProcess interface
206 std::vector<G4HadronicProcess*> process;
207 std::vector<G4HadronicInteraction*> model;
208 std::vector<G4String> modelName;
209 std::vector<PD> particle;
210 std::vector<G4int> wasPrinted;
211
212 std::multimap<PD,HP> p_map;
213 std::multimap<HP,HI> m_map;
214
215 // hadronic processes following only G4VProcess interface
216 std::vector<G4VProcess*> extraProcess;
217 std::multimap<PD,G4VProcess*> ep_map;
218
219 // counters and options
220 G4int n_proc;
221 G4int n_model;
222 G4int n_part;
223 G4int n_extra;
224
225 G4int verbose;
226 G4bool buildTableStart;
227
228 // cash
229 HP currentProcess;
230 PD currentParticle;
231
232 G4DynamicParticle localDP;
233
[1315]234 G4HadronicEPTestMessenger* theEPTestMessenger;
[966]235};
236
237
238#endif
239
Note: See TracBrowser for help on using the repository browser.