Skip to main content

Fast spin spherical transforms

Project description

SSHT Python Documentation

This guide is intended to explain the python interface of SSHT. For a description of the workings of SSHT see here. The python package also offers an interface to some of the functionality from ducc0, including forward and inverse transforms.

pyssht.forward

flm = pyssht.forward(f, int L, Spin=0, Method='MW', Reality=False)

Performs the forward spherical harmonic transform.

Inputs

  • f the signal on the sphere, numpy.ndarray type complex or real, ndim 2. NB different for 'MW_pole' sampling.
  • L the band limit of the signal, non-zero positive integer
  • Spin the spin of the signal, non-negative integer (default = 0)
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MW_pole' [McEwen & Wiaux sampling with the south pole as a separate double.]
    3. 'MWSS' [McEwen & Wiaux symmetric sampling]
    4. 'DH' [Driscoll & Healy sampling]
    5. 'GL' [Gauss-Legendre sampling]
  • Reality determines if the signal is real or complex, Boolean (default = False)
  • backend the backend that runs the transforms:
    1. 'SSHT' this package
    2. 'ducc' interface to ducc0. "MW_pole" is not available in this backend.
  • nthreads: number of threads when calling into the 'ducc' backend. Ignored otherwise.

Output

flm the spherical harmonic transform of f, 1D numpy.ndarray type complex

Note on 'MW_pole' sampling

This is the same as the 'MW' sampling however the south pole is expressed as a double only not a vector. Therefore the size of the array is one smaller on the (\theta) direction. f is now a tuple containing either:

  • (f_array, f_sp, phi_sp) if complex transform
  • (f_array, f_sp) if real transform

pyssht.inverse

f = pyssht.inverse(np.ndarray[ double complex, ndim=1, mode="c"] f_lm not None, L, Spin=0, Method='MW', Reality=False)

Performs the inverse spherical harmonic transform.

Inputs

  • flm the spherical harmonic transform of f, numpy.ndarray type complex, ndim 1
  • L the band limit of the signal, non-zero positive integer
  • Spin the spin of the signal, non-negative integer (default = 0)
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MW_pole' [McEwen & Wiaux sampling with the south pole as a separate double.]
    3. 'MWSS' [McEwen & Wiaux symmetric sampling]
    4. 'DH' [Driscoll & Healy sampling]
    5. 'GL' [Gauss-Legendre sampling]
  • Reality determines if the signal is real or complex, Boolean (default = False)
  • backend the backend that runs the transforms:
    1. 'SSHT' this package
    2. 'ducc' interface to ducc0. "MW_pole" is not available in this backend.
  • nthreads: number of threads when calling into the 'ducc' backend. Ignored otherwise.

Output

f the signal on the sphere, 2D numpy.ndarray type complex or real. NB different for 'MW_pole' sampling.

Note on 'MW_pole' sampling

This is the same as the 'MW' sampling however the south pole is expressed as a double only not a vector. Therefore the size of the array is one smaller on the (\theta) direction. f is now a tuple containing either:

  • (f_array, f_sp, phi_sp) if complex transform
  • (f_array, f_sp) if real transform

pyssht.forward_adjoint

f = pyssht.forward_adjoint(np.ndarray[ double complex, ndim=1, mode="c"] f_lm not None, L, Spin=0, Method='MW', Reality=False)

Performs the adjoint of the forward spherical harmonic transform.

Inputs

  • flm the spherical harmonic transform of f, numpy.ndarray type complex, ndim 1
  • L the band limit of the signal, non-zero positive integer
  • Spin the spin of the signal, non-negative integer (default = 0)
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MWSS' [McEwen & Wiaux symmetric sampling]
  • Reality determines if the signal is real or complex, Boolean (default = False)
  • backend the backend that runs the transforms:
    1. 'SSHT' this package
    2. 'ducc' interface to ducc0.
  • nthreads: number of threads when calling into the 'ducc' backend. Ignored otherwise.

Output

f the signal on the sphere, 2D numpy.ndarray type complex or real.

pyssht.inverse_adjoint

flm = pyssht.inverse_adjoint(f, int L, Spin=0, Method='MW', Reality=False)

Performs the adjoint of the inverse spherical harmonic transform.

Inputs

  • f the signal on the sphere, numpy.ndarray type complex or real, ndim 2.
  • L the band limit of the signal, non-zero positive integer
  • Spin the spin of the signal, non-negative integer (default = 0)
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MWSS' [McEwen & Wiaux symmetric sampling]
  • Reality determines if the signal is real or complex, Boolean (default = False)
  • backend the backend that runs the transforms:
    1. 'SSHT' this package
    2. 'ducc' interface to ducc0.
  • nthreads: number of threads when calling into the 'ducc' backend. Ignored otherwise.

Output

flm the spherical harmonic transform of f, 1D numpy.ndarray type complex

pyssht.elm2ind

index = pyssht.elm2ind( int el, int m)

Computes the index in the flm array of a particular harmonic coefficient (\ell ) and (m).

Inputs

  • el the scale parameter of the spherical harmonic coefficients, integer from (0) to (L-1), where (L) is the band limit.
  • em the azimuthal parameter, integer from -el to el.

Output

Index of the coefficient in flm array, integer

pyssht.ind2elm

(el, em) = pyssht.ind2elm(int ind)

Computes harmonic coefficient (\ell ) and (m) from the index in the flm array.

Inputs

  • ind index of the flm array

Output

Tuple containing (el, em)

  • el the scale parameter of the spherical harmonic coefficients, integer from (0) to (L-1), where (L) is the band limit.
  • em the azimuthal parameter, integer from -el to el.

pyssht.theta_to_index

p = pyssht.theta_to_index(double theta, int L, str Method="MW")

Outputs the (\theta) index (the first) in the 2 dimensional array used to store spherical images. The index returned is that of the closest (\theta) sample smaller then the angle given on input.

Inputs

  • theta the angle (\theta)
  • L the band limit of the signal, non-zero positive integer
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MW_pole' [McEwen & Wiaux sampling with the south pole as a separate double.]
    3. 'MWSS' [McEwen & Wiaux symmetric sampling]
    4. 'DH' [Driscoll & Healy sampling]
    5. 'GL' [Gauss-Legendre sampling]

