Displacement_strain_planet package

Submodules

Module contents

Displacement_strain_planet

Displacement_strain_planet provides several functions and example scripts for generating crustal thickness, displacement, gravity, lateral density variations, stress, and strain maps on a planet given a set of input constraints such as from observed gravity and topography data.

These functions solve the Banerdt (1986) thin shell model under different assumptions. Various improvements have been made to the model including the possibility to account for finite-amplitude correction and filtering (Wieczorek & Phillips, 1998), lateral density variations at any arbitrary depth and within the surface or moho relief (Wieczorek et al., 2013), and density difference between the surface topography and crust (Broquet & Wieczorek, 2019).

We note that some of these functions relies heavily on the pyshtools package.

Thin_shell_matrix

Solve for the Banerdt et al. (1986) system of 5 equations with the possibility to account for finite-amplitude corrections and lateral density variations with the surface topography or moho relief.

Thin_shell_matrix_nmax

Solve the Banerdt (1986) system of 5 equations with finite-amplitude correction and accounting for the potential presence of density variations within the surface or moho reliefs.

DownContFilter

Compute the downward minimum-amplitude or -curvature filter of Wieczorek & Phillips, (1998).

corr_nmax_drho

Calculate the difference in gravitational exterior to relief referenced to a spherical interface (with or without laterally varying density) between the mass-sheet case and when using the finite amplitude algorithm of Wieczorek & Phillips (1998).

SH_deriv

Compute on the spherical harmonic derivatives (first and second order) at a given single colatiude/longitude location.

SH_deriv_store

Compute and store or load spherical harmonic derivatives (first and second order) over the whole sphere or given a set of colatiudes/longitudes.

Displacement_strains

Computes the Banerdt (1986) equations to determine strains from displacements with a correction to the theta_phi term.

Displacement_strains_shtools

Computes the Banerdt (1986) equations to determine strains and stresses from the displacements. This function uses SHTOOLS to derive the spherical harmonic gradients.

Principal_strainstress_angle

Calculate principal strains, stresses, and their principal angles.

Strainstress_from_principal

Calculate strains or stresses, from their principal values.

Plt_tecto_Mars

Plot the Knampeyer et al. (2006) dataset of extensional and compressional tectonic features on Mars.

Displacement_strain_planet.Displacement_strains(A_lm, w_lm, E, v, R, Te, lmax, depth=0, colat_min=0, colat_max=180, lon_min=0, lon_max=360, grid='DH', lmaxgrid=None, Y_lm_d1_t=None, Y_lm_d1_p=None, Y_lm_d2_t=None, Y_lm_d2_p=None, Y_lm_d2_tp=None, y_lm=None, path=None, quiet=True)

Computes the Banerdt (1986) equations to determine strains and stresses from the displacements.

Returns:

  • stress_theta (array, size(2*lmax+2,2*(2*lmax+2))) – Array with the stress field with respect to colatitude. This is equation A12 from Banerdt (1986).

  • stress_phi (array, size(2,lmax+1,lmax+1)) – Array with the stress field with respect to longitude. This is equation A13 from Banerdt (1986).

  • stress_theta_phi (array, size(2,lmax+1,lmax+1)) – Array with the stress field with respect to colatitude and longitude. This is equation A14 from Banerdt (1986).

  • eps_theta (array, size(2,lmax+1,lmax+1)) – Array with the elongation with respect to colatitude. This is equation A16 from Banerdt (1986).

  • eps_phi (array, size(2,lmax+1,lmax+1)) – Array with the elongation with respect to longitude. This is equation A17 from Banerdt (1986).

  • omega (array, size(2,lmax+1,lmax+1)) – Array with the shearing deformation. This is equation A18 from Banerdt (1986).

  • kappa_theta (array, size(2,lmax+1,lmax+1)) – Array with the bending deformation with respect to colatitude. This is equation A19 from Banerdt (1986).

  • kappa_phi (array, size(2,lmax+1,lmax+1)) – Array with the bending deformation with respect to longitude. This is equation A20 from Banerdt (1986).

  • tau (array, size(2,lmax+1,lmax+1)) – Array with the twisting deformation. This is equation A21 from Banerdt (1986).

  • tot_theta (array, size(2,lmax+1,lmax+1)) – Array with the total deformation with respect to colatitude.

  • tot_phi (array, size(2,lmax+1,lmax+1)) – Array with the total deformation with respect to longitude.

  • tot_thetaphi (array, size(2,lmax+1,lmax+1)) – Array with the total deformation with respect to colatitude and longitude.

