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

#include <XrayTelPhysicsList.hh>

Inheritance diagram for XrayTelPhysicsList:
Collaboration diagram for XrayTelPhysicsList:

Public Member Functions

 XrayTelPhysicsList ()
 
 ~XrayTelPhysicsList ()
 
void SetCutForGamma (G4double)
 
void SetCutForElectron (G4double)
 
G4double GetCutForGamma () const
 
G4double GetCutForElectron () const
 
- 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 ConstructMesons ()
 
void ConstructBaryons ()
 
void ConstructIons ()
 
void ConstructAllShortLiveds ()
 
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 57 of file XrayTelPhysicsList.hh.

Constructor & Destructor Documentation

XrayTelPhysicsList::XrayTelPhysicsList ( )

Definition at line 70 of file XrayTelPhysicsList.cc.

71 {
72  // Default cut values
73  defaultCutValue = 2.0*mm;
74  cutForGamma = 1.0*m;
75  cutForElectron = 1.0*m;
76 
77  SetVerboseLevel(1);
78 }
static constexpr double mm
Definition: G4SIunits.hh:115
static constexpr double m
Definition: G4SIunits.hh:129
void SetVerboseLevel(G4int value)

Here is the call graph for this function:

XrayTelPhysicsList::~XrayTelPhysicsList ( )

Definition at line 80 of file XrayTelPhysicsList.cc.

81 {}

Member Function Documentation

void XrayTelPhysicsList::ConstructAllShortLiveds ( )
protected

Definition at line 143 of file XrayTelPhysicsList.cc.

144 {
145  //Short-lived
146  G4ShortLivedConstructor slConstructor;
147  slConstructor.ConstructParticle();
148 }

Here is the call graph for this function:

Here is the caller graph for this function:

void XrayTelPhysicsList::ConstructBaryons ( )
protected

Definition at line 129 of file XrayTelPhysicsList.cc.

130 {
131  // barions
132  G4BaryonConstructor bConstructor;
133  bConstructor.ConstructParticle();
134 }
static void ConstructParticle()

Here is the call graph for this function:

Here is the caller graph for this function:

void XrayTelPhysicsList::ConstructBosons ( )
protected

Definition at line 96 of file XrayTelPhysicsList.cc.

97 {
98  // pseudo-particles
101 
102  // gamma
104 
105  // optical photon
107 }
static G4Geantino * GeantinoDefinition()
Definition: G4Geantino.cc:82
static G4ChargedGeantino * ChargedGeantinoDefinition()
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 XrayTelPhysicsList::ConstructEM ( )
protected

Definition at line 182 of file XrayTelPhysicsList.cc.

183 {
185  particleIterator->reset();
186 
187  while( (*particleIterator)() )
188  {
189  G4ParticleDefinition* particle = particleIterator->value();
190  G4ProcessManager* pmanager = particle->GetProcessManager();
191  G4String particleName = particle->GetParticleName();
192 
193  if (particleName == "gamma")
194  {
195  //gamma
196  pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
197  pmanager->AddDiscreteProcess(new G4ComptonScattering());
198  pmanager->AddDiscreteProcess(new G4GammaConversion());
199  }
200  else if (particleName == "e-")
201  {
202  //electron
203  pmanager->AddProcess(new G4eMultipleScattering(),-1, 1,1);
204  pmanager->AddProcess(new G4eIonisation(), -1, 2,2);
205  pmanager->AddProcess(new G4eBremsstrahlung(), -1, 3,3);
206 
207  } else if (particleName == "e+") {
208 
209  //positron
210  pmanager->AddProcess(new G4eMultipleScattering(),-1, 1,1);
211  pmanager->AddProcess(new G4eIonisation(), -1, 2,2);
212  pmanager->AddProcess(new G4eBremsstrahlung(), -1, 3,3);
213  pmanager->AddProcess(new G4eplusAnnihilation(), 0,-1,4);
214 
215  } else if( particleName == "mu-" ||
216  particleName == "mu+" ) {
217 
218  //muon
219  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
220  pmanager->AddProcess(new G4MuIonisation, -1, 2, 2);
221  pmanager->AddProcess(new G4MuBremsstrahlung, -1, 3, 3);
222  pmanager->AddProcess(new G4MuPairProduction, -1, 4, 4);
223 
224  } else if( particleName == "pi-" ||
225  particleName == "pi+" ) {
226 
227  //pions
228  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
229  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
230  pmanager->AddProcess(new G4hBremsstrahlung, -1, 3, 3);
231  pmanager->AddProcess(new G4hPairProduction, -1, 4, 4);
232 
233  } else if( particleName == "proton" ) {
234 
235  //proton
236  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
237  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
238  pmanager->AddProcess(new G4hBremsstrahlung, -1, 3, 3);
239  pmanager->AddProcess(new G4hPairProduction, -1, 4, 4);
240 
241  } else if( particleName == "alpha" ||
242  particleName == "He3" ||
243  particleName == "GenericIon" ) {
244 
245  //Ions
246  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
247  pmanager->AddProcess(new G4ionIonisation, -1, 2, 2);
248 
249  } else if ((!particle->IsShortLived()) &&
250  (particle->GetPDGCharge() != 0.0) &&
251  (particle->GetParticleName() != "chargedgeantino")) {
252 
253  //all others charged particles except geantino
254  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
255  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
256 
257  }
258  }
259 }
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 XrayTelPhysicsList::ConstructGeneral ( )
protected

