avni.tools.bases module#

avni.tools.bases.eval_vbspl(radius: Union[list, tuple, numpy.ndarray], knots: Union[str, list, tuple, float, numpy.int64, numpy.ndarray, bool])[source]#
Evaluate the cubic spline knot with second derivative as 0 at end points.

In contrast to eval_splrem(), this function distributes the spline knots unevenly at the specified depths or radii.

Parameters
radiustp.Union[list,tuple,np.ndarray]

A radius value or array of radii (or alternatively depths) queried in km

knotstp.Union[str,list,tuple,float,np.int64,np.ndarray,bool]

A numpy array or list of radii (or depths) of spline knots in km

Returns
vercof, dvercof: float or np.ndarray

value of the polynomial coefficients at each depth and derivative. Both arrays have size (Nradius, Nsplines).

Authors

Raj Moulik (moulik@caa.columbia.edu)

Last Modified

2023.02.16 5.00

avni.tools.bases.eval_splrem(radius: Union[list, tuple, numpy.ndarray], radius_range: Union[list, tuple, numpy.ndarray], nsplines: int)[source]#
Evaluate the cubic spline knot with second derivative as 0 at end points.

In contrast to eval_vbspl(), this function distributes the spline knots evenly within the radius range.

Parameters
radiustp.Union[list,tuple,np.ndarray]

A radius value or array of radii (or alternatively depths) queried in km

radius_rangetp.Union[list,tuple,np.ndarray]

Limits of the radius (or depths) limits of the region

nsplinesint

number of splines within the range

Returns
vercof, dvercof: float or np.ndarray

value of the polynomial coefficients at each depth and derivative. Both arrays have size (Nradius, Nsplines).

Authors

Raj Moulik (moulik@caa.columbia.edu)

Last Modified

2023.02.16 5.00

avni.tools.bases.eval_polynomial(radius: Union[list, tuple, numpy.ndarray], radius_range: Union[list, tuple, numpy.ndarray], rnorm: float, types: list = ['CONSTANT', 'LINEAR'])[source]#

Evaluate a set of polynomial functions within a range of radii.

Parameters
radiustp.Union[list,tuple,np.ndarray]

A radius value or array of radii queried in km

radius_rangetp.Union[list,tuple,np.ndarray]

Limits of the radius limits of the region

rnormfloat

normalization for radius, usually the radius of the planet

typeslist, optional

polynomial coefficients to be used for calculation. Options are : TOP,BOTTOM, CONSTANT, LINEAR, QUADRATIC, CUBIC, by default [‘CONSTANT’,’LINEAR’]

Returns
vercof, dvercof: float or np.ndarray

value of the polynomial coefficients and derivative at each depth size (Nradius)

Authors

Raj Moulik (moulik@caa.columbia.edu)

Last Modified

2023.02.16 5.00

avni.tools.bases.eval_splcon(latitude: Union[list, tuple, numpy.ndarray], longitude: Union[list, tuple, numpy.ndarray], xlaspl: numpy.ndarray, xlospl: numpy.ndarray, xraspl: numpy.ndarray)[source]#

Evaluate spherical splines at a set of locations.

Parameters
latitudetp.Union[list,tuple,np.ndarray]

Latitudes of locations queried

longitudetp.Union[list,tuple,np.ndarray]

Longitudes of locations queried

xlasplnp.ndarray

Latitude locations of splines

xlosplnp.ndarray

Longitude locations of splines

xrasplnp.ndarray

Radius of splines

Returns
horcofnp.ndarray

Value of the horizontal coefficents at each location. Size of numpy array is [len(latitude) X ncoefhor]

avni.tools.bases.splcon(lat: float, lon: float, ncoefhor: int, xlaspl: numpy.ndarray, xlospl: numpy.ndarray, xraspl: numpy.ndarray)[source]#

Evaluate spherical splines at a given location. Modified from the Fortran code splcon.f that is based on Wang and Dahlen (1995) [WD95].

Parameters
latfloat

latitude query

lonfloat

longitude query

ncoefhorint

number of splines

xlasplnp.ndarray

spline latitudes

xlosplnp.ndarray

spline longitudes

xrasplnp.ndarray

spline radii

Returns
ncon,icon,con

number of splines, spline index, and value at each location

Authors

Raj Moulik (moulik@caa.columbia.edu)

Last Modified

2023.02.16 5.00

avni.tools.bases.eval_ylm(latitude: Union[list, tuple, numpy.ndarray], longitude: Union[list, tuple, numpy.ndarray], lmaxhor: int, weights: Union[None, list, tuple, numpy.ndarray] = None, grid: bool = False, norm: str = 'ylm')[source]#

Evaluate spherical harmonics with a specific normalization.

Parameters
latitudetp.Union[list,tuple,np.ndarray]

Latitudes of locations queried

longitudetp.Union[list,tuple,np.ndarray]

Longitudes of locations queried

lmaxhorint

Maximum spherical harmonic degree

weightstp.Union[None,list,tuple,np.ndarray], optional

Weights to multiply the bases with i.e. coefficients, by default None

gridbool, optional

Create a grid based on a combination of latitudes and longitudes, by default False

normstr, optional

Spherical harmonic normalization, by default ‘ylm’

Returns
horcof

Value of the horizontal coefficents at each location. Size of numpy array is [len(latitude) X ((lmaxhor+1)^2)]

Authors

Raj Moulik (moulik@caa.columbia.edu)

Last Modified

2023.02.16 5.00

avni.tools.bases.eval_pixel(latitude: Union[list, tuple, numpy.ndarray], longitude: Union[list, tuple, numpy.ndarray], xlapix: Union[list, tuple, numpy.ndarray], xlopix: Union[list, tuple, numpy.ndarray], xsipix: Union[list, tuple, numpy.ndarray])[source]#

Evaluate pixel values at a set of locations.

Parameters
latitudetp.Union[list,tuple,np.ndarray]

Latitudes of locations queried

longitudetp.Union[list,tuple,np.ndarray]

Longitudes of locations queried

xlapixtp.Union[list,tuple,np.ndarray]

Pixel center latitudes

xlopixtp.Union[list,tuple,np.ndarray]

Pixel center longitudes

xsipixtp.Union[list,tuple,np.ndarray]

Pixel sizes in degrees

Returns
horcof

Value of the horizontal coefficents at each location. Size of numpy array is [len(latitude) X len(xsipix)]

Authors

Raj Moulik (moulik@caa.columbia.edu)

Last Modified

2023.02.16 5.00