Geant4.10
 All Data Structures 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: G4EmPenelopePhysics.cc 79157 2014-02-19 15:35:01Z gcosmo $
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 "G4UrbanMscModel.hh"
93 #include "G4WentzelVIModel.hh"
94 #include "G4CoulombScattering.hh"
96 //
97 
98 #include "G4LossTableManager.hh"
99 #include "G4VAtomDeexcitation.hh"
100 #include "G4UAtomicDeexcitation.hh"
101 #include "G4EmProcessOptions.hh"
102 
103 // particles
104 
105 #include "G4Gamma.hh"
106 #include "G4Electron.hh"
107 #include "G4Positron.hh"
108 #include "G4MuonPlus.hh"
109 #include "G4MuonMinus.hh"
110 #include "G4PionPlus.hh"
111 #include "G4PionMinus.hh"
112 #include "G4KaonPlus.hh"
113 #include "G4KaonMinus.hh"
114 #include "G4Proton.hh"
115 #include "G4AntiProton.hh"
116 #include "G4Deuteron.hh"
117 #include "G4Triton.hh"
118 #include "G4He3.hh"
119 #include "G4Alpha.hh"
120 #include "G4GenericIon.hh"
121 
122 //
123 #include "G4PhysicsListHelper.hh"
124 #include "G4BuilderType.hh"
125 
126 // factory
128 //
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132 
134  : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
135 {
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141 
143  : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
144 {
147 }
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150 
152 {}
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
157 {
158 // gamma
159  G4Gamma::Gamma();
160 
161 // leptons
166 
167 // mesons
172 
173 // baryons
176 
177 // ions
180  G4He3::He3();
181  G4Alpha::Alpha();
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186 
188 {
189  if(verbose > 1) {
190  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
191  }
193 
194  // muon & hadron bremsstrahlung and pair production
203 
204  // muon & hadron multiple scattering
206  mumsc->AddEmModel(0, new G4WentzelVIModel());
208  //pimsc->AddEmModel(0, new G4WentzelVIModel());
210  //kmsc->AddEmModel(0, new G4WentzelVIModel());
212  //pmsc->AddEmModel(0, new G4WentzelVIModel());
213  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
214 
215  // high energy limit for e+- scattering models
216  G4double highEnergyLimit = 100*MeV;
217 
218  // nuclear stopping
219  G4NuclearStopping* ionnuc = new G4NuclearStopping();
220  G4NuclearStopping* pnuc = new G4NuclearStopping();
221 
222  // Add Penelope EM Processes
223  aParticleIterator->reset();
224 
225  while( (*aParticleIterator)() ){
226 
227  G4ParticleDefinition* particle = aParticleIterator->value();
228  G4String particleName = particle->GetParticleName();
229 
230  //Applicability range for Penelope models
231  //for higher energies, the Standard models are used
232  G4double PenelopeHighEnergyLimit = 1.0*GeV;
233 
234  if (particleName == "gamma") {
235 
236  //Photo-electric effect
237  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
238  G4PenelopePhotoElectricModel* thePEPenelopeModel = new
240  thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
241  thePhotoElectricEffect->SetEmModel(thePEPenelopeModel, 1);
242  ph->RegisterProcess(thePhotoElectricEffect, particle);
243 
244  //Compton scattering
245  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
246  G4PenelopeComptonModel* theComptonPenelopeModel =
248  theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
249  theComptonScattering->SetEmModel(theComptonPenelopeModel, 1);
250  ph->RegisterProcess(theComptonScattering, particle);
251 
252  //Gamma conversion
253  G4GammaConversion* theGammaConversion = new G4GammaConversion();
254  G4PenelopeGammaConversionModel* theGCPenelopeModel =
256  theGammaConversion->SetEmModel(theGCPenelopeModel,1);
257  ph->RegisterProcess(theGammaConversion, particle);
258 
259  //Rayleigh scattering
260  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
261  G4PenelopeRayleighModel* theRayleighPenelopeModel =
263  //theRayleighPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
264  theRayleigh->SetEmModel(theRayleighPenelopeModel, 1);
265  ph->RegisterProcess(theRayleigh, particle);
266 
267  } else if (particleName == "e-") {
268 
269  // multiple scattering
272  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
273  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
274  msc1->SetHighEnergyLimit(highEnergyLimit);
275  msc2->SetLowEnergyLimit(highEnergyLimit);
276  msc->SetRangeFactor(0.01);
277  msc->AddEmModel(0, msc1);
278  msc->AddEmModel(0, msc2);
279 
282  ss->SetEmModel(ssm, 1);
283  ss->SetMinKinEnergy(highEnergyLimit);
284  ssm->SetLowEnergyLimit(highEnergyLimit);
285  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
286 
287  //Ionisation
288  G4eIonisation* eIoni = new G4eIonisation();
289  G4PenelopeIonisationModel* theIoniPenelope =
291  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
292  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
293  eIoni->SetStepFunction(0.2, 100*um); //
294 
295  //Bremsstrahlung
296  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
297  G4PenelopeBremsstrahlungModel* theBremPenelope = new
299  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
300  eBrem->AddEmModel(0,theBremPenelope);
301 
302  // register processes
303  ph->RegisterProcess(msc, particle);
304  ph->RegisterProcess(eIoni, particle);
305  ph->RegisterProcess(eBrem, particle);
306  ph->RegisterProcess(ss, particle);
307 
308  } else if (particleName == "e+") {
309 
310  // multiple scattering
313  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
314  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
315  msc1->SetHighEnergyLimit(highEnergyLimit);
316  msc2->SetLowEnergyLimit(highEnergyLimit);
317  msc->SetRangeFactor(0.01);
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  G4PenelopeIonisationModel* theIoniPenelope =
332  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
333  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
334  eIoni->SetStepFunction(0.2, 100*um); //
335 
336  //Bremsstrahlung
337  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
338  G4PenelopeBremsstrahlungModel* theBremPenelope = new
340  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
341  eBrem->AddEmModel(0,theBremPenelope);
342 
343  //Annihilation
345  G4PenelopeAnnihilationModel* theAnnPenelope = new
347  theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
348  eAnni->AddEmModel(0,theAnnPenelope);
349 
350  // register processes
351  ph->RegisterProcess(msc, particle);
352  ph->RegisterProcess(eIoni, particle);
353  ph->RegisterProcess(eBrem, particle);
354  ph->RegisterProcess(eAnni, particle);
355  ph->RegisterProcess(ss, particle);
356 
357  } else if (particleName == "mu+" ||
358  particleName == "mu-" ) {
359 
360  // Identical to G4EmStandardPhysics_option3
361 
362  G4MuIonisation* muIoni = new G4MuIonisation();
363  muIoni->SetStepFunction(0.2, 50*um);
364 
365  ph->RegisterProcess(mumsc, particle);
366  ph->RegisterProcess(muIoni, particle);
367  ph->RegisterProcess(mub, particle);
368  ph->RegisterProcess(mup, particle);
369  ph->RegisterProcess(new G4CoulombScattering(), particle);
370 
371  } else if (particleName == "alpha" ||
372  particleName == "He3" ) {
373 
374  // Identical to G4EmStandardPhysics_option3
375 
377  G4ionIonisation* ionIoni = new G4ionIonisation();
378  ionIoni->SetStepFunction(0.1, 10*um);
379 
380  ph->RegisterProcess(msc, particle);
381  ph->RegisterProcess(ionIoni, particle);
382  ph->RegisterProcess(ionnuc, particle);
383 
384  } else if (particleName == "GenericIon") {
385 
386  // Identical to G4EmStandardPhysics_option3
387 
388  G4ionIonisation* ionIoni = new G4ionIonisation();
389  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
390  ionIoni->SetStepFunction(0.1, 1*um);
391 
392  ph->RegisterProcess(hmsc, particle);
393  ph->RegisterProcess(ionIoni, particle);
394  ph->RegisterProcess(ionnuc, particle);
395 
396  } else if (particleName == "pi+" ||
397  particleName == "pi-" ) {
398 
399  G4hIonisation* hIoni = new G4hIonisation();
400  hIoni->SetStepFunction(0.2, 50*um);
401 
402  ph->RegisterProcess(pimsc, particle);
403  ph->RegisterProcess(hIoni, particle);
404  ph->RegisterProcess(pib, particle);
405  ph->RegisterProcess(pip, particle);
406 
407  } else if (particleName == "kaon+" ||
408  particleName == "kaon-" ) {
409 
410  G4hIonisation* hIoni = new G4hIonisation();
411  hIoni->SetStepFunction(0.2, 50*um);
412 
413  ph->RegisterProcess(kmsc, particle);
414  ph->RegisterProcess(hIoni, particle);
415  ph->RegisterProcess(kb, particle);
416  ph->RegisterProcess(kp, particle);
417 
418  } else if (particleName == "proton" ||
419  particleName == "anti_proton") {
420 
421  G4hIonisation* hIoni = new G4hIonisation();
422  hIoni->SetStepFunction(0.2, 50*um);
423 
424  ph->RegisterProcess(pmsc, particle);
425  ph->RegisterProcess(hIoni, particle);
426  ph->RegisterProcess(pb, particle);
427  ph->RegisterProcess(pp, particle);
428  ph->RegisterProcess(pnuc, particle);
429 
430  } else if (particleName == "B+" ||
431  particleName == "B-" ||
432  particleName == "D+" ||
433  particleName == "D-" ||
434  particleName == "Ds+" ||
435  particleName == "Ds-" ||
436  particleName == "anti_He3" ||
437  particleName == "anti_alpha" ||
438  particleName == "anti_deuteron" ||
439  particleName == "anti_lambda_c+" ||
440  particleName == "anti_omega-" ||
441  particleName == "anti_sigma_c+" ||
442  particleName == "anti_sigma_c++" ||
443  particleName == "anti_sigma+" ||
444  particleName == "anti_sigma-" ||
445  particleName == "anti_triton" ||
446  particleName == "anti_xi_c+" ||
447  particleName == "anti_xi-" ||
448  particleName == "deuteron" ||
449  particleName == "lambda_c+" ||
450  particleName == "omega-" ||
451  particleName == "sigma_c+" ||
452  particleName == "sigma_c++" ||
453  particleName == "sigma+" ||
454  particleName == "sigma-" ||
455  particleName == "tau+" ||
456  particleName == "tau-" ||
457  particleName == "triton" ||
458  particleName == "xi_c+" ||
459  particleName == "xi-" ) {
460 
461  // Identical to G4EmStandardPhysics_option3
462 
463  ph->RegisterProcess(hmsc, particle);
464  ph->RegisterProcess(new G4hIonisation(), particle);
465  ph->RegisterProcess(pnuc, particle);
466  }
467  }
468 
469  // Em options
470  //
471  G4EmProcessOptions opt;
472  opt.SetVerbose(verbose);
473 
474  // Multiple Coulomb scattering
475  //
476  //opt.SetMscStepLimitation(fUseDistanceToBoundary);
477  //opt.SetMscRangeFactor(0.02);
478 
479  // Physics tables
480  //
481 
482  opt.SetMinEnergy(100*eV);
483  opt.SetMaxEnergy(10*TeV);
484  opt.SetDEDXBinning(220);
485  opt.SetLambdaBinning(220);
486 
487  // Nuclear stopping
488  pnuc->SetMaxKinEnergy(MeV);
489 
490  //opt.SetSplineFlag(true);
491  opt.SetPolarAngleLimit(CLHEP::pi);
492 
493  // Ionization
494  //
495  //opt.SetSubCutoff(true);
496 
497 
498  // Deexcitation
499  //
500  G4VAtomDeexcitation* deexcitation = new G4UAtomicDeexcitation();
502  deexcitation->SetFluo(true);
503 }
504 
505 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
static G4LossTableManager * Instance()
G4EmPenelopePhysics(G4int ver=1)
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
void SetMinEnergy(G4double val)
void SetStepFunction(G4double v1, G4double v2)
virtual void ConstructProcess()
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
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 G4Positron * Positron()
Definition: G4Positron.cc:94
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmPenelopePhysics)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
virtual void ConstructParticle()
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)