cmc C++ API Reference

struct beta_dist

Local definition of beta distribution.

The struct defines a pdf function, ie. the beta distribtuion not normed, though. It provides a rvs method to draw random samples from the beta distribution using gsl_rng

The parameters used for the beta distribibution (alpha and beta) are defined.

Functions are passed to I3MetropolisHastings and need to be static.

Public Functions

SET_LOGGER ("I3CascadeSimulation")

Public Static Functions

static inline double pdf(const double &x)

Beta Function \( x^{\alpha-1} (1-x)^{\beta-1}\)

Parameters:

x

static inline double rvs(const gsl_rng *random)

Returns a random variate drawn from the beta distribution

Uses the gsl_ran_beta method

Parameters:

random – gsl_rng instance

Public Static Attributes

static double alpha
static double beta
class I3CascadeDevelopment

This abstract class defines the interface of a cascade development description.

Cascade longitudinal base class for parametrization and simulation classes

copyright (C) 2007 the icecube collaboration

Rcs

I3CascadeDevelopment.h 32874 2007-06-04 18:22:36Z bvoigt

An implementation of this class would define a Simulate method that calculates the length and the longitudinal energy deposit profile of a cascade of given energy.

Version

Rcs

32874

Date

Rcs

2007-06-04 12:22:36 -0600 (Mon, 04 Jun 2007)

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by:

Rcs

bvoigt

After this method has been invoked, the length and longitudinal profile are accessible via the corresponding methods of this class.

Implementations can be either a pure parametrization or a simulation of the shower development.

Author

Bernhard Voigt

Subclassed by I3CascadeParametrization, I3CascadeSimulation

Public Functions

inline virtual ~I3CascadeDevelopment()

Virtual destructor

virtual void Simulate(I3Particle::ParticleType type, double energy) = 0

Simulates the cascade development

This method must be called in order to get the length and/or longitudinal profile.

Todo:

Add shower type parameter (EM or Hardonic)

Parameters:

energy – energy of the shower

double GetLength()

Returns the lenght of a shower units of radiation length

I3CascadeDevelopment::Simulate must be called previously in order to get the shower length

Returns:

the length of the shower in units of radiation length

double GetEnergyDeposit(int a, int b)

Returns the energy loss of the shower between a and b (units are radiation length)

I3CascadeDevelopment::Simualte must be called beforehand.

Parameters:
  • a – start position (units are radiation length)

  • b – end position (units are radiation length)

Returns:

fractional energy loss between a and b (excluding b)

inline std::vector<double> GetEnergyLossProfile()

Returns the energy loss profile filled in the simulate method.

Public Static Attributes

static const double RADIATION_LENGTH = 0.358 / 0.9216

radiation length in ice in meter

Implementation I3CascadeDevelopment class.

copyright (C) 2007 the icecube collaboration

$Id$

Version

$Revision$

Date

$LastChangedDate$

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by: $LastChangedBy$

Protected Attributes

std::vector<double> energyLossProfile_

pointer to the vector holding the energy loss profile at every radiation length filled in Simulate method.

Private Functions

SET_LOGGER ("I3CascadeDevelopment")
class I3CascadeMCService : public I3PropagatorService

Cascade Monte Carlo Service.

This module provides the functionality to replace cascade-like in the MCTree with with a list of sub-cascades to simulate the cascade longitudinal development.

Public Functions

inline ~I3CascadeMCService()
I3CascadeMCService(I3RandomServicePtr r)

Constructor

virtual std::vector<I3Particle> Propagate(I3Particle &p, DiagnosticMapPtr frame, I3FramePtr)
virtual void SetRandomNumberGenerator(I3RandomServicePtr random)
void Simulate(I3Particle &c, std::vector<I3Particle> &d)
void SetEnergyThresholdMuons(double t)
void SetMaxMuons(int i)
void SetThresholdSplit(double t)
void SetSegmentMaxEnergy(double t)
void SetEnergyThresholdSimulation(double t)
void SetStepWidth(int w)

Private Functions

SET_LOGGER ("I3CascadeMCService")

Private Members

I3RandomServicePtr random_

random service provided by the frame

double energyThresholdMuons_

