MuonGun C++ API Reference¶
-
class BasicSurfaceScalingFunction : public I3MuonGun::SurfaceScalingFunction¶
- #include <EnergyDependentSurfaceInjector.h>
Public Functions
-
BasicSurfaceScalingFunction()¶
-
virtual ~BasicSurfaceScalingFunction()¶
-
virtual SamplingSurfacePtr GetSurface(double energy) const override¶
Propose a target surface for the given energy.
-
virtual bool operator==(const SurfaceScalingFunction&) const override¶
Compare for equality.
-
void SetCapScaling(double energyScale, double scale, double offset, double power)¶
-
void SetSideScaling(double energyScale, double scale, double offset, double power)¶
-
void SetRadiusBounds(double rmin, double rmax)¶
-
void SetZBounds(double zmin, double zmax)¶
-
void SetCenterMin(double x, double y)¶
-
void SetCenterMax(double x, double y)¶
Private Functions
-
double GetMargin(double logenergy, double scale, double offset, double power) const¶
Private Members
Friends
- friend class icecube::serialization::access
-
BasicSurfaceScalingFunction()¶
-
class BMSSEnergyDistribution : public I3MuonGun::EnergyDistribution¶
- #include <EnergyDistribution.h>
Public Functions
-
BMSSEnergyDistribution()¶
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity, double radius, log_value log_energy) const¶
-
virtual std::vector<std::pair<double, double>> Generate(I3RandomService &rng, double depth, double cos_theta, unsigned multiplicity, unsigned samples) const¶
Sample samples (radius, energy) pairs.
-
virtual double GetMaxRadius() const override¶
-
virtual bool operator==(const EnergyDistribution&) const override¶
-
OffsetPowerLaw GetSpectrum(double depth, double cos_theta, unsigned m, double r) const¶
Private Members
-
double beta_¶
-
double g0_¶
-
double g1_¶
-
double e0a_¶
-
double e0b_¶
-
double e1a_¶
-
double e1b_¶
-
double a0_¶
-
double a1_¶
-
double b0a_¶
-
double b0b_¶
-
double b1a_¶
-
double b1b_¶
-
double q0_¶
-
double q1_¶
-
double c0a_¶
-
double c0b_¶
-
double c1_¶
-
double d0a_¶
-
double d0b_¶
-
double d1a_¶
-
double d1b_¶
Friends
- friend class icecube::serialization::access
-
BMSSEnergyDistribution()¶
-
class BMSSFlux : public I3MuonGun::Flux¶
- #include <Flux.h>
Total flux according to Becherini et al.
Public Functions
-
BMSSFlux()¶
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity) const¶
Private Members
-
double k0a_¶
-
double k0b_¶
-
double k1a_¶
-
double k1b_¶
-
double v0a_¶
-
double v0b_¶
-
double v0c_¶
-
double v1a_¶
-
double v1b_¶
Friends
- friend class icecube::serialization::access
-
BMSSFlux()¶
-
class BMSSRadialDistribution : public I3MuonGun::RadialDistribution¶
- #include <RadialDistribution.h>
Radial distribution according to Becherini et al.
Public Functions
-
BMSSRadialDistribution()¶
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity, double radius) const¶
Calculate the logarithm of the probability of obtaining the given radial offset.
See also
operator()
-
virtual double Generate(I3RandomService &rng, double depth, double cos_theta, unsigned multiplicity) const¶
Draw a sample from the distribution of radii.
- Parameters:
rng – [in] a random number generator
depth – [in] vertical depth in km
cos_theta – [in] cosine of zenith angle
multiplicity – [in] number of muons in the bundle
- Returns:
a radius in m
-
virtual bool operator==(const RadialDistribution&) const override¶
Private Functions
-
double GetMeanRadius(double, double, unsigned) const¶
-
double GetShapeParameter(double, double, unsigned) const¶
-
double GetGenerationProbability(double, double, double) const¶
Private Members
-
double rho0a_¶
-
double rho0b_¶
-
double rho1_¶
-
double theta0_¶
-
double f_¶
-
double alpha0a_¶
-
double alpha0b_¶
-
double alpha1a_¶
-
double alpha1b_¶
-
double rmax_¶
Friends
- friend class icecube::serialization::access
-
BMSSRadialDistribution()¶
-
struct BundleEntry¶
- #include <Generator.h>
The radial offset and energy of each muon in a bundle
Public Functions
-
inline BundleEntry(double r = 0., double e = 0.)¶
-
inline bool operator<(const BundleEntry &other) const¶
-
inline bool operator==(const BundleEntry &other) const¶
-
inline BundleEntry(double r = 0., double e = 0.)¶
-
struct BundleModel¶
- #include <WeightCalculator.h>
Convenience struct to hold the components of a muon bundle flux model.
Public Functions
-
inline BundleModel(FluxPtr f, RadialDistributionPtr r, EnergyDistributionPtr e)¶
-
inline BundleModel(FluxPtr f, RadialDistributionPtr r, EnergyDistributionPtr e)¶
-
struct Checkpoint¶
- #include <Track.h>
A point at which the absolute energy of the particle is known.
Public Functions
-
inline Checkpoint(double l, double e = 0., size_t o = 0)¶
-
inline Checkpoint(double l, double e = 0., size_t o = 0)¶
-
class CompactTrack¶
- #include <CompactTrack.h>
A compressed representation of a muon in a bundle.
Public Functions
-
inline CompactTrack()¶
-
CompactTrack(const I3Particle &track)¶
-
inline void SetRadius(double r)¶
-
inline double GetRadius() const¶
-
inline void SetEnergy(double r)¶
-
inline double GetEnergy() const¶
-
inline void SetTime(double r)¶
-
inline double GetTime() const¶
-
inline void SetType(I3Particle::ParticleType t)¶
-
inline I3Particle::ParticleType GetType() const¶
-
bool operator==(const CompactTrack&) const¶
Friends
- friend class icecube::serialization::access
-
inline CompactTrack()¶
-
class ConstantSurfaceScalingFunction : public I3MuonGun::SurfaceScalingFunction¶
- #include <EnergyDependentSurfaceInjector.h>
Public Functions
-
ConstantSurfaceScalingFunction(SamplingSurfacePtr surface)¶
-
virtual ~ConstantSurfaceScalingFunction()¶
-
virtual SamplingSurfacePtr GetSurface(double energy) const override¶
Propose a target surface for the given energy.
-
virtual bool operator==(const SurfaceScalingFunction&) const override¶
Compare for equality.
Private Functions
-
ConstantSurfaceScalingFunction()¶
Private Members
-
SamplingSurfacePtr surface_¶
Friends
- friend class icecube::serialization::access
-
ConstantSurfaceScalingFunction(SamplingSurfacePtr surface)¶
-
class CORSIKAGenerationProbability : public I3MuonGun::GenerationProbability¶
- #include <CORSIKAGenerationProbability.h>
A parametrization of the muon yield from direct air-shower simulation.
Public Functions
-
CORSIKAGenerationProbability(SamplingSurfacePtr, FluxPtr, RadialDistributionPtr, EnergyDistributionPtr)¶
-
inline SamplingSurfaceConstPtr GetSurface() const¶
-
inline FluxConstPtr GetFlux() const¶
-
inline RadialDistributionConstPtr GetRadialDistribution() const¶
-
inline EnergyDistributionConstPtr GetEnergyDistribution() const¶
-
inline virtual SamplingSurfaceConstPtr GetInjectionSurface() const override¶
Get the injection surface this generation scheme uses.
If generation schemes are combined, the injection surfaces must be identical!
-
virtual GenerationProbabilityPtr Clone() const¶
Copy self into a shared pointer
-
virtual bool IsCompatible(GenerationProbabilityConstPtr) const override¶
Compare to another GenerationProbability.
- Returns:
true if the argument is identical to *this to within a scale factor, false otherwise.
Protected Functions
-
virtual double GetLogGenerationProbability(const I3Particle &axis, const BundleConfiguration &bundle) const¶
Calculate the differential probability per event that the given configuration was generated.
For single muons, this is \( \log(dP/dE) [\log(1/GeV)]\), for bundles \( \log(d^2P/dEdr) [\log(1/GeV m)]\)
- Parameters:
axis – [in] the bundle axis
bundle – [in] the radial offset and energy of each muon in the bundle
-
CORSIKAGenerationProbability(SamplingSurfacePtr, FluxPtr, RadialDistributionPtr, EnergyDistributionPtr)¶
-
struct cosine¶
- #include <histogram.h>
Bin edges linear in \( \cos{\theta} \).
-
class Crust¶
- #include <MuonPropagator.h>
A set of nested media layers.
Public Functions
Add an inner layer
-
I3Particle Ingest(const I3Particle &p)¶
Propagate a muon to the outer boundary of the innermost layer
Private Members
-
boost::shared_ptr<MuonPropagator> defaultPropagator_¶
-
std::vector<I3Surfaces::SurfacePtr> boundaries_¶
-
std::vector<boost::shared_ptr<MuonPropagator>> propagators_¶
-
class Cylinder : public I3MuonGun::detail::UprightSurface<I3Surfaces::detail::CylinderBase<SamplingSurface>>¶
- #include <Cylinder.h>
A cylinder aligned with the z axis.
Public Functions
-
inline Cylinder(double length, double radius, I3Position center = I3Position(0, 0, 0))¶
-
bool operator==(const SamplingSurface&) const¶
-
inline virtual double GetLength() const¶
Private Types
-
typedef I3Surfaces::detail::CylinderBase<SamplingSurface> CylinderBase¶
-
typedef detail::UprightSurface<CylinderBase> Base¶
Friends
- friend class icecube::serialization::access
-
inline Cylinder(double length, double radius, I3Position center = I3Position(0, 0, 0))¶
-
class EnergyDependentSurfaceInjector : public I3MuonGun::StaticSurfaceInjector¶
- #include <EnergyDependentSurfaceInjector.h>
A rejection-sampling Generator with energy-dependent sampling surface.
EnergyDependentSurfaceInjector samples bundle impact points, angles, multiplicities, and radial distributions at their natural frequencies, but scales the sampling surface based on the highest-energy muon in the bundle: dim, low-energy muons are targeted only at a small inner surface, while the surface scales up to full size for potentially bright muons. This technique can be used to efficiently simulate background for an event selection that requires a thick veto for dim events (where the rates are also highest) but becomes more accepting for bright events.
Public Functions
-
virtual double GetLogGenerationProbability(const I3Particle &axis, const BundleConfiguration &bundle) const override¶
Calculate the differential probability per event that the given configuration was generated.
For single muons, this is \( \log(dP/dE) [\log(1/GeV)]\), for bundles \( \log(d^2P/dEdr) [\log(1/GeV m)]\)
- Parameters:
axis – [in] the bundle axis
bundle – [in] the radial offset and energy of each muon in the bundle
-
virtual GenerationProbabilityPtr Clone() const override¶
Copy self into a shared pointer
-
virtual bool IsCompatible(GenerationProbabilityConstPtr) const override¶
Compare to another GenerationProbability.
- Returns:
true if the argument is identical to *this to within a scale factor, false otherwise.
-
virtual void Generate(I3RandomService &rng, I3MCTree &tree, BundleConfiguration &bundle) const¶
Generate a muon bundle.
- Parameters:
rng – [in] A random number generator
tree – [out] An I3MCTree to fill the generated bundle in to. The bundle axis should be used as the primary, with its type set to I3Particle::unknown
bundle – [out] the radial offset and energy of each muon in the bundle
-
inline SurfaceScalingFunctionPtr GetScaling() const¶
-
inline void SetScaling(SurfaceScalingFunctionPtr &f)¶
-
SamplingSurfacePtr GetTargetSurface(double energy) const¶
Scale the sampling cylinder to a size appropriate for the given maximum muon energy by the given energy. This is not necessarily the fastest.
-
double GetTotalRate(SamplingSurfaceConstPtr surface) const¶
Integrate the flux to get the total rate on the surface. This is not necessarily the fastest.
Private Members
-
SurfaceScalingFunctionPtr scalingFunction_¶
Friends
- friend class icecube::serialization::access
-
virtual double GetLogGenerationProbability(const I3Particle &axis, const BundleConfiguration &bundle) const override¶
-
class EnergyDistribution¶
- #include <EnergyDistribution.h>
Normalized distribution of energies within a bundle.
The energy distribution \( dP/dE \,\, [GeV^{-1}] \) is parameterized in terms of vertical depth \( [km] \), cosine of the zenith angle, bundle multiplicity, perpendicular distance from the bundle axis \( [m] \), and energy \( [GeV] \).
Single muons are assumed to be colinear with the bundle axis.
Subclassed by I3MuonGun::BMSSEnergyDistribution, I3MuonGun::SplineEnergyDistribution
Public Types
-
typedef double result_type¶
Public Functions
-
inline EnergyDistribution()¶
-
virtual ~EnergyDistribution()¶
-
double operator()(double depth, double cos_theta, unsigned multiplicity, double radius, double energy) const¶
-
double Integrate(double depth, double cos_theta, unsigned multiplicity, double r_min, double r_max, double e_min, double e_max) const¶
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity, double radius, log_value log_energy) const = 0¶
-
virtual std::vector<std::pair<double, double>> Generate(I3RandomService &rng, double depth, double cos_theta, unsigned multiplicity, unsigned samples) const = 0¶
Sample samples (radius, energy) pairs.
-
inline double GetMax() const¶
-
inline double GetMin() const¶
-
inline void SetMax(double v)¶
-
inline void SetMin(double v)¶
-
virtual double GetMaxRadius() const = 0¶
-
virtual bool operator==(const EnergyDistribution&) const = 0¶
- template<typename Archive> void serialize (Archive &ar __attribute__((unused)), unsigned version __attribute__((unused)))
Friends
- friend class icecube::serialization::access
-
typedef double result_type¶
-
template<typename Signature>
class EnsembleSampler¶ - #include <EnsembleSampler.h>
An affine invariant Markov chain Monte Carlo (MCMC) sampler.
Goodman & Weare, Ensemble Samplers With Affine Invariance Comm. App. Math. Comp. Sci., Vol. 5 (2010), No. 1, 65–80
This implementation is a simplified C++ port of emcee’s EnsembleSampler: http://dan.iel.fm/emcee/current/
Public Functions
-
inline EnsembleSampler(boost::function<Signature> log_posterior, const std::vector<array_type> &initial_ensemble)¶
-
inline void Reset()¶
-
inline const std::vector<sample> &Sample(I3RandomService &rng)¶
-
inline double GetAcceptanceRate() const¶
Private Functions
-
inline double LogProbability(const array_type &point)¶
-
inline void ProposeStretch(I3RandomService &rng, unsigned pos, unsigned offset)¶
-
inline double Stretch(I3RandomService &rng) const¶
-
inline EnsembleSampler(boost::function<Signature> log_posterior, const std::vector<array_type> &initial_ensemble)¶
-
class ExtrudedPolygon : public I3MuonGun::detail::UprightSurface<I3Surfaces::ExtrudedPolygonBase<SamplingSurface>>¶
- #include <ExtrudedPolygon.h>
Public Functions
-
virtual ~ExtrudedPolygon()¶
-
inline ExtrudedPolygon(const std::vector<I3Position> &points, double padding = 0.)¶
-
inline virtual bool operator==(const SamplingSurface&) const¶
Protected Functions
-
inline virtual double GetTopArea() const¶
-
inline virtual double GetSideArea() const¶
-
inline virtual double GetLength() const¶
Private Types
-
typedef I3Surfaces::ExtrudedPolygonBase<SamplingSurface> ExtrudedPolygonBase¶
-
typedef detail::UprightSurface<ExtrudedPolygonBase> Base¶
Private Functions
-
inline ExtrudedPolygon()¶
Friends
- friend class icecube::serialization::access
-
virtual ~ExtrudedPolygon()¶
-
class Floodlight : public I3MuonGun::Generator¶
- #include <Floodlight.h>
Public Functions
-
virtual void Generate(I3RandomService &rng, I3MCTree &tree, BundleConfiguration &bundle) const override¶
Generate a muon bundle.
- Parameters:
rng – [in] A random number generator
tree – [out] An I3MCTree to fill the generated bundle in to. The bundle axis should be used as the primary, with its type set to I3Particle::unknown
bundle – [out] the radial offset and energy of each muon in the bundle
-
virtual GenerationProbabilityPtr Clone() const override¶
Copy self into a shared pointer
-
virtual bool IsCompatible(GenerationProbabilityConstPtr) const override¶
Compare to another GenerationProbability.
- Returns:
true if the argument is identical to *this to within a scale factor, false otherwise.
-
virtual double GetLogGenerationProbability(const I3Particle &axis, const BundleConfiguration &bundle) const override¶
Calculate the differential probability per event that the given configuration was generated.
For single muons, this is \( \log(dP/dE) [\log(1/GeV)]\), for bundles \( \log(d^2P/dEdr) [\log(1/GeV m)]\)
- Parameters:
axis – [in] the bundle axis
bundle – [in] the radial offset and energy of each muon in the bundle
-
inline virtual SamplingSurfaceConstPtr GetInjectionSurface() const override¶
Get the injection surface this generation scheme uses.
If generation schemes are combined, the injection surfaces must be identical!
Private Functions
-
inline Floodlight()¶
Private Members
-
SamplingSurfacePtr surface_¶
-
boost::shared_ptr<OffsetPowerLaw> energyGenerator_¶
-
double log_acceptance_¶
Friends
- friend class icecube::serialization::access
-
virtual void Generate(I3RandomService &rng, I3MCTree &tree, BundleConfiguration &bundle) const override¶
-
class Flux¶
- #include <Flux.h>
A parameterization of the total muon-bundle flux.
The total flux (in units of \( [ m^{-2} sr^{-2} s^{-1} ]\)) is parameterized in terms of vertical depth \( km \), the cosine of the zenith angle, and bundle multiplicity.
Subclassed by I3MuonGun::BMSSFlux, I3MuonGun::SplineFlux
Public Types
-
typedef double result_type¶
Public Functions
-
Flux()¶
-
virtual ~Flux()¶
-
double operator()(double depth, double cos_theta, unsigned multiplicity) const¶
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity) const = 0¶
-
inline unsigned GetMaxMultiplicity() const¶
-
inline unsigned GetMinMultiplicity() const¶
-
inline void SetMaxMultiplicity(unsigned m)¶
-
inline void SetMinMultiplicity(unsigned m)¶
Friends
- friend class icecube::serialization::access
-
typedef double result_type¶
-
class general : public I3MuonGun::histogram::binning::scheme¶
- #include <histogram.h>
A non-equispaced binning scheme.
In the general case, the proper bin can be found in logarithmic time by binary search
-
class GenerationProbability : public I3FrameObject¶
- #include <Generator.h>
A muon bundle generation scheme.
GenerationProbability represents the normalization required for WeightCalculator
Subclassed by I3MuonGun::CORSIKAGenerationProbability, I3MuonGun::GenerationProbabilityCollection, I3MuonGun::Generator
Public Functions
-
inline GenerationProbability()¶
-
virtual ~GenerationProbability()¶
-
inline void SetTotalEvents(double n)¶
-
inline double GetTotalEvents() const¶
-
double GetLogGeneratedEvents(const I3Particle &axis, const BundleConfiguration &bundle) const¶
Calculate the logarithm of the differential number of events that should have been generated by the represented scheme.
- Parameters:
axis – [in] the bundle axis
bundle – [in] the radial offset and energy of each muon in the bundle
- Returns:
The logarithm of a number whose units depend on the number of muons in the bundle. For single muons, this is \( dN/dA d\Omega dP/dE \) and has units of \( [ 1/GeV m^2 sr ] \). For multi-muon bundles it contains an additional doubly differential probability in radius and energy for each muon in the bundle.
-
double GetGeneratedEvents(const I3Particle &axis, const BundleConfiguration &bundle) const¶
Call GetLogGeneratedEvents and exponential the result.
-
virtual SamplingSurfaceConstPtr GetInjectionSurface() const = 0¶
Get the injection surface this generation scheme uses.
If generation schemes are combined, the injection surfaces must be identical!
-
virtual GenerationProbabilityPtr Clone() const = 0¶
Copy self into a shared pointer
-
virtual bool IsCompatible(GenerationProbabilityConstPtr) const = 0¶
Compare to another GenerationProbability.
- Returns:
true if the argument is identical to *this to within a scale factor, false otherwise.
Protected Functions
-
virtual double GetLogGenerationProbability(const I3Particle &axis, const BundleConfiguration &bundle) const = 0¶
Calculate the differential probability per event that the given configuration was generated.
For single muons, this is \( \log(dP/dE) [\log(1/GeV)]\), for bundles \( \log(d^2P/dEdr) [\log(1/GeV m)]\)
- Parameters:
axis – [in] the bundle axis
bundle – [in] the radial offset and energy of each muon in the bundle
Private Members
-
double numEvents_¶
The total number of events that should be generated.
Friends
- friend class icecube::serialization::access
-
inline GenerationProbability()¶
-
class GenerationProbabilityCollection : public I3MuonGun::GenerationProbability, public std::vector<GenerationProbabilityPtr>¶
- #include <Generator.h>
A collection of muon bundle generation schemes.
Public Functions
-
GenerationProbabilityCollection(GenerationProbabilityPtr, GenerationProbabilityPtr)¶
-
void push_back(const GenerationProbabilityPtr&)¶
-
virtual GenerationProbabilityPtr Clone() const override¶
Copy self into a shared pointer
-
virtual SamplingSurfaceConstPtr GetInjectionSurface() const override¶
Get the injection surface this generation scheme uses.
If generation schemes are combined, the injection surfaces must be identical!
-
virtual bool IsCompatible(GenerationProbabilityConstPtr) const override¶
Compare to another GenerationProbability.
- Returns:
true if the argument is identical to *this to within a scale factor, false otherwise.
Protected Functions
-
virtual double GetLogGenerationProbability(const I3Particle &axis, const BundleConfiguration &bundle) const override¶
Calculate the total probability that the given configuration was generated by any of the distributions in the colleciton.
Private Functions
-
inline GenerationProbabilityCollection()¶
Friends
- friend class icecube::serialization::access
-
GenerationProbabilityCollection(GenerationProbabilityPtr, GenerationProbabilityPtr)¶
-
class Generator : public I3MuonGun::GenerationProbability¶
- #include <Generator.h>
A muon bundle generator.
Generators know both how to draw muon bundles from some distribution and how to calculate the probability that they would have drawn some arbitrary bundle from that same distribution.
Subclassed by I3MuonGun::Floodlight, I3MuonGun::NaturalRateInjector, I3MuonGun::StaticSurfaceInjector
Public Functions
-
virtual ~Generator()¶
-
virtual void Generate(I3RandomService &rng, I3MCTree &tree, BundleConfiguration &bundle) const = 0¶
Generate a muon bundle.
- Parameters:
rng – [in] A random number generator
tree – [out] An I3MCTree to fill the generated bundle in to. The bundle axis should be used as the primary, with its type set to I3Particle::unknown
bundle – [out] the radial offset and energy of each muon in the bundle
Public Static Functions
-
static I3Particle CreateParallelTrack(double radius, double azimuth, const I3Surfaces::Surface &surface, const I3Particle &axis)¶
Place a muon at a given radial offset and rotation with respect to the shower axis.
- Parameters:
radius – [in] perpendicular distance from the shower axis
azimuth – [in] rotation about the shower axis
surface – [in] place the muon on this surface, with timing adjusted so that it remains in the shower plane
axis – [in] the shower axis
Friends
- friend class icecube::serialization::access
-
virtual ~Generator()¶
-
template<size_t N, typename T = double>
class histogram : public I3MuonGun::histogram::histogram_base¶ - #include <histogram.h>
An N-dimensional histogram.
- Template Parameters:
N – number of dimensions
T – type weight stored in the histogram
Public Types
Public Functions
-
inline histogram(const boost::array<std::vector<double>, N> &edges)¶
Construct with non-uniform bins in all dimensions
Construct with uniform bins in all dimensions
Construct with a mix of uniform and non-uniform bins in different dimensions
-
inline void fill(const boost::array<double, N> &values, T weight = 1)¶
Add a sample to the histogram
- Parameters:
values – [in] value of the sample in each dimension
weight – [in] weight to assign to the sample
-
inline virtual const char *raw_bincontent() const¶
Get a pointer to the memory backing the bin contents
-
inline virtual const char *raw_squaredweights() const¶
Get a pointer to the memory backing the sum of squared weights
-
inline virtual size_type ndim() const¶
Get the dimensionality of the histogram
Private Functions
-
inline void make_datacube()¶
-
class histogram_base¶
- #include <histogram.h>
A base class for histograms of arbitrary type and dimensionality.
Subclassed by I3MuonGun::histogram::histogram< N, T >
Public Types
-
typedef boost::multi_array_types::index index¶
-
typedef boost::multi_array_types::size_type size_type¶
-
typedef boost::multi_array_types::index index¶
-
struct identity¶
- #include <histogram.h>
Trivial linear mapping.
-
struct log10¶
- #include <histogram.h>
Bin edges linear in \( \log_{10}{x} \).
-
struct log_value¶
- #include <EnergyDistribution.h>
A tag type to enforce units.
Public Members
-
double value_¶
-
double value_¶
-
struct LossSum¶
- #include <Track.h>
The sum of stochastic energy losses since the last checkpoint.
Public Functions
-
inline LossSum(double l, double e = 0.)¶
-
inline LossSum(double l, double e = 0.)¶
-
class MuonBundleConverter : public I3ConverterImplementation<I3MCTree>¶
- #include <WeightCalculator.h>
Public Functions
-
MuonBundleConverter(size_t maxMultiplicity = 25, SamplingSurfaceConstPtr surface = SamplingSurfaceConstPtr())¶
Private Functions
- SET_LOGGER ("I3MuonGun::WeightCalculator")
-
MuonBundleConverter(size_t maxMultiplicity = 25, SamplingSurfaceConstPtr surface = SamplingSurfaceConstPtr())¶
-
class MuonGunBackgroundService : public I3GeneratorService¶
- #include <MuonGunBackgroundService.h>
Service Generator for MuonGun.
-
class Muonitron : public I3Module¶
- #include <Muonitron.h>
a utility module for recording the state of a muon bundle as it propagates to depth
Public Static Functions
-
static double GetOverburden(double z, double d = OriginDepth, double r = SurfaceRadius)¶
-
static double GetSurfaceZenith(double z, double d = OriginDepth, double r = SurfaceRadius)¶
-
static double GetGeocentricZenith(double z, double d = OriginDepth, double r = SurfaceRadius)¶
-
static double GetDetectorZenith(double z, double d = OriginDepth, double r = SurfaceRadius)¶
-
static I3Direction RotateToZenith(const I3Direction&, const I3Direction&)¶
-
static I3Position RotateToZenith(const I3Direction&, const I3Position&)¶
-
static I3Particle RotateToZenith(const I3Particle&, const I3Particle&)¶
-
static I3Position Impact(const I3Particle&)¶
Private Functions
-
bool PropagateTrack(I3Particle &target, double depth)¶
-
static double GetOverburden(double z, double d = OriginDepth, double r = SurfaceRadius)¶
-
class MuonPropagator¶
- #include <MuonPropagator.h>
A simple muon energy-loss calculator.
This hides the nasty details of PROPOSAL (a C++ translation of MMC)
Public Functions
-
MuonPropagator(const std::string &medium, double ecut = -1, double vcut = -1, double rho = 1.0)¶
- Parameters:
medium – [in] The name of the medium, e.g. “ice”
ecut – [in] Absolute energy above which an energy loss is considered stochastic \( [MeV] \)
vcut – [in] Proportion of the current muon energy above which an energy loss is considered stochastic
rho – [in] Density adjustment factor for the medium
-
~MuonPropagator()¶
- Parameters:
p – [in] Muon to propagate
distance – [in] Maximum distance to propagate
- Returns:
an I3Particle representing the muon at the end of propagation. If the muon stopped before the given distance, the energy will be set to zero. The length of the output I3Particle is the distance traveled to reach its current position from the position given as input.
-
double GetStochasticRate(double energy, double fraction, I3Particle::ParticleType type = I3Particle::MuMinus) const¶
Differential stochastic rate: d^2N/dv/dx [1/m]
-
double GetTotalStochasticRate(double energy, I3Particle::ParticleType type = I3Particle::MuMinus) const¶
total stochastic rate: dN/dx [1/m]
Public Static Functions
-
static void SetSeed(int seed)¶
Set the (global) state of the random number generator used in the implementation.
-
static std::string GetName(const I3Particle &p)¶
Get the internal MMC name associated with a particle type
Private Functions
-
I3Particle::ParticleType GenerateI3Type(const PROPOSAL::DynamicData &secondary)¶
Get the internal MMC name associated with a particle type
-
I3Particle GenerateI3Particle(const PROPOSAL::DynamicData &proposal_particle)¶
Private Members
-
PROPOSAL::Propagator *propagator_¶
-
MuonPropagator(const std::string &medium, double ecut = -1, double vcut = -1, double rho = 1.0)¶
-
class NaturalRateInjector : public I3MuonGun::Generator¶
- #include <NaturalRateInjector.h>
A natural-rate Generator.
NaturalRateInjector samples bundle impact points, angles, multiplicities, and radius/energy distributions at their natural frequencies on a fixed surface using a brain-dead acceptance/rejection technique. This is to MUPAGE, except that the flux parameterizations can be swapped out and the events can be combined with biased generation and reweighted to a different parameterization after the fact.
Public Functions
-
NaturalRateInjector()¶
-
NaturalRateInjector(SamplingSurfacePtr surface, FluxPtr flux, EnergyDistributionPtr edist)¶
-
virtual void Generate(I3RandomService &rng, I3MCTree &tree, BundleConfiguration &bundle) const¶
Generate a muon bundle.
- Parameters:
rng – [in] A random number generator
tree – [out] An I3MCTree to fill the generated bundle in to. The bundle axis should be used as the primary, with its type set to I3Particle::unknown
bundle – [out] the radial offset and energy of each muon in the bundle
-
virtual GenerationProbabilityPtr Clone() const¶
Copy self into a shared pointer
-
virtual bool IsCompatible(GenerationProbabilityConstPtr) const¶
Compare to another GenerationProbability.
- Returns:
true if the argument is identical to *this to within a scale factor, false otherwise.
-
virtual double GetLogGenerationProbability(const I3Particle &axis, const BundleConfiguration &bundle) const¶
Calculate the differential probability per event that the given configuration was generated.
For single muons, this is \( \log(dP/dE) [\log(1/GeV)]\), for bundles \( \log(d^2P/dEdr) [\log(1/GeV m)]\)
- Parameters:
axis – [in] the bundle axis
bundle – [in] the radial offset and energy of each muon in the bundle
-
inline virtual SamplingSurfaceConstPtr GetInjectionSurface() const¶
Get the injection surface this generation scheme uses.
If generation schemes are combined, the injection surfaces must be identical!
-
void SetSurface(SamplingSurfacePtr p)¶
-
inline SamplingSurfacePtr GetSurface()¶
-
void SetFlux(FluxPtr p)¶
-
inline FluxPtr GetFlux()¶
-
void SetEnergyDistribution(EnergyDistributionPtr e)¶
-
inline EnergyDistributionPtr GetEnergyDistribution()¶
-
double GetTotalRate() const¶
Integrate the configured flux over the sampling surface, summing over all allowed multiplicities.
- Returns:
a rate in units of \( [s^{-1}] \)
Protected Functions
-
void GenerateAxis(I3RandomService &rng, std::pair<I3Particle, unsigned> &axis) const¶
Draw a sample from the distribution of shower impact points
The shower axis and multiplcity are filled into the axis
-
void FillMCTree(I3RandomService &rng, const std::pair<I3Particle, unsigned> &axis, I3MCTree&, BundleConfiguration&) const¶
Distribute the given number of muons in the transverse plane and draw an energy for each
Protected Attributes
-
SamplingSurfacePtr surface_¶
-
FluxPtr flux_¶
-
EnergyDistributionPtr energyDistribution_¶
-
mutable double totalRate_¶
Friends
- friend class icecube::serialization::access
-
NaturalRateInjector()¶
-
class NeutrinoBinner¶
- #include <TrackBinner.h>
Public Types
-
class OffsetPowerLaw¶
- #include <EnergyDistribution.h>
An approximate form for the underground muon energy spectrum.
The deep underground muon spectrum very nearly follows a power law with a break at low energies caused by pile-up of minimum-ionizing muons. While the details of the spectrum depend on depth, zenith angle, distance from the bundle axis, etc., this approximate form can still be used to efficiently generate plausible energies that can later be weighted to reflect their frequency under a more precise model.
Public Types
-
typedef double result_type¶
Public Functions
-
OffsetPowerLaw()¶
-
OffsetPowerLaw(double gamma, double offset, double emin, double emax)¶
Create an offset power-law energy distribution of the form \( dP/dE_{\mu} \propto (E_{\mu} + b)^{-\gamma} \)
- Parameters:
gamma – [in] The power law index (positive)
offset – [in] The offset in \( (E_{\mu} + b) \)
emin – [in] minimum generated energy
emax – [in] maximum generated energy
-
double operator()(double energy) const¶
Calculate the probability that the given energy was generated
-
double GetLog(double energy) const¶
-
double GetLog(EnergyDistribution::log_value log_energy) const¶
-
double Generate(I3RandomService &rng) const¶
Draw an energy from the distribution
-
double InverseSurvivalFunction(double p) const¶
-
inline const double GetMin() const¶
-
inline const double GetMax() const¶
-
bool operator==(const OffsetPowerLaw &other) const¶
Private Members
-
double gamma_¶
-
double offset_¶
-
double emin_¶
-
double emax_¶
-
double nmin_¶
-
double nmax_¶
-
double norm_¶
-
double lognorm_¶
Friends
- friend class icecube::serialization::access
-
typedef double result_type¶
-
template<int N>
struct power¶ - #include <histogram.h>
Bin edges linear in \( x^N \).
- Template Parameters:
N – the exponent of the power law
-
class RadialDistribution¶
- #include <RadialDistribution.h>
The distribution of distance between a muon and the bundle axis.
Subclassed by I3MuonGun::BMSSRadialDistribution, I3MuonGun::SplineRadialDistribution
Public Functions
-
virtual ~RadialDistribution()¶
-
double operator()(double depth, double cos_theta, unsigned multiplicity, double radius) const¶
Calculate the probability of obtaining the given radial offset.
- Parameters:
depth – [in] vertical depth in km
cos_theta – [in] cosine of zenith angle
multiplicity – [in] number of muons in the bundle
radius – [in] distance to bundle axis
- Returns:
a properly normalized propability \( dP/dr \,\, [m^{-1}] \)
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity, double radius) const = 0¶
Calculate the logarithm of the probability of obtaining the given radial offset.
See also
operator()
-
virtual double Generate(I3RandomService &rng, double depth, double cos_theta, unsigned multiplicity) const = 0¶
Draw a sample from the distribution of radii.
- Parameters:
rng – [in] a random number generator
depth – [in] vertical depth in km
cos_theta – [in] cosine of zenith angle
multiplicity – [in] number of muons in the bundle
- Returns:
a radius in m
-
virtual bool operator==(const RadialDistribution&) const = 0¶
- template<typename Archive> void serialize (Archive &ar __attribute__((unused)), unsigned __attribute__((unused)))
Friends
- friend class icecube::serialization::access
-
virtual ~RadialDistribution()¶
-
struct sample¶
- #include <EnsembleSampler.h>
Public Functions
-
inline sample(const array_type &p, double logprob)¶
-
inline sample(const array_type &p, double logprob)¶
-
class SamplingSurface : public SamplingSurface¶
- #include <SamplingSurface.h>
A surface upon which muon bundles can be generated.
SamplingSurface knows how to calculate its projected area and integrate a flux over its surface. It is assumed to be azimuthally symmetric, but its projected area may vary with zenith angle.
Public Functions
-
virtual ~SamplingSurface()¶
-
virtual double GetAcceptance(double cosMin = 0, double cosMax = 1) const = 0¶
Get the integral of area*solid_angle over the given cos(theta) range
-
virtual double GetMinDepth() const = 0¶
Get the minimum vertical depth the surface occupies
-
virtual double IntegrateFlux(boost::function<double(double, double)> flux, double cosMin = 0, double cosMax = 1) const = 0¶
Integrate a flux (defined in terms of dN/dOmega(depth [km], cos(theta))) over the outer surface
-
virtual bool operator==(const SamplingSurface&) const = 0¶
Friends
- friend class icecube::serialization::access
-
virtual ~SamplingSurface()¶
-
class scheme¶
- #include <histogram.h>
Base class for binning schemes.
Subclassed by I3MuonGun::histogram::binning::general, I3MuonGun::histogram::binning::uniform< Transformation >
-
class SplineEnergyDistribution : public I3MuonGun::EnergyDistribution¶
- #include <EnergyDistribution.h>
Energy distribution fit to a tensor-product B-spline surface.
The surface is fit to \( d \log{P} / d\log(E) \), which is nearly polynomial and thus easier to represent with low-order splines.
Public Functions
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity, double radius, log_value log_energy) const¶
-
virtual std::vector<std::pair<double, double>> Generate(I3RandomService &rng, double depth, double cos_theta, unsigned multiplicity, unsigned samples) const¶
Sample samples (radius, energy) pairs.
-
virtual double GetMaxRadius() const override¶
-
virtual bool operator==(const EnergyDistribution&) const override¶
Private Functions
-
inline SplineEnergyDistribution()¶
Friends
- friend class icecube::serialization::access
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity, double radius, log_value log_energy) const¶
-
class SplineFlux : public I3MuonGun::Flux¶
- #include <Flux.h>
Total flux, fit to a tensor-product B-spline surface.
Public Functions
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity) const¶
Private Functions
-
inline SplineFlux()¶
Friends
- friend class icecube::serialization::access
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity) const¶
-
class SplineRadialDistribution : public I3MuonGun::RadialDistribution¶
- #include <RadialDistribution.h>
Radial distribution fit to a tensor-product B-spline surface.
The surface is fit to \( d \log{P} / d{r^2} \) to remove the factor of differential area implicit in \( dP / dr \)
Public Functions
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity, double radius) const¶
Calculate the logarithm of the probability of obtaining the given radial offset.
See also
operator()
-
virtual double Generate(I3RandomService &rng, double depth, double cos_theta, unsigned multiplicity) const¶
Draw a sample from the distribution of radii.
- Parameters:
rng – [in] a random number generator
depth – [in] vertical depth in km
cos_theta – [in] cosine of zenith angle
multiplicity – [in] number of muons in the bundle
- Returns:
a radius in m
-
virtual bool operator==(const RadialDistribution&) const override¶
Private Functions
-
inline SplineRadialDistribution()¶
Private Members
-
SplineTable spline_¶
Friends
- friend class icecube::serialization::access
-
virtual double GetLog(double depth, double cos_theta, unsigned multiplicity, double radius) const¶
-
class SplineTable : public I3FrameObject¶
- #include <SplineTable.h>
An encapsulated, serializeable interface to spline tables.
Public Functions
-
SplineTable(const std::string &path)¶
Read a spline table from a FITS file on disk.
- Parameters:
path – [in] File system path to the FITS file
- Throws:
std::runtime_error – If the file does not exist or is corrupt
-
inline SplineTable()¶
Default constructor; to be used only in serialization.
-
SplineTable(const SplineTable&) = default¶
-
SplineTable(SplineTable&&) = default¶
-
virtual ~SplineTable() = default¶
-
SplineTable &operator=(const SplineTable&) = default¶
-
SplineTable &operator=(SplineTable&&) = default¶
-
bool operator==(const SplineTable&) const¶
Deep comparison.
-
inline bool operator!=(const SplineTable &other) const¶
-
int Eval(const std::vector<double> &coordinates, double &result) const¶
Evaluate the spline surface at the given coordinates.
- Parameters:
coordinates – [in] Coordinates at which to evaluate
result – [out] Where to store the result
- Returns:
0 on success
-
inline unsigned GetNDim() const¶
- Returns:
Number of dimensions of the spline surface
-
std::pair<double, double> GetExtents(uint32_t dim) const¶
Table extents.
- Parameters:
dim – [in] The dimension to query
- Throws:
std::out_of_range – If
dim
is smaller than zero or larger or equal than the number of dimensions- Returns:
The smallest and largest coordinates along the given axis where the spline surface has full support
Private Types
-
using splinetable_t = photospline::splinetable<>¶
Private Functions
-
I3_SERIALIZATION_SPLIT_MEMBER()¶
Friends
- friend class icecube::serialization::access
-
SplineTable(const std::string &path)¶
-
class StaticSurfaceInjector : public I3MuonGun::Generator¶
- #include <StaticSurfaceInjector.h>
A simple rejection-sampling Generator.
StaticSurfaceInjector samples bundle impact points, angles, multiplicities, and radial distributions at their natural frequencies on a fixed surface using a brain-dead acceptance/rejection technique. Energies, on the other hand, are sampled from an OffsetPowerLaw, which both boosts efficiency and allows a measure of control over the generated energy distribution.
Subclassed by I3MuonGun::EnergyDependentSurfaceInjector
Public Functions
-
StaticSurfaceInjector()¶
-
virtual void Generate(I3RandomService &rng, I3MCTree &tree, BundleConfiguration &bundle) const override¶
Generate a muon bundle.
- Parameters:
rng – [in] A random number generator
tree – [out] An I3MCTree to fill the generated bundle in to. The bundle axis should be used as the primary, with its type set to I3Particle::unknown
bundle – [out] the radial offset and energy of each muon in the bundle
-
virtual GenerationProbabilityPtr Clone() const override¶
Copy self into a shared pointer
-
virtual bool IsCompatible(GenerationProbabilityConstPtr) const override¶
Compare to another GenerationProbability.
- Returns:
true if the argument is identical to *this to within a scale factor, false otherwise.
-
virtual double GetLogGenerationProbability(const I3Particle &axis, const BundleConfiguration &bundle) const override¶
Calculate the differential probability per event that the given configuration was generated.
For single muons, this is \( \log(dP/dE) [\log(1/GeV)]\), for bundles \( \log(d^2P/dEdr) [\log(1/GeV m)]\)
- Parameters:
axis – [in] the bundle axis
bundle – [in] the radial offset and energy of each muon in the bundle
-
inline virtual SamplingSurfaceConstPtr GetInjectionSurface() const override¶
Get the injection surface this generation scheme uses.
If generation schemes are combined, the injection surfaces must be identical!
-
void SetSurface(SamplingSurfacePtr p)¶
-
inline SamplingSurfacePtr GetSurface()¶
-
void SetFlux(FluxPtr p)¶
-
inline FluxPtr GetFlux()¶
-
inline void SetRadialDistribution(RadialDistributionPtr r)¶
-
inline RadialDistributionPtr GetRadialDistribution()¶
-
inline boost::shared_ptr<OffsetPowerLaw> GetEnergyDistribution()¶
-
double GetTotalRate() const¶
Integrate the configured flux over the sampling surface, summing over all allowed multiplicities.
- Returns:
a rate in units of \( [s^{-1}] \)
Protected Functions
-
void GenerateAxis(I3RandomService &rng, std::pair<I3Particle, unsigned> &axis) const¶
Draw a sample from the distribution of shower impact points
The shower axis and multiplcity are filled into the axis
-
void FillMCTree(I3RandomService &rng, const std::pair<I3Particle, unsigned> &axis, I3MCTree&, BundleConfiguration&) const¶
Distribute the given number of muons in the transverse plane and draw an energy for each
-
void CalculateMaxFlux()¶
-
double GetZenithNorm() const¶
Get the normalization term for relative weighting of zenith angles and multiplicities by integrating the flux at the top of the cylinder over all zenith angles and summing over all allowed multiplicities.
Protected Attributes
-
SamplingSurfacePtr surface_¶
-
FluxPtr flux_¶
-
boost::shared_ptr<OffsetPowerLaw> energyGenerator_¶
-
RadialDistributionPtr radialDistribution_¶
-
double maxFlux_¶
-
mutable double totalRate_¶
-
mutable double zenithNorm_¶
Friends
- friend class icecube::serialization::access
-
StaticSurfaceInjector()¶
-
class SurfaceScalingFunction¶
- #include <EnergyDependentSurfaceInjector.h>
Subclassed by I3MuonGun::BasicSurfaceScalingFunction, I3MuonGun::ConstantSurfaceScalingFunction
Public Functions
-
virtual ~SurfaceScalingFunction()¶
-
virtual SamplingSurfacePtr GetSurface(double energy) const = 0¶
Propose a target surface for the given energy.
-
virtual bool operator==(const SurfaceScalingFunction&) const = 0¶
Compare for equality.
- template<typename Archive> void serialize (Archive &ar __attribute__((unused)), unsigned version __attribute__((unused)))
Friends
- friend class icecube::serialization::access
-
virtual ~SurfaceScalingFunction()¶
-
class Track : public I3Particle¶
- #include <Track.h>
A particle that moves in a straight line, and loses energy as it does so.
Public Functions
-
inline Track()¶
-
Track(const I3MMCTrack &mmctrack, const I3MCTree::sibling_const_iterator &secondaries_begin, const I3MCTree::sibling_const_iterator &secondaries_end)¶
-
double GetEnergy(double length) const¶
Get the energy of the particle at the given down-track distance, assuming that the “continuous” contribution to the energy loss is constant between measurement points.
- Parameters:
length – [in] distance from the track origin
- Returns:
an energy if 0 <= length < range; otherwise 0
-
I3Position GetPos(double length) const¶
-
double GetTime(double length) const¶
-
inline std::vector<Checkpoint> GetCheckpoints() const¶
Public Static Functions
-
inline Track()¶
-
class TrackBinner¶
- #include <TrackBinner.h>
A utility class for filling muon bundles into histograms.
-
struct TrackBundle : public I3FrameObject, public std::map<double, std::vector<CompactTrack>>¶
- #include <CompactTrack.h>
The state of a muon bundle at a set of vertical depths.
Public Functions
-
virtual ~TrackBundle()¶
Friends
- friend class icecube::serialization::access
-
virtual ~TrackBundle()¶
-
template<typename Transformation = identity>
class uniform : public I3MuonGun::histogram::binning::scheme¶ - #include <histogram.h>
An equispaced binning scheme.
In this optimal case the bin edges are uniform under some transformation between set limits and the bin index can be found in constant time
- Template Parameters:
Transformation – the transformation that makes the bin edges equispaced
Public Functions
-
inline uniform(double low, double high, size_t nsteps)¶
-
inline virtual ~uniform()¶
-
inline virtual long index(double value) const¶
Get the bin index of the given value
Public Static Functions
-
static inline boost::shared_ptr<uniform<Transformation>> create(double low, double high, size_t nsteps)¶
-
template<typename Base>
class UprightSurface : public Base¶ - #include <UprightSurface.h>
A surface consisting only of vertical and horizontal faces.
Public Functions
-
inline double GetMinDepth() const¶
-
inline double IntegrateFlux(boost::function<double(double, double)> flux, double cosMin = 0, double cosMax = 1) const¶
Protected Functions
-
virtual double GetSideArea() const = 0¶
-
virtual double GetTopArea() const = 0¶
-
virtual double GetLength() const = 0¶
-
inline double GetDifferentialTopArea(double coszen) const¶
-
inline double GetDifferentialSideArea(double coszen) const¶
-
inline UprightSurface()¶
Friends
- friend class icecube::serialization::access
-
inline double GetMinDepth() const¶
-
class WeightCalculator¶
- #include <WeightCalculator.h>
Utility class to calculate weights for muon bundles.
Subclassed by I3MuonGun::WeightCalculatorModule
Public Functions
-
inline WeightCalculator(const BundleModel &m, GenerationProbabilityPtr g)¶
- Parameters:
m – [in] Target model of the muon flux as a function of angle, depth, multiplicity, radius, and energy
g – [in] Generation scheme according by which the events were generated.
-
double GetWeight(const I3Particle &axis, const BundleConfiguration &bundle) const¶
Calculate a weight corresponding to the frequency with which the given bundle configuration appears in the flux model.
- Parameters:
axis – [in] The bundle axis
bundle – [in] The radial offset and energy of each muon in the bundle
- Returns:
a weight in units of \( [s^{-1}] \)
-
inline SamplingSurfaceConstPtr GetSurface()¶
-
inline void SetSurface(SamplingSurfacePtr s)¶
Protected Functions
-
inline WeightCalculator()¶
-
inline WeightCalculator(const BundleModel &m, GenerationProbabilityPtr g)¶
-
class WeightCalculatorModule : public I3Module, protected I3MuonGun::WeightCalculator¶
Interface between WeightCalculator and IceTray.
WeightCalculatorModule handles the details of extracting energies and radial offsets of muons from an I3MCTree and MMCTrackList.
-
namespace [anonymous]¶
-
namespace I3MuonGun¶
Typedefs
-
typedef boost::bimap<I3Particle::ParticleType, std::string> bimap_ParticleType¶
-
typedef std::list<BundleEntry> BundleConfiguration¶
Functions
-
I3_POINTER_TYPEDEFS(TrackBundle)¶
-
I3_FORWARD_DECLARATION(SamplingSurface)¶
-
I3_FORWARD_DECLARATION(RadialDistribution)¶
-
I3_FORWARD_DECLARATION(EnergyDistribution)¶
-
I3_FORWARD_DECLARATION(OffsetPowerLaw)¶
-
I3_POINTER_TYPEDEFS(SurfaceScalingFunction)¶
-
I3_POINTER_TYPEDEFS(ConstantSurfaceScalingFunction)¶
-
I3_POINTER_TYPEDEFS(BasicSurfaceScalingFunction)¶
-
template<typename T, size_t N>
std::ostream &operator<<(std::ostream &o, const boost::array<T, N> &vals)¶
-
GenerationProbabilityPtr operator*=(GenerationProbabilityPtr p, double n)¶
-
GenerationProbabilityPtr operator*(GenerationProbabilityPtr p, double n)¶
-
GenerationProbabilityPtr operator*(double n, GenerationProbabilityPtr p)¶
Scale the distribution by the given number of events
-
GenerationProbabilityPtr operator+(GenerationProbabilityPtr p1, GenerationProbabilityPtr p2)¶
Combine the distributions to form a GenerationProbabilityCollection (or append to it if any of the arguments is already one)
-
double GetDepth(double z)¶
Convert an IceCube z-coordinate [m] to a vertical depth [km]
-
double Integrate(boost::function<double(double)> f, double low, double high, double epsabs, double epsrel, size_t limit)¶
-
inline std::string GetMMCName(I3Particle::ParticleType pt)¶
-
std::vector<I3Particle> GetMuonsAtSurface(I3FramePtr frame, I3Surfaces::SurfaceConstPtr surface)¶
-
I3_POINTER_TYPEDEFS(EnergyDistribution)¶
-
I3_POINTER_TYPEDEFS(ExtrudedPolygon)¶
-
I3_FORWARD_DECLARATION(GenerationProbability)¶
-
I3_POINTER_TYPEDEFS(GenerationProbability)¶
-
template<typename Signature>
double Integrate(boost::function<Signature> f, typename detail::traits<Signature>::array_type low, typename detail::traits<Signature>::array_type high, double epsabs = 1.49e-6, double epsrel = 1.49e-6, size_t limit = 0)¶
-
I3_POINTER_TYPEDEFS(RadialDistribution)¶
-
I3_POINTER_TYPEDEFS(BMSSRadialDistribution)¶
-
I3_POINTER_TYPEDEFS(SamplingSurface)¶
Variables
-
static const bimap_ParticleType I3_PROPOSAL_ParticleType_bimap¶
-
typedef boost::bimap<I3Particle::ParticleType, std::string> bimap_ParticleType¶
-
namespace [anonymous]¶
-
namespace [anonymous]¶
-
namespace [anonymous]¶
-
namespace [anonymous]¶
-
namespace [anonymous]¶
-
namespace detail¶
-
namespace histogram¶
-
namespace binning¶
-
namespace I3Surfaces¶
-
namespace std
STL namespace.
- file CompactTrack.cxx
- #include <MuonGun/CompactTrack.h>#include <icetray/serialization.h>#include <icetray/I3FrameObject.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file CompactTrack.h
- #include <dataclasses/physics/I3Particle.h>#include <dataclasses/I3Vector.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file CORSIKAGenerationProbability.cxx
- #include <MuonGun/I3MuonGun.h>#include <MuonGun/CORSIKAGenerationProbability.h>#include <dataclasses/physics/I3Particle.h>#include <MuonGun/Flux.h>#include <MuonGun/RadialDistribution.h>#include <MuonGun/EnergyDistribution.h>#include <boost/foreach.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file CORSIKAGenerationProbability.h
- #include <MuonGun/Generator.h>#include <MuonGun/SamplingSurface.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file cubature.h
Adaptive multidimensional integration of a vector of integrands.
Portions (see comments) based on HIntLib (also distributed under the GNU GPL, v2 or later), Copyright (c) 2002-2005 Rudolf Schuerer. (https://github.com/JohannesBuchner/HIntLib)
- Author
Steven G. Johnson stevenj@math.mit.edu
- Copyright
Copyright (c) 2005-2009 Steven G. Johnson SPDX-License-Identifier: GPL-2.0-or-later
Portions (see comments) based on GNU GSL (also distributed under the GNU GPL, v2 or later), Copyright (c) 1996-2000 Brian Gough. (http://www.gnu.org/software/gsl/)
This program 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 2 of the License, or (at your option) any later version.
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.
You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
Typedefs
-
typedef void (*integrand)(unsigned ndim, const double *x, void*, unsigned fdim, double *fval)¶
-
typedef void (*integrand_v)(unsigned ndim, unsigned npt, const double *x, void*, unsigned fdim, double *fval)¶
Functions
-
int adapt_integrate(unsigned fdim, integrand f, void *fdata, unsigned dim, const double *xmin, const double *xmax, unsigned maxEval, double reqAbsError, double reqRelError, double *val, double *err)¶
-
int adapt_integrate_v(unsigned fdim, integrand_v f, void *fdata, unsigned dim, const double *xmin, const double *xmax, unsigned maxEval, double reqAbsError, double reqRelError, double *val, double *err)¶
- file Cylinder.cxx
- #include <MuonGun/Cylinder.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file Cylinder.h
- #include <MuonGun/SamplingSurface.h>#include <phys-services/surfaces/detail/CylinderBase.h>#include <MuonGun/UprightSurface.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
- I3_CLASS_VERSION (I3MuonGun::Cylinder, 1)
- file EnergyDependentSurfaceInjector.cxx
- #include <MuonGun/EnergyDependentSurfaceInjector.h>#include <MuonGun/I3MuonGun.h>#include <MuonGun/SamplingSurface.h>#include <MuonGun/Cylinder.h>#include <MuonGun/Flux.h>#include <MuonGun/EnergyDistribution.h>#include <MuonGun/RadialDistribution.h>#include <dataclasses/I3Constants.h>#include <dataclasses/physics/I3Particle.h>#include <dataclasses/physics/I3MCTreeUtils.h>#include <phys-services/I3RandomService.h>#include <boost/bind.hpp>#include <boost/foreach.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
-
I3_SERIALIZABLE(I3MuonGun::SurfaceScalingFunction)¶
-
I3_SERIALIZABLE(I3MuonGun::ConstantSurfaceScalingFunction)¶
-
I3_SERIALIZABLE(I3MuonGun::BasicSurfaceScalingFunction)¶
-
I3_SERIALIZABLE(I3MuonGun::EnergyDependentSurfaceInjector)¶
- file EnergyDependentSurfaceInjector.h
- #include <MuonGun/Generator.h>#include <MuonGun/StaticSurfaceInjector.h>#include <boost/function.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
- I3_CLASS_VERSION (I3MuonGun::SurfaceScalingFunction, 0)
- I3_CLASS_VERSION (I3MuonGun::ConstantSurfaceScalingFunction, 0)
- I3_CLASS_VERSION (I3MuonGun::BasicSurfaceScalingFunction, 0)
- I3_CLASS_VERSION (I3MuonGun::EnergyDependentSurfaceInjector, 0)
- file EnergyDistribution.cxx
- #include <MuonGun/EnergyDistribution.h>#include <MuonGun/RadialDistribution.h>#include <phys-services/I3RandomService.h>#include <MuonGun/EnsembleSampler.h>#include <boost/bind.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
-
I3_SERIALIZABLE(I3MuonGun::EnergyDistribution)¶
-
I3_SERIALIZABLE(I3MuonGun::SplineEnergyDistribution)¶
-
I3_SERIALIZABLE(I3MuonGun::BMSSEnergyDistribution)¶
-
I3_SERIALIZABLE(I3MuonGun::OffsetPowerLaw)¶
- file EnergyDistribution.h
- #include <icetray/I3Units.h>#include <icetray/serialization.h>#include <icetray/I3PointerTypedefs.h>#include <MuonGun/SplineTable.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
- I3_CLASS_VERSION (I3MuonGun::EnergyDistribution, 0)
- I3_CLASS_VERSION (I3MuonGun::SplineEnergyDistribution, 0)
- I3_CLASS_VERSION (I3MuonGun::BMSSEnergyDistribution, 0)
- I3_CLASS_VERSION (I3MuonGun::OffsetPowerLaw, 0)
- file EnsembleSampler.h
- #include <MuonGun/I3MuonGun.h>#include <icetray/I3Logging.h>#include <vector>#include <cmath>#include <boost/foreach.hpp>#include <phys-services/I3RandomService.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen jakob.van.santen@desy.de
- file ExtrudedPolygon.cxx
- #include <MuonGun/ExtrudedPolygon.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen jakob.van.santen@desy.de
Functions
-
I3_SERIALIZABLE(I3MuonGun::ExtrudedPolygon)¶
- file ExtrudedPolygon.h
- #include <MuonGun/SamplingSurface.h>#include <phys-services/surfaces/detail/ExtrudedPolygonBase.h>#include <MuonGun/UprightSurface.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen jakob.van.santen@desy.de
Functions
- I3_CLASS_VERSION (I3MuonGun::ExtrudedPolygon, 0)
- file Floodlight.cxx
- #include <MuonGun/Floodlight.h>#include <MuonGun/I3MuonGun.h>#include <MuonGun/EnergyDependentSurfaceInjector.h>#include <MuonGun/Cylinder.h>#include <dataclasses/I3Position.h>#include <dataclasses/I3Direction.h>#include <dataclasses/physics/I3Particle.h>#include <dataclasses/physics/I3MCTreeUtils.h>#include <boost/foreach.hpp>
Functions
-
I3_SERIALIZABLE(I3MuonGun::Floodlight)¶
-
I3_SERIALIZABLE(I3MuonGun::Floodlight)¶
- file Floodlight.h
- #include <MuonGun/Generator.h>#include <MuonGun/SamplingSurface.h>#include <MuonGun/EnergyDistribution.h>
Functions
- I3_CLASS_VERSION (I3MuonGun::Floodlight, 1)
- file Flux.cxx
- #include <MuonGun/Flux.h>#include <icetray/I3Units.h>#include <icetray/I3Logging.h>#include <limits>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file Flux.h
- #include <MuonGun/SplineTable.h>#include <icetray/I3PointerTypedefs.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file Generator.cxx
- #include <MuonGun/Generator.h>#include <MuonGun/SamplingSurface.h>#include <icetray/I3Module.h>#include <dataclasses/physics/I3MCTreeUtils.h>#include <dataclasses/I3Double.h>#include <phys-services/I3RandomService.h>#include <boost/make_shared.hpp>#include <boost/foreach.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
-
I3_SERIALIZABLE(I3MuonGun::GenerationProbability)¶
-
I3_SERIALIZABLE(I3MuonGun::GenerationProbabilityCollection)¶
-
I3_MODULE(I3MuonGun::GeneratorModule)¶
- file Generator.h
- #include <list>#include <vector>#include <boost/make_shared.hpp>#include <icetray/I3PointerTypedefs.h>#include <icetray/I3FrameObject.h>#include <icetray/serialization.h>#include <I3/hash_map.h>#include <dataclasses/physics/detail/I3MCTree_fwd.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
- I3_CLASS_VERSION (I3MuonGun::Generator, 0)
- I3_CLASS_VERSION (I3MuonGun::GenerationProbability, 0)
- I3_CLASS_VERSION (I3MuonGun::GenerationProbabilityCollection, 0)
- file histogram.cxx
- #include <MuonGun/histogram.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file histogram.h
- #include <boost/function.hpp>#include <boost/make_shared.hpp>#include <boost/variant.hpp>#include <boost/multi_array.hpp>#include <boost/array.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file I3MuonGun.cxx
- #include <MuonGun/I3MuonGun.h>#include <dataclasses/physics/I3Particle.h>#include <dataclasses/I3Constants.h>#include <phys-services/I3RandomService.h>#include <icetray/I3Logging.h>#include <icetray/I3Units.h>#include <gsl/gsl_integration.h>#include <gsl/gsl_errno.h>#include <cubature/cubature.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Variables
-
gsl_error_handler_t *handler_¶
- file I3MuonGun.h
- #include <icetray/I3PointerTypedefs.h>#include <boost/function.hpp>#include <boost/array.hpp>#include <cubature/cubature.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file MuonGunBackgroundService.cxx
- #include <MuonGun/MuonGunBackgroundService.h>#include <MuonGun/SamplingSurface.h>#include <dataclasses/physics/I3MCTreeUtils.h>#include <phys-services/I3RandomService.h>#include <boost/make_shared.hpp>#include <boost/foreach.hpp>
- Rcs
Generator.cxx 137064 2015-08-31 18:24:47Z jvansanten
- Author
Jakob van Santen vansanten@wisc.edu
- Rcs
137064
- Rcs
2015-08-31 13:24:47 -0500 (Mon, 31 Aug 2015)
- file MuonGunBackgroundService.h
- #include <sim-services/I3GeneratorService.h>#include <MuonGun/Generator.h>#include <dataclasses/physics/I3MCTreeUtils.h>#include <phys-services/I3RandomService.h>
Functions
-
I3_POINTER_TYPEDEFS(MuonGunBackgroundService)¶
-
I3_POINTER_TYPEDEFS(MuonGunBackgroundService)¶
- file Muonitron.cxx
- #include <MuonGun/Muonitron.h>#include <MuonGun/CompactTrack.h>#include <dataclasses/physics/I3MCTree.h>#include <dataclasses/physics/I3MCTreeUtils.h>#include <phys-services/I3Calculator.h>#include <simclasses/I3MMCTrack.h>#include <boost/foreach.hpp>#include <boost/make_shared.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file Muonitron.h
- #include <icetray/I3Module.h>#include <dataclasses/I3Constants.h>#include <dataclasses/physics/I3Particle.h>#include <MuonGun/MuonPropagator.h>#include <phys-services/surfaces/Sphere.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file MuonPropagator.cxx
- #include <boost/assign.hpp>#include <boost/foreach.hpp>#include <boost/bimap.hpp>#include “MuonGun/MuonPropagator.h”
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file MuonPropagator.h
- #include “MuonGun/I3MuonGun.h”#include “PROPOSAL/PROPOSAL.h”#include “dataclasses/physics/I3Particle.h”#include “icetray/I3Units.h”#include “phys-services/surfaces/Surface.h”
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file NaturalRateInjector.cxx
- #include <MuonGun/I3MuonGun.h>#include <MuonGun/NaturalRateInjector.h>#include <MuonGun/Cylinder.h>#include <dataclasses/I3Constants.h>#include <boost/bind.hpp>#include <boost/make_shared.hpp>#include <boost/foreach.hpp>#include <icetray/I3Module.h>#include <dataclasses/I3Double.h>#include <dataclasses/physics/I3MCTree.h>#include <dataclasses/physics/I3MCTreeUtils.h>#include <phys-services/I3RandomService.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
-
I3_SERIALIZABLE(I3MuonGun::NaturalRateInjector)¶
- file NaturalRateInjector.h
- #include <MuonGun/Generator.h>#include <MuonGun/SamplingSurface.h>#include <MuonGun/Flux.h>#include <MuonGun/RadialDistribution.h>#include <MuonGun/EnergyDistribution.h>#include <dataclasses/physics/I3Particle.h>#include <boost/tuple/tuple.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
- I3_CLASS_VERSION (I3MuonGun::NaturalRateInjector, 0)
- file RadialDistribution.cxx
- #include <MuonGun/RadialDistribution.h>#include <phys-services/I3RandomService.h>#include <icetray/I3Units.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
-
I3_SERIALIZABLE(I3MuonGun::RadialDistribution)¶
-
I3_SERIALIZABLE(I3MuonGun::SplineRadialDistribution)¶
-
I3_SERIALIZABLE(I3MuonGun::BMSSRadialDistribution)¶
- file RadialDistribution.h
- #include <MuonGun/SplineTable.h>#include <icetray/I3PointerTypedefs.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
- I3_CLASS_VERSION (I3MuonGun::RadialDistribution, 0)
- I3_CLASS_VERSION (I3MuonGun::SplineRadialDistribution, 0)
- file SamplingSurface.cxx
- #include <MuonGun/SamplingSurface.h>#include <icetray/I3Logging.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
-
I3_SERIALIZABLE(I3MuonGun::SamplingSurface)¶
- file SamplingSurface.h
- #include <phys-services/surfaces/SamplingSurface.h>#include <boost/function.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
- I3_CLASS_VERSION (I3MuonGun::SamplingSurface, 0)
- file SplineTable.cxx
- #include “MuonGun/SplineTable.h”#include <cerrno>#include <cstdlib>#include <stdexcept>#include <serialization/binary_object.hpp>#include <icetray/I3Logging.h>
Implementation of spline table interface
Copyright (c) 2013 the IceCube Collaboration $Id$
- Date
$Date$
- Author
jvansanten
Functions
-
I3_SPLIT_SERIALIZABLE(I3MuonGun::SplineTable)¶
- file SplineTable.h
- #include <memory>#include <string>#include <vector>#include <utility>#include <photospline/splinetable.h>#include <icetray/I3FrameObject.h>#include <icetray/serialization.h>
Definition of spline table interface
Copyright (c) 2013 the IceCube Collaboration $Id$
- Date
$Date$
- Author
jvansanten
Functions
- I3_CLASS_VERSION (I3MuonGun::SplineTable, 0)
- file StaticSurfaceInjector.cxx
- #include <MuonGun/I3MuonGun.h>#include <MuonGun/StaticSurfaceInjector.h>#include <MuonGun/Cylinder.h>#include <dataclasses/I3Constants.h>#include <boost/bind.hpp>#include <boost/make_shared.hpp>#include <boost/foreach.hpp>#include <icetray/I3Module.h>#include <dataclasses/I3Double.h>#include <dataclasses/physics/I3MCTree.h>#include <dataclasses/physics/I3MCTreeUtils.h>#include <phys-services/I3RandomService.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
-
I3_SERIALIZABLE(I3MuonGun::StaticSurfaceInjector)¶
- file StaticSurfaceInjector.h
- #include <MuonGun/Generator.h>#include <MuonGun/SamplingSurface.h>#include <MuonGun/Flux.h>#include <MuonGun/RadialDistribution.h>#include <MuonGun/EnergyDistribution.h>#include <dataclasses/physics/I3Particle.h>#include <boost/tuple/tuple.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
- I3_CLASS_VERSION (I3MuonGun::StaticSurfaceInjector, 1)
- file Track.cxx
- #include <MuonGun/Track.h>#include <simclasses/I3MMCTrack.h>#include <boost/foreach.hpp>#include <boost/next_prior.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file Track.h
- #include <dataclasses/physics/I3Particle.h>#include <dataclasses/physics/I3MCTree.h>#include <boost/tuple/tuple.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Typedefs
-
typedef I3Vector<I3MMCTrack> I3MMCTrackList¶
Functions
-
I3_FORWARD_DECLARATION(I3MMCTrack)¶
- file TrackBinner.cxx
- #include <MuonGun/TrackBinner.h>#include <boost/assign.hpp>#include <boost/foreach.hpp>#include <icetray/I3Units.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file TrackBinner.h
- #include <MuonGun/histogram.h>#include <MuonGun/CompactTrack.h>#include <dataclasses/physics/I3MCTree.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- file UprightSurface.h
- #include <MuonGun/I3MuonGun.h>#include <boost/bind.hpp>
- file WeightCalculator.cxx
- #include <MuonGun/WeightCalculator.h>#include <MuonGun/Generator.h>#include <MuonGun/SamplingSurface.h>#include <MuonGun/Cylinder.h>#include <MuonGun/Flux.h>#include <MuonGun/RadialDistribution.h>#include <MuonGun/EnergyDistribution.h>#include <MuonGun/I3MuonGun.h>#include <MuonGun/Track.h>#include <boost/foreach.hpp>#include <icetray/I3Module.h>#include <dataclasses/physics/I3MCTreeUtils.h>#include <dataclasses/I3Double.h>#include <simclasses/I3MMCTrack.h>#include <phys-services/I3Calculator.h>#include <boost/make_shared.hpp>#include <boost/numeric/ublas/vector.hpp>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
Functions
-
I3_MODULE(I3MuonGun::WeightCalculatorModule)¶
- file WeightCalculator.h
- #include <MuonGun/Generator.h>#include <icetray/I3PointerTypedefs.h>#include <icetray/I3Frame.h>#include <dataclasses/physics/I3Particle.h>#include <dataclasses/physics/I3MCTree.h>#include <tableio/I3Converter.h>
$Id$
$Revision$ $Date$
- Author
Jakob van Santen vansanten@wisc.edu
- dir cubature
- dir icetray
- dir MuonGun
- dir MuonGun
- dir MuonGun
- dir private
- dir public