CommonVariables C++ API Reference

struct ArithMeanIterativeCalculator

Public Functions

inline ArithMeanIterativeCalculator()
inline double AsDouble()

Returns the calculated arithmetic mean value as a double value. The result of this method is only valid after the call of the FinalizeCalculation method.

inline void DoNextIteration(double value)

Does the next iteration step for the arithmetic mean value calculation.

inline void FinalizeCalculation()

Finalizes the arithmetic mean value calculation.

inline void Reset()

Resets the ArithMeanIterativeCalculator, in order to be able to use the same ArithMeanIterativeCalculator object more than once.

Public Members

uint64_t n_
double sum_
double arithMean_
struct COGIterativeCalculator

Public Functions

inline COGIterativeCalculator()
inline I3Position AsI3Position()

Returns the calculated COG as an I3Position object. The result of this method is only valid after the call of the FinalizeCalculation method, but this condition is not enforced.

inline void DoNextIteration(double omX, double omY, double omZ, double pulseCharge)

Does the next iteration step. This method is supposed to get called within the I3RecoPulseSeries loop for a particular pulse.

inline void FinalizeCalculation()

Finalizes the COG calculation. This method is supposed to get called after the I3RecoPulseSeriesMap loop.

inline void Reset()

Resets the COGIterativeCalculator, in order to be able to use the same COGIterativeCalculator object more than once.

Public Members

double x_
double y_
double z_
double pulseChargeSum_
struct DirectHitsClassData_t

A structure to store temporary data about the calculation of a particular class of direct hits.

Public Functions

inline DirectHitsClassData_t(const I3DirectHitsDefinition &dhDefinition_)

Public Members

I3DirectHitsDefinition dhDefinition

The definition of the class of the direct hits.

hm::details::NHitStringsIterativeCalculator nDirStringsIterCalc

The generated values for a particular direct hits time window.

uint32_t nDirDoms
uint64_t nDirPulses
double qDirPulses
hm::details::NHitStringsIterativeCalculator nEarlyStringsIterCalc
uint32_t nEarlyDoms
uint64_t nEarlyPulses
double qEarlyPulses
hm::details::NHitStringsIterativeCalculator nLateStringsIterCalc
uint32_t nLateDoms
uint64_t nLatePulses
double qLatePulses
bool currentDomGotEarlyPulses

The flags if the current processed DOM has early, direct, or late pulses for this particular direct hits time window.

bool currentDomGotDirectPulses
bool currentDomGotLatePulses
std::vector<double> dirHitDistancesAlongTrack

The vector with all direct hit distances along the track of the direct hit DOMs.

class I3DirectHitsDefinition

This class describes a definition of a class of direct hits.

Public Functions

I3DirectHitsDefinition()
I3DirectHitsDefinition(std::string name, double time0, double time1)
I3DirectHitsDefinition(const I3DirectHitsDefinition &dhd)

copy constructor

inline virtual ~I3DirectHitsDefinition()
inline std::string GetName() const
inline const I3DirectHitsTimeWindow &GetTimeWindow() const
inline I3DirectHitsTimeWindow &GetTimeWindow()
inline void SetName(std::string name)
inline void SetTimeWindow(const I3DirectHitsTimeWindow &timeWindow)
bool operator==(const I3DirectHitsDefinition &rhs) const

Protected Attributes

std::string name_

the name of the direct hits class

I3DirectHitsTimeWindow timeWindow_

the time window definition of this direct hits class

class I3DirectHitsTimeWindow

This class describes a time window to calculate the direct hits which fall inside this time window. NOTE: This class should be in principle a general I3TimeWindow class, but at the time this project has been written no such class did exist.

Public Functions

I3DirectHitsTimeWindow()
I3DirectHitsTimeWindow(double time0, double time1)
I3DirectHitsTimeWindow(const I3DirectHitsTimeWindow &tw)
inline double GetTime0() const
inline double GetTime1() const
inline void SetTime0(double time)
inline void SetTime1(double time)
inline bool IsTimeAfterTimeWindow(double time) const
inline bool IsTimeBeforeTimeWindow(double time) const
inline bool IsTimeInsideTimeWindow(double time) const
bool operator==(const I3DirectHitsTimeWindow &rhs) const