Output

Int p of corresponding to the angle (\theta).

pyssht.phi_to_index

q = pyssht.phi_to_index(double phi, int L, str Method="MW")

Outputs the (\phi) index (the second) in the 2 dimensional array used to store spherical images. The index returned is that of the closest (\phi) sample smaller then the angle given on input.

Inputs

  • phi the angle (\phi)
  • L the band limit of the signal, non-zero positive integer
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MW_pole' [McEwen & Wiaux sampling with the south pole as a separate double.]
    3. 'MWSS' [McEwen & Wiaux symmetric sampling]
    4. 'DH' [Driscoll & Healy sampling]
    5. 'GL' [Gauss-Legendre sampling]

Output

Int q of corresponding to the angle (\phi).

pyssht.sample_length

n = pyssht.sample_length(int L, Method = 'MW')

Outputs a size of the array used for storing the data on the sphere for different sampling schemes.

Inputs

  • L the band limit of the signal, non-zero positive integer
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MW_pole' [McEwen & Wiaux sampling with the south pole as a separate double.]
    3. 'MWSS' [McEwen & Wiaux symmetric sampling]
    4. 'DH' [Driscoll & Healy sampling]
    5. 'GL' [Gauss-Legendre sampling]

Output

Int n equal to the collapsed 1 dimensional size of the 2 dimensional array used to store data on the sphere.

pyssht.sample_shape

(n_theta, n_phi) = pyssht.sample_shape(int L, Method='MW')

Outputs a tuple with the shape of the array used for storing the data on the sphere for different sampling schemes.

Inputs

  • L the band limit of the signal, non-zero positive integer
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MW_pole' [McEwen & Wiaux sampling with the south pole as a separate double.]
    3. 'MWSS' [McEwen & Wiaux symmetric sampling]
    4. 'DH' [Driscoll & Healy sampling]
    5. 'GL' [Gauss-Legendre sampling]

Output

Tuple containing (n_theta, n_phi)

  • n_theta the number of samples in the (\theta) direction, integer
  • n_phi the number of samples in the (\phi) direction, integer

pyssht.sample_positions

(thetas, phis) = pyssht.sample_positions(int L, Method = 'MW', Grid=False)

Computes the positions on the sphere of the samples.

Inputs

  • L the band limit of the signal, non-zero positive integer
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MW_pole' [McEwen & Wiaux sampling with the south pole as a separate double.]
    3. 'MWSS' [McEwen & Wiaux symmetric sampling]
    4. 'DH' [Driscoll & Healy sampling]
    5. 'GL' [Gauss-Legendre sampling]
  • Grid describes if the output is a vector of the sample positions or a 2D array the same shape as the signal on the sphere, default False

Outputs

Tuple containing (thetas, phis)

  • thetas positions of the samples in the (\theta) direction
  • phis positions of the samples in the (\theta) direction

pyssht.s2_to_cart

(x, y, z) = pyssht.s2_to_cart(theta, phi)

Computes the (x), (y), and (z) coordinates from (\theta) and (\phi) on the sphere.

Inputs

  • theta (\theta) values, type numpy.ndarray
  • phi (\phi) values, type numpy.ndarray

Output

Tuple containing (x, y, z)

  • x the (x) coordinate of each point, type numpy.ndarray
  • y the (y) coordinate of each point, type numpy.ndarray
  • z the (z) coordinate of each point, type numpy.ndarray

pyssht.cart_to_s2

thetas, phis = cart_to_s2(x, y, z)

Computes the (\theta) and (\phi) on the coordinates on the sphere from (x), (y), and (z) coordinates.

Inputs

  • x (x) values, type numpy.ndarray
  • y (y) values, type numpy.ndarray
  • z (z) values, type numpy.ndarray

Output

Tuple containing (theta, phi)

  • theta the (\theta) coordinate of each point, type numpy.ndarray
  • phi the (\phi) coordinate of each point, type numpy.ndarray

pyssht.spherical_to_cart

(x, y, z) = pyssht.spherical_to_cart(r, theta, phi)

Computes the (x), (y), and (z) coordinates from the spherical coordinates (r), (\theta) and (\phi).

Inputs

  • r (r) values, type numpy.ndarray
  • theta (\theta) values, type numpy.ndarray
  • phi (\phi) values, type numpy.ndarray

Output

Tuple containing (x, y, z)

  • x the (x) coordinate of each point, type numpy.ndarray
  • y the (y) coordinate of each point, type numpy.ndarray
  • z the (z) coordinate of each point, type numpy.ndarray

pyssht.cart_to_spherical

r, theta, phi = pyssht.cart_to_spherical(x, y, z)

Computes the (\r), (\theta) and (\phi) on the spherical coordinates from (x), (y), and (z) coordinates.

Inputs

  • x (x) values, type numpy.ndarray
  • y (y) values, type numpy.ndarray
  • z (z) values, type numpy.ndarray

Output

Tuple containing (r, theta, phi)

  • r the (r) coordinate of each point, type numpy.ndarray
  • theta the (\theta) coordinate of each point, type numpy.ndarray
  • phi the (\phi) coordinate of each point, type numpy.ndarray

pyssht.theta_phi_to_ra_dec

(dec, ra) = pyssht.theta_phi_to_ra_dec(theta, phi, Degrees=False)

Computes the Right Assension and declination from an array of (\theta) and (\phi) values.

Inputs

  • theta (\theta) values, type numpy.ndarray
  • phi (\phi) values, type numpy.ndarray
  • Degrees defines if the output is in degrees or radians

Output

Tuple containing (dec, ra)

  • dec the declination angle, numpy.ndarray
  • ra the Right Assension, numpy.ndarray

pyssht.ra_dec_to_theta_phi

(theta, phi) = pyssht.ra_dec_to_theta_phi(ra, dec, Degrees=False)

Computes the (\theta) and (\phi) values from an array of Right Assension and declination values.

Inputs

  • dec the declination angle, numpy.ndarray
  • ra the Right Assension, numpy.ndarray
  • Degrees defines if the input is in degrees or radians, if degrees they are converted