hadronic particles with an energy exceeding this threshold will generate muons.

int maxMuons_

maximal number of generated muons.

double energyThresholdSplit_

particles with an energy exceeding this threshold will be split

double energyThresholdSimulation_

particles with an energy exceeding this threshold will be split and the the longitudinal development is simulated

int stepWidth_

distance between individual sub-cascades in units of radiation length

boost::shared_ptr<I3CascadeSplit> splitter_

cascade splitting class

boost::shared_ptr<I3CascadeMuonSplit> muonSplitter_

cascade moun splitting class

class I3CascadeMuonSplit : public I3CascadeMCCommon::Sampler

This class defines methods to generate muons from a cascade.

Generation of muons for a realistic treatment of the longitudinal development of cascasdes

copyright (C) 2007 the icecube collaboration

Rcs

I3CascadeSplit.h 33096 2007-06-11 15:09:40Z bvoigt

There is only one method of interest to the user which is the split method. It takes a I3Particle of a cascade typ, generates a set of muons, reduces the energy of the cascade and returns the corresponding I3Particles for the muons.

Version

Rcs

33096

Date

Rcs

2007-06-11 17:09:40 +0200 (Mon, 11 Jun 2007)

Author

Sebastian Panknin panknin@physik.hu-berlin.de Last changed by:

Rcs

bvoigt

Author

Sebastian Panknin

Public Functions

I3CascadeMuonSplit(I3RandomServicePtr random)

Creates a I3CascadeMuonSplit object

Parameters:

random – service pointer

virtual ~I3CascadeMuonSplit()

Default destructor

std::vector<I3Particle> GenerateMuons(I3Particle &cascade)

Generate Muons for a given I3Particle.

For a (hadronic) cascade-like particle muons according to the cascade energy are generated. The number and energies of the muons is taken from a parameterization. The muons have the same origian and the same direction as the orgriginal cascade. The are positivly charged (I3Particle::MuPlus). The energy of the cascade is reduced by the energy of the muons.

Parameters:

cascade – I3Particle

Returns:

vector of I3Particles

void SetEnergyThresholdMuons(double threshold)

set the energyThresholdMuons

void SetMaxMuons(int max)

set the maxMuons number

Public Static Attributes

static const double DEFAULT_ENERGY_THRESHOLD = 1 * I3Units::GeV
static const int DEFAULT_MAX_MUONS = 10

Private Functions

SET_LOGGER ("I3CascadeMuonSplit")
I3CascadeMuonSplit(const I3CascadeMuonSplit &cascadeMuonSplit)

Private Members

double energyThresholdMuons_

energy treshold for muons to consider

int maxMuons_

maximal number of muons to consider

class I3CascadeParametrization : public I3CascadeDevelopment, public I3CascadeMCCommon::Sampler

This class implements the longitudinal shower profile of electromagnetic cascades for energies below the LPM Threshold (roughly 1PeV in Ice)

Cascade longitudinal energy profile parametrization

copyright (C) 2007 the icecube collaboration

Rcs

I3CascadeParametrization.h 32874 2007-06-04 18:22:36Z bvoigt

The parametrization is take from the Particle Data Booklet, Section 26.5 The parameters a,b of the dE/dt function described their have been determined by Christoher Wiebusch using GEANT 3. They are taken from his thesis.

Version

Rcs

32874

Date

Rcs

2007-06-04 12:22:36 -0600 (Mon, 04 Jun 2007)

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by:

Rcs

bvoigt

Author

Bernhard Voigt

Public Functions

I3CascadeParametrization()

Default constructor

inline ~I3CascadeParametrization()

Default destructor

virtual void Simulate(I3Particle::ParticleType type, double energy)

Evaluates the longitudinal shower profile parametrization for an EM shower of given energy

The parametrization is taken from PDG Booklet Section (Particle Passage through matter)

\(\frac{dE}{dt} = E\ b\ \frac{(bt)^{a-1}\ e^{-bt}}{\Gamma(a)}\)

The parameters used are a, b which have been determined by Christopher Wiebusch and quoted in his thesis.

The energy deposit is calculated in steps of one radiation length. It is internaly stored and can be accessed using the I3CascadeDevelopment::GetEnergyDeposit method.

