Geant4  10.02.p03
Em10PhysicsList Class Reference

#include <Em10PhysicsList.hh>

Inheritance diagram for Em10PhysicsList:
Collaboration diagram for Em10PhysicsList:

Public Member Functions

 Em10PhysicsList (Em10DetectorConstruction *)
 
 ~Em10PhysicsList ()
 
void ConstructParticle ()
 
void ConstructProcess ()
 
void SetCuts ()
 
void SetGammaCut (G4double)
 
void SetElectronCut (G4double)
 
void SetRegGammaCut (G4double cut)
 
void SetRegElectronCut (G4double cut)
 
void SetRegPositronCut (G4double cut)
 
void SetRadiatorCuts ()
 
void SetDetectorCuts ()
 
void SetMaxStep (G4double)
 
void SetXTRModel (G4String m)
 
- 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 Member Functions

void ConstructBosons ()
 
void ConstructLeptons ()
 
void ConstructMesons ()
 
void ConstructBarions ()
 
void AddParameterisation ()
 
void ConstructGeneral ()
 
void ConstructEM ()
 

Private Attributes

G4double MaxChargedStep
 
G4ForwardXrayTRfForwardXrayTR
 
Em10StepCuttheeminusStepCut
 
Em10StepCuttheeplusStepCut
 
G4double cutForGamma
 
G4double cutForElectron
 
G4double cutForPositron
 
Em10DetectorConstructionpDet
 
Em10PhysicsListMessengerphysicsListMessenger
 
G4ProductionCutsfRadiatorCuts
 
G4ProductionCutsfDetectorCuts
 
G4double fElectronCut
 
G4double fGammaCut
 
G4double fPositronCut
 
G4String fXTRModel
 

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 Em10PhysicsList.hh.

Constructor & Destructor Documentation

◆ Em10PhysicsList()

Em10PhysicsList::Em10PhysicsList ( Em10DetectorConstruction p)

Definition at line 94 of file Em10PhysicsList.cc.

98  fRadiatorCuts(0),fDetectorCuts(0),fXTRModel("transpM")
99 {
100  pDet = p;
101 
102  // world cuts
103 
104  defaultCutValue = 1.*mm;
108 
109  // Region cuts
110 
114 
115  SetVerboseLevel(1);
117 }
G4ProductionCuts * fDetectorCuts
G4double MaxChargedStep
G4double cutForElectron
G4ProductionCuts * fRadiatorCuts
void SetVerboseLevel(G4int value)
Em10DetectorConstruction * pDet
Em10StepCut * theeminusStepCut
Em10StepCut * theeplusStepCut
Em10PhysicsListMessenger * physicsListMessenger
G4double cutForPositron
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:

◆ ~Em10PhysicsList()

Em10PhysicsList::~Em10PhysicsList ( )

Definition at line 119 of file Em10PhysicsList.cc.

120 {
121  delete physicsListMessenger;
122 }
Em10PhysicsListMessenger * physicsListMessenger

Member Function Documentation

◆ AddParameterisation()

void Em10PhysicsList::AddParameterisation ( )
private

◆ ConstructBarions()

void Em10PhysicsList::ConstructBarions ( )
private

Definition at line 179 of file Em10PhysicsList.cc.

180 {
181 // barions
182 
185 }
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 Em10PhysicsList::ConstructBosons ( )
private

Definition at line 141 of file Em10PhysicsList.cc.

142 {
143  // gamma
145 }
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 Em10PhysicsList::ConstructEM ( )
private

Definition at line 198 of file Em10PhysicsList.cc.