Output

Tuple containing (theta, phi])

  • theta (\theta) values, type numpy.ndarray
  • phi (\phi) values, type numpy.ndarray

pyssht.make_rotation_matrix

rot_much = pyssht.make_rotation_matrix(list rot)

Computes the 3 by 3 rotation matrix from the Euler angles given on input

Inputs

  • rot List of length 3. Each element are the Euler angles [alpha, beta, gamma]

Output

3 by 3 rot_matrix the rotation matrix type ndarray dtype float

pyssht.rot_cart

x_p, y_p, z_p = pyssht.rot_cart(x, y, z, list rot)

Computes the rotations of the cartesian coordinates given a set of Euler angles. The inputs can be any shape ndarrays. For speed if the arrays are 1 or 2 dimensional it is recommended to use pyssht.rot_cart_1D or pyssht.rot_cart_2D.

Inputs

  • x (x) values, type numpy.ndarray
  • y (y) values, type numpy.ndarray
  • z (z) values, type numpy.ndarray
  • rot List of length 3. Each element are the Euler angles [alpha, beta, gamma]

Output

Tuple containing (x_p, y_p, z_p) the rotated coordinates the same shape and type as the inputs.

pyssht.rot_cart_1d and pyssht.rot_cart_2d

(x_p, y_p, z_p) = pyssht.rot_cart_1d(np.ndarray[np.float64_t, ndim=1] x, np.ndarray[np.float64_t, ndim=1] y, np.ndarray[np.float64_t, ndim=1] z, list rot)
(x_p, y_p, z_p) = pyssht.rot_cart_2d(np.ndarray[np.float64_t, ndim=2] x, np.ndarray[np.float64_t, ndim=2] y, np.ndarray[np.float64_t, ndim=2] z, list rot)

Computes the rotations of the cartesian coordinates given a set of Euler angles. The inputs can be any shape ndarrays. Same as pyssht.rot_cart except optimised for arrays that are 1 or 2 dimensional.

Inputs

  • x (x) values, type numpy.ndarray, dtype float, ndim 1 or 2
  • y (y) values, type numpy.ndarray, dtype float, ndim 1 or 2
  • z (z) values, type numpy.ndarray, dtype float, ndim 1 or 2
  • rot List of length 3. Each element are the Euler angles [alpha, beta, gamma]

Output

Tuple containing (x_p, y_p, z_p) the rotated coordinates the same shape and type as the inputs.

pyssht.plot_sphere

pyssht.plot_sphere(
    f, L, Method='MW', Close=True, Parametric=False,
    Parametric_Scaling=[0.0,0.5], Output_File=None,
    Show=True, Color_Bar=True, Units=None, Color_Range=None,
    Axis=True
)

Plots data on to a sphere. It is really slow and not very good!

Inputs

  • f the signal on the sphere, numpy.ndarray type complex or real, ndim 2. NB different for 'MW_pole' sampling.
  • L the band limit of the signal, non-zero positive integer
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MW_pole' [McEwen & Wiaux sampling with the south pole as a separate double.]
    3. 'MWSS' [McEwen & Wiaux symmetric sampling]
    4. 'DH' [Driscoll & Healy sampling]
    5. 'GL' [Gauss-Legendre sampling]
  • Close if true the full sphere is plotted (without a gap after the last (\phi) position), default True
  • Parametric the radius of the object at a certain point is defined by the function (not just the color), default False
  • Parametric_Saling used if Parametric=True, defines the radius of the shape at a particular angle r = norm(f)*Parametric_Saling[1] + Parametric_Saling[0], default [0.0,0.5]
  • Output_File if set saves the plot to a file of that name
  • Show if True shows you the plot, default False
  • Color_Bar if True shows the color bar, default True
  • Units is set puts a label on the color bar
  • Color_Range if set saturates the color bar in that range, else the function min and max is used
  • Axis if True shows the 3d axis, default True

Outputs

None

pyssht.mollweide_projection

f_plot, mask = pyssht.mollweide_projection(
    f, L, resolution=500, rot=None,
    zoom_region=[np.sqrt(2.0)*2,np.sqrt(2.0)],
    Method="MW"
)

Creates an ndarray of the mollweide projection of a spherical image and a mask array. This is useful for plotting results, not to be used for analysis on the plane. Elements in the signal f that are NaNs are marked in the mask. This allows one to plot these regions the color of their choice.

Here is an example of using the function to plot real spherical data.

f_plot, mask = pyssht.mollweide_projection(f, L, Method="MW") # make projection
plt.figure() # start figure
imgplot = plt.imshow(f_real_plot,interpolation='nearest') # plot the projected image
plt.colorbar(imgplot,fraction=0.025, pad=0.04) # plot color bar (these extra keywords make the bar a reasonable size)
plt.imshow(mask_real, interpolation='nearest', cmap=cm.gray, vmin=-1., vmax=1.) # plot the NaN regions in grey
plt.gca().set_aspect("equal") # ensures the region is the correct proportions
plt.axis('off') # removes axis (looks better)
plt.show()

Inputs

  • f the signal on the sphere, numpy.ndarray type complex or real, ndim 2.
  • L the band limit of the signal, non-zero positive integer
  • resolution size of the projected image, default 500
  • rot If the image should be rotated before projecting, None. rot should be a list of length 1 or 3. If 1 then the image is rotated around the (z) axis by that amount. If 3 then the image is rotated by the Euler angles given in the list.
  • zoom_region the region of the sphere to be plotted, default [np.sqrt(2.0)*2,np.sqrt(2.0)] is the full sphere.
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MWSS' [McEwen & Wiaux symmetric sampling]
    3. 'DH' [Driscoll & Healy sampling]
    4. 'GL' [Gauss-Legendre sampling]

Outputs

If the input is real:

Tuple containing:

  • f_plot the projection of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask the projection of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.

If the input is complex:

Tuple containing:

  • f_plot_real the projection of the real part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_real the projection of the masked regions (NaNs in input f.real) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.
  • f_plot_imag the projection of the imaginary part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_imag the projection of the masked regions (NaNs in input f.imag) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.

