Geant4-11
G4ElNucleusSFcs.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#include <iostream>
29
30#include "G4SystemOfUnits.hh"
32#include "G4HadTmpUtil.hh"
33#include "G4ElNucleusSFcs.hh"
35
36// factory
38//
40
41
43{
45}
46
48{
49 if( fCHIPScs != nullptr ) delete fCHIPScs;
50}
51
52void
53G4ElNucleusSFcs::CrossSectionDescription( std::ostream& outFile) const
54{
55 outFile << "G4ElNucleusSFcs provides the inelastic\n"
56 << "cross section for e- and e+ interactions with nuclei according to the structure function approach."
57 << "all energies.\n";
58}
59
61{
62
63 G4double eKin = aPart->GetKineticEnergy();
64 G4double thTkin = ThresholdEnergy();
65
66 return ( Z > 0 && Z < 120 && eKin > thTkin );
67}
68
69
71 const G4Element* , const G4Material* )
72{
73 G4double xsc(0.), ratio(1.);
74
75 xsc = fCHIPScs->GetElementCrossSection( aPart, ZZ, nullptr ); //, const G4Material*);
76
77 ratio = GetRatio(ZZ,AA);
78
79 if( ratio > 0. ) xsc /= ratio;
80
81 return xsc;
82}
83
85//
86// Calculate ratio CHIPS-CS/SF-CS, averaged between 200 MeV and 100 GeV
87
89{
90 G4double ratio(1.);
91
92 if ( Z == 1 && A == 1 ) return 1.51;
93 else if ( Z == 1 && A == 2 ) return 0.33;
94 else if ( Z == 1 && A == 3 ) return 0.27;
95 else if ( Z == 2 && A == 4 ) return 1.81;
96 else if ( Z == 6 && A == 12 ) return 2.26;
97 else if ( Z == 7 && A == 14 ) return 2.47;
98 else if ( Z == 8 && A == 16 ) return 2.61;
99 else if ( Z == 13 && A == 27 ) return 2.57;
100 else if ( Z == 14 && A == 28 ) return 2.49;
101 else if ( Z == 18 && A == 40 ) return 2.72;
102 else if ( Z == 22 && A == 48 ) return 2.71;
103 else if ( Z == 26 && A == 56 ) return 2.79;
104 else if ( Z == 29 && A == 64 ) return 2.78;
105 else if ( Z == 32 && A == 73 ) return 2.87;
106 else if ( Z == 42 && A == 96 ) return 3.02;
107 else if ( Z == 46 && A == 106 ) return 3.02;
108 else if ( Z == 47 && A == 108 ) return 2.99;
109 else if ( Z == 48 && A == 112 ) return 3.00;
110 else if ( Z == 74 && A == 184 ) return 3.44;
111 else if ( Z == 79 && A == 200 ) return 3.49;
112 else if ( Z == 82 && A == 207 ) return 3.48;
113 else if ( Z == 92 && A == 238 ) return 3.88;
114 else
115 {
116 G4int it(0), iMax(19);
117 G4double zz = G4double(Z);
118
119 for ( it = 0; it < iMax; ++it ) if ( zz <= fZZ[it] ) break;
120
121 if ( it == 0 ) return fRR[0];
122 else if( it == iMax ) return fRR[iMax-1];
123 else
124 {
125 G4double x1 = fZZ[it-1];
126 G4double x2 = fZZ[it];
127 G4double y1 = fRR[it-1];
128 G4double y2 = fRR[it];
129
130 if( x1 >= x2 ) return fRR[it];
131 else
132 {
133 G4double angle = (y2-y1)/(x2-x1);
134 ratio = y1 + ( zz - x1 )*angle;
135 }
136 }
137 }
138 return ratio;
139}
140
142
144{
145 G4double thTkin = 134.9766*CLHEP::MeV; // Pi0 threshold for the nucleon
146
147 thTkin *= 1.3; // nucleon motion
148
149 return thTkin;
150}
151
153
155 {
156 2., 6., 7., 8., 13., 14., 18., 22., 26., 29.,
157 32., 42., 46., 47., 48., 74., 79., 82., 92.
158 };
159
161
163 {
164 1.81, 2.26, 2.47, 2.61, 2.57, 2.49, 2.72, 2.71, 2.79, 2.78,
165 2.87, 3.02, 3.02, 2.99, 3., 3.44, 3.49, 3.48, 2.88
166 };
G4_DECLARE_XS_FACTORY(G4ElNucleusSFcs)
static const G4double angle[DIMMOTT]
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4double GetKineticEnergy() const
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual void CrossSectionDescription(std::ostream &) const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
G4double ThresholdEnergy()
virtual ~G4ElNucleusSFcs()
static const G4double fRR[19]
G4double GetRatio(G4int Z, G4int A)
G4ElectroNuclearCrossSection * fCHIPScs
static const G4double fZZ[19]
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat)
static constexpr double MeV
static const G4double AA[5]
Definition: paraMaker.cc:44