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

#include <PurgMagPhysicsList.hh>

Inheritance diagram for PurgMagPhysicsList:
Collaboration diagram for PurgMagPhysicsList:

Public Member Functions

 PurgMagPhysicsList ()
 
 ~PurgMagPhysicsList ()
 
void SetGammaCut (G4double)
 
void SetElectronCut (G4double)
 
void SetPositronCut (G4double)
 
void SetProtonCut (G4double)
 
void SetGammaLowLimit (G4double)
 
void SetElectronLowLimit (G4double)
 
void SetPositronLowLimit (G4double)
 
void SetProtonLowLimit (G4double)
 
void SetGEPLowLimit (G4double)
 
void SetGELowLimit (G4double)
 
- 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

void ConstructParticle ()
 
void ConstructProcess ()
 
void SetCuts ()
 
void ConstructBosons ()
 
void ConstructLeptons ()
 
void ConstructBarions ()
 
void ConstructGeneral ()
 
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 ()
 
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 PurgMagPhysicsList.hh.

Constructor & Destructor Documentation

PurgMagPhysicsList::PurgMagPhysicsList ( )

Definition at line 54 of file PurgMagPhysicsList.cc.

55 {
57  cutForGamma = defaultCutValue;
58  cutForElectron = defaultCutValue;
59  cutForPositron = defaultCutValue;
60  cutForProton = defaultCutValue;
61 
62  SetVerboseLevel(1);
63 }
void SetVerboseLevel(G4int value)
static constexpr double micrometer
Definition: G4SIunits.hh:100

Here is the call graph for this function:

PurgMagPhysicsList::~PurgMagPhysicsList ( )

Definition at line 67 of file PurgMagPhysicsList.cc.

68 {}

Member Function Documentation

void PurgMagPhysicsList::ConstructBarions ( )
protected

Definition at line 106 of file PurgMagPhysicsList.cc.

107 {
108  // barions
111 }
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
static G4AntiProton * AntiProtonDefinition()
Definition: G4AntiProton.cc:88

Here is the call graph for this function:

Here is the caller graph for this function:

void PurgMagPhysicsList::ConstructBosons ( )
protected

Definition at line 86 of file PurgMagPhysicsList.cc.

87 {
88 
89  // gamma
91 
92  // optical photon
94 }
static G4OpticalPhoton * OpticalPhotonDefinition()
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81

Here is the call graph for this function:

Here is the caller graph for this function:

void PurgMagPhysicsList::ConstructEM ( )
protected

Definition at line 137 of file PurgMagPhysicsList.cc.

138 {
140  particleIterator->reset();
141  while( (*particleIterator)() )
142  {
143  G4ParticleDefinition* particle = particleIterator->value();
144  G4ProcessManager* pmanager = particle->GetProcessManager();
145  G4String particleName = particle->GetParticleName();
146 
147 
148  if (particleName == "gamma") {
149  //gamma
150  pmanager->AddDiscreteProcess(new G4GammaConversion);
153 
154  } else if (particleName == "e-") {
155  //electron
156  pmanager->AddProcess(new G4eBremsstrahlung, -1,-1,3);
157  pmanager->AddProcess(new G4eIonisation, -1, 2,2);
158  pmanager->AddProcess(new G4eMultipleScattering, -1, 1,1);
159 
160  } else if (particleName == "e+") {
161  //positron
162  pmanager->AddProcess(new G4eBremsstrahlung, -1,-1,3);
163  pmanager->AddProcess(new G4eIonisation, -1, 2,2);
164  pmanager->AddProcess(new G4eMultipleScattering, -1, 1,1);
165  pmanager->AddProcess(new G4eplusAnnihilation, 0,-1,4);
166 
167  } else if (particleName == "proton") {
168  //proton
169  pmanager->AddProcess(new G4hMultipleScattering, -1,1,1);
170 
171  } else if (particleName == "anti_proton") {
172  //antiproton
173  pmanager->AddProcess(new G4hMultipleScattering, -1,1,1);
174  }
175  }
176 }
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 PurgMagPhysicsList::ConstructGeneral ( )
protected

Definition at line 180 of file PurgMagPhysicsList.cc.

181 { }

Here is the caller graph for this function:

void PurgMagPhysicsList::ConstructLeptons ( )
protected

Definition at line 97 of file PurgMagPhysicsList.cc.

