Geant4  10.02.p03
F03PhysicsList Class Reference

#include <F03PhysicsList.hh>

Inheritance diagram for F03PhysicsList:
Collaboration diagram for F03PhysicsList:

Public Member Functions

 F03PhysicsList (F03DetectorConstruction *)
 
virtual ~F03PhysicsList ()
 
void SetGammaCut (G4double)
 
void SetElectronCut (G4double)
 
void SetMaxStep (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 ()
 

Public Attributes

G4double fMaxChargedStep
 

Protected Member Functions

virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
void SetCuts ()
 
void ConstructBosons ()
 
void ConstructLeptons ()
 
void ConstructMesons ()
 
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
 

Private Attributes

G4double fCutForGamma
 
G4double fCutForElectron
 
F03DetectorConstructionfDet
 
F03PhysicsListMessengerfPhysicsListMessenger
 

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 46 of file F03PhysicsList.hh.

Constructor & Destructor Documentation

◆ F03PhysicsList()

F03PhysicsList::F03PhysicsList ( F03DetectorConstruction p)

Definition at line 53 of file F03PhysicsList.cc.

55 {
56  fDet = p;
57 
58  defaultCutValue = 1.000*mm;
59 
62 
63  SetVerboseLevel(1);
65 }
G4double fCutForElectron
G4double fCutForGamma
void SetVerboseLevel(G4int value)
F03PhysicsListMessenger * fPhysicsListMessenger
G4double fMaxChargedStep
F03DetectorConstruction * fDet
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:

◆ ~F03PhysicsList()

F03PhysicsList::~F03PhysicsList ( )
virtual

Definition at line 69 of file F03PhysicsList.cc.

70 {
71  delete fPhysicsListMessenger;
72 }
F03PhysicsListMessenger * fPhysicsListMessenger

Member Function Documentation

◆ ConstructBarions()

void F03PhysicsList::ConstructBarions ( )
protected

Definition at line 137 of file F03PhysicsList.cc.

138 {
139  // barions
140 
143 }
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:

◆ ConstructBosons()

void F03PhysicsList::ConstructBosons ( )
protected

Definition at line 93 of file F03PhysicsList.cc.

94 {
95  // gamma
96 
98 
99  // charged geantino
100 
102 
103 }
static G4ChargedGeantino * ChargedGeantinoDefinition()
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructEM()

void F03PhysicsList::ConstructEM ( )
protected

Definition at line 177 of file F03PhysicsList.cc.

178 {
180  particleIterator->reset();
181 
182  while( (*particleIterator)() )
183  {
184  G4ParticleDefinition* particle = particleIterator->value();
185  G4ProcessManager* pmanager = particle->GetProcessManager();
186  G4String particleName = particle->GetParticleName();
187 
188  if (particleName == "gamma")
189  {
190  // Construct processes for gamma
191 
195 
196  pmanager->AddDiscreteProcess(fPhotoElectricEffect);
197  pmanager->AddDiscreteProcess(fComptonScattering);
198 
199  pmanager->AddDiscreteProcess(fGammaConversion);
200  }
201  else if (particleName == "e-")
202  {
203  // Construct processes for electron
204 
205  G4eIonisation* feminusIonisation = new G4eIonisation();
206  G4eBremsstrahlung* feminusBremsstrahlung = new G4eBremsstrahlung();
207  F03StepCut* feminusStepCut = new F03StepCut();
208  feminusStepCut->SetMaxStep(fMaxChargedStep);
209 
210  pmanager->AddProcess(feminusIonisation,-1,2,2);
211  pmanager->AddProcess(feminusBremsstrahlung,-1,-1,3);
212  pmanager->AddProcess(feminusStepCut,-1,-1,4);
213  }
214  else if (particleName == "e+")
215  {
216  // Construct processes for positron
217 
218  G4eIonisation* feplusIonisation = new G4eIonisation();
219  G4eBremsstrahlung* feplusBremsstrahlung = new G4eBremsstrahlung();
220  F03StepCut* feplusStepCut = new F03StepCut();
221  feplusStepCut->SetMaxStep(fMaxChargedStep);
222 
223  pmanager->AddProcess(feplusIonisation,-1,2,2);
224  pmanager->AddProcess(feplusBremsstrahlung,-1,-1,3);
225  pmanager->AddProcess(feplusStepCut,-1,-1,5);
226  }
227  else if( particleName == "mu+" || particleName == "mu-" )
228  {
229  // Construct processes for muon+
230 
231  F03StepCut* muonStepCut = new F03StepCut();
232  muonStepCut->SetMaxStep(fMaxChargedStep);
233  G4MuIonisation* muIonisation = new G4MuIonisation();
234 
235  pmanager->AddProcess(new G4MuMultipleScattering(),-1,1,1);
236  pmanager->AddProcess(muIonisation,-1,2,2);
237  pmanager->AddProcess(new G4MuBremsstrahlung(),-1,-1,3);
238  pmanager->AddProcess(new G4MuPairProduction(),-1,-1,4);
239  pmanager->AddProcess(muonStepCut,-1,-1,3);
240  }
241  else if ( particleName == "proton"
242  || particleName == "antiproton"
243  || particleName == "pi+"
244  || particleName == "pi-"
245  || particleName == "kaon+"
246  || particleName == "kaon-"
247  )
248  {
249  F03StepCut* theHadronStepCut = new F03StepCut();
250  theHadronStepCut->SetMaxStep(10*mm);
251 
252  G4hIonisation* thehIonisation = new G4hIonisation();
253  G4hMultipleScattering* thehMultipleScattering =
254  new G4hMultipleScattering();
255 
256  pmanager->AddProcess(thehMultipleScattering,-1,1,1);
257  pmanager->AddProcess(thehIonisation,-1,2,2);
258  pmanager->AddProcess(theHadronStepCut,-1,-1,3);
259  }
260  }
261 }
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ProcessManager * GetProcessManager() const
void SetMaxStep(G4double)
Definition: F03StepCut.cc:66
const G4String & GetParticleName() const
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
G4double fMaxChargedStep
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructGeneral()

void F03PhysicsList::ConstructGeneral ( )
protected

Definition at line 267 of file F03PhysicsList.cc.

268 {
269  // Add Decay Process
270 
271  G4Decay* theDecayProcess = new G4Decay();
273  particleIterator->reset();
274 
275  while( (*particleIterator)() )
276  {
277  G4ParticleDefinition* particle = particleIterator->value();
278  G4ProcessManager* pmanager = particle->GetProcessManager();
279 
280  if (theDecayProcess->IsApplicable(*particle))
281  {
282  pmanager ->AddProcess(theDecayProcess);
283 
284  // set ordering for PostStepDoIt and AtRestDoIt
285 
286  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
287  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
288  }
289  }
290 }
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:

◆ ConstructLeptons()

void F03PhysicsList::ConstructLeptons ( )
protected

Definition at line 107 of file F03PhysicsList.cc.

108 {
109  // leptons
110 
115 
120 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
static G4NeutrinoE * NeutrinoEDefinition()
Definition: G4NeutrinoE.cc:80
static G4AntiNeutrinoMu * AntiNeutrinoMuDefinition()
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
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:

◆ ConstructMesons()

void F03PhysicsList::ConstructMesons ( )
protected

Definition at line 124 of file F03PhysicsList.cc.

125 {
126  // mesons
127 
133 }
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
static G4PionZero * PionZeroDefinition()
Definition: G4PionZero.cc:103
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructParticle()

void F03PhysicsList::ConstructParticle ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 76 of file F03PhysicsList.cc.

77 {
78  // In this method, static member functions should be called
79  // for all particles which you want to use.
80  // This ensures that objects of these particle types will be
81  // created in the program.
82 
87 
89 }
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
void ConstructBosons()
Here is the call graph for this function:

◆ ConstructProcess()

void F03PhysicsList::ConstructProcess ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 147 of file F03PhysicsList.cc.

148 {
150 
151  ConstructEM();
153 }
Here is the call graph for this function:

◆ SetCuts()

void F03PhysicsList::SetCuts ( )
protectedvirtual

Reimplemented from G4VUserPhysicsList.

Definition at line 294 of file F03PhysicsList.cc.

295 {
296  G4Timer theTimer;
297  theTimer.Start();
298  if (verboseLevel >0)
299  {
300  G4cout << "F03PhysicsList::SetCuts:";
301  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
302  }
303  // set cut values for gamma at first and for e- second and next for e+,
304  // because some processes for e+/e- need cut values for gamma
305  SetCutValue(fCutForGamma,"gamma");
306 
309 
311 
312  theTimer.Stop();
313  G4cout.precision(6);
314  G4cout << G4endl;
315  G4cout << "total time(SetCuts)=" << theTimer.GetUserElapsed()
316  << " s " <<G4endl;
317 
318 }
G4double fCutForElectron
void SetCutValue(G4double aCut, const G4String &pname)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4double fCutForGamma
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
G4double GetUserElapsed() const
Definition: G4Timer.cc:130
void Stop()
#define G4endl
Definition: G4ios.hh:61
void Start()
Here is the call graph for this function:

◆ SetElectronCut()

void F03PhysicsList::SetElectronCut ( G4double  val)

Definition at line 329 of file F03PhysicsList.cc.

330 {
331  fCutForElectron = val;
332 }
G4double fCutForElectron
Here is the caller graph for this function:

◆ SetGammaCut()

void F03PhysicsList::SetGammaCut ( G4double  val)

Definition at line 322 of file F03PhysicsList.cc.

323 {
324  fCutForGamma = val;
325 }
G4double fCutForGamma
Here is the caller graph for this function:

◆ SetMaxStep()

void F03PhysicsList::SetMaxStep ( G4double  step)

Definition at line 336 of file F03PhysicsList.cc.

337 {
338  fMaxChargedStep = step;
339  G4cout << " MaxChargedStep=" << fMaxChargedStep << G4endl;
340  G4cout << G4endl;
341 }
G4GLOB_DLL std::ostream G4cout
G4double fMaxChargedStep
#define G4endl
Definition: G4ios.hh:61
Here is the caller graph for this function:

Member Data Documentation

◆ fCutForElectron

G4double F03PhysicsList::fCutForElectron
private

Definition at line 87 of file F03PhysicsList.hh.

◆ fCutForGamma

G4double F03PhysicsList::fCutForGamma
private

Definition at line 86 of file F03PhysicsList.hh.

◆ fDet

F03DetectorConstruction* F03PhysicsList::fDet
private

Definition at line 89 of file F03PhysicsList.hh.

◆ fMaxChargedStep

G4double F03PhysicsList::fMaxChargedStep

Definition at line 82 of file F03PhysicsList.hh.

◆ fPhysicsListMessenger

F03PhysicsListMessenger* F03PhysicsList::fPhysicsListMessenger
private

Definition at line 90 of file F03PhysicsList.hh.


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