Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 ()
 

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 51 of file OpNovicePhysicsList.hh.

Constructor & Destructor Documentation

OpNovicePhysicsList::OpNovicePhysicsList ( )
OpNovicePhysicsList::~OpNovicePhysicsList ( )
virtual

Definition at line 80 of file OpNovicePhysicsList.cc.

80 { delete fMessenger; }

Member Function Documentation

void OpNovicePhysicsList::ConstructDecay ( )

Definition at line 123 of file OpNovicePhysicsList.cc.

124 {
125  // Add Decay Process
126  G4Decay* theDecayProcess = new G4Decay();
128  particleIterator->reset();
129  while( (*particleIterator)() ){
130  G4ParticleDefinition* particle = particleIterator->value();
131  G4ProcessManager* pmanager = particle->GetProcessManager();
132  if (theDecayProcess->IsApplicable(*particle)) {
133  pmanager ->AddProcess(theDecayProcess);
134  // set ordering for PostStepDoIt and AtRestDoIt
135  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
136  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
137  }
138  }
139 }
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
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 OpNovicePhysicsList::ConstructEM ( )

Definition at line 163 of file OpNovicePhysicsList.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void OpNovicePhysicsList::ConstructOp ( )

Definition at line 218 of file OpNovicePhysicsList.cc.

219 {
220  fCerenkovProcess = new G4Cerenkov("Cerenkov");
221  fCerenkovProcess->SetMaxNumPhotonsPerStep(fMaxNumPhotonStep);
222  fCerenkovProcess->SetMaxBetaChangePerStep(10.0);
223  fCerenkovProcess->SetTrackSecondariesFirst(true);
224  fScintillationProcess = new G4Scintillation("Scintillation");
225  fScintillationProcess->SetScintillationYieldFactor(1.);
226  fScintillationProcess->SetTrackSecondariesFirst(true);
227  fAbsorptionProcess = new G4OpAbsorption();
228  fRayleighScatteringProcess = new G4OpRayleigh();
229  fMieHGScatteringProcess = new G4OpMieHG();
230  fBoundaryProcess = new G4OpBoundaryProcess();
231 
232  fCerenkovProcess->SetVerboseLevel(fVerboseLevel);
233  fScintillationProcess->SetVerboseLevel(fVerboseLevel);
234  fAbsorptionProcess->SetVerboseLevel(fVerboseLevel);
235  fRayleighScatteringProcess->SetVerboseLevel(fVerboseLevel);
236  fMieHGScatteringProcess->SetVerboseLevel(fVerboseLevel);
237  fBoundaryProcess->SetVerboseLevel(fVerboseLevel);
238 
239  // Use Birks Correction in the Scintillation process
241  {
242  G4EmSaturation* emSaturation =
244  fScintillationProcess->AddSaturation(emSaturation);
245  }
246 
248  particleIterator->reset();
249  while( (*particleIterator)() ){
250  G4ParticleDefinition* particle = particleIterator->value();
251  G4ProcessManager* pmanager = particle->GetProcessManager();
252  G4String particleName = particle->GetParticleName();
253  if (fCerenkovProcess->IsApplicable(*particle)) {
254  pmanager->AddProcess(fCerenkovProcess);
255  pmanager->SetProcessOrdering(fCerenkovProcess,idxPostStep);
256  }
257  if (fScintillationProcess->IsApplicable(*particle)) {
258  pmanager->AddProcess(fScintillationProcess);
259  pmanager->SetProcessOrderingToLast(fScintillationProcess, idxAtRest);
260  pmanager->SetProcessOrderingToLast(fScintillationProcess, idxPostStep);
261  }
262  if (particleName == "opticalphoton") {
263  G4cout << " AddDiscreteProcess to OpticalPhoton " << G4endl;
264  pmanager->AddDiscreteProcess(fAbsorptionProcess);
265  pmanager->AddDiscreteProcess(fRayleighScatteringProcess);
266  pmanager->AddDiscreteProcess(fMieHGScatteringProcess);
267  pmanager->AddDiscreteProcess(fBoundaryProcess);
268  }
269  }
270 }
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:150
static G4LossTableManager * Instance()
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:145
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
const G4String & GetParticleName() const
void SetScintillationYieldFactor(const G4double yieldfactor)
G4GLOB_DLL std::ostream G4cout
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
G4EmSaturation * EmSaturation()
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:155
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
void SetTrackSecondariesFirst(const G4bool state)
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
G4ProcessManager * GetProcessManager() const
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:133
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
#define G4endl
Definition: G4ios.hh:61
G4bool IsMasterThread()
Definition: G4Threading.cc:146
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
void AddSaturation(G4EmSaturation *sat)
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 OpNovicePhysicsList::ConstructParticle ( void  )
virtual

Implements G4VUserPhysicsList.

Definition at line 84 of file OpNovicePhysicsList.cc.

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

Here is the call graph for this function:

void OpNovicePhysicsList::ConstructProcess ( void  )
virtual

Implements G4VUserPhysicsList.

Definition at line 109 of file OpNovicePhysicsList.cc.

Here is the call graph for this function:

void OpNovicePhysicsList::SetCuts ( )
virtual

Reimplemented from G4VUserPhysicsList.

Definition at line 297 of file OpNovicePhysicsList.cc.

298 {
299  // " G4VUserPhysicsList::SetCutsWithDefault" method sets
300  // the default cut value for all particle types
301  //
303 
305 }
void DumpCutValuesTable(G4int flag=1)

Here is the call graph for this function:

void OpNovicePhysicsList::SetNbOfPhotonsCerenkov ( G4int  MaxNumber)

Definition at line 288 of file OpNovicePhysicsList.cc.

289 {
290  fMaxNumPhotonStep = MaxNumber;
291 
292  fCerenkovProcess->SetMaxNumPhotonsPerStep(fMaxNumPhotonStep);
293 }
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:155

Here is the call graph for this function:

void OpNovicePhysicsList::SetVerbose ( G4int  verbose)

Definition at line 274 of file OpNovicePhysicsList.cc.

275 {
276  fVerboseLevel = verbose;
277 
278  fCerenkovProcess->SetVerboseLevel(fVerboseLevel);
279  fScintillationProcess->SetVerboseLevel(fVerboseLevel);
280  fAbsorptionProcess->SetVerboseLevel(fVerboseLevel);
281  fRayleighScatteringProcess->SetVerboseLevel(fVerboseLevel);
282  fMieHGScatteringProcess->SetVerboseLevel(fVerboseLevel);
283  fBoundaryProcess->SetVerboseLevel(fVerboseLevel);
284 }
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437

Here is the call graph for this function:

Here is the caller graph for this function:


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