]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant4/TG4PhysicsConstructorHadron.h
Bugfix in AliPoints2Memory
[u/mrichter/AliRoot.git] / TGeant4 / TG4PhysicsConstructorHadron.h
1 // $Id$
2 // Category: physics
3 //
4 // Author: I. Hrivnacova
5 //
6 // Class TG4PhysicsConstructorHadron
7 // ---------------------------------
8 // Constructor of hadron physics.
9 // According to ExN04HadronPhysics.hh, GEANT4 tag Name: geant4-03-02
10
11 #ifndef TG4_PHYSICS_CONSTRUCTOR_HADRON_H
12 #define TG4_PHYSICS_CONSTRUCTOR_HADRON_H
13
14 #include "TG4VPhysicsConstructor.h"
15
16 #include <g4std/vector>
17 #include <globals.hh>
18
19 #include <G4MultipleScattering.hh>
20 #include <G4hIonisation.hh>
21
22 #include <G4HadronElasticProcess.hh>
23 #include <G4HadronFissionProcess.hh>
24 #include <G4HadronCaptureProcess.hh>
25
26 #include <G4PionPlusInelasticProcess.hh>
27 #include <G4PionMinusInelasticProcess.hh>
28 #include <G4KaonPlusInelasticProcess.hh>
29 #include <G4KaonZeroSInelasticProcess.hh>
30 #include <G4KaonZeroLInelasticProcess.hh>
31 #include <G4KaonMinusInelasticProcess.hh>
32 #include <G4ProtonInelasticProcess.hh>
33 #include <G4AntiProtonInelasticProcess.hh>
34 #include <G4NeutronInelasticProcess.hh>
35 #include <G4AntiNeutronInelasticProcess.hh>
36 #include <G4LambdaInelasticProcess.hh>
37 #include <G4AntiLambdaInelasticProcess.hh>
38 #include <G4SigmaPlusInelasticProcess.hh>
39 #include <G4SigmaMinusInelasticProcess.hh>
40 #include <G4AntiSigmaPlusInelasticProcess.hh>
41 #include <G4AntiSigmaMinusInelasticProcess.hh>
42 #include <G4XiZeroInelasticProcess.hh>
43 #include <G4XiMinusInelasticProcess.hh>
44 #include <G4AntiXiZeroInelasticProcess.hh>
45 #include <G4AntiXiMinusInelasticProcess.hh>
46 #include <G4DeuteronInelasticProcess.hh>
47 #include <G4TritonInelasticProcess.hh>
48 #include <G4AlphaInelasticProcess.hh>
49 #include <G4OmegaMinusInelasticProcess.hh>
50 #include <G4AntiOmegaMinusInelasticProcess.hh>
51
52 // Low-energy Models
53 #include <G4LElastic.hh>   
54 #include <G4LFission.hh>
55 #include <G4LCapture.hh>
56
57 #include <G4LEPionPlusInelastic.hh>
58 #include <G4LEPionMinusInelastic.hh>
59 #include <G4LEKaonPlusInelastic.hh>
60 #include <G4LEKaonZeroSInelastic.hh>
61 #include <G4LEKaonZeroLInelastic.hh>
62 #include <G4LEKaonMinusInelastic.hh>
63 #include <G4LEProtonInelastic.hh>
64 #include <G4LEAntiProtonInelastic.hh>
65 #include <G4LENeutronInelastic.hh>
66 #include <G4LEAntiNeutronInelastic.hh>
67 #include <G4LELambdaInelastic.hh>
68 #include <G4LEAntiLambdaInelastic.hh>
69 #include <G4LESigmaPlusInelastic.hh>
70 #include <G4LESigmaMinusInelastic.hh>
71 #include <G4LEAntiSigmaPlusInelastic.hh>
72 #include <G4LEAntiSigmaMinusInelastic.hh>
73 #include <G4LEXiZeroInelastic.hh>
74 #include <G4LEXiMinusInelastic.hh>
75 #include <G4LEAntiXiZeroInelastic.hh>
76 #include <G4LEAntiXiMinusInelastic.hh>
77 #include <G4LEDeuteronInelastic.hh>
78 #include <G4LETritonInelastic.hh>
79 #include <G4LEAlphaInelastic.hh>
80 #include <G4LEOmegaMinusInelastic.hh>
81 #include <G4LEAntiOmegaMinusInelastic.hh>
82
83 // High-energy Models
84
85 #include <G4HEPionPlusInelastic.hh>
86 #include <G4HEPionMinusInelastic.hh>
87 #include <G4HEKaonPlusInelastic.hh>
88 #include <G4HEKaonZeroInelastic.hh>
89 #include <G4HEKaonZeroInelastic.hh>
90 #include <G4HEKaonMinusInelastic.hh>
91 #include <G4HEProtonInelastic.hh>
92 #include <G4HEAntiProtonInelastic.hh>
93 #include <G4HENeutronInelastic.hh>
94 #include <G4HEAntiNeutronInelastic.hh>
95 #include <G4HELambdaInelastic.hh>
96 #include <G4HEAntiLambdaInelastic.hh>
97 #include <G4HESigmaPlusInelastic.hh>
98 #include <G4HESigmaMinusInelastic.hh>
99 #include <G4HEAntiSigmaPlusInelastic.hh>
100 #include <G4HEAntiSigmaMinusInelastic.hh>
101 #include <G4HEXiZeroInelastic.hh>
102 #include <G4HEXiMinusInelastic.hh>
103 #include <G4HEAntiXiZeroInelastic.hh>
104 #include <G4HEAntiXiMinusInelastic.hh>
105 #include <G4HEOmegaMinusInelastic.hh>
106 #include <G4HEAntiOmegaMinusInelastic.hh>
107
108 // Stopping processes
109 #include <G4AntiProtonAnnihilationAtRest.hh>
110 #include <G4AntiNeutronAnnihilationAtRest.hh>
111
112 #ifdef TRIUMF_STOP_PIMINUS
113 #include <G4PionMinusAbsorptionAtRest.hh>
114 #else
115 #include <G4PiMinusAbsorptionAtRest.hh>
116 #endif
117 #ifdef TRIUMF_STOP_KMINUS
118 #include <G4KaonMinusAbsorption.hh>
119 #else
120 #include <G4KaonMinusAbsorptionAtRest.hh>
121 #endif
122
123 class TG4PhysicsConstructorHadron: public TG4VPhysicsConstructor
124 {
125   typedef G4std::vector<G4VProcess*>  ProcessVector;
126
127   public:
128     TG4PhysicsConstructorHadron(const G4String& name = "Hadron");
129     TG4PhysicsConstructorHadron(G4int verboseLevel, 
130                                 G4bool setEM, G4bool setHadron,
131                                 const G4String& name = "Hadron");
132     // --> protected
133     // TG4PhysicsConstructorHadron(const TG4PhysicsConstructorHadron& right);
134     virtual ~TG4PhysicsConstructorHadron();
135
136   protected:
137     TG4PhysicsConstructorHadron(const TG4PhysicsConstructorHadron& right);
138     
139     // operators
140     TG4PhysicsConstructorHadron& operator=(
141                                 const TG4PhysicsConstructorHadron& right);
142   
143     // methods
144           // construct particle and physics
145     virtual void ConstructParticle();
146     virtual void ConstructProcess();
147
148     // data members
149          // Elastic Process
150     G4HadronElasticProcess  fElasticProcess;            //hadron elastic process
151     G4LElastic*             fElasticModel;              //elastic model
152   
153          // Pi + 
154     G4PionPlusInelasticProcess fPionPlusInelastic;      //pi+ inel process
155     G4LEPionPlusInelastic*     fLEPionPlusModel;        //pi+ LE inel model
156     G4HEPionPlusInelastic*     fHEPionPlusModel;        //pi+ HE inel model
157     G4MultipleScattering       fPionPlusMult;           //pi+ msc
158     G4hIonisation              fPionPlusIonisation;     //pi+ ionisation
159
160          // Pi -
161     G4PionMinusInelasticProcess  fPionMinusInelastic;   //pi- inel process
162     G4LEPionMinusInelastic*      fLEPionMinusModel;     //pi- LE inel model
163     G4HEPionMinusInelastic*      fHEPionMinusModel;     //pi- HE inel model
164     G4MultipleScattering         fPionMinusMult;        //pi- msc
165     G4hIonisation                fPionMinusIonisation;  //pi- ionisation
166 #ifdef TRIUMF_STOP_PIMINUS
167     G4PionMinusAbsorptionAtRest  fPionMinusAbsorption;  //pi- absorption
168 #else
169     G4PiMinusAbsorptionAtRest    fPionMinusAbsorption;  //pi- absorption
170 #endif
171
172          // K + 
173     G4KaonPlusInelasticProcess  fKaonPlusInelastic;     //kaon+ inel process
174     G4LEKaonPlusInelastic*      fLEKaonPlusModel;       //kaon+ LE inel model
175     G4HEKaonPlusInelastic*      fHEKaonPlusModel;       //kaon+ HE inel model
176     G4MultipleScattering        fKaonPlusMult;          //kaon+ msc
177     G4hIonisation               fKaonPlusIonisation;    //kaon+ ionisation
178
179          // K -
180     G4KaonMinusInelasticProcess  fKaonMinusInelastic;   //kaon- inel process
181     G4LEKaonMinusInelastic*      fLEKaonMinusModel;     //kaon- LE inel model
182     G4HEKaonMinusInelastic*      fHEKaonMinusModel;     //kaon- HE inel model
183     G4MultipleScattering         fKaonMinusMult;        //kaon- msc
184     G4hIonisation                fKaonMinusIonisation;  //kaon- ionisation
185 #ifdef TRIUMF_STOP_KMINUS
186     G4KaonMinusAbsorption        fKaonMinusAbsorption;  //kaon- absorption
187 #else
188     G4PiMinusAbsorptionAtRest    fKaonMinusAbsorption;  //kaon- absorption
189 #endif
190
191         // K0L
192     G4KaonZeroLInelasticProcess  fKaonZeroLInelastic;   //kaon0 inel process
193     G4LEKaonZeroLInelastic*      fLEKaonZeroLModel;     //kaon0 LE inel model
194     G4HEKaonZeroInelastic*       fHEKaonZeroLModel;     //kaon0 HE inel model
195
196         // K0S
197     G4KaonZeroSInelasticProcess  fKaonZeroSInelastic;   //kaon0S inel process
198     G4LEKaonZeroSInelastic*      fLEKaonZeroSModel;     //kaon0S LE inel model
199     G4HEKaonZeroInelastic*       fHEKaonZeroSModel;     //kaon0S HE inel mode
200
201         // Proton
202     G4ProtonInelasticProcess  fProtonInelastic;         //p inel process
203     G4LEProtonInelastic*      fLEProtonModel;           //p LE inel model
204     G4HEProtonInelastic*      fHEProtonModel;           //p HE inel model
205     G4MultipleScattering      fProtonMult;              //p msc
206     G4hIonisation             fProtonIonisation;        //p ionisation
207  
208         // anti-proton
209     G4AntiProtonInelasticProcess    fAntiProtonInelastic; //p_bar inel process
210     G4LEAntiProtonInelastic*        fLEAntiProtonModel;   //p_bar LE inel model
211     G4HEAntiProtonInelastic*        fHEAntiProtonModel;   //p_bar HE inel model
212     G4MultipleScattering            fAntiProtonMult;      //p_bar msc
213     G4hIonisation                   fAntiProtonIonisation;//p_bar ionisation
214     G4AntiProtonAnnihilationAtRest  fAntiProtonAnnihilation;//p_bar annihilation
215     
216        // neutron
217     G4NeutronInelasticProcess  fNeutronInelastic;       //n inel process
218     G4LENeutronInelastic*      fLENeutronModel;         //n LE inel model
219     G4HENeutronInelastic*      fHENeutronModel;         //n HE inel model
220     G4HadronFissionProcess     fNeutronFission;         //n fission
221     G4LFission*                fNeutronFissionModel;    //n fission model
222     G4HadronCaptureProcess     fNeutronCapture;         //n capture
223     G4LCapture*                fNeutronCaptureModel;    //n capture model
224
225        // anti-neutron
226     G4AntiNeutronInelasticProcess    fAntiNeutronInelastic;//n_bar inel process
227     G4LEAntiNeutronInelastic*        fLEAntiNeutronModel;  //n_bar LE inel model
228     G4HEAntiNeutronInelastic*        fHEAntiNeutronModel;  //n_bar HE inel model
229     G4AntiNeutronAnnihilationAtRest  fAntiNeutronAnnihilation;//n_bar ionisation
230      
231        // Lambda
232     G4LambdaInelasticProcess  fLambdaInelastic;         //lambda inel process
233     G4LELambdaInelastic*      fLELambdaModel;           //lambda LE inel model
234     G4HELambdaInelastic*      fHELambdaModel;           //lambda HE inel model
235   
236        // AntiLambda
237     G4AntiLambdaInelasticProcess  fAntiLambdaInelastic; //lambda_bar inel process
238     G4LEAntiLambdaInelastic*      fLEAntiLambdaModel;   //lambda_bar LE inel model
239     G4HEAntiLambdaInelastic*      fHEAntiLambdaModel;   //lambda_bar HE inel model
240   
241        // SigmaMinus
242     G4SigmaMinusInelasticProcess  fSigmaMinusInelastic; //sigma- inel process
243     G4LESigmaMinusInelastic*      fLESigmaMinusModel;   //sigma- LE inel model
244     G4HESigmaMinusInelastic*      fHESigmaMinusModel;   //sigma- HE inel model
245     G4MultipleScattering          fSigmaMinusMult;      //sigma- msc
246     G4hIonisation                 fSigmaMinusIonisation;//sigma- ionisation
247   
248        // AntiSigmaMinus
249     G4AntiSigmaMinusInelasticProcess  fAntiSigmaMinusInelastic; //sigma-_bar inel process
250     G4LEAntiSigmaMinusInelastic*      fLEAntiSigmaMinusModel;   //sigma-_bar LE inel model
251     G4HEAntiSigmaMinusInelastic*      fHEAntiSigmaMinusModel;   //sigma-_bar HE inel model
252     G4MultipleScattering              fAntiSigmaMinusMult;      //sigma-_bar msc
253     G4hIonisation                     fAntiSigmaMinusIonisation;//sigma-_bar ionisation
254    
255        // SigmaPlus
256     G4SigmaPlusInelasticProcess  fSigmaPlusInelastic;   //sigma+ inel process
257     G4LESigmaPlusInelastic*      fLESigmaPlusModel;     //sigma+ LE inel model
258     G4HESigmaPlusInelastic*      fHESigmaPlusModel;     //sigma+ HE inel model
259     G4MultipleScattering         fSigmaPlusMult;        //sigma+ msc
260     G4hIonisation                fSigmaPlusIonisation;  //sigma+ ionisation
261   
262        // AntiSigmaPlus
263     G4AntiSigmaPlusInelasticProcess  fAntiSigmaPlusInelastic;  //sigma+_bar inel process
264     G4LEAntiSigmaPlusInelastic*      fLEAntiSigmaPlusModel;    //sigma+_bar LE inel model
265     G4HEAntiSigmaPlusInelastic*      fHEAntiSigmaPlusModel;    //sigma+_bar HE inel model
266     G4MultipleScattering             fAntiSigmaPlusMult;       //sigma+_bar msc
267     G4hIonisation                    fAntiSigmaPlusIonisation; //sigma+_bar ionisation
268   
269       // XiZero
270     G4XiZeroInelasticProcess  fXiZeroInelastic;        //xi0 inel process
271     G4LEXiZeroInelastic*      fLEXiZeroModel;          //xi0 LE inel model
272     G4HEXiZeroInelastic*      fHEXiZeroModel;          //xi0 HE inel model
273   
274       // AntiXiZero
275     G4AntiXiZeroInelasticProcess  fAntiXiZeroInelastic;//xi0_bar inel process
276     G4LEAntiXiZeroInelastic*      fLEAntiXiZeroModel;  //xi0_bar LE inel model
277     G4HEAntiXiZeroInelastic*      fHEAntiXiZeroModel;  //xi0_bar HE inel model
278   
279       // XiMinus
280     G4XiMinusInelasticProcess  fXiMinusInelastic;      //xi- inel process
281     G4LEXiMinusInelastic*      fLEXiMinusModel;        //xi- LE inel model
282     G4HEXiMinusInelastic*      fHEXiMinusModel;        //xi- HE inel model
283     G4MultipleScattering       fXiMinusMult;           //xi- msc
284     G4hIonisation              fXiMinusIonisation;     //xi- ionisation
285
286       // AntiXiMinus
287     G4AntiXiMinusInelasticProcess  fAntiXiMinusInelastic; //xi-_bar inel process
288     G4LEAntiXiMinusInelastic*      fLEAntiXiMinusModel;   //xi-_bar LE inel model
289     G4HEAntiXiMinusInelastic*      fHEAntiXiMinusModel;   //xi-_bar HE inel model
290     G4MultipleScattering           fAntiXiMinusMult;      //xi-_bar msc
291     G4hIonisation                  fAntiXiMinusIonisation;//xi-_bar ionisation
292   
293       // OmegaMinus
294     G4OmegaMinusInelasticProcess  fOmegaMinusInelastic;   //omega- inel process
295     G4LEOmegaMinusInelastic*      fLEOmegaMinusModel;     //omega- LE inel model
296     G4HEOmegaMinusInelastic*      fHEOmegaMinusModel;     //omega- HE inel model
297     G4MultipleScattering          fOmegaMinusMult;        //omega- msc
298     G4hIonisation                 fOmegaMinusIonisation;  //omega- ionisation
299    
300       // AntiOmegaMinus
301     G4AntiOmegaMinusInelasticProcess  fAntiOmegaMinusInelastic; //omega-_bar inel process
302     G4LEAntiOmegaMinusInelastic*      fLEAntiOmegaMinusModel;   //omega-_bar LE inel model
303     G4HEAntiOmegaMinusInelastic*      fHEAntiOmegaMinusModel;   //omega-_bar HE inel model
304     G4MultipleScattering              fAntiOmegaMinusMult;      //omega-_bar msc
305     G4hIonisation                     fAntiOmegaMinusIonisation;//omega-_bar ionisation
306     
307       // Other
308     ProcessVector  fOtherProcesses; //other process
309     
310
311   private:
312     // methods
313
314     // EM processes
315     void ConstructEMProcessForPionPlus();
316     void ConstructEMProcessForPionMinus();
317     void ConstructEMProcessForKaonPlus();
318     void ConstructEMProcessForKaonMinus();
319     void ConstructEMProcessForProton();
320     void ConstructEMProcessForAntiProton();
321     void ConstructEMProcessForSigmaMinus();
322     void ConstructEMProcessForAntiSigmaMinus();
323     void ConstructEMProcessForSigmaPlus();
324     void ConstructEMProcessForAntiSigmaPlus();
325     void ConstructEMProcessForXiMinus();
326     void ConstructEMProcessForAntiXiMinus();
327     void ConstructEMProcessForOmegaMinus();
328     void ConstructEMProcessForAntiOmegaMinus();
329     void ConstructEMProcessForOther();
330
331     // hadron processes
332     void ConstructHadProcessForPionPlus();
333     void ConstructHadProcessForPionMinus();
334     void ConstructHadProcessForKaonPlus();
335     void ConstructHadProcessForKaonMinus();
336     void ConstructHadProcessForKaonZeroLong();
337     void ConstructHadProcessForKaonZeroShort();
338     void ConstructHadProcessForProton();
339     void ConstructHadProcessForAntiProton();
340     void ConstructHadProcessForNeutron();
341     void ConstructHadProcessForAntiNeutron();
342     void ConstructHadProcessForLambda();
343     void ConstructHadProcessForAntiLambda();
344     void ConstructHadProcessForSigmaMinus();
345     void ConstructHadProcessForAntiSigmaMinus();
346     void ConstructHadProcessForSigmaPlus();
347     void ConstructHadProcessForAntiSigmaPlus();
348     void ConstructHadProcessForXiMinus();
349     void ConstructHadProcessForAntiXiMinus();
350     void ConstructHadProcessForXiZero();
351     void ConstructHadProcessForAntiXiZero();
352     void ConstructHadProcessForOmegaMinus();
353     void ConstructHadProcessForAntiOmegaMinus();
354     void ConstructHadProcessForOther();
355
356     // data members
357     G4bool  fSetEM;    //if true - EM processes are constructed
358     G4bool  fSetHadron;//if true - hadron processes are constructed
359 };
360
361 #endif //TG4_PHYSICS_CONSTRUCTOR_HADRON_H
362