Geant4_10
G4EmStandardPhysics_option1.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_option1.cc 73153 2013-08-20 07:17:42Z gcosmo $
27 //
28 //---------------------------------------------------------------------------
29 //
30 // ClassName: G4EmStandardPhysics_option1
31 //
32 // Author: V.Ivanchenko 09.11.2005
33 //
34 // Modified:
35 // 19.12.2005 V.Ivanchenko StepFunction for electrons (1.0, 1.0*mm)
36 // 21.06.2006 V.Ivanchenko use recent msc with step limitation off
37 // 23.06.2006 V.Ivanchenko set dRoverRange = 0.8 for e- and e+
38 // 13.11.2006 V.Ivanchenko use G4hMultipleScattering
39 // 23.11.2006 V.Ivanchenko remove mscStepLimit option and improve cout
40 // 13.02.2007 V.Ivanchenko set skin=0.0
41 // 15.05.2007 V.Ivanchenko rename to _option1
42 // 21.04.2008 V.Ivanchenko add long-lived D and B mesons
43 // 19.11.2010 V.Ivanchenko added WentzelVI model for muons;
44 // disable ApplyCut option for compatibility with 9.2
45 //
46 //----------------------------------------------------------------------------
47 //
48 
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 
59 #include "G4eMultipleScattering.hh"
61 #include "G4hMultipleScattering.hh"
63 #include "G4CoulombScattering.hh"
64 #include "G4WentzelVIModel.hh"
65 #include "G4UrbanMscModel.hh"
66 
67 #include "G4eIonisation.hh"
68 #include "G4eBremsstrahlung.hh"
69 #include "G4eplusAnnihilation.hh"
70 #include "G4UAtomicDeexcitation.hh"
71 
72 #include "G4MuIonisation.hh"
73 #include "G4MuBremsstrahlung.hh"
74 #include "G4MuPairProduction.hh"
75 #include "G4hBremsstrahlung.hh"
76 #include "G4hPairProduction.hh"
77 
82 
83 #include "G4hIonisation.hh"
84 #include "G4ionIonisation.hh"
85 #include "G4alphaIonisation.hh"
86 
87 #include "G4Gamma.hh"
88 #include "G4Electron.hh"
89 #include "G4Positron.hh"
90 #include "G4MuonPlus.hh"
91 #include "G4MuonMinus.hh"
92 #include "G4PionPlus.hh"
93 #include "G4PionMinus.hh"
94 #include "G4KaonPlus.hh"
95 #include "G4KaonMinus.hh"
96 #include "G4Proton.hh"
97 #include "G4AntiProton.hh"
98 #include "G4Deuteron.hh"
99 #include "G4Triton.hh"
100 #include "G4He3.hh"
101 #include "G4Alpha.hh"
102 #include "G4GenericIon.hh"
103 
104 #include "G4PhysicsListHelper.hh"
105 #include "G4BuilderType.hh"
106 
107 // factory
109 //
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113 
115  : G4VPhysicsConstructor("G4EmStandard_opt1"), verbose(ver)
116 {
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122 
124  : G4VPhysicsConstructor("G4EmStandard_opt1"), verbose(ver)
125 {
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131 
133 {}
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136 
138 {
139 // gamma
140  G4Gamma::Gamma();
141 
142 // leptons
147 
148 // mesons
153 
154 // barions
157 
158 // ions
161  G4He3::He3();
162  G4Alpha::Alpha();
164 }
165 
166 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
167 
169 {
170  if(verbose > 1) {
171  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
172  }
174 
175  // muon & hadron bremsstrahlung and pair production
184 
185  // muon & hadron multiple scattering
187  mumsc->AddEmModel(0, new G4WentzelVIModel());
189  pimsc->AddEmModel(0, new G4WentzelVIModel());
191  kmsc->AddEmModel(0, new G4WentzelVIModel());
193  pmsc->AddEmModel(0, new G4WentzelVIModel());
194  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
195 
196  // high energy limit for e+- scattering models and bremsstrahlung
197  G4double highEnergyLimit = 100*MeV;
198 
199  // Add standard EM Processes
200  aParticleIterator->reset();
201  while( (*aParticleIterator)() ){
202  G4ParticleDefinition* particle = aParticleIterator->value();
203  G4String particleName = particle->GetParticleName();
204 
205  if (particleName == "gamma") {
206 
207  ph->RegisterProcess(new G4PhotoElectricEffect(), particle);
208  ph->RegisterProcess(new G4ComptonScattering(), particle);
209  ph->RegisterProcess(new G4GammaConversion(), particle);
210 
211  } else if (particleName == "e-") {
212 
213  G4eIonisation* eioni = new G4eIonisation();
214  eioni->SetStepFunction(0.8, 1.0*mm);
215 
218  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
219  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
220  msc1->SetHighEnergyLimit(highEnergyLimit);
221  msc2->SetLowEnergyLimit(highEnergyLimit);
222  msc->AddEmModel(0, msc1);
223  msc->AddEmModel(0, msc2);
224 
227  ss->SetEmModel(ssm, 1);
228  ss->SetMinKinEnergy(highEnergyLimit);
229  ssm->SetLowEnergyLimit(highEnergyLimit);
230  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
231 
232  ph->RegisterProcess(msc, particle);
233  ph->RegisterProcess(eioni, particle);
234  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
235  ph->RegisterProcess(ss, particle);
236 
237  } else if (particleName == "e+") {
238 
239  G4eIonisation* eioni = new G4eIonisation();
240  eioni->SetStepFunction(0.8, 1.0*mm);
241 
244  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
245  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
246  msc1->SetHighEnergyLimit(highEnergyLimit);
247  msc2->SetLowEnergyLimit(highEnergyLimit);
248  msc->AddEmModel(0, msc1);
249  msc->AddEmModel(0, msc2);
250 
253  ss->SetEmModel(ssm, 1);
254  ss->SetMinKinEnergy(highEnergyLimit);
255  ssm->SetLowEnergyLimit(highEnergyLimit);
256  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
257 
258  ph->RegisterProcess(msc, particle);
259  ph->RegisterProcess(eioni, particle);
260  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
261  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
262  ph->RegisterProcess(ss, particle);
263 
264  } else if (particleName == "mu+" ||
265  particleName == "mu-" ) {
266 
267  ph->RegisterProcess(mumsc, particle);
268  ph->RegisterProcess(new G4MuIonisation(), particle);
269  ph->RegisterProcess(mub, particle);
270  ph->RegisterProcess(mup, particle);
271  ph->RegisterProcess(new G4CoulombScattering(), particle);
272 
273  } else if (particleName == "alpha" ||
274  particleName == "He3" ) {
275 
276  //ph->RegisterProcess(hmsc, particle);
277  ph->RegisterProcess(new G4hMultipleScattering(), particle);
278  ph->RegisterProcess(new G4ionIonisation(), particle);
279 
280  } else if (particleName == "GenericIon") {
281 
282  ph->RegisterProcess(hmsc, particle);
283  ph->RegisterProcess(new G4ionIonisation(), particle);
284 
285  } else if (particleName == "pi+" ||
286  particleName == "pi-" ) {
287 
288  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
289  ph->RegisterProcess(pimsc, particle);
290  ph->RegisterProcess(new G4hIonisation(), particle);
291  ph->RegisterProcess(pib, particle);
292  ph->RegisterProcess(pip, particle);
293 
294  } else if (particleName == "kaon+" ||
295  particleName == "kaon-" ) {
296 
297  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
298  ph->RegisterProcess(kmsc, particle);
299  ph->RegisterProcess(new G4hIonisation(), particle);
300  ph->RegisterProcess(kb, particle);
301  ph->RegisterProcess(kp, particle);
302 
303  // } else if (particleName == "proton" ) {
304  } else if (particleName == "proton" ||
305  particleName == "anti_proton") {
306 
307  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
308  ph->RegisterProcess(pmsc, particle);
309  ph->RegisterProcess(new G4hIonisation(), particle);
310  ph->RegisterProcess(pb, particle);
311  ph->RegisterProcess(pp, particle);
312 
313  } else if (particleName == "B+" ||
314  particleName == "B-" ||
315  particleName == "D+" ||
316  particleName == "D-" ||
317  particleName == "Ds+" ||
318  particleName == "Ds-" ||
319  particleName == "anti_He3" ||
320  particleName == "anti_alpha" ||
321  particleName == "anti_deuteron" ||
322  particleName == "anti_lambda_c+" ||
323  particleName == "anti_omega-" ||
324  particleName == "anti_sigma_c+" ||
325  particleName == "anti_sigma_c++" ||
326  particleName == "anti_sigma+" ||
327  particleName == "anti_sigma-" ||
328  particleName == "anti_triton" ||
329  particleName == "anti_xi_c+" ||
330  particleName == "anti_xi-" ||
331  particleName == "deuteron" ||
332  particleName == "lambda_c+" ||
333  particleName == "omega-" ||
334  particleName == "sigma_c+" ||
335  particleName == "sigma_c++" ||
336  particleName == "sigma+" ||
337  particleName == "sigma-" ||
338  particleName == "tau+" ||
339  particleName == "tau-" ||
340  particleName == "triton" ||
341  particleName == "xi_c+" ||
342  particleName == "xi-" ) {
343 
344  ph->RegisterProcess(hmsc, particle);
345  ph->RegisterProcess(new G4hIonisation(), particle);
346  }
347  }
348  G4EmProcessOptions opt;
349  opt.SetVerbose(verbose);
350  opt.SetPolarAngleLimit(CLHEP::pi);
351  opt.SetApplyCuts(true);
352 
353  // Deexcitation
354  //
357 }
358 
359 //....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)
static G4PhysicsListHelper * GetPhysicsListHelper()
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 SetApplyCuts(G4bool val)
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)