Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmLowEPPhysics Class Reference

#include <G4EmLowEPPhysics.hh>

Inheritance diagram for G4EmLowEPPhysics:
Collaboration diagram for G4EmLowEPPhysics:

Public Member Functions

 G4EmLowEPPhysics (G4int ver=1, 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
 

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

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(1*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)
static constexpr double TeV
Definition: G4SIunits.hh:218
void SetNumberOfBinsPerDecade(G4int val)
static constexpr double eV
Definition: G4SIunits.hh:215
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:

G4EmLowEPPhysics::~G4EmLowEPPhysics ( )
virtual

Definition at line 143 of file G4EmLowEPPhysics.cc.

144 {}

Member Function Documentation

void G4EmLowEPPhysics::ConstructParticle ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 148 of file G4EmLowEPPhysics.cc.

149 {
150  // gamma
151  G4Gamma::Gamma();
152 
153  // leptons
158 
159  // mesons
164 
165  // baryons
168 
169  // ions
172  G4He3::He3();
173  G4Alpha::Alpha();
175 }
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 G4EmLowEPPhysics::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 179 of file G4EmLowEPPhysics.cc.

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

Here is the call graph for this function:


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