199 {
200 
201  // G4cout<<"fMinElectronEnergy = "<<fMinElectronEnergy/keV<<" keV"<<G4endl;
202  // G4cout<<"fMinGammaEnergy = "<<fMinGammaEnergy/keV<<" keV"<<G4endl;
203  G4cout<<"XTR model = "<<fXTRModel<<G4endl;
204 
205  const G4RegionStore* theRegionStore = G4RegionStore::GetInstance();
206  G4Region* gas = theRegionStore->GetRegion("XTRdEdxDetector");
207 
208  G4VXTRenergyLoss* processXTR = 0;
209 
210  if(fXTRModel == "gammaR" )
211  {
212  // G4GammaXTRadiator*
213  processXTR = new G4GammaXTRadiator(pDet->GetLogicalRadiator(),
214  100.,
215  100.,
217  pDet->GetGasMaterial(),
218  pDet->GetFoilThick(),
219  pDet->GetGasThick(),
220  pDet->GetFoilNumber(),
221  "GammaXTRadiator");
222  }
223  else if(fXTRModel == "gammaM" )
224  {
225  // G4XTRGammaRadModel*
226  processXTR = new G4XTRGammaRadModel(pDet->GetLogicalRadiator(),
227  100.,
228  100.,
230  pDet->GetGasMaterial(),
231  pDet->GetFoilThick(),
232  pDet->GetGasThick(),
233  pDet->GetFoilNumber(),
234  "GammaXTRadiator");
235  }
236  else if(fXTRModel == "strawR" )
237  {
238  // G4StrawTubeXTRadiator*
239  processXTR = new G4StrawTubeXTRadiator(pDet->GetLogicalRadiator(),
241  pDet->GetGasMaterial(),
242  0.53, // pDet->GetFoilThick(),
243  3.14159, // pDet->GetGasThick(),
245  true,
246  "strawXTRadiator");
247  }
248  else if(fXTRModel == "regR" )
249  {
250  // G4RegularXTRadiator*
251  processXTR = new G4RegularXTRadiator(pDet->GetLogicalRadiator(),
253  pDet->GetGasMaterial(),
254  pDet->GetFoilThick(),
255  pDet->GetGasThick(),
256  pDet->GetFoilNumber(),
257  "RegularXTRadiator");
258  }
259  else if(fXTRModel == "transpR" )
260  {
261  // G4TransparentRegXTRadiator*
264  pDet->GetGasMaterial(),
265  pDet->GetFoilThick(),
266  pDet->GetGasThick(),
267  pDet->GetFoilNumber(),
268  "RegularXTRadiator");
269  }
270  else if(fXTRModel == "regM" )
271  {
272  // G4XTRRegularRadModel*
273  processXTR = new G4XTRRegularRadModel(pDet->GetLogicalRadiator(),
275  pDet->GetGasMaterial(),
276  pDet->GetFoilThick(),
277  pDet->GetGasThick(),
278  pDet->GetFoilNumber(),
279  "RegularXTRadiator");
280 
281  }
282  else if(fXTRModel == "transpM" )
283  {
284  // G4XTRTransparentRegRadModel*
285  // processXTR = new G4XTRTransparentRegRadModel(pDet->GetLogicalRadiator(),
286  processXTR = new Em10XTRTransparentRegRadModel(pDet->GetLogicalRadiator(),
288  pDet->GetGasMaterial(),
289  pDet->GetFoilThick(),
290  pDet->GetGasThick(),
291  pDet->GetFoilNumber(),
292  "RegularXTRadiator");
293  }
294  else
295  {
296  G4Exception("Invalid XTR model name", "InvalidSetup",
297  FatalException, "XTR model name is out of the name list");
298  }
299  // processXTR->SetCompton(true);
300  processXTR->SetVerboseLevel(1);
301 
303  particleIterator->reset();
304 
305  while( (*particleIterator)() )
306  {
307  G4ParticleDefinition* particle = particleIterator->value();
308  G4ProcessManager* pmanager = particle->GetProcessManager();
309  G4String particleName = particle->GetParticleName();
310 
311  if (particleName == "gamma")
312  {
313  // Construct processes for gamma
314 
317  pmanager->AddDiscreteProcess(new G4GammaConversion);
318 
319  }
320  else if (particleName == "e-")
321  {
322  // Construct processes for electron
325  G4eIonisation* eioni = new G4eIonisation();
326  G4PAIModel* pai = new G4PAIModel(particle,"PAIModel");
327  eioni->AddEmModel(0,pai,pai,gas);
328 
329  pmanager->AddProcess(new G4eMultipleScattering,-1,1,1);
330  pmanager->AddProcess(eioni,-1,2,2);
331  pmanager->AddProcess(new G4eBremsstrahlung,-1,3,3);
332  pmanager->AddDiscreteProcess(processXTR);
335 
336  }
337  else if (particleName == "e+")
338  {
339  // Construct processes for positron
340 
343  G4eIonisation* eioni = new G4eIonisation();
344  G4PAIModel* pai = new G4PAIModel(particle,"PAIModel");
345  eioni->AddEmModel(0,pai,pai,gas);
346 
347  pmanager->AddProcess(new G4eMultipleScattering,-1,1,1);
348  pmanager->AddProcess(eioni,-1,2,2);
349  pmanager->AddProcess(new G4eBremsstrahlung,-1,3,3);
350  pmanager->AddProcess(new G4eplusAnnihilation,0,-1,4);
351  pmanager->AddDiscreteProcess(processXTR);
354 
355  }
356  else if( particleName == "mu+" ||
357  particleName == "mu-" )
358  {
359  // Construct processes for muon+
360 
361  Em10StepCut* muonStepCut = new Em10StepCut();
362  muonStepCut->SetMaxStep(MaxChargedStep) ;
363 
364  G4MuIonisation* muioni = new G4MuIonisation() ;
365 
366  G4PAIModel* pai = new G4PAIModel(particle,"PAIModel");
367  muioni->AddEmModel(0,pai,pai,gas);
368 
369  pmanager->AddProcess(new G4MuMultipleScattering(),-1,1,1);
370  pmanager->AddProcess(muioni,-1,2,2);
371  pmanager->AddProcess(new G4MuBremsstrahlung(),-1,3,3);
372  pmanager->AddProcess(new G4MuPairProduction(),-1,4,4);
373  pmanager->AddProcess( muonStepCut,-1,-1,5);
374 
375  }
376  else if (
377  particleName == "proton"
378  || particleName == "antiproton"
379  || particleName == "pi+"
380  || particleName == "pi-"
381  || particleName == "kaon+"
382  || particleName == "kaon-"
383  )
384  {
385  Em10StepCut* thehadronStepCut = new Em10StepCut();
386  thehadronStepCut->SetMaxStep(MaxChargedStep) ;
387 
388  G4hIonisation* thehIonisation = new G4hIonisation();
389  G4PAIModel* pai = new G4PAIModel(particle,"PAIModel");
390  thehIonisation->AddEmModel(0,pai,pai,gas);
391 
392  pmanager->AddProcess(new G4hMultipleScattering,-1,1,1);
393  pmanager->AddProcess(thehIonisation,-1,2,2);
394  pmanager->AddProcess( thehadronStepCut,-1,-1,3);
395 
396  }
397  }
398  G4EmProcessOptions opt;
399  opt.SetApplyCuts(true);
400 }
void SetMaxStep(G4double)
Definition: Em10StepCut.cc:56
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4double MaxChargedStep
G4ProcessManager * GetProcessManager() const
static G4RegionStore * GetInstance()
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
Em10DetectorConstruction * pDet
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Em10StepCut * theeminusStepCut
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
Em10StepCut * theeplusStepCut
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
#define G4endl
Definition: G4ios.hh:61
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4LogicalVolume * GetLogicalRadiator()
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
void SetApplyCuts(G4bool val)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructGeneral()

