Geant4  10.01.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 86233 2014-11-07 17:21:03Z 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 
118 // factory
120 //
122 
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125 
127  : G4VPhysicsConstructor("G4EmLivermorePolarizedPhysics"), verbose(ver)
128 {
130  param->SetVerbose(verbose);
131  param->SetMinEnergy(100*eV);
132  param->SetMaxEnergy(10*TeV);
133  param->SetNumberOfBinsPerDecade(20);
136 }
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139 
141  : G4VPhysicsConstructor("G4EmLivermorePolarizedPhysics"), verbose(ver)
142 {
144  param->SetVerbose(verbose);
145  param->SetMinEnergy(100*eV);
146  param->SetMaxEnergy(10*TeV);
147  param->SetNumberOfBinsPerDecade(20);
150 }
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153 
155 {}
156 
157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
158 
160 {
161 // gamma
162  G4Gamma::Gamma();
163 
164 // leptons
169 
170 // mesons
175 
176 // baryons
179 
180 // ions
183  G4He3::He3();
184  G4Alpha::Alpha();
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189 
191 {
192  if(verbose > 1) {
193  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
194  }
196 
197  // muon & hadron bremsstrahlung and pair production
206 
207  // muon & hadron multiple scattering
209  mumsc->AddEmModel(0, new G4WentzelVIModel());
211  pimsc->AddEmModel(0, new G4WentzelVIModel());
213  kmsc->AddEmModel(0, new G4WentzelVIModel());
215  pmsc->AddEmModel(0, new G4WentzelVIModel());
216  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
217 
218  // high energy limit for e+- scattering models
219  G4double highEnergyLimit = 100*MeV;
220 
221  // nuclear stopping
222  G4NuclearStopping* pnuc = new G4NuclearStopping();
223 
224  // Add Livermore EM Processes
225  aParticleIterator->reset();
226 
227  while( (*aParticleIterator)() ){
228 
229  G4ParticleDefinition* particle = aParticleIterator->value();
230  G4String particleName = particle->GetParticleName();
231 
232  //Applicability range for Livermore models
233  //for higher energies, the Standard models are used
234  G4double LivermoreHighEnergyLimit = GeV;
235 
236  if (particleName == "gamma") {
237 
238  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
240  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
241  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
242  ph->RegisterProcess(thePhotoElectricEffect, particle);
243 
244  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
245  G4LivermorePolarizedComptonModel* theLivermoreComptonModel = new G4LivermorePolarizedComptonModel();
246  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
247  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
248  ph->RegisterProcess(theComptonScattering, particle);
249 
250  G4GammaConversion* theGammaConversion = new G4GammaConversion();
252  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
253  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
254  ph->RegisterProcess(theGammaConversion, particle);
255 
256  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
258  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
259  theRayleigh->AddEmModel(0, theRayleighModel);
260  ph->RegisterProcess(theRayleigh, particle);
261 
262  } else if (particleName == "e-") {
263 
264  // multiple scattering
267  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
268  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
269  msc1->SetHighEnergyLimit(highEnergyLimit);
270  msc2->SetLowEnergyLimit(highEnergyLimit);
271  msc->AddEmModel(0, msc1);
272  msc->AddEmModel(0, msc2);
273 
276  ss->SetEmModel(ssm, 1);
277  ss->SetMinKinEnergy(highEnergyLimit);
278  ssm->SetLowEnergyLimit(highEnergyLimit);
279  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
280  ph->RegisterProcess(msc, particle);
281  ph->RegisterProcess(ss, particle);
282 
283  // Ionisation
284  G4eIonisation* eIoni = new G4eIonisation();
285  G4LivermoreIonisationModel* theIoniLivermore = new
287  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
288  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
289  eIoni->SetStepFunction(0.2, 100*um); //
290  ph->RegisterProcess(eIoni, particle);
291 
292  // Bremsstrahlung from standard
293  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
294  ph->RegisterProcess(eBrem, particle);
295 
296  } else if (particleName == "e+") {
297 
298  // multiple scattering
301  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
302  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
303  msc1->SetHighEnergyLimit(highEnergyLimit);
304  msc2->SetLowEnergyLimit(highEnergyLimit);
305  msc->AddEmModel(0, msc1);
306  msc->AddEmModel(0, msc2);
307 
310  ss->SetEmModel(ssm, 1);
311  ss->SetMinKinEnergy(highEnergyLimit);
312  ssm->SetLowEnergyLimit(highEnergyLimit);
313  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
314 
315  // Ionisation
316  G4eIonisation* eIoni = new G4eIonisation();
317  eIoni->SetStepFunction(0.2, 100*um);
318 
319  ph->RegisterProcess(msc, particle);
320  ph->RegisterProcess(eIoni, particle);
321  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
322  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
323  ph->RegisterProcess(ss, particle);
324 
325  } else if (particleName == "mu+" ||
326  particleName == "mu-" ) {
327 
328  G4MuIonisation* muIoni = new G4MuIonisation();
329  muIoni->SetStepFunction(0.2, 50*um);
330 
331  ph->RegisterProcess(mumsc, particle);
332  ph->RegisterProcess(muIoni, particle);
333  ph->RegisterProcess(mub, particle);
334  ph->RegisterProcess(mup, particle);
335  ph->RegisterProcess(new G4CoulombScattering(), particle);
336 
337  } else if (particleName == "alpha" ||
338  particleName == "He3" ) {
339 
340  // Identical to G4EmStandardPhysics_option3
341 
343  G4ionIonisation* ionIoni = new G4ionIonisation();
344  ionIoni->SetStepFunction(0.1, 10*um);
345 
346  ph->RegisterProcess(msc, particle);
347  ph->RegisterProcess(ionIoni, particle);
348  ph->RegisterProcess(pnuc, particle);
349 
350  } else if (particleName == "GenericIon") {
351 
352  // Identical to G4EmStandardPhysics_option3
353 
354  G4ionIonisation* ionIoni = new G4ionIonisation();
355  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
356  ionIoni->SetStepFunction(0.1, 1*um);
357 
358  ph->RegisterProcess(hmsc, particle);
359  ph->RegisterProcess(ionIoni, particle);
360  ph->RegisterProcess(pnuc, particle);
361 
362  } else if (particleName == "pi+" ||
363  particleName == "pi-" ) {
364 
365  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
366  G4hIonisation* hIoni = new G4hIonisation();
367  hIoni->SetStepFunction(0.2, 50*um);
368 
369  ph->RegisterProcess(pimsc, particle);
370  ph->RegisterProcess(hIoni, particle);
371  ph->RegisterProcess(pib, particle);
372  ph->RegisterProcess(pip, particle);
373  ph->RegisterProcess(new G4CoulombScattering(), particle);
374 
375  } else if (particleName == "kaon+" ||
376  particleName == "kaon-" ) {
377 
378  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
379  G4hIonisation* hIoni = new G4hIonisation();
380  hIoni->SetStepFunction(0.2, 50*um);
381 
382  ph->RegisterProcess(kmsc, particle);
383  ph->RegisterProcess(hIoni, particle);
384  ph->RegisterProcess(kb, particle);
385  ph->RegisterProcess(kp, particle);
386  ph->RegisterProcess(new G4CoulombScattering(), particle);
387 
388  } else if (particleName == "proton" ||
389  particleName == "anti_proton") {
390 
391  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
392  G4hIonisation* hIoni = new G4hIonisation();
393  hIoni->SetStepFunction(0.2, 50*um);
394 
395  ph->RegisterProcess(pmsc, particle);
396  ph->RegisterProcess(hIoni, particle);
397  ph->RegisterProcess(pb, particle);
398  ph->RegisterProcess(pp, particle);
399  ph->RegisterProcess(pnuc, particle);
400  ph->RegisterProcess(new G4CoulombScattering(), particle);
401 
402  } else if (particleName == "B+" ||
403  particleName == "B-" ||
404  particleName == "D+" ||
405  particleName == "D-" ||
406  particleName == "Ds+" ||
407  particleName == "Ds-" ||
408  particleName == "anti_He3" ||
409  particleName == "anti_alpha" ||
410  particleName == "anti_deuteron" ||
411  particleName == "anti_lambda_c+" ||
412  particleName == "anti_omega-" ||
413  particleName == "anti_sigma_c+" ||
414  particleName == "anti_sigma_c++" ||
415  particleName == "anti_sigma+" ||
416  particleName == "anti_sigma-" ||
417  particleName == "anti_triton" ||
418  particleName == "anti_xi_c+" ||
419  particleName == "anti_xi-" ||
420  particleName == "deuteron" ||
421  particleName == "lambda_c+" ||
422  particleName == "omega-" ||
423  particleName == "sigma_c+" ||
424  particleName == "sigma_c++" ||
425  particleName == "sigma+" ||
426  particleName == "sigma-" ||
427  particleName == "tau+" ||
428  particleName == "tau-" ||
429  particleName == "triton" ||
430  particleName == "xi_c+" ||
431  particleName == "xi-" ) {
432 
433  // Identical to G4EmStandardPhysics_option3
434 
435  ph->RegisterProcess(hmsc, particle);
436  ph->RegisterProcess(new G4hIonisation(), particle);
437  ph->RegisterProcess(pnuc, particle);
438  }
439  }
440 
441  // Nuclear stopping
442  pnuc->SetMaxKinEnergy(MeV);
443 
444  // Deexcitation
445  //
448  de->SetFluo(true);
449 }
450 
451 //....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:193
void SetVerbose(G4int val)
static G4LossTableManager * Instance()
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
void SetStepFunction(G4double v1, G4double v2)
int G4int
Definition: G4Types.hh:78
void SetMaxEnergy(G4double val)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:707
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)
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:196
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 SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:728
static const double eV
Definition: G4SIunits.hh:194
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)
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:197
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:714
static G4He3 * He3()
Definition: G4He3.cc:94
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetStepLimitType(G4MscStepLimitType val)