API Reference

Simulation.hpp

Some useful functions to perform simulations.

Functions

template<typename Vector3, typename Value>
std::vector<Event<Vector3, Value>> simulate_cross_sphere(Particle<Vector3, Value> &p, MediumBall<Vector3, Value> &ball, Value Tcut, const Vector3 &r0, Value R, RandomNumber<Value> &rn, Value weight0 = 1.0, unsigned long t_long_track = 1000000)

Run a simulation and collect events that cross a sphere.

Parameters:
  • Tcut – cutoff energy

  • r0 – center of the sphere

  • R – radius

template<typename Vector3, typename Value>
std::vector<std::vector<Event<Vector3, Value>>> simulate_cross_depth(Particle<Vector3, Value> &p, MediumBall<Vector3, Value> &ball, Value Tcut, std::vector<Value> depth, RandomNumber<Value> &rn, Value weight0 = 1.0, unsigned long t_long_track = 1000000)

Run a simulation and collect crossing events on several depths.

Parameters:
  • Tcut – cutoff energy

  • depth – collect events on these depth

template<typename Vector3, typename Value>
std::vector<Event<Vector3, Value>> simulate_track(Particle<Vector3, Value> &p, Medium<Vector3, Value> &m, Value Tcut, RandomNumber<Value> &rn, Value weight0 = 1.0, unsigned long t_long_track = 1000000)

Run a simulation and return the complete trajectory.

Parameters:

Tcut – cutoff energy

template<typename Vector3, typename Value>
std::vector<std::size_t> analyse_track(Particle<Vector3, Value> &p, MediumBall<Vector3, Value> &ball, Value Tcut, Value depth, RandomNumber<Value> &rn)

Initialization.hpp

Initial condition of dark matter.

Functions

template<typename Value>
Value sample_halo_speed(RandomNumber<Value> &rn, Value vmin = 0.0, Value vmax = constants::vesc, Value vesc = constants::vesc, Value v0 = constants::v0)

Sample 1D halo speed from the isotropic SHM.

Parameters:
  • rn – random number gen

  • vmin – lower bound

  • vmax – upper bound

  • vesc – escape speed

  • v0 – most probable speed

template<typename Vector3, typename Value>
Vector3 random_isotropic_vector(RandomNumber<Value> &rn, Value len = 1)

Sample a random vector isotropically with length len (default 1).

template<typename Vector3, typename Value>
Vector3 sample_halo_velocity_earth_frame(Vector3 v_earth, RandomNumber<Value> &rn, Value vmin = 0.0, Value vmax = -1.0, Value vesc = constants::vesc, Value v0 = constants::v0)

Sample 3D halo velocity in the Earth frame from the SHM.

The sampling range can be set to (vmin, vmax). If vmax <= 0, it will be set to v_earth.norm() + vesc.

Parameters:
  • v_earth – earth velocity

  • rn – random number generator

  • vmin – speed lower bound

  • vmax – speed upper bound

  • vesc – halo escape speed

  • v0 – halo most probable speed

template<typename Vector3, typename Value>
Vector3 random_incident_position(RandomNumber<Value> &rn, const Vector3 &vi, Value r_earth = constants::rEarth, Value R = 1.1 * constants::rEarth)

Sample an incident position according to the given direction.

template<typename Vector3, typename Value>
void init_halo(Particle<Vector3, Value> &p, Value t, Value vmin, const Vector3 &v_earth, RandomNumber<Value> &rn, Value R, Value vmax = -1.0, Value vesc = constants::vesc, Value v0 = constants::v0, Value radius = constants::rEarth)

Set halo DM random initial condition.

vmin and vmax are in the Earth frame.

Parameters:
  • t – initial time

  • vmin – speed lower bound

  • v_earth – velocity of the Earth

  • rn – random number generator

  • R – initial disk distance

  • vmax – speed upper bound

  • vesc – halo escape velocity

  • v0 – halo most probable speed

template<typename Value>
inline void check_non_positive(Value T)
template<typename Vector3, typename Value>
void init_Tbl(Particle<Vector3, Value> &p, Value T, Value b, Value l, RandomNumber<Value> &rn, Value t = 0.0, Value R = 1.1 * constants::rEarth, Value radius = constants::rEarth)

Set the initial condition of a particle with energy T and direction in the galactic coordinate (b, l).

