Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DMXPhysicsList Class Reference

#include <DMXPhysicsList.hh>

Inheritance diagram for DMXPhysicsList:
Collaboration diagram for DMXPhysicsList:

Public Member Functions

 DMXPhysicsList ()
 
 ~DMXPhysicsList ()
 
virtual void SetCuts ()
 
- Public Member Functions inherited from G4VUserPhysicsList
 G4VUserPhysicsList ()
 
virtual ~G4VUserPhysicsList ()
 
 G4VUserPhysicsList (const G4VUserPhysicsList &)
 
G4VUserPhysicsListoperator= (const G4VUserPhysicsList &)
 
void Construct ()
 
void UseCoupledTransportation (G4bool vl=true)
 
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 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=0)
 
void SetParticleCuts (G4double cut, const G4String &particleName, G4Region *region=0)
 
void SetCutsForRegion (G4double aCut, const G4String &rname)
 
void ResetCuts ()
 obsolete methods More...
 
void SetApplyCuts (G4bool value, const G4String &name)
 
G4bool GetApplyCuts (const G4String &name) const
 
void RemoveProcessManager ()
 
void AddProcessManager (G4ParticleDefinition *newParticle, G4ProcessManager *newManager=0)
 
void CheckParticleList ()
 
void DisableCheckParticleList ()
 
G4int GetInstanceID () const
 
void InitializeWorker ()
 

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 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 G4VUserPhysicsList
static const G4VUPLManagerGetSubInstanceManager ()
 
- Protected Attributes inherited from G4VUserPhysicsList
G4ParticleTabletheParticleTable
 
G4int verboseLevel
 
G4double defaultCutValue
 
G4bool isSetDefaultCutValue
 
G4ProductionCutsTablefCutsTable
 
G4bool fRetrievePhysicsTable
 
G4bool fStoredInAscii
 
G4bool fIsCheckedForRetrievePhysicsTable
 
G4bool fIsRestoredCutValues
 
G4String directoryPhysicsTable
 
G4bool fDisableCheckParticleList
 
G4int g4vuplInstanceID
 
- Static Protected Attributes inherited from G4VUserPhysicsList
static G4RUN_DLL G4VUPLManager subInstanceManager
 

Detailed Description

Definition at line 56 of file DMXPhysicsList.hh.

Constructor & Destructor Documentation

DMXPhysicsList::DMXPhysicsList ( )

Definition at line 72 of file DMXPhysicsList.cc.

73 {
74 
75  defaultCutValue = 1.0*micrometer; //
76  cutForGamma = defaultCutValue;
77  cutForElectron = 1.0*nanometer;
78  cutForPositron = defaultCutValue;
79 
80  VerboseLevel = 1;
81  OpVerbLevel = 0;
82 
83  SetVerboseLevel(VerboseLevel);
84 }
static constexpr double nanometer
Definition: G4SIunits.hh:101
void SetVerboseLevel(G4int value)
static constexpr double micrometer
Definition: G4SIunits.hh:100

Here is the call graph for this function:

DMXPhysicsList::~DMXPhysicsList ( )

Definition at line 88 of file DMXPhysicsList.cc.

89 {;}

Member Function Documentation

void DMXPhysicsList::AddTransportation ( )
protectedvirtual

Definition at line 196 of file DMXPhysicsList.cc.

196  {
197 
199 
201  particleIterator->reset();
202  while( (*particleIterator)() ){
203  G4ParticleDefinition* particle = particleIterator->value();
204  G4ProcessManager* pmanager = particle->GetProcessManager();
205  G4String particleName = particle->GetParticleName();
206  // time cuts for ONLY neutrons:
207  if(particleName == "neutron")
208  pmanager->AddDiscreteProcess(new DMXMaxTimeCuts());
209  // Energy cuts to kill charged (embedded in method) particles:
210  pmanager->AddDiscreteProcess(new DMXMinEkineCuts());
211 
212  // Step limit applied to all particles:
213  pmanager->AddProcess(new G4StepLimiter, -1,-1,1);
214 
215  }
216 }
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
const G4String & GetParticleName() const
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
G4ProcessManager * GetProcessManager() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const

Here is the call graph for this function:

Here is the caller graph for this function:

void DMXPhysicsList::ConstructEM ( )
protectedvirtual

Definition at line 273 of file DMXPhysicsList.cc.

273  {
274 
275  //set a finer grid of the physic tables in order to improve precision
276  //former LowEnergy models have 200 bins up to 100 GeV
278  param->SetMaxEnergy(100*GeV);
279  param->SetNumberOfBinsPerDecade(20);
281  param->SetFluo(true);
282  param->SetPixe(true);
283  param->SetAuger(true);
286  if(!ad) {
288  }
289 
291  particleIterator->reset();
292  while( (*particleIterator)() ){
293  G4ParticleDefinition* particle = particleIterator->value();
294  G4ProcessManager* pmanager = particle->GetProcessManager();
295  G4String particleName = particle->GetParticleName();
296  G4String particleType = particle->GetParticleType();
297  G4double charge = particle->GetPDGCharge();
298 
299  if (particleName == "gamma")
300  {
301  //gamma
302  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
303  pmanager->AddDiscreteProcess(theRayleigh);
304 
305  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
306  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
307  pmanager->AddDiscreteProcess(thePhotoElectricEffect);
308 
309  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
310  theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
311  pmanager->AddDiscreteProcess(theComptonScattering);
312 
313  G4GammaConversion* theGammaConversion = new G4GammaConversion();
314  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
315  pmanager->AddDiscreteProcess(theGammaConversion);
316 
317  }
318  else if (particleName == "e-")
319  {
320  //electron
321  // process ordering: AddProcess(name, at rest, along step, post step)
322  // Multiple scattering
325  pmanager->AddProcess(msc,-1, 1, -1);
326 
327  // Ionisation
328  G4eIonisation* eIonisation = new G4eIonisation();
329  eIonisation->SetEmModel(new G4LivermoreIonisationModel());
330  eIonisation->SetStepFunction(0.2, 100*um); //improved precision in tracking
331  pmanager->AddProcess(eIonisation,-1, 2, 2);
332 
333  // Bremsstrahlung
334  G4eBremsstrahlung* eBremsstrahlung = new G4eBremsstrahlung();
335  eBremsstrahlung->SetEmModel(new G4LivermoreBremsstrahlungModel());
336  pmanager->AddProcess(eBremsstrahlung, -1,-3, 3);
337  }
338  else if (particleName == "e+")
339  {
340  //positron
343  pmanager->AddProcess(msc,-1, 1, 1);
344 
345  // Ionisation
346  G4eIonisation* eIonisation = new G4eIonisation();
347  eIonisation->SetStepFunction(0.2, 100*um); //
348  pmanager->AddProcess(eIonisation, -1, 2, 2);
349 
350  //Bremsstrahlung (use default, no low-energy available)
351  pmanager->AddProcess(new G4eBremsstrahlung(), -1,-1, 3);
352 
353  //Annihilation
354  pmanager->AddProcess(new G4eplusAnnihilation(),0,-1, 4);
355  }
356  else if( particleName == "mu+" ||
357  particleName == "mu-" )
358  {
359  //muon
360  pmanager->AddProcess(new G4eMultipleScattering, -1, 1, 1);
361  pmanager->AddProcess(new G4MuIonisation(), -1, 2, 2);
362  pmanager->AddProcess(new G4MuBremsstrahlung(), -1,-1, 3);
363  pmanager->AddProcess(new G4MuPairProduction(), -1,-1, 4);
364  if( particleName == "mu-" )
365  pmanager->AddProcess(new G4MuonMinusCapture(), 0,-1,-1);
366  }
367  else if (particleName == "proton" ||
368  particleName == "pi+" ||
369  particleName == "pi-")
370  {
371  //multiple scattering
372  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
373 
374  //ionisation
375  G4hIonisation* hIonisation = new G4hIonisation();
376  hIonisation->SetStepFunction(0.2, 50*um);
377  pmanager->AddProcess(hIonisation, -1, 2, 2);
378 
379  //bremmstrahlung
380  pmanager->AddProcess(new G4hBremsstrahlung, -1,-3, 3);
381  }
382  else if(particleName == "alpha" ||
383  particleName == "deuteron" ||
384  particleName == "triton" ||
385  particleName == "He3")
386  {
387  //multiple scattering
388  pmanager->AddProcess(new G4hMultipleScattering,-1,1,1);
389 
390  //ionisation
391  G4ionIonisation* ionIoni = new G4ionIonisation();
392  ionIoni->SetStepFunction(0.1, 20*um);
393  pmanager->AddProcess(ionIoni, -1, 2, 2);
394  }
395  else if (particleName == "GenericIon")
396  {
397  // OBJECT may be dynamically created as either a GenericIon or nucleus
398  // G4Nucleus exists and therefore has particle type nucleus
399  // genericIon:
400 
401  //multiple scattering
402  pmanager->AddProcess(new G4hMultipleScattering,-1,1,1);
403 
404  //ionisation
405  G4ionIonisation* ionIoni = new G4ionIonisation();
406  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
407  ionIoni->SetStepFunction(0.1, 20*um);
408  pmanager->AddProcess(ionIoni, -1, 2, 2);
409  }
410 
411  else if ((!particle->IsShortLived()) &&
412  (charge != 0.0) &&
413  (particle->GetParticleName() != "chargedgeantino"))
414  {
415  //all others charged particles except geantino
416  G4hMultipleScattering* aMultipleScattering = new G4hMultipleScattering();
417  G4hIonisation* ahadronIon = new G4hIonisation();
418 
419  //multiple scattering
420  pmanager->AddProcess(aMultipleScattering,-1,1,1);
421 
422  //ionisation
423  pmanager->AddProcess(ahadronIon, -1,2,2);
424  }
425 
426  }
427 }
static G4LossTableManager * Instance()
void SetMscStepLimitType(G4MscStepLimitType val)
void SetAuger(G4bool val)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void SetMaxEnergy(G4double val)
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
const G4String & GetParticleName() const
void SetEmModel(G4VEmModel *, G4int index=1)
static constexpr double um
Definition: G4SIunits.hh:113
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetNumberOfBinsPerDecade(G4int val)
const G4String & GetParticleType() const
void SetPixe(G4bool val)
G4ProcessManager * GetProcessManager() const
static G4EmParameters * Instance()
void SetEmModel(G4VEmModel *, G4int index=1)
static constexpr double GeV
Definition: G4SIunits.hh:217
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetFluo(G4bool val)
void SetStepLimitType(G4MscStepLimitType val)

