Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
PhysListEmStandard_GS Class Reference

#include <PhysListEmStandard_GS.hh>

Inheritance diagram for PhysListEmStandard_GS:
G4VPhysicsConstructor G4VPhysicsConstructor

Public Member Functions

 PhysListEmStandard_GS (const G4String &name, DetectorConstruction *det)
 
 ~PhysListEmStandard_GS ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
 PhysListEmStandard_GS (const G4String &name)
 
 ~PhysListEmStandard_GS ()
 
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)
 
- 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 include/PhysListEmStandard_GS.hh.

Constructor & Destructor Documentation

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

Definition at line 63 of file src/PhysListEmStandard_GS.cc.

65 : G4VPhysicsConstructor(name), fDetector(det)
66 {}
G4VPhysicsConstructor(const G4String &="")
PhysListEmStandard_GS::~PhysListEmStandard_GS ( )

Definition at line 70 of file src/PhysListEmStandard_GS.cc.

71 {}
PhysListEmStandard_GS::PhysListEmStandard_GS ( const G4String name)

Definition at line 61 of file /src/PhysListEmStandard_GS.cc.

63 {}
G4VPhysicsConstructor(const G4String &="")
PhysListEmStandard_GS::~PhysListEmStandard_GS ( )

Member Function Documentation

virtual void PhysListEmStandard_GS::ConstructParticle ( void  )
inlinevirtual

Implements G4VPhysicsConstructor.

Definition at line 50 of file /include/PhysListEmStandard_GS.hh.

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

Implements G4VPhysicsConstructor.

Definition at line 52 of file include/PhysListEmStandard_GS.hh.

52 {};
virtual void PhysListEmStandard_GS::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

void PhysListEmStandard_GS::ConstructProcess ( void  )
virtual

pmanager->AddProcess(new G4eBremsstrahlung, -1, 3, 3);

pmanager->AddProcess(new G4eBremsstrahlung, -1, 3, 3);

Implements G4VPhysicsConstructor.

Definition at line 75 of file src/PhysListEmStandard_GS.cc.

References G4ProcessManager::AddDiscreteProcess(), G4VMultipleScattering::AddEmModel(), G4ProcessManager::AddProcess(), aParticleIterator, python.hepunit::eV, fUseDistanceToBoundary, G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetProcessManager(), python.hepunit::GeV, G4EmProcessOptions::SetBuildCSDARange(), MyKleinNishinaCompton::SetCSFactor(), G4EmProcessOptions::SetDEDXBinning(), G4EmProcessOptions::SetDEDXBinningForCSDARange(), G4VEmProcess::SetEmModel(), G4VEnergyLossProcess::SetEmModel(), G4EmProcessOptions::SetLambdaBinning(), G4EmProcessOptions::SetMaxEnergy(), G4EmProcessOptions::SetMaxEnergyForCSDARange(), G4EmProcessOptions::SetMinEnergy(), G4EmProcessOptions::SetMscStepLimitation(), and G4EmProcessOptions::SetStepFunction().

76 {
77  // Add standard EM Processes
78  //
79 
80  aParticleIterator->reset();
81  while( (*aParticleIterator)() ){
82  G4ParticleDefinition* particle = aParticleIterator->value();
83  G4ProcessManager* pmanager = particle->GetProcessManager();
84  G4String particleName = particle->GetParticleName();
85 
86  if (particleName == "gamma") {
87  // gamma
88 
90  MyKleinNishinaCompton* comptonModel = new MyKleinNishinaCompton(fDetector);
91  comptonModel->SetCSFactor(1000.);
92  compton->SetEmModel(comptonModel );
93 
95  pmanager->AddDiscreteProcess(compton);
97 
98  } else if (particleName == "e-") {
99  //electron
100 
103 
104  G4eIonisation* eIoni = new G4eIonisation();
105  eIoni->SetEmModel(new MyMollerBhabhaModel);
106 
107  pmanager->AddProcess(eMsc, -1, 1, 1);
108  pmanager->AddProcess(eIoni, -1, 2, 2);
109 /// pmanager->AddProcess(new G4eBremsstrahlung, -1, 3, 3);
110 
111  } else if (particleName == "e+") {
112  //positron
113 
116 
117  G4eIonisation* pIoni = new G4eIonisation();
118  pIoni->SetEmModel(new MyMollerBhabhaModel);
119 
120  pmanager->AddProcess(pMsc, -1, 1, 1);
121  pmanager->AddProcess(pIoni, -1, 2, 2);
122 /// pmanager->AddProcess(new G4eBremsstrahlung, -1, 3, 3);
123  pmanager->AddProcess(new G4eplusAnnihilation, 0,-1, 3);
124 
125  } else if( particleName == "proton" ) {
126  //proton
127  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
128  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
129  }
130  }
131 
132  // Em options
133  //
134  // Main options and setting parameters are shown here.
135  // Several of them have default values.
136  //
137  G4EmProcessOptions emOptions;
138 
139  //physics tables
140  //
141  emOptions.SetMinEnergy(100*eV); //default
142  emOptions.SetMaxEnergy(10*GeV); //default
143  emOptions.SetDEDXBinning(8*20); //default=8*7
144  emOptions.SetLambdaBinning(8*20); //default=8*7
145 
146  //multiple coulomb scattering
147  //
148  emOptions.SetMscStepLimitation(fUseDistanceToBoundary); //default=fUseSafety
149 
150  //energy loss
151  //
152  emOptions.SetStepFunction(0.2, 10*um); //default=(0.2, 1*mm)
153 
154  //build CSDA range
155  //
156  emOptions.SetBuildCSDARange(true); //default=false
157  emOptions.SetMaxEnergyForCSDARange(10*GeV);
158  emOptions.SetDEDXBinningForCSDARange(8*7); //default=8*7
159 }
void SetCSFactor(G4double factor)
void SetMinEnergy(G4double val)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void SetStepFunction(G4double v1, G4double v2)
void SetDEDXBinningForCSDARange(G4int val)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
void SetMaxEnergyForCSDARange(G4double val)
void SetDEDXBinning(G4int val)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetLambdaBinning(G4int val)
#define aParticleIterator
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetMaxEnergy(G4double val)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetBuildCSDARange(G4bool val)
void SetMscStepLimitation(G4MscStepLimitType val)

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