Geant4  10.01.p01
DMXPhysicsList.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 // --------------------------------------------------------------
28 // GEANT 4 - Underground Dark Matter Detector Advanced Example
29 //
30 // For information related to this code contact: Alex Howard
31 // e-mail: alexander.howard@cern.ch
32 // --------------------------------------------------------------
33 // Comments
34 //
35 // Underground Advanced
36 // by A. Howard and H. Araujo
37 // (27th November 2001)
38 //
39 // PhysicsList program
40 //
41 // Modified:
42 //
43 // 14-02-03 Fix bugs in msc and hIon instanciation + cut per region
44 //
45 // 05-02-05 AH - changes to G4Decay - added is not short lived protection
46 // and redefined particles to allow non-static creation
47 // i.e. changed construction to G4MesonConstructor, G4BaryonConstructor
48 //
49 // 23-10-09 LP - migrated EM physics from the LowEnergy processes (not supported) to
50 // the new G4Livermore model implementation. Results unchanged.
51 //
52 // --------------------------------------------------------------
53 
54 #include <iomanip>
55 
56 #include "DMXPhysicsList.hh"
57 
58 #include "globals.hh"
59 #include "G4SystemOfUnits.hh"
60 #include "G4ProcessManager.hh"
61 #include "G4ProcessVector.hh"
62 
63 #include "G4ParticleDefinition.hh"
64 #include "G4ParticleWithCuts.hh"
65 #include "G4ParticleTypes.hh"
66 #include "G4ParticleTable.hh"
67 
68 #include "G4ios.hh"
69 #include "G4UserLimits.hh"
70 
71 // Constructor /////////////////////////////////////////////////////////////
73 {
74 
75  defaultCutValue = 1.0*micrometer; //
79 
80  VerboseLevel = 1;
81  OpVerbLevel = 0;
82 
84 }
85 
86 
87 // Destructor //////////////////////////////////////////////////////////////
89 {;}
90 
91 
92 // Construct Particles /////////////////////////////////////////////////////
94 {
95 
96  // In this method, static member functions should be called
97  // for all particles which you want to use.
98  // This ensures that objects of these particle types will be
99  // created in the program.
100 
105 
106 }
107 
108 
109 // construct Bosons://///////////////////////////////////////////////////
111 {
112  // pseudo-particles
115 
116  // gamma
118 
119  //OpticalPhotons
121 
122 }
123 
124 
125 // construct Leptons://///////////////////////////////////////////////////
127 {
128  // leptons
133 
138 }
139 
140 
141 #include "G4MesonConstructor.hh"
142 #include "G4BaryonConstructor.hh"
143 #include "G4IonConstructor.hh"
144 
145 // construct Hadrons://///////////////////////////////////////////////////
147 {
148  // mesons
149  G4MesonConstructor mConstructor;
150  mConstructor.ConstructParticle();
151 
152  // baryons
153  G4BaryonConstructor bConstructor;
154  bConstructor.ConstructParticle();
155 
156  // ions
157  G4IonConstructor iConstructor;
158  iConstructor.ConstructParticle();
159 
160 }
161 
163 
164 // construct Shortliveds://///////////////////////////////////////////////////
166 {
167  G4ShortLivedConstructor slConstructor;
168  slConstructor.ConstructParticle();
169 }
170 
171 
172 
173 
174 // Construct Processes //////////////////////////////////////////////////////
176 {
177 
179 
180  ConstructEM();
181 
182  ConstructOp();
183 
184  ConstructHad();
185 
187 
188 }
189 
190 
191 // Transportation ///////////////////////////////////////////////////////////
192 #include "DMXMaxTimeCuts.hh"
193 #include "DMXMinEkineCuts.hh"
194 #include "G4StepLimiter.hh"
195 
197 
199 
200  theParticleIterator->reset();
201  while( (*theParticleIterator)() ){
202  G4ParticleDefinition* particle = theParticleIterator->value();
203  G4ProcessManager* pmanager = particle->GetProcessManager();
204  G4String particleName = particle->GetParticleName();
205  // time cuts for ONLY neutrons:
206  if(particleName == "neutron")
207  pmanager->AddDiscreteProcess(new DMXMaxTimeCuts());
208  // Energy cuts to kill charged (embedded in method) particles:
209  pmanager->AddDiscreteProcess(new DMXMinEkineCuts());
210 
211  // Step limit applied to all particles:
212  pmanager->AddProcess(new G4StepLimiter, -1,-1,1);
213 
214  }
215 }
216 
217 
218 // Electromagnetic Processes ////////////////////////////////////////////////
219 // all charged particles
220 
221 // gamma
222 #include "G4PhotoElectricEffect.hh"
224 
225 #include "G4ComptonScattering.hh"
227 
228 #include "G4GammaConversion.hh"
230 
231 #include "G4RayleighScattering.hh"
233 
234 
235 // e-
236 #include "G4eMultipleScattering.hh"
237 
238 #include "G4eIonisation.hh"
240 
241 #include "G4eBremsstrahlung.hh"
243 
244 
245 // e+
246 #include "G4eIonisation.hh"
247 #include "G4eBremsstrahlung.hh"
248 #include "G4eplusAnnihilation.hh"
249 
250 
251 // alpha and GenericIon and deuterons, triton, He3:
252 
253 //muon:
254 #include "G4MuIonisation.hh"
255 #include "G4MuBremsstrahlung.hh"
256 #include "G4MuPairProduction.hh"
257 #include "G4MuonMinusCapture.hh"
258 
259 //OTHERS:
260 #include "G4hIonisation.hh"
261 #include "G4hMultipleScattering.hh"
262 #include "G4hBremsstrahlung.hh"
263 #include "G4ionIonisation.hh"
265 
266 //em process options to allow msc step-limitation to be switched off
267 #include "G4EmProcessOptions.hh"
268 
270 
271  //set a finer grid of the physic tables in order to improve precision
272  //former LowEnergy models have 200 bins up to 100 GeV
273  G4EmProcessOptions opt;
274  opt.SetMaxEnergy(100*GeV);
275  opt.SetDEDXBinning(200);
276  opt.SetLambdaBinning(200);
277 
278  theParticleIterator->reset();
279  while( (*theParticleIterator)() ){
280  G4ParticleDefinition* particle = theParticleIterator->value();
281  G4ProcessManager* pmanager = particle->GetProcessManager();
282  G4String particleName = particle->GetParticleName();
283  G4String particleType = particle->GetParticleType();
284  G4double charge = particle->GetPDGCharge();
285 
286  if (particleName == "gamma")
287  {
288  //gamma
289  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
290  theRayleigh->SetEmModel(new G4LivermoreRayleighModel()); //not strictly necessary
291  pmanager->AddDiscreteProcess(theRayleigh);
292 
293  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
294  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
295  pmanager->AddDiscreteProcess(thePhotoElectricEffect);
296 
297  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
298  theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
299  pmanager->AddDiscreteProcess(theComptonScattering);
300 
301  G4GammaConversion* theGammaConversion = new G4GammaConversion();
302  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
303  pmanager->AddDiscreteProcess(theGammaConversion);
304 
305  }
306  else if (particleName == "e-")
307  {
308  //electron
309  // process ordering: AddProcess(name, at rest, along step, post step)
310  // Multiple scattering
313  pmanager->AddProcess(msc,-1, 1, 1);
314 
315  // Ionisation
316  G4eIonisation* eIonisation = new G4eIonisation();
317  eIonisation->SetEmModel(new G4LivermoreIonisationModel());
318  eIonisation->SetStepFunction(0.2, 100*um); //improved precision in tracking
319  pmanager->AddProcess(eIonisation,-1, 2, 2);
320 
321  // Bremsstrahlung
322  G4eBremsstrahlung* eBremsstrahlung = new G4eBremsstrahlung();
323  eBremsstrahlung->SetEmModel(new G4LivermoreBremsstrahlungModel());
324  pmanager->AddProcess(eBremsstrahlung, -1,-3, 3);
325  }
326  else if (particleName == "e+")
327  {
328  //positron
331  pmanager->AddProcess(msc,-1, 1, 1);
332 
333  // Ionisation
334  G4eIonisation* eIonisation = new G4eIonisation();
335  eIonisation->SetStepFunction(0.2, 100*um); //
336  pmanager->AddProcess(eIonisation, -1, 2, 2);
337 
338  //Bremsstrahlung (use default, no low-energy available)
339  pmanager->AddProcess(new G4eBremsstrahlung(), -1,-1, 3);
340 
341  //Annihilation
342  pmanager->AddProcess(new G4eplusAnnihilation(),0,-1, 4);
343  }
344  else if( particleName == "mu+" ||
345  particleName == "mu-" )
346  {
347  //muon
348  pmanager->AddProcess(new G4eMultipleScattering, -1, 1, 1);
349  pmanager->AddProcess(new G4MuIonisation(), -1, 2, 2);
350  pmanager->AddProcess(new G4MuBremsstrahlung(), -1,-1, 3);
351  pmanager->AddProcess(new G4MuPairProduction(), -1,-1, 4);
352  if( particleName == "mu-" )
353  pmanager->AddProcess(new G4MuonMinusCapture(), 0,-1,-1);
354  }
355  else if (particleName == "proton" ||
356  particleName == "pi+" ||
357  particleName == "pi-")
358  {
359  //multiple scattering
360  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
361 
362  //ionisation
363  G4hIonisation* hIonisation = new G4hIonisation();
364  hIonisation->SetStepFunction(0.2, 50*um);
365  pmanager->AddProcess(hIonisation, -1, 2, 2);
366 
367  //bremmstrahlung
368  pmanager->AddProcess(new G4hBremsstrahlung, -1,-3, 3);
369  }
370  else if(particleName == "alpha" ||
371  particleName == "deuteron" ||
372  particleName == "triton" ||
373  particleName == "He3")
374  {
375  //multiple scattering
376  pmanager->AddProcess(new G4hMultipleScattering,-1,1,1);
377 
378  //ionisation
379  G4ionIonisation* ionIoni = new G4ionIonisation();
380  ionIoni->SetStepFunction(0.1, 20*um);
381  pmanager->AddProcess(ionIoni, -1, 2, 2);
382  }
383  else if (particleName == "GenericIon")
384  {
385  // OBJECT may be dynamically created as either a GenericIon or nucleus
386  // G4Nucleus exists and therefore has particle type nucleus
387  // genericIon:
388 
389  //multiple scattering
390  pmanager->AddProcess(new G4hMultipleScattering,-1,1,1);
391 
392  //ionisation
393  G4ionIonisation* ionIoni = new G4ionIonisation();
394  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
395  ionIoni->SetStepFunction(0.1, 20*um);
396  pmanager->AddProcess(ionIoni, -1, 2, 2);
397  }
398 
399  else if ((!particle->IsShortLived()) &&
400  (charge != 0.0) &&
401  (particle->GetParticleName() != "chargedgeantino"))
402  {
403  //all others charged particles except geantino
404  G4hMultipleScattering* aMultipleScattering = new G4hMultipleScattering();
405  G4hIonisation* ahadronIon = new G4hIonisation();
406 
407  //multiple scattering
408  pmanager->AddProcess(aMultipleScattering,-1,1,1);
409 
410  //ionisation
411  pmanager->AddProcess(ahadronIon, -1,2,2);
412  }
413 
414  }
415 
416  // turn off msc step-limitation - especially as electron cut 1nm
418 
419  // switch on fluorescence, PIXE and Auger:
420  opt.SetFluo(true);
421  opt.SetPIXE(true);
422  opt.SetAuger(true);
423 
424 }
425 
426 
427 // Optical Processes ////////////////////////////////////////////////////////
428 #include "G4Scintillation.hh"
429 #include "G4OpAbsorption.hh"
430 //#include "G4OpRayleigh.hh"
431 #include "G4OpBoundaryProcess.hh"
432 
434 {
435  // default scintillation process
436  G4Scintillation* theScintProcessDef = new G4Scintillation("Scintillation");
437  // theScintProcessDef->DumpPhysicsTable();
438  theScintProcessDef->SetTrackSecondariesFirst(true);
439  theScintProcessDef->SetScintillationYieldFactor(1.0); //
440  theScintProcessDef->SetScintillationExcitationRatio(0.0); //
441  theScintProcessDef->SetVerboseLevel(OpVerbLevel);
442 
443  // scintillation process for alpha:
444  G4Scintillation* theScintProcessAlpha = new G4Scintillation("Scintillation");
445  // theScintProcessNuc->DumpPhysicsTable();
446  theScintProcessAlpha->SetTrackSecondariesFirst(true);
447  theScintProcessAlpha->SetScintillationYieldFactor(1.1);
448  theScintProcessAlpha->SetScintillationExcitationRatio(1.0);
449  theScintProcessAlpha->SetVerboseLevel(OpVerbLevel);
450 
451  // scintillation process for heavy nuclei
452  G4Scintillation* theScintProcessNuc = new G4Scintillation("Scintillation");
453  // theScintProcessNuc->DumpPhysicsTable();
454  theScintProcessNuc->SetTrackSecondariesFirst(true);
455  theScintProcessNuc->SetScintillationYieldFactor(0.2);
456  theScintProcessNuc->SetScintillationExcitationRatio(1.0);
457  theScintProcessNuc->SetVerboseLevel(OpVerbLevel);
458 
459  // optical processes
460  G4OpAbsorption* theAbsorptionProcess = new G4OpAbsorption();
461  // G4OpRayleigh* theRayleighScatteringProcess = new G4OpRayleigh();
462  G4OpBoundaryProcess* theBoundaryProcess = new G4OpBoundaryProcess();
463  // theAbsorptionProcess->DumpPhysicsTable();
464  // theRayleighScatteringProcess->DumpPhysicsTable();
465  theAbsorptionProcess->SetVerboseLevel(OpVerbLevel);
466  // theRayleighScatteringProcess->SetVerboseLevel(OpVerbLevel);
467  theBoundaryProcess->SetVerboseLevel(OpVerbLevel);
468 
469  theParticleIterator->reset();
470  while( (*theParticleIterator)() )
471  {
472  G4ParticleDefinition* particle = theParticleIterator->value();
473  G4ProcessManager* pmanager = particle->GetProcessManager();
474  G4String particleName = particle->GetParticleName();
475  if (theScintProcessDef->IsApplicable(*particle)) {
476  // if(particle->GetPDGMass() > 5.0*GeV)
477  if(particle->GetParticleName() == "GenericIon") {
478  pmanager->AddProcess(theScintProcessNuc); // AtRestDiscrete
479  pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxAtRest);
480  pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxPostStep);
481  }
482  else if(particle->GetParticleName() == "alpha") {
483  pmanager->AddProcess(theScintProcessAlpha);
484  pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxAtRest);
485  pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxPostStep);
486  }
487  else {
488  pmanager->AddProcess(theScintProcessDef);
489  pmanager->SetProcessOrderingToLast(theScintProcessDef,idxAtRest);
490  pmanager->SetProcessOrderingToLast(theScintProcessDef,idxPostStep);
491  }
492  }
493 
494  if (particleName == "opticalphoton") {
495  pmanager->AddDiscreteProcess(theAbsorptionProcess);
496  // pmanager->AddDiscreteProcess(theRayleighScatteringProcess);
497  pmanager->AddDiscreteProcess(theBoundaryProcess);
498  }
499  }
500 }
501 
502 
503 // Hadronic processes ////////////////////////////////////////////////////////
504 
505 // Elastic processes:
506 #include "G4HadronElasticProcess.hh"
507 #include "G4ChipsElasticModel.hh"
508 #include "G4ElasticHadrNucleusHE.hh"
509 
510 // Inelastic processes:
524 
525 // High energy FTFP model and Bertini cascade
526 #include "G4FTFModel.hh"
528 #include "G4ExcitedStringDecay.hh"
529 #include "G4PreCompoundModel.hh"
531 #include "G4TheoFSGenerator.hh"
532 #include "G4CascadeInterface.hh"
533 
534 // Cross sections
535 #include "G4VCrossSectionDataSet.hh"
537 
538 #include "G4CrossSectionElastic.hh"
539 #include "G4BGGPionElasticXS.hh"
540 #include "G4AntiNuclElastic.hh"
541 
544 #include "G4CrossSectionPairGG.hh"
548 
549 #include "G4HadronElastic.hh"
550 #include "G4HadronCaptureProcess.hh"
551 
552 // Neutron high-precision models: <20 MeV
553 #include "G4NeutronHPElastic.hh"
554 #include "G4NeutronHPElasticData.hh"
555 #include "G4NeutronHPCapture.hh"
556 #include "G4NeutronHPCaptureData.hh"
557 #include "G4NeutronHPInelastic.hh"
559 
560 // Stopping processes
564 
565 
566 
568 {
569  //Elastic models
570  const G4double elastic_elimitPi = 1.0*GeV;
571 
572  G4HadronElastic* elastic_lhep0 = new G4HadronElastic();
573  G4HadronElastic* elastic_lhep1 = new G4HadronElastic();
574  elastic_lhep1->SetMaxEnergy( elastic_elimitPi );
575  G4ChipsElasticModel* elastic_chip = new G4ChipsElasticModel();
576  G4ElasticHadrNucleusHE* elastic_he = new G4ElasticHadrNucleusHE();
577  elastic_he->SetMinEnergy( elastic_elimitPi );
578 
579 
580  // Inelastic scattering
581  const G4double theFTFMin0 = 0.0*GeV;
582  const G4double theFTFMin1 = 4.0*GeV;
583  const G4double theFTFMax = 100.0*TeV;
584  const G4double theBERTMin0 = 0.0*GeV;
585  const G4double theBERTMin1 = 19.0*MeV;
586  const G4double theBERTMax = 5.0*GeV;
587  const G4double theHPMin = 0.0*GeV;
588  const G4double theHPMax = 20.0*MeV;
589 
590  G4FTFModel * theStringModel = new G4FTFModel;
592  theStringModel->SetFragmentationModel( theStringDecay );
593  G4PreCompoundModel * thePreEquilib = new G4PreCompoundModel( new G4ExcitationHandler );
594  G4GeneratorPrecompoundInterface * theCascade = new G4GeneratorPrecompoundInterface( thePreEquilib );
595 
596  G4TheoFSGenerator * theFTFModel0 = new G4TheoFSGenerator( "FTFP" );
597  theFTFModel0->SetHighEnergyGenerator( theStringModel );
598  theFTFModel0->SetTransport( theCascade );
599  theFTFModel0->SetMinEnergy( theFTFMin0 );
600  theFTFModel0->SetMaxEnergy( theFTFMax );
601 
602  G4TheoFSGenerator * theFTFModel1 = new G4TheoFSGenerator( "FTFP" );
603  theFTFModel1->SetHighEnergyGenerator( theStringModel );
604  theFTFModel1->SetTransport( theCascade );
605  theFTFModel1->SetMinEnergy( theFTFMin1 );
606  theFTFModel1->SetMaxEnergy( theFTFMax );
607 
608  G4CascadeInterface * theBERTModel0 = new G4CascadeInterface;
609  theBERTModel0->SetMinEnergy( theBERTMin0 );
610  theBERTModel0->SetMaxEnergy( theBERTMax );
611 
612  G4CascadeInterface * theBERTModel1 = new G4CascadeInterface;
613  theBERTModel1->SetMinEnergy( theBERTMin1 );
614  theBERTModel1->SetMaxEnergy( theBERTMax );
615 
619  GetCrossSectionDataSet(G4GGNuclNuclCrossSection::Default_Name());
620 
621 
622  theParticleIterator->reset();
623  while ((*theParticleIterator)())
624  {
625  G4ParticleDefinition* particle = theParticleIterator->value();
626  G4ProcessManager* pmanager = particle->GetProcessManager();
627  G4String particleName = particle->GetParticleName();
628 
629  if (particleName == "pi+")
630  {
631  // Elastic scattering
632  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
633  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
634  theElasticProcess->RegisterMe( elastic_lhep1 );
635  theElasticProcess->RegisterMe( elastic_he );
636  pmanager->AddDiscreteProcess( theElasticProcess );
637  //Inelastic scattering
638  G4PionPlusInelasticProcess* theInelasticProcess =
639  new G4PionPlusInelasticProcess("inelastic");
640  theInelasticProcess->AddDataSet( thePiData );
641  theInelasticProcess->RegisterMe( theFTFModel1 );
642  theInelasticProcess->RegisterMe( theBERTModel0 );
643  pmanager->AddDiscreteProcess( theInelasticProcess );
644  }
645 
646  else if (particleName == "pi-")
647  {
648  // Elastic scattering
649  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
650  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
651  theElasticProcess->RegisterMe( elastic_lhep1 );
652  theElasticProcess->RegisterMe( elastic_he );
653  pmanager->AddDiscreteProcess( theElasticProcess );
654  //Inelastic scattering
655  G4PionMinusInelasticProcess* theInelasticProcess =
656  new G4PionMinusInelasticProcess("inelastic");
657  theInelasticProcess->AddDataSet( thePiData );
658  theInelasticProcess->RegisterMe( theFTFModel1 );
659  theInelasticProcess->RegisterMe( theBERTModel0 );
660  pmanager->AddDiscreteProcess( theInelasticProcess );
661  //Absorption
663  }
664 
665  else if (particleName == "kaon+")
666  {
667  // Elastic scattering
668  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
669  theElasticProcess->RegisterMe( elastic_lhep0 );
670  pmanager->AddDiscreteProcess( theElasticProcess );
671  // Inelastic scattering
672  G4KaonPlusInelasticProcess* theInelasticProcess =
673  new G4KaonPlusInelasticProcess("inelastic");
674  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
675  GetCrossSectionDataSet(G4ChipsKaonPlusInelasticXS::Default_Name()));
676  theInelasticProcess->RegisterMe( theFTFModel1 );
677  theInelasticProcess->RegisterMe( theBERTModel0 );
678  pmanager->AddDiscreteProcess( theInelasticProcess );
679  }
680 
681  else if (particleName == "kaon0S")
682  {
683  // Elastic scattering
684  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
685  theElasticProcess->RegisterMe( elastic_lhep0 );
686  pmanager->AddDiscreteProcess( theElasticProcess );
687  // Inelastic scattering
688  G4KaonZeroSInelasticProcess* theInelasticProcess =
689  new G4KaonZeroSInelasticProcess("inelastic");
690  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
691  GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
692  theInelasticProcess->RegisterMe( theFTFModel1 );
693  theInelasticProcess->RegisterMe( theBERTModel0 );
694  pmanager->AddDiscreteProcess( theInelasticProcess );
695  }
696 
697  else if (particleName == "kaon0L")
698  {
699  // Elastic scattering
700  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
701  theElasticProcess->RegisterMe( elastic_lhep0 );
702  pmanager->AddDiscreteProcess( theElasticProcess );
703  // Inelastic scattering
704  G4KaonZeroLInelasticProcess* theInelasticProcess =
705  new G4KaonZeroLInelasticProcess("inelastic");
706  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
707  GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
708  theInelasticProcess->RegisterMe( theFTFModel1 );
709  theInelasticProcess->RegisterMe( theBERTModel0 );
710  pmanager->AddDiscreteProcess( theInelasticProcess );
711  }
712 
713  else if (particleName == "kaon-")
714  {
715  // Elastic scattering
716  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
717  theElasticProcess->RegisterMe( elastic_lhep0 );
718  pmanager->AddDiscreteProcess( theElasticProcess );
719  // Inelastic scattering
720  G4KaonMinusInelasticProcess* theInelasticProcess =
721  new G4KaonMinusInelasticProcess("inelastic");
722  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
723  GetCrossSectionDataSet(G4ChipsKaonMinusInelasticXS::Default_Name()));
724  theInelasticProcess->RegisterMe( theFTFModel1 );
725  theInelasticProcess->RegisterMe( theBERTModel0 );
726  pmanager->AddDiscreteProcess( theInelasticProcess );
728  }
729 
730  else if (particleName == "proton")
731  {
732  // Elastic scattering
733  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
735  GetCrossSectionDataSet(G4ChipsProtonElasticXS::Default_Name()));
736  theElasticProcess->RegisterMe( elastic_chip );
737  pmanager->AddDiscreteProcess( theElasticProcess );
738  // Inelastic scattering
739  G4ProtonInelasticProcess* theInelasticProcess =
740  new G4ProtonInelasticProcess("inelastic");
741  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Proton::Proton() ) );
742  theInelasticProcess->RegisterMe( theFTFModel1 );
743  theInelasticProcess->RegisterMe( theBERTModel0 );
744  pmanager->AddDiscreteProcess( theInelasticProcess );
745  }
746  else if (particleName == "anti_proton")
747  {
748  // Elastic scattering
749  const G4double elastic_elimitAntiNuc = 100.0*MeV;
750  G4AntiNuclElastic* elastic_anuc = new G4AntiNuclElastic();
751  elastic_anuc->SetMinEnergy( elastic_elimitAntiNuc );
752  G4CrossSectionElastic* elastic_anucxs = new G4CrossSectionElastic( elastic_anuc->GetComponentCrossSection() );
753  G4HadronElastic* elastic_lhep2 = new G4HadronElastic();
754  elastic_lhep2->SetMaxEnergy( elastic_elimitAntiNuc );
755  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
756  theElasticProcess->AddDataSet( elastic_anucxs );
757  theElasticProcess->RegisterMe( elastic_lhep2 );
758  theElasticProcess->RegisterMe( elastic_anuc );
759  pmanager->AddDiscreteProcess( theElasticProcess );
760  // Inelastic scattering
761  G4AntiProtonInelasticProcess* theInelasticProcess =
762  new G4AntiProtonInelasticProcess("inelastic");
763  theInelasticProcess->AddDataSet( theAntiNucleonData );
764  theInelasticProcess->RegisterMe( theFTFModel0 );
765  pmanager->AddDiscreteProcess( theInelasticProcess );
766  // Absorption
768  }
769 
770  else if (particleName == "neutron") {
771  // elastic scattering
772  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
773  theElasticProcess->AddDataSet(G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsNeutronElasticXS::Default_Name()));
774  G4HadronElastic* elastic_neutronChipsModel = new G4ChipsElasticModel();
775  elastic_neutronChipsModel->SetMinEnergy( 19.0*MeV );
776  theElasticProcess->RegisterMe( elastic_neutronChipsModel );
777  G4NeutronHPElastic * theElasticNeutronHP = new G4NeutronHPElastic;
778  theElasticNeutronHP->SetMinEnergy( theHPMin );
779  theElasticNeutronHP->SetMaxEnergy( theHPMax );
780  theElasticProcess->RegisterMe( theElasticNeutronHP );
781  theElasticProcess->AddDataSet( new G4NeutronHPElasticData );
782  pmanager->AddDiscreteProcess( theElasticProcess );
783  // inelastic scattering
784  G4NeutronInelasticProcess* theInelasticProcess =
785  new G4NeutronInelasticProcess("inelastic");
786  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Neutron::Neutron() ) );
787  theInelasticProcess->RegisterMe( theFTFModel1 );
788  theInelasticProcess->RegisterMe( theBERTModel1 );
789  G4NeutronHPInelastic * theNeutronInelasticHPModel = new G4NeutronHPInelastic;
790  theNeutronInelasticHPModel->SetMinEnergy( theHPMin );
791  theNeutronInelasticHPModel->SetMaxEnergy( theHPMax );
792  theInelasticProcess->RegisterMe( theNeutronInelasticHPModel );
793  theInelasticProcess->AddDataSet( new G4NeutronHPInelasticData );
794  pmanager->AddDiscreteProcess(theInelasticProcess);
795  // capture
796  G4HadronCaptureProcess* theCaptureProcess =
798  G4NeutronHPCapture * theLENeutronCaptureModel = new G4NeutronHPCapture;
799  theLENeutronCaptureModel->SetMinEnergy(theHPMin);
800  theLENeutronCaptureModel->SetMaxEnergy(theHPMax);
801  theCaptureProcess->RegisterMe(theLENeutronCaptureModel);
802  theCaptureProcess->AddDataSet( new G4NeutronHPCaptureData);
803  pmanager->AddDiscreteProcess(theCaptureProcess);
804 
805  }
806  else if (particleName == "anti_neutron")
807  {
808  // Elastic scattering
809  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
810  theElasticProcess->RegisterMe( elastic_lhep0 );
811  pmanager->AddDiscreteProcess( theElasticProcess );
812  // Inelastic scattering (include annihilation on-fly)
813  G4AntiNeutronInelasticProcess* theInelasticProcess =
814  new G4AntiNeutronInelasticProcess("inelastic");
815  theInelasticProcess->AddDataSet( theAntiNucleonData );
816  theInelasticProcess->RegisterMe( theFTFModel0 );
817  pmanager->AddDiscreteProcess( theInelasticProcess );
818  }
819 
820  else if (particleName == "deuteron")
821  {
822  // Elastic scattering
823  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
824  theElasticProcess->RegisterMe( elastic_lhep0 );
825  pmanager->AddDiscreteProcess( theElasticProcess );
826  // Inelastic scattering
827  G4DeuteronInelasticProcess* theInelasticProcess =
828  new G4DeuteronInelasticProcess("inelastic");
829  theInelasticProcess->AddDataSet( theGGNuclNuclData );
830  theInelasticProcess->RegisterMe( theFTFModel1 );
831  theInelasticProcess->RegisterMe( theBERTModel0 );
832  pmanager->AddDiscreteProcess( theInelasticProcess );
833  }
834 
835  else if (particleName == "triton")
836  {
837  // Elastic scattering
838  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
839  theElasticProcess->RegisterMe( elastic_lhep0 );
840  pmanager->AddDiscreteProcess( theElasticProcess );
841  // Inelastic scattering
842  G4TritonInelasticProcess* theInelasticProcess =
843  new G4TritonInelasticProcess("inelastic");
844  theInelasticProcess->AddDataSet( theGGNuclNuclData );
845  theInelasticProcess->RegisterMe( theFTFModel1 );
846  theInelasticProcess->RegisterMe( theBERTModel0 );
847  pmanager->AddDiscreteProcess( theInelasticProcess );
848  }
849  else if (particleName == "alpha")
850  {
851  // Elastic scattering
852  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
853  theElasticProcess->RegisterMe( elastic_lhep0 );
854  pmanager->AddDiscreteProcess( theElasticProcess );
855  // Inelastic scattering
856  G4AlphaInelasticProcess* theInelasticProcess =
857  new G4AlphaInelasticProcess("inelastic");
858  theInelasticProcess->AddDataSet( theGGNuclNuclData );
859  theInelasticProcess->RegisterMe( theFTFModel1 );
860  theInelasticProcess->RegisterMe( theBERTModel0 );
861  pmanager->AddDiscreteProcess( theInelasticProcess );
862  }
863 
864  }
865 }
866 
867 
868 // Decays ///////////////////////////////////////////////////////////////////
869 #include "G4Decay.hh"
870 #include "G4RadioactiveDecay.hh"
871 #include "G4IonTable.hh"
872 #include "G4Ions.hh"
873 
875 
876  // Add Decay Process
877  G4Decay* theDecayProcess = new G4Decay();
878  theParticleIterator->reset();
879  while( (*theParticleIterator)() )
880  {
881  G4ParticleDefinition* particle = theParticleIterator->value();
882  G4ProcessManager* pmanager = particle->GetProcessManager();
883 
884  if (theDecayProcess->IsApplicable(*particle) && !particle->IsShortLived())
885  {
886  pmanager ->AddProcess(theDecayProcess);
887  // set ordering for PostStepDoIt and AtRestDoIt
888  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
889  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
890  }
891  }
892 
893  // Declare radioactive decay to the GenericIon in the IonTable.
894  const G4IonTable *theIonTable =
896  G4RadioactiveDecay *theRadioactiveDecay = new G4RadioactiveDecay();
897 
898  for (G4int i=0; i<theIonTable->Entries(); i++)
899  {
900  G4String particleName = theIonTable->GetParticle(i)->GetParticleName();
901  G4String particleType = theIonTable->GetParticle(i)->GetParticleType();
902 
903  if (particleName == "GenericIon")
904  {
905  G4ProcessManager* pmanager =
906  theIonTable->GetParticle(i)->GetProcessManager();
907  pmanager->SetVerboseLevel(VerboseLevel);
908  pmanager ->AddProcess(theRadioactiveDecay);
909  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxPostStep);
910  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxAtRest);
911  }
912  }
913 }
914 
915 // Cuts /////////////////////////////////////////////////////////////////////
917 {
918 
919  if (verboseLevel >1)
920  G4cout << "DMXPhysicsList::SetCuts:";
921 
922  if (verboseLevel>0){
923  G4cout << "DMXPhysicsList::SetCuts:";
924  G4cout << "CutLength : "
925  << G4BestUnit(defaultCutValue,"Length") << G4endl;
926  }
927 
928  //special for low energy physics
929  G4double lowlimit=250*eV;
931 
932  // set cut values for gamma at first and for e- second and next for e+,
933  // because some processes for e+/e- need cut values for gamma
934  SetCutValue(cutForGamma, "gamma");
937 
939 }
940 
static G4Geantino * GeantinoDefinition()
Definition: G4Geantino.cc:82
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static const double MeV
Definition: G4SIunits.hh:193
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
void ConstructMyLeptons()
void SetPIXE(G4bool val)
G4int Entries() const
Definition: G4IonTable.cc:1602
void SetCutValue(G4double aCut, const G4String &pname)
void SetEnergyRange(G4double lowedge, G4double highedge)
virtual void ConstructGeneral()
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
void SetFragmentationModel(G4VStringFragmentation *aModel)
static const double nanometer
Definition: G4SIunits.hh:91
void SetStepFunction(G4double v1, G4double v2)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void ConstructMyHadrons()
void SetVerboseLevel(G4int value)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static void ConstructParticle()
static void ConstructParticle()
static const char * Default_Name()
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
G4ComponentAntiNuclNuclearXS * GetComponentCrossSection()
const G4String & GetParticleName() const
void SetFluo(G4bool val)
virtual void ConstructProcess()
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
static void ConstructParticle()
void SetAuger(G4bool val)
void RegisterMe(G4HadronicInteraction *a)
void SetDEDXBinning(G4int val)
virtual void AddTransportation()
void SetScintillationYieldFactor(const G4double yieldfactor)
void SetMinEnergy(G4double anEnergy)
void SetEmModel(G4VEmModel *, G4int index=1)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
G4IonTable * GetIonTable() const
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
void SetLambdaBinning(G4int val)
void SetVerboseLevel(G4int value)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
static G4NeutrinoE * NeutrinoEDefinition()
Definition: G4NeutrinoE.cc:80
G4double cutForGamma
static const char * Default_Name()
static G4AntiNeutrinoMu * AntiNeutrinoMuDefinition()
static G4CrossSectionDataSetRegistry * Instance()
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double cutForElectron
static const double GeV
Definition: G4SIunits.hh:196
const G4String & GetParticleType() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
void SetMaxEnergy(G4double val)
virtual void ConstructHad()
void SetScintillationExcitationRatio(const G4double ratio)
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
static const double micrometer
Definition: G4SIunits.hh:90
void ConstructMyBosons()
static G4ProductionCutsTable * GetProductionCutsTable()
static const double eV
Definition: G4SIunits.hh:194
virtual void SetCuts()
static G4ParticleTable * GetParticleTable()
void SetTrackSecondariesFirst(const G4bool state)
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
void ConstructMyShortLiveds()
void SetEmModel(G4VEmModel *, G4int index=1)
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
virtual void ConstructOp()
void SetMaxEnergy(const G4double anEnergy)
static G4ChargedGeantino * ChargedGeantinoDefinition()
static const char * Default_Name()
#define G4endl
Definition: G4ios.hh:61
void SetTransport(G4VIntraNuclearTransportModel *const value)
static const double TeV
Definition: G4SIunits.hh:197
void SetMscStepLimitation(G4MscStepLimitType val)
G4int AddRestProcess(G4VProcess *aProcess, G4int ord=ordDefault)
static G4AntiNeutrinoE * AntiNeutrinoEDefinition()
static G4OpticalPhoton * OpticalPhotonDefinition()
double G4double
Definition: G4Types.hh:76
G4double cutForPositron
G4double GetPDGCharge() const
static G4NeutrinoMu * NeutrinoMuDefinition()
Definition: G4NeutrinoMu.cc:80
#define theParticleIterator
virtual void ConstructParticle()
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
virtual void ConstructEM()
void SetStepLimitType(G4MscStepLimitType val)
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
G4ParticleDefinition * GetParticle(G4int index) const
Definition: G4IonTable.cc:1556