.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "content/tutorials/02-linear_inversion/plot_inv_2_inversion_irls.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code. .. rst-class:: sphx-glr-example-title .. _sphx_glr_content_tutorials_02-linear_inversion_plot_inv_2_inversion_irls.py: Sparse Inversion with Iteratively Re-Weighted Least-Squares =========================================================== Least-squares inversion produces smooth models which may not be an accurate representation of the true model. Here we demonstrate the basics of inverting for sparse and/or blocky models. Here, we used the iteratively reweighted least-squares approach. For this tutorial, we focus on the following: - Defining the forward problem - Defining the inverse problem (data misfit, regularization, optimization) - Defining the paramters for the IRLS algorithm - Specifying directives for the inversion - Recovering a set of model parameters which explains the observations .. GENERATED FROM PYTHON SOURCE LINES 18-37 .. code-block:: Python import numpy as np import matplotlib.pyplot as plt from discretize import TensorMesh from simpeg import ( simulation, maps, data_misfit, directives, optimization, regularization, inverse_problem, inversion, ) # sphinx_gallery_thumbnail_number = 3 .. GENERATED FROM PYTHON SOURCE LINES 38-44 Defining the Model and Mapping ------------------------------ Here we generate a synthetic model and a mappig which goes from the model space to the row space of our linear operator. .. GENERATED FROM PYTHON SOURCE LINES 44-65 .. code-block:: Python nParam = 100 # Number of model paramters # A 1D mesh is used to define the row-space of the linear operator. mesh = TensorMesh([nParam]) # Creating the true model true_model = np.zeros(mesh.nC) true_model[mesh.cell_centers_x > 0.3] = 1.0 true_model[mesh.cell_centers_x > 0.45] = -0.5 true_model[mesh.cell_centers_x > 0.6] = 0 # Mapping from the model space to the row space of the linear operator model_map = maps.IdentityMap(mesh) # Plotting the true model fig = plt.figure(figsize=(8, 5)) ax = fig.add_subplot(111) ax.plot(mesh.cell_centers_x, true_model, "b-") ax.set_ylim([-2, 2]) .. image-sg:: /content/tutorials/02-linear_inversion/images/sphx_glr_plot_inv_2_inversion_irls_001.png :alt: plot inv 2 inversion irls :srcset: /content/tutorials/02-linear_inversion/images/sphx_glr_plot_inv_2_inversion_irls_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none (-2.0, 2.0) .. GENERATED FROM PYTHON SOURCE LINES 66-73 Defining the Linear Operator ---------------------------- Here we define the linear operator with dimensions (nData, nParam). In practive, you may have a problem-specific linear operator which you would like to construct or load here. .. GENERATED FROM PYTHON SOURCE LINES 73-104 .. code-block:: Python # Number of data observations (rows) nData = 20 # Create the linear operator for the tutorial. The columns of the linear operator # represents a set of decaying and oscillating functions. jk = np.linspace(1.0, 60.0, nData) p = -0.25 q = 0.25 def g(k): return np.exp(p * jk[k] * mesh.cell_centers_x) * np.cos( np.pi * q * jk[k] * mesh.cell_centers_x ) G = np.empty((nData, nParam)) for i in range(nData): G[i, :] = g(i) # Plot the columns of G fig = plt.figure(figsize=(8, 5)) ax = fig.add_subplot(111) for i in range(G.shape[0]): ax.plot(G[i, :]) ax.set_title("Columns of matrix G") .. image-sg:: /content/tutorials/02-linear_inversion/images/sphx_glr_plot_inv_2_inversion_irls_002.png :alt: Columns of matrix G :srcset: /content/tutorials/02-linear_inversion/images/sphx_glr_plot_inv_2_inversion_irls_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none Text(0.5, 1.0, 'Columns of matrix G') .. GENERATED FROM PYTHON SOURCE LINES 105-111 Defining the Simulation ----------------------- The simulation defines the relationship between the model parameters and predicted data. .. GENERATED FROM PYTHON SOURCE LINES 111-115 .. code-block:: Python sim = simulation.LinearSimulation(mesh, G=G, model_map=model_map) .. GENERATED FROM PYTHON SOURCE LINES 116-122 Predict Synthetic Data ---------------------- Here, we use the true model to create synthetic data which we will subsequently invert. .. GENERATED FROM PYTHON SOURCE LINES 122-130 .. code-block:: Python # Standard deviation of Gaussian noise being added std = 0.02 np.random.seed(1) # Create a SimPEG data object data_obj = sim.make_synthetic_data(true_model, noise_floor=std, add_noise=True) .. GENERATED FROM PYTHON SOURCE LINES 131-140 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 .. GENERATED FROM PYTHON SOURCE LINES 140-164 .. code-block:: Python # 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(simulation=sim, data=data_obj) # Define the regularization (model objective function). Here, 'p' defines the # the norm of the smallness term and 'q' defines the norm of the smoothness # term. reg = regularization.Sparse(mesh, mapping=model_map) reg.reference_model = np.zeros(nParam) p = 0.0 q = 0.0 reg.norms = [p, q] # Define how the optimization problem is solved. opt = optimization.ProjectedGNCG( maxIter=100, lower=-2.0, upper=2.0, maxIterLS=20, maxIterCG=30, tolCG=1e-4 ) # Here we define the inverse problem that is to be solved inv_prob = inverse_problem.BaseInvProblem(dmis, reg, opt) .. GENERATED FROM PYTHON SOURCE LINES 165-172 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. .. GENERATED FROM PYTHON SOURCE LINES 172-193 .. code-block:: Python # Add sensitivity weights but don't update at each beta sensitivity_weights = directives.UpdateSensitivityWeights(every_iteration=False) # Reach target misfit for L2 solution, then use IRLS until model stops changing. IRLS = directives.Update_IRLS(max_irls_iterations=40, minGNiter=1, f_min_change=1e-4) # 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) # Update the preconditionner update_Jacobi = directives.UpdatePreconditioner() # Save output at each iteration saveDict = directives.SaveOutputEveryIteration(save_txt=False) # Define the directives as a list directives_list = [sensitivity_weights, IRLS, starting_beta, update_Jacobi, saveDict] .. rst-class:: sphx-glr-script-out .. code-block:: none /home/vsts/work/1/s/tutorials/02-linear_inversion/plot_inv_2_inversion_irls.py:177: DeprecationWarning: Update_IRLS has been deprecated, please use InversionDirective. It will be removed in version 0.24.0 of SimPEG. .. GENERATED FROM PYTHON SOURCE LINES 194-200 Setting a Starting Model and 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. .. GENERATED FROM PYTHON SOURCE LINES 200-210 .. code-block:: Python # Here we combine the inverse problem and the set of directives inv = inversion.BaseInversion(inv_prob, directives_list) # Starting model starting_model = 1e-4 * np.ones(nParam) # Run inversion recovered_model = inv.run(starting_model) .. rst-class:: sphx-glr-script-out .. code-block:: none Running inversion with SimPEG v0.22.2.dev13+g048ef809f 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 1.74e+06 3.71e+03 1.02e-09 3.71e+03 2.00e+01 0 1 8.68e+05 1.92e+03 3.67e-04 2.24e+03 1.92e+01 0 2 4.34e+05 1.32e+03 8.61e-04 1.70e+03 1.87e+01 0 Skip BFGS 3 2.17e+05 7.79e+02 1.76e-03 1.16e+03 1.79e+01 0 Skip BFGS 4 1.08e+05 3.85e+02 3.03e-03 7.14e+02 1.65e+01 0 Skip BFGS 5 5.42e+04 1.62e+02 4.46e-03 4.04e+02 1.44e+01 0 Skip BFGS 6 2.71e+04 6.06e+01 5.74e-03 2.16e+02 1.28e+01 0 Skip BFGS 7 1.36e+04 2.22e+01 6.70e-03 1.13e+02 1.02e+01 0 Skip BFGS Reached starting chifact with l2-norm regularization: Start IRLS steps... irls_threshold 1.2962714242090152 8 6.78e+03 9.62e+00 9.96e-03 7.72e+01 4.15e+00 0 Skip BFGS 9 1.53e+04 7.95e+00 1.15e-02 1.83e+02 1.77e+01 0 10 1.26e+04 2.23e+01 1.11e-02 1.62e+02 7.20e+00 0 11 1.02e+04 2.32e+01 1.18e-02 1.43e+02 3.65e+00 0 Skip BFGS 12 8.25e+03 2.29e+01 1.24e-02 1.25e+02 2.93e+00 0 Skip BFGS 13 6.68e+03 2.30e+01 1.23e-02 1.05e+02 1.65e+01 0 Skip BFGS 14 6.68e+03 2.04e+01 1.21e-02 1.01e+02 8.49e+00 0 15 6.68e+03 2.12e+01 1.15e-02 9.82e+01 9.46e+00 0 16 6.68e+03 2.19e+01 1.08e-02 9.42e+01 1.02e+01 0 17 5.51e+03 2.22e+01 1.00e-02 7.74e+01 9.18e+00 0 18 5.51e+03 1.92e+01 9.37e-03 7.08e+01 9.14e+00 0 19 5.51e+03 1.81e+01 8.52e-03 6.51e+01 9.77e+00 0 Skip BFGS 20 8.68e+03 1.74e+01 7.66e-03 8.39e+01 1.61e+01 0 Skip BFGS 21 8.68e+03 2.17e+01 6.54e-03 7.84e+01 1.23e+01 0 22 7.20e+03 2.20e+01 5.72e-03 6.32e+01 8.39e+00 0 23 7.20e+03 2.01e+01 5.03e-03 5.63e+01 1.15e+01 0 24 7.20e+03 1.96e+01 4.33e-03 5.07e+01 1.13e+01 0 25 7.20e+03 1.91e+01 3.72e-03 4.59e+01 1.15e+01 0 26 7.20e+03 1.83e+01 3.19e-03 4.13e+01 1.17e+01 0 27 1.15e+04 1.66e+01 2.68e-03 4.75e+01 1.48e+01 0 28 1.85e+04 1.66e+01 2.08e-03 5.50e+01 1.53e+01 0 29 2.92e+04 1.73e+01 1.63e-03 6.50e+01 1.54e+01 0 30 2.92e+04 1.94e+01 1.29e-03 5.71e+01 1.27e+01 0 31 2.92e+04 1.83e+01 1.05e-03 4.88e+01 1.27e+01 0 32 4.76e+04 1.58e+01 8.37e-04 5.57e+01 1.40e+01 0 33 7.60e+04 1.67e+01 6.39e-04 6.53e+01 1.49e+01 0 34 7.60e+04 1.88e+01 5.10e-04 5.75e+01 1.19e+01 0 35 1.19e+05 1.78e+01 4.24e-04 6.81e+01 1.48e+01 0 36 1.19e+05 2.00e+01 3.43e-04 6.07e+01 1.21e+01 0 37 1.19e+05 1.92e+01 2.86e-04 5.32e+01 1.21e+01 0 38 1.86e+05 1.77e+01 2.41e-04 6.25e+01 1.51e+01 0 39 1.86e+05 1.93e+01 1.97e-04 5.59e+01 1.22e+01 0 40 1.86e+05 1.85e+01 1.65e-04 4.91e+01 1.27e+01 0 41 2.94e+05 1.71e+01 1.44e-04 5.94e+01 1.53e+01 0 42 2.94e+05 1.87e+01 1.14e-04 5.24e+01 1.26e+01 0 43 4.61e+05 1.77e+01 9.84e-05 6.30e+01 1.54e+01 0 44 4.61e+05 1.94e+01 7.91e-05 5.59e+01 1.23e+01 0 45 4.61e+05 1.84e+01 6.55e-05 4.86e+01 1.25e+01 0 46 7.34e+05 1.69e+01 5.73e-05 5.89e+01 1.53e+01 0 47 7.34e+05 1.83e+01 4.56e-05 5.17e+01 1.24e+01 0 Reach maximum number of IRLS cycles: 40 ------------------------- STOP! ------------------------- 1 : |fc-fOld| = 0.0000e+00 <= tolF*(1+|f0|) = 3.7086e+02 0 : |xc-x_last| = 1.0960e-01 <= tolX*(1+|x0|) = 1.0010e-01 0 : |proj(x-g)-x| = 1.2366e+01 <= tolG = 1.0000e-01 0 : |proj(x-g)-x| = 1.2366e+01 <= 1e3*eps = 1.0000e-02 0 : maxIter = 100 <= iter = 48 ------------------------- DONE! ------------------------- .. GENERATED FROM PYTHON SOURCE LINES 211-214 Plotting Results ---------------- .. GENERATED FROM PYTHON SOURCE LINES 214-251 .. code-block:: Python fig, ax = plt.subplots(1, 2, figsize=(12 * 1.2, 4 * 1.2)) # True versus recovered model ax[0].plot(mesh.cell_centers_x, true_model, "k-") ax[0].plot(mesh.cell_centers_x, inv_prob.l2model, "b-") ax[0].plot(mesh.cell_centers_x, recovered_model, "r-") ax[0].legend(("True Model", "Recovered L2 Model", "Recovered Sparse Model")) ax[0].set_ylim([-2, 2]) # Observed versus predicted data ax[1].plot(data_obj.dobs, "k-") ax[1].plot(inv_prob.dpred, "ko") ax[1].legend(("Observed Data", "Predicted Data")) # Plot convergence fig = plt.figure(figsize=(9, 5)) ax = fig.add_axes([0.2, 0.1, 0.7, 0.85]) ax.plot(saveDict.phi_d, "k", lw=2) twin = ax.twinx() twin.plot(saveDict.phi_m, "k--", lw=2) ax.plot(np.r_[IRLS.iterStart, IRLS.iterStart], np.r_[0, np.max(saveDict.phi_d)], "k:") ax.text( IRLS.iterStart, 0.0, "IRLS Start", va="bottom", ha="center", rotation="vertical", size=12, bbox={"facecolor": "white"}, ) ax.set_ylabel(r"$\phi_d$", size=16, rotation=0) ax.set_xlabel("Iterations", size=14) twin.set_ylabel(r"$\phi_m$", size=16, rotation=0) .. rst-class:: sphx-glr-horizontal * .. image-sg:: /content/tutorials/02-linear_inversion/images/sphx_glr_plot_inv_2_inversion_irls_003.png :alt: plot inv 2 inversion irls :srcset: /content/tutorials/02-linear_inversion/images/sphx_glr_plot_inv_2_inversion_irls_003.png :class: sphx-glr-multi-img * .. image-sg:: /content/tutorials/02-linear_inversion/images/sphx_glr_plot_inv_2_inversion_irls_004.png :alt: plot inv 2 inversion irls :srcset: /content/tutorials/02-linear_inversion/images/sphx_glr_plot_inv_2_inversion_irls_004.png :class: sphx-glr-multi-img .. rst-class:: sphx-glr-script-out .. code-block:: none Text(865.1527777777777, 0.5, '$\\phi_m$') .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 35.105 seconds) **Estimated memory usage:** 234 MB .. _sphx_glr_download_content_tutorials_02-linear_inversion_plot_inv_2_inversion_irls.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_inv_2_inversion_irls.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_inv_2_inversion_irls.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: plot_inv_2_inversion_irls.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_