C++ API Reference

class ANFlux : public nuflux::FluxFunction

Public Functions

ANFlux(const std::string &fluxName)
virtual double getFlux(ParticleType type, double energy, double cosZenith) const

Computes the expected flux for neutrinos of the given type, energy, and zenith angle.

virtual double getMinEnergy() const
virtual double getMaxEnergy() const

Public Static Functions

static boost::shared_ptr<FluxFunction> makeFlux(const std::string &fluxName)

Private Members

double emin_
double emax_
std::map<ParticleType, boost::shared_ptr<Evaluator>> fluxes_
class component

Public Functions

double getFlux(double energy, double cosZenith) const
void adjustPionKaonRatio(double pionAdjust, double kaonAdjust)

Private Functions

double inclination(double cosZenith) const

A correction for the mean interaction length of mesons at production altitude The form of this correction is derived in a report by D. Chirkin.

double transitionEnergy(double cosZenith) const

Private Members

double lowECoeffs[15]

These are the coefficients for the degree 5 polynomial in cos(zenith) and log10(energy) used for low energies

double lowECutoffs[4]

These are the maximum energies for which the low energy polynomial is to be used.

progenitorComponent directComponents[2]

These are (usually?) the contributions form pions and kaons.

progenitorComponent modifiedComponents[3]

Like the ‘direct’ components, but multiplied by the ‘modification’.

double highEIndex

Base powerlaw index at high energies.

double norm

Overall normalization factor.

double inclinationCoeffs[5]

The coefficients for the inclination correction to the cos(zenith)

double normTweaks[6]

Slight modifications to the normalization of the high energy calculation These correspond to the predefined cos(zenith) ranges: [0,.1), [.1,.2), [.2,.4), [.4,.8), [.8,1)

double modificationNorm

The normalization for the ‘modification’ term.

double modificationScale

The scale log10(energy) for the ‘modification’ term.

double pionAdjust

the factor by which to boost the pion decay contribution

double kaonAdjust

the factor by which to boost the kaon decay contribution

double pionKaonNorm

the additional normalization factor necessary to maintain the same total flux over all angles and energies after adjusting the pion/kaon ratio

Friends

friend std::istream &operator>>(std::istream&, component&)
class component

Public Functions

double getFlux(double energy, double cosZenith) const

Public Members

double eMin

Private Members

double coeffs[15]

These are the coefficients for the degree 5 polynomial in cos(zenith) and log10(energy)

Friends

friend std::istream &operator>>(std::istream&, component&)
class CubicSpline

Public Functions

inline CubicSpline()
inline CubicSpline(const std::vector<double> &x, const std::vector<double> &y)
inline double operator()(double x)
inline double derivative(double x) const

Private Functions

inline const splineSegment &findSegment(double x) const
Pre:

x is within the domain of the spline

Private Members

std::vector<splineSegment> data
double xlast
double mlast
double blast
struct dumbHistogram

A really simple histogram. Must be filled from left to right.

Public Functions

inline explicit dumbHistogram(double firstEdge)

Construct a histogram initialized with the lower edge of its first bin.

inline void push_back(double edge, double value)

Append the upper edge of the next bin and the value contained in that bin.

Pre:

edge > all edges previously added

inline double binWidth(unsigned int i) const
Returns:

the width of bin i

inline double value(unsigned int i) const
Returns:

the value in bin i

inline std::pair<std::vector<double>, std::vector<double>> accumulate() const

Compute the integral of the histogram

Returns:

a pair of vectors, which contain coordinates (bin edges) and accumulated values, respectively

Private Members

std::vector<double> binEdges
std::vector<double> values
class Evaluator

Public Functions

Evaluator(const std::string&)
~Evaluator()
inline size_t GetNDim() const
std::pair<double, double> GetExtents(size_t dim) const
double operator()(double energy, double cosZenith) const

Private Functions

Evaluator(const Evaluator&)
Evaluator &operator=(const Evaluator&)

Private Members

photospline::splinetable data
struct FluxFactory

Public Functions

inline FluxFactory(boost::shared_ptr<FluxFunction> (*factoryFn)(const std::string&), const std::string &deprecationReason)
inline boost::shared_ptr<FluxFunction> operator()(const std::string &name) const

Public Members

boost::shared_ptr<FluxFunction> (*factoryFn)(const std::string&)
std::string deprecationReason
class FluxFunction
#include <FluxFunction.h>

