Core Classes

group Core

Core classes used to build CRPropa.

Typedefs

typedef Grid<double> Grid1d
typedef Grid<float> Grid1f
typedef Grid<Vector3f> Grid3f
typedef Grid<Vector3d> Grid3d
typedef Vector3<double> Vector3d
typedef Vector3<float> Vector3f

Functions

std::string getDataPath(std::string filename)
std::string getInstallPrefix()
inline int digit(const int &value, const int &d)
template<typename T>
T clip(const T &x, const T &lower, const T &upper)
double interpolate(double x, const std::vector<double> &X, const std::vector<double> &Y)
double interpolate2d(double x, double y, const std::vector<double> &X, const std::vector<double> &Y, const std::vector<double> &Z)
double interpolateEquidistant(double x, double lo, double hi, const std::vector<double> &Y)
size_t closestIndex(double x, const std::vector<double> &X)
Vector3f meanFieldVector(ref_ptr<Grid3f> grid)

Evaluate the mean vector of all grid points.

Parameters

grid – a vector grid (Grid3f)

Returns

The mean of all grid points along each axis

double meanFieldStrength(ref_ptr<Grid1f> grid)

Evaluate the mean of all grid points.

Parameters

grid – a scalar grid (Grid1f)

Returns

The mean of all grid points

double meanFieldStrength(ref_ptr<Grid3f> grid)

Evaluate the mean of all grid points.

Parameters

grid – a vector grid (Grid3f)

Returns

The mean of all grid points

double rmsFieldStrength(ref_ptr<Grid1f> grid)

Evaluate the RMS of all grid points.

Parameters

grid – a scalar grid (Grid1f)

Returns

The total RMS of all grid points.

double rmsFieldStrength(ref_ptr<Grid3f> grid)

Evaluate the RMS of all grid points.

Parameters

grid – a vector grid (Grid3f)

Returns

The total RMS of all grid points.

std::array<float, 3> rmsFieldStrengthPerAxis(ref_ptr<Grid3f> grid)

Evaluate the RMS of all grid points per axis.

Parameters

grid – a vector grid (Grid3f)

Returns

An array of length 3 with the RMS field along each axis.

void scaleGrid(ref_ptr<Grid1f> grid, double a)

Multiply all grid values by a given factor.

Parameters
  • grid – a scalar grid (Grid1f)

  • a – scaling factor that will be used to multiply all points in grid

void scaleGrid(ref_ptr<Grid3f> grid, double a)

Multiply all grid values by a given factor.

Parameters
  • grid – a vector grid (Grid3f)

  • a – scaling factor that will be used to multiply all points in grid

void fromMagneticField(ref_ptr<Grid3f> grid, ref_ptr<MagneticField> field)

Fill vector grid from provided magnetic field.

Parameters
  • grid – a vector grid (Grid3f)

  • field – the magnetic field

void fromMagneticFieldStrength(ref_ptr<Grid1f> grid, ref_ptr<MagneticField> field)

Fill scalar grid from provided magnetic field.

Parameters
  • grid – a scalar grid (Grid1f)

  • field – the magnetic field

void loadGrid(ref_ptr<Grid3f> grid, std::string filename, double conversion = 1)

Load a Grid3f from a binary file with single precision.

Parameters
  • grid – a vector grid (Grid3f)

  • filename – name of input file

  • conversion – multiply every point in grid by a conversion factor

void loadGrid(ref_ptr<Grid1f> grid, std::string filename, double conversion = 1)

Load a Grid1f from a binary file with single precision.

Parameters
  • grid – a scalar grid (Grid1f)

  • filename – name of input file

  • conversion – multiply every point in grid by a conversion factor

void dumpGrid(ref_ptr<Grid3f> grid, std::string filename, double conversion = 1)

Dump a Grid3f to a binary file.

Parameters
  • grid – a vector grid (Grid3f)

  • filename – name of input file

  • conversion – multiply every point in grid by a conversion factor

void dumpGrid(ref_ptr<Grid1f> grid, std::string filename, double conversion = 1)

Dump a Grid1f to a binary file with single precision.

Parameters
  • grid – a scalar grid (Grid1f)

  • filename – name of input file

  • conversion – multiply every point in grid by a conversion factor

void loadGridFromTxt(ref_ptr<Grid3f> grid, std::string filename, double conversion = 1)

Load a Grid3f grid from a plain text file.

Parameters
  • grid – a vector grid (Grid3f) to which the points will be loaded

  • filename – name of input file

  • conversion – multiply every point in grid by a conversion factor

void loadGridFromTxt(ref_ptr<Grid1f> grid, std::string filename, double conversion = 1)

Load a Grid1f from a plain text file

Parameters
  • grid – a scalar grid (Grid1f) to which the points will be loaded

  • filename – name of input file

  • conversion – multiply every point in grid by a conversion factor

