Geant4  10.02.p03
G4EmStandardPhysics_option4 Class Reference

#include <G4EmStandardPhysics_option4.hh>

Inheritance diagram for G4EmStandardPhysics_option4:
Collaboration diagram for G4EmStandardPhysics_option4:

Public Member Functions

 G4EmStandardPhysics_option4 (G4int ver=1)
 
 G4EmStandardPhysics_option4 (G4int ver, const G4String &name)
 
virtual ~G4EmStandardPhysics_option4 ()
 
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
 

Private Attributes

G4int verbose
 

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 53 of file G4EmStandardPhysics_option4.hh.

Constructor & Destructor Documentation

◆ G4EmStandardPhysics_option4() [1/2]

G4EmStandardPhysics_option4::G4EmStandardPhysics_option4 ( G4int  ver = 1)

Definition at line 121 of file G4EmStandardPhysics_option4.cc.

122  : G4VPhysicsConstructor("G4EmStandard_opt4"), verbose(ver)
123 {
125  param->SetDefaults();
126  param->SetVerbose(verbose);
127  param->SetMinEnergy(100*eV);
128  param->SetMaxEnergy(10*TeV);
129  param->SetLowestElectronEnergy(100*eV);
130  param->SetNumberOfBinsPerDecade(20);
132  param->SetMscRangeFactor(0.02);
134  //param->SetLatDisplacementBeyondSafety(true);
135  param->SetFluo(true);
137 }
void SetVerbose(G4int val)
void SetLowestElectronEnergy(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMaxEnergy(G4double val)
void SetMscRangeFactor(G4double val)
void SetNumberOfBinsPerDecade(G4int val)
static const double eV
Definition: G4SIunits.hh:212
void SetMinEnergy(G4double val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
void ActivateAngularGeneratorForIonisation(G4bool val)
static const double TeV
Definition: G4SIunits.hh:215
void SetFluo(G4bool val)
Here is the call graph for this function:

◆ G4EmStandardPhysics_option4() [2/2]

G4EmStandardPhysics_option4::G4EmStandardPhysics_option4 ( G4int  ver,
const G4String name 
)

Definition at line 141 of file G4EmStandardPhysics_option4.cc.

143  : G4VPhysicsConstructor("G4EmStandard_opt4"), verbose(ver)
144 {
146  param->SetDefaults();
147  param->SetVerbose(verbose);
148  param->SetMinEnergy(100*eV);
149  param->SetMaxEnergy(10*TeV);
150  param->SetLowestElectronEnergy(100*eV);
151  param->SetNumberOfBinsPerDecade(20);
153  param->SetMscRangeFactor(0.02);
155  // param->SetLatDisplacementBeyondSafety(true);
156  param->SetFluo(true);
158 }
void SetVerbose(G4int val)
void SetLowestElectronEnergy(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMaxEnergy(G4double val)
void SetMscRangeFactor(G4double val)
void SetNumberOfBinsPerDecade(G4int val)
static const double eV
Definition: G4SIunits.hh:212
void SetMinEnergy(G4double val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
void ActivateAngularGeneratorForIonisation(G4bool val)
static const double TeV
Definition: G4SIunits.hh:215
void SetFluo(G4bool val)
Here is the call graph for this function:

◆ ~G4EmStandardPhysics_option4()

G4EmStandardPhysics_option4::~G4EmStandardPhysics_option4 ( )
virtual

Definition at line 162 of file G4EmStandardPhysics_option4.cc.

163 {}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option4::ConstructParticle ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 167 of file G4EmStandardPhysics_option4.cc.

168 {
169  // gamma
170  G4Gamma::Gamma();
171 
172  // leptons
177 
178  // mesons
183 
184  // barions
187 
188  // ions
191  G4He3::He3();
192  G4Alpha::Alpha();
194 
195  // dna
196  G4EmModelActivator mact;
197  mact.ConstructParticle();
198 }
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4He3 * He3()
Definition: G4He3.cc:94
Here is the call graph for this function:

◆ ConstructProcess()

void G4EmStandardPhysics_option4::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 202 of file G4EmStandardPhysics_option4.cc.

203 {
204  if(verbose > 1) {
205  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
206  }
208 
209  // muon & hadron bremsstrahlung and pair production
218 
219  // muon & hadron multiple scattering
221  mumsc->AddEmModel(0, new G4WentzelVIModel());
223 
225  pimsc->AddEmModel(0, new G4WentzelVIModel());
227 
229  kmsc->AddEmModel(0, new G4WentzelVIModel());
231 
233  pmsc->AddEmModel(0, new G4WentzelVIModel());
235 
236  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
237 
238  // energy limits for e+- scattering models
239  G4double highEnergyLimit = 100*MeV;
240  // energy limits for e+- ionisation models
241  G4double penEnergyLimit = 1*MeV;
242 
243  // nuclear stopping
244  G4NuclearStopping* pnuc = new G4NuclearStopping();
245 
246  // Add standard EM Processes
247  auto myParticleIterator=GetParticleIterator();
248  myParticleIterator->reset();
249  while( (*myParticleIterator)() ){
250  G4ParticleDefinition* particle = myParticleIterator->value();
251  G4String particleName = particle->GetParticleName();
252 
253  if (particleName == "gamma") {
254 
255  // Photoelectric
257  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
258  pe->SetEmModel(theLivermorePEModel,1);
259  ph->RegisterProcess(pe, particle);
260 
261  // Compton scattering
263  cs->SetEmModel(new G4KleinNishinaModel(),1);
264  G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
265  theLowEPComptonModel->SetHighEnergyLimit(20*MeV);
266  cs->AddEmModel(0, theLowEPComptonModel);
267  ph->RegisterProcess(cs, particle);
268 
269  // Gamma conversion
271  G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel();
272  thePenelopeGCModel->SetHighEnergyLimit(1*GeV);
273  gc->SetEmModel(thePenelopeGCModel,1);
274  ph->RegisterProcess(gc, particle);
275 
276  // Rayleigh scattering
277  ph->RegisterProcess(new G4RayleighScattering(), particle);
278 
279  } else if (particleName == "e-") {
280 
281  // multiple scattering
283  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
284  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
285  msc1->SetHighEnergyLimit(highEnergyLimit);
286  msc2->SetLowEnergyLimit(highEnergyLimit);
287  msc->AddEmModel(0, msc1);
288  msc->AddEmModel(0, msc2);
289 
292  ss->SetEmModel(ssm, 1);
293  ss->SetMinKinEnergy(highEnergyLimit);
294  ssm->SetLowEnergyLimit(highEnergyLimit);
295  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
296 
297  // ionisation
298  G4eIonisation* eIoni = new G4eIonisation();
299  eIoni->SetStepFunction(0.2, 100*um);
301  pen->SetHighEnergyLimit(penEnergyLimit);
302  eIoni->AddEmModel(0, pen, new G4UniversalFluctuation());
303 
304  // bremsstrahlung
305  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
310  brem->SetEmModel(br1,1);
311  brem->SetEmModel(br2,2);
312  br2->SetLowEnergyLimit(GeV);
313 
314  // register processes
315  ph->RegisterProcess(msc, particle);
316  ph->RegisterProcess(eIoni, particle);
317  ph->RegisterProcess(brem, particle);
318  ph->RegisterProcess(ss, particle);
319 
320  } else if (particleName == "e+") {
321 
322  // multiple scattering
324  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
325  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
326  msc1->SetHighEnergyLimit(highEnergyLimit);
327  msc2->SetLowEnergyLimit(highEnergyLimit);
328  msc->AddEmModel(0, msc1);
329  msc->AddEmModel(0, msc2);
330 
333  ss->SetEmModel(ssm, 1);
334  ss->SetMinKinEnergy(highEnergyLimit);
335  ssm->SetLowEnergyLimit(highEnergyLimit);
336  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
337 
338  // ionisation
339  G4eIonisation* eIoni = new G4eIonisation();
340  eIoni->SetStepFunction(0.2, 100*um);
342  pen->SetHighEnergyLimit(penEnergyLimit);
343  eIoni->AddEmModel(0, pen, new G4UniversalFluctuation());
344 
345  // bremsstrahlung
346  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
351  brem->SetEmModel(br1,1);
352  brem->SetEmModel(br2,2);
353  br2->SetLowEnergyLimit(GeV);
354 
355  // register processes
356  ph->RegisterProcess(msc, particle);
357  ph->RegisterProcess(eIoni, particle);
358  ph->RegisterProcess(brem, particle);
359  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
360  ph->RegisterProcess(ss, particle);
361 
362  } else if (particleName == "mu+" ||
363  particleName == "mu-" ) {
364 
365  G4MuIonisation* muIoni = new G4MuIonisation();
366  muIoni->SetStepFunction(0.2, 50*um);
367 
368  ph->RegisterProcess(mumsc, particle);
369  ph->RegisterProcess(muIoni, particle);
370  ph->RegisterProcess(mub, particle);
371  ph->RegisterProcess(mup, particle);
372  ph->RegisterProcess(muss, particle);
373 
374  } else if (particleName == "alpha" ||
375  particleName == "He3") {
376 
378  G4ionIonisation* ionIoni = new G4ionIonisation();
379  ionIoni->SetStepFunction(0.1, 10*um);
380 
381  ph->RegisterProcess(msc, particle);
382  ph->RegisterProcess(ionIoni, particle);
383  ph->RegisterProcess(pnuc, particle);
384 
385  } else if (particleName == "GenericIon") {
386 
387  G4ionIonisation* ionIoni = new G4ionIonisation();
388  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
389  ionIoni->SetStepFunction(0.1, 1*um);
390 
391  ph->RegisterProcess(hmsc, particle);
392  ph->RegisterProcess(ionIoni, particle);
393  ph->RegisterProcess(pnuc, particle);
394 
395  } else if (particleName == "pi+" ||
396  particleName == "pi-" ) {
397 
398  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
399  G4hIonisation* hIoni = new G4hIonisation();
400  hIoni->SetStepFunction(0.2, 50*um);
401 
402  ph->RegisterProcess(pimsc, particle);
403  ph->RegisterProcess(hIoni, particle);
404  ph->RegisterProcess(pib, particle);
405  ph->RegisterProcess(pip, particle);
406  ph->RegisterProcess(piss, particle);
407 
408  } else if (particleName == "kaon+" ||
409  particleName == "kaon-" ) {
410 
411  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
412  G4hIonisation* hIoni = new G4hIonisation();
413  hIoni->SetStepFunction(0.2, 50*um);
414 
415  ph->RegisterProcess(kmsc, particle);
416  ph->RegisterProcess(hIoni, particle);
417  ph->RegisterProcess(kb, particle);
418  ph->RegisterProcess(kp, particle);
419  ph->RegisterProcess(kss, particle);
420 
421  } else if (particleName == "proton" ||
422  particleName == "anti_proton") {
423 
424  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
425  G4hIonisation* hIoni = new G4hIonisation();
426  hIoni->SetStepFunction(0.1, 20*um);
427 
428  ph->RegisterProcess(pmsc, particle);
429  ph->RegisterProcess(hIoni, particle);
430  ph->RegisterProcess(pb, particle);
431  ph->RegisterProcess(pp, particle);
432  ph->RegisterProcess(pss, particle);
433  ph->RegisterProcess(pnuc, particle);
434 
435  } else if (particleName == "B+" ||
436  particleName == "B-" ||
437  particleName == "D+" ||
438  particleName == "D-" ||
439  particleName == "Ds+" ||
440  particleName == "Ds-" ||
441  particleName == "anti_He3" ||
442  particleName == "anti_alpha" ||
443  particleName == "anti_deuteron" ||
444  particleName == "anti_lambda_c+" ||
445  particleName == "anti_omega-" ||
446  particleName == "anti_sigma_c+" ||
447  particleName == "anti_sigma_c++" ||
448  particleName == "anti_sigma+" ||
449  particleName == "anti_sigma-" ||
450  particleName == "anti_triton" ||
451  particleName == "anti_xi_c+" ||
452  particleName == "anti_xi-" ||
453  particleName == "deuteron" ||
454  particleName == "lambda_c+" ||
455  particleName == "omega-" ||
456  particleName == "sigma_c+" ||
457  particleName == "sigma_c++" ||
458  particleName == "sigma+" ||
459  particleName == "sigma-" ||
460  particleName == "tau+" ||
461  particleName == "tau-" ||
462  particleName == "triton" ||
463  particleName == "xi_c+" ||
464  particleName == "xi-" ) {
465 
466  ph->RegisterProcess(hmsc, particle);
467  ph->RegisterProcess(new G4hIonisation(), particle);
468  ph->RegisterProcess(pnuc, particle);
469  }
470  }
471 
472  // Nuclear stopping
473  pnuc->SetMaxKinEnergy(MeV);
474 
475  // Deexcitation
478 
479  G4EmModelActivator mact;
480  mact.ConstructProcess();
481 }
static const double MeV
Definition: G4SIunits.hh:211
static G4LossTableManager * Instance()
void SetStepFunction(G4double v1, G4double v2)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
void SetEmModel(G4VEmModel *, G4int index=1)
const G4String & GetParticleName() const
const G4String & GetPhysicsName() const
G4GLOB_DLL std::ostream G4cout
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static const double GeV
Definition: G4SIunits.hh:214
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
static const double um
Definition: G4SIunits.hh:112
#define G4endl
Definition: G4ios.hh:61
void SetMinKinEnergy(G4double e)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
void SetAtomDeexcitation(G4VAtomDeexcitation *)
Here is the call graph for this function:

Member Data Documentation

◆ verbose

G4int G4EmStandardPhysics_option4::verbose
private

Definition at line 67 of file G4EmStandardPhysics_option4.hh.


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