OscProb
OscProb::PMNS_NSI Class Reference

Implementation of oscillations of neutrinos in matter in a three-neutrino framework with NSI. More...

#include <PMNS_NSI.h>

Inheritance diagram for OscProb::PMNS_NSI:
OscProb::PMNS_Fast OscProb::PMNS_Base OscProb::PMNS_SNSI

Public Member Functions

 PMNS_NSI ()
 Constructor. More...
 
virtual ~PMNS_NSI ()
 Destructor. More...
 
virtual void SetEps (int flvi, int flvj, double val, double phase)
 Set any given NSI parameter. More...
 
virtual complexD GetEps (int flvi, int flvj)
 Get any given NSI parameter. More...
 
void SetNSI (double eps_ee, double eps_emu, double eps_etau, double eps_mumu, double eps_mutau, double eps_tautau, double delta_emu=0, double delta_etau=0, double delta_mutau=0)
 Set the NSI parameters all at once. More...
 
virtual void SetEps_ee (double a)
 Set eps_ee parameter. More...
 
virtual void SetEps_mumu (double a)
 Set eps_mumu parameter. More...
 
virtual void SetEps_tautau (double a)
 Set eps_tautau parameter. More...
 
virtual void SetEps_emu (double a, double phi)
 Set eps_emu parameter. More...
 
virtual void SetEps_etau (double a, double phi)
 Set eps_etau parameter. More...
 
virtual void SetEps_mutau (double a, double phi)
 Set eps_mutau parameter. More...
 
virtual void SetElecCoup (double e)
 Set electron coupling. More...
 
virtual void SetUpCoup (double u)
 Set u-quark couling. More...
 
virtual void SetDownCoup (double d)
 Set d-quark couling. More...
 
virtual void SetFermCoup (double e, double u, double d)
 Set all fermion couplings. More...
 
virtual void SetCoupByIndex (double c, int i)
 Set a given fermion coupling. More...
 
virtual double GetElecCoup ()
 Get electron coupling. More...
 
virtual double GetUpCoup ()
 Get u-quark couling. More...
 
virtual double GetDownCoup ()
 Get d-quark couling. More...
 
virtual double GetZoACoup ()
 Get effective Z/A coupling. 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 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 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< NuPathGetPath ()
 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 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 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

complexD fEps [3][3]
 Stores each NSI parameter. More...
 
double fNSIcoup [3]
 Relative NSI couplings. 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< NuPathfNuPaths
 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< EigenPointfMixCache
 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...
 

Detailed Description

This class expands the PMNS_Fast class including a general matter potential matrix describing Non-Standard Interactions (NSI).

The matter potential is parametrized by dimensionless quantities epsilon which quantify the intensity of the NSI with respect to the standard matter effects from coherent forward scattering with electrons.

Reference: https://doi.org/10.1088/0034-4885/76/4/044201

See also
PMNS_Fast
Author
jcoelho@apc.in2p3.fr

Definition at line 28 of file PMNS_NSI.h.

Constructor & Destructor Documentation

◆ PMNS_NSI()

PMNS_NSI::PMNS_NSI ( )

Constructor.

See also
PMNS_Base::PMNS_Base

This class is restricted to 3 neutrino flavours.

Definition at line 25 of file PMNS_NSI.cxx.

25 : PMNS_Fast(), fEps()
26{
27 SetStdPath();
28 SetNSI(0., 0., 0., 0., 0., 0., 0., 0., 0.);
29 SetFermCoup(1, 0, 0);
30}
virtual void SetStdPath()
Set standard neutrino path.
Definition: PMNS_Base.cxx:205
PMNS_Fast()
Constructor.
Definition: PMNS_Fast.cxx:21
virtual void SetFermCoup(double e, double u, double d)
Set all fermion couplings.
Definition: PMNS_NSI.cxx:312
complexD fEps[3][3]
Stores each NSI parameter.
Definition: PMNS_NSI.h:76
void SetNSI(double eps_ee, double eps_emu, double eps_etau, double eps_mumu, double eps_mutau, double eps_tautau, double delta_emu=0, double delta_etau=0, double delta_mutau=0)
Set the NSI parameters all at once.
Definition: PMNS_NSI.cxx:55

References SetFermCoup(), SetNSI(), and OscProb::PMNS_Base::SetStdPath().

◆ ~PMNS_NSI()

PMNS_NSI::~PMNS_NSI ( )
virtual

Nothing to clean.

Definition at line 36 of file PMNS_NSI.cxx.

36{}

Member Function Documentation

◆ AddPath() [1/2]

void PMNS_Base::AddPath ( double  length,
double  density,
double  zoa = 0.5,
int  layer = 0 
)
virtualinherited

Add a path to the sequence defining attributes directly.

Parameters
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.

318{
319 AddPath(NuPath(length, density, zoa, layer));
320}
virtual void AddPath(NuPath p)
Add a path to the sequence.
Definition: PMNS_Base.cxx:307
A struct representing a neutrino path segment.
Definition: NuPath.h:34

References OscProb::PMNS_Base::AddPath().

◆ AddPath() [2/2]

void PMNS_Base::AddPath ( NuPath  p)
virtualinherited

Add a path to the sequence.

Parameters
p- A neutrino path segment

Definition at line 307 of file PMNS_Base.cxx.

307{ fNuPaths.push_back(p); }
std::vector< NuPath > fNuPaths
Vector of neutrino paths.
Definition: PMNS_Base.h:295

References OscProb::PMNS_Base::fNuPaths.

Referenced by OscProb::PMNS_Base::AddPath(), OscProb::PMNS_Base::SetAtt(), OscProb::PMNS_Base::SetPath(), and SetTestPath().

◆ AvgProb() [1/2]

double PMNS_Base::AvgProb ( int  flvi,
int  flvf,
double  E,
double  dE = 0 
)
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.
Parameters
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
Returns
Average neutrino oscillation probability

Definition at line 1500 of file PMNS_Base.cxx.

1501{
1502 ResetToFlavour(flvi);
1503
1504 return AvgProb(fNuState, flvf, E, dE);
1505}
vectorC fNuState
The neutrino current state.
Definition: PMNS_Base.h:283
virtual double AvgProb(vectorC nu_in, int flvf, double E, double dE=0)
Compute the average probability over a bin of energy.
Definition: PMNS_Base.cxx:1568
virtual void ResetToFlavour(int flv)
Reset neutrino state to pure flavour flv.
Definition: PMNS_Base.cxx:1034

References OscProb::PMNS_Base::AvgProb(), OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::ResetToFlavour().

◆ AvgProb() [2/2]

double PMNS_Base::AvgProb ( vectorC  nu_in,
int  flvf,
double  E,
double  dE = 0 
)
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.
Parameters
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
Returns
Average neutrino oscillation probability

Definition at line 1568 of file PMNS_Base.cxx.

1569{
1570 // Do nothing if energy is not positive
1571 if (E <= 0) return 0;
1572
1573 if (fNuPaths.empty()) return 0;
1574
1575 // Don't average zero width
1576 if (dE <= 0) return Prob(nu_in, flvf, E);
1577
1578 vectorD LoEbin = ConvertEtoLoE(E, dE);
1579
1580 // Compute average in LoE
1581 return AvgProbLoE(nu_in, flvf, LoEbin[0], LoEbin[1]);
1582}
virtual vectorD ConvertEtoLoE(double E, double dE)
Definition: PMNS_Base.cxx:1516
virtual double AvgProbLoE(vectorC nu_in, int flvf, double LoE, double dLoE=0)
Compute the average probability over a bin of L/E.
Definition: PMNS_Base.cxx:1643
virtual double Prob(vectorC nu_in, int flvf)
Compute the probability of nu_in going to flvf.
Definition: PMNS_Base.cxx:1114
std::vector< double > vectorD
Definition: Definitions.h:18

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().

◆ AvgProbLoE() [1/2]

double PMNS_Base::AvgProbLoE ( int  flvi,
int  flvf,
double  LoE,
double  dLoE = 0 
)
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.
Parameters
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
Returns
Average neutrino oscillation probability

Definition at line 1610 of file PMNS_Base.cxx.

1611{
1612 ResetToFlavour(flvi);
1613
1614 return AvgProbLoE(fNuState, flvf, LoE, dLoE);
1615}

References OscProb::PMNS_Base::AvgProbLoE(), OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::ResetToFlavour().

◆ AvgProbLoE() [2/2]

double PMNS_Base::AvgProbLoE ( vectorC  nu_in,
int  flvf,
double  LoE,
double  dLoE = 0 
)
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.
Parameters
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
Returns
Average neutrino oscillation probability

Definition at line 1643 of file PMNS_Base.cxx.

1644{
1645 // Do nothing if L/E is not positive
1646 if (LoE <= 0) return 0;
1647
1648 if (fNuPaths.empty()) return 0;
1649
1650 // Make sure fPath is set
1651 // Use average if multiple paths
1653
1654 // Set the energy at bin center
1655 SetEnergy(fPath.length / LoE);
1656
1657 // Don't average zero width
1658 if (dLoE <= 0) return Prob(nu_in, flvf);
1659
1660 // Get sample points for this bin
1661 vectorD samples = GetSamplePoints(LoE, dLoE);
1662
1663 // Variables to fill sample
1664 // probabilities and weights
1665 double sumw = 0;
1666 double prob = 0;
1667 double length = fPath.length;
1668
1669 // Loop over all sample points
1670 for (int j = 0; j < int(samples.size()); j++) {
1671 // Set (L/E)^-2 weights
1672 double w = 1. / pow(samples[j], 2);
1673
1674 // Add weighted probability
1675 prob += w * Prob(nu_in, flvf, length / samples[j]);
1676
1677 // Increment sum of weights
1678 sumw += w;
1679 }
1680
1681 // Return weighted average of probabilities
1682 return prob / sumw;
1683}
NuPath fPath
Current neutrino path.
Definition: PMNS_Base.h:296
virtual void SetEnergy(double E)
Set the neutrino energy in GeV.
Definition: PMNS_Base.cxx:226
virtual void SetCurPath(NuPath p)
Set the path currently in use by the class.
Definition: PMNS_Base.cxx:274
virtual vectorD GetSamplePoints(double LoE, double dLoE)
Compute the sample points for a bin of L/E with width dLoE.
Definition: PMNS_Base.cxx:1985
NuPath AvgPath(NuPath &p1, NuPath &p2)
Get the average of two paths.
Definition: NuPath.cxx:27
double length
The length of the path segment in km.
Definition: NuPath.h:78

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().

◆ AvgProbMatrix()

matrixD PMNS_Base::AvgProbMatrix ( int  nflvi,
int  nflvf,
double  E,
double  dE = 0 
)
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.

Parameters
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
Returns
Average neutrino oscillation probabilities

Definition at line 1861 of file PMNS_Base.cxx.

1862{
1863 matrixD probs(nflvi, vectorD(nflvf, 0));
1864
1865 // Do nothing if energy is not positive
1866 if (E <= 0) return probs;
1867
1868 if (fNuPaths.empty()) return probs;
1869
1870 // Don't average zero width
1871 if (dE <= 0) return ProbMatrix(nflvi, nflvf, E);
1872
1873 vectorD LoEbin = ConvertEtoLoE(E, dE);
1874
1875 // Compute average in LoE
1876 return AvgProbMatrixLoE(nflvi, nflvf, LoEbin[0], LoEbin[1]);
1877}
virtual matrixD AvgProbMatrixLoE(int nflvi, int nflvf, double LoE, double dLoE=0)
Compute the average probability matrix over a bin of L/E.
Definition: PMNS_Base.cxx:1900
virtual matrixD ProbMatrix(int nflvi, int nflvf)
Compute the probability matrix.
Definition: PMNS_Base.cxx:1387
std::vector< vectorD > matrixD
Definition: Definitions.h:19

