Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmLowEPPhysics.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: G4EmLowEPPhysics.cc 73153 2013-08-20 07:17:42Z gcosmo $
27 
28 #include "G4EmLowEPPhysics.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"
39 #include "G4LowEPComptonModel.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 "G4UrbanMscModel.hh"
87 #include "G4WentzelVIModel.hh"
89 #include "G4CoulombScattering.hh"
90 
91 // interfaces
92 #include "G4LossTableManager.hh"
93 #include "G4EmProcessOptions.hh"
94 #include "G4UAtomicDeexcitation.hh"
95 
96 // particles
97 
98 #include "G4Gamma.hh"
99 #include "G4Electron.hh"
100 #include "G4Positron.hh"
101 #include "G4MuonPlus.hh"
102 #include "G4MuonMinus.hh"
103 #include "G4PionPlus.hh"
104 #include "G4PionMinus.hh"
105 #include "G4KaonPlus.hh"
106 #include "G4KaonMinus.hh"
107 #include "G4Proton.hh"
108 #include "G4AntiProton.hh"
109 #include "G4Deuteron.hh"
110 #include "G4Triton.hh"
111 #include "G4He3.hh"
112 #include "G4Alpha.hh"
113 #include "G4GenericIon.hh"
114 
115 //
116 #include "G4PhysicsListHelper.hh"
117 #include "G4BuilderType.hh"
118 
119 // factory
121 //
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125 
127  : G4VPhysicsConstructor("G4EmLowEPPhysics"), verbose(ver)
128 {
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134 
136  : G4VPhysicsConstructor("G4EmLowEPPhysics"), verbose(ver)
137 {
140 }
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143 
145 {}
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
150 {
151 // gamma
152  G4Gamma::Gamma();
153 
154 // leptons
159 
160 // mesons
165 
166 // baryons
169 
170 // ions
173  G4He3::He3();
174  G4Alpha::Alpha();
176 }
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179 
181 {
183 
184  // muon & hadron bremsstrahlung and pair production
193 
194  // muon & hadron multiple scattering
196  mumsc->AddEmModel(0, new G4WentzelVIModel());
201 
202  // Add Livermore EM Processes
203  aParticleIterator->reset();
204 
205  while( (*aParticleIterator)() ){
206 
207  G4ParticleDefinition* particle = aParticleIterator->value();
208  G4String particleName = particle->GetParticleName();
209 
210  if(verbose > 1)
211  G4cout << "### " << GetPhysicsName() << " instantiates for "
212  << particleName << G4endl;
213 
214  //Applicability range for Livermore models
215  //for higher energies, the Standard models are used
216  G4double LivermoreHighEnergyLimit = GeV;
217 
218  if (particleName == "gamma") {
219 
220  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
221  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
223  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
224  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
225  ph->RegisterProcess(thePhotoElectricEffect, particle);
226 
227  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
228  G4LowEPComptonModel* theLowEPComptonModel =
229  new G4LowEPComptonModel();
230  theLowEPComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
231  theComptonScattering->AddEmModel(0, theLowEPComptonModel);
232  ph->RegisterProcess(theComptonScattering, particle);
233 
234  G4GammaConversion* theGammaConversion = new G4GammaConversion();
235  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
237  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
238  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
239  ph->RegisterProcess(theGammaConversion, particle);
240 
241  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
242  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
243  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
244  theRayleigh->AddEmModel(0, theRayleighModel);
245  ph->RegisterProcess(theRayleigh, particle);
246 
247  } else if (particleName == "e-") {
248 
250  //msc->AddEmModel(0, new G4UrbanMscModel());
251  //msc->AddEmModel(0, new G4GoudsmitSaundersonMscModel());
253  ph->RegisterProcess(msc, particle);
254 
255  // Ionisation
256  G4eIonisation* eIoni = new G4eIonisation();
257  G4LivermoreIonisationModel* theIoniLivermore = new
259  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
260  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
261  eIoni->SetStepFunction(0.2, 100*um); //
262  ph->RegisterProcess(eIoni, particle);
263 
264  // Bremsstrahlung
265  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
266  G4LivermoreBremsstrahlungModel* theBremLivermore = new
268  theBremLivermore->SetHighEnergyLimit(25*MeV);
269  theBremLivermore->SetAngularDistribution(new G4Generator2BS());
270  eBrem->AddEmModel(0, theBremLivermore);
271  ph->RegisterProcess(eBrem, particle);
272 
273  } else if (particleName == "e+") {
274 
275  // Identical to G4EmStandardPhysics_option3
276 
278  // msc->AddEmModel(0, new G4UrbanMscModel());
279  //msc->AddEmModel(0, new G4GoudsmitSaundersonMscModel());
281  G4eIonisation* eIoni = new G4eIonisation();
282  eIoni->SetStepFunction(0.2, 100*um);
283 
284  ph->RegisterProcess(msc, particle);
285  ph->RegisterProcess(eIoni, particle);
286  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
287  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
288 
289  } else if (particleName == "mu+" ||
290  particleName == "mu-" ) {
291 
292  G4MuIonisation* muIoni = new G4MuIonisation();
293  muIoni->SetStepFunction(0.2, 50*um);
294 
295  ph->RegisterProcess(mumsc, particle);
296  ph->RegisterProcess(muIoni, particle);
297  ph->RegisterProcess(mub, particle);
298  ph->RegisterProcess(mup, particle);
299  ph->RegisterProcess(new G4CoulombScattering(), particle);
300 
301  } else if (particleName == "alpha" ||
302  particleName == "He3" ) {
303 
304  // Identical to G4EmStandardPhysics_option3
305 
306  G4ionIonisation* ionIoni = new G4ionIonisation();
307  ionIoni->SetStepFunction(0.1, 10*um);
308 
309  ph->RegisterProcess(hmsc, particle);
310  ph->RegisterProcess(ionIoni, particle);
311  ph->RegisterProcess(new G4NuclearStopping(), particle);
312 
313  } else if (particleName == "GenericIon") {
314 
315  // Identical to G4EmStandardPhysics_option3
316 
317  G4ionIonisation* ionIoni = new G4ionIonisation();
318  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
319  ionIoni->SetStepFunction(0.1, 1*um);
320 
321  ph->RegisterProcess(hmsc, particle);
322  ph->RegisterProcess(ionIoni, particle);
323  ph->RegisterProcess(new G4NuclearStopping(), particle);
324 
325  } else if (particleName == "pi+" ||
326  particleName == "pi-" ) {
327 
328  G4hIonisation* hIoni = new G4hIonisation();
329  hIoni->SetStepFunction(0.2, 50*um);
330 
331  ph->RegisterProcess(pimsc, particle);
332  ph->RegisterProcess(hIoni, particle);
333  ph->RegisterProcess(pib, particle);
334  ph->RegisterProcess(pip, particle);
335 
336  } else if (particleName == "kaon+" ||
337  particleName == "kaon-" ) {
338 
339  G4hIonisation* hIoni = new G4hIonisation();
340  hIoni->SetStepFunction(0.2, 50*um);
341 
342  ph->RegisterProcess(kmsc, particle);
343  ph->RegisterProcess(hIoni, particle);
344  ph->RegisterProcess(kb, particle);
345  ph->RegisterProcess(kp, particle);
346 
347  } else if (particleName == "proton" ||
348  particleName == "anti_proton") {
349 
350  G4hIonisation* hIoni = new G4hIonisation();
351  hIoni->SetStepFunction(0.2, 50*um);
352 
353  ph->RegisterProcess(pmsc, particle);
354  ph->RegisterProcess(hIoni, particle);
355  ph->RegisterProcess(pb, particle);
356  ph->RegisterProcess(pp, particle);
357 
358  } else if (particleName == "B+" ||
359  particleName == "B-" ||
360  particleName == "D+" ||
361  particleName == "D-" ||
362  particleName == "Ds+" ||
363  particleName == "Ds-" ||
364  particleName == "anti_He3" ||
365  particleName == "anti_alpha" ||
366  particleName == "anti_deuteron" ||
367  particleName == "anti_lambda_c+" ||
368  particleName == "anti_omega-" ||
369  particleName == "anti_sigma_c+" ||
370  particleName == "anti_sigma_c++" ||
371  particleName == "anti_sigma+" ||
372  particleName == "anti_sigma-" ||
373  particleName == "anti_triton" ||
374  particleName == "anti_xi_c+" ||
375  particleName == "anti_xi-" ||
376  particleName == "deuteron" ||
377  particleName == "lambda_c+" ||
378  particleName == "omega-" ||
379  particleName == "sigma_c+" ||
380  particleName == "sigma_c++" ||
381  particleName == "sigma+" ||
382  particleName == "sigma-" ||
383  particleName == "tau+" ||
384  particleName == "tau-" ||
385  particleName == "triton" ||
386  particleName == "xi_c+" ||
387  particleName == "xi-" ) {
388 
389  // Identical to G4EmStandardPhysics_option3
390 
391  ph->RegisterProcess(hmsc, particle);
392  ph->RegisterProcess(new G4hIonisation(), particle);
393 
394  }
395  }
396 
397  // Em options
398  //
399  G4EmProcessOptions opt;
400  opt.SetVerbose(verbose);
401 
402  // Multiple Coulomb scattering
403  //
404  opt.SetPolarAngleLimit(CLHEP::pi);
405 
406  // Physics tables
407  //
408 
409  opt.SetMinEnergy(100*eV);
410  opt.SetMaxEnergy(10*TeV);
411  opt.SetDEDXBinning(220);
412  opt.SetLambdaBinning(220);
413 
414  // Ionization
415  //
416  //opt.SetSubCutoff(true);
417 
418  // Deexcitation
419  //
422  de->SetFluo(true);
423 }
424 
425 //....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()
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
virtual void ConstructParticle()
void SetMinEnergy(G4double val)
void SetStepFunction(G4double v1, G4double v2)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
void SetDEDXBinning(G4int val)
G4GLOB_DLL std::ostream G4cout
void SetLambdaBinning(G4int val)
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmLowEPPhysics)
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
virtual void ConstructProcess()
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
virtual ~G4EmLowEPPhysics()
void SetMaxEnergy(G4double val)
static G4Positron * Positron()
Definition: G4Positron.cc:94
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
G4EmLowEPPhysics(G4int ver=1)
double G4double
Definition: G4Types.hh:76
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)