Todo:

add shower type parameter (EM or Hardonic)

Parameters:

energy – energy of the shower

double dE_dt(double energy, double t)

Calculates the energy loss function dE/dt

Parameters:
  • energy – of the shower

  • t – is the depth in radiation length units

inline void SetThreshold(double threshold)

Sets the threshold for the parametrization

Calculation of energy deposit is quit, when the deposit in one radiation length step is less than this threshold, default is 1 TeV

Parameters:

threshold

inline double GetThreshold()

Returns the threshold of energy deposit calculation

Public Static Attributes

static const double DEFAULT_THRESHOLD = 1 * I3Units::TeV

threshold for energy deposit calculation

Implementation I3CascadeParametrization class.

copyright (C) 2007 the icecube collaboration

$Id$

Version

$Revision$

Date

$LastChangedDate$

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by: $LastChangedBy$

Private Functions

SET_LOGGER ("I3CascadeParametrization")

Private Members

double threshold_

Threshold for calculation energy deposit.

Calculation of energy deposit is quit, when the deposit in one radiation length step is less than this threshold

Private Static Attributes

static const double A_0

parameter for ‘a’ parameter calculation, no units.

From Wiebusch thesis

static const double A_1
static const double B

paramter for dE/dt function

From Wiebusch thesis

class I3CascadeSimulation : public I3CascadeDevelopment, public I3CascadeMCCommon::Sampler

This class implements a Cascade development simulation including the LPM effect.

Description of the simulation can be found here: http://wiki.icecube.wisc.edu/index.php/Cascade_Simulation

Cross section formulars used are defined in I3CascadeSimulationCrossSections

Author

Bernhard Voigt

Public Functions

I3CascadeSimulation(I3RandomServicePtr random)

Constructor.

Constructor with random service given

same as above, except random is not created

Parameters:

random – random number service from icetray frame

~I3CascadeSimulation()

Destructor.

Frees the gsl objects

Destructor - frees all gsl struct that have been initialized

virtual void SetRandomNumberGenerator(I3RandomServicePtr)
virtual void Simulate(I3Particle::ParticleType type, double energy)

Simulates a cascade of the given energy

Parameters:

energy

inline void SetThreshold(double thresholdEnergy)

Sets the energy threshold down to which particles are tracked (Default is 50 TeV)

This shouldn’t be lower than 1 GeV, since this is the limit for the Bremsstrahlung interaction

Parameters:

thresholdEnergy

inline double GetThreshold()

Returns the threshold energy down to which particles are tracked

Returns:

threshold energy

Public Static Attributes

static const double DEFAULT_THRESHOLD = 50 * I3Units::TeV

default energy threshold down to which particles are tracked

Private Functions

SET_LOGGER ("I3CascadeSimulation")
I3CascadeSimulation(const I3CascadeSimulation&)
double PairProductionMeanFreePath(double energy)

Cacluates the pair productio mean free path for the given energy

Parameters:

energy

double BremsstrahlungMeanFreePath(double energy)

Cacluates the bremsstrahlung radiation length (cut off energy is defined as 1 GeV)

Parameters:

energy

double SampleBremsstrahlungCrossSection(double energy)

Samples the bremsstrahlung cross section to get a fractional energy of a secondary particle produced in an interaction with a incident particle of the given energy

Parameters:

energy

double SamplePairProductionCrossSection(double energy)

Samples the pair production cross section to get a fractional energy of a secondary particle produced in an interaction with a incident particle of the given energy

Parameters:

energy

inline void Init()

Inits the simulation

Interpolation of the mean free paths is performed The metropolis hastings random samplers are initialized

Calls the interpolation routine and the mh-sampler initialization

void InitMetropolisHastings()

Inits the metropolis hastings random samplers

Creates different metropolis hasting samplers for pair production and bremsstrahlung cross sections

The parameters of the proposal distribution (beta distribution) are set individualy in the sample methods below.

There are different instances of the sampler since the sampler keeps the old state and they shouldn’t mix accross the different cross sections for high and low energies, as well as for the different proccesses.

void BurnIn()

Draw and discard samples from the MH samplers to make them independent of initial conditions