References OscProb::PMNS_Base::AvgProbMatrixLoE(), OscProb::PMNS_Base::ConvertEtoLoE(), OscProb::PMNS_Base::fNuPaths, and OscProb::PMNS_Base::ProbMatrix().

◆ AvgProbMatrixLoE()

matrixD PMNS_Base::AvgProbMatrixLoE ( int  nflvi,
int  nflvf,
double  LoE,
double  dLoE = 0 
)
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.

Parameters
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
Returns
Average neutrino oscillation probabilities

Definition at line 1900 of file PMNS_Base.cxx.

1902{
1903 matrixD probs(nflvi, vectorD(nflvf, 0));
1904
1905 // Do nothing if L/E is not positive
1906 if (LoE <= 0) return probs;
1907
1908 if (fNuPaths.empty()) return probs;
1909
1910 // Make sure fPath is set
1911 // Use average if multiple paths
1913
1914 // Set the energy at bin center
1915 SetEnergy(fPath.length / LoE);
1916
1917 // Don't average zero width
1918 if (dLoE <= 0) return ProbMatrix(nflvi, nflvf);
1919
1920 // Get sample points for this bin
1921 vectorD samples = GetSamplePoints(LoE, dLoE);
1922
1923 // Variables to fill sample
1924 // probabilities and weights
1925 double sumw = 0;
1926 double length = fPath.length;
1927
1928 // Loop over all sample points
1929 for (int j = 0; j < int(samples.size()); j++) {
1930 // Set (L/E)^-2 weights
1931 double w = 1. / pow(samples[j], 2);
1932
1933 matrixD sample_probs = ProbMatrix(nflvi, nflvf, length / samples[j]);
1934
1935 for (int flvi = 0; flvi < nflvi; flvi++) {
1936 for (int flvf = 0; flvf < nflvf; flvf++) {
1937 // Add weighted probability
1938 probs[flvi][flvf] += w * sample_probs[flvi][flvf];
1939 }
1940 }
1941 // Increment sum of weights
1942 sumw += w;
1943 }
1944
1945 for (int flvi = 0; flvi < nflvi; flvi++) {
1946 for (int flvf = 0; flvf < nflvf; flvf++) {
1947 // Divide by total sampling weight
1948 probs[flvi][flvf] /= sumw;
1949 }
1950 }
1951
1952 // Return weighted average of probabilities
1953 return probs;
1954}

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().

◆ AvgProbVector() [1/2]

vectorD PMNS_Base::AvgProbVector ( int  flvi,
double  E,
double  dE = 0 
)
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.

Parameters
flvi- The neutrino starting flavour.
E- The neutrino energy in the bin center in GeV
dE- The energy bin width in GeV
Returns
Average neutrino oscillation probabilities

Definition at line 1729 of file PMNS_Base.cxx.

1730{
1731 ResetToFlavour(flvi);
1732 return AvgProbVector(fNuState, E, dE);
1733}
virtual vectorD AvgProbVector(vectorC nu_in, double E, double dE=0)
Definition: PMNS_Base.cxx:1753

References OscProb::PMNS_Base::AvgProbVector(), OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::ResetToFlavour().

◆ AvgProbVector() [2/2]

vectorD PMNS_Base::AvgProbVector ( vectorC  nu_in,
double  E,
double  dE = 0 
)
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.

Parameters
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
Returns
Average neutrino oscillation probabilities

Definition at line 1753 of file PMNS_Base.cxx.

1754{
1755 vectorD probs(fNumNus, 0);
1756
1757 // Do nothing if energy is not positive
1758 if (E <= 0) return probs;
1759
1760 if (fNuPaths.empty()) return probs;
1761
1762 // Don't average zero width
1763 if (dE <= 0) return ProbVector(nu_in, E);
1764
1765 vectorD LoEbin = ConvertEtoLoE(E, dE);
1766
1767 // Compute average in LoE
1768 return AvgProbVectorLoE(nu_in, LoEbin[0], LoEbin[1]);
1769}
virtual vectorD AvgProbVectorLoE(vectorC nu_in, double LoE, double dLoE=0)
Compute the average probability vector over a bin of L/E.
Definition: PMNS_Base.cxx:1791
int fNumNus
Number of neutrino flavours.
Definition: PMNS_Base.h:277
virtual vectorD ProbVector(vectorC nu_in)
Definition: PMNS_Base.cxx:1250

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().

◆ AvgProbVectorLoE() [1/2]

vectorD PMNS_Base::AvgProbVectorLoE ( int  flvi,
double  LoE,
double  dLoE = 0 
)
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.

Parameters
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
Returns
Average neutrino oscillation probabilities

Definition at line 1705 of file PMNS_Base.cxx.

1706{
1707 ResetToFlavour(flvi);
1708 return AvgProbVectorLoE(fNuState, LoE, dLoE);
1709}

References OscProb::PMNS_Base::AvgProbVectorLoE(), OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::ResetToFlavour().

◆ AvgProbVectorLoE() [2/2]

vectorD PMNS_Base::AvgProbVectorLoE ( vectorC  nu_in,
double  LoE,
double  dLoE = 0 
)
virtualinherited

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.

Parameters
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
Returns
Average neutrino oscillation probabilities

Definition at line 1791 of file PMNS_Base.cxx.

1792{
1793 vectorD probs(fNumNus, 0);
1794
1795 // Do nothing if L/E is not positive
1796 if (LoE <= 0) return probs;
1797
1798 if (fNuPaths.empty()) return probs;
1799
1800 // Make sure fPath is set
1801 // Use average if multiple paths
1803
1804 // Set the energy at bin center
1805 SetEnergy(fPath.length / LoE);
1806
1807 // Don't average zero width
1808 if (dLoE <= 0) return ProbVector(nu_in);
1809
1810 // Get sample points for this bin
1811 vectorD samples = GetSamplePoints(LoE, dLoE);
1812
1813 // Variables to fill sample
1814 // probabilities and weights
1815 double sumw = 0;
1816 double length = fPath.length;
1817
1818 // Loop over all sample points
1819 for (int j = 0; j < int(samples.size()); j++) {
1820 // Set (L/E)^-2 weights
1821 double w = 1. / pow(samples[j], 2);
1822
1823 vectorD sample_probs = ProbVector(nu_in, length / samples[j]);
1824
1825 for (int i = 0; i < fNumNus; i++) {
1826 // Add weighted probability
1827 probs[i] += w * sample_probs[i];
1828 }
1829 // Increment sum of weights
1830 sumw += w;
1831 }
1832
1833 for (int i = 0; i < fNumNus; i++) {
1834 // Divide by total sampling weight
1835 probs[i] /= sumw;
1836 }
1837
1838 // Return weighted average of probabilities
1839 return probs;
1840}

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().

◆ BuildHms()

void PMNS_Base::BuildHms ( )
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.

956{
957 // Check if anything changed
958 if (fBuiltHms) return;
959
960 // Tag to recompute eigensystem
961 fGotES = false;
962
963 for (int j = 0; j < fNumNus; j++) {
964 // Set mass splitting
965 fHms[j][j] = fDm[j];
966 // Reset off-diagonal elements
967 for (int i = 0; i < j; i++) { fHms[i][j] = 0; }
968 // Rotate j neutrinos
969 for (int i = 0; i < j; i++) { RotateH(i, j, fHms); }
970 }
971
972 ClearCache();
973
974 // Tag as built
975 fBuiltHms = true;
976}
virtual void RotateH(int i, int j, matrixC &Ham)
Rotate the Hamiltonian by theta_ij and delta_ij.
Definition: PMNS_Base.cxx:822
matrixC fHms
matrix H*2E in eV^2
Definition: PMNS_Base.h:284
bool fGotES
Tag to avoid recalculating eigensystem.
Definition: PMNS_Base.h:299
bool fBuiltHms
Tag to avoid rebuilding Hms.
Definition: PMNS_Base.h:298
vectorD fDm
m^2_i - m^2_1 in vacuum
Definition: PMNS_Base.h:279
virtual void ClearCache()
Clear the cache.
Definition: PMNS_Base.cxx:111

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().

◆ ClearCache()

void PMNS_Base::ClearCache ( )
virtualinherited

Clear the cache

Definition at line 111 of file PMNS_Base.cxx.

112{
113 fMixCache.clear();
114
115 // Set some better hash table parameters
116 fMixCache.max_load_factor(0.25);
117 fMixCache.reserve(512);
118}
std::unordered_set< EigenPoint > fMixCache
Caching set of eigensystems.
Definition: PMNS_Base.h:307

References OscProb::PMNS_Base::fMixCache.

Referenced by OscProb::PMNS_Base::BuildHms(), OscProb::PMNS_Base::PMNS_Base(), OscProb::PMNS_NUNM::SetAlpha(), OscProb::PMNS_LIV::SetaT(), SetCoupByIndex(), OscProb::PMNS_LIV::SetcT(), SetEps(), and OscProb::PMNS_NUNM::SetFracVnc().

◆ ClearPath()

void PMNS_Base::ClearPath ( )
virtualinherited

Clear the path vector.

Definition at line 287 of file PMNS_Base.cxx.

287{ fNuPaths.clear(); }

References OscProb::PMNS_Base::fNuPaths.

Referenced by OscProb::PMNS_Base::SetAtt(), and OscProb::PMNS_Base::SetPath().

◆ ConvertEtoLoE()

vectorD PMNS_Base::ConvertEtoLoE ( double  E,
double  dE 
)
protectedvirtualinherited

Convert a bin of energy into a bin of L/E

Parameters
E- energy bin center in GeV
dE- energy bin width in GeV
Returns
The L/E bin center and width in km/GeV

Definition at line 1516 of file PMNS_Base.cxx.

1517{
1518 // Make sure fPath is set
1519 // Use average if multiple paths
1521
1522 // Define L/E variables
1523 vectorD LoEbin(2);
1524
1525 // Set a minimum energy
1526 double minE = 0.1 * E;
1527
1528 // Transform range to L/E
1529 // Full range if low edge > minE
1530 if (E - dE / 2 > minE) {
1531 LoEbin[0] =
1532 0.5 * (fPath.length / (E - dE / 2) + fPath.length / (E + dE / 2));
1533 LoEbin[1] = fPath.length / (E - dE / 2) - fPath.length / (E + dE / 2);
1534 }
1535 // Else start at minE
1536 else {
1537 LoEbin[0] = 0.5 * (fPath.length / minE + fPath.length / (E + dE / 2));
1538 LoEbin[1] = fPath.length / minE - fPath.length / (E + dE / 2);
1539 }
1540
1541 return LoEbin;
1542}

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().

◆ FillCache()

void PMNS_Base::FillCache ( )
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.

158{
159 if (fUseCache) {
160 if (fMixCache.size() > fMaxCache) { fMixCache.erase(fMixCache.begin()); }
162 for (int i = 0; i < fNumNus; i++) {
163 fProbe.fEval[i] = fEval[i];
164 for (int j = 0; j < fNumNus; j++) { fProbe.fEvec[i][j] = fEvec[i][j]; }
165 }
166 fMixCache.insert(fProbe);
167 }
168}
bool fIsNuBar
Anti-neutrino flag.
Definition: PMNS_Base.h:293
double fEnergy
Neutrino energy.
Definition: PMNS_Base.h:292
int fMaxCache
Maximum cache size.
Definition: PMNS_Base.h:303
matrixC fEvec
Eigenvectors of the Hamiltonian.
Definition: PMNS_Base.h:290
vectorD fEval
Eigenvalues of the Hamiltonian.
Definition: PMNS_Base.h:289
EigenPoint fProbe
EigenpPoint to try.
Definition: PMNS_Base.h:308
bool fUseCache
Flag for whether to use caching.
Definition: PMNS_Base.h:301
vectorD fEval
Eigenvalues to be cached.
Definition: EigenPoint.h:38
void SetVars(double e=0, NuPath p=NuPath(0, 0), bool n=false)
Set eigensystem parameters.
Definition: EigenPoint.cxx:39
matrixC fEvec
Eigenvectors to be cached.
Definition: EigenPoint.h:39

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().

