Geant4  10.03.p03
 All Classes 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 "G4ecpssrBaseKxsModel.hh"
43 #include "G4ecpssrBaseLixsModel.hh"
44 
48 
50  :G4VhShellCrossSection(nam),totalCS(0.0)
51 {
52  ecpssrShellMi = nullptr;
53  if (nam == "ECPSSR_Analytical")
54  {
55  ecpssrShellK = new G4ecpssrBaseKxsModel();
56  ecpssrShellLi = new G4ecpssrBaseLixsModel();
57  }
58  else if (nam == "ECPSSR_FormFactor")
59  {
60  ecpssrShellK = new G4ecpssrFormFactorKxsModel();
61  ecpssrShellLi = new G4ecpssrFormFactorLixsModel();
62  ecpssrShellMi = new G4ecpssrFormFactorMixsModel();
63  }
64  else
65  {
66  G4cout << "G4teoCrossSection::G4teoCrossSection: ERROR "
67  << " in cross section name ECPSSR_Analytical is used"
68  << G4endl;
69  ecpssrShellK = new G4ecpssrBaseKxsModel();
70  ecpssrShellLi = new G4ecpssrBaseLixsModel();
71  }
72 }
73 
75 {
76  delete ecpssrShellK;
77  delete ecpssrShellLi;
78  delete ecpssrShellMi;
79 }
80 
82  G4double incidentEnergy,
83  G4double mass,
84  G4double,
85  const G4Material*)
86 {
87  std::vector<G4double> crossSections;
88 
89  crossSections.push_back( ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy) );
90 
91  crossSections.push_back( ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy) );
92  crossSections.push_back( ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy) );
93  crossSections.push_back( ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy) );
94 
95  if (ecpssrShellMi) {
96 
97  crossSections.push_back( ecpssrShellMi->CalculateM1CrossSection(Z, mass, incidentEnergy) );
98  crossSections.push_back( ecpssrShellMi->CalculateM2CrossSection(Z, mass, incidentEnergy) );
99  crossSections.push_back( ecpssrShellMi->CalculateM3CrossSection(Z, mass, incidentEnergy) );
100  crossSections.push_back( ecpssrShellMi->CalculateM4CrossSection(Z, mass, incidentEnergy) );
101  crossSections.push_back( ecpssrShellMi->CalculateM5CrossSection(Z, mass, incidentEnergy) );
102 
103  }
104 
105 
106  return crossSections;
107 }
108 
110  G4double incidentEnergy,
111  G4double mass,
112  const G4Material*)
113 {
114  G4double res = 0.0;
115  if(shell > 3 && !ecpssrShellMi) {
116  return res;
117  }
118 
119  else if(shell > 8) {
120  return res;
121  }
122 
123  else if(fKShell == shell)
124  {
125  res = ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy);
126  }
127 
128  else if(fL1Shell == shell)
129  {
130  res = ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy);
131  }
132 
133  else if(fL2Shell == shell)
134  {
135  res = ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy);
136  }
137 
138  else if(fL3Shell == shell)
139  {
140  res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
141  }
142 
143  else if(fM1Shell == shell)
144  {
145  res = ecpssrShellMi->CalculateM1CrossSection(Z, mass, incidentEnergy);
146  }
147 
148  else if(fM2Shell == shell)
149  {
150  res = ecpssrShellMi->CalculateM2CrossSection(Z, mass, incidentEnergy);
151  }
152 
153  else if(fM3Shell == shell)
154  {
155  res = ecpssrShellMi->CalculateM3CrossSection(Z, mass, incidentEnergy);
156  }
157 
158  else if(fM4Shell == shell)
159  {
160  res = ecpssrShellMi->CalculateM4CrossSection(Z, mass, incidentEnergy);
161  }
162 
163  else if(fM5Shell == shell)
164  {
165  res = ecpssrShellMi->CalculateM5CrossSection(Z, mass, incidentEnergy);
166  }
167  return res;
168 }
169 
171  G4double incidentEnergy,
172  G4double mass,
173  G4double deltaEnergy,
174  const G4Material*)
175 {
176  std::vector<G4double> crossSections =
177  GetCrossSection(Z, incidentEnergy, mass, deltaEnergy);
178 
179  for (size_t i=0; i<crossSections.size(); i++ ) {
180 
181  if (totalCS) {
182  crossSections[i] = crossSections[i]/totalCS;
183  }
184 
185  }
186  return crossSections;
187 }
188 
190 
191  totalCS = val;
192  // G4cout << "totalXS set to: " << val / barn << " barns" << G4endl;
193 }
194 
195 
196 
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