Parameters:
  • A_lm (array, float, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the poloidal term of the tangential displacement.

  • w_lm (array, float, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the upward displacement.

  • E (float) – Young’s modulus.

  • v (float) – Poisson’s ratio.

  • R (float) – Mean radius of the planet.

  • Te (float) – Elastic thickness of the lithosphere.

  • lmax (int) – Maximum spherical harmonic degree for computations.

  • depth (float, optional, default = 0) – The depth at which stresses are estimated.

  • colat_min (float, optional, default = 0) – Minimum colatitude for grid computation of strains and stresses.

  • colat_max (float, optional, default = 180) – Maximum colatitude for grid computation of strains and stresses.

  • lon_min (float, optional, default = 0) – Minimum longitude for grid computation of strains and stresses.

  • lon_max (float, optional, default = 360) – Maximum longitude for grid computation of strains and stresses.

  • grid (string, optional, default = 'DH') – Either ‘DH’ or ‘GLQ’ for Driscoll and Healy grids or Gauss-Legendre Quadrature grids following the convention of SHTOOLs.

  • lmaxgrid (int, optional, default = None) – The maximum spherical harmonic degree resolvable by the grid. If None, this parameter is set to lmax. When grid==’GLQ’, the gridshape is (lmaxgrid+1, 2*lmaxgrid+1) and (2*lmaxgrid+2, 2*(2*lmaxgrid+2)) when grid==’DH’.

  • Y_lm_d1_t (array, float, size(2,lmax+1,lmax+1), optional, default = None) – Array with the first derivative of Legendre polynomials with respect to colatitude.

  • Y_lm_d1_p (array, float, size(2,lmax+1,lmax+1), optional, default = None) – Array with the first derivative of Legendre polynomials with respect to longitude.

  • Y_lm_d2_t (array, float, size(2,lmax+1,lmax+1), optional, default = None) – Array with the second derivative of Legendre polynomials with respect to colatitude.

  • Y_lm_d2_p (array, float, size(2,lmax+1,lmax+1), optional, default = None) – Array with the second derivative of Legendre polynomials with respect to longitude.

  • Y_lm_d2_tp (array, float, size(2,lmax+1,lmax+1), optional, default = None) – Array with the first derivative of Legendre polynomials with respect to colatitude and longitude.

  • y_lm (array, float, size(2,lmax+1,lmax+1), optional, default = None) – Array of spherical harmonic functions.

  • path (string, optional, default = None) – path where to find the stored Legendre polynomials.

  • quiet (bool, optional, default = True) – If True, suppress printing output.

Displacement_strain_planet.Displacement_strains_shtools(A_lm, w_lm, E, v, R, Te, lmax, depth=0, lmaxgrid=None, quiet=True)

Computes the Banerdt (1986) equations to determine strains and stresses from the displacements. This function uses SHTOOLS to derive the spherical harmonic gradients. This does not support GLQ grids.

Returns:

  • stress_theta (array, size(2*lmax+2,2*(2*lmax+2))) – Array with the stress field with respect to colatitude. This is equation A12 from Banerdt (1986).

  • stress_phi (array, size(2,lmax+1,lmax+1)) – Array with the stress field with respect to longitude. This is equation A13 from Banerdt (1986).

  • stress_theta_phi (array, size(2,lmax+1,lmax+1)) – Array with the stress field with respect to colatitude and longitude. This is equation A14 from Banerdt (1986).

  • eps_theta (array, size(2,lmax+1,lmax+1)) – Array with the elongation with respect to colatitude. This is equation A16 from Banerdt (1986).

  • eps_phi (array, size(2,lmax+1,lmax+1)) – Array with the elongation with respect to longitude. This is equation A17 from Banerdt (1986).

  • omega (array, size(2,lmax+1,lmax+1)) – Array with the shearing deformation. This is equation A18 from Banerdt (1986).

  • kappa_theta (array, size(2,lmax+1,lmax+1)) – Array with the bending deformation with respect to colatitude. This is equation A19 from Banerdt (1986).

  • kappa_phi (array, size(2,lmax+1,lmax+1)) – Array with the bending deformation with respect to longitude. This is equation A20 from Banerdt (1986).

  • tau (array, size(2,lmax+1,lmax+1)) – Array with the twisting deformation. This is equation A21 from Banerdt (1986).

  • tot_theta (array, size(2,lmax+1,lmax+1)) – Array with the total deformation with respect to colatitude.

  • tot_phi (array, size(2,lmax+1,lmax+1)) – Array with the total deformation with respect to longitude.

  • tot_thetaphi (array, size(2,lmax+1,lmax+1)) – Array with the total deformation with respect to colatitude and longitude.

Parameters:
  • A_lm (array, float, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the poloidal term of the tangential displacement.

  • w_lm (array, float, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the upward displacement.

  • E (float) – Young’s modulus.

  • v (float) – Poisson’s ratio.

  • R (float) – Mean radius of the planet.

  • Te (float) – Elastic thickness of the lithosphere.

  • lmax (int) – Maximum spherical harmonic degree for computations.

  • depth (float, optional, default = 0) – The depth at which stresses are estimated.

  • lmaxgrid (int, optional, default = None) – The maximum spherical harmonic degree resolvable by the grid. If None, this parameter is set to lmax. When grid==’GLQ’, the gridshape is (lmaxgrid+1, 2*lmaxgrid+1) and (2*lmaxgrid+2, 2*(2*lmaxgrid+2)) when grid==’DH’. Quadrature grids following the convention of SHTOOLs. If None, the grid is set to ‘GLQ’.

  • quiet (bool, optional, default = True) – If True, suppress printing output.

Displacement_strain_planet.DownContFilter(l, half, R_ref, D_relief, type='Mc', quiet=False)

Compute the downward minimum-amplitude or -curvature filter of Wieczorek & Phillips, (1998).

Returns:

Value of the filter at degrees l

Return type:

float

Parameters:
  • l (array) – Array of spherical harmonic degrees.

  • half (int) – The spherical harmonic degree where the filter is equal to 0.5.

  • R_ref (float) – The reference radius of the gravitational field.

  • D_relief (float) – The radius of the surface to downward continue to.

  • type (string, optional, default = "Mc") – Filter type, minimum amplitude (“Ma”) of curvature (“Mc”)

  • quiet (bool, optional, default = True) – if True, prints a warning when D_relief > R_ref.

Displacement_strain_planet.Plt_tecto_Mars(path, compression=False, extension=True, ax=None, compression_col='k', extension_col='purple', lw=1, legend_show=True, legend_loc='upper left')

Plot the Knampeyer et al. (2006) dataset of extensional and compressional tectonic features on Mars.

Parameters:
  • path (string) – path for the location of the Knameyer et al (2006) dataset.

  • compression (bool, optional, default = False) – If True, plot compressive tectonic features.

  • extension (bool, optional, default = True) – If True, plot extensive tectonic features.

  • ax (array of object, optional, default = None) – Matplotlib axes.

  • compression_col (string, optional, default = "k") – Color of compressive tectonic features.

  • extension_col (string, optional, default = "purple") – Color of extensive tectonic features.

  • lw (int, optional, default = 1) – Linewidth for the tectonic features

  • legend_show (bool, optional, default = True) – If True, add a legend to the plot.

  • legend_loc (string, optional, default = "upper left") – Determine the legend position.

Displacement_strain_planet.Principal_strainstress_angle(s_theta, s_phi, s_theta_phi)

Calculate principal strains, stresses, and their principal angles.

Returns:

  • min_strain (array, size same as input arrays) – Array with the minimum principal horizontal strain or stress.

  • max_strain (array, size same as input arrays) – Array with the maximum principal horizontal strain or stress.

  • sum_strain (array, size same as input arrays) – Array with the sum of the principal horizontal strain or stress.

  • principal_angle (array, size same as input arrays) – Array with the principal strain or stress direction in degrees.

Parameters:
  • s_theta (array, float, size(nlat, nlon)) – Array of the colatitude component of the stress or strain field.

  • s_phi (array, float, size(nlat, nlon)) – Array of the longitude component of the stress or strain field.

  • s_theta_phi (array, float, size(nlat, nlon)) – Array of the colatitude and longitude component of the stress or strain field.

Displacement_strain_planet.SH_deriv(theta, phi, lmax)

Compute spherical harmonic derivatives at a given location (first and second order).

Returns:

  • Y_lm_d1_theta_a (array, size(2,lmax+1,lmax+1)) – Array with the first derivative of Legendre polynomials with respect to colatitude.

  • Y_lm_d1_phi_a (array, size(2,lmax+1,lmax+1)) – Array with the first derivative of Legendre polynomials with respect to longitude.

  • Y_lm_d2_theta_a (array, size(2,lmax+1,lmax+1)) – Array with the second derivative of Legendre polynomials with respect to colatitude.

  • Y_lm_d2_phi_a (array, size(2,lmax+1,lmax+1)) – Array with the second derivative of Legendre polynomials with respect to longitude.

  • Y_lm_d2_thetaphi_a (array, size(2,lmax+1,lmax+1)) – Array with the first derivative of Legendre polynomials with respect to colatitude and longitude.

  • y_lm (array, size(2,lmax+1,lmax+1)) – Array of spherical harmonic functions.

Parameters:
  • theta (float) – Colatitude in radian.

  • phi (float) – Longitude in radian.

  • lmax (int) – Maximum spherical harmonic degree to compute for the derivatives.

Displacement_strain_planet.SH_deriv_store(lmax, path, colat_min=0, colat_max=180, lon_min=0, lon_max=360, grid='DH', dtype=<class 'numpy.float64'>, lmaxgrid=None, save=True, compressed=False, quiet=True)

Compute and store or load spherical harmonic derivatives (first and second order) over the entire sphere or given a set of colatiudes/longitudes bounds. The spherical harmonic degree and order correspond to the index l*(l+1)/2+m. This routine supports both Driscoll and Healy (DH) and Gauss-Legendre Quadrature (GLQ) grids.

Returns:

  • Y_lm_d1_theta_a (array, size(2,(lmax+1)*(lmax+2)/2)) – Array with the first derivative of Legendre polynomials with respect to colatitude.

  • Y_lm_d1_phi_a (array, size(2,(lmax+1)*(lmax+2)/2)) – Array with the first derivative of Legendre polynomials with respect to longitude.

  • Y_lm_d2_theta_a (array, size(2,(lmax+1)*(lmax+2)/2)) – Array with the second derivative of Legendre polynomials with respect to colatitude.

  • Y_lm_d2_phi_a (array, size(2,(lmax+1)*(lmax+2)/2)) – Array with the second derivative of Legendre polynomials with respect to longitude.

  • Y_lm_d2_thetaphi_a (array, size(2,(lmax+1)*(lmax+2)/2)) – Array with the first derivative of Legendre polynomials with respect to colatitude and longitude.

  • y_lm_save (array, size(2,(lmax+1)*(lmax+2)/2)) – Array of spherical harmonic functions.

Parameters:
  • lmax (int) – Maximum spherical harmonic degree to compute for the derivatives.

  • path (string) – Path to store or load spherical harmonic derivatives.

  • colat_min (float, optional, default = 0) – Minimum colatitude for grid computation of SH derivatives.

  • colat_max (float, optional, default = 180) – Maximum colatitude for grid computation of SH derivatives.

  • lon_min (float, optional, default = 0) – Minimum longitude for grid computation of SH derivatives.

  • lon_max (float, optional, default = 360) – Maximum longitude for grid computation of SH derivatives.

  • grid (string, optional, default = 'DH') – Either ‘DH’ or ‘GLQ’ for Driscoll and Healy grids or Gauss-Legendre Quadrature grids following the convention of SHTOOLs.

  • dtype (data-type, optional, default = numpy.float64) – The desired data-type for the arrays (default is that of numpy). This can help reducing the size of the stored array.

  • lmaxgrid (int, optional, default = None) – The maximum spherical harmonic degree resolvable by the grid. If None, this parameter is set to lmax. The gridshape is (2*lmaxgrid+2, 2*(2*lmaxgrid+2)), DH2 grid. If None, the grid is set to ‘GLQ’.

  • save (bool, optional, default = True) – If True, save the data at the given path location.

  • compressed (bool, optional, default = False) – If True, the data is saved in compressed .npz format instead of npy, which decreases the file size by about a factor 2. This is recommended when lmax > 75.

  • quiet (bool, optional, default = True) – If True, suppress printing output.

Displacement_strain_planet.Strainstress_from_principal(min_strain, max_strain, sum_strain, principal_angle)

Calculate strains or stresses, from their principal values.

Returns:

  • s_theta (array, float, size same as input arrays) – Array of the colatitude component of the stress or strain field.

  • s_phi (array, float, size same as input arrays) – Array of the longitude component of the stress or strain field.

  • s_theta_phi (array, float, size same as input arrays) – Array of the colatitude and longitude component of the stress or strain field.

Parameters:
  • min_strain (array, size(nlat, nlon)) – Array with the minimum principal horizontal strain or stress.

  • max_strain (array, size(nlat, nlon)) – Array with the maximum principal horizontal strain or stress.

  • sum_strain (array, size(nlat, nlon)) – Array with the sum of the principal horizontal strain or stress.

  • principal_angle (array, size(nlat, nlon)) – Array with the principal strain or stress direction in degrees.

Displacement_strain_planet.Thin_shell_matrix(g0, R, c, Te, rhom, rhoc, rhol, rhobar, lmax, E, v, base_drho=50000.0, top_drho=0, filter_in=None, filter=None, filter_half=50, H_lm=None, drhom_lm=None, dc_lm=None, w_lm=None, omega_lm=None, q_lm=None, G_lm=None, Gc_lm=None, add_equation=None, add_arrays=None, quiet=False, remove_equation=None, w_corr=None, wdc_corr=None, H_corr=None, drho_omega_corr=None, drho_q_corr=None, COM=True, lambdify_func=None, first_inv=True)

Solve for the Banerdt et al. (1986) system of equations with the possibility to account for finite-amplitude corrections and lateral density variations with the surface topography or moho relief.

Returns:

  • w_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the upward displacement.

  • A_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the poloidal term of the tangential displacement.

  • moho_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the moho-relief.

  • dc_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the crustal root variations.

  • drhom_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the lateral density variations.

  • omega_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the tangential load potential.

  • q_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the net load on the lithosphere.

  • Gc_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the geoid at the moho depth.

  • G_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the geoid at the surface.

  • H_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the surface topography.

  • lambdify_func (array, size(2,lmax+1,lmax+1)) – Array with the lambda functions (size lmax+1) of all components. Lambda functions can be used to re-calculate the same problem with different inputs very fast.

Parameters:
  • g0 (float) – Gravitational attraction at the surface.

  • R (float) – Mean radius of the planet.

  • c (float) – Average crustal thickness.

  • Te (float) – Elastic thickness of the lithosphere.

  • rhom (float) – Density of the mantle.

  • rhoc (float) – Density of the crust.

  • rhol (float) – Density of the surface topography.

  • rhobar (float) – Mean density of the planet.

  • lmax (int) – Maximum spherical harmonic degree of calculations.

  • E (float) – Young’s modulus.

  • v (float) – Poisson’s ratio.

  • base_drho (float, optional, default = 50e3) – Lower depth for the of the density contrast.

  • top_drho (float, optional, default = 0) – Upper depth for the of the density contrast.

  • filter_in (array, size(lmax+1), optional, default = None.) – Array with the input filter to use.

  • filter (string, optional, default = None) – If ‘Ma’ or ‘Mc’, apply minimum-amplitude or minimum-curvature filtering. If None, no filtering.

  • filter_half (int, default = 50) – Spherical harmonic degree at which the filter equals 0.5.

  • H_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the surface topography.

  • drhom_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the lateral density variations.

  • dc_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the crustal root variations.

  • w_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the upward displacement.

  • omega_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the tangential load potential.

  • q_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the net load on the lithosphere.

  • G_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the geoid at the surface.

  • Gc_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the geoid at the moho depth.

  • add_equation (string, optional, default = None) – Equation to be added to the system. This must include at least one of the 8 parameters aboves.

  • add_arrays (array size(N, 2,lmax+1,lmax+1), optional, default = None) – N arrays of spherical harmonics to be added in ‘add_equation’, which are written ‘add_array1’ ‘add_array2’ etc. Order is important.

  • quiet (bool, optional, default = False) – if True, print various outputs.

  • remove_equation (string, optional, default = None) – String of the equation to be removed. This must be either ‘G_lm’, ‘Gc_lm’, ‘w_lm’, ‘omega_lm’, or ‘q_lm’.

  • w_corr (array size(2,lmax+1,lmax+1), optional, default = None) – Array with spherical harmonic coefficients for finite-amplitude and or lateral density variations corrections of the w_lm relief.

  • wdc_corr (array size(2,lmax+1,lmax+1), optional, default = None) – Array with spherical harmonic coefficients for finite-amplitude and or lateral density variations corrections of the moho_lm relief.

  • H_corr (array size(2,lmax+1,lmax+1), optional, default = None) – Array with spherical harmonic coefficients for finite-amplitude and or lateral density variations corrections of the H_lm relief.

  • drho_omega_corr (array size(2,lmax+1,lmax+1), optional, default = None) – Array with spherical harmonic coefficients for lateral lateral density variations corrections for omega_lm.

  • drho_q_corr (array size(2,lmax+1,lmax+1), optional, default = None) – Array with spherical harmonic coefficients for lateral lateral density variations corrections for q_lm.

  • COM (bool, optional, default = True) – if True, force the model to be in a center-of-mass frame by setting the degree-1 geoid terms to zero.

  • lambdify_func (array size(lmax+1), optional, default = None) – Use the lambidfy functions (i.e. design of the inversion matrix, without the specific inputs) of another run.

  • first_inv (bool, optional, default = True) – If True, the code assumes that this is the first time doing the inversion in this setup, and will store the lambdify results in ‘lambdify_func’

Displacement_strain_planet.Thin_shell_matrix_nmax(g0, R, c, Te, rhom, rhoc, rhol, lmax, E, v, mass, filter_in=None, filter=None, filter_half=50, nmax=5, lmaxgrid=None, H_lm=None, drhom_lm=None, dc_lm=None, w_lm=None, omega_lm=None, q_lm=None, G_lm=None, Gc_lm=None, C_lm=None, add_equation=None, add_arrays=None, quiet=True, remove_equation=None, COM=True, base_drho=50000.0, top_drho=0, delta_max=5, iter_max=250, delta_out=500000.0, iterate=True)

Solve the Banerdt (1986) system of 5 equations with finite-amplitude correction and accounting for the potential presence of density variations within the surface or moho reliefs.

Returns:

  • w_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the upward displacement.

  • A_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the poloidal term of the tangential displacement.

  • moho_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the moho-relief.

  • dc_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the crustal root variations.

  • drhom_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the lateral density variations.

  • omega_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the tangential load potential.

  • q_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the net load on the lithosphere.

  • Gc_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the geoid at the moho depth.

  • G_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the geoid at the surface.

  • H_lm (array, size(2,lmax+1,lmax+1)) – Array with the spherical harmonic coefficients of the surface topography.

  • lambdify_func (array, size(2,lmax+1,lmax+1)) – Array with the lambda functions (i.e., the design of the inversion matrix without the inputs) of all components. Lambda functions can be used to re-calculate the same problem with different inputs very fast.

Parameters:
  • g0 (float) – Gravitational attraction at the surface.

  • R (float) – Mean radius of the planet.

  • c (float) – Average crustal thickness.

  • Te (float) – Elastic thickness of the lithosphere.

  • rhom (float) – Density of the mantle.

  • rhoc (float) – Density of the crust.

  • rhol (float) – Density of the surface topography.

  • lmax (int) – Maximum spherical harmonic degree of calculations.

  • E (float) – Young’s modulus.

  • v (float) – Poisson’s ratio.

  • mass (float) – Mass of the planet.

  • filter_in (array, size(lmax+1), optional, default = None.) – Array with the input filter to use.

  • filter (string, optional, default = None) – If ‘Ma’ or ‘Mc’, apply minimum-amplitude or minimum-curvature filtering. If None, no filtering.

  • filter_half (int, optional, default = 50) – Spherical harmonic degree at which the filter equals 0.5.

  • nmax (int, optional, default = 5) – Maximum order of the finite-amplitude correction.

  • lmaxgrid (int, optional, default = None) – If None, this parameter is set to 3*lmax. Resolution of the input grid for the finite-amplitude correction routines. 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.

  • H_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the surface topography.

  • drhom_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the lateral density variations.

  • dc_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the crustal root variations.

  • w_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the upward displacement.

  • omega_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the tangential load potential.

  • q_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the net load on the lithosphere.

  • G_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the geoid at the surface.

  • Gc_lm (array, size(2,lmax+1,lmax+1), optional, default = None) – Array with the spherical harmonic coefficients of the geoid at the moho depth.

  • add_equation (string, optional, default = None) – Equation to be added to the system. This must include at least one of the 8 parameters aboves.

  • add_arrays (array size(N, 2,lmax+1,lmax+1), optional, default = None) – N arrays of spherical harmonics to be added in ‘add_equation’, which are written ‘add_array1’ ‘add_array2’ etc. Order is important.

  • quiet (bool, optional, default = False) – if True, print various outputs.

  • COM (bool, optional, default = True) – if True, force the model to be in a center-of-mass frame by setting the degree-1 geoid terms to zero.

  • remove_equation (string, optional, default = None) – String of the equation to be removed. This must be either ‘G_lm’, ‘Gc_lm’, ‘w_lm’, ‘omega_lm’, or ‘q_lm’.

  • base_drho (float, optional, default = 50e3) – Lower depth for the of the density contrast.

  • top_drho (float, optional, default = 0) – Upper depth for the of the density contrast.

  • delta_max (float, optional, default = 5) – The algorithm will continue to iterate until the maximum difference in relief (or density contrast) between solutions is less than this value (in meters or kg m-3).

  • iter_max (int, optional, default = 250) – Maximum number of iterations before the algorithm stops.

  • delta_out (float, optional, default = 500e3) – If the delta is larger than this value, the algorithm stops and prints that it is not converging.

  • iterate (bool, optional, default = True) – if False, solve the system without any corrections.

Displacement_strain_planet.corr_nmax_drho(dr_lm, shape_grid, rho_grid, lmax, mass, nmax, drho, R, density_var=False)

Calculate the gravitational difference (with or without laterally varying density) between the mass-sheet case and when using the finite amplitude algorithm of Wieczorek & Phillips (1998).

Returns:

Array with the spherical harmonic coefficients of the difference between the mass-sheet and finite-ampltiude geoid.

Return type:

array, size of input dr_lm

Parameters:
  • dr_lm (array, size (2,lmax+1,lmax+1)) – Array with spherical harmonic coefficients of the relief.

  • shape_grid (array, size (2,2*(lmax+1),2*2(lmax+1))) – Array with a grid of the relief.

  • rho_grid (array, size (2,2*(lmax+1),2*2(lmax+1))) – Array with a grid of the lateral density contrast.

  • lmax (int) – Maximum spherical harmonic degree to compute for the derivatives.

  • mass (float) – Mass of the planet.

  • nmax (int) – Order of the finite-amplitude correction.

  • drho (float) – Mean density contrast.

  • R (float) – Mean radius of the planet.

  • density_var (bool, optional, default = False) – If True, correct for density variations.