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

#include <G4ErrorPhysicsList.hh>

+ Inheritance diagram for G4ErrorPhysicsList:

Public Member Functions

 G4ErrorPhysicsList ()
 
virtual ~G4ErrorPhysicsList ()
 
- Public Member Functions inherited from G4VUserPhysicsList
 G4VUserPhysicsList ()
 
virtual ~G4VUserPhysicsList ()
 
 G4VUserPhysicsList (const G4VUserPhysicsList &)
 
G4VUserPhysicsListoperator= (const G4VUserPhysicsList &)
 
virtual void ConstructParticle ()=0
 
void Construct ()
 
virtual void ConstructProcess ()=0
 
void UseCoupledTransportation (G4bool vl=true)
 
virtual void SetCuts ()
 
void SetDefaultCutValue (G4double newCutValue)
 
G4double GetDefaultCutValue () const
 
void BuildPhysicsTable ()
 
void PreparePhysicsTable (G4ParticleDefinition *)
 
void BuildPhysicsTable (G4ParticleDefinition *)
 
G4bool StorePhysicsTable (const G4String &directory=".")
 
G4bool IsPhysicsTableRetrieved () const
 
G4bool IsStoredInAscii () const
 
const G4StringGetPhysicsTableDirectory () const
 
void SetPhysicsTableRetrieved (const G4String &directory="")
 
void SetStoredInAscii ()
 
void ResetPhysicsTableRetrieved ()
 
void ResetStoredInAscii ()
 
void DumpList () const
 
void DumpCutValuesTable (G4int flag=1)
 
void DumpCutValuesTableIfRequested ()
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
void 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
 
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 ()
 

Protected Member Functions

virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
virtual void SetCuts ()
 
virtual void ConstructEM ()
 
- 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 ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VUserPhysicsList
G4ParticleTabletheParticleTable
 
G4ParticleTable::G4PTblDicIteratortheParticleIterator
 
G4UserPhysicsListMessengertheMessenger
 
G4int verboseLevel
 
G4double defaultCutValue
 
G4bool isSetDefaultCutValue
 
G4ProductionCutsTablefCutsTable
 
G4bool fRetrievePhysicsTable
 
G4bool fStoredInAscii
 
G4bool fIsCheckedForRetrievePhysicsTable
 
G4bool fIsRestoredCutValues
 
G4String directoryPhysicsTable
 
G4int fDisplayThreshold
 
G4bool fIsPhysicsTableBuilt
 
G4bool fDisableCheckParticleList
 
G4PhysicsListHelperthePLHelper
 

Detailed Description

Definition at line 49 of file G4ErrorPhysicsList.hh.

Constructor & Destructor Documentation

◆ G4ErrorPhysicsList()

G4ErrorPhysicsList::G4ErrorPhysicsList ( )

Definition at line 72 of file G4ErrorPhysicsList.cc.

73{
74 defaultCutValue = 1.0E+9*cm; // set big step so that AlongStep computes all the energy
75}

◆ ~G4ErrorPhysicsList()

G4ErrorPhysicsList::~G4ErrorPhysicsList ( )
virtual

Definition at line 79 of file G4ErrorPhysicsList.cc.

80{
81}

Member Function Documentation

◆ ConstructEM()

void G4ErrorPhysicsList::ConstructEM ( )
protectedvirtual

Definition at line 169 of file G4ErrorPhysicsList.cc.