The interface for all neutrino fluxes.

Subclassed by nuflux::ANFlux, nuflux::IntegralPreservingFlux, nuflux::LegacyConventionalFlux, nuflux::LegacyPromptFlux, nuflux::SimpleSplineFlux, nuflux::SplineFlux2

Public Functions

inline FluxFunction(std::string name)
inline virtual ~FluxFunction()
inline std::string getName() const
virtual double getFlux(ParticleType type, double energy, double cosZenith) const = 0

Computes the expected flux for neutrinos of the given type, energy, and zenith angle.

inline virtual double getFlux(ParticleType type, double energy, double azimuth, double cosZenith) const

Computes the expected flux for neutrinos of the given type, energy, azimuth, and zenith angle If the derived classes do not have a specific azimuth dependence then just use the normal flux call.

inline virtual double getMinEnergy() const
inline virtual double getMaxEnergy() const
inline virtual double readExtents(ParticleType type) const

Protected Attributes

std::string name
struct FluxRegisterererer

Public Functions

inline FluxRegisterererer(const std::string &name, boost::shared_ptr<FluxFunction> (*factoryFn)(const std::string&))
inline FluxRegisterererer(const std::string &name, boost::shared_ptr<FluxFunction> (*factoryFn)(const std::string&), const std::string &deprecationReason)
class IntegralPreservingFlux : public nuflux::FluxFunction

Public Functions

IntegralPreservingFlux(const std::string &fluxName)
virtual ~IntegralPreservingFlux()
virtual double getMinEnergy() const
virtual double getMaxEnergy() const
virtual double getFlux(ParticleType type, double energy, double cosZenith) const

Computes the expected flux for neutrinos of the given type, energy, and zenith angle.

virtual double getFlux(ParticleType type, double energy, double azimuth, double cosZenith) const

Computes the expected flux for neutrinos of the given type, energy, azimuth, and zenith angle If the derived classes do not have a specific azimuth dependence then just use the normal flux call.

Public Members

bool SetSubGeV

Public Static Functions

static boost::shared_ptr<FluxFunction> makeFlux(const std::string &fluxName)

Private Functions

void loadTables(const std::string &fluxName, ParticleType type)
double evaluate2D(ParticleType type, double energy, double cosZenith) const
double InterpolateAzimuth(ParticleType type, double energy, double azimuth, double cosZenith) const
double InterpolateCZFlux(ParticleType type, double energy, double azimuth, double cosZenith) const

Private Members

std::map<ParticleType, std::map<double, CubicSpline>> energySplines2D
std::map<ParticleType, std::map<std::pair<double, double>, CubicSpline>> energySplines3D
class kneeFunction

Public Functions

kneeFunction()
double operator()(double energy) const

Private Members

double c0
double c1
double x0
double e0

Friends

friend std::istream &operator>>(std::istream&, kneeFunction&)
struct KneeRegisterererer

Public Functions

inline KneeRegisterererer(const std::string baseModel, const std::string name)
class KneeReweightable

Subclassed by nuflux::LegacyConventionalFlux, nuflux::LegacyPromptFlux

Public Functions

inline KneeReweightable(std::string kneeName = "none")
virtual void setKneeReweightingModel(std::string reweightModel) = 0
inline std::string getKneeReweightingModel() const
inline virtual ~KneeReweightable()

Protected Attributes

std::string kneeCorrectionName
class kneeSpline

Public Functions

kneeSpline()

Constructs a kneeSpline whose value is 1 for all energies.

double operator()(double energy) const

Computes the correction factor for this knee model, relative to the base model

Parameters:

energy – the neutrino energy in GeV

Private Members

std::vector<double> x
std::vector<double> y
std::vector<double> b
std::vector<double> c
std::vector<double> d

Friends

friend std::istream &operator>>(std::istream&, kneeSpline&)
class LegacyConventionalFlux : public nuflux::FluxFunction, public nuflux::KneeReweightable, public nuflux::PionKaonAdjustable

Public Functions

LegacyConventionalFlux(const std::string &fluxName)
virtual double getFlux(ParticleType type, double energy, double cosZenith) const

Computes the expected flux for neutrinos of the given type, energy, and zenith angle.

virtual double getMinEnergy() const
virtual double getMaxEnergy() const
virtual void setRelativePionContribution(double adjust)

