Geant4-11
G4MuNeutrinoNucleusTotXsc.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 "G4NistManager.hh"
36#include "G4Material.hh"
37#include "G4Element.hh"
38#include "G4Isotope.hh"
39#include "G4ElementVector.hh"
40
41#include "G4MuonMinus.hh"
42#include "G4MuonPlus.hh"
43
44using namespace std;
45using namespace CLHEP;
46
48 : G4VCrossSectionDataSet("NuMuNuclTotXsc")
49{
50 fCofXsc = 1.e-38*cm2/GeV;
51
52 // G4cout<<"fCofXsc = "<<fCofXsc*GeV/cm2<<" cm2/GeV"<<G4endl;
53
54 // PDG2016: sin^2 theta Weinberg
55
56 fSin2tW = 0.23129; // 0.2312;
57
58 // 9 <-> 6, 5/9 or 5/6 ?
59
60 fCofS = 5.*fSin2tW*fSin2tW/9.;
61 fCofL = 1. - fSin2tW + fCofS;
62
63 // G4cout<<"fCosL = "<<fCofL<<", fCofS = "<<fCofS<<G4endl;
64
65 fCutEnergy = 0.; // default value
66
67 fBiasingFactor = 1.; // default as physics
68
69 fIndex = 50;
70
71 fTotXsc = 0.;
72 fCcTotRatio = 0.75; // from nc/cc~0.33 ratio
73 fCcFactor = fNcFactor = 1.;
74 fQEratio = 0.5; // mean in the 1 GeV range
75
78}
79
81{}
82
84
85G4bool
87{
88 G4bool result = false;
89 G4String pName = aPart->GetDefinition()->GetParticleName();
90
91 if( pName == "nu_mu" || pName == "anti_nu_mu" )
92 {
93 result = true;
94 }
95 return result;
96}
97
99
101 G4int Z, const G4Material* mat )
102{
103 G4int Zi(0);
104 size_t i(0), j(0);
105 const G4ElementVector* theElementVector = mat->GetElementVector();
106
107 for ( i = 0; i < theElementVector->size(); ++i )
108 {
109 Zi = (*theElementVector)[i]->GetZasInt();
110 if( Zi == Z ) break;
111 }
112 const G4Element* elm = (*theElementVector)[i];
113 size_t nIso = elm->GetNumberOfIsotopes();
114 G4double fact = 0.0;
115 G4double xsec = 0.0;
116 const G4Isotope* iso = nullptr;
117 const G4IsotopeVector* isoVector = elm->GetIsotopeVector();
118 const G4double* abundVector = elm->GetRelativeAbundanceVector();
119
120 for (j = 0; j<nIso; ++j)
121 {
122 iso = (*isoVector)[j];
123 G4int A = iso->GetN();
124
125 if( abundVector[j] > 0.0 && IsIsoApplicable(part, Z, A, elm, mat) )
126 {
127 fact += abundVector[j];
128 xsec += abundVector[j]*GetIsoCrossSection( part, Z, A, iso, elm, mat);
129 }
130 }
131 if( fact > 0.0) { xsec /= fact; }
132 return xsec;
133}
134
136//
137//
138
140 const G4Isotope*, const G4Element*, const G4Material* )
141{
142 fCcFactor = fNcFactor = 1.;
143 fCcTotRatio = 0.25;
144
145 G4double ccnuXsc, ccanuXsc, ncXsc, totXsc(0.);
146
147 G4double energy = aPart->GetTotalEnergy();
148 G4String pName = aPart->GetDefinition()->GetParticleName();
149
150 G4int index = GetEnergyIndex(energy);
151
152 if( index >= fIndex )
153 {
155 G4double s2 = 2.*energy*pm+pm*pm;
156 G4double aa = 1.;
157 G4double bb = 1.085;
158 G4double mw = 80.385*GeV;
159 fCcFactor = bb/(1.+ aa*s2/mw/mw);
160
161 G4double mz = 91.1876*GeV;
162 fNcFactor = bb/(1.+ aa*s2/mz/mz);
163 }
164 ccnuXsc = GetNuMuTotCsXsc(index, energy, Z, A);
165 ccnuXsc *= fCcFactor;
166 ccanuXsc = GetANuMuTotCsXsc(index, energy, Z, A);
167 ccanuXsc *= fCcFactor;
168
169 if( pName == "nu_mu")
170 {
171 ncXsc = fCofL*ccnuXsc + fCofS*ccanuXsc;
172 ncXsc *= fNcFactor/fCcFactor;
173 totXsc = ccnuXsc + ncXsc;
174 if( totXsc > 0.) fCcTotRatio = ccnuXsc/totXsc;
175 }
176 else if( pName == "anti_nu_mu")
177 {
178 ncXsc = fCofL*ccanuXsc + fCofS*ccnuXsc;
179 ncXsc *= fNcFactor/fCcFactor;
180 totXsc = ccanuXsc + ncXsc;
181 if( totXsc > 0.) fCcTotRatio = ccanuXsc/totXsc;
182 }
183 else return totXsc;
184
185 totXsc *= fCofXsc;
186 totXsc *= energy;
187 // totXsc *= A; // incoherent sum over all isotope nucleons
188
189 totXsc *= fBiasingFactor; // biasing up, if set >1
190
191 fTotXsc = totXsc;
192
193 return totXsc;
194}
195
197//
198// Return index of nu/anu energy array corresponding to the neutrino energy
199
201{
202 G4int i, eIndex = 0;
203
204 for( i = 0; i < fIndex; i++)
205 {
206 if( energy <= fNuMuEnergy[i]*GeV )
207 {
208 eIndex = i;
209 break;
210 }
211 }
212 if( i >= fIndex ) eIndex = i;
213 // G4cout<<"eIndex = "<<eIndex<<G4endl;
214 return eIndex;
215}
216
218//
219// nu_mu xsc for index-1, index linear over energy
220
222{
223 G4double xsc(0.), qexsc(0.), inxsc(0.);
224 G4int nn = aa - zz;
225 if(nn < 1) nn = 0;
226
227 if( index <= 0 || energy < theMuonMinus->GetPDGMass() ) xsc = aa*fNuMuInXsc[0] + nn*fNuMuQeXsc[0];
228 else if (index >= fIndex) xsc = aa*fNuMuInXsc[fIndex-1] + nn*fNuMuQeXsc[fIndex-1];
229 else
230 {
231 G4double x1 = fNuMuEnergy[index-1]*GeV;
232 G4double x2 = fNuMuEnergy[index]*GeV;
233 G4double y1 = fNuMuInXsc[index-1];
234 G4double y2 = fNuMuInXsc[index];
235 G4double z1 = fNuMuQeXsc[index-1];
236 G4double z2 = fNuMuQeXsc[index];
237
238 if(x1 >= x2) return aa*fNuMuInXsc[index] + nn*fNuMuQeXsc[index];
239 else
240 {
241 G4double angle = (y2-y1)/(x2-x1);
242 inxsc = y1 + (energy-x1)*angle;
243 angle = (z2-z1)/(x2-x1);
244 qexsc = z1 + (energy-x1)*angle;
245 qexsc *= nn;
246 xsc = inxsc*aa + qexsc;
247
248 if( xsc > 0.) fQEratio = qexsc/xsc;
249 }
250 }
251 return xsc;
252}
253
255//
256// anu_mu xsc for index-1, index linear over energy
257
259{
260 G4double xsc(0.), qexsc(0.), inxsc(0.);
261
262 if( index <= 0 || energy < theMuonPlus->GetPDGMass() ) xsc = aa*fANuMuInXsc[0] + zz*fANuMuQeXsc[0];
263 else if (index >= fIndex) xsc = aa*fANuMuInXsc[fIndex-1] + zz*fANuMuQeXsc[fIndex-1];
264 else
265 {
266 G4double x1 = fNuMuEnergy[index-1]*GeV;
267 G4double x2 = fNuMuEnergy[index]*GeV;
268 G4double y1 = fANuMuInXsc[index-1];
269 G4double y2 = fANuMuInXsc[index];
270 G4double z1 = fANuMuQeXsc[index-1];
271 G4double z2 = fANuMuQeXsc[index];
272
273 if( x1 >= x2 ) return aa*fANuMuInXsc[index] + zz*fANuMuQeXsc[index];
274 else
275 {
276 G4double angle = (y2-y1)/(x2-x1);
277 inxsc = y1 + (energy-x1)*angle;
278
279 angle = (z2-z1)/(x2-x1);
280 qexsc = z1 + (energy-x1)*angle;
281 qexsc *= zz;
282 xsc = inxsc*aa + qexsc;
283
284 if( xsc > 0.) fQEratio = qexsc/xsc;
285 }
286 }
287 return xsc;
288}
289
291//
292// return fNuMuTotXsc[index] if the index is in the array range
293
295{
296 if( index >= 0 && index < fIndex) return fNuMuInXsc[index] + fNuMuInXsc[index];
297 else
298 {
299 G4cout<<"Improper index of fNuMuTotXsc array"<<G4endl;
300 return 0.;
301 }
302}
303
305//
306// return fANuMuTotXsc[index] if the index is in the array range
307
309{
310 if( index >= 0 && index < fIndex) return fANuMuInXsc[index] + fANuMuQeXsc[index];
311 else
312 {
313 G4cout<<"Improper index of fANuMuTotXsc array"<<G4endl;
314 return 0.;
315 }
316}
317
318
320//
321// E_nu in GeV, ( Eth = 111.603 MeV, EthW = 330.994 MeV)
322
324{
325 0.12, 0.141136, 0.165996, 0.195233, 0.229621,
326 0.270066, 0.317634, 0.373581, 0.439382, 0.516773,
327 0.607795, 0.714849, 0.84076, 0.988848, 1.16302,
328 1.36787, 1.6088, 1.89217, 2.22545, 2.61743,
329 3.07845, 3.62068, 4.25841, 5.00847, 5.89065,
330 6.9282, 8.14851, 9.58376, 11.2718, 13.2572,
331 15.5922, 18.3386, 21.5687, 25.3677, 29.8359,
332 35.0911, 41.2719, 48.5413, 57.0912, 67.147,
333 78.974, 92.8842, 109.244, 128.486, 151.117,
334 177.735, 209.04, 245.86, 289.164, 340.097 };
335
337//
338// XS/E arrays in 10^-38cm2/GeV
339
341{
342 0, 0, 0, 0, 0,
343 0, 0, 0.0166853, 0.0649693, 0.132346,
344 0.209102, 0.286795, 0.3595, 0.423961, 0.479009,
345 0.524797, 0.562165, 0.592225, 0.61612, 0.63491,
346 0.649524, 0.660751, 0.669245, 0.675546, 0.680092,
347 0.683247, 0.685307, 0.686521, 0.687093, 0.687184,
348 0.686919, 0.686384, 0.685631, 0.684689, 0.68357,
349 0.682275, 0.680806, 0.67917, 0.677376, 0.675442,
350 0.673387, 0.671229, 0.668985, 0.666665, 0.664272,
351 0.661804, 0.65925, 0.656593, 0.65381, 0.650871 };
352
354{
355 0.20787, 0.411055, 0.570762, 0.705379, 0.814702,
356 0.89543, 0.944299, 0.959743, 0.942906, 0.897917,
357 0.831331, 0.750948, 0.66443, 0.578191, 0.496828,
358 0.423071, 0.358103, 0.302016, 0.254241, 0.213889,
359 0.179971, 0.151527, 0.12769, 0.107706, 0.0909373,
360 0.0768491, 0.0649975, 0.0550143, 0.0465948, 0.0394861,
361 0.0334782, 0.0283964, 0.0240945, 0.0204506, 0.0173623,
362 0.0147437, 0.0125223, 0.0106374, 0.00903737, 0.00767892,
363 0.00652531, 0.00554547, 0.0047131, 0.0040059, 0.003405,
364 0.00289436, 0.00246039, 0.00209155, 0.00177804, 0.00151152 };
365
366
367
369
371{
372 0, 0, 0, 0, 0,
373 0, 0, 0.00437363, 0.0161485, 0.0333162,
374 0.0557621, 0.0814548, 0.108838, 0.136598, 0.163526,
375 0.188908, 0.212041, 0.232727, 0.250872, 0.26631,
376 0.279467, 0.290341, 0.299177, 0.306299, 0.311864,
377 0.316108, 0.319378, 0.321892, 0.323583, 0.324909,
378 0.325841, 0.326568, 0.327111, 0.327623, 0.32798,
379 0.328412, 0.328704, 0.328988, 0.329326, 0.329559,
380 0.329791, 0.330051, 0.330327, 0.33057, 0.330834,
381 0.331115, 0.331416, 0.331678, 0.33192, 0.332124 };
382
384
386{
387 0.0770264, 0.138754, 0.177006, 0.202417, 0.21804,
388 0.225742, 0.227151, 0.223805, 0.21709, 0.208137,
389 0.197763, 0.186496, 0.174651, 0.162429, 0.14999,
390 0.137498, 0.125127, 0.113057, 0.101455, 0.0904642,
391 0.0801914, 0.0707075, 0.0620483, 0.0542192, 0.0472011,
392 0.0409571, 0.0354377, 0.0305862, 0.0263422, 0.0226451,
393 0.0194358, 0.0166585, 0.0142613, 0.0121968, 0.0104221,
394 0.00889912, 0.00759389, 0.00647662, 0.00552119, 0.00470487,
395 0.00400791, 0.00341322, 0.00290607, 0.00247377, 0.0021054,
396 0.00179162, 0.00152441, 0.00129691, 0.00110323, 0.000938345 };
397
399
std::vector< const G4Element * > G4ElementVector
std::vector< G4Isotope * > G4IsotopeVector
static constexpr double cm2
Definition: G4SIunits.hh:100
static constexpr double GeV
Definition: G4SIunits.hh:203
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]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:163
G4int GetN() const
Definition: G4Isotope.hh:93
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
static const G4double fNuMuInXsc[50]
G4double GetANuMuTotCsArray(G4int index)
G4double GetNuMuTotCsXsc(G4int index, G4double energy, G4int Z, G4int A)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *aPart, G4int Z, G4int A, const G4Isotope *, const G4Element *, const G4Material *)
G4double GetANuMuTotCsXsc(G4int index, G4double energy, G4int Z, G4int A)
G4int GetEnergyIndex(G4double energy)
static const G4double fANuMuQeXsc[50]
virtual G4double GetElementCrossSection(const G4DynamicParticle *dynPart, G4int Z, const G4Material *mat)
static const G4double fNuMuEnergy[50]
static const G4double fANuMuInXsc[50]
G4double GetNuMuTotCsArray(G4int index)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *)
static const G4double fNuMuQeXsc[50]
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
const G4String & GetParticleName() const
Definition: DoubConv.h:17
G4double energy(const ThreeVector &p, const G4double m)
float proton_mass_c2
Definition: hepunit.py:274