Sparse Norm Inversion of Gravity Anomaly Data#

Here we invert gravity anomaly data to recover a density contrast model. We formulate the inverse problem as an iteratively re-weighted least-squares (IRLS) optimization problem. For this tutorial, we focus on the following:

  • Defining the survey from xyz formatted data

  • Generating a mesh based on survey geometry

  • Including surface topography

  • Defining the inverse problem (data misfit, regularization, optimization)

  • Specifying directives for the inversion

  • Setting sparse and blocky norms

  • Plotting the recovered model and data misfit

Although we consider gravity anomaly data in this tutorial, the same approach can be used to invert gradiometry and other types of geophysical data.

Import modules#

import os
import numpy as np
import matplotlib as mpl
import matplotlib.pyplot as plt
import tarfile

from discretize import TensorMesh
from discretize.utils import active_from_xyz
from SimPEG.utils import plot2Ddata, model_builder
from SimPEG.potential_fields import gravity
from SimPEG import (
    maps,
    data,
    data_misfit,
    inverse_problem,
    regularization,
    optimization,
    directives,
    inversion,
    utils,
)

# sphinx_gallery_thumbnail_number = 3

Define File Names#

File paths for assets we are loading. To set up the inversion, we require topography and field observations. The true model defined on the whole mesh is loaded to compare with the inversion result. These files are stored as a tar-file on our google cloud bucket: “https://storage.googleapis.com/simpeg/doc-assets/gravity.tar.gz

# storage bucket where we have the data
data_source = "https://storage.googleapis.com/simpeg/doc-assets/gravity.tar.gz"

# download the data
downloaded_data = utils.download(data_source, overwrite=True)

# unzip the tarfile
tar = tarfile.open(downloaded_data, "r")
tar.extractall()
tar.close()

# path to the directory containing our data
dir_path = downloaded_data.split(".")[0] + os.path.sep

# files to work with
topo_filename = dir_path + "gravity_topo.txt"
data_filename = dir_path + "gravity_data.obs"
model_filename = dir_path + "true_model.txt"
overwriting /home/vsts/work/1/s/tutorials/03-gravity/gravity.tar.gz
Downloading https://storage.googleapis.com/simpeg/doc-assets/gravity.tar.gz
   saved to: /home/vsts/work/1/s/tutorials/03-gravity/gravity.tar.gz
Download completed!

Load Data and Plot#

Here we load and plot synthetic gravity anomaly data. Topography is generally defined as an (N, 3) array. Gravity data is generally defined with 4 columns: x, y, z and data.

# Load topography
xyz_topo = np.loadtxt(str(topo_filename))

# Load field data
dobs = np.loadtxt(str(data_filename))

# Define receiver locations and observed data
receiver_locations = dobs[:, 0:3]
dobs = dobs[:, -1]

# Plot
mpl.rcParams.update({"font.size": 12})
fig = plt.figure(figsize=(7, 5))

ax1 = fig.add_axes([0.1, 0.1, 0.73, 0.85])
plot2Ddata(receiver_locations, dobs, ax=ax1, contourOpts={"cmap": "bwr"})
ax1.set_title("Gravity Anomaly")
ax1.set_xlabel("x (m)")
ax1.set_ylabel("y (m)")

ax2 = fig.add_axes([0.8, 0.1, 0.03, 0.85])
norm = mpl.colors.Normalize(vmin=-np.max(np.abs(dobs)), vmax=np.max(np.abs(dobs)))
cbar = mpl.colorbar.ColorbarBase(
    ax2, norm=norm, orientation="vertical", cmap=mpl.cm.bwr, format="%.1e"
)
cbar.set_label("$mgal$", rotation=270, labelpad=15, size=12)

plt.show()
Gravity Anomaly

Assign Uncertainties#

Inversion with SimPEG requires that we define standard deviation on our data. This represents our estimate of the noise in our data. For gravity inversion, a constant floor value is generally applied to all data. For this tutorial, the standard deviation on each datum will be 1% of the maximum observed gravity anomaly value.