boosts the contribution from pion decays

Parameters:

adjust – the factor by which to increase the pion contribution

virtual void setRelativeKaonContribution(double adjust)

boosts the contribution from kaon decays

Parameters:

adjust – the factor by which to increase the kaon contribution

virtual void setKneeReweightingModel(std::string reweightModel)

Public Static Functions

static boost::shared_ptr<FluxFunction> makeFlux(const std::string &fluxName)

Private Members

std::map<ParticleType, component> components
kneeSpline kneeCorrection

Friends

friend std::istream &operator>>(std::istream&, component&)
friend component readConvComponent(const std::string &path)
class LegacyPromptFlux : public nuflux::FluxFunction, public nuflux::KneeReweightable

Public Functions

virtual double getMinEnergy() const
virtual double getMaxEnergy() const
LegacyPromptFlux(const std::string &fluxName)
virtual double getFlux(ParticleType type, double energy, double cosZenith) const

Computes the expected flux for neutrinos of the given type, energy, and zenith angle.

virtual void setKneeReweightingModel(std::string reweightModel)

Public Static Functions

static boost::shared_ptr<FluxFunction> makeFlux(const std::string &fluxName)

Private Members

std::map<ParticleType, component> components
kneeFunction kneeCorrection

Friends

friend std::istream &operator>>(std::istream&, component&)
friend component readPromptComponent(const std::string &path)
friend std::istream &operator>>(std::istream&, kneeFunction&)
class PionKaonAdjustable

Subclassed by nuflux::LegacyConventionalFlux

Public Functions

inline PionKaonAdjustable(double pionAdjust = 1.0, double kaonAdjust = 1.0)
virtual void setRelativePionContribution(double adjust) = 0

boosts the contribution from pion decays

Parameters:

adjust – the factor by which to increase the pion contribution

inline double getRelativePionContribution() const
virtual void setRelativeKaonContribution(double adjust) = 0

boosts the contribution from kaon decays

Parameters:

adjust – the factor by which to increase the kaon contribution

inline double getRelativeKaonContribution() const
inline virtual ~PionKaonAdjustable()

Protected Attributes

double pionAdjust
double kaonAdjust
struct progenitorComponent

A component for an equation like Gaisser 6.7 or 7.5.

Public Functions

double operator()(double energy, double inclination) const

Public Members

double A
double B
double E_c
bool kaonic
class SimpleSplineFlux : public nuflux::FluxFunction

Public Functions

SimpleSplineFlux(const std::string &fluxName)
virtual double getFlux(ParticleType type, double energy, double cosZenith) const

Computes the expected flux for neutrinos of the given type, energy, and zenith angle.

virtual double getMinEnergy() const
virtual double getMaxEnergy() const
virtual double readExtents(ParticleType type) const

Public Static Functions

static boost::shared_ptr<FluxFunction> makeFlux(const std::string &fluxName)

Private Members

boost::shared_ptr<photospline::splinetable<>> spline
class SplineFlux2 : public nuflux::FluxFunction

Public Functions

SplineFlux2(const std::string &fluxName)
virtual double getFlux(ParticleType type, double energy, double cosZenith) const

Computes the expected flux for neutrinos of the given type, energy, and zenith angle.

virtual double readExtents(ParticleType type) const
inline bool PathExist(const std::string &name)
virtual double getMinEnergy() const
virtual double getMaxEnergy() const

Public Static Functions

static boost::shared_ptr<FluxFunction> makeFlux(const std::string &fluxName)

Private Members

std::map<ParticleType, boost::shared_ptr<photospline::splinetable<>>> components
struct splineSegment

Public Functions

inline double operator()(double x) const
inline double derivative(double x) const

Public Members

double x
double a
double b
double c
double d
namespace nuflux

Enums

enum ParticleType

Values:

enumerator NuE
enumerator NuEBar
enumerator NuMu
enumerator NuMuBar
enumerator NuTau
enumerator NuTauBar

Functions

bool isNeutrino(ParticleType pdgid)

Determines if particle type is a neutrino.

boost::shared_ptr<FluxFunction> makeFlux(const std::string &fluxName)

Constructs a FluxFunction object for the named flux \fluxName the name of the flux being requested

Throws:

std::invalid_argument – if fluxName does not correspond to a known flux

std::vector<std::string> availableFluxes()

Returns a list of the names of the supported flux models (which may be used when calling makeFlux)