void Em10PhysicsList::ConstructGeneral ( )
private

Definition at line 404 of file Em10PhysicsList.cc.

405 {
406  // Add Decay Process
407 
408  G4Decay* theDecayProcess = new G4Decay();
410  particleIterator->reset();
411 
412  while( (*particleIterator)() )
413  {
414  G4ParticleDefinition* particle = particleIterator->value();
415  G4ProcessManager* pmanager = particle->GetProcessManager();
416 
417  if (theDecayProcess->IsApplicable(*particle))
418  {
419  pmanager ->AddProcess(theDecayProcess);
420 
421  // set ordering for PostStepDoIt and AtRestDoIt
422 
423  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
424  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
425  }
426  }
427 }
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 Em10PhysicsList::ConstructLeptons ( )
private

Definition at line 149 of file Em10PhysicsList.cc.

150 {
151  // leptons
152 
157 
162 }
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 Em10PhysicsList::ConstructMesons ( )
private

Definition at line 166 of file Em10PhysicsList.cc.

167 {
168  // mesons
169 
175 }
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 Em10PhysicsList::ConstructParticle ( void  )
virtual

Reimplemented from G4VModularPhysicsList.

Definition at line 126 of file Em10PhysicsList.cc.

127 {
128  // In this method, static member functions should be called
129  // for all particles which you want to use.
130  // This ensures that objects of these particle types will be
131  // created in the program.
132 
133  ConstructBosons();
135  ConstructMesons();
137 }
Here is the call graph for this function:

◆ ConstructProcess()

void Em10PhysicsList::ConstructProcess ( void  )
virtual

Reimplemented from G4VModularPhysicsList.

Definition at line 189 of file Em10PhysicsList.cc.

190 {
192  ConstructEM();
194 }
Here is the call graph for this function:

◆ SetCuts()

void Em10PhysicsList::SetCuts ( )
virtual

Reimplemented from G4VUserPhysicsList.

Definition at line 431 of file Em10PhysicsList.cc.

432 {
433  // set cut values for gamma at first and for e- second and next for e+,
434  // because some processes for e+/e- need cut values for gamma
435 
436  SetCutValue(cutForGamma, "gamma", "DefaultRegionForTheWorld");
437  SetCutValue(cutForElectron, "e-", "DefaultRegionForTheWorld");
438  SetCutValue(cutForPositron, "e+", "DefaultRegionForTheWorld");
439 
440  if (verboseLevel > 0)
441  {
442  G4cout << "Em10PhysicsList::SetCuts:";
443  G4cout << "CutLength for e-, e+ and gamma is: "
444  << G4BestUnit(defaultCutValue,"Length") << G4endl;
445  }
446 
448 
449  G4Region* region = (G4RegionStore::GetInstance())->GetRegion("XTRradiator");
451  G4cout << "Radiator cuts are set" << G4endl;
452 
454  region = (G4RegionStore::GetInstance())->GetRegion("XTRdEdxDetector");
456  G4cout << "Detector cuts are set" << G4endl;
457 
458  if (verboseLevel > 1) DumpCutValuesTable();
459 }
G4ProductionCuts * fDetectorCuts
void SetCutValue(G4double aCut, const G4String &pname)
G4double cutForElectron
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4ProductionCuts * fRadiatorCuts
static G4RegionStore * GetInstance()
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void SetProductionCuts(G4ProductionCuts *cut)
G4double cutForPositron
Here is the call graph for this function:

◆ SetDetectorCuts()

void Em10PhysicsList::SetDetectorCuts ( )

Definition at line 501 of file Em10PhysicsList.cc.

502 {
504 
508 
509  G4cout<<"Detector gamma cut = "<<fGammaCut/mm<<" mm"<<G4endl;
510  G4cout<<"Detector electron cut = "<<fElectronCut/mm<<" mm"<<G4endl;
511  G4cout<<"Detector positron cut = "<<fPositronCut/mm<<" mm"<<G4endl;
512 
513 }
G4ProductionCuts * fDetectorCuts
void SetProductionCut(G4double cut, G4int index=-1)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetElectronCut()

void Em10PhysicsList::SetElectronCut ( G4double  val)

Definition at line 470 of file Em10PhysicsList.cc.

471 {
472  cutForElectron = val;
473 }
G4double cutForElectron
Here is the caller graph for this function:

◆ SetGammaCut()

void Em10PhysicsList::SetGammaCut ( G4double  val)

Definition at line 463 of file Em10PhysicsList.cc.

464 {
465  cutForGamma = val;
466 }
Here is the caller graph for this function:

◆ SetMaxStep()

void Em10PhysicsList::SetMaxStep ( G4double  step)

Definition at line 477 of file Em10PhysicsList.cc.

478 {
479  MaxChargedStep = step ;
480  G4cout << " MaxChargedStep=" << MaxChargedStep << G4endl;
481  G4cout << G4endl;
482 }
G4double MaxChargedStep
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the caller graph for this function:

◆ SetRadiatorCuts()

void Em10PhysicsList::SetRadiatorCuts ( )

Definition at line 486 of file Em10PhysicsList.cc.

487 {
489 
493 
494  G4cout<<"Radiator gamma cut = "<<fGammaCut/mm<<" mm"<<G4endl;
495  G4cout<<"Radiator electron cut = "<<fElectronCut/mm<<" mm"<<G4endl;
496  G4cout<<"Radiator positron cut = "<<fPositronCut/mm<<" mm"<<G4endl;
497 }
void SetProductionCut(G4double cut, G4int index=-1)
G4ProductionCuts * fRadiatorCuts
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetRegElectronCut()

void Em10PhysicsList::SetRegElectronCut ( G4double  cut)
inline

Definition at line 79 of file Em10PhysicsList.hh.

79 {fElectronCut = cut;};
Here is the caller graph for this function:

◆ SetRegGammaCut()

void Em10PhysicsList::SetRegGammaCut ( G4double  cut)
inline

Definition at line 78 of file Em10PhysicsList.hh.

78 {fGammaCut = cut;};
Here is the caller graph for this function:

◆ SetRegPositronCut()

void Em10PhysicsList::SetRegPositronCut ( G4double  cut)
inline

Definition at line 80 of file Em10PhysicsList.hh.

80 {fPositronCut = cut;};
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetXTRModel()

void Em10PhysicsList::SetXTRModel ( G4String  m)
inline

Definition at line 86 of file Em10PhysicsList.hh.

G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static const double m
Definition: G4SIunits.hh:128
Here is the caller graph for this function:

Member Data Documentation

◆ cutForElectron

G4double Em10PhysicsList::cutForElectron
private

Definition at line 98 of file Em10PhysicsList.hh.

◆ cutForGamma

G4double Em10PhysicsList::cutForGamma
private

Definition at line 97 of file Em10PhysicsList.hh.

◆ cutForPositron

G4double Em10PhysicsList::cutForPositron
private

Definition at line 98 of file Em10PhysicsList.hh.

◆ fDetectorCuts

G4ProductionCuts* Em10PhysicsList::fDetectorCuts
private

Definition at line 105 of file Em10PhysicsList.hh.

◆ fElectronCut

G4double Em10PhysicsList::fElectronCut
private

Definition at line 106 of file Em10PhysicsList.hh.

◆ fForwardXrayTR

G4ForwardXrayTR* Em10PhysicsList::fForwardXrayTR
private

Definition at line 92 of file Em10PhysicsList.hh.

◆ fGammaCut

G4double Em10PhysicsList::fGammaCut
private

Definition at line 106 of file Em10PhysicsList.hh.

◆ fPositronCut

G4double Em10PhysicsList::fPositronCut
private

Definition at line 106 of file Em10PhysicsList.hh.

◆ fRadiatorCuts

G4ProductionCuts* Em10PhysicsList::fRadiatorCuts
private

Definition at line 104 of file Em10PhysicsList.hh.

◆ fXTRModel

G4String Em10PhysicsList::fXTRModel
private

Definition at line 107 of file Em10PhysicsList.hh.

◆ MaxChargedStep

G4double Em10PhysicsList::MaxChargedStep
private

Definition at line 86 of file Em10PhysicsList.hh.

◆ pDet

Em10DetectorConstruction* Em10PhysicsList::pDet
private

Definition at line 100 of file Em10PhysicsList.hh.

◆ physicsListMessenger

Em10PhysicsListMessenger* Em10PhysicsList::physicsListMessenger
private

Definition at line 102 of file Em10PhysicsList.hh.

◆ theeminusStepCut

Em10StepCut* Em10PhysicsList::theeminusStepCut
private

Definition at line 94 of file Em10PhysicsList.hh.

◆ theeplusStepCut

Em10StepCut* Em10PhysicsList::theeplusStepCut
private

Definition at line 95 of file Em10PhysicsList.hh.


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