Geant4  10.02.p03
OpNovicePhysicsList Class Reference

#include <OpNovicePhysicsList.hh>

Inheritance diagram for OpNovicePhysicsList:
Collaboration diagram for OpNovicePhysicsList:

Public Member Functions

 OpNovicePhysicsList ()
 
virtual ~OpNovicePhysicsList ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
virtual void SetCuts ()
 
void ConstructDecay ()
 
void ConstructEM ()
 
void ConstructOp ()
 
void SetVerbose (G4int)
 
void SetNbOfPhotonsCerenkov (G4int)
 
- 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 ()
 

Private Attributes

OpNovicePhysicsListMessengerfMessenger
 

Static Private Attributes

static G4ThreadLocal G4int fVerboseLevel = 1
 
static G4ThreadLocal G4int fMaxNumPhotonStep = 20
 
static G4ThreadLocal G4CerenkovfCerenkovProcess = 0
 
static G4ThreadLocal G4ScintillationfScintillationProcess = 0
 
static G4ThreadLocal G4OpAbsorptionfAbsorptionProcess = 0
 
static G4ThreadLocal G4OpRayleighfRayleighScatteringProcess = 0
 
static G4ThreadLocal G4OpMieHGfMieHGScatteringProcess = 0
 
static G4ThreadLocal G4OpBoundaryProcessfBoundaryProcess = 0
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VUserPhysicsList
static const G4VUPLManagerGetSubInstanceManager ()
 
- 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
 
- 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 48 of file OpNovicePhysicsList.hh.

Constructor & Destructor Documentation

◆ OpNovicePhysicsList()

OpNovicePhysicsList::OpNovicePhysicsList ( )

◆ ~OpNovicePhysicsList()

OpNovicePhysicsList::~OpNovicePhysicsList ( )
virtual

Definition at line 77 of file OpNovicePhysicsList.cc.

77 { delete fMessenger; }
OpNovicePhysicsListMessenger * fMessenger

Member Function Documentation

◆ ConstructDecay()

void OpNovicePhysicsList::ConstructDecay ( )

Definition at line 120 of file OpNovicePhysicsList.cc.

121 {
122  // Add Decay Process
123  G4Decay* theDecayProcess = new G4Decay();
125  particleIterator->reset();
126  while( (*particleIterator)() ){
127  G4ParticleDefinition* particle = particleIterator->value();
128  G4ProcessManager* pmanager = particle->GetProcessManager();
129  if (theDecayProcess->IsApplicable(*particle)) {
130  pmanager ->AddProcess(theDecayProcess);
131  // set ordering for PostStepDoIt and AtRestDoIt
132  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
133  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
134  }
135  }
136 }
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:

◆ ConstructEM()

void OpNovicePhysicsList::ConstructEM ( )

Definition at line 160 of file OpNovicePhysicsList.cc.

161 {
163  particleIterator->reset();
164  while( (*particleIterator)() ){
165  G4ParticleDefinition* particle = particleIterator->value();
166  G4ProcessManager* pmanager = particle->GetProcessManager();
167  G4String particleName = particle->GetParticleName();
168 
169  if (particleName == "gamma") {
170  // gamma
171  // Construct processes for gamma
172  pmanager->AddDiscreteProcess(new G4GammaConversion());
173  pmanager->AddDiscreteProcess(new G4ComptonScattering());
174  pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
175 
176  } else if (particleName == "e-") {
177  //electron
178  // Construct processes for electron
179  pmanager->AddProcess(new G4eMultipleScattering(),-1, 1, 1);
180  pmanager->AddProcess(new G4eIonisation(), -1, 2, 2);
181  pmanager->AddProcess(new G4eBremsstrahlung(), -1, 3, 3);
182 
183  } else if (particleName == "e+") {
184  //positron
185  // Construct processes for positron
186  pmanager->AddProcess(new G4eMultipleScattering(),-1, 1, 1);
187  pmanager->AddProcess(new G4eIonisation(), -1, 2, 2);
188  pmanager->AddProcess(new G4eBremsstrahlung(), -1, 3, 3);
189  pmanager->AddProcess(new G4eplusAnnihilation(), 0,-1, 4);
190 
191  } else if( particleName == "mu+" ||
192  particleName == "mu-" ) {
193  //muon
194  // Construct processes for muon
195  pmanager->AddProcess(new G4MuMultipleScattering(),-1, 1, 1);
196  pmanager->AddProcess(new G4MuIonisation(), -1, 2, 2);
197  pmanager->AddProcess(new G4MuBremsstrahlung(), -1, 3, 3);
198  pmanager->AddProcess(new G4MuPairProduction(), -1, 4, 4);
199 
200  } else {
201  if ((particle->GetPDGCharge() != 0.0) &&
202  (particle->GetParticleName() != "chargedgeantino") &&
203  !particle->IsShortLived()) {
204  // all others charged particles except geantino
205  pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
206  pmanager->AddProcess(new G4hIonisation(), -1,2,2);
207  }
208  }
209  }
210 }
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:

◆ ConstructOp()

void OpNovicePhysicsList::ConstructOp ( )

Definition at line 215 of file OpNovicePhysicsList.cc.

216 {
217  fCerenkovProcess = new G4Cerenkov("Cerenkov");
221  fScintillationProcess = new G4Scintillation("Scintillation");
228 
235 
236  // Use Birks Correction in the Scintillation process
238  {
239  G4EmSaturation* emSaturation =
241  fScintillationProcess->AddSaturation(emSaturation);
242  }
243 
245  particleIterator->reset();
246  while( (*particleIterator)() ){
247  G4ParticleDefinition* particle = particleIterator->value();
248  G4ProcessManager* pmanager = particle->GetProcessManager();
249  G4String particleName = particle->GetParticleName();
250  if (fCerenkovProcess->IsApplicable(*particle)) {
251  pmanager->AddProcess(fCerenkovProcess);
253  }
254  if (fScintillationProcess->IsApplicable(*particle)) {
258  }
259  if (particleName == "opticalphoton") {
260  G4cout << " AddDiscreteProcess to OpticalPhoton " << G4endl;
265  }
266  }
267 }
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:151
static G4LossTableManager * Instance()
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:146
static G4ThreadLocal G4Cerenkov * fCerenkovProcess
static G4ThreadLocal G4int fMaxNumPhotonStep
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
G4ProcessManager * GetProcessManager() const
static G4ThreadLocal G4OpRayleigh * fRayleighScatteringProcess
void SetScintillationYieldFactor(const G4double yieldfactor)
const G4String & GetParticleName() const
static G4ThreadLocal G4OpBoundaryProcess * fBoundaryProcess
G4GLOB_DLL std::ostream G4cout
void AddSaturation(G4EmSaturation *)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
G4EmSaturation * EmSaturation()
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:156
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
static G4ThreadLocal G4OpMieHG * fMieHGScatteringProcess
static G4ThreadLocal G4OpAbsorption * fAbsorptionProcess
void SetTrackSecondariesFirst(const G4bool state)
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
#define G4endl
Definition: G4ios.hh:61
G4bool IsMasterThread()
Definition: G4Threading.cc:136
static G4ThreadLocal G4int fVerboseLevel
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
static G4ThreadLocal G4Scintillation * fScintillationProcess
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Definition: G4Cerenkov.cc:135
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructParticle()

void OpNovicePhysicsList::ConstructParticle ( void  )
virtual

Implements G4VUserPhysicsList.

Definition at line 81 of file OpNovicePhysicsList.cc.

82 {
83  // In this method, static member functions should be called
84  // for all particles which you want to use.
85  // This ensures that objects of these particle types will be
86  // created in the program.
87 
88  G4BosonConstructor bConstructor;
89  bConstructor.ConstructParticle();
90 
91  G4LeptonConstructor lConstructor;
92  lConstructor.ConstructParticle();
93 
94  G4MesonConstructor mConstructor;
95  mConstructor.ConstructParticle();
96 
97  G4BaryonConstructor rConstructor;
98  rConstructor.ConstructParticle();
99 
100  G4IonConstructor iConstructor;
101  iConstructor.ConstructParticle();
102 }
static void ConstructParticle()
static void ConstructParticle()
static void ConstructParticle()
static void ConstructParticle()
static void ConstructParticle()
Here is the call graph for this function:

