Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmPenelopePhysics.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 "G4EmPenelopePhysics.hh"
29 #include "G4ParticleDefinition.hh"
30 #include "G4SystemOfUnits.hh"
31 
32 // *** Processes and models
33 
34 // gamma
35 
36 #include "G4PhotoElectricEffect.hh"
38 
39 #include "G4ComptonScattering.hh"
41 
42 #include "G4GammaConversion.hh"
44 
45 #include "G4RayleighScattering.hh"
47 
48 // e- and e+
49 
50 #include "G4eMultipleScattering.hh"
52 
53 #include "G4eIonisation.hh"
55 
56 #include "G4eBremsstrahlung.hh"
58 
59 // e+ only
60 
61 #include "G4eplusAnnihilation.hh"
63 
64 // mu
65 
67 #include "G4MuIonisation.hh"
68 #include "G4MuBremsstrahlung.hh"
69 #include "G4MuPairProduction.hh"
70 
75 
76 // hadrons
77 
78 #include "G4hMultipleScattering.hh"
79 #include "G4MscStepLimitType.hh"
80 
81 #include "G4hBremsstrahlung.hh"
82 #include "G4hPairProduction.hh"
83 
84 #include "G4hIonisation.hh"
85 #include "G4ionIonisation.hh"
86 #include "G4alphaIonisation.hh"
88 #include "G4NuclearStopping.hh"
89 
90 // msc models
91 #include "G4UrbanMscModel93.hh"
92 #include "G4UrbanMscModel95.hh"
94 #include "G4WentzelVIModel.hh"
95 #include "G4CoulombScattering.hh"
97 //
98 
99 #include "G4LossTableManager.hh"
100 #include "G4VAtomDeexcitation.hh"
101 #include "G4UAtomicDeexcitation.hh"
102 #include "G4EmProcessOptions.hh"
103 
104 // particles
105 
106 #include "G4Gamma.hh"
107 #include "G4Electron.hh"
108 #include "G4Positron.hh"
109 #include "G4MuonPlus.hh"
110 #include "G4MuonMinus.hh"
111 #include "G4PionPlus.hh"
112 #include "G4PionMinus.hh"
113 #include "G4KaonPlus.hh"
114 #include "G4KaonMinus.hh"
115 #include "G4Proton.hh"
116 #include "G4AntiProton.hh"
117 #include "G4Deuteron.hh"
118 #include "G4Triton.hh"
119 #include "G4He3.hh"
120 #include "G4Alpha.hh"
121 #include "G4GenericIon.hh"
122 
123 //
124 #include "G4PhysicsListHelper.hh"
125 #include "G4BuilderType.hh"
126 
127 // factory
129 //
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133 
135  : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
136 {
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142 
144  : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
145 {
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151 
153 {}
154 
155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
156 
158 {
159 // gamma
160  G4Gamma::Gamma();
161 
162 // leptons
167 
168 // mesons
173 
174 // baryons
177 
178 // ions
181  G4He3::He3();
182  G4Alpha::Alpha();
184 }
185 
186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
187 
189 {
191 
192  // muon & hadron bremsstrahlung and pair production
201 
202  // muon & hadron multiple scattering
204  mumsc->AddEmModel(0, new G4WentzelVIModel());
206  pimsc->AddEmModel(0, new G4WentzelVIModel());
208  kmsc->AddEmModel(0, new G4WentzelVIModel());
210  pmsc->AddEmModel(0, new G4WentzelVIModel());
211  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
212 
213  // high energy limit for e+- scattering models
214  G4double highEnergyLimit = 100*MeV;
215 
216  // nuclear stopping
217  G4NuclearStopping* ionnuc = new G4NuclearStopping();
218  G4NuclearStopping* pnuc = new G4NuclearStopping();
219 
220  // Add Penelope EM Processes
222 
223  while( (*theParticleIterator)() ){
224 
226  G4String particleName = particle->GetParticleName();
227 
228  if(verbose > 1)
229  G4cout << "### " << GetPhysicsName() << " instantiates for "
230  << particleName << G4endl;
231 
232  //Applicability range for Penelope models
233  //for higher energies, the Standard models are used
234  G4double PenelopeHighEnergyLimit = 1.0*GeV;
235 
236  if (particleName == "gamma") {
237 
238  //Photo-electric effect
239  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
240  G4PenelopePhotoElectricModel* thePEPenelopeModel = new
242  thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
243  thePhotoElectricEffect->SetEmModel(thePEPenelopeModel, 1);
244  ph->RegisterProcess(thePhotoElectricEffect, particle);
245 
246  //Compton scattering
247  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
248  G4PenelopeComptonModel* theComptonPenelopeModel =
250  theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
251  theComptonScattering->SetEmModel(theComptonPenelopeModel, 1);
252  ph->RegisterProcess(theComptonScattering, particle);
253 
254  //Gamma conversion
255  G4GammaConversion* theGammaConversion = new G4GammaConversion();
256  G4PenelopeGammaConversionModel* theGCPenelopeModel =
258  theGammaConversion->SetEmModel(theGCPenelopeModel,1);
259  ph->RegisterProcess(theGammaConversion, particle);
260 
261  //Rayleigh scattering
262  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
263  G4PenelopeRayleighModel* theRayleighPenelopeModel =
265  //theRayleighPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
266  theRayleigh->SetEmModel(theRayleighPenelopeModel, 1);
267  ph->RegisterProcess(theRayleigh, particle);
268 
269  } else if (particleName == "e-") {
270 
271  // multiple scattering
274  G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
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  ph->RegisterProcess(msc, particle);
288  ph->RegisterProcess(ss, particle);
289 
290  //Ionisation
291  G4eIonisation* eIoni = new G4eIonisation();
292  G4PenelopeIonisationModel* theIoniPenelope =
294  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
295  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
296  eIoni->SetStepFunction(0.2, 100*um); //
297  ph->RegisterProcess(eIoni, particle);
298 
299  //Bremsstrahlung
300  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
301  G4PenelopeBremsstrahlungModel* theBremPenelope = new
303  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
304  eBrem->AddEmModel(0,theBremPenelope);
305  ph->RegisterProcess(eBrem, particle);
306 
307  } else if (particleName == "e+") {
308 
309  // multiple scattering
312  G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
313  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
314  msc1->SetHighEnergyLimit(highEnergyLimit);
315  msc2->SetLowEnergyLimit(highEnergyLimit);
316  msc->AddEmModel(0, msc1);
317  msc->AddEmModel(0, msc2);
318 
321  ss->SetEmModel(ssm, 1);
322  ss->SetMinKinEnergy(highEnergyLimit);
323  ssm->SetLowEnergyLimit(highEnergyLimit);
324  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
325  ph->RegisterProcess(msc, particle);
326  ph->RegisterProcess(ss, particle);
327 
328  //Ionisation
329  G4eIonisation* eIoni = new G4eIonisation();
330  G4PenelopeIonisationModel* theIoniPenelope =
332  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
333  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
334  eIoni->SetStepFunction(0.2, 100*um); //
335  ph->RegisterProcess(eIoni, particle);
336 
337  //Bremsstrahlung
338  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
339  G4PenelopeBremsstrahlungModel* theBremPenelope = new
341  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
342  eBrem->AddEmModel(0,theBremPenelope);
343  ph->RegisterProcess(eBrem, particle);
344 
345  //Annihilation
347  G4PenelopeAnnihilationModel* theAnnPenelope = new
349  theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
350  eAnni->AddEmModel(0,theAnnPenelope);
351  ph->RegisterProcess(eAnni, particle);
352 
353  } else if (particleName == "mu+" ||
354  particleName == "mu-" ) {
355 
356  // Identical to G4EmStandardPhysics_option3
357 
358  G4MuIonisation* muIoni = new G4MuIonisation();
359  muIoni->SetStepFunction(0.2, 50*um);
360 
361  ph->RegisterProcess(mumsc, particle);
362  ph->RegisterProcess(muIoni, particle);
363  ph->RegisterProcess(mub, particle);
364  ph->RegisterProcess(mup, particle);
365  ph->RegisterProcess(new G4CoulombScattering(), particle);
366 
367  } else if (particleName == "alpha" ||
368  particleName == "He3" ) {
369 
370  // Identical to G4EmStandardPhysics_option3
371 
373  G4ionIonisation* ionIoni = new G4ionIonisation();
374  ionIoni->SetStepFunction(0.1, 10*um);
375 
376  ph->RegisterProcess(msc, particle);
377  ph->RegisterProcess(ionIoni, particle);
378  ph->RegisterProcess(ionnuc, particle);
379 
380  } else if (particleName == "GenericIon") {
381 
382  // Identical to G4EmStandardPhysics_option3
383 
384  G4ionIonisation* ionIoni = new G4ionIonisation();
385  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
386  ionIoni->SetStepFunction(0.1, 1*um);
387 
388  ph->RegisterProcess(hmsc, particle);
389  ph->RegisterProcess(ionIoni, particle);
390  ph->RegisterProcess(ionnuc, particle);
391 
392  } else if (particleName == "pi+" ||
393  particleName == "pi-" ) {
394 
395  G4hIonisation* hIoni = new G4hIonisation();
396  hIoni->SetStepFunction(0.2, 50*um);
397 
398  ph->RegisterProcess(pimsc, particle);
399  ph->RegisterProcess(hIoni, particle);
400  ph->RegisterProcess(pib, particle);
401  ph->RegisterProcess(pip, particle);
402 
403  } else if (particleName == "kaon+" ||
404  particleName == "kaon-" ) {
405 
406  G4hIonisation* hIoni = new G4hIonisation();
407  hIoni->SetStepFunction(0.2, 50*um);
408 
409  ph->RegisterProcess(kmsc, particle);
410  ph->RegisterProcess(hIoni, particle);
411  ph->RegisterProcess(kb, particle);
412  ph->RegisterProcess(kp, particle);
413 
414  } else if (particleName == "proton" ||
415  particleName == "anti_proton") {
416 
417  G4hIonisation* hIoni = new G4hIonisation();
418  hIoni->SetStepFunction(0.2, 50*um);
419 
420  ph->RegisterProcess(pmsc, particle);
421  ph->RegisterProcess(hIoni, particle);
422  ph->RegisterProcess(pb, particle);
423  ph->RegisterProcess(pp, particle);
424  ph->RegisterProcess(pnuc, particle);
425 
426  } else if (particleName == "B+" ||
427  particleName == "B-" ||
428  particleName == "D+" ||
429  particleName == "D-" ||
430  particleName == "Ds+" ||
431  particleName == "Ds-" ||
432  particleName == "anti_He3" ||
433  particleName == "anti_alpha" ||
434  particleName == "anti_deuteron" ||
435  particleName == "anti_lambda_c+" ||
436  particleName == "anti_omega-" ||
437  particleName == "anti_sigma_c+" ||
438  particleName == "anti_sigma_c++" ||
439  particleName == "anti_sigma+" ||
440  particleName == "anti_sigma-" ||
441  particleName == "anti_triton" ||
442  particleName == "anti_xi_c+" ||
443  particleName == "anti_xi-" ||
444  particleName == "deuteron" ||
445  particleName == "lambda_c+" ||
446  particleName == "omega-" ||
447  particleName == "sigma_c+" ||
448  particleName == "sigma_c++" ||
449  particleName == "sigma+" ||
450  particleName == "sigma-" ||
451  particleName == "tau+" ||
452  particleName == "tau-" ||
453  particleName == "triton" ||
454  particleName == "xi_c+" ||
455  particleName == "xi-" ) {
456 
457  // Identical to G4EmStandardPhysics_option3
458 
459  ph->RegisterProcess(hmsc, particle);
460  ph->RegisterProcess(new G4hIonisation(), particle);
461  ph->RegisterProcess(pnuc, particle);
462  }
463  }
464 
465  // Em options
466  //
467  G4EmProcessOptions opt;
468  opt.SetVerbose(verbose);
469 
470  // Multiple Coulomb scattering
471  //
472  //opt.SetMscStepLimitation(fUseDistanceToBoundary);
473  //opt.SetMscRangeFactor(0.02);
474 
475  // Physics tables
476  //
477 
478  opt.SetMinEnergy(100*eV);
479  opt.SetMaxEnergy(10*TeV);
480  opt.SetDEDXBinning(220);
481  opt.SetLambdaBinning(220);
482 
483  // Nuclear stopping
484  pnuc->SetMaxKinEnergy(MeV);
485 
486  //opt.SetSplineFlag(true);
487  opt.SetPolarAngleLimit(CLHEP::pi);
488 
489  // Ionization
490  //
491  //opt.SetSubCutoff(true);
492 
493 
494  // Deexcitation
495  //
496  G4VAtomDeexcitation* deexcitation = new G4UAtomicDeexcitation();
498  deexcitation->SetFluo(true);
499 }
500 
501 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......