Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4teoCrossSection.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 //
28 //
29 // History:
30 // -----------
31 // 21 Apr 2009 ALF 1st implementation
32 // 29 Apr 2009 ALF Updated Desing for Integration
33 // 15 Mar 2011 ALF introduced the usage of G4AtomicShellEnumerator
34 // 20 Oct 2011 ALF updated to take into account ECPSSR Form Factor
35 // 09 Mar 2012 LP update methods
36 // 09 Mar 2012 ALF update for M-shells Simulation
37 //
38 
39 #include "globals.hh"
40 #include "G4teoCrossSection.hh"
41 #include "G4Proton.hh"
42 //#include "G4Alpha.hh"
43 #include "G4ecpssrBaseKxsModel.hh"
44 #include "G4ecpssrBaseLixsModel.hh"
45 
49 
51  :G4VhShellCrossSection(nam),totalCS(0.0),ecpssrShellK(0),ecpssrShellLi(0), ecpssrShellMi(0)
52 {
53 
54  if (nam == "Analytical")
55  {
56  ecpssrShellK = new G4ecpssrBaseKxsModel();
57  ecpssrShellLi = new G4ecpssrBaseLixsModel();
58  }
59  else if (nam == "ECPSSR_FormFactor")
60  {
61  ecpssrShellK = new G4ecpssrFormFactorKxsModel();
62  ecpssrShellLi = new G4ecpssrFormFactorLixsModel();
63  ecpssrShellMi = new G4ecpssrFormFactorMixsModel();
64  }
65  else { G4cout << "ERROR" << G4endl;}
66 
67 
68 }
69 
71 {
72  delete ecpssrShellK;
73  delete ecpssrShellLi;
74  if (ecpssrShellMi) {delete ecpssrShellMi;}
75 }
76 
77 std::vector<G4double> G4teoCrossSection::GetCrossSection(G4int Z,
78  G4double incidentEnergy,
79  G4double mass,
80  G4double,
81  const G4Material*)
82 {
83  std::vector<G4double> crossSections;
84 
85  crossSections.push_back( ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy) );
86 
87  crossSections.push_back( ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy) );
88  crossSections.push_back( ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy) );
89  crossSections.push_back( ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy) );
90 
91  if (ecpssrShellMi) {
92 
93  crossSections.push_back( ecpssrShellMi->CalculateM1CrossSection(Z, mass, incidentEnergy) );
94  crossSections.push_back( ecpssrShellMi->CalculateM2CrossSection(Z, mass, incidentEnergy) );
95  crossSections.push_back( ecpssrShellMi->CalculateM3CrossSection(Z, mass, incidentEnergy) );
96  crossSections.push_back( ecpssrShellMi->CalculateM4CrossSection(Z, mass, incidentEnergy) );
97  crossSections.push_back( ecpssrShellMi->CalculateM5CrossSection(Z, mass, incidentEnergy) );
98 
99  }
100 
101 
102  return crossSections;
103 }
104 
106  G4double incidentEnergy,
107  G4double mass,
108  const G4Material*)
109 {
110  G4double res = 0.0;
111  if(shell > 3 && !ecpssrShellMi) {
112  return res;
113  }
114 
115  else if(shell > 8) {
116  return res;
117  }
118 
119  else if(fKShell == shell)
120  {
121  res = ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy);
122  }
123 
124  else if(fL1Shell == shell)
125  {
126  res = ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy);
127  }
128 
129  else if(fL2Shell == shell)
130  {
131  res = ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy);
132  }
133 
134  else if(fL3Shell == shell)
135  {
136  res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
137  }
138 
139  else if(fM1Shell == shell)
140  {
141  res = ecpssrShellMi->CalculateM1CrossSection(Z, mass, incidentEnergy);
142  }
143 
144  else if(fM2Shell == shell)
145  {
146  res = ecpssrShellMi->CalculateM2CrossSection(Z, mass, incidentEnergy);
147  }
148 
149  else if(fM3Shell == shell)
150  {
151  res = ecpssrShellMi->CalculateM3CrossSection(Z, mass, incidentEnergy);
152  }
153 
154  else if(fM4Shell == shell)
155  {
156  res = ecpssrShellMi->CalculateM4CrossSection(Z, mass, incidentEnergy);
157  }
158 
159  else if(fM5Shell == shell)
160  {
161  res = ecpssrShellMi->CalculateM5CrossSection(Z, mass, incidentEnergy);
162  }
163  return res;
164 }
165 
166 std::vector<G4double> G4teoCrossSection::Probabilities(G4int Z,
167  G4double incidentEnergy,
168  G4double mass,
169  G4double deltaEnergy,
170  const G4Material*)
171 {
172  std::vector<G4double> crossSections =
173  GetCrossSection(Z, incidentEnergy, mass, deltaEnergy);
174 
175  for (size_t i=0; i<crossSections.size(); i++ ) {
176 
177  if (totalCS) {
178  crossSections[i] = crossSections[i]/totalCS;
179  }
180 
181  }
182  return crossSections;
183 }
184 
186 
187  totalCS = val;
188  // G4cout << "totalXS set to: " << val / barn << " barns" << G4endl;
189 }
190 
191 
192 
G4teoCrossSection(const G4String &name)
virtual G4double CalculateM1CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
void SetTotalCS(G4double)
virtual G4double CalculateM3CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateCrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
std::vector< G4double > GetCrossSection(G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=0)
virtual G4double CalculateM2CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateL2CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateM4CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateM5CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)
virtual G4double CalculateL1CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
#define G4endl
Definition: G4ios.hh:61
virtual G4double CalculateL3CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
double G4double
Definition: G4Types.hh:76
std::vector< G4double > Probabilities(G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=0)
virtual ~G4teoCrossSection()
G4AtomicShellEnumerator