Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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, 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
 

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::G4EmStandardPhysics_option4 ( G4int  ver = 1,
const G4String name = "" 
)
explicit

Definition at line 122 of file G4EmStandardPhysics_option4.cc.

124  : G4VPhysicsConstructor("G4EmStandard_opt4"), verbose(ver)
125 {
127  param->SetDefaults();
128  param->SetVerbose(verbose);
129  param->SetMinEnergy(100*eV);
130  param->SetMaxEnergy(10*TeV);
131  param->SetLowestElectronEnergy(100*eV);
132  param->SetNumberOfBinsPerDecade(20);
134  param->SetMscRangeFactor(0.02);
136  param->SetMuHadLateralDisplacement(true);
137  // param->SetLatDisplacementBeyondSafety(true);
138  param->SetFluo(true);
140 }
void SetVerbose(G4int val)
void SetLowestElectronEnergy(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMaxEnergy(G4double val)
static constexpr double TeV
Definition: G4SIunits.hh:218
void SetMscRangeFactor(G4double val)
void SetNumberOfBinsPerDecade(G4int val)
static constexpr double eV
Definition: G4SIunits.hh:215
void SetMuHadLateralDisplacement(G4bool val)
void SetMinEnergy(G4double val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetFluo(G4bool val)

Here is the call graph for this function:

G4EmStandardPhysics_option4::~G4EmStandardPhysics_option4 ( )
virtual

Definition at line 144 of file G4EmStandardPhysics_option4.cc.

145 {}

Member Function Documentation

void G4EmStandardPhysics_option4::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 149 of file G4EmStandardPhysics_option4.cc.

150 {
151  // gamma
152  G4Gamma::Gamma();
153 
154  // leptons
159 
160  // mesons
165 
166  // barions
169 
170  // ions
173  G4He3::He3();
174  G4Alpha::Alpha();
176 }
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:

void G4EmStandardPhysics_option4::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 180 of file G4EmStandardPhysics_option4.cc.

181 {
182  if(verbose > 1) {
183  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
184  }
186 
187  // muon & hadron bremsstrahlung and pair production
197 
198  // muon & hadron multiple scattering
200  mumsc->AddEmModel(0, new G4WentzelVIModel());
202 
204  pimsc->AddEmModel(0, new G4WentzelVIModel());
206 
208  kmsc->AddEmModel(0, new G4WentzelVIModel());
210 
212  pmsc->AddEmModel(0, new G4WentzelVIModel());
214 
215  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
216 
217  // energy limits for e+- scattering models
218  G4double highEnergyLimit = 100*MeV;
219  // energy limits for e+- ionisation models
220  G4double penEnergyLimit = 1*MeV;
221 
222  // nuclear stopping
223  G4NuclearStopping* pnuc = new G4NuclearStopping();
224 
225  // Add standard EM Processes
226  auto myParticleIterator=GetParticleIterator();
227  myParticleIterator->reset();
228  while( (*myParticleIterator)() ){
229  G4ParticleDefinition* particle = myParticleIterator->value();
230  G4String particleName = particle->GetParticleName();
231 
232  if (particleName == "gamma") {
233 
234  // Photoelectric
236  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
237  pe->SetEmModel(theLivermorePEModel,1);
238  ph->RegisterProcess(pe, particle);
239 
240  // Compton scattering
242  cs->SetEmModel(new G4KleinNishinaModel(),1);
243  G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
244  theLowEPComptonModel->SetHighEnergyLimit(20*MeV);
245  cs->AddEmModel(0, theLowEPComptonModel);
246  ph->RegisterProcess(cs, particle);
247 
248  // Gamma conversion
250  G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel();
251  thePenelopeGCModel->SetHighEnergyLimit(1*GeV);
252  gc->SetEmModel(thePenelopeGCModel,1);
253  ph->RegisterProcess(gc, particle);
254 
255  // Rayleigh scattering
256  ph->RegisterProcess(new G4RayleighScattering(), particle);
257 
258  } else if (particleName == "e-") {
259 
260  // multiple scattering
262  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
263  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
264  msc1->SetHighEnergyLimit(highEnergyLimit);
265  msc2->SetLowEnergyLimit(highEnergyLimit);
266  msc->AddEmModel(0, msc1);
267  msc->AddEmModel(0, msc2);
268 
271  ss->SetEmModel(ssm, 1);
272  ss->SetMinKinEnergy(highEnergyLimit);
273  ssm->SetLowEnergyLimit(highEnergyLimit);
274  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
275 
276  // ionisation
277  G4eIonisation* eIoni = new G4eIonisation();
278  eIoni->SetStepFunction(0.2, 100*um);
280  pen->SetHighEnergyLimit(penEnergyLimit);
281  eIoni->AddEmModel(0, pen, new G4UniversalFluctuation());
282 
283  // bremsstrahlung
284  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
289  brem->SetEmModel(br1,1);
290  brem->SetEmModel(br2,2);
291  br2->SetLowEnergyLimit(GeV);
292 
293  // register processes
294  ph->RegisterProcess(msc, particle);
295  ph->RegisterProcess(eIoni, particle);
296  ph->RegisterProcess(brem, particle);
297  ph->RegisterProcess(ee, particle);
298  ph->RegisterProcess(ss, particle);
299 
300  } else if (particleName == "e+") {
301 
302  // multiple scattering
304  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
305  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
306  msc1->SetHighEnergyLimit(highEnergyLimit);
307  msc2->SetLowEnergyLimit(highEnergyLimit);
308  msc->AddEmModel(0, msc1);
309  msc->AddEmModel(0, msc2);
310 
313  ss->SetEmModel(ssm, 1);
314  ss->SetMinKinEnergy(highEnergyLimit);
315  ssm->SetLowEnergyLimit(highEnergyLimit);
316  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
317 
318  // ionisation
319  G4eIonisation* eIoni = new G4eIonisation();
320  eIoni->SetStepFunction(0.2, 100*um);
322  pen->SetHighEnergyLimit(penEnergyLimit);
323  eIoni->AddEmModel(0, pen, new G4UniversalFluctuation());
324 
325  // bremsstrahlung
326  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
331  brem->SetEmModel(br1,1);
332  brem->SetEmModel(br2,2);
333  br2->SetLowEnergyLimit(GeV);
334 
335  // register processes
336  ph->RegisterProcess(msc, particle);
337  ph->RegisterProcess(eIoni, particle);
338  ph->RegisterProcess(brem, particle);
339  ph->RegisterProcess(ee, particle);
340  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
341  ph->RegisterProcess(ss, particle);
342 
343  } else if (particleName == "mu+" ||
344  particleName == "mu-" ) {
345 
346  G4MuIonisation* muIoni = new G4MuIonisation();
347  muIoni->SetStepFunction(0.2, 50*um);
348 
349  ph->RegisterProcess(mumsc, particle);
350  ph->RegisterProcess(muIoni, particle);
351  ph->RegisterProcess(mub, particle);
352  ph->RegisterProcess(mup, particle);
353  ph->RegisterProcess(muss, particle);
354 
355  } else if (particleName == "alpha" ||
356  particleName == "He3") {
357 
359  G4ionIonisation* ionIoni = new G4ionIonisation();
360  ionIoni->SetStepFunction(0.1, 10*um);
361 
362  ph->RegisterProcess(msc, particle);
363  ph->RegisterProcess(ionIoni, particle);
364  ph->RegisterProcess(pnuc, particle);
365 
366  } else if (particleName == "GenericIon") {
367 
368  G4ionIonisation* ionIoni = new G4ionIonisation();
369  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
370  ionIoni->SetStepFunction(0.1, 1*um);
371 
372  ph->RegisterProcess(hmsc, particle);
373  ph->RegisterProcess(ionIoni, particle);
374  ph->RegisterProcess(pnuc, particle);
375 
376  } else if (particleName == "pi+" ||
377  particleName == "pi-" ) {
378 
379  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
380  G4hIonisation* hIoni = new G4hIonisation();
381  hIoni->SetStepFunction(0.2, 50*um);
382 
383  ph->RegisterProcess(pimsc, particle);
384  ph->RegisterProcess(hIoni, particle);
385  ph->RegisterProcess(pib, particle);
386  ph->RegisterProcess(pip, particle);
387  ph->RegisterProcess(piss, particle);
388 
389  } else if (particleName == "kaon+" ||
390  particleName == "kaon-" ) {
391 
392  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
393  G4hIonisation* hIoni = new G4hIonisation();
394  hIoni->SetStepFunction(0.2, 50*um);
395 
396  ph->RegisterProcess(kmsc, particle);
397  ph->RegisterProcess(hIoni, particle);
398  ph->RegisterProcess(kb, particle);
399  ph->RegisterProcess(kp, particle);
400  ph->RegisterProcess(kss, particle);
401 
402  } else if (particleName == "proton" ||
403  particleName == "anti_proton") {
404 
405  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
406  G4hIonisation* hIoni = new G4hIonisation();
407  hIoni->SetStepFunction(0.1, 20*um);
408 
409  ph->RegisterProcess(pmsc, particle);
410  ph->RegisterProcess(hIoni, particle);
411  ph->RegisterProcess(pb, particle);
412  ph->RegisterProcess(pp, particle);
413  ph->RegisterProcess(pss, particle);
414  ph->RegisterProcess(pnuc, particle);
415 
416  } else if (particleName == "B+" ||
417  particleName == "B-" ||
418  particleName == "D+" ||
419  particleName == "D-" ||
420  particleName == "Ds+" ||
421  particleName == "Ds-" ||
422  particleName == "anti_He3" ||
423  particleName == "anti_alpha" ||
424  particleName == "anti_deuteron" ||
425  particleName == "anti_lambda_c+" ||
426  particleName == "anti_omega-" ||
427  particleName == "anti_sigma_c+" ||
428  particleName == "anti_sigma_c++" ||
429  particleName == "anti_sigma+" ||
430  particleName == "anti_sigma-" ||
431  particleName == "anti_triton" ||
432  particleName == "anti_xi_c+" ||
433  particleName == "anti_xi-" ||
434  particleName == "deuteron" ||
435  particleName == "lambda_c+" ||
436  particleName == "omega-" ||
437  particleName == "sigma_c+" ||
438  particleName == "sigma_c++" ||
439  particleName == "sigma+" ||
440  particleName == "sigma-" ||
441  particleName == "tau+" ||
442  particleName == "tau-" ||
443  particleName == "triton" ||
444  particleName == "xi_c+" ||
445  particleName == "xi-" ) {
446 
447  ph->RegisterProcess(hmsc, particle);
448  ph->RegisterProcess(new G4hIonisation(), particle);
449  ph->RegisterProcess(pnuc, particle);
450  }
451  }
452 
453  // Nuclear stopping
454  pnuc->SetMaxKinEnergy(MeV);
455 
456  // Deexcitation
459 
461 }
static G4LossTableManager * Instance()
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
static constexpr double um
Definition: G4SIunits.hh:113
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
const G4String & GetPhysicsName() const
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:745
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetMaxKinEnergy(G4double e)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:623
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
static constexpr double GeV
Definition: G4SIunits.hh:217
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=nullptr)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
void SetAtomDeexcitation(G4VAtomDeexcitation *)

Here is the call graph for this function:


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