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

#include <G4EmStandardPhysics_option3.hh>

+ Inheritance diagram for G4EmStandardPhysics_option3:

Public Member Functions

 G4EmStandardPhysics_option3 (G4int ver=1, const G4String &name="")
 
 ~G4EmStandardPhysics_option3 () 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 51 of file G4EmStandardPhysics_option3.hh.

Constructor & Destructor Documentation

◆ G4EmStandardPhysics_option3()

G4EmStandardPhysics_option3::G4EmStandardPhysics_option3 ( G4int ver = 1,
const G4String & name = "" )
explicit

Definition at line 101 of file G4EmStandardPhysics_option3.cc.

103 : G4VPhysicsConstructor("G4EmStandard_opt3")
104{
105 SetVerboseLevel(ver);
106 G4EmParameters* param = G4EmParameters::Instance();
107 param->SetDefaults();
108 param->SetVerbose(ver);
109 param->SetGeneralProcessActive(true);
110 param->SetMinEnergy(10*CLHEP::eV);
111 param->SetLowestElectronEnergy(100*CLHEP::eV);
112 param->SetNumberOfBinsPerDecade(20);
114 param->SetUseMottCorrection(true);
115 param->SetStepFunction(0.2, 100*CLHEP::um);
116 param->SetStepFunctionMuHad(0.2, 50*CLHEP::um);
117 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
118 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
120 param->SetMuHadLateralDisplacement(true);
121 param->SetLateralDisplacementAlg96(true);
122 param->SetUseICRU90Data(true);
124 param->SetFluo(true);
125 param->SetMaxNIELEnergy(1*CLHEP::MeV);
128}
@ bElectromagnetic
@ fUrbanFluctuation
@ fAllisonPositronium
@ fUseDistanceToBoundary
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetGeneralProcessActive(G4bool val)
void SetLateralDisplacementAlg96(G4bool val)
void SetPositronAtRestModelType(G4PositronAtRestModelType val)
void SetMuHadLateralDisplacement(G4bool val)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetStepFunction(G4double v1, G4double v2)
void SetFluo(G4bool val)
void SetFluctuationType(G4EmFluctuationType val)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetVerbose(G4int val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)
G4VPhysicsConstructor(const G4String &="")
void SetVerboseLevel(G4int value)

◆ ~G4EmStandardPhysics_option3()

G4EmStandardPhysics_option3::~G4EmStandardPhysics_option3 ( )
overridedefault

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option3::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 136 of file G4EmStandardPhysics_option3.cc.

137{
138 // minimal set of particles for EM physics
140}
static void ConstructMinimalEmSet()

◆ ConstructProcess()

void G4EmStandardPhysics_option3::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 144 of file G4EmStandardPhysics_option3.cc.

