MAGPIE: Monte cArlo weiGhted PIxel rEmapping

https://badge.fury.io/py/magpie-pkg.svg https://anaconda.org/knaidoo29/magpie-pkg/badges/version.svg Documentation Status https://circleci.com/gh/knaidoo29/magpie/tree/master.svg?style=svg https://codecov.io/gh/knaidoo29/magpie/branch/master/graph/badge.svg?token=P7H8FAJT43 https://anaconda.org/knaidoo29/magpie-pkg/badges/license.svg

Author

Krishna Naidoo

Version

0.2.2

Repository

https://github.com/knaidoo29/magpie

Documentation

https://magpie-doc.readthedocs.io/

Warning

MAGPIE is currently in development. Functions and classes may change. Use with caution.

Contents

Introduction

MAGPIE is a python module for remapping bins (in 1D), pixels (in 2D) and cells (in 3D) into different coordinate systems. The package will enable data to be remapped from cartesian to polar coordinates, cartesian to spherical polar coordinates and will enable rotations of these systems. When carrying out coordinate transformations we typically take the center of the pixel or cell and apply the transform without considering the surface area of the pixel (or volume of the cell). In MAGPIE this is taken into account by applying two remapping schemes. The first weights pixels to a new coordinate grid using monte-carlo integration. The second uses a higher-resolution mesh (denser than the new coordinate grid by some integer factor along each dimension) which is rebinned to the target coordinate grid. In both cases we sample the surface area or volume of the new coordinate pixels to accurately remap the data. The monte-carlo method is more accurate but scales poorly to 3D. For 2D this scheme will work very well even for moderately large datasets. The higher-resolution mesh method, while less accurate, is very fast and should be used for large data sets and all 3D transformations. In 1D these are computed exactly without requiring the approximate schemes above.

Note

MPI functionality is currently unavailable but the plan will be to implement this using mpi4py and an additional library MPIutils which will handle all of the MPI enabled functions.

Tutorials and API

Coords

Conventions

All angular coordinates are defined in radians.

Polar Coordinates

Polar coordinates are related to 2D cartesian coordinates by the relations

\[\begin{split}& x = r\cos(\phi) + x_{0},\\ & y = r\sin(\phi) + y_{0},\end{split}\]

where \(x_{0}\) and \(y_{0}\) define the center of the polar coordinate grid. The polar coordinates \(r\) lies in the range \([0, \infty)\) while \(\phi\) lies in the range \([0, 2\pi)\). The conversion from 2D cartesian coordinates to polar is given by

\[\begin{split}& r = \sqrt{(x-x_{0})^{2} + (y-y_{0})^{2}},\\ & \phi = \arctan \left(\frac{y-y_{0}}{x-x_{0}}\right).\end{split}\]

Note

This means \(\phi=0\) lies along the x-axis i.e. y=0.

Spherical Polar Coordinates

API

Conversion
magpie.coords.cart2polar
magpie.coords.cart2polar(x, y[, center=[0., 0.]])

Returns the polar coordinates for a given set of cartesian coordinates.

Parameters
xarray

x coordinates.

yarray

y coordinates.

centerlist

Center point of polar coordinate grid.

Returns
rarray

Radial coordinate.

phiarray

Phi coordinate.

magpie.coords.cart2sphere
magpie.coords.cart2sphere(x, y, z[, center=[0., 0., 0.]])

Return polar coordinates for a given set of cartesian coordinates.

Parameters
xarray

x coordinate.

yarray

y coordinate.

centerlist

Center point of polar coordinate grid.

Returns
rarray

Radial coordinates.

phiarray

Phi coordinates [0, 2pi].

thetaarray

Theta coordinates [0, pi].

magpie.coords.ortho2cart
magpie.coords.ortho2cart(x, y[, r= 1., fill_value=np.nan])

Orthographic projection to cartesian coordinates.

Parameters
xarray

X value in the orthographic projection.

yarray

Y value in the orthographic projection.

rfloat

Radius of the sphere.

fill_valuefloat

Fill values outside the sphere with this value.

