Geant4  10.03
LBE.icc
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 //
29 // For information related to this code contact: Alex Howard
30 // e-mail: alexander.howard@cern.ch
31 // --------------------------------------------------------------
32 // Comments
33 //
34 // Underground Advanced
35 //
36 // This physics list is taken from the underground_physics example with small
37 // modifications. It is an example of a "flat" physics list with no dependence
38 // on builders. The physics covered would be suitable for a low background
39 // experiment including the neutron_hp package
40 //
41 //
42 //
43 // PhysicsList program
44 //
45 // Modified:
46 //
47 // 14-02-03 Fix bugs in msc and hIon instanciation + cut per region
48 // 16-08-10 Remove inclusion of obsolete class of G4ParticleWithCuts
49 // 20-10-10 Migrate LowEnergy process to Livermore models, LP
50 // 28-03-13 Replace LEP/HEP with FTFP+BERT (A.R.)
51 // --------------------------------------------------------------
52 
53 #include <iomanip>
54 
55 #include "globals.hh"
56 #include "G4ios.hh"
57 #include "G4ProcessManager.hh"
58 #include "G4ProcessVector.hh"
59 
60 #include "G4ParticleTypes.hh"
61 #include "G4ParticleTable.hh"
62 #include "G4ProductionCutsTable.hh"
63 
64 #include "G4UserLimits.hh"
65 #include "G4DataQuestionaire.hh"
66 #include "G4WarnPLStatus.hh"
67 
68 // Builder for all stopping processes
69 #include "G4StoppingPhysics.hh"
70 
71 // Constructor /////////////////////////////////////////////////////////////
72 template<class T> TLBE<T>::TLBE(G4int ver) :T()
73 {
74 
75  G4DataQuestionaire it(photon, lowenergy, neutron, radioactive);
76  G4cout << "You are using the simulation engine: LBE 5.3"<<G4endl;
77  G4cout <<G4endl<<G4endl;
78  this->defaultCutValue = 1.0*CLHEP::micrometer; //
79  cutForGamma = this->defaultCutValue;
80 // cutForElectron = 1.0*CLHEP::nanometer;
81  cutForElectron = 1.0*CLHEP::micrometer;
82  cutForPositron = this->defaultCutValue;
83  //not used:
84  // cutForProton = this->defaultCutValue;
85  // cutForAlpha = 1.0*CLHEP::nanometer;
86  // cutForGenericIon = 1.0*CLHEP::nanometer;
87 
88  stoppingPhysics = new G4StoppingPhysics;
89 
90  VerboseLevel = ver;
91  OpVerbLevel = 0;
92 
93  this->SetVerboseLevel(VerboseLevel);
94 }
95 
96 
97 // Destructor //////////////////////////////////////////////////////////////
98 template<class T> TLBE<T>::~TLBE()
99 {
100  delete stoppingPhysics;
101 }
102 
103 
104 // Construct Particles /////////////////////////////////////////////////////
105  template<class T> void TLBE<T>::ConstructParticle()
106 {
107 
108  // In this method, static member functions should be called
109  // for all particles which you want to use.
110  // This ensures that objects of these particle types will be
111  // created in the program.
112 
113  ConstructMyBosons();
114  ConstructMyLeptons();
115  ConstructMyMesons();
116  ConstructMyBaryons();
117  ConstructMyIons();
118  ConstructMyShortLiveds();
119  stoppingPhysics->ConstructParticle(); // Anything not included above
120 }
121 
122 
123 // construct Bosons://///////////////////////////////////////////////////
124  template<class T> void TLBE<T>::ConstructMyBosons()
125 {
126  // pseudo-particles
127  G4Geantino::GeantinoDefinition();
128  G4ChargedGeantino::ChargedGeantinoDefinition();
129 
130  // gamma
131  G4Gamma::GammaDefinition();
132 
133  //OpticalPhotons
134  G4OpticalPhoton::OpticalPhotonDefinition();
135 }
136 
137 
138 // construct Leptons://///////////////////////////////////////////////////
139  template<class T> void TLBE<T>::ConstructMyLeptons()
140 {
141  // leptons
142  G4Electron::ElectronDefinition();
143  G4Positron::PositronDefinition();
144  G4MuonPlus::MuonPlusDefinition();
145  G4MuonMinus::MuonMinusDefinition();
146 
147  G4NeutrinoE::NeutrinoEDefinition();
148  G4AntiNeutrinoE::AntiNeutrinoEDefinition();
149  G4NeutrinoMu::NeutrinoMuDefinition();
150  G4AntiNeutrinoMu::AntiNeutrinoMuDefinition();
151 }
152 
153 #include "G4MesonConstructor.hh"
154 #include "G4BaryonConstructor.hh"
155 #include "G4IonConstructor.hh"
156 
157 
158 // construct Mesons://///////////////////////////////////////////////////
159  template<class T> void TLBE<T>::ConstructMyMesons()
160 {
161  // mesons
162  G4MesonConstructor mConstructor;
163  mConstructor.ConstructParticle();
164 
165 }
166 
167 
168 // construct Baryons://///////////////////////////////////////////////////
169  template<class T> void TLBE<T>::ConstructMyBaryons()
170 {
171  // baryons
172  G4BaryonConstructor bConstructor;
173  bConstructor.ConstructParticle();
174 
175 }
176 
177 
178 // construct Ions://///////////////////////////////////////////////////
179  template<class T> void TLBE<T>::ConstructMyIons()
180 {
181  // ions
182  G4IonConstructor iConstructor;
183  iConstructor.ConstructParticle();
184 
185 }
186 
187 // construct Shortliveds://///////////////////////////////////////////////////
188  template<class T> void TLBE<T>::ConstructMyShortLiveds()
189 {
190  // ShortLiveds
191  ;
192 }
193 
194 
195 
196 
197 // Construct Processes //////////////////////////////////////////////////////
198  template<class T> void TLBE<T>::ConstructProcess()
199 {
200  AddTransportation();
201  ConstructEM();
202  ConstructOp();
203  ConstructHad();
204  ConstructGeneral();
205 }
206 
207 
208 // Transportation ///////////////////////////////////////////////////////////
209 #include "G4MaxTimeCuts.hh"
210 #include "G4MinEkineCuts.hh"
211 
212  template<class T> void TLBE<T>::AddTransportation() {
213 
214  G4VUserPhysicsList::AddTransportation();
215 
216  auto myParticleIterator=G4ParticleTable::GetParticleTable()->GetIterator();
217  myParticleIterator->reset();
218  while( (*(myParticleIterator))() ){
219  G4ParticleDefinition* particle = myParticleIterator->value();
220  G4ProcessManager* pmanager = particle->GetProcessManager();
221  G4String particleName = particle->GetParticleName();
222  // time cuts for ONLY neutrons:
223  if(particleName == "neutron")
224  pmanager->AddDiscreteProcess(new G4MaxTimeCuts());
225  // Energy cuts to kill charged (embedded in method) particles:
226  pmanager->AddDiscreteProcess(new G4MinEkineCuts());
227  }
228 }
229 
230 
231 // Electromagnetic Processes ////////////////////////////////////////////////
232 // all charged particles
233 
234 #include "G4eMultipleScattering.hh"
235 #include "G4MuMultipleScattering.hh"
236 #include "G4hMultipleScattering.hh"
237 
238 // gamma. Use Livermore models
239 #include "G4PhotoElectricEffect.hh"
240 #include "G4LivermorePhotoElectricModel.hh"
241 #include "G4ComptonScattering.hh"
242 #include "G4LivermoreComptonModel.hh"
243 #include "G4GammaConversion.hh"
244 #include "G4LivermoreGammaConversionModel.hh"
245 #include "G4RayleighScattering.hh"
246 #include "G4LivermoreRayleighModel.hh"
247 
248 
249 // e-
250 #include "G4eMultipleScattering.hh"
251 #include "G4UniversalFluctuation.hh"
252 #include "G4UrbanMscModel.hh"
253 
254 #include "G4eIonisation.hh"
255 #include "G4LivermoreIonisationModel.hh"
256 
257 #include "G4eBremsstrahlung.hh"
258 #include "G4LivermoreBremsstrahlungModel.hh"
259 
260 // e+
261 #include "G4eplusAnnihilation.hh"
262 
263 
264 // alpha and GenericIon and deuterons, triton, He3:
265 #include "G4ionIonisation.hh"
266 #include "G4hIonisation.hh"
267 #include "G4hBremsstrahlung.hh"
268 //
269 #include "G4IonParametrisedLossModel.hh"
270 #include "G4NuclearStopping.hh"
271 #include "G4EnergyLossTables.hh"
272 
273 //muon:
274 #include "G4MuIonisation.hh"
275 #include "G4MuBremsstrahlung.hh"
276 #include "G4MuPairProduction.hh"
277 #include "G4MuonMinusCapture.hh"
278 
279 //OTHERS:
280 //#include "G4hIonisation.hh" // standard hadron ionisation
281 
282  template<class T> void TLBE<T>::ConstructEM() {
283 
284  // models & processes:
285  // Use Livermore models up to 20 MeV, and standard
286  // models for higher energy
287  G4double LivermoreHighEnergyLimit = 20*CLHEP::MeV;
288  //
289  auto myParticleIterator=G4ParticleTable::GetParticleTable()->GetIterator();
290  myParticleIterator->reset();
291  while( (*(myParticleIterator))() ){
292  G4ParticleDefinition* particle = myParticleIterator->value();
293  G4ProcessManager* pmanager = particle->GetProcessManager();
294  G4String particleName = particle->GetParticleName();
295  G4String particleType = particle->GetParticleType();
296  G4double charge = particle->GetPDGCharge();
297 
298  if (particleName == "gamma")
299  {
300  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
301  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
302  new G4LivermorePhotoElectricModel();
303  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
304  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
305  pmanager->AddDiscreteProcess(thePhotoElectricEffect);
306 
307  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
308  G4LivermoreComptonModel* theLivermoreComptonModel =
309  new G4LivermoreComptonModel();
310  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
311  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
312  pmanager->AddDiscreteProcess(theComptonScattering);
313 
314  G4GammaConversion* theGammaConversion = new G4GammaConversion();
315  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
316  new G4LivermoreGammaConversionModel();
317  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
318  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
319  pmanager->AddDiscreteProcess(theGammaConversion);
320 
321  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
322  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
323  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
324  theRayleigh->AddEmModel(0, theRayleighModel);
325  pmanager->AddDiscreteProcess(theRayleigh);
326 
327  }
328  else if (particleName == "e-")
329  {
330  //electron
331  // process ordering: AddProcess(name, at rest, along step, post step)
332  // -1 = not implemented, then ordering
333  G4eMultipleScattering* msc = new G4eMultipleScattering();
334  //msc->AddEmModel(0, new G4UrbanMscModel());
335  msc->SetStepLimitType(fUseDistanceToBoundary);
336  pmanager->AddProcess(msc, -1, 1, 1);
337 
338  // Ionisation
339  G4eIonisation* eIoni = new G4eIonisation();
340  G4LivermoreIonisationModel* theIoniLivermore = new
341  G4LivermoreIonisationModel();
342  theIoniLivermore->SetHighEnergyLimit(1*CLHEP::MeV);
343  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
344  eIoni->SetStepFunction(0.2, 100*CLHEP::um); //
345  pmanager->AddProcess(eIoni, -1, 2, 2);
346 
347  // Bremsstrahlung
348  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
349  G4LivermoreBremsstrahlungModel* theBremLivermore = new
350  G4LivermoreBremsstrahlungModel();
351  theBremLivermore->SetHighEnergyLimit(LivermoreHighEnergyLimit);
352  eBrem->AddEmModel(0, theBremLivermore);
353  pmanager->AddProcess(eBrem, -1,-3, 3);
354  }
355  else if (particleName == "e+")
356  {
357  //positron
358  G4eMultipleScattering* msc = new G4eMultipleScattering();
359  //msc->AddEmModel(0, new G4UrbanMscModel());
360  msc->SetStepLimitType(fUseDistanceToBoundary);
361  pmanager->AddProcess(msc, -1, 1, 1);
362  G4eIonisation* eIoni = new G4eIonisation();
363  eIoni->SetStepFunction(0.2, 100*CLHEP::um);
364  pmanager->AddProcess(eIoni, -1, 2, 2);
365  pmanager->AddProcess(new G4eBremsstrahlung, -1,-3, 3);
366  pmanager->AddProcess(new G4eplusAnnihilation,0,-1, 4);
367  }
368  else if( particleName == "mu+" ||
369  particleName == "mu-" )
370  {
371  //muon
372  G4MuMultipleScattering* aMultipleScattering = new G4MuMultipleScattering();
373  pmanager->AddProcess(aMultipleScattering, -1, 1, 1);
374  pmanager->AddProcess(new G4MuIonisation(), -1, 2, 2);
375  pmanager->AddProcess(new G4MuBremsstrahlung(), -1,-1, 3);
376  pmanager->AddProcess(new G4MuPairProduction(), -1,-1, 4);
377  if( particleName == "mu-" )
378  pmanager->AddProcess(new G4MuonMinusCapture(), 0,-1,-1);
379  }
380  else if (particleName == "GenericIon")
381  {
382  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
383  G4ionIonisation* ionIoni = new G4ionIonisation();
384  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
385  ionIoni->SetStepFunction(0.1, 10*CLHEP::um);
386  pmanager->AddProcess(ionIoni, -1, 2, 2);
387  pmanager->AddProcess(new G4NuclearStopping(), -1, 3,-1);
388  }
389  else if (particleName == "alpha" || particleName == "He3")
390  {
391  //MSC, ion-Ionisation, Nuclear Stopping
392  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
393 
394  G4ionIonisation* ionIoni = new G4ionIonisation();
395  ionIoni->SetStepFunction(0.1, 20*CLHEP::um);
396  pmanager->AddProcess(ionIoni, -1, 2, 2);
397  pmanager->AddProcess(new G4NuclearStopping(), -1, 3,-1);
398  }
399  else if (particleName == "proton" ||
400  particleName == "deuteron" ||
401  particleName == "triton" ||
402  particleName == "pi+" ||
403  particleName == "pi-" ||
404  particleName == "kaon+" ||
405  particleName == "kaon-")
406  {
407  //MSC, h-ionisation, bremsstrahlung
408  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
409  G4hIonisation* hIoni = new G4hIonisation();
410  hIoni->SetStepFunction(0.2, 50*CLHEP::um);
411  pmanager->AddProcess(hIoni, -1, 2, 2);
412  pmanager->AddProcess(new G4hBremsstrahlung, -1,-3, 3);
413  }
414  else if ((!particle->IsShortLived()) &&
415  (charge != 0.0) &&
416  (particle->GetParticleName() != "chargedgeantino"))
417  {
418  //all others charged particles except geantino
419  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
420  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
421  }
422 
423  }
424 }
425 
426 
427 // Optical Processes ////////////////////////////////////////////////////////
428 #include "G4Scintillation.hh"
429 #include "G4OpAbsorption.hh"
430 //#include "G4OpRayleigh.hh"
431 #include "G4OpBoundaryProcess.hh"
432 
433  template<class T> void TLBE<T>::ConstructOp()
434 {
435  // default scintillation process
436  //Coverity report: check that the process is actually used, if not must delete
437  G4bool theScintProcessDefNeverUsed = true;
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  G4bool theScintProcessAlphaNeverUsed = true;
447  G4Scintillation* theScintProcessAlpha = new G4Scintillation("Scintillation");
448  // theScintProcessNuc->DumpPhysicsTable();
449  theScintProcessAlpha->SetTrackSecondariesFirst(true);
450  theScintProcessAlpha->SetScintillationYieldFactor(1.1);
451  theScintProcessAlpha->SetScintillationExcitationRatio(1.0);
452  theScintProcessAlpha->SetVerboseLevel(OpVerbLevel);
453 
454  // scintillation process for heavy nuclei
455  G4bool theScintProcessNucNeverUsed = true;
456  G4Scintillation* theScintProcessNuc = new G4Scintillation("Scintillation");
457  // theScintProcessNuc->DumpPhysicsTable();
458  theScintProcessNuc->SetTrackSecondariesFirst(true);
459  theScintProcessNuc->SetScintillationYieldFactor(0.2);
460  theScintProcessNuc->SetScintillationExcitationRatio(1.0);
461  theScintProcessNuc->SetVerboseLevel(OpVerbLevel);
462 
463  // optical processes
464  G4bool theAbsorptionProcessNeverUsed = true;
465  G4OpAbsorption* theAbsorptionProcess = new G4OpAbsorption();
466  // G4OpRayleigh* theRayleighScatteringProcess = new G4OpRayleigh();
467  G4bool theBoundaryProcessNeverUsed = true;
468  G4OpBoundaryProcess* theBoundaryProcess = new G4OpBoundaryProcess();
469  // theAbsorptionProcess->DumpPhysicsTable();
470  // theRayleighScatteringProcess->DumpPhysicsTable();
471  theAbsorptionProcess->SetVerboseLevel(OpVerbLevel);
472  // theRayleighScatteringProcess->SetVerboseLevel(OpVerbLevel);
473  theBoundaryProcess->SetVerboseLevel(OpVerbLevel);
474 
475  auto myParticleIterator=G4ParticleTable::GetParticleTable()->GetIterator();
476  myParticleIterator->reset();
477  while( (*(myParticleIterator))() )
478  {
479  G4ParticleDefinition* particle = myParticleIterator->value();
480  G4ProcessManager* pmanager = particle->GetProcessManager();
481  G4String particleName = particle->GetParticleName();
482  if (theScintProcessDef->IsApplicable(*particle)) {
483  // if(particle->GetPDGMass() > 5.0*CLHEP::GeV)
484  if(particle->GetParticleName() == "GenericIon") {
485  pmanager->AddProcess(theScintProcessNuc); // AtRestDiscrete
486  pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxAtRest);
487  pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxPostStep);
488  theScintProcessNucNeverUsed = false;
489  }
490  else if(particle->GetParticleName() == "alpha") {
491  pmanager->AddProcess(theScintProcessAlpha);
492  pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxAtRest);
493  pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxPostStep);
494  theScintProcessAlphaNeverUsed = false;
495  }
496  else {
497  pmanager->AddProcess(theScintProcessDef);
498  pmanager->SetProcessOrderingToLast(theScintProcessDef,idxAtRest);
499  pmanager->SetProcessOrderingToLast(theScintProcessDef,idxPostStep);
500  theScintProcessDefNeverUsed = false;
501  }
502  }
503 
504  if (particleName == "opticalphoton") {
505  pmanager->AddDiscreteProcess(theAbsorptionProcess);
506  theAbsorptionProcessNeverUsed = false;
507  // pmanager->AddDiscreteProcess(theRayleighScatteringProcess);
508  theBoundaryProcessNeverUsed = false;
509  pmanager->AddDiscreteProcess(theBoundaryProcess);
510  }
511  }
512  if ( theScintProcessDefNeverUsed ) delete theScintProcessDef;
513  if ( theScintProcessAlphaNeverUsed ) delete theScintProcessAlpha;
514  if ( theScintProcessNucNeverUsed ) delete theScintProcessNuc;
515  if ( theBoundaryProcessNeverUsed ) delete theBoundaryProcess;
516  if ( theAbsorptionProcessNeverUsed ) delete theAbsorptionProcess;
517 }
518 
519 
520 // Hadronic processes ////////////////////////////////////////////////////////
521 
522 // Elastic processes:
523 #include "G4HadronElasticProcess.hh"
524 #include "G4HadronCaptureProcess.hh"
525 #include "G4HadronElastic.hh"
526 #include "G4ChipsElasticModel.hh"
527 #include "G4ElasticHadrNucleusHE.hh"
528 #include "G4AntiNuclElastic.hh"
529 #include "G4BGGPionElasticXS.hh"
530 #include "G4CrossSectionDataSetRegistry.hh"
531 #include "G4ChipsProtonElasticXS.hh"
532 #include "G4ChipsNeutronElasticXS.hh"
533 #include "G4ComponentAntiNuclNuclearXS.hh"
534 #include "G4CrossSectionElastic.hh"
535 
536 // Inelastic processes:
537 #include "G4PionPlusInelasticProcess.hh"
538 #include "G4PionMinusInelasticProcess.hh"
539 #include "G4KaonPlusInelasticProcess.hh"
540 #include "G4KaonZeroSInelasticProcess.hh"
541 #include "G4KaonZeroLInelasticProcess.hh"
542 #include "G4KaonMinusInelasticProcess.hh"
543 #include "G4ProtonInelasticProcess.hh"
544 #include "G4AntiProtonInelasticProcess.hh"
545 #include "G4NeutronInelasticProcess.hh"
546 #include "G4AntiNeutronInelasticProcess.hh"
547 #include "G4DeuteronInelasticProcess.hh"
548 #include "G4TritonInelasticProcess.hh"
549 #include "G4AlphaInelasticProcess.hh"
550 
551 // FTFP + BERT model
552 #include "G4TheoFSGenerator.hh"
553 #include "G4ExcitationHandler.hh"
554 #include "G4PreCompoundModel.hh"
555 #include "G4GeneratorPrecompoundInterface.hh"
556 #include "G4FTFModel.hh"
557 #include "G4LundStringFragmentation.hh"
558 #include "G4ExcitedStringDecay.hh"
559 #include "G4CascadeInterface.hh"
560 #include "G4CrossSectionInelastic.hh"
561 #include "G4PiNuclearCrossSection.hh"
562 #include "G4CrossSectionPairGG.hh"
563 #include "G4ChipsKaonMinusInelasticXS.hh"
564 #include "G4ChipsKaonPlusInelasticXS.hh"
565 #include "G4ChipsKaonZeroInelasticXS.hh"
566 #include "G4CrossSectionDataSetRegistry.hh"
567 #include "G4BGGNucleonInelasticXS.hh"
568 #include "G4ComponentAntiNuclNuclearXS.hh"
569 #include "G4ComponentGGNuclNuclXsc.hh"
570 
571 // Neutron high-precision models: <20 MeV
572 #include "G4ParticleHPElastic.hh"
573 #include "G4ParticleHPElasticData.hh"
574 #include "G4ParticleHPCapture.hh"
575 #include "G4ParticleHPCaptureData.hh"
576 #include "G4ParticleHPInelastic.hh"
577 #include "G4ParticleHPInelasticData.hh"
578 #include "G4NeutronCaptureXS.hh"
579 #include "G4NeutronRadCapture.hh"
580 
581 // Binary light ion cascade for alpha, deuteron and triton
582 #include "G4BinaryLightIonReaction.hh"
583 
584 // ConstructHad()
585 // Makes discrete physics processes for the hadrons, at present limited
586 // to those particles with GHEISHA interactions (INTRC > 0).
587 // The processes are: Elastic scattering and Inelastic scattering.
588 // F.W.Jones 09-JUL-1998
589  template<class T> void TLBE<T>::ConstructHad()
590 {
591  // Elastic scattering
592  const G4double elastic_elimitPi = 1.0*CLHEP::GeV;
593 
594  G4HadronElastic* elastic_lhep0 = new G4HadronElastic();
595  G4HadronElastic* elastic_lhep1 = new G4HadronElastic();
596  elastic_lhep1->SetMaxEnergy( elastic_elimitPi );
597 
598  G4ChipsElasticModel* elastic_chip = new G4ChipsElasticModel();
599 
600  G4ElasticHadrNucleusHE* elastic_he = new G4ElasticHadrNucleusHE();
601  elastic_he->SetMinEnergy( elastic_elimitPi );
602 
603  // Inelastic scattering
604  const G4double theFTFMin0 = 0.0*CLHEP::GeV;
605  const G4double theFTFMin1 = 4.0*CLHEP::GeV;
606  const G4double theFTFMax = 100.0*CLHEP::TeV;
607  const G4double theBERTMin0 = 0.0*CLHEP::GeV;
608  const G4double theBERTMin1 = 19.0*CLHEP::MeV;
609  const G4double theBERTMax = 5.0*CLHEP::GeV;
610  const G4double theHPMin = 0.0*CLHEP::GeV;
611  const G4double theHPMax = 20.0*CLHEP::MeV;
612  const G4double theIonBCMin = 0.0*CLHEP::GeV;
613  const G4double theIonBCMax = 5.0*CLHEP::GeV;
614 
615 
616  G4FTFModel * theStringModel = new G4FTFModel;
617  G4ExcitedStringDecay * theStringDecay = new G4ExcitedStringDecay( new G4LundStringFragmentation );
618  theStringModel->SetFragmentationModel( theStringDecay );
619  G4PreCompoundModel * thePreEquilib = new G4PreCompoundModel( new G4ExcitationHandler );
620  G4GeneratorPrecompoundInterface * theCascade = new G4GeneratorPrecompoundInterface( thePreEquilib );
621 
622  G4TheoFSGenerator * theFTFModel0 = new G4TheoFSGenerator( "FTFP" );
623  theFTFModel0->SetHighEnergyGenerator( theStringModel );
624  theFTFModel0->SetTransport( theCascade );
625  theFTFModel0->SetMinEnergy( theFTFMin0 );
626  theFTFModel0->SetMaxEnergy( theFTFMax );
627 
628  G4TheoFSGenerator * theFTFModel1 = new G4TheoFSGenerator( "FTFP" );
629  theFTFModel1->SetHighEnergyGenerator( theStringModel );
630  theFTFModel1->SetTransport( theCascade );
631  theFTFModel1->SetMinEnergy( theFTFMin1 );
632  theFTFModel1->SetMaxEnergy( theFTFMax );
633 
634  G4CascadeInterface * theBERTModel0 = new G4CascadeInterface;
635  theBERTModel0->SetMinEnergy( theBERTMin0 );
636  theBERTModel0->SetMaxEnergy( theBERTMax );
637 
638  G4CascadeInterface * theBERTModel1 = new G4CascadeInterface;
639  theBERTModel1->SetMinEnergy( theBERTMin1 );
640  theBERTModel1->SetMaxEnergy( theBERTMax );
641 
642  // Binary Cascade
643  G4BinaryLightIonReaction * theIonBC = new G4BinaryLightIonReaction( thePreEquilib );
644  theIonBC->SetMinEnergy( theIonBCMin );
645  theIonBC->SetMaxEnergy( theIonBCMax );
646 
647  G4VCrossSectionDataSet * thePiData = new G4CrossSectionPairGG(
648  (G4PiNuclearCrossSection*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4PiNuclearCrossSection::Default_Name()), 91*CLHEP::GeV );
649  G4VCrossSectionDataSet * theAntiNucleonData = new G4CrossSectionInelastic( new G4ComponentAntiNuclNuclearXS );
650  G4ComponentGGNuclNuclXsc * ggNuclNuclXsec = new G4ComponentGGNuclNuclXsc();
651  G4VCrossSectionDataSet * theGGNuclNuclData = new G4CrossSectionInelastic(ggNuclNuclXsec);
652 
653  auto myParticleIterator=G4ParticleTable::GetParticleTable()->GetIterator();
654  myParticleIterator->reset();
655  while ((*(myParticleIterator))())
656  {
657  G4ParticleDefinition* particle = myParticleIterator->value();
658  G4ProcessManager* pmanager = particle->GetProcessManager();
659  G4String particleName = particle->GetParticleName();
660 
661  if (particleName == "pi+")
662  {
663  // Elastic scattering
664  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
665  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
666  theElasticProcess->RegisterMe( elastic_lhep1 );
667  theElasticProcess->RegisterMe( elastic_he );
668  pmanager->AddDiscreteProcess( theElasticProcess );
669  // Inelastic scattering
670  G4PionPlusInelasticProcess* theInelasticProcess = new G4PionPlusInelasticProcess("inelastic");
671  theInelasticProcess->AddDataSet( thePiData );
672  theInelasticProcess->RegisterMe( theFTFModel1 );
673  theInelasticProcess->RegisterMe( theBERTModel0 );
674  pmanager->AddDiscreteProcess( theInelasticProcess );
675  }
676 
677  else if (particleName == "pi-")
678  {
679  // Elastic scattering
680  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
681  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
682  theElasticProcess->RegisterMe( elastic_lhep1 );
683  theElasticProcess->RegisterMe( elastic_he );
684  pmanager->AddDiscreteProcess( theElasticProcess );
685  // Inelastic scattering
686  G4PionMinusInelasticProcess* theInelasticProcess = new G4PionMinusInelasticProcess("inelastic");
687  theInelasticProcess->AddDataSet( thePiData );
688  theInelasticProcess->RegisterMe( theFTFModel1 );
689  theInelasticProcess->RegisterMe( theBERTModel0 );
690  pmanager->AddDiscreteProcess( theInelasticProcess );
691  }
692 
693  else if (particleName == "kaon+")
694  {
695  // Elastic scattering
696  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
697  theElasticProcess->RegisterMe( elastic_lhep0 );
698  pmanager->AddDiscreteProcess( theElasticProcess );
699  // Inelastic scattering
700  G4KaonPlusInelasticProcess* theInelasticProcess = new G4KaonPlusInelasticProcess("inelastic");
701  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonPlusInelasticXS::Default_Name()));
702  theInelasticProcess->RegisterMe( theFTFModel1 );
703  theInelasticProcess->RegisterMe( theBERTModel0 );
704  pmanager->AddDiscreteProcess( theInelasticProcess );
705  }
706 
707  else if (particleName == "kaon0S")
708  {
709  // Elastic scattering
710  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
711  theElasticProcess->RegisterMe( elastic_lhep0 );
712  pmanager->AddDiscreteProcess( theElasticProcess );
713  // Inelastic scattering
714  G4KaonZeroSInelasticProcess* theInelasticProcess = new G4KaonZeroSInelasticProcess("inelastic");
715  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
716  theInelasticProcess->RegisterMe( theFTFModel1 );
717  theInelasticProcess->RegisterMe( theBERTModel0 );
718  pmanager->AddDiscreteProcess( theInelasticProcess );
719  }
720 
721  else if (particleName == "kaon0L")
722  {
723  // Elastic scattering
724  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
725  theElasticProcess->RegisterMe( elastic_lhep0 );
726  pmanager->AddDiscreteProcess( theElasticProcess );
727  // Inelastic scattering
728  G4KaonZeroLInelasticProcess* theInelasticProcess = new G4KaonZeroLInelasticProcess("inelastic");
729  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
730  theInelasticProcess->RegisterMe( theFTFModel1 );
731  theInelasticProcess->RegisterMe( theBERTModel0 );
732  pmanager->AddDiscreteProcess( theInelasticProcess );
733  }
734 
735  else if (particleName == "kaon-")
736  {
737  // Elastic scattering
738  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
739  theElasticProcess->RegisterMe( elastic_lhep0 );
740  pmanager->AddDiscreteProcess( theElasticProcess );
741  // Inelastic scattering
742  G4KaonMinusInelasticProcess* theInelasticProcess = new G4KaonMinusInelasticProcess("inelastic");
743  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonMinusInelasticXS::Default_Name()));
744  theInelasticProcess->RegisterMe( theFTFModel1 );
745  theInelasticProcess->RegisterMe( theBERTModel0 );
746  pmanager->AddDiscreteProcess( theInelasticProcess );
747  }
748 
749  else if (particleName == "proton")
750  {
751  // Elastic scattering
752  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
753  theElasticProcess->AddDataSet(G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsProtonElasticXS::Default_Name()));
754  theElasticProcess->RegisterMe( elastic_chip );
755  pmanager->AddDiscreteProcess( theElasticProcess );
756  // Inelastic scattering
757  G4ProtonInelasticProcess* theInelasticProcess = new G4ProtonInelasticProcess("inelastic");
758  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Proton::Proton() ) );
759  theInelasticProcess->RegisterMe( theFTFModel1 );
760  theInelasticProcess->RegisterMe( theBERTModel0 );
761  pmanager->AddDiscreteProcess( theInelasticProcess );
762  }
763 
764  else if (particleName == "anti_proton")
765  {
766  // Elastic scattering
767  const G4double elastic_elimitAntiNuc = 100.0*CLHEP::MeV;
768  G4AntiNuclElastic* elastic_anuc = new G4AntiNuclElastic();
769  elastic_anuc->SetMinEnergy( elastic_elimitAntiNuc );
770  G4CrossSectionElastic* elastic_anucxs = new G4CrossSectionElastic( elastic_anuc->GetComponentCrossSection() );
771  G4HadronElastic* elastic_lhep2 = new G4HadronElastic();
772  elastic_lhep2->SetMaxEnergy( elastic_elimitAntiNuc );
773  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
774  theElasticProcess->AddDataSet( elastic_anucxs );
775  theElasticProcess->RegisterMe( elastic_lhep2 );
776  theElasticProcess->RegisterMe( elastic_anuc );
777  pmanager->AddDiscreteProcess( theElasticProcess );
778  // Inelastic scattering
779  G4AntiProtonInelasticProcess* theInelasticProcess = new G4AntiProtonInelasticProcess("inelastic");
780  theInelasticProcess->AddDataSet( theAntiNucleonData );
781  theInelasticProcess->RegisterMe( theFTFModel0 );
782  pmanager->AddDiscreteProcess( theInelasticProcess );
783  }
784 
785  else if (particleName == "neutron") {
786  // elastic scattering
787  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
788  theElasticProcess->AddDataSet(G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsNeutronElasticXS::Default_Name()));
789  G4HadronElastic* elastic_neutronChipsModel = new G4ChipsElasticModel();
790  elastic_neutronChipsModel->SetMinEnergy( 19.0*CLHEP::MeV );
791  theElasticProcess->RegisterMe( elastic_neutronChipsModel );
792  G4ParticleHPElastic * theElasticNeutronHP = new G4ParticleHPElastic;
793  theElasticNeutronHP->SetMinEnergy( theHPMin );
794  theElasticNeutronHP->SetMaxEnergy( theHPMax );
795  theElasticProcess->RegisterMe( theElasticNeutronHP );
796  theElasticProcess->AddDataSet( new G4ParticleHPElasticData );
797  pmanager->AddDiscreteProcess( theElasticProcess );
798  // inelastic scattering
799  G4NeutronInelasticProcess* theInelasticProcess = new G4NeutronInelasticProcess("inelastic");
800  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Neutron::Neutron() ) );
801  theInelasticProcess->RegisterMe( theFTFModel1 );
802  theInelasticProcess->RegisterMe( theBERTModel1 );
803  G4ParticleHPInelastic * theNeutronInelasticHPModel = new G4ParticleHPInelastic;
804  theNeutronInelasticHPModel->SetMinEnergy( theHPMin );
805  theNeutronInelasticHPModel->SetMaxEnergy( theHPMax );
806  theInelasticProcess->RegisterMe( theNeutronInelasticHPModel );
807  theInelasticProcess->AddDataSet( new G4ParticleHPInelasticData );
808  pmanager->AddDiscreteProcess(theInelasticProcess);
809  // capture
810  G4HadronCaptureProcess* theCaptureProcess = new G4HadronCaptureProcess;
811  G4ParticleHPCapture * theNeutronCaptureHPModel = new G4ParticleHPCapture;
812  theNeutronCaptureHPModel->SetMinEnergy( theHPMin );
813  theNeutronCaptureHPModel->SetMaxEnergy( theHPMax );
814  G4NeutronRadCapture* theNeutronRadCapture = new G4NeutronRadCapture();
815  theNeutronRadCapture->SetMinEnergy(theHPMax*0.99);
816  theCaptureProcess->RegisterMe( theNeutronCaptureHPModel );
817  theCaptureProcess->RegisterMe( theNeutronRadCapture);
818  theCaptureProcess->AddDataSet( new G4ParticleHPCaptureData );
819  theCaptureProcess->AddDataSet((G4NeutronCaptureXS*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4NeutronCaptureXS::Default_Name()));
820  pmanager->AddDiscreteProcess(theCaptureProcess);
821  }
822  else if (particleName == "anti_neutron")
823  {
824  // Elastic scattering
825  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
826  theElasticProcess->RegisterMe( elastic_lhep0 );
827  pmanager->AddDiscreteProcess( theElasticProcess );
828  // Inelastic scattering
829  G4AntiNeutronInelasticProcess* theInelasticProcess = new G4AntiNeutronInelasticProcess("inelastic");
830  theInelasticProcess->AddDataSet( theAntiNucleonData );
831  theInelasticProcess->RegisterMe( theFTFModel0 );
832  pmanager->AddDiscreteProcess( theInelasticProcess );
833  }
834 
835  else if (particleName == "deuteron")
836  {
837  // Elastic scattering
838  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
839  theElasticProcess->RegisterMe( elastic_lhep0 );
840  pmanager->AddDiscreteProcess( theElasticProcess );
841  // Inelastic scattering
842  G4DeuteronInelasticProcess* theInelasticProcess = new G4DeuteronInelasticProcess("inelastic");
843  theInelasticProcess->AddDataSet( theGGNuclNuclData );
844  theInelasticProcess->RegisterMe( theFTFModel1 );
845  theInelasticProcess->RegisterMe( theIonBC );
846  pmanager->AddDiscreteProcess( theInelasticProcess );
847  }
848 
849  else if (particleName == "triton")
850  {
851  // Elastic scattering
852  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
853  theElasticProcess->RegisterMe( elastic_lhep0 );
854  pmanager->AddDiscreteProcess( theElasticProcess );
855  // Inelastic scattering
856  G4TritonInelasticProcess* theInelasticProcess = new G4TritonInelasticProcess("inelastic");
857  theInelasticProcess->AddDataSet( theGGNuclNuclData );
858  theInelasticProcess->RegisterMe( theFTFModel1 );
859  theInelasticProcess->RegisterMe( theIonBC );
860  pmanager->AddDiscreteProcess( theInelasticProcess );
861  }
862 
863  else if (particleName == "alpha")
864  {
865  // Elastic scattering
866  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
867  theElasticProcess->RegisterMe( elastic_lhep0 );
868  pmanager->AddDiscreteProcess( theElasticProcess );
869  // Inelastic scattering
870  G4AlphaInelasticProcess* theInelasticProcess = new G4AlphaInelasticProcess("inelastic");
871  theInelasticProcess->AddDataSet( theGGNuclNuclData );
872  theInelasticProcess->RegisterMe( theFTFModel1 );
873  theInelasticProcess->RegisterMe( theIonBC );
874  pmanager->AddDiscreteProcess( theInelasticProcess );
875  }
876  } // while ((*(myParticleIterator))())
877 
878  // Add stopping processes with builder
879  stoppingPhysics->ConstructProcess();
880 }
881 
882 
883 // Decays ///////////////////////////////////////////////////////////////////
884 #include "G4Decay.hh"
885 #include "G4RadioactiveDecay.hh"
886 #include "G4IonTable.hh"
887 #include "G4Ions.hh"
888 
889  template<class T> void TLBE<T>::ConstructGeneral() {
890 
891  // Add Decay Process
892  G4Decay* theDecayProcess = new G4Decay();
893  G4bool theDecayProcessNeverUsed = true; //Check if theDecayProcess will be used
894  auto myParticleIterator=G4ParticleTable::GetParticleTable()->GetIterator();
895  myParticleIterator->reset();
896  while( (*(myParticleIterator))() )
897  {
898  G4ParticleDefinition* particle = myParticleIterator->value();
899  G4ProcessManager* pmanager = particle->GetProcessManager();
900 
901  if (theDecayProcess->IsApplicable(*particle) && !particle->IsShortLived())
902  {
903  theDecayProcessNeverUsed = false;
904  pmanager ->AddProcess(theDecayProcess);
905  // set ordering for PostStepDoIt and AtRestDoIt
906  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
907  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
908  }
909  }
910 
911  // Declare radioactive decay to the GenericIon in the IonTable.
912  const G4IonTable *theIonTable =
913  G4ParticleTable::GetParticleTable()->GetIonTable();
914  G4RadioactiveDecay *theRadioactiveDecay = new G4RadioactiveDecay();
915 
916  for (G4int i=0; i<theIonTable->Entries(); i++)
917  {
918  G4String particleName = theIonTable->GetParticle(i)->GetParticleName();
919  G4String particleType = theIonTable->GetParticle(i)->GetParticleType();
920 
921  if (particleName == "GenericIon")
922  {
923  G4ProcessManager* pmanager =
924  theIonTable->GetParticle(i)->GetProcessManager();
925  pmanager->SetVerboseLevel(VerboseLevel);
926  pmanager ->AddProcess(theRadioactiveDecay);
927  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxPostStep);
928  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxAtRest);
929  }
930  }
931  //If we actually never used the process, delete it
932  //From Coverity report
933  if ( theDecayProcessNeverUsed ) delete theDecayProcess;
934 }
935 
936 // Cuts /////////////////////////////////////////////////////////////////////
937 template<class T> void TLBE<T>::SetCuts()
938 {
939 
940  if (this->verboseLevel >1)
941  G4cout << "LBE::SetCuts:";
942 
943  if (this->verboseLevel>0){
944  G4cout << "LBE::SetCuts:";
945  G4cout << "CutLength : "
946  << G4BestUnit(this->defaultCutValue,"Length") << G4endl;
947  }
948 
949  //special for low energy physics
950  G4double lowlimit=250*CLHEP::eV;
951  G4ProductionCutsTable * aPCTable = G4ProductionCutsTable::GetProductionCutsTable();
952  aPCTable->SetEnergyRange(lowlimit,100*CLHEP::GeV);
953 
954  // set cut values for gamma at first and for e- second and next for e+,
955  // because some processes for e+/e- need cut values for gamma
956  this->SetCutValue(cutForGamma, "gamma");
957  this->SetCutValue(cutForElectron, "e-");
958  this->SetCutValue(cutForPositron, "e+");
959 
960  // this->SetCutValue(cutForProton, "proton");
961  // this->SetCutValue(cutForProton, "anti_proton");
962  // this->SetCutValue(cutForAlpha, "alpha");
963  // this->SetCutValue(cutForGenericIon, "GenericIon");
964 
965  // this->SetCutValueForOthers(this->defaultCutValue);
966 
967  if (this->verboseLevel>0) this->DumpCutValuesTable();
968 }
969