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