Returns
zarray

Returns the z value of the cartesian coordinates.

magpie.coords.polar2cart
magpie.coords.polar2cart(r, p[, center=[0., 0.]])

Return cartesian coordinates for a given set of polar coordinates.

Parameters
rarray

Radial coordinate.

phiarray

Phi coordinate.

centerlist

Center point of polar coordinate grid.

Returns
xarray

x coordinate.

yarray

y coordinate.

magpie.coords.sphere2cart
magpie.coords.sphere2cart(r, phi, theta[, center=[0., 0., 0.]])

Converts spherical polar coordinates into cartesian coordinates.

Parameters
rarray

Radial distance.

phiarray

Longitudinal coordinates (radians = [0, 2pi]).

thetaarray

Latitude coordinates (radians = [0, pi]).

centerlist

Center point of spherical polar coordinate grid.

Returns
x, y, zarray

Euclidean coordinates.

Rotations
Utility
magpie.coords.usphere_area
magpie.coords.usphere_area(phi_min, phi_max, theta_min, theta_max)

Returns the area for a ‘square’ segment of a unit sphere given in spherical coordinates phi, theta.

Parameters
phi_minfloat

Minimum latitudinal coordinate in radians (where phi lies [0, 2pi]).

phi_maxfloat

Maximum latitudinal coordinate in radians (where phi lies [0, 2pi]).

theta_minfloat

Minimum longitudinal coordinate in radians (where theta lies [0, pi]).

theta_maxfloat

Maximum longitudinal coordinate in radians (where theta lies [0, pi]).

Returns
Areafloat

Area in square radians.

magpie.coords.sphere2lonlat
magpie.coords.sphere2lonlat(theta)

