Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 68750 2013-04-05 10:19:04Z 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 "G4EmProcessOptions.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 {
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134 
136  : G4VPhysicsConstructor("G4EmLivermorePolarizedPhysics"), verbose(ver)
137 {
140 }
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143 
145 {}
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
150 {
151 // gamma
152  G4Gamma::Gamma();
153 
154 // leptons
159 
160 // mesons
165 
166 // baryons
169 
170 // ions
173  G4He3::He3();
174  G4Alpha::Alpha();
176 }
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179 
181 {
182  if(verbose > 1) {
183  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
184  }
186 
187  // muon & hadron bremsstrahlung and pair production
196 
197  // muon & hadron multiple scattering
199  mumsc->AddEmModel(0, new G4WentzelVIModel());
201  pimsc->AddEmModel(0, new G4WentzelVIModel());
203  kmsc->AddEmModel(0, new G4WentzelVIModel());
205  pmsc->AddEmModel(0, new G4WentzelVIModel());
206  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
207 
208  // high energy limit for e+- scattering models
209  G4double highEnergyLimit = 100*MeV;
210 
211  // nuclear stopping
212  G4NuclearStopping* ionnuc = new G4NuclearStopping();
213  G4NuclearStopping* pnuc = new G4NuclearStopping();
214 
215  // Add Livermore EM Processes
216  aParticleIterator->reset();
217 
218  while( (*aParticleIterator)() ){
219 
220  G4ParticleDefinition* particle = aParticleIterator->value();
221  G4String particleName = particle->GetParticleName();
222 
223  //Applicability range for Livermore models
224  //for higher energies, the Standard models are used
225  G4double LivermoreHighEnergyLimit = GeV;
226 
227  if (particleName == "gamma") {
228 
229  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
231  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
232  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
233  ph->RegisterProcess(thePhotoElectricEffect, particle);
234 
235  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
236  G4LivermorePolarizedComptonModel* theLivermoreComptonModel = new G4LivermorePolarizedComptonModel();
237  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
238  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
239  ph->RegisterProcess(theComptonScattering, particle);
240 
241  G4GammaConversion* theGammaConversion = new G4GammaConversion();
243  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
244  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
245  ph->RegisterProcess(theGammaConversion, particle);
246 
247  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
249  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
250  theRayleigh->AddEmModel(0, theRayleighModel);
251  ph->RegisterProcess(theRayleigh, particle);
252 
253  } else if (particleName == "e-") {
254 
255  // multiple scattering
258  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
259  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
260  msc1->SetHighEnergyLimit(highEnergyLimit);
261  msc2->SetLowEnergyLimit(highEnergyLimit);
262  msc->AddEmModel(0, msc1);
263  msc->AddEmModel(0, msc2);
264 
267  ss->SetEmModel(ssm, 1);
268  ss->SetMinKinEnergy(highEnergyLimit);
269  ssm->SetLowEnergyLimit(highEnergyLimit);
270  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
271  ph->RegisterProcess(msc, particle);
272  ph->RegisterProcess(ss, particle);
273 
274  // Ionisation
275  G4eIonisation* eIoni = new G4eIonisation();
276  G4LivermoreIonisationModel* theIoniLivermore = new
278  theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
279  eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
280  eIoni->SetStepFunction(0.2, 100*um); //
281  ph->RegisterProcess(eIoni, particle);
282 
283  // Bremsstrahlung from standard
284  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
285  ph->RegisterProcess(eBrem, particle);
286 
287  } else if (particleName == "e+") {
288 
289  // multiple scattering
292  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
293  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
294  msc1->SetHighEnergyLimit(highEnergyLimit);
295  msc2->SetLowEnergyLimit(highEnergyLimit);
296  msc->AddEmModel(0, msc1);
297  msc->AddEmModel(0, msc2);
298 
301  ss->SetEmModel(ssm, 1);
302  ss->SetMinKinEnergy(highEnergyLimit);
303  ssm->SetLowEnergyLimit(highEnergyLimit);
304  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
305 
306  // Ionisation
307  G4eIonisation* eIoni = new G4eIonisation();
308  eIoni->SetStepFunction(0.2, 100*um);
309 
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 
331  // Identical to G4EmStandardPhysics_option3
332 
334  G4ionIonisation* ionIoni = new G4ionIonisation();
335  ionIoni->SetStepFunction(0.1, 10*um);
336 
337  ph->RegisterProcess(msc, particle);
338  ph->RegisterProcess(ionIoni, particle);
339  ph->RegisterProcess(ionnuc, particle);
340 
341  } else if (particleName == "GenericIon") {
342 
343  // Identical to G4EmStandardPhysics_option3
344 
345  G4ionIonisation* ionIoni = new G4ionIonisation();
346  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
347  ionIoni->SetStepFunction(0.1, 1*um);
348 
349  ph->RegisterProcess(hmsc, particle);
350  ph->RegisterProcess(ionIoni, particle);
351  ph->RegisterProcess(ionnuc, particle);
352 
353  } else if (particleName == "pi+" ||
354  particleName == "pi-" ) {
355 
356  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
357  G4hIonisation* hIoni = new G4hIonisation();
358  hIoni->SetStepFunction(0.2, 50*um);
359 
360  ph->RegisterProcess(pimsc, particle);
361  ph->RegisterProcess(hIoni, particle);
362  ph->RegisterProcess(pib, particle);
363  ph->RegisterProcess(pip, particle);
364 
365  } else if (particleName == "kaon+" ||
366  particleName == "kaon-" ) {
367 
368  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
369  G4hIonisation* hIoni = new G4hIonisation();
370  hIoni->SetStepFunction(0.2, 50*um);
371 
372  ph->RegisterProcess(kmsc, particle);
373  ph->RegisterProcess(hIoni, particle);
374  ph->RegisterProcess(kb, particle);
375  ph->RegisterProcess(kp, particle);
376 
377  } else if (particleName == "proton" ||
378  particleName == "anti_proton") {
379 
380  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
381  G4hIonisation* hIoni = new G4hIonisation();
382  hIoni->SetStepFunction(0.2, 50*um);
383 
384  ph->RegisterProcess(pmsc, particle);
385  ph->RegisterProcess(hIoni, particle);
386  ph->RegisterProcess(pb, particle);
387  ph->RegisterProcess(pp, particle);
388  ph->RegisterProcess(pnuc, particle);
389 
390  } else if (particleName == "B+" ||
391  particleName == "B-" ||
392  particleName == "D+" ||
393  particleName == "D-" ||
394  particleName == "Ds+" ||
395  particleName == "Ds-" ||
396  particleName == "anti_He3" ||
397  particleName == "anti_alpha" ||
398  particleName == "anti_deuteron" ||
399  particleName == "anti_lambda_c+" ||
400  particleName == "anti_omega-" ||
401  particleName == "anti_sigma_c+" ||
402  particleName == "anti_sigma_c++" ||
403  particleName == "anti_sigma+" ||
404  particleName == "anti_sigma-" ||
405  particleName == "anti_triton" ||
406  particleName == "anti_xi_c+" ||
407  particleName == "anti_xi-" ||
408  particleName == "deuteron" ||
409  particleName == "lambda_c+" ||
410  particleName == "omega-" ||
411  particleName == "sigma_c+" ||
412  particleName == "sigma_c++" ||
413  particleName == "sigma+" ||
414  particleName == "sigma-" ||
415  particleName == "tau+" ||
416  particleName == "tau-" ||
417  particleName == "triton" ||
418  particleName == "xi_c+" ||
419  particleName == "xi-" ) {
420 
421  // Identical to G4EmStandardPhysics_option3
422 
423  ph->RegisterProcess(hmsc, particle);
424  ph->RegisterProcess(new G4hIonisation(), particle);
425  ph->RegisterProcess(pnuc, particle);
426  }
427  }
428 
429  // Em options
430  //
431  G4EmProcessOptions opt;
432  opt.SetVerbose(verbose);
433 
434  // Multiple Coulomb scattering
435  //
436  opt.SetPolarAngleLimit(CLHEP::pi);
437 
438  // Physics tables
439  //
440 
441  opt.SetMinEnergy(100*eV);
442  opt.SetMaxEnergy(10*TeV);
443  opt.SetDEDXBinning(220);
444  opt.SetLambdaBinning(220);
445 
446  // Nuclear stopping
447  pnuc->SetMaxKinEnergy(MeV);
448 
449  // Ionization
450  //
451  //opt.SetSubCutoff(true);
452 
453  // Deexcitation
454  //
457  de->SetFluo(true);
458 }
459 
460 //....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 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)
void SetVerbose(G4int val, const G4String &name="all", G4bool worker=false)
void SetPolarAngleLimit(G4double val)