Geant4  10.02.p01
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 92821 2015-09-17 15:23:49Z 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 #include "G4EmModelActivator.hh"
125 
126 // factory
128 //
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132 
134  : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
135 {
137  param->SetDefaults();
138  param->SetVerbose(verbose);
139  param->SetMinEnergy(100*eV);
140  param->SetMaxEnergy(10*TeV);
141  param->SetLowestElectronEnergy(100*eV);
142  param->SetNumberOfBinsPerDecade(20);
143  param->SetMscRangeFactor(0.02);
145  param->SetFluo(true);
146  param->SetPIXEElectronCrossSectionModel("Penelope");
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151 
153  : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
154 {
156  param->SetDefaults();
157  param->SetVerbose(verbose);
158  param->SetMinEnergy(100*eV);
159  param->SetMaxEnergy(10*TeV);
160  param->SetLowestElectronEnergy(100*eV);
161  param->SetNumberOfBinsPerDecade(20);
162  param->SetMscRangeFactor(0.02);
164  param->SetFluo(true);
165  param->SetPIXEElectronCrossSectionModel("Penelope");
167 }
168 
169 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
170 
172 {}
173 
174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
175 
177 {
178  // gamma
179  G4Gamma::Gamma();
180 
181  // leptons
186 
187  // mesons
192 
193  // baryons
196 
197  // ions
200  G4He3::He3();
201  G4Alpha::Alpha();
203 
204  // dna
205  G4EmModelActivator mact;
206  mact.ConstructParticle();
207 }
208 
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210 
212 {
213  if(verbose > 1) {
214  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
215  }
217 
218  // muon & hadron bremsstrahlung and pair production
227 
228  // muon & hadron multiple scattering
230  mumsc->AddEmModel(0, new G4WentzelVIModel());
231  //G4MuMultipleScattering* pimsc = new G4MuMultipleScattering();
232  //pimsc->AddEmModel(0, new G4WentzelVIModel());
233  //G4MuMultipleScattering* kmsc = new G4MuMultipleScattering();
234  //kmsc->AddEmModel(0, new G4WentzelVIModel());
235  //G4MuMultipleScattering* pmsc = new G4MuMultipleScattering();
236  //pmsc->AddEmModel(0, new G4WentzelVIModel());
237  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
238 
239  // high energy limit for e+- scattering models
240  G4double highEnergyLimit = 100*MeV;
241 
242  // nuclear stopping
243  G4NuclearStopping* pnuc = new G4NuclearStopping();
244 
245  // Add Penelope EM Processes
246  aParticleIterator->reset();
247 
248  while( (*aParticleIterator)() ){
249 
250  G4ParticleDefinition* particle = aParticleIterator->value();
251  G4String particleName = particle->GetParticleName();
252 
253  //Applicability range for Penelope models
254  //for higher energies, the Standard models are used
255  G4double PenelopeHighEnergyLimit = 1.0*GeV;
256 
257  if (particleName == "gamma") {
258 
259  //Photo-electric effect
260  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
261  G4PenelopePhotoElectricModel* thePEPenelopeModel = new
263  thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
264  thePhotoElectricEffect->SetEmModel(thePEPenelopeModel, 1);
265  ph->RegisterProcess(thePhotoElectricEffect, particle);
266 
267  //Compton scattering
268  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
269  G4PenelopeComptonModel* theComptonPenelopeModel =
271  theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
272  theComptonScattering->SetEmModel(theComptonPenelopeModel, 1);
273  ph->RegisterProcess(theComptonScattering, particle);
274 
275  //Gamma conversion
276  G4GammaConversion* theGammaConversion = new G4GammaConversion();
277  G4PenelopeGammaConversionModel* theGCPenelopeModel =
279  theGammaConversion->SetEmModel(theGCPenelopeModel,1);
280  ph->RegisterProcess(theGammaConversion, particle);
281 
282  //Rayleigh scattering
283  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
284  G4PenelopeRayleighModel* theRayleighPenelopeModel =
286  //theRayleighPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
287  theRayleigh->SetEmModel(theRayleighPenelopeModel, 1);
288  ph->RegisterProcess(theRayleigh, particle);
289 
290  } else if (particleName == "e-") {
291 
292  // multiple scattering
294  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
295  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
296  msc1->SetHighEnergyLimit(highEnergyLimit);
297  msc2->SetLowEnergyLimit(highEnergyLimit);
298  msc->AddEmModel(0, msc1);
299  msc->AddEmModel(0, msc2);
300 
303  ss->SetEmModel(ssm, 1);
304  ss->SetMinKinEnergy(highEnergyLimit);
305  ssm->SetLowEnergyLimit(highEnergyLimit);
306  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
307 
308  //Ionisation
309  G4eIonisation* eIoni = new G4eIonisation();
310  G4PenelopeIonisationModel* theIoniPenelope =
312  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
313  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
314  eIoni->SetStepFunction(0.2, 100*um); //
315 
316  //Bremsstrahlung
317  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
318  G4PenelopeBremsstrahlungModel* theBremPenelope = new
320  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
321  eBrem->AddEmModel(0,theBremPenelope);
322 
323  // register processes
324  ph->RegisterProcess(msc, particle);
325  ph->RegisterProcess(eIoni, particle);
326  ph->RegisterProcess(eBrem, particle);
327  ph->RegisterProcess(ss, particle);
328 
329  } else if (particleName == "e+") {
330 
331  // multiple scattering
333  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
334  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
335  msc1->SetHighEnergyLimit(highEnergyLimit);
336  msc2->SetLowEnergyLimit(highEnergyLimit);
337  msc->AddEmModel(0, msc1);
338  msc->AddEmModel(0, msc2);
339 
342  ss->SetEmModel(ssm, 1);
343  ss->SetMinKinEnergy(highEnergyLimit);
344  ssm->SetLowEnergyLimit(highEnergyLimit);
345  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
346 
347  //Ionisation
348  G4eIonisation* eIoni = new G4eIonisation();
349  G4PenelopeIonisationModel* theIoniPenelope =
351  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
352  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
353  eIoni->SetStepFunction(0.2, 100*um); //
354 
355  //Bremsstrahlung
356  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
357  G4PenelopeBremsstrahlungModel* theBremPenelope = new
359  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
360  eBrem->AddEmModel(0,theBremPenelope);
361 
362  //Annihilation
364  G4PenelopeAnnihilationModel* theAnnPenelope = new
366  theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
367  eAnni->AddEmModel(0,theAnnPenelope);
368 
369  // register processes
370  ph->RegisterProcess(msc, particle);
371  ph->RegisterProcess(eIoni, particle);
372  ph->RegisterProcess(eBrem, particle);
373  ph->RegisterProcess(eAnni, particle);
374  ph->RegisterProcess(ss, particle);
375 
376  } else if (particleName == "mu+" ||
377  particleName == "mu-" ) {
378 
379  G4MuIonisation* muIoni = new G4MuIonisation();
380  muIoni->SetStepFunction(0.2, 50*um);
381 
382  ph->RegisterProcess(mumsc, particle);
383  ph->RegisterProcess(muIoni, particle);
384  ph->RegisterProcess(mub, particle);
385  ph->RegisterProcess(mup, particle);
386  ph->RegisterProcess(new G4CoulombScattering(), particle);
387 
388  } else if (particleName == "alpha" ||
389  particleName == "He3" ) {
390 
392  G4ionIonisation* ionIoni = new G4ionIonisation();
393  ionIoni->SetStepFunction(0.1, 10*um);
394 
395  ph->RegisterProcess(msc, particle);
396  ph->RegisterProcess(ionIoni, particle);
397  ph->RegisterProcess(pnuc, particle);
398 
399  } else if (particleName == "GenericIon") {
400 
401  G4ionIonisation* ionIoni = new G4ionIonisation();
402  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
403  ionIoni->SetStepFunction(0.1, 1*um);
404 
405  ph->RegisterProcess(hmsc, particle);
406  ph->RegisterProcess(ionIoni, particle);
407  ph->RegisterProcess(pnuc, particle);
408 
409  } else if (particleName == "pi+" ||
410  particleName == "pi-" ) {
411 
413  G4hIonisation* hIoni = new G4hIonisation();
414  hIoni->SetStepFunction(0.2, 50*um);
415 
416  ph->RegisterProcess(pimsc, particle);
417  ph->RegisterProcess(hIoni, particle);
418  ph->RegisterProcess(pib, particle);
419  ph->RegisterProcess(pip, particle);
420 
421  } else if (particleName == "kaon+" ||
422  particleName == "kaon-" ) {
423 
425  G4hIonisation* hIoni = new G4hIonisation();
426  hIoni->SetStepFunction(0.2, 50*um);
427 
428  ph->RegisterProcess(kmsc, particle);
429  ph->RegisterProcess(hIoni, particle);
430  ph->RegisterProcess(kb, particle);
431  ph->RegisterProcess(kp, particle);
432 
433  } else if (particleName == "proton" ||
434  particleName == "anti_proton") {
435 
437  G4hIonisation* hIoni = new G4hIonisation();
438  hIoni->SetStepFunction(0.2, 50*um);
439 
440  ph->RegisterProcess(pmsc, particle);
441  ph->RegisterProcess(hIoni, particle);
442  ph->RegisterProcess(pb, particle);
443  ph->RegisterProcess(pp, particle);
444  ph->RegisterProcess(pnuc, particle);
445 
446  } else if (particleName == "B+" ||
447  particleName == "B-" ||
448  particleName == "D+" ||
449  particleName == "D-" ||
450  particleName == "Ds+" ||
451  particleName == "Ds-" ||
452  particleName == "anti_He3" ||
453  particleName == "anti_alpha" ||
454  particleName == "anti_deuteron" ||
455  particleName == "anti_lambda_c+" ||
456  particleName == "anti_omega-" ||
457  particleName == "anti_sigma_c+" ||
458  particleName == "anti_sigma_c++" ||
459  particleName == "anti_sigma+" ||
460  particleName == "anti_sigma-" ||
461  particleName == "anti_triton" ||
462  particleName == "anti_xi_c+" ||
463  particleName == "anti_xi-" ||
464  particleName == "deuteron" ||
465  particleName == "lambda_c+" ||
466  particleName == "omega-" ||
467  particleName == "sigma_c+" ||
468  particleName == "sigma_c++" ||
469  particleName == "sigma+" ||
470  particleName == "sigma-" ||
471  particleName == "tau+" ||
472  particleName == "tau-" ||
473  particleName == "triton" ||
474  particleName == "xi_c+" ||
475  particleName == "xi-" ) {
476 
477  ph->RegisterProcess(hmsc, particle);
478  ph->RegisterProcess(new G4hIonisation(), particle);
479  }
480  }
481 
482  // Nuclear stopping
483  pnuc->SetMaxKinEnergy(MeV);
484 
485  // Deexcitation
486  //
487  G4VAtomDeexcitation* deexcitation = new G4UAtomicDeexcitation();
489 
490  G4EmModelActivator mact;
491  mact.ConstructProcess();
492 }
493 
494 //....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:211
void SetVerbose(G4int val)
static G4LossTableManager * Instance()
G4EmPenelopePhysics(G4int ver=1)
void SetLowestElectronEnergy(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
void SetStepFunction(G4double v1, G4double v2)
virtual void ConstructProcess()
void SetPIXEElectronCrossSectionModel(const G4String &)
int G4int
Definition: G4Types.hh:78
void SetMaxEnergy(G4double val)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
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
void SetMscRangeFactor(G4double val)
#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:214
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:746
static const double eV
Definition: G4SIunits.hh:212
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 const double um
Definition: G4SIunits.hh:112
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:215
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static G4He3 * He3()
Definition: G4He3.cc:94
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetFluo(G4bool val)