Displacement_strain_planet.B1986_nmax module

Functions for calculating the Banerdt (1986) system of equations.

Displacement_strain_planet.B1986_nmax.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.B1986_nmax.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.B1986_nmax.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.B1986_nmax.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.

Displacement_strain_planet.B1986_nmax.test_symb(str_symb, arr, constraint_test, not_constraint, arr_symb)

This function return None or the input array depending on the input constraints in Thin_shell_matrix.

Returns:

Input array or None

Return type:

array, size of input arr or None

Parameters:
  • str_symb (string) – String of the investigated symbol.

  • arr (array, size (2,lmax+1,lmax+1)) – Array with spherical harmonic coefficients of the input array.

  • constraint_test (list of strings, size variable) – List of input constraints (i.e., ‘G_lm’, ‘drhom_lm’…).

  • not_constraint (list of strings, size variable) – List of strings that are not input constraints (i.e., ‘Gc_lm’).

  • arr_symb (list of sympy symbols) – List of all sympy symbols.