pyssht.equatorial_projection

f_proj_real, mask_real, (f_proj_imag, mask_imag)\
            = pyssht.equatorial_projection(f, int L, int resolution=500,\
             rot=None, list zoom_region=[-1,-1], str Method="MW", \
             str Projection="MERCATOR", int Spin=0)

Creates ndarrays of the projections of a spherical image and a mask array. This is useful for plotting results and performing analysis on the plane. All the spherical samples that fall in one planar pixel is averaged, if no samples fall in a pixel then the pixel is assigned the value of the closest spherical sample. Elements in the signal f that are NaNs are marked in the mask. This allows one to plot these regions the color of their choice.

There are two projections supported.

  1. The Mercator projection often used in maps.
  2. The Sinusoidal projection a simple equal area projection.

Here is an example of using the function to plot real spherical data using the Mercator projection.

f_proj, mask \
        = pyssht.equatorial_projection(f, L, resolution=500, Method="MW", \
        Projection="MERCATOR")
plt.figure() # start figure
imgplot = plt.imshow(f_proj,interpolation='nearest')# plot the projected image (north part)
plt.colorbar(imgplot) # plot color bar
plt.imshow(mask, interpolation='nearest', cmap=cm.gray, vmin=-1., vmax=1.) # plot the NaN regions in grey
plt.axis('off') # removes axis

plt.show()

Inputs

  • f the signal on the sphere, numpy.ndarray type complex or real, ndim 2.
  • L the band limit of the signal, non-zero positive integer
  • resolution size of the projected image, default 500
  • rot If the image should be rotated before projecting, None. rot should be a list of length 1 or 3. If 1 then the image is rotated around the (z) axis by that amount. If 3 then the image is rotated by the Euler angles given in the list.
  • zoom_region the region of the sphere to be plotted in radians. The first element is the angle left and right of the centre, default is np.pi for both projections. The second element is up and down of the equator, default is np.pi/2 for the Sinusoidal projection and 7*np.pi/16 for the Mercator projection.
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MWSS' [McEwen & Wiaux symmetric sampling]
    3. 'DH' [Driscoll & Healy sampling]
    4. 'GL' [Gauss-Legendre sampling]
  • Projection string describing which of the projections to use. Use "MERCATOR" for the Mercator projection and "SINE" for the Sinusoidal projection, default is "MERCATOR"
  • Spin the spin of the signal. If the signal has non-zero spin then on projection the signal must be rotated to account for the changing direction of the definition of the signal. By setting this to a non-zero integer will ensure this rotation is performed.

Outputs

If the input is real:

Tuple containing:

  • f_plot the projection of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask the projection of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.

If the input is complex:

Tuple containing:

  • f_plot_real the projection of the real part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_real the projection of the real part of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.
  • f_plot_imag the projection of the imaginary part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_imag the projection of the imaginary part of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.

pyssht.polar_projection

f_proj_north_real, mask_north_real, f_proj_south_real, mask_south_real,\
(f_proj_north_imag, mask_north_imag, f_proj_south_imag, mask_south_imag\ )
        = pyssht.polar_projection(f, int L, int resolution=500, rot=None,\
        float zoom_region=-1, str Method="MW", str Projection="OP",int Spin=0):

Creates an two ndarrays of the polar projection of a spherical image and a mask array. This is useful for plotting results and performing analysis on the plane. All the spherical samples that fall in one planar pixel is averaged, if no samples fall in a pixel then the pixel is assigned the value of the closest spherical sample. Elements in the signal f that are NaNs are marked in the mask. This allows one to plot these regions the color of their choice.

All the projections are centred around a pole. There are three projections supported.

  1. The Gnomic projection, defined by drawing a line from the centre of the circle trough the sphere on to the plane.
  2. The Stereographic projection, defined by drawing a line starting at the opposite pole through the sphere to the plane.
  3. The Orthographic, defined by a vertical projection on the plane.

Here is an example of using the function to plot real spherical data.

f_proj_north, mask_north, f_proj_south, mask_south \
        = pyssht.polar_projection(f, L, resolution=500, Method="MW", \
        Projection="OP")
plt.figure() # start figure
imgplot = plt.imshow(f_proj_north,interpolation='nearest')# plot the projected image (north part)
plt.colorbar(imgplot) # plot color bar
plt.imshow(mask_north, interpolation='nearest', cmap=cm.gray, vmin=-1., vmax=1.) # plot the NaN regions in grey
plt.axis('off') # removes axis (looks better)

plt.figure()
imgplot = plt.imshow(f_proj_south,interpolation='nearest')# plot the projected image (south part)
plt.colorbar(imgplot)
plt.imshow(mask_south, interpolation='nearest', cmap=cm.gray, vmin=-1., vmax=1.)
plt.title("orthographic projection south")
plt.axis('off')
plt.show()

Inputs

  • f the signal on the sphere, numpy.ndarray type complex or real, ndim 2.
  • L the band limit of the signal, non-zero positive integer
  • resolution size of the projected image, default 500
  • rot If the image should be rotated before projecting, default None. rot should be a list of length 1 or 3. If 1 then the image is rotated around the (z) axis by that amount. If 3 then the image is rotated by the Euler angles given in the list.
  • zoom_region the region of the sphere to be plotted in radians, default np.pi/2 is the full half sphere for the orthographic and stereographic projections and np.pi/4 for the gnomic projection as the equator is at infinity in this projection.
  • Method the sampling scheme used, string:
    1. 'MW' [McEwen & Wiaux sampling (default)]
    2. 'MWSS' [McEwen & Wiaux symmetric sampling]
    3. 'DH' [Driscoll & Healy sampling]
    4. 'GL' [Gauss-Legendre sampling]
  • Projection string describing which of the projections to use. Use "GP" for the Gnomic projection, "SP" for the Stereographic projection and "OP" for the Orthographic projection, default is "OP"
  • Spin the spin of the signal. If the signal has non-zero spin then on projection the signal must be rotated to account for the changing direction of the definition of the signal. By setting this to a non-zero integer will ensure this rotation is performed.

Outputs

If the input is real:

Tuple containing:

  • f_plot_north the projection of the north part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_north the projection of the north part of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.
  • f_plot_south the projection of the south part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_south the projection of the south part of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.

If the input is complex:

Tuple containing:

  • f_plot_north_real the projection of the north part of the real part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_north_real the projection of the north part of the real part of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.
  • f_plot_south_real the projection of the south part of the real part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_south_real the projection of the south part of the real part of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.
  • f_plot_north_imag the projection of the north part of the imaginary part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_north_imag the projection of the north part of the imaginary part of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.
  • f_plot_south_imag the projection of the south part of the imaginary part of the image as a 2 dimensional ndarray of type float. Masked regions and regions not in the sphere are NaNs to make them clear when plotted
  • mask_south_imag the projection of the south part of the imaginary part of the masked regions (NaNs in input f) as a 2 dimensional ndarray of type float. Masked regions have a value 0.0 and regions not in the sphere are NaNs to make them clear when plotted.

pyssht.dl_beta_recurse

dl = pyssht.dl_beta_recurse(np.ndarray[ double, ndim=2, mode="c"] dl not None,\
            double beta, int L, int el, \
            np.ndarray[ double, ndim=1, mode="c"] sqrt_tbl not None,\
            np.ndarray[ double, ndim=1, mode="c"] signs not None)

Compute the el-th plane of the Wigner small-d functions (from the (el-1)-th plane) using Risbo's method.

Inputs

  • dl the Wigner plane for all m and n, indexed dl[m][n] of size (2L-1)(2*L-1)
  • beta angle to calculate Wigner D matrix at, type double
  • L the band limit of the signal, non-zero positive integer
  • el is the current harmonic degree (i.e. dl input should already be computed for el-1, and dl output will be computed for el)
  • sqrt_tbl precomputed square-roots from (0) to (2*(L-1)+1)
  • signs precomputed ((-1)^m) signs from (m=0) to (L)

Outputs

Numpy ndarray dl, type float_ of the Wigner plane for all m and n

pyssht.dln_beta_recurse

dl = pyssht.dln_beta_recurse(np.ndarray[ double, ndim=1, mode="c"] dl not None,\
            np.ndarray[ double, ndim=1, mode="c"] dlm1 not None, double beta,\
            int L, int el, int n, np.ndarray[ double, ndim=1, mode="c"] sqrt_tbl not None,\
            np.ndarray[ double, ndim=1, mode="c"] signs not None)

Compute the el-th line of the Wigner small-d functions for given n (from the (el-1)-th and (el-2)-th lines) using 3-term recursion of Kostelec.

Inputs

  • dl the Wigner line for el for non-negative m and given n of size L
  • dlm1 is the line for el-1 and dlp1 is the line computed for el+1
  • beta angle to calculate Wigner D matrix at, type double
  • L the band limit of the signal, non-zero positive integer
  • el el is the current harmonic degree
  • n the third index in Wigner D matrices
  • sqrt_tbl precomputed square-roots from (0) to (2*(L-1)+1)
  • signs precomputed ((-1)^m) signs from (m=0) to (L)

Outputs

Numpy ndarray dl, type float_ the Wigner line for el for non-negative m and given n of size

pyssht.generate_dl

dl_array = pyssht.generate_dl(double beta, int L)

Generates the small Wigner D matrices up to a given band limit for a given (\beta)

Inputs

  • beta angle to calculate Wigner D matrix at, type double
  • L the band limit of the signal, non-zero positive integer

Outputs

Numpy ndarray dl_array, type float_ of the small Wigner D matrices

pyssht.rotate_flms

flm_rotated = pyssht.rotate_flms(
                np.ndarray[ double complex, ndim=1, mode="c"] f_lm not None,\
                double alpha, double beta, double gamma, int L, dl_array=None,\
                M=None, Axisymmetric=False, Keep_dl=False)

Function to rotate a set of spherical harmonic coefficients by the set of Euler angles (\alpha, \beta, \gamma ) using the (z,y,z) convention.

Inputs

  • flm the spherical harmonic transform of f, numpy.ndarray type complex, ndim 1
  • alpha rotation angle (\alpha), type double
  • beta rotation angle (\beta), type double
  • gamma rotation angle (\gamma), type double
  • L the band limit of the signal, non-zero positive integer
  • dl_array if set should be the precomputed small Wigner D matrix for angle (\beta) and harmonic band limit L. If not set this is calculated in the function. (This parameter is ignored when using the ducc backend.)
  • M if set is the azimuthal band limit of the function to be rotated, default M=L.
  • Axisymmetric set if the function is axisymmetric and axisymmetric harmonic coefficients are parsed.
  • Keep_dl if set the output is changed to allow one to keep the computed dl_array. (This parameter is ignored when using the ducc backend.)
  • backend the backend that runs the transforms:
    1. 'SSHT' this package
    2. 'ducc' interface to ducc0

Output

If Keep_dl is not set the output is the rotated set of spherical harmonic coefficients. If it is the output is a tuple (flm_rotated, dl_array), ie the rotated harmonic coefficients and the small Wigner D matrix computed for that band limit and (\alpha) value.

pyssht.guassian_smoothing

fs_lm = pyssht.guassian_smoothing(np.ndarray[ double complex, ndim=1, mode="c"] f_lm not None, int L, sigma_in=None, bl_in = None)

Smooths a set of harmonic coefficients either with a precomputed smoothing kernel bl or with a Gaussian given on input.

Inputs

  • f_lm the spherical harmonic transform of f, numpy.ndarray type complex, ndim 1
  • L the band limit of the signal, non-zero positive integer
  • sigma_in the input sigma of the Gaussian to smooth the signal with, default None
  • bl_in the smoothing kernel to smooth the signal with, default None

Output

fs_lm the smoothed harmonic coefficients.

pyssht.create_ylm

ylm = pyssht.create_ylm(thetas, phis, int L, int Spin=0, str recursion='Kostelec')

Computes spherical harmonic functions for all el and all 0<=|m|<= el using various recursions.