Definition at line 263 of file XrayTelPhysicsList.cc.

264 {
265  G4Decay* theDecayProcess = new G4Decay();
267  particleIterator->reset();
268  while( (*particleIterator)() ){
269  G4ParticleDefinition* particle = particleIterator->value();
270  G4ProcessManager* pmanager = particle->GetProcessManager();
271  if (theDecayProcess->IsApplicable(*particle)) {
272  pmanager ->AddProcess(theDecayProcess);
273  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
274  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
275  }
276  }
277 }
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 XrayTelPhysicsList::ConstructIons ( )
protected

Definition at line 136 of file XrayTelPhysicsList.cc.

137 {
138  // ions
139  G4IonConstructor iConstructor;
140  iConstructor.ConstructParticle();
141 }
static void ConstructParticle()

Here is the call graph for this function:

Here is the caller graph for this function:

void XrayTelPhysicsList::ConstructLeptons ( )
protected

Definition at line 108 of file XrayTelPhysicsList.cc.

109 {
110  // leptons
113 
118 
121 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4MuonMinus * Definition()
Definition: G4MuonMinus.cc:52
static G4NeutrinoE * NeutrinoEDefinition()
Definition: G4NeutrinoE.cc:80
static G4MuonPlus * Definition()
Definition: G4MuonPlus.cc:52
static G4AntiNeutrinoMu * AntiNeutrinoMuDefinition()
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
static G4AntiNeutrinoE * AntiNeutrinoEDefinition()
static G4NeutrinoMu * NeutrinoMuDefinition()
Definition: G4NeutrinoMu.cc:80

Here is the call graph for this function:

Here is the caller graph for this function:

void XrayTelPhysicsList::ConstructMesons ( )
protected

Definition at line 122 of file XrayTelPhysicsList.cc.

123 {
124  // mesons
125  G4MesonConstructor mConstructor;
126  mConstructor.ConstructParticle();
127 }
static void ConstructParticle()

Here is the call graph for this function:

Here is the caller graph for this function:

void XrayTelPhysicsList::ConstructParticle ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 83 of file XrayTelPhysicsList.cc.

84 {
85  // Here are constructed all particles
90  ConstructIons();
92 }

Here is the call graph for this function:

void XrayTelPhysicsList::ConstructProcess ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 151 of file XrayTelPhysicsList.cc.

152 {
153  // Transportation, electromagnetic and general processes
155  ConstructEM();
157 }

Here is the call graph for this function:

G4double XrayTelPhysicsList::GetCutForElectron ( ) const

Definition at line 316 of file XrayTelPhysicsList.cc.

317 {
318  return cutForElectron;
319 }
G4double XrayTelPhysicsList::GetCutForGamma ( ) const

Definition at line 311 of file XrayTelPhysicsList.cc.

312 {
313  return cutForGamma;
314 }
void XrayTelPhysicsList::SetCutForElectron ( G4double  cut)

Definition at line 305 of file XrayTelPhysicsList.cc.

306 {
307  ResetCuts();
308  cutForElectron = cut;
309 }
void ResetCuts()
obsolete methods

Here is the call graph for this function:

void XrayTelPhysicsList::SetCutForGamma ( G4double  cut)

Definition at line 299 of file XrayTelPhysicsList.cc.

300 {
301  ResetCuts();
302  cutForGamma = cut;
303 }
void ResetCuts()
obsolete methods

Here is the call graph for this function:

void XrayTelPhysicsList::SetCuts ( )
protectedvirtual

Reimplemented from G4VUserPhysicsList.

Definition at line 279 of file XrayTelPhysicsList.cc.

280 {
281  // defaultCutValue you have typed in is used
282 
283  if (verboseLevel >1){
284  G4cout << "XrayTelPhysicsList::SetCuts:" << G4endl;
285  }
286 
287  // set cut values for gamma at first and for e- second
288  SetCutValue(cutForGamma, "gamma");
289  SetCutValue(cutForElectron, "e-");
290  SetCutValue(cutForElectron, "e+");
291 
292  // SetCutValueForOthers(defaultCutValue);
293 
294  if (verboseLevel >1) {
296  }
297 }
void SetCutValue(G4double aCut, const G4String &pname)
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:


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