◆ ConstructProcess()

void OpNovicePhysicsList::ConstructProcess ( void  )
virtual

Implements G4VUserPhysicsList.

Definition at line 106 of file OpNovicePhysicsList.cc.

Here is the call graph for this function:

◆ SetCuts()

void OpNovicePhysicsList::SetCuts ( )
virtual

Reimplemented from G4VUserPhysicsList.

Definition at line 294 of file OpNovicePhysicsList.cc.

295 {
296  // " G4VUserPhysicsList::SetCutsWithDefault" method sets
297  // the default cut value for all particle types
298  //
300 
302 }
void DumpCutValuesTable(G4int flag=1)
Here is the call graph for this function:

◆ SetNbOfPhotonsCerenkov()

void OpNovicePhysicsList::SetNbOfPhotonsCerenkov ( G4int  MaxNumber)

Definition at line 285 of file OpNovicePhysicsList.cc.

286 {
287  fMaxNumPhotonStep = MaxNumber;
288 
290 }
static G4ThreadLocal G4Cerenkov * fCerenkovProcess
static G4ThreadLocal G4int fMaxNumPhotonStep
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:156
Here is the call graph for this function:

◆ SetVerbose()

void OpNovicePhysicsList::SetVerbose ( G4int  verbose)

Definition at line 271 of file OpNovicePhysicsList.cc.

272 {
273  fVerboseLevel = verbose;
274 
281 }
static G4ThreadLocal G4Cerenkov * fCerenkovProcess
static G4ThreadLocal G4OpRayleigh * fRayleighScatteringProcess
static G4ThreadLocal G4OpBoundaryProcess * fBoundaryProcess
static G4ThreadLocal G4OpMieHG * fMieHGScatteringProcess
static G4ThreadLocal G4OpAbsorption * fAbsorptionProcess
static G4ThreadLocal G4int fVerboseLevel
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
static G4ThreadLocal G4Scintillation * fScintillationProcess
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fAbsorptionProcess

G4ThreadLocal G4OpAbsorption * OpNovicePhysicsList::fAbsorptionProcess = 0
staticprivate

Definition at line 80 of file OpNovicePhysicsList.hh.

◆ fBoundaryProcess

G4ThreadLocal G4OpBoundaryProcess * OpNovicePhysicsList::fBoundaryProcess = 0
staticprivate

Definition at line 83 of file OpNovicePhysicsList.hh.

◆ fCerenkovProcess

G4ThreadLocal G4Cerenkov * OpNovicePhysicsList::fCerenkovProcess = 0
staticprivate

Definition at line 78 of file OpNovicePhysicsList.hh.

◆ fMaxNumPhotonStep

G4ThreadLocal G4int OpNovicePhysicsList::fMaxNumPhotonStep = 20
staticprivate

Definition at line 76 of file OpNovicePhysicsList.hh.

◆ fMessenger

OpNovicePhysicsListMessenger* OpNovicePhysicsList::fMessenger
private

Definition at line 73 of file OpNovicePhysicsList.hh.

◆ fMieHGScatteringProcess

G4ThreadLocal G4OpMieHG * OpNovicePhysicsList::fMieHGScatteringProcess = 0
staticprivate

Definition at line 82 of file OpNovicePhysicsList.hh.

◆ fRayleighScatteringProcess

G4ThreadLocal G4OpRayleigh * OpNovicePhysicsList::fRayleighScatteringProcess = 0
staticprivate

Definition at line 81 of file OpNovicePhysicsList.hh.

◆ fScintillationProcess

G4ThreadLocal G4Scintillation * OpNovicePhysicsList::fScintillationProcess = 0
staticprivate

Definition at line 79 of file OpNovicePhysicsList.hh.

◆ fVerboseLevel

G4ThreadLocal G4int OpNovicePhysicsList::fVerboseLevel = 1
staticprivate

Definition at line 75 of file OpNovicePhysicsList.hh.


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