void InterpolPairProductionMeanFreePath(double lower = 2., double upper = 13.0, int points = 100)

Calcuation and interpolation of the pair production mean free path in an energy range specified by lower and upper (in log10) The interpolation used the given number of points. Default energy range is from 100 GeV to 10 EeV

Uses the integration and interpolation routines from gsl

Parameters:
  • lower

  • upper

  • points – integrate and interpolate pair production cross section to get the mean free path

void InterpolBremsstrahlungMeanFreePath(double lower = 2., double upper = 13.0, int points = 100)

Calcuation and interpolation of the bremsstrahlung radiation length in an energy range specified by lower and upper (in log10) The interpolation used the given number of points. Default energy range is from 100 GeV to 10 EeV

Uses the integration and interpolation routines from gsl

The integration of the cross section is done down to a lower edge of I3CascadeSimulation::BREMSSTRAHLUNGCUTOFF this is due to the infrared limit, where the cross section raises to infinite values.

Parameters:
  • lower

  • upper

  • points – integrate and interpolate bremsstrahlung cross section to get mean free path

double MeanFreePath(Particle particle)

Returns the mean free path of a particle

Parameters:

particle

void Propagate(Particle particle)

Propagates a particle

Draws the next interaction point Calls the interaction method to produce secondaries

Parameters:

particle

void Interaction(Particle particle)

Creates secondaries produced in an interaction of the given particle

Parameters:

particle

void BremsstrahlungInteraction(Particle electron)

Produces a secondary gamma

void PairProductionInteraction(Particle photon)

Produces a electron, positron pair

void Deposit(Particle particle)

Deposits the particle

Calculates the energy loss profile of the given particle and adds it to the total energy loss profile of this simulation

Private Members

std::stack<Particle> particles_

stack holding particles that have to be propagated

unsigned int particlesCreated_

counter how many particles have been created

unsigned int particlesDeleted_

counter how many particles have been deleted

beta_dist betaDist_

struct holding beta distribution parameters, provides function evaluation and random numbers

I3MetropolisHastings lowBremsSampler_

random number sampler for low energy bremsstrahlung

I3MetropolisHastings highBremsSampler_

random number sampler for high energy bremsstrahlung

I3MetropolisHastings lowPairSampler_

random number sampler for low energy pair production

I3MetropolisHastings highPairSampler_

random number sampler for high energy pair production

I3CascadeParametrization parametrization_

parametrization of energy loss profile for low energy particles

gsl_interp_accel *interpolAccBrems_

gsl spline accelerator object used for spline evaluation (bremsstrahlung)

gsl_spline *splineBrems_

gsl spline - interpolation of Bresmstrahlung mean free path

gsl_interp_accel *interpolAccPair_

gsl spline accelerator object used for spline evaluation (pair production)

gsl_spline *splinePair_

gsl spline - interpolation of Pairproduction mean free path

double threshold_

threshold down to which particles are tracked, default is 50 TeV

bool perEventBurnIn_

repeat burn-in for each event

This ensures that the events are truly independent and can be reproduced given the state of the random number generator.

Private Static Attributes

static const double BREMSSTRAHLUNGCUTOFF = 1 * I3Units::GeV

lowest energy of bremsstrahlung photons

class I3CascadeSimulationCrossSections

This class implements Bremsstrahlung and Pairproduction cross sections including suppresion due to the LPM effect.

Bremsstrahlung and Pairproduction cross sections including LPM suppresion

copyright (C) 2007 the icecube collaboration

Rcs

I3CascadeSimulationCrossSections.h 48260 2008-08-18 13:10:15Z bvoigt

The cross section are take from Spencer Klein’s review: “Suppression of bremsstrahlung and pair production due to environmental factors” http://arxiv.org/pdf/hep-ph/9802442

Version

Rcs

48260

Date

Rcs

2008-08-18 07:10:15 -0600 (Mon, 18 Aug 2008)

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by:

Rcs

bvoigt

Author

Bernhard Voigt

Public Static Functions

static double BremsstrahlungCrossSection(const double &energy, const double &y)

Caclulates the bremsstrahlung differential cross section for a given energy and given energy fraction of the secondary particle

