Magnetic Fields¶
The gamale class provides functionality to deal with galactic magnetic field lenses. The lenses are based on the healpy framework and maps extragalactic directions to observed directions, depending on the respective magnetic field model and the rigidity of the particle.
Module for handling galactic magnetic field lenses. The lenses can be created with the lens-factory: https://git.rwth-aachen.de/astro/lens-factory
-
class
gamale.
Lens
(cfname=None)[source]¶ Bases:
object
Galactic magnetic field lens class. Lenses can be created with the lens-factory tool (requires CRPropa): https://git.rwth-aachen.de/astro/lens-factory
- We use the following conventions:
the lens maps directions at the galactic border (pointing outwards back to the source) to observed directions on Earth (pointing outwards)
the Galactic coordinate system is used
spherical coordinates are avoided
for each logarithmic energy bin there is a lens part represented by a matrix
energies are given in log10(energy[eV])
the matrices (lens_parts) are in compressed sparse column format (scipy.sparse.csc)
- for each matrix M_ij
the row number i indexes the observed direction
the column number j the direction at the Galactic edge
indices are HEALPix pixel in ring scheme.
-
check_lens_part
(lp)[source]¶ Perform sanity checks and set HEALpix nside parameter.
- Parameters
lp – lens part as scipy.sparse matrix (alternative: float for log10(R/V))
-
get_lens_part
(log10e, z=1, cache=True, force=False)[source]¶ Return the matrix corresponding to a given energy log10e [log_10(energy[eV])] and charge number Z
- Parameters
log10e – energy in units log_10(energy / eV) of the lens part
z – charge number z of the lens part
cache – Caches all the loaded lens parts (increases speed, but may consume a lot of memory!)
force – Forces to take the closest available bin, even if not directly covered
- Returns
the specified lens part as scipy.sparse matrix
-
load
(cfname)[source]¶ Load and configure the lens from a config file (columns: filename minR maxR …) in order of ascending rigidity. For conventions see: https://git.rwth-aachen.de/astro/lens-factory
- Parameters
cfname – path where the config filename can be found
-
gamale.
extragalactic_vector
(mat, i)[source]¶ Return the HEALpix vector of extragalactic directions for a given matrix and observed pixel i.
- Parameters
mat – lens part as scipy.sparse matrix
i – pixel of the observed direction
- Returns
extragalactic distribution for observed direction i (HEALpix vector of size npix)
-
gamale.
flux_map
(mat, observed_map=None)[source]¶ Computes the flux (transparency) of the galactic magnetic field outside the Galaxy
- Parameters
mat – lens part, scipy sparse matrix with shape (npix, npix)
observed_map – HEALPix map of assumed observed map, if None uniform distribution
- Returns
flux map as HEALpix vector of size npix
-
gamale.
load_lens_part
(fname, transpose=False)[source]¶ Load a lens part from the given filename (should be .npz or .mldat).
- Parameters
fname – file name to save the lens part (either .npz or .mldat)
transpose – transposes the matrix (for use with NY convention)
- Returns
lens part as scipy.sparse matrix
-
gamale.
mat2nside
(mat)[source]¶ Calculate nside from a given lenspart matrice.
- Parameters
mat – lens part as scipy.sparse matrix
- Returns
Healpy nside of the lens part
-
gamale.
mean_deflection
(mat, skymap=False)[source]¶ Computes the mean deflection in radian of the given matrix.
- Parameters
mat – lens part, scipy sparse matrix with shape (npix, npix)
skymap – if not False: returns entire skymap of size npix
- Returns
mean deflection in radians
-
gamale.
observed_vector
(mat, j)[source]¶ Return the HEALpix vector of observed directions for a given matrix and extragalactic pixel j.
- Parameters
mat – lens part as scipy.sparse matrix
j – pixel of the extragalactic direction
- Returns
observed distribution for extragalactic direction j (HEALpix vector of size npix)