template<typename Vector3, typename Value>
void init_Tbl_vertical(Particle<Vector3, Value> &p, Value T, Value b, Value l, Value t = 0.0, Value R = constants::rEarth)
template<typename Vector3, typename Value>
void init_T_point_source(Particle<Vector3, Value> &p, Value T, const Vector3 &ep, RandomNumber<Value> &rn, Value t = 0.0, Value R = 1.1 * constants::rEarth, Value radius = constants::rEarth)
template<typename Vector3, typename Value>
void init_T_isotropic(Particle<Vector3, Value> &p, Value T, RandomNumber<Value> &rn, Value t = 0.0, Value R = 1.1 * constants::rEarth, Value radius = constants::rEarth)

Set the initial condition of a particle isotropically with energy T.

template<typename Vector3, typename Value>
void init_rT_isotropic(Particle<Vector3, Value> &p, Value r, Value T, RandomNumber<Value> &rn, Value t = 0.0)

Set the initial condition of a particle isotropically on a surface with energy T.

template<typename Vector3, typename Value>
void inject(Particle<Vector3, Value> &p, Value radius = constants::rEarth)

Place the particle on a surface.

IO.hpp

Functions

template<typename T>
constexpr hid_t hid_type(T var = T())

Determine HDF5 hid_t type.

template<typename T>
void write_attr(hid_t &dset, const char *attr_name, T attr_value)
template<typename T>
void write_attr(hid_t &group, const char *dset_name, const char *attr_name, T attr_value)
template<typename D, typename A>
void write_scalar(hid_t &group, const char *dset_name, D dset_value, const char *attr_name, A attr_value)
template<typename D, typename A>
void write_scalar(const char *filename, const char *dset_name, D dset_value, const char *attr_name, A attr_value)
template<typename D, typename A>
void write_scalar(const char *filename, const char *groupname, const char *dset_name, D dset_value, const char *attr_name, A attr_value)
template<typename D>
void write_scalar(const char *filename, const char *dset_name, D value)
template<typename C>
void append_array(hid_t &group, const char *dset_name, const C &new_data)

Extend at the 0th dimension.

template<typename C>
void write_array(hid_t &group, const char *dset_name, const C &data, const std::vector<hsize_t> &dims, hid_t &type, hsize_t chunk_size)
template<typename TStore = float, typename Vector3, typename Value>
void h5dump_events_galactic(const char *filename, const char *groupname, const std::vector<Event<Vector3, Value>> &events, Value R, Value m, hsize_t chunk_size = 10000)

Dump events into a hdf5 file using galactic coordinates b and l.

This function is used to store large data. The stored data type is the same as R and m. Please refer to hdfgroup to determine the chunk size, the default is 10000.

Parameters:
  • R – radius [km]

  • m – mass [GeV]

template<typename TStore = float, typename Vector3, typename Value>
void h5dump_events_galactic(const std::string &filename, const std::string &groupname, const std::vector<Event<Vector3, Value>> &events, Value R, Value m, hsize_t chunk_size = 10000)
Parameters:
  • R – radius [km]

  • m – mass [GeV]

template<typename arr>
void h5dump_array(const char *filename, const char *groupname, const char *dset_name, const arr &v, const std::vector<hsize_t> &dims, hsize_t chunk_size)
template<typename arr>
void h5dump_array(const std::string &filename, const std::string &groupname, const std::string &dset_name, const arr &v, const std::vector<hsize_t> &dims, hsize_t chunk_size)
template<typename TStore = float, typename Vector3, typename Value>
void h5dump_events(const char *filename, const char *groupname, const std::vector<Event<Vector3, Value>> &events, hsize_t chunk_size = 0)

Dump events into a hdf5 file using Cartesian coordinates.

This function is used to dump small data at once. Complete but also redundant information of events is stored. The stored data type is Value. Please refer to hdfgroup to determine the chunk size, the default is not chunked (chunk_size = 0).