static double PairProductionCrossSection(const double &energy, const double &y)

Caclulates the pair production differential cross section for a given energy and given energy fraction of the secondary particle

Private Functions

SET_LOGGER ("I3CascadeSimulationCrossSections")

Private Static Functions

static double MeanZ(const double &exp)

Calculates the mean Z of a component of different elements

Parameters:

exp – is the power the sum of components can be raised before the mean is calcualted (a+b+c)**power/3

static double Xi(const double &s)

Function used in the cross section calculation

static double Phi(const double &s)

Function used in the cross section calculation

static double Psi(const double &s)

Function used in the cross section calculation

static double G(const double &s)

Function used in the cross section calculation

Private Static Attributes

static const int NUMBER_OF_COMPONENTS = 3

number of elements of the material

Implementation I3CascadeSimulationCrossSections static methods.

All formulars are taken from S.Klein’s review. There’s not much documentation since it’s just copied from the paper and there’s no logic, just formulars.

copyright (C) 2007 the icecube collaboration

$Id$

Version

$Revision$

Date

$LastChangedDate$

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by: $LastChangedBy$

static const int Z[] = {1, 1, 8}

atomic charge of material, filled in the constructor with values for wather this is (1, 1, 8) (H_2 O), set in the implemenatin file

static const double X0 = 36.0

radiation length in material (here water/ice)

static const double ALPHA = 1 / 137.

fine structur constant

static const double ELECTRON_RADIUS = 2.817e-13

electron radius

static const double N = 6.022e23

Avogadro’s number.

static const double DENSITY = 0.917

Density of material (here ice)

static const double A = 18.0153

atomic weight of water molecule

static const double Z_A_RATIO = .55509

ratio of Z/A, the atomic charge over atomic weight (here for water)

static const double ENERGY_THRESHOLD = 7.7e3 * 36.0 * I3Units::GeV

energy threshold for which the LPM suppression sets in

static const double NUMBER_OF_NUKLEONS = 1 / Z_A_RATIO * N * DENSITY / A * 4 * ALPHA * pow(ELECTRON_RADIUS, 2.0) * MeanZ(2.0) * log(184.0 / (MeanZ(1.0 / 3.0)))

number of nucleons in the material seen by an electron/photon in the interaction. Taken from the reference paper and particle data booklet The value is set in the implementation file

class I3CascadeSplit : public I3CascadeMCCommon::Sampler
#include <I3CascadeSplit.h>

This class defines methods to split a cascade into sub-showers.

Cascade splitting for a realistic treatment of the longitudinal development of cascasdes

copyright (C) 2007 the icecube collaboration

Rcs

I3CascadeSplit.h 77043 2011-06-21 19:43:30Z olivas

There is only one method of interest to the user which is the split method. It takes a I3Particle of a cascade typ, splits it into a set of sub-cascades and returns the corresponding I3Particles.

Version

Rcs

77043

Date

Rcs

2011-06-21 13:43:30 -0600 (Tue, 21 Jun 2011)

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by:

Rcs

olivas

The parametrizations used for the splitting are defined in a I3CascadeParametrization object which is passed to the constructor of this class. Update docu!!!

Author

Bernhard Voigt

Public Functions

I3CascadeSplit(I3RandomServicePtr random)

Creates a I3CascadeSplit object

Parameters:

random – service pointer

virtual void SetRandomNumberGenerator(I3RandomServicePtr r)
virtual ~I3CascadeSplit()

Default destructor

std::vector<I3Particle> SplitCascade(I3Particle &cascade)

Splits a given I3Particle into multiple I3Partitcle.

A cascade-like particle is split into multiple sub-cascades to simulate the longitudinal development of the shower. The length and shower profile are obtained from the I3CascadeDevelopment object passed to the constructor of this object. A vector of I3Particles is returned, with individual particle information, each of a EM cascade type (I3Particle::EMinus).

Parameters:

cascade – I3Particle

Returns:

vector of I3Particles

inline void SetStepWidth(int stepWidth)

Set the step width between sub-cascades

Parameters:

stepWidth – distance between sub-cascades in units of radiation length

inline int GetStepWidth()

Returns the step width, distance between sub-cascades

inline void EnableSimulation(bool state)

