avni.tools.harmonics module#
- avni.tools.harmonics.getdepthsfolder(folder: str = '.', extension: str = '.epix', delimiter: str = '.') list [source]#
Get list of depths from filenames to iterate through
- Parameters
- folderstr, optional
Folder to search, by default ‘.’
- extensionstr, optional
File extension to search, by default ‘.epix’
- delimiterstr, optional
delimiter in file names, by default ‘.’
- Returns
- list
List of depths
- Authors
Raj Moulik (moulik@caa.columbia.edu)
- Last Modified
2023.02.16 5.00
- avni.tools.harmonics.rdswpsh(filename: str) Tuple[numpy.ndarray, dict, list] [source]#
Code to get spherical harmonic coefficients in the ylm normalization.
- Parameters
- filenamestr
Input file with spherical harmonic coefficients
- Returns
- tp.Tuple[np.ndarray, dict, list]
First element is a numpy array with coefficients in the ylm normalization.
Second element are metadata from input fields if specified.
Third element are all other comments except lines containing metadata.
- Raises
- IOError
File not found in the file system
- Authors
Raj Moulik (moulik@caa.columbia.edu)
- Last Modified
2023.02.16 5.00
- avni.tools.harmonics.swp_to_epix(infile: str, grid: int = 1, lmax: Union[None, int] = None) Tuple[numpy.ndarray, dict, list] [source]#
Convert spherical harmonics coefficients to extended pixel (.epix) format.
- Parameters
- infilestr
Input spherical harmonic coefficient file
- gridint, optional
Grid size for pixels, by default 1
- lmaxtp.Union[None,int], optional
Maximum angular degree, by default None which results in the maximum specified on file
- Returns
- tp.Tuple[np.ndarray, dict, list]
First element is an array containing (latitude, longitude, pixel_size, value).
Second element are metadata from input fields if specified.
Third element are all other comments except lines containing metadata.
- Authors
Raj Moulik (moulik@caa.columbia.edu)
- Last Modified
2023.02.16 5.00
- avni.tools.harmonics.wrswpsh(filename: str, shmatrix: numpy.ndarray, metadata: dict = {'FORMAT': '0'}, comments: Optional[list] = None, lmax: Union[None, int] = None)[source]#
Write spherical harmonic coefficients from the ylm normalization to a file.
- Parameters
- filenamestr
Output file name
- shmatrixnp.ndarray
Numpy array with coefficients in the ylm normalization
- metadata_type_, optional
Metadata from input fields if specified., by default {‘FORMAT’:’0’}
- commentslist, optional
All other comments except lines containing metadata., by default None
- lmaxtp.Union[None,int], optional
Maximum angular degree, by default None which results in the maximum specified on file
- :Authors:
Raj Moulik (moulik@caa.columbia.edu)
- :Last Modified:
2023.02.16 5.00
- avni.tools.harmonics.get_coefficients(shmatrix: numpy.ndarray, lmin: Union[None, int] = None, lmax: Union[None, int] = None) numpy.ndarray [source]#
Read the spherical harmonic coefficients into a named numpy array
- Parameters
- shmatrixnp.ndarray
Numpy array with coefficients in the ylm normalization
- lmintp.Union[None,int], optional
Minimum angular degree, by default None which results in the minimum specified on file
- lmaxtp.Union[None,int], optional
Maximum angular degree, by default None which results in the maximum specified on file
- Returns
- ——-
- np.ndarray
Named numpy array of spherical harmonic coefficients
- :Authors:
Raj Moulik (moulik@caa.columbia.edu)
- :Last Modified:
2023.02.16 5.00
- avni.tools.harmonics.swp_to_xarray(shmatrix: numpy.ndarray, grid: int = 10, lmax: Union[None, int] = None) xarray.core.dataarray.DataArray [source]#
Convert from spherical harmonic coefficients in ylm normalization to a multi-dimensional pixel grid.
- Parameters
- shmatrixnp.ndarray
Numpy array with coefficients in the ylm normalization
- gridint, optional
Grid size in degrees for pixels, by default 10
- lmaxtp.Union[None,int], optional
Maximum angular degree, by default None which results in the maximum specified on file
- Returns
- xr.DataArray
A multi-dimensional xarray DataArray containing the pixel grid
- Authors
Raj Moulik (moulik@caa.columbia.edu)
- Last Modified
2023.02.16 5.00
- avni.tools.harmonics.convert_to_swp(data: Union[numpy.ndarray, xarray.core.dataarray.DataArray], lmax: int) numpy.ndarray [source]#
Convert multi-dimensional pixel grid to spherical harmonic coefficients in ylm normalization.
- Parameters
- datatp.Union[np.ndarray,xr.DataArray]
A multi-dimensional xarray DataArray or named numpy array containing the pixel grid
- lmaxint
Maximum angular degree to evaluate coefficients to
- Returns
- np.ndarray
Numpy array with coefficients in the ylm normalization
- Authors
Raj Moulik (moulik@caa.columbia.edu)
- Last Modified
2023.02.16 5.00
- avni.tools.harmonics.calcshpar2(shmatrix: numpy.ndarray, lmax: Union[None, int] = None) Tuple[float, float, float, numpy.ndarray] [source]#
This calculates the mean, roughness, RMS and power of spherical harmonic coefficients in ylm normalization
- Parameters
- shmatrixnp.ndarray
Numpy array with coefficients in the ylm normalization
- lmaxtp.Union[None,int], optional
Maximum angular degree, by default None which results in the maximum specified on file
- Returns
- tp.Tuple[float,float,float,np.ndarray]
First element the globally averaged value on the unit sphere.
Second and third elements are RMS and roughness values, respectively.
Fourth element is a numpy array containing the power at each degree.
- Authors
Raj Moulik (moulik@caa.columbia.edu)
- Last Modified
2023.02.16 5.00
- avni.tools.harmonics.swp_correlation(shmatrix1: numpy.ndarray, shmatrix2: numpy.ndarray, lmax: Union[None, int] = None) Tuple[float, float, numpy.ndarray, numpy.ndarray] [source]#
Calculate RMS and correlation between two sets of spherical harmonic coefficients.
- Parameters
- shmatrix1np.ndarray
First numpy array with coefficients in the ylm normalization
- shmatrix2np.ndarray
Second numpy array with coefficients in the ylm normalization
- lmaxtp.Union[None,int], optional
Maximum angular degree, by default None which results in the maximum specified on file
- Returns
- tp.Tuple[float,float,np.ndarray,np.ndarray]
First and second elements are the RMS values for the two sets of coefficients.
Third element is the numpy array containing the correlation at each spherical harmonic degree.
Fourth element is a numpy array containing cumulative correlation up to each spherical harmonic degree.
- Authors
Raj Moulik (moulik@caa.columbia.edu)
- Last Modified
2023.02.16 5.00