gulliver C++ API Reference

class FastLogSum
#include <FastLogSum.h>

Accelerate evaluation of sums of logs, avoiding numerical limits.

In log-likelihood based reconstruction, we often need to compute sums of logs: L = sum_i log(x_i) Mathematically that is equivalent to the log of a product: L = log( prod_i x_i ) In principle the latter would be much faster to compute, because there is only one log computation involved. However, the practical difficulty here is that one can have many small terms x_i, which, when multiplied, might easily exit the available range for (double) floating point numbers. This class circumvents this difficulty by rescaling the product back to the [1.0,2.0) ranger after every x_i multiplication. Thanks to the way ieee754 double precision floats are laid out in memory, this is a very cheap operation. Using the ieee754_double “union” it is possible to access the bit patterns in a platform independent way (thanks, Troy!).

The speed gain from using this class depends on the length of the sums (the more terms the better), the platform (hardware, OS, compiler) and the relative fraction of the CPU time that your program spends on computing (sums of) logs. E.g. if the computation of x_i itself takes a long time, then the gain from using FastLogSum is negligible, but if x_i is a relatively simple expression or function that is cheaper than a log, then you may expect speedups of a factor 2 to 3 for logsums of a hundred x_i values.

Numerical security: normally when you evaluate an expression with many terms/factors, the result will be NAN if only one of the inputs is equal to NAN. In FastLogSum this is not the case, if one of the x_i is bad (NAN) then the result is still a normal number. The Add(x_i) method checks that each input x_i is not NAN, INF or negative. The FastAdd(x_i) leaves out those checks and is 20-40% faster. We use FastAdd for processing because we think the computation of the x_i is safe enough. If you are debugging/developing likelihood code you could consider using Add instead.

Public Functions

inline FastLogSum()

constructor

inline ~FastLogSum()

destructor

inline void Reset()

set the product back to 1 (so the sum_i log(x_i) = 0).

inline void Add(double x_i)

Add a term/factor x_i (from log(prod_i x_i) == sum_i log(x_i)) The x_i input is checked against negativity, NAN and INF.

We assume that the double precision floating point variables are stored in memory according to the IEEE 754 convention, that is 1 sign bit, 11 exponent bits and 52 fraction bits. After multiplying with x_i we take out the exponent (powers of 2) from the product and add that to the powersum_ datamember.

inline void FastAdd(double x_i)

Add a term/factor x_i (from log(prod_i x_i) == sum_i log(x_i)) The x_i input is NOT checked against negativity/NAN/INF. This saves 20-50% more computing time on the logs. Use this only if you are absolutely confident that the x_i are positive and finite.

We assume that the double precision floating point variables are stored in memory according to the IEEE 754 convention, that is 1 sign bit, 11 exponent bits and 52 fraction bits. After multiplying with x_i we take out the exponent (powers of 2) from the product and add that to the powersum_ datamember.

inline double GetLogSum()

Compute the sum of the logs of the x_i collected sofar. If any of the x_i was NAN or INF or not positive, then the product will be NAN.

Private Members

bool allpos_

a flag which is true if all x_i were positive

bool infnan_
ieee754_double prod_
int powersum_
class I3EventHypothesis

copyright (C) 2006 the icecube collaboration $Id$

Reconstructions in IceTray should always be stored in an object of type I3Particle, and any additional variables which cannot be accommodated in such an object (or in its vector of children, for composite solutions) should be stored in a custom frame object with “fit parameters”. IceTray does not provide a method to unambiguously link these two, so that’s why we have this hyperparticle in Gulliver.

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

This object is for use in Gulliver services and modules, not for storage in the I3Frame; particle and nonstd are stored separately.

Note that the “nonstd” part is intended solely for fittable variables related to the physics of the kind of event you are trying to reconstruct; for instance, for a double bang maybe you would store the times and energies of the first and second bang. You would not store e.g. the number of the number of hits associated to the two bangs.

Todo:

Actually, why should we not allow unfittable parameters in nonstd?

Public Functions

inline I3EventHypothesis()

default constructor

inline I3EventHypothesis(const I3Particle &p)

partial copy constructor

Parameters:

p[in] This particle is copied into this hypothesis

inline I3EventHypothesis(I3ParticlePtr p, I3FrameObjectPtr nstd)

partial copy constructor

Parameters:
  • p[in] This pointer is copied, no deep/shallow copy is made.

  • nstd[in] This pointer is copied, no deep/shallow copy is made.

inline ~I3EventHypothesis()

destructor

Public Members

I3ParticlePtr particle

the standard physics variables

I3FrameObjectPtr nonstd

object holding any non-standard physics variables (may be NULL)

class I3EventLogLikelihoodBase

Base class for event log-likelihood functions (service)

You only need to know about this class if you would like to implement a (wrapper for a) likelihood function, to be used by a Gulliver-based reconstruction module. You do not need to know about it when you are implementing any other Gulliver components, such as a parametrizations and minimizers.

If you write a Gulliver-based reconstruction module then (in principle) you only need to know that this class exists and that you need to provide (a pointer to) an object of the I3MinimizerBase class to the I3Gulliver core of your module. This pure base class is the base class for event log-likelihood functions.

If you implement a new I3EventLogLikelihoodBase subclass, then please also make a service factory which creates & configures an object of this new class.

See also

I3Gulliver

See also

I3GulliverBase

Subclassed by I3EventLogLikelihoodCombiner

Public Functions

inline I3EventLogLikelihoodBase()

base class constructor (idle)

inline virtual ~I3EventLogLikelihoodBase()

base class destructor (idle)

virtual void SetGeometry(const I3Geometry &f) = 0

This should be called on each geometry frame The implementation should get the geometry from the frame and store it in a format which is convenient and efficient for calculating the likelihood

virtual void SetEvent(const I3Frame &f) = 0

This will or should be called before doing a new fit. The implementation should get the relevant event data (hits, pulses, waveforms, geometry) from the frame and store it in a format which is convenient and efficient for calculating the likelihood.

virtual double GetLogLikelihood(const I3EventHypothesis &ehypo) = 0

Get +log(likelihood) for a particular emission hypothesis

The event data should have been provided earlier via the SetEvent method.

If you implement a likelihood for some non-standard topology (e.g. double muons) then you should check that the nonstd of ehypo indeed contains an object of the class relevant the particular topology for which your likelihood function is intended.

For likelihood functions for vanilla hypothesis (muon, cascade,…) implementations should of course check that the particle has indeed the appropriate type/shape.

Todo:

Maybe add a bool TestHypothesisType(const I3EventHypothesis&) method, so that this type checking needs to be done only once per event instead of for each minimizer iteration.

inline virtual double GetLogLikelihoodWithGradient(const I3EventHypothesis &ehypo, I3EventHypothesis &gradient, double weight = 1)

Same as GetLogLikelihood, plus gradient calculation.

The log(likelihood) is computed just like GetLogLikelihood does, but also the derivatives w.r.t. any requested variables are computed.

The “gradient” argument serves both to specify which derivatives should be computed as well as to return the results. E.g. if derivatives w.r.t. x, y and z should be computed, then the calling code should make sure that all (floating point) datamembers of gradient.particle are set to 0. except for x, y and z. The GetLogLikelihoodWithGradient computes the x, y and z derivatives and stores them in gradient.particle.

