Geant4  10.02.p01
G4empCrossSection.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: G4empCrossSection.cc 92765 2015-09-15 15:20:47Z gcosmo $
27 //
28 //
29 //
30 // History:
31 // -----------
32 // 29 Apr 2009 ALF 1st implementation
33 // 15 Mar 2011 ALF introduced the usage of G4AtomicShellEnumerator
34 // 09 Mar 2012 LP updated methods
35 //
36 
37 
38 #include "globals.hh"
39 #include "G4empCrossSection.hh"
40 #include "G4Proton.hh"
41 
42 
44  :G4VhShellCrossSection(nam),totalCS(0.0)
45 {
46  if (nam == "Empirical")
47  {
48  paulShellK = new G4PaulKxsModel();
50  flag=0;
51  }
52  else
53  {
54  G4cout << "G4empCrossSection::G4empCrossSection: "
55  << "ERROR in G4empCrossSection name; Paul+Orlic is selected."
56  << G4endl;
57  paulShellK = new G4PaulKxsModel();
59  flag=0;
60  }
61 
62 }
63 
65 {
66  delete paulShellK;
67  delete orlicShellLi;
68 }
69 
70 std::vector<G4double> G4empCrossSection::GetCrossSection(G4int Z,
71  G4double incidentEnergy,
72  G4double mass,
73  G4double,
74  const G4Material*)
75 {
76  std::vector<G4double> crossSections;
78 
79  crossSections.push_back( paulShellK->CalculateKCrossSection(Z, mass, incidentEnergy) );
80 
81  // this check should be done in the Orlic class, that can handle only protons;
82  // however this would lead up tp three checks of the mass, while here we have only one
83  // moreover, at the present time,this class handles explicitly Paul and Orlic models,
84  // so it can hadle the responsibility of this check too
85 
86  if (mass == aProton->GetPDGMass()) {
87 
88  if (flag==0)
89  {
90  crossSections.push_back( orlicShellLi->CalculateL1CrossSection(Z, incidentEnergy) );
91  crossSections.push_back( orlicShellLi->CalculateL2CrossSection(Z, incidentEnergy) );
92  crossSections.push_back( orlicShellLi->CalculateL3CrossSection(Z, incidentEnergy) );
93  }
94 
95  }
96 
97  else {
98  crossSections.push_back( 0. );
99  crossSections.push_back( 0. );
100  crossSections.push_back( 0. );
101  }
102  return crossSections;
103 
104 }
105 
107  G4double incidentEnergy,
108  G4double mass,
109  const G4Material*)
110 {
111  G4double res = 0.0;
113 
114  if(fKShell == shell) {
115  res = paulShellK->CalculateKCrossSection(Z, mass, incidentEnergy);
116  }
117  // this check should be done in the Orlic class, that can handle only protons;
118  // however this would lead up tp three checks of the mass, while here we have only one
119  // moreover, at the present time,this class handles explicitly Paul and Orlic models,
120  // so it can hadle the responsibility of this check too
121 
122  else if (mass == aProton->GetPDGMass()) {
123 
124  if(fL1Shell == shell) {
125  if (flag==0) res = orlicShellLi->CalculateL1CrossSection(Z, incidentEnergy);
126  }
127  else if(fL2Shell == shell) {
128  if (flag==0) res = orlicShellLi->CalculateL2CrossSection(Z, incidentEnergy);
129  }
130  else if(fL3Shell == shell) {
131  if (flag==0) res = orlicShellLi->CalculateL3CrossSection(Z, incidentEnergy);
132  }
133  }
134  return res;
135 }
136 
137 std::vector<G4double> G4empCrossSection::Probabilities(G4int Z,
138  G4double incidentEnergy,
139  G4double mass,
140  G4double deltaEnergy,
141  const G4Material* mat)
142 {
143 
144  std::vector<G4double> crossSections = GetCrossSection(Z, incidentEnergy, mass, deltaEnergy,mat);
145 
146  for (size_t i=0; i<crossSections.size(); i++ ) {
147 
148  if (totalCS) {
149  crossSections[i] = crossSections[i]/totalCS;
150  }
151 
152  }
153 
154  return crossSections;
155 
156 }
157 
158 
160 
161  totalCS = val;
162 
163 }
164 
165 
166 
std::vector< G4double > GetCrossSection(G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy, const G4Material *mat)
G4double CalculateKCrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)
G4OrlicLiXsModel * orlicShellLi
virtual ~G4empCrossSection()
int G4int
Definition: G4Types.hh:78
G4empCrossSection(const G4String &nam="")
G4GLOB_DLL std::ostream G4cout
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double CalculateL3CrossSection(G4int zTarget, G4double energyIncident)
G4double GetPDGMass() const
void SetTotalCS(G4double)
G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)
std::vector< G4double > Probabilities(G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy, const G4Material *mat)
#define G4endl
Definition: G4ios.hh:61
G4PaulKxsModel * paulShellK
double G4double
Definition: G4Types.hh:76
G4double CalculateL2CrossSection(G4int zTarget, G4double energyIncident)
G4double CalculateL1CrossSection(G4int zTarget, G4double energyIncident)
G4AtomicShellEnumerator