4 // Author: I. Hrivnacova
6 // Class TG4PhysicsConstructorHadron
7 // ---------------------------------
8 // Constructor of hadron physics.
9 // According to ExN04HadronPhysics.hh, GEANT4 tag Name: geant4-03-02
11 #ifndef TG4_PHYSICS_CONSTRUCTOR_HADRON_H
12 #define TG4_PHYSICS_CONSTRUCTOR_HADRON_H
14 #include <G4VPhysicsConstructor.hh>
15 #include <g4std/vector>
18 #include <G4MultipleScattering.hh>
19 #include <G4hIonisation.hh>
21 #include <G4HadronElasticProcess.hh>
22 #include <G4HadronFissionProcess.hh>
23 #include <G4HadronCaptureProcess.hh>
25 #include <G4PionPlusInelasticProcess.hh>
26 #include <G4PionMinusInelasticProcess.hh>
27 #include <G4KaonPlusInelasticProcess.hh>
28 #include <G4KaonZeroSInelasticProcess.hh>
29 #include <G4KaonZeroLInelasticProcess.hh>
30 #include <G4KaonMinusInelasticProcess.hh>
31 #include <G4ProtonInelasticProcess.hh>
32 #include <G4AntiProtonInelasticProcess.hh>
33 #include <G4NeutronInelasticProcess.hh>
34 #include <G4AntiNeutronInelasticProcess.hh>
35 #include <G4LambdaInelasticProcess.hh>
36 #include <G4AntiLambdaInelasticProcess.hh>
37 #include <G4SigmaPlusInelasticProcess.hh>
38 #include <G4SigmaMinusInelasticProcess.hh>
39 #include <G4AntiSigmaPlusInelasticProcess.hh>
40 #include <G4AntiSigmaMinusInelasticProcess.hh>
41 #include <G4XiZeroInelasticProcess.hh>
42 #include <G4XiMinusInelasticProcess.hh>
43 #include <G4AntiXiZeroInelasticProcess.hh>
44 #include <G4AntiXiMinusInelasticProcess.hh>
45 #include <G4DeuteronInelasticProcess.hh>
46 #include <G4TritonInelasticProcess.hh>
47 #include <G4AlphaInelasticProcess.hh>
48 #include <G4OmegaMinusInelasticProcess.hh>
49 #include <G4AntiOmegaMinusInelasticProcess.hh>
52 #include <G4LElastic.hh>
53 #include <G4LFission.hh>
54 #include <G4LCapture.hh>
56 #include <G4LEPionPlusInelastic.hh>
57 #include <G4LEPionMinusInelastic.hh>
58 #include <G4LEKaonPlusInelastic.hh>
59 #include <G4LEKaonZeroSInelastic.hh>
60 #include <G4LEKaonZeroLInelastic.hh>
61 #include <G4LEKaonMinusInelastic.hh>
62 #include <G4LEProtonInelastic.hh>
63 #include <G4LEAntiProtonInelastic.hh>
64 #include <G4LENeutronInelastic.hh>
65 #include <G4LEAntiNeutronInelastic.hh>
66 #include <G4LELambdaInelastic.hh>
67 #include <G4LEAntiLambdaInelastic.hh>
68 #include <G4LESigmaPlusInelastic.hh>
69 #include <G4LESigmaMinusInelastic.hh>
70 #include <G4LEAntiSigmaPlusInelastic.hh>
71 #include <G4LEAntiSigmaMinusInelastic.hh>
72 #include <G4LEXiZeroInelastic.hh>
73 #include <G4LEXiMinusInelastic.hh>
74 #include <G4LEAntiXiZeroInelastic.hh>
75 #include <G4LEAntiXiMinusInelastic.hh>
76 #include <G4LEDeuteronInelastic.hh>
77 #include <G4LETritonInelastic.hh>
78 #include <G4LEAlphaInelastic.hh>
79 #include <G4LEOmegaMinusInelastic.hh>
80 #include <G4LEAntiOmegaMinusInelastic.hh>
84 #include <G4HEPionPlusInelastic.hh>
85 #include <G4HEPionMinusInelastic.hh>
86 #include <G4HEKaonPlusInelastic.hh>
87 #include <G4HEKaonZeroInelastic.hh>
88 #include <G4HEKaonZeroInelastic.hh>
89 #include <G4HEKaonMinusInelastic.hh>
90 #include <G4HEProtonInelastic.hh>
91 #include <G4HEAntiProtonInelastic.hh>
92 #include <G4HENeutronInelastic.hh>
93 #include <G4HEAntiNeutronInelastic.hh>
94 #include <G4HELambdaInelastic.hh>
95 #include <G4HEAntiLambdaInelastic.hh>
96 #include <G4HESigmaPlusInelastic.hh>
97 #include <G4HESigmaMinusInelastic.hh>
98 #include <G4HEAntiSigmaPlusInelastic.hh>
99 #include <G4HEAntiSigmaMinusInelastic.hh>
100 #include <G4HEXiZeroInelastic.hh>
101 #include <G4HEXiMinusInelastic.hh>
102 #include <G4HEAntiXiZeroInelastic.hh>
103 #include <G4HEAntiXiMinusInelastic.hh>
104 #include <G4HEOmegaMinusInelastic.hh>
105 #include <G4HEAntiOmegaMinusInelastic.hh>
107 // Stopping processes
108 #include <G4AntiProtonAnnihilationAtRest.hh>
109 #include <G4AntiNeutronAnnihilationAtRest.hh>
111 #ifdef TRIUMF_STOP_PIMINUS
112 #include <G4PionMinusAbsorptionAtRest.hh>
114 #include <G4PiMinusAbsorptionAtRest.hh>
116 #ifdef TRIUMF_STOP_KMINUS
117 #include <G4KaonMinusAbsorption.hh>
119 #include <G4KaonMinusAbsorptionAtRest.hh>
122 class TG4PhysicsConstructorHadron: public G4VPhysicsConstructor
124 typedef G4std::vector<G4VProcess*> ProcessVector;
127 TG4PhysicsConstructorHadron(const G4String& name = "Hadron");
129 // TG4PhysicsConstructorHadron(const TG4PhysicsConstructorHadron& right);
130 virtual ~TG4PhysicsConstructorHadron();
133 TG4PhysicsConstructorHadron(const TG4PhysicsConstructorHadron& right);
136 TG4PhysicsConstructorHadron& operator=(
137 const TG4PhysicsConstructorHadron& right);
140 // construct particle and physics
141 virtual void ConstructParticle();
142 virtual void ConstructProcess();
146 G4HadronElasticProcess fElasticProcess; //hadron elastic process
147 G4LElastic* fElasticModel; //elastic model
150 G4PionPlusInelasticProcess fPionPlusInelastic; //pi+ inel process
151 G4LEPionPlusInelastic* fLEPionPlusModel; //pi+ LE inel model
152 G4HEPionPlusInelastic* fHEPionPlusModel; //pi+ HE inel model
153 G4MultipleScattering fPionPlusMult; //pi+ msc
154 G4hIonisation fPionPlusIonisation; //pi+ ionisation
157 G4PionMinusInelasticProcess fPionMinusInelastic; //pi- inel process
158 G4LEPionMinusInelastic* fLEPionMinusModel; //pi- LE inel model
159 G4HEPionMinusInelastic* fHEPionMinusModel; //pi- HE inel model
160 G4MultipleScattering fPionMinusMult; //pi- msc
161 G4hIonisation fPionMinusIonisation; //pi- ionisation
162 #ifdef TRIUMF_STOP_PIMINUS
163 G4PionMinusAbsorptionAtRest fPionMinusAbsorption; //pi- absorption
165 G4PiMinusAbsorptionAtRest fPionMinusAbsorption; //pi- absorption
169 G4KaonPlusInelasticProcess fKaonPlusInelastic; //kaon+ inel process
170 G4LEKaonPlusInelastic* fLEKaonPlusModel; //kaon+ LE inel model
171 G4HEKaonPlusInelastic* fHEKaonPlusModel; //kaon+ HE inel model
172 G4MultipleScattering fKaonPlusMult; //kaon+ msc
173 G4hIonisation fKaonPlusIonisation; //kaon+ ionisation
176 G4KaonMinusInelasticProcess fKaonMinusInelastic; //kaon- inel process
177 G4LEKaonMinusInelastic* fLEKaonMinusModel; //kaon- LE inel model
178 G4HEKaonMinusInelastic* fHEKaonMinusModel; //kaon- HE inel model
179 G4MultipleScattering fKaonMinusMult; //kaon- msc
180 G4hIonisation fKaonMinusIonisation; //kaon- ionisation
181 #ifdef TRIUMF_STOP_KMINUS
182 G4KaonMinusAbsorption fKaonMinusAbsorption; //kaon- absorption
184 G4PiMinusAbsorptionAtRest fKaonMinusAbsorption; //kaon- absorption
188 G4KaonZeroLInelasticProcess fKaonZeroLInelastic; //kaon0 inel process
189 G4LEKaonZeroLInelastic* fLEKaonZeroLModel; //kaon0 LE inel model
190 G4HEKaonZeroInelastic* fHEKaonZeroLModel; //kaon0 HE inel model
193 G4KaonZeroSInelasticProcess fKaonZeroSInelastic; //kaon0S inel process
194 G4LEKaonZeroSInelastic* fLEKaonZeroSModel; //kaon0S LE inel model
195 G4HEKaonZeroInelastic* fHEKaonZeroSModel; //kaon0S HE inel mode
198 G4ProtonInelasticProcess fProtonInelastic; //p inel process
199 G4LEProtonInelastic* fLEProtonModel; //p LE inel model
200 G4HEProtonInelastic* fHEProtonModel; //p HE inel model
201 G4MultipleScattering fProtonMult; //p msc
202 G4hIonisation fProtonIonisation; //p ionisation
205 G4AntiProtonInelasticProcess fAntiProtonInelastic; //p_bar inel process
206 G4LEAntiProtonInelastic* fLEAntiProtonModel; //p_bar LE inel model
207 G4HEAntiProtonInelastic* fHEAntiProtonModel; //p_bar HE inel model
208 G4MultipleScattering fAntiProtonMult; //p_bar msc
209 G4hIonisation fAntiProtonIonisation;//p_bar ionisation
210 G4AntiProtonAnnihilationAtRest fAntiProtonAnnihilation;//p_bar annihilation
213 G4NeutronInelasticProcess fNeutronInelastic; //n inel process
214 G4LENeutronInelastic* fLENeutronModel; //n LE inel model
215 G4HENeutronInelastic* fHENeutronModel; //n HE inel model
216 G4HadronFissionProcess fNeutronFission; //n fission
217 G4LFission* fNeutronFissionModel; //n fission model
218 G4HadronCaptureProcess fNeutronCapture; //n capture
219 G4LCapture* fNeutronCaptureModel; //n capture model
222 G4AntiNeutronInelasticProcess fAntiNeutronInelastic;//n_bar inel process
223 G4LEAntiNeutronInelastic* fLEAntiNeutronModel; //n_bar LE inel model
224 G4HEAntiNeutronInelastic* fHEAntiNeutronModel; //n_bar HE inel model
225 G4AntiNeutronAnnihilationAtRest fAntiNeutronAnnihilation;//n_bar ionisation
228 G4LambdaInelasticProcess fLambdaInelastic; //lambda inel process
229 G4LELambdaInelastic* fLELambdaModel; //lambda LE inel model
230 G4HELambdaInelastic* fHELambdaModel; //lambda HE inel model
233 G4AntiLambdaInelasticProcess fAntiLambdaInelastic; //lambda_bar inel process
234 G4LEAntiLambdaInelastic* fLEAntiLambdaModel; //lambda_bar LE inel model
235 G4HEAntiLambdaInelastic* fHEAntiLambdaModel; //lambda_bar HE inel model
238 G4SigmaMinusInelasticProcess fSigmaMinusInelastic; //sigma- inel process
239 G4LESigmaMinusInelastic* fLESigmaMinusModel; //sigma- LE inel model
240 G4HESigmaMinusInelastic* fHESigmaMinusModel; //sigma- HE inel model
241 G4MultipleScattering fSigmaMinusMult; //sigma- msc
242 G4hIonisation fSigmaMinusIonisation;//sigma- ionisation
245 G4AntiSigmaMinusInelasticProcess fAntiSigmaMinusInelastic; //sigma-_bar inel process
246 G4LEAntiSigmaMinusInelastic* fLEAntiSigmaMinusModel; //sigma-_bar LE inel model
247 G4HEAntiSigmaMinusInelastic* fHEAntiSigmaMinusModel; //sigma-_bar HE inel model
248 G4MultipleScattering fAntiSigmaMinusMult; //sigma-_bar msc
249 G4hIonisation fAntiSigmaMinusIonisation;//sigma-_bar ionisation
252 G4SigmaPlusInelasticProcess fSigmaPlusInelastic; //sigma+ inel process
253 G4LESigmaPlusInelastic* fLESigmaPlusModel; //sigma+ LE inel model
254 G4HESigmaPlusInelastic* fHESigmaPlusModel; //sigma+ HE inel model
255 G4MultipleScattering fSigmaPlusMult; //sigma+ msc
256 G4hIonisation fSigmaPlusIonisation; //sigma+ ionisation
259 G4AntiSigmaPlusInelasticProcess fAntiSigmaPlusInelastic; //sigma+_bar inel process
260 G4LEAntiSigmaPlusInelastic* fLEAntiSigmaPlusModel; //sigma+_bar LE inel model
261 G4HEAntiSigmaPlusInelastic* fHEAntiSigmaPlusModel; //sigma+_bar HE inel model
262 G4MultipleScattering fAntiSigmaPlusMult; //sigma+_bar msc
263 G4hIonisation fAntiSigmaPlusIonisation; //sigma+_bar ionisation
266 G4XiZeroInelasticProcess fXiZeroInelastic; //xi0 inel process
267 G4LEXiZeroInelastic* fLEXiZeroModel; //xi0 LE inel model
268 G4HEXiZeroInelastic* fHEXiZeroModel; //xi0 HE inel model
271 G4AntiXiZeroInelasticProcess fAntiXiZeroInelastic;//xi0_bar inel process
272 G4LEAntiXiZeroInelastic* fLEAntiXiZeroModel; //xi0_bar LE inel model
273 G4HEAntiXiZeroInelastic* fHEAntiXiZeroModel; //xi0_bar HE inel model
276 G4XiMinusInelasticProcess fXiMinusInelastic; //xi- inel process
277 G4LEXiMinusInelastic* fLEXiMinusModel; //xi- LE inel model
278 G4HEXiMinusInelastic* fHEXiMinusModel; //xi- HE inel model
279 G4MultipleScattering fXiMinusMult; //xi- msc
280 G4hIonisation fXiMinusIonisation; //xi- ionisation
283 G4AntiXiMinusInelasticProcess fAntiXiMinusInelastic; //xi-_bar inel process
284 G4LEAntiXiMinusInelastic* fLEAntiXiMinusModel; //xi-_bar LE inel model
285 G4HEAntiXiMinusInelastic* fHEAntiXiMinusModel; //xi-_bar HE inel model
286 G4MultipleScattering fAntiXiMinusMult; //xi-_bar msc
287 G4hIonisation fAntiXiMinusIonisation;//xi-_bar ionisation
290 G4OmegaMinusInelasticProcess fOmegaMinusInelastic; //omega- inel process
291 G4LEOmegaMinusInelastic* fLEOmegaMinusModel; //omega- LE inel model
292 G4HEOmegaMinusInelastic* fHEOmegaMinusModel; //omega- HE inel model
293 G4MultipleScattering fOmegaMinusMult; //omega- msc
294 G4hIonisation fOmegaMinusIonisation; //omega- ionisation
297 G4AntiOmegaMinusInelasticProcess fAntiOmegaMinusInelastic; //omega-_bar inel process
298 G4LEAntiOmegaMinusInelastic* fLEAntiOmegaMinusModel; //omega-_bar LE inel model
299 G4HEAntiOmegaMinusInelastic* fHEAntiOmegaMinusModel; //omega-_bar HE inel model
300 G4MultipleScattering fAntiOmegaMinusMult; //omega-_bar msc
301 G4hIonisation fAntiOmegaMinusIonisation;//omega-_bar ionisation
304 ProcessVector fOtherProcesses; //other process
309 void ConstructProcessForPionPlus();
310 void ConstructProcessForPionMinus();
311 void ConstructProcessForKaonPlus();
312 void ConstructProcessForKaonMinus();
313 void ConstructProcessForKaonZeroLong();
314 void ConstructProcessForKaonZeroShort();
315 void ConstructProcessForProton();
316 void ConstructProcessForAntiProton();
317 void ConstructProcessForNeutron();
318 void ConstructProcessForAntiNeutron();
319 void ConstructProcessForLambda();
320 void ConstructProcessForAntiLambda();
321 void ConstructProcessForSigmaMinus();
322 void ConstructProcessForAntiSigmaMinus();
323 void ConstructProcessForSigmaPlus();
324 void ConstructProcessForAntiSigmaPlus();
325 void ConstructProcessForXiMinus();
326 void ConstructProcessForAntiXiMinus();
327 void ConstructProcessForXiZero();
328 void ConstructProcessForAntiXiZero();
329 void ConstructProcessForOmegaMinus();
330 void ConstructProcessForAntiOmegaMinus();
331 void ConstructProcessForOther();
334 #endif //TG4_PHYSICS_CONSTRUCTOR_HADRON_H