◆ GetAngle()

double PMNS_Base::GetAngle ( int  i,
int  j 
)
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.

Parameters
i,j- the indices of theta_ij

Definition at line 570 of file PMNS_Base.cxx.

571{
572 if (i > j) {
573 cerr << "WARNING: First argument should be smaller than second argument"
574 << endl
575 << " Setting reverse order (Theta" << j << i << "). " << endl;
576 int temp = i;
577 i = j;
578 j = temp;
579 }
580 if (i < 1 || i > fNumNus - 1 || j < 2 || j > fNumNus) {
581 cerr << "ERROR: Theta" << i << j << " not valid for " << fNumNus
582 << " neutrinos. Returning zero." << endl;
583 return 0;
584 }
585
586 return fTheta[i - 1][j - 1];
587}
matrixD fTheta
theta[i][j] mixing angle
Definition: PMNS_Base.h:280

References OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::fTheta.

◆ GetDelta()

double PMNS_Base::GetDelta ( int  i,
int  j 
)
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.

Parameters
i,j- the indices of delta_ij

Definition at line 638 of file PMNS_Base.cxx.

639{
640 if (i > j) {
641 cerr << "WARNING: First argument should be smaller than second argument"
642 << endl
643 << " Setting reverse order (Delta" << j << i << "). " << endl;
644 int temp = i;
645 i = j;
646 j = temp;
647 }
648 if (i < 1 || i > fNumNus - 1 || j < 2 || j > fNumNus) {
649 cerr << "ERROR: Delta" << i << j << " not valid for " << fNumNus
650 << " neutrinos. Returning zero." << endl;
651 return 0;
652 }
653 if (i + 1 == j) {
654 cerr << "WARNING: Rotation " << i << j << " is real. Returning zero."
655 << endl;
656 return 0;
657 }
658
659 return fDelta[i - 1][j - 1];
660}
matrixD fDelta
delta[i][j] CP violating phase
Definition: PMNS_Base.h:281

References OscProb::PMNS_Base::fDelta, and OscProb::PMNS_Base::fNumNus.

◆ GetDm()

double PMNS_Base::GetDm ( int  j)
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.

Parameters
j- the index of dm_j1

Definition at line 696 of file PMNS_Base.cxx.

697{
698 if (j < 2 || j > fNumNus) {
699 cerr << "ERROR: Dm" << j << "1 not valid for " << fNumNus
700 << " neutrinos. Returning zero." << endl;
701 return 0;
702 }
703
704 return fDm[j - 1];
705}

References OscProb::PMNS_Base::fDm, and OscProb::PMNS_Base::fNumNus.

◆ GetDmEff()

double PMNS_Base::GetDmEff ( int  j)
virtualinherited

Get the effective mass-splitting dm_j1 in matter in eV^2

Requires that j>1. Will notify you if input is wrong.

Parameters
j- the index of dm_j1

Definition at line 732 of file PMNS_Base.cxx.

733{
734 if (j < 2 || j > fNumNus) {
735 cerr << "ERROR: Dm_" << j << "1 not valid for " << fNumNus
736 << " neutrinos. Returning zero." << endl;
737 return 0;
738 }
739
740 // Solve the Hamiltonian to update eigenvalues
741 SolveHam();
742
743 // Sort eigenvalues in same order as vacuum Dm^2
744 vectorI dm_idx = GetSortedIndices(fDm);
745 vectorD dm_idx_double(dm_idx.begin(), dm_idx.end());
746 dm_idx = GetSortedIndices(dm_idx_double);
748
749 // Return difference in eigenvalues * 2E
750 return (fEval[ev_idx[dm_idx[j - 1]]] - fEval[ev_idx[dm_idx[0]]]) * 2 *
752}
virtual void SolveHam()=0
static const double kGeV2eV
GeV to eV.
Definition: PMNS_Base.h:217
virtual std::vector< int > GetSortedIndices(const vectorD x)
Get indices that sort a vector.
Definition: PMNS_Base.cxx:715
std::vector< int > vectorI
Definition: Definitions.h:16

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().

◆ GetDownCoup()

double PMNS_NSI::GetDownCoup ( )
virtual

Get the NSI relative coupling to d-quarks.

Definition at line 335 of file PMNS_NSI.cxx.

335{ return fNSIcoup[2]; }
double fNSIcoup[3]
Relative NSI couplings.
Definition: PMNS_NSI.h:78

References fNSIcoup.

◆ GetElecCoup()

double PMNS_NSI::GetElecCoup ( )
virtual

Get the NSI relative coupling to electrons.

Definition at line 323 of file PMNS_NSI.cxx.

323{ return fNSIcoup[0]; }

References fNSIcoup.

◆ GetEnergy()

double PMNS_Base::GetEnergy ( )
virtualinherited

Get the neutrino energy in GeV.

Definition at line 255 of file PMNS_Base.cxx.

255{ return fEnergy; }

References OscProb::PMNS_Base::fEnergy.

◆ GetEps()

complexD PMNS_NSI::GetEps ( int  flvi,
int  flvj 
)
virtual

Get any given NSI parameter.

Flavours are:

  • 0 = nue, 1 = numu, 2 = nutau

Requires that flvi < flvj. Will notify you if input is wrong. If flvi > flvj, will assume reverse order and swap flvi and flvj.

Parameters
flvi- The first flavour index
flvj- The second flavour index

Definition at line 129 of file PMNS_NSI.cxx.

130{
131 if (flvi > flvj) {
132 cerr << "WARNING: First argument should be smaller or equal to second "
133 "argument"
134 << endl
135 << "Setting reverse order (Eps_" << flvj << flvi << "). " << endl;
136 int temp = flvi;
137 flvi = flvj;
138 flvj = temp;
139 }
140 if (flvi < 0 || flvi > 2 || flvj < flvi || flvj > 2) {
141 cerr << "WARNING: Eps_" << flvi << flvj << " not valid for " << fNumNus
142 << " neutrinos. Returning 0." << endl;
143 return zero;
144 }
145
146 return fEps[flvi][flvj];
147}
static const complexD zero
zero in complex
Definition: PMNS_Base.h:211

References fEps, OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::zero.

◆ GetIsNuBar()

bool PMNS_Base::GetIsNuBar ( )
virtualinherited

Get the anti-neutrino flag.

Definition at line 261 of file PMNS_Base.cxx.

261{ return fIsNuBar; }

References OscProb::PMNS_Base::fIsNuBar.

◆ GetMassEigenstate()

vectorC PMNS_Base::GetMassEigenstate ( int  mi)
virtualinherited

Get the neutrino mass eigenstate in vacuum

States are:

  0 = m_1, 1 = m_2, 2 = m_3, etc.
Parameters
mi- the mass eigenstate index
Returns
The mass eigenstate

Definition at line 795 of file PMNS_Base.cxx.

796{
797 vectorC oldState = fNuState;
798
799 ResetToFlavour(mi);
800
801 for (int j = 0; j < fNumNus; j++) {
802 for (int i = 0; i < j; i++) { RotateState(i, j); }
803 }
804
805 vectorC newState = fNuState;
806 fNuState = oldState;
807
808 return newState;
809}
virtual void RotateState(int i, int j)
Rotate the neutrino state by theta_ij and delta_ij.
Definition: PMNS_Base.cxx:760
std::vector< complexD > vectorC
Definition: Definitions.h:22

References OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fNuState, OscProb::PMNS_Base::ResetToFlavour(), and OscProb::PMNS_Base::RotateState().

◆ GetPath()

vector< NuPath > PMNS_Base::GetPath ( )
virtualinherited

Get the vector of neutrino paths.

Definition at line 300 of file PMNS_Base.cxx.

300{ return fNuPaths; }

References OscProb::PMNS_Base::fNuPaths.

◆ GetProbVector()

vectorD PMNS_Base::GetProbVector ( )
protectedvirtualinherited

Return vector of probabilities from final state

Get the vector of probabilities for current state

Returns
Neutrino oscillation probabilities

Definition at line 1233 of file PMNS_Base.cxx.

1234{
1235 vectorD probs(fNumNus);
1236
1237 for (int i = 0; i < probs.size(); i++) { probs[i] = P(i); }
1238
1239 return probs;
1240}
virtual double P(int flv)
Return the probability of final state in flavour flv.
Definition: PMNS_Base.cxx:1058

References OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::P().

Referenced by OscProb::PMNS_Base::ProbVector().

◆ GetSamplePoints()

vectorD PMNS_Base::GetSamplePoints ( double  LoE,
double  dLoE 
)
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.

Parameters
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.

1986{
1987 // Solve Hamiltonian to get eigenvalues
1988 SolveHam();
1989
1990 // Define conversion factor [km/GeV -> 1/(4 eV^2)]
1991 const double k1267 = kKm2eV / (4 * kGeV2eV);
1992
1993 // Get list of all effective Dm^2
1994 vectorD effDm;
1995
1996 for (int i = 0; i < fNumNus - 1; i++) {
1997 for (int j = i + 1; j < fNumNus; j++) {
1998 effDm.push_back(2 * kGeV2eV * fEnergy * fabs(fEval[j] - fEval[i]));
1999 }
2000 }
2001
2002 int numDm = effDm.size();
2003
2004 // Sort the effective Dm^2 list
2005 sort(effDm.begin(), effDm.end());
2006
2007 // Set a number of sub-divisions to achieve "good" accuracy
2008 // This needs to be studied better
2009 int n_div = ceil(200 * pow(dLoE / LoE, 0.8) / sqrt(fAvgProbPrec / 1e-4));
2010 // int n_div = 1;
2011
2012 // A vector to store sample points
2013 vectorD allSamples;
2014
2015 // Loop over sub-divisions
2016 for (int k = 0; k < n_div; k++) {
2017 // Define sub-division center and width
2018 double bctr = LoE - dLoE / 2 + (k + 0.5) * dLoE / n_div;
2019 double bwdt = dLoE / n_div;
2020
2021 // Make a vector of L/E sample values
2022 // Initialized in the sub-division center
2023 vectorD samples;
2024 samples.push_back(bctr);
2025
2026 // Loop over all Dm^2 to average each frequency
2027 // This will recursively sample points in smaller
2028 // bins so that all relevant frequencies are used
2029 for (int i = 0; i < numDm; i++) {
2030 // Copy the list of sample L/E values
2031 vectorD prev = samples;
2032
2033 // Redefine bin width to lie within full sub-division
2034 double Width =
2035 2 * min(prev[0] - (bctr - bwdt / 2), (bctr + bwdt / 2) - prev[0]);
2036
2037 // Compute oscillation argument sorted from lowest to highest
2038 const double arg = k1267 * effDm[i] * Width;
2039
2040 // Skip small oscillation values.
2041 // If it's the last one, lower the tolerance
2042 if (i < numDm - 1) {
2043 if (arg < 0.9) continue;
2044 }
2045 else {
2046 if (arg < 0.1) continue;
2047 }
2048
2049 // Reset samples to redefine them
2050 samples.clear();
2051
2052 // Loop over previous samples
2053 for (int j = 0; j < int(prev.size()); j++) {
2054 // Compute new sample points around old samples
2055 // This is based on a oscillatory quadrature rule
2056 double sample = (1 / sqrt(3)) * (Width / 2);
2057 if (arg != 0) sample = acos(sin(arg) / arg) / arg * (Width / 2);
2058
2059 // Add samples above and below center
2060 samples.push_back(prev[j] - sample);
2061 samples.push_back(prev[j] + sample);
2062 }
2063
2064 } // End of loop over Dm^2
2065
2066 // Add sub-division samples to the end of allSamples vector
2067 allSamples.insert(allSamples.end(), samples.begin(), samples.end());
2068
2069 } // End of loop over sub-divisions
2070
2071 // Return all sample points
2072 return allSamples;
2073}
static const double kKm2eV
km to eV^-1
Definition: PMNS_Base.h:215
double fAvgProbPrec
AvgProb precision.
Definition: PMNS_Base.h:305

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().

