OscProb
|
Base class implementing general functions for computing neutrino oscillations. More...
#include <PMNS_Base.h>
Public Member Functions | |
PMNS_Base (int numNus=3) | |
Constructor. More... | |
virtual | ~PMNS_Base () |
Destructor. More... | |
virtual double | Prob (vectorC nu_in, int flvf) |
Compute the probability of nu_in going to flvf. More... | |
virtual double | Prob (vectorC nu_in, int flvf, double E) |
virtual double | Prob (vectorC nu_in, int flvf, double E, double L) |
virtual double | Prob (int flvi, int flvf) |
Compute the probability of flvi going to flvf. More... | |
virtual double | Prob (int flvi, int flvf, double E) |
virtual double | Prob (int flvi, int flvf, double E, double L) |
virtual vectorD | ProbVector (vectorC nu_in) |
virtual vectorD | ProbVector (vectorC nu_in, double E) |
flavours for energy E More... | |
virtual vectorD | ProbVector (vectorC nu_in, double E, double L) |
virtual vectorD | ProbVector (int flvi) |
virtual vectorD | ProbVector (int flvi, double E) |
virtual vectorD | ProbVector (int flvi, double E, double L) |
virtual matrixD | ProbMatrix (int nflvi, int nflvf) |
Compute the probability matrix. More... | |
virtual matrixD | ProbMatrix (int nflvi, int nflvf, double E) |
Compute the probability matrix for energy E. More... | |
virtual matrixD | ProbMatrix (int nflvi, int nflvf, double E, double L) |
virtual double | AvgProb (vectorC nu_in, int flvf, double E, double dE=0) |
Compute the average probability over a bin of energy. More... | |
virtual double | AvgProbLoE (vectorC nu_in, int flvf, double LoE, double dLoE=0) |
Compute the average probability over a bin of L/E. More... | |
virtual double | AvgProb (int flvi, int flvf, double E, double dE=0) |
Compute the average probability over a bin of energy. More... | |
virtual double | AvgProbLoE (int flvi, int flvf, double LoE, double dLoE=0) |
Compute the average probability over a bin of L/E. More... | |
virtual vectorD | AvgProbVector (vectorC nu_in, double E, double dE=0) |
virtual vectorD | AvgProbVectorLoE (vectorC nu_in, double LoE, double dLoE=0) |
Compute the average probability vector over a bin of L/E. More... | |
virtual vectorD | AvgProbVector (int flvi, double E, double dE=0) |
virtual vectorD | AvgProbVectorLoE (int flvi, double LoE, double dLoE=0) |
Compute the average probability vector over a bin of L/E. More... | |
virtual matrixD | AvgProbMatrix (int nflvi, int nflvf, double E, double dE=0) |
virtual matrixD | AvgProbMatrixLoE (int nflvi, int nflvf, double LoE, double dLoE=0) |
Compute the average probability matrix over a bin of L/E. More... | |
virtual vectorC | GetMassEigenstate (int mi) |
Get a neutrino mass eigenstate. More... | |
virtual void | SetAngle (int i, int j, double th) |
Set the mixing angle theta_ij. More... | |
virtual void | SetDelta (int i, int j, double delta) |
Set the CP phase delta_ij. More... | |
virtual void | SetDm (int j, double dm) |
Set the mass-splitting dm_j1 in eV^2. More... | |
virtual double | GetAngle (int i, int j) |
Get the mixing angle theta_ij. More... | |
virtual double | GetDelta (int i, int j) |
Get the CP phase delta_ij. More... | |
virtual double | GetDm (int j) |
Get the mass-splitting dm_j1 in eV^2. More... | |
virtual double | GetDmEff (int j) |
Get the effective mass-splitting dm_j1 in eV^2. More... | |
virtual void | SetStdPars () |
Set PDG 3-flavor parameters. More... | |
virtual void | SetEnergy (double E) |
Set the neutrino energy in GeV. More... | |
virtual void | SetIsNuBar (bool isNuBar) |
Set the anti-neutrino flag. More... | |
virtual double | GetEnergy () |
Get the neutrino energy in GeV. More... | |
virtual bool | GetIsNuBar () |
Get the anti-neutrino flag. More... | |
virtual void | SetPath (NuPath p) |
Set a single path. More... | |
virtual void | SetPath (double length, double density, double zoa=0.5, int layer=0) |
Set a single path. More... | |
virtual void | SetPath (std::vector< NuPath > paths) |
Set a path sequence. More... | |
virtual void | AddPath (NuPath p) |
Add a path to the sequence. More... | |
virtual void | AddPath (double length, double density, double zoa=0.5, int layer=0) |
Add a path to the sequence. More... | |
virtual void | ClearPath () |
Clear the path vector. More... | |
virtual void | SetLength (double L) |
Set a single path lentgh in km. More... | |
virtual void | SetDensity (double rho) |
Set single path density in g/cm^3. More... | |
virtual void | SetZoA (double zoa) |
Set Z/A value for single path. More... | |
virtual void | SetLength (vectorD L) |
Set multiple path lengths. More... | |
virtual void | SetDensity (vectorD rho) |
Set multiple path densities. More... | |
virtual void | SetZoA (vectorD zoa) |
Set multiple path Z/A values. More... | |
virtual void | SetLayers (std::vector< int > lay) |
Set multiple path layer indices. More... | |
virtual void | SetStdPath () |
Set standard neutrino path. More... | |
virtual std::vector< NuPath > | GetPath () |
Get the neutrino path sequence. More... | |
virtual vectorD | GetSamplePoints (double LoE, double dLoE) |
Compute the sample points for a bin of L/E with width dLoE. More... | |
virtual void | SetUseCache (bool u=true) |
Set caching on/off. More... | |
virtual void | ClearCache () |
Clear the cache. More... | |
virtual void | SetMaxCache (int mc=1e6) |
Set max cache size. More... | |
virtual void | SetAvgProbPrec (double prec) |
Set the AvgProb precision. More... | |
Protected Member Functions | |
virtual void | InitializeVectors () |
virtual bool | TryCache () |
Try to find a cached eigensystem. More... | |
virtual void | FillCache () |
Cache the current eigensystem. More... | |
virtual void | SetCurPath (NuPath p) |
Set the path currently in use by the class. More... | |
virtual void | SetAtt (double att, int idx) |
Set one of the path attributes. More... | |
virtual void | SetAtt (vectorD att, int idx) |
Set all values of a path attribute. More... | |
virtual void | RotateH (int i, int j, matrixC &Ham) |
Rotate the Hamiltonian by theta_ij and delta_ij. More... | |
virtual void | RotateState (int i, int j) |
Rotate the neutrino state by theta_ij and delta_ij. More... | |
virtual void | BuildHms () |
Build the matrix of masses squared. More... | |
virtual void | SolveHam ()=0 |
virtual void | ResetToFlavour (int flv) |
Reset neutrino state to pure flavour flv. More... | |
virtual void | SetPureState (vectorC nu_in) |
Set the initial state from a pure state. More... | |
virtual void | PropagatePath (NuPath p) |
Propagate neutrino through a single path. More... | |
virtual void | Propagate () |
Propagate neutrino through full path. More... | |
virtual double | P (int flv) |
Return the probability of final state in flavour flv. More... | |
virtual vectorD | GetProbVector () |
virtual std::vector< int > | GetSortedIndices (const vectorD x) |
Get indices that sort a vector. More... | |
virtual vectorD | ConvertEtoLoE (double E, double dE) |
Protected Attributes | |
int | fNumNus |
Number of neutrino flavours. More... | |
vectorD | fDm |
m^2_i - m^2_1 in vacuum More... | |
matrixD | fTheta |
theta[i][j] mixing angle More... | |
matrixD | fDelta |
delta[i][j] CP violating phase More... | |
vectorC | fNuState |
The neutrino current state. More... | |
matrixC | fHms |
matrix H*2E in eV^2 More... | |
vectorC | fPhases |
Buffer for oscillation phases. More... | |
vectorC | fBuffer |
Buffer for neutrino state tranformations. More... | |
vectorD | fEval |
Eigenvalues of the Hamiltonian. More... | |
matrixC | fEvec |
Eigenvectors of the Hamiltonian. More... | |
double | fEnergy |
Neutrino energy. More... | |
bool | fIsNuBar |
Anti-neutrino flag. More... | |
std::vector< NuPath > | fNuPaths |
Vector of neutrino paths. More... | |
NuPath | fPath |
Current neutrino path. More... | |
bool | fBuiltHms |
Tag to avoid rebuilding Hms. More... | |
bool | fGotES |
Tag to avoid recalculating eigensystem. More... | |
bool | fUseCache |
Flag for whether to use caching. More... | |
double | fCachePrec |
Precision of cache matching. More... | |
int | fMaxCache |
Maximum cache size. More... | |
double | fAvgProbPrec |
AvgProb precision. More... | |
std::unordered_set< EigenPoint > | fMixCache |
Caching set of eigensystems. More... | |
EigenPoint | fProbe |
EigenpPoint to try. More... | |
Static Protected Attributes | |
static const complexD | zero |
zero in complex More... | |
static const complexD | one |
one in complex More... | |
static const double | kKm2eV = 1.0 / 1.973269788e-10 |
km to eV^-1 More... | |
static const double | kK2 |
mol/GeV^2/cm^3 to eV More... | |
static const double | kGeV2eV = 1.0e+09 |
GeV to eV. More... | |
static const double | kNA = 6.022140857e23 |
Avogadro constant. More... | |
static const double | kGf = 1.1663787e-05 |
G_F in units of GeV^-2. More... | |
This is an abstract class implementing the general functions needed for setting up an oscillation calculator. The method for solving the eigensystem for the Hamiltonian must be defined in the derived classes.
Definition at line 26 of file PMNS_Base.h.
PMNS_Base::PMNS_Base | ( | int | numNus = 3 | ) |
Constructor.
Sets the number of neutrinos and initializes attributes
Default starts with a 2 GeV muon neutrino.
Path is set to the default 1000 km in crust density.
Oscillation parameters are from PDG for NH by default.
numNus | - the number of neutrino flavours |
Definition at line 47 of file PMNS_Base.cxx.
References ClearCache(), fNumNus, InitializeVectors(), ResetToFlavour(), SetAvgProbPrec(), SetEnergy(), SetIsNuBar(), SetStdPars(), SetStdPath(), and SetUseCache().
|
virtual |
|
virtual |
Add a path to the sequence defining attributes directly.
length | - The length of the path segment in km |
density | - The density of the path segment in g/cm^3 |
zoa | - The effective Z/A of the path segment |
layer | - An index to identify the layer type (e.g. earth inner core) |
Definition at line 317 of file PMNS_Base.cxx.
References AddPath().
|
virtual |
Add a path to the sequence.
p | - A neutrino path segment |
Definition at line 307 of file PMNS_Base.cxx.
References fNuPaths.
Referenced by AddPath(), SetAtt(), SetPath(), and SetTestPath().
|
virtual |
Compute the average probability of flvi going to flvf over a bin of energy E with width dE.
This gets transformed into L/E, since the oscillation terms have arguments linear in L/E and not E.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller energy ranges.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flvi | - The neutrino starting flavour. |
flvf | - The neutrino final flavour. |
E | - The neutrino energy in the bin center in GeV |
dE | - The energy bin width in GeV |
Definition at line 1500 of file PMNS_Base.cxx.
References AvgProb(), fNuState, and ResetToFlavour().
|
virtual |
Compute the average probability of nu_in going to flvf over a bin of energy E with width dE.
This gets transformed into L/E, since the oscillation terms have arguments linear in L/E and not E.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller energy ranges.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
nu_in | - The neutrino initial state in flavour. |
flvf | - The neutrino final flavour. |
E | - The neutrino energy in the bin center in GeV |
dE | - The energy bin width in GeV |
Definition at line 1568 of file PMNS_Base.cxx.
References AvgProbLoE(), ConvertEtoLoE(), fNuPaths, and Prob().
Referenced by AvgProb(), CheckProb(), and SaveTestFile().
|
virtual |
Compute the average probability of flvi going to flvf over a bin of L/E with width dLoE.
The probabilities are weighted by (L/E)^-2 so that event density is flat in energy. This avoids giving too much weight to low energies. Better approximations would be achieved if we used an interpolated event density.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller L/E ranges.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flvi | - The neutrino starting flavour. |
flvf | - The neutrino final flavour. |
LoE | - The neutrino L/E value in the bin center in km/GeV |
dLoE | - The L/E bin width in km/GeV |
Definition at line 1610 of file PMNS_Base.cxx.
References AvgProbLoE(), fNuState, and ResetToFlavour().
|
virtual |
Compute the average probability of nu_in going to flvf over a bin of L/E with width dLoE.
The probabilities are weighted by (L/E)^-2 so that event density is flat in energy. This avoids giving too much weight to low energies. Better approximations would be achieved if we used an interpolated event density.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller L/E ranges.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
nu_in | - The neutrino intial state in flavour basis. |
flvf | - The neutrino final flavour. |
LoE | - The neutrino L/E value in the bin center in km/GeV |
dLoE | - The L/E bin width in km/GeV |
Definition at line 1643 of file PMNS_Base.cxx.
References OscProb::AvgPath(), fNuPaths, fPath, GetSamplePoints(), OscProb::NuPath::length, Prob(), SetCurPath(), and SetEnergy().
Referenced by AvgProb(), and AvgProbLoE().
|
virtual |
Compute the average probability matrix over a bin of energy
Compute the average probability matrix for nflvi and nflvf over a bin of energy E with width dE.
This gets transformed into L/E, since the oscillation terms have arguments linear in L/E and not E.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller energy ranges.
nflvi | - The number of initial flavours in the matrix. |
nflvf | - The number of final flavours in the matrix. |
E | - The neutrino energy in the bin center in GeV |
dE | - The energy bin width in GeV |
Definition at line 1861 of file PMNS_Base.cxx.
References AvgProbMatrixLoE(), ConvertEtoLoE(), fNuPaths, and ProbMatrix().
|
virtual |
Compute the average probability matrix for nflvi and nflvf over a bin of L/E with width dLoE.
The probabilities are weighted by (L/E)^-2 so that event density is flat in energy. This avoids giving too much weight to low energies. Better approximations would be achieved if we used an interpolated event density.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller L/E ranges.
nflvi | - The number of initial flavours in the matrix. |
nflvf | - The number of final flavours in the matrix. |
LoE | - The neutrino L/E value in the bin center in km/GeV |
dLoE | - The L/E bin width in km/GeV |
Definition at line 1900 of file PMNS_Base.cxx.
References OscProb::AvgPath(), fNuPaths, fPath, GetSamplePoints(), OscProb::NuPath::length, ProbMatrix(), SetCurPath(), and SetEnergy().
Referenced by AvgProbMatrix().
|
virtual |
Compute the average probability vector over a bin of energy
Compute the average probability of nu_in going to all flavours over a bin of energy E with width dE.
This gets transformed into L/E, since the oscillation terms have arguments linear in L/E and not E.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller energy ranges.
flvi | - The neutrino starting flavour. |
E | - The neutrino energy in the bin center in GeV |
dE | - The energy bin width in GeV |
Definition at line 1729 of file PMNS_Base.cxx.
References AvgProbVector(), fNuState, and ResetToFlavour().
Compute the average probability vector over a bin of energy
Compute the average probability of nu_in going to all flavours over a bin of energy E with width dE.
This gets transformed into L/E, since the oscillation terms have arguments linear in L/E and not E.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller energy ranges.
nu_in | - The neutrino initial state in flavour. |
E | - The neutrino energy in the bin center in GeV |
dE | - The energy bin width in GeV |
Definition at line 1753 of file PMNS_Base.cxx.
References AvgProbVectorLoE(), ConvertEtoLoE(), fNumNus, fNuPaths, and ProbVector().
Referenced by AvgProbVector().
|
virtual |
Compute the average probability of flvi going to all flavours over a bin of L/E with width dLoE.
The probabilities are weighted by (L/E)^-2 so that event density is flat in energy. This avoids giving too much weight to low energies. Better approximations would be achieved if we used an interpolated event density.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller L/E ranges.
flvi | - The neutrino starting flavour. |
LoE | - The neutrino L/E value in the bin center in km/GeV |
dLoE | - The L/E bin width in km/GeV |
Definition at line 1705 of file PMNS_Base.cxx.
References AvgProbVectorLoE(), fNuState, and ResetToFlavour().
Compute the average probability of nu_in going to all flavours over a bin of L/E with width dLoE.
The probabilities are weighted by (L/E)^-2 so that event density is flat in energy. This avoids giving too much weight to low energies. Better approximations would be achieved if we used an interpolated event density.
This function works best for single paths. In multiple paths the accuracy may be somewhat worse. If needed, average over smaller L/E ranges.
nu_in | - The neutrino intial state in flavour basis. |
LoE | - The neutrino L/E value in the bin center in km/GeV |
dLoE | - The L/E bin width in km/GeV |
Definition at line 1791 of file PMNS_Base.cxx.
References OscProb::AvgPath(), fNumNus, fNuPaths, fPath, GetSamplePoints(), OscProb::NuPath::length, ProbVector(), SetCurPath(), and SetEnergy().
Referenced by AvgProbVector(), and AvgProbVectorLoE().
|
protectedvirtual |
Build Hms = H*2E, where H is the Hamiltonian in vacuum on flavour basis and E is the neutrino energy in eV. Hms is effectively the matrix of masses squared.
This is a hermitian matrix, so only the upper triangular part needs to be filled
The construction of the Hamiltonian avoids computing terms that are simply zero. This has a big impact in the computation time.
Reimplemented in OscProb::PMNS_Decay, and OscProb::PMNS_SNSI.
Definition at line 955 of file PMNS_Base.cxx.
References ClearCache(), fBuiltHms, fDm, fGotES, fHms, fNumNus, and RotateH().
Referenced by OscProb::PMNS_Fast::SolveHam(), and OscProb::PMNS_Sterile::SolveHam().
|
virtual |
Clear the cache
Definition at line 111 of file PMNS_Base.cxx.
References fMixCache.
Referenced by BuildHms(), PMNS_Base(), OscProb::PMNS_NUNM::SetAlpha(), OscProb::PMNS_LIV::SetaT(), OscProb::PMNS_NSI::SetCoupByIndex(), OscProb::PMNS_LIV::SetcT(), OscProb::PMNS_NSI::SetEps(), and OscProb::PMNS_NUNM::SetFracVnc().
|
virtual |
|
protectedvirtual |
Convert a bin of energy into a bin of L/E
E | - energy bin center in GeV |
dE | - energy bin width in GeV |
Definition at line 1516 of file PMNS_Base.cxx.
References OscProb::AvgPath(), fNuPaths, fPath, OscProb::NuPath::length, and SetCurPath().
Referenced by AvgProb(), AvgProbMatrix(), and AvgProbVector().
|
protectedvirtual |
If using caching, save the eigensystem in memory
Reimplemented in OscProb::PMNS_LIV, and OscProb::PMNS_SNSI.
Definition at line 157 of file PMNS_Base.cxx.
References fEnergy, OscProb::EigenPoint::fEval, fEval, OscProb::EigenPoint::fEvec, fEvec, fIsNuBar, fMaxCache, fMixCache, fNumNus, fPath, fProbe, fUseCache, and OscProb::EigenPoint::SetVars().
Referenced by OscProb::PMNS_Fast::SolveHam(), and OscProb::PMNS_Sterile::SolveHam().
|
virtual |
Get the mixing angle theta_ij in radians.
Requires that i<j. Will notify you if input is wrong. If i>j, will assume reverse order and swap i and j.
i,j | - the indices of theta_ij |
Definition at line 570 of file PMNS_Base.cxx.
|
virtual |
Get the CP phase delta_ij in radians.
Requires that i+1<j. Will notify you if input is wrong. If i>j, will assume reverse order and swap i and j.
i,j | - the indices of delta_ij |
Definition at line 638 of file PMNS_Base.cxx.
|
virtual |
Get the mass-splitting dm_j1 = (m_j^2 - m_1^2) in eV^2
Requires that j>1. Will notify you if input is wrong.
j | - the index of dm_j1 |
Definition at line 696 of file PMNS_Base.cxx.
|
virtual |
Get the effective mass-splitting dm_j1 in matter in eV^2
Requires that j>1. Will notify you if input is wrong.
j | - the index of dm_j1 |
Definition at line 732 of file PMNS_Base.cxx.
References fDm, fEnergy, fEval, fNumNus, GetSortedIndices(), kGeV2eV, and SolveHam().
|
virtual |
|
virtual |
|
virtual |
Get the neutrino mass eigenstate in vacuum
States are:
0 = m_1, 1 = m_2, 2 = m_3, etc.
mi | - the mass eigenstate index |
Definition at line 795 of file PMNS_Base.cxx.
References fNumNus, fNuState, ResetToFlavour(), and RotateState().
|
virtual |
Get the vector of neutrino paths.
Definition at line 300 of file PMNS_Base.cxx.
References fNuPaths.
|
protectedvirtual |
Return vector of probabilities from final state
Get the vector of probabilities for current state
Definition at line 1233 of file PMNS_Base.cxx.
Referenced by ProbVector().
|
virtual |
Compute the sample points for a bin of L/E with width dLoE
This is used for averaging the probability over a bin of L/E. It should be a private function, but I'm keeping it public for now for debugging purposes. The number of sample points seems too high for most purposes. The number of subdivisions needs to be optimized.
LoE | - The neutrino L/E value in the bin center in km/GeV |
dLoE | - The L/E bin width in km/GeV |
Definition at line 1985 of file PMNS_Base.cxx.
References fAvgProbPrec, fEnergy, fEval, fNumNus, kGeV2eV, kKm2eV, and SolveHam().
Referenced by AvgProbLoE(), AvgProbMatrixLoE(), and AvgProbVectorLoE().
Get the indices of the sorted x vector
x | - input vector |
Definition at line 715 of file PMNS_Base.cxx.
Referenced by GetDmEff().
|
protectedvirtual |
Initialize all member vectors with zeros
Set vector sizes and initialize elements to zero.
Definition at line 79 of file PMNS_Base.cxx.
Referenced by PMNS_Base().
|
protectedvirtual |
Compute oscillation probability of flavour flv from current state
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flv | - The neutrino final flavour. |
Reimplemented in OscProb::PMNS_Deco.
Definition at line 1058 of file PMNS_Base.cxx.
References fNumNus, and fNuState.
Referenced by GetProbVector(), and Prob().
|
virtual |
Compute the probability of flvi going to flvf.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flvi | - The neutrino starting flavour. |
flvf | - The neutrino final flavour. |
Definition at line 1091 of file PMNS_Base.cxx.
References P(), Propagate(), and ResetToFlavour().
|
virtual |
Compute the probability of flvi going to flvf for energy E
Compute the probability of flvi going to flvf for a given energy in GeV.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flvi | - The neutrino starting flavour. |
flvf | - The neutrino final flavour. |
E | - The neutrino energy in GeV |
Definition at line 1160 of file PMNS_Base.cxx.
References Prob(), and SetEnergy().
|
virtual |
Compute the probability of flvi going to flvf for energy E and distance L
Compute the probability of flvi going to flvf for a given energy in GeV and distance in km in a single path.
If the path sequence is not a single path, a new single path will be created and the previous sequence will be lost.
Don't use this if you want to propagate over multiple path segments.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flvi | - The neutrino starting flavour. |
flvf | - The neutrino final flavour. |
E | - The neutrino energy in GeV |
L | - The neutrino path length in km |
Definition at line 1219 of file PMNS_Base.cxx.
References Prob(), SetEnergy(), and SetLength().
|
virtual |
Compute the probability of nu_in going to flvf.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
nu_in | - The neutrino initial state in flavour basis. |
flvf | - The neutrino final flavour. |
Definition at line 1114 of file PMNS_Base.cxx.
References P(), Propagate(), and SetPureState().
Referenced by AvgProb(), AvgProbLoE(), and Prob().
|
virtual |
Compute the probability of nu_in going to flvf for energy E
Compute the probability of nu_in going to flvf for a given energy in GeV.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
nu_in | - The neutrino initial state in flavour basis. |
flvf | - The neutrino final flavour. |
E | - The neutrino energy in GeV |
Definition at line 1138 of file PMNS_Base.cxx.
References Prob(), and SetEnergy().
|
virtual |
Compute the probability of nu_in going to flvf for energy E and distance L
Compute the probability of nu_in going to flvf for a given energy in GeV and distance in km in a single path.
If the path sequence is not a single path, a new single path will be created and the previous sequence will be lost.
Don't use this if you want to propagate over multiple path segments.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
nu_in | - The neutrino initial state in flavour basis. |
flvf | - The neutrino final flavour. |
E | - The neutrino energy in GeV |
L | - The neutrino path length in km |
Definition at line 1189 of file PMNS_Base.cxx.
References Prob(), SetEnergy(), and SetLength().
|
virtual |
Compute the probability matrix for the first nflvi and nflvf states.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
nflvi | - The number of initial flavours in the matrix. |
nflvf | - The number of final flavours in the matrix. |
Reimplemented in OscProb::PMNS_Deco, OscProb::PMNS_Deco, and OscProb::PMNS_NUNM.
Definition at line 1387 of file PMNS_Base.cxx.
References fNumNus, fNuPaths, fNuState, PropagatePath(), and ResetToFlavour().
Referenced by AvgProbMatrix(), AvgProbMatrixLoE(), and ProbMatrix().
|
virtual |
Compute the probability matrix for the first nflvi and nflvf states for a given energy in GeV.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
nflvi | - The number of initial flavours in the matrix. |
nflvf | - The number of final flavours in the matrix. |
E | - The neutrino energy in GeV |
Reimplemented in OscProb::PMNS_Deco.
Definition at line 1439 of file PMNS_Base.cxx.
References ProbMatrix(), and SetEnergy().
|
virtual |
Compute the probability matrix for energy E and distance L
Compute the probability matrix for the first nflvi and nflvf states for a given energy in GeV and distance in km in a single path.
If the path sequence is not a single path, a new single path will be created and the previous sequence will be lost.
Don't use this if you want to propagate over multiple path segments.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
nflvi | - The number of initial flavours in the matrix. |
nflvf | - The number of final flavours in the matrix. |
E | - The neutrino energy in GeV |
L | - The neutrino path length in km |
Reimplemented in OscProb::PMNS_Deco.
Definition at line 1468 of file PMNS_Base.cxx.
References ProbMatrix(), SetEnergy(), and SetLength().
|
virtual |
Compute the probabilities of flvi going to all flavours
Compute the probability of flvi going to all flavours.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flvi | - The neutrino starting flavour. |
Definition at line 1272 of file PMNS_Base.cxx.
References GetProbVector(), Propagate(), and ResetToFlavour().
|
virtual |
Compute the probabilities of flvi going to all flavours for energy E
Compute the probability of flvi going to all flavours for a given energy in GeV.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flvi | - The neutrino starting flavour. |
E | - The neutrino energy in GeV |
Definition at line 1313 of file PMNS_Base.cxx.
References ProbVector(), and SetEnergy().
|
virtual |
Compute the probabilities of flvi going to all flavours for energy E and distance L
Compute the probability of flvi going to all flavours for a given energy in GeV and distance in km in a single path.
If the path sequence is not a single path, a new single path will be created and the previous sequence will be lost.
Don't use this if you want to propagate over multiple path segments.
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flvi | - The neutrino starting flavour. |
E | - The neutrino energy in GeV |
L | - The neutrino path length in km |
Definition at line 1365 of file PMNS_Base.cxx.
References ProbVector(), SetEnergy(), and SetLength().
Compute the probabilities of nu_in going to all flavours
Compute the probability of nu_in going to all flavours.
nu_in | - The neutrino initial state in flavour basis. |
Definition at line 1250 of file PMNS_Base.cxx.
References GetProbVector(), Propagate(), and SetPureState().
Referenced by AvgProbVector(), AvgProbVectorLoE(), and ProbVector().
Compute the probabilities of nu_in going to all
Compute the probability of nu_in going to all flavours for a given energy in GeV.
nu_in | - The neutrino initial state in flavour basis. |
E | - The neutrino energy in GeV |
Definition at line 1291 of file PMNS_Base.cxx.
References ProbVector(), and SetEnergy().
Compute the probabilities of nu_in going to all flavours for energy E and distance L
Compute the probability of nu_in going to all flavours for a given energy in GeV and distance in km in a single path.
If the path sequence is not a single path, a new single path will be created and the previous sequence will be lost.
Don't use this if you want to propagate over multiple path segments.
nu_in | - The neutrino initial state in flavour basis. |
E | - The neutrino energy in GeV |
L | - The neutrino path length in km |
Definition at line 1336 of file PMNS_Base.cxx.
References ProbVector(), SetEnergy(), and SetLength().
|
protectedvirtual |
Propagate neutrino state through full path
Reimplemented in OscProb::PMNS_NUNM.
Definition at line 1018 of file PMNS_Base.cxx.
References fNuPaths, and PropagatePath().
Referenced by Prob(), ProbVector(), and OscProb::PMNS_NUNM::Propagate().
|
protectedvirtual |
Propagate the current neutrino state through a given path
p | - A neutrino path segment |
Reimplemented in OscProb::PMNS_Decay, OscProb::PMNS_Deco, OscProb::PMNS_Iter, and OscProb::PMNS_NUNM.
Definition at line 983 of file PMNS_Base.cxx.
References fBuffer, fEval, fEvec, fNumNus, fNuState, fPhases, kKm2eV, OscProb::NuPath::length, SetCurPath(), and SolveHam().
Referenced by ProbMatrix(), Propagate(), OscProb::PMNS_Iter::PropagatePath(), and OscProb::PMNS_NUNM::PropagatePath().
|
protectedvirtual |
Reset the neutrino state back to a pure flavour where it starts
Flavours are:
0 = nue, 1 = numu, 2 = nutau 3 = sterile_1, 4 = sterile_2, etc.
flv | - The neutrino starting flavour. |
Reimplemented in OscProb::PMNS_Deco.
Definition at line 1034 of file PMNS_Base.cxx.
References fNumNus, fNuState, one, and zero.
Referenced by AvgProb(), AvgProbLoE(), AvgProbVector(), AvgProbVectorLoE(), GetMassEigenstate(), PMNS_Base(), Prob(), ProbMatrix(), OscProb::PMNS_NUNM::ProbMatrix(), ProbVector(), and OscProb::PMNS_Deco::ResetToFlavour().
|
protectedvirtual |
Rotate the Hamiltonian by the angle theta_ij and phase delta_ij.
The rotations assume all off-diagonal elements with i > j are zero. This is correct if the order of rotations is chosen appropriately and it speeds up computation by skipping null terms
i,j | - the indices of the rotation ij |
Ham | - the Hamiltonian to be rotated |
Definition at line 822 of file PMNS_Base.cxx.
References fDelta, and fTheta.
Referenced by BuildHms(), OscProb::PMNS_Decay::BuildHms(), and OscProb::PMNS_SNSI::BuildHms().
|
protectedvirtual |
Rotate the neutrino state by the angle theta_ij and phase delta_ij.
i,j | - the indices of the rotation ij |
Definition at line 760 of file PMNS_Base.cxx.
References fDelta, fNuState, and fTheta.
Referenced by GetMassEigenstate().
|
virtual |
Set the mixing angle theta_ij in radians.
Requires that i<j. Will notify you if input is wrong. If i>j, will assume reverse order and swap i and j.
This will check if value is changing to keep track of whether the hamiltonian needs to be rebuilt.
i,j | - the indices of theta_ij |
th | - the value of theta_ij |
Definition at line 539 of file PMNS_Base.cxx.
References fBuiltHms, fNumNus, and fTheta.
Referenced by GetSterile(), OscProb::PMNS_Decay::SetMix(), OscProb::PMNS_Fast::SetMix(), SetNominalPars(), and SetStdPars().
|
protectedvirtual |
Set some single path attribute.
An auxiliary function to set individual attributes in a single path.
If the path sequence is not a single path, a new single path will be created and the previous sequence will be lost. Use with care.
att | - The value of the attribute |
idx | - The index of the attribute (0,1,2,3) = (L, Rho, Z/A, Layer) |
Definition at line 364 of file PMNS_Base.cxx.
References fNuPaths, and SetStdPath().
Referenced by SetDensity(), SetLayers(), SetLength(), and SetZoA().
|
protectedvirtual |
Set all values of a path attribute.
An auxiliary function to set individual attributes in a path sequence.
If the path sequence is of a different size, a new path sequence will be created and the previous sequence will be lost. Use with care.
att | - The values of the attribute |
idx | - The index of the attribute (0,1,2,3) = (L, Rho, Z/A, Layer) |
Definition at line 427 of file PMNS_Base.cxx.
References AddPath(), ClearPath(), OscProb::NuPath::density, fNuPaths, OscProb::NuPath::layer, OscProb::NuPath::length, SetStdPath(), and OscProb::NuPath::zoa.
|
virtual |
Set the precision for the AvgProb method
prec | - AvgProb precision |
Definition at line 1962 of file PMNS_Base.cxx.
References fAvgProbPrec.
Referenced by PMNS_Base().
|
protectedvirtual |
Set the path currentlyin use by the class.
This will be used to know what path to propagate through next.
It will also check if values are changing to keep track of whether the eigensystem needs to be recomputed.
p | - A neutrino path segment |
Definition at line 274 of file PMNS_Base.cxx.
References OscProb::NuPath::density, fGotES, fPath, and OscProb::NuPath::zoa.
Referenced by AvgProbLoE(), AvgProbMatrixLoE(), AvgProbVectorLoE(), ConvertEtoLoE(), PropagatePath(), OscProb::PMNS_Decay::PropagatePath(), and OscProb::PMNS_Deco::PropagatePath().
|
virtual |
Set the CP phase delta_ij in radians.
Requires that i+1<j. Will notify you if input is wrong. If i>j, will assume reverse order and swap i and j.
This will check if value is changing to keep track of whether the hamiltonian needs to be rebuilt.
i,j | - the indices of delta_ij |
delta | - the value of delta_ij |
Definition at line 602 of file PMNS_Base.cxx.
References fBuiltHms, fDelta, and fNumNus.
Referenced by OscProb::PMNS_Decay::SetMix(), OscProb::PMNS_Fast::SetMix(), and SetNominalPars().
|
virtual |
Set single path density.
If the path sequence is not a single path, a new single path will be created and the previous sequence will be lost. Use with care.
rho | - The density of the path segment in g/cm^3 |
Definition at line 402 of file PMNS_Base.cxx.
References SetAtt().
|
virtual |
Set multiple path densities.
If the path sequence is of a different size, a new path sequence will be created and the previous sequence will be lost. Use with care.
rho | - The densities of the path segments in g/cm^3 |
Definition at line 497 of file PMNS_Base.cxx.
References SetAtt().
|
virtual |
Set the mass-splitting dm_j1 = (m_j^2 - m_1^2) in eV^2
Requires that j>1. Will notify you if input is wrong.
This will check if value is changing to keep track of whether the hamiltonian needs to be rebuilt.
j | - the index of dm_j1 |
dm | - the value of dm_j1 |
Definition at line 674 of file PMNS_Base.cxx.
References fBuiltHms, fDm, and fNumNus.
Referenced by GetSterile(), OscProb::PMNS_Decay::SetDeltaMsqrs(), OscProb::PMNS_Fast::SetDeltaMsqrs(), SetNominalPars(), and SetStdPars().
|
virtual |
Set neutrino energy in GeV.
This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.
E | - The neutrino energy in GeV |
Definition at line 226 of file PMNS_Base.cxx.
References fEnergy, and fGotES.
Referenced by AvgProbLoE(), AvgProbMatrixLoE(), AvgProbVectorLoE(), PMNS_Base(), Prob(), ProbMatrix(), and ProbVector().
|
virtual |
Set anti-neutrino flag.
This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.
isNuBar | - Set to true for anti-neutrino and false for neutrino. |
Reimplemented in OscProb::PMNS_Decay, and OscProb::PMNS_Iter.
Definition at line 243 of file PMNS_Base.cxx.
References fGotES, and fIsNuBar.
Referenced by CheckProb(), PMNS_Base(), and SaveTestFile().
|
virtual |
Set multiple path layer indices.
If the path sequence is of a different size, a new path sequence will be created and the previous sequence will be lost. Use with care.
lay | - Indices to identify the layer types (e.g. earth inner core) |
Definition at line 519 of file PMNS_Base.cxx.
References SetAtt().
|
virtual |
Set the length for a single path.
If the path sequence is not a single path, a new single path will be created and the previous sequence will be lost. Use with care.
L | - The length of the path segment in km |
Definition at line 391 of file PMNS_Base.cxx.
References SetAtt().
Referenced by Prob(), ProbMatrix(), and ProbVector().
|
virtual |
Set multiple path lengths.
If the path sequence is of a different size, a new path sequence will be created and the previous sequence will be lost. Use with care.
L | - The lengths of the path segments in km |
Definition at line 486 of file PMNS_Base.cxx.
References SetAtt().
|
virtual |
Set maximum number of cached eigensystems. Finding eigensystems can become slow and take up memory. This protects the cache from becoming too large.
mc | - Max cache size (default: 1e6) |
Definition at line 128 of file PMNS_Base.cxx.
References fMaxCache.
|
virtual |
Set a single path defining attributes directly.
This destroys the current path sequence and creates a new first path.
length | - The length of the path segment in km |
density | - The density of the path segment in g/cm^3 |
zoa | - The effective Z/A of the path segment |
layer | - An index to identify the layer type (e.g. earth inner core) |
Definition at line 347 of file PMNS_Base.cxx.
References SetPath().
|
virtual |
Set a single path.
This destroys the current path sequence and creates a new first path.
p | - A neutrino path segment |
Definition at line 330 of file PMNS_Base.cxx.
References AddPath(), and ClearPath().
Referenced by SetPath(), SetStdPath(), and SetTestPath().
|
virtual |
Set vector of neutrino paths.
paths | - A sequence of neutrino paths |
Definition at line 294 of file PMNS_Base.cxx.
References fNuPaths.
|
protectedvirtual |
Set the initial state from a pure state
nu_in | - The neutrino initial state in flavour basis. |
Reimplemented in OscProb::PMNS_Deco.
Definition at line 1070 of file PMNS_Base.cxx.
References fNumNus, and fNuState.
Referenced by Prob(), and ProbVector().
|
virtual |
Set standard oscillation parameters from PDG 2015.
For two neutrinos, Dm is set to the muon disappearance effective mass-splitting and mixing angle.
Definition at line 177 of file PMNS_Base.cxx.
References fNumNus, SetAngle(), and SetDm().
Referenced by PMNS_Base().
|
virtual |
Set standard single path.
Length is 1000 km, so ~2 GeV peak energy.
Density is approximate from CRUST2.0 (~2.8 g/cm^3). Z/A is set to a round 0.5.
Definition at line 205 of file PMNS_Base.cxx.
References OscProb::NuPath::density, OscProb::NuPath::layer, OscProb::NuPath::length, SetPath(), and OscProb::NuPath::zoa.
Referenced by PMNS_Base(), OscProb::PMNS_Deco::PMNS_Deco(), OscProb::PMNS_LIV::PMNS_LIV(), OscProb::PMNS_NSI::PMNS_NSI(), OscProb::PMNS_NUNM::PMNS_NUNM(), and SetAtt().
|
virtual |
Turn on/off caching of eigensystems. This can save a lot of CPU time by avoiding recomputing eigensystems if we've already seen them recently. Especially useful when running over multiple earth layers and even more if multiple baselines will be computed, e.g. for atmospheric neutrinos.
u | - flag to set caching on (default: true) |
Definition at line 105 of file PMNS_Base.cxx.
References fUseCache.
Referenced by PMNS_Base().
|
virtual |
Set single path Z/A.
If the path sequence is not a single path, a new single path will be created and the previous sequence will be lost. Use with care.
zoa | - The effective Z/A of the path segment |
Definition at line 413 of file PMNS_Base.cxx.
References SetAtt().
|
virtual |
Set multiple path Z/A values.
If the path sequence is of a different size, a new path sequence will be created and the previous sequence will be lost. Use with care.
zoa | - The effective Z/A of the path segments |
Definition at line 508 of file PMNS_Base.cxx.
References SetAtt().
|
protectedpure virtual |
Solve the full Hamiltonian
Not implemented in base class. Solve the full Hamiltonian for eigenvectors and eigenvalues
Implemented in OscProb::PMNS_Decay, OscProb::PMNS_Fast, OscProb::PMNS_Iter, and OscProb::PMNS_Sterile.
Referenced by GetDmEff(), GetSamplePoints(), and PropagatePath().
|
protectedvirtual |
Try to find a cached version of this eigensystem.
Definition at line 134 of file PMNS_Base.cxx.
References fEnergy, fEval, fEvec, fIsNuBar, fMixCache, fNumNus, fPath, fProbe, fUseCache, and OscProb::EigenPoint::SetVars().
Referenced by OscProb::PMNS_Fast::SolveHam(), and OscProb::PMNS_Sterile::SolveHam().
|
protected |
Definition at line 305 of file PMNS_Base.h.
Referenced by GetSamplePoints(), and SetAvgProbPrec().
|
protected |
Definition at line 287 of file PMNS_Base.h.
Referenced by PropagatePath(), and OscProb::PMNS_Decay::PropagatePath().
|
protected |
Definition at line 298 of file PMNS_Base.h.
Referenced by BuildHms(), OscProb::PMNS_Decay::BuildHms(), OscProb::PMNS_SNSI::BuildHms(), OscProb::PMNS_Decay::SetAlpha2(), OscProb::PMNS_Decay::SetAlpha3(), SetAngle(), SetDelta(), SetDm(), OscProb::PMNS_Decay::SetIsNuBar(), OscProb::PMNS_Iter::SetIsNuBar(), OscProb::PMNS_SNSI::SetLowestMass(), and OscProb::PMNS_Iter::SolveHam().
|
protected |
Definition at line 302 of file PMNS_Base.h.
|
protected |
Definition at line 281 of file PMNS_Base.h.
Referenced by GetDelta(), RotateH(), RotateState(), SetDelta(), and OscProb::PMNS_Fast::SetVacuumEigensystem().
|
protected |
Definition at line 279 of file PMNS_Base.h.
Referenced by BuildHms(), OscProb::PMNS_Decay::BuildHms(), OscProb::PMNS_SNSI::BuildHms(), GetDm(), GetDmEff(), OscProb::PMNS_Deco::ProbMatrix(), OscProb::PMNS_Deco::PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), SetDm(), and OscProb::PMNS_Fast::SetVacuumEigensystem().
|
protected |
Definition at line 292 of file PMNS_Base.h.
Referenced by FillCache(), GetDmEff(), GetEnergy(), GetSamplePoints(), OscProb::PMNS_Deco::PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), SetEnergy(), OscProb::PMNS_Fast::SetVacuumEigensystem(), OscProb::PMNS_Iter::SolveHam(), TryCache(), OscProb::PMNS_Decay::UpdateHam(), OscProb::PMNS_Fast::UpdateHam(), OscProb::PMNS_LIV::UpdateHam(), OscProb::PMNS_NSI::UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), OscProb::PMNS_SNSI::UpdateHam(), and OscProb::PMNS_Sterile::UpdateHam().
|
protected |
Definition at line 289 of file PMNS_Base.h.
Referenced by FillCache(), GetDmEff(), GetSamplePoints(), PropagatePath(), OscProb::PMNS_Deco::PropagatePath(), OscProb::PMNS_Fast::SetVacuumEigensystem(), OscProb::PMNS_Sterile::SolveEigenSystem(), OscProb::PMNS_Decay::SolveHam(), OscProb::PMNS_Fast::SolveHam(), OscProb::PMNS_Iter::SolveHam(), and TryCache().
|
protected |
Definition at line 290 of file PMNS_Base.h.
Referenced by FillCache(), PropagatePath(), OscProb::PMNS_Deco::RotateState(), OscProb::PMNS_Fast::SetVacuumEigensystem(), OscProb::PMNS_Sterile::SolveEigenSystem(), OscProb::PMNS_Fast::SolveHam(), and TryCache().
|
protected |
Definition at line 299 of file PMNS_Base.h.
Referenced by BuildHms(), OscProb::PMNS_Decay::BuildHms(), OscProb::PMNS_SNSI::BuildHms(), OscProb::PMNS_NUNM::SetAlpha(), OscProb::PMNS_LIV::SetaT(), OscProb::PMNS_NSI::SetCoupByIndex(), OscProb::PMNS_LIV::SetcT(), SetCurPath(), SetEnergy(), OscProb::PMNS_NSI::SetEps(), OscProb::PMNS_NUNM::SetFracVnc(), SetIsNuBar(), OscProb::PMNS_Decay::SetIsNuBar(), OscProb::PMNS_Decay::SolveHam(), OscProb::PMNS_Fast::SolveHam(), OscProb::PMNS_Iter::SolveHam(), and OscProb::PMNS_Sterile::SolveHam().
|
protected |
Definition at line 284 of file PMNS_Base.h.
Referenced by BuildHms(), OscProb::PMNS_Decay::BuildHms(), OscProb::PMNS_SNSI::BuildHms(), OscProb::PMNS_Decay::UpdateHam(), OscProb::PMNS_Fast::UpdateHam(), OscProb::PMNS_LIV::UpdateHam(), OscProb::PMNS_NSI::UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), OscProb::PMNS_SNSI::UpdateHam(), and OscProb::PMNS_Sterile::UpdateHam().
|
protected |
Definition at line 293 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Decay::BuildHms(), FillCache(), GetIsNuBar(), OscProb::PMNS_Iter::SetExpVL(), SetIsNuBar(), OscProb::PMNS_Decay::SetIsNuBar(), OscProb::PMNS_Iter::SetIsNuBar(), OscProb::PMNS_Fast::SetVacuumEigensystem(), TryCache(), OscProb::PMNS_Decay::UpdateHam(), OscProb::PMNS_Fast::UpdateHam(), OscProb::PMNS_LIV::UpdateHam(), OscProb::PMNS_NSI::UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), OscProb::PMNS_SNSI::UpdateHam(), and OscProb::PMNS_Sterile::UpdateHam().
|
protected |
Definition at line 303 of file PMNS_Base.h.
Referenced by FillCache(), and SetMaxCache().
|
protected |
Definition at line 307 of file PMNS_Base.h.
Referenced by ClearCache(), FillCache(), and TryCache().
|
protected |
Definition at line 277 of file PMNS_Base.h.
Referenced by OscProb::PMNS_NUNM::ApplyAlphaDagger(), AvgProbVector(), AvgProbVectorLoE(), BuildHms(), OscProb::PMNS_Decay::BuildHms(), OscProb::PMNS_SNSI::BuildHms(), FillCache(), OscProb::PMNS_NUNM::GetAlpha(), GetAngle(), OscProb::PMNS_LIV::GetaT(), OscProb::PMNS_LIV::GetcT(), GetDelta(), GetDm(), GetDmEff(), OscProb::PMNS_NSI::GetEps(), OscProb::PMNS_Deco::GetGamma(), GetMassEigenstate(), GetProbVector(), GetSamplePoints(), P(), OscProb::PMNS_Deco::P(), PMNS_Base(), OscProb::PMNS_Decay::PMNS_Decay(), ProbMatrix(), OscProb::PMNS_NUNM::ProbMatrix(), OscProb::PMNS_Deco::ProbMatrix(), PropagatePath(), OscProb::PMNS_Decay::PropagatePath(), OscProb::PMNS_Deco::PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), ResetToFlavour(), OscProb::PMNS_Deco::ResetToFlavour(), OscProb::PMNS_Deco::RotateState(), OscProb::PMNS_NUNM::SetAlpha(), SetAngle(), OscProb::PMNS_LIV::SetaT(), OscProb::PMNS_LIV::SetcT(), SetDelta(), SetDm(), OscProb::PMNS_NSI::SetEps(), OscProb::PMNS_Deco::SetGamma(), SetPureState(), OscProb::PMNS_Deco::SetPureState(), SetStdPars(), OscProb::PMNS_Sterile::SolveEigenSystem(), OscProb::PMNS_Fast::SolveHam(), OscProb::PMNS_Iter::SolveHam(), OscProb::PMNS_Sterile::SolveHam(), TryCache(), OscProb::PMNS_Decay::UpdateHam(), OscProb::PMNS_Fast::UpdateHam(), OscProb::PMNS_LIV::UpdateHam(), OscProb::PMNS_NSI::UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), OscProb::PMNS_SNSI::UpdateHam(), and OscProb::PMNS_Sterile::UpdateHam().
|
protected |
Definition at line 295 of file PMNS_Base.h.
Referenced by AddPath(), AvgProb(), AvgProbLoE(), AvgProbMatrix(), AvgProbMatrixLoE(), AvgProbVector(), AvgProbVectorLoE(), ClearPath(), ConvertEtoLoE(), GetPath(), ProbMatrix(), OscProb::PMNS_Deco::ProbMatrix(), OscProb::PMNS_NUNM::ProbMatrix(), Propagate(), SetAtt(), and SetPath().
|
protected |
Definition at line 283 of file PMNS_Base.h.
Referenced by AvgProb(), AvgProbLoE(), AvgProbVector(), AvgProbVectorLoE(), GetMassEigenstate(), P(), ProbMatrix(), OscProb::PMNS_NUNM::ProbMatrix(), OscProb::PMNS_NUNM::Propagate(), PropagatePath(), OscProb::PMNS_Decay::PropagatePath(), OscProb::PMNS_Iter::PropMatter(), ResetToFlavour(), RotateState(), and SetPureState().
|
protected |
Definition at line 296 of file PMNS_Base.h.
Referenced by AvgProbLoE(), AvgProbMatrixLoE(), AvgProbVectorLoE(), ConvertEtoLoE(), FillCache(), OscProb::PMNS_NSI::GetZoACoup(), SetCurPath(), OscProb::PMNS_Fast::SolveHam(), TryCache(), OscProb::PMNS_Decay::UpdateHam(), OscProb::PMNS_Fast::UpdateHam(), OscProb::PMNS_LIV::UpdateHam(), OscProb::PMNS_NSI::UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), OscProb::PMNS_SNSI::UpdateHam(), and OscProb::PMNS_Sterile::UpdateHam().
|
protected |
Definition at line 286 of file PMNS_Base.h.
Referenced by PropagatePath().
|
protected |
Definition at line 308 of file PMNS_Base.h.
Referenced by FillCache(), and TryCache().
|
protected |
Definition at line 280 of file PMNS_Base.h.
Referenced by GetAngle(), OscProb::PMNS_Deco::ProbMatrix(), RotateH(), RotateState(), SetAngle(), and OscProb::PMNS_Fast::SetVacuumEigensystem().
|
protected |
Definition at line 301 of file PMNS_Base.h.
Referenced by FillCache(), SetUseCache(), and TryCache().
|
staticprotected |
Definition at line 217 of file PMNS_Base.h.
Referenced by GetDmEff(), GetSamplePoints(), OscProb::PMNS_Deco::PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), OscProb::PMNS_Fast::SetVacuumEigensystem(), OscProb::PMNS_Decay::UpdateHam(), OscProb::PMNS_Fast::UpdateHam(), OscProb::PMNS_LIV::UpdateHam(), OscProb::PMNS_NSI::UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), OscProb::PMNS_SNSI::UpdateHam(), and OscProb::PMNS_Sterile::UpdateHam().
|
staticprotected |
Definition at line 220 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Iter::SetExpVL(), OscProb::PMNS_Decay::UpdateHam(), OscProb::PMNS_Fast::UpdateHam(), OscProb::PMNS_LIV::UpdateHam(), OscProb::PMNS_NSI::UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), OscProb::PMNS_SNSI::UpdateHam(), and OscProb::PMNS_Sterile::UpdateHam().
|
staticprotected |
Definition at line 216 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Iter::SetExpVL(), OscProb::PMNS_Decay::UpdateHam(), OscProb::PMNS_Fast::UpdateHam(), OscProb::PMNS_LIV::UpdateHam(), OscProb::PMNS_NSI::UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), OscProb::PMNS_SNSI::UpdateHam(), and OscProb::PMNS_Sterile::UpdateHam().
|
staticprotected |
Definition at line 215 of file PMNS_Base.h.
Referenced by GetSamplePoints(), PropagatePath(), OscProb::PMNS_Decay::PropagatePath(), OscProb::PMNS_Deco::PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), and OscProb::PMNS_Iter::SetExpVL().
|
staticprotected |
Definition at line 218 of file PMNS_Base.h.
|
staticprotected |
Definition at line 212 of file PMNS_Base.h.
Referenced by ResetToFlavour(), and OscProb::PMNS_Deco::ResetToFlavour().
|
staticprotected |
Definition at line 211 of file PMNS_Base.h.
Referenced by OscProb::PMNS_NUNM::GetAlpha(), OscProb::PMNS_LIV::GetaT(), OscProb::PMNS_LIV::GetcT(), OscProb::PMNS_NSI::GetEps(), OscProb::PMNS_Decay::PMNS_Decay(), ResetToFlavour(), and OscProb::PMNS_Deco::ResetToFlavour().