template<typename Value>
void h5load_array(const char *filename, const char *groupname, const char *dset_name, std::vector<Value> &v, std::vector<hsize_t> &dims, const hid_t &fapl = H5P_DEFAULT)
template<typename Value>
Value h5load_scalar(const char *filename, const char *groupname, const char *dset_name, const hid_t &fapl = H5P_DEFAULT)
template<typename Value>
void merge_dataset(const char *tofilename, const char *fromfilename, const char *groupname, const char *dset_name, hsize_t chunk_size = 10000)
template<typename Value>
void merge_dataset(const std::string &tofilename, const std::string &fromfilename, const std::string &groupname, const std::string &dset_name)
inline hid_t make_group(hid_t loc_id, const std::string &gname)
template<typename TStore = float>
hid_t make_dset(hid_t loc_id, const std::string &dsetname, std::vector<hsize_t> cur_dims = {}, std::vector<hsize_t> max_dims = {}, std::vector<hsize_t> chunk = {})
template<typename TStore = float, typename C>
void write_dset(hid_t ds, const C &data, std::vector<hsize_t> offset, std::vector<hsize_t> count, hid_t xferpl = H5P_DEFAULT)
template<typename TStore = float>
inline void write_dset(hid_t ds, TStore data, hid_t xferpl = H5P_DEFAULT)
template<typename TStore = float, typename Value>
std::map<std::string, hid_t> make_dsets(hid_t fid, const std::vector<std::string> &depth_groups, const std::vector<std::string> &ring_groups, const std::vector<Value> &depth, const std::vector<Value> &ring_bins, Value radius, std::size_t sample_sz, std::size_t chunk_size, bool store_time = false)
inline void close_dsets(std::map<std::string, hid_t> &dsets)
template<typename TStore = float, typename Vector3, typename Value>
void h5dump_events_galactic(std::map<std::string, hid_t> dsets, const std::string &groupname, const std::vector<Event<Vector3, Value>> &events, std::size_t offset, bool store_time = false)
template<typename T>
std::vector<T> read_dset(hid_t ds, const std::vector<hsize_t> &offset = {}, const std::vector<hsize_t> &count = {}, hid_t xferpl = H5P_DEFAULT)
template<typename T>
std::vector<T> read_dset(hid_t fd, const std::string &dsetname, const std::vector<hsize_t> &offset = {}, const std::vector<hsize_t> &count = {}, hid_t xferpl = H5P_DEFAULT)
template<typename T>
std::vector<T> read_dset(const std::string &filename, const std::string &dsetname, const std::vector<hsize_t> &offset = {}, const std::vector<hsize_t> &count = {}, unsigned int flag = H5F_ACC_RDONLY, hid_t fapl = H5P_DEFAULT, hid_t xferpl = H5P_DEFAULT)

Read dataset from a hdf5 file.

Return a one dimensional std::vector of type T flattened from the data read.

The offset and count parameters are used to specify the part of the data set to be read, if they are empty, read the entire dataset.

Parameters:
  • filename – hdf5 filename

  • dsetname – dataset name (full path)

  • offset – offset

  • count – count

  • flag – flag pass to H5Fopen

  • fapl – fapl pass to H5Fopen

  • xferpl – xferpl pass to H5Dread

template<typename Vector3, typename Value>
void dump_track_csv(const std::vector<Event<Vector3, Value>> &events, const std::string &filename)

Injector.hpp

Injectors for injecting dark matter particle.

namespace darkprop
template<typename Vector3, typename Value>
class FluxInjector : public darkprop::Injector<Vector3, Value>

Public Functions

inline explicit FluxInjector(int seed = -1, Value t_radius = constants::rEarth)
inline FluxInjector(const std::vector<Value> &t_T, const std::vector<Value> &t_flux, Value t_radius = constants::rEarth, int seed = -1)
inline void init(const std::vector<Value> &t_T, const std::vector<Value> &t_flux)
inline ~FluxInjector()
inline virtual Value next(Particle<Vector3, Value> &p) override

Public Members

Value m_point_source
Value m_Tmin
Value m_Tmax
Value m_norm

[cm^-2 sr^-1]

GSLInterpolator<Value> spectrum
template<typename Vector3, typename Value>
class HaloInjector : public darkprop::Injector<Vector3, Value>

Public Functions

inline explicit HaloInjector(Value t_vmin = 0, Value t_vmax = -1, Vector3 t_vearth = {0, 0, 240 * units::km / units::sec}, Value t_vesc = constants::vesc, Value t_v0 = constants::v0, Value t_radius = constants::rEarth, Value t_t0 = 0, int t_seed = -1)
inline void set_norm()
inline virtual Value next(Particle<Vector3, Value> &p) override

Public Members

Value m_vmin
Value m_vmax
Vector3 m_vearth
Value m_vesc
Value m_v0
Value m_norm
template<typename Vector3, typename Value>
class Injector
#include <Injector.hpp>

