Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
LBE Class Reference

#include <LBE.hh>

+ Inheritance diagram for LBE:

Public Member Functions

 LBE (G4int ver=1)
 
virtual ~LBE ()
 
 LBE (const LBE &)=delete
 
LBEoperator= (const LBE &right)=delete
 
virtual void SetCuts ()
 
- Public Member Functions inherited from G4VModularPhysicsList
 G4VModularPhysicsList ()
 
virtual ~G4VModularPhysicsList ()
 
virtual void ConstructParticle () override
 
virtual void ConstructProcess () override
 
void RegisterPhysics (G4VPhysicsConstructor *)
 
const G4VPhysicsConstructorGetPhysics (G4int index) const
 
const G4VPhysicsConstructorGetPhysics (const G4String &name) const
 
const G4VPhysicsConstructorGetPhysicsWithType (G4int physics_type) const
 
void ReplacePhysics (G4VPhysicsConstructor *)
 
void RemovePhysics (G4VPhysicsConstructor *)
 
void RemovePhysics (G4int type)
 
void RemovePhysics (const G4String &name)
 
G4int GetInstanceID () const
 
virtual void TerminateWorker () override
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
- Public Member Functions inherited from G4VUserPhysicsList
 G4VUserPhysicsList ()
 
virtual ~G4VUserPhysicsList ()
 
 G4VUserPhysicsList (const G4VUserPhysicsList &)
 
G4VUserPhysicsListoperator= (const G4VUserPhysicsList &)
 
virtual void ConstructParticle ()=0
 
void Construct ()
 
virtual void ConstructProcess ()=0
 
virtual void SetCuts ()
 
void SetDefaultCutValue (G4double newCutValue)
 
G4double GetDefaultCutValue () const
 
void BuildPhysicsTable ()
 
void PreparePhysicsTable (G4ParticleDefinition *)
 
void BuildPhysicsTable (G4ParticleDefinition *)
 
G4bool StorePhysicsTable (const G4String &directory=".")
 
G4bool IsPhysicsTableRetrieved () const
 
G4bool IsStoredInAscii () const
 
const G4StringGetPhysicsTableDirectory () const
 
void SetPhysicsTableRetrieved (const G4String &directory="")
 
void SetStoredInAscii ()
 
void ResetPhysicsTableRetrieved ()
 
void ResetStoredInAscii ()
 
void DumpList () const
 
void DumpCutValuesTable (G4int flag=1)
 
void DumpCutValuesTableIfRequested ()
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
void UseCoupledTransportation (G4bool vl=true)
 
void SetCutsWithDefault ()
 
void SetCutValue (G4double aCut, const G4String &pname)
 
G4double GetCutValue (const G4String &pname) const
 
void SetCutValue (G4double aCut, const G4String &pname, const G4String &rname)
 
void SetParticleCuts (G4double cut, G4ParticleDefinition *particle, G4Region *region=nullptr)
 
void SetParticleCuts (G4double cut, const G4String &particleName, G4Region *region=nullptr)
 
void SetCutsForRegion (G4double aCut, const G4String &rname)
 
void SetApplyCuts (G4bool value, const G4String &name)
 
G4bool GetApplyCuts (const G4String &name) const
 
void RemoveProcessManager ()
 
void RemoveTrackingManager ()
 
void AddProcessManager (G4ParticleDefinition *newParticle, G4ProcessManager *newManager=nullptr)
 
void CheckParticleList ()
 
void DisableCheckParticleList ()
 
G4int GetInstanceID () const
 
virtual void InitializeWorker ()
 
virtual void TerminateWorker ()
 

Protected Member Functions

virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
virtual void ConstructGeneral ()
 
virtual void ConstructEM ()
 
virtual void ConstructHad ()
 
virtual void ConstructOp ()
 
virtual void AddTransportation ()
 
- Protected Member Functions inherited from G4VModularPhysicsList
 G4VModularPhysicsList (const G4VModularPhysicsList &)
 
G4VModularPhysicsListoperator= (const G4VModularPhysicsList &)
 
- Protected Member Functions inherited from G4VUserPhysicsList
void AddTransportation ()
 
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
void BuildIntegralPhysicsTable (G4VProcess *, G4ParticleDefinition *)
 