Note that for likelihoods which use the ehypo.nonstd datamember, also gradients w.r.t. those nonstd variables should be computed (when requested).

Minor wrinkle: likelihoods should add gradients to the already-existing fields in the gradient hypothesis with the given weight.

If this functionality is implemented, then make sure to also overload the HasGradient method to return true.

inline virtual bool HasGradient()

If GetLogLikelihoodWithGradient is implemented, then this method should be overloaded to return “true”.

inline virtual bool HasHessian()
inline virtual bool HasHVP()
inline virtual double GetLogLikelihood(const I3EventHypothesis &ehypo, I3EventHypothesis *gradient, bool maximize_extra_dimensions, double weight = 1)

Some likelihood functions may be able to calculate gradients and/or analytically maximize out extra dimensions. Overload this function if yours is one of the fancy few.

inline virtual double GetLogLikelihoodHess(const I3EventHypothesis &ehypo, I3EventHypothesis &gradient, I3EventHypothesis &hessian)
inline virtual void GetLogLikelihoodHVP(const std::vector<double> &pars, const std::vector<double> &arb_vector, I3EventHypothesis &hvp_result)
virtual unsigned int GetMultiplicity() = 0

Get the multiplicity of the current input event (e.g. number of good hits). The number of degrees of freedom for the fit will be calculated as: multiplicity minus number of fit parameters.

inline virtual I3FrameObjectPtr GetDiagnostics(const I3EventHypothesis &fitresult)

Get any extra information that might be interesting for the analysis. This is very specific to the particulars of the likelihood implmentation. Gulliver modules will get this frame object through the fit result pointer from the I3Gulliver object, and if it is indeed defined (pointer not NULL) then they should store it into the frame.

virtual const std::string GetName() const = 0

tell your name

class I3EventLogLikelihoodCombiner : public I3ServiceBase, public I3EventLogLikelihoodBase

Auxiliary class for combining event log-likelihood functions (service)

This class can be used to combine two or more likelihood services which are Added and configured separately to the tray. The most obvious (and maybe only) application is the Bayesian likelihood reconstrucion, in which a weight function (which only depends on the event hypothesis, not on the event data) is combined with an arbitrary other likelihood, e.g. convoluted Pandel.

Configuration parameters:

  • InputLogLikelihoods List of names of log-likelihood services

  • RelativeWeights List of relative weights to apply

  • Multiplicity Multiplicity calculation mode

For the “multiplicity calculation mode” there are several options: you can give the name of one of the likelihood services, then the multiplicity as calculated by that service will also be used for the combined likelihood. Alternatively, you can specify “Max” (default), “Min” or “Sum”; with those options, the combiner lets all likelihood services compute the multiplicity and then defines the multiplicity of the combined likelihood as the maximum, minimum or sum of those mulitplicities, respectively.

See also

I3Gulliver

See also

I3GulliverBase

See also

I3ServiceBase

Public Functions

I3EventLogLikelihoodCombiner(const I3Context&)

constructor

virtual ~I3EventLogLikelihoodCombiner()

destructor

virtual void Configure()

This method is inherited from I3ServiceBase. The factory for this service is constructed using the I3SingleServiceFactory or I3MultiServiceFactory.

virtual void SetGeometry(const I3Geometry &geo)

This should be called for every geometry frame The likelihood combiner will not use the geometry itself, it will pass the geometry to the to-be-combined likelihood services.

virtual void SetEvent(const I3Frame &f)

This will or should be called before doing a new fit. The likelihood combiner will not use the frame itself, it will pass the frame on to the to-be-combined likelihood services.

This will or should be called before doing a new fit.

virtual double GetLogLikelihood(const I3EventHypothesis &ehypo)

Get +log(likelihood) for a particular emission hypothesis

This returns the weighted sum of the log-likelihood values returned by each of the configured likelihood services.

virtual bool HasGradient()

If GetLogLikelihoodWithGradient is implemented, then this method should be overloaded to return “true”.

virtual double GetLogLikelihoodWithGradient(const I3EventHypothesis &ehypo, I3EventHypothesis&, double weight)

Same as GetLogLikelihood, plus gradient calculation.

The log(likelihood) is computed just like GetLogLikelihood does, but also the derivatives w.r.t. any requested variables are computed.

The “gradient” argument serves both to specify which derivatives should be computed as well as to return the results. E.g. if derivatives w.r.t. x, y and z should be computed, then the calling code should make sure that all (floating point) datamembers of gradient.particle are set to 0. except for x, y and z. The GetLogLikelihoodWithGradient computes the x, y and z derivatives and stores them in gradient.particle.

Note that for likelihoods which use the ehypo.nonstd datamember, also gradients w.r.t. those nonstd variables should be computed (when requested).

Minor wrinkle: likelihoods should add gradients to the already-existing fields in the gradient hypothesis with the given weight.

If this functionality is implemented, then make sure to also overload the HasGradient method to return true.

virtual unsigned int GetMultiplicity()

Get the multiplicity of the current input event (e.g. number of good hits). The number of degrees of freedom for the fit will be calculated as: multiplicity minus number of fit parameters.

Need to think how to define this.

virtual I3FrameObjectPtr GetDiagnostics(const I3EventHypothesis&)

Get any extra information that might be interesting for the analysis. This is very specific to the particulars of the likelihood implmentation. Gulliver modules will get this frame object through the fit result pointer from the I3Gulliver object, and if it is indeed defined (pointer not NULL) then they should store it into the frame.

inline virtual const std::string GetName() const

tell your name

Private Types

enum MultiModes

classifier type: how to compute multiplicity of an event

Values:

enumerator MMode_Max
enumerator MMode_Min

max of multiplicities reported by llh services

enumerator MMode_Sum

min of multiplicities reported by llh services

enumerator MMode_Favorite

sum of multiplicities reported by llh services

multiplicity reported by one particular llh service

Private Functions

SET_LOGGER ("I3EventLogLikelihoodCombiner")

Private Members

std::vector<std::string> llhNamesVector_

the names of likelihoods to combine

std::vector<I3EventLogLikelihoodBasePtr> logLikelihoods_

the pointers to the likelihood services to combine

std::vector<double> weightsVector_

the respective weights of likelihoods to combine

enum I3EventLogLikelihoodCombiner::MultiModes multiMode_
std::string multiModeString_

option string for multiplicity mode

class I3FitParameterInitSpecs

specs for a fitting parameter (as seen by minimizer)

copyright (C) 2004 the icecube collaboration $Id$

Simple auxiliary struct/class to I3GulliverBase and I3Parametrization

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

See also

I3GulliverBase

See also

I3Parametrization

Public Functions

inline I3FitParameterInitSpecs(std::string name)
inline ~I3FitParameterInitSpecs()
inline bool operator==(const I3FitParameterInitSpecs &o)

Public Members

std::string name_

parameter name (some minimizers want that)

double initval_

initval start value

double stepsize_

stepsize step size (expected order of magnitude for initial variation)

double minval_

minimum value of parameter (NAN if unlimited)

double maxval_

maximum value of parameter (NAN if unlimited)

class I3Gulliver : public I3GulliverBase
#include <I3Gulliver.h>

Fitter core of the Gulliver utility.

You need to know about this class if you would like to implement a new Gulliver-based reconstruction module. You do not need to know about it when you are implementing Gulliver components, such as a likelihood functions, parametrizations and minimizers.

