Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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$
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 "G4UrbanMscModel95.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 {
185 
186  // muon & hadron bremsstrahlung and pair production
195 
196  // muon & hadron multiple scattering
198  mumsc->AddEmModel(0, new G4WentzelVIModel());
200  pimsc->AddEmModel(0, new G4WentzelVIModel());
202  kmsc->AddEmModel(0, new G4WentzelVIModel());
204  pmsc->AddEmModel(0, new G4WentzelVIModel());
205  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
206 
207  // high energy limit for e+- scattering models
208  G4double highEnergyLimit = 100*MeV;
209 
210  // nuclear stopping
211  G4NuclearStopping* ionnuc = new G4NuclearStopping();
212  G4NuclearStopping* pnuc = new G4NuclearStopping();
213 
214  // Add Livermore EM Processes
216 
217  while( (*theParticleIterator)() ){
218 
220  G4String particleName = particle->GetParticleName();
221 
222  if(verbose > 1)
223  G4cout << "### " << GetPhysicsName() << " instantiates for "
224  << particleName << G4endl;
225 
226  //Applicability range for Livermore models
227  //for higher energies, the Standard models are used
228  G4double LivermoreHighEnergyLimit = GeV;
229 
230  if (particleName == "gamma") {
231 
232  // Photoelectric effect - define low-energy model
233  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
234  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
236  thePhotoElectricEffect->SetEmModel(theLivermorePhotoElectricModel);
237  ph->RegisterProcess(thePhotoElectricEffect, particle);
238 
239  // Compton scattering - define low-energy model
240  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
241  G4LivermoreComptonModel* theLivermoreComptonModel =
243  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
244  theComptonScattering->SetEmModel(theLivermoreComptonModel, 1);
245  ph->RegisterProcess(theComptonScattering, particle);
246 
247  // gamma conversion - define low-energy model
248  G4GammaConversion* theGammaConversion = new G4GammaConversion();
249  G4VEmModel* theLivermoreGammaConversionModel =
251  theGammaConversion->SetEmModel(theLivermoreGammaConversionModel, 1);
252  ph->RegisterProcess(theGammaConversion, particle);
253 
254  // default Rayleigh scattering is Livermore
255  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
256  ph->RegisterProcess(theRayleigh, particle);
257 
258  } else if (particleName == "e-") {
259 
260  // multiple scattering
263  G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
264  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
265  msc1->SetHighEnergyLimit(highEnergyLimit);
266  msc2->SetLowEnergyLimit(highEnergyLimit);
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  ph->RegisterProcess(msc, particle);
277  ph->RegisterProcess(ss, particle);
278 
279  // Ionisation - Livermore should be used only for low energies
280  G4eIonisation* eIoni = new G4eIonisation();
281  G4LivermoreIonisationModel* theIoniLivermore = new
283  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
284  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
285  eIoni->SetStepFunction(0.2, 100*um); //
286  ph->RegisterProcess(eIoni, particle);
287 
288  // Bremsstrahlung
289  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
290  G4VEmModel* theBremLivermore = new G4LivermoreBremsstrahlungModel();
291  theBremLivermore->SetHighEnergyLimit(1*GeV);
292  //theBremLivermore->SetAngularDistribution(new G4Generator2BS());
293  eBrem->SetEmModel(theBremLivermore,1);
294 
295  ph->RegisterProcess(eBrem, particle);
296 
297  } else if (particleName == "e+") {
298 
299  // Identical to G4EmStandardPhysics_option3
300 
301  // multiple scattering
304  G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
305  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
306  msc1->SetHighEnergyLimit(highEnergyLimit);
307  msc2->SetLowEnergyLimit(highEnergyLimit);
308  msc->AddEmModel(0, msc1);
309  msc->AddEmModel(0, msc2);
310 
313  ss->SetEmModel(ssm, 1);
314  ss->SetMinKinEnergy(highEnergyLimit);
315  ssm->SetLowEnergyLimit(highEnergyLimit);
316  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
317 
318  G4eIonisation* eIoni = new G4eIonisation();
319  eIoni->SetStepFunction(0.2, 100*um);
320 
321  ph->RegisterProcess(msc, particle);
322  ph->RegisterProcess(eIoni, particle);
323  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
324  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
325  ph->RegisterProcess(ss, particle);
326 
327  } else if (particleName == "mu+" ||
328  particleName == "mu-" ) {
329 
330  G4MuIonisation* muIoni = new G4MuIonisation();
331  muIoni->SetStepFunction(0.2, 50*um);
332 
333  ph->RegisterProcess(mumsc, particle);
334  ph->RegisterProcess(muIoni, particle);
335  ph->RegisterProcess(mub, particle);
336  ph->RegisterProcess(mup, particle);
337  ph->RegisterProcess(new G4CoulombScattering(), particle);
338 
339  } else if (particleName == "alpha" ||
340  particleName == "He3" ) {
341 
342  // Identical to G4EmStandardPhysics_option3
343 
345  G4ionIonisation* ionIoni = new G4ionIonisation();
346  ionIoni->SetStepFunction(0.1, 10*um);
347 
348  ph->RegisterProcess(msc, particle);
349  ph->RegisterProcess(ionIoni, particle);
350  ph->RegisterProcess(ionnuc, particle);
351 
352  } else if (particleName == "GenericIon") {
353 
354  // Identical to G4EmStandardPhysics_option3
355 
356  G4ionIonisation* ionIoni = new G4ionIonisation();
357  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
358  ionIoni->SetStepFunction(0.1, 1*um);
359 
360  ph->RegisterProcess(hmsc, particle);
361  ph->RegisterProcess(ionIoni, particle);
362  ph->RegisterProcess(ionnuc, particle);
363 
364  } else if (particleName == "pi+" ||
365  particleName == "pi-" ) {
366 
367  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
368  G4hIonisation* hIoni = new G4hIonisation();
369  hIoni->SetStepFunction(0.2, 50*um);
370 
371  ph->RegisterProcess(pimsc, particle);
372  ph->RegisterProcess(hIoni, particle);
373  ph->RegisterProcess(pib, particle);
374  ph->RegisterProcess(pip, particle);
375 
376  } else if (particleName == "kaon+" ||
377  particleName == "kaon-" ) {
378 
379  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
380  G4hIonisation* hIoni = new G4hIonisation();
381  hIoni->SetStepFunction(0.2, 50*um);
382 
383  ph->RegisterProcess(kmsc, particle);
384  ph->RegisterProcess(hIoni, particle);
385  ph->RegisterProcess(kb, particle);
386  ph->RegisterProcess(kp, particle);
387 
388  } else if (particleName == "proton" ||
389  particleName == "anti_proton") {
390 
391  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
392  G4hIonisation* hIoni = new G4hIonisation();
393  hIoni->SetStepFunction(0.2, 50*um);
394 
395  ph->RegisterProcess(pmsc, particle);
396  ph->RegisterProcess(hIoni, particle);
397  ph->RegisterProcess(pb, particle);
398  ph->RegisterProcess(pp, particle);
399  ph->RegisterProcess(pnuc, particle);
400 
401  } else if (particleName == "B+" ||
402  particleName == "B-" ||
403  particleName == "D+" ||
404  particleName == "D-" ||
405  particleName == "Ds+" ||
406  particleName == "Ds-" ||
407  particleName == "anti_He3" ||
408  particleName == "anti_alpha" ||
409  particleName == "anti_deuteron" ||
410  particleName == "anti_lambda_c+" ||
411  particleName == "anti_omega-" ||
412  particleName == "anti_sigma_c+" ||
413  particleName == "anti_sigma_c++" ||
414  particleName == "anti_sigma+" ||
415  particleName == "anti_sigma-" ||
416  particleName == "anti_triton" ||
417  particleName == "anti_xi_c+" ||
418  particleName == "anti_xi-" ||
419  particleName == "deuteron" ||
420  particleName == "lambda_c+" ||
421  particleName == "omega-" ||
422  particleName == "sigma_c+" ||
423  particleName == "sigma_c++" ||
424  particleName == "sigma+" ||
425  particleName == "sigma-" ||
426  particleName == "tau+" ||
427  particleName == "tau-" ||
428  particleName == "triton" ||
429  particleName == "xi_c+" ||
430  particleName == "xi-" ) {
431 
432  // Identical to G4EmStandardPhysics_option3
433 
434  ph->RegisterProcess(hmsc, particle);
435  ph->RegisterProcess(new G4hIonisation(), particle);
436  ph->RegisterProcess(pnuc, particle);
437  }
438  }
439 
440  // Em options
441  //
442  G4EmProcessOptions opt;
443  opt.SetVerbose(verbose);
444 
445  // Multiple Coulomb scattering
446  //
447  opt.SetPolarAngleLimit(CLHEP::pi);
448 
449  // Physics tables
450  //
451 
452  opt.SetMinEnergy(100*eV);
453  opt.SetMaxEnergy(10*TeV);
454  opt.SetDEDXBinning(220);
455  opt.SetLambdaBinning(220);
456 
457  // Nuclear stopping
458  pnuc->SetMaxKinEnergy(MeV);
459 
460  // Ionization
461  //
462  //opt.SetSubCutoff(true);
463 
464  // Deexcitation
465  //
468  de->SetFluo(true);
469 }
470 
471 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......