virtual void RetrievePhysicsTable (G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
void InitializeProcessManager ()
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VModularPhysicsList
static const G4VMPLManagerGetSubInstanceManager ()
 
- Static Public Member Functions inherited from G4VUserPhysicsList
static const G4VUPLManagerGetSubInstanceManager ()
 
- Protected Types inherited from G4VModularPhysicsList
using G4PhysConstVector = G4VMPLData::G4PhysConstVectorData
 
- Protected Attributes inherited from G4VModularPhysicsList
G4int verboseLevel = 0
 
G4int g4vmplInstanceID = 0
 
- Protected Attributes inherited from G4VUserPhysicsList
G4ParticleTabletheParticleTable = nullptr
 
G4int verboseLevel = 1
 
G4double defaultCutValue = 1.0
 
G4bool isSetDefaultCutValue = false
 
G4ProductionCutsTablefCutsTable = nullptr
 
G4bool fRetrievePhysicsTable = false
 
G4bool fStoredInAscii = true
 
G4bool fIsCheckedForRetrievePhysicsTable = false
 
G4bool fIsRestoredCutValues = false
 
G4String directoryPhysicsTable = "."
 
G4bool fDisableCheckParticleList = false
 
G4int g4vuplInstanceID = 0
 
- Static Protected Attributes inherited from G4VModularPhysicsList
static G4RUN_DLL G4VMPLManager G4VMPLsubInstanceManager
 
- Static Protected Attributes inherited from G4VUserPhysicsList
static G4RUN_DLL G4VUPLManager subInstanceManager
 

Detailed Description

Definition at line 58 of file LBE.hh.

Constructor & Destructor Documentation

◆ LBE() [1/2]

LBE::LBE ( G4int  ver = 1)

Definition at line 77 of file LBE.cc.

78{
79 if(ver > 0) {
80 G4cout << "You are using the simulation engine: LBE"<<G4endl;
81 G4cout <<G4endl;
82 }
83 defaultCutValue = 1.0*CLHEP::micrometer; //
84 cutForGamma = defaultCutValue;
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}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetVerboseLevel(G4int value)

◆ ~LBE()

LBE::~LBE ( )
virtual

Definition at line 102 of file LBE.cc.

103{
104 delete stoppingPhysics;
105}

◆ LBE() [2/2]

LBE::LBE ( const LBE )
delete

Member Function Documentation

◆ AddTransportation()

void LBE::AddTransportation ( )
protectedvirtual

Definition at line 217 of file LBE.cc.

217 {
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}
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
void reset(G4bool ifSkipIon=true)
G4PTblDicIterator * GetIterator() const
static G4ParticleTable * GetParticleTable()
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)

Referenced by ConstructProcess().

◆ ConstructEM()

void LBE::ConstructEM ( )
protectedvirtual

Definition at line 287 of file LBE.cc.

287 {
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}
@ fUseDistanceToBoundary
double G4double
Definition: G4Types.hh:83
const G4String & GetParticleType() const
G4double GetPDGCharge() const
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetStepFunction(G4double v1, G4double v2)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetStepLimitType(G4MscStepLimitType val)

Referenced by ConstructProcess().

◆ ConstructGeneral()

void LBE::ConstructGeneral ( )
protectedvirtual

Definition at line 886 of file LBE.cc.

