Geant4  10.02.p03
B03PhysicsList Class Reference

#include <B03PhysicsList.hh>

Inheritance diagram for B03PhysicsList:
Collaboration diagram for B03PhysicsList:

Public Member Functions

 B03PhysicsList ()
 
virtual ~B03PhysicsList ()
 
void AddParallelWorldName (G4String &pname)
 
void AddBiasing (G4GeometrySampler *mgs, G4String &pname)
 
- 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 SetCuts ()
 
virtual void ConstructGeneral ()
 
virtual void ConstructEM ()
 
virtual void ConstructHad ()
 
virtual void ConstructLeptHad ()
 
void AddScoringProcess ()
 
void AddBiasingProcess ()
 
void ConstructAllBosons ()
 
void ConstructAllLeptons ()
 
void ConstructAllMesons ()
 
void ConstructAllBaryons ()
 
void ConstructAllIons ()
 
void ConstructAllShortLiveds ()
 
- 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
 

Private Attributes

std::vector< G4StringfParaWorldName
 
G4String fBiasWorldName
 
G4GeometrySampler * fGeomSampler
 

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 41 of file B03PhysicsList.hh.

Constructor & Destructor Documentation

◆ B03PhysicsList()

B03PhysicsList::B03PhysicsList ( )

Definition at line 56 of file B03PhysicsList.cc.

57 {
58  fParaWorldName.clear();
59  SetVerboseLevel(1);
60 }
void SetVerboseLevel(G4int value)
std::vector< G4String > fParaWorldName
Here is the call graph for this function:

◆ ~B03PhysicsList()

B03PhysicsList::~B03PhysicsList ( )
virtual

Definition at line 64 of file B03PhysicsList.cc.

65 {
66  fParaWorldName.clear();
67 }
std::vector< G4String > fParaWorldName

Member Function Documentation

◆ AddBiasing()

void B03PhysicsList::AddBiasing ( G4GeometrySampler *  mgs,
G4String pname 
)
inline

Definition at line 51 of file B03PhysicsList.hh.

string pname
Definition: eplot.py:33
G4String fBiasWorldName
G4GeometrySampler * fGeomSampler
Here is the call graph for this function:

◆ AddBiasingProcess()

void B03PhysicsList::AddBiasingProcess ( )
protected

Definition at line 674 of file B03PhysicsList.cc.

674  {
675 
676 
677  G4cout << " Preparing Importance Sampling with GhostWorld "
678  << fBiasWorldName << G4endl;
679  fGeomSampler->SetParallel(true); // parallelworld
681  fGeomSampler->SetWorld(iStore->GetParallelWorldVolumePointer());
682  // fGeomSampler->PrepareImportanceSampling(G4IStore::
683  // GetInstance(fBiasWorldName), 0);
684  static G4bool first = true;
685  if(first) {
686  fGeomSampler->PrepareImportanceSampling(iStore, 0);
687 
688  fGeomSampler->Configure();
689  G4cout << " GeomSampler Configured!!! " << G4endl;
690  first = false;
691  }
692 
693 #ifdef G4MULTITHREADED
694  fGeomSampler->AddProcess();
695 #else
696  G4cout << " Running in singlethreaded mode!!! " << G4endl;
697 #endif
698 
699 }
virtual const G4VPhysicalVolume * GetParallelWorldVolumePointer() const
Definition: G4IStore.cc:93
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
static G4IStore * GetInstance()
Definition: G4IStore.cc:210
#define G4endl
Definition: G4ios.hh:61
G4String fBiasWorldName
G4GeometrySampler * fGeomSampler
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AddParallelWorldName()

void B03PhysicsList::AddParallelWorldName ( G4String pname)
inline

Definition at line 48 of file B03PhysicsList.hh.

49  {fParaWorldName.push_back(pname);}
std::vector< G4String > fParaWorldName
Here is the caller graph for this function:

◆ AddScoringProcess()

void B03PhysicsList::AddScoringProcess ( )
protected

Definition at line 646 of file B03PhysicsList.cc.

