Geant4  10.02.p03
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 
201  particleIterator->reset();
202  while( (*particleIterator)() ){
203  G4ParticleDefinition* particle = particleIterator->value();
204  G4ProcessManager* pmanager = particle->GetProcessManager();
205  G4String particleName = particle->GetParticleName();
206  // time cuts for ONLY neutrons:
207  if(particleName == "neutron")
208  pmanager->AddDiscreteProcess(new DMXMaxTimeCuts());
209  // Energy cuts to kill charged (embedded in method) particles:
210  pmanager->AddDiscreteProcess(new DMXMinEkineCuts());
211 
212  // Step limit applied to all particles:
213  pmanager->AddProcess(new G4StepLimiter, -1,-1,1);
214 
215  }
216 }
217 
218 
219 // Electromagnetic Processes ////////////////////////////////////////////////
220 // all charged particles
221 
222 // gamma
223 #include "G4PhotoElectricEffect.hh"
224 #include "G4LivermorePhotoElectricModel.hh"
225 
226 #include "G4ComptonScattering.hh"
228 
229 #include "G4GammaConversion.hh"
231 
232 #include "G4RayleighScattering.hh"
234 
235 
236 // e-
237 #include "G4eMultipleScattering.hh"
238 
239 #include "G4eIonisation.hh"
241 
242 #include "G4eBremsstrahlung.hh"
244 
245 
246 // e+
247 #include "G4eIonisation.hh"
248 #include "G4eBremsstrahlung.hh"
249 #include "G4eplusAnnihilation.hh"
250 
251 
252 // alpha and GenericIon and deuterons, triton, He3:
253 
254 //muon:
255 #include "G4MuIonisation.hh"
256 #include "G4MuBremsstrahlung.hh"
257 #include "G4MuPairProduction.hh"
258 #include "G4MuonMinusCapture.hh"
259 
260 //OTHERS:
261 #include "G4hIonisation.hh"
262 #include "G4hMultipleScattering.hh"
263 #include "G4hBremsstrahlung.hh"
264 #include "G4ionIonisation.hh"
266 
267 //em process options to allow msc step-limitation to be switched off
268 #include "G4EmProcessOptions.hh"
269 
271 
272  //set a finer grid of the physic tables in order to improve precision
273  //former LowEnergy models have 200 bins up to 100 GeV
274  G4EmProcessOptions opt;
275  opt.SetMaxEnergy(100*GeV);
276  opt.SetDEDXBinning(200);
277  opt.SetLambdaBinning(200);
278 
280  particleIterator->reset();
281  while( (*particleIterator)() ){
282  G4ParticleDefinition* particle = particleIterator->value();
283  G4ProcessManager* pmanager = particle->GetProcessManager();
284  G4String particleName = particle->GetParticleName();
285  G4String particleType = particle->GetParticleType();
286  G4double charge = particle->GetPDGCharge();
287 
288  if (particleName == "gamma")
289  {
290  //gamma
291  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
292  theRayleigh->SetEmModel(new G4LivermoreRayleighModel()); //not strictly necessary
293  pmanager->AddDiscreteProcess(theRayleigh);
294 
295  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
296  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
297  pmanager->AddDiscreteProcess(thePhotoElectricEffect);
298 
299  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
300  theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
301  pmanager->AddDiscreteProcess(theComptonScattering);
302 
303  G4GammaConversion* theGammaConversion = new G4GammaConversion();
304  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
305  pmanager->AddDiscreteProcess(theGammaConversion);
306 
307  }
308  else if (particleName == "e-")
309  {
310  //electron
311  // process ordering: AddProcess(name, at rest, along step, post step)
312  // Multiple scattering
315  pmanager->AddProcess(msc,-1, 1, 1);
316 
317  // Ionisation
318  G4eIonisation* eIonisation = new G4eIonisation();
319  eIonisation->SetEmModel(new G4LivermoreIonisationModel());
320  eIonisation->SetStepFunction(0.2, 100*um); //improved precision in tracking
321  pmanager->AddProcess(eIonisation,-1, 2, 2);
322 
323  // Bremsstrahlung
324  G4eBremsstrahlung* eBremsstrahlung = new G4eBremsstrahlung();
325  eBremsstrahlung->SetEmModel(new G4LivermoreBremsstrahlungModel());
326  pmanager->AddProcess(eBremsstrahlung, -1,-3, 3);
327  }
328  else if (particleName == "e+")
329  {
330  //positron
333  pmanager->AddProcess(msc,-1, 1, 1);
334 
335  // Ionisation
336  G4eIonisation* eIonisation = new G4eIonisation();
337  eIonisation->SetStepFunction(0.2, 100*um); //
338  pmanager->AddProcess(eIonisation, -1, 2, 2);
339 
340  //Bremsstrahlung (use default, no low-energy available)
341  pmanager->AddProcess(new G4eBremsstrahlung(), -1,-1, 3);
342 
343  //Annihilation
344  pmanager->AddProcess(new G4eplusAnnihilation(),0,-1, 4);
345  }
346  else if( particleName == "mu+" ||
347  particleName == "mu-" )
348  {
349  //muon
350  pmanager->AddProcess(new G4eMultipleScattering, -1, 1, 1);
351  pmanager->AddProcess(new G4MuIonisation(), -1, 2, 2);
352  pmanager->AddProcess(new G4MuBremsstrahlung(), -1,-1, 3);
353  pmanager->AddProcess(new G4MuPairProduction(), -1,-1, 4);
354  if( particleName == "mu-" )
355  pmanager->AddProcess(new G4MuonMinusCapture(), 0,-1,-1);
356  }
357  else if (particleName == "proton" ||
358  particleName == "pi+" ||
359  particleName == "pi-")
360  {
361  //multiple scattering
362  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
363 
364  //ionisation
365  G4hIonisation* hIonisation = new G4hIonisation();
366  hIonisation->SetStepFunction(0.2, 50*um);
367  pmanager->AddProcess(hIonisation, -1, 2, 2);
368 
369  //bremmstrahlung
370  pmanager->AddProcess(new G4hBremsstrahlung, -1,-3, 3);
371  }
372  else if(particleName == "alpha" ||
373  particleName == "deuteron" ||
374  particleName == "triton" ||
375  particleName == "He3")
376  {
377  //multiple scattering
378  pmanager->AddProcess(new G4hMultipleScattering,-1,1,1);
379 
380  //ionisation
381  G4ionIonisation* ionIoni = new G4ionIonisation();
382  ionIoni->SetStepFunction(0.1, 20*um);
383  pmanager->AddProcess(ionIoni, -1, 2, 2);
384  }
385  else if (particleName == "GenericIon")
386  {
387  // OBJECT may be dynamically created as either a GenericIon or nucleus
388  // G4Nucleus exists and therefore has particle type nucleus
389  // genericIon:
390 
391  //multiple scattering
392  pmanager->AddProcess(new G4hMultipleScattering,-1,1,1);
393 
394  //ionisation
395  G4ionIonisation* ionIoni = new G4ionIonisation();
396  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
397  ionIoni->SetStepFunction(0.1, 20*um);
398  pmanager->AddProcess(ionIoni, -1, 2, 2);
399  }
400 
401  else if ((!particle->IsShortLived()) &&
402  (charge != 0.0) &&
403  (particle->GetParticleName() != "chargedgeantino"))
404  {
405  //all others charged particles except geantino
406  G4hMultipleScattering* aMultipleScattering = new G4hMultipleScattering();
407  G4hIonisation* ahadronIon = new G4hIonisation();
408 
409  //multiple scattering
410  pmanager->AddProcess(aMultipleScattering,-1,1,1);
411 
412  //ionisation
413  pmanager->AddProcess(ahadronIon, -1,2,2);
414  }
415 
416  }
417 
418  // turn off msc step-limitation - especially as electron cut 1nm
420 
421  // switch on fluorescence, PIXE and Auger:
422  opt.SetFluo(true);
423  opt.SetPIXE(true);
424  opt.SetAuger(true);
425 
426 }
427 
428 
429 // Optical Processes ////////////////////////////////////////////////////////
430 #include "G4Scintillation.hh"
431 #include "G4OpAbsorption.hh"
432 //#include "G4OpRayleigh.hh"
433 #include "G4OpBoundaryProcess.hh"
434 
436 {
437  // default scintillation process
438  G4Scintillation* theScintProcessDef = new G4Scintillation("Scintillation");
439  // theScintProcessDef->DumpPhysicsTable();
440  theScintProcessDef->SetTrackSecondariesFirst(true);
441  theScintProcessDef->SetScintillationYieldFactor(1.0); //
442  theScintProcessDef->SetScintillationExcitationRatio(0.0); //
443  theScintProcessDef->SetVerboseLevel(OpVerbLevel);
444 
445  // scintillation process for alpha:
446  G4Scintillation* theScintProcessAlpha = new G4Scintillation("Scintillation");
447  // theScintProcessNuc->DumpPhysicsTable();
448  theScintProcessAlpha->SetTrackSecondariesFirst(true);
449  theScintProcessAlpha->SetScintillationYieldFactor(1.1);
450  theScintProcessAlpha->SetScintillationExcitationRatio(1.0);
451  theScintProcessAlpha->SetVerboseLevel(OpVerbLevel);
452 
453  // scintillation process for heavy nuclei
454  G4Scintillation* theScintProcessNuc = new G4Scintillation("Scintillation");
455  // theScintProcessNuc->DumpPhysicsTable();
456  theScintProcessNuc->SetTrackSecondariesFirst(true);
457  theScintProcessNuc->SetScintillationYieldFactor(0.2);
458  theScintProcessNuc->SetScintillationExcitationRatio(1.0);
459  theScintProcessNuc->SetVerboseLevel(OpVerbLevel);
460 
461  // optical processes
462  G4OpAbsorption* theAbsorptionProcess = new G4OpAbsorption();
463  // G4OpRayleigh* theRayleighScatteringProcess = new G4OpRayleigh();
464  G4OpBoundaryProcess* theBoundaryProcess = new G4OpBoundaryProcess();
465  // theAbsorptionProcess->DumpPhysicsTable();
466  // theRayleighScatteringProcess->DumpPhysicsTable();
467  theAbsorptionProcess->SetVerboseLevel(OpVerbLevel);
468  // theRayleighScatteringProcess->SetVerboseLevel(OpVerbLevel);
469  theBoundaryProcess->SetVerboseLevel(OpVerbLevel);
470 
472  particleIterator->reset();
473  while( (*particleIterator)() )
474  {
475  G4ParticleDefinition* particle = particleIterator->value();
476  G4ProcessManager* pmanager = particle->GetProcessManager();
477  G4String particleName = particle->GetParticleName();
478  if (theScintProcessDef->IsApplicable(*particle)) {
479  // if(particle->GetPDGMass() > 5.0*GeV)
480  if(particle->GetParticleName() == "GenericIon") {
481  pmanager->AddProcess(theScintProcessNuc); // AtRestDiscrete
482  pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxAtRest);
483  pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxPostStep);
484  }
485  else if(particle->GetParticleName() == "alpha") {
486  pmanager->AddProcess(theScintProcessAlpha);
487  pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxAtRest);
488  pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxPostStep);
489  }
490  else {
491  pmanager->AddProcess(theScintProcessDef);
492  pmanager->SetProcessOrderingToLast(theScintProcessDef,idxAtRest);
493  pmanager->SetProcessOrderingToLast(theScintProcessDef,idxPostStep);
494  }
495  }
496 
497  if (particleName == "opticalphoton") {
498  pmanager->AddDiscreteProcess(theAbsorptionProcess);
499  // pmanager->AddDiscreteProcess(theRayleighScatteringProcess);
500  pmanager->AddDiscreteProcess(theBoundaryProcess);
501  }
502  }
503 }
504 
505 
506 // Hadronic processes ////////////////////////////////////////////////////////
507 
508 // Elastic processes:
509 #include "G4HadronElasticProcess.hh"
510 #include "G4ChipsElasticModel.hh"
511 #include "G4ElasticHadrNucleusHE.hh"
512 
513 // Inelastic processes:
527 
528 // High energy FTFP model and Bertini cascade
529 #include "G4FTFModel.hh"
531 #include "G4ExcitedStringDecay.hh"
532 #include "G4PreCompoundModel.hh"
534 #include "G4TheoFSGenerator.hh"
535 #include "G4CascadeInterface.hh"
536 
537 // Cross sections
538 #include "G4VCrossSectionDataSet.hh"
540 
541 #include "G4CrossSectionElastic.hh"
542 #include "G4BGGPionElasticXS.hh"
543 #include "G4AntiNuclElastic.hh"
544 
547 #include "G4CrossSectionPairGG.hh"
551 
552 #include "G4HadronElastic.hh"
553 #include "G4HadronCaptureProcess.hh"
554 
555 // Neutron high-precision models: <20 MeV
556 #include "G4ParticleHPElastic.hh"
558 #include "G4ParticleHPCapture.hh"
560 #include "G4ParticleHPInelastic.hh"
562 
563 // Stopping processes
567 
568 
569 
571 {
572  //Elastic models
573  const G4double elastic_elimitPi = 1.0*GeV;
574 
575  G4HadronElastic* elastic_lhep0 = new G4HadronElastic();
576  G4HadronElastic* elastic_lhep1 = new G4HadronElastic();
577  elastic_lhep1->SetMaxEnergy( elastic_elimitPi );
578  G4ChipsElasticModel* elastic_chip = new G4ChipsElasticModel();
579  G4ElasticHadrNucleusHE* elastic_he = new G4ElasticHadrNucleusHE();
580  elastic_he->SetMinEnergy( elastic_elimitPi );
581 
582 
583  // Inelastic scattering
584  const G4double theFTFMin0 = 0.0*GeV;
585  const G4double theFTFMin1 = 4.0*GeV;
586  const G4double theFTFMax = 100.0*TeV;
587  const G4double theBERTMin0 = 0.0*GeV;
588  const G4double theBERTMin1 = 19.0*MeV;
589  const G4double theBERTMax = 5.0*GeV;
590  const G4double theHPMin = 0.0*GeV;
591  const G4double theHPMax = 20.0*MeV;
592 
593  G4FTFModel * theStringModel = new G4FTFModel;
595  theStringModel->SetFragmentationModel( theStringDecay );
596  G4PreCompoundModel * thePreEquilib = new G4PreCompoundModel( new G4ExcitationHandler );
597  G4GeneratorPrecompoundInterface * theCascade = new G4GeneratorPrecompoundInterface( thePreEquilib );
598 
599  G4TheoFSGenerator * theFTFModel0 = new G4TheoFSGenerator( "FTFP" );
600  theFTFModel0->SetHighEnergyGenerator( theStringModel );
601  theFTFModel0->SetTransport( theCascade );
602  theFTFModel0->SetMinEnergy( theFTFMin0 );
603  theFTFModel0->SetMaxEnergy( theFTFMax );
604 
605  G4TheoFSGenerator * theFTFModel1 = new G4TheoFSGenerator( "FTFP" );
606  theFTFModel1->SetHighEnergyGenerator( theStringModel );
607  theFTFModel1->SetTransport( theCascade );
608  theFTFModel1->SetMinEnergy( theFTFMin1 );
609  theFTFModel1->SetMaxEnergy( theFTFMax );
610 
611  G4CascadeInterface * theBERTModel0 = new G4CascadeInterface;
612  theBERTModel0->SetMinEnergy( theBERTMin0 );
613  theBERTModel0->SetMaxEnergy( theBERTMax );
614 
615  G4CascadeInterface * theBERTModel1 = new G4CascadeInterface;
616  theBERTModel1->SetMinEnergy( theBERTMin1 );
617  theBERTModel1->SetMaxEnergy( theBERTMax );
618 
621  G4ComponentGGNuclNuclXsc * ggNuclNuclXsec = new G4ComponentGGNuclNuclXsc();
622  G4VCrossSectionDataSet * theGGNuclNuclData = new G4CrossSectionInelastic(ggNuclNuclXsec);
623 
625  particleIterator->reset();
626  while ((*particleIterator)())
627  {
628  G4ParticleDefinition* particle = particleIterator->value();
629  G4ProcessManager* pmanager = particle->GetProcessManager();
630  G4String particleName = particle->GetParticleName();
631 
632  if (particleName == "pi+")
633  {
634  // Elastic scattering
635  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
636  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
637  theElasticProcess->RegisterMe( elastic_lhep1 );
638  theElasticProcess->RegisterMe( elastic_he );
639  pmanager->AddDiscreteProcess( theElasticProcess );
640  //Inelastic scattering
641  G4PionPlusInelasticProcess* theInelasticProcess =
642  new G4PionPlusInelasticProcess("inelastic");
643  theInelasticProcess->AddDataSet( thePiData );
644  theInelasticProcess->RegisterMe( theFTFModel1 );
645  theInelasticProcess->RegisterMe( theBERTModel0 );
646  pmanager->AddDiscreteProcess( theInelasticProcess );
647  }
648 
649  else if (particleName == "pi-")
650  {
651  // Elastic scattering
652  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
653  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
654  theElasticProcess->RegisterMe( elastic_lhep1 );
655  theElasticProcess->RegisterMe( elastic_he );
656  pmanager->AddDiscreteProcess( theElasticProcess );
657  //Inelastic scattering
658  G4PionMinusInelasticProcess* theInelasticProcess =
659  new G4PionMinusInelasticProcess("inelastic");
660  theInelasticProcess->AddDataSet( thePiData );
661  theInelasticProcess->RegisterMe( theFTFModel1 );
662  theInelasticProcess->RegisterMe( theBERTModel0 );
663  pmanager->AddDiscreteProcess( theInelasticProcess );
664  //Absorption
666  }
667 
668  else if (particleName == "kaon+")
669  {
670  // Elastic scattering
671  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
672  theElasticProcess->RegisterMe( elastic_lhep0 );
673  pmanager->AddDiscreteProcess( theElasticProcess );
674  // Inelastic scattering
675  G4KaonPlusInelasticProcess* theInelasticProcess =
676  new G4KaonPlusInelasticProcess("inelastic");
677  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
678  GetCrossSectionDataSet(G4ChipsKaonPlusInelasticXS::Default_Name()));
679  theInelasticProcess->RegisterMe( theFTFModel1 );
680  theInelasticProcess->RegisterMe( theBERTModel0 );
681  pmanager->AddDiscreteProcess( theInelasticProcess );
682  }
683 
684  else if (particleName == "kaon0S")
685  {
686  // Elastic scattering
687  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
688  theElasticProcess->RegisterMe( elastic_lhep0 );
689  pmanager->AddDiscreteProcess( theElasticProcess );
690  // Inelastic scattering
691  G4KaonZeroSInelasticProcess* theInelasticProcess =
692  new G4KaonZeroSInelasticProcess("inelastic");
693  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
694  GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
695  theInelasticProcess->RegisterMe( theFTFModel1 );
696  theInelasticProcess->RegisterMe( theBERTModel0 );
697  pmanager->AddDiscreteProcess( theInelasticProcess );
698  }
699 
700  else if (particleName == "kaon0L")
701  {
702  // Elastic scattering
703  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
704  theElasticProcess->RegisterMe( elastic_lhep0 );
705  pmanager->AddDiscreteProcess( theElasticProcess );
706  // Inelastic scattering
707  G4KaonZeroLInelasticProcess* theInelasticProcess =
708  new G4KaonZeroLInelasticProcess("inelastic");
709  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
710  GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
711  theInelasticProcess->RegisterMe( theFTFModel1 );
712  theInelasticProcess->RegisterMe( theBERTModel0 );
713  pmanager->AddDiscreteProcess( theInelasticProcess );
714  }
715 
716  else if (particleName == "kaon-")
717  {
718  // Elastic scattering
719  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
720  theElasticProcess->RegisterMe( elastic_lhep0 );
721  pmanager->AddDiscreteProcess( theElasticProcess );
722  // Inelastic scattering
723  G4KaonMinusInelasticProcess* theInelasticProcess =
724  new G4KaonMinusInelasticProcess("inelastic");
725  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
726  GetCrossSectionDataSet(G4ChipsKaonMinusInelasticXS::Default_Name()));
727  theInelasticProcess->RegisterMe( theFTFModel1 );
728  theInelasticProcess->RegisterMe( theBERTModel0 );
729  pmanager->AddDiscreteProcess( theInelasticProcess );
731  }
732 
733  else if (particleName == "proton")
734  {
735  // Elastic scattering
736  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
738  GetCrossSectionDataSet(G4ChipsProtonElasticXS::Default_Name()));
739  theElasticProcess->RegisterMe( elastic_chip );
740  pmanager->AddDiscreteProcess( theElasticProcess );
741  // Inelastic scattering
742  G4ProtonInelasticProcess* theInelasticProcess =
743  new G4ProtonInelasticProcess("inelastic");
744  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Proton::Proton() ) );
745  theInelasticProcess->RegisterMe( theFTFModel1 );
746  theInelasticProcess->RegisterMe( theBERTModel0 );
747  pmanager->AddDiscreteProcess( theInelasticProcess );
748  }
749  else if (particleName == "anti_proton")
750  {
751  // Elastic scattering
752  const G4double elastic_elimitAntiNuc = 100.0*MeV;
753  G4AntiNuclElastic* elastic_anuc = new G4AntiNuclElastic();
754  elastic_anuc->SetMinEnergy( elastic_elimitAntiNuc );
755  G4CrossSectionElastic* elastic_anucxs = new G4CrossSectionElastic( elastic_anuc->GetComponentCrossSection() );
756  G4HadronElastic* elastic_lhep2 = new G4HadronElastic();
757  elastic_lhep2->SetMaxEnergy( elastic_elimitAntiNuc );
758  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
759  theElasticProcess->AddDataSet( elastic_anucxs );
760  theElasticProcess->RegisterMe( elastic_lhep2 );
761  theElasticProcess->RegisterMe( elastic_anuc );
762  pmanager->AddDiscreteProcess( theElasticProcess );
763  // Inelastic scattering
764  G4AntiProtonInelasticProcess* theInelasticProcess =
765  new G4AntiProtonInelasticProcess("inelastic");
766  theInelasticProcess->AddDataSet( theAntiNucleonData );
767  theInelasticProcess->RegisterMe( theFTFModel0 );
768  pmanager->AddDiscreteProcess( theInelasticProcess );
769  // Absorption
771  }
772 
773  else if (particleName == "neutron") {
774  // elastic scattering
775  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
776  theElasticProcess->AddDataSet(G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsNeutronElasticXS::Default_Name()));
777  G4HadronElastic* elastic_neutronChipsModel = new G4ChipsElasticModel();
778  elastic_neutronChipsModel->SetMinEnergy( 19.0*MeV );
779  theElasticProcess->RegisterMe( elastic_neutronChipsModel );
780  G4ParticleHPElastic * theElasticNeutronHP = new G4ParticleHPElastic;
781  theElasticNeutronHP->SetMinEnergy( theHPMin );
782  theElasticNeutronHP->SetMaxEnergy( theHPMax );
783  theElasticProcess->RegisterMe( theElasticNeutronHP );
784  theElasticProcess->AddDataSet( new G4ParticleHPElasticData );
785  pmanager->AddDiscreteProcess( theElasticProcess );
786  // inelastic scattering
787  G4NeutronInelasticProcess* theInelasticProcess =
788  new G4NeutronInelasticProcess("inelastic");
789  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Neutron::Neutron() ) );
790  theInelasticProcess->RegisterMe( theFTFModel1 );
791  theInelasticProcess->RegisterMe( theBERTModel1 );
792  G4ParticleHPInelastic * theNeutronInelasticHPModel = new G4ParticleHPInelastic;
793  theNeutronInelasticHPModel->SetMinEnergy( theHPMin );
794  theNeutronInelasticHPModel->SetMaxEnergy( theHPMax );
795  theInelasticProcess->RegisterMe( theNeutronInelasticHPModel );
796  theInelasticProcess->AddDataSet( new G4ParticleHPInelasticData );
797  pmanager->AddDiscreteProcess(theInelasticProcess);
798  // capture
799  G4HadronCaptureProcess* theCaptureProcess =
801  G4ParticleHPCapture * theLENeutronCaptureModel = new G4ParticleHPCapture;
802  theLENeutronCaptureModel->SetMinEnergy(theHPMin);
803  theLENeutronCaptureModel->SetMaxEnergy(theHPMax);
804  theCaptureProcess->RegisterMe(theLENeutronCaptureModel);
805  theCaptureProcess->AddDataSet( new G4ParticleHPCaptureData);
806  pmanager->AddDiscreteProcess(theCaptureProcess);
807 
808  }
809  else if (particleName == "anti_neutron")
810  {
811  // Elastic scattering
812  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
813  theElasticProcess->RegisterMe( elastic_lhep0 );
814  pmanager->AddDiscreteProcess( theElasticProcess );
815  // Inelastic scattering (include annihilation on-fly)
816  G4AntiNeutronInelasticProcess* theInelasticProcess =
817  new G4AntiNeutronInelasticProcess("inelastic");
818  theInelasticProcess->AddDataSet( theAntiNucleonData );
819  theInelasticProcess->RegisterMe( theFTFModel0 );
820  pmanager->AddDiscreteProcess( theInelasticProcess );
821  }
822 
823  else if (particleName == "deuteron")
824  {
825  // Elastic scattering
826  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
827  theElasticProcess->RegisterMe( elastic_lhep0 );
828  pmanager->AddDiscreteProcess( theElasticProcess );
829  // Inelastic scattering
830  G4DeuteronInelasticProcess* theInelasticProcess =
831  new G4DeuteronInelasticProcess("inelastic");
832  theInelasticProcess->AddDataSet( theGGNuclNuclData );
833  theInelasticProcess->RegisterMe( theFTFModel1 );
834  theInelasticProcess->RegisterMe( theBERTModel0 );
835  pmanager->AddDiscreteProcess( theInelasticProcess );
836  }
837 
838  else if (particleName == "triton")
839  {
840  // Elastic scattering
841  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
842  theElasticProcess->RegisterMe( elastic_lhep0 );
843  pmanager->AddDiscreteProcess( theElasticProcess );
844  // Inelastic scattering
845  G4TritonInelasticProcess* theInelasticProcess =
846  new G4TritonInelasticProcess("inelastic");
847  theInelasticProcess->AddDataSet( theGGNuclNuclData );
848  theInelasticProcess->RegisterMe( theFTFModel1 );
849  theInelasticProcess->RegisterMe( theBERTModel0 );
850  pmanager->AddDiscreteProcess( theInelasticProcess );
851  }
852  else if (particleName == "alpha")
853  {
854  // Elastic scattering
855  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
856  theElasticProcess->RegisterMe( elastic_lhep0 );
857  pmanager->AddDiscreteProcess( theElasticProcess );
858  // Inelastic scattering
859  G4AlphaInelasticProcess* theInelasticProcess =
860  new G4AlphaInelasticProcess("inelastic");
861  theInelasticProcess->AddDataSet( theGGNuclNuclData );
862  theInelasticProcess->RegisterMe( theFTFModel1 );
863  theInelasticProcess->RegisterMe( theBERTModel0 );
864  pmanager->AddDiscreteProcess( theInelasticProcess );
865  }
866 
867  }
868 }
869 
870 
871 // Decays ///////////////////////////////////////////////////////////////////
872 #include "G4Decay.hh"
873 #include "G4RadioactiveDecay.hh"
874 #include "G4IonTable.hh"
875 #include "G4Ions.hh"
876 
878 
879  // Add Decay Process
880  G4Decay* theDecayProcess = new G4Decay();
882  particleIterator->reset();
883  while( (*particleIterator)() )
884  {
885  G4ParticleDefinition* particle = particleIterator->value();
886  G4ProcessManager* pmanager = particle->GetProcessManager();
887 
888  if (theDecayProcess->IsApplicable(*particle) && !particle->IsShortLived())
889  {
890  pmanager ->AddProcess(theDecayProcess);
891  // set ordering for PostStepDoIt and AtRestDoIt
892  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
893  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
894  }
895  }
896 
897  // Declare radioactive decay to the GenericIon in the IonTable.
898  const G4IonTable *theIonTable =
900  G4RadioactiveDecay *theRadioactiveDecay = new G4RadioactiveDecay();
901 
902  for (G4int i=0; i<theIonTable->Entries(); i++)
903  {
904  G4String particleName = theIonTable->GetParticle(i)->GetParticleName();
905  G4String particleType = theIonTable->GetParticle(i)->GetParticleType();
906 
907  if (particleName == "GenericIon")
908  {
909  G4ProcessManager* pmanager =
910  theIonTable->GetParticle(i)->GetProcessManager();
911  pmanager->SetVerboseLevel(VerboseLevel);
912  pmanager ->AddProcess(theRadioactiveDecay);
913  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxPostStep);
914  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxAtRest);
915  }
916  }
917 }
918 
919 // Cuts /////////////////////////////////////////////////////////////////////
921 {
922 
923  if (verboseLevel >1)
924  G4cout << "DMXPhysicsList::SetCuts:";
925 
926  if (verboseLevel>0){
927  G4cout << "DMXPhysicsList::SetCuts:";
928  G4cout << "CutLength : "
929  << G4BestUnit(defaultCutValue,"Length") << G4endl;
930  }
931 
932  //special for low energy physics
933  G4double lowlimit=250*eV;
935 
936  // set cut values for gamma at first and for e- second and next for e+,
937  // because some processes for e+/e- need cut values for gamma
938  SetCutValue(cutForGamma, "gamma");
941 
943 }
944 
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)
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)
const G4String & GetParticleType() const
G4ProcessManager * GetProcessManager() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static void ConstructParticle()
G4ParticleDefinition * GetParticle(G4int index) const
Definition: G4IonTable.cc:1586
static void ConstructParticle()
int G4int
Definition: G4Types.hh:78
G4ComponentAntiNuclNuclearXS * GetComponentCrossSection()
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)
const G4String & GetParticleName() const
G4int Entries() const
Definition: G4IonTable.cc:1632
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
void SetLambdaBinning(G4int val)
G4IonTable * GetIonTable() const
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
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)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
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)
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
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
static G4NeutrinoMu * NeutrinoMuDefinition()
Definition: G4NeutrinoMu.cc:80
G4double GetPDGCharge() const
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