886 {
887
888 // Add Decay Process
889 G4Decay* theDecayProcess = new G4Decay();
890 G4bool theDecayProcessNeverUsed = true; //Check if theDecayProcess will be used
891 auto myParticleIterator=G4ParticleTable::GetParticleTable()->GetIterator();
892 myParticleIterator->reset();
893 while( (*(myParticleIterator))() )
894 {
895 G4ParticleDefinition* particle = myParticleIterator->value();
896 G4ProcessManager* pmanager = particle->GetProcessManager();
897
898 if (theDecayProcess->IsApplicable(*particle) && !particle->IsShortLived())
899 {
900 theDecayProcessNeverUsed = false;
901 pmanager ->AddProcess(theDecayProcess);
902 // set ordering for PostStepDoIt and AtRestDoIt
903 pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
904 pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
905 }
906 }
907
908 // Declare radioactive decay to the GenericIon in the IonTable.
909 const G4IonTable *theIonTable =
911 G4RadioactiveDecay* theRadioactiveDecay = new G4RadioactiveDecay();
912
913 //Fix for activation of RadioactiveDecay, based on G4RadioactiveDecayPhysics
915 param->SetAugerCascade(true);
916 param->AddPhysics("world","G4RadioactiveDecay");
917
919 deex->SetStoreAllLevels(true);
920 deex->SetMaxLifeTime(G4NuclideTable::GetInstance()->GetThresholdOfHalfLife()
921 /std::log(2.));
922
925 if(!ad) {
926 ad = new G4UAtomicDeexcitation();
927 man->SetAtomDeexcitation(ad);
929 }
930
931 for (G4int i=0; i<theIonTable->Entries(); i++)
932 {
933 G4String particleName = theIonTable->GetParticle(i)->GetParticleName();
934 G4String particleType = theIonTable->GetParticle(i)->GetParticleType();
935
936 if (particleName == "GenericIon")
937 {
938 G4ProcessManager* pmanager =
939 theIonTable->GetParticle(i)->GetProcessManager();
940 pmanager->SetVerboseLevel(VerboseLevel);
941 pmanager ->AddProcess(theRadioactiveDecay);
942 pmanager ->SetProcessOrdering(theRadioactiveDecay, idxPostStep);
943 pmanager ->SetProcessOrdering(theRadioactiveDecay, idxAtRest);
944 }
945 }
946 //If we actually never used the process, delete it
947 //From Coverity report
948 if ( theDecayProcessNeverUsed ) delete theDecayProcess;
949}
@ idxPostStep
@ idxAtRest
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
Definition: G4Decay.cc:88
static G4EmParameters * Instance()
void SetAugerCascade(G4bool val)
void AddPhysics(const G4String &region, const G4String &type)
G4ParticleDefinition * GetParticle(G4int index) const
Definition: G4IonTable.cc:1905
G4int Entries() const
Definition: G4IonTable.cc:1962
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
static G4NuclideTable * GetInstance()
G4IonTable * GetIonTable() const
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
void SetVerboseLevel(G4int value)

Referenced by ConstructProcess().

◆ ConstructHad()

void LBE::ConstructHad ( )
protectedvirtual

Definition at line 579 of file LBE.cc.

