Geant4_10
G4EmStandardPhysics_option2.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: G4EmStandardPhysics_option2.cc 75169 2013-10-29 09:21:54Z gcosmo $
27 //
28 //---------------------------------------------------------------------------
29 //
30 // ClassName: G4EmStandardPhysics_option2
31 //
32 // Author: V.Ivanchenko 09.11.2005
33 //
34 // Modified:
35 // 19.12.2005 V.Ivanchenko rename 71 -> 72
36 // 15.06.2006 V.Ivanchenko use this class as a constructor of fast EM physics
37 // 13.11.2006 V.Ivanchenko use G4hMultipleScattering
38 // 14.11.2006 V.Ivanchenko use sub-cutoff option for all particles
39 // 13.02.2007 V.Ivanchenko use default msc
40 // 15.05.2007 V.Ivanchenko rename to _option2
41 // 13.03.2008 V.Ivanchenko use G4eMultipleScattering
42 // 21.04.2008 V.Ivanchenko add long-lived D and B mesons; use spline
43 // 28.05.2008 V.Ivanchenko linLossLimit=0.01; added hBrem and hPairProd processes
44 //
45 //----------------------------------------------------------------------------
46 //
47 
49 
50 #include "G4SystemOfUnits.hh"
51 #include "G4ParticleDefinition.hh"
52 #include "G4LossTableManager.hh"
53 #include "G4EmProcessOptions.hh"
54 
55 #include "G4ComptonScattering.hh"
56 #include "G4GammaConversion.hh"
57 #include "G4PhotoElectricEffect.hh"
58 #include "G4PEEffectFluoModel.hh"
59 #include "G4KleinNishinaModel.hh"
60 #include "G4RayleighScattering.hh"
61 
62 #include "G4eMultipleScattering.hh"
63 #include "G4hMultipleScattering.hh"
65 #include "G4CoulombScattering.hh"
67 #include "G4UrbanMscModel.hh"
68 #include "G4WentzelVIModel.hh"
69 
70 #include "G4eIonisation.hh"
71 #include "G4eBremsstrahlung.hh"
73 #include "G4eplusAnnihilation.hh"
74 #include "G4Generator2BS.hh"
75 #include "G4SeltzerBergerModel.hh"
76 #include "G4UAtomicDeexcitation.hh"
77 
78 #include "G4MuIonisation.hh"
79 #include "G4MuBremsstrahlung.hh"
80 #include "G4MuPairProduction.hh"
81 #include "G4hBremsstrahlung.hh"
82 #include "G4hPairProduction.hh"
83 
88 
89 #include "G4hIonisation.hh"
90 #include "G4ionIonisation.hh"
91 #include "G4alphaIonisation.hh"
92 
93 #include "G4Gamma.hh"
94 #include "G4Electron.hh"
95 #include "G4Positron.hh"
96 #include "G4MuonPlus.hh"
97 #include "G4MuonMinus.hh"
98 #include "G4PionPlus.hh"
99 #include "G4PionMinus.hh"
100 #include "G4KaonPlus.hh"
101 #include "G4KaonMinus.hh"
102 #include "G4Proton.hh"
103 #include "G4AntiProton.hh"
104 #include "G4Deuteron.hh"
105 #include "G4Triton.hh"
106 #include "G4He3.hh"
107 #include "G4Alpha.hh"
108 #include "G4GenericIon.hh"
109 
110 #include "G4PhysicsListHelper.hh"
111 #include "G4BuilderType.hh"
112 
113 // factory
115 //
117 
118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119 
121  : G4VPhysicsConstructor("G4EmStandard_opt2"), verbose(ver)
122 {
125 }
126 
127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128 
130  : G4VPhysicsConstructor("G4EmStandard_opt2"), verbose(ver)
131 {
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137 
139 {}
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142 
144 {
145 // gamma
146  G4Gamma::Gamma();
147 
148 // leptons
153 
154 // mesons
159 
160 // barions
163 
164 // ions
167  G4He3::He3();
168  G4Alpha::Alpha();
170 }
171 
172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
173 
175 {
176  if(verbose > 1) {
177  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
178  }
180 
181  // muon & hadron bremsstrahlung and pair production
190 
191  // muon & hadron multiple scattering
193  mumsc->AddEmModel(0, new G4WentzelVIModel());
195  pmsc->AddEmModel(0, new G4WentzelVIModel());
197  pimsc->AddEmModel(0, new G4WentzelVIModel());
199  kmsc->AddEmModel(0, new G4WentzelVIModel());
200  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
201 
202  // high energy limit for e+- scattering models and bremsstrahlung
203  G4double highEnergyLimit = 100*MeV;
204 
205  // Add standard EM Processes
206  aParticleIterator->reset();
207  while( (*aParticleIterator)() ){
208  G4ParticleDefinition* particle = aParticleIterator->value();
209  G4String particleName = particle->GetParticleName();
210 
211  if (particleName == "gamma") {
212 
213  ph->RegisterProcess(new G4PhotoElectricEffect(), particle);
214  ph->RegisterProcess(new G4ComptonScattering(), particle);
215  ph->RegisterProcess(new G4GammaConversion(), particle);
216  ph->RegisterProcess(new G4RayleighScattering(), particle);
217 
218  } else if (particleName == "e-") {
219 
220  G4eIonisation* eioni = new G4eIonisation();
221  eioni->SetStepFunction(0.8, 1.0*mm);
222 
225  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
226  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
227  msc1->SetHighEnergyLimit(highEnergyLimit);
228  msc2->SetLowEnergyLimit(highEnergyLimit);
229  msc->AddEmModel(0, msc1);
230  msc->AddEmModel(0, msc2);
231 
234  ss->SetEmModel(ssm, 1);
235  ss->SetMinKinEnergy(highEnergyLimit);
236  ssm->SetLowEnergyLimit(highEnergyLimit);
237  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
238 
239  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
244  brem->SetEmModel(br1,1);
245  brem->SetEmModel(br2,2);
246  br2->SetLowEnergyLimit(GeV);
247 
248  ph->RegisterProcess(msc, particle);
249  ph->RegisterProcess(eioni, particle);
250  ph->RegisterProcess(brem, particle);
251  ph->RegisterProcess(ss, particle);
252 
253  } else if (particleName == "e+") {
254 
255  G4eIonisation* eioni = new G4eIonisation();
256  eioni->SetStepFunction(0.8, 1.0*mm);
257 
260  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
261  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
262  msc1->SetHighEnergyLimit(highEnergyLimit);
263  msc2->SetLowEnergyLimit(highEnergyLimit);
264  msc->AddEmModel(0, msc1);
265  msc->AddEmModel(0, msc2);
266 
269  ss->SetEmModel(ssm, 1);
270  ss->SetMinKinEnergy(highEnergyLimit);
271  ssm->SetLowEnergyLimit(highEnergyLimit);
272  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
273 
274  ph->RegisterProcess(msc, particle);
275  ph->RegisterProcess(eioni, particle);
276  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
277  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
278  ph->RegisterProcess(ss, particle);
279 
280  } else if (particleName == "mu+" ||
281  particleName == "mu-" ) {
282 
283  ph->RegisterProcess(mumsc, particle);
284  ph->RegisterProcess(new G4MuIonisation(), particle);
285  ph->RegisterProcess(mub, particle);
286  ph->RegisterProcess(mup, particle);
287  ph->RegisterProcess(new G4CoulombScattering(), particle);
288 
289  } else if (particleName == "alpha" ||
290  particleName == "He3") {
291 
292  //ph->RegisterProcess(hmsc, particle);
293  ph->RegisterProcess(new G4hMultipleScattering(), particle);
294  ph->RegisterProcess(new G4ionIonisation(), particle);
295 
296  } else if (particleName == "GenericIon") {
297 
298  G4ionIonisation* ionIoni = new G4ionIonisation();
299  //ionIoni->SetEmModel(new G4IonParametrisedLossModel());
300  //ionIoni->SetStepFunction(0.1, 20*um);
301 
302  ph->RegisterProcess(hmsc, particle);
303  ph->RegisterProcess(ionIoni, particle);
304 
305  } else if (particleName == "pi+" ||
306  particleName == "pi-" ) {
307 
308  ph->RegisterProcess(pimsc, particle);
309  ph->RegisterProcess(new G4hIonisation(), particle);
310  ph->RegisterProcess(pib, particle);
311  ph->RegisterProcess(pip, particle);
312 
313  } else if (particleName == "kaon+" ||
314  particleName == "kaon-" ) {
315 
316  ph->RegisterProcess(kmsc, particle);
317  ph->RegisterProcess(new G4hIonisation(), particle);
318  ph->RegisterProcess(kb, particle);
319  ph->RegisterProcess(kp, particle);
320 
321  } else if (particleName == "proton" ||
322  particleName == "anti_proton") {
323 
324  ph->RegisterProcess(pmsc, particle);
325  ph->RegisterProcess(new G4hIonisation(), particle);
326  ph->RegisterProcess(pb, particle);
327  ph->RegisterProcess(pp, particle);
328 
329  } else if (particleName == "B+" ||
330  particleName == "B-" ||
331  particleName == "D+" ||
332  particleName == "D-" ||
333  particleName == "Ds+" ||
334  particleName == "Ds-" ||
335  particleName == "anti_He3" ||
336  particleName == "anti_alpha" ||
337  particleName == "anti_deuteron" ||
338  particleName == "anti_lambda_c+" ||
339  particleName == "anti_omega-" ||
340  particleName == "anti_sigma_c+" ||
341  particleName == "anti_sigma_c++" ||
342  particleName == "anti_sigma+" ||
343  particleName == "anti_sigma-" ||
344  particleName == "anti_triton" ||
345  particleName == "anti_xi_c+" ||
346  particleName == "anti_xi-" ||
347  particleName == "deuteron" ||
348  particleName == "lambda_c+" ||
349  particleName == "omega-" ||
350  particleName == "sigma_c+" ||
351  particleName == "sigma_c++" ||
352  particleName == "sigma+" ||
353  particleName == "sigma-" ||
354  particleName == "tau+" ||
355  particleName == "tau-" ||
356  particleName == "triton" ||
357  particleName == "xi_c+" ||
358  particleName == "xi-" ) {
359 
360  ph->RegisterProcess(hmsc, particle);
361  ph->RegisterProcess(new G4hIonisation(), particle);
362  }
363  }
364 
365  // Em options
366  //
367  G4EmProcessOptions opt;
368  opt.SetVerbose(verbose);
369 
370  // Scattering options
371  //
372  opt.SetPolarAngleLimit(CLHEP::pi);
373 
374  // Ionization
375  //
376  //opt.SetSubCutoff(true);
377 
378  // Deexcitation
379  //
382 
383 }
384 
385 //....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
void SetStepFunction(G4double v1, G4double v2)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
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)
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:93
const G4String & GetPhysicsName() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:704
static G4Positron * Positron()
Definition: G4Positron.cc:94
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
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
static G4He3 * He3()
Definition: G4He3.cc:94
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetStepLimitType(G4MscStepLimitType val)
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
void SetVerbose(G4int val, const G4String &name="all", G4bool worker=false)
void SetPolarAngleLimit(G4double val)