◆ GetSortedIndices()

vectorI PMNS_Base::GetSortedIndices ( const vectorD  x)
protectedvirtualinherited

Get the indices of the sorted x vector

Parameters
x- input vector
Returns
The vector of sorted indices

Definition at line 715 of file PMNS_Base.cxx.

716{
717 vectorI idx(x.size(), 0);
718 for (int i = 0; i < x.size(); i++) idx[i] = i;
719 sort(idx.begin(), idx.end(), IdxCompare(x));
720
721 return idx;
722}
An index sorting comparator.
Definition: PMNS_Base.h:312

Referenced by OscProb::PMNS_Base::GetDmEff().

◆ GetUpCoup()

double PMNS_NSI::GetUpCoup ( )
virtual

Get the NSI relative coupling to u-quarks.

Definition at line 329 of file PMNS_NSI.cxx.

329{ return fNSIcoup[1]; }

References fNSIcoup.

◆ GetZoACoup()

double PMNS_NSI::GetZoACoup ( )
virtual

Get the effective NSI Z/A dependence.

Definition at line 341 of file PMNS_NSI.cxx.

342{
343 return fNSIcoup[0] * fPath.zoa // electrons: Z
344 + fNSIcoup[1] * (1 + fPath.zoa) // u-quarks: A + Z
345 + fNSIcoup[2] * (2 - fPath.zoa); // d-quarks: 2A - Z
346}
double zoa
The effective Z/A value of the path segment.
Definition: NuPath.h:80

References fNSIcoup, OscProb::PMNS_Base::fPath, and OscProb::NuPath::zoa.

Referenced by UpdateHam(), and OscProb::PMNS_SNSI::UpdateHam().

◆ InitializeVectors()

void PMNS_Base::InitializeVectors ( )
protectedvirtualinherited

Initialize all member vectors with zeros

Set vector sizes and initialize elements to zero.

Definition at line 79 of file PMNS_Base.cxx.

80{
81 fDm = vectorD(fNumNus, 0);
84
87
90
91 fEval = vectorD(fNumNus, 0);
93}
vectorC fBuffer
Buffer for neutrino state tranformations.
Definition: PMNS_Base.h:287
vectorC fPhases
Buffer for oscillation phases.
Definition: PMNS_Base.h:286
std::vector< vectorC > matrixC
Definition: Definitions.h:23

Referenced by OscProb::PMNS_Base::PMNS_Base().

◆ P()

double PMNS_Base::P ( int  flv)
protectedvirtualinherited

Compute oscillation probability of flavour flv from current state

Flavours are:

  0 = nue, 1 = numu, 2 = nutau
  3 = sterile_1, 4 = sterile_2, etc.
Parameters
flv- The neutrino final flavour.
Returns
Neutrino oscillation probability

Reimplemented in OscProb::PMNS_Deco.

Definition at line 1058 of file PMNS_Base.cxx.

1059{
1060 assert(flv >= 0 && flv < fNumNus);
1061 return norm(fNuState[flv]);
1062}

References OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::fNuState.

Referenced by OscProb::PMNS_Base::GetProbVector(), and OscProb::PMNS_Base::Prob().

◆ Prob() [1/6]

double PMNS_Base::Prob ( int  flvi,
int  flvf 
)
virtualinherited

Compute the probability of flvi going to flvf.

Flavours are:

  0 = nue, 1 = numu, 2 = nutau
  3 = sterile_1, 4 = sterile_2, etc.
Parameters
flvi- The neutrino starting flavour.
flvf- The neutrino final flavour.
Returns
Neutrino oscillation probability

Definition at line 1091 of file PMNS_Base.cxx.

1092{
1093 ResetToFlavour(flvi);
1094
1095 Propagate();
1096
1097 return P(flvf);
1098}
virtual void Propagate()
Propagate neutrino through full path.
Definition: PMNS_Base.cxx:1018

References OscProb::PMNS_Base::P(), OscProb::PMNS_Base::Propagate(), and OscProb::PMNS_Base::ResetToFlavour().

◆ Prob() [2/6]

double PMNS_Base::Prob ( int  flvi,
int  flvf,
double  E 
)
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.
Parameters
flvi- The neutrino starting flavour.
flvf- The neutrino final flavour.
E- The neutrino energy in GeV
Returns
Neutrino oscillation probability

Definition at line 1160 of file PMNS_Base.cxx.

1161{
1162 SetEnergy(E);
1163
1164 return Prob(flvi, flvf);
1165}

References OscProb::PMNS_Base::Prob(), and OscProb::PMNS_Base::SetEnergy().

◆ Prob() [3/6]

double PMNS_Base::Prob ( int  flvi,
int  flvf,
double  E,
double  L 
)
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.
Parameters
flvi- The neutrino starting flavour.
flvf- The neutrino final flavour.
E- The neutrino energy in GeV
L- The neutrino path length in km
Returns
Neutrino oscillation probability

Definition at line 1219 of file PMNS_Base.cxx.

1220{
1221 SetEnergy(E);
1222 SetLength(L);
1223
1224 return Prob(flvi, flvf);
1225}
virtual void SetLength(double L)
Set a single path lentgh in km.
Definition: PMNS_Base.cxx:391

References OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::SetEnergy(), and OscProb::PMNS_Base::SetLength().

◆ Prob() [4/6]

double PMNS_Base::Prob ( vectorC  nu_in,
int  flvf 
)
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.
Parameters
nu_in- The neutrino initial state in flavour basis.
flvf- The neutrino final flavour.
Returns
Neutrino oscillation probability

Definition at line 1114 of file PMNS_Base.cxx.

1115{
1116 SetPureState(nu_in);
1117
1118 Propagate();
1119
1120 return P(flvf);
1121}
virtual void SetPureState(vectorC nu_in)
Set the initial state from a pure state.
Definition: PMNS_Base.cxx:1070

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().

◆ Prob() [5/6]

double PMNS_Base::Prob ( vectorC  nu_in,
int  flvf,
double  E 
)
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.
Parameters
nu_in- The neutrino initial state in flavour basis.
flvf- The neutrino final flavour.
E- The neutrino energy in GeV
Returns
Neutrino oscillation probability

Definition at line 1138 of file PMNS_Base.cxx.

1139{
1140 SetEnergy(E);
1141
1142 return Prob(nu_in, flvf);
1143}

References OscProb::PMNS_Base::Prob(), and OscProb::PMNS_Base::SetEnergy().

◆ Prob() [6/6]

double PMNS_Base::Prob ( vectorC  nu_in,
int  flvf,
double  E,
double  L 
)
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.
Parameters
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
Returns
Neutrino oscillation probability

Definition at line 1189 of file PMNS_Base.cxx.

1190{
1191 SetEnergy(E);
1192 SetLength(L);
1193
1194 return Prob(nu_in, flvf);
1195}

References OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::SetEnergy(), and OscProb::PMNS_Base::SetLength().

◆ ProbMatrix() [1/3]

matrixD PMNS_Base::ProbMatrix ( int  nflvi,
int  nflvf 
)
virtualinherited

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.
Parameters
nflvi- The number of initial flavours in the matrix.
nflvf- The number of final flavours in the matrix.
Returns
Neutrino oscillation probabilities

Reimplemented in OscProb::PMNS_Deco, OscProb::PMNS_Deco, and OscProb::PMNS_NUNM.

Definition at line 1387 of file PMNS_Base.cxx.

1388{
1389 assert(nflvi <= fNumNus && nflvi >= 0);
1390 assert(nflvf <= fNumNus && nflvf >= 0);
1391
1392 // Output probabilities
1393 matrixD probs(nflvi, vectorD(nflvf));
1394
1395 // List of states
1396 matrixC allstates(nflvi, vectorC(fNumNus));
1397
1398 // Reset all initial states
1399 for (int i = 0; i < nflvi; i++) {
1400 ResetToFlavour(i);
1401 allstates[i] = fNuState;
1402 }
1403
1404 // Propagate all states in parallel
1405 for (int i = 0; i < int(fNuPaths.size()); i++) {
1406 for (int flvi = 0; flvi < nflvi; flvi++) {
1407 fNuState = allstates[flvi];
1409 allstates[flvi] = fNuState;
1410 }
1411 }
1412
1413 // Get all probabilities
1414 for (int flvi = 0; flvi < nflvi; flvi++) {
1415 for (int flvj = 0; flvj < nflvf; flvj++) {
1416 probs[flvi][flvj] = norm(allstates[flvi][flvj]);
1417 }
1418 }
1419
1420 return probs;
1421}
virtual void PropagatePath(NuPath p)
Propagate neutrino through a single path.
Definition: PMNS_Base.cxx:983

References OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fNuPaths, OscProb::PMNS_Base::fNuState, OscProb::PMNS_Base::PropagatePath(), and OscProb::PMNS_Base::ResetToFlavour().

Referenced by OscProb::PMNS_Base::AvgProbMatrix(), OscProb::PMNS_Base::AvgProbMatrixLoE(), and OscProb::PMNS_Base::ProbMatrix().

◆ ProbMatrix() [2/3]

matrixD PMNS_Base::ProbMatrix ( int  nflvi,
int  nflvf,
double  E 
)
virtualinherited

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.
Parameters
nflvi- The number of initial flavours in the matrix.
nflvf- The number of final flavours in the matrix.
E- The neutrino energy in GeV
Returns
Neutrino oscillation probabilities

Reimplemented in OscProb::PMNS_Deco.

Definition at line 1439 of file PMNS_Base.cxx.

1440{
1441 SetEnergy(E);
1442
1443 return ProbMatrix(nflvi, nflvf);
1444}

References OscProb::PMNS_Base::ProbMatrix(), and OscProb::PMNS_Base::SetEnergy().

◆ ProbMatrix() [3/3]

matrixD PMNS_Base::ProbMatrix ( int  nflvi,
int  nflvf,
double  E,
double  L 
)
virtualinherited

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.
Parameters
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
Returns
Neutrino oscillation probabilities

Reimplemented in OscProb::PMNS_Deco.

Definition at line 1468 of file PMNS_Base.cxx.

1469{
1470 SetEnergy(E);
1471 SetLength(L);
1472
1473 return ProbMatrix(nflvi, nflvf);
1474}

References OscProb::PMNS_Base::ProbMatrix(), OscProb::PMNS_Base::SetEnergy(), and OscProb::PMNS_Base::SetLength().

◆ ProbVector() [1/6]

vectorD PMNS_Base::ProbVector ( int  flvi)
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.
Parameters
flvi- The neutrino starting flavour.
Returns
Neutrino oscillation probabilities

Definition at line 1272 of file PMNS_Base.cxx.

1273{
1274 ResetToFlavour(flvi);
1275
1276 Propagate();
1277
1278 return GetProbVector();
1279}
virtual vectorD GetProbVector()
Definition: PMNS_Base.cxx:1233

References OscProb::PMNS_Base::GetProbVector(), OscProb::PMNS_Base::Propagate(), and OscProb::PMNS_Base::ResetToFlavour().

◆ ProbVector() [2/6]