Inputs

  • thetas positions of the samples in the (\theta) direction
  • phis positions of the samples in the (\phi) direction
  • L the band limit of the signal, non-zero positive integer
  • Spin the spin of the signal, non-negative integer (default = 0)
  • recursion the recursion scheme used, string:
    1. 'Kostelec' [3-term recursion, e.g. Kostelec (default)]
    2. 'Risbo' [Risbo recursion]
    3. 'NumericalRecipes' [Numerical Recipes]

Output

ylm the spherical harmonics indexed ylm[ind][theta][phi], where ind = pyssht.elm2ind(el, m)

Project details


Download files

Download the file for your platform. If you're not sure which to choose, learn more about installing packages.

Source Distribution

pyssht-1.5.3.tar.gz (2.4 MB view details)

Uploaded Source

Built Distributions

pyssht-1.5.3-cp313-cp313-manylinux_2_17_x86_64.manylinux2014_x86_64.whl (1.3 MB view details)

Uploaded CPython 3.13manylinux: glibc 2.17+ x86-64

pyssht-1.5.3-cp313-cp313-manylinux_2_17_aarch64.manylinux2014_aarch64.whl (1.2 MB view details)

Uploaded CPython 3.13manylinux: glibc 2.17+ ARM64

pyssht-1.5.3-cp313-cp313-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl (1.1 MB view details)

Uploaded CPython 3.13manylinux: glibc 2.17+ i686manylinux: glibc 2.5+ i686

pyssht-1.5.3-cp313-cp313-macosx_11_0_arm64.whl (959.5 kB view details)

Uploaded CPython 3.13macOS 11.0+ ARM64

pyssht-1.5.3-cp313-cp313-macosx_10_13_x86_64.whl (1.1 MB view details)

Uploaded CPython 3.13macOS 10.13+ x86-64

pyssht-1.5.3-cp312-cp312-manylinux_2_17_x86_64.manylinux2014_x86_64.whl (1.2 MB view details)

Uploaded CPython 3.12manylinux: glibc 2.17+ x86-64

pyssht-1.5.3-cp312-cp312-manylinux_2_17_aarch64.manylinux2014_aarch64.whl (1.2 MB view details)

Uploaded CPython 3.12manylinux: glibc 2.17+ ARM64

pyssht-1.5.3-cp312-cp312-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl (1.1 MB view details)

Uploaded CPython 3.12manylinux: glibc 2.17+ i686manylinux: glibc 2.5+ i686

pyssht-1.5.3-cp312-cp312-macosx_11_0_arm64.whl (964.2 kB view details)

Uploaded CPython 3.12macOS 11.0+ ARM64

pyssht-1.5.3-cp312-cp312-macosx_10_13_x86_64.whl (1.1 MB view details)

Uploaded CPython 3.12macOS 10.13+ x86-64

pyssht-1.5.3-cp311-cp311-manylinux_2_17_x86_64.manylinux2014_x86_64.whl (1.3 MB view details)

Uploaded CPython 3.11manylinux: glibc 2.17+ x86-64

pyssht-1.5.3-cp311-cp311-manylinux_2_17_aarch64.manylinux2014_aarch64.whl (1.2 MB view details)

Uploaded CPython 3.11manylinux: glibc 2.17+ ARM64

pyssht-1.5.3-cp311-cp311-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl (1.2 MB view details)

Uploaded CPython 3.11manylinux: glibc 2.17+ i686manylinux: glibc 2.5+ i686

pyssht-1.5.3-cp311-cp311-macosx_11_0_arm64.whl (967.9 kB view details)

Uploaded CPython 3.11macOS 11.0+ ARM64

pyssht-1.5.3-cp311-cp311-macosx_10_9_x86_64.whl (1.1 MB view details)

Uploaded CPython 3.11macOS 10.9+ x86-64

pyssht-1.5.3-cp310-cp310-manylinux_2_17_x86_64.manylinux2014_x86_64.whl (1.3 MB view details)

Uploaded CPython 3.10manylinux: glibc 2.17+ x86-64

pyssht-1.5.3-cp310-cp310-manylinux_2_17_aarch64.manylinux2014_aarch64.whl (1.2 MB view details)

Uploaded CPython 3.10manylinux: glibc 2.17+ ARM64

pyssht-1.5.3-cp310-cp310-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl (1.2 MB view details)

Uploaded CPython 3.10manylinux: glibc 2.17+ i686manylinux: glibc 2.5+ i686

pyssht-1.5.3-cp310-cp310-macosx_11_0_arm64.whl (966.9 kB view details)

Uploaded CPython 3.10macOS 11.0+ ARM64

pyssht-1.5.3-cp310-cp310-macosx_10_9_x86_64.whl (1.1 MB view details)

Uploaded CPython 3.10macOS 10.9+ x86-64

pyssht-1.5.3-cp39-cp39-manylinux_2_17_x86_64.manylinux2014_x86_64.whl (1.3 MB view details)

Uploaded CPython 3.9manylinux: glibc 2.17+ x86-64

pyssht-1.5.3-cp39-cp39-manylinux_2_17_aarch64.manylinux2014_aarch64.whl (1.2 MB view details)

Uploaded CPython 3.9manylinux: glibc 2.17+ ARM64

pyssht-1.5.3-cp39-cp39-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl (1.2 MB view details)

Uploaded CPython 3.9manylinux: glibc 2.17+ i686manylinux: glibc 2.5+ i686

pyssht-1.5.3-cp39-cp39-macosx_11_0_arm64.whl (967.2 kB view details)

Uploaded CPython 3.9macOS 11.0+ ARM64

pyssht-1.5.3-cp39-cp39-macosx_10_9_x86_64.whl (1.1 MB view details)

Uploaded CPython 3.9macOS 10.9+ x86-64

File details

Details for the file pyssht-1.5.3.tar.gz.

File metadata

  • Download URL: pyssht-1.5.3.tar.gz
  • Upload date:
  • Size: 2.4 MB
  • Tags: Source
  • Uploaded using Trusted Publishing? Yes
  • Uploaded via: uv/0.6.5

File hashes