580{
581 // Elastic scattering
582 G4HadronElastic* elastic_lhep0 = new G4HadronElastic();
583 G4ChipsElasticModel* elastic_chip = new G4ChipsElasticModel();
585
586 const G4double elastic_elimitAntiNuc = 100.0*CLHEP::MeV;
587 G4AntiNuclElastic* elastic_anuc = new G4AntiNuclElastic();
588 elastic_anuc->SetMinEnergy( elastic_elimitAntiNuc );
589 G4CrossSectionElastic* elastic_anucxs = new G4CrossSectionElastic( elastic_anuc->GetComponentCrossSection() );
590 G4HadronElastic* elastic_lhep2 = new G4HadronElastic();
591 elastic_lhep2->SetMaxEnergy( elastic_elimitAntiNuc );
592
593 // Inelastic scattering
594 const G4double theFTFMin0 = 0.0*CLHEP::GeV;
595 const G4double theFTFMin1 = 4.0*CLHEP::GeV;
597 const G4double theBERTMin0 = 0.0*CLHEP::GeV;
598 const G4double theBERTMin1 = 19.0*CLHEP::MeV;
599 const G4double theBERTMax = 5.0*CLHEP::GeV;
600 const G4double theHPMin = 0.0*CLHEP::GeV;
601 const G4double theHPMax = 20.0*CLHEP::MeV;
602 const G4double theIonBCMin = 0.0*CLHEP::GeV;
603 const G4double theIonBCMax = 5.0*CLHEP::GeV;
604
605
606 G4FTFModel * theStringModel = new G4FTFModel;
608 theStringModel->SetFragmentationModel( theStringDecay );
609 G4PreCompoundModel * thePreEquilib = new G4PreCompoundModel( new G4ExcitationHandler );
610 G4GeneratorPrecompoundInterface * theCascade = new G4GeneratorPrecompoundInterface( thePreEquilib );
611
612 G4TheoFSGenerator * theFTFModel0 = new G4TheoFSGenerator( "FTFP" );
613 theFTFModel0->SetHighEnergyGenerator( theStringModel );
614 theFTFModel0->SetTransport( theCascade );
615 theFTFModel0->SetMinEnergy( theFTFMin0 );
616 theFTFModel0->SetMaxEnergy( theFTFMax );
617
618 G4TheoFSGenerator * theFTFModel1 = new G4TheoFSGenerator( "FTFP" );
619 theFTFModel1->SetHighEnergyGenerator( theStringModel );
620 theFTFModel1->SetTransport( theCascade );
621 theFTFModel1->SetMinEnergy( theFTFMin1 );
622 theFTFModel1->SetMaxEnergy( theFTFMax );
623
624 G4CascadeInterface * theBERTModel0 = new G4CascadeInterface;
625 theBERTModel0->SetMinEnergy( theBERTMin0 );
626 theBERTModel0->SetMaxEnergy( theBERTMax );
627
628 G4CascadeInterface * theBERTModel1 = new G4CascadeInterface;
629 theBERTModel1->SetMinEnergy( theBERTMin1 );
630 theBERTModel1->SetMaxEnergy( theBERTMax );
631
632 // Binary Cascade
633 G4BinaryLightIonReaction * theIonBC = new G4BinaryLightIonReaction( thePreEquilib );
634 theIonBC->SetMinEnergy( theIonBCMin );
635 theIonBC->SetMaxEnergy( theIonBCMax );
636
638 G4ComponentGGNuclNuclXsc * ggNuclNuclXsec = new G4ComponentGGNuclNuclXsc();
639 G4VCrossSectionDataSet * theGGNuclNuclData = new G4CrossSectionInelastic(ggNuclNuclXsec);
640
641 auto myParticleIterator=G4ParticleTable::GetParticleTable()->GetIterator();
642 myParticleIterator->reset();
643 while ((*(myParticleIterator))())
644 {
645 G4ParticleDefinition* particle = myParticleIterator->value();
646 G4ProcessManager* pmanager = particle->GetProcessManager();
647 G4String particleName = particle->GetParticleName();
648
649 if (particleName == "pi+")
650 {
651 // Elastic scattering
652 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
653 theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
654 theElasticProcess->RegisterMe( elastic_he );
655 pmanager->AddDiscreteProcess( theElasticProcess );
656 // Inelastic scattering
657 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4PionPlus::Definition() );
658 theInelasticProcess->AddDataSet( new G4BGGPionInelasticXS( G4PionPlus::Definition() ) );
659 theInelasticProcess->RegisterMe( theFTFModel1 );
660 theInelasticProcess->RegisterMe( theBERTModel0 );
661 pmanager->AddDiscreteProcess( theInelasticProcess );
662 }
663
664 else if (particleName == "pi-")
665 {
666 // Elastic scattering
667 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
668 theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
669 theElasticProcess->RegisterMe( elastic_he );
670 pmanager->AddDiscreteProcess( theElasticProcess );
671 // Inelastic scattering
672 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4PionMinus::Definition() );
673 theInelasticProcess->AddDataSet( new G4BGGPionInelasticXS( G4PionMinus::Definition() ) );
674 theInelasticProcess->RegisterMe( theFTFModel1 );
675 theInelasticProcess->RegisterMe( theBERTModel0 );
676 pmanager->AddDiscreteProcess( theInelasticProcess );
677 }
678
679 else if (particleName == "kaon+")
680 {
681 // Elastic scattering
682 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
683 theElasticProcess->RegisterMe( elastic_lhep0 );
684 theElasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonPlusElasticXS::Default_Name()));
685 pmanager->AddDiscreteProcess( theElasticProcess );
686 // Inelastic scattering
687 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4KaonPlus::Definition() );
688
689 theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonPlusInelasticXS::Default_Name()));
690 theInelasticProcess->RegisterMe( theFTFModel1 );
691 theInelasticProcess->RegisterMe( theBERTModel0 );
692 pmanager->AddDiscreteProcess( theInelasticProcess );
693 }
694
695 else if (particleName == "kaon0S")
696 {
697 // Elastic scattering
698 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
699 theElasticProcess->RegisterMe( elastic_lhep0 );
700 theElasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonZeroElasticXS::Default_Name()));
701 pmanager->AddDiscreteProcess( theElasticProcess );
702 // Inelastic scattering
703 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4KaonZeroShort::Definition() );
704 theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
705 theInelasticProcess->RegisterMe( theFTFModel1 );
706 theInelasticProcess->RegisterMe( theBERTModel0 );
707 pmanager->AddDiscreteProcess( theInelasticProcess );
708 }
709
710 else if (particleName == "kaon0L")
711 {
712 // Elastic scattering
713 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
714 theElasticProcess->RegisterMe( elastic_lhep0 );
715 theElasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonZeroElasticXS::Default_Name()));
716 pmanager->AddDiscreteProcess( theElasticProcess );
717 // Inelastic scattering
718 //G4KaonZeroLInelasticProcess* theInelasticProcess = new G4KaonZeroLInelasticProcess("inelastic");
719 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4KaonZeroLong::Definition() );
720 theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
721 theInelasticProcess->RegisterMe( theFTFModel1 );
722 theInelasticProcess->RegisterMe( theBERTModel0 );
723 pmanager->AddDiscreteProcess( theInelasticProcess );
724 }
725
726 else if (particleName == "kaon-")
727 {
728 // Elastic scattering
729 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
730 theElasticProcess->RegisterMe( elastic_lhep0 );
732 pmanager->AddDiscreteProcess( theElasticProcess );
733 // Inelastic scattering
734 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4KaonMinus::Definition() );
735 theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsKaonMinusInelasticXS::Default_Name()));
736 theInelasticProcess->RegisterMe( theFTFModel1 );
737 theInelasticProcess->RegisterMe( theBERTModel0 );
738 pmanager->AddDiscreteProcess( theInelasticProcess );
739 }
740
741 else if (particleName == "proton")
742 {
743 // Elastic scattering
744 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
745 theElasticProcess->AddDataSet(G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsProtonElasticXS::Default_Name()));
746 theElasticProcess->AddDataSet( new G4BGGNucleonElasticXS( G4Proton::Proton() ) );
747 theElasticProcess->RegisterMe( elastic_chip );
748 pmanager->AddDiscreteProcess( theElasticProcess );
749 // Inelastic scattering
750 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4Proton::Definition() );
751 theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Proton::Proton() ) );
752 theInelasticProcess->RegisterMe( theFTFModel1 );
753 theInelasticProcess->RegisterMe( theBERTModel0 );
754 pmanager->AddDiscreteProcess( theInelasticProcess );
755 }
756
757 else if (particleName == "anti_proton")
758 {
759 // Elastic scattering
760 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
761 theElasticProcess->AddDataSet( elastic_anucxs );
762 theElasticProcess->RegisterMe( elastic_lhep2 );
763 theElasticProcess->RegisterMe( elastic_anuc );
764 pmanager->AddDiscreteProcess( theElasticProcess );
765 // Inelastic scattering
766 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4AntiProton::Definition() );
767 theInelasticProcess->AddDataSet( theAntiNucleonData );
768 theInelasticProcess->RegisterMe( theFTFModel0 );
769 pmanager->AddDiscreteProcess( theInelasticProcess );
770 }
771
772 else if (particleName == "neutron") {
773 // elastic scattering
774 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
776 G4HadronElastic* elastic_neutronChipsModel = new G4ChipsElasticModel();
777 elastic_neutronChipsModel->SetMinEnergy( 19.0*CLHEP::MeV );
778 theElasticProcess->RegisterMe( elastic_neutronChipsModel );
779 G4ParticleHPElastic * theElasticNeutronHP = new G4ParticleHPElastic;
780 theElasticNeutronHP->SetMinEnergy( theHPMin );
781 theElasticNeutronHP->SetMaxEnergy( theHPMax );
782 theElasticProcess->RegisterMe( theElasticNeutronHP );
783 theElasticProcess->AddDataSet( new G4ParticleHPElasticData );
784 pmanager->AddDiscreteProcess( theElasticProcess );
785 // inelastic scattering
786 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4Neutron::Definition() );
787 theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Neutron::Neutron() ) );
788 theInelasticProcess->RegisterMe( theFTFModel1 );
789 theInelasticProcess->RegisterMe( theBERTModel1 );
790 G4ParticleHPInelastic * theNeutronInelasticHPModel = new G4ParticleHPInelastic;
791 theNeutronInelasticHPModel->SetMinEnergy( theHPMin );
792 theNeutronInelasticHPModel->SetMaxEnergy( theHPMax );
793 theInelasticProcess->RegisterMe( theNeutronInelasticHPModel );
794 theInelasticProcess->AddDataSet( new G4ParticleHPInelasticData );
795 pmanager->AddDiscreteProcess(theInelasticProcess);
796 // capture
797 G4NeutronCaptureProcess* theCaptureProcess = new G4NeutronCaptureProcess;
798 G4ParticleHPCapture * theNeutronCaptureHPModel = new G4ParticleHPCapture;
799 theNeutronCaptureHPModel->SetMinEnergy( theHPMin );
800 theNeutronCaptureHPModel->SetMaxEnergy( theHPMax );
801 G4NeutronRadCapture* theNeutronRadCapture = new G4NeutronRadCapture();
802 theNeutronRadCapture->SetMinEnergy(theHPMax*0.99);
803 theCaptureProcess->RegisterMe( theNeutronCaptureHPModel );
804 theCaptureProcess->RegisterMe( theNeutronRadCapture);
805 theCaptureProcess->AddDataSet( new G4ParticleHPCaptureData );
807 pmanager->AddDiscreteProcess(theCaptureProcess);
808 }
809 else if (particleName == "anti_neutron")
810 {
811 // Elastic scattering
812 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
813 theElasticProcess->AddDataSet( elastic_anucxs );
814 theElasticProcess->RegisterMe( elastic_lhep2 );
815 theElasticProcess->RegisterMe( elastic_anuc );
816 pmanager->AddDiscreteProcess( theElasticProcess );
817 // Inelastic scattering
818 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4AntiNeutron::Definition() );
819 theInelasticProcess->AddDataSet( theAntiNucleonData );
820 theInelasticProcess->RegisterMe( theFTFModel0 );
821 pmanager->AddDiscreteProcess( theInelasticProcess );
822 }
823
824 else if (particleName == "deuteron")
825 {
826 // Elastic scattering
827 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
828 theElasticProcess->RegisterMe( elastic_lhep0 );
829 theElasticProcess->AddDataSet( theGGNuclNuclData );
830 pmanager->AddDiscreteProcess( theElasticProcess );
831 // Inelastic scattering
832 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4Deuteron::Definition() );
833 theInelasticProcess->AddDataSet( theGGNuclNuclData );
834 theInelasticProcess->RegisterMe( theFTFModel1 );
835 theInelasticProcess->RegisterMe( theIonBC );
836 pmanager->AddDiscreteProcess( theInelasticProcess );
837 }
838
839 else if (particleName == "triton")
840 {
841 // Elastic scattering
842 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
843 theElasticProcess->RegisterMe( elastic_lhep0 );
844 theElasticProcess->AddDataSet( theGGNuclNuclData );
845 pmanager->AddDiscreteProcess( theElasticProcess );
846 // Inelastic scattering
847 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4Triton::Definition() );
848 theInelasticProcess->AddDataSet( theGGNuclNuclData );
849 theInelasticProcess->RegisterMe( theFTFModel1 );
850 theInelasticProcess->RegisterMe( theIonBC );
851 pmanager->AddDiscreteProcess( theInelasticProcess );
852 }
853
854 else if (particleName == "alpha")
855 {
856 // Elastic scattering
857 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
858 theElasticProcess->RegisterMe( elastic_lhep0 );
859 theElasticProcess->AddDataSet( theGGNuclNuclData );
860 pmanager->AddDiscreteProcess( theElasticProcess );
861 // Inelastic scattering
862 G4HadronInelasticProcess* theInelasticProcess = new G4HadronInelasticProcess( "inelastic", G4Alpha::Definition() );
863 theInelasticProcess->AddDataSet( theGGNuclNuclData );
864 theInelasticProcess->RegisterMe( theFTFModel1 );
865 theInelasticProcess->RegisterMe( theIonBC );
866 pmanager->AddDiscreteProcess( theInelasticProcess );
867 }
868 } // while ((*(myParticleIterator))())
869
870 // Add stopping processes with builder
871 stoppingPhysics->ConstructProcess();
872}
static G4Alpha * Definition()
Definition: G4Alpha.cc:48
static G4AntiNeutron * Definition()
G4ComponentAntiNuclNuclearXS * GetComponentCrossSection()
static G4AntiProton * Definition()
Definition: G4AntiProton.cc:50
static const char * Default_Name()
static const char * Default_Name()
static const char * Default_Name()
static const char * Default_Name()
static G4CrossSectionDataSetRegistry * Instance()
static G4Deuteron * Definition()
Definition: G4Deuteron.cc:49
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
G4double GetMaxEnergy() const
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void RegisterMe(G4HadronicInteraction *a)
static G4KaonMinus * Definition()
Definition: G4KaonMinus.cc:53
static G4KaonPlus * Definition()
Definition: G4KaonPlus.cc:53
static G4KaonZeroLong * Definition()
static G4KaonZeroShort * Definition()
static const char * Default_Name()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4Neutron * Definition()
Definition: G4Neutron.cc:53
static G4PionMinus * Definition()
Definition: G4PionMinus.cc:51
static G4PionPlus * Definition()
Definition: G4PionPlus.cc:51
static G4Proton * Definition()
Definition: G4Proton.cc:48
static G4Proton * Proton()
Definition: G4Proton.cc:92
virtual void ConstructProcess() override
void SetTransport(G4VIntraNuclearTransportModel *const value)
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
static G4Triton * Definition()
Definition: G4Triton.cc:49
void SetFragmentationModel(G4VStringFragmentation *aModel)