170{
171
172 G4ErrorEnergyLoss* eLossProcess = new G4ErrorEnergyLoss;
175 new G4ErrorMessenger( stepLengthLimitProcess, magFieldLimitProcess, eLossProcess );
176
178 while( (*theParticleIterator)() ){
180 G4ProcessManager* pmanager = particle->GetProcessManager();
181 G4String particleName = particle->GetParticleName();
182
183 if (particleName == "gamma") {
184 // gamma
185 pmanager->AddDiscreteProcess(new G4GammaConversion());
186 pmanager->AddDiscreteProcess(new G4ComptonScattering());
188
189 // } else if (particleName == "e-" || particleName == "e+"
190 // || particleName == "mu+" || particleName == "mu-" ) {
191 }else if (!particle->IsShortLived() && particle->GetPDGCharge() != 0 ) {
192
193 pmanager->AddContinuousProcess(eLossProcess,1);
194 pmanager->AddDiscreteProcess( stepLengthLimitProcess, 2 );
195 pmanager->AddDiscreteProcess( magFieldLimitProcess, 3 );
196
197 /* } else if ((!particle->IsShortLived()) &&
198 (particle->GetPDGCharge() != 0.0) &&
199 (particle->GetParticleName() != "chargedgeantino")) {
200 // all others charged particles except geantino
201 // G4VProcess* aMultipleScattering = new G4MultipleScattering();
202 G4VProcess* anIonisation = new G4hIonisation();
203 ////G4VProcess* theUserCuts = new G4UserSpecialCuts();
204
205 //
206 // add processes
207 pmanager->AddProcess(anIonisation);
208 // pmanager->AddProcess(aMultipleScattering);
209 ////pmanager->AddProcess(theUserCuts);
210
211 //
212 // set ordering for AlongStepDoIt
213 // pmanager->SetProcessOrdering(aMultipleScattering, idxAlongStep,1);
214 pmanager->SetProcessOrdering(anIonisation, idxAlongStep,1);
215
216 //
217 // set ordering for PostStepDoIt
218 // pmanager->SetProcessOrdering(aMultipleScattering, idxPostStep,1);
219 pmanager->SetProcessOrdering(anIonisation, idxPostStep,1);
220 ////pmanager->SetProcessOrdering(theUserCuts, idxPostStep,2);
221 */
222 }
223 }
224}
G4ProcessManager * GetProcessManager() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4int AddContinuousProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ParticleTable::G4PTblDicIterator * theParticleIterator

Referenced by ConstructProcess().

◆ ConstructParticle()

void G4ErrorPhysicsList::ConstructParticle ( )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 85 of file G4ErrorPhysicsList.cc.

86{
87// In this method, static member functions should be called
88 // for all particles which you want to use.
89 // This ensures that objects of these particle types will be
90 // created in the program.
91 // gamma
93 // e+/-
96 // mu+/-
99
100 // pi+/-
103
104 // proton
106
107}
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88

◆ ConstructProcess()

void G4ErrorPhysicsList::ConstructProcess ( )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 111 of file G4ErrorPhysicsList.cc.

112{
113 G4Transportation* theTransportationProcess= new G4Transportation();
114
115#ifdef G4VERBOSE
116 if (verboseLevel >= 4){
117 G4cout << "G4VUserPhysicsList::ConstructProcess() "<< G4endl;
118 }
119#endif
120
121 // loop over all particles in G4ParticleTable
123 while( (*theParticleIterator)() ){
125 G4ProcessManager* pmanager = particle->GetProcessManager();
126 if (!particle->IsShortLived()) {
127 G4cout << particle << "G4ErrorPhysicsList:: particle process manager " << particle->GetParticleName() << " = " << particle->GetProcessManager() << G4endl;
128 // Add transportation process for all particles other than "shortlived"
129 if ( pmanager == 0) {
130 // Error !! no process manager
131 G4String particleName = particle->GetParticleName();
132 G4Exception("G4ErrorPhysicsList::ConstructProcess","No process manager",
133 RunMustBeAborted, particleName );
134 } else {
135 // add transportation with ordering = ( -1, "first", "first" )
136 pmanager ->AddProcess(theTransportationProcess);
137 pmanager ->SetProcessOrderingToFirst(theTransportationProcess, idxAlongStep);
138 pmanager ->SetProcessOrderingToFirst(theTransportationProcess, idxPostStep);
139 }
140 } else {
141 // shortlived particle case
142 }
143 }
144
145 ConstructEM();
146}
@ RunMustBeAborted
@ idxPostStep
@ idxAlongStep
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
virtual void ConstructEM()
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetProcessOrderingToFirst(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ SetCuts()

void G4ErrorPhysicsList::SetCuts ( )
protectedvirtual

Reimplemented from G4VUserPhysicsList.

Definition at line 228 of file G4ErrorPhysicsList.cc.

229{
230 // " G4VUserPhysicsList::SetCutsWithDefault" method sets
231 // the default cut value or all particle types
233 // if (verboseLevel>0)
234 // DumpCutValuesTable();
235}

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