Simulate high energy cascades, rather than using the enrgy loss profile parameterisation

Parameters:

state

inline bool IsEnableSimulation()

Returns true when the simulation for high energy particles is enabled

inline void SetSimulationThreshold(double threshold)

Simulate cascades with energies above this threshold

Parameters:

threshold – energy threshold above which cascades are simulated

inline double GetSimulationThreshold()

Returns the threshold above which cascades are simulated

inline void SetSegmentMaxEnergy(double threshold)

Ensure that all output cascades are split up to be no larger than this

Parameters:

threshold – energy threshold above which ouptut cascades are split

inline double GetSegmentMaxEnergy()

Returns the threshold above which ouptut cascades are split

Public Static Attributes

static const double DEFAULT_SIMULATION_THRESHOLD = 1 * I3Units::PeV
static const int DEFAULT_STEP_WIDTH = 3
static const double DEFAULT_MAX_ENERGY = 1 * I3Units::PeV

Private Functions

SET_LOGGER ("I3CascadeSplit")
I3CascadeSplit(const I3CascadeSplit &cascadeSplit)

Copy constructor

Parameters:

cascadeSplitI3CascadeSplit

I3Particle CreateParticle(I3Particle &particle, double energy, int steps, int length)

Creates a I3Particle from a given Particle with given energy, length in a distance of steps * radiation length away from the original I3Particle position in the direction of the original I3Particle

Parameters:
  • particle

  • energy

  • steps

  • length

Returns:

newly created I3Particle

Private Members

I3CascadeParametrization cascadeParametrization_

PSI_CascadeDevelopmentParametrization object providing a parametrization of the longitudinal shower profile and shower length.

I3CascadeSimulation cascadeSimulation_

PSI_CascadeDevelopment object to obtain the longitudinal shower profile and length from a simulation.

bool enableSimulation_

switch to enable the cascade simulation

double simulationThreshold_

threshold energy for cascade simulation, if simulation is enable showers above this energy will be simulated

int stepWidth_

distance between sub-cascades in units of radiation length

double segmentMaxEnergy_

maximum energy to output in a single cascade

Output cascades which would exceed this energy will be further split into subcascades of this energy, all at the same position.

class I3MetropolisHastings : public I3CascadeMCCommon::Sampler

This class implements a Metropolis Hastings sampler.

A metropolis hasting sampler for bremsstrahlung and pair production cross sections

copyright (C) 2007 the icecube collaboration

$Id$

The proposal distribution used is a beta distribution, which matches almost the shape of the bremsstrahlung and pair production cross sections, using the right parameters. See I3CascadeSimulation for the usage.

Version

$Revision$

Date

$LastChangedDate$

Author

Bernhard Voigt bernhard.voigt@desy.de $LastChangedBy$ It uses a pdf and a proposal distribution to efficiently produce random samples from the pdf

For details on the algorithm, take a look at : http://en.wikipedia.org/wiki/Metropolis-Hastings_algorithm

Public Functions

inline I3MetropolisHastings()

Default Constructor

I3MetropolisHastings(I3RandomServicePtr rng, double (*pdf)(const double&, const double&), double (*proposal)(const gsl_rng*), double (*proposalDistribution)(const double&), std::string name)

Constructor

The proposal gsl random instance as well as the proposalDistribution are implemented in a struct beta_dist in the header file I3CascadeSimulation.h

Implementation of Metropolis-Hastings sampler

copyright (C) 2007 the icecube collaboration

$Id$

Version

$Revision$

Date

$LastChangedDate$

Author

Bernhard Voigt bernhard.voigt@desy.de $LastChangedBy$

Parameters:
  • rng – a gsl_random instance

  • pdf – a function with two parameters, ie. the energy and fractional energy of the cross section

  • proposal – a gsl random instance sampling from the proposal distribution which is a beta distribution

  • proposalDistribution – a beta function in x

  • name – name of this instance for logging purposes

inline unsigned long GetCounter()

Get the total sample loop count

inline void ResetCounter()

Reset the total sample loop count

inline unsigned long GetYield()

Get the number of returned samples

inline void ResetYield()

Reset the number of returned samples

void BurnIn(double energy, double start, unsigned int n = 100)