Here is the call graph for this function:

Here is the caller graph for this function:

void DMXPhysicsList::ConstructGeneral ( )
protectedvirtual

Definition at line 878 of file DMXPhysicsList.cc.

878  {
879 
880  // Add Decay Process
881  G4Decay* theDecayProcess = new G4Decay();
883  particleIterator->reset();
884  while( (*particleIterator)() )
885  {
886  G4ParticleDefinition* particle = particleIterator->value();
887  G4ProcessManager* pmanager = particle->GetProcessManager();
888 
889  if (theDecayProcess->IsApplicable(*particle) && !particle->IsShortLived())
890  {
891  pmanager ->AddProcess(theDecayProcess);
892  // set ordering for PostStepDoIt and AtRestDoIt
893  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
894  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
895  }
896  }
897 
898  // Declare radioactive decay to the GenericIon in the IonTable.
899  const G4IonTable *theIonTable =
901  G4RadioactiveDecay *theRadioactiveDecay = new G4RadioactiveDecay();
902 
903  for (G4int i=0; i<theIonTable->Entries(); i++)
904  {
905  G4String particleName = theIonTable->GetParticle(i)->GetParticleName();
906  G4String particleType = theIonTable->GetParticle(i)->GetParticleType();
907 
908  if (particleName == "GenericIon")
909  {
910  G4ProcessManager* pmanager =
911  theIonTable->GetParticle(i)->GetProcessManager();
912  pmanager->SetVerboseLevel(VerboseLevel);
913  pmanager ->AddProcess(theRadioactiveDecay);
914  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxPostStep);
915  pmanager ->SetProcessOrdering(theRadioactiveDecay, idxAtRest);
916  }
917  }
918 }
G4int Entries() const
Definition: G4IonTable.cc:1689
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
void SetVerboseLevel(G4int value)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
const G4String & GetParticleType() const
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
static G4ParticleTable * GetParticleTable()
G4ProcessManager * GetProcessManager() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
G4ParticleDefinition * GetParticle(G4int index) const
Definition: G4IonTable.cc:1643

