Geant4  10.02.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 "G4ParticleHPElastic.hh"
555 #include "G4ParticleHPCapture.hh"
557 #include "G4ParticleHPInelastic.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 
618  G4ComponentGGNuclNuclXsc * ggNuclNuclXsec = new G4ComponentGGNuclNuclXsc();
619  G4VCrossSectionDataSet * theGGNuclNuclData = new G4CrossSectionInelastic(ggNuclNuclXsec);
620 
621  theParticleIterator->reset();
622  while ((*theParticleIterator)())
623  {
624  G4ParticleDefinition* particle = theParticleIterator->value();
625  G4ProcessManager* pmanager = particle->GetProcessManager();
626  G4String particleName = particle->GetParticleName();
627 
628  if (particleName == "pi+")
629  {
630  // Elastic scattering
631  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
632  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
633  theElasticProcess->RegisterMe( elastic_lhep1 );
634  theElasticProcess->RegisterMe( elastic_he );
635  pmanager->AddDiscreteProcess( theElasticProcess );
636  //Inelastic scattering
637  G4PionPlusInelasticProcess* theInelasticProcess =
638  new G4PionPlusInelasticProcess("inelastic");
639  theInelasticProcess->AddDataSet( thePiData );
640  theInelasticProcess->RegisterMe( theFTFModel1 );
641  theInelasticProcess->RegisterMe( theBERTModel0 );
642  pmanager->AddDiscreteProcess( theInelasticProcess );
643  }
644 
645  else if (particleName == "pi-")
646  {
647  // Elastic scattering
648  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
649  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
650  theElasticProcess->RegisterMe( elastic_lhep1 );
651  theElasticProcess->RegisterMe( elastic_he );
652  pmanager->AddDiscreteProcess( theElasticProcess );
653  //Inelastic scattering
654  G4PionMinusInelasticProcess* theInelasticProcess =
655  new G4PionMinusInelasticProcess("inelastic");
656  theInelasticProcess->AddDataSet( thePiData );
657  theInelasticProcess->RegisterMe( theFTFModel1 );
658  theInelasticProcess->RegisterMe( theBERTModel0 );
659  pmanager->AddDiscreteProcess( theInelasticProcess );
660  //Absorption
662  }
663 
664  else if (particleName == "kaon+")
665  {
666  // Elastic scattering
667  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
668  theElasticProcess->RegisterMe( elastic_lhep0 );
669  pmanager->AddDiscreteProcess( theElasticProcess );
670  // Inelastic scattering
671  G4KaonPlusInelasticProcess* theInelasticProcess =
672  new G4KaonPlusInelasticProcess("inelastic");
673  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
674  GetCrossSectionDataSet(G4ChipsKaonPlusInelasticXS::Default_Name()));
675  theInelasticProcess->RegisterMe( theFTFModel1 );
676  theInelasticProcess->RegisterMe( theBERTModel0 );
677  pmanager->AddDiscreteProcess( theInelasticProcess );
678  }
679 
680  else if (particleName == "kaon0S")
681  {
682  // Elastic scattering
683  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
684  theElasticProcess->RegisterMe( elastic_lhep0 );
685  pmanager->AddDiscreteProcess( theElasticProcess );
686  // Inelastic scattering
687  G4KaonZeroSInelasticProcess* theInelasticProcess =
688  new G4KaonZeroSInelasticProcess("inelastic");
689  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
690  GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
691  theInelasticProcess->RegisterMe( theFTFModel1 );
692  theInelasticProcess->RegisterMe( theBERTModel0 );
693  pmanager->AddDiscreteProcess( theInelasticProcess );
694  }
695 
696  else if (particleName == "kaon0L")
697  {
698  // Elastic scattering
699  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
700  theElasticProcess->RegisterMe( elastic_lhep0 );
701  pmanager->AddDiscreteProcess( theElasticProcess );
702  // Inelastic scattering
703  G4KaonZeroLInelasticProcess* theInelasticProcess =
704  new G4KaonZeroLInelasticProcess("inelastic");
705  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
706  GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
707  theInelasticProcess->RegisterMe( theFTFModel1 );
708  theInelasticProcess->RegisterMe( theBERTModel0 );
709  pmanager->AddDiscreteProcess( theInelasticProcess );
710  }
711 
712  else if (particleName == "kaon-")
713  {
714  // Elastic scattering
715  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
716  theElasticProcess->RegisterMe( elastic_lhep0 );
717  pmanager->AddDiscreteProcess( theElasticProcess );
718  // Inelastic scattering
719  G4KaonMinusInelasticProcess* theInelasticProcess =
720  new G4KaonMinusInelasticProcess("inelastic");
721  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
722  GetCrossSectionDataSet(G4ChipsKaonMinusInelasticXS::Default_Name()));
723  theInelasticProcess->RegisterMe( theFTFModel1 );
724  theInelasticProcess->RegisterMe( theBERTModel0 );
725  pmanager->AddDiscreteProcess( theInelasticProcess );
727  }
728 
729  else if (particleName == "proton")
730  {
731  // Elastic scattering
732  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
734  GetCrossSectionDataSet(G4ChipsProtonElasticXS::Default_Name()));
735  theElasticProcess->RegisterMe( elastic_chip );
736  pmanager->AddDiscreteProcess( theElasticProcess );
737  // Inelastic scattering
738  G4ProtonInelasticProcess* theInelasticProcess =
739  new G4ProtonInelasticProcess("inelastic");
740  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Proton::Proton() ) );
741  theInelasticProcess->RegisterMe( theFTFModel1 );
742  theInelasticProcess->RegisterMe( theBERTModel0 );
743  pmanager->AddDiscreteProcess( theInelasticProcess );
744  }
745  else if (particleName == "anti_proton")
746  {
747  // Elastic scattering
748  const G4double elastic_elimitAntiNuc = 100.0*MeV;
749  G4AntiNuclElastic* elastic_anuc = new G4AntiNuclElastic();
750  elastic_anuc->SetMinEnergy( elastic_elimitAntiNuc );
751  G4CrossSectionElastic* elastic_anucxs = new G4CrossSectionElastic( elastic_anuc->GetComponentCrossSection() );
752  G4HadronElastic* elastic_lhep2 = new G4HadronElastic();
753  elastic_lhep2->SetMaxEnergy( elastic_elimitAntiNuc );
754  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
755  theElasticProcess->AddDataSet( elastic_anucxs );
756  theElasticProcess->RegisterMe( elastic_lhep2 );
757  theElasticProcess->RegisterMe( elastic_anuc );
758  pmanager->AddDiscreteProcess( theElasticProcess );
759  // Inelastic scattering
760  G4AntiProtonInelasticProcess* theInelasticProcess =
761  new G4AntiProtonInelasticProcess("inelastic");
762  theInelasticProcess->AddDataSet( theAntiNucleonData );
763  theInelasticProcess->RegisterMe( theFTFModel0 );
764  pmanager->AddDiscreteProcess( theInelasticProcess );
765  // Absorption
767  }
768 
769  else if (particleName == "neutron") {
770  // elastic scattering
771  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
772  theElasticProcess->AddDataSet(G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsNeutronElasticXS::Default_Name()));
773  G4HadronElastic* elastic_neutronChipsModel = new G4ChipsElasticModel();
774  elastic_neutronChipsModel->SetMinEnergy( 19.0*MeV );
775  theElasticProcess->RegisterMe( elastic_neutronChipsModel );
776  G4ParticleHPElastic * theElasticNeutronHP = new G4ParticleHPElastic;
777  theElasticNeutronHP->SetMinEnergy( theHPMin );
778  theElasticNeutronHP->SetMaxEnergy( theHPMax );
779  theElasticProcess->RegisterMe( theElasticNeutronHP );
780  theElasticProcess->AddDataSet( new G4ParticleHPElasticData );
781  pmanager->AddDiscreteProcess( theElasticProcess );
782  // inelastic scattering
783  G4NeutronInelasticProcess* theInelasticProcess =
784  new G4NeutronInelasticProcess("inelastic");
785  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Neutron::Neutron() ) );
786  theInelasticProcess->RegisterMe( theFTFModel1 );
787  theInelasticProcess->RegisterMe( theBERTModel1 );
788  G4ParticleHPInelastic * theNeutronInelasticHPModel = new G4ParticleHPInelastic;
789  theNeutronInelasticHPModel->SetMinEnergy( theHPMin );
790  theNeutronInelasticHPModel->SetMaxEnergy( theHPMax );
791  theInelasticProcess->RegisterMe( theNeutronInelasticHPModel );
792  theInelasticProcess->AddDataSet( new G4ParticleHPInelasticData );
793  pmanager->AddDiscreteProcess(theInelasticProcess);
794  // capture
795  G4HadronCaptureProcess* theCaptureProcess =
797  G4ParticleHPCapture * theLENeutronCaptureModel = new G4ParticleHPCapture;
798  theLENeutronCaptureModel->SetMinEnergy(theHPMin);
799  theLENeutronCaptureModel->SetMaxEnergy(theHPMax);
800  theCaptureProcess->RegisterMe(theLENeutronCaptureModel);
801  theCaptureProcess->AddDataSet( new G4ParticleHPCaptureData);
802  pmanager->AddDiscreteProcess(theCaptureProcess);
803 
804  }
805  else if (particleName == "anti_neutron")
806  {
807  // Elastic scattering
808  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
809  theElasticProcess->RegisterMe( elastic_lhep0 );
810  pmanager->AddDiscreteProcess( theElasticProcess );
811  // Inelastic scattering (include annihilation on-fly)
812  G4AntiNeutronInelasticProcess* theInelasticProcess =
813  new G4AntiNeutronInelasticProcess("inelastic");
814  theInelasticProcess->AddDataSet( theAntiNucleonData );
815  theInelasticProcess->RegisterMe( theFTFModel0 );
816  pmanager->AddDiscreteProcess( theInelasticProcess );
817  }
818 
819  else if (particleName == "deuteron")
820  {
821  // Elastic scattering
822  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
823  theElasticProcess->RegisterMe( elastic_lhep0 );
824  pmanager->AddDiscreteProcess( theElasticProcess );
825  // Inelastic scattering
826  G4DeuteronInelasticProcess* theInelasticProcess =
827  new G4DeuteronInelasticProcess("inelastic");
828  theInelasticProcess->AddDataSet( theGGNuclNuclData );
829  theInelasticProcess->RegisterMe( theFTFModel1 );
830  theInelasticProcess->RegisterMe( theBERTModel0 );
831  pmanager->AddDiscreteProcess( theInelasticProcess );
832  }
833 
834  else if (particleName == "triton")
835  {
836  // Elastic scattering
837  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
838  theElasticProcess->RegisterMe( elastic_lhep0 );
839  pmanager->AddDiscreteProcess( theElasticProcess );
840  // Inelastic scattering
841  G4TritonInelasticProcess* theInelasticProcess =
842  new G4TritonInelasticProcess("inelastic");
843  theInelasticProcess->AddDataSet( theGGNuclNuclData );
844  theInelasticProcess->RegisterMe( theFTFModel1 );
845  theInelasticProcess->RegisterMe( theBERTModel0 );
846  pmanager->AddDiscreteProcess( theInelasticProcess );
847  }
848  else if (particleName == "alpha")
849  {
850  // Elastic scattering
851  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
852  theElasticProcess->RegisterMe( elastic_lhep0 );
853  pmanager->AddDiscreteProcess( theElasticProcess );
854  // Inelastic scattering
855  G4AlphaInelasticProcess* theInelasticProcess =
856  new G4AlphaInelasticProcess("inelastic");
857  theInelasticProcess->AddDataSet( theGGNuclNuclData );
858  theInelasticProcess->RegisterMe( theFTFModel1 );
859  theInelasticProcess->RegisterMe( theBERTModel0 );
860  pmanager->AddDiscreteProcess( theInelasticProcess );
861  }
862 
863  }
864 }
865 
866 
867 // Decays ///////////////////////////////////////////////////////////////////
868 #include "G4Decay.hh"
869 #include "G4RadioactiveDecay.hh"
870 #include "G4IonTable.hh"
871 #include "G4Ions.hh"
872 
874 
875  // Add Decay Process
876  G4Decay* theDecayProcess = new G4Decay();
877  theParticleIterator->reset();
878  while( (*theParticleIterator)() )
879  {
880  G4ParticleDefinition* particle = theParticleIterator->value();
881  G4ProcessManager* pmanager = particle->GetProcessManager();
882 
883  if (theDecayProcess->IsApplicable(*particle) && !particle->IsShortLived())
884  {
885  pmanager ->AddProcess(theDecayProcess);
886  // set ordering for PostStepDoIt and AtRestDoIt
887  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
888  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
889  }
890  }
891 
892  // Declare radioactive decay to the GenericIon in the IonTable.
893  const G4IonTable *theIonTable =
895  G4RadioactiveDecay *theRadioactiveDecay = new G4RadioactiveDecay();
896 
897  for (G4int i=0; i<theIonTable->Entries(); i++)
898  {
899  G4String particleName = theIonTable->GetParticle(i)->GetParticleName();
900  G4String particleType = theIonTable->GetParticle(i)->GetParticleType();
901 
902  if (particleName == "GenericIon")
903  {
904  G4ProcessManager* pmanager =
905  theIonTable->GetParticle(i)->GetProcessManager();
906  pmanager->SetVerboseLevel(VerboseLevel);
907  pmanager ->AddProcess(theRadioactiveDecay);
908  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxPostStep);
909  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxAtRest);
910  }
911  }
912 }
913 
914 // Cuts /////////////////////////////////////////////////////////////////////
916 {
917 
918  if (verboseLevel >1)
919  G4cout << "DMXPhysicsList::SetCuts:";
920 
921  if (verboseLevel>0){
922  G4cout << "DMXPhysicsList::SetCuts:";
923  G4cout << "CutLength : "
924  << G4BestUnit(defaultCutValue,"Length") << G4endl;
925  }
926 
927  //special for low energy physics
928  G4double lowlimit=250*eV;
930 
931  // set cut values for gamma at first and for e- second and next for e+,
932  // because some processes for e+/e- need cut values for gamma
933  SetCutValue(cutForGamma, "gamma");
936 
938 }
939 
static G4Geantino * GeantinoDefinition()
Definition: G4Geantino.cc:82
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static const double MeV
Definition: G4SIunits.hh:211
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
void ConstructMyLeptons()
void SetPIXE(G4bool val)
G4int Entries() const
Definition: G4IonTable.cc:1632
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:100
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()
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:214
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:99
void ConstructMyBosons()
static G4ProductionCutsTable * GetProductionCutsTable()
static const double eV
Definition: G4SIunits.hh:212
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()
static const double um
Definition: G4SIunits.hh:112
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:215
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:1586