Defining the Survey#

Here, we define survey that will be used for this tutorial. Gravity surveys are simple to create. The user only needs an (N, 3) array to define the xyz locations of the observation locations. From this, the user can define the receivers and the source field.

# Define the receivers. The data consist of vertical gravity anomaly measurements.
# The set of receivers must be defined as a list.
receiver_list = gravity.receivers.Point(receiver_locations, components="gz")

receiver_list = [receiver_list]

# Define the source field
source_field = gravity.sources.SourceField(receiver_list=receiver_list)

# Define the survey
survey = gravity.survey.Survey(source_field)

Defining the Data#

Here is where we define the data that are inverted. The data are defined by the survey, the observation values and the standard deviation.

data_object = data.Data(survey, dobs=dobs, standard_deviation=uncertainties)

Defining a Tensor Mesh#

Here, we create the tensor mesh that will be used to invert gravity anomaly data. If desired, we could define an OcTree mesh.

dh = 5.0
hx = [(dh, 5, -1.3), (dh, 40), (dh, 5, 1.3)]
hy = [(dh, 5, -1.3), (dh, 40), (dh, 5, 1.3)]
hz = [(dh, 5, -1.3), (dh, 15)]
mesh = TensorMesh([hx, hy, hz], "CCN")

Starting/Reference Model and Mapping on Tensor Mesh#

Here, we create starting and/or reference models for the inversion as well as the mapping from the model space to the active cells. Starting and reference models can be a constant background value or contain a-priori structures.

# Find the indices of the active cells in forward model (ones below surface)
ind_active = active_from_xyz(mesh, xyz_topo)

# Define mapping from model to active cells
nC = int(ind_active.sum())
model_map = maps.IdentityMap(nP=nC)  # model consists of a value for each active cell

# Define and plot starting model
starting_model = np.zeros(nC)

Define the Physics#

Here, we define the physics of the gravity problem by using the simulation class.

simulation = gravity.simulation.Simulation3DIntegral(
    survey=survey, mesh=mesh, rhoMap=model_map, ind_active=ind_active
)

Define the Inverse Problem#

The inverse problem is defined by 3 things:

  1. Data Misfit: a measure of how well our recovered model explains the field data

  2. Regularization: constraints placed on the recovered model and a priori information

  3. Optimization: the numerical approach used to solve the inverse problem

# Define the data misfit. Here the data misfit is the L2 norm of the weighted
# residual between the observed data and the data predicted for a given model.
# Within the data misfit, the residual between predicted and observed data are
# normalized by the data's standard deviation.
dmis = data_misfit.L2DataMisfit(data=data_object, simulation=simulation)
dmis.W = utils.sdiag(1 / uncertainties)

# Define the regularization (model objective function).
reg = regularization.Sparse(mesh, active_cells=ind_active, mapping=model_map)
reg.norms = [0, 2, 2, 2]

# Define how the optimization problem is solved. Here we will use a projected
# Gauss-Newton approach that employs the conjugate gradient solver.
opt = optimization.ProjectedGNCG(
    maxIter=100, lower=-1.0, upper=1.0, maxIterLS=20, maxIterCG=10, tolCG=1e-3
)

# Here we define the inverse problem that is to be solved
inv_prob = inverse_problem.BaseInvProblem(dmis, reg, opt)

Define Inversion Directives#

Here we define any directiveas that are carried out during the inversion. This includes the cooling schedule for the trade-off parameter (beta), stopping criteria for the inversion and saving inversion results at each iteration.

# Defining a starting value for the trade-off parameter (beta) between the data
# misfit and the regularization.
starting_beta = directives.BetaEstimate_ByEig(beta0_ratio=1e0)

# Defines the directives for the IRLS regularization. This includes setting
# the cooling schedule for the trade-off parameter.
update_IRLS = directives.Update_IRLS(
    f_min_change=1e-4,
    max_irls_iterations=30,
    coolEpsFact=1.5,
    beta_tol=1e-2,
)

