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

#include <G4EmStandardPhysics_option4.hh>

+ Inheritance diagram for G4EmStandardPhysics_option4:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmStandardPhysics_option4()

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

Definition at line 108 of file G4EmStandardPhysics_option4.cc.

110 : G4VPhysicsConstructor("G4EmStandard_opt4")
111{
112 SetVerboseLevel(ver);
113 G4EmParameters* param = G4EmParameters::Instance();
114 param->SetDefaults();
115 param->SetVerbose(ver);
116 param->SetGeneralProcessActive(true);
117 param->SetMinEnergy(100*CLHEP::eV);
118 param->SetLowestElectronEnergy(100*CLHEP::eV);
119 param->SetNumberOfBinsPerDecade(20);
121 param->SetStepFunction(0.2, 10*CLHEP::um);
122 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
123 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
124 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
125 param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
126 param->SetMscStepLimitType(fUseSafetyPlus); // for e-/e+ msc gs
127 param->SetMscSkin(3); // error-free stepping for e-/e+ msc gs
128 param->SetMscRangeFactor(0.08); // error-free stepping for e-/e+ msc gs
129 param->SetMuHadLateralDisplacement(true);
130 param->SetFluo(true);
131 param->SetUseICRU90Data(true);
132 param->Set3GammaAnnihilationOnFly(true);
134 param->SetMaxNIELEnergy(1*CLHEP::MeV);
137}
@ bElectromagnetic
@ fUrbanFluctuation
@ fAllisonPositronium
@ fUseSafetyPlus
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 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 SetMscSkin(G4double val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)
void Set3GammaAnnihilationOnFly(G4bool v)
G4VPhysicsConstructor(const G4String &="")
void SetVerboseLevel(G4int value)

◆ ~G4EmStandardPhysics_option4()

G4EmStandardPhysics_option4::~G4EmStandardPhysics_option4 ( )
overridedefault

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option4::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 145 of file G4EmStandardPhysics_option4.cc.

146{
147 // minimal set of particles for EM physics
149}
static void ConstructMinimalEmSet()

◆ ConstructProcess()

void G4EmStandardPhysics_option4::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 153 of file G4EmStandardPhysics_option4.cc.