Here is the call graph for this function:

Here is the caller graph for this function:

void DMXPhysicsList::ConstructHad ( )
protectedvirtual

Definition at line 571 of file DMXPhysicsList.cc.

572 {
573  //Elastic models
574  const G4double elastic_elimitPi = 1.0*GeV;
575 
576  G4HadronElastic* elastic_lhep0 = new G4HadronElastic();
577  G4HadronElastic* elastic_lhep1 = new G4HadronElastic();
578  elastic_lhep1->SetMaxEnergy( elastic_elimitPi );
579  G4ChipsElasticModel* elastic_chip = new G4ChipsElasticModel();
580  G4ElasticHadrNucleusHE* elastic_he = new G4ElasticHadrNucleusHE();
581  elastic_he->SetMinEnergy( elastic_elimitPi );
582 
583 
584  // Inelastic scattering
585  const G4double theFTFMin0 = 0.0*GeV;
586  const G4double theFTFMin1 = 4.0*GeV;
587  const G4double theFTFMax = 100.0*TeV;
588  const G4double theBERTMin0 = 0.0*GeV;
589  const G4double theBERTMin1 = 19.0*MeV;
590  const G4double theBERTMax = 5.0*GeV;
591  const G4double theHPMin = 0.0*GeV;
592  const G4double theHPMax = 20.0*MeV;
593 
594  G4FTFModel * theStringModel = new G4FTFModel;
596  theStringModel->SetFragmentationModel( theStringDecay );
597  G4PreCompoundModel * thePreEquilib = new G4PreCompoundModel( new G4ExcitationHandler );
598  G4GeneratorPrecompoundInterface * theCascade = new G4GeneratorPrecompoundInterface( thePreEquilib );
599 
600  G4TheoFSGenerator * theFTFModel0 = new G4TheoFSGenerator( "FTFP" );
601  theFTFModel0->SetHighEnergyGenerator( theStringModel );
602  theFTFModel0->SetTransport( theCascade );
603  theFTFModel0->SetMinEnergy( theFTFMin0 );
604  theFTFModel0->SetMaxEnergy( theFTFMax );
605 
606  G4TheoFSGenerator * theFTFModel1 = new G4TheoFSGenerator( "FTFP" );
607  theFTFModel1->SetHighEnergyGenerator( theStringModel );
608  theFTFModel1->SetTransport( theCascade );
609  theFTFModel1->SetMinEnergy( theFTFMin1 );
610  theFTFModel1->SetMaxEnergy( theFTFMax );
611 
612  G4CascadeInterface * theBERTModel0 = new G4CascadeInterface;
613  theBERTModel0->SetMinEnergy( theBERTMin0 );
614  theBERTModel0->SetMaxEnergy( theBERTMax );
615 
616  G4CascadeInterface * theBERTModel1 = new G4CascadeInterface;
617  theBERTModel1->SetMinEnergy( theBERTMin1 );
618  theBERTModel1->SetMaxEnergy( theBERTMax );
619 
622  G4ComponentGGNuclNuclXsc * ggNuclNuclXsec = new G4ComponentGGNuclNuclXsc();
623  G4VCrossSectionDataSet * theGGNuclNuclData = new G4CrossSectionInelastic(ggNuclNuclXsec);
624 
626  particleIterator->reset();
627  while ((*particleIterator)())
628  {
629  G4ParticleDefinition* particle = particleIterator->value();
630  G4ProcessManager* pmanager = particle->GetProcessManager();
631  G4String particleName = particle->GetParticleName();
632 
633  if (particleName == "pi+")
634  {
635  // Elastic scattering
636  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
637  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
638  theElasticProcess->RegisterMe( elastic_lhep1 );
639  theElasticProcess->RegisterMe( elastic_he );
640  pmanager->AddDiscreteProcess( theElasticProcess );
641  //Inelastic scattering
642  G4PionPlusInelasticProcess* theInelasticProcess =
643  new G4PionPlusInelasticProcess("inelastic");
644  theInelasticProcess->AddDataSet( thePiData );
645  theInelasticProcess->RegisterMe( theFTFModel1 );
646  theInelasticProcess->RegisterMe( theBERTModel0 );
647  pmanager->AddDiscreteProcess( theInelasticProcess );
648  }
649 
650  else if (particleName == "pi-")
651  {
652  // Elastic scattering
653  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
654  theElasticProcess->AddDataSet( new G4BGGPionElasticXS( particle ) );
655  theElasticProcess->RegisterMe( elastic_lhep1 );
656  theElasticProcess->RegisterMe( elastic_he );
657  pmanager->AddDiscreteProcess( theElasticProcess );
658  //Inelastic scattering
659  G4PionMinusInelasticProcess* theInelasticProcess =
660  new G4PionMinusInelasticProcess("inelastic");
661  theInelasticProcess->AddDataSet( thePiData );
662  theInelasticProcess->RegisterMe( theFTFModel1 );
663  theInelasticProcess->RegisterMe( theBERTModel0 );
664  pmanager->AddDiscreteProcess( theInelasticProcess );
665  //Absorption
667  }
668 
669  else if (particleName == "kaon+")
670  {
671  // Elastic scattering
672  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
673  theElasticProcess->RegisterMe( elastic_lhep0 );
674  pmanager->AddDiscreteProcess( theElasticProcess );
675  // Inelastic scattering
676  G4KaonPlusInelasticProcess* theInelasticProcess =
677  new G4KaonPlusInelasticProcess("inelastic");
678  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
679  GetCrossSectionDataSet(G4ChipsKaonPlusInelasticXS::Default_Name()));
680  theInelasticProcess->RegisterMe( theFTFModel1 );
681  theInelasticProcess->RegisterMe( theBERTModel0 );
682  pmanager->AddDiscreteProcess( theInelasticProcess );
683  }
684 
685  else if (particleName == "kaon0S")
686  {
687  // Elastic scattering
688  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
689  theElasticProcess->RegisterMe( elastic_lhep0 );
690  pmanager->AddDiscreteProcess( theElasticProcess );
691  // Inelastic scattering
692  G4KaonZeroSInelasticProcess* theInelasticProcess =
693  new G4KaonZeroSInelasticProcess("inelastic");
694  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
695  GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
696  theInelasticProcess->RegisterMe( theFTFModel1 );
697  theInelasticProcess->RegisterMe( theBERTModel0 );
698  pmanager->AddDiscreteProcess( theInelasticProcess );
699  }
700 
701  else if (particleName == "kaon0L")
702  {
703  // Elastic scattering
704  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
705  theElasticProcess->RegisterMe( elastic_lhep0 );
706  pmanager->AddDiscreteProcess( theElasticProcess );
707  // Inelastic scattering
708  G4KaonZeroLInelasticProcess* theInelasticProcess =
709  new G4KaonZeroLInelasticProcess("inelastic");
710  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
711  GetCrossSectionDataSet(G4ChipsKaonZeroInelasticXS::Default_Name()));
712  theInelasticProcess->RegisterMe( theFTFModel1 );
713  theInelasticProcess->RegisterMe( theBERTModel0 );
714  pmanager->AddDiscreteProcess( theInelasticProcess );
715  }
716 
717  else if (particleName == "kaon-")
718  {
719  // Elastic scattering
720  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
721  theElasticProcess->RegisterMe( elastic_lhep0 );
722  pmanager->AddDiscreteProcess( theElasticProcess );
723  // Inelastic scattering
724  G4KaonMinusInelasticProcess* theInelasticProcess =
725  new G4KaonMinusInelasticProcess("inelastic");
726  theInelasticProcess->AddDataSet( G4CrossSectionDataSetRegistry::Instance()->
727  GetCrossSectionDataSet(G4ChipsKaonMinusInelasticXS::Default_Name()));
728  theInelasticProcess->RegisterMe( theFTFModel1 );
729  theInelasticProcess->RegisterMe( theBERTModel0 );
730  pmanager->AddDiscreteProcess( theInelasticProcess );
732  }
733 
734  else if (particleName == "proton")
735  {
736  // Elastic scattering
737  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
739  GetCrossSectionDataSet(G4ChipsProtonElasticXS::Default_Name()));
740  theElasticProcess->RegisterMe( elastic_chip );
741  pmanager->AddDiscreteProcess( theElasticProcess );
742  // Inelastic scattering
743  G4ProtonInelasticProcess* theInelasticProcess =
744  new G4ProtonInelasticProcess("inelastic");
745  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Proton::Proton() ) );
746  theInelasticProcess->RegisterMe( theFTFModel1 );
747  theInelasticProcess->RegisterMe( theBERTModel0 );
748  pmanager->AddDiscreteProcess( theInelasticProcess );
749  }
750  else if (particleName == "anti_proton")
751  {
752  // Elastic scattering
753  const G4double elastic_elimitAntiNuc = 100.0*MeV;
754  G4AntiNuclElastic* elastic_anuc = new G4AntiNuclElastic();
755  elastic_anuc->SetMinEnergy( elastic_elimitAntiNuc );
756  G4CrossSectionElastic* elastic_anucxs = new G4CrossSectionElastic( elastic_anuc->GetComponentCrossSection() );
757  G4HadronElastic* elastic_lhep2 = new G4HadronElastic();
758  elastic_lhep2->SetMaxEnergy( elastic_elimitAntiNuc );
759  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
760  theElasticProcess->AddDataSet( elastic_anucxs );
761  theElasticProcess->RegisterMe( elastic_lhep2 );
762  theElasticProcess->RegisterMe( elastic_anuc );
763  pmanager->AddDiscreteProcess( theElasticProcess );
764  // Inelastic scattering
765  G4AntiProtonInelasticProcess* theInelasticProcess =
766  new G4AntiProtonInelasticProcess("inelastic");
767  theInelasticProcess->AddDataSet( theAntiNucleonData );
768  theInelasticProcess->RegisterMe( theFTFModel0 );
769  pmanager->AddDiscreteProcess( theInelasticProcess );
770  // Absorption
772  }
773 
774  else if (particleName == "neutron") {
775  // elastic scattering
776  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
777  theElasticProcess->AddDataSet(G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsNeutronElasticXS::Default_Name()));
778  G4HadronElastic* elastic_neutronChipsModel = new G4ChipsElasticModel();
779  elastic_neutronChipsModel->SetMinEnergy( 19.0*MeV );
780  theElasticProcess->RegisterMe( elastic_neutronChipsModel );
781  G4ParticleHPElastic * theElasticNeutronHP = new G4ParticleHPElastic;
782  theElasticNeutronHP->SetMinEnergy( theHPMin );
783  theElasticNeutronHP->SetMaxEnergy( theHPMax );
784  theElasticProcess->RegisterMe( theElasticNeutronHP );
785  theElasticProcess->AddDataSet( new G4ParticleHPElasticData );
786  pmanager->AddDiscreteProcess( theElasticProcess );
787  // inelastic scattering
788  G4NeutronInelasticProcess* theInelasticProcess =
789  new G4NeutronInelasticProcess("inelastic");
790  theInelasticProcess->AddDataSet( new G4BGGNucleonInelasticXS( G4Neutron::Neutron() ) );
791  theInelasticProcess->RegisterMe( theFTFModel1 );
792  theInelasticProcess->RegisterMe( theBERTModel1 );
793  G4ParticleHPInelastic * theNeutronInelasticHPModel = new G4ParticleHPInelastic;
794  theNeutronInelasticHPModel->SetMinEnergy( theHPMin );
795  theNeutronInelasticHPModel->SetMaxEnergy( theHPMax );
796  theInelasticProcess->RegisterMe( theNeutronInelasticHPModel );
797  theInelasticProcess->AddDataSet( new G4ParticleHPInelasticData );
798  pmanager->AddDiscreteProcess(theInelasticProcess);
799  // capture
800  G4HadronCaptureProcess* theCaptureProcess =
802  G4ParticleHPCapture * theLENeutronCaptureModel = new G4ParticleHPCapture;
803  theLENeutronCaptureModel->SetMinEnergy(theHPMin);
804  theLENeutronCaptureModel->SetMaxEnergy(theHPMax);
805  theCaptureProcess->RegisterMe(theLENeutronCaptureModel);
806  theCaptureProcess->AddDataSet( new G4ParticleHPCaptureData);
807  pmanager->AddDiscreteProcess(theCaptureProcess);
808 
809  }
810  else if (particleName == "anti_neutron")
811  {
812  // Elastic scattering
813  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
814  theElasticProcess->RegisterMe( elastic_lhep0 );
815  pmanager->AddDiscreteProcess( theElasticProcess );
816  // Inelastic scattering (include annihilation on-fly)
817  G4AntiNeutronInelasticProcess* theInelasticProcess =
818  new G4AntiNeutronInelasticProcess("inelastic");
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  pmanager->AddDiscreteProcess( theElasticProcess );
830  // Inelastic scattering
831  G4DeuteronInelasticProcess* theInelasticProcess =
832  new G4DeuteronInelasticProcess("inelastic");
833  theInelasticProcess->AddDataSet( theGGNuclNuclData );
834  theInelasticProcess->RegisterMe( theFTFModel1 );
835  theInelasticProcess->RegisterMe( theBERTModel0 );
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  pmanager->AddDiscreteProcess( theElasticProcess );
845  // Inelastic scattering
846  G4TritonInelasticProcess* theInelasticProcess =
847  new G4TritonInelasticProcess("inelastic");
848  theInelasticProcess->AddDataSet( theGGNuclNuclData );
849  theInelasticProcess->RegisterMe( theFTFModel1 );
850  theInelasticProcess->RegisterMe( theBERTModel0 );
851  pmanager->AddDiscreteProcess( theInelasticProcess );
852  }
853  else if (particleName == "alpha")
854  {
855  // Elastic scattering
856  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
857  theElasticProcess->RegisterMe( elastic_lhep0 );
858  pmanager->AddDiscreteProcess( theElasticProcess );
859  // Inelastic scattering
860  G4AlphaInelasticProcess* theInelasticProcess =
861  new G4AlphaInelasticProcess("inelastic");
862  theInelasticProcess->AddDataSet( theGGNuclNuclData );
863  theInelasticProcess->RegisterMe( theFTFModel1 );
864  theInelasticProcess->RegisterMe( theBERTModel0 );
865  pmanager->AddDiscreteProcess( theInelasticProcess );
866  }
867 
868  }
869 }
void SetFragmentationModel(G4VStringFragmentation *aModel)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ComponentAntiNuclNuclearXS * GetComponentCrossSection()
const G4String & GetParticleName() const
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
void RegisterMe(G4HadronicInteraction *a)
static constexpr double TeV
Definition: G4SIunits.hh:218
void SetMinEnergy(G4double anEnergy)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
static const char * Default_Name()
static G4CrossSectionDataSetRegistry * Instance()
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4ProcessManager * GetProcessManager() const
static constexpr double GeV
Definition: G4SIunits.hh:217
void SetMaxEnergy(const G4double anEnergy)
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
static const char * Default_Name()
void SetTransport(G4VIntraNuclearTransportModel *const value)
static constexpr double MeV
Definition: G4SIunits.hh:214
G4int AddRestProcess(G4VProcess *aProcess, G4int ord=ordDefault)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void DMXPhysicsList::ConstructOp ( )
protectedvirtual

