Geant4  10.01.p02
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 86233 2014-11-07 17:21:03Z 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 
102 // particles
103 
104 #include "G4Gamma.hh"
105 #include "G4Electron.hh"
106 #include "G4Positron.hh"
107 #include "G4MuonPlus.hh"
108 #include "G4MuonMinus.hh"
109 #include "G4PionPlus.hh"
110 #include "G4PionMinus.hh"
111 #include "G4KaonPlus.hh"
112 #include "G4KaonMinus.hh"
113 #include "G4Proton.hh"
114 #include "G4AntiProton.hh"
115 #include "G4Deuteron.hh"
116 #include "G4Triton.hh"
117 #include "G4He3.hh"
118 #include "G4Alpha.hh"
119 #include "G4GenericIon.hh"
120 
121 //
122 #include "G4PhysicsListHelper.hh"
123 #include "G4BuilderType.hh"
124 
125 // factory
127 //
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131 
133  : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
134 {
136  param->SetVerbose(verbose);
137  param->SetMinEnergy(100*eV);
138  param->SetMaxEnergy(10*TeV);
139  param->SetNumberOfBinsPerDecade(20);
141 }
142 
143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144 
146  : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
147 {
149  param->SetVerbose(verbose);
150  param->SetMinEnergy(100*eV);
151  param->SetMaxEnergy(10*TeV);
152  param->SetNumberOfBinsPerDecade(20);
154 }
155 
156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
157 
159 {}
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 
164 {
165 // gamma
166  G4Gamma::Gamma();
167 
168 // leptons
173 
174 // mesons
179 
180 // baryons
183 
184 // ions
187  G4He3::He3();
188  G4Alpha::Alpha();
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193 
195 {
196  if(verbose > 1) {
197  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
198  }
200 
201  // muon & hadron bremsstrahlung and pair production
210 
211  // muon & hadron multiple scattering
213  mumsc->AddEmModel(0, new G4WentzelVIModel());
214  //G4MuMultipleScattering* pimsc = new G4MuMultipleScattering();
215  //pimsc->AddEmModel(0, new G4WentzelVIModel());
216  //G4MuMultipleScattering* kmsc = new G4MuMultipleScattering();
217  //kmsc->AddEmModel(0, new G4WentzelVIModel());
218  //G4MuMultipleScattering* pmsc = new G4MuMultipleScattering();
219  //pmsc->AddEmModel(0, new G4WentzelVIModel());
220  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
221 
222  // high energy limit for e+- scattering models
223  G4double highEnergyLimit = 100*MeV;
224 
225  // nuclear stopping
226  G4NuclearStopping* pnuc = new G4NuclearStopping();
227 
228  // Add Penelope EM Processes
229  aParticleIterator->reset();
230 
231  while( (*aParticleIterator)() ){
232 
233  G4ParticleDefinition* particle = aParticleIterator->value();
234  G4String particleName = particle->GetParticleName();
235 
236  //Applicability range for Penelope models
237  //for higher energies, the Standard models are used
238  G4double PenelopeHighEnergyLimit = 1.0*GeV;
239 
240  if (particleName == "gamma") {
241 
242  //Photo-electric effect
243  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
244  G4PenelopePhotoElectricModel* thePEPenelopeModel = new
246  thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
247  thePhotoElectricEffect->SetEmModel(thePEPenelopeModel, 1);
248  ph->RegisterProcess(thePhotoElectricEffect, particle);
249 
250  //Compton scattering
251  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
252  G4PenelopeComptonModel* theComptonPenelopeModel =
254  theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
255  theComptonScattering->SetEmModel(theComptonPenelopeModel, 1);
256  ph->RegisterProcess(theComptonScattering, particle);
257 
258  //Gamma conversion
259  G4GammaConversion* theGammaConversion = new G4GammaConversion();
260  G4PenelopeGammaConversionModel* theGCPenelopeModel =
262  theGammaConversion->SetEmModel(theGCPenelopeModel,1);
263  ph->RegisterProcess(theGammaConversion, particle);
264 
265  //Rayleigh scattering
266  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
267  G4PenelopeRayleighModel* theRayleighPenelopeModel =
269  //theRayleighPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
270  theRayleigh->SetEmModel(theRayleighPenelopeModel, 1);
271  ph->RegisterProcess(theRayleigh, particle);
272 
273  } else if (particleName == "e-") {
274 
275  // multiple scattering
278  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
279  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
280  msc1->SetHighEnergyLimit(highEnergyLimit);
281  msc2->SetLowEnergyLimit(highEnergyLimit);
282  msc->SetRangeFactor(0.01);
283  msc->AddEmModel(0, msc1);
284  msc->AddEmModel(0, msc2);
285 
288  ss->SetEmModel(ssm, 1);
289  ss->SetMinKinEnergy(highEnergyLimit);
290  ssm->SetLowEnergyLimit(highEnergyLimit);
291  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
292 
293  //Ionisation
294  G4eIonisation* eIoni = new G4eIonisation();
295  G4PenelopeIonisationModel* theIoniPenelope =
297  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
298  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
299  eIoni->SetStepFunction(0.2, 100*um); //
300 
301  //Bremsstrahlung
302  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
303  G4PenelopeBremsstrahlungModel* theBremPenelope = new
305  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
306  eBrem->AddEmModel(0,theBremPenelope);
307 
308  // register processes
309  ph->RegisterProcess(msc, particle);
310  ph->RegisterProcess(eIoni, particle);
311  ph->RegisterProcess(eBrem, particle);
312  ph->RegisterProcess(ss, particle);
313 
314  } else if (particleName == "e+") {
315 
316  // multiple scattering
319  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
320  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
321  msc1->SetHighEnergyLimit(highEnergyLimit);
322  msc2->SetLowEnergyLimit(highEnergyLimit);
323  msc->SetRangeFactor(0.01);
324  msc->AddEmModel(0, msc1);
325  msc->AddEmModel(0, msc2);
326 
329  ss->SetEmModel(ssm, 1);
330  ss->SetMinKinEnergy(highEnergyLimit);
331  ssm->SetLowEnergyLimit(highEnergyLimit);
332  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
333 
334  //Ionisation
335  G4eIonisation* eIoni = new G4eIonisation();
336  G4PenelopeIonisationModel* theIoniPenelope =
338  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
339  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
340  eIoni->SetStepFunction(0.2, 100*um); //
341 
342  //Bremsstrahlung
343  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
344  G4PenelopeBremsstrahlungModel* theBremPenelope = new
346  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
347  eBrem->AddEmModel(0,theBremPenelope);
348 
349  //Annihilation
351  G4PenelopeAnnihilationModel* theAnnPenelope = new
353  theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
354  eAnni->AddEmModel(0,theAnnPenelope);
355 
356  // register processes
357  ph->RegisterProcess(msc, particle);
358  ph->RegisterProcess(eIoni, particle);
359  ph->RegisterProcess(eBrem, particle);
360  ph->RegisterProcess(eAnni, particle);
361  ph->RegisterProcess(ss, particle);
362 
363  } else if (particleName == "mu+" ||
364  particleName == "mu-" ) {
365 
366  G4MuIonisation* muIoni = new G4MuIonisation();
367  muIoni->SetStepFunction(0.2, 50*um);
368 
369  ph->RegisterProcess(mumsc, particle);
370  ph->RegisterProcess(muIoni, particle);
371  ph->RegisterProcess(mub, particle);
372  ph->RegisterProcess(mup, particle);
373  ph->RegisterProcess(new G4CoulombScattering(), particle);
374 
375  } else if (particleName == "alpha" ||
376  particleName == "He3" ) {
377 
379  G4ionIonisation* ionIoni = new G4ionIonisation();
380  ionIoni->SetStepFunction(0.1, 10*um);
381 
382  ph->RegisterProcess(msc, particle);
383  ph->RegisterProcess(ionIoni, particle);
384  ph->RegisterProcess(pnuc, particle);
385 
386  } else if (particleName == "GenericIon") {
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(pnuc, particle);
395 
396  } else if (particleName == "pi+" ||
397  particleName == "pi-" ) {
398 
400  G4hIonisation* hIoni = new G4hIonisation();
401  hIoni->SetStepFunction(0.2, 50*um);
402 
403  ph->RegisterProcess(pimsc, particle);
404  ph->RegisterProcess(hIoni, particle);
405  ph->RegisterProcess(pib, particle);
406  ph->RegisterProcess(pip, particle);
407 
408  } else if (particleName == "kaon+" ||
409  particleName == "kaon-" ) {
410 
412  G4hIonisation* hIoni = new G4hIonisation();
413  hIoni->SetStepFunction(0.2, 50*um);
414 
415  ph->RegisterProcess(kmsc, particle);
416  ph->RegisterProcess(hIoni, particle);
417  ph->RegisterProcess(kb, particle);
418  ph->RegisterProcess(kp, particle);
419 
420  } else if (particleName == "proton" ||
421  particleName == "anti_proton") {
422 
424  G4hIonisation* hIoni = new G4hIonisation();
425  hIoni->SetStepFunction(0.2, 50*um);
426 
427  ph->RegisterProcess(pmsc, particle);
428  ph->RegisterProcess(hIoni, particle);
429  ph->RegisterProcess(pb, particle);
430  ph->RegisterProcess(pp, particle);
431  ph->RegisterProcess(pnuc, particle);
432 
433  } else if (particleName == "B+" ||
434  particleName == "B-" ||
435  particleName == "D+" ||
436  particleName == "D-" ||
437  particleName == "Ds+" ||
438  particleName == "Ds-" ||
439  particleName == "anti_He3" ||
440  particleName == "anti_alpha" ||
441  particleName == "anti_deuteron" ||
442  particleName == "anti_lambda_c+" ||
443  particleName == "anti_omega-" ||
444  particleName == "anti_sigma_c+" ||
445  particleName == "anti_sigma_c++" ||
446  particleName == "anti_sigma+" ||
447  particleName == "anti_sigma-" ||
448  particleName == "anti_triton" ||
449  particleName == "anti_xi_c+" ||
450  particleName == "anti_xi-" ||
451  particleName == "deuteron" ||
452  particleName == "lambda_c+" ||
453  particleName == "omega-" ||
454  particleName == "sigma_c+" ||
455  particleName == "sigma_c++" ||
456  particleName == "sigma+" ||
457  particleName == "sigma-" ||
458  particleName == "tau+" ||
459  particleName == "tau-" ||
460  particleName == "triton" ||
461  particleName == "xi_c+" ||
462  particleName == "xi-" ) {
463 
464  ph->RegisterProcess(hmsc, particle);
465  ph->RegisterProcess(new G4hIonisation(), particle);
466  }
467  }
468 
469  // Nuclear stopping
470  pnuc->SetMaxKinEnergy(MeV);
471 
472  // Deexcitation
473  //
474  G4VAtomDeexcitation* deexcitation = new G4UAtomicDeexcitation();
476  deexcitation->SetFluo(true);
477 }
478 
479 //....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 const double MeV
Definition: G4SIunits.hh:193
void SetVerbose(G4int val)
static G4LossTableManager * Instance()
G4EmPenelopePhysics(G4int ver=1)
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
void SetStepFunction(G4double v1, G4double v2)
virtual void ConstructProcess()
int G4int
Definition: G4Types.hh:78
void SetMaxEnergy(G4double val)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:707
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
#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: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 SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:728
static const double eV
Definition: G4SIunits.hh:194
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetMinEnergy(G4double val)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4EmParameters * Instance()
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 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:714
void SetRangeFactor(G4double val)
static G4He3 * He3()
Definition: G4He3.cc:94
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetStepLimitType(G4MscStepLimitType val)