Protected Attributes

double time0_
double time1_
struct NHitStringsIterativeCalculator

Public Functions

inline NHitStringsIterativeCalculator()
inline uint32_t AsUInt32()

Returns the calculated number of hit strings as an uint32_t value. The result of this method is only valid after the call of the FinalizeCalculation method.

inline void DoNextIteration(int string)

Does the next iteration step for the number of hit strings calculation.

inline void FinalizeCalculation()

Finalizes the number of hit strings calculation.

inline void Reset()

Resets the NHitStringsIterativeCalculator, in order to be able to use the same NHitStringsIterativeCalculator object more than once.

Public Members

std::set<int> hitStrings_
uint32_t nHitStrings_
struct S1Variables_t

Public Functions

inline S1Variables_t(const I3Position &cog_, double cogZSigma_, double minPulseTime_, double maxPulseTime_, double qMaxDoms_, double qTotPulses_, double zMin_, double zMax_, double zMean_, double zSigma_, double zTravel_)
inline S1Variables_t(const S1Variables_t &rhs)

Public Members

I3Position cog
double cogZSigma
double minPulseTime
double maxPulseTime
double qMaxDoms
double qTotPulses
double zMin
double zMax
double zMean
double zSigma
double zTravel
struct SigmaIterativeCalculator

Public Functions

inline SigmaIterativeCalculator()
inline double AsDouble()

Returns the calculated sigma value as a double value. The result of this method is only valid after the call of the FinalizeCalculation method, but this condition is not enforced.

inline void DoNextIteration(double value, double weight)

Does the next iteration step for the sigma calculation.

inline void FinalizeCalculation()

Finalizes the sigma value calculation:

\( sigma = sqrt(variance) = sqrt( n_eq/(n_eq - 1) * (E[x^2] - E[x]^2)) \)

where \( n_eq = (\sum w_i)^2 / (\sum w_i^2) \) to account for the unbiased mean.

inline void Reset()

Resets the SigmaIterativeCalculator, in order to be able to use the same SigmaIterativeCalculator object more than once.

Public Members

uint64_t n_
double weightSum_
double weightSquaredSum_
double weightTimesValueSum_
double weightTimesValueSquaredSum_
double sigma_
struct TrackHitsSeparationLengthIterativeCalculator

Public Functions

inline TrackHitsSeparationLengthIterativeCalculator(const I3Particle &particle)
inline double AsDouble()

Returns the calculated TrackHitsSeparationLength variable as a double value. The result of this method is only valid after the call of the FinalizeCalculation method, but this condition is not enforced.

inline void DoNextIteration(double omX, double omY, double omZ, double hitTime, double hitCharge)

Does the next iteration step. This method is supposed to get called within the I3RecoPulseSeriesMap loop for a particular DOM. Usually the hitTime is the time of the first pulse, and hitCharge is the integrated charge of the OM.

inline void FinalizeCalculation()

Finalizes the TrackHitsSeparationLength calculation. This method is supposed to get called after the I3RecoPulseSeriesMap loop.

inline void Reset(const I3Particle &particle)

Resets the TrackHitsSeparationLengthIterativeCalculator, in order to be able to use the same TrackHitsSeparationLengthIterativeCalculator object more than once.

Public Members

I3Particle particle_
std::vector<double> omXs_
std::vector<double> omYs_
std::vector<double> omZs_
std::vector<double> hitTimes_
std::vector<double> hitCharges_
double trackHitsSeparationLength_
struct WeightedMeanIterativeCalculator

Public Functions

inline WeightedMeanIterativeCalculator()
inline double AsDouble()

Returns the calculated weighted mean value as a double value. The result of this method is only valid after the call of the FinalizeCalculation method.

inline void DoNextIteration(double value, double weight)

Does the next iteration step for the weighted mean value calculation.

inline void FinalizeCalculation()

Finalizes the weighted mean value calculation.

inline void Reset()

Resets the WeightedMeanIterativeCalculator, in order to be able to use the same WeightedMeanIterativeCalculator object more than once.

Public Members

double weightTimesValueSum_
double weightSum_
double weightedMean_
struct ZTravelIterativeCalculator

Public Functions

inline ZTravelIterativeCalculator()
inline double AsDouble()

Returns the calculated ZTravel variable as a double value. The result of this method is only valid after the call of the FinalizeCalculation method, but this condition is not enforced.

inline void DoNextIteration(double pulseTime, double omZ)

Does the next iteration step. This method is supposed to get called within the I3RecoPulseSeriesMap loop for a particular DOM.

inline void FinalizeCalculation()

Finalizes the ZTravel calculation. This method is supposed to get called after the I3RecoPulseSeriesMap loop.

inline void Reset()

Resets the ZTravelIterativeCalculator, in order to be able to use the same ZTravelIterativeCalculator object more than once.

Public Members

std::vector<ZTravelPulse_t> zTravelPulseVec_
double zTravel_

Public Static Functions

static inline bool compare_ZTravelPulse_pulseTime(ZTravelPulse_t lhs, ZTravelPulse_t rhs)

Compares two ZTravelPulse_t objects after the time of the pulse.

struct ZTravelPulse_t

A structure for sorting the pulses after time.

Public Functions

inline ZTravelPulse_t(double pulseTime_, double omZ_)

Public Members

double pulseTime
double omZ
namespace common_variables

Author

Anna Pollmann anna.pollmann@uni-wuppertal.de

namespace details
namespace direct_hits

Functions

I3DirectHitsDefinitionSeries GetDefaultDefinitions()

Gets an I3DirectHitsDefinitionSeries object with the default direct hits definitions (e.g. time windows) A, B, C, D:

  • time window “A”: [-15ns; + 15ns]

  • time window “B”: [-15ns; + 25ns]

  • time window “C”: [-15ns; + 75ns]

  • time window “D”: [-15ns; +125ns]

I3DirectHitsValuesMapPtr CalculateDirectHits(const I3DirectHitsDefinitionSeries &dhDefinitions, const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap, const I3Particle &particle)

