Geant4  10.02.p03
G4EmLivermorePhysics Class Reference

#include <G4EmLivermorePhysics.hh>

Inheritance diagram for G4EmLivermorePhysics:
Collaboration diagram for G4EmLivermorePhysics:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmLivermorePhysics() [1/2]

G4EmLivermorePhysics::G4EmLivermorePhysics ( G4int  ver = 1)

Definition at line 129 of file G4EmLivermorePhysics.cc.

130  : G4VPhysicsConstructor("G4EmLivermorePhysics"), verbose(ver)
131 {
133  param->SetDefaults();
134  param->SetVerbose(verbose);
135  param->SetMinEnergy(100*eV);
136  param->SetMaxEnergy(10*TeV);
137  param->SetLowestElectronEnergy(100*eV);
138  param->SetNumberOfBinsPerDecade(20);
140  param->SetMscRangeFactor(0.02);
142  param->SetFluo(true);
144 }
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:

◆ G4EmLivermorePhysics() [2/2]

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

Definition at line 148 of file G4EmLivermorePhysics.cc.

149  : G4VPhysicsConstructor("G4EmLivermorePhysics"), verbose(ver)
150 {
152  param->SetDefaults();
153  param->SetVerbose(verbose);
154  param->SetMinEnergy(100*eV);
155  param->SetMaxEnergy(10*TeV);
156  param->SetLowestElectronEnergy(100*eV);
157  param->SetNumberOfBinsPerDecade(20);
159  param->SetMscRangeFactor(0.02);
161  param->SetFluo(true);
163 }
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:

◆ ~G4EmLivermorePhysics()

G4EmLivermorePhysics::~G4EmLivermorePhysics ( )
virtual

Definition at line 167 of file G4EmLivermorePhysics.cc.

168 {}

Member Function Documentation

◆ ConstructParticle()

void G4EmLivermorePhysics::ConstructParticle ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 172 of file G4EmLivermorePhysics.cc.

173 {
174  // gamma
175  G4Gamma::Gamma();
176 
177  // leptons
182 
183  // mesons
188 
189  // baryons
192 
193  // ions
196  G4He3::He3();
197  G4Alpha::Alpha();
199 
200  // dna
201  G4EmModelActivator mact;
202  mact.ConstructParticle();
203 }
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 G4EmLivermorePhysics::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 207 of file G4EmLivermorePhysics.cc.

208 {
209  if(verbose > 1) {
210  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
211  }
213 
214  // muon & hadron bremsstrahlung and pair production
223 
224  // muon & hadron multiple scattering
226  mumsc->AddEmModel(0, new G4WentzelVIModel());
227  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
228  //pimsc->AddEmModel(0, new G4WentzelVIModel());
229  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
230  //kmsc->AddEmModel(0, new G4WentzelVIModel());
231  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
232  //pmsc->AddEmModel(0, new G4WentzelVIModel());
233  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
234 
235  // high energy limit for e+- scattering models
236  G4double highEnergyLimit = 100*MeV;
237 
238  // nuclear stopping
239  G4NuclearStopping* pnuc = new G4NuclearStopping();
240 
241  // Add Livermore EM Processes
242  auto myParticleIterator=GetParticleIterator();
243  myParticleIterator->reset();
244 
245  while( (*myParticleIterator)() ){
246 
247  G4ParticleDefinition* particle = myParticleIterator->value();
248  G4String particleName = particle->GetParticleName();
249 
250  if (particleName == "gamma") {
251 
252  // photoelectric effect - Livermore model only
253  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
254  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel(), 1);
255  ph->RegisterProcess(thePhotoElectricEffect, particle);
256 
257  // Compton scattering - Livermore model only
258  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
259  theComptonScattering->SetEmModel(new G4LivermoreComptonModel(),1);
260  ph->RegisterProcess(theComptonScattering, particle);
261 
262  // gamma conversion - Livermore model below 80 GeV
263  G4GammaConversion* theGammaConversion = new G4GammaConversion();
264  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel(),1);
265  ph->RegisterProcess(theGammaConversion, particle);
266 
267  // default Rayleigh scattering is Livermore
268  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
269  ph->RegisterProcess(theRayleigh, particle);
270 
271  } else if (particleName == "e-") {
272 
273  // multiple scattering
275  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
276  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
277  msc1->SetHighEnergyLimit(highEnergyLimit);
278  msc2->SetLowEnergyLimit(highEnergyLimit);
279  msc->AddEmModel(0, msc1);
280  msc->AddEmModel(0, msc2);
281 
284  ss->SetEmModel(ssm, 1);
285  ss->SetMinKinEnergy(highEnergyLimit);
286  ssm->SetLowEnergyLimit(highEnergyLimit);
287  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
288 
289  // Ionisation - Livermore should be used only for low energies
290  G4eIonisation* eIoni = new G4eIonisation();
291  G4LivermoreIonisationModel* theIoniLivermore = new
293  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
294  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
295  eIoni->SetStepFunction(0.2, 100*um); //
296 
297  // Bremsstrahlung
298  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
299  G4VEmModel* theBremLivermore = new G4LivermoreBremsstrahlungModel();
300  theBremLivermore->SetHighEnergyLimit(1*GeV);
301  theBremLivermore->SetAngularDistribution(new G4Generator2BS());
302  eBrem->SetEmModel(theBremLivermore,1);
303 
304  // register processes
305  ph->RegisterProcess(msc, particle);
306  ph->RegisterProcess(eIoni, particle);
307  ph->RegisterProcess(eBrem, particle);
308  ph->RegisterProcess(ss, particle);
309 
310  } else if (particleName == "e+") {
311 
312  // multiple scattering
314  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
315  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
316  msc1->SetHighEnergyLimit(highEnergyLimit);
317  msc2->SetLowEnergyLimit(highEnergyLimit);
318  msc->AddEmModel(0, msc1);
319  msc->AddEmModel(0, msc2);
320 
323  ss->SetEmModel(ssm, 1);
324  ss->SetMinKinEnergy(highEnergyLimit);
325  ssm->SetLowEnergyLimit(highEnergyLimit);
326  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
327 
328  // ionisation
329  G4eIonisation* eIoni = new G4eIonisation();
330  eIoni->SetStepFunction(0.2, 100*um);
331 
332  // register processes
333  ph->RegisterProcess(msc, particle);
334  ph->RegisterProcess(eIoni, particle);
335  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
336  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
337  ph->RegisterProcess(ss, particle);
338 
339  } else if (particleName == "mu+" ||
340  particleName == "mu-" ) {
341 
342  G4MuIonisation* muIoni = new G4MuIonisation();
343  muIoni->SetStepFunction(0.2, 50*um);
344 
345  ph->RegisterProcess(mumsc, particle);
346  ph->RegisterProcess(muIoni, particle);
347  ph->RegisterProcess(mub, particle);
348  ph->RegisterProcess(mup, particle);
349  ph->RegisterProcess(new G4CoulombScattering(), particle);
350 
351  } else if (particleName == "alpha" ||
352  particleName == "He3" ) {
353 
355  G4ionIonisation* ionIoni = new G4ionIonisation();
356  ionIoni->SetStepFunction(0.1, 10*um);
357 
358  ph->RegisterProcess(msc, particle);
359  ph->RegisterProcess(ionIoni, particle);
360  ph->RegisterProcess(pnuc, particle);
361 
362  } else if (particleName == "GenericIon") {
363 
364  G4ionIonisation* ionIoni = new G4ionIonisation();
365  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
366  ionIoni->SetStepFunction(0.1, 1*um);
367 
368  ph->RegisterProcess(hmsc, particle);
369  ph->RegisterProcess(ionIoni, particle);
370  ph->RegisterProcess(pnuc, particle);
371 
372  } else if (particleName == "pi+" ||
373  particleName == "pi-" ) {
374 
376  G4hIonisation* hIoni = new G4hIonisation();
377  hIoni->SetStepFunction(0.2, 50*um);
378 
379  ph->RegisterProcess(pimsc, particle);
380  ph->RegisterProcess(hIoni, particle);
381  ph->RegisterProcess(pib, particle);
382  ph->RegisterProcess(pip, particle);
383 
384  } else if (particleName == "kaon+" ||
385  particleName == "kaon-" ) {
386 
388  G4hIonisation* hIoni = new G4hIonisation();
389  hIoni->SetStepFunction(0.2, 50*um);
390 
391  ph->RegisterProcess(kmsc, particle);
392  ph->RegisterProcess(hIoni, particle);
393  ph->RegisterProcess(kb, particle);
394  ph->RegisterProcess(kp, particle);
395 
396  } else if (particleName == "proton" ||
397  particleName == "anti_proton") {
398 
400  G4hIonisation* hIoni = new G4hIonisation();
401  hIoni->SetStepFunction(0.2, 50*um);
402 
403  ph->RegisterProcess(pmsc, particle);
404  ph->RegisterProcess(hIoni, particle);
405  ph->RegisterProcess(pb, particle);
406  ph->RegisterProcess(pp, particle);
407  ph->RegisterProcess(pnuc, particle);
408 
409  } else if (particleName == "B+" ||
410  particleName == "B-" ||
411  particleName == "D+" ||
412  particleName == "D-" ||
413  particleName == "Ds+" ||
414  particleName == "Ds-" ||
415  particleName == "anti_He3" ||
416  particleName == "anti_alpha" ||
417  particleName == "anti_deuteron" ||
418  particleName == "anti_lambda_c+" ||
419  particleName == "anti_omega-" ||
420  particleName == "anti_sigma_c+" ||
421  particleName == "anti_sigma_c++" ||
422  particleName == "anti_sigma+" ||
423  particleName == "anti_sigma-" ||
424  particleName == "anti_triton" ||
425  particleName == "anti_xi_c+" ||
426  particleName == "anti_xi-" ||
427  particleName == "deuteron" ||
428  particleName == "lambda_c+" ||
429  particleName == "omega-" ||
430  particleName == "sigma_c+" ||
431  particleName == "sigma_c++" ||
432  particleName == "sigma+" ||
433  particleName == "sigma-" ||
434  particleName == "tau+" ||
435  particleName == "tau-" ||
436  particleName == "triton" ||
437  particleName == "xi_c+" ||
438  particleName == "xi-" ) {
439 
440  ph->RegisterProcess(hmsc, particle);
441  ph->RegisterProcess(new G4hIonisation(), particle);
442  ph->RegisterProcess(pnuc, particle);
443  }
444  }
445 
446  // Nuclear stopping
447  pnuc->SetMaxKinEnergy(MeV);
448 
449  // Deexcitation
450  //
453 
454  G4EmModelActivator mact;
455  mact.ConstructProcess();
456 }
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 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 G4EmLivermorePhysics::verbose
private

Definition at line 51 of file G4EmLivermorePhysics.hh.


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