Definition at line 436 of file DMXPhysicsList.cc.

437 {
438  // default scintillation process
439  G4Scintillation* theScintProcessDef = new G4Scintillation("Scintillation");
440  // theScintProcessDef->DumpPhysicsTable();
441  theScintProcessDef->SetTrackSecondariesFirst(true);
442  theScintProcessDef->SetScintillationYieldFactor(1.0); //
443  theScintProcessDef->SetScintillationExcitationRatio(0.0); //
444  theScintProcessDef->SetVerboseLevel(OpVerbLevel);
445 
446  // scintillation process for alpha:
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  G4Scintillation* theScintProcessNuc = new G4Scintillation("Scintillation");
456  // theScintProcessNuc->DumpPhysicsTable();
457  theScintProcessNuc->SetTrackSecondariesFirst(true);
458  theScintProcessNuc->SetScintillationYieldFactor(0.2);
459  theScintProcessNuc->SetScintillationExcitationRatio(1.0);
460  theScintProcessNuc->SetVerboseLevel(OpVerbLevel);
461 
462  // optical processes
463  G4OpAbsorption* theAbsorptionProcess = new G4OpAbsorption();
464  // G4OpRayleigh* theRayleighScatteringProcess = new G4OpRayleigh();
465  G4OpBoundaryProcess* theBoundaryProcess = new G4OpBoundaryProcess();
466  // theAbsorptionProcess->DumpPhysicsTable();
467  // theRayleighScatteringProcess->DumpPhysicsTable();
468  theAbsorptionProcess->SetVerboseLevel(OpVerbLevel);
469  // theRayleighScatteringProcess->SetVerboseLevel(OpVerbLevel);
470  theBoundaryProcess->SetVerboseLevel(OpVerbLevel);
471 
473  particleIterator->reset();
474  while( (*particleIterator)() )
475  {
476  G4ParticleDefinition* particle = particleIterator->value();
477  G4ProcessManager* pmanager = particle->GetProcessManager();
478  G4String particleName = particle->GetParticleName();
479  if (theScintProcessDef->IsApplicable(*particle)) {
480  // if(particle->GetPDGMass() > 5.0*GeV)
481  if(particle->GetParticleName() == "GenericIon") {
482  pmanager->AddProcess(theScintProcessNuc); // AtRestDiscrete
483  pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxAtRest);
484  pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxPostStep);
485  }
486  else if(particle->GetParticleName() == "alpha") {
487  pmanager->AddProcess(theScintProcessAlpha);
488  pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxAtRest);
489  pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxPostStep);
490  }
491  else {
492  pmanager->AddProcess(theScintProcessDef);
493  pmanager->SetProcessOrderingToLast(theScintProcessDef,idxAtRest);
494  pmanager->SetProcessOrderingToLast(theScintProcessDef,idxPostStep);
495  }
496  }
497 
498  if (particleName == "opticalphoton") {
499  pmanager->AddDiscreteProcess(theAbsorptionProcess);
500  // pmanager->AddDiscreteProcess(theRayleighScatteringProcess);
501  pmanager->AddDiscreteProcess(theBoundaryProcess);
502  }
503  }
504 }
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
const G4String & GetParticleName() const
void SetScintillationYieldFactor(const G4double yieldfactor)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetScintillationExcitationRatio(const G4double ratio)
void SetTrackSecondariesFirst(const G4bool state)
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
G4ProcessManager * GetProcessManager() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437

