Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MuElecInelastic.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 //
27 // G4MuElecInelastic.cc, 2011/08/29 A.Valentin, M. Raine
28 //
29 // Based on the following publications
30 //
31 // - Inelastic cross-sections of low energy electrons in silicon
32 // for the simulation of heavy ion tracks with theGeant4-DNA toolkit,
33 // NSS Conf. Record 2010, pp. 80-85.
34 // - Geant4 physics processes for microdosimetry simulation:
35 // very low energy electromagnetic models for electrons in Si,
36 // NIM B, vol. 288, pp. 66 - 73, 2012.
37 // - Geant4 physics processes for microdosimetry simulation:
38 // very low energy electromagnetic models for protons and
39 // heavy ions in Si, NIM B, vol. 287, pp. 124 - 129, 2012.
40 //
41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42 
43 
44 #include "G4MuElecInelastic.hh"
45 #include "G4SystemOfUnits.hh"
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
49 using namespace std;
50 
52  G4ProcessType type):G4VEmProcess (processName, type),
53  isInitialised(false)
54 {
56 }
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
61 {}
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64 
66 {
67  return (&p == G4Electron::Electron() ||
68  &p == G4Proton::Proton() ||
69  (p.GetPDGCharge() != 0.0 && !p.IsShortLived() && p.GetParticleType() == "nucleus"));
70 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73 
75 {
76  if(!isInitialised)
77  {
78  isInitialised = true;
79  SetBuildTableFlag(false);
81 
82  if(name == "e-")
83  {
85  Model()->SetLowEnergyLimit(16.7*eV);
87 
88  AddEmModel(1, Model());
89  }
90 
91  else if(name == "proton")
92  {
94  Model()->SetLowEnergyLimit(50.*keV);
96 
97  AddEmModel(1, Model());
98  }
99 
100  else
101  {
103  Model()->SetLowEnergyLimit(50.*keV);
104  Model()->SetHighEnergyLimit(100.*GeV);
105 
106  AddEmModel(1, Model());
107  }
108  }
109 }
110 
111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112 
114 {
115  if (Model(2))
116  {
117  G4cout
118  << " Total cross sections computed from "
119  << Model(1)->GetName()
120  << " and "
121  << Model(2)->GetName()
122  << " models"
123  << G4endl;
124  }
125  else
126  {
127  G4cout
128  << " Total cross sections computed from "
129  << Model()->GetName()
130  << G4endl;
131  }
132 }
133 
134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......