Converts the spherical coordinates theta to the longitude and latitude convention (where theta lies [-pi/2., pi/2.].

Parameters
thetaarray

Latitude given in the range [0., pi] where theta = 0 at the north pole.

Returns
Latitudearray

Latitude given in the range [-pi/2, pi/2].

magpie.coords.lonlat2sphere
magpie.coords.lonlat2sphere(latitude)

Converts from latitude to spherical coordinate convention.

Parameters
Latitudearray

Latitude given in the range [-pi/2, pi/2].

Returns
thetaarray

Latitude given in the range [0., pi] where theta = 0 at the north pole.

Grids

Theory

Tutorial

API

Cartesian
magpie.grids.get_xedges
magpie.grids.get_xedges(length, ngrid[, xmin=0.])

Returns the edges of a uniform grid along one axis.

Parameters
lengthfloat

Length of the grid.

ngridint

Number of grid points.

xminfloat, optional

Minimum value of the grid.

Returns
xedgesarray

The edges of a uniform grid along one axis.

magpie.grids.xedges2mid
magpie.grids.xedges2mid(xedges)

Returns the mid-point of a uniform grid from the grid edges.

Parameters
xedgesarray

The edges of a uniform grid along one axis.

Returns
xmidarray

The mid-point of a uniform grid.

magpie.grids.xmid2edges
magpie.grids.xmid2edges(xmid)

Returns the mid-point of a uniform grid from the grid edges.

Parameters
xmidarray

The mid-point of a uniform grid.

Returns
xedgesarray

The edges of a uniform grid along one axis.

magpie.grids.grid1d
magpie.grids.grid1d(length, ngrid[, xmin=0., return_edges=False])

Returns the edges of a uniform grid along one axis.

Parameters
Lengthfloat

Length of the grid.

ngridint

Minimum value of the grid.

xminfloat, optional

Minimum value of the grid.

return_edgesbool, optional

If True then the edges of the uniform grid are also supplied.

Returns
xmidarray

The mid-point of a uniform grid.

xedgesarray, optional

The edges of a uniform grid along one axis.

magpie.grids.grid2d
magpie.grids.grid2d(lengths, ngrids[, mins=[None, None], return1d=False])

Returns the mid-point of a uniform grid.

Parameters
lengthsfloat or list[float]

Length of the grid.

ngridsint or list[int]

Number of grid points.

minslist[float], optional

Minimum values of the grid along each axis.

return1dbool, optional

Returns 1d mid-points.

Returns
x2d, y2d2darray

The uniform 2d grid.

xmid, ymidarray, optional

The mid-point of the uniform grid.

magpie.grids.grid3d
magpie.grids.grid3d(lengths, ngrids[, mins=[None, None, None], return1d=False])

Returns the mid-point of a uniform grid.

Parameters
lengthsfloat or list[float]

Length of the grid.

ngridsint or list[int]

Number of grid points.

minslist[float], optional

Minimum values of the grid along each axis.

return1dbool, optional

Returns 1d mid-points.

Returns
x3d, y3d, z3d3darray

The uniform 3d grid.

xmid, ymid, zmidarray, optional

The mid-point of the uniform grid.

Polar
magpie.grids.polargrid
magpie.grids.polargrid(nr, nphi[, rmin=0., rmax=1., phimin=0., phimax=2*np.pi, returns1d=False])

Constructs a polar coordinate grid.

Parameters
nrint

Number of divisions along the r-axis.

nphiint

Number of divisions along the p-axis.

rminfloat, optional

Minimum radial value, default=0.

rmaxfloat, optional

Maximum radial value, default=1.

phiminfloat, optional

Minimum phi value, default=0.

phimaxfloat, optional

Maximum phi value, default=2pi.

return1dbool, optional

Returns 1d mid-points.

Returns
r2d, p2d2darray

Radial and phi grid points in 2darray.

rmid, pmidarray, optional

The mid-point of the polar grid.

magpie.grids.polarEA_grid
magpie.grids.polarEA_grid(nr[, rmax=1., base_nphi=4])

Constructs an equal area polar grid.

Parameters
nrint

Number of divisions along the r-axis.

rmaxfloat, optional

Maximum radial value, default=1.

base_nphiint, optional

Number of pixels around r=0, default=4.

Returns
r, parray

Radial and phi mid points for each pixel in the equal area polar grid.

magpie.grids.polarEA_area
magpie.grids.polarEA_area(nr[, rmax=1., base_nphi=4])

Returns the area for an equal area polar grid.

Parameters
nrint

Number of divisions along the r-axis.

rmaxfloat, optional

Maximum radial value, default=1.

base_nphiint, optional

Number of pixels around r=0, default=4.

Returns
areafloat

The area of pixels in the equal area polar grid.

magpie.grids.polarEA_npix
magpie.grids.polarEA_npix(nr[, rmax=1., base_nphi=4])

Returns the number of pixels in an equal area polar grid.

Parameters
nrint

Number of divisions along the r-axis.

base_nphiint, optional

Number of pixels around r=0, default=4.

Returns
npixint

Number of pixels in an equal area polar grid.

Pixel

Theory

Tutorial

API

Pixel Shapes
magpie.pixel.get_square
magpie.pixel.get_square(xmin, xmax, ymin, ymax[, steps=40])

Returns the coordinates of the perimeter of a rectangle.

Parameters
xminfloat

Minimum x-value.

xmaxfloat

Maximum x-value.

yminfloat

Minimum y-value.

ymaxfloat

Maximum y-value.

stepsint, optional

Number of divisions across each line segment, default = 4.

Returns
x, yarray

Perimeter coordinates of a rectangle.

magpie.pixel.get_box
magpie.pixel.get_box(xmin, xmax, ymin, ymax, zmin, zmax[, steps=120, return_nearest=False, center=[0., 0., 0.]])

Returns 3D coordinates for a box.

Parameters
xminfloat

Minimum x value.

xmaxfloat

Maximum x value.

yminfloat

Minimum y value.

ymaxfloat

Maximum y value.

zminfloat

Minimum z value.

zmaxfloat

Maximum z value.

stepsint, optional

Number of divisions in each line element.

return_nearestbool, optional

If True only the lines for the nearest or visible faces from a defined center will be outputed.

centerlist, optional

Coordinate center used for return_nearest=True

Returns
x, y, zarray

Coordinates of the box edges.

magpie.pixel.get_arc
magpie.pixel.get_arc(radius[, phimin=0., phimax=2.*np.pi, center=[0., 0.], steps=10])

Returns the coordinates of the arc of a circle, default settings will return the coordinates along a complete circle.

Parameters
radiusfloat

Radius of the circle.

phiminfloat, optional

Minimum phi value, default=0.

phimaxfloat, optional

Maximum phi value, default=2pi.

centerlist, optional

Central x and y coordinates for the arc.

stepsint, optional

The number of points along the arc curve.

Returns
x, yarray

Coordinates along the arc of a circle.

magpie.pixel.get_disc
magpie.pixel.get_disc([rmin=0., rmax=1., phimin=0., phimax=2.*np.pi, center=[0., 0.], steps=40])

Returns the coordinates of a disc segment.

Parameters
rminfloat, optional

Minimum radial value, default=0.

rmaxfloat, optional

Maximum radial value, default=1.

phiminfloat, optional

Minimum phi value, default=0.

phimaxfloat, optional

Maximum phi value, default=2pi.

centerlist, optional

Central x and y coordinates for the arc.

stepsint, optional

The number of points along the arc curve.

Returns
x, yarray

Coordinates along the arc of a circle.

magpie.pixel.healpix_boundary_bot
healpix_boundary_bot(p, nside[, steps=20, reverse=False])

Returns the bottom side of a healpix pixel in healpix x and y coordinates.

Parameters
pint

Healpix pixel index.

nsideint

Healpix Nside.

stepsint, optional

Number of steps for the top function.

reversebool, optional

Reverse the order so the output has descending x.

Returns
xarray

X-coordinates of the bottom side of the healpix pixel.

yarray

Y-coordinates of the bottom side of the healpix pixel.

magpie.pixel.healpix_boundary_top
healpix_boundary_top(p, nside[, steps=20, reverse=False])

Returns the top side of a healpix pixel in healpix x and y coordinates.

Parameters
pint

Healpix pixel index.

nsideint

Healpix Nside.

stepsint, optional

Number of steps for the top function.

reversebool, optional

Reverse the order so the output has descending x.

Returns
xarray

X-coordinates of the top side of the healpix pixel.

yarray

Y-coordinates of the top side of the healpix pixel.

magpie.pixel.healpix_boundary
healpix_boundary(p, nside[, steps=40, reverse=False])

Returns the boundary of a healpix pixel in healpix x and y coordinates, default given in clockwise directions.

Parameters
pint

Healpix pixel index.

nsideint

Healpix Nside.

stepsint, optional

Number of steps.

reversebool, optional

Reverse the order so the output is anti-clockwise.

Returns
xarray

X-coordinates of the boundaries of the healpix pixel.

yarray

Y-coordinates of the boundaries of the healpix pixel.

Pixel Indices

Randoms

Theory

Tutorial

API

Cartesian
magpie.randoms.randoms_1d
magpie.randoms.randoms_1d(size[, xmin=0., xmax=1.])

Returns uniform randoms in 1D.

Parameters
sizeint

Number of randoms to produce.

xminfloat, optional

Minimum value.

xmaxfloat, optional

Maximum value.

Returns
xrandsarray

Uniform randoms in 1D.

magpie.randoms.randoms_2d
magpie.randoms.randoms_2d(size[, mins=[0., 0.], maxs=[1., 1.]])

Returns uniform randoms in 2D.

Parameters
sizeint

Number of randoms to produce.

minsfloat, optional

Minimum values.

maxsfloat, optional

Maximum values.

Returns
xrands, yrandsarray

Uniform randoms in 2D.

magpie.randoms.randoms_3d
magpie.randoms.randoms_3d(size[, mins=[0., 0., 0.], maxs=[1., 1., 1.]])

Returns uniform randoms in 3D.

Parameters
sizeint

Number of randoms to produce.

minsfloat, optional

Minimum values.

maxsfloat, optional

Maximum values.

Returns
xrands, yrands, zrandsarray

Uniform randoms in 3D.

Polar
magpie.randoms.randoms_polar
magpie.randoms.randoms_polar(size[, r_min=0., r_max=1., phi_min=0., phi_max=2.*np.pi])

Generates randoms for polar coordinates. Default will produce randoms within a unit circle. This can be specified to a ring segment, i.e. with inner radius r_min and outer radius r_max and specifying the angle of the ring segment.

Parameters
sizeint

Number of randoms.

r_minfloat

Minimum r.

r_maxfloat

Maximum r.

phi_minfloat

Minimum phi.

phi_maxfloat

Maximum phi.

Returns
rarray

Random radial coordinates.

phiarray

Random phi coordinates.

Unit Sphere
magpie.randoms.randoms_usphere
magpie.randoms.randoms_usphere(size[, phi_min=0., phi_max=2.*np.pi, theta_min=0., theta_max=np.pi])

Random points on the unit sphere or more generally across the surface of a sphere. The default will give randoms on the full unit sphere.

Coordinate convention:
  • phi lies in the range [0, 2pi]

  • theta lies in the rang [0, pi].

Parameters
sizeint

Number of randoms to generate.

phi_minfloat

Minimum longitude in radians.

phi_maxfloat

Maximum longitude in radians.

theta_minfloat

Minimum latitude in radians.

theta_maxfloat

Maximum longitude in radians.

Returns
phiarray

Random phi.

thetaarray

Random theta.

magpie.randoms.randoms_healpix_pixel
magpie.randoms.randoms_healpix_pixel(size, pix, nside)

Returns roughly size number of randoms inside a HEALPix pixel.

Parameters
sizeint

Average number of randoms per pixel.

pint

Pixel identifier for healpix map.

nsideint

Nside of the healpix map.

Returns
phiarray

Random phi (latitude angle) in radians.

thetaarray

Random theta (longitude angle) in radians.

Spherical
magpie.randoms.randoms_sphere_r
magpie.randoms.randoms_sphere_r(size[, r_min=0., r_max=1.])

Random radial points for a segment of a sphere (default will give randoms within a unit sphere).

Parameters
sizeint

Number of randoms to generate.

r_minfloat

Minimum radius.

r_maxfloat

Maximum radius.

Returns
rarray

Random r.

magpie.randoms.randoms_sphere
magpie.randoms.randoms_sphere(size[, r_min=0., r_max=1., phi_min=0., phi_max=2*np.pi, theta_min=0., theta_max=np.pi])

Random points inside a sphere (default will give randoms within a unit sphere). You can specify the inner and outer radii to get randoms in a shell and the region on the sky.

Coordinate convention:
  • phi lies in the range [0, 2pi]

  • theta lies in the rang [0, pi].

Parameters
sizeint

Number of randoms to generate.

r_minfloat

Minimum radius.

r_maxfloat

Maximum radius.

phi_minfloat

Minimum longitude in radians.

phi_maxfloat

Maximum longitude in radians.

theta_minfloat

Minimum latitude in radians.

theta_maxfloat

Maximum longitude in radians.

Returns
rarray

Random r.

phiarray

Random phi.

thetaarray

Random theta.

Sample PDF/CDF
magpie.randoms.pdf2cdf
magpie.randoms.pdf2cdf(xmid, pdf[, return_normpdf=True])

Calculates the CDF from a given PDF.

Parameters
xmidarray

Linearly spaced x-values given at the middle of a bin of length dx.

pdfarray

Probabilty distribution function.

return_normpdfbool, optional

Normalise PDF is also outputed.

Returns
xarray

X-coordinates.

cdfarray

Cumulative distribution function with extreme points set 0 and 1.

normpdfarray, optional

Normalised PDF.

magpie.randoms.randoms_cdf
magpie.randoms.randoms_cdf(x, cdf, size[, kind='cubic'])

Generates randoms from a given cumulative distribution function.

Parameters
xarray

X-coordinates.

cdfarray

Cumulative distribution function, extreme points must be 0 and 1 i.e. cdf[0] = 0 and cdf[-1] = 1.

sizeint

Size of the random sample.

kindstr, optional

Scipy CDF interpolation kind.

Returns
randsarray

Randoms drawn from sample CDF.

magpie.randoms.randoms_pdf
magpie.randoms.randoms_pdf(x, pdf, size[, kind='cubic'])

Generates randoms from a given probability distribution function by first calculating a CDF.

Parameters
xmidarray

Linearly spaced x-values given at the middle of a bin of length dx.

pdfarray

Probabilty distribution function.

sizeint

Size of the random sample.

kindstr, optional

Scipy CDF interpolation kind.

Returns
randsarray

Randoms drawn from sample PDF.

Subsampling
magpie.randoms.shuffle
magpie.randoms.shuffle(sample)

Shuffles the ordering of a sample.

Parameters
samplearray

Input sample data.

magpie.randoms.random_draw
magpie.randoms.random_draw(sample, size)

Draws a random sample from an input function, the algorithm ensures there can be no repeats.

Parameters
samplearray

Input sample data.

sizeint

Size of the random draws.

Returns
randsamparray

Random subsample.

magpie.randoms.random_prob_draw
magpie.randoms.random_prob_draw(sample, prob, size)

Probabilistic draw from an input sample.

Parameters
samplearray

Input sample data.

probarray, optional

The probability assigned to each sample.

sizeint, optional

Size of the probabilistic draw.

Returns
randsamparray

Random subsample.

magpie.randoms.stochastic_integer_weights
magpie.randoms.stochastic_integer_weights(weights)

Returns stochastic integer weights for an input weight. This is useful for point processes that require integer weights, where a non-integer weight can be achieved by superposition of many realisations.

Parameters
weightsarray

Input weights.

Returns
weights_SIarray

Stochastic integer weights.

magpie.randoms.stochastic_binary_weights
magpie.randoms.stochastic_binary_weights(weights)

Returns stochastic binary integer weights for an input weight. This is useful for point processes that require binary integer weights, where a non-integer weight can be achieved by superposition of many realisations.

Parameters
weightsarray

Input weights.

Returns
weights_SBarray

Stochastic binary weights.

Remap

Remapping data from one coordinate system mesh to another can play a critical role in data science. Converting one dimensional data, image-like data in two dimensions and three dimensional mesh data into another coordinate systems can be critical in enabling certain analysis pipelines to be calculated and on a more fundamental level enable better intuitive understanding of certain phenomena. Typically we assume that data provided on a grid can be defined by the cells (bins in 1D and pixels in 2D) center and map those central positions into a new coordinate system. However, this fails to take into account the pixel’s volume (length in 1D and area in 2D) and can lead to some strange artifacts (such as empty pixels in the new coordinate mesh). A possible solution to this problem is to use interpolation however while interpolation works well for theoretically generated data, for measured data this makes the implicit assumption that data measured and binned in a cell is equivalent to the value at the cells center. In general this assumption is not hugely problematic but can lead to some issues particularly when the cells are integrated values across the cell. In magpie we take a different approach, we consider each pixels volume and map to another coordinate system mesh by taking the volume-weighted mean of overlapping pixels in the original coordinate system. This is achieved through two methods, the first uses monte carlo (MC) integration to determine the weights but for large data sets this can be slow so a second scheme approximates these weights by using a higher resolution (HR) mesh.

Theory

Remapping

Remapping data \(d\) (defined in coordinate system \(X\)) to a new coordinate system \(Y\) is calculated by taking the weighted mean,

\[\tilde{d}(Y_{j}) = \sum_{i=1}^{N}w(X_{i}, Y_{j})d(X_{i}),\]

where \(\tilde{d}\) is \(d\) mapped onto coordinate system \(Y\) and \(i\) and \(j\) indicate pixels in \(X\) and \(Y\) respectively. The weights \(w(X, Y)\) are given by

\[w(X, Y) = \frac{A(X\cap Y)}{A(Y)},\]

where \(A(Y)\) is the area (length in 1D and volume in 3D) of pixel Y and \(A(X\cap Y)\) is the area of the overlapping region for pixels X and Y. For variance we can use propagation of errors to determine the variance on the remapped coordinate system although note we assume covariant terms are zero. The remapped variance is given by

\[\tilde{\Delta d}^{2}(Y_{j}) = \sum_{i=1}^{N}U(X_{i}, Y_{j})w^{2}(X_{i}, Y_{j})\Delta d^{2}(X_{i}).\]

Note, the variance for a subset of \(X\) called \(s\) is given by the relation

\[\Delta d^{2} (s) = \frac{A(s)}{A(X)}\Delta d^{2}(X),\]

this assumes that the variance observed in a bin is actually the standard error on the mean for that bin. The above relation shows that dividing by the area gives a normalised variance,

\[\delta d^{2}(X) = \frac{1}{A(X)}\Delta d^{2}(X),\]

and it follows for \(\tilde{\delta d}(Y)\equiv\delta d(X)\) for constant variance that \(U(X_{i}, Y_{j})=A(X_{i})/A(Y_{j})\) which means,

\[\tilde{\Delta d}^{2}(Y_{j}) = \sum_{i=1}^{N}\frac{A(X_{i})}{A(Y_{j})}w^{2}(X_{i}, Y_{j})\Delta d^{2}(X_{i}).\]
Weight Calculation
Monte-Carlo Method

Weights are computed by monte carlo integration – random points are distributed equally in each pixel \(Y\) and the weights are simply the normalised counts which fall into pixels in \(X\).

Higher-Resolution Method

A higher resolution mesh, determined by some input integer factor \(N_{H}\), is constructed. This means for the new coordinate system \(Y\) there are \(N_{H}\) HR pixels in each of the coordinate dimensions. The weights are simply given by

\[w(X_{i}, Y_{j}) = \frac{1}{A(Y_{j})}\sum_{k=1}^{N_{H}^{D}} A(y_{k})S(y_{k}, X_{i})\]

where \(D\) is the dimension of the coordinate system, \(y\) indicates pixels in the higher resolution mesh and \(S(y_{k}, X_{i})=1\) if the pixel center for \(y\) is within \(X_{i}\) and \(0\) otherwise.

Tutorial

API

API

Coords

Grids

Randoms

Remap

Dependencies

MAGPIE is being developed in Python 3.9 but should work on all versions >=3.4. MAGPIE is written mostly in python but with some heavy computation carried out in Fortran. Compiling the Fortran source code will require the availability of a fortran compiler such as gfortran or ifort.

The following Python modules are required:

For testing you will require nose or pytest .

Installation

MAGPIE can be installed in a variety of ways – using pip, conda or by directly cloning the repository. If you are having trouble installing MAGPIE or running MAGPIE we recommend using the conda install as this will setup the environment.

  1. Using pip:

    pip install magpie-pkg
    
  2. Using conda:

    conda install -c knaidoo29 magpie-pkg
    
  3. By cloning the github repository:

    git clone https://github.com/knaidoo29/magpie.git
    cd magpie
    python setup.py build
    python setup.py install
    

Once this is done you should be able to call MAGPIE from python:

import magpie

Support

If you have any issues with the code or want to suggest ways to improve it please open a new issue (here) or (if you don’t have a github account) email krishna.naidoo.11@ucl.ac.uk.

Version History

  • Version 0.2:
    • Restructured layout and the incorporation of documentation and unit testing for eventual first release.

    • Randoms in a variety of coordinate systems: cartesian, polar, spherical and from an input PDF/CDF and subsampling and stochastic weights.

    • Cartesian coordinate remapping in 1D, 2D and 3D.

  • Version 0.1:
    • Coordinate transformations between cartesian, polar and spherical polar coordinates.

    • Rebinning in 1D (computed exactly), in 2D and 3D via monte-carlo weighted remapping and a higher-resolution mesh.

    • Randoms in cartesian, polar and spherical polar coordinates.

    • Rotation transformations.

    • Polar coordinate utilities and integration for polar grid to radial profiles.

    • Plotting routine for orthographic projection of unit sphere data.