forked from joaoabcoelho/OscProb
-
Notifications
You must be signed in to change notification settings - Fork 0
/
PMNS_Sterile.cxx
135 lines (111 loc) · 4.34 KB
/
PMNS_Sterile.cxx
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
////////////////////////////////////////////////////////////////////////
//
// Implementation of oscillations of neutrinos in matter in a
// n-neutrino framework.
//
//......................................................................
//
// Throughout I have taken:
// - L to be the neutrino flight distance in km
// - E to be the neutrino energy in GeV
// - Dmij to be the differences between the mass-squares in eV^2
// - Rho to be the matter density in g/cm^3
// - theta_ij and delta_ij to be in radians
//
////////////////////////////////////////////////////////////////////////
#include "PMNS_Sterile.h"
#include <iostream>
#include <cassert>
#include <stdlib.h>
#include <gsl/gsl_complex_math.h>
using namespace std;
using namespace OscProb;
//......................................................................
///
/// Constructor. \sa PMNS_Base::PMNS_Base
///
/// This class can implement any number of neutrino flavours.
///
/// @param numNus - the number of neutrino flavours
///
PMNS_Sterile::PMNS_Sterile(int numNus) : PMNS_Base(numNus),
fEvalGSL(0), fEvecGSL(0), H_GSL(0), W_GSL(0)
{
// Allocate memory for the GSL objects
fEvalGSL = gsl_vector_alloc(fNumNus);
fEvecGSL = gsl_matrix_complex_alloc(fNumNus, fNumNus);
H_GSL = gsl_matrix_complex_alloc(fNumNus, fNumNus);
W_GSL = gsl_eigen_hermv_alloc(fNumNus);
}
//......................................................................
///
/// Destructor.
///
/// Clear GSL objects.
///
PMNS_Sterile::~PMNS_Sterile(){
// Free memory from GSL objects
gsl_matrix_complex_free(H_GSL); H_GSL = 0;
gsl_eigen_hermv_free(W_GSL); W_GSL = 0;
gsl_matrix_complex_free(fEvecGSL); fEvecGSL = 0;
gsl_vector_free(fEvalGSL); fEvalGSL = 0;
}
//......................................................................
///
/// Build and solve the full Hamiltonian for eigenvectors and eigenvalues.
///
/// This uses GSL to solve the eigensystem.
///
/// Here we divide the mass squared matrix Hms by the 2E
/// to obtain the vacuum Hamiltonian in eV. Then, the matter
/// potential is added to the electron and sterile components.
///
/// The matter potential is implemented as in arXiv:hep-ph/0606054v3 (eq. 3.15).
///
void PMNS_Sterile::SolveHam()
{
// Build Hamiltonian
BuildHms();
// Check if anything changed
if(fGotES) return;
// Try caching if activated
if(TryCache()) return;
double rho = fPath.density;
double zoa = fPath.zoa;
double lv = 2 * kGeV2eV*fEnergy; // 2E in eV
double kr2GNe = kK2*M_SQRT2*kGf * rho*zoa; // Electron matter potential in eV
double kr2GNn = kK2*M_SQRT2*kGf * rho*(1-zoa)/2; // Neutron matter potential in eV
// Finish building Hamiltonian in matter with dimension of eV
for(size_t i=0;i<size_t(fNumNus);i++){
complexD buf = fHms[i][i]/lv;
*gsl_matrix_complex_ptr(H_GSL, i, i) = gsl_complex_rect(real(buf),0);
for(size_t j=i+1;j<size_t(fNumNus);j++){
buf = fHms[i][j]/lv;
if(!fIsNuBar) *gsl_matrix_complex_ptr(H_GSL, j, i) = gsl_complex_rect(real(buf),-imag(buf));
else *gsl_matrix_complex_ptr(H_GSL, j, i) = gsl_complex_rect(real(buf), imag(buf));
}
if(i>2){
// Subtract NC coherent forward scattering from sterile neutrinos. See arXiv:hep-ph/0606054v3, eq. 3.15, for example.
if(!fIsNuBar) *gsl_matrix_complex_ptr(H_GSL, i, i) = gsl_complex_add_real(gsl_matrix_complex_get(H_GSL,i,i) , kr2GNn);
else *gsl_matrix_complex_ptr(H_GSL, i, i) = gsl_complex_add_real(gsl_matrix_complex_get(H_GSL,i,i) , -kr2GNn);;
}
}
// Add nue CC coherent forward scattering from sterile neutrinos.
if(!fIsNuBar) *gsl_matrix_complex_ptr(H_GSL, 0, 0) = gsl_complex_add_real(gsl_matrix_complex_get(H_GSL,0,0) , kr2GNe);
else *gsl_matrix_complex_ptr(H_GSL, 0, 0) = gsl_complex_add_real(gsl_matrix_complex_get(H_GSL,0,0) , -kr2GNe);
// Solve Hamiltonian for eigensystem
gsl_eigen_hermv(H_GSL, fEvalGSL, fEvecGSL, W_GSL);
// Fill fEval and fEvec vectors from GSL objects
for(int i=0;i<fNumNus;i++){
fEval[i] = gsl_vector_get(fEvalGSL,i);
for(int j=0;j<fNumNus;j++){
gsl_complex buf = gsl_matrix_complex_get(fEvecGSL,i,j);
fEvec[i][j] = complexD( GSL_REAL(buf), GSL_IMAG(buf) );
}
}
fGotES = true;
// Fill cache if activated
FillCache();
}
////////////////////////////////////////////////////////////////////////