rpdf C++ API Reference¶
-
class FastConvolutedPandel : public rpdf::PhotoElectronProbability¶
- #include <pandel.h>
Contains functions for calculating the convoluted Pandel function and it’s survival function using approximations of these functions which were optimized to be as fast as possible.
Intended to be used with I3RecoLLH parameter PEProb=”GaussConvoluted”
Public Functions
-
inline FastConvolutedPandel(const double jitter, const rpdf::IceModel &ice_model)¶
Creates a FastConvolutedPandel object for calculating probabilities from the convoluted pandel function
- Parameters:
jitter – the width of the Gaussian which is convoluted with the pandel function
ice_model – a struct containing the description of optical properties of the ice
-
virtual double pdf(const double t_res, const double eff_distance) const¶
Calculate the convoluted pandel function using the approximation by van Eindhoven et al, by breaking it into 5 regions and evaluating a different approximation in each region.
- Parameters:
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
- Returns:
the probability density of a photon arriving at t=t_res
-
virtual double sf(const double t_res, const double eff_distance) const¶
Calculate complementary cumulative distribution function of the convoluted pandel function using an approximate integral by D. Chirkin.
The approximation replaces the Gaussian with a box of the same first and second moment. An additional term is added to account for the behavior at the singularity at t=0. This approximation is very accurate. See https://icecube.wisc.edu/~dima/work/WISC/cpdf/a.pdf for more information
- Parameters:
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
- Returns:
the probability of a photon arriving at t>=t_res
-
inline FastConvolutedPandel(const double jitter, const rpdf::IceModel &ice_model)¶
-
struct I3HitCache¶
Structure for caching the information needed to calculate the likelihood. This data is the same regardless of the hypothesis so it is computed once per event and cached so it doesn’t have to be calculated for each iteration
Public Members
-
double total_npe¶
Total number of Photoelectrons observed by the DOM in the event.
-
double first_pulse_time¶
The time of the first pulsse.
-
I3Position pos¶
position of the DOM
-
I3Direction om_ori¶
direction of the PMT axis (for Gen-1, this’ll be downward)
-
double total_npe¶
-
class I3RecoLLH : public I3EventLogLikelihoodBase¶
- #include <I3RecoLLH.h>
A gulliver likelihood service for reconstructions which use the Pandel function.
It calculates the likelihood of a muon track hypothesis using an analytic description of light propagation through the ice. This analytic description is referred to as the Pandel function. This implementation assumes all ice in the detector has uniform optical properties.
Public Functions
-
I3RecoLLH(const std::string &input_readout, const std::string &likelihood, const std::string &peprob, const double jitter, const double noise, const rpdf::IceModel &ice, const bool gen1fast = true, const double bundleradius = SMALLEST_RADIUS)¶
Create a new I3RecoLLH object.
- Parameters:
input_readout – The name of the I3RecoPulseSeriesMap to read from the frame
likelihood – The name of the DOM likelihood algorithm to use. Options are “GaussConvoluted” or “UnconvolutedPandel”
peprob – The name of the photoelectron probability calculation to use Options are “SPE1st” or “MPE”
jitter – the width of the Gaussian which is convoluted with the pandel function
noise – The frequency of noise hits to use in the reconstruction
ice – An IceModel struct containing the description of optical properties of the ice
gen1fast – The backward-compatible “fast” version for which all PMT’s are assumed downward and no bundle radius
bundleradius – The radius of the “truncated” part of the bundle-hypothesis truncated cone (default = a very small number)
-
virtual ~I3RecoLLH()¶
-
void SetGeometry(const I3Geometry &geo)¶
This is called when the reader gets a new geometry. It saves the geometry for SetEvent() to save the location of each DOM in the hit cache
- Parameters:
geo – the Geometry to use to calculate likelihoods
-
void SetEvent(const I3Frame &frame)¶
Called when a new event occurs, this reads the new event and the pulse map in the frame. and calls SetHitCache()
- Parameters:
frame – the frame to extract the I3RecoPulseSeriesMap from
-
void SetPulseMap(const I3RecoPulseSeriesMap &pulse_map)¶
Calculates the hit cache from the geometry and given pulse map
- Parameters:
pulse_map – map containing the pulses for calculating the likelihood
-
double GetLogLikelihood(const I3EventHypothesis &event_hypothesis)¶
Calculate the log likelihood of an event hypothesis for the current event
- Parameters:
event_hypothesis – the gulliver event hypothesis (which contains an I3Particle) with which to calculate the likelihood
- Returns:
the likelihood of the hypothesis
-
unsigned int GetMultiplicity()¶
- Returns:
the multiplicity of the event in question: the number of hit DOMs
-
const std::string GetName() const¶
- Returns:
a string which describes this particular instance of I3RecoLLH
- SET_LOGGER ("I3RecoLLH")
Private Members
-
const std::string peprob_¶
The name of the Photoelectron probability to use: FastConvoluted or Unconvoluted.
-
const double jitter_¶
The time scale of the DOM jitter to consider when reconstructing.
-
const double noise_¶
The frequency of noise hits to use in the reconstruction.
-
const bool gen1fast_¶
Run in Gen-1-style “fast mode” (default = true), with no bundle radius or non-vertical PMT orientations?
-
double bundle_radius_¶
The bundle radius (default = a very small number)
-
rpdf::DOMLikelihoodFunction dom_likelihood_func_¶
Instance of the DOM likelihood function object.
-
std::shared_ptr<rpdf::PhotoElectronProbability> pe_prob_¶
Instance of the Photoelectron probability function object.
-
I3GeometryConstPtr geoptr_¶
A pointer to the geometry to store.
-
std::vector<I3HitCache> hit_cache_¶
a list of the pertinent information for each hit in an event stored in a vector for fast access
-
I3RecoLLH(const std::string &input_readout, const std::string &likelihood, const std::string &peprob, const double jitter, const double noise, const rpdf::IceModel &ice, const bool gen1fast = true, const double bundleradius = SMALLEST_RADIUS)¶
-
class I3RecoLLHFactory : public I3ServiceFactory¶
- #include <I3RecoLLHFactory.h>
A I3Service factory to install an I3RecoLLH into the context
Public Functions
-
virtual ~I3RecoLLHFactory()¶
Private Functions
- SET_LOGGER ("I3RecoLLHFactory")
Private Members
-
double jitter_¶
The time scale of the DOM jitter to consider when reconstructing.
-
double noiseProb_¶
The frequency of noise hits to use in the reconstruction.
-
std::string peprob_¶
The name of the Photoelectron probability to use: FastConvoluted or Unconvoluted.
-
I3EventLogLikelihoodBasePtr llh_¶
pointer to the instance of the object installed in the context
-
bool gen1fast_¶
Run in Gen-1-style fast mode (no bundle radius, or non-vertical PMT’s)
-
double bundle_radius_¶
For truncated-cone (Devil’s Tower) style bundles: the radius of the bundle.
-
virtual ~I3RecoLLHFactory()¶
-
struct IceModel¶
- #include <geometry.h>
This structure holds all the relevant parameters to describe the optical properties of ice used for reconstruction
Public Functions
-
inline IceModel(const double a, const double t, const double s, const double p1, const double p0c0, const double p0c1, const double p0c2)¶
IceModel is mostly a dumb struct, however the constructor computes rho from the other parameters
- Parameters:
a – absorption_length - The distance a photon can go before being absorbed on average
t – tau_scale - The pandel tau parameter, represents the non-linear behavior of multiple photon scatters
s – scattering_length - The distance a photon can go on average before being absorbed
p1 – P1 - The coefficient in front of the distance when calculating effective distance
p0c0 – P0_CS0 - The constant coefficient when calculating effective distance
p0c1 – P0_CS1 - The coefficient in front of \(\cos\eta\) when calculating effective distance
p0c2 – P0_CS2 - The coefficient in front of \(\cos^2\eta\) when calculating effective distance
Public Members
-
double absorption_length¶
The distance a photon can go before being absorbed on average.
-
double tau_scale¶
The pandel tau parameter, represents the non-linear behavior of multiple photon scatters.
-
double scattering_length¶
The distance a photon can go on average before being absorbed.
-
double P1¶
The coefficient in front of the distance when calculating effective distance.
-
double P0_CS0¶
The constant coefficient when calculating effective distance.
-
double P0_CS1¶
The coefficient in front of \(\cos\eta\) when calculating effective distance.
-
double P0_CS2¶
The coefficient in front of \(\cos^2\eta\) when calculating effective distance.
-
double rho¶
rho is a combination of the other optical properties, it is the scale parameter in the Pandel function
-
inline IceModel(const double a, const double t, const double s, const double p1, const double p0c0, const double p0c1, const double p0c2)¶
-
struct MPEfunc¶
- #include <pandel.h>
Public Functions
-
double operator()(const PhotoElectronProbability &p, const double t_res, const double eff_distance, const double Npe) const¶
Boost.Function for calculating the MPE likelihood for an individual DOM.
MPE is the likelihood of observing the first photon at t=t_res given that the DOM saw Npe photoelectrons total. This is intended for use with I3RecoLLH parameter DOMLikelihood=MPE.
- Parameters:
p – the class representing the photoelectron probability calculation
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
Npe – the number of photoelectrons observed by this DOM in this event
- Returns:
the MPE likelihood for this DOM
-
double operator()(const PhotoElectronProbability &p, const double t_res, const double eff_distance, const double Npe) const¶
-
class PhotoElectronProbability¶
- #include <pandel.h>
Base Class for Photoelectron Probabilities. I3RecoLLH uses a pluggable system to pick which PE probability to use. Those object inherent from this class.
Subclassed by rpdf::FastConvolutedPandel, rpdf::UnconvolutedPandel
Public Functions
-
inline virtual ~PhotoElectronProbability()¶
-
virtual double pdf(const double t_res, const double eff_distance) const = 0¶
virtual function for the probability density of a photon arriving at t=t_res at a distance d from the track
- Parameters:
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
- Returns:
the probability density of a photon arriving at t=t_res
-
virtual double sf(const double t_res, const double eff_distance) const = 0¶
virtual function for complementary cumulative distribution function: the probability of observing a hit at t_res or later a distance eff_distance away from the track
- Parameters:
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
- Returns:
the probability of a photon arriving at t>=t_res
-
inline virtual ~PhotoElectronProbability()¶
-
struct SPEfunc¶
- #include <pandel.h>
Public Functions
-
double operator()(const PhotoElectronProbability &p, const double t_res, const double eff_distance, const double Npe) const¶
Boost.Function for calculating the SPE1st likelihood for an individual DOM.
It takes the time of the first from the hit series and simply returns the PDF for that hit. Intended for use with I3RecoLLH parameter DOMLikelihood=SPE1st.
- Parameters:
p – the class representing the photoelectron probability calculation
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
Npe – the number of photoelectrons observed by this DOM in this event
- Returns:
the SPE1st Likelihood for this DOM
-
double operator()(const PhotoElectronProbability &p, const double t_res, const double eff_distance, const double Npe) const¶
-
class UnconvolutedPandel : public rpdf::PhotoElectronProbability¶
- #include <pandel.h>
contains functions for calculating the pandel function without any convolution. Intended for use with I3RecoLLH with parameter PEProb=”UnconvolutedPandel”
Public Functions
-
inline UnconvolutedPandel(const rpdf::IceModel &ice_model)¶
Construct an UnconvolutedPandel object
- Parameters:
ice_model – ice struct containing the description of the ice
-
virtual double pdf(const double t_res, const double eff_distance) const¶
Calculate the unconvouted pandel function using boost’s implementation of the gamma distribution.
- Parameters:
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
- Returns:
the probability density of a photon arriving at t=t_res
-
virtual double sf(const double t_res, const double eff_distance) const¶
calculate complementary cumulative distribution function of the unconvouted pandel function using gsl’s implementation of the gamma distribution
- Parameters:
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
- Returns:
the probability of a photon arriving at t>=t_res
-
inline UnconvolutedPandel(const rpdf::IceModel &ice_model)¶
-
namespace [anonymous]
-
namespace rpdf¶
Typedefs
-
typedef boost::function<double(const PhotoElectronProbability &p, const double t_res, const double d_eff, const double Npe)> DOMLikelihoodFunction¶
function declaration for the DOM Likelihoods
Functions
-
double effective_distance(const double distance, const double cos_eta, const IceModel &ice_model)¶
Calculate the effective distance between an emitter and a DOM due to scattering. This is a small correction for DOMs orientated away from the source.
- Parameters:
distance – the actual distance traveled by the photon
cos_eta – the cosine of the angle between the photon’s trajectory and the direction the DOM is pointing
ice_model – the object holding the description of the ice
- Returns:
the equivalent distance that would have experienced the same amount of scattering if the DOM was orientated toward the photon path
-
double cherenkov_time(const double d_track, const double d_approach)¶
Calculate the time Cherenkov photon arrives at a position d_track along the track, to a DOM located d_approach from the track.
- Parameters:
d_track – The distance along the track from the vertex to the point of closest approach
d_approach – The distance of closest approach between the DOM and track
- Returns:
the time at which the Cherenkov photon arrives at the DOM
-
double planewave_time(const double d_track)¶
Calculate the time a planewave wavefront arrives at a position d_track along the track. Useful for “Devil’s Tower” style reconstructions.
- Parameters:
d_track – The distance along the track from the vertex to the point of closest approach
- Returns:
the time at which the planewave arrives at the DOM
-
std::pair<double, double> muon_geometry(const I3Position &om, const I3Particle &track, const IceModel &ice_model)¶
This calculates the two geometrical parameters needed by the pandel function to calculate the PDF.
Because the calculations relating to these function are very tightly coupled they are both calculated in a single function
- Parameters:
om – the position of the optical module being hit
track – the muon track hypothesis
ice_model – the model used to describe the optical properties of the ice
- Returns:
t_hit the difference in time between the hit and the expected time of a direct photon
- Returns:
effective_distance the effective distance between the track hypothesis and the OM including scattering needed to hit an OM facing the other way
-
std::pair<double, double> muon_bundle_geometry(const I3Position &om, const I3Particle &track, const double &bundle_radius, const IceModel &ice_model, const I3Direction &omori = I3Direction(0, 0, -1))¶
Similar to the above, but for muon bundles, under the “Devil’s Tower” model (11/2023, K. Rawlins)
The “Devil’s Tower” model, for muon bundles, in which the Cherenkov cone is “truncated”: modeled as a planewave inside some radius (the “bundle_radius”), and a Cherenkov cone outside. New 11/2023 (K. Rawlins):
In theory, such a model will reduce to the same thing as the “muon_geometry” function above if the bundle_radius is zero. So in theory, only one function is needed; there is some duplicate code here. But since this version is probably a little slower — and new/experimental — I’ll keep it seperate for now, and leave the original one undisturbed. If we like it, maybe later we can merge them into one function.
This version of the function also takes the PMT orientation as an argument, for computing d_eff for any PMT, including those who might not be oriented vertically downward (for instance, mDOM’s, etc.)
-
double pandel_pdf(const double t_res, const double eff_distance, const IceModel &ice_model)¶
Pure function implementation of the Pandel PDF: the probability of observing a hit at t_res a distance eff_distance away from the track
- Parameters:
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
ice_model – a struct containing the description of optical properties of the ice
- Returns:
the probability density of a photon arriving at t=t_res
-
double pandel_sample(const double eff_distance, const IceModel &ice_model, I3RandomService &rng)¶
Sample a photon arrival time from the Pandel PDF.
- Parameters:
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
ice_model – a struct containing the description of optical properties of the ice
rng – an I3RandomService to use to sample random numbers
- Returns:
a time residual relative to the direct Cherenkov time at which the photon hits a DOM (see t_res above).
-
double pandel_sf(const double t_res, const double eff_distance, const IceModel &ice)¶
Pure function implementation of the Pandel complementary cumulative distribution function: the probability of observing a hit at t_res or later a distance eff_distance away from the track
- Parameters:
t_res – the time residual of the of a photon which hits a DOM
eff_distance – the distance of the DOM from the hypothesis track accounting for extra scatters caused by the PMT looking away from the track
ice – An IceModel struct containing the description of optical properties of the ice
- Returns:
the probability of a photon arriving at t>=t_res
-
double gslConvoluted1F1Diff(const double xi, const double eta)¶
-
double gslConvolutedU(const double xi, const double eta)¶
-
double fastConvolutedHyperg(const double xi, const double eta)¶
Compute the hypergeometric portion of the Gaussian convoluted Pandel PDF.
1F1(0.5*xi,0.5,0.5*eta*eta)/gamma(0.5*(xi+1)) - sqrt2*eta*1F1(0.5*(xi+1),1.5,0.5*eta*eta)/gamma(0.5*xi) for 0.5*eta*eta < 100 and xi < 5 using the hypergeometric power series with the denominator precomputed. This is >10x faster than GSL.
Variables
- IceModel H0 (98., 596.0, 36.93, 0.9045, 4.249,-6.629, 5.430)
AMANDA ice model 0.
- IceModel H1 (98., 578.0, 35.04, 0.8682, 3.615,-5.081, 5.015)
AMANDA ice model 1.
- IceModel H2 (98., 556.7, 33.29, 0.8395, 3.094,-3.946, 4.636)
AMANDA ice model 2, this is what is reported in the AMANDA muon reconstruction paper. It is the default ice model in I3RecoLLH.
- IceModel H3 (98., 542.2, 31.72, 0.8106, 2.683,-3.090, 4.436)
AMANDA ice model 3.
- IceModel H4 (98., 526.3, 29.94, 0.7790, 2.139,-0.9614, 4.020)
AMANDA ice model 4.
-
typedef boost::function<double(const PhotoElectronProbability &p, const double t_res, const double d_eff, const double Npe)> DOMLikelihoodFunction¶
-
namespace constants¶
Variables
-
const double C_VACUUM = 2.997924e-1¶
Speed of light in vacuum, less precision than what is found in I3Constants because this is what ipdf did.
-
const double N_ICE_G = I3Constants::n_ice_group¶
Group index of refraction of ice.
-
const double SIN_CHERENKOV = sin(I3Constants::theta_cherenkov)¶
the sine of the Cherenkov angle
-
const double TAN_CHERENKOV = tan(I3Constants::theta_cherenkov)¶
the tangent of the Cherenkov angle
-
const double EFF_TAN_CHERENKOV = N_ICE_G / SIN_CHERENKOV - 1 / TAN_CHERENKOV¶
a grouping of constants that is always the same in the calculation. It is called the effective tangent of the Cherenkov angle because It is that if the group and phase angles are the same
-
const double C_VACUUM = 2.997924e-1¶
-
namespace std
STL namespace.
- file geometry.cxx
- #include “dataclasses/physics/I3Particle.h”#include “rpdf/geometry.h”
This file contains functions pertaining to geometrical calculations needed for reconstructing muons.
- Copyright
(C) 2018 The Icecube Collaboration
- Author
Kevin Meagher
- Date
January 2018
- file geometry.h
- #include <utility>#include “dataclasses/I3Constants.h”
functions for calculating geometrical parameters for muon reconstructions
- Copyright
(C) 2018 The Icecube Collaboration
- Author
Kevin Meagher
- Date
January 2018
- file I3RecoLLH.cxx
- #include “gulliver/I3EventHypothesis.h”#include “rpdf/I3RecoLLH.h”
Provides a Gulliver likelihood service for muon based reconstructions.
- Copyright
(C) 2018 The Icecube Collaboration
- Author
Kevin Meagher
- Date
January 2018
- file I3RecoLLH.h
- #include “dataclasses/geometry/I3Geometry.h”#include “gulliver/I3EventLogLikelihoodBase.h”#include “rpdf/pandel.h”#include “rpdf/geometry.h”
Provides a Gulliver likelihood service for muon based reconstructions.
(c) 2018 the IceCube Collaboration
- Author
Kevin Meagher
- Date
January 2018
Defines
-
SMALLEST_RADIUS¶
- file I3RecoLLHFactory.cxx
- #include “rpdf/I3RecoLLH.h”#include “rpdf/I3RecoLLHFactory.h”
- file I3RecoLLHFactory.h
- #include “icetray/I3ServiceFactory.h”#include “gulliver/I3EventLogLikelihoodBase.h”#include “rpdf/geometry.h”
Provides a Gulliver likelihood service for muon based reconstructions.
(c) 2018 the IceCube Collaboration
- Author
Kevin Meagher
- Date
January 2018
- file pandel.cxx
- #include <gsl/gsl_math.h>#include <gsl/gsl_rng.h>#include <gsl/gsl_randist.h>#include <gsl/gsl_sf_gamma.h>#include <gsl/gsl_sf_hyperg.h>#include <gsl/gsl_sf_erf.h>#include <gsl/gsl_errno.h>#include <boost/math/distributions/gamma.hpp>#include “dataclasses/physics/I3Particle.h”#include “rpdf/geometry.h”#include “rpdf/pandel.h”
functions for calculating photoelectron probabilities for muon reconstructions
- Copyright
(C) 2018 The Icecube Collaboration
- Author
Kevin Meagher
- Date
January 2018
- file pandel.h
- #include <boost/function.hpp>#include “dataclasses/physics/I3RecoPulse.h”#include “phys-services/I3RandomService.h”#include “rpdf/geometry.h”
functions for calculating photoelectron probabilities for muon reconstructions
- Copyright
(C) 2018 The Icecube Collaboration
- Author
Kevin Meagher
- Date
January 2018
- dir icetray
- dir private
- dir public
- dir rpdf
- dir rpdf
- dir rpdf