# Options for outputting recovered models and predicted data for each beta.
save_iteration = directives.SaveOutputEveryIteration(save_txt=False)

# Updating the preconditionner if it is model dependent.
update_jacobi = directives.UpdatePreconditioner()

# Add sensitivity weights
sensitivity_weights = directives.UpdateSensitivityWeights(everyIter=False)

# The directives are defined as a list.
directives_list = [
    update_IRLS,
    sensitivity_weights,
    starting_beta,
    save_iteration,
    update_jacobi,
]
/home/vsts/work/1/s/tutorials/03-gravity/plot_inv_1b_gravity_anomaly_irls.py:277: UserWarning:

'everyIter' property is deprecated and will be removed in SimPEG 0.20.0.Please use 'every_iteration'.

Running the Inversion#

To define the inversion object, we need to define the inversion problem and the set of directives. We can then run the inversion.

# Here we combine the inverse problem and the set of directives
inv = inversion.BaseInversion(inv_prob, directives_list)

# Run inversion
recovered_model = inv.run(starting_model)
SimPEG.InvProblem will set Regularization.reference_model to m0.
SimPEG.InvProblem will set Regularization.reference_model to m0.
SimPEG.InvProblem will set Regularization.reference_model to m0.
SimPEG.InvProblem will set Regularization.reference_model to m0.

                    SimPEG.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
                    ***Done using the default solver Pardiso and no solver_opts.***

model has any nan: 0
=============================== Projected GNCG ===============================
  #     beta     phi_d     phi_m       f      |proj(x-g)-x|  LS    Comment
-----------------------------------------------------------------------------
x0 has any nan: 0
   0  3.21e+03  1.33e+05  0.00e+00  1.33e+05    2.19e+02      0
   1  1.60e+03  9.43e+03  5.73e+00  1.86e+04    2.13e+02      0
   2  8.01e+02  4.10e+03  8.03e+00  1.05e+04    2.10e+02      0   Skip BFGS
   3  4.01e+02  1.63e+03  1.02e+01  5.70e+03    2.02e+02      0   Skip BFGS
   4  2.00e+02  6.20e+02  1.19e+01  3.00e+03    1.86e+02      0   Skip BFGS
   5  1.00e+02  2.36e+02  1.32e+01  1.56e+03    1.72e+02      0   Skip BFGS