Referenced by ConstructProcess().

◆ ConstructOp()

void LBE::ConstructOp ( )
protectedvirtual

Definition at line 438 of file LBE.cc.

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->SetVerboseLevel(OpVerbLevel);
447
448 // scintillation process for alpha:
449 G4bool theScintProcessAlphaNeverUsed = true;
450 G4Scintillation* theScintProcessAlpha = new G4Scintillation("Scintillation");
451 // theScintProcessNuc->DumpPhysicsTable();
452 theScintProcessAlpha->SetTrackSecondariesFirst(true);
453 theScintProcessAlpha->SetVerboseLevel(OpVerbLevel);
454
455 // scintillation process for heavy nuclei
456 G4bool theScintProcessNucNeverUsed = true;
457 G4Scintillation* theScintProcessNuc = new G4Scintillation("Scintillation");
458 // theScintProcessNuc->DumpPhysicsTable();
459 theScintProcessNuc->SetTrackSecondariesFirst(true);
460 theScintProcessNuc->SetVerboseLevel(OpVerbLevel);
461
462 // optical processes
463 G4bool theAbsorptionProcessNeverUsed = true;
464 G4OpAbsorption* theAbsorptionProcess = new G4OpAbsorption();
465 // G4OpRayleigh* theRayleighScatteringProcess = new G4OpRayleigh();
466 G4bool theBoundaryProcessNeverUsed = true;
467 G4OpBoundaryProcess* theBoundaryProcess = new G4OpBoundaryProcess();
468 // theAbsorptionProcess->DumpPhysicsTable();
469 // theRayleighScatteringProcess->DumpPhysicsTable();
470 theAbsorptionProcess->SetVerboseLevel(OpVerbLevel);
471 // theRayleighScatteringProcess->SetVerboseLevel(OpVerbLevel);
472 theBoundaryProcess->SetVerboseLevel(OpVerbLevel);
473
474 auto myParticleIterator=G4ParticleTable::GetParticleTable()->GetIterator();
475 myParticleIterator->reset();
476 while( (*(myParticleIterator))() )
477 {
478 G4ParticleDefinition* particle = myParticleIterator->value();
479 G4ProcessManager* pmanager = particle->GetProcessManager();
480 G4String particleName = particle->GetParticleName();
481 if (theScintProcessDef->IsApplicable(*particle)) {
482 // if(particle->GetPDGMass() > 5.0*CLHEP::GeV)
483 if(particle->GetParticleName() == "GenericIon") {
484 pmanager->AddProcess(theScintProcessNuc); // AtRestDiscrete
485 pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxAtRest);
486 pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxPostStep);
487 theScintProcessNucNeverUsed = false;
488 }
489 else if(particle->GetParticleName() == "alpha") {
490 pmanager->AddProcess(theScintProcessAlpha);
491 pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxAtRest);
492 pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxPostStep);
493 theScintProcessAlphaNeverUsed = false;
494 }
495 else {
496 pmanager->AddProcess(theScintProcessDef);
497 pmanager->SetProcessOrderingToLast(theScintProcessDef,idxAtRest);
498 pmanager->SetProcessOrderingToLast(theScintProcessDef,idxPostStep);
499 theScintProcessDefNeverUsed = false;
500 }
501 }
502
503 if (particleName == "opticalphoton") {
504 pmanager->AddDiscreteProcess(theAbsorptionProcess);
505 theAbsorptionProcessNeverUsed = false;
506 // pmanager->AddDiscreteProcess(theRayleighScatteringProcess);
507 theBoundaryProcessNeverUsed = false;
508 pmanager->AddDiscreteProcess(theBoundaryProcess);
509 }
510 }
511 if ( theScintProcessDefNeverUsed ) delete theScintProcessDef;
512 if ( theScintProcessAlphaNeverUsed ) delete theScintProcessAlpha;
513 if ( theScintProcessNucNeverUsed ) delete theScintProcessNuc;
514 if ( theBoundaryProcessNeverUsed ) delete theBoundaryProcess;
515 if ( theAbsorptionProcessNeverUsed ) delete theAbsorptionProcess;
516}
void SetVerboseLevel(G4int)
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
void SetTrackSecondariesFirst(const G4bool state)
void SetVerboseLevel(G4int)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override

