Tesseroids with variable density

Tesseroids with variable density#

The harmonica.tesseroid_gravity is capable of computing the gravitational effects of tesseroids whose density is defined through a continuous function of the radial coordinate. This is achieved by the application of the method introduced in [Soler2021].

To do so we need to define a regular Python function for the density, which should have a single argument (the radius coordinate) and return the density of the tesseroids at that radial coordinate. In addition, we need to decorate the density function with numba.jit(nopython=True) or numba.njit for short.

On this example we will show how we can compute the gravitational effect of four tesseroids whose densities are given by a custom linear density function.

tesseroid variable density
[[15.65236342 15.84941723 16.049421   ... 17.73452066 17.50786618
  17.28361491]
 [15.9363258  16.1466073  16.36042117 ... 18.17283119 17.92655366
  17.68346615]
 [16.2278511  16.45233251 16.68102124 ... 18.63085735 18.36295406
  18.0991969 ]
 ...
 [15.52701487 15.82032918 16.1231775  ... 18.50416206 18.12429333
  17.75648614]
 [15.24953447 15.52751852 15.81396405 ... 18.0637015  17.70888034
  17.36442893]
 [14.97982069 15.24343382 15.51455753 ... 17.64113169 17.30937979
  16.98653028]]

import boule as bl
import pygmt
import verde as vd
from numba import njit

import harmonica as hm

# Use the WGS84 ellipsoid to obtain the mean Earth radius which we'll use to
# reference the tesseroid
ellipsoid = bl.WGS84
mean_radius = ellipsoid.mean_radius

# Define tesseroid with top surface at the mean Earth radius, a thickness of
# 10km and a linear density function
tesseroids = (
    [-70, -60, -40, -30, mean_radius - 3e3, mean_radius],
    [-70, -60, -30, -20, mean_radius - 5e3, mean_radius],
    [-60, -50, -40, -30, mean_radius - 7e3, mean_radius],
    [-60, -50, -30, -20, mean_radius - 10e3, mean_radius],
)

# Define a linear density function. We should use the jit decorator so Numba
# can run the forward model efficiently.


@njit
def density(radius):
    """Linear density function"""
    top = mean_radius
    bottom = mean_radius - 10e3
    density_top = 2670
    density_bottom = 3000
    slope = (density_top - density_bottom) / (top - bottom)
    return slope * (radius - bottom) + density_bottom


# Define computation points on a regular grid at 100km above the mean Earth
# radius
coordinates = vd.grid_coordinates(
    region=[-80, -40, -50, -10],
    shape=(80, 80),
    extra_coords=100e3 + ellipsoid.mean_radius,
)

# Compute the radial component of the acceleration
gravity = hm.tesseroid_gravity(coordinates, tesseroids, density, field="g_z")
print(gravity)
grid = vd.make_xarray_grid(
    coordinates, gravity, data_names="gravity", extra_coords_names="extra"
)

# Plot the gravitational field
fig = pygmt.Figure()

title = "Downward component of gravitational acceleration"

with pygmt.config(FONT_TITLE="16p"):
    fig.grdimage(
        region=[-80, -40, -50, -10],
        projection="M-60/-30/10c",
        grid=grid.gravity,
        frame=["a", f"+t{title}"],
        cmap="viridis",
    )

fig.colorbar(cmap=True, frame=["a200f50", "x+lmGal"])

fig.coast(shorelines="1p,black")

fig.show()

Total running time of the script: (0 minutes 3.348 seconds)

Gallery generated by Sphinx-Gallery