Geant4  10.02.p03
F04PhysicsList Class Reference

#include <F04PhysicsList.hh>

Inheritance diagram for F04PhysicsList:
Collaboration diagram for F04PhysicsList:

Public Member Functions

 F04PhysicsList (G4String)
 
virtual ~F04PhysicsList ()
 
virtual void SetCuts ()
 
void SetCutForGamma (G4double)
 
void SetCutForElectron (G4double)
 
void SetCutForPositron (G4double)
 
void SetStepMax (G4double)
 
F04StepMaxGetStepMaxProcess ()
 
void AddStepMax ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
- Public Member Functions inherited from G4VModularPhysicsList
 G4VModularPhysicsList ()
 
virtual ~G4VModularPhysicsList ()
 
void RegisterPhysics (G4VPhysicsConstructor *)
 
const G4VPhysicsConstructorGetPhysics (G4int index) const
 
const G4VPhysicsConstructorGetPhysics (const G4String &name) const
 
const G4VPhysicsConstructorGetPhysicsWithType (G4int physics_type) const
 
void ReplacePhysics (G4VPhysicsConstructor *)
 
void RemovePhysics (G4VPhysicsConstructor *)
 
void RemovePhysics (G4int type)
 
void RemovePhysics (const G4String &name)
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4int GetInstanceID () 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 ()
 

Private Attributes

G4double fCutForGamma
 
G4double fCutForElectron
 
G4double fCutForPositron
 
G4double fMaxChargedStep
 
F04PhysicsListMessengerfMessenger
 

Static Private Attributes

static G4ThreadLocal F04StepMaxfStepMaxProcess = 0
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VModularPhysicsList
static const G4VMPLManagerGetSubInstanceManager ()
 
- Static Public Member Functions inherited from G4VUserPhysicsList
static const G4VUPLManagerGetSubInstanceManager ()
 
- Protected Types inherited from G4VModularPhysicsList
typedef G4VMPLData::G4PhysConstVectorData G4PhysConstVector
 
- Protected Member Functions inherited from G4VModularPhysicsList
 G4VModularPhysicsList (const G4VModularPhysicsList &)
 
G4VModularPhysicsListoperator= (const G4VModularPhysicsList &)
 
- 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 G4VModularPhysicsList
G4int verboseLevel
 
G4int g4vmplInstanceID
 
- 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 G4VModularPhysicsList
static G4RUN_DLL G4VMPLManager G4VMPLsubInstanceManager
 
- Static Protected Attributes inherited from G4VUserPhysicsList
static G4RUN_DLL G4VUPLManager subInstanceManager
 

Detailed Description

Definition at line 43 of file F04PhysicsList.hh.

Constructor & Destructor Documentation

◆ F04PhysicsList()

F04PhysicsList::F04PhysicsList ( G4String  physName)

Definition at line 70 of file F04PhysicsList.cc.

71 {
73 
74  defaultCutValue = 1.*mm;
78 
80 
81  SetVerboseLevel(1);
82 
83 // G4PhysListFactory factory;
84  G4VModularPhysicsList* phys = NULL;
85  if (physName == "QGSP_BERT") {
86  phys = new QGSP_BERT;
87  } else {
88  phys = new FTFP_BERT;
89  }
90 
91 // if (factory.IsReferencePhysList(physName))
92 // phys =factory.GetReferencePhysList(physName);
93 
94  // Physics List is defined via environment variable PHYSLIST
95 // if (!phys) phys = factory.ReferencePhysList();
96 
97  if (!phys) G4Exception("WLSPhysicsList::WLSPhysicsList","InvalidSetup",
98  FatalException,"PhysicsList does not exist");
99 
100  for (G4int i = 0; ; ++i) {
101  G4VPhysicsConstructor* elem =
102  const_cast<G4VPhysicsConstructor*> (phys->GetPhysics(i));
103  if (elem == NULL) break;
104  G4cout << "RegisterPhysics: " << elem->GetPhysicsName() << G4endl;
105  RegisterPhysics(elem);
106  }
107 
110 
112 }
TQGSP_BERT< G4VModularPhysicsList > QGSP_BERT
Definition: QGSP_BERT.hh:63
void RegisterPhysics(G4VPhysicsConstructor *)
static G4LossTableManager * Instance()
G4double fCutForGamma
G4double fCutForPositron
Provide control of the physics list and cut parameters.
int G4int
Definition: G4Types.hh:78
const G4VPhysicsConstructor * GetPhysics(G4int index) const
G4double fCutForElectron
const G4String & GetPhysicsName() const
G4GLOB_DLL std::ostream G4cout
F04PhysicsListMessenger * fMessenger
void SetVerboseLevel(G4int value)
G4double fMaxChargedStep
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
TFTFP_BERT< G4VModularPhysicsList > FTFP_BERT
Definition: FTFP_BERT.hh:63
#define G4endl
Definition: G4ios.hh:61
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:

◆ ~F04PhysicsList()

F04PhysicsList::~F04PhysicsList ( )
virtual

Definition at line 116 of file F04PhysicsList.cc.

117 {
118  delete fMessenger;
119 
120  //delete fStepMaxProcess;
121 }
F04PhysicsListMessenger * fMessenger

Member Function Documentation

◆ AddStepMax()

void F04PhysicsList::AddStepMax ( )

Definition at line 308 of file F04PhysicsList.cc.

309 {
310  // Step limitation seen as a process
311 
313  particleIterator->reset();
314  while ((*particleIterator)()){
315  G4ParticleDefinition* particle = particleIterator->value();
316  G4ProcessManager* pmanager = particle->GetProcessManager();
317 
318  if (fStepMaxProcess->IsApplicable(*particle) && !particle->IsShortLived())
319  {
320  if (pmanager) pmanager ->AddDiscreteProcess(fStepMaxProcess);
321  }
322  }
323 }
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: F04StepMax.cc:57
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ProcessManager * GetProcessManager() const
static G4ThreadLocal F04StepMax * fStepMaxProcess
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:

◆ ConstructParticle()

void F04PhysicsList::ConstructParticle ( void  )
virtual

Reimplemented from G4VModularPhysicsList.

Definition at line 125 of file F04PhysicsList.cc.

126 {
128 
130 
131  G4DecayTable* muonPlusDecayTable = new G4DecayTable();
132  muonPlusDecayTable -> Insert(new
133  G4MuonDecayChannelWithSpin("mu+",0.986));
134  muonPlusDecayTable -> Insert(new
136  G4MuonPlus::MuonPlusDefinition() -> SetDecayTable(muonPlusDecayTable);
137 
138  G4DecayTable* muonMinusDecayTable = new G4DecayTable();
139  muonMinusDecayTable -> Insert(new
140  G4MuonDecayChannelWithSpin("mu-",0.986));
141  muonMinusDecayTable -> Insert(new
143  G4MuonMinus::MuonMinusDefinition() -> SetDecayTable(muonMinusDecayTable);
144 }
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
void Insert(const PVNodeID *pvPath, size_t pathLength, G4int index, Node *node)
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
Here is the call graph for this function:

◆ ConstructProcess()

void F04PhysicsList::ConstructProcess ( void  )
virtual

Reimplemented from G4VModularPhysicsList.

Definition at line 148 of file F04PhysicsList.cc.

149 {
151 
152  fStepMaxProcess = new F04StepMax();
154 
155  G4DecayWithSpin* decayWithSpin = new G4DecayWithSpin();
156 
158 
159  G4VProcess* decay;
160  decay = processTable->FindProcess("Decay",G4MuonPlus::MuonPlus());
161 
162  G4ProcessManager* pmanager;
163  pmanager = G4MuonPlus::MuonPlus()->GetProcessManager();
164 
165  if (pmanager) {
166  if (decay) pmanager->RemoveProcess(decay);
167  pmanager->AddProcess(decayWithSpin);
168  // set ordering for PostStepDoIt and AtRestDoIt
169  pmanager ->SetProcessOrdering(decayWithSpin, idxPostStep);
170  pmanager ->SetProcessOrdering(decayWithSpin, idxAtRest);
171  }
172 
173  decay = processTable->FindProcess("Decay",G4MuonMinus::MuonMinus());
174 
176 
177  if (pmanager) {
178  if (decay) pmanager->RemoveProcess(decay);
179  pmanager->AddProcess(decayWithSpin);
180  // set ordering for PostStepDoIt and AtRestDoIt
181  pmanager ->SetProcessOrdering(decayWithSpin, idxPostStep);
182  pmanager ->SetProcessOrdering(decayWithSpin, idxAtRest);
183  }
184 
185  G4PionDecayMakeSpin* poldecay = new G4PionDecayMakeSpin();
186 
187  decay = processTable->FindProcess("Decay",G4PionPlus::PionPlus());
188 
189  pmanager = G4PionPlus::PionPlus()->GetProcessManager();
190 
191  if (pmanager) {
192  if (decay) pmanager->RemoveProcess(decay);
193  pmanager->AddProcess(poldecay);
194  // set ordering for PostStepDoIt and AtRestDoIt
195  pmanager ->SetProcessOrdering(poldecay, idxPostStep);
196  pmanager ->SetProcessOrdering(poldecay, idxAtRest);
197  }
198 
199  decay = processTable->FindProcess("Decay",G4PionMinus::PionMinus());
200 
202 
203  if (pmanager) {
204  if (decay) pmanager->RemoveProcess(decay);
205  pmanager->AddProcess(poldecay);
206  // set ordering for PostStepDoIt and AtRestDoIt
207  pmanager ->SetProcessOrdering(poldecay, idxPostStep);
208  pmanager ->SetProcessOrdering(poldecay, idxAtRest);
209  }
210 
211  AddStepMax();
212 }
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
G4VProcess * FindProcess(const G4String &processName, const G4String &particleName) const
G4ProcessManager * GetProcessManager() const
static G4ThreadLocal F04StepMax * fStepMaxProcess
void Register(T *inst)
Definition: G4AutoDelete.hh:65
ParticleList decay(Cluster *const c)
Carries out a cluster decay.
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
G4VProcess * RemoveProcess(G4VProcess *aProcess)
static G4ProcessTable * GetProcessTable()
Here is the call graph for this function:

◆ GetStepMaxProcess()

F04StepMax * F04PhysicsList::GetStepMaxProcess ( )

Definition at line 301 of file F04PhysicsList.cc.

302 {
303  return fStepMaxProcess;
304 }
static G4ThreadLocal F04StepMax * fStepMaxProcess

◆ SetCutForElectron()

void F04PhysicsList::SetCutForElectron ( G4double  cut)

Definition at line 277 of file F04PhysicsList.cc.

278 {
279  fCutForElectron = cut;
281 }
G4double fCutForElectron
void SetParticleCuts(G4double cut, G4ParticleDefinition *particle, G4Region *region=0)
static G4Electron * Electron()
Definition: G4Electron.cc:94
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetCutForGamma()

void F04PhysicsList::SetCutForGamma ( G4double  cut)

Definition at line 269 of file F04PhysicsList.cc.

270 {
271  fCutForGamma = cut;
273 }
G4double fCutForGamma
void SetParticleCuts(G4double cut, G4ParticleDefinition *particle, G4Region *region=0)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetCutForPositron()

void F04PhysicsList::SetCutForPositron ( G4double  cut)

Definition at line 285 of file F04PhysicsList.cc.

286 {
287  fCutForPositron = cut;
289 }
G4double fCutForPositron
void SetParticleCuts(G4double cut, G4ParticleDefinition *particle, G4Region *region=0)
static G4Positron * Positron()
Definition: G4Positron.cc:94
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetCuts()

void F04PhysicsList::SetCuts ( )
virtual

Reimplemented from G4VUserPhysicsList.

Definition at line 250 of file F04PhysicsList.cc.

251 {
252  if (verboseLevel >0) {
253  G4cout << "F04PhysicsList::SetCuts:";
254  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length")
255  << G4endl;
256  }
257 
258  // set cut values for gamma at first and for e- second and next for e+,
259  // because some processes for e+/e- need cut values for gamma
260  SetCutValue(fCutForGamma, "gamma");
263 
265 }
G4double fCutForGamma
G4double fCutForPositron
void SetCutValue(G4double aCut, const G4String &pname)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4double fCutForElectron
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ SetStepMax()

void F04PhysicsList::SetStepMax ( G4double  step)

Definition at line 293 of file F04PhysicsList.cc.

294 {
295  fMaxChargedStep = step ;
297 }
void SetStepMax(G4double)
Definition: F04StepMax.cc:64
static G4ThreadLocal F04StepMax * fStepMaxProcess
G4double fMaxChargedStep
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fCutForElectron

G4double F04PhysicsList::fCutForElectron
private

Definition at line 71 of file F04PhysicsList.hh.

◆ fCutForGamma

G4double F04PhysicsList::fCutForGamma
private

Definition at line 70 of file F04PhysicsList.hh.

◆ fCutForPositron

G4double F04PhysicsList::fCutForPositron
private

Definition at line 72 of file F04PhysicsList.hh.

◆ fMaxChargedStep

G4double F04PhysicsList::fMaxChargedStep
private

Definition at line 74 of file F04PhysicsList.hh.

◆ fMessenger

F04PhysicsListMessenger* F04PhysicsList::fMessenger
private

Definition at line 77 of file F04PhysicsList.hh.

◆ fStepMaxProcess

G4ThreadLocal F04StepMax * F04PhysicsList::fStepMaxProcess = 0
staticprivate

Definition at line 75 of file F04PhysicsList.hh.


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