Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmStandardPhysicsSS Class Reference

#include <G4EmStandardPhysicsSS.hh>

Inheritance diagram for G4EmStandardPhysicsSS:

Public Member Functions

 G4EmStandardPhysicsSS (G4int ver=1)
 ~G4EmStandardPhysicsSS () override
void ConstructParticle () override
void ConstructProcess () override
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
G4int GetInstanceID () const
virtual void TerminateWorker ()
void SetVerboseLevel (G4int value)
G4int GetVerboseLevel () const

Additional Inherited Members

Static Public Member Functions inherited from G4VPhysicsConstructor
static const G4VPCManagerGetSubInstanceManager ()
Protected Types inherited from G4VPhysicsConstructor
using PhysicsBuilder_V = G4VPCData::PhysicsBuilders_V
Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
PhysicsBuilder_V GetBuilders () const
void AddBuilder (G4PhysicsBuilderInterface *bld)
Protected Attributes inherited from G4VPhysicsConstructor
G4int verboseLevel = 0
G4String namePhysics = ""
G4int typePhysics = 0
G4ParticleTabletheParticleTable = nullptr
G4int g4vpcInstanceID = 0
Static Protected Attributes inherited from G4VPhysicsConstructor
static G4RUN_DLL G4VPCManager subInstanceManager

Detailed Description

Definition at line 50 of file G4EmStandardPhysicsSS.hh.

Constructor & Destructor Documentation

◆ G4EmStandardPhysicsSS()

G4EmStandardPhysicsSS::G4EmStandardPhysicsSS ( G4int ver = 1)
explicit

Definition at line 99 of file G4EmStandardPhysicsSS.cc.

100 : G4VPhysicsConstructor("G4EmStandardSS")
101{
102 SetVerboseLevel(ver);
103 G4EmParameters* param = G4EmParameters::Instance();
104 param->SetDefaults();
105 param->SetVerbose(ver);
106 param->SetLowestElectronEnergy(10*CLHEP::eV);
107 param->SetMscThetaLimit(0.0);
108 param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
109 param->SetAuger(true);
110 param->SetPixe(true);
112}
@ bElectromagnetic
void SetLowestElectronEnergy(G4double val)
static G4EmParameters * Instance()
void SetMscThetaLimit(G4double val)
void SetVerbose(G4int val)
void SetPixe(G4bool val)
void SetAuger(G4bool val)
void SetUseMottCorrection(G4bool val)
G4VPhysicsConstructor(const G4String &="")
void SetVerboseLevel(G4int value)

◆ ~G4EmStandardPhysicsSS()

G4EmStandardPhysicsSS::~G4EmStandardPhysicsSS ( )
override

Definition at line 116 of file G4EmStandardPhysicsSS.cc.

117{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysicsSS::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 121 of file G4EmStandardPhysicsSS.cc.

122{
123 // minimal set of particles for EM physics
125}
static void ConstructMinimalEmSet()

◆ ConstructProcess()

void G4EmStandardPhysicsSS::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 129 of file G4EmStandardPhysicsSS.cc.

130{
131 if(verboseLevel > 1) {
132 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
133 }
135
136 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
137 G4EmParameters* param = G4EmParameters::Instance();
138
139 // processes used by several particles
140 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
141
142 // Add gamma EM Processes
143 G4ParticleDefinition* particle = G4Gamma::Gamma();
144
145 // Photoelectric
146 G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
147 G4VEmModel* peModel = new G4LivermorePhotoElectricModel();
148 pe->SetEmModel(peModel);
149 if(param->EnablePolarisation()) {
150 peModel->SetAngularDistribution(new G4PhotoElectricAngularGeneratorPolarized());
151 }
152
153 // Compton scattering
154 G4ComptonScattering* cs = new G4ComptonScattering;
155 cs->SetEmModel(new G4KleinNishinaModel());
156
157 // Gamma conversion
158 G4GammaConversion* gc = new G4GammaConversion();
159 G4VEmModel* conv = new G4BetheHeitler5DModel();
160 gc->SetEmModel(conv);
161
162 // default Rayleigh scattering is Livermore
163 G4RayleighScattering* rl = new G4RayleighScattering();
164 if(param->EnablePolarisation()) {
165 rl->SetEmModel(new G4LivermorePolarizedRayleighModel());
166 }
167
168 if(param->GeneralProcessActive()) {
169 G4GammaGeneralProcess* sp = new G4GammaGeneralProcess();
170 sp->AddEmProcess(pe);
171 sp->AddEmProcess(cs);
172 sp->AddEmProcess(gc);
173 sp->AddEmProcess(rl);
175 ph->RegisterProcess(sp, particle);
176 } else {
177 ph->RegisterProcess(pe, particle);
178 ph->RegisterProcess(cs, particle);
179 ph->RegisterProcess(gc, particle);
180 ph->RegisterProcess(rl, particle);
181 }
182 // e-
183 particle = G4Electron::Electron();
184
185 G4VEmModel* ss = nullptr;
186 if(param->UseMottCorrection()) {
187 ss = new G4eDPWACoulombScatteringModel();
188 } else {
189 ss = new G4eCoulombScatteringModel(false);
190 }
191 ph->RegisterProcess(new G4eIonisation(), particle);
192 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
193
194 G4ePairProduction* ee = new G4ePairProduction();
195 ph->RegisterProcess(ee, particle);
197
198 // e+
199 particle = G4Positron::Positron();
200
201 if(param->UseMottCorrection()) {
202 ss = new G4eDPWACoulombScatteringModel();
203 } else {
204 ss = new G4eCoulombScatteringModel(false);
205 }
206 ph->RegisterProcess(new G4eIonisation(), particle);
207 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
208 ph->RegisterProcess(ee, particle);
209 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
211
212 // generic ion
213 particle = G4GenericIon::GenericIon();
214 G4ionIonisation* ionIoni = new G4ionIonisation();
215 auto fluc = new G4IonFluctuations();
216 ionIoni->SetFluctModel(fluc);
217 ionIoni->SetEmModel(new G4LindhardSorensenIonModel());
218 ph->RegisterProcess(ionIoni, particle);
219 ph->RegisterProcess(new G4CoulombScattering(false), particle);
220
221 // muons, hadrons, ions
223
224 // extra configuration
225 G4EmModelActivator mact(GetPhysicsName());
226}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition G4Electron.cc:91
static void ConstructChargedSS(G4hMultipleScattering *hmsc)
static void ConstructElectronSSProcess(G4VEmModel *ss, G4ParticleDefinition *particle)
static void PrepareEMPhysics()
G4bool EnablePolarisation() const
G4bool UseMottCorrection() const
G4bool GeneralProcessActive() const
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
static G4LossTableManager * Instance()
void SetGammaGeneralProcess(G4VEmProcess *)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition G4Positron.cc:90
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const

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