Geant4_10
G4EmStandardPhysics_option4.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_option4.cc 73153 2013-08-20 07:17:42Z gcosmo $
27 //
28 //---------------------------------------------------------------------------
29 //
30 // ClassName: G4EmStandardPhysics_option4
31 //
32 // Author: V.Ivanchenko 28.09.2012 from Option3 physics constructor
33 //
34 // Modified:
35 //
36 //----------------------------------------------------------------------------
37 //
38 
40 
41 #include "G4SystemOfUnits.hh"
42 #include "G4ParticleDefinition.hh"
43 #include "G4LossTableManager.hh"
44 #include "G4EmProcessOptions.hh"
45 
46 #include "G4ComptonScattering.hh"
47 #include "G4GammaConversion.hh"
48 #include "G4PhotoElectricEffect.hh"
49 #include "G4RayleighScattering.hh"
50 #include "G4PEEffectFluoModel.hh"
51 #include "G4KleinNishinaModel.hh"
52 #include "G4LowEPComptonModel.hh"
55 
56 #include "G4eMultipleScattering.hh"
58 #include "G4hMultipleScattering.hh"
59 #include "G4MscStepLimitType.hh"
60 #include "G4UrbanMscModel.hh"
61 #include "G4DummyModel.hh"
62 #include "G4WentzelVIModel.hh"
63 #include "G4CoulombScattering.hh"
65 
66 #include "G4eIonisation.hh"
67 #include "G4eBremsstrahlung.hh"
68 #include "G4Generator2BS.hh"
69 #include "G4Generator2BN.hh"
70 #include "G4SeltzerBergerModel.hh"
73 
74 #include "G4eplusAnnihilation.hh"
75 #include "G4UAtomicDeexcitation.hh"
76 
77 #include "G4MuIonisation.hh"
78 #include "G4MuBremsstrahlung.hh"
79 #include "G4MuPairProduction.hh"
80 #include "G4hBremsstrahlung.hh"
81 #include "G4hPairProduction.hh"
82 
87 
88 #include "G4hIonisation.hh"
89 #include "G4ionIonisation.hh"
91 #include "G4NuclearStopping.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_opt4"), verbose(ver)
122 {
125 }
126 
127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128 
130  : G4VPhysicsConstructor("G4EmStandard_opt3"), 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  pimsc->AddEmModel(0, new G4WentzelVIModel());
197  kmsc->AddEmModel(0, new G4WentzelVIModel());
199  pmsc->AddEmModel(0, new G4WentzelVIModel());
200  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
201 
202  // high energy limit for e+- scattering models
203  G4double highEnergyLimit = 100*MeV;
204 
205  // nuclear stopping
206  G4NuclearStopping* ionnuc = new G4NuclearStopping();
207  G4NuclearStopping* pnuc = new G4NuclearStopping();
208 
209  // Add standard EM Processes
210  aParticleIterator->reset();
211  while( (*aParticleIterator)() ){
212  G4ParticleDefinition* particle = aParticleIterator->value();
213  G4String particleName = particle->GetParticleName();
214 
215  if (particleName == "gamma") {
216 
217  // Compton scattering
219  cs->SetEmModel(new G4KleinNishinaModel(),1);
220  G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
221  theLowEPComptonModel->SetHighEnergyLimit(20*MeV);
222  cs->AddEmModel(0, theLowEPComptonModel);
223  ph->RegisterProcess(cs, particle);
224 
225  // Photoelectric
227  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
228  theLivermorePEModel->SetHighEnergyLimit(10*GeV);
229  pe->SetEmModel(theLivermorePEModel,1);
230  ph->RegisterProcess(pe, particle);
231 
232  // Gamma conversion
234  G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel();
235  thePenelopeGCModel->SetHighEnergyLimit(1*GeV);
236  gc->SetEmModel(thePenelopeGCModel,1);
237  ph->RegisterProcess(gc, particle);
238 
239  // Rayleigh scattering
240  ph->RegisterProcess(new G4RayleighScattering(), particle);
241 
242  } else if (particleName == "e-") {
243 
244  // multiple scattering
247  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
248  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
249  msc1->SetHighEnergyLimit(highEnergyLimit);
250  msc2->SetLowEnergyLimit(highEnergyLimit);
251  msc->AddEmModel(0, msc1);
252  msc->AddEmModel(0, msc2);
253 
256  ss->SetEmModel(ssm, 1);
257  ss->SetMinKinEnergy(highEnergyLimit);
258  ssm->SetLowEnergyLimit(highEnergyLimit);
259  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
260 
261  // ionisation
262  G4eIonisation* eIoni = new G4eIonisation();
263  eIoni->SetStepFunction(0.2, 100*um);
264  G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
265  theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
266  eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
267 
268  // bremsstrahlung
269  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
270 
271  ph->RegisterProcess(msc, particle);
272  ph->RegisterProcess(eIoni, particle);
273  ph->RegisterProcess(eBrem, particle);
274  ph->RegisterProcess(ss, particle);
275 
276  } else if (particleName == "e+") {
277 
278  // multiple scattering
281  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
282  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
283  msc1->SetHighEnergyLimit(highEnergyLimit);
284  msc2->SetLowEnergyLimit(highEnergyLimit);
285  msc->AddEmModel(0, msc1);
286  msc->AddEmModel(0, msc2);
287 
290  ss->SetEmModel(ssm, 1);
291  ss->SetMinKinEnergy(highEnergyLimit);
292  ssm->SetLowEnergyLimit(highEnergyLimit);
293  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
294 
295  // ionisation
296  G4eIonisation* eIoni = new G4eIonisation();
297  eIoni->SetStepFunction(0.2, 100*um);
298  G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
299  theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
300  eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
301 
302  // bremsstrahlung
303  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
304 
305  ph->RegisterProcess(msc, particle);
306  ph->RegisterProcess(eIoni, particle);
307  ph->RegisterProcess(eBrem, particle);
308  ph->RegisterProcess(ss, particle);
309 
310  // annihilation at rest and in flight
311  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
312 
313  } else if (particleName == "mu+" ||
314  particleName == "mu-" ) {
315 
316  G4MuIonisation* muIoni = new G4MuIonisation();
317  muIoni->SetStepFunction(0.2, 50*um);
318 
319  ph->RegisterProcess(mumsc, particle);
320  ph->RegisterProcess(muIoni, particle);
321  ph->RegisterProcess(mub, particle);
322  ph->RegisterProcess(mup, particle);
323  ph->RegisterProcess(new G4CoulombScattering(), particle);
324 
325  } else if (particleName == "alpha" ||
326  particleName == "He3") {
327 
329  G4ionIonisation* ionIoni = new G4ionIonisation();
330  ionIoni->SetStepFunction(0.1, 10*um);
331 
332  ph->RegisterProcess(msc, particle);
333  ph->RegisterProcess(ionIoni, particle);
334  ph->RegisterProcess(ionnuc, particle);
335 
336  } else if (particleName == "GenericIon") {
337 
338  G4ionIonisation* ionIoni = new G4ionIonisation();
339  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
340  ionIoni->SetStepFunction(0.1, 1*um);
341 
342  ph->RegisterProcess(hmsc, particle);
343  ph->RegisterProcess(ionIoni, particle);
344  ph->RegisterProcess(ionnuc, particle);
345 
346  } else if (particleName == "pi+" ||
347  particleName == "pi-" ) {
348 
349  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
350  G4hIonisation* hIoni = new G4hIonisation();
351  hIoni->SetStepFunction(0.2, 50*um);
352 
353  ph->RegisterProcess(pimsc, particle);
354  ph->RegisterProcess(hIoni, particle);
355  ph->RegisterProcess(pib, particle);
356  ph->RegisterProcess(pip, particle);
357 
358  } else if (particleName == "kaon+" ||
359  particleName == "kaon-" ) {
360 
361  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
362  G4hIonisation* hIoni = new G4hIonisation();
363  hIoni->SetStepFunction(0.2, 50*um);
364 
365  ph->RegisterProcess(kmsc, particle);
366  ph->RegisterProcess(hIoni, particle);
367  ph->RegisterProcess(kb, particle);
368  ph->RegisterProcess(kp, particle);
369 
370  } else if (particleName == "proton" ||
371  particleName == "anti_proton") {
372 
373  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
374  G4hIonisation* hIoni = new G4hIonisation();
375  hIoni->SetStepFunction(0.2, 50*um);
376 
377  ph->RegisterProcess(pmsc, particle);
378  ph->RegisterProcess(hIoni, particle);
379  ph->RegisterProcess(pb, particle);
380  ph->RegisterProcess(pp, particle);
381  ph->RegisterProcess(pnuc, particle);
382 
383  } else if (particleName == "B+" ||
384  particleName == "B-" ||
385  particleName == "D+" ||
386  particleName == "D-" ||
387  particleName == "Ds+" ||
388  particleName == "Ds-" ||
389  particleName == "anti_He3" ||
390  particleName == "anti_alpha" ||
391  particleName == "anti_deuteron" ||
392  particleName == "anti_lambda_c+" ||
393  particleName == "anti_omega-" ||
394  particleName == "anti_sigma_c+" ||
395  particleName == "anti_sigma_c++" ||
396  particleName == "anti_sigma+" ||
397  particleName == "anti_sigma-" ||
398  particleName == "anti_triton" ||
399  particleName == "anti_xi_c+" ||
400  particleName == "anti_xi-" ||
401  particleName == "deuteron" ||
402  particleName == "lambda_c+" ||
403  particleName == "omega-" ||
404  particleName == "sigma_c+" ||
405  particleName == "sigma_c++" ||
406  particleName == "sigma+" ||
407  particleName == "sigma-" ||
408  particleName == "tau+" ||
409  particleName == "tau-" ||
410  particleName == "triton" ||
411  particleName == "xi_c+" ||
412  particleName == "xi-" ) {
413 
414  ph->RegisterProcess(hmsc, particle);
415  ph->RegisterProcess(new G4hIonisation(), particle);
416  ph->RegisterProcess(pnuc, particle);
417  }
418  }
419 
420  // Em options
421  //
422  G4EmProcessOptions opt;
423  opt.SetVerbose(verbose);
424 
425  // Multiple Coulomb scattering
426  //
427  opt.SetPolarAngleLimit(CLHEP::pi);
428 
429  // Physics tables
430  //
431  opt.SetMinEnergy(100*eV);
432  opt.SetMaxEnergy(10*TeV);
433  opt.SetDEDXBinning(220);
434  opt.SetLambdaBinning(220);
435 
436  // Nuclear stopping
437  pnuc->SetMaxKinEnergy(MeV);
438 
439  // Ionization
440  //
441  //opt.SetSubCutoff(true);
442 
443  // Deexcitation
446  de->SetFluo(true);
447 }
448 
449 //....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 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)
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
void SetLambdaBinning(G4int val)
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
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
void SetMaxEnergy(G4double val)
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:704
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
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)