Reached starting chifact with l2-norm regularization: Start IRLS steps...
irls_threshold 0.054708593865006566
   6  5.01e+01  9.51e+01  2.72e+01  1.46e+03    1.81e+02      0   Skip BFGS
   7  1.28e+02  4.66e+01  3.74e+01  4.82e+03    2.12e+02      0
   8  8.31e+01  2.63e+02  3.70e+01  3.34e+03    1.95e+02      0
   9  7.00e+01  1.54e+02  3.80e+01  2.82e+03    2.09e+02      0
  10  5.99e+01  1.50e+02  3.43e+01  2.21e+03    1.98e+02      0
  11  9.59e+01  1.20e+02  2.97e+01  2.96e+03    2.12e+02      0
  12  7.33e+01  1.85e+02  2.50e+01  2.02e+03    2.03e+02      0
  13  1.14e+02  1.30e+02  2.17e+01  2.60e+03    2.15e+02      0   Skip BFGS
  14  8.92e+01  1.77e+02  1.79e+01  1.78e+03    2.00e+02      0
  15  1.40e+02  1.28e+02  1.60e+01  2.36e+03    2.16e+02      0   Skip BFGS
  16  1.11e+02  1.72e+02  1.37e+01  1.69e+03    2.09e+02      0
  17  1.69e+02  1.38e+02  1.27e+01  2.28e+03    2.17e+02      0   Skip BFGS
  18  1.27e+02  1.91e+02  1.13e+01  1.63e+03    2.05e+02      0
  19  1.10e+02  1.49e+02  1.08e+01  1.34e+03    2.14e+02      0   Skip BFGS
  20  1.71e+02  1.29e+02  1.05e+01  1.92e+03    2.18e+02      0   Skip BFGS
  21  1.34e+02  1.76e+02  9.63e+00  1.47e+03    2.11e+02      0
  22  2.04e+02  1.40e+02  9.42e+00  2.06e+03    2.18e+02      0   Skip BFGS
  23  1.60e+02  1.77e+02  8.85e+00  1.59e+03    2.12e+02      0
  24  2.43e+02  1.38e+02  8.73e+00  2.26e+03    2.18e+02      0   Skip BFGS
  25  1.91e+02  1.75e+02  8.31e+00  1.77e+03    2.13e+02      0
  26  2.91e+02  1.38e+02  8.21e+00  2.53e+03    2.18e+02      0   Skip BFGS
  27  2.29e+02  1.77e+02  7.88e+00  1.98e+03    2.13e+02      0
  28  3.45e+02  1.42e+02  7.89e+00  2.86e+03    2.18e+02      0
  29  2.64e+02  1.86e+02  7.58e+00  2.18e+03    2.11e+02      0
  30  3.98e+02  1.41e+02  7.56e+00  3.15e+03    2.18e+02      0   Skip BFGS
  31  3.00e+02  1.90e+02  7.30e+00  2.38e+03    2.10e+02      0
  32  3.00e+02  1.46e+02  7.34e+00  2.35e+03    2.14e+02      0
  33  4.53e+02  1.43e+02  7.29e+00  3.44e+03    2.18e+02      0   Skip BFGS
  34  3.34e+02  1.99e+02  7.03e+00  2.55e+03    2.08e+02      0
  35  2.89e+02  1.47e+02  7.09e+00  2.20e+03    2.14e+02      0
Reach maximum number of IRLS cycles: 30
------------------------- STOP! -------------------------
1 : |fc-fOld| = 0.0000e+00 <= tolF*(1+|f0|) = 1.3311e+04
1 : |xc-x_last| = 4.6490e-02 <= tolX*(1+|x0|) = 1.0000e-01
0 : |proj(x-g)-x|    = 2.1446e+02 <= tolG          = 1.0000e-01
0 : |proj(x-g)-x|    = 2.1446e+02 <= 1e3*eps       = 1.0000e-02
0 : maxIter   =     100    <= iter          =     36
------------------------- DONE! -------------------------

Recreate True Model#

# Define density contrast values for each unit in g/cc
background_density = 0.0
block_density = -0.2
sphere_density = 0.2

# Define model. Models in SimPEG are vector arrays.
true_model = background_density * np.ones(nC)

# You could find the indicies of specific cells within the model and change their
# value to add structures.
ind_block = (
    (mesh.gridCC[ind_active, 0] > -50.0)
    & (mesh.gridCC[ind_active, 0] < -20.0)
    & (mesh.gridCC[ind_active, 1] > -15.0)
    & (mesh.gridCC[ind_active, 1] < 15.0)
    & (mesh.gridCC[ind_active, 2] > -50.0)
    & (mesh.gridCC[ind_active, 2] < -30.0)
)
true_model[ind_block] = block_density

# You can also use SimPEG utilities to add structures to the model more concisely
ind_sphere = model_builder.getIndicesSphere(np.r_[35.0, 0.0, -40.0], 15.0, mesh.gridCC)
ind_sphere = ind_sphere[ind_active]
true_model[ind_sphere] = sphere_density

Plotting True Model and Recovered Model#

# Plot True Model
fig = plt.figure(figsize=(9, 4))
plotting_map = maps.InjectActiveCells(mesh, ind_active, np.nan)

