Geant4  10.02.p03
G4EmLowEPPhysics Class Reference

#include <G4EmLowEPPhysics.hh>

Inheritance diagram for G4EmLowEPPhysics:
Collaboration diagram for G4EmLowEPPhysics:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmLowEPPhysics() [1/2]

G4EmLowEPPhysics::G4EmLowEPPhysics ( G4int  ver = 1)

Definition at line 126 of file G4EmLowEPPhysics.cc.

127  : G4VPhysicsConstructor("G4EmLowEPPhysics"), verbose(ver)
128 {
130  param->SetDefaults();
131  param->SetVerbose(verbose);
132  param->SetMinEnergy(100*eV);
133  param->SetMaxEnergy(10*TeV);
134  param->SetLowestElectronEnergy(100*eV);
135  param->SetNumberOfBinsPerDecade(20);
137  param->SetFluo(true);
139 }
void SetVerbose(G4int val)
void SetLowestElectronEnergy(G4double val)
void SetMaxEnergy(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:

◆ G4EmLowEPPhysics() [2/2]

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

Definition at line 143 of file G4EmLowEPPhysics.cc.

144  : G4VPhysicsConstructor("G4EmLowEPPhysics"), verbose(ver)
145 {
147  param->SetDefaults();
148  param->SetVerbose(verbose);
149  param->SetMinEnergy(100*eV);
150  param->SetMaxEnergy(10*TeV);
151  param->SetLowestElectronEnergy(100*eV);
152  param->SetNumberOfBinsPerDecade(20);
154  param->SetFluo(true);
156 }
void SetVerbose(G4int val)
void SetLowestElectronEnergy(G4double val)
void SetMaxEnergy(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:

◆ ~G4EmLowEPPhysics()

G4EmLowEPPhysics::~G4EmLowEPPhysics ( )
virtual

Definition at line 160 of file G4EmLowEPPhysics.cc.

161 {}

Member Function Documentation

◆ ConstructParticle()

void G4EmLowEPPhysics::ConstructParticle ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 165 of file G4EmLowEPPhysics.cc.

166 {
167  // gamma
168  G4Gamma::Gamma();
169 
170  // leptons
175 
176  // mesons
181 
182  // baryons
185 
186  // ions
189  G4He3::He3();
190  G4Alpha::Alpha();
192 
193  // dna
194  G4EmModelActivator mact;
195  mact.ConstructParticle();
196 }
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 G4EmLowEPPhysics::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 200 of file G4EmLowEPPhysics.cc.

201 {
203 
204  // muon & hadron bremsstrahlung and pair production
213 
214  // muon & hadron multiple scattering
216  mumsc->SetEmModel(new G4LowEWentzelVIModel());
218  hmsc->SetEmModel(new G4LowEWentzelVIModel());
220  pmsc->SetEmModel(new G4LowEWentzelVIModel());
222  pimsc->SetEmModel(new G4LowEWentzelVIModel());
224  kmsc->SetEmModel(new G4LowEWentzelVIModel());
225 
226  // nuclear stopping
227  G4NuclearStopping* ionnuc = new G4NuclearStopping();
228  G4NuclearStopping* pnuc = new G4NuclearStopping();
229 
230  // Add Livermore EM Processes
231  auto myParticleIterator=GetParticleIterator();
232  myParticleIterator->reset();
233 
234  while( (*myParticleIterator)() ){
235 
236  G4ParticleDefinition* particle = myParticleIterator->value();
237  G4String particleName = particle->GetParticleName();
238 
239  if(verbose > 1)
240  G4cout << "### " << GetPhysicsName() << " instantiates for "
241  << particleName << G4endl;
242 
243  if (particleName == "gamma") {
244 
245  // Photoelectric effect - Livermore model only
246  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
247  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel(), 1);
248  ph->RegisterProcess(thePhotoElectricEffect, particle);
249 
250  // Compton scattering - Livermore model above 20 MeV, Monarsh's model below
251  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
252  theComptonScattering->SetEmModel(new G4LivermoreComptonModel(),1);
253  G4LowEPComptonModel* theLowEPComptonModel =
254  new G4LowEPComptonModel();
255  theLowEPComptonModel->SetHighEnergyLimit(20*MeV);
256  theComptonScattering->AddEmModel(0, theLowEPComptonModel);
257  ph->RegisterProcess(theComptonScattering, particle);
258 
259  // gamma conversion - Livermore model below 80 GeV
260  G4GammaConversion* theGammaConversion = new G4GammaConversion();
261  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel(),1);
262  ph->RegisterProcess(theGammaConversion, particle);
263 
264  // default Rayleigh scattering is Livermore
265  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
266  ph->RegisterProcess(theRayleigh, particle);
267 
268  } else if (particleName == "e-") {
269 
270  // multiple scattering
272  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
273 
274  // Ionisation - Livermore should be used only for low energies
275  G4eIonisation* eIoni = new G4eIonisation();
276  G4LivermoreIonisationModel* theIoniLivermore = new
278  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
279  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
280  eIoni->SetStepFunction(0.2, 100*um); //
281 
282  // Bremsstrahlung
283  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
284  G4VEmModel* theBrem = new G4SeltzerBergerModel();
285  theBrem->SetHighEnergyLimit(1*GeV);
286  theBrem->SetAngularDistribution(new G4Generator2BS());
287  eBrem->SetEmModel(theBrem, 1);
288 
289  // register processes
290  ph->RegisterProcess(msc, particle);
291  ph->RegisterProcess(eIoni, particle);
292  ph->RegisterProcess(eBrem, particle);
293 
294  } else if (particleName == "e+") {
295 
296  // multiple scattering
298  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
299 
300  // Standard ionisation
301  G4eIonisation* eIoni = new G4eIonisation();
302  eIoni->SetStepFunction(0.2, 100*um);
303 
304  // Bremsstrahlung
305  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
306  G4VEmModel* theBrem = new G4SeltzerBergerModel();
307  theBrem->SetHighEnergyLimit(1*GeV);
308  theBrem->SetAngularDistribution(new G4Generator2BS());
309  eBrem->SetEmModel(theBrem, 1);
310 
311  // register processes
312  ph->RegisterProcess(msc, particle);
313  ph->RegisterProcess(eIoni, particle);
314  ph->RegisterProcess(eBrem, particle);
315  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
316 
317  } else if (particleName == "mu+" ||
318  particleName == "mu-" ) {
319 
320  G4MuIonisation* muIoni = new G4MuIonisation();
321  muIoni->SetStepFunction(0.2, 50*um);
322 
323  ph->RegisterProcess(mumsc, particle);
324  ph->RegisterProcess(muIoni, particle);
325  ph->RegisterProcess(mub, particle);
326  ph->RegisterProcess(mup, particle);
327  //ph->RegisterProcess(new G4CoulombScattering(), particle);
328 
329  } else if (particleName == "alpha" ||
330  particleName == "He3" ) {
331 
332  G4ionIonisation* ionIoni = new G4ionIonisation();
333  ionIoni->SetStepFunction(0.1, 10*um);
334 
335  ph->RegisterProcess(hmsc, particle);
336  ph->RegisterProcess(ionIoni, particle);
337  ph->RegisterProcess(ionnuc, particle);
338 
339  } else if (particleName == "GenericIon") {
340 
341  G4ionIonisation* ionIoni = new G4ionIonisation();
342  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
343  ionIoni->SetStepFunction(0.1, 1*um);
344 
345  ph->RegisterProcess(hmsc, particle);
346  ph->RegisterProcess(ionIoni, particle);
347  ph->RegisterProcess(ionnuc, particle);
348 
349  } else if (particleName == "pi+" ||
350  particleName == "pi-" ) {
351 
352  G4hIonisation* hIoni = new G4hIonisation();
353  hIoni->SetStepFunction(0.2, 50*um);
354 
355  ph->RegisterProcess(pimsc, particle);
356  ph->RegisterProcess(hIoni, particle);
357  ph->RegisterProcess(pib, particle);
358  ph->RegisterProcess(pip, particle);
359 
360  } else if (particleName == "kaon+" ||
361  particleName == "kaon-" ) {
362 
363  G4hIonisation* hIoni = new G4hIonisation();
364  hIoni->SetStepFunction(0.2, 50*um);
365 
366  ph->RegisterProcess(kmsc, particle);
367  ph->RegisterProcess(hIoni, particle);
368  ph->RegisterProcess(kb, particle);
369  ph->RegisterProcess(kp, particle);
370 
371  } else if (particleName == "proton" ||
372  particleName == "anti_proton") {
373 
374  G4hIonisation* hIoni = new G4hIonisation();
375  hIoni->SetStepFunction(0.2, 50*um);
376 
377  ph->RegisterProcess(pmsc, particle);
378  ph->RegisterProcess(hIoni, particle);
379  ph->RegisterProcess(pb, particle);
380  ph->RegisterProcess(pp, particle);
381  ph->RegisterProcess(pnuc, particle);
382 
383  } else if (particleName == "B+" ||
384  particleName == "B-" ||
385  particleName == "D+" ||
386  particleName == "D-" ||
387  particleName == "Ds+" ||
388  particleName == "Ds-" ||
389  particleName == "anti_He3" ||
390  particleName == "anti_alpha" ||
391  particleName == "anti_deuteron" ||
392  particleName == "anti_lambda_c+" ||
393  particleName == "anti_omega-" ||
394  particleName == "anti_sigma_c+" ||
395  particleName == "anti_sigma_c++" ||
396  particleName == "anti_sigma+" ||
397  particleName == "anti_sigma-" ||
398  particleName == "anti_triton" ||
399  particleName == "anti_xi_c+" ||
400  particleName == "anti_xi-" ||
401  particleName == "deuteron" ||
402  particleName == "lambda_c+" ||
403  particleName == "omega-" ||
404  particleName == "sigma_c+" ||
405  particleName == "sigma_c++" ||
406  particleName == "sigma+" ||
407  particleName == "sigma-" ||
408  particleName == "tau+" ||
409  particleName == "tau-" ||
410  particleName == "triton" ||
411  particleName == "xi_c+" ||
412  particleName == "xi-" ) {
413 
414  ph->RegisterProcess(hmsc, particle);
415  ph->RegisterProcess(new G4hIonisation(), particle);
416  }
417  }
418 
419  // Nuclear stopping
420  pnuc->SetMaxKinEnergy(MeV);
421 
422  // Deexcitation
423  //
426 
427  G4EmModelActivator mact;
428  mact.ConstructProcess();
429 }
static const double MeV
Definition: G4SIunits.hh:211
static G4LossTableManager * Instance()
void SetStepFunction(G4double v1, G4double v2)
void SetEmModel(G4VMscModel *, G4int index=1)
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 SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, 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
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetAtomDeexcitation(G4VAtomDeexcitation *)
Here is the call graph for this function:

Member Data Documentation

◆ verbose

G4int G4EmLowEPPhysics::verbose
private

Definition at line 50 of file G4EmLowEPPhysics.hh.


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