vectorD PMNS_Base::ProbVector ( int  flvi,
double  E 
)
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.
Parameters
flvi- The neutrino starting flavour.
E- The neutrino energy in GeV
Returns
Neutrino oscillation probability

Definition at line 1313 of file PMNS_Base.cxx.

1314{
1315 SetEnergy(E);
1316
1317 return ProbVector(flvi);
1318}

References OscProb::PMNS_Base::ProbVector(), and OscProb::PMNS_Base::SetEnergy().

◆ ProbVector() [3/6]

vectorD PMNS_Base::ProbVector ( int  flvi,
double  E,
double  L 
)
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.
Parameters
flvi- The neutrino starting flavour.
E- The neutrino energy in GeV
L- The neutrino path length in km
Returns
Neutrino oscillation probability

Definition at line 1365 of file PMNS_Base.cxx.

1366{
1367 SetEnergy(E);
1368 SetLength(L);
1369
1370 return ProbVector(flvi);
1371}

References OscProb::PMNS_Base::ProbVector(), OscProb::PMNS_Base::SetEnergy(), and OscProb::PMNS_Base::SetLength().

◆ ProbVector() [4/6]

vectorD PMNS_Base::ProbVector ( vectorC  nu_in)
virtualinherited

Compute the probabilities of nu_in going to all flavours

Compute the probability of nu_in going to all flavours.

Parameters
nu_in- The neutrino initial state in flavour basis.
Returns
Neutrino oscillation probabilities

Definition at line 1250 of file PMNS_Base.cxx.

1251{
1252 SetPureState(nu_in);
1253
1254 Propagate();
1255
1256 return GetProbVector();
1257}

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().

◆ ProbVector() [5/6]

vectorD PMNS_Base::ProbVector ( vectorC  nu_in,
double  E 
)
virtualinherited

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.

Parameters
nu_in- The neutrino initial state in flavour basis.
E- The neutrino energy in GeV
Returns
Neutrino oscillation probabilities

Definition at line 1291 of file PMNS_Base.cxx.

1292{
1293 SetEnergy(E);
1294
1295 return ProbVector(nu_in);
1296}

References OscProb::PMNS_Base::ProbVector(), and OscProb::PMNS_Base::SetEnergy().

◆ ProbVector() [6/6]

vectorD PMNS_Base::ProbVector ( vectorC  nu_in,
double  E,
double  L 
)
virtualinherited

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.

Parameters
nu_in- The neutrino initial state in flavour basis.
E- The neutrino energy in GeV
L- The neutrino path length in km
Returns
Neutrino oscillation probabilities

Definition at line 1336 of file PMNS_Base.cxx.

1337{
1338 SetEnergy(E);
1339 SetLength(L);
1340
1341 return ProbVector(nu_in);
1342}

References OscProb::PMNS_Base::ProbVector(), OscProb::PMNS_Base::SetEnergy(), and OscProb::PMNS_Base::SetLength().

◆ Propagate()

void PMNS_Base::Propagate ( )
protectedvirtualinherited

Propagate neutrino state through full path

Reimplemented in OscProb::PMNS_NUNM.

Definition at line 1018 of file PMNS_Base.cxx.

1019{
1020 for (int i = 0; i < int(fNuPaths.size()); i++) { PropagatePath(fNuPaths[i]); }
1021}

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().

◆ PropagatePath()

void PMNS_Base::PropagatePath ( NuPath  p)
protectedvirtualinherited

Propagate the current neutrino state through a given path

Parameters
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.

984{
985 // Set the neutrino path
986 SetCurPath(p);
987
988 // Solve for eigensystem
989 SolveHam();
990
991 double LengthIneV = kKm2eV * p.length;
992 for (int i = 0; i < fNumNus; i++) {
993 double arg = fEval[i] * LengthIneV;
994 fPhases[i] = complexD(cos(arg), -sin(arg));
995 }
996
997 for (int i = 0; i < fNumNus; i++) {
998 fBuffer[i] = 0;
999 for (int j = 0; j < fNumNus; j++) {
1000 fBuffer[i] += conj(fEvec[j][i]) * fNuState[j];
1001 }
1002 fBuffer[i] *= fPhases[i];
1003 }
1004
1005 // Propagate neutrino state
1006 for (int i = 0; i < fNumNus; i++) {
1007 fNuState[i] = 0;
1008 for (int j = 0; j < fNumNus; j++) {
1009 fNuState[i] += fEvec[i][j] * fBuffer[j];
1010 }
1011 }
1012}
std::complex< double > complexD
Definition: Definitions.h:21

References OscProb::PMNS_Base::fBuffer, OscProb::PMNS_Base::fEval, OscProb::PMNS_Base::fEvec, OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fNuState, OscProb::PMNS_Base::fPhases, OscProb::PMNS_Base::kKm2eV, OscProb::NuPath::length, OscProb::PMNS_Base::SetCurPath(), and OscProb::PMNS_Base::SolveHam().

Referenced by OscProb::PMNS_Base::ProbMatrix(), OscProb::PMNS_Base::Propagate(), OscProb::PMNS_Iter::PropagatePath(), and OscProb::PMNS_NUNM::PropagatePath().

◆ ResetToFlavour()

void PMNS_Base::ResetToFlavour ( int  flv)
protectedvirtualinherited

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.
Parameters
flv- The neutrino starting flavour.

Reimplemented in OscProb::PMNS_Deco.

Definition at line 1034 of file PMNS_Base.cxx.

1035{
1036 assert(flv >= 0 && flv < fNumNus);
1037 for (int i = 0; i < fNumNus; ++i) {
1038 if (i == flv)
1039 fNuState[i] = one;
1040 else
1041 fNuState[i] = zero;
1042 }
1043}
static const complexD one
one in complex
Definition: PMNS_Base.h:212

References OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fNuState, OscProb::PMNS_Base::one, and OscProb::PMNS_Base::zero.

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::PMNS_Base(), OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::ProbMatrix(), OscProb::PMNS_NUNM::ProbMatrix(), OscProb::PMNS_Base::ProbVector(), and OscProb::PMNS_Deco::ResetToFlavour().

◆ RotateH()

void PMNS_Base::RotateH ( int  i,
int  j,
matrixC Ham 
)
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

Parameters
i,j- the indices of the rotation ij
Ham- the Hamiltonian to be rotated

Definition at line 822 of file PMNS_Base.cxx.

823{
824 // Do nothing if angle is zero
825 if (fTheta[i][j] == 0) return;
826
827 double fSinBuffer = sin(fTheta[i][j]);
828 double fCosBuffer = cos(fTheta[i][j]);
829
830 double fHmsBufferD;
831 complexD fHmsBufferC;
832
833 // With Delta
834 if (i + 1 < j) {
835 complexD fExpBuffer = complexD(cos(fDelta[i][j]), -sin(fDelta[i][j]));
836
837 // General case
838 if (i > 0) {
839 // Top columns
840 for (int k = 0; k < i; k++) {
841 fHmsBufferC = Ham[k][i];
842
843 Ham[k][i] *= fCosBuffer;
844 Ham[k][i] += Ham[k][j] * fSinBuffer * conj(fExpBuffer);
845
846 Ham[k][j] *= fCosBuffer;
847 Ham[k][j] -= fHmsBufferC * fSinBuffer * fExpBuffer;
848 }
849
850 // Middle row and column
851 for (int k = i + 1; k < j; k++) {
852 fHmsBufferC = Ham[k][j];
853
854 Ham[k][j] *= fCosBuffer;
855 Ham[k][j] -= conj(Ham[i][k]) * fSinBuffer * fExpBuffer;
856
857 Ham[i][k] *= fCosBuffer;
858 Ham[i][k] += fSinBuffer * fExpBuffer * conj(fHmsBufferC);
859 }
860
861 // Nodes ij
862 fHmsBufferC = Ham[i][i];
863 fHmsBufferD = real(Ham[j][j]);
864
865 Ham[i][i] *= fCosBuffer * fCosBuffer;
866 Ham[i][i] +=
867 2 * fSinBuffer * fCosBuffer * real(Ham[i][j] * conj(fExpBuffer));
868 Ham[i][i] += fSinBuffer * Ham[j][j] * fSinBuffer;
869
870 Ham[j][j] *= fCosBuffer * fCosBuffer;
871 Ham[j][j] += fSinBuffer * fHmsBufferC * fSinBuffer;
872 Ham[j][j] -=
873 2 * fSinBuffer * fCosBuffer * real(Ham[i][j] * conj(fExpBuffer));
874
875 Ham[i][j] -= 2 * fSinBuffer * real(Ham[i][j] * conj(fExpBuffer)) *
876 fSinBuffer * fExpBuffer;
877 Ham[i][j] +=
878 fSinBuffer * fCosBuffer * (fHmsBufferD - fHmsBufferC) * fExpBuffer;
879 }
880 // First rotation on j (No top columns)
881 else {
882 // Middle rows and columns
883 for (int k = i + 1; k < j; k++) {
884 Ham[k][j] = -conj(Ham[i][k]) * fSinBuffer * fExpBuffer;
885
886 Ham[i][k] *= fCosBuffer;
887 }
888
889 // Nodes ij
890 fHmsBufferD = real(Ham[i][i]);
891
892 Ham[i][j] =
893 fSinBuffer * fCosBuffer * (Ham[j][j] - fHmsBufferD) * fExpBuffer;
894
895 Ham[i][i] *= fCosBuffer * fCosBuffer;
896 Ham[i][i] += fSinBuffer * Ham[j][j] * fSinBuffer;
897
898 Ham[j][j] *= fCosBuffer * fCosBuffer;
899 Ham[j][j] += fSinBuffer * fHmsBufferD * fSinBuffer;
900 }
901 }
902 // Without Delta (No middle rows or columns: j = i+1)
903 else {
904 // General case
905 if (i > 0) {
906 // Top columns
907 for (int k = 0; k < i; k++) {
908 fHmsBufferC = Ham[k][i];
909
910 Ham[k][i] *= fCosBuffer;
911 Ham[k][i] += Ham[k][j] * fSinBuffer;
912
913 Ham[k][j] *= fCosBuffer;
914 Ham[k][j] -= fHmsBufferC * fSinBuffer;
915 }
916
917 // Nodes ij
918 fHmsBufferC = Ham[i][i];
919 fHmsBufferD = real(Ham[j][j]);
920
921 Ham[i][i] *= fCosBuffer * fCosBuffer;
922 Ham[i][i] += 2 * fSinBuffer * fCosBuffer * real(Ham[i][j]);
923 Ham[i][i] += fSinBuffer * Ham[j][j] * fSinBuffer;
924
925 Ham[j][j] *= fCosBuffer * fCosBuffer;
926 Ham[j][j] += fSinBuffer * fHmsBufferC * fSinBuffer;
927 Ham[j][j] -= 2 * fSinBuffer * fCosBuffer * real(Ham[i][j]);
928
929 Ham[i][j] -= 2 * fSinBuffer * real(Ham[i][j]) * fSinBuffer;
930 Ham[i][j] += fSinBuffer * fCosBuffer * (fHmsBufferD - fHmsBufferC);
931 }
932 // First rotation (theta12)
933 else {
934 Ham[i][j] = fSinBuffer * fCosBuffer * Ham[j][j];
935
936 Ham[i][i] = fSinBuffer * Ham[j][j] * fSinBuffer;
937
938 Ham[j][j] *= fCosBuffer * fCosBuffer;
939 }
940 }
941}

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().

◆ RotateState()

void PMNS_Base::RotateState ( int  i,
int  j 
)
protectedvirtualinherited

Rotate the neutrino state by the angle theta_ij and phase delta_ij.

Parameters
i,j- the indices of the rotation ij

Definition at line 760 of file PMNS_Base.cxx.