Hashes for pyssht-1.5.3.tar.gz
Algorithm Hash digest
SHA256 3690c74186db05efb2b12006a510ff8bb71cc07f42946d1cba1ce8e998a7b130
MD5 052b35594673ae7e7c19b46f65dbfc04
BLAKE2b-256 b00103837219c8fca0e045a921be2c1ed7da30c8d16bc7b22bd6cae97fff3f5d

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp313-cp313-manylinux_2_17_x86_64.manylinux2014_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp313-cp313-manylinux_2_17_x86_64.manylinux2014_x86_64.whl
Algorithm Hash digest
SHA256 5f1c8cd45158b68760d97e6b97120cedf2e03a3be4fb9e4ad1739af838e8253e
MD5 c4ed34d7fd625e1c38ccc9ab8bcf93f7
BLAKE2b-256 7307cf30fe1473a80c26d3aeac951142c50212f1a3faf65391e88516ed2f00ea

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp313-cp313-manylinux_2_17_aarch64.manylinux2014_aarch64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp313-cp313-manylinux_2_17_aarch64.manylinux2014_aarch64.whl
Algorithm Hash digest
SHA256 4acca572b0c2189962763adeac93c817b99d69394f6b05e4e12f894d3f642351
MD5 c7f225f5bf94e37663b6cca57cb6af87
BLAKE2b-256 5220d9daf6d2346745f1c47c46a354b15cb02a398eaaa34c125fb3958baaed98

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp313-cp313-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp313-cp313-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl
Algorithm Hash digest
SHA256 5b3124908a3fbaffc943150289c46397a6567f236f917f816ded1efeddf138c4
MD5 4ca9009325d0fa291065ae09ba584bc5
BLAKE2b-256 36008a9fc6081bb7e7f6096c5c61e9b1b7efa4546c52b00fe9408bb9928b8068

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp313-cp313-macosx_11_0_arm64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp313-cp313-macosx_11_0_arm64.whl
Algorithm Hash digest
SHA256 c571a3eb8b9e4e3f1e9cacbf0a2daff71a7439cace56ec5bd64cb2a79eadfbef
MD5 559c64f994edf75cbedda919cf973582
BLAKE2b-256 7c26f09322744916812948da5c48def4cdfdf1d2f21587aa71f0e5e232c2647b

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp313-cp313-macosx_10_13_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp313-cp313-macosx_10_13_x86_64.whl
Algorithm Hash digest
SHA256 148e9c952d2d50a34cc22bc6f9ba76061feab53d2510a6e75d5bd9a2cf64f958
MD5 64a45c9e45643e6051b860028e0fdc9b
BLAKE2b-256 02c4bb8ec6274bdbc70e8c03d7b6c8ccd65555da2d65e5a17b11aa766629ea90

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp312-cp312-manylinux_2_17_x86_64.manylinux2014_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp312-cp312-manylinux_2_17_x86_64.manylinux2014_x86_64.whl
Algorithm Hash digest
SHA256 e5883b4cb301ca8323b5f2241b0e5017adbc1e44f2987cefbdfc368b32172119
MD5 2ce74cc0c4fe73984558a3fc976f8ddd
BLAKE2b-256 8d121eb646dcdb9e36b91ca0e54e2f01326f7c74432d5a955bd3fd776a5829e4

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp312-cp312-manylinux_2_17_aarch64.manylinux2014_aarch64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp312-cp312-manylinux_2_17_aarch64.manylinux2014_aarch64.whl
Algorithm Hash digest
SHA256 8a47a9e586f875057bbef4d08b3a4d7e753477cb8fbf63cf3449c585dd3d9646
MD5 9f287702b87f5b821611e6c69407ec32
BLAKE2b-256 bb44332d3ccd6e0e93028cb254d1f696f5bd82df932ab20e3daee3ed5c6f9a2d

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp312-cp312-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp312-cp312-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl
Algorithm Hash digest
SHA256 4a05f4da404035dfb3be1f7f62f0ca28d34638965147cf86a48bd3631358101f
MD5 dc6707c7ab9a5e3b2f6f9560065ad193
BLAKE2b-256 7063db1a90b5542464f327ad33332d94c113cffb7c863ad3650f1eeb53e08844

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp312-cp312-macosx_11_0_arm64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp312-cp312-macosx_11_0_arm64.whl
Algorithm Hash digest
SHA256 a9d25f529ca411936a2f994b2f62bcbc6671015b17612730cbbc6313ff9830d2
MD5 85394533a02359b6a6f72ed82ed1d748
BLAKE2b-256 d417e3573f44af6e33131444f31d5d85e992ccba7ffd96789bbe860c9a45f55c

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp312-cp312-macosx_10_13_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp312-cp312-macosx_10_13_x86_64.whl
Algorithm Hash digest
SHA256 b0b52b7001efb7479c2352baf26fc3977dabbe9e022da064f70be468298704a8
MD5 b4d217cdd1b68b0652298638380f30c1
BLAKE2b-256 408daf3a6bcd5e21e4484b557b3db3fe014e213167cc584b5b0016563cd39740

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp311-cp311-manylinux_2_17_x86_64.manylinux2014_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp311-cp311-manylinux_2_17_x86_64.manylinux2014_x86_64.whl
Algorithm Hash digest
SHA256 a4ebab8c12d03c9e877e614986c0390dc915b88cbd2baa9afb066ff5a465f41a
MD5 5e83c8dfb02810cbe4055f0d2b294962
BLAKE2b-256 7efc85f7ebed76589e14697b0fe4e7cfac6388ae631bc99dae0d75cfa75ee7b1

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp311-cp311-manylinux_2_17_aarch64.manylinux2014_aarch64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp311-cp311-manylinux_2_17_aarch64.manylinux2014_aarch64.whl
Algorithm Hash digest
SHA256 77c70c1682c70cb9d4f07a75d0cf0baf9ba60a7167fdece338f660acf067c23b
MD5 ba52c6d07ed4baf184b18b55520baed7
BLAKE2b-256 3e8f47eb3f68513caa844451e8873153ea9b23e69391c6fe021e36d8111ef8f3

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp311-cp311-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp311-cp311-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl
Algorithm Hash digest
SHA256 31d4f66b22041d324e0b710e51613aa1888f9371ea613c73f03615194c8b4942
MD5 cb47802b59ff6ca84256891ac53a2814
BLAKE2b-256 a928e77fef40c9ba0849b45d90e2ebc2ddaadbc5ec7caa455a5fff1f490ac822

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp311-cp311-macosx_11_0_arm64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp311-cp311-macosx_11_0_arm64.whl
Algorithm Hash digest
SHA256 777e7d2d69a2fd20e6b12a7bd8f3c54e23a97d7720448b8e504796867b511a02
MD5 09709ed82bd054fab608bafcf65dc9a8
BLAKE2b-256 39cdac2d001877f66c0d08fcae05acd4162545544cf41a04f0f69f26dc275be7

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp311-cp311-macosx_10_9_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp311-cp311-macosx_10_9_x86_64.whl
Algorithm Hash digest
SHA256 c5f6b7b2019894599d5638334d0736d54b9ba5f1058a6646ebb7a9dbd1a99940
MD5 ad4da9bbe5a6d49de5bc8136b3dd895b
BLAKE2b-256 b0028db3ae81e0444d7f164a3aa370820e2c8d012ad3f375ce29f87a98c228cf

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp310-cp310-manylinux_2_17_x86_64.manylinux2014_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp310-cp310-manylinux_2_17_x86_64.manylinux2014_x86_64.whl
Algorithm Hash digest
SHA256 4f39dbdd335f52d057c090282d9a096c45a978827f662a70cdb1c902d585daca
MD5 c19687e879f793fc77e44a355163636a
BLAKE2b-256 4165e42e84a81ba2ea44bda75d76bbaab49a6f998b8783b89396ab4e8fd67dcf

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp310-cp310-manylinux_2_17_aarch64.manylinux2014_aarch64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp310-cp310-manylinux_2_17_aarch64.manylinux2014_aarch64.whl
Algorithm Hash digest
SHA256 017d22c9b96a2503020bdfaa8579766d671080bee143def809cfa15bc15ca4f1
MD5 f86404e241e4ef5028817771d5755dec
BLAKE2b-256 b26038323291fed7658de705f77b1d2ddf206d8983ff53f8e394ce56ea2fd40c

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp310-cp310-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp310-cp310-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl
Algorithm Hash digest
SHA256 c3e84e6b2f1a540f58e6537e50343218e7b7c403f228e2ef03f359aa35603109
MD5 44bc986b97f7be95d24e8a52606b3be1
BLAKE2b-256 bae895bcb825b1f73972eda7faa90a4f29d3f882a8776a9488a3333155f34c45

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp310-cp310-macosx_11_0_arm64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp310-cp310-macosx_11_0_arm64.whl
Algorithm Hash digest
SHA256 ee1a0a838eb3a4b8481981bf02a8d1c1f6338a8a3e034438da670d88b654c09f
MD5 40ab91020a63055487feedc805a48e4f
BLAKE2b-256 127c3f7b2818a485cde16c532fb8f403eae2160760040efd3ff1c7e6b236f782

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp310-cp310-macosx_10_9_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp310-cp310-macosx_10_9_x86_64.whl
Algorithm Hash digest
SHA256 808804d85b2fa889413e92008d0f460069841a61ce5b97cd70661cc03c216fb5
MD5 a6eab36ea74fbca23431d3814d648f62
BLAKE2b-256 6ef18e199bc55cce6da46a90d98eb4021fdced90129c136bc06bcc0f78184b44

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp39-cp39-manylinux_2_17_x86_64.manylinux2014_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp39-cp39-manylinux_2_17_x86_64.manylinux2014_x86_64.whl
Algorithm Hash digest
SHA256 92a78b413995ffd92a30d75f2a3b4e38931785ecf699ab49f5ed71350ba84a38
MD5 982c0e013a82d7d976cd4a613f8cae24
BLAKE2b-256 2d1d14ab725d8427550016f1ae1a53de34d6dfd4770f259f4f8cbc6738c6361d

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp39-cp39-manylinux_2_17_aarch64.manylinux2014_aarch64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp39-cp39-manylinux_2_17_aarch64.manylinux2014_aarch64.whl
Algorithm Hash digest
SHA256 036dc3f6fce165dbf2f8f5fb0ba4fc0b87f818cc9be822e02bf4af0545b3486f
MD5 9702e531990ae1565126c3baf3a515c1
BLAKE2b-256 49e94a2496ab1c009039364bbe34bbb263dcf4dba853cff46d90f672b0eeebc5

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp39-cp39-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp39-cp39-manylinux_2_5_i686.manylinux1_i686.manylinux_2_17_i686.manylinux2014_i686.whl
Algorithm Hash digest
SHA256 8a7c728a28179f7e0ea89e2063988513485bcaa2c9d07b9cba3891d083ef5ec3
MD5 fd3829e83921ac52361d4f820f2cadc5
BLAKE2b-256 cc42fdd955de43f6083b12409f386189e50aaf3e79ec7ac21e42800c4d97ae1c

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp39-cp39-macosx_11_0_arm64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp39-cp39-macosx_11_0_arm64.whl
Algorithm Hash digest
SHA256 74a39e77959ed3bb23e5fc667f00cb54fd06a4d1850179722d608e35a4f937f9
MD5 8bf0d654bb2132edda3d86848af6d5d3
BLAKE2b-256 f689d7f5988c0821242dbbf48046806a1e992032b12da8e827ce5874df077cfc

See more details on using hashes here.

File details

Details for the file pyssht-1.5.3-cp39-cp39-macosx_10_9_x86_64.whl.

File metadata

File hashes

Hashes for pyssht-1.5.3-cp39-cp39-macosx_10_9_x86_64.whl
Algorithm Hash digest
SHA256 4b0fdc45627db024fbe1f484a89c2a5f65c4e5b242c3b58439cb61471025d684
MD5 f98f5d495bb8188c56e5c635d0759e97
BLAKE2b-256 593f3dfa8e0067ba0de71e3ea70eda3523d5726a57cf9df060fb4d3e44748f5a

See more details on using hashes here.

Supported by

AWS Cloud computing and Security Sponsor Datadog Monitoring Fastly CDN Google Download Analytics Pingdom Monitoring Sentry Error logging StatusPage Status page