This class connects the physics-aware log-likelihood function, the physics-ignorant minimizer and the parametrization of the physics track together. Typically your module will have a I3Gulliver datamember, during Configuration you retrieve services for minimization, parametrizaton and likelihood calculation and you plug them into this Gulliver object. Then in the Physics method you call the SetEvent(const I3Frame&) method of the likelihood service, followed by one or more calls to I3Gulliver::Fit(I3LogLikelihoodFitPtr) for various seeds (stored in an I3EventHypothesis), from which you may derive the reconstruction result.

The “particle” and “nonstd” members of the I3LogLikelihoodFitPtr which you provide to the Fit method should have been initilialized with a seed track. After the Fit call they contain the physics part of the fit result, while the “fitparams” datamember contains the reconstruction diagnostics such as the (reduced) likelihood, number of degrees of freedom and the number of times that the minimizer evaluated the likelihoodfunction.

If you’d like to store additional fit diagnostics (e.g. the paraboloid results), then can e.g. derive a new fit parameter class from I3LogLikelihoodFitParams to add additional variables and fill those after the Fit call.

See also

I3EventHypothesis, I3GulliverBase, I3SimpleFitter

Public Functions

I3Gulliver(I3ParametrizationBasePtr par, I3EventLogLikelihoodBasePtr llh, I3MinimizerBasePtr m, std::string name)

constructor

Parameters:
  • par[in] Translator between physics object and array of doubles

  • llh[in] Log-Likelihood calculator (physics model)

  • m[in] minimizer (physics ignorant)

  • name[in]

~I3Gulliver()

destructor

inline void ChangeParametrization(I3ParametrizationBasePtr newpar)

set/change parametrization If you write a reconstruction module in which you’d like to do the same fit with different parametrizations, then you use this method to swap the parametrization (keeping the event llh service and the minimizer unchanged).

Parameters:

newpar[in] shared pointer to alternative parametrization service

inline void ChangeFunction(I3EventLogLikelihoodBasePtr newllh)

set/change log-likelihood function If you write a reconstruction module in which you’d like to do the same fit with different likelihood functions, then you use this method to swap the eventllh service (keeping the parametrization and the minimizer unchanged).

Parameters:

newllh[in] shared pointer to alternative eventllh service

inline void ChangeMinimizer(I3MinimizerBasePtr m)

set/change minimizer If you write a reconstruction module in which you’d like to do the same fit with different minimizers, then you use this method to swap the minizer service (keeping the parametrization and the likelihood unchanged).

Parameters:

m[in] shared pointer to alternative minimizer service

void UseGradients(bool)

Set up gradient support, checking that both the parametrization and likelihood also support that. If not, then this throws a log_fatal().

void UseHessians(bool)
void UseHVP(bool)
virtual double operator()(const std::vector<double> &par)

Get the negative log-likelihood function value (via parameter array) This implements the virtual method/operator inherited from I3GulliverBase. This method is called by the minimizer to evaluate the function varue for a particular set of parameter values.

The parameter values are converted into an event hypothesis (using the parametrization service), from which the event loglikelihood service will compute the function value (negative log-likelihood).

Parameters:

par[in] array of parameters

Returns:

-log(likelihood)

virtual double operator()(const std::vector<double> &par, std::vector<double> &grad)

Same as above, but also compute gradient w.r.t. parameters.

virtual double operator()(const std::vector<double> &par, std::vector<double> &grad, I3Matrix &hessian)
virtual void operator()(const std::vector<double> &par, const std::vector<double> &arb_vector, std::vector<double> &hvp_result)
inline double GetFunctionValue(const I3EventHypothesis &p)

