Surface wave and receiver function - joint inversion (field data)#

Open In Colab

If you are running this notebook locally, make sure you’ve followed steps here to set up the environment. (This environment.yml file specifies a list of packages required to run the notebooks)


What we do in this notebook#

Here we extend the example where CoFI has been used for the inversion of Rayleigh wave phase velocities for a 1D layered earth to a joint inversion where we also account for receiver functions.

In this example, we extend from the synthetic data notebook and apply the joint inversion on a field data presented in the Computer Programs in Seismology (see: https://www.eas.slu.edu/eqc/eqc_cps/TUTORIAL/STRUCT/index.html).

Learning outcomes

  • A demonstration of CoFI’s ability to switch between parameter estimation and ensemble methods.

  • An application of CoFI for a joint inversion, here of Rayleigh wave pahse velocity and receiver function data, to a field dataset

Utilities preparation#

In this example, we look at a joint inversion problem of surface wave and receiver function.

We use pysurf96 for computing the forward step of surface wave, and use pyhk for receiver function calculations.

# -------------------------------------------------------- #
#                                                          #
#     Uncomment below to set up environment on "colab"     #
#                                                          #
# -------------------------------------------------------- #

# !pip install -U cofi
# !pip install git+https://github.com/inlab-geo/pysurf96.git
# !git clone https://github.com/inlab-geo/cofi-examples.git
# %cd cofi-examples/examples/sw_rf_joint
# If this notebook is run locally pysurf96 needs to be installed separately by uncommenting the following line,
# that is by removing the # and the white space between it and the exclammation mark.
# !pip install -U cofi git+https://github.com/inlab-geo/pysurf96.git@ctypes
import glob
import numpy as np
import scipy
import matplotlib.pyplot as plt
import matplotlib.gridspec as gridspec
import arviz as az

import pysurf96
import bayesbay
import cofi
import pyhk
np.seterr(all="ignore");
{'divide': 'warn', 'over': 'warn', 'under': 'ignore', 'invalid': 'warn'}

Model vector

# display theory on the 1D model parameterisation
from IPython.display import display, Markdown

with open("../../theory/misc_1d_model_parameterisation.md", "r") as f:
    content = f.read()

display(Markdown(content))
<IPython.core.display.Markdown object>
# layercake model utilities
def form_layercake_model(thicknesses, vs):
    model = np.zeros((len(vs)*2-1,))
    model[1::2] = thicknesses
    model[::2] = vs
    return model

def split_layercake_model(model):
    thicknesses = model[1::2]
    vs = model[::2]
    return thicknesses, vs
# voronoi model utilities
def form_voronoi_model(voronoi_sites, vs):
    return np.hstack((vs, voronoi_sites))

def split_voronoi_model(model):
    voronoi_sites = model[len(model)//2:]
    vs = model[:len(model)//2]
    return voronoi_sites, vs
def voronoi_to_layercake(voronoi_vector: np.ndarray) -> np.ndarray:
    n_layers = len(voronoi_vector) // 2
    velocities = voronoi_vector[:n_layers]
    voronoi_sites = voronoi_vector[n_layers:]
    depths = (voronoi_sites[:-1] + voronoi_sites[1:]) / 2
    thicknesses = depths - np.insert(depths[:-1], 0, 0)
    layercake_vector = np.zeros((2*n_layers-1,))
    layercake_vector[::2] = velocities
    layercake_vector[1::2] = thicknesses
    return layercake_vector

def layercake_to_voronoi(layercake_vector: np.ndarray, first_voronoi_site: float = 0.0) -> np.ndarray:
    n_layers = len(layercake_vector) // 2 + 1
    thicknesses = layercake_vector[1::2]
    velocities = layercake_vector[::2]
    depths = np.cumsum(thicknesses)
    voronoi_sites = np.zeros((n_layers,))
    for i in range(1,len(voronoi_sites)):
        voronoi_sites[i] = 2 * depths[i-1] - voronoi_sites[i-1]
    voronoi_vector = np.hstack((velocities, voronoi_sites))
    return voronoi_vector

Interfacing to pysurf96

# display theory on the using the forward solver
with open("../../theory/geo_surface_wave_dispersion2.md", "r") as f:
    content = f.read()

display(Markdown(content))
<IPython.core.display.Markdown object>
VP_VS = 1.77
RHO_VP_K = 0.32
RHO_VP_B = 0.77
# forward through pysurf96
def forward_sw(model, periods):
    thicknesses, vs = split_layercake_model(model)
    thicknesses = np.append(thicknesses, 10)
    vp = vs * VP_VS
    rho = RHO_VP_K * vp + RHO_VP_B
    return pysurf96.surf96(
        thicknesses,
        vp,
        vs,
        rho,
        periods,
        wave="rayleigh",
        mode=1,
        velocity="phase",
        flat_earth=False,
    )
def forward_sw_interp(model, periods):
    pysurf_periods = np.logspace(
        np.log(np.min(periods)),
        np.log(np.max(periods+1)),
        60,
        base=np.e,
    )
    pysurf_dpred = forward_sw(model, pysurf_periods)
    interp_func = scipy.interpolate.interp1d(pysurf_periods,
                                             pysurf_dpred,
                                             fill_value="extrapolate")
    dpred = interp_func(periods)
    return dpred
t_shift = 5
t_duration = 25
t_sampling_interval = 0.1
gauss = 1
ray_param_s_km = 0.07
# forward through rf in pyhf
def forward_rf(
    model,
    t_shift=t_shift,
    t_duration=t_duration,
    t_sampling_interval=t_sampling_interval,
    gauss=gauss,
    ray_param_s_km=ray_param_s_km,
    return_times=False
):
    h, vs = split_layercake_model(model)
    data = pyhk.rfcalc(ps=0, thik=h, beta=vs, kapa=np.ones((len(vs),))*VP_VS, p=ray_param_s_km,
                      duration=t_duration, dt=t_sampling_interval, shft=t_shift, gauss=gauss)
    if return_times:
        times = np.arange(len(data)) * t_sampling_interval - t_shift
        return data, times
    else:
        return data
t_shift = 5
t_duration = 70
t_sampling_interval = 0.5

def forward_rf_interp(model, gauss, ray_param):
    dpred, times = forward_rf(model, t_shift, t_duration, t_sampling_interval,
                       gauss=gauss, ray_param_s_km=ray_param, return_times=True)
    interp_func = scipy.interpolate.interp1d(times, dpred, fill_value="extrapolate")
    return interp_func(rf_field_times)

Numerical Jacobian

def jacobian(model, data_length, fwd, fwd_kwargs, relative_step=0.01):
    jac = np.zeros((data_length, len(model)))
    original_dpred = fwd(model, **fwd_kwargs)
    for i in range(len(model)):
        perturbed_model = model.copy()
        step = relative_step * model[i]
        perturbed_model[i] += step
        perturbed_dpred = fwd(perturbed_model, **fwd_kwargs)
        derivative = (perturbed_dpred - original_dpred) / step
        jac[:, i] = derivative
    return jac

Plotting

def plot_model(model, ax=None, title="model", **kwargs):
    # process data
    thicknesses = np.append(model[1::2], max(model[1::2]))
    velocities = model[::2]
    y = np.insert(np.cumsum(thicknesses), 0, 0)
    x = np.insert(velocities, 0, velocities[0])

    # plot depth profile
    if ax is None:
        _, ax = plt.subplots()
    plotting_style = {
        "linewidth": kwargs.pop("linewidth", kwargs.pop("lw", 0.5)),
        "alpha": 0.2,
        "color": kwargs.pop("color", kwargs.pop("c", "blue")),
    }
    plotting_style.update(kwargs)
    ax.step(x, y, where="post", **plotting_style)
    if ax.get_ylim()[0] < ax.get_ylim()[1]:
        ax.invert_yaxis()
    ax.set_xlabel("Vs (km/s)")
    ax.set_ylabel("Depth (km)")
    ax.set_title(title)
    return ax
def plot_data(x, y, ax=None, scatter=False, xlabel=None, ylabel=None,
              title="surface wave data", **kwargs):
    if ax is None:
        _, ax = plt.subplots()
    plotting_style = {
        "linewidth": kwargs.pop("linewidth", kwargs.pop("lw", 1)),
        "alpha": 1,
        "color": kwargs.pop("color", kwargs.pop("c", "blue")),
    }
    plotting_style.update(**kwargs)
    if scatter:
        ax.scatter(x, y, **plotting_style)
    else:
        ax.plot(x, y, **plotting_style)
    ax.set_xlabel(xlabel)
    ax.set_ylabel(ylabel)
    ax.set_title(title)
    return ax
def plot_sw_data(rayleigh_phase_velocities, periods, ax=None, scatter=False,
                 title="surface wave data", **kwargs):
    return plot_data(x=periods,
                     y=rayleigh_phase_velocities,
                     ax=ax,
                     scatter=scatter,
                     title=title,
                     xlabel="Periods (s)",
                     ylabel="Rayleigh phase velocities (km/s)",
                     **kwargs)

def plot_rf_data(rf_data, rf_times, ax=None, scatter=False,
                 title="receiver function data", **kwargs):
    return plot_data(x=rf_times,
                     y=rf_data,
                     ax=ax,
                     scatter=scatter,
                     title=title,
                     xlabel="Times (s)",
                     ylabel="Receiver function amplitudes",
                     **kwargs)

Read data#

Rayleigh observations

file_surf_data = "../../data/sw_rf_joint/data/SURF/nnall.dsp"

with open(file_surf_data, "r") as file:
    lines = file.readlines()
    surf_data = []
    for line in lines:
        row = line.strip().split()
        if "C" in row:
            surf_data.append([float(e) for e in row[5:8]])

field_d = np.array(surf_data)
field_d_periods = field_d[:,0]
rayleigh_field_d_obs = field_d[:,1]
ax = plot_sw_data(rayleigh_field_d_obs, field_d_periods, color="brown", s=5, scatter=True,
             label="d_obs")
ax.legend();
<matplotlib.legend.Legend object at 0x7efdbd2a0830>

Receiver function observations

files_rftn_data = glob.glob("../../data/sw_rf_joint/data/RFTN/rf_*.txt")

all_gauss = []
all_ray_params = []
all_rf_field_dobs = []
rf_field_times = None

for file in files_rftn_data:
    gauss, ray_param_s_km = file.split("_")[-2:]
    all_gauss.append(float(gauss))
    all_ray_params.append(float(ray_param_s_km[:-4]))
    dataset = np.loadtxt(file)
    if rf_field_times is None:
        rf_field_times = dataset[:, 0]
    all_rf_field_dobs.append(dataset[:, 1])
def plot_rf_field_data(all_rf_data, all_gauss, all_ray_params, rf_field_times,
                       axes=None, color="darkblue", label=None, **kwargs):
    if axes is None:
        _, axes = plt.subplots(13, 3, figsize=(10, 12))
    for i, ax in enumerate(axes.flat):
        ax.set_xlim(-5, 20)
        gauss = all_gauss[i]
        ray_param = all_ray_params[i]
        rf_dobs = all_rf_data[i]
        title=f"ray (s/km) = {ray_param}, gauss = {gauss}"
        plot_rf_data(rf_dobs, rf_field_times, ax=ax, color=color,
                     title=title, label=label, **kwargs)
        ax.set_ylabel("Amplitude")

    for ax in axes[:-1, :].flat:
        ax.set_xlabel("")
        ax.tick_params(labelbottom=False)
    for ax in axes[:, 1:].flat:
        ax.set_ylabel("")
    plt.tight_layout()
    return axes
plot_rf_field_data(all_rf_field_dobs, all_gauss, all_ray_params, rf_field_times);
array([[<Axes: title={'center': 'ray (s/km) = 0.0716, gauss = 1.0'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.0751, gauss = 1.0'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0705, gauss = 1.0'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0713, gauss = 2.5'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.0738, gauss = 1.0'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0658, gauss = 2.5'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0698, gauss = 1.0'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.0665, gauss = 1.0'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0698, gauss = 2.5'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0746, gauss = 1.0'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.076, gauss = 1.0'}>,
        <Axes: title={'center': 'ray (s/km) = 0.069, gauss = 2.5'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0687, gauss = 2.5'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.0698, gauss = 2.5'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0746, gauss = 1.0'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0705, gauss = 2.5'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.069, gauss = 1.0'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0724, gauss = 1.0'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0658, gauss = 1.0'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.0713, gauss = 1.0'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0724, gauss = 2.5'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0687, gauss = 1.0'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.0665, gauss = 2.5'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0716, gauss = 2.5'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0787, gauss = 1.0'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.07, gauss = 1.0'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0732, gauss = 1.0'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0724, gauss = 1.0'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.0704, gauss = 1.0'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0724, gauss = 2.5'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0698, gauss = 1.0'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.0746, gauss = 2.5'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0739, gauss = 2.5'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0739, gauss = 1.0'}, ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.0738, gauss = 2.5'}>,
        <Axes: title={'center': 'ray (s/km) = 0.0716, gauss = 1.0'}>],
       [<Axes: title={'center': 'ray (s/km) = 0.0751, gauss = 2.5'}, xlabel='Times (s)', ylabel='Amplitude'>,
        <Axes: title={'center': 'ray (s/km) = 0.076, gauss = 2.5'}, xlabel='Times (s)'>,
        <Axes: title={'center': 'ray (s/km) = 0.07, gauss = 2.5'}, xlabel='Times (s)'>]],
      dtype=object)
# stacking the data into gauss=1.0 and gauss=2.5 groups
all_rf_field_dobs_1_0 = []
all_rf_field_dobs_2_5 = []
ray_params_1_0 = []
ray_params_2_5 = []
for gauss, ray_param, dobs in zip(all_gauss, all_ray_params, all_rf_field_dobs):
    if gauss == 1:
        all_rf_field_dobs_1_0.append(dobs)
        ray_params_1_0.append(ray_param)
    else:
        all_rf_field_dobs_2_5.append(dobs)
        ray_params_2_5.append(ray_param)

rf_field_dobs_1_0 = np.mean(all_rf_field_dobs_1_0, axis=0)
rf_field_dobs_2_5 = np.mean(all_rf_field_dobs_2_5, axis=0)
ray_param_1_0 = np.mean(ray_params_1_0)
ray_param_2_5 = np.mean(ray_params_2_5)
_, axes = plt.subplots(2, 1, figsize=(5,4))

axes[0].set_xlim(-5, 20)
axes[1].set_xlim(-5, 20)

plot_rf_data(
    rf_field_dobs_1_0,
    rf_field_times,
    title=f"ray (s/km) = {ray_param_1_0}, gauss = 1.0",
    ax=axes[0]
)
axes[0].set_ylabel("Amplitude")
plot_rf_data(
    rf_field_dobs_2_5,
    rf_field_times,
    title=f"ray (s/km) = {ray_param_2_5}, gauss = 2.5",
    ax=axes[1]
)
axes[1].set_ylabel("Amplitude")
plt.tight_layout()

Reference good model

file_end_mod = "../../data/sw_rf_joint/data/JOINT/end.mod"

with open(file_end_mod, "r") as file:
    lines = file.readlines()
    ref_good_model = []
    for line in lines[12:]:
        row = line.strip().split()
        ref_good_model.append([float(row[0]), float(row[2])])

ref_good_model = np.array(ref_good_model)
ref_good_model = form_layercake_model(ref_good_model[:-1,0], ref_good_model[:,1])
_, ax = plt.subplots(figsize=(4,6))
plot_model(ref_good_model, ax=ax, alpha=1);
<Axes: title={'center': 'model'}, xlabel='Vs (km/s)', ylabel='Depth (km)'>

Optimisation#

n_dims = 29

init_thicknesses = np.ones((n_dims//2,)) * 10
init_vs = np.ones((n_dims//2+1,)) * 4.0
init_model = form_layercake_model(init_thicknesses, init_vs)
my_reg = cofi.utils.QuadraticReg(
    weighting_matrix="damping",
    model_shape=(n_dims,),
    reference_model=init_model
)
def my_objective(model, fwd_funcs, d_obs_list, lamda=1.0):
    data_misfit = 0
    for (fwd, fwd_kwargs), d_obs in zip(fwd_funcs, d_obs_list):
        d_pred = fwd(model, **fwd_kwargs)
        data_misfit += np.sum((d_obs - d_pred) ** 2)
    reg = my_reg(model)
    return data_misfit + lamda * reg

def my_objective_gradient(model, fwd_funcs, d_obs_list, lamda=1.0):
    data_misfit_grad = 0
    for (fwd, fwd_kwargs), d_obs in zip(fwd_funcs, d_obs_list):
        d_pred = fwd(model, **fwd_kwargs)
        jac = jacobian(model, len(d_obs), fwd, fwd_kwargs)
        data_misfit_grad += -2 * jac.T @ (d_obs - d_pred)
    reg_grad = my_reg.gradient(model)
    return data_misfit_grad + lamda * reg_grad

def my_objective_hessian(model, fwd_funcs, d_obs_list, lamda=1.0):
    data_misfit_hess = 0
    for (fwd, fwd_kwargs), d_obs in zip(fwd_funcs, d_obs_list):
        jac = jacobian(model, len(d_obs), fwd, fwd_kwargs)
        data_misfit_hess += 2 * jac.T @ jac
    reg_hess = my_reg.hessian(model)
    return data_misfit_hess + lamda * reg_hess
fwd_funcs = [
    (forward_sw_interp, {"periods": field_d_periods}),
    (forward_rf_interp, {"gauss": 1, "ray_param": ray_param_1_0}),
    (forward_rf_interp, {"gauss": 2.5, "ray_param": ray_param_2_5})
]

d_obs_list = [rayleigh_field_d_obs, rf_field_dobs_1_0, rf_field_dobs_2_5]

Optimisation with no damping#

lamda = 0

kwargs = {
    "fwd_funcs": fwd_funcs,
    "d_obs_list": d_obs_list,
    "lamda": lamda
}
joint_field_problem_no_reg = cofi.BaseProblem()
joint_field_problem_no_reg.set_objective(my_objective, kwargs=kwargs)
joint_field_problem_no_reg.set_gradient(my_objective_gradient, kwargs=kwargs)
joint_field_problem_no_reg.set_hessian(my_objective_hessian, kwargs=kwargs)
joint_field_problem_no_reg.set_initial_model(init_model)

Define ``InversionOptions``

inv_options_optimiser = cofi.InversionOptions()
inv_options_optimiser.set_tool("scipy.optimize.minimize")
inv_options_optimiser.set_params(method="trust-exact")

Define ``Inversion`` and run

inv_optimiser_no_reg_field = cofi.Inversion(joint_field_problem_no_reg, inv_options_optimiser)
inv_res_optimiser_field_no_reg = inv_optimiser_no_reg_field.run()
/home/kitc/actions-runner/_work/cofi/cofi/src/cofi/tools/_scipy_opt_min.py:120: RuntimeWarning: Method trust-exact does not use Hessian-vector product information (hessp).
  return minimize(

Plot results

_, (ax1, ax2) = plt.subplots(1, 2, figsize=(10, 4), gridspec_kw={'width_ratios': [1, 2.5]})

ax1.set_ylim(100)

plot_model(inv_res_optimiser_field_no_reg.model, ax=ax1, color="green", alpha=1,
           label="model inverted from field data")
plot_model(ref_good_model, ax=ax1, color="red", alpha=1,
           label="reference good model")
plot_model(init_model, ax=ax1, alpha=1, lw=1.5, color="purple",
           label="initial model for damped solution")

field_d_periods_logspace = np.logspace(
    np.log(np.min(field_d_periods)),
    np.log(np.max(field_d_periods+1)),
    60,
    base=np.e,
)

plot_sw_data(rayleigh_field_d_obs, field_d_periods, ax=ax2, color="orange", s=5, scatter=True,
             label="d_obs")
plot_sw_data(forward_sw_interp(inv_res_optimiser_field_no_reg.model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2, color="green",
             label="d_pred from inverted model")
plot_sw_data(forward_sw_interp(ref_good_model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2, color="red",
             label="d_pred from reference good model")
plot_sw_data(forward_sw_interp(init_model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2,
             alpha=1, lw=1.5, linestyle="--", color="purple",
             label="d_pred from initial model")

ax1.legend(loc="lower center", bbox_to_anchor=(0.5, -0.4))
ax2.legend(loc="lower center", bbox_to_anchor=(0.5, -0.46));
<matplotlib.legend.Legend object at 0x7efd9b8c9950>
all_rf_dpred = []
all_rf_dpred_init_m = []

for gauss, ray_param in zip(all_gauss, all_ray_params):
    dpred = forward_rf_interp(inv_res_optimiser_field_no_reg.model, gauss, ray_param)
    all_rf_dpred.append(dpred)
    dpred_init_m = forward_rf_interp(init_model, gauss, ray_param)
    all_rf_dpred_init_m.append(dpred_init_m)

axes = plot_rf_field_data(all_rf_dpred, all_gauss, all_ray_params, rf_field_times,
                          color="darkblue", linestyle="dashed",
                          label="d_pred from inverted model")
plot_rf_field_data(all_rf_dpred_init_m, all_gauss, all_ray_params, rf_field_times,
                   axes=axes, color="gray",
                   label="d_pred from starting model")
plot_rf_field_data(all_rf_field_dobs, all_gauss, all_ray_params, rf_field_times,
                   axes=axes, color="brown", linestyle="dashed",
                   label="d_obs")

axes[-1,-1].legend(loc="lower center", bbox_to_anchor=(0.5, -2.5));
<matplotlib.legend.Legend object at 0x7efd9ab8e210>

Optimal damping#

lambdas = np.logspace(-6, 6, 15)

my_lcurve_problems = []
for lamb in lambdas:
    my_problem = cofi.BaseProblem()
    kwargs = {
        "fwd_funcs": fwd_funcs,
        "d_obs_list": d_obs_list,
        "lamda": lamb
    }
    my_problem.set_objective(my_objective, kwargs=kwargs)
    my_problem.set_gradient(my_objective_gradient, kwargs=kwargs)
    my_problem.set_hessian(my_objective_hessian, kwargs=kwargs)
    my_problem.set_initial_model(init_model)
    my_lcurve_problems.append(my_problem)

def my_callback(inv_result, i):
    m = inv_result.model
    res_norm = 0
    for (fwd, fwd_kwargs), d_obs in zip(fwd_funcs, d_obs_list):
        d_pred = fwd(m, **fwd_kwargs)
        res_norm += np.sum((d_obs - d_pred) ** 2)
    reg_norm = np.sqrt(my_reg(m))
    print(f"Finished inversion with lambda={lambdas[i]}: {res_norm}, {reg_norm}")
    return res_norm, reg_norm

my_inversion_pool = cofi.utils.InversionPool(
    list_of_inv_problems=my_lcurve_problems,
    list_of_inv_options=inv_options_optimiser,
    callback=my_callback,
    parallel=False
)
all_res, all_cb_returns = my_inversion_pool.run()

l_curve_points = list(zip(*all_cb_returns))
Finished inversion with lambda=1e-06: 4.850218722067227, 4.937498526465971
Finished inversion with lambda=7.196856730011514e-06: 4.849774255571649, 4.912326656546407
Finished inversion with lambda=5.1794746792312125e-05: 4.853709052508744, 4.867685127850994
Finished inversion with lambda=0.0003727593720314938: 4.852300649703785, 3.3034374883142714
Finished inversion with lambda=0.0026826957952797246: 4.824152114403086, 3.772815609355825
Finished inversion with lambda=0.019306977288832496: 4.88136388170419, 1.9370955776488272
Finished inversion with lambda=0.1389495494373136: 4.930653800789861, 1.7211066298716369
Finished inversion with lambda=1.0: 5.3616057330004345, 1.4439460966126472
Finished inversion with lambda=7.196856730011514: 9.114989453128079, 1.017225839248105
Finished inversion with lambda=51.79474679231202: 26.641067987253134, 0.3682224494151776
Finished inversion with lambda=372.7593720314938: 39.773055639246586, 0.06361431709711397
Finished inversion with lambda=2682.6957952797275: 42.41683448109172, 0.009154583594688922
Finished inversion with lambda=19306.977288832455: 42.80554641294274, 0.0012794364910603618
Finished inversion with lambda=138949.5494373136: 42.86001723232693, 0.00017792903609310876
Finished inversion with lambda=1000000.0: 42.86754960946024, 2.4726095628418436e-05
# print all the lambdas
lambdas
array([1.00000000e-06, 7.19685673e-06, 5.17947468e-05, 3.72759372e-04,
       2.68269580e-03, 1.93069773e-02, 1.38949549e-01, 1.00000000e+00,
       7.19685673e+00, 5.17947468e+01, 3.72759372e+02, 2.68269580e+03,
       1.93069773e+04, 1.38949549e+05, 1.00000000e+06])

Plot L-curve

# plot the L-curve
res_norm, reg_norm = l_curve_points
plt.plot(reg_norm, res_norm, '.-')
plt.xlabel(r'Norm of regularization term $||Wm||_2$')
plt.ylabel(r'Norm of residual $||g(m)-d||_2$')
for i in range(0, len(lambdas), 2):
    plt.annotate(f'{lambdas[i]:.1e}', (reg_norm[i], res_norm[i]), fontsize=8)

Optimisation with damping#

lamda = 1

kwargs = {
    "fwd_funcs": fwd_funcs,
    "d_obs_list": d_obs_list,
    "lamda": lamda
}
joint_field_problem = cofi.BaseProblem()
joint_field_problem.set_objective(my_objective, kwargs=kwargs)
joint_field_problem.set_gradient(my_objective_gradient, kwargs=kwargs)
joint_field_problem.set_hessian(my_objective_hessian, kwargs=kwargs)
joint_field_problem.set_initial_model(init_model)

Define ``Inversion`` and run

inv_optimiser_field = cofi.Inversion(joint_field_problem, inv_options_optimiser)
inv_res_optimiser_field = inv_optimiser_field.run()

Plot results

_, (ax1, ax2) = plt.subplots(1, 2, figsize=(10, 4), gridspec_kw={'width_ratios': [1, 2.5]})

ax1.set_ylim(100)

plot_model(inv_res_optimiser_field.model, ax=ax1, color="green", alpha=1,
           label="model inverted from field data")
plot_model(ref_good_model, ax=ax1, color="red", alpha=1,
           label="reference good model")
plot_model(init_model, ax=ax1, alpha=1, lw=1.5, color="purple",
           label="initial model for damped solution")

field_d_periods_logspace = np.logspace(
    np.log(np.min(field_d_periods)),
    np.log(np.max(field_d_periods+1)),
    60,
    base=np.e,
)

plot_sw_data(rayleigh_field_d_obs, field_d_periods, ax=ax2, color="orange", s=5, scatter=True,
             label="d_obs")
plot_sw_data(forward_sw_interp(inv_res_optimiser_field.model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2, color="green",
             label="d_pred from inverted model")
plot_sw_data(forward_sw_interp(ref_good_model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2, color="red",
             label="d_pred from reference good model")
plot_sw_data(forward_sw_interp(init_model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2,
             alpha=1, lw=1.5, linestyle="--", color="purple",
             label="d_pred from initial model")

ax1.legend(loc="lower center", bbox_to_anchor=(0.5, -0.4))
ax2.legend(loc="lower center", bbox_to_anchor=(0.5, -0.46));
<matplotlib.legend.Legend object at 0x7efdbd317890>
all_rf_dpred = []
all_rf_dpred_init_m = []

for gauss, ray_param in zip(all_gauss, all_ray_params):
    dpred = forward_rf_interp(inv_res_optimiser_field.model, gauss, ray_param)
    all_rf_dpred.append(dpred)
    dpred_init_m = forward_rf_interp(init_model, gauss, ray_param)
    all_rf_dpred_init_m.append(dpred_init_m)

axes = plot_rf_field_data(all_rf_dpred, all_gauss, all_ray_params, rf_field_times,
                          color="darkblue", linestyle="dashed",
                          label="d_pred from inverted model")
plot_rf_field_data(all_rf_dpred_init_m, all_gauss, all_ray_params, rf_field_times,
                   axes=axes, color="gray",
                   label="d_pred from starting model")
plot_rf_field_data(all_rf_field_dobs, all_gauss, all_ray_params, rf_field_times,
                   axes=axes, color="brown", linestyle="dashed",
                   label="d_obs")

axes[-1,-1].legend(loc="lower center", bbox_to_anchor=(0.5, -2.5));
<matplotlib.legend.Legend object at 0x7efd9aaaa0d0>

Fixed-dimensional sampling#

Prepare ``BaseProblem`` for fixed-dimensional sampling

thick_min = 3
thick_max = 10
vs_min = 2
vs_max = 5

def my_log_prior(model):
    thicknesses, vs = split_layercake_model(model)
    thicknesses_out_of_bounds = (thicknesses < thick_min) | (thicknesses > thick_max)
    vs_out_of_bounds = (vs < vs_min) | (vs > vs_max)
    if np.any(thicknesses_out_of_bounds) or np.any(vs_out_of_bounds):
        return float("-inf")
    log_prior = - np.log(thick_max - thick_min) * len(thicknesses) \
                - np.log(vs_max - vs_min) * len(vs)
    return log_prior
# estimate the data noise
rayleigh_dpred = forward_sw_interp(ref_good_model, field_d_periods)
rayleigh_std = np.std(rayleigh_dpred - rayleigh_field_d_obs)
rf_dpred_1_0 = forward_rf_interp(ref_good_model, 1, ray_param_1_0)
rf_1_0_std = np.std(rf_dpred_1_0 - rf_field_dobs_1_0)
rf_dpred_2_5 = forward_rf_interp(ref_good_model, 2.5, ray_param_2_5)
rf_2_5_std = np.std(rf_dpred_2_5 - rf_field_dobs_2_5)
# inverse data covariance matrix
Cdinv_rayleigh = np.eye(len(rayleigh_field_d_obs)) / (rayleigh_std**2)
Cdinv_rf_1_0 = np.eye(len(rf_field_dobs_1_0)) / (rf_1_0_std**2)
Cdinv_rf_2_5 = np.eye(len(rf_field_dobs_2_5)) / (rf_2_5_std**2)

Cdinv_list = [Cdinv_rayleigh, Cdinv_rf_1_0, Cdinv_rf_2_5]

def my_log_likelihood(
    model,
    fwd_funcs=fwd_funcs,
    d_obs_list=d_obs_list,
    Cdinv_list=Cdinv_list
):
    log_like_sum = 0
    for (fwd, fwd_kwargs), d_obs, Cdinv in zip(fwd_funcs, d_obs_list, Cdinv_list):
        try:
            d_pred = fwd(model, **fwd_kwargs)
        except:
            return float("-inf")
        residual = d_obs - d_pred
        log_like_sum += -0.5 * residual @ (Cdinv @ residual).T
    return log_like_sum
n_walkers = 60

my_walkers_start = np.ones((n_walkers, n_dims)) * inv_res_optimiser_field.model
for i in range(n_walkers):
    my_walkers_start[i,:] += np.random.normal(0, 0.5, n_dims)
joint_field_problem.set_log_prior(my_log_prior)
joint_field_problem.set_log_likelihood(my_log_likelihood)

Define ``InversionOptions``

inv_options_fixed_d_sampling = cofi.InversionOptions()
inv_options_fixed_d_sampling.set_tool("emcee")
inv_options_fixed_d_sampling.set_params(
    nwalkers=n_walkers,
    nsteps=20_000,
    initial_state=my_walkers_start,
    skip_initial_state_check=True,
    progress=True
)

Define ``Inversion`` and run

inv_fixed_d_sampler_field = cofi.Inversion(joint_field_problem, inv_options_fixed_d_sampling)
inv_res_fixed_d_sampler_field = inv_fixed_d_sampler_field.run()
  0%|          | 0/20000 [00:00<?, ?it/s]
  0%|          | 90/20000 [00:00<00:22, 893.48it/s]
  1%|          | 182/20000 [00:00<00:21, 905.11it/s]
  1%|▏         | 273/20000 [00:00<00:22, 888.49it/s]
  2%|▏         | 363/20000 [00:00<00:22, 891.76it/s]
  2%|▏         | 454/20000 [00:00<00:21, 897.91it/s]
  3%|▎         | 544/20000 [00:00<00:21, 891.25it/s]
  3%|▎         | 636/20000 [00:00<00:21, 900.15it/s]
  4%|▎         | 727/20000 [00:00<00:21, 901.43it/s]
  4%|▍         | 818/20000 [00:00<00:21, 897.97it/s]
  5%|▍         | 908/20000 [00:01<00:21, 898.26it/s]
  5%|▍         | 999/20000 [00:01<00:21, 901.36it/s]
  5%|▌         | 1090/20000 [00:01<00:21, 898.95it/s]
  6%|▌         | 1180/20000 [00:01<00:20, 898.84it/s]
  6%|▋         | 1273/20000 [00:01<00:20, 905.79it/s]
  7%|▋         | 1364/20000 [00:01<00:20, 902.94it/s]
  7%|▋         | 1455/20000 [00:01<00:20, 897.61it/s]
  8%|▊         | 1548/20000 [00:01<00:20, 905.29it/s]
  8%|▊         | 1639/20000 [00:01<00:20, 900.80it/s]
  9%|▊         | 1731/20000 [00:01<00:20, 906.19it/s]
  9%|▉         | 1822/20000 [00:02<00:20, 905.98it/s]
 10%|▉         | 1913/20000 [00:02<00:19, 906.85it/s]
 10%|█         | 2004/20000 [00:02<00:19, 902.45it/s]
 10%|█         | 2095/20000 [00:02<00:19, 897.22it/s]
 11%|█         | 2187/20000 [00:02<00:19, 902.01it/s]
 11%|█▏        | 2278/20000 [00:02<00:19, 901.60it/s]
 12%|█▏        | 2369/20000 [00:02<00:19, 899.41it/s]
 12%|█▏        | 2460/20000 [00:02<00:19, 900.37it/s]
 13%|█▎        | 2552/20000 [00:02<00:19, 903.28it/s]
 13%|█▎        | 2644/20000 [00:02<00:19, 906.55it/s]
 14%|█▎        | 2735/20000 [00:03<00:19, 901.13it/s]
 14%|█▍        | 2829/20000 [00:03<00:18, 910.01it/s]
 15%|█▍        | 2923/20000 [00:03<00:18, 916.80it/s]
 15%|█▌        | 3016/20000 [00:03<00:18, 920.52it/s]
 16%|█▌        | 3109/20000 [00:03<00:18, 920.03it/s]
 16%|█▌        | 3202/20000 [00:03<00:18, 915.55it/s]
 16%|█▋        | 3294/20000 [00:03<00:18, 907.42it/s]
 17%|█▋        | 3387/20000 [00:03<00:18, 911.89it/s]
 17%|█▋        | 3479/20000 [00:03<00:18, 912.52it/s]
 18%|█▊        | 3573/20000 [00:03<00:17, 918.29it/s]
 18%|█▊        | 3666/20000 [00:04<00:17, 920.16it/s]
 19%|█▉        | 3759/20000 [00:04<00:17, 921.95it/s]
 19%|█▉        | 3852/20000 [00:04<00:17, 917.69it/s]
 20%|█▉        | 3945/20000 [00:04<00:17, 919.06it/s]
 20%|██        | 4037/20000 [00:04<00:17, 916.37it/s]
 21%|██        | 4129/20000 [00:04<00:17, 903.88it/s]
 21%|██        | 4220/20000 [00:04<00:17, 903.24it/s]
 22%|██▏       | 4311/20000 [00:04<00:17, 898.58it/s]
 22%|██▏       | 4401/20000 [00:04<00:17, 889.26it/s]
 22%|██▏       | 4492/20000 [00:04<00:17, 894.82it/s]
 23%|██▎       | 4583/20000 [00:05<00:17, 898.08it/s]
 23%|██▎       | 4674/20000 [00:05<00:17, 900.25it/s]
 24%|██▍       | 4766/20000 [00:05<00:16, 904.07it/s]
 24%|██▍       | 4859/20000 [00:05<00:16, 911.53it/s]
 25%|██▍       | 4953/20000 [00:05<00:16, 919.66it/s]
 25%|██▌       | 5046/20000 [00:05<00:16, 921.93it/s]
 26%|██▌       | 5139/20000 [00:05<00:16, 918.53it/s]
 26%|██▌       | 5231/20000 [00:05<00:16, 918.08it/s]
 27%|██▋       | 5324/20000 [00:05<00:15, 921.18it/s]
 27%|██▋       | 5417/20000 [00:05<00:15, 919.46it/s]
 28%|██▊       | 5509/20000 [00:06<00:15, 916.77it/s]
 28%|██▊       | 5603/20000 [00:06<00:15, 922.14it/s]
 28%|██▊       | 5696/20000 [00:06<00:15, 921.30it/s]
 29%|██▉       | 5789/20000 [00:06<00:15, 919.93it/s]
 29%|██▉       | 5882/20000 [00:06<00:15, 921.07it/s]
 30%|██▉       | 5975/20000 [00:06<00:15, 916.89it/s]
 30%|███       | 6069/20000 [00:06<00:15, 921.80it/s]
 31%|███       | 6162/20000 [00:06<00:15, 922.43it/s]
 31%|███▏      | 6256/20000 [00:06<00:14, 925.62it/s]
 32%|███▏      | 6349/20000 [00:06<00:14, 923.91it/s]
 32%|███▏      | 6442/20000 [00:07<00:14, 922.25it/s]
 33%|███▎      | 6535/20000 [00:07<00:14, 920.13it/s]
 33%|███▎      | 6629/20000 [00:07<00:14, 924.46it/s]
 34%|███▎      | 6724/20000 [00:07<00:14, 929.11it/s]
 34%|███▍      | 6819/20000 [00:07<00:14, 932.44it/s]
 35%|███▍      | 6913/20000 [00:07<00:14, 917.86it/s]
 35%|███▌      | 7006/20000 [00:07<00:14, 921.41it/s]
 35%|███▌      | 7099/20000 [00:07<00:14, 917.04it/s]
 36%|███▌      | 7193/20000 [00:07<00:13, 921.20it/s]
 36%|███▋      | 7286/20000 [00:08<00:13, 915.54it/s]
 37%|███▋      | 7378/20000 [00:08<00:13, 916.35it/s]
 37%|███▋      | 7471/20000 [00:08<00:13, 918.97it/s]
 38%|███▊      | 7563/20000 [00:08<00:13, 917.29it/s]
 38%|███▊      | 7655/20000 [00:08<00:13, 916.05it/s]
 39%|███▊      | 7747/20000 [00:08<00:13, 916.34it/s]
 39%|███▉      | 7839/20000 [00:08<00:13, 912.97it/s]
 40%|███▉      | 7934/20000 [00:08<00:13, 921.39it/s]
 40%|████      | 8027/20000 [00:08<00:13, 912.23it/s]
 41%|████      | 8119/20000 [00:08<00:13, 913.13it/s]
 41%|████      | 8212/20000 [00:09<00:12, 917.52it/s]
 42%|████▏     | 8305/20000 [00:09<00:12, 918.77it/s]
 42%|████▏     | 8399/20000 [00:09<00:12, 922.83it/s]
 42%|████▏     | 8492/20000 [00:09<00:12, 921.43it/s]
 43%|████▎     | 8585/20000 [00:09<00:12, 922.02it/s]
 43%|████▎     | 8678/20000 [00:09<00:12, 921.04it/s]
 44%|████▍     | 8771/20000 [00:09<00:12, 914.55it/s]
 44%|████▍     | 8863/20000 [00:09<00:12, 914.99it/s]
 45%|████▍     | 8955/20000 [00:09<00:12, 904.31it/s]
 45%|████▌     | 9046/20000 [00:09<00:12, 901.71it/s]
 46%|████▌     | 9138/20000 [00:10<00:11, 906.25it/s]
 46%|████▌     | 9229/20000 [00:10<00:11, 904.32it/s]
 47%|████▋     | 9324/20000 [00:10<00:11, 916.75it/s]
 47%|████▋     | 9416/20000 [00:10<00:11, 917.16it/s]
 48%|████▊     | 9508/20000 [00:10<00:11, 914.37it/s]
 48%|████▊     | 9601/20000 [00:10<00:11, 918.66it/s]
 48%|████▊     | 9694/20000 [00:10<00:11, 920.40it/s]
 49%|████▉     | 9788/20000 [00:10<00:11, 924.08it/s]
 49%|████▉     | 9881/20000 [00:10<00:10, 921.40it/s]
 50%|████▉     | 9974/20000 [00:10<00:10, 923.11it/s]
 50%|█████     | 10067/20000 [00:11<00:10, 922.85it/s]
 51%|█████     | 10160/20000 [00:11<00:10, 920.47it/s]
 51%|█████▏    | 10253/20000 [00:11<00:10, 915.43it/s]
 52%|█████▏    | 10346/20000 [00:11<00:10, 917.15it/s]
 52%|█████▏    | 10440/20000 [00:11<00:10, 922.13it/s]
 53%|█████▎    | 10533/20000 [00:11<00:10, 921.37it/s]
 53%|█████▎    | 10627/20000 [00:11<00:10, 924.12it/s]
 54%|█████▎    | 10720/20000 [00:11<00:10, 923.27it/s]
 54%|█████▍    | 10813/20000 [00:11<00:09, 925.22it/s]
 55%|█████▍    | 10906/20000 [00:11<00:09, 917.00it/s]
 55%|█████▍    | 10998/20000 [00:12<00:09, 915.74it/s]
 55%|█████▌    | 11092/20000 [00:12<00:09, 920.24it/s]
 56%|█████▌    | 11185/20000 [00:12<00:09, 921.94it/s]
 56%|█████▋    | 11278/20000 [00:12<00:09, 919.04it/s]
 57%|█████▋    | 11371/20000 [00:12<00:09, 920.56it/s]
 57%|█████▋    | 11464/20000 [00:12<00:09, 920.66it/s]
 58%|█████▊    | 11557/20000 [00:12<00:09, 912.10it/s]
 58%|█████▊    | 11649/20000 [00:12<00:09, 910.48it/s]
 59%|█████▊    | 11741/20000 [00:12<00:09, 911.57it/s]
 59%|█████▉    | 11833/20000 [00:12<00:08, 912.37it/s]
 60%|█████▉    | 11926/20000 [00:13<00:08, 917.53it/s]
 60%|██████    | 12018/20000 [00:13<00:08, 911.37it/s]
 61%|██████    | 12110/20000 [00:13<00:08, 908.16it/s]
 61%|██████    | 12201/20000 [00:13<00:08, 903.21it/s]
 61%|██████▏   | 12293/20000 [00:13<00:08, 906.01it/s]
 62%|██████▏   | 12384/20000 [00:13<00:08, 904.37it/s]
 62%|██████▏   | 12475/20000 [00:13<00:08, 899.51it/s]
 63%|██████▎   | 12565/20000 [00:13<00:08, 898.02it/s]
 63%|██████▎   | 12655/20000 [00:13<00:08, 895.97it/s]
 64%|██████▎   | 12746/20000 [00:13<00:08, 899.33it/s]
 64%|██████▍   | 12836/20000 [00:14<00:08, 894.80it/s]
 65%|██████▍   | 12928/20000 [00:14<00:07, 899.97it/s]
 65%|██████▌   | 13020/20000 [00:14<00:07, 905.68it/s]
 66%|██████▌   | 13113/20000 [00:14<00:07, 910.86it/s]
 66%|██████▌   | 13207/20000 [00:14<00:07, 919.23it/s]
 67%|██████▋   | 13301/20000 [00:14<00:07, 923.03it/s]
 67%|██████▋   | 13394/20000 [00:14<00:07, 921.34it/s]
 67%|██████▋   | 13487/20000 [00:14<00:07, 917.56it/s]
 68%|██████▊   | 13579/20000 [00:14<00:07, 909.61it/s]
 68%|██████▊   | 13670/20000 [00:14<00:07, 903.45it/s]
 69%|██████▉   | 13761/20000 [00:15<00:06, 899.03it/s]
 69%|██████▉   | 13851/20000 [00:15<00:06, 898.77it/s]
 70%|██████▉   | 13943/20000 [00:15<00:06, 904.33it/s]
 70%|███████   | 14038/20000 [00:15<00:06, 915.33it/s]
 71%|███████   | 14132/20000 [00:15<00:06, 921.98it/s]
 71%|███████   | 14226/20000 [00:15<00:06, 925.77it/s]
 72%|███████▏  | 14320/20000 [00:15<00:06, 927.86it/s]
 72%|███████▏  | 14415/20000 [00:15<00:05, 933.33it/s]
 73%|███████▎  | 14509/20000 [00:15<00:05, 932.08it/s]
 73%|███████▎  | 14603/20000 [00:15<00:05, 928.58it/s]
 73%|███████▎  | 14696/20000 [00:16<00:05, 927.69it/s]
 74%|███████▍  | 14789/20000 [00:16<00:05, 927.70it/s]
 74%|███████▍  | 14883/20000 [00:16<00:05, 931.00it/s]
 75%|███████▍  | 14977/20000 [00:16<00:05, 930.23it/s]
 75%|███████▌  | 15071/20000 [00:16<00:05, 931.50it/s]
 76%|███████▌  | 15165/20000 [00:16<00:05, 932.12it/s]
 76%|███████▋  | 15259/20000 [00:16<00:05, 928.50it/s]
 77%|███████▋  | 15352/20000 [00:16<00:05, 924.13it/s]
 77%|███████▋  | 15445/20000 [00:16<00:04, 921.85it/s]
 78%|███████▊  | 15539/20000 [00:17<00:04, 927.01it/s]
 78%|███████▊  | 15632/20000 [00:17<00:04, 926.62it/s]
 79%|███████▊  | 15725/20000 [00:17<00:04, 919.29it/s]
 79%|███████▉  | 15819/20000 [00:17<00:04, 923.05it/s]
 80%|███████▉  | 15913/20000 [00:17<00:04, 927.03it/s]
 80%|████████  | 16006/20000 [00:17<00:04, 923.86it/s]
 80%|████████  | 16099/20000 [00:17<00:04, 924.95it/s]
 81%|████████  | 16194/20000 [00:17<00:04, 930.43it/s]
 81%|████████▏ | 16288/20000 [00:17<00:03, 930.99it/s]
 82%|████████▏ | 16382/20000 [00:17<00:03, 928.50it/s]
 82%|████████▏ | 16476/20000 [00:18<00:03, 931.32it/s]
 83%|████████▎ | 16570/20000 [00:18<00:03, 930.23it/s]
 83%|████████▎ | 16664/20000 [00:18<00:03, 928.25it/s]
 84%|████████▍ | 16759/20000 [00:18<00:03, 932.74it/s]
 84%|████████▍ | 16853/20000 [00:18<00:03, 930.14it/s]
 85%|████████▍ | 16947/20000 [00:18<00:03, 930.05it/s]
 85%|████████▌ | 17041/20000 [00:18<00:03, 931.65it/s]
 86%|████████▌ | 17135/20000 [00:18<00:03, 929.06it/s]
 86%|████████▌ | 17228/20000 [00:18<00:02, 926.10it/s]
 87%|████████▋ | 17321/20000 [00:18<00:02, 925.82it/s]
 87%|████████▋ | 17414/20000 [00:19<00:02, 926.99it/s]
 88%|████████▊ | 17508/20000 [00:19<00:02, 928.07it/s]
 88%|████████▊ | 17601/20000 [00:19<00:02, 928.21it/s]
 88%|████████▊ | 17695/20000 [00:19<00:02, 931.32it/s]
 89%|████████▉ | 17789/20000 [00:19<00:02, 923.01it/s]
 89%|████████▉ | 17882/20000 [00:19<00:02, 923.37it/s]
 90%|████████▉ | 17976/20000 [00:19<00:02, 925.51it/s]
 90%|█████████ | 18070/20000 [00:19<00:02, 928.64it/s]
 91%|█████████ | 18163/20000 [00:19<00:01, 919.52it/s]
 91%|█████████▏| 18256/20000 [00:19<00:01, 921.12it/s]
 92%|█████████▏| 18349/20000 [00:20<00:01, 906.32it/s]
 92%|█████████▏| 18440/20000 [00:20<00:01, 905.67it/s]
 93%|█████████▎| 18531/20000 [00:20<00:01, 890.79it/s]
 93%|█████████▎| 18621/20000 [00:20<00:01, 892.72it/s]
 94%|█████████▎| 18711/20000 [00:20<00:01, 894.68it/s]
 94%|█████████▍| 18804/20000 [00:20<00:01, 902.62it/s]
 94%|█████████▍| 18898/20000 [00:20<00:01, 911.10it/s]
 95%|█████████▍| 18991/20000 [00:20<00:01, 914.96it/s]
 95%|█████████▌| 19084/20000 [00:20<00:00, 918.49it/s]
 96%|█████████▌| 19176/20000 [00:20<00:00, 917.73it/s]
 96%|█████████▋| 19270/20000 [00:21<00:00, 923.56it/s]
 97%|█████████▋| 19363/20000 [00:21<00:00, 923.64it/s]
 97%|█████████▋| 19456/20000 [00:21<00:00, 921.25it/s]
 98%|█████████▊| 19549/20000 [00:21<00:00, 923.06it/s]
 98%|█████████▊| 19643/20000 [00:21<00:00, 926.61it/s]
 99%|█████████▊| 19736/20000 [00:21<00:00, 927.20it/s]
 99%|█████████▉| 19829/20000 [00:21<00:00, 921.73it/s]
100%|█████████▉| 19923/20000 [00:21<00:00, 926.85it/s]
100%|██████████| 20000/20000 [00:21<00:00, 915.60it/s]
labels_v = [f"v{i}" for i in range(n_dims//2+1)]
labels_t = [f"t{i}" for i in range(n_dims//2)]
labels = [0] * n_dims
labels[::2] = labels_v
labels[1::2] = labels_t

sampler = inv_res_fixed_d_sampler_field.sampler
az_idata = az.from_emcee(sampler, var_names=labels)
az_idata.get("posterior")
<xarray.DataTree 'posterior'>
Group: /posterior
    Dimensions:  (draw: 20000, chain: 60)
    Coordinates:
      * draw     (draw) int64 160kB 0 1 2 3 4 5 ... 19995 19996 19997 19998 19999
      * chain    (chain) int64 480B 0 1 2 3 4 5 6 7 8 ... 51 52 53 54 55 56 57 58 59
    Data variables: (12/29)
        v0       (draw, chain) float64 10MB 3.722 4.249 3.103 ... 2.724 2.518 3.722
        t0       (draw, chain) float64 10MB 9.36 9.94 10.2 ... 9.557 10.26 10.22
        v1       (draw, chain) float64 10MB 3.12 3.433 3.743 ... 3.604 3.116 3.337
        t1       (draw, chain) float64 10MB 9.929 9.494 10.1 ... 9.333 10.16 10.14
        v2       (draw, chain) float64 10MB 3.981 3.585 4.127 ... 3.782 3.857 3.779
        t2       (draw, chain) float64 10MB 10.78 10.78 10.25 ... 9.283 9.734 10.23
        ...       ...
        t11      (draw, chain) float64 10MB 9.927 10.28 9.149 ... 9.369 10.45 9.77
        v12      (draw, chain) float64 10MB 3.534 4.445 3.58 ... 4.661 4.572 3.706
        t12      (draw, chain) float64 10MB 9.684 9.16 9.754 ... 9.279 10.24 9.336
        v13      (draw, chain) float64 10MB 4.49 4.171 4.702 ... 3.372 3.446 3.827
        t13      (draw, chain) float64 10MB 10.51 10.44 9.992 ... 10.31 9.546 9.99
        v14      (draw, chain) float64 10MB 5.419 4.819 5.116 ... 5.24 5.438 4.642
    Attributes:
        created_at:                 2026-05-25T05:48:00.251361+00:00
        creation_library:           ArviZ
        creation_library_version:   1.1.0
        creation_library_language:  Python
        inference_library:          emcee
        inference_library_version:  3.1.6
        sample_dims:                ['draw', 'chain']


flat_samples = sampler.get_chain(discard=10_000, thin=500, flat=True)
rand_indices = np.random.randint(len(flat_samples), size=100)

gs = gridspec.GridSpec(2, 3, width_ratios=[1, 2, 2])
fig = plt.figure(figsize=(12, 5))
ax1 = fig.add_subplot(gs[:, 0])
ax2 = fig.add_subplot(gs[:, 1])
ax3 = fig.add_subplot(gs[0, 2])
ax4 = fig.add_subplot(gs[1, 2])

ax1.set_ylim(100)
ax3.set_xlim(-5, 20)
ax3.set_ylim(-0.2, 0.6)
ax4.set_xlim(-5, 20)
ax4.set_ylim(-0.4, 1.2)

# plot samples and data predictions from samples
for idx in rand_indices:
    sample = flat_samples[idx]
    plot_model(sample, ax=ax1, alpha=0.2, lw=0.5, color="gray")
    plot_sw_data(forward_sw_interp(sample, field_d_periods_logspace),
                 field_d_periods_logspace,
                 ax=ax2, alpha=0.2, lw=0.5, color="gray")
    plot_rf_data(forward_rf_interp(sample, 1, ray_param_1_0), rf_field_times,
                 ax=ax3, alpha=0.2, lw=0.5, color="gray")
    plot_rf_data(forward_rf_interp(sample, 2.5, ray_param_2_5), rf_field_times,
                 ax=ax4, alpha=0.2, lw=0.5, color="gray")

# add labels to samples
sample_0 = flat_samples[rand_indices[0]]
plot_model(sample_0, ax=ax1, alpha=0.5, lw=0.5, color="gray", label="samples")
plot_sw_data(forward_sw_interp(sample_0, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2,
             alpha=0.5, lw=0.5, color="gray", label="rayleigh_dpred from samples")
plot_rf_data(forward_rf_interp(sample_0, 1, ray_param_1_0), rf_field_times, ax=ax3,
             alpha=0.5, lw=0.5, color="gray", label="rf_dpred from samples")
plot_rf_data(forward_rf_interp(sample_0, 2.5, ray_param_2_5), rf_field_times, ax=ax4,
             alpha=0.5, lw=0.5, color="gray", label="rf_dpred from samples")

# plot reference good model and data observations
plot_model(ref_good_model, ax=ax1, alpha=1, color="r", label="reference good model")
plot_sw_data(rayleigh_field_d_obs, field_d_periods, ax=ax2, scatter=True, color="r", s=4,
          label="rayleigh_dobs")
plot_rf_data(rf_field_dobs_1_0, rf_field_times, ax=ax3, scatter=True, color="r", s=2,
          label="rf_dobs")
plot_rf_data(rf_field_dobs_2_5, rf_field_times, ax=ax4, scatter=True, color="r", s=2,
          label="rf_dobs")

# plot damped optimisation result
plot_model(inv_res_optimiser_field.model, ax=ax1, alpha=1, color="green",
           label="damped solution")
plot_sw_data(forward_sw_interp(inv_res_optimiser_field.model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2, color="green",
             label="rayleigh_dpred from damped solution")
plot_rf_data(forward_rf_interp(inv_res_optimiser_field.model, 1, ray_param_1_0),
             rf_field_times, ax=ax3, color="green",
             label="rf_dpred from damped solution")
plot_rf_data(forward_rf_interp(inv_res_optimiser_field.model, 2.5, ray_param_2_5),
             rf_field_times, ax=ax4, color="green",
             label="rf_dpred from damped solution")

# plot initial model for dampied optimsiation
plot_model(init_model, ax=ax1, alpha=1, color="purple",
           label="initial model for damped solution")
plot_sw_data(forward_sw_interp(init_model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2, color="purple",
             label="rayleigh_dpred from initial model for damped solution")
plot_rf_data(forward_rf_interp(init_model, 1, ray_param_1_0), rf_field_times,
             ax=ax3, color="purple",
             label="rf_dpred from initial model for damped solution")
plot_rf_data(forward_rf_interp(init_model, 2.5, ray_param_2_5), rf_field_times,
             ax=ax4, color="purple",
             label="rf_dpred from initial model for damped solution")

ax1.legend(loc="upper center", bbox_to_anchor=(0.5, -0.18))
ax2.legend(loc="upper center", bbox_to_anchor=(0.5, -0.18))
ax4.legend(loc="upper center", bbox_to_anchor=(0.5, -0.6))
ax1.grid()
ax2.grid()
ax3.grid()
ax4.grid()
ax3.set_title("receiver function data (ray=1.0)")
ax4.set_title("receiver function data (ray=2.5)")

plt.tight_layout()

Trans-dimensional sampling#

Prepare utilities for trans-d sampling

def fwd_for_bayesbay(state, fwd_func, **kwargs):
    vs = state["voronoi"]["vs"]
    voronoi_sites = state["voronoi"]["discretization"]
    depths = (voronoi_sites[:-1] + voronoi_sites[1:]) / 2
    thicknesses = depths - np.insert(depths[:-1], 0, 0)
    model = form_layercake_model(thicknesses, vs)
    return fwd_func(model, **kwargs)
targets = [
    bayesbay.Target("rayleigh", rayleigh_field_d_obs, covariance_mat_inv=1/rayleigh_std**2),
    bayesbay.Target("rf_1_0", rf_field_dobs_1_0, covariance_mat_inv=1/rf_1_0_std**2),
    bayesbay.Target("rf_2_5", rf_field_dobs_2_5, covariance_mat_inv=1/rf_2_5_std**2)
]
forward_funcs = [
    (fwd_for_bayesbay, {"fwd_func": forward_sw_interp, "periods": field_d_periods}),
    (fwd_for_bayesbay, {"fwd_func": forward_rf_interp, "gauss": 1, "ray_param": ray_param_1_0}),
    (fwd_for_bayesbay, {"fwd_func": forward_rf_interp, "gauss": 2.5, "ray_param": ray_param_2_5})
]

my_log_likelihood = bayesbay.LogLikelihood(targets, forward_funcs)
/home/kitc/actions-runner/_work/cofi/cofi/docs/source/examples/scripts_field_data/surface_wave_receiver_function_joint.py:1111: DeprecationWarning: The 'Target' class has been moved to the 'likelihood' module. Please use 'from bayesbay.likelihood import Target' instead.
  bayesbay.Target("rayleigh", rayleigh_field_d_obs, covariance_mat_inv=1/rayleigh_std**2),
/home/kitc/actions-runner/_work/cofi/cofi/docs/source/examples/scripts_field_data/surface_wave_receiver_function_joint.py:1112: DeprecationWarning: The 'Target' class has been moved to the 'likelihood' module. Please use 'from bayesbay.likelihood import Target' instead.
  bayesbay.Target("rf_1_0", rf_field_dobs_1_0, covariance_mat_inv=1/rf_1_0_std**2),
/home/kitc/actions-runner/_work/cofi/cofi/docs/source/examples/scripts_field_data/surface_wave_receiver_function_joint.py:1113: DeprecationWarning: The 'Target' class has been moved to the 'likelihood' module. Please use 'from bayesbay.likelihood import Target' instead.
  bayesbay.Target("rf_2_5", rf_field_dobs_2_5, covariance_mat_inv=1/rf_2_5_std**2)
/home/kitc/actions-runner/_work/cofi/cofi/docs/source/examples/scripts_field_data/surface_wave_receiver_function_joint.py:1121: DeprecationWarning: The 'LogLikelihood' class has been moved to the 'likelihood' module. Please use 'from bayesbay.likelihood import LogLikelihood' instead.
  my_log_likelihood = bayesbay.LogLikelihood(targets, forward_funcs)
param_vs = bayesbay.prior.UniformPrior(
    name="vs",
    vmin=[2.7, 3.2, 3.75],
    vmax=[4, 4.75, 5],
    position=[0, 40, 80],
    perturb_std=0.15
)

def param_vs_initialize(param, positions):
    vmin, vmax = param.get_vmin_vmax(positions)
    sorted_vals = np.sort(np.random.uniform(vmin, vmax, positions.size))
    for i in range (len(sorted_vals)):
        val = sorted_vals[i]
        vmin_i = vmin if np.isscalar(vmin) else vmin[i]
        vmax_i = vmax if np.isscalar(vmax) else vmax[i]
        if val < vmin_i or val > vmax_i:
            if val > vmax_i: sorted_vals[i] = vmax_i
            if val < vmin_i: sorted_vals[i] = vmin_i
    return sorted_vals

param_vs.set_custom_initialize(param_vs_initialize)
parameterization = bayesbay.parameterization.Parameterization(
    bayesbay.discretization.Voronoi1D(
        name="voronoi",
        vmin=0,
        vmax=150,
        perturb_std=10,
        n_dimensions=None,
        n_dimensions_min=4,
        n_dimensions_max=15,
        parameters=[param_vs],
    )
)
my_perturbation_funcs = parameterization.perturbation_funcs
n_chains=12
walkers_start = []
for i in range(n_chains):
    walkers_start.append(parameterization.initialize())

Define ``InversionOptions``

inv_options_trans_d = cofi.InversionOptions()
inv_options_trans_d.set_tool("bayesbay")
inv_options_trans_d.set_params(
    walkers_starting_states=walkers_start,
    perturbation_funcs=my_perturbation_funcs,
    log_like_ratio_func=my_log_likelihood,
    n_chains=n_chains,
    n_iterations=2_000,
    burnin_iterations=1_000,
    verbose=False,
    save_every=200,
)

Define ``Inversion`` and run

inv_trans_d_sampler_field = cofi.Inversion(joint_field_problem, inv_options_trans_d)
inv_res_trans_d_sampler_field = inv_trans_d_sampler_field.run()
saved_states = inv_res_trans_d_sampler_field.models
samples_voronoi = saved_states["voronoi.discretization"]
samples_vs = saved_states["voronoi.vs"]
interp_depths = np.arange(150, dtype=float)
rand_indices = np.random.randint(len(samples_voronoi), size=100)

gs = gridspec.GridSpec(2, 3, width_ratios=[1, 2, 2])
fig = plt.figure(figsize=(12, 5))
ax1 = fig.add_subplot(gs[:, 0])
ax2 = fig.add_subplot(gs[:, 1])
ax3 = fig.add_subplot(gs[0, 2])
ax4 = fig.add_subplot(gs[1, 2])

ax1.set_ylim(100)
ax3.set_xlim(-5, 20)
ax3.set_ylim(-0.2, 0.6)
ax4.set_xlim(-5, 20)
ax4.set_ylim(-0.4, 1.2)

# plot samples and data predictions from samples
for idx in rand_indices:
    sample_voronoi = form_voronoi_model(samples_voronoi[idx], samples_vs[idx])
    sample = voronoi_to_layercake(sample_voronoi)
    plot_model(sample, ax=ax1, alpha=0.2, lw=0.5, color="gray")
    plot_sw_data(forward_sw_interp(sample, field_d_periods_logspace),
                 field_d_periods_logspace,
                 ax=ax2, alpha=0.2, lw=0.5, color="gray")
    plot_rf_data(forward_rf_interp(sample, 1, ray_param_1_0), rf_field_times,
                 ax=ax3, alpha=0.2, lw=0.5, color="gray")
    plot_rf_data(forward_rf_interp(sample, 2.5, ray_param_2_5), rf_field_times,
                 ax=ax4, alpha=0.2, lw=0.5, color="gray")

# add labels to samples
sample_0_voronoi = form_voronoi_model(samples_voronoi[0], samples_vs[0])
sample_0 = voronoi_to_layercake(sample_0_voronoi)
plot_model(sample_0, ax=ax1, alpha=0.5, lw=0.5, color="gray", label="samples")
plot_sw_data(forward_sw_interp(sample_0, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2,
             alpha=0.5, lw=0.5, color="gray", label="rayleigh_dpred from samples")
plot_rf_data(forward_rf_interp(sample_0, 1, ray_param_1_0), rf_field_times, ax=ax3,
             alpha=0.5, lw=0.5, color="gray", label="rf_dpred from samples")
plot_rf_data(forward_rf_interp(sample_0, 2.5, ray_param_2_5), rf_field_times, ax=ax4,
             alpha=0.5, lw=0.5, color="gray", label="rf_dpred from samples")

# plot reference good model and data observations
plot_model(ref_good_model, ax=ax1, alpha=1, color="r", label="reference good model")
plot_sw_data(rayleigh_field_d_obs, field_d_periods, ax=ax2, scatter=True, color="r", s=4,
          label="rayleigh_dobs")
plot_rf_data(rf_field_dobs_1_0, rf_field_times, ax=ax3, scatter=True, color="r", s=2,
          label="rf_dobs")
plot_rf_data(rf_field_dobs_2_5, rf_field_times, ax=ax4, scatter=True, color="r", s=2,
          label="rf_dobs")

# plot damped optimisation result
plot_model(inv_res_optimiser_field.model, ax=ax1, alpha=1, color="green",
           label="damped solution")
plot_sw_data(forward_sw_interp(inv_res_optimiser_field.model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2, color="green",
             label="rayleigh_dpred from damped solution")
plot_rf_data(forward_rf_interp(inv_res_optimiser_field.model, 1, ray_param_1_0),
             rf_field_times, ax=ax3, color="green",
             label="rf_dpred from damped solution")
plot_rf_data(forward_rf_interp(inv_res_optimiser_field.model, 2.5, ray_param_2_5),
             rf_field_times, ax=ax4, color="green",
             label="rf_dpred from damped solution")

# plot initial model for dampied optimsiation
plot_model(init_model, ax=ax1, alpha=1, color="purple",
           label="initial model for damped solution")
plot_sw_data(forward_sw_interp(init_model, field_d_periods_logspace),
             field_d_periods_logspace, ax=ax2, color="purple",
             label="rayleigh_dpred from initial model for damped solution")
plot_rf_data(forward_rf_interp(init_model, 1, ray_param_1_0), rf_field_times,
             ax=ax3, color="purple",
             label="rf_dpred from initial model for damped solution")
plot_rf_data(forward_rf_interp(init_model, 2.5, ray_param_2_5), rf_field_times,
             ax=ax4, color="purple",
             label="rf_dpred from initial model for damped solution")

ax1.legend(loc="upper center", bbox_to_anchor=(0.5, -0.18))
ax2.legend(loc="upper center", bbox_to_anchor=(0.5, -0.18))
ax4.legend(loc="upper center", bbox_to_anchor=(0.5, -0.6))
ax1.grid()
ax2.grid()
ax3.grid()
ax4.grid()
ax3.set_title("receiver function data (ray=1.0)")
ax4.set_title("receiver function data (ray=2.5)")

plt.tight_layout()

Watermark#

watermark_list = ["cofi", "numpy", "matplotlib", "scipy", "bayesbay"]
for pkg in watermark_list:
    pkg_var = __import__(pkg)
    print(pkg, getattr(pkg_var, "__version__"))
cofi 0.2.11+71.gb28b5b0
numpy 2.2.6
matplotlib 3.10.9
scipy 1.17.1
bayesbay 0.3.10

sphinx_gallery_thumbnail_number = -1

Total running time of the script: (2 minutes 15.392 seconds)

Gallery generated by Sphinx-Gallery