Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmLivermorePhysics.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: G4EmLivermorePhysics.cc 99938 2016-10-12 08:06:52Z gcosmo $
27 
28 #include "G4EmLivermorePhysics.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"
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"
91 
92 // interfaces
93 #include "G4LossTableManager.hh"
94 #include "G4UAtomicDeexcitation.hh"
95 #include "G4EmParameters.hh"
96 
97 // particles
98 
99 #include "G4Gamma.hh"
100 #include "G4Electron.hh"
101 #include "G4Positron.hh"
102 #include "G4MuonPlus.hh"
103 #include "G4MuonMinus.hh"
104 #include "G4PionPlus.hh"
105 #include "G4PionMinus.hh"
106 #include "G4KaonPlus.hh"
107 #include "G4KaonMinus.hh"
108 #include "G4Proton.hh"
109 #include "G4AntiProton.hh"
110 #include "G4Deuteron.hh"
111 #include "G4Triton.hh"
112 #include "G4He3.hh"
113 #include "G4Alpha.hh"
114 #include "G4GenericIon.hh"
115 
116 //
117 #include "G4PhysicsListHelper.hh"
118 #include "G4BuilderType.hh"
119 #include "G4EmModelActivator.hh"
120 
121 // factory
123 //
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127 
129  : G4VPhysicsConstructor("G4EmLivermore"), verbose(ver)
130 {
132  param->SetDefaults();
133  param->SetVerbose(verbose);
134  param->SetMinEnergy(100*eV);
135  param->SetMaxEnergy(1*TeV);
136  param->SetLowestElectronEnergy(100*eV);
137  param->SetNumberOfBinsPerDecade(20);
139  param->SetMscRangeFactor(0.02);
140  param->SetMuHadLateralDisplacement(true);
142  param->SetFluo(true);
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147 
149 {}
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152 
154 {
155  // gamma
156  G4Gamma::Gamma();
157 
158  // leptons
163 
164  // mesons
169 
170  // baryons
173 
174  // ions
177  G4He3::He3();
178  G4Alpha::Alpha();
180 }
181 
182 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183 
185 {
186  if(verbose > 1) {
187  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
188  }
190 
191  // muon & hadron bremsstrahlung and pair production
200 
201  // muon & hadron multiple scattering
203  mumsc->AddEmModel(0, new G4WentzelVIModel());
204  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
205  //pimsc->AddEmModel(0, new G4WentzelVIModel());
206  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
207  //kmsc->AddEmModel(0, new G4WentzelVIModel());
208  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
209  //pmsc->AddEmModel(0, new G4WentzelVIModel());
210  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
211 
212  // high energy limit for e+- scattering models
213  G4double highEnergyLimit = 100*MeV;
214 
215  // nuclear stopping
216  G4NuclearStopping* pnuc = new G4NuclearStopping();
217 
218  // Add Livermore EM Processes
219  auto myParticleIterator=GetParticleIterator();
220  myParticleIterator->reset();
221 
222  while( (*myParticleIterator)() ){
223 
224  G4ParticleDefinition* particle = myParticleIterator->value();
225  G4String particleName = particle->GetParticleName();
226 
227  if (particleName == "gamma") {
228 
229  // photoelectric effect - Livermore model only
230  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
231  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel(), 1);
232  ph->RegisterProcess(thePhotoElectricEffect, particle);
233 
234  // Compton scattering - Livermore model only
235  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
236  theComptonScattering->SetEmModel(new G4LivermoreComptonModel(),1);
237  ph->RegisterProcess(theComptonScattering, particle);
238 
239  // gamma conversion - Livermore model below 80 GeV
240  G4GammaConversion* theGammaConversion = new G4GammaConversion();
241  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel(),1);
242  ph->RegisterProcess(theGammaConversion, particle);
243 
244  // default Rayleigh scattering is Livermore
245  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
246  ph->RegisterProcess(theRayleigh, particle);
247 
248  } else if (particleName == "e-") {
249 
250  // multiple scattering
252  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
253  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
254  msc1->SetHighEnergyLimit(highEnergyLimit);
255  msc2->SetLowEnergyLimit(highEnergyLimit);
256  msc->AddEmModel(0, msc1);
257  msc->AddEmModel(0, msc2);
258 
261  ss->SetEmModel(ssm, 1);
262  ss->SetMinKinEnergy(highEnergyLimit);
263  ssm->SetLowEnergyLimit(highEnergyLimit);
264  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
265 
266  // Ionisation - Livermore should be used only for low energies
267  G4eIonisation* eIoni = new G4eIonisation();
268  G4LivermoreIonisationModel* theIoniLivermore = new
270  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
271  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
272  eIoni->SetStepFunction(0.2, 100*um); //
273 
274  // Bremsstrahlung
275  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
276  G4VEmModel* theBremLivermore = new G4LivermoreBremsstrahlungModel();
277  theBremLivermore->SetHighEnergyLimit(1*GeV);
278  theBremLivermore->SetAngularDistribution(new G4Generator2BS());
279  eBrem->SetEmModel(theBremLivermore,1);
280 
281  // register processes
282  ph->RegisterProcess(msc, particle);
283  ph->RegisterProcess(eIoni, particle);
284  ph->RegisterProcess(eBrem, particle);
285  ph->RegisterProcess(ss, particle);
286 
287  } else if (particleName == "e+") {
288 
289  // multiple scattering
291  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
292  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
293  msc1->SetHighEnergyLimit(highEnergyLimit);
294  msc2->SetLowEnergyLimit(highEnergyLimit);
295  msc->AddEmModel(0, msc1);
296  msc->AddEmModel(0, msc2);
297 
300  ss->SetEmModel(ssm, 1);
301  ss->SetMinKinEnergy(highEnergyLimit);
302  ssm->SetLowEnergyLimit(highEnergyLimit);
303  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
304 
305  // ionisation
306  G4eIonisation* eIoni = new G4eIonisation();
307  eIoni->SetStepFunction(0.2, 100*um);
308 
309  // register processes
310  ph->RegisterProcess(msc, particle);
311  ph->RegisterProcess(eIoni, particle);
312  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
313  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
314  ph->RegisterProcess(ss, particle);
315 
316  } else if (particleName == "mu+" ||
317  particleName == "mu-" ) {
318 
319  G4MuIonisation* muIoni = new G4MuIonisation();
320  muIoni->SetStepFunction(0.2, 50*um);
321 
322  ph->RegisterProcess(mumsc, particle);
323  ph->RegisterProcess(muIoni, particle);
324  ph->RegisterProcess(mub, particle);
325  ph->RegisterProcess(mup, particle);
326  ph->RegisterProcess(new G4CoulombScattering(), particle);
327 
328  } else if (particleName == "alpha" ||
329  particleName == "He3" ) {
330 
332  G4ionIonisation* ionIoni = new G4ionIonisation();
333  ionIoni->SetStepFunction(0.1, 10*um);
334 
335  ph->RegisterProcess(msc, particle);
336  ph->RegisterProcess(ionIoni, particle);
337  ph->RegisterProcess(pnuc, particle);
338 
339  } else if (particleName == "GenericIon") {
340 
341  G4ionIonisation* ionIoni = new G4ionIonisation();
342  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
343  ionIoni->SetStepFunction(0.1, 1*um);
344 
345  ph->RegisterProcess(hmsc, particle);
346  ph->RegisterProcess(ionIoni, particle);
347  ph->RegisterProcess(pnuc, particle);
348 
349  } else if (particleName == "pi+" ||
350  particleName == "pi-" ) {
351 
353  G4hIonisation* hIoni = new G4hIonisation();
354  hIoni->SetStepFunction(0.2, 50*um);
355 
356  ph->RegisterProcess(pimsc, particle);
357  ph->RegisterProcess(hIoni, particle);
358  ph->RegisterProcess(pib, particle);
359  ph->RegisterProcess(pip, particle);
360 
361  } else if (particleName == "kaon+" ||
362  particleName == "kaon-" ) {
363 
365  G4hIonisation* hIoni = new G4hIonisation();
366  hIoni->SetStepFunction(0.2, 50*um);
367 
368  ph->RegisterProcess(kmsc, particle);
369  ph->RegisterProcess(hIoni, particle);
370  ph->RegisterProcess(kb, particle);
371  ph->RegisterProcess(kp, particle);
372 
373  } else if (particleName == "proton" ||
374  particleName == "anti_proton") {
375 
377  G4hIonisation* hIoni = new G4hIonisation();
378  hIoni->SetStepFunction(0.2, 50*um);
379 
380  ph->RegisterProcess(pmsc, particle);
381  ph->RegisterProcess(hIoni, particle);
382  ph->RegisterProcess(pb, particle);
383  ph->RegisterProcess(pp, particle);
384  ph->RegisterProcess(pnuc, particle);
385 
386  } else if (particleName == "B+" ||
387  particleName == "B-" ||
388  particleName == "D+" ||
389  particleName == "D-" ||
390  particleName == "Ds+" ||
391  particleName == "Ds-" ||
392  particleName == "anti_He3" ||
393  particleName == "anti_alpha" ||
394  particleName == "anti_deuteron" ||
395  particleName == "anti_lambda_c+" ||
396  particleName == "anti_omega-" ||
397  particleName == "anti_sigma_c+" ||
398  particleName == "anti_sigma_c++" ||
399  particleName == "anti_sigma+" ||
400  particleName == "anti_sigma-" ||
401  particleName == "anti_triton" ||
402  particleName == "anti_xi_c+" ||
403  particleName == "anti_xi-" ||
404  particleName == "deuteron" ||
405  particleName == "lambda_c+" ||
406  particleName == "omega-" ||
407  particleName == "sigma_c+" ||
408  particleName == "sigma_c++" ||
409  particleName == "sigma+" ||
410  particleName == "sigma-" ||
411  particleName == "tau+" ||
412  particleName == "tau-" ||
413  particleName == "triton" ||
414  particleName == "xi_c+" ||
415  particleName == "xi-" ) {
416 
417  ph->RegisterProcess(hmsc, particle);
418  ph->RegisterProcess(new G4hIonisation(), particle);
419  ph->RegisterProcess(pnuc, particle);
420  }
421  }
422 
423  // Nuclear stopping
424  pnuc->SetMaxKinEnergy(MeV);
425 
426  // Deexcitation
427  //
430 
432 }
433 
434 //....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
void SetVerbose(G4int val)
virtual void ConstructParticle()
static G4LossTableManager * Instance()
void SetLowestElectronEnergy(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
int G4int
Definition: G4Types.hh:78
void SetMaxEnergy(G4double val)
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static constexpr double TeV
Definition: G4SIunits.hh:218
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
static constexpr double um
Definition: G4SIunits.hh:113
G4EmLivermorePhysics(G4int ver=1, const G4String &name="")
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
void SetMscRangeFactor(G4double val)
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 constexpr double eV
Definition: G4SIunits.hh:215
const G4String & GetPhysicsName() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetMaxKinEnergy(G4double e)
void SetMuHadLateralDisplacement(G4bool val)
void SetMinEnergy(G4double val)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:626
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4EmParameters * Instance()
void SetEmModel(G4VEmModel *, G4int index=1)
static constexpr double GeV
Definition: G4SIunits.hh:217
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=nullptr)
void ActivateAngularGeneratorForIonisation(G4bool val)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
static G4He3 * He3()
Definition: G4He3.cc:94
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetFluo(G4bool val)
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)