OscProb
|
Implementation of oscillations of neutrinos in matter in a three-neutrino framework with decoherence. More...
#include <PMNS_Deco.h>
Public Member Functions | |
PMNS_Deco () | |
Constructor. More... | |
virtual | ~PMNS_Deco () |
Destructor. More... | |
virtual void | SetGamma (int j, double val) |
Set any given decoherence parameter. More... | |
virtual void | SetGamma32 (double val) |
Set the parameter. More... | |
virtual void | SetDecoAngle (double th) |
Set the decoherence angle. More... | |
virtual void | SetPower (double n) |
Set the power index. More... | |
virtual double | GetGamma (int i, int j) |
Get any given decoherence parameter. More... | |
virtual double | GetDecoAngle () |
Get the decoherence angle. More... | |
virtual double | GetPower () |
Get the power index. More... | |
virtual matrixD | ProbMatrix (int nflvi, int nflvf) |
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. More... | |
virtual matrixD | ProbMatrix (int nflvi, int nflvf, double E, double L) |
Compute the probability matrix. More... | |
virtual void | SetMix (double th12, double th23, double th13, double deltacp) |
Set the all mixing parameters at once. More... | |
virtual void | SetDeltaMsqrs (double dm21, double dm32) |
Set both mass-splittings at once. 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 double | AvgProb (vectorC nu_in, int flvf, double E, double dE=0) |
Compute the average probability over a bin of energy. 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 (vectorC nu_in, int flvf, double LoE, double dLoE=0) |
Compute the average probability over a bin of L/E. 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 | AvgProbVector (int flvi, 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 | 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 | SetLength (vectorD L) |
Set multiple path lengths. More... | |
virtual void | SetDensity (double rho) |
Set single path density in g/cm^3. More... | |
virtual void | SetDensity (vectorD rho) |
Set multiple path densities. More... | |
virtual void | SetZoA (double zoa) |
Set Z/A value for single path. 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 | ResetToFlavour (int flv) |
Reset neutrino state to pure flavour flv. More... | |
virtual void | SetPureState (vectorC nu_in) |
Set the density matrix from a pure state. More... | |
virtual void | PropagatePath (NuPath p) |
Propagate neutrino through a single path. More... | |
virtual double | P (int flv) |
Return the probability of final state in flavour flv. More... | |
virtual void | RotateState (bool to_mass) |
Rotate rho to/from mass basis. More... | |
virtual void | UpdateHam () |
Build the full Hamiltonian. More... | |
virtual void | SolveHam () |
Solve the full Hamiltonian for eigenvectors and eigenvalues. More... | |
virtual void | SetVacuumEigensystem () |
Set the eigensystem to the analytic solution of the vacuum Hamiltonian. More... | |
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 | Propagate () |
Propagate neutrino through full path. 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 | |
double | fGamma [3] |
Stores each decoherence parameter. More... | |
double | fPower |
Stores the power index parameter. More... | |
matrixC | fRho |
The neutrino density matrix state. More... | |
matrixC | fMBuffer |
Some memory buffer for matrix operations. More... | |
complexD | fHam [3][3] |
The full hamiltonian. More... | |
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 class expands the PMNS_Fast class including a effects from decoherence in an increasing entropy and energy conserving model.
The model assumes a power law energy dependence of the decoherence parameters and that decoherence occurs in the effective mass basis.
Reference: https://doi.org/10.1140/epjc/s10052-013-2434-6
Definition at line 27 of file PMNS_Deco.h.
PMNS_Deco::PMNS_Deco | ( | ) |
Constructor.
This class is restricted to 3 neutrino flavours.
Definition at line 26 of file PMNS_Deco.cxx.
References SetDecoAngle(), SetGamma(), SetPower(), and OscProb::PMNS_Base::SetStdPath().
|
virtual |
|
virtualinherited |
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 OscProb::PMNS_Base::AddPath().
|
virtualinherited |
Add a path to the sequence.
p | - A neutrino path segment |
Definition at line 307 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fNuPaths.
Referenced by OscProb::PMNS_Base::AddPath(), OscProb::PMNS_Base::SetAtt(), OscProb::PMNS_Base::SetPath(), and SetTestPath().
|
virtualinherited |
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 OscProb::PMNS_Base::AvgProb(), OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::ResetToFlavour().
|
virtualinherited |
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 OscProb::PMNS_Base::AvgProbLoE(), OscProb::PMNS_Base::ConvertEtoLoE(), OscProb::PMNS_Base::fNuPaths, and OscProb::PMNS_Base::Prob().
Referenced by OscProb::PMNS_Base::AvgProb(), CheckProb(), and SaveTestFile().
|
virtualinherited |
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 OscProb::PMNS_Base::AvgProbLoE(), OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::ResetToFlavour().
|
virtualinherited |
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(), OscProb::PMNS_Base::fNuPaths, OscProb::PMNS_Base::fPath, OscProb::PMNS_Base::GetSamplePoints(), OscProb::NuPath::length, OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::SetCurPath(), and OscProb::PMNS_Base::SetEnergy().
Referenced by OscProb::PMNS_Base::AvgProb(), and OscProb::PMNS_Base::AvgProbLoE().
|
virtualinherited |
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 OscProb::PMNS_Base::AvgProbMatrixLoE(), OscProb::PMNS_Base::ConvertEtoLoE(), OscProb::PMNS_Base::fNuPaths, and OscProb::PMNS_Base::ProbMatrix().
|
virtualinherited |
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(), OscProb::PMNS_Base::fNuPaths, OscProb::PMNS_Base::fPath, OscProb::PMNS_Base::GetSamplePoints(), OscProb::NuPath::length, OscProb::PMNS_Base::ProbMatrix(), OscProb::PMNS_Base::SetCurPath(), and OscProb::PMNS_Base::SetEnergy().
Referenced by OscProb::PMNS_Base::AvgProbMatrix().
|
virtualinherited |
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 OscProb::PMNS_Base::AvgProbVector(), OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::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 OscProb::PMNS_Base::AvgProbVectorLoE(), OscProb::PMNS_Base::ConvertEtoLoE(), OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fNuPaths, and OscProb::PMNS_Base::ProbVector().
Referenced by OscProb::PMNS_Base::AvgProbVector().
|
virtualinherited |
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 OscProb::PMNS_Base::AvgProbVectorLoE(), OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::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(), OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fNuPaths, OscProb::PMNS_Base::fPath, OscProb::PMNS_Base::GetSamplePoints(), OscProb::NuPath::length, OscProb::PMNS_Base::ProbVector(), OscProb::PMNS_Base::SetCurPath(), and OscProb::PMNS_Base::SetEnergy().
Referenced by OscProb::PMNS_Base::AvgProbVector(), and OscProb::PMNS_Base::AvgProbVectorLoE().
|
protectedvirtualinherited |
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 OscProb::PMNS_Base::ClearCache(), OscProb::PMNS_Base::fBuiltHms, OscProb::PMNS_Base::fDm, OscProb::PMNS_Base::fGotES, OscProb::PMNS_Base::fHms, OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::RotateH().
Referenced by OscProb::PMNS_Fast::SolveHam(), and OscProb::PMNS_Sterile::SolveHam().
|
virtualinherited |
Clear the cache
Definition at line 111 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fMixCache.
Referenced by OscProb::PMNS_Base::BuildHms(), OscProb::PMNS_Base::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().
|
virtualinherited |
Clear the path vector.
Definition at line 287 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fNuPaths.
Referenced by OscProb::PMNS_Base::SetAtt(), and OscProb::PMNS_Base::SetPath().
|
protectedvirtualinherited |
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(), OscProb::PMNS_Base::fNuPaths, OscProb::PMNS_Base::fPath, OscProb::NuPath::length, and OscProb::PMNS_Base::SetCurPath().
Referenced by OscProb::PMNS_Base::AvgProb(), OscProb::PMNS_Base::AvgProbMatrix(), and OscProb::PMNS_Base::AvgProbVector().
|
protectedvirtualinherited |
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 OscProb::PMNS_Base::fEnergy, OscProb::EigenPoint::fEval, OscProb::PMNS_Base::fEval, OscProb::EigenPoint::fEvec, OscProb::PMNS_Base::fEvec, OscProb::PMNS_Base::fIsNuBar, OscProb::PMNS_Base::fMaxCache, OscProb::PMNS_Base::fMixCache, OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fPath, OscProb::PMNS_Base::fProbe, OscProb::PMNS_Base::fUseCache, and OscProb::EigenPoint::SetVars().
Referenced by OscProb::PMNS_Fast::SolveHam(), and OscProb::PMNS_Sterile::SolveHam().
|
virtualinherited |
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.
References OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::fTheta.
|
virtual |
|
virtualinherited |
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.
References OscProb::PMNS_Base::fDelta, and OscProb::PMNS_Base::fNumNus.
|
virtualinherited |
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.
References OscProb::PMNS_Base::fDm, and OscProb::PMNS_Base::fNumNus.
|
virtualinherited |
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 OscProb::PMNS_Base::fDm, OscProb::PMNS_Base::fEnergy, OscProb::PMNS_Base::fEval, OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::GetSortedIndices(), OscProb::PMNS_Base::kGeV2eV, and OscProb::PMNS_Base::SolveHam().
|
virtualinherited |
Get the neutrino energy in GeV.
Definition at line 255 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fEnergy.
|
virtual |
Get any given decoherence parameter.
Requires that i > j. Will notify you if input is wrong. If i < j, will assume reverse order and swap i and j.
i | - The first mass index |
j | - The second mass index |
Definition at line 185 of file PMNS_Deco.cxx.
References fGamma, and OscProb::PMNS_Base::fNumNus.
Referenced by PropagatePath(), and SetGamma32().
|
virtualinherited |
Get the anti-neutrino flag.
Definition at line 261 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fIsNuBar.
|
virtualinherited |
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 OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fNuState, OscProb::PMNS_Base::ResetToFlavour(), and OscProb::PMNS_Base::RotateState().
|
virtualinherited |
Get the vector of neutrino paths.
Definition at line 300 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fNuPaths.
|
virtual |
Get the power index for energy dependence.
Definition at line 173 of file PMNS_Deco.cxx.
References fPower.
|
protectedvirtualinherited |
Return vector of probabilities from final state
Get the vector of probabilities for current state
Definition at line 1233 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::P().
Referenced by OscProb::PMNS_Base::ProbVector().
|
virtualinherited |
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 OscProb::PMNS_Base::fAvgProbPrec, OscProb::PMNS_Base::fEnergy, OscProb::PMNS_Base::fEval, OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::kGeV2eV, OscProb::PMNS_Base::kKm2eV, and OscProb::PMNS_Base::SolveHam().
Referenced by OscProb::PMNS_Base::AvgProbLoE(), OscProb::PMNS_Base::AvgProbMatrixLoE(), and OscProb::PMNS_Base::AvgProbVectorLoE().
Get the indices of the sorted x vector
x | - input vector |
Definition at line 715 of file PMNS_Base.cxx.
Referenced by OscProb::PMNS_Base::GetDmEff().
|
protectedvirtualinherited |
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 OscProb::PMNS_Base::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 from OscProb::PMNS_Base.
Definition at line 384 of file PMNS_Deco.cxx.
References OscProb::PMNS_Base::fNumNus, and fRho.
|
virtualinherited |
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 OscProb::PMNS_Base::P(), OscProb::PMNS_Base::Propagate(), and OscProb::PMNS_Base::ResetToFlavour().
|
virtualinherited |
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 OscProb::PMNS_Base::Prob(), and OscProb::PMNS_Base::SetEnergy().
|
virtualinherited |
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 OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::SetEnergy(), and OscProb::PMNS_Base::SetLength().
|
virtualinherited |
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 OscProb::PMNS_Base::P(), OscProb::PMNS_Base::Propagate(), and OscProb::PMNS_Base::SetPureState().
Referenced by OscProb::PMNS_Base::AvgProb(), OscProb::PMNS_Base::AvgProbLoE(), and OscProb::PMNS_Base::Prob().
|
virtualinherited |
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 OscProb::PMNS_Base::Prob(), and OscProb::PMNS_Base::SetEnergy().
|
virtualinherited |
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 OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::SetEnergy(), and OscProb::PMNS_Base::SetLength().
|
virtual |
Reimplemented from OscProb::PMNS_Base.
Definition at line 78 of file PMNS_Base.cxx.
|
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 from OscProb::PMNS_Base.
Definition at line 421 of file PMNS_Deco.cxx.
References OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fNuPaths, fRho, PropagatePath(), and ResetToFlavour().
|
virtual |
Reimplemented from OscProb::PMNS_Base.
Definition at line 80 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fDm, OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::fTheta.
|
virtual |
Reimplemented from OscProb::PMNS_Base.
Definition at line 83 of file PMNS_Base.cxx.
|
virtualinherited |
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 OscProb::PMNS_Base::GetProbVector(), OscProb::PMNS_Base::Propagate(), and OscProb::PMNS_Base::ResetToFlavour().
|
virtualinherited |
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 OscProb::PMNS_Base::ProbVector(), and OscProb::PMNS_Base::SetEnergy().
|
virtualinherited |
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 OscProb::PMNS_Base::ProbVector(), OscProb::PMNS_Base::SetEnergy(), and OscProb::PMNS_Base::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 OscProb::PMNS_Base::GetProbVector(), OscProb::PMNS_Base::Propagate(), and OscProb::PMNS_Base::SetPureState().
Referenced by OscProb::PMNS_Base::AvgProbVector(), OscProb::PMNS_Base::AvgProbVectorLoE(), and OscProb::PMNS_Base::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 OscProb::PMNS_Base::ProbVector(), and OscProb::PMNS_Base::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 OscProb::PMNS_Base::ProbVector(), OscProb::PMNS_Base::SetEnergy(), and OscProb::PMNS_Base::SetLength().
|
protectedvirtualinherited |
Propagate neutrino state through full path
Reimplemented in OscProb::PMNS_NUNM.
Definition at line 1018 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fNuPaths, and OscProb::PMNS_Base::PropagatePath().
Referenced by OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::ProbVector(), and OscProb::PMNS_NUNM::Propagate().
|
protectedvirtual |
Propagate the current neutrino state through a given path
p | - A neutrino path segment |
Reimplemented from OscProb::PMNS_Base.
Definition at line 300 of file PMNS_Deco.cxx.
References OscProb::PMNS_Base::fDm, OscProb::PMNS_Base::fEnergy, OscProb::PMNS_Base::fEval, OscProb::PMNS_Base::fNumNus, fPower, fRho, GetGamma(), OscProb::PMNS_Base::kGeV2eV, OscProb::PMNS_Base::kKm2eV, OscProb::NuPath::length, RotateState(), OscProb::PMNS_Base::SetCurPath(), OscProb::PMNS_Fast::SolveHam(), and sort3().
Referenced by ProbMatrix().
|
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 from OscProb::PMNS_Base.
Definition at line 356 of file PMNS_Deco.cxx.
References OscProb::PMNS_Base::fNumNus, fRho, OscProb::PMNS_Base::one, OscProb::PMNS_Base::ResetToFlavour(), and OscProb::PMNS_Base::zero.
Referenced by ProbMatrix().
|
protectedvirtualinherited |
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 OscProb::PMNS_Base::fDelta, and OscProb::PMNS_Base::fTheta.
Referenced by OscProb::PMNS_Base::BuildHms(), OscProb::PMNS_Decay::BuildHms(), and OscProb::PMNS_SNSI::BuildHms().
|
protectedvirtual |
Rotate the density matrix to or from the mass basis
to_mass | - true if to mass basis |
Definition at line 220 of file PMNS_Deco.cxx.
References OscProb::PMNS_Base::fEvec, fMBuffer, OscProb::PMNS_Base::fNumNus, and fRho.
Referenced by PropagatePath().
|
protectedvirtualinherited |
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 OscProb::PMNS_Base::fDelta, OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::fTheta.
Referenced by OscProb::PMNS_Base::GetMassEigenstate().
|
virtualinherited |
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 OscProb::PMNS_Base::fBuiltHms, OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::fTheta.
Referenced by GetSterile(), OscProb::PMNS_Decay::SetMix(), OscProb::PMNS_Fast::SetMix(), SetNominalPars(), and OscProb::PMNS_Base::SetStdPars().
|
protectedvirtualinherited |
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 OscProb::PMNS_Base::fNuPaths, and OscProb::PMNS_Base::SetStdPath().
Referenced by OscProb::PMNS_Base::SetDensity(), OscProb::PMNS_Base::SetLayers(), OscProb::PMNS_Base::SetLength(), and OscProb::PMNS_Base::SetZoA().
|
protectedvirtualinherited |
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 OscProb::PMNS_Base::AddPath(), OscProb::PMNS_Base::ClearPath(), OscProb::NuPath::density, OscProb::PMNS_Base::fNuPaths, OscProb::NuPath::layer, OscProb::NuPath::length, OscProb::PMNS_Base::SetStdPath(), and OscProb::NuPath::zoa.
|
virtualinherited |
Set the precision for the AvgProb method
prec | - AvgProb precision |
Definition at line 1962 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fAvgProbPrec.
Referenced by OscProb::PMNS_Base::PMNS_Base().
|
protectedvirtualinherited |
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, OscProb::PMNS_Base::fGotES, OscProb::PMNS_Base::fPath, and OscProb::NuPath::zoa.
Referenced by OscProb::PMNS_Base::AvgProbLoE(), OscProb::PMNS_Base::AvgProbMatrixLoE(), OscProb::PMNS_Base::AvgProbVectorLoE(), OscProb::PMNS_Base::ConvertEtoLoE(), OscProb::PMNS_Base::PropagatePath(), OscProb::PMNS_Decay::PropagatePath(), and PropagatePath().
|
virtual |
Set the decoherence angle. This will define the relationship:
th | - decoherence angle |
Definition at line 148 of file PMNS_Deco.cxx.
References fGamma.
Referenced by PMNS_Deco().
|
virtualinherited |
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 OscProb::PMNS_Base::fBuiltHms, OscProb::PMNS_Base::fDelta, and OscProb::PMNS_Base::fNumNus.
Referenced by OscProb::PMNS_Decay::SetMix(), OscProb::PMNS_Fast::SetMix(), and SetNominalPars().
|
virtualinherited |
Set both mass-splittings at once.
These are Dm_21 and Dm_32 in eV^2.
The corresponding Dm_31 is set in the class attributes.
dm21 | - The solar mass-splitting Dm_21 |
dm32 | - The atmospheric mass-splitting Dm_32 |
Definition at line 55 of file PMNS_Fast.cxx.
References OscProb::PMNS_Base::SetDm().
|
virtualinherited |
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 OscProb::PMNS_Base::SetAtt().
|
virtualinherited |
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 OscProb::PMNS_Base::SetAtt().
|
virtualinherited |
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 OscProb::PMNS_Base::fBuiltHms, OscProb::PMNS_Base::fDm, and OscProb::PMNS_Base::fNumNus.
Referenced by GetSterile(), OscProb::PMNS_Decay::SetDeltaMsqrs(), OscProb::PMNS_Fast::SetDeltaMsqrs(), SetNominalPars(), and OscProb::PMNS_Base::SetStdPars().
|
virtualinherited |
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 OscProb::PMNS_Base::fEnergy, and OscProb::PMNS_Base::fGotES.
Referenced by OscProb::PMNS_Base::AvgProbLoE(), OscProb::PMNS_Base::AvgProbMatrixLoE(), OscProb::PMNS_Base::AvgProbVectorLoE(), OscProb::PMNS_Base::PMNS_Base(), OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::ProbMatrix(), and OscProb::PMNS_Base::ProbVector().
|
virtual |
Set the decoherence parameter .
Requires that j = 2 or 3. Will notify you if input is wrong.
j | - The first mass index |
val | - The absolute value of the parameter |
Definition at line 51 of file PMNS_Deco.cxx.
References fGamma, and OscProb::PMNS_Base::fNumNus.
Referenced by GetDeco(), and PMNS_Deco().
|
virtual |
Set the decoherence parameter .
In practice, this sets internally via the formula:
IMPORTANT: Note this needs to be used AFTER defining and .
val | - The absolute value of the parameter |
Definition at line 84 of file PMNS_Deco.cxx.
References fGamma, and GetGamma().
|
virtualinherited |
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 OscProb::PMNS_Base::fGotES, and OscProb::PMNS_Base::fIsNuBar.
Referenced by CheckProb(), OscProb::PMNS_Base::PMNS_Base(), and SaveTestFile().
|
virtualinherited |
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 OscProb::PMNS_Base::SetAtt().
|
virtualinherited |
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 OscProb::PMNS_Base::SetAtt().
Referenced by OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::ProbMatrix(), and OscProb::PMNS_Base::ProbVector().
|
virtualinherited |
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 OscProb::PMNS_Base::SetAtt().
|
virtualinherited |
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 OscProb::PMNS_Base::fMaxCache.
|
virtualinherited |
Set all mixing parameters at once.
th12 | - The value of the mixing angle theta_12 |
th23 | - The value of the mixing angle theta_23 |
th13 | - The value of the mixing angle theta_13 |
deltacp | - The value of the CP phase delta_13 |
Definition at line 37 of file PMNS_Fast.cxx.
References OscProb::PMNS_Base::SetAngle(), and OscProb::PMNS_Base::SetDelta().
|
virtualinherited |
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 OscProb::PMNS_Base::SetPath().
|
virtualinherited |
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 OscProb::PMNS_Base::AddPath(), and OscProb::PMNS_Base::ClearPath().
Referenced by OscProb::PMNS_Base::SetPath(), OscProb::PMNS_Base::SetStdPath(), and SetTestPath().
|
virtualinherited |
Set vector of neutrino paths.
paths | - A sequence of neutrino paths |
Definition at line 294 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fNuPaths.
|
virtual |
Set the power index of the decoherence energy dependence. This will multiply the Gammas such that the final parameters are:
n | - power index |
Definition at line 161 of file PMNS_Deco.cxx.
References fPower.
Referenced by PMNS_Deco().
|
protectedvirtual |
Set the density matrix from a pure state
nu_in | - The neutrino initial state in flavour basis. |
Reimplemented from OscProb::PMNS_Base.
Definition at line 396 of file PMNS_Deco.cxx.
References OscProb::PMNS_Base::fNumNus, and fRho.
|
virtualinherited |
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 OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::SetAngle(), and OscProb::PMNS_Base::SetDm().
Referenced by OscProb::PMNS_Base::PMNS_Base().
|
virtualinherited |
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, OscProb::PMNS_Base::SetPath(), and OscProb::NuPath::zoa.
Referenced by OscProb::PMNS_Base::PMNS_Base(), PMNS_Deco(), OscProb::PMNS_LIV::PMNS_LIV(), OscProb::PMNS_NSI::PMNS_NSI(), OscProb::PMNS_NUNM::PMNS_NUNM(), and OscProb::PMNS_Base::SetAtt().
|
virtualinherited |
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 OscProb::PMNS_Base::fUseCache.
Referenced by OscProb::PMNS_Base::PMNS_Base().
|
protectedvirtualinherited |
Set the eigensystem to the analytic solution in vacuum.
We know the vacuum eigensystem, so just write it explicitly
Definition at line 143 of file PMNS_Fast.cxx.
References OscProb::PMNS_Base::fDelta, OscProb::PMNS_Base::fDm, OscProb::PMNS_Base::fEnergy, OscProb::PMNS_Base::fEval, OscProb::PMNS_Base::fEvec, OscProb::PMNS_Base::fIsNuBar, OscProb::PMNS_Base::fTheta, and OscProb::PMNS_Base::kGeV2eV.
Referenced by OscProb::PMNS_Fast::SolveHam(), and OscProb::PMNS_Iter::SolveHam().
|
virtualinherited |
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 OscProb::PMNS_Base::SetAtt().
|
virtualinherited |
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 OscProb::PMNS_Base::SetAtt().
|
protectedvirtualinherited |
Solve the full Hamiltonian for eigenvectors and eigenvalues.
This is using a method from the GLoBES software available at http://www.mpi-hd.mpg.de/personalhomes/globes/3x3/
We should cite them accordingly
Implements OscProb::PMNS_Base.
Reimplemented in OscProb::PMNS_Iter.
Definition at line 100 of file PMNS_Fast.cxx.
References OscProb::PMNS_Base::BuildHms(), OscProb::NuPath::density, OscProb::PMNS_Base::fEval, OscProb::PMNS_Base::fEvec, OscProb::PMNS_Base::fGotES, OscProb::PMNS_Fast::fHam, OscProb::PMNS_Base::FillCache(), OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fPath, OscProb::PMNS_Fast::SetVacuumEigensystem(), OscProb::PMNS_Base::TryCache(), OscProb::PMNS_Fast::UpdateHam(), and zheevh3().
Referenced by PropagatePath().
|
protectedvirtualinherited |
Try to find a cached version of this eigensystem.
Definition at line 134 of file PMNS_Base.cxx.
References OscProb::PMNS_Base::fEnergy, OscProb::PMNS_Base::fEval, OscProb::PMNS_Base::fEvec, OscProb::PMNS_Base::fIsNuBar, OscProb::PMNS_Base::fMixCache, OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fPath, OscProb::PMNS_Base::fProbe, OscProb::PMNS_Base::fUseCache, and OscProb::EigenPoint::SetVars().
Referenced by OscProb::PMNS_Fast::SolveHam(), and OscProb::PMNS_Sterile::SolveHam().
|
protectedvirtualinherited |
Build the full Hamiltonian in matter.
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 component.
Reimplemented in OscProb::PMNS_LIV, OscProb::PMNS_NSI, OscProb::PMNS_NUNM, and OscProb::PMNS_SNSI.
Definition at line 69 of file PMNS_Fast.cxx.
References OscProb::NuPath::density, OscProb::PMNS_Base::fEnergy, OscProb::PMNS_Fast::fHam, OscProb::PMNS_Base::fHms, OscProb::PMNS_Base::fIsNuBar, OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fPath, OscProb::PMNS_Base::kGeV2eV, OscProb::PMNS_Base::kGf, OscProb::PMNS_Base::kK2, and OscProb::NuPath::zoa.
Referenced by OscProb::PMNS_Fast::SolveHam().
|
protectedinherited |
Definition at line 305 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::GetSamplePoints(), and OscProb::PMNS_Base::SetAvgProbPrec().
|
protectedinherited |
Definition at line 287 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::PropagatePath(), and OscProb::PMNS_Decay::PropagatePath().
|
protectedinherited |
Definition at line 298 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::BuildHms(), OscProb::PMNS_Decay::BuildHms(), OscProb::PMNS_SNSI::BuildHms(), OscProb::PMNS_Decay::SetAlpha2(), OscProb::PMNS_Decay::SetAlpha3(), OscProb::PMNS_Base::SetAngle(), OscProb::PMNS_Base::SetDelta(), OscProb::PMNS_Base::SetDm(), OscProb::PMNS_Decay::SetIsNuBar(), OscProb::PMNS_Iter::SetIsNuBar(), OscProb::PMNS_SNSI::SetLowestMass(), and OscProb::PMNS_Iter::SolveHam().
|
protectedinherited |
Definition at line 302 of file PMNS_Base.h.
|
protectedinherited |
Definition at line 281 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::GetDelta(), OscProb::PMNS_Base::RotateH(), OscProb::PMNS_Base::RotateState(), OscProb::PMNS_Base::SetDelta(), and OscProb::PMNS_Fast::SetVacuumEigensystem().
|
protectedinherited |
Definition at line 279 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::BuildHms(), OscProb::PMNS_Decay::BuildHms(), OscProb::PMNS_SNSI::BuildHms(), OscProb::PMNS_Base::GetDm(), OscProb::PMNS_Base::GetDmEff(), ProbMatrix(), PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), OscProb::PMNS_Base::SetDm(), and OscProb::PMNS_Fast::SetVacuumEigensystem().
|
protectedinherited |
Definition at line 292 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::FillCache(), OscProb::PMNS_Base::GetDmEff(), OscProb::PMNS_Base::GetEnergy(), OscProb::PMNS_Base::GetSamplePoints(), PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), OscProb::PMNS_Base::SetEnergy(), OscProb::PMNS_Fast::SetVacuumEigensystem(), OscProb::PMNS_Iter::SolveHam(), OscProb::PMNS_Base::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().
|
protectedinherited |
Definition at line 289 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::FillCache(), OscProb::PMNS_Base::GetDmEff(), OscProb::PMNS_Base::GetSamplePoints(), OscProb::PMNS_Base::PropagatePath(), PropagatePath(), OscProb::PMNS_Fast::SetVacuumEigensystem(), OscProb::PMNS_Sterile::SolveEigenSystem(), OscProb::PMNS_Decay::SolveHam(), OscProb::PMNS_Fast::SolveHam(), OscProb::PMNS_Iter::SolveHam(), and OscProb::PMNS_Base::TryCache().
|
protectedinherited |
Definition at line 290 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::FillCache(), OscProb::PMNS_Base::PropagatePath(), RotateState(), OscProb::PMNS_Fast::SetVacuumEigensystem(), OscProb::PMNS_Sterile::SolveEigenSystem(), OscProb::PMNS_Fast::SolveHam(), and OscProb::PMNS_Base::TryCache().
|
protected |
Definition at line 72 of file PMNS_Deco.h.
Referenced by GetDecoAngle(), GetGamma(), SetDecoAngle(), SetGamma(), and SetGamma32().
|
protectedinherited |
Definition at line 299 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::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(), OscProb::PMNS_Base::SetCurPath(), OscProb::PMNS_Base::SetEnergy(), OscProb::PMNS_NSI::SetEps(), OscProb::PMNS_NUNM::SetFracVnc(), OscProb::PMNS_Base::SetIsNuBar(), OscProb::PMNS_Decay::SetIsNuBar(), OscProb::PMNS_Decay::SolveHam(), OscProb::PMNS_Fast::SolveHam(), OscProb::PMNS_Iter::SolveHam(), and OscProb::PMNS_Sterile::SolveHam().
|
protectedinherited |
Definition at line 62 of file PMNS_Fast.h.
Referenced by OscProb::PMNS_Fast::SolveHam(), OscProb::PMNS_Fast::UpdateHam(), OscProb::PMNS_LIV::UpdateHam(), OscProb::PMNS_NSI::UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), and OscProb::PMNS_SNSI::UpdateHam().
|
protectedinherited |
Definition at line 284 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::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().
|
protectedinherited |
Definition at line 293 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Decay::BuildHms(), OscProb::PMNS_Base::FillCache(), OscProb::PMNS_Base::GetIsNuBar(), OscProb::PMNS_Iter::SetExpVL(), OscProb::PMNS_Base::SetIsNuBar(), OscProb::PMNS_Decay::SetIsNuBar(), OscProb::PMNS_Iter::SetIsNuBar(), OscProb::PMNS_Fast::SetVacuumEigensystem(), OscProb::PMNS_Base::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().
|
protectedinherited |
Definition at line 303 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::FillCache(), and OscProb::PMNS_Base::SetMaxCache().
|
protected |
Definition at line 77 of file PMNS_Deco.h.
Referenced by RotateState().
|
protectedinherited |
Definition at line 307 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::ClearCache(), OscProb::PMNS_Base::FillCache(), and OscProb::PMNS_Base::TryCache().
|
protectedinherited |
Definition at line 277 of file PMNS_Base.h.
Referenced by OscProb::PMNS_NUNM::ApplyAlphaDagger(), OscProb::PMNS_Base::AvgProbVector(), OscProb::PMNS_Base::AvgProbVectorLoE(), OscProb::PMNS_Base::BuildHms(), OscProb::PMNS_Decay::BuildHms(), OscProb::PMNS_SNSI::BuildHms(), OscProb::PMNS_Base::FillCache(), OscProb::PMNS_NUNM::GetAlpha(), OscProb::PMNS_Base::GetAngle(), OscProb::PMNS_LIV::GetaT(), OscProb::PMNS_LIV::GetcT(), OscProb::PMNS_Base::GetDelta(), OscProb::PMNS_Base::GetDm(), OscProb::PMNS_Base::GetDmEff(), OscProb::PMNS_NSI::GetEps(), GetGamma(), OscProb::PMNS_Base::GetMassEigenstate(), OscProb::PMNS_Base::GetProbVector(), OscProb::PMNS_Base::GetSamplePoints(), OscProb::PMNS_Base::P(), P(), OscProb::PMNS_Base::PMNS_Base(), OscProb::PMNS_Decay::PMNS_Decay(), OscProb::PMNS_Base::ProbMatrix(), OscProb::PMNS_NUNM::ProbMatrix(), ProbMatrix(), OscProb::PMNS_Base::PropagatePath(), OscProb::PMNS_Decay::PropagatePath(), PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), OscProb::PMNS_Base::ResetToFlavour(), ResetToFlavour(), RotateState(), OscProb::PMNS_NUNM::SetAlpha(), OscProb::PMNS_Base::SetAngle(), OscProb::PMNS_LIV::SetaT(), OscProb::PMNS_LIV::SetcT(), OscProb::PMNS_Base::SetDelta(), OscProb::PMNS_Base::SetDm(), OscProb::PMNS_NSI::SetEps(), SetGamma(), OscProb::PMNS_Base::SetPureState(), SetPureState(), OscProb::PMNS_Base::SetStdPars(), OscProb::PMNS_Sterile::SolveEigenSystem(), OscProb::PMNS_Fast::SolveHam(), OscProb::PMNS_Iter::SolveHam(), OscProb::PMNS_Sterile::SolveHam(), OscProb::PMNS_Base::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().
|
protectedinherited |
Definition at line 295 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::AddPath(), OscProb::PMNS_Base::AvgProb(), OscProb::PMNS_Base::AvgProbLoE(), OscProb::PMNS_Base::AvgProbMatrix(), OscProb::PMNS_Base::AvgProbMatrixLoE(), OscProb::PMNS_Base::AvgProbVector(), OscProb::PMNS_Base::AvgProbVectorLoE(), OscProb::PMNS_Base::ClearPath(), OscProb::PMNS_Base::ConvertEtoLoE(), OscProb::PMNS_Base::GetPath(), OscProb::PMNS_Base::ProbMatrix(), ProbMatrix(), OscProb::PMNS_NUNM::ProbMatrix(), OscProb::PMNS_Base::Propagate(), OscProb::PMNS_Base::SetAtt(), and OscProb::PMNS_Base::SetPath().
|
protectedinherited |
Definition at line 283 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::AvgProb(), OscProb::PMNS_Base::AvgProbLoE(), OscProb::PMNS_Base::AvgProbVector(), OscProb::PMNS_Base::AvgProbVectorLoE(), OscProb::PMNS_Base::GetMassEigenstate(), OscProb::PMNS_Base::P(), OscProb::PMNS_Base::ProbMatrix(), OscProb::PMNS_NUNM::ProbMatrix(), OscProb::PMNS_NUNM::Propagate(), OscProb::PMNS_Base::PropagatePath(), OscProb::PMNS_Decay::PropagatePath(), OscProb::PMNS_Iter::PropMatter(), OscProb::PMNS_Base::ResetToFlavour(), OscProb::PMNS_Base::RotateState(), and OscProb::PMNS_Base::SetPureState().
|
protectedinherited |
Definition at line 296 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::AvgProbLoE(), OscProb::PMNS_Base::AvgProbMatrixLoE(), OscProb::PMNS_Base::AvgProbVectorLoE(), OscProb::PMNS_Base::ConvertEtoLoE(), OscProb::PMNS_Base::FillCache(), OscProb::PMNS_NSI::GetZoACoup(), OscProb::PMNS_Base::SetCurPath(), OscProb::PMNS_Fast::SolveHam(), OscProb::PMNS_Base::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().
|
protectedinherited |
Definition at line 286 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::PropagatePath().
|
protected |
Definition at line 73 of file PMNS_Deco.h.
Referenced by GetPower(), PropagatePath(), and SetPower().
|
protectedinherited |
Definition at line 308 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::FillCache(), and OscProb::PMNS_Base::TryCache().
|
protected |
Definition at line 75 of file PMNS_Deco.h.
Referenced by P(), ProbMatrix(), PropagatePath(), ResetToFlavour(), RotateState(), and SetPureState().
|
protectedinherited |
Definition at line 280 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::GetAngle(), ProbMatrix(), OscProb::PMNS_Base::RotateH(), OscProb::PMNS_Base::RotateState(), OscProb::PMNS_Base::SetAngle(), and OscProb::PMNS_Fast::SetVacuumEigensystem().
|
protectedinherited |
Definition at line 301 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::FillCache(), OscProb::PMNS_Base::SetUseCache(), and OscProb::PMNS_Base::TryCache().
|
staticprotectedinherited |
Definition at line 217 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::GetDmEff(), OscProb::PMNS_Base::GetSamplePoints(), 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().
|
staticprotectedinherited |
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().
|
staticprotectedinherited |
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().
|
staticprotectedinherited |
Definition at line 215 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::GetSamplePoints(), OscProb::PMNS_Base::PropagatePath(), OscProb::PMNS_Decay::PropagatePath(), PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), and OscProb::PMNS_Iter::SetExpVL().
|
staticprotectedinherited |
Definition at line 218 of file PMNS_Base.h.
|
staticprotectedinherited |
Definition at line 212 of file PMNS_Base.h.
Referenced by OscProb::PMNS_Base::ResetToFlavour(), and ResetToFlavour().
|
staticprotectedinherited |
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(), OscProb::PMNS_Base::ResetToFlavour(), and ResetToFlavour().