Referenced by ConstructProcess().

◆ ConstructParticle()

void LBE::ConstructParticle ( )
protectedvirtual

Reimplemented from G4VModularPhysicsList.

Definition at line 109 of file LBE.cc.

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}
virtual void ConstructParticle() override

◆ ConstructProcess()

void LBE::ConstructProcess ( )
protectedvirtual

Reimplemented from G4VModularPhysicsList.

Definition at line 203 of file LBE.cc.

204{
206 ConstructEM();
207 ConstructOp();
208 ConstructHad();
210}
virtual void ConstructEM()
Definition: LBE.cc:287
virtual void AddTransportation()
Definition: LBE.cc:217
virtual void ConstructOp()
Definition: LBE.cc:438
virtual void ConstructHad()
Definition: LBE.cc:579
virtual void ConstructGeneral()
Definition: LBE.cc:886

◆ operator=()

LBE & LBE::operator= ( const LBE right)
delete

◆ SetCuts()

void LBE::SetCuts ( )
virtual

Reimplemented from G4VUserPhysicsList.

Definition at line 952 of file LBE.cc.

953{
954
955 if (verboseLevel >1)
956 G4cout << "LBE::SetCuts:";
957
958 if (verboseLevel>0){
959 G4cout << "LBE::SetCuts:";
960 G4cout << "CutLength : "
961 << G4BestUnit(defaultCutValue,"Length") << G4endl;
962 }
963
964 //special for low energy physics
965 G4double lowlimit=250*CLHEP::eV;
967 aPCTable->SetEnergyRange(lowlimit,100*CLHEP::GeV);
968
969 // set cut values for gamma at first and for e- second and next for e+,
970 // because some processes for e+/e- need cut values for gamma
971 SetCutValue(cutForGamma, "gamma");
972 SetCutValue(cutForElectron, "e-");
973 SetCutValue(cutForPositron, "e+");
974
975 // SetCutValue(cutForProton, "proton");
976 // SetCutValue(cutForProton, "anti_proton");
977 // SetCutValue(cutForAlpha, "alpha");
978 // SetCutValue(cutForGenericIon, "GenericIon");
979
980 // SetCutValueForOthers(defaultCutValue);
981
983}
#define G4BestUnit(a, b)
void SetEnergyRange(G4double lowedge, G4double highedge)
static G4ProductionCutsTable * GetProductionCutsTable()
void SetCutValue(G4double aCut, const G4String &pname)
void DumpCutValuesTable(G4int flag=1)

The documentation for this class was generated from the following files: