Geant4  10.00.p01
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 79157 2014-02-19 15:35:01Z 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 "G4EmProcessOptions.hh"
95 #include "G4UAtomicDeexcitation.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 
120 // factory
122 //
124 
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127 
129  : G4VPhysicsConstructor("G4EmLivermorePhysics"), verbose(ver)
130 {
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136 
138  : G4VPhysicsConstructor("G4EmLivermorePhysics"), verbose(ver)
139 {
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145 
147 {}
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150 
152 {
153 // gamma
154  G4Gamma::Gamma();
155 
156 // leptons
161 
162 // mesons
167 
168 // baryons
171 
172 // ions
175  G4He3::He3();
176  G4Alpha::Alpha();
178 }
179 
180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181 
183 {
184  if(verbose > 1) {
185  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
186  }
188 
189  // muon & hadron bremsstrahlung and pair production
198 
199  // muon & hadron multiple scattering
201  mumsc->AddEmModel(0, new G4WentzelVIModel());
203  //pimsc->AddEmModel(0, new G4WentzelVIModel());
205  //kmsc->AddEmModel(0, new G4WentzelVIModel());
207  //pmsc->AddEmModel(0, new G4WentzelVIModel());
208  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
209 
210  // high energy limit for e+- scattering models
211  G4double highEnergyLimit = 100*MeV;
212 
213  // nuclear stopping
214  G4NuclearStopping* ionnuc = new G4NuclearStopping();
215  G4NuclearStopping* pnuc = new G4NuclearStopping();
216 
217  // Add Livermore EM Processes
218  aParticleIterator->reset();
219 
220  while( (*aParticleIterator)() ){
221 
222  G4ParticleDefinition* particle = aParticleIterator->value();
223  G4String particleName = particle->GetParticleName();
224 
225  //Applicability range for Livermore models
226  //for higher energies, the Standard models are used
227  G4double LivermoreHighEnergyLimit = GeV;
228 
229  if (particleName == "gamma") {
230 
231  // Photoelectric effect - define low-energy model
232  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
233  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
235  thePhotoElectricEffect->SetEmModel(theLivermorePhotoElectricModel);
236  ph->RegisterProcess(thePhotoElectricEffect, particle);
237 
238  // Compton scattering - define low-energy model
239  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
240  G4LivermoreComptonModel* theLivermoreComptonModel =
242  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
243  theComptonScattering->SetEmModel(theLivermoreComptonModel, 1);
244  ph->RegisterProcess(theComptonScattering, particle);
245 
246  // gamma conversion - define low-energy model
247  G4GammaConversion* theGammaConversion = new G4GammaConversion();
248  G4VEmModel* theLivermoreGammaConversionModel =
250  theGammaConversion->SetEmModel(theLivermoreGammaConversionModel, 1);
251  ph->RegisterProcess(theGammaConversion, particle);
252 
253  // default Rayleigh scattering is Livermore
254  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
255  ph->RegisterProcess(theRayleigh, particle);
256 
257  } else if (particleName == "e-") {
258 
259  // multiple scattering
262  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
263  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
264  msc1->SetHighEnergyLimit(highEnergyLimit);
265  msc2->SetLowEnergyLimit(highEnergyLimit);
266  msc->SetRangeFactor(0.01);
267  msc->AddEmModel(0, msc1);
268  msc->AddEmModel(0, msc2);
269 
272  ss->SetEmModel(ssm, 1);
273  ss->SetMinKinEnergy(highEnergyLimit);
274  ssm->SetLowEnergyLimit(highEnergyLimit);
275  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
276 
277  // Ionisation - Livermore should be used only for low energies
278  G4eIonisation* eIoni = new G4eIonisation();
279  G4LivermoreIonisationModel* theIoniLivermore = new
281  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
282  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
283  eIoni->SetStepFunction(0.2, 100*um); //
284 
285  // Bremsstrahlung
286  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
287  G4VEmModel* theBremLivermore = new G4LivermoreBremsstrahlungModel();
288  theBremLivermore->SetHighEnergyLimit(1*GeV);
289  //theBremLivermore->SetAngularDistribution(new G4Generator2BS());
290  eBrem->SetEmModel(theBremLivermore,1);
291 
292  // register processes
293  ph->RegisterProcess(msc, particle);
294  ph->RegisterProcess(eIoni, particle);
295  ph->RegisterProcess(eBrem, particle);
296  ph->RegisterProcess(ss, particle);
297 
298  } else if (particleName == "e+") {
299 
300  // Identical to G4EmStandardPhysics_option3
301 
302  // multiple scattering
305  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
306  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
307  msc1->SetHighEnergyLimit(highEnergyLimit);
308  msc2->SetLowEnergyLimit(highEnergyLimit);
309  msc->SetRangeFactor(0.01);
310  msc->AddEmModel(0, msc1);
311  msc->AddEmModel(0, msc2);
312 
315  ss->SetEmModel(ssm, 1);
316  ss->SetMinKinEnergy(highEnergyLimit);
317  ssm->SetLowEnergyLimit(highEnergyLimit);
318  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
319 
320  G4eIonisation* eIoni = new G4eIonisation();
321  eIoni->SetStepFunction(0.2, 100*um);
322 
323  // register processes
324  ph->RegisterProcess(msc, particle);
325  ph->RegisterProcess(eIoni, particle);
326  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
327  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
328  ph->RegisterProcess(ss, particle);
329 
330  } else if (particleName == "mu+" ||
331  particleName == "mu-" ) {
332 
333  G4MuIonisation* muIoni = new G4MuIonisation();
334  muIoni->SetStepFunction(0.2, 50*um);
335 
336  ph->RegisterProcess(mumsc, particle);
337  ph->RegisterProcess(muIoni, particle);
338  ph->RegisterProcess(mub, particle);
339  ph->RegisterProcess(mup, particle);
340  ph->RegisterProcess(new G4CoulombScattering(), particle);
341 
342  } else if (particleName == "alpha" ||
343  particleName == "He3" ) {
344 
345  // Identical to G4EmStandardPhysics_option3
346 
348  G4ionIonisation* ionIoni = new G4ionIonisation();
349  ionIoni->SetStepFunction(0.1, 10*um);
350 
351  ph->RegisterProcess(msc, particle);
352  ph->RegisterProcess(ionIoni, particle);
353  ph->RegisterProcess(ionnuc, particle);
354 
355  } else if (particleName == "GenericIon") {
356 
357  // Identical to G4EmStandardPhysics_option3
358 
359  G4ionIonisation* ionIoni = new G4ionIonisation();
360  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
361  ionIoni->SetStepFunction(0.1, 1*um);
362 
363  ph->RegisterProcess(hmsc, particle);
364  ph->RegisterProcess(ionIoni, particle);
365  ph->RegisterProcess(ionnuc, particle);
366 
367  } else if (particleName == "pi+" ||
368  particleName == "pi-" ) {
369 
370  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
371  G4hIonisation* hIoni = new G4hIonisation();
372  hIoni->SetStepFunction(0.2, 50*um);
373 
374  ph->RegisterProcess(pimsc, particle);
375  ph->RegisterProcess(hIoni, particle);
376  ph->RegisterProcess(pib, particle);
377  ph->RegisterProcess(pip, particle);
378 
379  } else if (particleName == "kaon+" ||
380  particleName == "kaon-" ) {
381 
382  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
383  G4hIonisation* hIoni = new G4hIonisation();
384  hIoni->SetStepFunction(0.2, 50*um);
385 
386  ph->RegisterProcess(kmsc, particle);
387  ph->RegisterProcess(hIoni, particle);
388  ph->RegisterProcess(kb, particle);
389  ph->RegisterProcess(kp, particle);
390 
391  } else if (particleName == "proton" ||
392  particleName == "anti_proton") {
393 
394  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
395  G4hIonisation* hIoni = new G4hIonisation();
396  hIoni->SetStepFunction(0.2, 50*um);
397 
398  ph->RegisterProcess(pmsc, particle);
399  ph->RegisterProcess(hIoni, particle);
400  ph->RegisterProcess(pb, particle);
401  ph->RegisterProcess(pp, particle);
402  ph->RegisterProcess(pnuc, particle);
403 
404  } else if (particleName == "B+" ||
405  particleName == "B-" ||
406  particleName == "D+" ||
407  particleName == "D-" ||
408  particleName == "Ds+" ||
409  particleName == "Ds-" ||
410  particleName == "anti_He3" ||
411  particleName == "anti_alpha" ||
412  particleName == "anti_deuteron" ||
413  particleName == "anti_lambda_c+" ||
414  particleName == "anti_omega-" ||
415  particleName == "anti_sigma_c+" ||
416  particleName == "anti_sigma_c++" ||
417  particleName == "anti_sigma+" ||
418  particleName == "anti_sigma-" ||
419  particleName == "anti_triton" ||
420  particleName == "anti_xi_c+" ||
421  particleName == "anti_xi-" ||
422  particleName == "deuteron" ||
423  particleName == "lambda_c+" ||
424  particleName == "omega-" ||
425  particleName == "sigma_c+" ||
426  particleName == "sigma_c++" ||
427  particleName == "sigma+" ||
428  particleName == "sigma-" ||
429  particleName == "tau+" ||
430  particleName == "tau-" ||
431  particleName == "triton" ||
432  particleName == "xi_c+" ||
433  particleName == "xi-" ) {
434 
435  // Identical to G4EmStandardPhysics_option3
436 
437  ph->RegisterProcess(hmsc, particle);
438  ph->RegisterProcess(new G4hIonisation(), particle);
439  ph->RegisterProcess(pnuc, particle);
440  }
441  }
442 
443  // Em options
444  //
445  G4EmProcessOptions opt;
446  opt.SetVerbose(verbose);
447 
448  // Multiple Coulomb scattering
449  //
451 
452  // Physics tables
453  //
454 
455  opt.SetMinEnergy(100*eV);
456  opt.SetMaxEnergy(10*TeV);
457  opt.SetDEDXBinning(220);
458  opt.SetLambdaBinning(220);
459 
460  // Nuclear stopping
461  pnuc->SetMaxKinEnergy(MeV);
462 
463  // Ionization
464  //
465  //opt.SetSubCutoff(true);
466 
467  // Deexcitation
468  //
471  de->SetFluo(true);
472 }
473 
474 //....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:193
virtual void ConstructParticle()
static G4LossTableManager * Instance()
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
const G4double pi
void SetMinEnergy(G4double val)
void SetStepFunction(G4double v1, G4double v2)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
void SetDEDXBinning(G4int val)
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
void SetLambdaBinning(G4int val)
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
#define aParticleIterator
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
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:196
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 SetMaxEnergy(G4double val)
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:704
static const double eV
Definition: G4SIunits.hh:194
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
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:197
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
void SetRangeFactor(G4double val)
static G4He3 * He3()
Definition: G4He3.cc:94
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetStepLimitType(G4MscStepLimitType val)
void SetVerbose(G4int val, const G4String &name="all", G4bool worker=false)
void SetPolarAngleLimit(G4double val)