761{
762 // Do nothing if angle is zero
763 if (fTheta[i][j] == 0) return;
764
765 double sij = sin(fTheta[i][j]);
766 double cij = cos(fTheta[i][j]);
767
768 complexD buffer;
769
770 if (i + 1 == j || fDelta[i][j] == 0) {
771 buffer = cij * fNuState[i] + sij * fNuState[j];
772 fNuState[j] = cij * fNuState[j] - sij * fNuState[i];
773 }
774 else {
775 complexD eij = complexD(cos(fDelta[i][j]), -sin(fDelta[i][j]));
776 buffer = cij * fNuState[i] + sij * eij * fNuState[j];
777 fNuState[j] = cij * fNuState[j] - sij * conj(eij) * fNuState[i];
778 }
779
780 fNuState[i] = buffer;
781}

References OscProb::PMNS_Base::fDelta, OscProb::PMNS_Base::fNuState, and OscProb::PMNS_Base::fTheta.

Referenced by OscProb::PMNS_Base::GetMassEigenstate().

◆ SetAngle()

void PMNS_Base::SetAngle ( int  i,
int  j,
double  th 
)
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.

Parameters
i,j- the indices of theta_ij
th- the value of theta_ij

Definition at line 539 of file PMNS_Base.cxx.

540{
541 if (i > j) {
542 cerr << "WARNING: First argument should be smaller than second argument"
543 << endl
544 << " Setting reverse order (Theta" << j << i << "). " << endl;
545 int temp = i;
546 i = j;
547 j = temp;
548 }
549 if (i < 1 || i > fNumNus - 1 || j < 2 || j > fNumNus) {
550 cerr << "ERROR: Theta" << i << j << " not valid for " << fNumNus
551 << " neutrinos. Doing nothing." << endl;
552 return;
553 }
554
555 // Check if value is actually changing
556 fBuiltHms *= (fTheta[i - 1][j - 1] == th);
557
558 fTheta[i - 1][j - 1] = th;
559}

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().

◆ SetAtt() [1/2]

void PMNS_Base::SetAtt ( double  att,
int  idx 
)
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.

Parameters
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.

365{
366 if (fNuPaths.size() != 1) {
367 cerr << "WARNING: Clearing path vector and starting new single path."
368 << endl
369 << "To avoid possible issues, use the SetPath function." << endl;
370
371 SetStdPath();
372 }
373
374 switch (idx) {
375 case 0: fNuPaths[0].length = att; break;
376 case 1: fNuPaths[0].density = att; break;
377 case 2: fNuPaths[0].zoa = att; break;
378 case 3: fNuPaths[0].layer = att; break;
379 }
380}

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().

◆ SetAtt() [2/2]

void PMNS_Base::SetAtt ( vectorD  att,
int  idx 
)
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.

Parameters
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.

428{
429 // Get the sizes of the attribute and
430 // path sequence vectors
431 int nA = att.size();
432 int nP = fNuPaths.size();
433
434 // If the vector sizes are equal, update this attribute
435 if (nA == nP) {
436 for (int i = 0; i < nP; i++) {
437 switch (idx) {
438 case 0: fNuPaths[i].length = att[i]; break;
439 case 1: fNuPaths[i].density = att[i]; break;
440 case 2: fNuPaths[i].zoa = att[i]; break;
441 case 3: fNuPaths[i].layer = att[i]; break;
442 }
443 }
444 }
445 // If the vector sizes differ, create a new path sequence
446 // and set value for this attribute. Other attributes will
447 // be taken from default single path.
448 else {
449 cerr << "WARNING: New vector size. Starting new path vector." << endl
450 << "To avoid possible issues, use the SetPath function." << endl;
451
452 // Start a new standard path just
453 // to set default values
454 SetStdPath();
455
456 // Create a path segment with default values
457 NuPath p = fNuPaths[0];
458
459 // Clear the path sequence
460 ClearPath();
461
462 // Set this particular attribute's value
463 // and add the path segment to the sequence
464 for (int i = 0; i < nA; i++) {
465 switch (idx) {
466 case 0: p.length = att[i]; break;
467 case 1: p.density = att[i]; break;
468 case 2: p.zoa = att[i]; break;
469 case 3: p.layer = att[i]; break;
470 }
471
472 AddPath(p);
473 }
474 }
475}
virtual void ClearPath()
Clear the path vector.
Definition: PMNS_Base.cxx:287
int layer
An index to identify the matter type.
Definition: NuPath.h:81
double density
The density of the path segment in g/cm^3.
Definition: NuPath.h:79

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.

◆ SetAvgProbPrec()

void PMNS_Base::SetAvgProbPrec ( double  prec)
virtualinherited

Set the precision for the AvgProb method

Parameters
prec- AvgProb precision

Definition at line 1962 of file PMNS_Base.cxx.

1963{
1964 if (prec < 1e-8) {
1965 cerr << "WARNING: Cannot set AvgProb precision lower that 1e-8."
1966 << "Setting to 1e-8." << endl;
1967 prec = 1e-8;
1968 }
1969 fAvgProbPrec = prec;
1970}

References OscProb::PMNS_Base::fAvgProbPrec.

Referenced by OscProb::PMNS_Base::PMNS_Base().

◆ SetCoupByIndex()

void PMNS_NSI::SetCoupByIndex ( double  c,
int  i 
)
virtual

Set the NSI relative coupling to a given fermion. This factor represents what fraction of the fermion are seen by the interaction.

Parameters
c- fermion coupling strength
i- fermion index (0,1,2) -> (e,u,d)

Definition at line 261 of file PMNS_NSI.cxx.

262{
263 bool isSame = (fNSIcoup[i] == c);
264
265 if (!isSame) ClearCache();
266
267 fGotES *= isSame;
268
269 fNSIcoup[i] = c;
270}

References OscProb::PMNS_Base::ClearCache(), OscProb::PMNS_Base::fGotES, and fNSIcoup.

Referenced by SetDownCoup(), SetElecCoup(), and SetUpCoup().

◆ SetCurPath()

void PMNS_Base::SetCurPath ( NuPath  p)
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.

Parameters
p- A neutrino path segment

Definition at line 274 of file PMNS_Base.cxx.

275{
276 // Check if relevant value are actually changing
277 fGotES *= (fPath.density == p.density);
278 fGotES *= (fPath.zoa == p.zoa);
279
280 fPath = p;
281}

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 OscProb::PMNS_Deco::PropagatePath().

◆ SetDelta()

void PMNS_Base::SetDelta ( int  i,
int  j,
double  delta 
)
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.

Parameters
i,j- the indices of delta_ij
delta- the value of delta_ij

Definition at line 602 of file PMNS_Base.cxx.

603{
604 if (i > j) {
605 cerr << "WARNING: First argument should be smaller than second argument"
606 << endl
607 << " Setting reverse order (Delta" << j << i << "). " << endl;
608 int temp = i;
609 i = j;
610 j = temp;
611 }
612 if (i < 1 || i > fNumNus - 1 || j < 2 || j > fNumNus) {
613 cerr << "ERROR: Delta" << i << j << " not valid for " << fNumNus
614 << " neutrinos. Doing nothing." << endl;
615 return;
616 }
617 if (i + 1 == j) {
618 cerr << "WARNING: Rotation " << i << j << " is real. Doing nothing."
619 << endl;
620 return;
621 }
622
623 // Check if value is actually changing
624 fBuiltHms *= (fDelta[i - 1][j - 1] == delta);
625
626 fDelta[i - 1][j - 1] = delta;
627}

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().

◆ SetDeltaMsqrs()

void PMNS_Fast::SetDeltaMsqrs ( double  dm21,
double  dm32 
)
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.

Parameters
dm21- The solar mass-splitting Dm_21
dm32- The atmospheric mass-splitting Dm_32

Definition at line 55 of file PMNS_Fast.cxx.

56{
57 SetDm(2, dm21);
58 SetDm(3, dm32 + dm21);
59}
virtual void SetDm(int j, double dm)
Set the mass-splitting dm_j1 in eV^2.
Definition: PMNS_Base.cxx:674

References OscProb::PMNS_Base::SetDm().

◆ SetDensity() [1/2]

void PMNS_Base::SetDensity ( double  rho)
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.

Parameters
rho- The density of the path segment in g/cm^3

Definition at line 402 of file PMNS_Base.cxx.

402{ SetAtt(rho, 1); }
virtual void SetAtt(double att, int idx)
Set one of the path attributes.
Definition: PMNS_Base.cxx:364

References OscProb::PMNS_Base::SetAtt().

◆ SetDensity() [2/2]

void PMNS_Base::SetDensity ( vectorD  rho)
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.

Parameters
rho- The densities of the path segments in g/cm^3

Definition at line 497 of file PMNS_Base.cxx.

497{ SetAtt(rho, 1); }

References OscProb::PMNS_Base::SetAtt().

◆ SetDm()

void PMNS_Base::SetDm ( int  j,
double  dm 
)
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.

Parameters
j- the index of dm_j1
dm- the value of dm_j1

Definition at line 674 of file PMNS_Base.cxx.

675{
676 if (j < 2 || j > fNumNus) {
677 cerr << "ERROR: Dm" << j << "1 not valid for " << fNumNus
678 << " neutrinos. Doing nothing." << endl;
679 return;
680 }
681
682 // Check if value is actually changing
683 fBuiltHms *= (fDm[j - 1] == dm);
684
685 fDm[j - 1] = dm;
686}

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().

◆ SetDownCoup()

void PMNS_NSI::SetDownCoup ( double  d)
virtual

Set the NSI relative coupling to d-quarks. This factor represents what fraction of the d-quarks are seen by the interaction.

Parameters
d- d-quark coupling strength

Definition at line 300 of file PMNS_NSI.cxx.

300{ SetCoupByIndex(d, 2); }
virtual void SetCoupByIndex(double c, int i)
Set a given fermion coupling.
Definition: PMNS_NSI.cxx:261

References SetCoupByIndex().

Referenced by SetFermCoup().

◆ SetElecCoup()

void PMNS_NSI::SetElecCoup ( double  e)
virtual

Set the NSI relative coupling to electrons. This factor represents what fraction of the electrons are seen by the interaction.

Parameters
e- electron coupling strength

Definition at line 280 of file PMNS_NSI.cxx.

280{ SetCoupByIndex(e, 0); }

References SetCoupByIndex().

Referenced by SetFermCoup().

◆ SetEnergy()

void PMNS_Base::SetEnergy ( double  E)
virtualinherited

Set neutrino energy in GeV.

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Parameters
E- The neutrino energy in GeV

Definition at line 226 of file PMNS_Base.cxx.

227{
228 // Check if value is actually changing
229 fGotES *= (fEnergy == E);
230
231 fEnergy = E;
232}

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().

◆ SetEps()

void PMNS_NSI::SetEps ( int  flvi,
int  flvj,
double  val,
double  phase 
)
virtual

Set any given NSI parameter.

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Flavours are:

  • 0 = nue, 1 = numu, 2 = nutau

Requires that flvi < flvj. Will notify you if input is wrong. If flvi > flvj, will assume reverse order and swap flvi and flvj.

Parameters
flvi- The first flavour index
flvj- The second flavour index
val- The absolute value of the parameter
phase- The complex phase of the parameter in radians

Definition at line 86 of file PMNS_NSI.cxx.

87{
88 if (flvi > flvj) {
89 cerr << "WARNING: First argument should be smaller or equal to second "
90 "argument"
91 << endl
92 << "Setting reverse order (Eps_" << flvj << flvi << "). " << endl;
93 int temp = flvi;
94 flvi = flvj;
95 flvj = temp;
96 }
97 if (flvi < 0 || flvi > 2 || flvj < flvi || flvj > 2) {
98 cerr << "WARNING: Eps_" << flvi << flvj << " not valid for " << fNumNus
99 << " neutrinos. Doing nothing." << endl;
100 return;
101 }
102
103 complexD h = val;
104
105 if (flvi != flvj) h *= complexD(cos(phase), sin(phase));
106
107 bool isSame = (fEps[flvi][flvj] == h);
108
109 if (!isSame) ClearCache();
110
111 fGotES *= isSame;
112
113 fEps[flvi][flvj] = h;
114}