The abstract base class for injecting DM.

Subclassed by darkprop::FluxInjector< Vector3, Value >, darkprop::HaloInjector< Vector3, Value >, darkprop::IntensityInjector< Vector3, Value >, darkprop::SampleInjector< Vector3, Value >, darkprop::SourceInjector< Vector3, Value >

Public Functions

inline virtual ~Injector()
virtual Value next(Particle<Vector3, Value>&) = 0
inline virtual Value operator()(Particle<Vector3, Value> &p)
inline void setSeed(int t_seed)
inline void sett0(Value t_0)

Public Members

Value m_radius
Value m_t0
RandomNumber<Value> m_rn
template<typename Vector3, typename Value>
class IntensityInjector : public darkprop::Injector<Vector3, Value>

Public Functions

inline explicit IntensityInjector(int seed = -1, Value t_radius = constants::rEarth)
inline IntensityInjector(const std::vector<Value> &t_b, const std::vector<Value> &t_l, const std::vector<Value> &t_T, const std::vector<Value> &t_I, int seed = -1, Value t_radius = constants::rEarth)
inline void init(const std::vector<Value> &t_b, const std::vector<Value> &t_l, const std::vector<Value> &t_T, const std::vector<Value> &t_I)
inline void setSigma(Value t_sigma)
inline void setBaseSigma(Value t_base_sigma)
inline Value getNorm()
inline Value intensity(Value b, Value l, Value T)
inline ~IntensityInjector()
inline virtual Value next(Particle<Vector3, Value> &p) override

Private Members

GridLinearInterpolator<3, std::vector<Value>> intensity_core
Value m_Tmin
Value m_Tmax
Value m_bmin
Value m_bmax
Value m_lmin
Value m_lmax
Value m_norm

[cm^-2 sr^-1]

Value m_sigma

[cm^2]

Value m_base_sigma
template<typename Vector3, typename Value>
class SampleInjector : public darkprop::Injector<Vector3, Value>

Public Functions

inline explicit SampleInjector(int seed = -1, Value t_radius = constants::rEarth)
inline explicit SampleInjector(const std::string &filename, unsigned flag = H5F_ACC_RDONLY, int seed = -1, Value t_radius = constants::rEarth)
SampleInjector(const SampleInjector&) = delete
SampleInjector &operator=(const SampleInjector&) = delete
inline ~SampleInjector()
inline void open(const std::string &filename, unsigned int flag = H5F_ACC_RDONLY)
inline void setBufferSize(hsize_t sz)
inline void setDataRange(hsize_t t_start, hsize_t t_end)
inline void setWeightMassIdx(hsize_t idx)
inline virtual Value next(Particle<Vector3, Value> &p) override

Public Members

Value Tmin = std::numeric_limits<Value>::max()
Value Tmax = 0

Private Functions

inline void loadData()
inline void closeDsets()

Private Members

hsize_t start
hsize_t end
hsize_t buf_idx
hsize_t buf_size
hsize_t mass_idx
std::string datafile
std::vector<Value> T
std::vector<Value> b
std::vector<Value> l
std::vector<Value> weight
hid_t fd
hid_t dxpl
hid_t fapl
hid_t ds_T
hid_t ds_bl
hid_t ds_w
template<typename Vector3 = Vector3d<double>, typename Value = double>
class SourceInjector : public darkprop::Injector<Vector3d<double>, double>

Public Functions

inline explicit SourceInjector(int t_seed = -1)
inline SourceInjector(const std::vector<Value> &t_xi_T, const std::vector<Value> &t_inv_cdf_T, const std::vector<Value> &t_xi_r, const std::vector<Value> &t_inv_cdf_r, const std::string &t_method = "linear", int t_seed = -1)
inline SourceInjector(const std::string &t_r_file, const std::string &t_T_file, const std::string &t_method = "linear", int seed = -1)
inline void interpolate_r(const std::vector<Value> &t_r, const std::vector<Value> &t_spectrum)
inline void interpolate_T(const std::vector<Value> &t_T, const std::vector<Value> &t_spectrum)
inline ~SourceInjector()
inline virtual Value next(Particle<Vector3, Value> &p) override

Private Members

std::string method
GSLInterpolator<Value> inverse_cdf_r
GSLInterpolator<Value> inverse_cdf_T
namespace numerical