Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 ()
 

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 ( G4String  physName)

Definition at line 73 of file F04PhysicsList.cc.

74 {
76 
77  defaultCutValue = 1.*mm;
78  fCutForGamma = defaultCutValue;
79  fCutForElectron = defaultCutValue;
80  fCutForPositron = defaultCutValue;
81 
82  fMessenger = new F04PhysicsListMessenger(this);
83 
84  SetVerboseLevel(1);
85 
86 // G4PhysListFactory factory;
87  G4VModularPhysicsList* phys = NULL;
88  if (physName == "QGSP_BERT") {
89  phys = new QGSP_BERT;
90  } else {
91  phys = new FTFP_BERT;
92  }
93 
94 // if (factory.IsReferencePhysList(physName))
95 // phys =factory.GetReferencePhysList(physName);
96 
97  // Physics List is defined via environment variable PHYSLIST
98 // if (!phys) phys = factory.ReferencePhysList();
99 
100  if (!phys) G4Exception("WLSPhysicsList::WLSPhysicsList","InvalidSetup",
101  FatalException,"PhysicsList does not exist");
102 
103  for (G4int i = 0; ; ++i) {
104  G4VPhysicsConstructor* elem =
105  const_cast<G4VPhysicsConstructor*> (phys->GetPhysics(i));
106  if (elem == NULL) break;
107  G4cout << "RegisterPhysics: " << elem->GetPhysicsName() << G4endl;
108  RegisterPhysics(elem);
109  }
110 
113 
114  fMaxChargedStep = DBL_MAX;
115 }
TQGSP_BERT< G4VModularPhysicsList > QGSP_BERT
Definition: QGSP_BERT.hh:63
void RegisterPhysics(G4VPhysicsConstructor *)
static G4LossTableManager * Instance()
static constexpr double mm
Definition: G4SIunits.hh:115
Provide control of the physics list and cut parameters.
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
const G4VPhysicsConstructor * GetPhysics(G4int index) const
const G4String & GetPhysicsName() const
void SetVerboseLevel(G4int value)
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

Here is the call graph for this function:

F04PhysicsList::~F04PhysicsList ( )
virtual

Definition at line 119 of file F04PhysicsList.cc.

120 {
121  delete fMessenger;
122 
123  //delete fStepMaxProcess;
124 }

Member Function Documentation

void F04PhysicsList::AddStepMax ( )

Definition at line 320 of file F04PhysicsList.cc.

321 {
322  // Step limitation seen as a process
323 
325  particleIterator->reset();
326  while ((*particleIterator)()){
327  G4ParticleDefinition* particle = particleIterator->value();
328  G4ProcessManager* pmanager = particle->GetProcessManager();
329 
330  if (fStepMaxProcess->IsApplicable(*particle) && !particle->IsShortLived())
331  {
332  if (pmanager) pmanager ->AddDiscreteProcess(fStepMaxProcess);
333  }
334  }
335 }
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: F04StepMax.cc:57
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=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 F04PhysicsList::ConstructParticle ( void  )
virtual

Reimplemented from G4VModularPhysicsList.

Definition at line 128 of file F04PhysicsList.cc.

129 {
131 
133 
134  G4DecayTable* muonPlusDecayTable = new G4DecayTable();
135  muonPlusDecayTable -> Insert(new
136  G4MuonDecayChannelWithSpin("mu+",0.986));
137  muonPlusDecayTable -> Insert(new
139  G4MuonPlus::MuonPlusDefinition() -> SetDecayTable(muonPlusDecayTable);
140 
141  G4DecayTable* muonMinusDecayTable = new G4DecayTable();
142  muonMinusDecayTable -> Insert(new
143  G4MuonDecayChannelWithSpin("mu-",0.986));
144  muonMinusDecayTable -> Insert(new
146  G4MuonMinus::MuonMinusDefinition() -> SetDecayTable(muonMinusDecayTable);
147 }
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:

void F04PhysicsList::ConstructProcess ( void  )
virtual

Reimplemented from G4VModularPhysicsList.

Definition at line 151 of file F04PhysicsList.cc.

152 {
154 
155  fStepMaxProcess = new F04StepMax();
156  G4AutoDelete::Register(fStepMaxProcess);
157 
158  G4DecayWithSpin* decayWithSpin = new G4DecayWithSpin();
159 
161 
162  G4VProcess* decay;
163  decay = processTable->FindProcess("Decay",G4MuonPlus::MuonPlus());
164 
165  G4ProcessManager* pmanager;
166  pmanager = G4MuonPlus::MuonPlus()->GetProcessManager();
167 
168  if (pmanager) {
169  if (decay) pmanager->RemoveProcess(decay);
170  pmanager->AddProcess(decayWithSpin);
171  // set ordering for PostStepDoIt and AtRestDoIt
172  pmanager ->SetProcessOrdering(decayWithSpin, idxPostStep);
173  pmanager ->SetProcessOrdering(decayWithSpin, idxAtRest);
174  }
175 
176  decay = processTable->FindProcess("Decay",G4MuonMinus::MuonMinus());
177 
179 
180  if (pmanager) {
181  if (decay) pmanager->RemoveProcess(decay);
182  pmanager->AddProcess(decayWithSpin);
183  // set ordering for PostStepDoIt and AtRestDoIt
184  pmanager ->SetProcessOrdering(decayWithSpin, idxPostStep);
185  pmanager ->SetProcessOrdering(decayWithSpin, idxAtRest);
186  }
187 
188  G4VProcess* process = processTable->
189  FindProcess("muMinusCaptureAtRest",G4MuonMinus::MuonMinus());
190 
191  if (pmanager) {
192  if (process) pmanager->RemoveProcess(process);
193  process = new G4MuonMinusCapture(new G4MuMinusCapturePrecompound());
194  pmanager->AddRestProcess(process);
195  }
196 
197  G4PionDecayMakeSpin* poldecay = new G4PionDecayMakeSpin();
198 
199  decay = processTable->FindProcess("Decay",G4PionPlus::PionPlus());
200 
201  pmanager = G4PionPlus::PionPlus()->GetProcessManager();
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  decay = processTable->FindProcess("Decay",G4PionMinus::PionMinus());
212 
214 
215  if (pmanager) {
216  if (decay) pmanager->RemoveProcess(decay);
217  pmanager->AddProcess(poldecay);
218  // set ordering for PostStepDoIt and AtRestDoIt
219  pmanager ->SetProcessOrdering(poldecay, idxPostStep);
220  pmanager ->SetProcessOrdering(poldecay, idxAtRest);
221  }
222 
223  AddStepMax();
224 }
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
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)
G4ProcessManager * GetProcessManager() const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
G4int AddRestProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4VProcess * RemoveProcess(G4VProcess *aProcess)
static G4ProcessTable * GetProcessTable()
G4VProcess * FindProcess(const G4String &processName, const G4String &particleName) const

Here is the call graph for this function:

F04StepMax * F04PhysicsList::GetStepMaxProcess ( )

Definition at line 313 of file F04PhysicsList.cc.

314 {
315  return fStepMaxProcess;
316 }
void F04PhysicsList::SetCutForElectron ( G4double  cut)

Definition at line 289 of file F04PhysicsList.cc.

290 {
291  fCutForElectron = cut;
292  SetParticleCuts(fCutForElectron, G4Electron::Electron());
293 }
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:

void F04PhysicsList::SetCutForGamma ( G4double  cut)

Definition at line 281 of file F04PhysicsList.cc.

282 {
283  fCutForGamma = cut;
284  SetParticleCuts(fCutForGamma, G4Gamma::Gamma());
285 }
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:

void F04PhysicsList::SetCutForPositron ( G4double  cut)

Definition at line 297 of file F04PhysicsList.cc.

298 {
299  fCutForPositron = cut;
300  SetParticleCuts(fCutForPositron, G4Positron::Positron());
301 }
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:

void F04PhysicsList::SetCuts ( )
virtual

Reimplemented from G4VUserPhysicsList.

Definition at line 262 of file F04PhysicsList.cc.

263 {
264  if (verboseLevel >0) {
265  G4cout << "F04PhysicsList::SetCuts:";
266  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length")
267  << G4endl;
268  }
269 
270  // set cut values for gamma at first and for e- second and next for e+,
271  // because some processes for e+/e- need cut values for gamma
272  SetCutValue(fCutForGamma, "gamma");
273  SetCutValue(fCutForElectron, "e-");
274  SetCutValue(fCutForPositron, "e+");
275 
277 }
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 F04PhysicsList::SetStepMax ( G4double  step)

Definition at line 305 of file F04PhysicsList.cc.

306 {
307  fMaxChargedStep = step ;
308  fStepMaxProcess->SetStepMax(fMaxChargedStep);
309 }
void SetStepMax(G4double)
Definition: F04StepMax.cc:64

Here is the call graph for this function:

Here is the caller graph for this function:


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