145{
146 if(verboseLevel > 1) {
147 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
148 }
150
151 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
152 G4EmParameters* param = G4EmParameters::Instance();
153
154 // processes used by several particles
155 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
156
157 // nuclear stopping is enabled if th eenergy limit above zero
158 G4double nielEnergyLimit = param->MaxNIELEnergy();
159 G4NuclearStopping* pnuc = nullptr;
160 if(nielEnergyLimit > 0.0) {
161 pnuc = new G4NuclearStopping();
162 pnuc->SetMaxKinEnergy(nielEnergyLimit);
163 }
164
165 // Add gamma EM Processes
166 G4ParticleDefinition* particle = G4Gamma::Gamma();
167
168 G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
169 G4VEmModel* peModel = new G4LivermorePhotoElectricModel();
170 pe->SetEmModel(peModel);
171 if(param->EnablePolarisation()) {
172 peModel->SetAngularDistribution(new G4PhotoElectricAngularGeneratorPolarized());
173 }
174
175 G4ComptonScattering* cs = new G4ComptonScattering();
176 cs->SetEmModel(new G4KleinNishinaModel());
177
178 G4GammaConversion* gc = new G4GammaConversion();
179 if(param->EnablePolarisation()) {
180 gc->SetEmModel(new G4BetheHeitler5DModel());
181 }
182
183 G4RayleighScattering* rl = new G4RayleighScattering();
184 if(param->EnablePolarisation()) {
185 rl->SetEmModel(new G4LivermorePolarizedRayleighModel());
186 }
187
188 if(G4EmParameters::Instance()->GeneralProcessActive()) {
189 G4GammaGeneralProcess* sp = new G4GammaGeneralProcess();
190 sp->AddEmProcess(pe);
191 sp->AddEmProcess(cs);
192 sp->AddEmProcess(gc);
193 sp->AddEmProcess(rl);
195 ph->RegisterProcess(sp, particle);
196 } else {
197 ph->RegisterProcess(pe, particle);
198 ph->RegisterProcess(cs, particle);
199 ph->RegisterProcess(gc, particle);
200 ph->RegisterProcess(rl, particle);
201 }
202
203 // e-
204 particle = G4Electron::Electron();
205
206 G4UrbanMscModel* msc1 = new G4UrbanMscModel();
207 G4EmBuilder::ConstructElectronMscProcess(msc1, nullptr, particle);
208
209 G4eIonisation* eIoni = new G4eIonisation();
210
211 G4eBremsstrahlung* brem = new G4eBremsstrahlung();
212 G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel();
213 G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel();
214 br1->SetAngularDistribution(new G4Generator2BS());
215 br2->SetAngularDistribution(new G4Generator2BS());
216 brem->SetEmModel(br1);
217 brem->SetEmModel(br2);
218 br2->SetLowEnergyLimit(CLHEP::GeV);
219
220 G4ePairProduction* ee = new G4ePairProduction();
221
222 ph->RegisterProcess(eIoni, particle);
223 ph->RegisterProcess(brem, particle);
224 ph->RegisterProcess(ee, particle);
225
226 // e+
227 particle = G4Positron::Positron();
228
229 msc1 = new G4UrbanMscModel();
230 G4EmBuilder::ConstructElectronMscProcess(msc1, nullptr, particle);
231
232 eIoni = new G4eIonisation();
233
234 brem = new G4eBremsstrahlung();
235 br1 = new G4SeltzerBergerModel();
236 br2 = new G4eBremsstrahlungRelModel();
237 br1->SetAngularDistribution(new G4Generator2BS());
238 br2->SetAngularDistribution(new G4Generator2BS());
239 brem->SetEmModel(br1);
240 brem->SetEmModel(br2);
241 br2->SetLowEnergyLimit(CLHEP::GeV);
242
243 ph->RegisterProcess(eIoni, particle);
244 ph->RegisterProcess(brem, particle);
245
246 // annihilation
247 auto anni = new G4eplusAnnihilation();
248 if (param->Use3GammaAnnihilationOnFly()) {
249 anni->SetEmModel(new G4eplusTo2or3GammaModel());
250 }
251 ph->RegisterProcess(ee, particle);
252 ph->RegisterProcess(anni, particle);
253
254 // generic ion
255 particle = G4GenericIon::GenericIon();
256 G4ionIonisation* ionIoni = new G4ionIonisation();
257 auto fluc = new G4IonFluctuations();
258 ionIoni->SetFluctModel(fluc);
259 ionIoni->SetEmModel(new G4LindhardSorensenIonModel());
260 ph->RegisterProcess(hmsc, particle);
261 ph->RegisterProcess(ionIoni, particle);
262 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
263
264 // muons, hadrons, ions
265 G4EmBuilder::ConstructCharged(hmsc, pnuc, false);
266
267 // extra configuration
268 G4EmModelActivator mact(GetPhysicsName());
269}
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition G4Electron.cc:91
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
static void ConstructElectronMscProcess(G4VMscModel *msc1, G4VMscModel *msc2, G4ParticleDefinition *particle)
static void PrepareEMPhysics()
G4bool EnablePolarisation() const
G4double MaxNIELEnergy() const
G4bool Use3GammaAnnihilationOnFly() 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 SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const

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