Samples n times starting at start position

Used to forget the initial state

Parameters:
  • energy – parameter of underlying pdf

  • start – starting position

  • n – number of samples

double Sample(double energy, double lower = 0, double higher = 1)

Samples the underlying pdf in the given range

Parameters:
  • energy – parameter to the underlying pdf

  • lower – edge of sample range

  • higher – edge of sample range

Returns:

random number

Private Functions

SET_LOGGER ("I3MetropolisHastings")

Private Members

double (*pdf_)(const double&, const double&)

pdf, ie a cross section with energy and fractional energy parameters

double (*proposal_)(const gsl_rng*)

a gls random instance drawing samples from the proposal distribution

double (*proposalDistribution_)(const double&)

a proposal distribution, in this case a beta distribution in x

unsigned long counter_

loop counter

unsigned long yield_

counts how many successful attemps have been made to draw a new random number

std::string name_

name of this instance, for logging purpose

double x_

last drawn sample

double y_

last value of pdf(x)

double propDistY_

last value of proposalDist(x)

struct Particle

Local definition of a Particle struct.

Cascade simulation based on pair production and bremsstrahlung interactions

copyright (C) 2007 the icecube collaboration

Rcs

I3CascadeSimulation.h 35916 2007-08-23 09:34:12Z bvoigt

The particle struct holds energy and postion informations of particles in the simulation.

Version

Rcs

35916

Date

Rcs

2007-08-23 03:34:12 -0600 (Thu, 23 Aug 2007)

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by:

Rcs

bvoigt

Public Functions

inline Particle(double energy, double x, short type)

Particle with energy, at position x (1-dimensional), of type type (1 is electron, 0 is photon)

Parameters:
  • energy

  • x

  • type

Public Members

double energy
double x
short type

Public Static Attributes

static const short ELECTRON = 1
static const short PHOTON = 0
class Sampler

a mix-in class for class that hold a random number generator

Subclassed by I3CascadeMuonSplit, I3CascadeParametrization, I3CascadeSimulation, I3CascadeSplit, I3MetropolisHastings

Public Functions

inline Sampler(I3RandomServicePtr rng = I3RandomServicePtr())
inline virtual void SetRandomNumberGenerator(I3RandomServicePtr rng)

Protected Attributes

I3RandomServicePtr random_
namespace I3CascadeMCCommon

Functions

double HadronEnergyCorrection(double energy, boost::shared_ptr<I3RandomService> rand)

Scales hadronic shower energy.

At high energies hadronic showers are becoming more electromagnetic like and hence the constant scaling with .8 applied in Photonics is not suitable. Here a energy dependent scaling factor is applied Taken from M. Kowalski’s internal note on the simulation of Hadr. Cascades

Parameters:
  • energy – energy of the hadron shower

  • rand

Returns:

rescaled energy

namespace std

STL namespace.

file I3CascadeDevelopment.cxx
#include <vector>
#include “icetray/I3Units.h”
#include “I3CascadeDevelopment.h
file I3CascadeDevelopment.h
#include “icetray/I3Logging.h”
#include “dataclasses/physics/I3Particle.h”
file I3CascadeMCCommon.cxx
#include “I3CascadeMCCommon.h
#include “icetray/I3Units.h”
#include “dataclasses/physics/I3Particle.h”
#include “phys-services/I3RandomService.h”
file I3CascadeMCCommon.h
#include “icetray/I3PointerTypedefs.h”
#include “dataclasses/physics/I3Particle.h”

Functions

I3_FORWARD_DECLARATION(I3RandomService)

This file is the header file of the I3CascadeMCModule class.

copyright (C) 2007 the icecube collaboration

Rcs

I3CascadeMCModule.h 43507 2008-03-20 09:32:11Z bvoigt

Version

Rcs

43507

Date

Rcs

2008-03-20 05:32:11 -0400 (Thu, 20 Mar 2008)

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by:

Rcs

bvoigt