References OscProb::PMNS_Base::ClearCache(), fEps, OscProb::PMNS_Base::fGotES, and OscProb::PMNS_Base::fNumNus.

Referenced by GetNSI(), GetSNSI(), SetEps_ee(), SetEps_emu(), SetEps_etau(), SetEps_mumu(), SetEps_mutau(), SetEps_tautau(), and SetNSI().

◆ SetEps_ee()

void PMNS_NSI::SetEps_ee ( double  a)
virtual

Set eps_ee parameter

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Parameters
a- The value of the real parameter eps_ee

Definition at line 158 of file PMNS_NSI.cxx.

158{ SetEps(0, 0, a, 0); }
virtual void SetEps(int flvi, int flvj, double val, double phase)
Set any given NSI parameter.
Definition: PMNS_NSI.cxx:86

References SetEps().

◆ SetEps_emu()

void PMNS_NSI::SetEps_emu ( double  a,
double  phi 
)
virtual

Set eps_emu parameter

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Parameters
a- The absolute value of the parameter eps_emu
phi- The phase of the complex parameter eps_emu in radians

Definition at line 192 of file PMNS_NSI.cxx.

192{ SetEps(0, 1, a, phi); }

References SetEps().

◆ SetEps_etau()

void PMNS_NSI::SetEps_etau ( double  a,
double  phi 
)
virtual

Set eps_etau parameter

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Parameters
a- The absolute value of the parameter eps_etau
phi- The phase of the complex parameter eps_etau in radians

Definition at line 204 of file PMNS_NSI.cxx.

204{ SetEps(0, 2, a, phi); }

References SetEps().

◆ SetEps_mumu()

void PMNS_NSI::SetEps_mumu ( double  a)
virtual

Set eps_mumu parameter

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Parameters
a- The value of the real parameter eps_mumu

Definition at line 169 of file PMNS_NSI.cxx.

169{ SetEps(1, 1, a, 0); }

References SetEps().

◆ SetEps_mutau()

void PMNS_NSI::SetEps_mutau ( double  a,
double  phi 
)
virtual

Set eps_mutau parameter

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Parameters
a- The absolute value of the parameter eps_mutau
phi- The phase of the complex parameter eps_mutau in radians

Definition at line 216 of file PMNS_NSI.cxx.

216{ SetEps(1, 2, a, phi); }

References SetEps().

◆ SetEps_tautau()

void PMNS_NSI::SetEps_tautau ( double  a)
virtual

Set eps_tautau parameter

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Parameters
a- The value of the real parameter eps_tautau

Definition at line 180 of file PMNS_NSI.cxx.

180{ SetEps(2, 2, a, 0); }

References SetEps().

◆ SetFermCoup()

void PMNS_NSI::SetFermCoup ( double  e,
double  u,
double  d 
)
virtual

Set the NSI relative couplings to each fermion. These factors represent what fraction of each fermion are seen by the interaction.

Parameters
e- electron coupling strength
u- u-quark coupling strength
d- d-quark coupling strength

Definition at line 312 of file PMNS_NSI.cxx.

313{
314 SetElecCoup(e); // electron coupling
315 SetUpCoup(u); // u-quark coupling
316 SetDownCoup(d); // d-quark coupling
317}
virtual void SetUpCoup(double u)
Set u-quark couling.
Definition: PMNS_NSI.cxx:290
virtual void SetElecCoup(double e)
Set electron coupling.
Definition: PMNS_NSI.cxx:280
virtual void SetDownCoup(double d)
Set d-quark couling.
Definition: PMNS_NSI.cxx:300

References SetDownCoup(), SetElecCoup(), and SetUpCoup().

Referenced by PMNS_NSI().

◆ SetIsNuBar()

void PMNS_Base::SetIsNuBar ( bool  isNuBar)
virtualinherited

Set anti-neutrino flag.

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Parameters
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.

244{
245 // Check if value is actually changing
246 fGotES *= (fIsNuBar == isNuBar);
247
248 fIsNuBar = isNuBar;
249}

References OscProb::PMNS_Base::fGotES, and OscProb::PMNS_Base::fIsNuBar.

Referenced by CheckProb(), OscProb::PMNS_Base::PMNS_Base(), and SaveTestFile().

◆ SetLayers()

void PMNS_Base::SetLayers ( std::vector< int >  lay)
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.

Parameters
lay- Indices to identify the layer types (e.g. earth inner core)

Definition at line 519 of file PMNS_Base.cxx.

520{
521 vectorD lay_double(lay.begin(), lay.end());
522
523 SetAtt(lay_double, 3);
524}

References OscProb::PMNS_Base::SetAtt().

◆ SetLength() [1/2]

void PMNS_Base::SetLength ( double  L)
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.

Parameters
L- The length of the path segment in km

Definition at line 391 of file PMNS_Base.cxx.

391{ SetAtt(L, 0); }

References OscProb::PMNS_Base::SetAtt().

Referenced by OscProb::PMNS_Base::Prob(), OscProb::PMNS_Base::ProbMatrix(), and OscProb::PMNS_Base::ProbVector().

◆ SetLength() [2/2]

void PMNS_Base::SetLength ( vectorD  L)
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.

Parameters
L- The lengths of the path segments in km

Definition at line 486 of file PMNS_Base.cxx.

486{ SetAtt(L, 0); }

References OscProb::PMNS_Base::SetAtt().

◆ SetMaxCache()

void PMNS_Base::SetMaxCache ( int  mc = 1e6)
virtualinherited

Set maximum number of cached eigensystems. Finding eigensystems can become slow and take up memory. This protects the cache from becoming too large.

Parameters
mc- Max cache size (default: 1e6)

Definition at line 128 of file PMNS_Base.cxx.

128{ fMaxCache = mc; }

References OscProb::PMNS_Base::fMaxCache.

◆ SetMix()

void PMNS_Fast::SetMix ( double  th12,
double  th23,
double  th13,
double  deltacp 
)
virtualinherited

Set all mixing parameters at once.

Parameters
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.

38{
39 SetAngle(1, 2, th12);
40 SetAngle(1, 3, th13);
41 SetAngle(2, 3, th23);
42 SetDelta(1, 3, deltacp);
43}
virtual void SetDelta(int i, int j, double delta)
Set the CP phase delta_ij.
Definition: PMNS_Base.cxx:602
virtual void SetAngle(int i, int j, double th)
Set the mixing angle theta_ij.
Definition: PMNS_Base.cxx:539

References OscProb::PMNS_Base::SetAngle(), and OscProb::PMNS_Base::SetDelta().

◆ SetNSI()

void PMNS_NSI::SetNSI ( double  eps_ee,
double  eps_emu,
double  eps_etau,
double  eps_mumu,
double  eps_mutau,
double  eps_tautau,
double  delta_emu = 0,
double  delta_etau = 0,
double  delta_mutau = 0 
)

Set all NSI parameters at once.

This will check if value is changing to keep track of whether the eigensystem needs to be recomputed.

Parameters
eps_ee- The real parameter eps_ee
eps_emu- The absolute value of the complex parameter eps_emu
eps_etau- The absolute value of the complex parameter eps_etau
eps_mumu- The real parameter eps_mumu
eps_mutau- The absolute value of the complex parameter eps_mutau
eps_tautau- The real parameter eps_tautau
delta_emu- The phase of the complex parameter eps_emu in radians
delta_etau- The phase of the complex parameter eps_etau in radians
delta_mutau- The phase of the complex parameter eps_mutau in radians

Definition at line 55 of file PMNS_NSI.cxx.

58{
59 SetEps(0, 0, eps_ee, 0);
60 SetEps(1, 1, eps_mumu, 0);
61 SetEps(2, 2, eps_tautau, 0);
62
63 SetEps(0, 1, eps_emu, delta_emu);
64 SetEps(0, 2, eps_etau, delta_etau);
65 SetEps(1, 2, eps_mutau, delta_mutau);
66}

References SetEps().

Referenced by PMNS_NSI().

◆ SetPath() [1/3]

void PMNS_Base::SetPath ( double  length,
double  density,
double  zoa = 0.5,
int  layer = 0 
)
virtualinherited

Set a single path defining attributes directly.

This destroys the current path sequence and creates a new first path.

Parameters
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.

348{
349 SetPath(NuPath(length, density, zoa, layer));
350}
virtual void SetPath(NuPath p)
Set a single path.
Definition: PMNS_Base.cxx:330

References OscProb::PMNS_Base::SetPath().

◆ SetPath() [2/3]

void PMNS_Base::SetPath ( NuPath  p)
virtualinherited

Set a single path.

This destroys the current path sequence and creates a new first path.

Parameters
p- A neutrino path segment

Definition at line 330 of file PMNS_Base.cxx.

331{
332 ClearPath();
333 AddPath(p);
334}

References OscProb::PMNS_Base::AddPath(), and OscProb::PMNS_Base::ClearPath().

Referenced by OscProb::PMNS_Base::SetPath(), OscProb::PMNS_Base::SetStdPath(), and SetTestPath().

◆ SetPath() [3/3]

void PMNS_Base::SetPath ( std::vector< NuPath paths)
virtualinherited

Set vector of neutrino paths.

Parameters
paths- A sequence of neutrino paths

Definition at line 294 of file PMNS_Base.cxx.

294{ fNuPaths = paths; }

References OscProb::PMNS_Base::fNuPaths.

◆ SetPureState()

void PMNS_Base::SetPureState ( vectorC  nu_in)
protectedvirtualinherited

Set the initial state from a pure state

Parameters
nu_in- The neutrino initial state in flavour basis.

Reimplemented in OscProb::PMNS_Deco.

Definition at line 1070 of file PMNS_Base.cxx.

1071{
1072 assert(nu_in.size() == fNumNus);
1073
1074 fNuState = nu_in;
1075}

References OscProb::PMNS_Base::fNumNus, and OscProb::PMNS_Base::fNuState.

Referenced by OscProb::PMNS_Base::Prob(), and OscProb::PMNS_Base::ProbVector().

◆ SetStdPars()

void PMNS_Base::SetStdPars ( )
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.

178{
179 if (fNumNus > 2) {
180 // PDG values for 3 neutrinos
181 // Also applicable for 3+N neutrinos
182 SetAngle(1, 2, asin(sqrt(0.304)));
183 SetAngle(1, 3, asin(sqrt(0.0219)));
184 SetAngle(2, 3, asin(sqrt(0.514)));
185 SetDm(2, 7.53e-5);
186 SetDm(3, 2.52e-3);
187 }
188 else if (fNumNus == 2) {
189 // Effective muon disappearance values
190 // for two-flavour approximation
191 SetAngle(1, 2, 0.788);
192 SetDm(2, 2.47e-3);
193 }
194}

References OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::SetAngle(), and OscProb::PMNS_Base::SetDm().

Referenced by OscProb::PMNS_Base::PMNS_Base().

◆ SetStdPath()

void PMNS_Base::SetStdPath ( )
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.

206{
207 NuPath p;
208
209 p.length = 1000; // 1000 km default
210 p.density = 2.8; // Crust density
211 p.zoa = 0.5; // Crust Z/A
212 p.layer = 0; // Single layer
213
214 SetPath(p);
215}

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(), OscProb::PMNS_Deco::PMNS_Deco(), OscProb::PMNS_LIV::PMNS_LIV(), PMNS_NSI(), OscProb::PMNS_NUNM::PMNS_NUNM(), and OscProb::PMNS_Base::SetAtt().