Get the negative log-likelihood function value for a particular event hypothesis, using the actual physics object (instead of the array of parameters, as in operator() (see above).

This can be useful for gulliver-based modules, to compute -log(L) for some specific hypothesis in the current event.

Parameters:

p[in] const ref to event hypothesis object

Returns:

-log(likelihood)

inline double GetFunctionValue(I3EventHypothesisConstPtr ptr)

Get the negative log-likelihood function value for a particular event hypothesis, using the actual physics object (instead of the array of parameters, as in operator() (see above).

Parameters:

ptr[in] shared pointer to event hypothesis object

Returns:

-log(likelihood)

void SetGeometry(const I3Geometry &geo)

Pass the geometry to the llh service.

Parameters:

geo – geometry

int SetEvent(const I3Frame &f)

Pass the event & geometry to the llh service.

Parameters:

f – current frame

Returns:

number of degrees of freedom

bool Fit(I3LogLikelihoodFitPtr fitptr)

do fit, starting from specified seed

Before you call this, either the SetEvent(f) or the Fit(f,fitptr) method should have been called for the current frame/event.

See also

I3Gulliver::GetMinimizerDetails() to get more diagnostics about fit result

Parameters:

fitptr – both input and output: it’s assumed to be initialized with some seed and will hold the result of the fit (if it converged)

Returns:

true if converged, false if fit failed

bool Fit(const I3Frame &f, I3LogLikelihoodFitPtr fitptr)

This convenience method is equivalent to running SetEvent(f) and Fit(fitptr), respectively.

inline I3MinimizerResult GetMinimizationDetails() const

get details about fit errors, number of iterations etc., of the latest fit

This is for debugging your recomodule. Also, if you’d like to store the fit errors, and the minimizer actually computes them, then you can retrieve them from this object.

Todo:

for simplicity I just return a copy now; maybe it should be a const ref or so.

void Trace()

Enable tracing, in order to debug the minimization process. This call also resets (reinitializes) the trace.

See also

GetTrace, NoTrace

void NoTrace()

Disable tracing

Discards any trace collected so far, and for upcoming fits no trace information will be stored.

See also

GetTrace, Trace

I3VectorDoublePtr GetTrace()

Return the current trace. This returns a (smart pointer to a) simple vector of doubles.

With npar fit parameters and nmini function calls, the trace has a length of (npar+1)*nmini. For every function call the npar parameter values and the llh value is stored, so in order to retrieve the function value of the 7th function call in a 5 parameter fit you would check the 41st (41=(7-1)*(5+1)+5) element of the vector.

Note that a new I3VectorDouble object will be made every time when the trace is reset using the Trace() call. So if you are interested in many fits performed with an I3Gulliver object, then it is not sufficient to retrieve the I3VectorDoublePtr only at the beginning and checking its contents many times. Before every fit you call Trace() to reinitialize the trace vector, and after each fit you retrieve it.

See also

Trace, NoTrace

Public Static Functions

static bool AnglesInRange(I3Particle &p, const std::string &name)

convenience method: set theta between 0 and pi, phi between 0 and 2pi (would be nice to have this as a I3Particle method) DJB Oct 31: return false if zenith/azimuth is pathologically unfixable Current definition of unfixable: fabs(angle)>1e6. It’s unfixable because for large floats, 2pi gets smaller than the floating point number precision.

Private Functions

SET_LOGGER ("I3Gulliver")
void NANParError(const std::vector<double> &par)

If the function call with NAN parameters occurs then this function is used to print the error message. (This should never happen, unless there is a critical error in the parametrization, minimizer or likelihood.)

void CheckGradientCompatibility()

Checks whether minimizer wants to use gradients, and if so, that the parametrization and likelihood also support that. If not, then this throws a log_fatal().

Checked at the start of every fit, after checking that the event multiplicity is OK.

void CheckHessianCompatibility()
void CheckHVPCompatibility()

Private Members

I3ParametrizationBasePtr parametrization_

smart pointer to the parametrization agent for a particle p

I3EventLogLikelihoodBasePtr eventllh_

smart pointer to the provider of the event log-likelihood for a particle p

I3MinimizerBasePtr minimizer_

smart pointer to the minimizer

I3MinimizerResult minimizerResult_

minimizer details of latest fit

unsigned int nFunctionEvaluations_

Counting how often the minimizer queries the likelihood function with proper parameter values, for the current fit. Will be reset to zero at the beginning of every Fit() call.

unsigned int nFunctionEvaluationsTotal_

Counting how often the minimizer queries the likelihood function with proper parameter values, for the lifetime of this I3Gulliver object (typically the duration of your processing script).

If the log-level for I3Gulliver is set to INFO or lower, then this number gets printed when the I3Gulliver object goes out of scope (end of run, typically).

unsigned int nNANfromMinimizer_

Count how often the minimizer queries the likelihood function with NAN parameter values, for the current fit. Will be reset to zero at the beginning of every Fit() call.

unsigned int nNANfromMinimizerTotal_

Count how often the minimizer queries the likelihood function with NAN parameter values, for the lifetime of this I3Gulliver object (typically the duration of your processing script).

If the log-level for I3Gulliver is set to INFO or lower, and there were indeed NAN parameter values, then this number gets printed when the I3Gulliver object goes out of scope (end of run, typically).

bool withGradient_

Fitting with or without gradient?

bool withHessian_
bool withHVP_
const std::string name_

Identifier/name to use in log messages. E.g. name of module.

I3VectorDoublePtr trace_

Trace information.

class I3GulliverBase
#include <I3GulliverBase.h>

Functor interface to likelihood function, for the minimizer algorithms.

copyright (C) 2004 the icecube collaboration $Id$

This base class is the medium through which the minimizer talks with the likelihood function and the “recomethod” aka “strategy” aka “global minimizer”.

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

You only need to know about this class if you would like to implement a (Gulliver wrapper of a) minimizer algorithm. You do not need to know about it if you are implementing other Gulliver components, such as a likelihood functions and parametrizations, or if you are writing a new Gulliver-based reconstruction module.

This virtual base class is the medium through which the minimizer evaluates the likelihood function. Its only implementation is the I3Gulliver class, which uses the parametrization object to convert a vector of doubles to a physics object (I3EventHypothesis) and evaluates the likelihood function with that. This base class does not depend on any physics-related classes.

Subclassed by I3Gulliver

Public Functions

inline I3GulliverBase()

constructor (idle)

inline virtual ~I3GulliverBase()

destructor (idle)

virtual double operator()(const std::vector<double> &par) = 0

Get function value

Parameters:

par[in] vector with fit parameter values

Returns:

function value

virtual double operator()(const std::vector<double> &par, std::vector<double> &grad) = 0

Get function value and its gradient

Parameters:
  • par[in] vector with fit parameter values

  • grad – with gradients w.r.t. parameters

Returns:

function value

virtual double operator()(const std::vector<double> &par, std::vector<double> &grad, I3Matrix &hessian) = 0
virtual void operator()(const std::vector<double> &par, const std::vector<double> &arb_vector, std::vector<double> &hvp_result) = 0
class I3LogLikelihoodFit

complete fit result for an arbitrary event type (“emission hypothesis”)

copyright (C) 2004 the icecube collaboration $Id$

You need to know about this class if you would like to implement a new Gulliver-based reconstruction module. You do not need to know about it when you are implementing Gulliver components, such as a likelihood functions, parametrizations and minimizers.

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

This class has just two datamembers, one just for the physics result, and one for the details and diagnostics of the likelihood fit.

Public Functions

inline I3LogLikelihoodFit()

dummy constructor

inline I3LogLikelihoodFit(I3EventHypothesisPtr ehypo, I3LogLikelihoodFitParamsPtr params)

handy constructor: initialize with seed pointer

Parameters:
  • ehypo[in] Pointer to event hypothesis (the pointer itself is used, no deep/shallow copy of the object is made)

  • params[in] Pointer to fit params (the pointer itself is used, no deep/shallow copy of the object is made)

inline I3LogLikelihoodFit(const I3EventHypothesis &ehypo, const I3LogLikelihoodFitParams &params)

handy constructor: initialize with seed reference

Parameters:
  • ehypo[in] Reference to event hypothesis (NOTE: only a shallow copy of the object is made)

  • params[in] Reference to fit params (full copy of the object is made)

inline I3LogLikelihoodFit(const I3EventHypothesis &ehypo)

handy constructor: initialize with seed reference

Parameters:

ehypo[in] Reference to event hypothesis (NOTE: only a shallow copy of the object is made)

inline I3LogLikelihoodFit(const I3Particle &p, const I3LogLikelihoodFitParams &params)

handy constructor for vanilla event types: initialize with seed reference

Parameters:
  • p[in] Reference to particle (full copy of the particle object is made)

  • params[in] Reference to fit params (full copy of the params is made)

Public Members

I3EventHypothesisPtr hypothesis_

the fit (to be)

I3LogLikelihoodFitParamsPtr fitparams_

object holding statistics and reconstruction info: likelihood etc.

I3FrameObjectPtr llhdiagnostics_

object holding specific diagnostics from the llh service

I3FrameObjectPtr minidiagnostics_

object holding specific diagnostics from the minimizer service

I3FrameObjectPtr paradiagnostics_

object holding specific diagnostics from the parametrization service

Public Static Attributes

static const std::string PARTICLE_SUFFIX = ""
static const std::string NONSTD_SUFFIX = "Params"
static const std::string FITPARAMS_SUFFIX = "FitParams"
static const std::string PARTICLEVECT_SUFFIX = "Vect"
static const std::string FITPARAMSVECT_SUFFIX = "FitParamsVect"
static const std::string NONSTDVECT_SUFFIX = "ParamsVect"
class I3LogLikelihoodFitParams : public I3FrameObject

Likelihood result and fitting diagnostics.

copyright (C) 2004 the icecube collaboration $Id$

This (base) class holds some generic parameters which are relevant for most or all loglikelihood reconstructions. The datamembers can be set by the Gulliver fit; for a particular sophisticated fitter (e.g. paraboloid) you might want to derive a fitparams from this, or just define a separate result class.

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

Todo:

Think of a way to include fitting errors

Maybe the minimizer algorithm should create this object

Public Functions

inline I3LogLikelihoodFitParams(double l = NAN, double r = NAN, int dof = -1, int mini = -1)

constructor

inline virtual ~I3LogLikelihoodFitParams()

destructor

std::ostream &Print(std::ostream&) const override
virtual void Reset()

set logl_ and rlogl_ to NAN, ndof_ and nmini_ to -1

Public Members

double logl_

-log(total event likelihood) (sometimes named “chi^2”).

double rlogl_

Reduced -log(event likelihood) (sometimes named “rchi^2”). Defined as logl divided by the number of degrees of freedom.

The number of degrees of freedom is usually the number of selected “hits” &#8212; or pulses, waveform fragments &#8212; minus the number of fit parameters.

int ndof_

the number of degrees of freedom for the fit zero if fit failed negative if this is for some reason ill-defined

int nmini_

the number of times that the likelihood function was evaluated by the minimizer useful diagnostic (variable for NN?), depends also on minimizer algorithm though

Protected Functions

template<class Archive>
void serialize(Archive &ar, unsigned version)

copyright (C) 2004 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

SET_LOGGER ("I3LogLikelihoodFitParams")

Friends

friend class icecube::serialization::access
class I3LogLikelihoodFitParamsConverter : public I3ConverterImplementation<I3LogLikelihoodFitParams>

copyright (C) 2010 The Icecube Collaboration

$Id$

Version

$Revision$

Date

$LastChangedDate$

Author

Eike Middell eike.middell@desy.de $LastChangedBy$

Private Functions

I3TableRowDescriptionPtr CreateDescription(const I3LogLikelihoodFitParams &params)

copyright (C) 2010 The Icecube Collaboration

$Id$

Version

$Revision$

Date

$LastChangedDate$

Author

Eike Middell eike.middell@desy.de $LastChangedBy$

size_t FillRows(const I3LogLikelihoodFitParams &params, I3TableRowPtr rows)
class I3LogLRatioFilter : public I3IcePick

An I3IcePick to select events with a maximum log(likelihood ratio). The ratio logL1-logL2 is calculated from two different likelihood reconstructions. Both reconstructions have to be applied before this modul and the resulting particle keys have to be contained in the frame.

Public Functions

explicit I3LogLRatioFilter(const I3Context &context)
virtual ~I3LogLRatioFilter()
void Configure()
bool SelectFrame(I3Frame &frame)

Private Functions

I3LogLRatioFilter(const I3LogLRatioFilter&)
I3LogLRatioFilter &operator=(const I3LogLRatioFilter&)
SET_LOGGER ("I3LogLRatioFilter")

Private Members

std::string particleKey1_
std::string particleKey2_
double maxLogLRatio_
class I3MinimizerBase
#include <I3MinimizerBase.h>

base class for minimizer algorithms, for generic loglh reconstruction

copyright (C) 2004 the icecube collaboration $Id$

This base class is the base class for the interfaces to minimizer algorithms used in log-likelihood based reconstruction.

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

If you would like to implement a new minimizer algorithm, you basically need to implement the “Minimizer(..)” method.

See also

I3GulliverBase

See also

I3Gulliver

If you write a Gulliver-based reconstruction module then you only need to know that this class exists and that you need to provide (a pointer to) an object of the I3MinimizerBase class to the I3Gulliver core of your module.

If you implement a new I3MinimizerBase subclass, then please also make a service factory which creates & configures an object of this new class.

Todo:

maybe let the minimizer tell how close it thinks it is to the minimum

Public Functions

inline I3MinimizerBase()
inline virtual ~I3MinimizerBase()
virtual void SetTolerance(double newtol) = 0
virtual double GetTolerance() const = 0
virtual unsigned int GetMaxIterations() const = 0

Get/Set max number of iterations What constitutes an “iteration” may differ per minimizer algorithm. It should be a measure for how hard the minimizer will try to find a minimum.

We do not actually store the number of iterations; rather the number of times that the minimizer evaluates the likelihood function.

See also

I3Gulliver

virtual void SetMaxIterations(unsigned int) = 0
virtual const std::string GetName() const = 0
virtual I3MinimizerResult Minimize(I3GulliverBase &gullifunctor, const std::vector<I3FitParameterInitSpecs> &parinit) = 0

Finds a minimum of a function, given initial parameter specs. The function simply takes a vector of doubles.

It is the caller’s responsibility to make sure that the number of parameters expected by the function is equal to the number of provided parameter specs. (This may be solved a bit more elegantly; but e.g. a template<int npar> won’t do here, because then the number of fit parameters can only be some const number in the caller’s code. Too restrictive.)

Todo:

is this a convenient enough interface?

Parameters:
  • gullifunctor – The interface to the likelihood function

  • parinit – Specs of number, names, initvalues & bounds of fit parameters

Returns:

an object with convergence info, fit parameters & errors

inline virtual bool UsesGradient()

For a minimizer implementation that uses the gradient w.r.t. the parameters, this method should return true. This method is called just before Minimize(), so it is possible to have a minimizer that only uses the gradient depending on e.g. the configuration by the user.

inline virtual bool UsesHessian()
inline virtual bool UsesHVP()
class I3MinimizerResult

Holds results of a minimization. Could be a struct, except that I think it’s convenient that the result vectors are already intialized to the right length and with a meaningful default result.

Public Functions

inline I3MinimizerResult(int npar)
inline virtual ~I3MinimizerResult()

Public Members

bool converged_

Success status of minimization.

double minval_

Value of the found minimum (NAN if not converged)

std::vector<double> par_

Parameter values which yielded the found minimum.

std::vector<double> err_

Uncertainties in parameter values which yielded the found minimum.

I3FrameObjectPtr diagnostics_

Extra diagnostics, specific to the minimizer algorithm, that could possibly be interesting for analysis.

Private Functions

I3MinimizerResult()
class I3ParametrizationBase

This base class defines the methods of parametrized particles which are relevant for the event loglikelihood “functor” classes (which mediate between the physics likelihood function and the non-phyisics minimizer).

The physics of the “emission hypothesis” is given in an I3EventHypothesis, which is just struct with an I3ParticlePtr with regular particle information and an I3FrameObjectPtr datamember to accommodate any non-standard fitvariables.

If you define a new parametrization, you only need to define the two methods which translate from physics variables to minimizer parameters and vice versa, and in the constructor you should initialize the parameters (par_) and parameter specs (parspecs_). Oh yes, you should also implement the trivial GetName() method.

If you (implementer of new parametrization) actually use the non-standard part of the I3EventHypothesis, then it is up to you to check that the “nonstd” I3FrameObjectPtr is indeed of the type for which your parametrization is intended.

Todo:

Maybe implement a Test() method which parametrizes back and forth and checks that the same result is obtained.

Public Functions

virtual const std::string GetName() const = 0

just return a name, useful for diagnostic blurbs

inline virtual bool InitChainRule(bool wantgradient, bool wanthessian)

Initilialize gradient capability.

If wantgradient is false, then no gradients are desired (until further notice) and any operations for applying the chain rule can be skipped. Return value is ignored in this case.

If wantgradient is true, then gradients are desired and any operations needed for applying the chain rule should be performed. In particular, the gradient_ data member should be initialized such that all datamembers of gradient_.particle (and gradient.nonstd_) that will be used in the chain rule calculation are set to zero.

If the parametrization subclass indeed provides the ApplyChainRule method (possibly depending on user configuration), then the return value should be true. Otherwise it should return false.

Note

I3Gulliver calls this init method for every minimizer function call.

inline virtual void SetEvent(const I3Frame &f)

Some parametrizations may want to access frame information at the start of an event. E.g. in order to align a parametrization to another fit or to the COG of the pulses.

inline virtual I3FrameObjectPtr GetDiagnostics(const I3EventHypothesis &fitresult)

Get any extra information that might be interesting for the analysis. This is very specific to the particulars of the parametrization. Gulliver modules will get this frame object through the fit result pointer from the I3Gulliver object, and if it is indeed defined (pointer not NULL) then they should store it into the frame.

inline virtual void PassCovariance(const boost::numeric::ublas::symmetric_matrix<double> &cov)

This is called by Gulliver after the minimization has finished to provide the covariance matrix. Support in Gulliver is tentative at this point, you might get only NaNs, or only the diagonal elements. For backward compatibility, the default implementation of this function does nothing. To do something useful with the passed covariance matrix, overwrite this function in your derived class.

inline I3ParametrizationBase(I3EventHypothesisPtr eh)

constructor: an already existing I3EventHypothesis object is parametrized.

inline virtual ~I3ParametrizationBase()

destructor

void SetHypothesisPtr(I3EventHypothesisPtr eh)

Tie the parametrization to a different fit object.

Parameters:

eh[in] pointer of new hypothesis (NULL pointers are forbidden) Note that the pointer itself is used, no copy is made. The pointee will be modified by the parametrization service.

inline I3EventHypothesisPtr GetGradientPtr()

Get pointer to current event hypothesis.

inline I3EventHypothesisPtr GetHessianPtr()

Get pointer to current event hypothesis.

inline I3EventHypothesisPtr GetHypothesisPtr()

Get pointer to current event hypothesis.

I3EventHypothesisConstPtr GetHypothesisPtr(const vector<double> &par)

this calculates the physics particle from the parameter vector and returns it

void GetGradient(vector<double> &grad)

gradient (in par) is computed and copied into grad.

void GetHessian(vector<double> &grad, I3Matrix &hessian)

gradient (in par) is computed and copied into grad.

inline unsigned int GetNPar() const

Get number of free parameters.

inline virtual bool CheckEnergyFit()
const std::vector<I3FitParameterInitSpecs> &GetParInitSpecs(I3EventHypothesisPtr eh = I3EventHypothesisPtr())

update the parameters according to the current event hypothesis then return a I3FitParameterInitSpecs object which contains besides the parameter values also the ranges, stepsizes etc. Note: the default input to this fill is supposed to create a null pointer

SET_LOGGER ("I3ParametrizationBase")

macro for log_message system

Protected Functions

virtual void UpdatePhysicsVariables() = 0

This should calculate physics datamembers of the particle (in I3Particle and maybe I3FrameObject) from the values in the parametrization (par_).

virtual void UpdateParameters() = 0

This should calculate the values in par_ from datamembers of the particle If relevant it should also update stepsizes.

inline virtual void ApplyChainRule()

This should calculate the values in par_gradient_, using the current gradient_ and hypothesis_.

inline virtual void ApplyChainRuleHessian()

Protected Attributes

I3EventHypothesisPtr hypothesis_

the physics variables

I3EventHypothesisPtr gradient_

the gradient in physics variables (dlogl/dhypothesis)

I3EventHypothesisPtr hessian_
std::vector<I3FitParameterInitSpecs> parspecs_

parametrization info for the minimizer: stepsizes etcetera.

std::vector<double> par_

the parametrization of the free physics variables

std::vector<double> par_gradient_

the gradient in parameters (dlogl/dpar)

boost::numeric::ublas::symmetric_matrix<double> par_hessian_

the hessian in parameters

Private Functions

I3ParametrizationBase()

Prohibit the default constructor by making it private. A parameterization should refer to an already existing object, so it makes no sense to have a constructor without particle argument.

class I3PDFBase
#include <I3PDFBase.h>

This class is an interface to the PDFs, providing the wf module with a single function to call regardless of the PDF.

Public Functions

inline I3PDFBase()
inline virtual ~I3PDFBase()
virtual void GetProbability(const I3Position &ompos, const I3Position &hypo_pos, const I3Direction &hypo_dir, double hypo_time_residual, double hypo_energy, double &expected_amplitude, double &time_probability) = 0

GetProbability is the virtual function overloaded in the derived classes of I3PDFBase. Interface to the LLH reconstruction is always the same, implementation of function depends on the PDF. Currently, Photorec and patched pandel are supported. Input arguments: I3Position ompos … position of current OM I3Position hypo_pos … position of hypothesis particle I3Direction hypo_dir … direction of hypothesis particle double hypo_time_residual … time residual between direct hit time of hypothsis particle and actual hit_time double hypo_energy … energy of hypothesis particle

Returns:

expected_amplitude and time_probability of the hypothesis

inline virtual void GetProbability(const I3Position &ompos, const I3Particle &track, double hypo_time_residual, double &expected_amplitude, double &time_probability)
class I3PLogLFilter : public I3IcePick
#include <I3PLogLFilter.h>

An I3IcePick to select events which have a particle with reduced log(likelihood) obtained by a likelihood fit smaller than the given value MaxPLogL. Note that a small reduced log(likelhood) indicates a good reconstruction quality.

Public Functions

explicit I3PLogLFilter(const I3Context &context)
virtual ~I3PLogLFilter()
void Configure()
bool SelectFrame(I3Frame &frame)

Private Functions

I3PLogLFilter(const I3PLogLFilter&)
I3PLogLFilter &operator=(const I3PLogLFilter&)
SET_LOGGER ("I3PLogLFilter")

Private Members

std::string pulseSeriesName_
std::string particleKey_
std::string cutValuesKey_
double subFromNCh_
double maxPLogL_
class I3RLogLFilter : public I3IcePick
#include <I3RLogLFilter.h>

An I3IcePick to select events which have a particle with reduced log(likelihood) obtained by a likelihood fit smaller than the given value MaxRLogL. Note that a small reduced log(likelhood) indicates a good reconstruction quality.

Public Functions

explicit I3RLogLFilter(const I3Context &context)
virtual ~I3RLogLFilter()
void Configure()
bool SelectFrame(I3Frame &frame)

Private Functions

I3RLogLFilter(const I3RLogLFilter&)
I3RLogLFilter &operator=(const I3RLogLFilter&)
SET_LOGGER ("I3RLogLFilter")

Private Members

std::string particleKey_
double maxRLogL_
class I3SeedServiceBase

Base class seed track collection & preparation services.

Seed services are supposed to collect the first guess tracks, and do any desired transformations on them.

See also

I3Gulliver

See also

I3GulliverBase

Public Functions

inline I3SeedServiceBase()

Base class constructor (idle)

inline virtual ~I3SeedServiceBase()

Base class destructor (idle)

virtual unsigned int SetEvent(const I3Frame &f) = 0

Provide event data:

  • purge old seed tracks

  • get the first guess tracks

  • fill in the blanks

  • do tweaks (time/space corrections)

The generated seeds should all be good, have OK fit status and position, direction and time datamembers should not be NAN.

Todo:

Think about no-NAN requirement for energy and length.

Parameters:

f[in] Frame with event data

Returns:

number of available seeds

virtual I3EventHypothesis GetSeed(unsigned int iseed) const = 0

Return a seed. May call log_fatal() if iseed is larger than the number of available seeds.

When you implement this method: make sure that the returned seed is a deep copy (not a shallow copy) of some internally stored seed, because the receiver will most likely change its values. (This method is declared const to enforce this requirement…)

inline virtual I3EventHypothesisPtr GetSeedPtr(unsigned int iseed) const
inline virtual I3EventHypothesis GetDummy() const

Get dummy seed; useful in case all first guesses fail.

In particular if the seeds normally have a non-NULL “nonstd” data member, then this should be initialized in this dummy seed as well. (I would find it more logical and safer to have reconstruction results always in some kind of sortable result container, which could be set to zero results in case of total failure. But that does not seem to be compatible with with the IceTray design: there is no possibility to store a vector of fit parameters.)

When you implement this method: make sure that the returned dummy seed is a deep copy (not a shallow copy) of some internally stored seed, because the receiver will likely change its values. (This method is declared const to enforce this requirement…)

inline virtual I3EventHypothesisPtr GetDummyPtr() const
inline virtual I3EventHypothesis GetCopy(const I3EventHypothesis &eh) const

Get a (deep) copy of an existing event hypothesis.

The caller (a log-likelihood fitter module) does not need to know what kind of track exactly it is. E.g. for a tau seed service, this method correctly copies over the tau-specific variables in the “nonstd” data member. (This is not so relevant for vanilla muons and cascades; you could also just copy the I3Particle information.) Call fatal if the hypothesis is of the wrong type.

The default implementation just makes a (deep) copy of the I3Particle data member of the I3EventHypothesis. The “nonstd” data member is not copied.

inline virtual I3EventHypothesisPtr GetCopyPtr(const I3EventHypothesis &eh) const
inline virtual void Tweak(I3EventHypothesis &eh) const

When a seed service generates a seed, it may perform some tweaks or transformations on it intended to facilitate the minimization. The generic example is the vertex correction for an infinite track: purely mathematically it is irrelevant where the vertex (base point) is chosen, one can shift the vertex over an arbitrary length L along the track, as long as the vertex time is adjusted by L/c. However, a minimizer will more easily navigate the phase space if the vertex is chosen well within the clouds of hit DOMs. It can also be good to do a correction to the vertex time, for instance for first guess algorithms which use the average of the hit times, one can improve by choosing a vertex time based on the distribution of time residuals.

Now, make sure you do those tweaks in this method. Because then gulliver modules which generate new seed solutions (e.g. iterative and paraboloid) can use it as well.

The default implementation does no tweaking.

inline virtual void FillInTheBlanks(I3EventHypothesis &eh) const

Many first guess algorithms do not provide the complete description for the kind of event we are trying to reconstruct. A simple example is the energy of a cascade or muon. The seed service might take care of filling in a sensible initial value for such a quantity.

For track types other than cascades or muons, there might be other less trivial preparations to do. E.g. for a “sugardaddy” tau track, the point where the tau decays into a muon should be set to some sensible value, in the “nonstd” data member of the event hypothesis.

The default implementation does no filling.

union ieee754_double
#include <ieee754.h>

Public Members

double d
unsigned int negative
unsigned int exponent
unsigned int mantissa0
unsigned int mantissa1
struct ieee754_double::[anonymous] ieee
unsigned int quiet_nan
struct ieee754_double::[anonymous] ieee_nan
union ieee754_float
#include <ieee754.h>

Public Members

float f
unsigned int negative
unsigned int exponent
unsigned int mantissa
struct ieee754_float::[anonymous] ieee
unsigned int quiet_nan
struct ieee754_float::[anonymous] ieee_nan
union ieee854_long_double
#include <ieee754.h>

Public Members

long double d
unsigned int negative
unsigned int exponent
unsigned int empty
unsigned int mantissa0
unsigned int mantissa1
struct ieee854_long_double::[anonymous] ieee
unsigned int one
unsigned int quiet_nan
struct ieee854_long_double::[anonymous] ieee_nan
namespace std

STL namespace.

file FastLogSum.h
#include <cmath>
#include <cassert>

copyright (C) 2007 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu (with advice from Troy Straszheim troy@resophonic.com and Martin Merck mmerck@icecube.wisc.edu)

Defines

CAST_DOUBLE_TO_IEEE754
file I3EventHypothesis.h
#include “icetray/IcetrayFwd.h”
#include “dataclasses/physics/I3Particle.h”
#include “icetray/I3FrameObject.h”

Functions

I3_POINTER_TYPEDEFS(I3EventHypothesis)
file I3EventLogLikelihoodBase.h
#include “boost/shared_ptr.hpp”
#include <string>
#include “icetray/IcetrayFwd.h”
#include “icetray/I3DefaultName.h”
#include “icetray/I3FrameObject.h”

copyright (C) 2004 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

Functions

I3_POINTER_TYPEDEFS(I3EventLogLikelihoodBase)
file I3EventLogLikelihoodCombiner.cxx
#include <algorithm>
#include <icetray/I3SingleServiceFactory.h>
#include <dataclasses/geometry/I3Geometry.h>
#include <boost/foreach.hpp>
#include <boost/python/list.hpp>

copyright (C) 2007 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

Typedefs

typedef I3SingleServiceFactory<I3EventLogLikelihoodCombiner, I3EventLogLikelihoodBase> I3EventLogLikelihoodCombinerFactory

Functions

I3_SERVICE_FACTORY(I3EventLogLikelihoodCombinerFactory)

Variables

static const char *inputllhs_optionname = "InputLogLikelihoods"
static const char *relweights_optionname = "RelativeWeights"
static const char *multiplicity_optionname = "Multiplicity"
file I3EventLogLikelihoodCombiner.h
#include <string>
#include <vector>
#include “icetray/IcetrayFwd.h”
#include “icetray/I3ServiceBase.h”

copyright (C) 2007 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

Functions

I3_POINTER_TYPEDEFS(I3EventLogLikelihoodCombiner)
file I3FitParameterInitSpecs.h
#include <string>
#include <cmath>
file I3Gulliver.cxx
#include “gulliver/I3Gulliver.h
#include <algorithm>
#include <functional>
#include <numeric>
#include <dataclasses/I3Matrix.h>
#include <cmath>
#include <boost/numeric/ublas/symmetric.hpp>

copyright (C) 2004 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

file I3Gulliver.h
#include <vector>
#include “icetray/IcetrayFwd.h”
#include “icetray/I3Logging.h”
#include <dataclasses/I3Matrix.h>

copyright (C) 2004 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

Functions

I3_POINTER_TYPEDEFS(I3Gulliver)
file I3GulliverBase.cxx
file I3GulliverBase.h
#include <vector>
#include “icetray/IcetrayFwd.h”
#include <dataclasses/I3Matrix.h>

Functions

I3_POINTER_TYPEDEFS(I3GulliverBase)
file I3LogLikelihoodFit.cxx

copyright (C) 2004 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

Functions

bool operator<(const I3LogLikelihoodFit &lhs, const I3LogLikelihoodFit &rhs)

“less than” comparison operator, for sorting solutions

bool operator>(const I3LogLikelihoodFit &lhs, const I3LogLikelihoodFit &rhs)

“greater than” comparison operator, for sorting solutions

file I3LogLikelihoodFit.h
#include <boost/shared_ptr.hpp>
#include <string>
#include “dataclasses/physics/I3Particle.h”
#include “icetray/I3FrameObject.h”

Functions

bool operator<(const I3LogLikelihoodFit &lhs, const I3LogLikelihoodFit &rhs)

“less than” comparison operator, for sorting solutions

bool operator>(const I3LogLikelihoodFit &lhs, const I3LogLikelihoodFit &rhs)

“greater than” comparison operator, for sorting solutions

I3_POINTER_TYPEDEFS(I3LogLikelihoodFit)
file I3LogLikelihoodFitParams.cxx
#include “icetray/serialization.h”

Functions

std::ostream &operator<<(std::ostream &oss, const I3LogLikelihoodFitParams &p)
I3_CLASS_VERSION (I3LogLikelihoodFitParams, 0)
I3_SERIALIZABLE(I3LogLikelihoodFitParams)
I3_SERIALIZABLE(I3LogLikelihoodFitParamsVect)
file I3LogLikelihoodFitParams.h
#include <boost/shared_ptr.hpp>
#include “dataclasses/Utility.h”
#include “icetray/I3FrameObject.h”
#include “icetray/I3Logging.h”
#include “dataclasses/I3Vector.h”
#include <cmath>

Typedefs

typedef I3Vector<I3LogLikelihoodFitParams> I3LogLikelihoodFitParamsVect

Functions

inline bool operator==(const I3LogLikelihoodFitParams &lhs, const I3LogLikelihoodFitParams &rhs)
std::ostream &operator<<(std::ostream &oss, const I3LogLikelihoodFitParams &d)
I3_POINTER_TYPEDEFS(I3LogLikelihoodFitParams)
I3_POINTER_TYPEDEFS(I3LogLikelihoodFitParamsVect)
file I3LogLikelihoodFitParamsConverter.cxx
file I3LogLikelihoodFitParamsConverter.h
#include “tableio/I3Converter.h”
file I3LogLRatioFilter.cxx
#include <icetray/I3Context.h>
#include <icetray/I3Frame.h>
#include <icetray/I3IcePickInstaller.h>
#include <interfaces/I3IcePickModule.h>
#include <interfaces/I3IceForkModule.h>
#include <dataclasses/physics/I3Particle.h>

copyright (C) 2007 the IceCube collaboration $Id:$

Version

$Revision$

Author

Anna Franckowiak

Date

Oct 22 2007

file I3LogLRatioFilter.h
#include <string>
#include <icetray/I3IcePick.h>
#include <icetray/I3Logging.h>

copyright (C) 2007 the IceCube collaboration $Id:$

Version

$Revision$

Author

Anna Franckowiak

Date

Oct 22 2007

file I3MinimizerBase.h
#include <vector>
#include <boost/shared_ptr.hpp>

Functions

I3_POINTER_TYPEDEFS(I3MinimizerBase)
file I3MinimizerResult.h
#include <cmath>
#include <vector>
#include “icetray/I3FrameObject.h”
file I3ParametrizationBase.cxx

copyright (C) 2007 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

file I3ParametrizationBase.h
#include <cassert>
#include <vector>
#include <boost/numeric/ublas/symmetric.hpp>
#include <boost/numeric/ublas/matrix.hpp>
#include <boost/numeric/ublas/io.hpp>
#include “icetray/IcetrayFwd.h”
#include <dataclasses/I3Matrix.h>

Functions

I3_POINTER_TYPEDEFS(I3ParametrizationBase)
file I3PDFBase.h
#include “dataclasses/I3Position.h”
#include “dataclasses/I3Direction.h”
#include “dataclasses/physics/I3Particle.h”

Functions

I3_POINTER_TYPEDEFS(I3PDFBase)
file I3PLogLFilter.cxx
#include <icetray/I3Context.h>
#include <icetray/I3Frame.h>
#include <icetray/I3IcePickInstaller.h>
#include <interfaces/I3IcePickModule.h>
#include <interfaces/I3IceForkModule.h>
#include <dataclasses/physics/I3Particle.h>
#include <dataclasses/physics/I3RecoPulse.h>
#include <phys-services/I3CutValues.h>

copyright (C) 2007 the IceCube collaboration $Id:$

Version

$Revision$

Author

Anna Franckowiak

Date

Oct 22 2007

Functions

I3_MODULE(I3IcePickModule<I3PLogLFilter>)
I3_MODULE(I3IceForkModule<I3PLogLFilter>)
I3_SERVICE_FACTORY(I3IcePickInstaller<I3PLogLFilter>)
file I3PLogLFilter.h
#include <string>
#include <icetray/I3IcePick.h>
#include <icetray/I3Logging.h>

copyright (C) 2007 the IceCube collaboration $Id:$

Version

$Revision$

Author

Anna Franckowiak

Date

Oct 22 2007

file I3RLogLFilter.cxx
#include <icetray/I3Context.h>
#include <icetray/I3Frame.h>
#include <icetray/I3IcePickInstaller.h>
#include <interfaces/I3IcePickModule.h>
#include <interfaces/I3IceForkModule.h>
#include <dataclasses/physics/I3Particle.h>

copyright (C) 2007 the IceCube collaboration $Id:$

Version

$Revision$

Author

Anna Franckowiak

Date

Oct 22 2007

Functions

I3_MODULE(I3IcePickModule<I3RLogLFilter>)
I3_MODULE(I3IceForkModule<I3RLogLFilter>)
I3_SERVICE_FACTORY(I3IcePickInstaller<I3RLogLFilter>)
file I3RLogLFilter.h
#include <string>
#include <icetray/I3IcePick.h>
#include <icetray/I3Logging.h>

copyright (C) 2007 the IceCube collaboration $Id:$

Version

$Revision$

Author

Anna Franckowiak

Date

Oct 22 2007

file I3SeedServiceBase.h
#include “boost/shared_ptr.hpp”
#include “icetray/IcetrayFwd.h”
#include “icetray/I3DefaultName.h”

copyright (C) 2004 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

Functions

I3_POINTER_TYPEDEFS(I3SeedServiceBase)
I3_DEFAULT_NAME(I3SeedServiceBase)
file ieee754.h

Defines

IEEE754_FLOAT_BIAS
IEEE754_DOUBLE_BIAS
IEEE854_LONG_DOUBLE_BIAS
file ordinal.cxx
#include <stdio.h>

copyright (C) 2007 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

Functions

const char *ordinal(unsigned int i)

utility function to pretty-print ordinal numbers So 1 returns “1st”, 5 returns “5th”, 122 returns “122nd” etc..

file ordinal.h

copyright (C) 2007 the icecube collaboration $Id$

Version

$Revision$

Date

$Date$

Author

David Boersma boersma@icecube.wisc.edu

Functions

const char *ordinal(unsigned int i)

utility function to pretty-print ordinal numbers So 1 returns “1st”, 5 returns “5th”, 122 returns “122nd” etc..

page todo

Class I3EventHypothesis

Actually, why should we not allow unfittable parameters in nonstd?

Member I3EventLogLikelihoodBase::GetLogLikelihood  (const I3EventHypothesis &ehypo)=0

Maybe add a bool TestHypothesisType(const I3EventHypothesis&) method, so that this type checking needs to be done only once per event instead of for each minimizer iteration.

Member I3Gulliver::GetMinimizationDetails  () const

for simplicity I just return a copy now; maybe it should be a const ref or so.

Class I3LogLikelihoodFitParams

Think of a way to include fitting errors

Maybe the minimizer algorithm should create this object

Class I3MinimizerBase

maybe let the minimizer tell how close it thinks it is to the minimum

Member I3MinimizerBase::Minimize  (I3GulliverBase &gullifunctor, const std::vector< I3FitParameterInitSpecs > &parinit)=0

is this a convenient enough interface?

Class I3ParametrizationBase

Maybe implement a Test() method which parametrizes back and forth and checks that the same result is obtained.

Member I3SeedServiceBase::SetEvent  (const I3Frame &f)=0

Think about no-NAN requirement for energy and length.

dir converter
dir gulliver
dir gulliver
dir gulliver
dir icetray
dir private
dir public
dir utilities