Calculates the direct hits for the given time windows (given through the I3DirectHitsDefinitionSeries object, the given I3Geometry, the given I3RecoPulseSeriesMap and the given I3Particle.

I3RecoPulseSeriesMapPtr GetDirectHitsPulseMap(const I3DirectHitsDefinition dhDefinition, const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap, const I3Particle &particle, bool allPulsesOfDirectDoms)

Returns an I3RecoPulseSeriesMapPtr with the direct hit DOMs for a given direct hits definition, e.g. direct hits time window, for the given I3Geometry, the given I3RecoPulseSeriesMap, and the given I3Particle.

If the allPulsesOfDirectDoms option (default is false) is set to true, the I3RecoPulseSeries object of each direct hit DOM will contain all pulses of the original I3RecoPulseSeries of that DOM. Otherwise only the direct pulses will be included.

I3VectorDoublePtr GetDirectPulsesTimeResiduals(const I3DirectHitsDefinition dhDefinition, const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap, const I3Particle &particle)

Returns an I3VectorDoublePtr containing the time residuals of all direct pulses for a given direct hits definition, e.g. direct hits time window, for the given I3Geometry, the given I3RecoPulseSeriesMap, and the given I3Particle.

I3DirectHitsValuesMapPtr CalculateDirectHits(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap, const I3Particle &particle)

Calculates the direct hits for the default direct hits time windows (defined through the GetDefaultDefinitions() function), the given I3Geometry, the given I3RecoPulseSeriesMap and the given I3Particle.

namespace details

Typedefs

typedef std::vector<DirectHitsClassData_t> DirectHitsClassDataVector_t
namespace hit_multiplicity

Functions

I3HitMultiplicityValuesPtr CalculateHitMultiplicity(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the hit multiplicity, e.g. NString, NChannel, NHit, and NOneHit for the given I3Geometry, and the given I3RecoPulseSeriesMap.

If the “considerOnlyFirstPulse” argument is set to “true” (the default value), only the first pulse in each of the I3RecoPulseSeries will be used for the calculation.

namespace details
namespace hit_statistics

Functions

I3HitStatisticsValuesPtr CalculateHitStatistics(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the hit statistics, e.g. COG, PulseTimespan, QMax, QTot, ZMean, ZSigma, and ZTravel for the given I3Geometry, and the given I3RecoPulseSeriesMap.

I3PositionPtr CalculateCOG(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the COG variable for the given I3Geometry, and the given I3RecoPulseSeriesMap.

I3DoublePtr CalculateCOGZSigma(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the sigma value of the z-component of the COG position for the given I3Geometry, and the given I3RecoPulseSeriesMap.

I3DoublePtr CalculateMinPulseTime(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the minimal time of all pulses, e.g. the time of the first pulse, for the given I3Geometry, and the given I3RecoPulseSeriesMap.

I3DoublePtr CalculateMaxPulseTime(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the maximal time of all pulses, e.g. the time of the last pulse, for the given I3Geometry, and the given I3RecoPulseSeriesMap.

I3DoublePtr CalculateQMaxDoms(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the QMaxDoms variable for the given I3Geometry, and the given I3RecoPulseSeriesMap.

I3DoublePtr CalculateQTotPulses(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the QTot variable for the given I3Geometry, and the given I3RecoPulseSeriesMap.

I3DoublePtr CalculateZMin(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the minimal OM_Z value of all hits.

I3DoublePtr CalculateZMax(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the maximal OM_Z value of all hits.

I3DoublePtr CalculateZMean(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the mean value of OM_Z.

I3DoublePtr CalculateZSigma(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the sigma (RMSD) value of the OM_Z values for the given I3Geometry, and the given I3RecoPulseSeriesMap.

Note: The mean of the OM_Z values will be calculated using the CalculateZMean utility function.

I3DoublePtr CalculateZTravel(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the ZTravel variable for the given I3Geometry, and the given I3RecoPulseSeriesMap.

namespace details

Functions

S1Variables_t CalculateS1Variables(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)
namespace time_characteristics

Functions

I3TimeCharacteristicsValuesPtr CalculateTimeCharacteristics(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap)

Calculates the track characteristics values, e.g. TrackPulseTimespan, HitsTrackLength, EmptyHitsTrackLength, TrackHitsSeperationLength, and TrackHitsDistributionSmoothness, for the given I3Geometry, the given I3RecoPulseSeriesMap, and the given track, given through the I3Particle object, for hits within the given cylinder radius around the track.

namespace details

Functions

bool sortPulsesForTime(I3RecoPulse i, I3RecoPulse j)
bool sortPairPulsesForTime(std::pair<I3RecoPulse, I3Position> i, std::pair<I3RecoPulse, I3Position> j)
bool sortTimes(double i, double j)

This is the core implementation of the calculation of the EmptyHitsTrackLength variable, done by the CalculateTiming function.

NOTE: The vector of hit distances along the track must be sorted by its values from lower to higher values!

Returns NAN if there are less then 2 hits.

namespace track_characteristics

Functions

I3TrackCharacteristicsValuesPtr CalculateTrackCharacteristics(const I3Geometry &geometry, const I3RecoPulseSeriesMap &pulsesMap, const I3Particle &particle, double trackCylinderRadius)

Calculates the track characteristics values, e.g. TrackPulseTimespan, HitsTrackLength, EmptyHitsTrackLength, TrackHitsSeperationLength, and TrackHitsDistributionSmoothness, for the given I3Geometry, the given I3RecoPulseSeriesMap, and the given track, given through the I3Particle object, for hits within the given cylinder radius around the track.

namespace details

Functions

double CalculateEmptyHitsTrackLength(const std::vector<double> &sortedHitDistancesAlongTrack)

This is the core implementation of the calculation of the EmptyHitsTrackLength variable, done by the CalculateTrackCharacteristics function.

NOTE: The vector of hit distances along the track must be sorted by its values from lower to higher values!

Returns NAN if there are less then 2 hits.

double CalculateTrackHitsDistributionSmoothness(const std::vector<double> &sortedHitDistancesAlongTrack)

This is the core implementation of the calculation of the TrackHitsDistributionSmoothness variable, done by the CalculateTrackCharacteristics function.

NOTE: The vector of hit distances along the track must be sorted by its values from lower to higher values!

Returns NAN if there are less then 3 hits.

namespace std

STL namespace.

file ArithMeanIterativeCalculator.h
#include <cmath>
#include <stdint.h>

This file contains the definition and implementation of the iterative calculator for an arithmetic mean value. The ArithMeanIterativeCalculator has been introduced in order to have only one place for the implementation of an arithmetic mean value, and to minimize the number of needed loops over an I3RecoPulseSeriesMap.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file direct_hits/calculator.cxx
#include <cmath>
#include <vector>
#include “phys-services/I3Calculator.h”

This file contains the implementation of the calculator functions that are defined inside the common_variables::direct_hits namespace.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file hit_multiplicity/calculator.cxx

This file contains the implementation of the calculator functions that are defined inside the common_variables::hit_multiplicity namespace.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file hit_statistics/calculator.cxx
#include <cmath>

This file contains the implementation of the calculator functions that are defined inside the common_variables::hit_statistics namespace.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file time_characteristics/calculator.cxx
#include <cmath>
#include <vector>
#include <algorithm>
#include <iostream>
#include “phys-services/I3Calculator.h”
file track_characteristics/calculator.cxx
#include <cmath>
#include <vector>
#include “phys-services/I3Calculator.h”

This file contains the implementation of the utility functions defined in the common_variables::track_characteristics namespace.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file direct_hits/calculator.h
#include “dataclasses/I3Vector.h”
#include “dataclasses/geometry/I3Geometry.h”
#include “dataclasses/physics/I3Particle.h”
#include “dataclasses/physics/I3RecoPulse.h”
#include “recclasses/I3DirectHitsValues.h”

This file contains the definition of the calculator functions that are defined inside the common_variables::direct_hits namespace.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file hit_multiplicity/calculator.h
#include “dataclasses/geometry/I3Geometry.h”
#include “dataclasses/physics/I3RecoPulse.h”
#include “recclasses/I3HitMultiplicityValues.h”

This file contains the definition of the calculator functions that are defined inside the common_variables::hit_multiplicity namespace.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file hit_statistics/calculator.h
#include “dataclasses/I3Double.h”
#include “dataclasses/I3Position.h”
#include “dataclasses/geometry/I3Geometry.h”
#include “dataclasses/physics/I3RecoPulse.h”
#include “recclasses/I3HitStatisticsValues.h”

This file contains the definition of the calculator functions that are defined inside the common_variables::hit_statistics namespace.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file time_characteristics/calculator.h
#include <vector>
#include “dataclasses/geometry/I3Geometry.h”
#include “dataclasses/physics/I3Particle.h”
#include “dataclasses/physics/I3RecoPulse.h”
#include “recclasses/I3TimeCharacteristicsValues.h”
file track_characteristics/calculator.h
#include <vector>
#include “dataclasses/geometry/I3Geometry.h”
#include “dataclasses/physics/I3Particle.h”
#include “dataclasses/physics/I3RecoPulse.h”
#include “recclasses/I3TrackCharacteristicsValues.h”

This file contains the definition of the utility functions defined in the common_variables::track_characteristics namespace.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file COGIterativeCalculator.h
#include <cmath>
#include “dataclasses/I3Position.h”

This file contains the definition and implementation of the iterative calculator for the COG position. The COGIterativeCalculator has been introduced in order to have only one place for the implementation of the COG, and to minimize the number of needed loops over an I3RecoPulseSeriesMap.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file I3DirectHitsDefinition.cxx
#include <iostream>
#include “icetray/I3Units.h”

This file contains the implementation of the I3DirectHitsDefinition class, which is an I3FrameObject holding the values for a particular class of direct hits.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

Functions

std::ostream &operator<<(std::ostream &oss, const I3DirectHitsDefinition &rhs)
file I3DirectHitsDefinition.h
#include <string>
#include <vector>
#include “icetray/I3Logging.h”
#include “icetray/I3PointerTypedefs.h”

This file contains the definition of the I3DirectHitsDefinition class, describing a particular class of direct hits.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

Typedefs

typedef std::vector<I3DirectHitsDefinition> I3DirectHitsDefinitionSeries

Functions

std::ostream &operator<<(std::ostream &oss, const I3DirectHitsDefinition &d)
I3_POINTER_TYPEDEFS(I3DirectHitsDefinition)
I3_POINTER_TYPEDEFS(I3DirectHitsDefinitionSeries)
file I3DirectHitsTimeWindow.cxx
#include <iostream>
#include “icetray/I3Units.h”

This file contains the implementation of the I3DirectHitsTimeWindow class, which is an I3FrameObject holding the values for a particular class of direct hits.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

Functions

std::ostream &operator<<(std::ostream &oss, const I3DirectHitsTimeWindow &rhs)
file I3DirectHitsTimeWindow.h
#include <iostream>
#include <string>
#include <vector>
#include “icetray/I3Logging.h”
#include “icetray/I3PointerTypedefs.h”

This file contains the definition of the I3DirectHitsTimeWindow class, describing a time window of a particular class of direct hits.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

Functions

std::ostream &operator<<(std::ostream &oss, const I3DirectHitsTimeWindow &rhs)
I3_POINTER_TYPEDEFS(I3DirectHitsTimeWindow)
file NHitStringsIterativeCalculator.h
#include <set>

This file contains the definition and implementation of the iterative calculator for the number of hit strings. The NHitStringsIterativeCalculator has been introduced in order to have only one place for the implementation of the counting of hit strings.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file SigmaIterativeCalculator.h
#include <cmath>

This file contains the definition and implementation of the iterative calculator for a sigma value (for a given base value). The SigmaIterativeCalculator has been introduced in order to have only one place for the implementation of a sigma value, and to minimize the number of needed loops over an I3RecoPulseSeriesMap.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file TrackHitsSeparationLengthIterativeCalculator.h
#include <algorithm>
#include <cmath>
#include <vector>
#include “dataclasses/physics/I3Particle.h”
#include “phys-services/I3Calculator.h”

This file contains the definition and implementation of the iterative calculator for the TrackHitsSeparationLength variable. The TrackHitsSeparationLengthIterativeCalculator has been introduced in order to have only one place for the implementation of the TrackHitsSeparationLength variable.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file WeightedMeanIterativeCalculator.h
#include <cmath>

This file contains the definition and implementation of the iterative calculator for a weighted mean value. The WeightedMeanIterativeCalculator has been introduced in order to have only one place for the implementation of a weighted mean value, and to minimize the number of needed loops over an I3RecoPulseSeriesMap.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

file ZTravelIterativeCalculator.h
#include <algorithm>
#include <cmath>

This file contains the definition and implementation of the iterative calculator for the ZTravel cut variable. The ZTravelIterativeCalculator has been introduced in order to have only one compact place for the implementation of the ZTravel variable.

$Id$

Copyright (C) 2012 Martin Wolf martin.wolf@icecube.wisc.edu and the IceCube Collaboration http://www.icecube.wisc.edu

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

Version

$Revision$

Date

$Date$

Author

Martin Wolf martin.wolf@icecube.wisc.edu This file is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/

dir CommonVariables
dir CommonVariables
dir CommonVariables
dir direct_hits
dir direct_hits
dir hit_multiplicity
dir hit_multiplicity
dir hit_statistics
dir hit_statistics
dir icetray
dir private
dir public
dir time_characteristics
dir time_characteristics
dir track_characteristics
dir track_characteristics