file I3CascadeMCService.cxx
#include “cmc/I3CascadeSplit.h
#include “cmc/I3CascadeMuonSplit.h
#include “I3CascadeMCCommon.h
#include “icetray/I3Units.h”
#include “dataclasses/physics/I3MCTreeUtils.h”
#include “sim-services/I3SimConstants.h”
#include <boost/foreach.hpp>
file I3CascadeMCService.h
#include “dataclasses/physics/I3Particle.h”
#include “dataclasses/physics/I3MCTree.h”
#include “phys-services/I3RandomService.h”
#include “sim-services/I3PropagatorService.h”
file I3CascadeMuonSplit.cxx
#include <cmath>
#include <vector>
#include <algorithm>
#include “I3CascadeMuonSplit.h
#include “dataclasses/I3Constants.h”
#include “icetray/I3Units.h”
#include “phys-services/I3GSLRandomService.h”
#include <boost/math/special_functions/relative_difference.hpp>
file I3CascadeMuonSplit.h
#include <cmath>
#include <vector>
#include “cmc/I3CascadeMCCommon.h
#include “dataclasses/physics/I3Particle.h”
#include “phys-services/I3RandomService.h”
file I3CascadeParametrization.cxx
#include <cmath>
#include <vector>
#include <numeric>
#include <gsl/gsl_sf_gamma.h>
#include “icetray/I3Units.h”
#include “sim-services/I3SimConstants.h”
#include “phys-services/I3RandomService.h”
#include “I3CascadeDevelopment.h
file I3CascadeParametrization.h
#include <vector>
#include “cmc/I3CascadeMCCommon.h
file I3CascadeSimulation.cxx
#include <algorithm>
#include <cmath>
#include <iostream>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_randist.h>
#include <gsl/gsl_spline.h>
#include <gsl/gsl_integration.h>
#include <gsl/gsl_histogram.h>
#include <climits>
#include “phys-services/I3RandomService.h”
#include “phys-services/I3GSLRandomService.h”
#include “icetray/I3Units.h”

Functions

double BremsstrahlungCrossSectionWrapper(double y, void *params)

Implementation I3CascadeSimulation class.

copyright (C) 2007 the icecube collaboration

$Id$

Version

$Revision$

Date

$LastChangedDate$

Author

Bernhard Voigt bernhard.voigt@desy.de Last changed by: $LastChangedBy$ Local definition of a wrapper function needed to call the gsl_integration routine

Parameters:
  • params – energy fraction

  • y – double value

double PairProductionCrossSectionWrapper(double y, void *params)

Local definition of a wrapper function needed to call the gsl_integration routine

Parameters:
  • params – energy fraction

  • y – double value

file I3CascadeSimulation.h
#include <vector>
#include <stack>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_randist.h>
#include <gsl/gsl_spline.h>
#include <gsl/gsl_integration.h>
#include “cmc/I3CascadeMCCommon.h
#include “phys-services/I3RandomService.h”
file I3CascadeSimulationCrossSections.cxx
#include <cmath>
#include “dataclasses/I3Constants.h”
file I3CascadeSimulationCrossSections.h
#include “icetray/I3Units.h”
#include “icetray/I3Logging.h”
file I3CascadeSplit.cxx
#include <cmath>
#include <vector>
#include <algorithm>
#include “I3CascadeDevelopment.h
#include “I3CascadeSplit.h
#include “dataclasses/I3Constants.h”
#include “icetray/I3Units.h”
file I3CascadeSplit.h
#include <cmath>
#include <vector>
#include “cmc/I3CascadeMCCommon.h
#include “dataclasses/physics/I3Particle.h”
#include “phys-services/I3RandomService.h”
file I3MetropolisHastings.cxx
#include “I3MetropolisHastings.h
#include “phys-services/I3RandomService.h”
#include “icetray/I3Logging.h”
file I3MetropolisHastings.h
#include “icetray/I3Logging.h”
#include “cmc/I3CascadeMCCommon.h
#include <gsl/gsl_rng.h>
#include <string>
page todo

Member I3CascadeDevelopment::Simulate  (I3Particle::ParticleType type, double energy)=0

Add shower type parameter (EM or Hardonic)

Member I3CascadeParametrization::Simulate  (I3Particle::ParticleType type, double energy)

add shower type parameter (EM or Hardonic)

dir cmc
dir cmc
dir cmc
dir icetray
dir private
dir public