std::vector<std::string> kneesForFlux(std::string fluxName)

Returns a list of the names of the supported cosmic ray knee reweighting models for the given base neutrino flux

Parameters:

fluxName – the base neutrino flux

void printModels()

Print an listing of all supported models.

std::string getVersion()

Return the current nuflux version.

LegacyConventionalFlux::component readConvComponent(const std::string &fname)
bool classifyCriticalEnergy(double E_c)
std::istream &operator>>(std::istream &is, LegacyConventionalFlux::component::progenitorComponent &c)
std::istream &operator>>(std::istream &is, LegacyConventionalFlux::component &c)
std::istream &operator>>(std::istream &is, LegacyConventionalFlux::kneeSpline &s)
LegacyPromptFlux::component readPromptComponent(const std::string &path)
std::istream &operator>>(std::istream &is, LegacyPromptFlux::component &c)
std::istream &operator>>(std::istream &is, LegacyPromptFlux::kneeFunction &kf)
namespace detail

Things which are useful or necessary internally for the library but which aren’t useful for users.

Functions

void registerNeutrinoFlux(const std::string &name, boost::shared_ptr<FluxFunction> (*factoryFn)(const std::string&))

Add a flux to the registry.

void registerNeutrinoFlux(const std::string &name, boost::shared_ptr<FluxFunction> (*factoryFn)(const std::string&), const std::string &deprecationReason)

Add a flux to the registry with a warning that it is deprecated.

std::string getDataPath(std::string fname)
void registerKneeModel(const std::string baseModel, const std::string name)
const char *getBaseDir()

Variables

std::map<std::string, FluxFactory> *registry = NULL
std::multimap<std::string, std::string> *kneeRegistry = NULL
file ANFlux.h
#include <photospline/splinetable.h>
#include <nuflux/nuflux.h>
file Fluxes.h
file FluxFunction.h
#include <cmath>
#include <string>
#include <stdexcept>
#include <boost/preprocessor/cat.hpp>
#include <boost/shared_ptr.hpp>
#include <nuflux/nuflux.h>

Defines

NUFLUX_ParticleType
file Interfaces.h
#include <string>
file IPLEFlux.h
#include <nuflux/nuflux.h>
file LegacyConventionalFlux.h
#include <nuflux/nuflux.h>
#include <iosfwd>
#include <map>
file LegacyPromptFlux.h
#include <nuflux/nuflux.h>
#include <iosfwd>
#include <map>
file nuflux.h
#include <vector>
#include <nuflux/Fluxes.h>

Functions

std::string getDataPath(std::string fname)
file SplineFlux.h
#include <photospline/splinetable.h>
#include <nuflux/nuflux.h>
file SplineFlux2.h
#include <photospline/splinetable.h>
#include <nuflux/nuflux.h>
file ANFlux.cpp
#include <boost/lexical_cast.hpp>
#include <boost/make_shared.hpp>
#include “nuflux/ANFlux.h
file detail.cpp
#include <string>
#include “config.h”
file FluxFunction.cpp
#include <nuflux/nuflux.h>
#include <iostream>
#include “config.h”
file IPLEFlux.cpp
#include <iostream>
#include <cmath>
#include <fstream>
#include <vector>
#include <boost/make_shared.hpp>
#include <boost/math/constants/constants.hpp>
#include <nuflux/IPLEFlux.h>
file LegacyConventionalFlux.cpp
#include <fstream>
#include <boost/lexical_cast.hpp>
#include <boost/make_shared.hpp>

Variables

const double MIN_ENERGY = 10
const double MAX_ENERGY = 1e9
file LegacyPromptFlux.cpp
#include <fstream>
#include <boost/lexical_cast.hpp>
#include <boost/make_shared.hpp>

Variables

const float eMax = 1e10
file SplineFlux.cpp
#include <fstream>
#include <iostream>
#include <boost/make_shared.hpp>
#include <boost/lexical_cast.hpp>
file SplineFlux2.cpp
#include <fstream>
#include <iostream>
#include <boost/make_shared.hpp>
#include <boost/lexical_cast.hpp>
dir /home/runner/work/nuflux/nuflux/src/include
dir /home/runner/work/nuflux/nuflux/src/library
dir /home/runner/work/nuflux/nuflux/src/include/nuflux
dir /home/runner/work/nuflux/nuflux/src