ax1 = fig.add_axes([0.1, 0.1, 0.73, 0.8])
mesh.plot_slice(
    plotting_map * true_model,
    normal="Y",
    ax=ax1,
    ind=int(mesh.shape_cells[1] / 2),
    grid=True,
    clim=(np.min(true_model), np.max(true_model)),
    pcolor_opts={"cmap": "viridis"},
)
ax1.set_title("Model slice at y = 0 m")


ax2 = fig.add_axes([0.85, 0.1, 0.05, 0.8])
norm = mpl.colors.Normalize(vmin=np.min(true_model), vmax=np.max(true_model))
cbar = mpl.colorbar.ColorbarBase(
    ax2, norm=norm, orientation="vertical", cmap=mpl.cm.viridis, format="%.1e"
)
cbar.set_label("$g/cm^3$", rotation=270, labelpad=15, size=12)

plt.show()

# Plot Recovered Model
fig = plt.figure(figsize=(9, 4))
plotting_map = maps.InjectActiveCells(mesh, ind_active, np.nan)

ax1 = fig.add_axes([0.1, 0.1, 0.73, 0.8])
mesh.plot_slice(
    plotting_map * recovered_model,
    normal="Y",
    ax=ax1,
    ind=int(mesh.shape_cells[1] / 2),
    grid=True,
    clim=(np.min(recovered_model), np.max(recovered_model)),
    pcolor_opts={"cmap": "viridis"},
)
ax1.set_title("Model slice at y = 0 m")

ax2 = fig.add_axes([0.85, 0.1, 0.05, 0.8])
norm = mpl.colors.Normalize(vmin=np.min(recovered_model), vmax=np.max(recovered_model))
cbar = mpl.colorbar.ColorbarBase(
    ax2, norm=norm, orientation="vertical", cmap=mpl.cm.viridis
)
cbar.set_label("$g/cm^3$", rotation=270, labelpad=15, size=12)

plt.show()
  • Model slice at y = 0 m
  • Model slice at y = 0 m

Plotting Predicted Data and Normalized Misfit#

# Predicted data with final recovered model
# SimPEG uses right handed coordinate where Z is positive upward.
# This causes gravity signals look "inconsistent" with density values in visualization.
dpred = inv_prob.dpred

# Observed data | Predicted data | Normalized data misfit
data_array = np.c_[dobs, dpred, (dobs - dpred) / uncertainties]

fig = plt.figure(figsize=(17, 4))
plot_title = ["Observed", "Predicted", "Normalized Misfit"]
plot_units = ["mgal", "mgal", ""]

ax1 = 3 * [None]
ax2 = 3 * [None]
norm = 3 * [None]
cbar = 3 * [None]
cplot = 3 * [None]
v_lim = [np.max(np.abs(dobs)), np.max(np.abs(dobs)), np.max(np.abs(data_array[:, 2]))]

for ii in range(0, 3):
    ax1[ii] = fig.add_axes([0.33 * ii + 0.03, 0.11, 0.23, 0.84])
    cplot[ii] = plot2Ddata(
        receiver_list[0].locations,
        data_array[:, ii],
        ax=ax1[ii],
        ncontour=30,
        clim=(-v_lim[ii], v_lim[ii]),
        contourOpts={"cmap": "bwr"},
    )
    ax1[ii].set_title(plot_title[ii])
    ax1[ii].set_xlabel("x (m)")
    ax1[ii].set_ylabel("y (m)")

    ax2[ii] = fig.add_axes([0.33 * ii + 0.25, 0.11, 0.01, 0.85])
    norm[ii] = mpl.colors.Normalize(vmin=-v_lim[ii], vmax=v_lim[ii])
    cbar[ii] = mpl.colorbar.ColorbarBase(
        ax2[ii], norm=norm[ii], orientation="vertical", cmap=mpl.cm.bwr
    )
    cbar[ii].set_label(plot_units[ii], rotation=270, labelpad=15, size=12)

plt.show()
Observed, Predicted, Normalized Misfit

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

Estimated memory usage: 8 MB

Gallery generated by Sphinx-Gallery