ddddr C++ API Reference¶
-
struct BundleBin¶
- #include <I3TrueMuonEnergy.h>
Struct that contains the energy of a muon bundle at a certain position along its track.
-
struct Checkpoint
- #include <MuonGunTrack.h>
A point at which the absolute energy of the particle is known.
Public Functions
-
inline Checkpoint(double l, double e = 0., size_t o = 0)
Public Members
-
double length
-
double energy
-
ptrdiff_t offset
-
inline Checkpoint(double l, double e = 0., size_t o = 0)
-
struct DOMDATA¶
- #include <I3MuonEnergyProfile.h>
Struct that holds all data of a certain DOM needed for the energy loss distribution.
-
class ExpoFcn : public NoGradientFcn<ExpoFcn>¶
- #include <ExpoFcn.h>
Minimizer function based on an exponential function (ExpoFunction).
Public Functions
-
inline ExpoFcn(const std::vector<double> &meas, const std::vector<double> &pos)¶
Constructs object for given set of measurements.
- Parameters:
meas – Measured values (y-values)
pos – Position of the values (x-values)
-
inline ~ExpoFcn()¶
-
double operator()(const std::vector<double>&)¶
Returns the log likelihood for given set of parameters.
-
double f(const std::vector<double>&, double) const¶
Returns function value for given set of parameters.
Not part of the Minuit2 standard, but here used as a convenience method so that the chi squared and the reduced log likelihood can be calculated easier in the I3MuonEnergy and I3TrueMuonEnergy modules.
- Parameters:
par – Vector of parameters
X – Slant depth
- Returns:
ExpoFunction(X) for given parameters.
-
inline void setErrorDef(double def)¶
-
inline ExpoFcn(const std::vector<double> &meas, const std::vector<double> &pos)¶
-
class ExpoFunction¶
- #include <ExpoFunction.h>
Class that allows to call an exponential function for a given set of parameters.
Used by ExpoFcn.h.
-
class I3MuonEnergy : public I3ConditionalModule¶
- #include <I3MuonEnergy.h>
This module estimates and fits the energy loss along a given track.
The module is the implementation of the Data-Derived Differential Deposition Reconstruction (DDDDR) algorithm, described in detail in the IceCube wiki at https://wiki.icecube.wisc.edu/index.php/IC79_Atmospheric_Muons/DDDDR. It takes any kind of reconstructed track as input and estimates the energy loss around the track up to a specified impact parameter, using a table of depth-dependent values for the light attenuation coefficient lambda. The estimated (binned) energy loss for the individual DOMs along the track can be saved optionally, together with the result of the fit to the energy loss.
The result of the fit is saved together with some properties of the energy loss distribution in an object of the type I3MuonEnergyParams. More about the track parameters at https://wiki.icecube.wisc.edu/index.php/Muon_Bundle_Suppression#Stochasticity.
The parameters I3MuonEnergyParams::outputDOMs_ and I3MuonEnergyParams::saveDomResults_ can be used to save the full energy loss distribution used for the fit. OutputDOMsName+Slant/Impact/Depth are vectors with the slant depth, distance to the track or depth of each individual DOM while OutputDOMsName+(e)dEdX contain the (errors on the) energy losses. The OutputDOMsNamebinned vectors contain the same information for the binned energy loss distribution. In this case, the depth corresponds to the bincenters of the slant depth bins.
Public Functions
-
~I3MuonEnergy()¶
-
void Configure()¶
-
void Geometry(I3FramePtr)¶
-
void Calibration(I3FramePtr)¶
-
void DetectorStatus(I3FramePtr)¶
-
void Physics(I3FramePtr)¶
-
void Finish()¶
- SET_LOGGER ("I3MuonEnergy")
Private Functions
-
I3ParticlePtr getSeedTrack(I3FramePtr&)¶
Determines the track for the energy loss estimation, either reconstructed track specified in I3MuonEnergy::trackName_ or true Monte Carlo track.
-
boost::shared_ptr<I3MuonEnergyProfile> getMuonEnergyProfile(I3FramePtr&, double, double, std::vector<DOMDATA>&)¶
Determines bins for the energy distribution through the DOMs along the track or the length of the Monte Carlo track and fills it with the estimated energy losses.
The I3TrueMuonEnergy module also determines the energy loss in bins along the track of the primary particle. This method makes sure to use the same binning when using the Monte Carlo track as input track, e.g. when doing an comparison of true and reconstructed energy loss.
-
std::vector<I3FitParameterInitSpecs> getFitParameterSpecs(FitFunction&)¶
Creates parameter for the chosen fit function.
The fit function is either an exponential function or Tom Feusels’ function for muon energy loss (s. TomFFunction.h). In the case of the latter, an additional parameter is needed, and parameters are fixed of specified so in the module parameters.
- Parameters:
fitFunction –
- Returns:
Vector of FitParameterSpecs
Saves results from fit into an object of type I3MuonEnergyParams.
In addition to the fit parameters, the chi^2 and the reduced log likelihood are determined, as well as peak energy and sigma, median and mean of the energy loss distribution.
- Parameters:
eprofile – The histogram with the energy losses
- Returns:
I3MuonEnergyParamsPtr containing the results and the bundle suppression parameters (s. https://wiki.icecube.wisc.edu/index.php/Muon_Bundle_Suppression#Stochasticity).
-
I3ParticlePtr getWeightedPrimary(I3FramePtr&)¶
Saves the binned energy profile and the individual DOMs considered for the profile in the frame.
Saves the results from the cascade energy reconstruction into an object of type I3MuonEnergyCascadeParams.
Private Members
-
std::string mcTreeName_¶
Only needed if I3MuonEnergy::useMonteCarloTrack_ is true.
-
std::string mmcName_¶
Name of the MMCTracklist frame object. Only needed if I3MuonEnergy::useMonteCarloTrack_ is true.
-
std::string outputPars_¶
Name for the I3MuonEnergyParams frame object.
-
std::string outputCascadePars_¶
Name for the I3MuonEnergyCascadeParams frame object.
-
double levelDist_¶
Coefficient resp. for flattening of the exponential function used to estimate the energy loss for a DOM.
-
FitFunction fitFunction_¶
Either exponential function or Tom Feusels’s function to fit the energy loss.
-
int fitFunctionInt_¶
-
bool fixB_¶
fix parameter b in fit with Tom Feusels’s function
-
bool fixGamma_¶
fix parameter gamma in fit with Tom Feusels’s function
-
bool saveDomResults_¶
Save energy loss, impact factor, depth and slant depth both per DOM and binned as vectors to the frame.
-
bool useMonteCarloTrack_¶
If set to true, I3MuonEnergy::trackName_ is ignored.
-
int slantBinNo_¶
Number of bins to divide the track into. Only needed if I3MuonEnergy::useMonteCarloTrack_ is true.
-
double purityAdjust_¶
Decrease in the atten. length at the top.
-
double maxImpact_¶
Max. distance between DOM and track.
-
double fitBinWidth_¶
Bin width of slant depth along the track. Only needed if module is used with reconstructed track as input.
-
int nDoF_¶
Number of degrees of freedom for the fit.
-
double tolerance_¶
-
unsigned int maxIterations_¶
-
I3GeometryConstPtr geometry_¶
-
I3OMGeoMapPtr omGeo_¶
-
I3DOMCalibrationMapPtr domCal_¶
-
I3VectorOMKeyConstPtr badDomList_¶
-
int tracksFitted_¶
-
int pulsemap_warned_¶
-
double surface_height_¶
-
~I3MuonEnergy()¶
-
class I3MuonEnergyCascadeParams : public I3FrameObject¶
- #include <I3MuonEnergyCascadeParams.h>
Class to save the results of a fit of the energy loss distribution to an I3Frame.
Public Members
-
double nDOMsCascade¶
Number of DOMs used to determine the cascade energy.
-
double cascade_energy¶
Reconstructed cascade energy.
-
double cascade_energy_sigma¶
Statistical error on the reconstructed cascade energy, calculated from the standard deviation of the distribution of the reconstructed energies for all DOMs used for the cascade energy reconstruction.
-
I3Position cascade_position¶
Position of the leading casascade on the track.
-
double cascade_slant_depth¶
Slant depth of the leading cascade.
Friends
- friend class icecube::serialization::access
-
double nDOMsCascade¶
-
class I3MuonEnergyParams : public I3FrameObject¶
- #include <I3MuonEnergyParams.h>
Class to save the results of a fit of the energy loss distribution to an I3Frame.
Public Members
-
double N¶
Normalization constant of Exp. or TomF function.
-
double N_err¶
-
double b¶
Slope of Exp. or TomF function.
-
double b_err¶
-
double gamma¶
Exponent for TomF function.
-
double gamma_err¶
-
double nDOMs¶
Number of DOMs used to determine energy distribution.
-
double rllh¶
Reduced log likelihood of the fit.
-
double chi2¶
Chi square of the fit.
-
double chi2ndof¶
Chi square per degree of freedom of the fit.
-
double peak_energy¶
Peak energy of the energy loss distribution.
-
double peak_sigma¶
Uncertainty of the peak energy.
-
double mean¶
Mean of the energy loss distribution.
-
double median¶
Median of the energy loss distribution.
-
double bin_width¶
Bin width of the energy loss distribution.
-
I3Particle::FitStatus status¶
Status of the fit.
Friends
- friend class icecube::serialization::access
-
friend std::ostream &operator<<(std::ostream&, const I3MuonEnergyParams&)¶
-
double N¶
-
class I3MuonEnergyProfile¶
- #include <I3MuonEnergyProfile.h>
A class used to create a binned energy loss distribution from the individual energy losses of all DOMs.
I3MuonEnergyProfile is a reimplementation of some basic features of ROOT’s TProfile class. Each bin holds the average value and statistical error of all data points in this bin. The binning is determined by the bin width and the min. and max. values for the slant. Instead of the limits of slant depth, data points of slant depth and energy loss can be given to the constructor to determine the limits automatically.
Public Functions
-
I3MuonEnergyProfile(std::vector<double>, std::vector<double>, double, const std::string errors = "default")¶
-
~I3MuonEnergyProfile()¶
-
int Fill(double, double)¶
Fill a pair of slant depth and energy loss in the histogram.
-
int FindBin(double)¶
Find the bin that contains the input energy loss.
-
double GetBinCenter(int)¶
Get position of given bin center.
-
double GetBinContent(int)¶
Get content of given bin.
-
double GetBinError(int)¶
Get error of given bin.
-
double GetBinWidth()¶
Get bin width.
-
int GetMaxBin()¶
Get bin with max. energy loss.
-
double GetMeanEnergyLoss()¶
Get the mean energy loss.
The mean is calculated using the sum and number of all entries, not as average of the bins.
- Returns:
mean energy loss
-
double GetMedianEnergyLoss()¶
Get median energy loss.
The median of the energy content of all non-empty bins.
- Returns:
median energy loss
-
int GetNBins()¶
-
int GetNNonZeroBins()¶
-
int GetNEntries()¶
-
int GetNEnergyLosses()¶
-
double GetXMin()¶
-
double GetXMax()¶
Higher edge of the last bin.
Private Members
-
std::string errors_¶
Type of error to be returned when calling I3MuonEnergyProfile::GetBinError.
-
double nBins_¶
Number of bins.
-
int nBinEntries_¶
Total number of bin entries.
-
double binWidth_¶
Bin width.
-
double sumWeightsy_¶
Total sum of all bin entries.
-
double xmax_¶
Upper limit for slant depth.
-
double xmin_¶
Lower limit for slant depth.
-
double nEnergyLosses_¶
Number of DOMs that had energy loss.
-
I3MuonEnergyProfile(std::vector<double>, std::vector<double>, double, const std::string errors = "default")¶
-
class I3TrueMuonEnergy : public I3ConditionalModule¶
- #include <I3TrueMuonEnergy.h>
This module calculates the total energy of a muon bundle by summing up the contribution of each muon for each bin along the track of the whole bundle.
It also performs a fit of an exponential function to either the distribution of the total energy or the energy loss.
Public Functions
-
~I3TrueMuonEnergy()¶
-
void Configure()¶
-
void Physics(I3FramePtr)¶
- SET_LOGGER ("I3TrueMuonEnergy")
Private Functions
Private Members
-
double slantBinWidth_¶
The bin width along the track.
-
int slantBinNo_¶
Number of bins to use.
-
bool energyLoss_¶
If true, the energy loss is used for the fit and saved to the frame.
-
bool saveEnergyLosses_¶
If true, save the energy losses as vectors in the frame.
-
double maxRadius_¶
Max radius of the cylinder around the detector. Muons outside of it will not be considered.
-
double tolerance_¶
-
unsigned int maxIterations_¶
-
~I3TrueMuonEnergy()¶
-
struct LossSum
- #include <MuonGunTrack.h>
The sum of stochastic energy losses since the last checkpoint.
Public Functions
-
inline LossSum(double l, double e = 0.)
Public Members
-
double length
-
double energy
-
inline LossSum(double l, double e = 0.)
-
template<typename Base>
class NoGradientFcn : public I3GulliverBase¶ - #include <NoGradientFcn.h>
-
class TomFFcn : public NoGradientFcn<TomFFcn>¶
- #include <TomFFcn.h>
Minimizer function based on an Tom Feusels’ function for muon bundle energy loss (TomFFunction.h).
Public Functions
-
inline TomFFcn(const std::vector<double> &meas, const std::vector<double> &pos)¶
Constructs object for given set of measurements.
- Parameters:
meas – Measured values (y-values)
pos – Position of the values (x-values)
-
inline ~TomFFcn()¶
-
double operator()(const std::vector<double>&)¶
Returns the log likelihood for given set of parameters.
-
double f(const std::vector<double>&, double) const¶
Returns function value for given set of parameters.
Not part of the Minuit2 standard, but here used as a convenience method so that the chi squared and the reduced log likelihood can be calculated easier in the I3MuonEnergy and I3TrueMuonEnergy modules.
- Parameters:
par – Vector of parameters
X – Slant depth
- Returns:
TomFFunction(X) for given parameters.
-
inline void setErrorDef(double def)¶
-
inline TomFFcn(const std::vector<double> &meas, const std::vector<double> &pos)¶
-
class TomFFunction¶
- #include <TomFFunction.h>
Class that allows to call an Tom Feusels’ function for a given set of parameters.
Used by TomFFcn.h.
Public Functions
-
inline TomFFunction(double N, double b, double gamma)¶
-
inline ~TomFFunction()¶
-
inline double N() const¶
-
inline double b() const¶
-
inline double gamma() const¶
-
inline double operator()(double x) const¶
Computes the TomF function.
See https://arxiv.org/abs/0912.4668 for reference. Further information at https://srl.utu.fi/AuxDOC/kocharov/ICRC2009/pdf/icrc0518.pdf,
and in “Measurement of cosmic ray composition and energy spectrum between
1 PeV and 1 EeV with IceTop and IceCube”, PhD thesis Tom Feusels 2013, UGent.
https://hdl.handle.net/1854/LU-4337238
-
inline TomFFunction(double N, double b, double gamma)¶
-
class Track : public I3Particle
- #include <MuonGunTrack.h>
A particle that moves in a straight line, and loses energy as it does so.
Public Functions
-
inline Track()
-
Track(const I3MMCTrack &mmctrack, const I3MCTree::sibling_const_iterator &secondaries_begin, const I3MCTree::sibling_const_iterator &secondaries_end)
-
double GetEnergy(double length) const
Get the energy of the particle at the given down-track distance, assuming that the “continuous” contribution to the energy loss is constant between measurement points.
- Parameters:
length – [in] distance from the track origin
- Returns:
an energy if 0 <= length < range; otherwise 0
-
I3Position GetPos(double length) const
-
double GetTime(double length) const
-
inline std::vector<Checkpoint> GetCheckpoints() const
Public Static Functions
-
static std::list<Track> Harvest(const I3MCTree&, const I3MMCTrackList&)
Extract energy losses from frame objects.
Find the stochastic energy losses associated with each I3MMCTrack in the I3MCTree, and store them together in a Track.
-
inline Track()
-
namespace I3MuonGun
Functions
-
I3_POINTER_TYPEDEFS(Track)
-
I3_POINTER_TYPEDEFS(Track)
-
namespace [anonymous]¶
-
namespace MuonEnergyFunctions¶
Contains set of static helper functions used in I3MuonEnergy and I3TrueMuonEnergy.
Functions
-
double eLoss(double, double, double, double, double domEff = 1., double levelDist = 25., double fscale = 0.02)¶
Energy loss measured by a given DOM.
- Parameters:
charge – Charge measured by DOM
impact – Impact parameter of DOM relative to track
lambda – Attenuation parameter lambda at depth of the DOM
domEff – DOM efficiency
levelDist – Flattening parameter of the energy loss function
fscale – Scaling factor
- Returns:
Energy loss, s. https://wiki.icecube.wisc.edu/index.php/IC79_Atmospheric_Muons/DDDDR#Fit_Results
-
double eLossCascade(double, double, double, double, double domEff = 1., double levelDist = 25., double fscale = 0.02)¶
Energy loss around a certain DOM.
-
double getLambda(double, std::stringstream&)¶
Read attenuation parameter lambda for a given depth from a text stream.
- Parameters:
depth – Depth in the ice
fin – String stream of the input file with the lambda table
- Returns:
Ice attenuation parameter lambda
-
double getLambda(double, std::vector<std::vector<double>>)¶
Read attenuation parameter lambda for a given depth from a table.
- Parameters:
depth – Depth in the ice
lambdaFile – Vector of table rows (vectors of depth and lamdba) representing lambda vs. depth
- Returns:
Ice attenuation parameter lambda
-
double getSlantDepth(I3Particle&, I3Position&, double top = 1950.)¶
Get slant depth of a particle at a given position.
Slant depth of the muon (bundle) as seen from the position of a DOM at a given position.
- Parameters:
track – Muon (bundle) track
po – Position of the DOM
top – Z position of the ice surface
- Returns:
Slant depth along the track at the closest approach position of the muon (bundle) relative to the position of the DOM
-
double getVerticalDepthZ(I3Particle&, double, double top = 1950.)¶
Get vertical depth corresponding to given slant depth along a particle track.
- Parameters:
track – Muon (bundle) track
slant – Sland depth along the track
top – Z position of the ice surface
- Returns:
Vertical depth in z
-
double ln_poisson(double, double, double min_exp = 1e-30)¶
-
I3ParticlePtr getWeightedPrimary(I3FramePtr&, std::string)¶
Method that finds the primary particle of a simulated event. Code is taken from the weighting module https://github.com/icecube/icetray/blob/main/weighting/python/__init__.py.
- Parameters:
frame – The frame containing the event
mcTreeName – The I3MCTree containing the primary particle
- Returns:
the weighted primary of the event
-
double eLoss(double, double, double, double, double domEff = 1., double levelDist = 25., double fscale = 0.02)¶
-
namespace std
STL namespace.
- file ExpoFcn.cxx
- #include “ddddr/ExpoFcn.h”#include “ddddr/ExpoFunction.h”#include “gulliver/I3GulliverBase.h”#include <cassert>
- file ExpoFcn.h
- #include “gulliver/I3GulliverBase.h”#include “NoGradientFcn.h”#include “MuonEnergyFunctions.h”#include <vector>
Declaration of ExpoFcn
- Author
Hans-Peter Bretz
- file ExpoFunction.h
- #include <math.h>
Implementation of ExpoFunction
- Author
Hans-Peter Bretz
- file I3MuonEnergy.cxx
- #include “limits.h”#include “ddddr/ExpoFcn.h”#include “ddddr/TomFFcn.h”#include “ddddr/I3MuonEnergy.h”#include “ddddr/MuonEnergyFunctions.h”#include “ddddr/ExpoFunction.h”#include “ddddr/TomFFunction.h”#include “dataclasses/physics/I3EventHeader.h”#include “dataclasses/geometry/I3Geometry.h”#include “dataclasses/geometry/I3TankGeo.h”#include “dataclasses/geometry/I3OMGeo.h”#include “dataclasses/status/I3DetectorStatus.h”#include “dataclasses/physics/I3RecoPulse.h”#include “phys-services/I3Calculator.h”#include “icetray/OMKey.h”#include “dataclasses/I3Map.h”#include “dataclasses/I3Vector.h”#include “dataclasses/physics/I3MCTree.h”#include “dataclasses/physics/I3MCTreeUtils.h”#include “simclasses/I3MMCTrack.h”
< Implementation of I3MuonEnergy (DDDDR)
- Date
$Date: 14/01/23
- Author
Hans-Peter Bretz (hbretz@icecube.wisc.edu)
- Author
Patrick Berghaus
Functions
-
bool EarlierThan(const DOMDATA &r1, const DOMDATA &r2)¶
Helper function that compares two DOMDATA structs by their position along the muon bundle track.
- Parameters:
r1 –
r2 –
- Returns:
True if first DOMDATA is earlier than second.
-
bool smallerThan(const DOMDATA &r1, const DOMDATA &r2)¶
Helper function that compares the energy loss measured by two different DOMs, given in DOMDATA structs.
- Parameters:
r1 –
r2 –
- Returns:
True if first DOMDATA contains energy loss smaller than the second one.
-
double geometricDistance(const DOMDATA r1, I3Position p1)¶
-
I3_MODULE(I3MuonEnergy)¶
Variables
-
static const double SURFACE_HEIGHT = 1948.07¶
-
static const bool FIXB = true¶
-
static const bool FIXGAMMA = true¶
-
static const double EXPO_NORM_INIT = 10¶
-
static const double EXPO_NORM_STEPS = 0.1¶
-
static const double EXPO_NORM_MINV = 0¶
-
static const double EXPO_NORM_MAXV = 15¶
-
static const double EXPO_EXP_INIT = 1e-3¶
-
static const double EXPO_EXP_STEPS = 1e-1¶
-
static const double EXPO_EXP_MINV = -1e-1¶
-
static const double EXPO_EXP_MAXV = 1e-1¶
-
static const double TOMF_NORM_INIT = 1e3¶
-
static const double TOMF_NORM_STEPS = 100¶
-
static const double TOMF_NORM_MINV = 1e-10¶
-
static const double TOMF_NORM_MAXV = 1e10¶
-
static const double TOMF_B_INIT = 3.3e-4¶
-
static const double TOMF_B_STEPS = 1e-3¶
-
static const double TOMF_B_MINV = 1e-4¶
-
static const double TOMF_B_MAXV = 1e-1¶
-
static const double TOMF_GAMMA_INIT = 1.757¶
-
static const double TOMF_GAMMA_STEPS = 1e-1¶
-
static const double TOMF_GAMMA_MINV = 1¶
-
static const double TOMF_GAMMA_MAXV = 4¶
-
static const int CASCADE_CONTAINMENT_DISTANCE = 150¶
- file I3MuonEnergy.h
- #include “icetray/I3ConditionalModule.h”#include “dataclasses/physics/I3Particle.h”#include “dataclasses/geometry/I3Geometry.h”#include “dataclasses/calibration/I3Calibration.h”#include “dataclasses/calibration/I3DOMCalibration.h”#include <icetray/serialization.h>#include <iostream>#include <fstream>#include <cstdlib>#include <icetray/I3Logging.h>#include <string>#include <map>#include <vector>#include <algorithm>#include <boost/function.hpp>#include “ddddr/I3MuonEnergyProfile.h”#include “ddddr/I3MuonEnergyParams.h”#include “ddddr/I3MuonEnergyCascadeParams.h”#include “gulliver/I3FitParameterInitSpecs.h”#include “gulliver/I3GulliverBase.h”#include “lilliput/minimizer/I3GSLSimplex.h”#include “ddddr/TomFFcn.h”#include “ddddr/ExpoFcn.h”
Implementation of I3MuonEnergy (DDDDR)
- Date
$Date: 14/01/23
- Author
Hans-Peter Bretz (hbretz@icecube.wisc.edu)
- Author
Patrick Berghaus
- file I3MuonEnergyCascadeParams.cxx
- #include <icetray/serialization.h>#include “ddddr/I3MuonEnergyCascadeParams.h”
Implementation of the I3MuonEnergyCascadeParams class.
Functions
-
I3_SERIALIZABLE(I3MuonEnergyCascadeParams)¶
-
I3_SERIALIZABLE(I3MuonEnergyCascadeParams)¶
- file I3MuonEnergyCascadeParams.h
- #include <icetray/I3FrameObject.h>#include “dataclasses/Utility.h”#include “dataclasses/physics/I3Particle.h”
Definition of the I3MuonEnergyCascadeParams class.
Functions
-
I3_POINTER_TYPEDEFS(I3MuonEnergyCascadeParams)¶
-
I3_POINTER_TYPEDEFS(I3MuonEnergyCascadeParams)¶
- file I3MuonEnergyParams.cxx
- #include <icetray/serialization.h>#include “ddddr/I3MuonEnergyParams.h”
Implementation of the I3MuonEnergyParams class.
Functions
-
std::ostream &operator<<(std::ostream &s, const I3MuonEnergyParams ¶ms)¶
-
I3_SERIALIZABLE(I3MuonEnergyParams)¶
-
std::ostream &operator<<(std::ostream &s, const I3MuonEnergyParams ¶ms)¶
- file I3MuonEnergyParams.h
- #include <icetray/I3FrameObject.h>#include “dataclasses/Utility.h”#include “dataclasses/physics/I3Particle.h”
Definition of the I3MuonEnergyParams class.
Functions
-
I3_POINTER_TYPEDEFS(I3MuonEnergyParams)¶
-
I3_POINTER_TYPEDEFS(I3MuonEnergyParams)¶
- file I3MuonEnergyProfile.cxx
- #include “ddddr/I3MuonEnergyProfile.h”
- file I3MuonEnergyProfile.h
- #include “dataclasses/I3Vector.h”#include <icetray/I3Logging.h>
Definition of I3MuonEnergyProfile (DDDDR)
- Date
$Date: 14/01/23
- Author
Hans-Peter Bretz (hbretz@icecube.wisc.edu)
- file I3TrueMuonEnergy.cxx
- #include “ddddr/I3TrueMuonEnergy.h”#include “ddddr/MuonEnergyFunctions.h”#include “dataclasses/physics/I3EventHeader.h”#include “dataclasses/geometry/I3Geometry.h”#include “dataclasses/geometry/I3TankGeo.h”#include “dataclasses/geometry/I3OMGeo.h”#include “dataclasses/status/I3DetectorStatus.h”#include “dataclasses/physics/I3RecoPulse.h”#include “phys-services/I3Calculator.h”#include “icetray/OMKey.h”#include “dataclasses/I3Map.h”#include “MuonGunTrack.h”#include <boost/foreach.hpp>
Functions
-
I3_MODULE(I3TrueMuonEnergy)¶
Variables
-
static const double SURFACE_HEIGHT = 1948.07
-
I3_MODULE(I3TrueMuonEnergy)¶
- file I3TrueMuonEnergy.h
- #include “icetray/I3ConditionalModule.h”#include “dataclasses/physics/I3Particle.h”#include “dataclasses/geometry/I3Geometry.h”#include “dataclasses/physics/I3MCTree.h”#include “dataclasses/physics/I3MCTreeUtils.h”#include “simclasses/I3MMCTrack.h”#include <icetray/I3Logging.h>#include <string>#include <map>#include <vector>#include <algorithm>#include “ddddr/I3MuonEnergyParams.h”#include “ddddr/I3MuonEnergyCascadeParams.h”#include “gulliver/I3FitParameterInitSpecs.h”#include “gulliver/I3GulliverBase.h”#include “lilliput/minimizer/I3GSLSimplex.h”#include “ddddr/ExpoFcn.h”
Implementation of I3MuonEnergy (DDDDR)
- Date
$Date: 14/01/23
- Author
Hans-Peter Bretz (hbretz@icecube.wisc.edu)
- Author
Patrick Berghaus
- file MuonEnergyFunctions.cxx
- #include “dataclasses/I3Map.h”#include “dataclasses/physics/I3RecoPulse.h”#include “ddddr/I3MuonEnergy.h”#include “ddddr/MuonEnergyFunctions.h”#include “phys-services/I3Calculator.h”#include “dataclasses/physics/I3MCTree.h”#include “dataclasses/physics/I3MCTreeUtils.h”#include <iostream>#include <fstream>
- file MuonEnergyFunctions.h
- #include <iostream>#include <fstream>#include “dataclasses/physics/I3Particle.h”#include “simclasses/I3MMCTrack.h”
- file MuonGunTrack.cxx
- #include <ddddr/MuonGunTrack.h>#include <simclasses/I3MMCTrack.h>#include <boost/foreach.hpp>
- Rcs
Track.cxx 128654 2015-02-04 18:34:51Z jvansanten
- Author
Jakob van Santen vansanten@wisc.edu
- Rcs
128654
- Rcs
2015-02-04 19:34:51 +0100 (Wed, 04 Feb 2015)
- file MuonGunTrack.h
- #include <dataclasses/physics/I3Particle.h>#include <dataclasses/physics/I3MCTree.h>#include <boost/tuple/tuple.hpp>
- Rcs
Track.h 128654 2015-02-04 18:34:51Z jvansanten
- Author
Jakob van Santen vansanten@wisc.edu
- Rcs
128654
- Rcs
2015-02-04 19:34:51 +0100 (Wed, 04 Feb 2015)
Typedefs
-
typedef I3Vector<I3MMCTrack> I3MMCTrackList
Functions
-
I3_FORWARD_DECLARATION(I3MMCTrack)
- file NoGradientFcn.h
- #include “gulliver/I3GulliverBase.h”#include <vector>
- file TomFFcn.cxx
- #include “ddddr/TomFFcn.h”#include “ddddr/TomFFunction.h”#include “ddddr/MuonEnergyFunctions.h”#include <cassert>
- file TomFFcn.h
- #include “gulliver/I3GulliverBase.h”#include “NoGradientFcn.h”#include <vector>
- file TomFFunction.h
- #include <math.h>
- dir ddddr
- dir ddddr
- dir ddddr
- dir icetray
- dir private
- dir public