Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4NeutrinoElectronNcXsc.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
30#include "G4SystemOfUnits.hh"
31#include "G4DynamicParticle.hh"
32#include "G4ParticleTable.hh"
33#include "G4IonTable.hh"
34#include "G4HadTmpUtil.hh"
35#include "G4Proton.hh"
36#include "G4NistManager.hh"
37
38using namespace std;
39using namespace CLHEP;
40
42 : G4VCrossSectionDataSet("NuElectronNcXsc")
43{
44 // PDG2016: Gf=1.1663787(6)e-5*(hc)^3/GeV^2
45 // fCofXsc = Gf*Gf*MeC2*2/pi
46
47 fCofXsc = 1.36044e-22;
48 fCofXsc *= hbarc*hbarc*electron_mass_c2;
49 fCofXsc /= halfpi;
50
51 // G4cout<<"hbarc = "<<hbarc/MeV/fermi<<" MeV*fermi"<<G4endl;
52
53 // PDG2016: sin^2 theta Weinberg
54
55 fSin2tW = 0.23129; // 0.2312;
56
57 fCutEnergy = 0.; // default value
58 fBiasingFactor = 1.;
59}
60
62{}
63
64//////////////////////////////////////////////////////
65
66G4bool
68{
69 G4bool result = false;
70 G4String pName = aPart->GetDefinition()->GetParticleName();
71 G4double minEnergy = 0., energy = aPart->GetTotalEnergy();
72 // Z *= 1;
73 if( fCutEnergy > 0. ) // min detected recoil electron energy
74 {
75 minEnergy = 0.5*(fCutEnergy+sqrt(fCutEnergy*(fCutEnergy+2.*electron_mass_c2)));
76 }
77 if( ( pName == "nu_e" || pName == "anti_nu_e" ||
78 pName == "nu_mu" || pName == "anti_nu_mu" ||
79 pName == "nu_tau" || pName == "anti_nu_tau" ) &&
80 energy > minEnergy )
81 {
82 result = true;
83 }
84 return result;
85}
86
87////////////////////////////////////////////////////
88
91 const G4Material*)
92{
93 G4double result = 0., cofL, cofR, cofL2, cofR2, cofLR;
94
95 G4double energy = aPart->GetTotalEnergy();
96 G4String pName = aPart->GetDefinition()->GetParticleName();
97
98 if( pName == "nu_e")
99 {
100 cofL = 0.5 + fSin2tW;
101 cofR = fSin2tW;
102 }
103 else if( pName == "anti_nu_e")
104 {
105 cofL = fSin2tW;
106 cofR = 0.5 + fSin2tW;
107 }
108 else if( pName == "nu_mu")
109 {
110 cofL = -0.5 + fSin2tW;
111 cofR = fSin2tW;
112 }
113 else if( pName == "anti_nu_mu")
114 {
115 cofL = fSin2tW;
116 cofR = -0.5 + fSin2tW;
117 }
118 else if( pName == "nu_tau") // vmg: nu_tau as nu_mu ???
119 {
120 cofL = -0.5 + fSin2tW;
121 cofR = fSin2tW;
122 }
123 else if( pName == "anti_nu_tau")
124 {
125 cofL = fSin2tW;
126 cofR = -0.5 + fSin2tW;
127 }
128 else
129 {
130 return result;
131 }
132 // if( energy <= electron_mass_c2 ) return result;
133
134 cofL2 = cofL*cofL;
135 cofR2 = cofR*cofR;
136 cofLR = cofL*cofR;
137
138 if( fCutEnergy > 0. )
139 {
140 G4double tM = 2.*energy*energy/(electron_mass_c2 + 2.*energy);
141 G4double tM2 = tM*tM;
142 G4double tM3 = tM*tM2;
143 G4double tC = fCutEnergy;
144 G4double tC2 = tC*tC;
145 G4double tC3 = tC*tC2;
146
147 result = (cofL2+cofR2)*(tM-tC);
148 result -= (cofR2+cofLR*0.5*electron_mass_c2/energy)*(tM2-tC2)/energy;
149 result += cofR2*(tM3-tC3)/energy/energy/3.;
150 }
151 else
152 {
153 G4double rtM = 2.*energy/(electron_mass_c2 + 2.*energy);
154 G4double rtM2 = rtM*rtM;
155 G4double rtM3 = rtM*rtM2;
156
157 result = (cofL2+cofR2)*rtM*energy;
158 result -= (cofR2*energy+cofLR*0.5*electron_mass_c2)*rtM2;
159 result += cofR2*rtM3*energy/3.;
160 }
161 // result = cofL*cofL + cofR*cofR/3.;
162 // G4cout<<"cofL2 + cofR2/3. = "<<result<<G4endl;
163 // result -= 0.5*cofL*cofR*electron_mass_c2/energy;
164
165 G4double aa = 1.;
166 G4double bb = 1.7;
167 G4double gw = 2.141*GeV;
168 G4double dd = 5000.;
169 G4double mw = 80.385*GeV;
170 G4double mz = 91.1876*GeV;
171
172 G4double emass = electron_mass_c2;
173 G4double totS = 2.*energy*emass + emass*emass;
174
175 if( energy > 50.*GeV )
176 {
177 result *= bb;
178 result /= 1.+ aa*totS/mz/mz;
179
180 if( pName == "anti_nu_e")
181 {
182 result *= 1. + dd*gw*gw*totS/( (totS-mw*mw)*(totS-mw*mw)+gw*gw*mw*mw );
183 }
184 }
185
186 result *= fCofXsc; //*energy;
187
188 result *= ZZ; // incoherent sum over all element electrons
189
190 result *= fBiasingFactor;
191
192 return result;
193}
194
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
const G4String & GetParticleName() const
Definition: DoubConv.h:17