Py_Admittance package
Submodules
Module contents
Py_Admittance
Py_Admittance provides functions and an example scripts for calculating global/localized admittances of gravity and topography. These can be used to invert for parameters such as the elastic thickness of the lithosphere and density of the crust and surface.
- ForwardGravity
Compute the theoretical gravity field given a set of input parameters.
- TransferTGz
Compute theoretical transfer functions for calculating the admittance and correlation given a set of input parameters.
- LocalAdmitCorr
Compute the localized admittance of two functions.
- GlobalAdmitCorr
Compute the global admittance of two functions.
- Py_Admittance.ForwardGravity(topo, G, mass, g0, Tc, Te, rhoc, rhom, lmax, ratio_L=0, alpha_L=1, depth_L=50000.0, E=100000000000.0, v=0.25, nmax=5, R_ref=None, rhol=None, lmaxgrid=None, option_deg1=None, option_deg2=None, return_corr=False)
Compute the theoretical gravity field given the input parameters. For more information, see Broquet & Wieczorek (2019).
- Returns:
grav_Th (array, dimension (2, lmax+1, lmax+1)) – Theoretical global gravity field in mGal.
Corr_Th (array, dimension (lmax+1)) – If return_corr is True, returns the theoretical correlation for surface/internal loads that are out of phase (alpha_L != 1).
- Parameters:
topo (array, dimension (2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the surface topography.
G (float) – Gravitational constant.
mass (float) – Mass of the planet.
g0 (float) – Gravitational attraction at the surface.
Tc (float) – Average crustal thickness.
Te (float) – Elastic thickness of the lithosphere.
rhoc (float) – Density of the crust.
rhom (float) – Density of the mantle.
lmax (int) – Maximum spherical harmonic degree for calculations.
ratio_L (float, optional, default = 0) – Ratio of the internal / surface load.
alpha_L (float, optional, default = 1) – Phase relationship for the internal / surface load. This parameter is experimental.
depth_L (float, optional, default = 50e3) – Depth of the internal load.
E (float, optional, default = 100e9) – Young’s modulus.
v (float, optional, default = 0.25) – Poisson’s ratio.
nmax (int, optional, default = 5) – Order of the finite-amplitude correction.
R_ref (float, optional, default = None) – Reference radius for gravity field calculations. If None, this parameter is set to the mean radius of the topography file.
rhol (float, optional, default = None) – Density of the surface topography. If None, this parameter is set to rhoc.
lmaxgrid (int, optional, default = None) – Resolution of the input grid for the finite-amplitude correction routines. If None, this parameter is set to 3*lmax. For accurate results, this parameter should be about 3 times lmax, though this should be verified for each application. Lowering this parameter significantly increases speed.
option_deg1 (string, optional, default = None) – How to treat degree-1 displacement. If set to “Zero”, the degree-1 displacement is zeroed out. If set to “Airy”, the degree-1 topography is assumed to be Airy compensated. If anything else, no special treatment is applied to degree-1.
option_deg2 (string/float, optional, default = None) – How to treat degree-2 topography. If set to “Zero”, the C20 topography is zeroed out. If set to “Flat”, the C20 topography is not considered as a load, but is added back for finite-amplitude calculations. If set to a float, only option_deg2 * topography is used as a load. If anything else, no special treatment is applied to degree-2.
return_corr (string, optional, default = False) – If set to True, return the theoretical global correlation.
- Py_Admittance.GlobalAdmitCorr(topo, grav, lmax=None)
Compute the global admittance and correlation functions from the input gravity and topography. For more information, see Broquet & Wieczorek (2019).
- Returns:
Admittance (array, dimension (lmax+1)) – Global admittance function in mGal/km.
Correlation (array, dimension (lmax+1)) – Global correlation function.
Admit_error (array, dimension (lmax+1)) – Global admittance uncertainty in mGal/km.
- Parameters:
topo (dimension (2, lmax+1, lmax+1)) – Spherical harmonic coefficients of the topography (km).
grav (dimension (2, lmax+1, lmax+1)) – Spherical harmonic coefficients of the gravity field (mGal).
lmax (int, optional, default = None) – Maximum degree at which the admittance and correlation are computed. If None, lmax = min(lmax_topo, lmax_grav). lmax must be <= min(lmax_topo, lmax_grav).
- Py_Admittance.LocalAdmitCorr(topo, grav, lat, lon, theta, lwin, lmax=None, quiet=True)
Compute the localized admittance and correlation functions from the input gravity and topography. For more information, see Broquet & Wieczorek (2019).
- Returns:
Admittance (array, dimension (lmax-lwin+1)) – Localized admittance function in mGal/km.
Correlation (array, dimension (lmax-lwin+1)) – Localized correlation function.
Admit_error (array, dimension (lmax-lwin+1)) – Localized admittance uncertainty in mGal/km.
- Parameters:
topo (dimension (2, lmax+1, lmax+1)) – Spherical harmonic coefficients of the topography (km).
grav (dimension (2, lmax+1, lmax+1)) – Spherical harmonic coefficients of the gravity field (mGal).
lat (float) – Central latitude (°) of the localization window.
lon (float) – Central longitude (°) of the localization window.
theta (float) – Angular radius (°) of the localization window.
lwin (int) – Bandwidth of the localization window.
lmax (int, optional, default = None) – Maximum degree at which the admittance and correlation are computed. If None, lmax = min(lmax_topo, lmax_grav). lmax must be <= min(lmax_topo, lmax_grav).
quiet (string, optional, default = True) – If False, the function will provide information regarding the spatio-spectral concentration of the localization window.
- Py_Admittance.TransferTGz(g0, R, Tc, Te, rhol, rhoc, rhom, rhobar, lmax, ratio_L=0, alpha_L=1, depth_L=50000.0, E=100000000000.0, v=0.25)
Compute theoretical transfer function given the input parameters, including the admittance and correlation. For more information, see Broquet & Wieczorek (2019).
- Returns:
T_s (array, dimension (lmax+1)) – Transfer function for flexure due to surface loads.
Qw_l (array, dimension (lmax+1)) – Transfer function for flexure due to surface/internal loads.
Qw_lz (array, dimension (lmax+1)) – Transfer function for flexure du to internal load.
Trsf_s (array, dimension (lmax+1)) – Transfer function for gravity from surface load.
Trsf_L (array, dimension (lmax+1)) – Transfer for gravity from internal load.
Correlation (array, dimension (lmax+1)) – Theoretical global correlation for out of phase surface/internal loads (alpha_L != 1).
- Parameters:
g0 (float) – Gravitational attraction at the surface.
R (float) – Mean radius of the planet.
Tc (float) – Average crustal thickness.
Te (float) – Elastic thickness of the lithosphere.
rhol (float) – Density of the surface topography.
rhoc (float) – Density of the crust.
rhom (float) – Density of the mantle.
lmax (int) – Maximum spherical harmonic degree for calculations.
ratio_L (float, optional, default = 0) – Ratio of the internal / surface load.
alpha_L (float, optional, default = 1) – Phase relationship for the internal / surface load. This parameter is experimental.
depth_L (float, optional, default = 50e3) – Depth of the internal load.
E (float, optional, default = 100e9) – Young’s modulus.
v (float, optional, default = 0.25) – Poisson’s ratio.