Geant4  10.02.p03
G4EmLivermorePolarizedPhysics.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: G4EmLivermorePolarizedPhysics.cc 102321 2017-01-23 09:51:59Z gcosmo $
27 
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 
57 // e+
58 #include "G4eplusAnnihilation.hh"
59 
60 // mu+-
62 #include "G4MuIonisation.hh"
63 #include "G4MuBremsstrahlung.hh"
64 #include "G4MuPairProduction.hh"
65 
70 
71 // hadrons
72 #include "G4hMultipleScattering.hh"
73 #include "G4MscStepLimitType.hh"
74 
75 #include "G4hBremsstrahlung.hh"
76 #include "G4hPairProduction.hh"
77 
78 #include "G4hIonisation.hh"
79 #include "G4ionIonisation.hh"
80 #include "G4alphaIonisation.hh"
82 #include "G4NuclearStopping.hh"
83 
84 // msc models
85 #include "G4UrbanMscModel.hh"
86 #include "G4WentzelVIModel.hh"
88 #include "G4CoulombScattering.hh"
90 
91 // interfaces
92 #include "G4LossTableManager.hh"
93 #include "G4EmParameters.hh"
94 #include "G4UAtomicDeexcitation.hh"
95 
96 // particles
97 #include "G4Gamma.hh"
98 #include "G4Electron.hh"
99 #include "G4Positron.hh"
100 #include "G4MuonPlus.hh"
101 #include "G4MuonMinus.hh"
102 #include "G4PionPlus.hh"
103 #include "G4PionMinus.hh"
104 #include "G4KaonPlus.hh"
105 #include "G4KaonMinus.hh"
106 #include "G4Proton.hh"
107 #include "G4AntiProton.hh"
108 #include "G4Deuteron.hh"
109 #include "G4Triton.hh"
110 #include "G4He3.hh"
111 #include "G4Alpha.hh"
112 #include "G4GenericIon.hh"
113 
114 //
115 #include "G4PhysicsListHelper.hh"
116 #include "G4BuilderType.hh"
117 #include "G4EmModelActivator.hh"
118 
119 // factory
121 //
123 
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126 
128  : G4VPhysicsConstructor("G4EmLivermorePolarizedPhysics"), verbose(ver)
129 {
131  param->SetDefaults();
132  param->SetVerbose(verbose);
133  param->SetMinEnergy(100*eV);
134  param->SetMaxEnergy(10*TeV);
135  param->SetLowestElectronEnergy(100*eV);
136  param->SetNumberOfBinsPerDecade(20);
138  param->SetFluo(true);
140 }
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143 
145  : G4VPhysicsConstructor("G4EmLivermorePolarizedPhysics"), verbose(ver)
146 {
148  param->SetDefaults();
149  param->SetVerbose(verbose);
150  param->SetMinEnergy(100*eV);
151  param->SetMaxEnergy(10*TeV);
152  param->SetNumberOfBinsPerDecade(20);
154  param->SetFluo(true);
156 }
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159 
161 {}
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164 
166 {
167  // gamma
168  G4Gamma::Gamma();
169 
170  // leptons
175 
176  // mesons
181 
182  // baryons
185 
186  // ions
189  G4He3::He3();
190  G4Alpha::Alpha();
192 
193  // dna
194  G4EmModelActivator mact;
195  mact.ConstructParticle();
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199 
201 {
202  if(verbose > 1) {
203  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
204  }
206 
207  // muon & hadron bremsstrahlung and pair production
216 
217  // muon & hadron multiple scattering
219  mumsc->AddEmModel(0, new G4WentzelVIModel());
221  pimsc->AddEmModel(0, new G4WentzelVIModel());
223  kmsc->AddEmModel(0, new G4WentzelVIModel());
225  pmsc->AddEmModel(0, new G4WentzelVIModel());
226  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
227 
228  // high energy limit for e+- scattering models
229  G4double highEnergyLimit = 100*MeV;
230 
231  // nuclear stopping
232  G4NuclearStopping* pnuc = new G4NuclearStopping();
233 
234  // Add Livermore EM Processes
235  auto myParticleIterator=GetParticleIterator();
236  myParticleIterator->reset();
237 
238  while( (*myParticleIterator)() ){
239 
240  G4ParticleDefinition* particle = myParticleIterator->value();
241  G4String particleName = particle->GetParticleName();
242 
243  //Applicability range for Livermore models
244  //for higher energies, the Standard models are used
245  G4double LivermoreHighEnergyLimit = GeV;
246 
247  if (particleName == "gamma") {
248 
249  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
251  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
252  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
253  ph->RegisterProcess(thePhotoElectricEffect, particle);
254 
255  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
256  G4LivermorePolarizedComptonModel* theLivermoreComptonModel = new G4LivermorePolarizedComptonModel();
257  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
258  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
259  ph->RegisterProcess(theComptonScattering, particle);
260 
261  G4GammaConversion* theGammaConversion = new G4GammaConversion();
263  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
264  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
265  ph->RegisterProcess(theGammaConversion, particle);
266 
267  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
269  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
270  theRayleigh->AddEmModel(0, theRayleighModel);
271  ph->RegisterProcess(theRayleigh, particle);
272 
273  } else if (particleName == "e-") {
274 
275  // multiple scattering
278  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
279  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
280  msc1->SetHighEnergyLimit(highEnergyLimit);
281  msc2->SetLowEnergyLimit(highEnergyLimit);
282  msc->AddEmModel(0, msc1);
283  msc->AddEmModel(0, msc2);
284 
287  ss->SetEmModel(ssm, 1);
288  ss->SetMinKinEnergy(highEnergyLimit);
289  ssm->SetLowEnergyLimit(highEnergyLimit);
290  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
291  ph->RegisterProcess(msc, particle);
292  ph->RegisterProcess(ss, particle);
293 
294  // Ionisation
295  G4eIonisation* eIoni = new G4eIonisation();
296  G4LivermoreIonisationModel* theIoniLivermore = new
298  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
299  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
300  eIoni->SetStepFunction(0.2, 100*um); //
301  ph->RegisterProcess(eIoni, particle);
302 
303  // Bremsstrahlung from standard
304  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
305  ph->RegisterProcess(eBrem, particle);
306 
307  } else if (particleName == "e+") {
308 
309  // multiple scattering
312  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
313  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
314  msc1->SetHighEnergyLimit(highEnergyLimit);
315  msc2->SetLowEnergyLimit(highEnergyLimit);
316  msc->AddEmModel(0, msc1);
317  msc->AddEmModel(0, msc2);
318 
321  ss->SetEmModel(ssm, 1);
322  ss->SetMinKinEnergy(highEnergyLimit);
323  ssm->SetLowEnergyLimit(highEnergyLimit);
324  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
325 
326  // Ionisation
327  G4eIonisation* eIoni = new G4eIonisation();
328  eIoni->SetStepFunction(0.2, 100*um);
329 
330  ph->RegisterProcess(msc, particle);
331  ph->RegisterProcess(eIoni, particle);
332  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
333  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
334  ph->RegisterProcess(ss, particle);
335 
336  } else if (particleName == "mu+" ||
337  particleName == "mu-" ) {
338 
339  G4MuIonisation* muIoni = new G4MuIonisation();
340  muIoni->SetStepFunction(0.2, 50*um);
341 
342  ph->RegisterProcess(mumsc, particle);
343  ph->RegisterProcess(muIoni, particle);
344  ph->RegisterProcess(mub, particle);
345  ph->RegisterProcess(mup, particle);
346  ph->RegisterProcess(new G4CoulombScattering(), particle);
347 
348  } else if (particleName == "alpha" ||
349  particleName == "He3" ) {
350 
351  // Identical to G4EmStandardPhysics_option3
352 
354  G4ionIonisation* ionIoni = new G4ionIonisation();
355  ionIoni->SetStepFunction(0.1, 10*um);
356 
357  ph->RegisterProcess(msc, particle);
358  ph->RegisterProcess(ionIoni, particle);
359  ph->RegisterProcess(pnuc, particle);
360 
361  } else if (particleName == "GenericIon") {
362 
363  // Identical to G4EmStandardPhysics_option3
364 
365  G4ionIonisation* ionIoni = new G4ionIonisation();
366  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
367  ionIoni->SetStepFunction(0.1, 1*um);
368 
369  ph->RegisterProcess(hmsc, particle);
370  ph->RegisterProcess(ionIoni, particle);
371  ph->RegisterProcess(pnuc, particle);
372 
373  } else if (particleName == "pi+" ||
374  particleName == "pi-" ) {
375 
376  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
377  G4hIonisation* hIoni = new G4hIonisation();
378  hIoni->SetStepFunction(0.2, 50*um);
379 
380  ph->RegisterProcess(pimsc, particle);
381  ph->RegisterProcess(hIoni, particle);
382  ph->RegisterProcess(pib, particle);
383  ph->RegisterProcess(pip, particle);
384  ph->RegisterProcess(new G4CoulombScattering(), particle);
385 
386  } else if (particleName == "kaon+" ||
387  particleName == "kaon-" ) {
388 
389  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
390  G4hIonisation* hIoni = new G4hIonisation();
391  hIoni->SetStepFunction(0.2, 50*um);
392 
393  ph->RegisterProcess(kmsc, particle);
394  ph->RegisterProcess(hIoni, particle);
395  ph->RegisterProcess(kb, particle);
396  ph->RegisterProcess(kp, particle);
397  ph->RegisterProcess(new G4CoulombScattering(), particle);
398 
399  } else if (particleName == "proton" ||
400  particleName == "anti_proton") {
401 
402  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
403  G4hIonisation* hIoni = new G4hIonisation();
404  hIoni->SetStepFunction(0.2, 50*um);
405 
406  ph->RegisterProcess(pmsc, particle);
407  ph->RegisterProcess(hIoni, particle);
408  ph->RegisterProcess(pb, particle);
409  ph->RegisterProcess(pp, particle);
410  ph->RegisterProcess(pnuc, particle);
411  ph->RegisterProcess(new G4CoulombScattering(), particle);
412 
413  } else if (particleName == "B+" ||
414  particleName == "B-" ||
415  particleName == "D+" ||
416  particleName == "D-" ||
417  particleName == "Ds+" ||
418  particleName == "Ds-" ||
419  particleName == "anti_He3" ||
420  particleName == "anti_alpha" ||
421  particleName == "anti_deuteron" ||
422  particleName == "anti_lambda_c+" ||
423  particleName == "anti_omega-" ||
424  particleName == "anti_sigma_c+" ||
425  particleName == "anti_sigma_c++" ||
426  particleName == "anti_sigma+" ||
427  particleName == "anti_sigma-" ||
428  particleName == "anti_triton" ||
429  particleName == "anti_xi_c+" ||
430  particleName == "anti_xi-" ||
431  particleName == "deuteron" ||
432  particleName == "lambda_c+" ||
433  particleName == "omega-" ||
434  particleName == "sigma_c+" ||
435  particleName == "sigma_c++" ||
436  particleName == "sigma+" ||
437  particleName == "sigma-" ||
438  particleName == "tau+" ||
439  particleName == "tau-" ||
440  particleName == "triton" ||
441  particleName == "xi_c+" ||
442  particleName == "xi-" ) {
443 
444  // Identical to G4EmStandardPhysics_option3
445 
446  ph->RegisterProcess(hmsc, particle);
447  ph->RegisterProcess(new G4hIonisation(), particle);
448  ph->RegisterProcess(pnuc, particle);
449  }
450  }
451 
452  // Nuclear stopping
453  pnuc->SetMaxKinEnergy(MeV);
454 
455  // Deexcitation
456  //
459 
460  G4EmModelActivator mact;
461  mact.ConstructProcess();
462 }
463 
464 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmLivermorePolarizedPhysics)
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
static const double MeV
Definition: G4SIunits.hh:211
void SetVerbose(G4int val)
static G4LossTableManager * Instance()
void SetLowestElectronEnergy(G4double val)
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
void SetStepFunction(G4double v1, G4double v2)
int G4int
Definition: G4Types.hh:78
void SetMaxEnergy(G4double val)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
void SetEmModel(G4VEmModel *, G4int index=1)
const G4String & GetParticleName() const
const G4String & GetPhysicsName() const
G4GLOB_DLL std::ostream G4cout
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
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 const double GeV
Definition: G4SIunits.hh:214
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
static const double eV
Definition: G4SIunits.hh:212
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetMinEnergy(G4double val)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4EmParameters * Instance()
void SetEmModel(G4VEmModel *, G4int index=1)
static const double um
Definition: G4SIunits.hh:112
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 const double TeV
Definition: G4SIunits.hh:215
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void SetMinKinEnergy(G4double e)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
static G4He3 * He3()
Definition: G4He3.cc:94
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetFluo(G4bool val)
void SetStepLimitType(G4MscStepLimitType val)