Geant4  10.02
G4EmLivermorePhysics.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4EmLivermorePhysics.cc 92821 2015-09-17 15:23:49Z gcosmo $
27 
28 #include "G4EmLivermorePhysics.hh"
29 #include "G4ParticleDefinition.hh"
30 #include "G4SystemOfUnits.hh"
31 
32 // *** Processes and models
33 
34 // gamma
35 #include "G4PhotoElectricEffect.hh"
37 
38 #include "G4ComptonScattering.hh"
40 
41 #include "G4GammaConversion.hh"
43 
44 #include "G4RayleighScattering.hh"
46 
47 // e+-
48 #include "G4eMultipleScattering.hh"
50 
51 #include "G4eIonisation.hh"
53 
54 #include "G4eBremsstrahlung.hh"
56 #include "G4Generator2BS.hh"
57 
58 // e+
59 #include "G4eplusAnnihilation.hh"
60 
61 // mu+-
63 #include "G4MuIonisation.hh"
64 #include "G4MuBremsstrahlung.hh"
65 #include "G4MuPairProduction.hh"
66 
71 
72 // hadrons
73 #include "G4hMultipleScattering.hh"
74 #include "G4MscStepLimitType.hh"
75 
76 #include "G4hBremsstrahlung.hh"
77 #include "G4hPairProduction.hh"
78 
79 #include "G4hIonisation.hh"
80 #include "G4ionIonisation.hh"
81 #include "G4alphaIonisation.hh"
83 #include "G4NuclearStopping.hh"
84 
85 // msc models
86 #include "G4UrbanMscModel.hh"
87 #include "G4WentzelVIModel.hh"
89 #include "G4CoulombScattering.hh"
91 
92 // interfaces
93 #include "G4LossTableManager.hh"
94 #include "G4UAtomicDeexcitation.hh"
95 #include "G4EmParameters.hh"
96 
97 // particles
98 
99 #include "G4Gamma.hh"
100 #include "G4Electron.hh"
101 #include "G4Positron.hh"
102 #include "G4MuonPlus.hh"
103 #include "G4MuonMinus.hh"
104 #include "G4PionPlus.hh"
105 #include "G4PionMinus.hh"
106 #include "G4KaonPlus.hh"
107 #include "G4KaonMinus.hh"
108 #include "G4Proton.hh"
109 #include "G4AntiProton.hh"
110 #include "G4Deuteron.hh"
111 #include "G4Triton.hh"
112 #include "G4He3.hh"
113 #include "G4Alpha.hh"
114 #include "G4GenericIon.hh"
115 
116 //
117 #include "G4PhysicsListHelper.hh"
118 #include "G4BuilderType.hh"
119 #include "G4EmModelActivator.hh"
120 
121 // factory
123 //
125 
126 
127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128 
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 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147 
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 }
164 
165 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
166 
168 {}
169 
170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171 
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 }
204 
205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
206 
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  aParticleIterator->reset();
243 
244  while( (*aParticleIterator)() ){
245 
246  G4ParticleDefinition* particle = aParticleIterator->value();
247  G4String particleName = particle->GetParticleName();
248 
249  if (particleName == "gamma") {
250 
251  // photoelectric effect - Livermore model only
252  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
253  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel(), 1);
254  ph->RegisterProcess(thePhotoElectricEffect, particle);
255 
256  // Compton scattering - Livermore model only
257  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
258  theComptonScattering->SetEmModel(new G4LivermoreComptonModel(),1);
259  ph->RegisterProcess(theComptonScattering, particle);
260 
261  // gamma conversion - Livermore model below 80 GeV
262  G4GammaConversion* theGammaConversion = new G4GammaConversion();
263  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel(),1);
264  ph->RegisterProcess(theGammaConversion, particle);
265 
266  // default Rayleigh scattering is Livermore
267  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
268  ph->RegisterProcess(theRayleigh, particle);
269 
270  } else if (particleName == "e-") {
271 
272  // multiple scattering
274  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
275  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
276  msc1->SetHighEnergyLimit(highEnergyLimit);
277  msc2->SetLowEnergyLimit(highEnergyLimit);
278  msc->AddEmModel(0, msc1);
279  msc->AddEmModel(0, msc2);
280 
283  ss->SetEmModel(ssm, 1);
284  ss->SetMinKinEnergy(highEnergyLimit);
285  ssm->SetLowEnergyLimit(highEnergyLimit);
286  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
287 
288  // Ionisation - Livermore should be used only for low energies
289  G4eIonisation* eIoni = new G4eIonisation();
290  G4LivermoreIonisationModel* theIoniLivermore = new
292  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
293  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
294  eIoni->SetStepFunction(0.2, 100*um); //
295 
296  // Bremsstrahlung
297  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
298  G4VEmModel* theBremLivermore = new G4LivermoreBremsstrahlungModel();
299  theBremLivermore->SetHighEnergyLimit(1*GeV);
300  theBremLivermore->SetAngularDistribution(new G4Generator2BS());
301  eBrem->SetEmModel(theBremLivermore,1);
302 
303  // register processes
304  ph->RegisterProcess(msc, particle);
305  ph->RegisterProcess(eIoni, particle);
306  ph->RegisterProcess(eBrem, particle);
307  ph->RegisterProcess(ss, particle);
308 
309  } else if (particleName == "e+") {
310 
311  // multiple scattering
313  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
314  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
315  msc1->SetHighEnergyLimit(highEnergyLimit);
316  msc2->SetLowEnergyLimit(highEnergyLimit);
317  msc->AddEmModel(0, msc1);
318  msc->AddEmModel(0, msc2);
319 
322  ss->SetEmModel(ssm, 1);
323  ss->SetMinKinEnergy(highEnergyLimit);
324  ssm->SetLowEnergyLimit(highEnergyLimit);
325  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
326 
327  // ionisation
328  G4eIonisation* eIoni = new G4eIonisation();
329  eIoni->SetStepFunction(0.2, 100*um);
330 
331  // register processes
332  ph->RegisterProcess(msc, particle);
333  ph->RegisterProcess(eIoni, particle);
334  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
335  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
336  ph->RegisterProcess(ss, particle);
337 
338  } else if (particleName == "mu+" ||
339  particleName == "mu-" ) {
340 
341  G4MuIonisation* muIoni = new G4MuIonisation();
342  muIoni->SetStepFunction(0.2, 50*um);
343 
344  ph->RegisterProcess(mumsc, particle);
345  ph->RegisterProcess(muIoni, particle);
346  ph->RegisterProcess(mub, particle);
347  ph->RegisterProcess(mup, particle);
348  ph->RegisterProcess(new G4CoulombScattering(), particle);
349 
350  } else if (particleName == "alpha" ||
351  particleName == "He3" ) {
352 
354  G4ionIonisation* ionIoni = new G4ionIonisation();
355  ionIoni->SetStepFunction(0.1, 10*um);
356 
357  ph->RegisterProcess(msc, particle);
358  ph->RegisterProcess(ionIoni, particle);
359  ph->RegisterProcess(pnuc, particle);
360 
361  } else if (particleName == "GenericIon") {
362 
363  G4ionIonisation* ionIoni = new G4ionIonisation();
364  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
365  ionIoni->SetStepFunction(0.1, 1*um);
366 
367  ph->RegisterProcess(hmsc, particle);
368  ph->RegisterProcess(ionIoni, particle);
369  ph->RegisterProcess(pnuc, particle);
370 
371  } else if (particleName == "pi+" ||
372  particleName == "pi-" ) {
373 
375  G4hIonisation* hIoni = new G4hIonisation();
376  hIoni->SetStepFunction(0.2, 50*um);
377 
378  ph->RegisterProcess(pimsc, particle);
379  ph->RegisterProcess(hIoni, particle);
380  ph->RegisterProcess(pib, particle);
381  ph->RegisterProcess(pip, particle);
382 
383  } else if (particleName == "kaon+" ||
384  particleName == "kaon-" ) {
385 
387  G4hIonisation* hIoni = new G4hIonisation();
388  hIoni->SetStepFunction(0.2, 50*um);
389 
390  ph->RegisterProcess(kmsc, particle);
391  ph->RegisterProcess(hIoni, particle);
392  ph->RegisterProcess(kb, particle);
393  ph->RegisterProcess(kp, particle);
394 
395  } else if (particleName == "proton" ||
396  particleName == "anti_proton") {
397 
399  G4hIonisation* hIoni = new G4hIonisation();
400  hIoni->SetStepFunction(0.2, 50*um);
401 
402  ph->RegisterProcess(pmsc, particle);
403  ph->RegisterProcess(hIoni, particle);
404  ph->RegisterProcess(pb, particle);
405  ph->RegisterProcess(pp, particle);
406  ph->RegisterProcess(pnuc, particle);
407 
408  } else if (particleName == "B+" ||
409  particleName == "B-" ||
410  particleName == "D+" ||
411  particleName == "D-" ||
412  particleName == "Ds+" ||
413  particleName == "Ds-" ||
414  particleName == "anti_He3" ||
415  particleName == "anti_alpha" ||
416  particleName == "anti_deuteron" ||
417  particleName == "anti_lambda_c+" ||
418  particleName == "anti_omega-" ||
419  particleName == "anti_sigma_c+" ||
420  particleName == "anti_sigma_c++" ||
421  particleName == "anti_sigma+" ||
422  particleName == "anti_sigma-" ||
423  particleName == "anti_triton" ||
424  particleName == "anti_xi_c+" ||
425  particleName == "anti_xi-" ||
426  particleName == "deuteron" ||
427  particleName == "lambda_c+" ||
428  particleName == "omega-" ||
429  particleName == "sigma_c+" ||
430  particleName == "sigma_c++" ||
431  particleName == "sigma+" ||
432  particleName == "sigma-" ||
433  particleName == "tau+" ||
434  particleName == "tau-" ||
435  particleName == "triton" ||
436  particleName == "xi_c+" ||
437  particleName == "xi-" ) {
438 
439  ph->RegisterProcess(hmsc, particle);
440  ph->RegisterProcess(new G4hIonisation(), particle);
441  ph->RegisterProcess(pnuc, particle);
442  }
443  }
444 
445  // Nuclear stopping
446  pnuc->SetMaxKinEnergy(MeV);
447 
448  // Deexcitation
449  //
452 
453  G4EmModelActivator mact;
454  mact.ConstructProcess();
455 }
456 
457 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmLivermorePhysics)
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
static const double MeV
Definition: G4SIunits.hh:211
void SetVerbose(G4int val)
virtual void ConstructParticle()
static G4LossTableManager * Instance()
void SetLowestElectronEnergy(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
void SetStepFunction(G4double v1, G4double v2)
int G4int
Definition: G4Types.hh:78
void SetMaxEnergy(G4double val)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
void SetMscRangeFactor(G4double val)
#define aParticleIterator
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
void SetNumberOfBinsPerDecade(G4int val)
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:93
static const double GeV
Definition: G4SIunits.hh:214
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
const G4String & GetPhysicsName() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
static const double eV
Definition: G4SIunits.hh:212
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetMaxKinEnergy(G4double e)
void SetMinEnergy(G4double val)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4EmParameters * Instance()
void SetEmModel(G4VEmModel *, G4int index=1)
static const double um
Definition: G4SIunits.hh:112
void ActivateAngularGeneratorForIonisation(G4bool val)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static const double TeV
Definition: G4SIunits.hh:215
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static G4He3 * He3()
Definition: G4He3.cc:94
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetFluo(G4bool val)