Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ecpssrFormFactorLixsModel.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// History:
27// -----------
28// 01 Oct 2011 A.M., S.I. - 1st implementation
29//
30// Class description
31// ----------------
32// Computation of K, L & M shell ECPSSR ionisation cross sections for protons and alphas
33// Based on the work of A. Taborda et al.
34// X-Ray Spectrom. 2011, 40, 127-134
35// ---------------------------------------------------------------------------------------
36
37#include <fstream>
38#include <iomanip>
39
41
42#include "globals.hh"
43#include "G4SystemOfUnits.hh"
44#include "G4ios.hh"
45
46#include "G4EMDataSet.hh"
47#include "G4LinInterpolation.hh"
48#include "G4Proton.hh"
49#include "G4Alpha.hh"
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
54{
55 interpolation = new G4LinInterpolation();
56
57 for (G4int i=11; i<93; i++)
58 {
59 protonL1DataSetMap[i] = new G4EMDataSet(i,interpolation);
60 protonL1DataSetMap[i]->LoadData("pixe/ecpssr/proton/l1-i01m001c01-");
61
62 protonL2DataSetMap[i] = new G4EMDataSet(i,interpolation);
63 protonL2DataSetMap[i]->LoadData("pixe/ecpssr/proton/l2-i01m001c01-");
64
65 protonL3DataSetMap[i] = new G4EMDataSet(i,interpolation);
66 protonL3DataSetMap[i]->LoadData("pixe/ecpssr/proton/l3-i01m001c01-");
67 }
68
69 for (G4int i=11; i<93; i++)
70 {
71 alphaL1DataSetMap[i] = new G4EMDataSet(i,interpolation);
72 alphaL1DataSetMap[i]->LoadData("pixe/ecpssr/alpha/l1-i02m004c02-");
73
74 alphaL2DataSetMap[i] = new G4EMDataSet(i,interpolation);
75 alphaL2DataSetMap[i]->LoadData("pixe/ecpssr/alpha/l2-i02m004c02-");
76
77 alphaL3DataSetMap[i] = new G4EMDataSet(i,interpolation);
78 alphaL3DataSetMap[i]->LoadData("pixe/ecpssr/alpha/l3-i02m004c02-");
79 }
80
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
86{
87 protonL1DataSetMap.clear();
88 alphaL1DataSetMap.clear();
89
90 protonL2DataSetMap.clear();
91 alphaL2DataSetMap.clear();
92
93 protonL3DataSetMap.clear();
94 alphaL3DataSetMap.clear();
95
96 delete interpolation;
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
102{
103 G4Proton* aProton = G4Proton::Proton();
104 G4Alpha* aAlpha = G4Alpha::Alpha();
105 G4double sigma = 0;
106
107 if (energyIncident > 0.1*MeV && energyIncident < 100.*MeV && zTarget < 93 && zTarget > 10) {
108
109 if (massIncident == aProton->GetPDGMass())
110 {
111 sigma = protonL1DataSetMap[zTarget]->FindValue(energyIncident/MeV);
112 if (sigma !=0 && energyIncident > protonL1DataSetMap[zTarget]->GetEnergies(0).back()*MeV) return 0.;
113 }
114 else if (massIncident == aAlpha->GetPDGMass())
115 {
116 sigma = alphaL1DataSetMap[zTarget]->FindValue(energyIncident/MeV);
117 if (sigma !=0 && energyIncident > alphaL1DataSetMap[zTarget]->GetEnergies(0).back()*MeV) return 0.;
118 }
119 else
120 {
121 sigma = 0.;
122 }
123 }
124
125 // sigma is in internal units: it has been converted from
126 // the input file in barns bt the EmDataset
127 return sigma;
128}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131
133{
134 G4Proton* aProton = G4Proton::Proton();
135 G4Alpha* aAlpha = G4Alpha::Alpha();
136 G4double sigma = 0;
137
138 if (energyIncident > 0.1*MeV && energyIncident < 100.*MeV && zTarget < 93 && zTarget > 10) {
139
140 if (massIncident == aProton->GetPDGMass())
141 {
142 sigma = protonL2DataSetMap[zTarget]->FindValue(energyIncident/MeV);
143 if (sigma !=0 && energyIncident > protonL2DataSetMap[zTarget]->GetEnergies(0).back()*MeV) return 0.;
144 }
145 else if (massIncident == aAlpha->GetPDGMass())
146 {
147 sigma = alphaL2DataSetMap[zTarget]->FindValue(energyIncident/MeV);
148 if (sigma !=0 && energyIncident > alphaL2DataSetMap[zTarget]->GetEnergies(0).back()*MeV) return 0.;
149 }
150 else
151 {
152 sigma = 0.;
153 }
154 }
155
156 // sigma is in internal units: it has been converted from
157 // the input file in barns bt the EmDataset
158 return sigma;
159}
160
161//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162
164{
165 G4Proton* aProton = G4Proton::Proton();
166 G4Alpha* aAlpha = G4Alpha::Alpha();
167 G4double sigma = 0;
168
169 if (energyIncident > 0.1*MeV && energyIncident < 100.*MeV && zTarget < 93 && zTarget > 10) {
170
171 if (massIncident == aProton->GetPDGMass())
172 {
173 sigma = protonL3DataSetMap[zTarget]->FindValue(energyIncident/MeV);
174 if (sigma !=0 && energyIncident > protonL3DataSetMap[zTarget]->GetEnergies(0).back()*MeV) return 0.;
175 }
176 else if (massIncident == aAlpha->GetPDGMass())
177 {
178 sigma = alphaL3DataSetMap[zTarget]->FindValue(energyIncident/MeV);
179 if (sigma !=0 && energyIncident > alphaL3DataSetMap[zTarget]->GetEnergies(0).back()*MeV) return 0.;
180 }
181 else
182 {
183 sigma = 0.;
184 }
185 }
186
187 // sigma is in internal units: it has been converted from
188 // the input file in barns bt the EmDataset
189 return sigma;
190}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4double CalculateL1CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)
G4double CalculateL2CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)
G4double CalculateL3CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)