154{
155 if(verboseLevel > 1) {
156 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
157 }
159
160 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
161 G4EmParameters* param = G4EmParameters::Instance();
162
163 // processes used by several particles
164 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
165
166 // nuclear stopping is enabled if the energy limit above zero
167 G4double nielEnergyLimit = param->MaxNIELEnergy();
168 G4NuclearStopping* pnuc = nullptr;
169 if(nielEnergyLimit > 0.0) {
170 pnuc = new G4NuclearStopping();
171 pnuc->SetMaxKinEnergy(nielEnergyLimit);
172 }
173
174 // high energy limit for e+- scattering models and bremsstrahlung
175 G4double highEnergyLimit = param->MscEnergyLimit();
176
177 // Add gamma EM Processes
178 G4ParticleDefinition* particle = G4Gamma::Gamma();
179 G4bool polar = param->EnablePolarisation();
180
181 // Photoelectric
182 G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
183 G4VEmModel* peModel = new G4LivermorePhotoElectricModel();
184 pe->SetEmModel(peModel);
185 if(polar) {
186 peModel->SetAngularDistribution(new G4PhotoElectricAngularGeneratorPolarized());
187 }
188
189 // Compton scattering
190 G4ComptonScattering* cs = new G4ComptonScattering;
191 cs->SetEmModel(new G4KleinNishinaModel());
192 G4VEmModel* cModel = nullptr;
193 if(polar) {
194 cModel = new G4LowEPPolarizedComptonModel();
195 } else {
196 cModel = new G4LowEPComptonModel();
197 }
198 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
199 cs->AddEmModel(0, cModel);
200
201 // Gamma conversion
202 G4GammaConversion* gc = new G4GammaConversion();
203 G4VEmModel* conv = new G4BetheHeitler5DModel();
204 gc->SetEmModel(conv);
205
206 // default Rayleigh scattering is Livermore
207 G4RayleighScattering* rl = new G4RayleighScattering();
208 if(polar) {
209 rl->SetEmModel(new G4LivermorePolarizedRayleighModel());
210 }
211
212 if(param->GeneralProcessActive()) {
213 G4GammaGeneralProcess* sp = new G4GammaGeneralProcess();
214 sp->AddEmProcess(pe);
215 sp->AddEmProcess(cs);
216 sp->AddEmProcess(gc);
217 sp->AddEmProcess(rl);
219 ph->RegisterProcess(sp, particle);
220 } else {
221 ph->RegisterProcess(pe, particle);
222 ph->RegisterProcess(cs, particle);
223 ph->RegisterProcess(gc, particle);
224 ph->RegisterProcess(rl, particle);
225 }
226
227 // e-
228 particle = G4Electron::Electron();
229
230 // e-/e+ msc gs with Mott-correction
231 // (Mott-correction is set through G4EmParameters)
232 G4GoudsmitSaundersonMscModel* msc1 = new G4GoudsmitSaundersonMscModel();
233 G4WentzelVIModel* msc2 = new G4WentzelVIModel();
234 msc1->SetHighEnergyLimit(highEnergyLimit);
235 msc2->SetLowEnergyLimit(highEnergyLimit);
236 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
237
238 G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
239 G4CoulombScattering* ss = new G4CoulombScattering();
240 ss->SetEmModel(ssm);
241 ss->SetMinKinEnergy(highEnergyLimit);
242 ssm->SetLowEnergyLimit(highEnergyLimit);
243 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
244
245 // ionisation
246 G4eIonisation* eioni = new G4eIonisation();
248 G4VEmModel* theIoniMod = new G4PenelopeIonisationModel();
249 theIoniMod->SetHighEnergyLimit(0.1*CLHEP::MeV);
250 eioni->AddEmModel(0, theIoniMod);
251
252 // bremsstrahlung
253 G4eBremsstrahlung* brem = new G4eBremsstrahlung();
254 G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel();
255 G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel();
256 br1->SetAngularDistribution(new G4Generator2BS());
257 br2->SetAngularDistribution(new G4Generator2BS());
258 brem->SetEmModel(br1);
259 brem->SetEmModel(br2);
260 br1->SetHighEnergyLimit(CLHEP::GeV);
261
262 G4ePairProduction* ee = new G4ePairProduction();
263
264 // register processes
265 ph->RegisterProcess(eioni, particle);
266 ph->RegisterProcess(brem, particle);
267 ph->RegisterProcess(ee, particle);
268 ph->RegisterProcess(ss, particle);
269
270 // e+
271 particle = G4Positron::Positron();
272
273 // e-/e+ msc gs with Mott-correction
274 // (Mott-correction is set through G4EmParameters)
275 msc1 = new G4GoudsmitSaundersonMscModel();
276 msc2 = new G4WentzelVIModel();
277 msc1->SetHighEnergyLimit(highEnergyLimit);
278 msc2->SetLowEnergyLimit(highEnergyLimit);
279 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
280
281 ssm = new G4eCoulombScatteringModel();
282 ss = new G4CoulombScattering();
283 ss->SetEmModel(ssm);
284 ss->SetMinKinEnergy(highEnergyLimit);
285 ssm->SetLowEnergyLimit(highEnergyLimit);
286 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
287
288 // ionisation
289 eioni = new G4eIonisation();
291 G4VEmModel* pen = new G4PenelopeIonisationModel();
292 pen->SetHighEnergyLimit(0.1*CLHEP::MeV);
293 eioni->AddEmModel(0, pen);
294
295 // bremsstrahlung
296 brem = new G4eBremsstrahlung();
297 br1 = new G4SeltzerBergerModel();
298 br2 = new G4eBremsstrahlungRelModel();
299 br1->SetAngularDistribution(new G4Generator2BS());
300 br2->SetAngularDistribution(new G4Generator2BS());
301 brem->SetEmModel(br1);
302 brem->SetEmModel(br2);
303 br1->SetHighEnergyLimit(CLHEP::GeV);
304
305 // annihilation
306 auto anni = new G4eplusAnnihilation();
307 if (param->Use3GammaAnnihilationOnFly()) {
308 anni->SetEmModel(new G4eplusTo2or3GammaModel());
309 }
310
311 // register processes
312 ph->RegisterProcess(eioni, particle);
313 ph->RegisterProcess(brem, particle);
314 ph->RegisterProcess(ee, particle);
315 ph->RegisterProcess(anni, particle);
316 ph->RegisterProcess(ss, particle);
317
318 // generic ion
319 particle = G4GenericIon::GenericIon();
320 G4ionIonisation* ionIoni = new G4ionIonisation();
321 ionIoni->SetEmModel(new G4LindhardSorensenIonModel());
322 ph->RegisterProcess(hmsc, particle);
323 ph->RegisterProcess(ionIoni, particle);
324 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
325
326 // muons, hadrons, ions
328
329 // extra configuration
330 G4EmModelActivator mact(GetPhysicsName());
331}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#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
G4double MscEnergyLimit() const
G4bool Use3GammaAnnihilationOnFly() const
G4bool GeneralProcessActive() const
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
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 SetHighEnergyLimit(G4double)
void SetActivationLowEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetMinKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const

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