2.5D DC Resistivity Inversion with Sparse Norms#

Here we invert a line of DC resistivity data to recover an electrical conductivity model. We formulate the inverse problem as a least-squares optimization problem. For this tutorial, we focus on the following:

  • Defining the survey

  • Generating a mesh based on survey geometry

  • Including surface topography

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

  • Applying sensitivity weighting

  • Plotting the recovered model and data misfit

Import modules#

import os
import numpy as np
import matplotlib as mpl
import matplotlib.pyplot as plt
from matplotlib.colors import LogNorm
import tarfile

from discretize import TreeMesh
from discretize.utils import mkvc, active_from_xyz

from simpeg.utils import model_builder
from simpeg import (
    maps,
    data_misfit,
    regularization,
    optimization,
    inverse_problem,
    inversion,
    directives,
    utils,
)
from simpeg.electromagnetics.static import resistivity as dc
from simpeg.electromagnetics.static.utils.static_utils import (
    plot_pseudosection,
    apparent_resistivity_from_voltage,
)
from simpeg.utils.io_utils.io_utils_electromagnetics import read_dcip2d_ubc

mpl.rcParams.update({"font.size": 16})
# sphinx_gallery_thumbnail_number = 3

Define File Names#

Here we provide the file paths to assets we need to run the inversion. The path to the true model conductivity and chargeability models are also provided for comparison with the inversion results. These files are stored as a tar-file on our google cloud bucket: “https://storage.googleapis.com/simpeg/doc-assets/dcr2d.tar.gz

# storage bucket where we have the data
data_source = "https://storage.googleapis.com/simpeg/doc-assets/dcr2d.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 + "topo_xyz.txt"
data_filename = dir_path + "dc_data.obs"
overwriting /home/vsts/work/1/s/tutorials/05-dcr/dcr2d.tar.gz
Downloading https://storage.googleapis.com/simpeg/doc-assets/dcr2d.tar.gz
   saved to: /home/vsts/work/1/s/tutorials/05-dcr/dcr2d.tar.gz
Download completed!

Load Data, Define Survey and Plot#

Here we load the observed data, define the DC and IP survey geometry and plot the data values using pseudo-sections. Warning: In the following example, the observations file is assumed to be sorted by sources

# Load data
topo_xyz = np.loadtxt(str(topo_filename))
dc_data = read_dcip2d_ubc(data_filename, "volt", "general")

Plot Observed Data in Pseudo-Section#

Here, we demonstrate how to plot 2D data in pseudo-section. First, we plot the actual data (voltages) in pseudo-section as a scatter plot. This allows us to visualize the pseudo-sensitivity locations for our survey. Next, we plot the data as apparent conductivities in pseudo-section with a filled contour plot.

# Plot voltages pseudo-section
fig = plt.figure(figsize=(12, 5))
ax1 = fig.add_axes([0.1, 0.15, 0.75, 0.78])
plot_pseudosection(
    dc_data,
    plot_type="scatter",
    ax=ax1,
    scale="log",
    cbar_label="V/A",
    scatter_opts={"cmap": mpl.cm.viridis},
)
ax1.set_title("Normalized Voltages")
plt.show()

# Get apparent conductivities from volts and survey geometry
apparent_conductivities = 1 / apparent_resistivity_from_voltage(
    dc_data.survey, dc_data.dobs
)

# Plot apparent conductivity pseudo-section
fig = plt.figure(figsize=(12, 5))
ax1 = fig.add_axes([0.1, 0.15, 0.75, 0.78])
plot_pseudosection(
    dc_data.survey,
    apparent_conductivities,
    plot_type="contourf",
    ax=ax1,
    scale="log",
    cbar_label="S/m",
    mask_topography=True,
    contourf_opts={"levels": 20, "cmap": mpl.cm.viridis},
)
ax1.set_title("Apparent Conductivity")
plt.show()
  • Normalized Voltages
  • Apparent Conductivity

Assign Uncertainties#

Inversion with SimPEG requires that we define the uncertainties on our data. This represents our estimate of the standard deviation of the noise in our data. For DC data, the uncertainties are 10% of the absolute value.

dc_data.standard_deviation = 0.05 * np.abs(dc_data.dobs)

Create Tree Mesh#

Here, we create the Tree mesh that will be used invert the DC data

dh = 4  # base cell width
dom_width_x = 3200.0  # domain width x
dom_width_z = 2400.0  # domain width z
nbcx = 2 ** int(np.round(np.log(dom_width_x / dh) / np.log(2.0)))  # num. base cells x
nbcz = 2 ** int(np.round(np.log(dom_width_z / dh) / np.log(2.0)))  # num. base cells z

# Define the base mesh
hx = [(dh, nbcx)]
hz = [(dh, nbcz)]
mesh = TreeMesh([hx, hz], x0="CN")

# Mesh refinement based on topography
mesh.refine_surface(
    topo_xyz[:, [0, 2]],
    padding_cells_by_level=[0, 0, 4, 4],
    finalize=False,
)

# Mesh refinement near transmitters and receivers. First we need to obtain the
# set of unique electrode locations.
electrode_locations = np.c_[
    dc_data.survey.locations_a,
    dc_data.survey.locations_b,
    dc_data.survey.locations_m,
    dc_data.survey.locations_n,
]

unique_locations = np.unique(
    np.reshape(electrode_locations, (4 * dc_data.survey.nD, 2)), axis=0
)

mesh.refine_points(unique_locations, padding_cells_by_level=[4, 4], finalize=False)

# Refine core mesh region
xp, zp = np.meshgrid([-600.0, 600.0], [-400.0, 0.0])
xyz = np.c_[mkvc(xp), mkvc(zp)]
mesh.refine_bounding_box(xyz, padding_cells_by_level=[0, 0, 2, 8], finalize=False)

mesh.finalize()

Project Surveys to Discretized Topography#

It is important that electrodes are not model as being in the air. Even if the electrodes are properly located along surface topography, they may lie above the discretized topography. This step is carried out to ensure all electrodes like on the discretized surface.

# Create 2D topography. Since our 3D topography only changes in the x direction,
# it is easy to define the 2D topography projected along the survey line. For
# arbitrary topography and for an arbitrary survey orientation, the user must
# define the 2D topography along the survey line.
topo_2d = np.unique(topo_xyz[:, [0, 2]], axis=0)

# Find cells that lie below surface topography
ind_active = active_from_xyz(mesh, topo_2d)

# Extract survey from data object
survey = dc_data.survey

# Shift electrodes to the surface of discretized topography
survey.drape_electrodes_on_topography(mesh, ind_active, option="top")

# Reset survey in data object
dc_data.survey = survey

Starting/Reference Model and Mapping on Tree Mesh#

Here, we would create starting and/or reference models for the DC 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. Here, the starting model is the natural log of 0.01 S/m.

# Define conductivity model in S/m (or resistivity model in Ohm m)
air_conductivity = np.log(1e-8)
background_conductivity = np.log(1e-2)

active_map = maps.InjectActiveCells(mesh, ind_active, np.exp(air_conductivity))
nC = int(ind_active.sum())

conductivity_map = active_map * maps.ExpMap()

# Define model
starting_conductivity_model = background_conductivity * np.ones(nC)

Define the Physics of the DC Simulation#

Here, we define the physics of the DC resistivity problem.

# Define the problem. Define the cells below topography and the mapping
simulation = dc.simulation_2d.Simulation2DNodal(
    mesh, survey=survey, sigmaMap=conductivity_map, storeJ=True
)

Define DC 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=dc_data, simulation=simulation)

# Define the regularization (model objective function). Here, 'p' defines the
# the norm of the smallness term, 'qx' defines the norm of the smoothness
# in x and 'qz' defines the norm of the smoothness in z.
regmap = maps.IdentityMap(nP=int(ind_active.sum()))

reg = regularization.Sparse(
    mesh,
    active_cells=ind_active,
    reference_model=starting_conductivity_model,
    mapping=regmap,
    gradient_type="total",
    alpha_s=0.01,
    alpha_x=1,
    alpha_y=1,
)

reg.reference_model_in_smooth = True  # Include reference model in smoothness

p = 0
qx = 1
qz = 1
reg.norms = [p, qx, qz]

# Define how the optimization problem is solved. Here we will use an inexact
# Gauss-Newton approach.
opt = optimization.InexactGaussNewton(maxIter=40)

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

Define DC Inversion Directives#

Here we define any directives 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.

# Apply and update sensitivity weighting as the model updates
update_sensitivity_weighting = directives.UpdateSensitivityWeights()

# Reach target misfit for L2 solution, then use IRLS until model stops changing.
update_IRLS = directives.UpdateIRLS(max_irls_iterations=25, chifact_start=1.0)

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

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

# Update preconditioner
update_jacobi = directives.UpdatePreconditioner()

directives_list = [
    update_sensitivity_weighting,
    update_IRLS,
    starting_beta,
    save_iteration,
    update_jacobi,
]

Running the DC 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
dc_inversion = inversion.BaseInversion(inv_prob, directiveList=directives_list)

# Run inversion
recovered_conductivity_model = dc_inversion.run(starting_conductivity_model)
Running inversion with SimPEG v0.23.0
/home/vsts/work/1/s/simpeg/simulation.py:197: DefaultSolverWarning:

Using the default solver: Pardiso.

If you would like to suppress this notification, add
warnings.filterwarnings('ignore', simpeg.utils.solver_utils.DefaultSolverWarning)
 to your script.


                        simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
                        ***Done using same Solver, and solver_opts as the Simulation2DNodal problem***

/usr/share/miniconda/envs/simpeg-test/lib/python3.10/site-packages/pymatsolver/solvers.py:415: FutureWarning:

In Future pymatsolver v0.4.0, passing a vector of shape (n, 1) to the solve method will return an array with shape (n, 1), instead of always returning a flattened array. This is to be consistent with numpy.linalg.solve broadcasting.

model has any nan: 0
============================ Inexact Gauss Newton ============================
  #     beta     phi_d     phi_m       f      |proj(x-g)-x|  LS    Comment
-----------------------------------------------------------------------------
x0 has any nan: 0
   0  1.30e+03  3.42e+04  0.00e+00  3.42e+04    2.90e+03      0
   1  6.50e+02  4.58e+03  3.13e+00  6.61e+03    4.47e+02      0
   2  3.25e+02  1.54e+03  5.65e+00  3.37e+03    1.41e+02      0   Skip BFGS
   3  1.62e+02  6.61e+02  7.35e+00  1.85e+03    8.59e+01      0   Skip BFGS
   4  8.12e+01  2.79e+02  8.94e+00  1.01e+03    4.76e+01      0   Skip BFGS
Reached starting chifact with l2-norm regularization: Start IRLS steps...
irls_threshold 2.1967395367573377
   5  8.12e+01  1.25e+02  1.53e+01  1.37e+03    3.87e+01      0   Skip BFGS
   6  4.66e+01  2.01e+02  1.65e+01  9.69e+02    3.06e+01      0
   7  6.12e+01  1.04e+02  2.03e+01  1.35e+03    4.64e+01      0
   8  3.60e+01  1.89e+02  2.18e+01  9.73e+02    3.21e+01      0
   9  4.40e+01  1.20e+02  2.56e+01  1.24e+03    5.41e+01      0
  10  2.55e+01  1.95e+02  2.58e+01  8.52e+02    3.27e+01      0
  11  3.19e+01  1.14e+02  2.67e+01  9.67e+02    3.87e+01      0
  12  4.00e+01  1.48e+02  2.55e+01  1.17e+03    5.25e+01      0
  13  2.28e+01  2.04e+02  2.44e+01  7.60e+02    4.22e+01      0
  14  2.94e+01  1.08e+02  2.58e+01  8.68e+02    4.78e+01      0
  15  3.79e+01  1.39e+02  2.46e+01  1.07e+03    6.30e+01      0
  16  2.22e+01  1.91e+02  2.33e+01  7.10e+02    4.50e+01      0
  17  2.83e+01  1.11e+02  2.46e+01  8.08e+02    5.16e+01      0
  18  3.60e+01  1.40e+02  2.32e+01  9.75e+02    6.73e+01      0
  19  2.15e+01  1.83e+02  2.16e+01  6.48e+02    5.01e+01      0
  20  2.65e+01  1.17e+02  2.24e+01  7.12e+02    5.12e+01      0
  21  3.28e+01  1.39e+02  2.10e+01  8.27e+02    6.74e+01      0
  22  2.00e+01  1.71e+02  1.97e+01  5.66e+02    4.96e+01      0
  23  2.43e+01  1.21e+02  2.03e+01  6.14e+02    4.83e+01      0
  24  2.95e+01  1.36e+02  1.87e+01  6.88e+02    6.25e+01      0
  25  1.86e+01  1.57e+02  1.74e+01  4.81e+02    5.06e+01      0
  26  2.31e+01  1.17e+02  1.76e+01  5.23e+02    5.03e+01      0
  27  2.86e+01  1.25e+02  1.64e+01  5.93e+02    6.61e+01      0
  28  3.54e+01  1.40e+02  1.53e+01  6.82e+02    8.53e+01      0
  29  2.22e+01  1.61e+02  1.41e+01  4.74e+02    5.20e+01      0
Reach maximum number of IRLS cycles: 25
------------------------- STOP! -------------------------
1 : |fc-fOld| = 0.0000e+00 <= tolF*(1+|f0|) = 3.4201e+03
1 : |xc-x_last| = 2.4280e+00 <= tolX*(1+|x0|) = 8.2347e+01
0 : |proj(x-g)-x|    = 5.2025e+01 <= tolG          = 1.0000e-01
0 : |proj(x-g)-x|    = 5.2025e+01 <= 1e3*eps       = 1.0000e-02
0 : maxIter   =      40    <= iter          =     30
------------------------- DONE! -------------------------

Recreate True Conductivity Model#

Plotting True and Recovered Conductivity Model#

# Get L2 and sparse recovered model in base 10
l2_conductivity = conductivity_map * inv_prob.l2model
l2_conductivity[~ind_active] = np.nan

recovered_conductivity = conductivity_map * recovered_conductivity_model
recovered_conductivity[~ind_active] = np.nan

# Plot True Model
norm = LogNorm(vmin=1e-3, vmax=1e-1)

fig = plt.figure(figsize=(9, 15))
ax1 = 3 * [None]
ax2 = 3 * [None]
title_str = [
    "True Conductivity Model",
    "Smooth Recovered Model",
    "Sparse Recovered Model",
]
plotting_model = [
    true_conductivity_model,
    l2_conductivity,
    recovered_conductivity,
]

for ii in range(0, 3):
    ax1[ii] = fig.add_axes([0.14, 0.75 - 0.3 * ii, 0.68, 0.2])
    mesh.plot_image(
        plotting_model[ii],
        ax=ax1[ii],
        grid=False,
        range_x=[-700, 700],
        range_y=[-600, 0],
        pcolor_opts={"norm": norm},
    )
    ax1[ii].set_xlim(-600, 600)
    ax1[ii].set_ylim(-600, 0)
    ax1[ii].set_title(title_str[ii])
    ax1[ii].set_xlabel("x (m)")
    ax1[ii].set_ylabel("z (m)")

    ax2[ii] = fig.add_axes([0.84, 0.75 - 0.3 * ii, 0.03, 0.2])
    cbar = mpl.colorbar.ColorbarBase(ax2[ii], norm=norm, orientation="vertical")
    cbar.set_label(r"$\sigma$ (S/m)", rotation=270, labelpad=15, size=12)

plt.show()
True Conductivity Model, Smooth Recovered Model, Sparse Recovered Model

Plotting Predicted DC Data and Misfit#

# Predicted data from recovered model
dpred = inv_prob.dpred
dobs = dc_data.dobs
std = dc_data.standard_deviation

# Plot
fig = plt.figure(figsize=(9, 13))
data_array = [np.abs(dobs), np.abs(dpred), (dobs - dpred) / std]
plot_title = ["Observed Voltage", "Predicted Voltage", "Normalized Misfit"]
plot_units = ["V/A", "V/A", ""]
scale = ["log", "log", "linear"]

ax1 = 3 * [None]
cax1 = 3 * [None]
cbar = 3 * [None]
cplot = 3 * [None]

for ii in range(0, 3):
    ax1[ii] = fig.add_axes([0.15, 0.72 - 0.33 * ii, 0.65, 0.21])
    cax1[ii] = fig.add_axes([0.81, 0.72 - 0.33 * ii, 0.03, 0.21])
    cplot[ii] = plot_pseudosection(
        survey,
        data_array[ii],
        "contourf",
        ax=ax1[ii],
        cax=cax1[ii],
        scale=scale[ii],
        cbar_label=plot_units[ii],
        mask_topography=True,
        contourf_opts={"levels": 25, "cmap": mpl.cm.viridis},
    )
    ax1[ii].set_title(plot_title[ii])

plt.show()
Observed Voltage, Predicted Voltage, Normalized Misfit

Total running time of the script: (17 minutes 38.330 seconds)

Estimated memory usage: 309 MB

Gallery generated by Sphinx-Gallery