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

#include <PhysListEmStandard_option3.hh>

Inheritance diagram for PhysListEmStandard_option3:
Collaboration diagram for PhysListEmStandard_option3:

Public Member Functions

 PhysListEmStandard_option3 (const G4String &name, DetectorConstruction *det)
 
 ~PhysListEmStandard_option3 ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
 PhysListEmStandard_option3 (const G4String &name)
 
 ~PhysListEmStandard_option3 ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
- Public Member Functions inherited from G4VPhysicsConstructor
 G4VPhysicsConstructor (const G4String &="")
 
 G4VPhysicsConstructor (const G4String &name, G4int physics_type)
 
virtual ~G4VPhysicsConstructor ()
 
void SetPhysicsName (const G4String &="")
 
const G4StringGetPhysicsName () const
 
void SetPhysicsType (G4int)
 
G4int GetPhysicsType () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4int GetInstanceID () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPhysicsConstructor
static const G4VPCManagerGetSubInstanceManager ()
 
- Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 
- Protected Attributes inherited from G4VPhysicsConstructor
G4int verboseLevel
 
G4String namePhysics
 
G4int typePhysics
 
G4ParticleTabletheParticleTable
 
G4int g4vpcInstanceID
 
- Static Protected Attributes inherited from G4VPhysicsConstructor
static G4RUN_DLL G4VPCManager subInstanceManager
 

Detailed Description

Definition at line 44 of file PhysListEmStandard_option3.hh.

Constructor & Destructor Documentation

PhysListEmStandard_option3::PhysListEmStandard_option3 ( const G4String name,
DetectorConstruction det 
)

Definition at line 69 of file PhysListEmStandard_option3.cc.

71 : G4VPhysicsConstructor(name), fDetector(det)
72 {
74  param->SetDefaults();
75  param->SetVerbose(1);
76  param->SetMinEnergy(100*eV);
77  param->SetMaxEnergy(10*GeV);
78  param->SetNumberOfBinsPerDecade(20);
79  param->SetLowestElectronEnergy(10*eV);
80  param->SetBuildCSDARange(true);
81  param->SetMaxEnergyForCSDARange(10*GeV);
84 }
void SetVerbose(G4int val)
void SetLowestElectronEnergy(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMaxEnergyForCSDARange(G4double val)
void SetMaxEnergy(G4double val)
void SetNumberOfBinsPerDecade(G4int val)
static constexpr double eV
Definition: G4SIunits.hh:215
void SetBuildCSDARange(G4bool val)
void SetMinEnergy(G4double val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
static constexpr double GeV
Definition: G4SIunits.hh:217

Here is the call graph for this function:

PhysListEmStandard_option3::~PhysListEmStandard_option3 ( )

Definition at line 88 of file PhysListEmStandard_option3.cc.

89 {}
PhysListEmStandard_option3::PhysListEmStandard_option3 ( const G4String name)

Definition at line 61 of file PhysListEmStandard_option3.cc.

63 {
65  param->SetDefaults();
66  param->SetVerbose(1);
67  param->SetMinEnergy(100*eV);
68  param->SetMaxEnergy(10*GeV);
69  param->SetNumberOfBinsPerDecade(20);
70  param->SetLowestElectronEnergy(10*eV);
71  param->SetBuildCSDARange(true);
72  param->SetMaxEnergyForCSDARange(10*GeV);
75 }
void SetVerbose(G4int val)
void SetLowestElectronEnergy(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMaxEnergyForCSDARange(G4double val)
void SetMaxEnergy(G4double val)
void SetNumberOfBinsPerDecade(G4int val)
static constexpr double eV
Definition: G4SIunits.hh:215
void SetBuildCSDARange(G4bool val)
void SetMinEnergy(G4double val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
static constexpr double GeV
Definition: G4SIunits.hh:217

Here is the call graph for this function:

PhysListEmStandard_option3::~PhysListEmStandard_option3 ( )

Member Function Documentation

virtual void PhysListEmStandard_option3::ConstructParticle ( void  )
inlinevirtual

Implements G4VPhysicsConstructor.

Definition at line 50 of file PhysListEmStandard_option3.hh.

50 {};
virtual void PhysListEmStandard_option3::ConstructParticle ( void  )
inlinevirtual

Implements G4VPhysicsConstructor.

Definition at line 53 of file PhysListEmStandard_option3.hh.

53 {};
virtual void PhysListEmStandard_option3::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

void PhysListEmStandard_option3::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 93 of file PhysListEmStandard_option3.cc.

94 {
95  // Add standard EM Processes
96  //
97 
99  particleIterator->reset();
100  while( (*particleIterator)() ){
101  G4ParticleDefinition* particle = particleIterator->value();
102  G4ProcessManager* pmanager = particle->GetProcessManager();
103  G4String particleName = particle->GetParticleName();
104 
105  if (particleName == "gamma") {
106  // gamma
107 
108  G4ComptonScattering* compton = new G4ComptonScattering();
109  MyKleinNishinaCompton* comptonModel =
110  new MyKleinNishinaCompton(fDetector);
111  comptonModel->SetCSFactor(1000.);
112  compton->SetEmModel(comptonModel );
113 
115  pmanager->AddDiscreteProcess(compton);
116  pmanager->AddDiscreteProcess(new G4GammaConversion);
117 
118  } else if (particleName == "e-") {
119  //electron
121 
122  G4eIonisation* eIoni = new G4eIonisation();
123  eIoni->SetEmModel(new MyMollerBhabhaModel);
124  eIoni->SetStepFunction(0.2, 10*um);
125 
126  pmanager->AddProcess(msc, -1, 1, 1);
127  pmanager->AddProcess(eIoni, -1, 2, 2);
128 
129  } else if (particleName == "e+") {
130  //positron
132 
133  G4eIonisation* pIoni = new G4eIonisation();
134  pIoni->SetEmModel(new MyMollerBhabhaModel);
135 
136  pmanager->AddProcess(msc, -1, 1, 1);
137  pmanager->AddProcess(pIoni, -1, 2, 2);
138  pmanager->AddProcess(new G4eplusAnnihilation, 0,-1, 3);
139 
140  } else if( particleName == "proton" ) {
141  //proton
142  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
143  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
144  }
145  }
146 }
void SetCSFactor(G4double factor)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
const G4String & GetParticleName() const
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetEmModel(G4VEmModel *, G4int index=1)
static constexpr double um
Definition: G4SIunits.hh:113
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
G4ProcessManager * GetProcessManager() const
void SetEmModel(G4VEmModel *, G4int index=1)
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65

Here is the call graph for this function:


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