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 | // $Id: RE01HadronPhysics.cc,v 1.2 2006/06/29 17:44:00 gunter Exp $ |
---|
28 | // GEANT4 tag $Name: $ |
---|
29 | // |
---|
30 | |
---|
31 | #include "RE01HadronPhysics.hh" |
---|
32 | #include "G4MesonConstructor.hh" |
---|
33 | #include "G4BaryonConstructor.hh" |
---|
34 | #include "G4ShortLivedConstructor.hh" |
---|
35 | |
---|
36 | // processes |
---|
37 | #include "G4ParticleTable.hh" |
---|
38 | #include "G4ProcessManager.hh" |
---|
39 | #include "G4MultipleScattering.hh" |
---|
40 | #include "G4hIonisation.hh" |
---|
41 | #include "G4HadronElasticProcess.hh" |
---|
42 | #include "G4PionPlusInelasticProcess.hh" |
---|
43 | #include "G4PionMinusInelasticProcess.hh" |
---|
44 | #include "G4KaonPlusInelasticProcess.hh" |
---|
45 | #include "G4KaonMinusInelasticProcess.hh" |
---|
46 | #include "G4KaonZeroLInelasticProcess.hh" |
---|
47 | #include "G4KaonZeroSInelasticProcess.hh" |
---|
48 | #include "G4ProtonInelasticProcess.hh" |
---|
49 | #include "G4NeutronInelasticProcess.hh" |
---|
50 | #include "G4HadronFissionProcess.hh" |
---|
51 | #include "G4HadronCaptureProcess.hh" |
---|
52 | #include "G4AntiProtonInelasticProcess.hh" |
---|
53 | #include "G4AntiNeutronInelasticProcess.hh" |
---|
54 | #include "G4LambdaInelasticProcess.hh" |
---|
55 | #include "G4AntiLambdaInelasticProcess.hh" |
---|
56 | #include "G4SigmaPlusInelasticProcess.hh" |
---|
57 | #include "G4SigmaMinusInelasticProcess.hh" |
---|
58 | #include "G4AntiSigmaPlusInelasticProcess.hh" |
---|
59 | #include "G4AntiSigmaMinusInelasticProcess.hh" |
---|
60 | #include "G4XiZeroInelasticProcess.hh" |
---|
61 | #include "G4XiMinusInelasticProcess.hh" |
---|
62 | #include "G4AntiXiZeroInelasticProcess.hh" |
---|
63 | #include "G4AntiXiMinusInelasticProcess.hh" |
---|
64 | #include "G4OmegaMinusInelasticProcess.hh" |
---|
65 | #include "G4AntiOmegaMinusInelasticProcess.hh" |
---|
66 | #include "G4AntiProtonAnnihilationAtRest.hh" |
---|
67 | #include "G4AntiNeutronAnnihilationAtRest.hh" |
---|
68 | #include "G4PionMinusAbsorptionAtRest.hh" |
---|
69 | #include "G4KaonMinusAbsorption.hh" |
---|
70 | |
---|
71 | // cross sections |
---|
72 | #include "G4PiNuclearCrossSection.hh" |
---|
73 | #include "G4ProtonInelasticCrossSection.hh" |
---|
74 | #include "G4NeutronInelasticCrossSection.hh" |
---|
75 | |
---|
76 | // models |
---|
77 | #include "G4LElastic.hh" |
---|
78 | #include "G4CascadeInterface.hh" |
---|
79 | #include "G4LEPionPlusInelastic.hh" |
---|
80 | #include "G4LEPionMinusInelastic.hh" |
---|
81 | #include "G4LEKaonPlusInelastic.hh" |
---|
82 | #include "G4LEKaonMinusInelastic.hh" |
---|
83 | #include "G4LEKaonZeroLInelastic.hh" |
---|
84 | #include "G4LEKaonZeroSInelastic.hh" |
---|
85 | #include "G4LEProtonInelastic.hh" |
---|
86 | #include "G4LENeutronInelastic.hh" |
---|
87 | #include "G4LFission.hh" |
---|
88 | #include "G4LCapture.hh" |
---|
89 | #include "G4LEAntiProtonInelastic.hh" |
---|
90 | #include "G4LEAntiNeutronInelastic.hh" |
---|
91 | #include "G4LELambdaInelastic.hh" |
---|
92 | #include "G4LEAntiLambdaInelastic.hh" |
---|
93 | #include "G4LESigmaPlusInelastic.hh" |
---|
94 | #include "G4LESigmaMinusInelastic.hh" |
---|
95 | #include "G4LEAntiSigmaPlusInelastic.hh" |
---|
96 | #include "G4LEAntiSigmaMinusInelastic.hh" |
---|
97 | #include "G4LEXiZeroInelastic.hh" |
---|
98 | #include "G4LEXiMinusInelastic.hh" |
---|
99 | #include "G4LEAntiXiZeroInelastic.hh" |
---|
100 | #include "G4LEAntiXiMinusInelastic.hh" |
---|
101 | #include "G4LEOmegaMinusInelastic.hh" |
---|
102 | #include "G4LEAntiOmegaMinusInelastic.hh" |
---|
103 | |
---|
104 | #include "G4HEAntiProtonInelastic.hh" |
---|
105 | #include "G4HEAntiNeutronInelastic.hh" |
---|
106 | #include "G4HELambdaInelastic.hh" |
---|
107 | #include "G4HEAntiLambdaInelastic.hh" |
---|
108 | #include "G4HESigmaPlusInelastic.hh" |
---|
109 | #include "G4HESigmaMinusInelastic.hh" |
---|
110 | #include "G4HEAntiSigmaPlusInelastic.hh" |
---|
111 | #include "G4HEAntiSigmaMinusInelastic.hh" |
---|
112 | #include "G4HEXiZeroInelastic.hh" |
---|
113 | #include "G4HEXiMinusInelastic.hh" |
---|
114 | #include "G4HEAntiXiZeroInelastic.hh" |
---|
115 | #include "G4HEAntiXiMinusInelastic.hh" |
---|
116 | #include "G4HEOmegaMinusInelastic.hh" |
---|
117 | #include "G4HEAntiOmegaMinusInelastic.hh" |
---|
118 | |
---|
119 | #include "G4TheoFSGenerator.hh" |
---|
120 | #include "G4GeneratorPrecompoundInterface.hh" |
---|
121 | #include "G4ExcitationHandler.hh" |
---|
122 | #include "G4PreCompoundModel.hh" |
---|
123 | #include "G4QGSMFragmentation.hh" |
---|
124 | #include "G4ExcitedStringDecay.hh" |
---|
125 | #include "G4QGSParticipants.hh" |
---|
126 | #include "G4QGSModel.hh" |
---|
127 | |
---|
128 | |
---|
129 | RE01HadronPhysics::RE01HadronPhysics(const G4String& name) |
---|
130 | :G4VPhysicsConstructor(name) |
---|
131 | {;} |
---|
132 | |
---|
133 | RE01HadronPhysics::~RE01HadronPhysics() |
---|
134 | {;} |
---|
135 | |
---|
136 | |
---|
137 | void RE01HadronPhysics::ConstructParticle() |
---|
138 | { |
---|
139 | // Construct all mesons |
---|
140 | G4MesonConstructor pMesonConstructor; |
---|
141 | pMesonConstructor.ConstructParticle(); |
---|
142 | |
---|
143 | // Construct all baryons |
---|
144 | G4BaryonConstructor pBaryonConstructor; |
---|
145 | pBaryonConstructor.ConstructParticle(); |
---|
146 | |
---|
147 | // Construct resonaces and quarks |
---|
148 | G4ShortLivedConstructor pShortLivedConstructor; |
---|
149 | pShortLivedConstructor.ConstructParticle(); |
---|
150 | } |
---|
151 | |
---|
152 | |
---|
153 | void RE01HadronPhysics::ConstructProcess() |
---|
154 | { |
---|
155 | // Hadronic Elastic Process and Model (the same for all hadrons) |
---|
156 | |
---|
157 | G4HadronElasticProcess* elasticProcess = new G4HadronElasticProcess(); |
---|
158 | G4LElastic* elasticModel = new G4LElastic(); |
---|
159 | elasticProcess->RegisterMe(elasticModel); |
---|
160 | |
---|
161 | // Hadronic inelastic models |
---|
162 | |
---|
163 | // Bertini cascade model: use for p,n,pi+,pi- between 0 and 9.9 GeV |
---|
164 | G4CascadeInterface* bertiniModel = new G4CascadeInterface(); |
---|
165 | bertiniModel->SetMaxEnergy(9.9*GeV); |
---|
166 | |
---|
167 | // Low energy parameterized models : use between 9.5 and 25 GeV |
---|
168 | G4double LEPUpperLimit = 25*GeV; |
---|
169 | G4double LEPpnpiLimit = 9.5*GeV; |
---|
170 | |
---|
171 | G4LEKaonZeroLInelastic* LEPk0LModel = new G4LEKaonZeroLInelastic(); |
---|
172 | LEPk0LModel->SetMaxEnergy(LEPUpperLimit); |
---|
173 | |
---|
174 | G4LEKaonZeroSInelastic* LEPk0SModel = new G4LEKaonZeroSInelastic(); |
---|
175 | LEPk0SModel->SetMaxEnergy(LEPUpperLimit); |
---|
176 | |
---|
177 | // Quark-Gluon String Model: use for p,n,pi+,pi- between 12 GeV and 100 TeV |
---|
178 | |
---|
179 | G4TheoFSGenerator* QGSPModel = new G4TheoFSGenerator(); |
---|
180 | G4GeneratorPrecompoundInterface* theCascade = |
---|
181 | new G4GeneratorPrecompoundInterface(); |
---|
182 | G4ExcitationHandler* exHandler = new G4ExcitationHandler(); |
---|
183 | G4PreCompoundModel* preCompound = new G4PreCompoundModel(exHandler); |
---|
184 | theCascade->SetDeExcitation(preCompound); |
---|
185 | QGSPModel->SetTransport(theCascade); |
---|
186 | G4QGSMFragmentation* frag = new G4QGSMFragmentation(); |
---|
187 | G4ExcitedStringDecay* stringDecay = new G4ExcitedStringDecay(frag); |
---|
188 | G4QGSModel<G4QGSParticipants>* stringModel = |
---|
189 | new G4QGSModel<G4QGSParticipants>(); |
---|
190 | stringModel->SetFragmentationModel(stringDecay); |
---|
191 | QGSPModel->SetHighEnergyGenerator(stringModel); |
---|
192 | QGSPModel->SetMinEnergy(12*GeV); |
---|
193 | QGSPModel->SetMaxEnergy(100*TeV); |
---|
194 | |
---|
195 | // |
---|
196 | G4ProcessManager * pManager = 0; |
---|
197 | |
---|
198 | /////////////////// |
---|
199 | // // |
---|
200 | // pi+ physics // |
---|
201 | // // |
---|
202 | /////////////////// |
---|
203 | |
---|
204 | pManager = G4PionPlus::PionPlus()->GetProcessManager(); |
---|
205 | |
---|
206 | // EM processes |
---|
207 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
208 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
209 | |
---|
210 | // hadron elastic |
---|
211 | pManager->AddDiscreteProcess(elasticProcess); |
---|
212 | |
---|
213 | // hadron inelastic |
---|
214 | G4PionPlusInelasticProcess* pipinelProc = new G4PionPlusInelasticProcess(); |
---|
215 | G4PiNuclearCrossSection* pion_XC = new G4PiNuclearCrossSection(); |
---|
216 | pipinelProc->AddDataSet(pion_XC); |
---|
217 | pipinelProc->RegisterMe(bertiniModel); |
---|
218 | |
---|
219 | G4LEPionPlusInelastic* LEPpipModel = new G4LEPionPlusInelastic(); |
---|
220 | LEPpipModel->SetMinEnergy(LEPpnpiLimit); |
---|
221 | LEPpipModel->SetMaxEnergy(LEPUpperLimit); |
---|
222 | pipinelProc->RegisterMe(LEPpipModel); |
---|
223 | |
---|
224 | pipinelProc->RegisterMe(QGSPModel); |
---|
225 | pManager->AddDiscreteProcess(pipinelProc); |
---|
226 | |
---|
227 | /////////////////// |
---|
228 | // // |
---|
229 | // pi- physics // |
---|
230 | // // |
---|
231 | /////////////////// |
---|
232 | |
---|
233 | pManager = G4PionMinus::PionMinus()->GetProcessManager(); |
---|
234 | |
---|
235 | // EM processes |
---|
236 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
237 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
238 | |
---|
239 | // hadron elastic |
---|
240 | pManager->AddDiscreteProcess(elasticProcess); |
---|
241 | |
---|
242 | // hadron inelastic |
---|
243 | G4PionMinusInelasticProcess* piminelProc = new G4PionMinusInelasticProcess(); |
---|
244 | piminelProc->AddDataSet(pion_XC); |
---|
245 | piminelProc->RegisterMe(bertiniModel); |
---|
246 | |
---|
247 | G4LEPionMinusInelastic* LEPpimModel = new G4LEPionMinusInelastic(); |
---|
248 | LEPpimModel->SetMinEnergy(LEPpnpiLimit); |
---|
249 | LEPpimModel->SetMaxEnergy(LEPUpperLimit); |
---|
250 | piminelProc->RegisterMe(LEPpimModel); |
---|
251 | |
---|
252 | piminelProc->RegisterMe(QGSPModel); |
---|
253 | pManager->AddDiscreteProcess(piminelProc); |
---|
254 | |
---|
255 | // pi- absorption at rest |
---|
256 | G4PionMinusAbsorptionAtRest* pimAbsorb = new G4PionMinusAbsorptionAtRest(); |
---|
257 | pManager->AddRestProcess(pimAbsorb); |
---|
258 | |
---|
259 | /////////////////// |
---|
260 | // // |
---|
261 | // K+ physics // |
---|
262 | // // |
---|
263 | /////////////////// |
---|
264 | |
---|
265 | pManager = G4KaonPlus::KaonPlus()->GetProcessManager(); |
---|
266 | |
---|
267 | // EM processes |
---|
268 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
269 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
270 | |
---|
271 | // hadron elastic |
---|
272 | pManager->AddDiscreteProcess(elasticProcess); |
---|
273 | |
---|
274 | // hadron inelastic |
---|
275 | G4KaonPlusInelasticProcess* kpinelProc = new G4KaonPlusInelasticProcess(); |
---|
276 | G4LEKaonPlusInelastic* LEPkpModel = new G4LEKaonPlusInelastic(); |
---|
277 | LEPkpModel->SetMaxEnergy(LEPUpperLimit); |
---|
278 | kpinelProc->RegisterMe(LEPkpModel); |
---|
279 | kpinelProc->RegisterMe(QGSPModel); |
---|
280 | pManager->AddDiscreteProcess(kpinelProc); |
---|
281 | |
---|
282 | /////////////////// |
---|
283 | // // |
---|
284 | // K- physics // |
---|
285 | // // |
---|
286 | /////////////////// |
---|
287 | |
---|
288 | pManager = G4KaonMinus::KaonMinus()->GetProcessManager(); |
---|
289 | |
---|
290 | // EM processes |
---|
291 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
292 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
293 | |
---|
294 | // hadron elastic |
---|
295 | pManager->AddDiscreteProcess(elasticProcess); |
---|
296 | |
---|
297 | // hadron inelastic |
---|
298 | G4KaonMinusInelasticProcess* kminelProc = new G4KaonMinusInelasticProcess(); |
---|
299 | G4LEKaonMinusInelastic* LEPkmModel = new G4LEKaonMinusInelastic(); |
---|
300 | LEPkmModel->SetMaxEnergy(LEPUpperLimit); |
---|
301 | kminelProc->RegisterMe(LEPkmModel); |
---|
302 | kminelProc->RegisterMe(QGSPModel); |
---|
303 | pManager->AddDiscreteProcess(kminelProc); |
---|
304 | |
---|
305 | // K- absorption at rest |
---|
306 | G4KaonMinusAbsorption* kmAbsorb = new G4KaonMinusAbsorption(); |
---|
307 | pManager->AddRestProcess(kmAbsorb); |
---|
308 | |
---|
309 | /////////////////// |
---|
310 | // // |
---|
311 | // K0L physics // |
---|
312 | // // |
---|
313 | /////////////////// |
---|
314 | |
---|
315 | pManager = G4KaonZeroLong::KaonZeroLong()->GetProcessManager(); |
---|
316 | |
---|
317 | // hadron elastic |
---|
318 | pManager->AddDiscreteProcess(elasticProcess); |
---|
319 | |
---|
320 | // hadron inelastic |
---|
321 | G4KaonZeroLInelasticProcess* k0LinelProc = new G4KaonZeroLInelasticProcess(); |
---|
322 | k0LinelProc->RegisterMe(LEPk0LModel); |
---|
323 | k0LinelProc->RegisterMe(QGSPModel); |
---|
324 | pManager->AddDiscreteProcess(k0LinelProc); |
---|
325 | |
---|
326 | /////////////////// |
---|
327 | // // |
---|
328 | // K0S physics // |
---|
329 | // // |
---|
330 | /////////////////// |
---|
331 | |
---|
332 | pManager = G4KaonZeroShort::KaonZeroShort()->GetProcessManager(); |
---|
333 | |
---|
334 | // hadron elastic |
---|
335 | pManager->AddDiscreteProcess(elasticProcess); |
---|
336 | |
---|
337 | // hadron inelastic |
---|
338 | G4KaonZeroSInelasticProcess* k0SinelProc = new G4KaonZeroSInelasticProcess(); |
---|
339 | k0SinelProc->RegisterMe(LEPk0SModel); |
---|
340 | k0SinelProc->RegisterMe(QGSPModel); |
---|
341 | pManager->AddDiscreteProcess(k0SinelProc); |
---|
342 | |
---|
343 | /////////////////// |
---|
344 | // // |
---|
345 | // Proton // |
---|
346 | // // |
---|
347 | /////////////////// |
---|
348 | |
---|
349 | pManager = G4Proton::Proton()->GetProcessManager(); |
---|
350 | |
---|
351 | // EM processes |
---|
352 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
353 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
354 | |
---|
355 | // hadron elastic |
---|
356 | pManager->AddDiscreteProcess(elasticProcess); |
---|
357 | |
---|
358 | // hadron inelastic |
---|
359 | G4ProtonInelasticProcess* pinelProc = new G4ProtonInelasticProcess(); |
---|
360 | G4ProtonInelasticCrossSection* proton_XC = |
---|
361 | new G4ProtonInelasticCrossSection(); |
---|
362 | pinelProc->AddDataSet(proton_XC); |
---|
363 | pinelProc->RegisterMe(bertiniModel); |
---|
364 | |
---|
365 | G4LEProtonInelastic* LEPpModel = new G4LEProtonInelastic(); |
---|
366 | LEPpModel->SetMinEnergy(LEPpnpiLimit); |
---|
367 | LEPpModel->SetMaxEnergy(LEPUpperLimit); |
---|
368 | pinelProc->RegisterMe(LEPpModel); |
---|
369 | |
---|
370 | pinelProc->RegisterMe(QGSPModel); |
---|
371 | pManager->AddDiscreteProcess(pinelProc); |
---|
372 | |
---|
373 | /////////////////// |
---|
374 | // // |
---|
375 | // Anti-Proton // |
---|
376 | // // |
---|
377 | /////////////////// |
---|
378 | |
---|
379 | pManager = G4AntiProton::AntiProton()->GetProcessManager(); |
---|
380 | |
---|
381 | // EM processes |
---|
382 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
383 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
384 | |
---|
385 | // hadron elastic |
---|
386 | pManager->AddDiscreteProcess(elasticProcess); |
---|
387 | |
---|
388 | // hadron inelastic |
---|
389 | G4AntiProtonInelasticProcess* apinelProc = |
---|
390 | new G4AntiProtonInelasticProcess(); |
---|
391 | G4LEAntiProtonInelastic* LEPapModel = new G4LEAntiProtonInelastic(); |
---|
392 | apinelProc->RegisterMe(LEPapModel); |
---|
393 | G4HEAntiProtonInelastic* HEPapModel = new G4HEAntiProtonInelastic(); |
---|
394 | apinelProc->RegisterMe(HEPapModel); |
---|
395 | pManager->AddDiscreteProcess(apinelProc); |
---|
396 | |
---|
397 | // anti-proton annihilation at rest |
---|
398 | G4AntiProtonAnnihilationAtRest* apAnnihil = |
---|
399 | new G4AntiProtonAnnihilationAtRest(); |
---|
400 | pManager->AddRestProcess(apAnnihil); |
---|
401 | |
---|
402 | /////////////////// |
---|
403 | // // |
---|
404 | // Neutron // |
---|
405 | // // |
---|
406 | /////////////////// |
---|
407 | |
---|
408 | pManager = G4Neutron::Neutron()->GetProcessManager(); |
---|
409 | |
---|
410 | // hadron elastic |
---|
411 | pManager->AddDiscreteProcess(elasticProcess); |
---|
412 | |
---|
413 | // hadron inelastic |
---|
414 | G4NeutronInelasticProcess* ninelProc = new G4NeutronInelasticProcess(); |
---|
415 | G4NeutronInelasticCrossSection* neutron_XC = |
---|
416 | new G4NeutronInelasticCrossSection(); |
---|
417 | ninelProc->AddDataSet(neutron_XC); |
---|
418 | ninelProc->RegisterMe(bertiniModel); |
---|
419 | |
---|
420 | G4LENeutronInelastic* LEPnModel = new G4LENeutronInelastic(); |
---|
421 | LEPnModel->SetMinEnergy(LEPpnpiLimit); |
---|
422 | LEPnModel->SetMaxEnergy(LEPUpperLimit); |
---|
423 | ninelProc->RegisterMe(LEPnModel); |
---|
424 | |
---|
425 | ninelProc->RegisterMe(QGSPModel); |
---|
426 | pManager->AddDiscreteProcess(ninelProc); |
---|
427 | |
---|
428 | // neutron-induced fission |
---|
429 | G4HadronFissionProcess* neutronFission = new G4HadronFissionProcess(); |
---|
430 | G4LFission* neutronFissionModel = new G4LFission(); |
---|
431 | neutronFissionModel->SetMinEnergy(0.); |
---|
432 | neutronFissionModel->SetMaxEnergy(20*TeV); |
---|
433 | neutronFission->RegisterMe(neutronFissionModel); |
---|
434 | pManager->AddDiscreteProcess(neutronFission); |
---|
435 | |
---|
436 | // neutron capture |
---|
437 | G4HadronCaptureProcess* neutronCapture = new G4HadronCaptureProcess(); |
---|
438 | G4LCapture* neutronCaptureModel = new G4LCapture(); |
---|
439 | neutronCaptureModel->SetMinEnergy(0.); |
---|
440 | neutronCaptureModel->SetMaxEnergy(20*TeV); |
---|
441 | neutronCapture->RegisterMe(neutronCaptureModel); |
---|
442 | pManager->AddDiscreteProcess(neutronCapture); |
---|
443 | |
---|
444 | /////////////////// |
---|
445 | // // |
---|
446 | // Anti-Neutron // |
---|
447 | // // |
---|
448 | /////////////////// |
---|
449 | |
---|
450 | pManager = G4AntiNeutron::AntiNeutron()->GetProcessManager(); |
---|
451 | |
---|
452 | // hadron elastic |
---|
453 | pManager->AddDiscreteProcess(elasticProcess); |
---|
454 | |
---|
455 | // hadron inelastic |
---|
456 | G4AntiNeutronInelasticProcess* aninelProc = |
---|
457 | new G4AntiNeutronInelasticProcess(); |
---|
458 | G4LEAntiNeutronInelastic* LEPanModel = new G4LEAntiNeutronInelastic(); |
---|
459 | aninelProc->RegisterMe(LEPanModel); |
---|
460 | G4HEAntiNeutronInelastic* HEPanModel = new G4HEAntiNeutronInelastic(); |
---|
461 | aninelProc->RegisterMe(HEPanModel); |
---|
462 | pManager->AddDiscreteProcess(aninelProc); |
---|
463 | |
---|
464 | // anti-neutron annihilation at rest |
---|
465 | G4AntiNeutronAnnihilationAtRest* anAnnihil = |
---|
466 | new G4AntiNeutronAnnihilationAtRest(); |
---|
467 | pManager->AddRestProcess(anAnnihil); |
---|
468 | |
---|
469 | /////////////////// |
---|
470 | // // |
---|
471 | // Lambda // |
---|
472 | // // |
---|
473 | /////////////////// |
---|
474 | |
---|
475 | pManager = G4Lambda::Lambda()->GetProcessManager(); |
---|
476 | |
---|
477 | // hadron elastic |
---|
478 | pManager->AddDiscreteProcess(elasticProcess); |
---|
479 | |
---|
480 | // hadron inelastic |
---|
481 | G4LambdaInelasticProcess* linelProc = |
---|
482 | new G4LambdaInelasticProcess(); |
---|
483 | G4LELambdaInelastic* LEPlModel = new G4LELambdaInelastic(); |
---|
484 | linelProc->RegisterMe(LEPlModel); |
---|
485 | G4HELambdaInelastic* HEPlModel = new G4HELambdaInelastic(); |
---|
486 | linelProc->RegisterMe(HEPlModel); |
---|
487 | pManager->AddDiscreteProcess(linelProc); |
---|
488 | |
---|
489 | /////////////////// |
---|
490 | // // |
---|
491 | // Anti-Lambda // |
---|
492 | // // |
---|
493 | /////////////////// |
---|
494 | |
---|
495 | pManager = G4AntiLambda::AntiLambda()->GetProcessManager(); |
---|
496 | |
---|
497 | // hadron elastic |
---|
498 | pManager->AddDiscreteProcess(elasticProcess); |
---|
499 | |
---|
500 | // hadron inelastic |
---|
501 | G4AntiLambdaInelasticProcess* alinelProc = |
---|
502 | new G4AntiLambdaInelasticProcess(); |
---|
503 | G4LEAntiLambdaInelastic* LEPalModel = new G4LEAntiLambdaInelastic(); |
---|
504 | alinelProc->RegisterMe(LEPalModel); |
---|
505 | G4HEAntiLambdaInelastic* HEPalModel = new G4HEAntiLambdaInelastic(); |
---|
506 | alinelProc->RegisterMe(HEPalModel); |
---|
507 | pManager->AddDiscreteProcess(alinelProc); |
---|
508 | |
---|
509 | /////////////////// |
---|
510 | // // |
---|
511 | // Sigma- // |
---|
512 | // // |
---|
513 | /////////////////// |
---|
514 | |
---|
515 | pManager = G4SigmaMinus::SigmaMinus()->GetProcessManager(); |
---|
516 | |
---|
517 | // EM processes |
---|
518 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
519 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
520 | |
---|
521 | // hadron elastic |
---|
522 | pManager->AddDiscreteProcess(elasticProcess); |
---|
523 | |
---|
524 | // hadron inelastic |
---|
525 | G4SigmaMinusInelasticProcess* sminelProc = |
---|
526 | new G4SigmaMinusInelasticProcess(); |
---|
527 | G4LESigmaMinusInelastic* LEPsmModel = new G4LESigmaMinusInelastic(); |
---|
528 | sminelProc->RegisterMe(LEPsmModel); |
---|
529 | G4HESigmaMinusInelastic* HEPsmModel = new G4HESigmaMinusInelastic(); |
---|
530 | sminelProc->RegisterMe(HEPsmModel); |
---|
531 | pManager->AddDiscreteProcess(sminelProc); |
---|
532 | |
---|
533 | /////////////////// |
---|
534 | // // |
---|
535 | // Anti-Sigma- // |
---|
536 | // // |
---|
537 | /////////////////// |
---|
538 | |
---|
539 | pManager = G4AntiSigmaMinus::AntiSigmaMinus()->GetProcessManager(); |
---|
540 | |
---|
541 | // EM processes |
---|
542 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
543 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
544 | |
---|
545 | // hadron elastic |
---|
546 | pManager->AddDiscreteProcess(elasticProcess); |
---|
547 | |
---|
548 | // hadron inelastic |
---|
549 | G4AntiSigmaMinusInelasticProcess* asminelProc = |
---|
550 | new G4AntiSigmaMinusInelasticProcess(); |
---|
551 | G4LEAntiSigmaMinusInelastic* LEPasmModel = |
---|
552 | new G4LEAntiSigmaMinusInelastic(); |
---|
553 | asminelProc->RegisterMe(LEPasmModel); |
---|
554 | G4HEAntiSigmaMinusInelastic* HEPasmModel = |
---|
555 | new G4HEAntiSigmaMinusInelastic(); |
---|
556 | asminelProc->RegisterMe(HEPasmModel); |
---|
557 | pManager->AddDiscreteProcess(asminelProc); |
---|
558 | |
---|
559 | /////////////////// |
---|
560 | // // |
---|
561 | // Sigma+ // |
---|
562 | // // |
---|
563 | /////////////////// |
---|
564 | |
---|
565 | pManager = G4SigmaPlus::SigmaPlus()->GetProcessManager(); |
---|
566 | |
---|
567 | // EM processes |
---|
568 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
569 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
570 | |
---|
571 | // hadron elastic |
---|
572 | pManager->AddDiscreteProcess(elasticProcess); |
---|
573 | |
---|
574 | // hadron inelastic |
---|
575 | G4SigmaPlusInelasticProcess* spinelProc = new G4SigmaPlusInelasticProcess(); |
---|
576 | G4LESigmaPlusInelastic* LEPspModel = new G4LESigmaPlusInelastic(); |
---|
577 | spinelProc->RegisterMe(LEPspModel); |
---|
578 | G4HESigmaPlusInelastic* HEPspModel = new G4HESigmaPlusInelastic(); |
---|
579 | spinelProc->RegisterMe(HEPspModel); |
---|
580 | pManager->AddDiscreteProcess(spinelProc); |
---|
581 | |
---|
582 | /////////////////// |
---|
583 | // // |
---|
584 | // Anti-Sigma+ // |
---|
585 | // // |
---|
586 | /////////////////// |
---|
587 | |
---|
588 | pManager = G4AntiSigmaPlus::AntiSigmaPlus()->GetProcessManager(); |
---|
589 | |
---|
590 | // EM processes |
---|
591 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
592 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
593 | |
---|
594 | // hadron elastic |
---|
595 | pManager->AddDiscreteProcess(elasticProcess); |
---|
596 | |
---|
597 | // hadron inelastic |
---|
598 | G4AntiSigmaPlusInelasticProcess* aspinelProc = |
---|
599 | new G4AntiSigmaPlusInelasticProcess(); |
---|
600 | G4LEAntiSigmaPlusInelastic* LEPaspModel = |
---|
601 | new G4LEAntiSigmaPlusInelastic(); |
---|
602 | aspinelProc->RegisterMe(LEPaspModel); |
---|
603 | G4HEAntiSigmaPlusInelastic* HEPaspModel = |
---|
604 | new G4HEAntiSigmaPlusInelastic(); |
---|
605 | aspinelProc->RegisterMe(HEPaspModel); |
---|
606 | pManager->AddDiscreteProcess(aspinelProc); |
---|
607 | |
---|
608 | /////////////////// |
---|
609 | // // |
---|
610 | // Xi- // |
---|
611 | // // |
---|
612 | /////////////////// |
---|
613 | |
---|
614 | pManager = G4XiMinus::XiMinus()->GetProcessManager(); |
---|
615 | |
---|
616 | // EM processes |
---|
617 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
618 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
619 | |
---|
620 | // hadron elastic |
---|
621 | pManager->AddDiscreteProcess(elasticProcess); |
---|
622 | |
---|
623 | // hadron inelastic |
---|
624 | G4XiMinusInelasticProcess* xminelProc = new G4XiMinusInelasticProcess(); |
---|
625 | G4LEXiMinusInelastic* LEPxmModel = new G4LEXiMinusInelastic(); |
---|
626 | xminelProc->RegisterMe(LEPxmModel); |
---|
627 | G4HEXiMinusInelastic* HEPxmModel = new G4HEXiMinusInelastic(); |
---|
628 | xminelProc->RegisterMe(HEPxmModel); |
---|
629 | pManager->AddDiscreteProcess(xminelProc); |
---|
630 | |
---|
631 | /////////////////// |
---|
632 | // // |
---|
633 | // Anti-Xi- // |
---|
634 | // // |
---|
635 | /////////////////// |
---|
636 | |
---|
637 | pManager = G4AntiXiMinus::AntiXiMinus()->GetProcessManager(); |
---|
638 | |
---|
639 | // EM processes |
---|
640 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
641 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
642 | |
---|
643 | // hadron elastic |
---|
644 | pManager->AddDiscreteProcess(elasticProcess); |
---|
645 | |
---|
646 | // hadron inelastic |
---|
647 | G4AntiXiMinusInelasticProcess* axminelProc = |
---|
648 | new G4AntiXiMinusInelasticProcess(); |
---|
649 | G4LEAntiXiMinusInelastic* LEPaxmModel = new G4LEAntiXiMinusInelastic(); |
---|
650 | axminelProc->RegisterMe(LEPaxmModel); |
---|
651 | G4HEAntiXiMinusInelastic* HEPaxmModel = new G4HEAntiXiMinusInelastic(); |
---|
652 | axminelProc->RegisterMe(HEPaxmModel); |
---|
653 | pManager->AddDiscreteProcess(axminelProc); |
---|
654 | |
---|
655 | /////////////////// |
---|
656 | // // |
---|
657 | // Xi0 // |
---|
658 | // // |
---|
659 | /////////////////// |
---|
660 | |
---|
661 | pManager = G4XiZero::XiZero()->GetProcessManager(); |
---|
662 | |
---|
663 | // hadron elastic |
---|
664 | pManager->AddDiscreteProcess(elasticProcess); |
---|
665 | |
---|
666 | // hadron inelastic |
---|
667 | G4XiZeroInelasticProcess* x0inelProc = new G4XiZeroInelasticProcess(); |
---|
668 | G4LEXiZeroInelastic* LEPx0Model = new G4LEXiZeroInelastic(); |
---|
669 | x0inelProc->RegisterMe(LEPx0Model); |
---|
670 | G4HEXiZeroInelastic* HEPx0Model = new G4HEXiZeroInelastic(); |
---|
671 | x0inelProc->RegisterMe(HEPx0Model); |
---|
672 | pManager->AddDiscreteProcess(x0inelProc); |
---|
673 | |
---|
674 | /////////////////// |
---|
675 | // // |
---|
676 | // Anti-Xi0 // |
---|
677 | // // |
---|
678 | /////////////////// |
---|
679 | |
---|
680 | pManager = G4AntiXiZero::AntiXiZero()->GetProcessManager(); |
---|
681 | |
---|
682 | // hadron elastic |
---|
683 | pManager->AddDiscreteProcess(elasticProcess); |
---|
684 | |
---|
685 | // hadron inelastic |
---|
686 | G4AntiXiZeroInelasticProcess* ax0inelProc = |
---|
687 | new G4AntiXiZeroInelasticProcess(); |
---|
688 | G4LEAntiXiZeroInelastic* LEPax0Model = new G4LEAntiXiZeroInelastic(); |
---|
689 | ax0inelProc->RegisterMe(LEPax0Model); |
---|
690 | G4HEAntiXiZeroInelastic* HEPax0Model = new G4HEAntiXiZeroInelastic(); |
---|
691 | ax0inelProc->RegisterMe(HEPax0Model); |
---|
692 | pManager->AddDiscreteProcess(ax0inelProc); |
---|
693 | |
---|
694 | /////////////////// |
---|
695 | // // |
---|
696 | // Omega- // |
---|
697 | // // |
---|
698 | /////////////////// |
---|
699 | |
---|
700 | pManager = G4OmegaMinus::OmegaMinus()->GetProcessManager(); |
---|
701 | |
---|
702 | // EM processes |
---|
703 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
704 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
705 | |
---|
706 | // hadron elastic |
---|
707 | pManager->AddDiscreteProcess(elasticProcess); |
---|
708 | |
---|
709 | // hadron inelastic |
---|
710 | G4OmegaMinusInelasticProcess* ominelProc = |
---|
711 | new G4OmegaMinusInelasticProcess(); |
---|
712 | G4LEOmegaMinusInelastic* LEPomModel = new G4LEOmegaMinusInelastic(); |
---|
713 | ominelProc->RegisterMe(LEPomModel); |
---|
714 | G4HEOmegaMinusInelastic* HEPomModel = new G4HEOmegaMinusInelastic(); |
---|
715 | ominelProc->RegisterMe(HEPomModel); |
---|
716 | pManager->AddDiscreteProcess(ominelProc); |
---|
717 | |
---|
718 | /////////////////// |
---|
719 | // // |
---|
720 | // Anti-Omega- // |
---|
721 | // // |
---|
722 | /////////////////// |
---|
723 | |
---|
724 | pManager = G4AntiOmegaMinus::AntiOmegaMinus()->GetProcessManager(); |
---|
725 | |
---|
726 | // EM processes |
---|
727 | pManager->AddProcess(new G4MultipleScattering(), -1, 1, 1); |
---|
728 | pManager->AddProcess(new G4hIonisation(), -1, 2, 2); |
---|
729 | |
---|
730 | // hadron elastic |
---|
731 | pManager->AddDiscreteProcess(elasticProcess); |
---|
732 | |
---|
733 | // hadron inelastic |
---|
734 | G4AntiOmegaMinusInelasticProcess* aominelProc = |
---|
735 | new G4AntiOmegaMinusInelasticProcess(); |
---|
736 | G4LEAntiOmegaMinusInelastic* LEPaomModel = |
---|
737 | new G4LEAntiOmegaMinusInelastic(); |
---|
738 | aominelProc->RegisterMe(LEPaomModel); |
---|
739 | G4HEAntiOmegaMinusInelastic* HEPaomModel = |
---|
740 | new G4HEAntiOmegaMinusInelastic(); |
---|
741 | aominelProc->RegisterMe(HEPaomModel); |
---|
742 | pManager->AddDiscreteProcess(aominelProc); |
---|
743 | |
---|
744 | } |
---|