SToGS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros
SToGS_HadronPhysicsList.cc
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 //
27 
29 
30 #include "globals.hh"
31 #include "G4ios.hh"
32 #include "G4Version.hh"
33 
34 #include <iomanip>
35 
36 // to deal with version dependant interfaces
37 #if G4VERSION_NUMBER < 1000
38 
39 //*******************************************************
40 //* GEANT4.9 VERSION *
41 //*******************************************************
42 
43 #include "G4LossTableManager.hh"
44 #include "G4ProcessManager.hh"
45 #include "G4ProcessVector.hh"
46 #include "G4PhysicsListHelper.hh"
47 
48 // PARTICLE DEFINITIONS
49 #include "G4ParticleDefinition.hh"
50 #include "G4ParticleWithCuts.hh"
51 #include "G4ParticleTypes.hh"
52 #include "G4ParticleTable.hh"
53 #include "G4BosonConstructor.hh"
54 #include "G4LeptonConstructor.hh"
55 #include "G4MesonConstructor.hh"
56 #include "G4BosonConstructor.hh"
57 #include "G4BaryonConstructor.hh"
58 #include "G4IonConstructor.hh"
59 #include "G4ShortLivedConstructor.hh"
60 #include "G4Neutron.hh"
61 
62 // PROCESSES FOR GAMMA
63 #include "G4PhotoElectricEffect.hh"
64 #include "G4LivermorePhotoElectricModel.hh"
65 #include "G4ComptonScattering.hh"
66 #include "G4LivermoreComptonModel.hh"
67 #include "G4GammaConversion.hh"
68 #include "G4LivermoreGammaConversionModel.hh"
69 #include "G4RayleighScattering.hh"
70 #include "G4LivermoreRayleighModel.hh"
71 
72 // PROCESSES FOR e-
73 #include "G4eMultipleScattering.hh"
74 #include "G4UniversalFluctuation.hh"
75 #include "G4eIonisation.hh"
76 #include "G4LivermoreIonisationModel.hh"
77 #include "G4eBremsstrahlung.hh"
78 #include "G4LivermoreBremsstrahlungModel.hh"
79 #include "G4eIonisation.hh"
80 #include "G4eBremsstrahlung.hh"
81 #include "G4eplusAnnihilation.hh"
82 #include "G4KleinNishinaModel.hh"
83 #include "G4eMultipleScattering.hh"
84 #include "G4MuMultipleScattering.hh"
85 #include "G4MuIonisation.hh"
86 #include "G4MuBremsstrahlung.hh"
87 #include "G4MuPairProduction.hh"
88 #include "G4hMultipleScattering.hh"
89 #include "G4hIonisation.hh"
90 #include "G4hBremsstrahlung.hh"
91 #include "G4hPairProduction.hh"
92 #include "G4ionIonisation.hh"
93 #include "G4IonParametrisedLossModel.hh"
94 #include "G4NuclearStopping.hh"
95 #include "G4SystemOfUnits.hh"
96 #include "G4HadronElasticPhysics.hh"
97 #include "G4HadronicProcess.hh"
98 #include "G4HadronElastic.hh"
99 #include "G4HadronElasticProcess.hh"
100 
101 //c-s
102 #include "G4TripathiCrossSection.hh"
103 #include "G4IonsShenCrossSection.hh"
104 #include "G4ProtonInelasticCrossSection.hh"
105 #include "G4NeutronInelasticCrossSection.hh"
106 
107 // RadioactiveDecay
108 #include "G4RadioactiveDecay.hh"
109 #include "G4GenericIon.hh"
110 
111 #include "G4UrbanMscModel96.hh"
112 
113 // hadronic
114 
115 #include "G4HadronElasticPhysicsHP.hh"
116 #include "HadronPhysicsFTFP_BERT_HP.hh"
117 
118 //HPNeutron
119 #include "G4NeutronHPThermalScatteringData.hh"
120 #include "G4NeutronHPThermalScattering.hh"
121 //
122 #include "G4HadronFissionProcess.hh"
123 #include "G4NeutronHPFissionData.hh"
124 #include "G4NeutronHPFission.hh"
125 //
126 #include "G4NeutronHPElastic.hh"
127 #include "G4NeutronHPElasticData.hh"
128 //
129 #include "G4HadronCaptureProcess.hh"
130 #include "G4NeutronHPCapture.hh"
131 #include "G4NeutronHPCaptureData.hh"
132 //
133 #include "G4NeutronInelasticProcess.hh"
134 #include "G4NeutronHPInelastic.hh"
135 #include "G4NeutronHPInelasticData.hh"
136 
137 
138 
139 /*
140 
141 
142 
143 // to deal with version dependant interfaces
144 #if G4VERSION_NUMBER < 1000
145 #if G4VERSION_NUMBER >= 960
146 #include "G4UrbanMscModel96.hh"
147 #define URBANMSCMODEL G4UrbanMscModel96
148 #else
149 #include "G4UrbanMscModel95.hh"
150 #define URBANMSCMODEL G4UrbanMscModel95
151 #endif
152 #else
153 #include "G4UrbanMscModel.hh"
154 #define URBANMSCMODEL G4UrbanMscModel
155 #endif
156 
157 // hadronic processes
158 #include "G4HadronElasticProcess.hh"
159 #include "G4HadronElasticPhysicsHP.hh"
160 #include "G4HadronElasticPhysics.hh"
161 #include "G4HadronicProcess.hh"
162 #include "G4HadronElastic.hh"
163 // to deal with version dependant interfaces
164 #if G4VERSION_NUMBER < 1000
165 #else
166 #include "FTFP_BERT_HP.hh"
167 #endif
168 */
169 
170 
171 
172 
173 SToGS::HadronPhysicsList::HadronPhysicsList(const G4String& name): G4VPhysicsConstructor(name)
174 {
175 }
176 
178 {
179  G4LossTableManager::Instance();
180  // defaultCutValue = 1.0*mm;
181 
182  G4cout << G4endl << " ------ INFO ------ You are working with an hadronic PhysicsList (transportation, electromagnetic and hadronic processes)" << G4endl;
183 }
184 
185 
187 {
188 
189  G4LeptonConstructor lepton;
190  lepton.ConstructParticle();
191 
192  G4BosonConstructor boson;
193  boson.ConstructParticle();
194 
195  G4MesonConstructor meson;
196  meson.ConstructParticle();
197 
198  G4BaryonConstructor baryon;
199  baryon.ConstructParticle();
200 
201  G4ShortLivedConstructor shortLived;
202  shortLived.ConstructParticle();
203 
204  G4IonConstructor ion;
205  ion.ConstructParticle();
206 
207 
208  G4cout << " ------ INFO ------ PARTICLES DEFINITION FINISHED" << G4endl;
209 }
210 
212 {
213 
214 
215 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
216 
217 
218 theParticleIterator->reset();
219 //
220 while( (*theParticleIterator)() ){
221  G4ParticleDefinition* particle = theParticleIterator->value();
222  G4ProcessManager* pmanager = particle->GetProcessManager();
223  G4String particleName = particle->GetParticleName();
224 
225 
226 
227 
228 
229  // define EM physics for gamma
230  if (particleName == "gamma") {
231 
232  // activate photoelectric effect
233  ph->RegisterProcess(new G4PhotoElectricEffect, particle);
234 
235  // activate Compton scattering
236  G4ComptonScattering* cs = new G4ComptonScattering;
237  cs->SetEmModel(new G4KleinNishinaModel());
238  ph->RegisterProcess(cs, particle);
239 
240  // activate gamma conversion
241  ph->RegisterProcess(new G4GammaConversion, particle);
242 
243  }
244 
245  // define EM physics for electron
246  else if (particleName == "e-") {
247 
248  // activate multiple scattering
249  G4eMultipleScattering* msc = new G4eMultipleScattering();
250  msc -> AddEmModel(0, new G4UrbanMscModel96());
251  ph->RegisterProcess(msc, particle);
252 
253  // activate ionisation
254  G4eIonisation* eIoni = new G4eIonisation();
255  eIoni->SetStepFunction(0.1, 100*um);
256  ph->RegisterProcess(eIoni, particle);
257 
258  // activate bremsstrahlung
259  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
260 
261  }
262 
263  // define EM physics for positron
264  else if (particleName == "e+") {
265 
266  // activate multiple scattering
267  G4eMultipleScattering* msc = new G4eMultipleScattering();
268  msc -> AddEmModel(0, new G4UrbanMscModel96());
269  ph->RegisterProcess(msc, particle);
270 
271  // activate ionisation
272  G4eIonisation* eIoni = new G4eIonisation();
273  eIoni->SetStepFunction(0.1, 100*um);
274  ph->RegisterProcess(eIoni, particle);
275 
276  // activate bremsstrahlung
277  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
278 
279  // activate annihilation
280  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
281 
282  }
283 
284  // define EM physics for muons
285  else if (particleName == "mu+" || particleName == "mu-") {
286 
287  // activate multiple scattering
288  ph->RegisterProcess(new G4MuMultipleScattering(), particle);
289 
290  // activate ionisation
291  G4MuIonisation* muIoni = new G4MuIonisation();
292  muIoni->SetStepFunction(0.1, 50*um);
293  ph->RegisterProcess(muIoni, particle);
294 
295  // activate bremsstrahlung
296  ph->RegisterProcess(new G4MuBremsstrahlung(), particle);
297 
298  // activate pair production
299  ph->RegisterProcess(new G4MuPairProduction(), particle);
300 
301  }
302 
303  // define EM physics for protons and pions
304  else if( particleName == "proton" ||
305  particleName == "pi-" ||
306  particleName == "pi+") {
307 
308  // activate multiple scattering
309  ph->RegisterProcess(new G4hMultipleScattering(), particle);
310 
311  // activate ionisation
312  G4hIonisation* hIoni = new G4hIonisation();
313  hIoni->SetStepFunction(0.1, 20*um);
314  ph->RegisterProcess(hIoni, particle);
315 
316  // activate bremsstrahlung
317  ph->RegisterProcess(new G4hBremsstrahlung(), particle);
318 
319  // activate pair production
320  ph->RegisterProcess(new G4hPairProduction(), particle);
321 
322  }
323 
324  // define EM physics for alpha and helium-3
325  else if( particleName == "alpha" ||
326  particleName == "He3") {
327 
328  // activate multiple scattering
329  ph->RegisterProcess(new G4hMultipleScattering(), particle);
330 
331  // activate ionisation
332  G4ionIonisation* ionIoni = new G4ionIonisation();
333  ionIoni->SetStepFunction(0.1, 1*um);
334  ph->RegisterProcess(ionIoni, particle);
335 
336  // activate nuclear stopping
337  ph->RegisterProcess(new G4NuclearStopping(), particle);
338 
339  }
340 
341  // define EM physics for all the other ions
342  else if( particleName == "GenericIon" ) {
343 
344  // activate multiple scattering
345  ph->RegisterProcess(new G4hMultipleScattering(), particle);
346 
347  // activate ionisation
348  G4ionIonisation* ionIoni = new G4ionIonisation();
349  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
350  ionIoni->SetStepFunction(0.1, 1*um);
351  ph->RegisterProcess(ionIoni, particle);
352 
353  // activate nuclear stopping
354  ph->RegisterProcess(new G4NuclearStopping(), particle);
355 
356  }
357 
358  // activate EM physics for all others charged particles except geantino
359  else if ((!particle->IsShortLived()) &&
360  (particle->GetPDGCharge() != 0.0) &&
361  (particle->GetParticleName() != "chargedgeantino")) {
362 
363  // activate multiple scattering
364  ph->RegisterProcess(new G4hMultipleScattering(), particle);
365 
366  // activate ionisation
367  ph->RegisterProcess(new G4hIonisation(), particle);
368  }
369  }
370 
371 
372 
373 
374  /* ***********************************************
375  * HADRONIC PHYSICS *
376  *********************************************** */
377 
378 
379  // define hadronic elastic processes
380  G4HadronElasticPhysicsHP* hadronicElastic = new G4HadronElasticPhysicsHP();
381  hadronicElastic->ConstructProcess();
382 
383  // define hadronic inelastic processes
384  HadronPhysicsFTFP_BERT_HP* hadronicPhysics = new HadronPhysicsFTFP_BERT_HP();
385  hadronicPhysics->ConstructProcess();
386 }
387 
389 {
390  /*
391  if ( verboseLevel > 0 ) {
392  G4cout << " HadronPhysicsList::SetCuts:";
393  G4cout << " CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
394  }
395 
396  // set cut values for gamma at first and for e- second and next for e+,
397  // because some processes for e+/e- need cut values for gamma
398  //
399  SetCutValue(defaultCutValue, "gamma");
400  SetCutValue(defaultCutValue, "e-");
401  SetCutValue(defaultCutValue, "e+");
402 
403  if (verboseLevel>0) DumpCutValuesTable();
404  */
405 }
406 
407 
408 
409 
410 #else
411 
412 //*******************************************************
413 //* GEANT4.10 VERSION *
414 //*******************************************************
415 
416 //#if G4VERSION_NUMBER >= 1000
417 #include "G4LossTableManager.hh"
418 #include "G4ProcessManager.hh"
419 #include "G4ProcessVector.hh"
420 #include "G4PhysicsListHelper.hh"
421 
422 // PARTICLE DEFINITIONS
423 #include "G4ParticleDefinition.hh"
424 #include "G4ParticleWithCuts.hh"
425 #include "G4ParticleTypes.hh"
426 #include "G4ParticleTable.hh"
427 #include "G4BosonConstructor.hh"
428 #include "G4LeptonConstructor.hh"
429 #include "G4MesonConstructor.hh"
430 #include "G4BosonConstructor.hh"
431 #include "G4BaryonConstructor.hh"
432 #include "G4IonConstructor.hh"
433 #include "G4ShortLivedConstructor.hh"
434 #include "G4Neutron.hh"
435 
436 // PROCESSES FOR GAMMA
437 #include "G4PhotoElectricEffect.hh"
438 #include "G4LivermorePhotoElectricModel.hh"
439 #include "G4ComptonScattering.hh"
440 #include "G4LivermoreComptonModel.hh"
441 #include "G4GammaConversion.hh"
442 #include "G4LivermoreGammaConversionModel.hh"
443 #include "G4RayleighScattering.hh"
444 #include "G4LivermoreRayleighModel.hh"
445 
446 // PROCESSES FOR e-
447 #include "G4eMultipleScattering.hh"
448 #include "G4UniversalFluctuation.hh"
449 #include "G4eIonisation.hh"
450 #include "G4LivermoreIonisationModel.hh"
451 #include "G4eBremsstrahlung.hh"
452 #include "G4LivermoreBremsstrahlungModel.hh"
453 #include "G4eIonisation.hh"
454 #include "G4eBremsstrahlung.hh"
455 #include "G4eplusAnnihilation.hh"
456 #include "G4KleinNishinaModel.hh"
457 #include "G4eMultipleScattering.hh"
458 #include "G4MuMultipleScattering.hh"
459 #include "G4MuIonisation.hh"
460 #include "G4MuBremsstrahlung.hh"
461 #include "G4MuPairProduction.hh"
462 #include "G4hMultipleScattering.hh"
463 #include "G4hIonisation.hh"
464 #include "G4hBremsstrahlung.hh"
465 #include "G4hPairProduction.hh"
466 #include "G4ionIonisation.hh"
467 #include "G4IonParametrisedLossModel.hh"
468 #include "G4NuclearStopping.hh"
469 #include "G4SystemOfUnits.hh"
470 #include "G4HadronElasticPhysics.hh"
471 #include "G4HadronicProcess.hh"
472 #include "G4HadronElastic.hh"
473 #include "G4HadronElasticProcess.hh"
474 
475 //c-s
476 #include "G4TripathiCrossSection.hh"
477 #include "G4IonsShenCrossSection.hh"
478 #include "G4ProtonInelasticCrossSection.hh"
479 #include "G4NeutronInelasticCrossSection.hh"
480 
481 // RadioactiveDecay
482 #include "G4RadioactiveDecay.hh"
483 #include "G4GenericIon.hh"
484 
485 #include "G4UrbanMscModel.hh"
486 #include "G4HadronElasticPhysicsHP.hh"
487 #include "QGSP_BERT_HP.hh"
488 
489 //HPNeutron
490 #include "G4NeutronHPThermalScatteringData.hh"
491 #include "G4NeutronHPThermalScattering.hh"
492 //
493 #include "G4HadronFissionProcess.hh"
494 #include "G4NeutronHPFissionData.hh"
495 #include "G4NeutronHPFission.hh"
496 //
497 #include "G4NeutronHPElastic.hh"
498 #include "G4NeutronHPElasticData.hh"
499 //
500 #include "G4HadronCaptureProcess.hh"
501 #include "G4NeutronHPCapture.hh"
502 #include "G4NeutronHPCaptureData.hh"
503 //
504 #include "G4NeutronInelasticProcess.hh"
505 #include "G4NeutronHPInelastic.hh"
506 #include "G4NeutronHPInelasticData.hh"
507 
508 //HPNeutron
509 #include "G4NeutronHPThermalScatteringData.hh"
510 #include "G4NeutronHPThermalScattering.hh"
511 //
512 #include "G4HadronFissionProcess.hh"
513 #include "G4NeutronHPFissionData.hh"
514 #include "G4NeutronHPFission.hh"
515 //
516 #include "G4NeutronHPElastic.hh"
517 #include "G4NeutronHPElasticData.hh"
518 //
519 #include "G4HadronCaptureProcess.hh"
520 #include "G4NeutronHPCapture.hh"
521 #include "G4NeutronHPCaptureData.hh"
522 //
523 #include "G4NeutronInelasticProcess.hh"
524 #include "G4NeutronHPInelastic.hh"
525 #include "G4NeutronHPInelasticData.hh"
526 
527 
528 
529 SToGS::HadronPhysicsList::HadronPhysicsList(const G4String& name): G4VPhysicsConstructor(name)
530 {
531 }
532 
534 {
535  G4LossTableManager::Instance();
536  // defaultCutValue = 1.0*CLHEP::mm;
537 
538  G4cout << G4endl << " ------ INFO ------ You are working with an hadronic PhysicsList (transportation, electromagnetic and hadronic processes)" << G4endl;
539 }
540 
541 
543 {
544 
545  G4LeptonConstructor lepton;
546  lepton.ConstructParticle();
547 
548  G4BosonConstructor boson;
549  boson.ConstructParticle();
550 
551  G4MesonConstructor meson;
552  meson.ConstructParticle();
553 
554  G4BaryonConstructor baryon;
555  baryon.ConstructParticle();
556 
557  G4ShortLivedConstructor shortLived;
558  shortLived.ConstructParticle();
559 
560  G4IonConstructor ion;
561  ion.ConstructParticle();
562 
563 
564  G4cout << " ------ INFO ------ PARTICLES DEFINITION FINISHED" << G4endl;
565 }
566 
568 {
569 
570 
571  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
572 
573  theParticleTable->GetIterator()->reset();
574 
575  //
576  while( (*theParticleTable->GetIterator())() ){
577  G4ParticleDefinition* particle = theParticleTable->GetIterator()->value();
578  G4ProcessManager* pmanager = particle->GetProcessManager();
579  G4String particleName = particle->GetParticleName();
580 
581  /*theParticleIterator->reset();
582  //
583  while( (*theParticleIterator)() ){
584  G4ParticleDefinition* particle = theParticleIterator->value();
585  G4ProcessManager* pmanager = particle->GetProcessManager();
586  G4String particleName = particle->GetParticleName();
587 
588 */
589 
590 
591 
592  // define EM physics for gamma
593  if (particleName == "gamma") {
594 
595  // activate photoelectric effect
596  ph->RegisterProcess(new G4PhotoElectricEffect, particle);
597 
598  // activate Compton scattering
599  G4ComptonScattering* cs = new G4ComptonScattering;
600  cs->SetEmModel(new G4KleinNishinaModel());
601  ph->RegisterProcess(cs, particle);
602 
603  // activate gamma conversion
604  ph->RegisterProcess(new G4GammaConversion, particle);
605 
606  }
607 
608 
609  // define EM physics for electron
610  else if (particleName == "e-") {
611 
612  // activate multiple scattering
613  G4eMultipleScattering* msc = new G4eMultipleScattering();
614  msc -> AddEmModel(0, new G4UrbanMscModel());
615  ph->RegisterProcess(msc, particle);
616 
617  // activate ionisation
618  G4eIonisation* eIoni = new G4eIonisation();
619  eIoni->SetStepFunction(0.1, 100*CLHEP::um);
620  ph->RegisterProcess(eIoni, particle);
621 
622  // activate bremsstrahlung
623  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
624 
625  }
626 
627  // define EM physics for positron
628  else if (particleName == "e+") {
629 
630  // activate multiple scattering
631  G4eMultipleScattering* msc = new G4eMultipleScattering();
632  msc -> AddEmModel(0, new G4UrbanMscModel());
633  ph->RegisterProcess(msc, particle);
634 
635  // activate ionisation
636  G4eIonisation* eIoni = new G4eIonisation();
637  eIoni->SetStepFunction(0.1, 100*CLHEP::um);
638  ph->RegisterProcess(eIoni, particle);
639 
640  // activate bremsstrahlung
641  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
642 
643  // activate annihilation
644  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
645 
646  }
647 
648 
649  // define EM physics for muons
650  else if (particleName == "mu+" || particleName == "mu-") {
651 
652  // activate multiple scattering
653  ph->RegisterProcess(new G4MuMultipleScattering(), particle);
654 
655  // activate ionisation
656  G4MuIonisation* muIoni = new G4MuIonisation();
657  muIoni->SetStepFunction(0.1, 50*CLHEP::um);
658  ph->RegisterProcess(muIoni, particle);
659 
660  // activate bremsstrahlung
661  ph->RegisterProcess(new G4MuBremsstrahlung(), particle);
662 
663  // activate pair production
664  ph->RegisterProcess(new G4MuPairProduction(), particle);
665 
666  }
667 
668 
669  // define EM physics for protons and pions
670  else if( particleName == "proton" ||
671  particleName == "pi-" ||
672  particleName == "pi+") {
673 
674  // activate multiple scattering
675  ph->RegisterProcess(new G4hMultipleScattering(), particle);
676 
677  // activate ionisation
678  G4hIonisation* hIoni = new G4hIonisation();
679  hIoni->SetStepFunction(0.1, 20*CLHEP::um);
680  ph->RegisterProcess(hIoni, particle);
681 
682  // activate bremsstrahlung
683  ph->RegisterProcess(new G4hBremsstrahlung(), particle);
684 
685  // activate pair production
686  ph->RegisterProcess(new G4hPairProduction(), particle);
687 
688  }
689 
690  // define EM physics for alpha and helium-3
691  else if( particleName == "alpha" ||
692  particleName == "He3") {
693 
694  // activate multiple scattering
695  ph->RegisterProcess(new G4hMultipleScattering(), particle);
696 
697  // activate ionisation
698  G4ionIonisation* ionIoni = new G4ionIonisation();
699  ionIoni->SetStepFunction(0.1, 1*CLHEP::um);
700  ph->RegisterProcess(ionIoni, particle);
701 
702  // activate nuclear stopping
703  ph->RegisterProcess(new G4NuclearStopping(), particle);
704 
705  }
706 
707  // define EM physics for all the other ions
708  else if( particleName == "GenericIon" ) {
709 
710  // activate multiple scattering
711  ph->RegisterProcess(new G4hMultipleScattering(), particle);
712 
713  // activate ionisation
714  G4ionIonisation* ionIoni = new G4ionIonisation();
715  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
716  ionIoni->SetStepFunction(0.1, 1*CLHEP::um);
717  ph->RegisterProcess(ionIoni, particle);
718 
719  // activate nuclear stopping
720  ph->RegisterProcess(new G4NuclearStopping(), particle);
721 
722  }
723 
724  // activate EM physics for all others charged particles except geantino
725  else if ((!particle->IsShortLived()) &&
726  (particle->GetPDGCharge() != 0.0) &&
727  (particle->GetParticleName() != "chargedgeantino")) {
728 
729  // activate multiple scattering
730  ph->RegisterProcess(new G4hMultipleScattering(), particle);
731 
732  // activate ionisation
733  ph->RegisterProcess(new G4hIonisation(), particle);
734  }
735 
736 
737 
738  }
739 
740  /* ***********************************************
741  * HADRONIC PHYSICS *
742  *********************************************** */
743 
744 
745  // define hadronic elastic processes
746  G4HadronElasticPhysicsHP* hadronicElastic = new G4HadronElasticPhysicsHP();
747  hadronicElastic->ConstructProcess();
748 
749  // define hadronic inelastic processes
750  QGSP_BERT_HP* hadronicPhysics = new QGSP_BERT_HP();
751  hadronicPhysics->ConstructProcess();
752 }
753 
755 {
756  /*
757  if ( verboseLevel > 0 ) {
758  G4cout << " HadronPhysicsList::SetCuts:";
759  G4cout << " CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
760  }
761 
762  // set cut values for gamma at first and for e- second and next for e+,
763  // because some processes for e+/e- need cut values for gamma
764  //
765  SetCutValue(defaultCutValue, "gamma");
766  SetCutValue(defaultCutValue, "e-");
767  SetCutValue(defaultCutValue, "e+");
768 
769  if (verboseLevel>0) DumpCutValuesTable();
770  */
771 }
772 
773 #endif
774 
775 
void ConstructProcess()
these methods Construct processes
void ConstructParticle()
these methods Construct particles
HadronPhysicsList(const G4String &name="SToGSHadron")