void dumpGridToTxt(ref_ptr<Grid3f> grid, std::string filename, double conversion = 1)

Dump a Grid3f to a plain text file.

Parameters
  • grid – a vector grid (Grid3f)

  • filename – name of output file

  • conversion – multiply every point in grid by a conversion factor

void dumpGridToTxt(ref_ptr<Grid1f> grid, std::string filename, double conversion = 1)

Dump a Grid1f to a plain text file.

Parameters
  • grid – a scalar grid (Grid1f)

  • filename – name of output file

  • conversion – multiply every point in grid by a conversion factor

std::vector<std::pair<int, float>> gridPowerSpectrum(ref_ptr<Grid3f> grid)

Calculate the omnidirectional power spectrum E(k) for a given turbulent field

Parameters

grid – a three-dimensional grid

Returns

Returns a vector of pairs (k_i, E(k_i)).

inline void intrusive_ptr_add_ref(Referenced *p)
inline void intrusive_ptr_release(Referenced *p)
template<class T>
inline void swap(ref_ptr<T> &rp1, ref_ptr<T> &rp2)
template<class T>
inline T *get_pointer(const ref_ptr<T> &rp)
template<class T, class Y>
inline ref_ptr<T> static_pointer_cast(const ref_ptr<Y> &rp)
template<class T, class Y>
inline ref_ptr<T> dynamic_pointer_cast(const ref_ptr<Y> &rp)
template<class T, class Y>
inline ref_ptr<T> const_pointer_cast(const ref_ptr<Y> &rp)
template<typename T>
inline std::ostream &operator<<(std::ostream &out, const Vector3<T> &v)
template<typename T>
inline std::istream &operator>>(std::istream &in, Vector3<T> &v)
template<typename T>
inline Vector3<T> operator*(T f, const Vector3<T> &v)
class Candidate : public Referenced
#include <include/crpropa/Candidate.h>

All information about the cosmic ray.

The Candidate is a passive object, that holds the information about the state of the cosmic ray and the simulation itself.

class Surface : public Referenced
#include <Geometry.h>

A geometrical surface.

Defines a surface. Can be queried if the candidate has crossed the surface in the last step.

Subclassed by ParaxialBox, Plane, Sphere

class Plane : public Surface
#include <Geometry.h>

A plane given by a point x0 and two axes v1 and v2 with normal n = v1.cross(v2) or the normal n. Note that distance is negative on one side of the plane and positive on the other, depending on the orientation of the normal vector.

class Sphere : public Surface
#include <Geometry.h>

A sphere around point _center with radius _radius.

class ParaxialBox : public Surface
#include <Geometry.h>

A box with perpendicular surfaces aligned to the x,y,z-axes.

class GridProperties : public Referenced
#include <Grid.h>

Combines parameters that uniquely define Grid class.

template<typename T>
class Grid : public Referenced
#include <Grid.h>

Template class for fields on a periodic grid with trilinear interpolation.

The grid spacing is constant with diffrent resolution along all three axes. Values are calculated by trilinear interpolation of the surrounding 8 grid points. The grid is periodically (default) or reflectively extended. The grid sample positions are at 1/2 * size/N, 3/2 * size/N … (2N-1)/2 * size/N.

class ParticleState
#include <ParticleState.h>

State of the particle: ID, energy, position, direction.

The ParticleState defines the state of an ultra-high energy cosmic ray, which is assumed to be traveling at the exact speed of light. The cosmic ray state is defined by particle ID, energy and position and direction vector. For faster lookup mass and charge of the particle are stored as members.

class Random
#include <Random.h>

Random number generator.

Mersenne Twister random number generator &#8212; a C++ class Random Based on code by Makoto Matsumoto, Takuji Nishimura, and Shawn Cokus Richard J. Wagner v1.0 15 May 2003 rjwagner@writeme.com

class Referenced
#include <Referenced.h>

Base class for reference counting.

A form of memory management is needed to prevent memory leaks when using MPC in Python via SWIG. This base class enables reference counting. Every reference increases the reference counter, every dereference decreases it. When the counter is decreased to 0, the object is deleted. Candidate, Module, MagneticField and Source inherit from this class

Subclassed by AdvectionField, Candidate, CylindricalProjectionMap, Density, EmissionMap, Grid< T >, GridProperties, MagneticField, Module, ObserverFeature, PhotonField, SourceFeature, SourceInterface, StepLengthModifier, Surface, TurbulenceSpectrum

template<class T>
class ref_ptr
#include <Referenced.h>

Referenced pointer.

template<typename T>
class Vector3
#include <Vector3.h>

Template class for 3-vectors of type float, double, …

Allows accessing and changing the elements x, y, z directly or through the corresponding get and set methods.

Angle definitions are phi [-pi, pi]: azimuthal angle in the x-y plane, 0 pointing in x-direction theta [0, pi]: zenith angle towards the z axis, 0 pointing in z-direction