646  {
647 
648  G4int npw = fParaWorldName.size();
649  for ( G4int i = 0; i < npw; i++){
650  G4String procName = "ParaWorldProc_"+fParaWorldName[i];
651  G4ParallelWorldProcess* theParallelWorldProcess
652  = new G4ParallelWorldProcess(procName);
653  theParallelWorldProcess->SetParallelWorld(fParaWorldName[i]);
654 
656  particleIterator->reset();
657  while( (*particleIterator)() ){
658  G4ParticleDefinition* particle = particleIterator->value();
659  G4ProcessManager* pmanager = particle->GetProcessManager();
660  pmanager->AddProcess(theParallelWorldProcess);
661  if(theParallelWorldProcess->IsAtRestRequired(particle))
662  {pmanager->SetProcessOrdering(theParallelWorldProcess, idxAtRest, 9900);}
663  pmanager->SetProcessOrderingToSecond(theParallelWorldProcess, idxAlongStep);
664  pmanager->SetProcessOrdering(theParallelWorldProcess, idxPostStep, 9900);
665  }
666  }
667 
668 }
void SetProcessOrderingToSecond(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
G4bool IsAtRestRequired(G4ParticleDefinition *)
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
void SetParallelWorld(G4String parallelWorldName)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
std::vector< G4String > fParaWorldName
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructAllBaryons()

void B03PhysicsList::ConstructAllBaryons ( )
protected

Definition at line 115 of file B03PhysicsList.cc.

116 {
117  // Construct all barions
118  G4BaryonConstructor pConstructor;
119  pConstructor.ConstructParticle();
120 }
static void ConstructParticle()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructAllBosons()

void B03PhysicsList::ConstructAllBosons ( )
protected

Definition at line 88 of file B03PhysicsList.cc.

89 {
90  // Construct all bosons
91  G4BosonConstructor pConstructor;
92  pConstructor.ConstructParticle();
93 }
static void ConstructParticle()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructAllIons()

void B03PhysicsList::ConstructAllIons ( )
protected

Definition at line 124 of file B03PhysicsList.cc.

125 {
126  // Construct light ions
127  G4IonConstructor pConstructor;
128  pConstructor.ConstructParticle();
129 }
static void ConstructParticle()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructAllLeptons()

void B03PhysicsList::ConstructAllLeptons ( )
protected

Definition at line 97 of file B03PhysicsList.cc.

98 {
99  // Construct all leptons
100  G4LeptonConstructor pConstructor;
101  pConstructor.ConstructParticle();
102 }
static void ConstructParticle()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructAllMesons()

void B03PhysicsList::ConstructAllMesons ( )
protected

Definition at line 106 of file B03PhysicsList.cc.

107 {
108  // Construct all mesons
109  G4MesonConstructor pConstructor;
110  pConstructor.ConstructParticle();
111 }
static void ConstructParticle()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructAllShortLiveds()

void B03PhysicsList::ConstructAllShortLiveds ( )
protected

Definition at line 133 of file B03PhysicsList.cc.

134 {
135  // Construct resonaces and quarks
136  G4ShortLivedConstructor pConstructor;
137  pConstructor.ConstructParticle();
138 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructEM()

void B03PhysicsList::ConstructEM ( )
protectedvirtual

Definition at line 173 of file B03PhysicsList.cc.

174 {
176  particleIterator->reset();
177  while( (*particleIterator)() ){
178  G4ParticleDefinition* particle = particleIterator->value();
179  G4ProcessManager* pmanager = particle->GetProcessManager();
180  G4String particleName = particle->GetParticleName();
181 
182  if (particleName == "gamma") {
183  // gamma
184  // Construct processes for gamma
185  pmanager->AddDiscreteProcess(new G4GammaConversion());
186  pmanager->AddDiscreteProcess(new G4ComptonScattering());
187  pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
188 
189  } else if (particleName == "e-") {
190  //electron
191  // Construct processes for electron
192  pmanager->AddProcess(new G4eMultipleScattering(),-1,1,1);
193  pmanager->AddProcess(new G4eIonisation(),-1,2,2);
194  pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
195 
196  } else if (particleName == "e+") {
197  //positron
198  // Construct processes for positron
199  pmanager->AddProcess(new G4eMultipleScattering(),-1,1,1);
200 
201  pmanager->AddProcess(new G4eIonisation(),-1,2,2);
202  pmanager->AddProcess(new G4eBremsstrahlung(),-1,-1,3);
203  pmanager->AddProcess(new G4eplusAnnihilation(),0,-1,4);
204 
205  } else if( particleName == "mu+" ||
206  particleName == "mu-" ) {
207  //muon
208  // Construct processes for muon+
209  pmanager->AddProcess(new G4MuMultipleScattering(),-1,1,1);
210  pmanager->AddProcess(new G4MuIonisation(),-1,2,2);
211  pmanager->AddProcess(new G4MuBremsstrahlung(),-1,-1,3);
212  pmanager->AddProcess(new G4MuPairProduction(),-1,-1,4);
213 
214  } else if( particleName == "GenericIon" ) {
215  pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
216  pmanager->AddProcess(new G4hIonisation(),-1,2,2);
217  } else {
218  if ((particle->GetPDGCharge() != 0.0) &&
219  (particle->GetParticleName() != "chargedgeantino")&&
220  (!particle->IsShortLived()) ) {
221  // all others charged particles except geantino
222  pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
223  pmanager->AddProcess(new G4hIonisation(),-1,2,2);
224  }
225  }
226  }
227 }
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4double GetPDGCharge() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructGeneral()

void B03PhysicsList::ConstructGeneral ( )
protectedvirtual

Definition at line 613 of file B03PhysicsList.cc.

614 {
615  G4Decay* theDecayProcess = new G4Decay();
617  particleIterator->reset();
618  while( (*particleIterator)() ){
619  G4ParticleDefinition* particle = particleIterator->value();
620  G4ProcessManager* pmanager = particle->GetProcessManager();
621  if (theDecayProcess->IsApplicable(*particle)) {
622  pmanager ->AddProcess(theDecayProcess);
623  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
624  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
625  }
626  }
627 }
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
G4ProcessManager * GetProcessManager() const
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructHad()

void B03PhysicsList::ConstructHad ( )
protectedvirtual

Definition at line 301 of file B03PhysicsList.cc.

302 {
303  // this will be the model class for high energies
304  G4TheoFSGenerator* theTheoModel = new G4TheoFSGenerator;
305  G4TheoFSGenerator* antiBHighEnergyModel = new G4TheoFSGenerator;
306 
307  // all models for treatment of thermal nucleus
308  G4Evaporation* theEvaporation = new G4Evaporation;
309  G4FermiBreakUp* theFermiBreakUp = new G4FermiBreakUp;
310  G4StatMF* theMF = new G4StatMF;
311 
312  // Evaporation logic
313  G4ExcitationHandler* theHandler = new G4ExcitationHandler;
314  theHandler->SetEvaporation(theEvaporation);
315  theHandler->SetFermiModel(theFermiBreakUp);
316  theHandler->SetMultiFragmentation(theMF);
317  theHandler->SetMaxAandZForFermiBreakUp(12, 6);
318  theHandler->SetMinEForMultiFrag(3*MeV);
319 
320  // Pre equilibrium stage
321  G4PreCompoundModel* thePreEquilib = new G4PreCompoundModel(theHandler);
322 
323  // a no-cascade generator-precompound interaface
324  G4GeneratorPrecompoundInterface* theCascade =
326  theCascade->SetDeExcitation(thePreEquilib);
327 
328  // Bertini cascade
329  G4CascadeInterface* bertini = new G4CascadeInterface;
330  bertini->SetMaxEnergy(22*MeV);
331 
332  // here come the high energy parts
333  G4VPartonStringModel* theStringModel;
334  theStringModel = new G4FTFModel;
335  theTheoModel->SetTransport(theCascade);
336  theTheoModel->SetHighEnergyGenerator(theStringModel);
337  theTheoModel->SetMinEnergy(19*GeV);
338  theTheoModel->SetMaxEnergy(100*TeV);
339 
340  G4VLongitudinalStringDecay* theFragmentation = new G4QGSMFragmentation;
341  G4ExcitedStringDecay* theStringDecay =
342  new G4ExcitedStringDecay(theFragmentation);
343  theStringModel->SetFragmentationModel(theStringDecay);
344 
345  // high energy model for anti-baryons
346  antiBHighEnergyModel = new G4TheoFSGenerator("ANTI-FTFP");
347  G4FTFModel* antiBStringModel = new G4FTFModel;
348  G4ExcitedStringDecay* stringDecay =
350  antiBStringModel->SetFragmentationModel(stringDecay);
351 
352  G4GeneratorPrecompoundInterface* antiBCascade =
354  G4PreCompoundModel* preEquilib =
356  antiBCascade->SetDeExcitation(preEquilib);
357 
358  antiBHighEnergyModel->SetTransport(antiBCascade);
359  antiBHighEnergyModel->SetHighEnergyGenerator(antiBStringModel);
360  antiBHighEnergyModel->SetMinEnergy(0.0);
361  antiBHighEnergyModel->SetMaxEnergy(20*TeV);
362 
363  // Light ion models
364  G4BinaryLightIonReaction* binaryCascade = new G4BinaryLightIonReaction;
365  binaryCascade->SetMinEnergy(0.0);
366  binaryCascade->SetMaxEnergy(110*MeV);
367 
368  G4QMDReaction* qmd = new G4QMDReaction;
369  qmd->SetMinEnergy(100*MeV);
370  qmd->SetMaxEnergy(10*GeV);
371 
375 
376  // Elastic process
377  G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
378  G4HadronElastic* theElasticModel = new G4HadronElastic;
379  theElasticProcess->RegisterMe(theElasticModel);
380 
382  particleIterator->reset();
383  while ((*particleIterator)()) {
384  G4ParticleDefinition* particle = particleIterator->value();
385  G4ProcessManager* pmanager = particle->GetProcessManager();
386  G4String particleName = particle->GetParticleName();
387 
388  if (particleName == "pi+") {
389  pmanager->AddDiscreteProcess(theElasticProcess);
390  G4PionPlusInelasticProcess* theInelasticProcess =
391  new G4PionPlusInelasticProcess("inelastic");
392  theInelasticProcess->RegisterMe(bertini);
393  theInelasticProcess->RegisterMe(theTheoModel);
394  pmanager->AddDiscreteProcess(theInelasticProcess);
395  } else if (particleName == "pi-") {
396  pmanager->AddDiscreteProcess(theElasticProcess);
397  G4PionMinusInelasticProcess* theInelasticProcess =
398  new G4PionMinusInelasticProcess("inelastic");
399  theInelasticProcess->RegisterMe(bertini);
400  theInelasticProcess->RegisterMe(theTheoModel);
401  pmanager->AddDiscreteProcess(theInelasticProcess);
402  } else if (particleName == "kaon+") {
403  pmanager->AddDiscreteProcess(theElasticProcess);
404  G4KaonPlusInelasticProcess* theInelasticProcess =
405  new G4KaonPlusInelasticProcess("inelastic");
406  theInelasticProcess->RegisterMe(bertini);
407  theInelasticProcess->RegisterMe(theTheoModel);
408  pmanager->AddDiscreteProcess(theInelasticProcess);
409  }
410  else if (particleName == "kaon0S") {
411  pmanager->AddDiscreteProcess(theElasticProcess);
412  G4KaonZeroSInelasticProcess* theInelasticProcess =
413  new G4KaonZeroSInelasticProcess("inelastic");
414  theInelasticProcess->RegisterMe(bertini);
415  theInelasticProcess->RegisterMe(theTheoModel);
416  pmanager->AddDiscreteProcess(theInelasticProcess);
417  }
418  else if (particleName == "kaon0L") {
419  pmanager->AddDiscreteProcess(theElasticProcess);
420  G4KaonZeroLInelasticProcess* theInelasticProcess =
421  new G4KaonZeroLInelasticProcess("inelastic");
422  theInelasticProcess->RegisterMe(bertini);
423  theInelasticProcess->RegisterMe(theTheoModel);
424  pmanager->AddDiscreteProcess(theInelasticProcess);
425  }
426  else if (particleName == "kaon-") {
427  pmanager->AddDiscreteProcess(theElasticProcess);
428  G4KaonMinusInelasticProcess* theInelasticProcess =
429  new G4KaonMinusInelasticProcess("inelastic");
430  theInelasticProcess->RegisterMe(bertini);
431  theInelasticProcess->RegisterMe(theTheoModel);
432  pmanager->AddDiscreteProcess(theInelasticProcess);
433  }
434  else if (particleName == "proton") {
435  pmanager->AddDiscreteProcess(theElasticProcess);
436  G4ProtonInelasticProcess* theInelasticProcess =
437  new G4ProtonInelasticProcess("inelastic");
438  theInelasticProcess->RegisterMe(bertini);
439  theInelasticProcess->RegisterMe(theTheoModel);
440  pmanager->AddDiscreteProcess(theInelasticProcess);
441  }
442  else if (particleName == "anti_proton") {
443  pmanager->AddDiscreteProcess(theElasticProcess);
444  G4AntiProtonInelasticProcess* theInelasticProcess =
445  new G4AntiProtonInelasticProcess("inelastic");
446  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
447  pmanager->AddDiscreteProcess(theInelasticProcess);
448 
449  } else if (particleName == "neutron") {
450  // elastic scattering
451  pmanager->AddDiscreteProcess(theElasticProcess);
452 
453  // inelastic scattering
454  G4NeutronInelasticProcess* theInelasticProcess =
455  new G4NeutronInelasticProcess("inelastic");
456  theInelasticProcess->RegisterMe(bertini);
457  theInelasticProcess->RegisterMe(theTheoModel);
458  pmanager->AddDiscreteProcess(theInelasticProcess);
459 
460  // fission
461  G4HadronFissionProcess* theFissionProcess = new G4HadronFissionProcess;
462  G4LFission* theFissionModel = new G4LFission;
463  theFissionProcess->RegisterMe(theFissionModel);
464  pmanager->AddDiscreteProcess(theFissionProcess);
465 
466  // capture
467  G4HadronCaptureProcess* theCaptureProcess = new G4HadronCaptureProcess;
468  G4NeutronRadCapture* theCaptureModel = new G4NeutronRadCapture;
469  theCaptureProcess->RegisterMe(theCaptureModel);
470  pmanager->AddDiscreteProcess(theCaptureProcess);
471 
472  } else if (particleName == "anti_neutron") {
473  pmanager->AddDiscreteProcess(theElasticProcess);
474  G4AntiNeutronInelasticProcess* theInelasticProcess =
475  new G4AntiNeutronInelasticProcess("inelastic");
476  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
477  pmanager->AddDiscreteProcess(theInelasticProcess);
478 
479  } else if (particleName == "lambda") {
480  pmanager->AddDiscreteProcess(theElasticProcess);
481  G4LambdaInelasticProcess* theInelasticProcess =
482  new G4LambdaInelasticProcess("inelastic");
483  theInelasticProcess->RegisterMe(bertini);
484  theInelasticProcess->RegisterMe(theTheoModel);
485  pmanager->AddDiscreteProcess(theInelasticProcess);
486  }
487  else if (particleName == "anti_lambda") {
488  pmanager->AddDiscreteProcess(theElasticProcess);
489  G4AntiLambdaInelasticProcess* theInelasticProcess =
490  new G4AntiLambdaInelasticProcess("inelastic");
491  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
492  pmanager->AddDiscreteProcess(theInelasticProcess);
493  }
494  else if (particleName == "sigma+") {
495  pmanager->AddDiscreteProcess(theElasticProcess);
496  G4SigmaPlusInelasticProcess* theInelasticProcess =
497  new G4SigmaPlusInelasticProcess("inelastic");
498  theInelasticProcess->RegisterMe(bertini);
499  theInelasticProcess->RegisterMe(theTheoModel);
500  pmanager->AddDiscreteProcess(theInelasticProcess);
501  }
502  else if (particleName == "sigma-") {
503  pmanager->AddDiscreteProcess(theElasticProcess);
504  G4SigmaMinusInelasticProcess* theInelasticProcess =
505  new G4SigmaMinusInelasticProcess("inelastic");
506  theInelasticProcess->RegisterMe(bertini);
507  theInelasticProcess->RegisterMe(theTheoModel);
508  pmanager->AddDiscreteProcess(theInelasticProcess);
509  }
510  else if (particleName == "anti_sigma+") {
511  pmanager->AddDiscreteProcess(theElasticProcess);
512  G4AntiSigmaPlusInelasticProcess* theInelasticProcess =
513  new G4AntiSigmaPlusInelasticProcess("inelastic");
514  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
515  pmanager->AddDiscreteProcess(theInelasticProcess);
516  }
517  else if (particleName == "anti_sigma-") {
518  pmanager->AddDiscreteProcess(theElasticProcess);
519  G4AntiSigmaMinusInelasticProcess* theInelasticProcess =
520  new G4AntiSigmaMinusInelasticProcess("inelastic");
521  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
522  pmanager->AddDiscreteProcess(theInelasticProcess);
523  }
524  else if (particleName == "xi0") {
525  pmanager->AddDiscreteProcess(theElasticProcess);
526  G4XiZeroInelasticProcess* theInelasticProcess =
527  new G4XiZeroInelasticProcess("inelastic");
528  theInelasticProcess->RegisterMe(bertini);
529  theInelasticProcess->RegisterMe(theTheoModel);
530  pmanager->AddDiscreteProcess(theInelasticProcess);
531  }
532  else if (particleName == "xi-") {
533  pmanager->AddDiscreteProcess(theElasticProcess);
534  G4XiMinusInelasticProcess* theInelasticProcess =
535  new G4XiMinusInelasticProcess("inelastic");
536  theInelasticProcess->RegisterMe(bertini);
537  theInelasticProcess->RegisterMe(theTheoModel);
538  pmanager->AddDiscreteProcess(theInelasticProcess);
539  }
540  else if (particleName == "anti_xi0") {
541  pmanager->AddDiscreteProcess(theElasticProcess);
542  G4AntiXiZeroInelasticProcess* theInelasticProcess =
543  new G4AntiXiZeroInelasticProcess("inelastic");
544  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
545  pmanager->AddDiscreteProcess(theInelasticProcess);
546  }
547  else if (particleName == "anti_xi-") {
548  pmanager->AddDiscreteProcess(theElasticProcess);
549  G4AntiXiMinusInelasticProcess* theInelasticProcess =
550  new G4AntiXiMinusInelasticProcess("inelastic");
551  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
552  pmanager->AddDiscreteProcess(theInelasticProcess);
553  }
554  else if (particleName == "deuteron") {
555  pmanager->AddDiscreteProcess(theElasticProcess);
556  G4DeuteronInelasticProcess* theInelasticProcess =
557  new G4DeuteronInelasticProcess("inelastic");
558  theInelasticProcess->RegisterMe(binaryCascade);
559  theInelasticProcess->RegisterMe(qmd);
560  theInelasticProcess->AddDataSet(shenXS);
561  theInelasticProcess->AddDataSet(tripXS);
562  theInelasticProcess->AddDataSet(tripLightXS);
563  pmanager->AddDiscreteProcess(theInelasticProcess);
564  }
565  else if (particleName == "triton") {
566  pmanager->AddDiscreteProcess(theElasticProcess);
567  G4TritonInelasticProcess* theInelasticProcess =
568  new G4TritonInelasticProcess("inelastic");
569  theInelasticProcess->RegisterMe(binaryCascade);
570  theInelasticProcess->RegisterMe(qmd);
571  theInelasticProcess->AddDataSet(shenXS);
572  theInelasticProcess->AddDataSet(tripXS);
573  theInelasticProcess->AddDataSet(tripLightXS);
574  pmanager->AddDiscreteProcess(theInelasticProcess);
575  }
576  else if (particleName == "alpha") {
577  pmanager->AddDiscreteProcess(theElasticProcess);
578  G4AlphaInelasticProcess* theInelasticProcess =
579  new G4AlphaInelasticProcess("inelastic");
580  theInelasticProcess->RegisterMe(binaryCascade);
581  theInelasticProcess->RegisterMe(qmd);
582  theInelasticProcess->AddDataSet(shenXS);
583  theInelasticProcess->AddDataSet(tripXS);
584  theInelasticProcess->AddDataSet(tripLightXS);
585  pmanager->AddDiscreteProcess(theInelasticProcess);
586 
587  } else if (particleName == "omega-") {
588  pmanager->AddDiscreteProcess(theElasticProcess);
589  G4OmegaMinusInelasticProcess* theInelasticProcess =
590  new G4OmegaMinusInelasticProcess("inelastic");
591  theInelasticProcess->RegisterMe(bertini);
592  theInelasticProcess->RegisterMe(theTheoModel);
593  pmanager->AddDiscreteProcess(theInelasticProcess);
594 
595  } else if (particleName == "anti_omega-") {
596  pmanager->AddDiscreteProcess(theElasticProcess);
597  G4AntiOmegaMinusInelasticProcess* theInelasticProcess =
598  new G4AntiOmegaMinusInelasticProcess("inelastic");
599  theInelasticProcess->RegisterMe(antiBHighEnergyModel);
600  pmanager->AddDiscreteProcess(theInelasticProcess);
601  }
602  }
603 }
static const double MeV
Definition: G4SIunits.hh:211
void SetFragmentationModel(G4VStringFragmentation *aModel)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void SetMinEForMultiFrag(G4double anE)
G4ProcessManager * GetProcessManager() const
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
void RegisterMe(G4HadronicInteraction *a)
void SetMinEnergy(G4double anEnergy)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
const G4String & GetParticleName() const
void SetFermiModel(G4VFermiBreakUp *ptr)
void SetMultiFragmentation(G4VMultiFragmentation *ptr)
static const double GeV
Definition: G4SIunits.hh:214
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetMaxAandZForFermiBreakUp(G4int anA, G4int aZ)
void SetEvaporation(G4VEvaporation *ptr)
void SetMaxEnergy(const G4double anEnergy)
void SetDeExcitation(G4VPreCompoundModel *ptr)
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
void SetTransport(G4VIntraNuclearTransportModel *const value)
static const double TeV
Definition: G4SIunits.hh:215
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructLeptHad()

void B03PhysicsList::ConstructLeptHad ( )
protectedvirtual

Definition at line 607 of file B03PhysicsList.cc.

608 {;}
Here is the caller graph for this function:

◆ ConstructParticle()

void B03PhysicsList::ConstructParticle ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 71 of file B03PhysicsList.cc.

72 {
73  // In this method, static member functions should be called
74  // for all particles which you want to use.
75  // This ensures that objects of these particle types will be
76  // created in the program.
77 
84 }
void ConstructAllBosons()
void ConstructAllLeptons()
void ConstructAllShortLiveds()
void ConstructAllBaryons()
void ConstructAllMesons()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructProcess()

void B03PhysicsList::ConstructProcess ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 142 of file B03PhysicsList.cc.

143 {
147  ConstructEM();
149  ConstructHad();
151 }
void AddScoringProcess()
void AddBiasingProcess()
virtual void ConstructGeneral()
virtual void ConstructLeptHad()
virtual void ConstructHad()
virtual void ConstructEM()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetCuts()

void B03PhysicsList::SetCuts ( )
protectedvirtual

Reimplemented from G4VUserPhysicsList.

Definition at line 631 of file B03PhysicsList.cc.

632 {
633  if (verboseLevel >0)
634  {
635  G4cout << "B03PhysicsList::SetCuts:";
636  G4cout << "CutLength : " << defaultCutValue/mm << " (mm)" << G4endl;
637  }
638  // "G4VUserPhysicsList::SetCutsWithDefault" method sets
639  // the default cut value for all particle types
641 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fBiasWorldName

G4String B03PhysicsList::fBiasWorldName
private

Definition at line 82 of file B03PhysicsList.hh.

◆ fGeomSampler

G4GeometrySampler* B03PhysicsList::fGeomSampler
private

Definition at line 83 of file B03PhysicsList.hh.

◆ fParaWorldName

std::vector<G4String> B03PhysicsList::fParaWorldName
private

Definition at line 81 of file B03PhysicsList.hh.


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