Here is the call graph for this function:

Here is the caller graph for this function:

void DMXPhysicsList::ConstructParticle ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 93 of file DMXPhysicsList.cc.

94 {
95 
96  // In this method, static member functions should be called
97  // for all particles which you want to use.
98  // This ensures that objects of these particle types will be
99  // created in the program.
100 
101  ConstructMyBosons();
102  ConstructMyLeptons();
103  ConstructMyHadrons();
104  ConstructMyShortLiveds();
105 
106 }
void DMXPhysicsList::ConstructProcess ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 175 of file DMXPhysicsList.cc.

176 {
177 
179 
180  ConstructEM();
181 
182  ConstructOp();
183 
184  ConstructHad();
185 
187 
188 }
virtual void ConstructGeneral()
virtual void AddTransportation()
virtual void ConstructHad()
virtual void ConstructOp()
virtual void ConstructEM()

Here is the call graph for this function:

void DMXPhysicsList::SetCuts ( )
virtual

Reimplemented from G4VUserPhysicsList.

Definition at line 921 of file DMXPhysicsList.cc.

922 {
923 
924  if (verboseLevel >1)
925  G4cout << "DMXPhysicsList::SetCuts:";
926 
927  if (verboseLevel>0){
928  G4cout << "DMXPhysicsList::SetCuts:";
929  G4cout << "CutLength : "
930  << G4BestUnit(defaultCutValue,"Length") << G4endl;
931  }
932 
933  //special for low energy physics
934  G4double lowlimit=250*eV;
936 
937  // set cut values for gamma at first and for e- second and next for e+,
938  // because some processes for e+/e- need cut values for gamma
939  SetCutValue(cutForGamma, "gamma");
940  SetCutValue(cutForElectron, "e-");
941  SetCutValue(cutForPositron, "e+");
942 
944 }
void SetCutValue(G4double aCut, const G4String &pname)
void SetEnergyRange(G4double lowedge, G4double highedge)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
static G4ProductionCutsTable * GetProductionCutsTable()
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


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