98 {
99  // leptons
102 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89

Here is the call graph for this function:

Here is the caller graph for this function:

void PurgMagPhysicsList::ConstructParticle ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 72 of file PurgMagPhysicsList.cc.

73 {
74  // In this method, static member functions should be called
75  // for all particles which you want to use.
76  // This ensures that objects of these particle types will be
77  // created in the program.
78 
82 }

Here is the call graph for this function:

void PurgMagPhysicsList::ConstructProcess ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 115 of file PurgMagPhysicsList.cc.

116 {
118  ConstructEM();
120 }

Here is the call graph for this function:

void PurgMagPhysicsList::SetCuts ( )
protectedvirtual

Reimplemented from G4VUserPhysicsList.

Definition at line 185 of file PurgMagPhysicsList.cc.

186 {
187  if (verboseLevel >0){
188  G4cout << "PurgMagPhysicsList::SetCuts:";
189  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
190  }
191 
192  // set cut values for gamma at first and for e- second and next for e+,
193  // because some processes for e+/e- need cut values for gamma
194  SetCutValue(cutForGamma, "gamma");
195  SetCutValue(cutForElectron, "e-");
196  SetCutValue(cutForPositron, "e+");
197 
198  // set cut values for proton and anti_proton before all other hadrons
199  // because some processes for hadrons need cut values for proton/anti_proton
200  SetCutValue(cutForProton, "proton");
201  SetCutValue(cutForProton, "anti_proton");
202 
203  // SetCutValueForOthers(defaultCutValue);
204 
206 }
void SetCutValue(G4double aCut, const G4String &pname)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void PurgMagPhysicsList::SetElectronCut ( G4double  val)

Definition at line 306 of file PurgMagPhysicsList.cc.

307 {
308  // ResetCuts();
309  cutForElectron = val;
310 }
void PurgMagPhysicsList::SetElectronLowLimit ( G4double  lowcut)

Definition at line 223 of file PurgMagPhysicsList.cc.

224 {
225  if (verboseLevel >0){
226 
227  G4cout << "PurgMagPhysicsList::SetCuts:";
228  G4cout << "Electron cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
229  }
230 
231  // G4Electron::SetEnergyRange(lowcut,1e5);
232  SetGELowLimit(lowcut);
233 }
void SetGELowLimit(G4double)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the call graph for this function:

void PurgMagPhysicsList::SetGammaCut ( G4double  val)

Definition at line 298 of file PurgMagPhysicsList.cc.

299 {
300  ResetCuts();
301  cutForGamma = val;
302 }
void ResetCuts()
obsolete methods

Here is the call graph for this function:

void PurgMagPhysicsList::SetGammaLowLimit ( G4double  lowcut)

Definition at line 210 of file PurgMagPhysicsList.cc.

211 {
212  if (verboseLevel >0){
213  G4cout << "PurgMagPhysicsList::SetCuts:";
214  G4cout << "Gamma cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
215  }
216 
217  // G4Gamma::SetEnergyRange(lowcut,1e5);
218  SetGELowLimit(lowcut);
219 }
void SetGELowLimit(G4double)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the call graph for this function:

void PurgMagPhysicsList::SetGELowLimit ( G4double  lowcut)

Definition at line 289 of file PurgMagPhysicsList.cc.

290 {
291  if (verboseLevel >0){
292  G4cout << "PurgMagPhysicsList::SetGELowLimit:";
293  G4cout << "Gamma and Electron cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
294  }
295 
297 }
void SetEnergyRange(G4double lowedge, G4double highedge)
G4GLOB_DLL std::ostream G4cout
static G4ProductionCutsTable * GetProductionCutsTable()
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the call graph for this function:

Here is the caller graph for this function:

void PurgMagPhysicsList::SetGEPLowLimit ( G4double  lowcut)

Definition at line 272 of file PurgMagPhysicsList.cc.

273 {
274  if (verboseLevel >0){
275  G4cout << "PurgMagPhysicsList::SetGEPLowLimit:";
276  G4cout << "Gamma and Electron cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
277  }
278 
279  // G4Gamma::SetEnergyRange(lowcut,1e5);
280  // G4Electron::SetEnergyRange(lowcut,1e5);
281  // G4Positron::SetEnergyRange(lowcut,1e5);
282  this->SetGELowLimit(lowcut);
283 
284  G4cerr << " SetGEPLowLimit : Uncertain whether setting Positron low limit " << G4endl;
285 }
void SetGELowLimit(G4double)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4GLOB_DLL std::ostream G4cerr

Here is the call graph for this function:

void PurgMagPhysicsList::SetPositronCut ( G4double  val)

Definition at line 314 of file PurgMagPhysicsList.cc.

315 {
316  // ResetCuts();
317  cutForPositron = val;
318 }
void PurgMagPhysicsList::SetPositronLowLimit ( G4double  lowcut)

Definition at line 237 of file PurgMagPhysicsList.cc.

238 {
239  if (verboseLevel >0){
240 
241  G4cout << "PurgMagPhysicsList::SetCuts:";
242  G4cout << "Positron cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
243  }
244 
245  G4cerr << "PurgMagPhysicsList::SetPositronLowLimit: Not currently able to set Positron LowLimit." << G4endl;
246  G4Exception("PurgMagPhysicsList::SetPositronLowLimit()","PurMag001",
247  FatalException,"Positron Low Limit: not implemented in PurgMagPhysicsList");
248  //
249  // G4Positron::SetEnergyRange(lowcut,1e5);
250 }
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4GLOB_DLL std::ostream G4cerr

Here is the call graph for this function:

void PurgMagPhysicsList::SetProtonCut ( G4double  val)

Definition at line 322 of file PurgMagPhysicsList.cc.

323 {
324  //ResetCuts();
325  cutForProton = val;
326 }
void PurgMagPhysicsList::SetProtonLowLimit ( G4double  lowcut)

Definition at line 254 of file PurgMagPhysicsList.cc.

255 {
256  if (verboseLevel >0){
257 
258  G4cout << "PurgMagPhysicsList::SetCuts:";
259  G4cout << "Proton cut in energy: " << lowcut*MeV << " (MeV)" << G4endl;
260  }
261 
262  G4cerr << "PurgMagPhysicsList::SetProtonLowLimit: Not currently able to set Proton LowLimit." << G4endl;
263  G4Exception("PurgMagPhysicsList::SetProtonLowLimit()","PurMag002",
264  FatalException,"Proton Low Limit: not implemented in PurgMagPhysicsList");
265  //
266  // G4Proton::SetEnergyRange(lowcut,1e5);
267  // G4AntiProton::SetEnergyRange(lowcut,1e5);
268 }
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4GLOB_DLL std::ostream G4cerr

Here is the call graph for this function:


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