◆ SetUpCoup()

void PMNS_NSI::SetUpCoup ( double  u)
virtual

Set the NSI relative coupling to u-quarks. This factor represents what fraction of the u-quarks are seen by the interaction.

Parameters
u- u-quark coupling strength

Definition at line 290 of file PMNS_NSI.cxx.

290{ SetCoupByIndex(u, 1); }

References SetCoupByIndex().

Referenced by SetFermCoup().

◆ SetUseCache()

void PMNS_Base::SetUseCache ( bool  u = true)
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.

Parameters
u- flag to set caching on (default: true)

Definition at line 105 of file PMNS_Base.cxx.

105{ fUseCache = u; }

References OscProb::PMNS_Base::fUseCache.

Referenced by OscProb::PMNS_Base::PMNS_Base().

◆ SetVacuumEigensystem()

void PMNS_Fast::SetVacuumEigensystem ( )
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.

144{
145 double s12, s23, s13, c12, c23, c13;
146 complexD idelta(0.0, fDelta[0][2]);
147 if (fIsNuBar) idelta = conj(idelta);
148
149 s12 = sin(fTheta[0][1]);
150 s23 = sin(fTheta[1][2]);
151 s13 = sin(fTheta[0][2]);
152 c12 = cos(fTheta[0][1]);
153 c23 = cos(fTheta[1][2]);
154 c13 = cos(fTheta[0][2]);
155
156 fEvec[0][0] = c12 * c13;
157 fEvec[0][1] = s12 * c13;
158 fEvec[0][2] = s13 * exp(-idelta);
159
160 fEvec[1][0] = -s12 * c23 - c12 * s23 * s13 * exp(idelta);
161 fEvec[1][1] = c12 * c23 - s12 * s23 * s13 * exp(idelta);
162 fEvec[1][2] = s23 * c13;
163
164 fEvec[2][0] = s12 * s23 - c12 * c23 * s13 * exp(idelta);
165 fEvec[2][1] = -c12 * s23 - s12 * c23 * s13 * exp(idelta);
166 fEvec[2][2] = c23 * c13;
167
168 fEval[0] = 0;
169 fEval[1] = fDm[1] / (2 * kGeV2eV * fEnergy);
170 fEval[2] = fDm[2] / (2 * kGeV2eV * fEnergy);
171}

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().

◆ SetZoA() [1/2]

void PMNS_Base::SetZoA ( double  zoa)
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.

Parameters
zoa- The effective Z/A of the path segment

Definition at line 413 of file PMNS_Base.cxx.

413{ SetAtt(zoa, 2); }

References OscProb::PMNS_Base::SetAtt().

◆ SetZoA() [2/2]

void PMNS_Base::SetZoA ( vectorD  zoa)
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.

Parameters
zoa- The effective Z/A of the path segments

Definition at line 508 of file PMNS_Base.cxx.

508{ SetAtt(zoa, 2); }

References OscProb::PMNS_Base::SetAtt().

◆ SolveHam()

void PMNS_Fast::SolveHam ( )
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.

101{
102 // Do vacuum oscillation in low density
103 if (fPath.density < 1.0e-6) {
105 return;
106 }
107
108 // Build Hamiltonian
109 BuildHms();
110
111 // Check if anything changed
112 if (fGotES) return;
113
114 // Try caching if activated
115 if (TryCache()) return;
116
117 UpdateHam();
118
119 double fEvalGLoBES[3];
120 complexD fEvecGLoBES[3][3];
121
122 // Solve Hamiltonian for eigensystem using the GLoBES method
123 zheevh3(fHam, fEvecGLoBES, fEvalGLoBES);
124
125 // Fill fEval and fEvec vectors from GLoBES arrays
126 for (int i = 0; i < fNumNus; i++) {
127 fEval[i] = fEvalGLoBES[i];
128 for (int j = 0; j < fNumNus; j++) { fEvec[i][j] = fEvecGLoBES[i][j]; }
129 }
130
131 fGotES = true;
132
133 // Fill cache if activated
134 FillCache();
135}
virtual void FillCache()
Cache the current eigensystem.
Definition: PMNS_Base.cxx:157
virtual bool TryCache()
Try to find a cached eigensystem.
Definition: PMNS_Base.cxx:134
virtual void BuildHms()
Build the matrix of masses squared.
Definition: PMNS_Base.cxx:955
virtual void UpdateHam()
Build the full Hamiltonian.
Definition: PMNS_Fast.cxx:69
virtual void SetVacuumEigensystem()
Set the eigensystem to the analytic solution of the vacuum Hamiltonian.
Definition: PMNS_Fast.cxx:143
complexD fHam[3][3]
The full hamiltonian.
Definition: PMNS_Fast.h:62
int zheevh3(std::complex< double > A[3][3], std::complex< double > Q[3][3], double w[3])
Definition: zheevh3.cxx:31

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 OscProb::PMNS_Deco::PropagatePath().

◆ TryCache()

bool PMNS_Base::TryCache ( )
protectedvirtualinherited

Try to find a cached version of this eigensystem.

Definition at line 134 of file PMNS_Base.cxx.

135{
136 if (fUseCache && !fMixCache.empty()) {
138
139 unordered_set<EigenPoint>::iterator it = fMixCache.find(fProbe);
140
141 if (it != fMixCache.end()) {
142 for (int i = 0; i < fNumNus; i++) {
143 fEval[i] = (*it).fEval[i] * (*it).fEnergy / fEnergy;
144 for (int j = 0; j < fNumNus; j++) { fEvec[i][j] = (*it).fEvec[i][j]; }
145 }
146 return true;
147 }
148 }
149
150 return false;
151}

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().

◆ UpdateHam()

void PMNS_NSI::UpdateHam ( )
protectedvirtual

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 each flavour pair.

Reimplemented from OscProb::PMNS_Fast.

Reimplemented in OscProb::PMNS_SNSI.

Definition at line 226 of file PMNS_NSI.cxx.

227{
228 double lv = 2 * kGeV2eV * fEnergy; // 2*E in eV
229
230 double kr2GNe = kK2 * M_SQRT2 * kGf * fPath.density;
231 double kr2GNnsi = kr2GNe;
232
233 kr2GNe *= fPath.zoa; // Std matter potential in eV
234 kr2GNnsi *= GetZoACoup(); // NSI matter potential in eV
235
236 // Finish build Hamiltonian in matter with dimension of eV
237 for (int i = 0; i < fNumNus; i++) {
238 for (int j = i; j < fNumNus; j++) {
239 if (!fIsNuBar)
240 fHam[i][j] = fHms[i][j] / lv + kr2GNnsi * fEps[i][j];
241 else
242 fHam[i][j] = conj(fHms[i][j] / lv - kr2GNnsi * fEps[i][j]);
243 }
244 }
245
246 if (!fIsNuBar)
247 fHam[0][0] += kr2GNe;
248 else
249 fHam[0][0] -= kr2GNe;
250}
static const double kK2
mol/GeV^2/cm^3 to eV
Definition: PMNS_Base.h:216
static const double kGf
G_F in units of GeV^-2.
Definition: PMNS_Base.h:220
virtual double GetZoACoup()
Get effective Z/A coupling.
Definition: PMNS_NSI.cxx:341

References OscProb::NuPath::density, OscProb::PMNS_Base::fEnergy, fEps, OscProb::PMNS_Fast::fHam, OscProb::PMNS_Base::fHms, OscProb::PMNS_Base::fIsNuBar, OscProb::PMNS_Base::fNumNus, OscProb::PMNS_Base::fPath, GetZoACoup(), OscProb::PMNS_Base::kGeV2eV, OscProb::PMNS_Base::kGf, OscProb::PMNS_Base::kK2, and OscProb::NuPath::zoa.

Member Data Documentation

◆ fAvgProbPrec

double OscProb::PMNS_Base::fAvgProbPrec
protectedinherited

◆ fBuffer

vectorC OscProb::PMNS_Base::fBuffer
protectedinherited

◆ fBuiltHms

◆ fCachePrec

double OscProb::PMNS_Base::fCachePrec
protectedinherited

Definition at line 302 of file PMNS_Base.h.

◆ fDelta

◆ fDm

◆ fEnergy

◆ fEps

complexD OscProb::PMNS_NSI::fEps[3][3]
protected

Definition at line 76 of file PMNS_NSI.h.

Referenced by GetEps(), SetEps(), UpdateHam(), and OscProb::PMNS_SNSI::UpdateHam().

◆ fEval

◆ fEvec

◆ fGotES

◆ fHam

◆ fHms

◆ fIsNuBar

◆ fMaxCache

int OscProb::PMNS_Base::fMaxCache
protectedinherited

Definition at line 303 of file PMNS_Base.h.

Referenced by OscProb::PMNS_Base::FillCache(), and OscProb::PMNS_Base::SetMaxCache().

◆ fMixCache

std::unordered_set<EigenPoint> OscProb::PMNS_Base::fMixCache
protectedinherited

◆ fNSIcoup

double OscProb::PMNS_NSI::fNSIcoup[3]
protected

Definition at line 78 of file PMNS_NSI.h.

Referenced by GetDownCoup(), GetElecCoup(), GetUpCoup(), GetZoACoup(), and SetCoupByIndex().

◆ fNumNus

int OscProb::PMNS_Base::fNumNus
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(), GetEps(), OscProb::PMNS_Deco::GetGamma(), OscProb::PMNS_Base::GetMassEigenstate(), OscProb::PMNS_Base::GetProbVector(), OscProb::PMNS_Base::GetSamplePoints(), OscProb::PMNS_Base::P(), OscProb::PMNS_Deco::P(), OscProb::PMNS_Base::PMNS_Base(), OscProb::PMNS_Decay::PMNS_Decay(), OscProb::PMNS_Base::ProbMatrix(), OscProb::PMNS_NUNM::ProbMatrix(), OscProb::PMNS_Deco::ProbMatrix(), OscProb::PMNS_Base::PropagatePath(), OscProb::PMNS_Decay::PropagatePath(), OscProb::PMNS_Deco::PropagatePath(), OscProb::PMNS_Iter::PropagatePath(), OscProb::PMNS_Base::ResetToFlavour(), OscProb::PMNS_Deco::ResetToFlavour(), OscProb::PMNS_Deco::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(), SetEps(), OscProb::PMNS_Deco::SetGamma(), OscProb::PMNS_Base::SetPureState(), OscProb::PMNS_Deco::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(), UpdateHam(), OscProb::PMNS_NUNM::UpdateHam(), OscProb::PMNS_SNSI::UpdateHam(), and OscProb::PMNS_Sterile::UpdateHam().

◆ fNuPaths

◆ fNuState

◆ fPath

◆ fPhases

vectorC OscProb::PMNS_Base::fPhases
protectedinherited

Definition at line 286 of file PMNS_Base.h.

Referenced by OscProb::PMNS_Base::PropagatePath().

◆ fProbe

EigenPoint OscProb::PMNS_Base::fProbe
protectedinherited

Definition at line 308 of file PMNS_Base.h.

Referenced by OscProb::PMNS_Base::FillCache(), and OscProb::PMNS_Base::TryCache().

◆ fTheta

◆ fUseCache

bool OscProb::PMNS_Base::fUseCache
protectedinherited

◆ kGeV2eV

◆ kGf

◆ kK2

const double PMNS_Base::kK2
staticprotectedinherited

◆ kKm2eV

◆ kNA

const double PMNS_Base::kNA = 6.022140857e23
staticprotectedinherited

Definition at line 218 of file PMNS_Base.h.

◆ one

const complexD PMNS_Base::one
staticprotectedinherited

◆ zero


The documentation for this class was generated from the following files: