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

#include <WLSPhysicsList.hh>

Inheritance diagram for WLSPhysicsList:
Collaboration diagram for WLSPhysicsList:

Public Member Functions

 WLSPhysicsList (G4String)
 
virtual ~WLSPhysicsList ()
 
void SetCuts ()
 
void SetCutForGamma (G4double)
 
void SetCutForElectron (G4double)
 
void SetCutForPositron (G4double)
 
void SetStepMax (G4double)
 
WLSStepMaxGetStepMaxProcess ()
 
void AddStepMax ()
 
void RemoveFromPhysicsList (const G4String &)
 Remove specific physics from physics list. More...
 
void ClearPhysics ()
 Make sure that the physics list is empty. More...
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
void SetAbsorption (G4bool)
 
void SetNbOfPhotonsCerenkov (G4int)
 
void SetVerbose (G4int)
 
- 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 47 of file WLSPhysicsList.hh.

Constructor & Destructor Documentation

WLSPhysicsList::WLSPhysicsList ( G4String  physName)

Definition at line 69 of file WLSPhysicsList.cc.

70 {
72 
73  defaultCutValue = 1.*mm;
74  fCutForGamma = defaultCutValue;
75  fCutForElectron = defaultCutValue;
76  fCutForPositron = defaultCutValue;
77 
78 // G4PhysListFactory factory;
79  G4VModularPhysicsList* phys = NULL;
80  if (physName == "QGSP_BERT_HP") {
81  phys = new QGSP_BERT_HP;
82  } else {
83  phys = new FTFP_BERT;
84  }
85 // if (factory.IsReferencePhysList(physName)) {
86 // phys = factory.GetReferencePhysList(physName);
87 // if(!phys)G4Exception("WLSPhysicsList::WLSPhysicsList","InvalidSetup",
88 // FatalException,"PhysicsList does not exist");
89  fMessenger = new WLSPhysicsListMessenger(this);
90 // }
91 
92  for (G4int i = 0; ; ++i) {
93  G4VPhysicsConstructor* elem =
94  const_cast<G4VPhysicsConstructor*> (phys->GetPhysics(i));
95  if (elem == NULL) break;
96  G4cout << "RegisterPhysics: " << elem->GetPhysicsName() << G4endl;
97  RegisterPhysics(elem);
98  }
99 
100  fAbsorptionOn = true;
101 
103  RegisterPhysics(fOpticalPhysics = new WLSOpticalPhysics(fAbsorptionOn));
104 
106 
107  fStepMaxProcess = new WLSStepMax();
108 }
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
TFTFP_BERT< G4VModularPhysicsList > FTFP_BERT
Definition: FTFP_BERT.hh:63
#define G4endl
Definition: G4ios.hh:61
TQGSP_BERT_HP< G4VModularPhysicsList > QGSP_BERT_HP
Definition: QGSP_BERT_HP.hh:62

Here is the call graph for this function:

WLSPhysicsList::~WLSPhysicsList ( )
virtual

Definition at line 112 of file WLSPhysicsList.cc.

113 {
114  delete fMessenger;
115 
116  delete fStepMaxProcess;
117 }

Member Function Documentation

void WLSPhysicsList::AddStepMax ( )

Definition at line 312 of file WLSPhysicsList.cc.

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

Make sure that the physics list is empty.

Definition at line 121 of file WLSPhysicsList.cc.

122 {
123  for (G4PhysConstVector::iterator p = fPhysicsVector->begin();
124  p != fPhysicsVector->end(); ++p) {
125  delete (*p);
126  }
127  fPhysicsVector->clear();
128 }
const char * p
Definition: xmltok.h:285

Here is the caller graph for this function:

void WLSPhysicsList::ConstructParticle ( void  )
virtual

Reimplemented from G4VModularPhysicsList.

Definition at line 132 of file WLSPhysicsList.cc.

133 {
135 
136  G4DecayTable* MuonPlusDecayTable = new G4DecayTable();
137  MuonPlusDecayTable -> Insert(new
138  G4MuonDecayChannelWithSpin("mu+",0.986));
139  MuonPlusDecayTable -> Insert(new
141  G4MuonPlus::MuonPlusDefinition() -> SetDecayTable(MuonPlusDecayTable);
142 
143  G4DecayTable* MuonMinusDecayTable = new G4DecayTable();
144  MuonMinusDecayTable -> Insert(new
145  G4MuonDecayChannelWithSpin("mu-",0.986));
146  MuonMinusDecayTable -> Insert(new
148  G4MuonMinus::MuonMinusDefinition() -> SetDecayTable(MuonMinusDecayTable);
149 
150 }
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 WLSPhysicsList::ConstructProcess ( void  )
virtual

Reimplemented from G4VModularPhysicsList.

Definition at line 154 of file WLSPhysicsList.cc.

155 {
157 
158  SetVerbose(0);
159 
160  G4DecayWithSpin* decayWithSpin = new G4DecayWithSpin();
161 
163 
164  G4VProcess* decay;
165  decay = processTable->FindProcess("Decay",G4MuonPlus::MuonPlus());
166 
167  G4ProcessManager* pManager;
168  pManager = G4MuonPlus::MuonPlus()->GetProcessManager();
169 
170  if (pManager) {
171  if (decay) pManager->RemoveProcess(decay);
172  pManager->AddProcess(decayWithSpin);
173  // set ordering for PostStepDoIt and AtRestDoIt
174  pManager ->SetProcessOrdering(decayWithSpin, idxPostStep);
175  pManager ->SetProcessOrdering(decayWithSpin, idxAtRest);
176  }
177 
178  decay = processTable->FindProcess("Decay",G4MuonMinus::MuonMinus());
179 
181 
182  if (pManager) {
183  if (decay) pManager->RemoveProcess(decay);
184  pManager->AddProcess(decayWithSpin);
185  // set ordering for PostStepDoIt and AtRestDoIt
186  pManager ->SetProcessOrdering(decayWithSpin, idxPostStep);
187  pManager ->SetProcessOrdering(decayWithSpin, idxAtRest);
188  }
189 
190  G4PionDecayMakeSpin* poldecay = new G4PionDecayMakeSpin();
191 
192  decay = processTable->FindProcess("Decay",G4PionPlus::PionPlus());
193 
194  pManager = G4PionPlus::PionPlus()->GetProcessManager();
195 
196  if (pManager) {
197  if (decay) pManager->RemoveProcess(decay);
198  pManager->AddProcess(poldecay);
199  // set ordering for PostStepDoIt and AtRestDoIt
200  pManager ->SetProcessOrdering(poldecay, idxPostStep);
201  pManager ->SetProcessOrdering(poldecay, idxAtRest);
202  }
203 
204  decay = processTable->FindProcess("Decay",G4PionMinus::PionMinus());
205 
207 
208  if (pManager) {
209  if (decay) pManager->RemoveProcess(decay);
210  pManager->AddProcess(poldecay);
211  // set ordering for PostStepDoIt and AtRestDoIt
212  pManager ->SetProcessOrdering(poldecay, idxPostStep);
213  pManager ->SetProcessOrdering(poldecay, idxAtRest);
214  }
215 
216  AddStepMax();
217 
218 }
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
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
G4VProcess * RemoveProcess(G4VProcess *aProcess)
static G4ProcessTable * GetProcessTable()
void SetVerbose(G4int)
G4VProcess * FindProcess(const G4String &processName, const G4String &particleName) const

Here is the call graph for this function:

WLSStepMax * WLSPhysicsList::GetStepMaxProcess ( )

Definition at line 305 of file WLSPhysicsList.cc.

306 {
307  return fStepMaxProcess;
308 }
void WLSPhysicsList::RemoveFromPhysicsList ( const G4String name)

Remove specific physics from physics list.

Definition at line 222 of file WLSPhysicsList.cc.

223 {
224  G4bool success = false;
225  for (G4PhysConstVector::iterator p = fPhysicsVector->begin();
226  p != fPhysicsVector->end(); ++p) {
227  G4VPhysicsConstructor* e = (*p);
228  if (e->GetPhysicsName() == name) {
229  fPhysicsVector->erase(p);
230  success = true;
231  break;
232  }
233  }
234  if (!success) {
235  G4ExceptionDescription message;
236  message << "PhysicsList::RemoveFromEMPhysicsList "<< name << "not found";
237  G4Exception("example WLSPhysicsList::RemoveFromPhysicsList()",
238  "ExamWLSPhysicsList01",FatalException,message);
239  }
240 }
const XML_Char * name
Definition: expat.h:151
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const char * p
Definition: xmltok.h:285
bool G4bool
Definition: G4Types.hh:79
const G4String & GetPhysicsName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Here is the caller graph for this function:

void WLSPhysicsList::SetAbsorption ( G4bool  toggle)

Definition at line 244 of file WLSPhysicsList.cc.

245 {
246  fAbsorptionOn = toggle;
247  RemoveFromPhysicsList("Optical");
248  fPhysicsVector->
249  push_back(fOpticalPhysics = new WLSOpticalPhysics(toggle));
250  fOpticalPhysics->ConstructProcess();
251 }
virtual void ConstructProcess()
void RemoveFromPhysicsList(const G4String &)
Remove specific physics from physics list.

Here is the call graph for this function:

Here is the caller graph for this function:

void WLSPhysicsList::SetCutForElectron ( G4double  cut)

Definition at line 282 of file WLSPhysicsList.cc.

283 {
284  fCutForElectron = cut;
285  SetParticleCuts(fCutForElectron, G4Electron::Electron());
286 }
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 WLSPhysicsList::SetCutForGamma ( G4double  cut)

Definition at line 274 of file WLSPhysicsList.cc.

275 {
276  fCutForGamma = cut;
277  SetParticleCuts(fCutForGamma, G4Gamma::Gamma());
278 }
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 WLSPhysicsList::SetCutForPositron ( G4double  cut)

Definition at line 290 of file WLSPhysicsList.cc.

291 {
292  fCutForPositron = cut;
293  SetParticleCuts(fCutForPositron, G4Positron::Positron());
294 }
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 WLSPhysicsList::SetCuts ( )
virtual

Reimplemented from G4VUserPhysicsList.

Definition at line 255 of file WLSPhysicsList.cc.

256 {
257  if (verboseLevel >0) {
258  G4cout << "WLSPhysicsList::SetCuts:";
259  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length")
260  << G4endl;
261  }
262 
263  // set cut values for gamma at first and for e- second and next for e+,
264  // because some processes for e+/e- need cut values for gamma
265  SetCutValue(fCutForGamma, "gamma");
266  SetCutValue(fCutForElectron, "e-");
267  SetCutValue(fCutForPositron, "e+");
268 
270 }
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 WLSPhysicsList::SetNbOfPhotonsCerenkov ( G4int  maxNumber)

Definition at line 331 of file WLSPhysicsList.cc.

332 {
333  fOpticalPhysics->SetNbOfPhotonsCerenkov(maxNumber);
334 }
void SetNbOfPhotonsCerenkov(G4int)

Here is the call graph for this function:

void WLSPhysicsList::SetStepMax ( G4double  step)

Definition at line 298 of file WLSPhysicsList.cc.

299 {
300  fStepMaxProcess->SetStepMax(step);
301 }
void SetStepMax(G4double)
Definition: WLSStepMax.cc:64

Here is the call graph for this function:

Here is the caller graph for this function:

void WLSPhysicsList::SetVerbose ( G4int  verbose)

Definition at line 338 of file WLSPhysicsList.cc.

339 {
340  fOpticalPhysics->GetCerenkovProcess()->SetVerboseLevel(verbose);
341  fOpticalPhysics->GetScintillationProcess()->SetVerboseLevel(verbose);
342  fOpticalPhysics->GetAbsorptionProcess()->SetVerboseLevel(verbose);
343  fOpticalPhysics->GetRayleighScatteringProcess()->SetVerboseLevel(verbose);
344  fOpticalPhysics->GetMieHGScatteringProcess()->SetVerboseLevel(verbose);
345  fOpticalPhysics->GetBoundaryProcess()->SetVerboseLevel(verbose);
346 }
G4Scintillation * GetScintillationProcess()
G4OpMieHG * GetMieHGScatteringProcess()
G4Cerenkov * GetCerenkovProcess()
G4OpAbsorption * GetAbsorptionProcess()
G4OpBoundaryProcess * GetBoundaryProcess()
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
G4OpRayleigh * GetRayleighScatteringProcess()

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: