.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "content/tutorials/07-fdem/plot_inv_1_em1dfm.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_07-fdem_plot_inv_1_em1dfm.py: 1D Inversion of for a Single Sounding ===================================== Here we use the module *simpeg.electromangetics.frequency_domain_1d* to invert frequency domain data and recover a 1D electrical conductivity model. In this tutorial, we focus on the following: - How to define sources and receivers from a survey file - How to define the survey - Sparse 1D inversion of with iteratively re-weighted least-squares For this tutorial, we will invert 1D frequency domain data for a single sounding. The end product is layered Earth model which explains the data. The survey consisted of a vertical magnetic dipole source located 30 m above the surface. The receiver measured the vertical component of the secondary field at a 10 m offset from the source in ppm. .. GENERATED FROM PYTHON SOURCE LINES 22-25 Import modules -------------- .. GENERATED FROM PYTHON SOURCE LINES 25-51 .. code-block:: Python import os import tarfile import numpy as np import matplotlib.pyplot as plt from discretize import TensorMesh import simpeg.electromagnetics.frequency_domain as fdem from simpeg.utils import mkvc, plot_1d_layer_model from simpeg import ( maps, data, data_misfit, inverse_problem, regularization, optimization, directives, inversion, utils, ) plt.rcParams.update({"font.size": 16, "lines.linewidth": 2, "lines.markersize": 8}) # sphinx_gallery_thumbnail_number = 2 .. GENERATED FROM PYTHON SOURCE LINES 52-60 Download Test Data File ----------------------- Here we provide the file path to the data we plan on inverting. The path to the data file is stored as a tar-file on our google cloud bucket: "https://storage.googleapis.com/simpeg/doc-assets/em1dfm.tar.gz" .. GENERATED FROM PYTHON SOURCE LINES 60-79 .. code-block:: Python # storage bucket where we have the data data_source = "https://storage.googleapis.com/simpeg/doc-assets/em1dfm.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 data_filename = dir_path + "em1dfm_data.txt" .. rst-class:: sphx-glr-script-out .. code-block:: none Downloading https://storage.googleapis.com/simpeg/doc-assets/em1dfm.tar.gz saved to: /home/vsts/work/1/s/tutorials/07-fdem/em1dfm.tar.gz Download completed! .. GENERATED FROM PYTHON SOURCE LINES 80-86 Load Data and Plot ------------------ Here we load and plot the 1D sounding data. In this case, we have the secondary field response in ppm for a set of frequencies. .. GENERATED FROM PYTHON SOURCE LINES 86-104 .. code-block:: Python # Load field data # dobs = np.loadtxt(str(data_filename)) dobs = np.loadtxt(str(data_filename), skiprows=1) # Define receiver locations and observed data frequencies = dobs[:, 0] dobs = mkvc(dobs[:, 1:].T) fig, ax = plt.subplots(1, 1, figsize=(7, 7)) ax.loglog(frequencies, np.abs(dobs[0::2]), "k-o", lw=3) ax.loglog(frequencies, np.abs(dobs[1::2]), "k:o", lw=3) ax.set_xlabel("Frequency (Hz)") ax.set_ylabel("|Hs/Hp| (ppm)") ax.set_title("Magnetic Field as a Function of Frequency") ax.legend(["Real", "Imaginary"]) .. image-sg:: /content/tutorials/07-fdem/images/sphx_glr_plot_inv_1_em1dfm_001.png :alt: Magnetic Field as a Function of Frequency :srcset: /content/tutorials/07-fdem/images/sphx_glr_plot_inv_1_em1dfm_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 105-113 Defining the Survey ------------------- Here we demonstrate a general way to define the receivers, sources and survey. The survey consisted of a vertical magnetic dipole source located 30 m above the surface. The receiver measured the vertical component of the secondary field at a 10 m offset from the source in ppm. .. GENERATED FROM PYTHON SOURCE LINES 113-157 .. code-block:: Python source_location = np.array([0.0, 0.0, 30.0]) moment = 1.0 receiver_location = np.array([10.0, 0.0, 30.0]) receiver_orientation = "z" data_type = "ppm" # Receiver list receiver_list = [] receiver_list.append( fdem.receivers.PointMagneticFieldSecondary( receiver_location, orientation=receiver_orientation, data_type=data_type, component="real", ) ) receiver_list.append( fdem.receivers.PointMagneticFieldSecondary( receiver_location, orientation=receiver_orientation, data_type=data_type, component="imag", ) ) # Define source list source_list = [] for freq in frequencies: source_list.append( fdem.sources.MagDipole( receiver_list=receiver_list, frequency=freq, location=source_location, orientation="z", moment=moment, ) ) # Survey survey = fdem.survey.Survey(source_list) .. GENERATED FROM PYTHON SOURCE LINES 158-164 Assign Uncertainties and Define the Data Object ----------------------------------------------- Here is where we define the data that are inverted. The data are defined by the survey, the observation values and the uncertainties. .. GENERATED FROM PYTHON SOURCE LINES 164-172 .. code-block:: Python # 5% of the absolute value uncertainties = 0.05 * np.abs(dobs) * np.ones(np.shape(dobs)) # Define the data object data_object = data.Data(survey, dobs=dobs, noise_floor=uncertainties) .. GENERATED FROM PYTHON SOURCE LINES 173-179 Defining a 1D Layered Earth (1D Tensor Mesh) -------------------------------------------- Here, we define the layer thicknesses for our 1D simulation. To do this, we use the TensorMesh class. .. GENERATED FROM PYTHON SOURCE LINES 179-187 .. code-block:: Python # Layer thicknesses inv_thicknesses = np.logspace(0, 1.5, 25) # Define a mesh for plotting and regularization. mesh = TensorMesh([(np.r_[inv_thicknesses, inv_thicknesses[-1]])], "0") .. GENERATED FROM PYTHON SOURCE LINES 188-200 Define a Starting and/or Reference Model and the Mapping -------------------------------------------------------- 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. Here, the starting model is log(0.1) S/m. Define log-conductivity values for each layer since our model is the log-conductivity. Don't make the values 0! Otherwise the gradient for the 1st iteration is zero and the inversion will not converge. .. GENERATED FROM PYTHON SOURCE LINES 200-208 .. code-block:: Python # Define model. A resistivity (Ohm meters) or conductivity (S/m) for each layer. starting_model = np.log(0.1 * np.ones(mesh.nC)) # Define mapping from model to active cells. model_mapping = maps.ExpMap() .. GENERATED FROM PYTHON SOURCE LINES 209-212 Define the Physics using a Simulation Object -------------------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 212-218 .. code-block:: Python simulation = fdem.Simulation1DLayered( survey=survey, thicknesses=inv_thicknesses, sigmaMap=model_mapping ) .. GENERATED FROM PYTHON SOURCE LINES 219-229 Define 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 229-253 .. 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. # The weighting is defined by the reciprocal of the uncertainties. dmis = data_misfit.L2DataMisfit(simulation=simulation, data=data_object) # Define the regularization (model objective function) reg_map = maps.IdentityMap(nP=mesh.nC) reg = regularization.Sparse(mesh, mapping=reg_map, alpha_s=0.025, alpha_x=1.0) # reference model reg.reference_model = starting_model # Define sparse and blocky norms p, q reg.norms = [0, 0] # Define how the optimization problem is solved. Here we will use an inexact # Gauss-Newton approach that employs the conjugate gradient solver. opt = optimization.ProjectedGNCG(maxIter=50, maxIterLS=20, maxIterCG=30, tolCG=1e-3) # Define the inverse problem inv_prob = inverse_problem.BaseInvProblem(dmis, reg, opt) .. GENERATED FROM PYTHON SOURCE LINES 254-261 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 261-290 .. code-block:: Python # 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) # Update the preconditionner update_Jacobi = directives.UpdatePreconditioner() # Options for outputting recovered models and predicted data for each beta. save_iteration = directives.SaveOutputEveryIteration(save_txt=False) # Directive for the IRLS update_IRLS = directives.UpdateIRLS(max_irls_iterations=30, irls_cooling_factor=1.5) # Updating the preconditionner if it is model dependent. update_jacobi = directives.UpdatePreconditioner() # Add sensitivity weights sensitivity_weights = directives.UpdateSensitivityWeights() # The directives are defined as a list. directives_list = [ sensitivity_weights, starting_beta, save_iteration, update_IRLS, update_jacobi, ] .. GENERATED FROM PYTHON SOURCE LINES 291-297 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 297-305 .. code-block:: Python # Here we combine the inverse problem and the set of directives inv = inversion.BaseInversion(inv_prob, directives_list) # Run the inversion recovered_model = inv.run(starting_model) .. rst-class:: sphx-glr-script-out .. code-block:: none 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 Simulation1DLayered problem*** 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 2.45e+02 2.81e+02 0.00e+00 2.81e+02 7.91e+01 0 1 1.23e+02 4.30e+01 2.19e-01 6.99e+01 3.68e+01 0 Reached starting chifact with l2-norm regularization: Start IRLS steps... irls_threshold 2.730655856743607 2 1.23e+02 7.75e+00 2.27e-01 3.56e+01 7.93e+00 0 3 1.65e+02 7.41e+00 2.16e-01 4.29e+01 1.23e+01 0 Skip BFGS 4 2.21e+02 1.09e+01 1.89e-01 5.26e+01 1.92e+01 0 5 1.22e+02 1.62e+01 1.67e-01 3.66e+01 5.71e+00 0 6 7.20e+01 1.39e+01 1.83e-01 2.70e+01 5.15e+00 0 7 4.56e+01 1.15e+01 1.79e-01 1.97e+01 4.68e+00 0 8 2.89e+01 9.06e+00 1.60e-01 1.37e+01 4.69e+00 0 9 4.14e+01 6.47e+00 1.51e-01 1.27e+01 7.70e+00 0 10 6.15e+01 5.98e+00 1.18e-01 1.33e+01 1.39e+01 0 11 8.85e+01 6.39e+00 8.52e-02 1.39e+01 2.12e+01 0 12 1.22e+02 7.02e+00 5.77e-02 1.41e+01 2.85e+01 0 13 1.63e+02 7.50e+00 3.71e-02 1.35e+01 3.27e+01 0 14 2.14e+02 7.72e+00 2.41e-02 1.29e+01 3.02e+01 0 15 2.80e+02 7.77e+00 1.60e-02 1.23e+01 2.60e+01 0 Skip BFGS 16 3.67e+02 7.75e+00 1.08e-02 1.17e+01 2.39e+01 0 Skip BFGS 17 4.83e+02 7.70e+00 7.33e-03 1.12e+01 2.33e+01 0 Skip BFGS 18 6.37e+02 7.66e+00 4.96e-03 1.08e+01 2.30e+01 0 Skip BFGS 19 8.43e+02 7.62e+00 3.35e-03 1.04e+01 2.29e+01 0 Skip BFGS 20 1.12e+03 7.58e+00 2.26e-03 1.01e+01 2.29e+01 0 Skip BFGS 21 1.49e+03 7.54e+00 1.51e-03 9.79e+00 2.28e+01 0 Skip BFGS 22 1.98e+03 7.51e+00 1.01e-03 9.52e+00 2.28e+01 0 23 2.65e+03 7.47e+00 6.76e-04 9.26e+00 2.28e+01 0 24 3.54e+03 7.44e+00 4.49e-04 9.02e+00 2.27e+01 0 25 4.76e+03 7.40e+00 2.96e-04 8.81e+00 2.27e+01 0 26 6.40e+03 7.35e+00 1.95e-04 8.60e+00 2.26e+01 0 Skip BFGS 27 8.65e+03 7.31e+00 1.28e-04 8.42e+00 2.24e+01 0 Skip BFGS 28 1.17e+04 7.27e+00 8.45e-05 8.26e+00 2.23e+01 0 Skip BFGS 29 1.59e+04 7.23e+00 5.58e-05 8.12e+00 2.22e+01 0 Skip BFGS 30 2.16e+04 7.20e+00 3.70e-05 8.00e+00 2.21e+01 0 Skip BFGS 31 2.95e+04 7.18e+00 2.45e-05 7.90e+00 2.20e+01 0 Skip BFGS Reach maximum number of IRLS cycles: 30 ------------------------- STOP! ------------------------- 1 : |fc-fOld| = 0.0000e+00 <= tolF*(1+|f0|) = 2.8196e+01 1 : |xc-x_last| = 5.5783e-02 <= tolX*(1+|x0|) = 1.2741e+00 0 : |proj(x-g)-x| = 2.2006e+01 <= tolG = 1.0000e-01 0 : |proj(x-g)-x| = 2.2006e+01 <= 1e3*eps = 1.0000e-02 0 : maxIter = 50 <= iter = 32 ------------------------- DONE! ------------------------- .. GENERATED FROM PYTHON SOURCE LINES 306-309 Plotting Results --------------------- .. GENERATED FROM PYTHON SOURCE LINES 309-366 .. code-block:: Python # Load the true model and layer thicknesses true_model = np.array([0.1, 1.0, 0.1]) true_layers = np.r_[20.0, 40.0, 160.0] # Extract Least-Squares model l2_model = inv_prob.l2model # Plot true model and recovered model fig = plt.figure(figsize=(8, 9)) x_min = np.min( np.r_[model_mapping * recovered_model, model_mapping * l2_model, true_model] ) x_max = np.max( np.r_[model_mapping * recovered_model, model_mapping * l2_model, true_model] ) ax1 = fig.add_axes([0.2, 0.15, 0.7, 0.7]) plot_1d_layer_model(true_layers, true_model, ax=ax1, show_layers=False, color="k") plot_1d_layer_model( mesh.h[0], model_mapping * l2_model, ax=ax1, show_layers=False, color="b" ) plot_1d_layer_model( mesh.h[0], model_mapping * recovered_model, ax=ax1, show_layers=False, color="r" ) ax1.set_xlim(0.01, 10) ax1.set_title("True and Recovered Models") ax1.legend(["True Model", "L2-Model", "Sparse Model"]) plt.gca().invert_yaxis() # Plot predicted and observed data dpred_l2 = simulation.dpred(l2_model) dpred_final = simulation.dpred(recovered_model) fig = plt.figure(figsize=(11, 6)) ax1 = fig.add_axes([0.2, 0.1, 0.6, 0.8]) ax1.loglog(frequencies, np.abs(dobs[0::2]), "k-o") ax1.loglog(frequencies, np.abs(dobs[1::2]), "k:o") ax1.loglog(frequencies, np.abs(dpred_l2[0::2]), "b-o") ax1.loglog(frequencies, np.abs(dpred_l2[1::2]), "b:o") ax1.loglog(frequencies, np.abs(dpred_final[0::2]), "r-o") ax1.loglog(frequencies, np.abs(dpred_final[1::2]), "r:o") ax1.set_xlabel("Frequencies (Hz)") ax1.set_ylabel("|Hs/Hp| (ppm)") ax1.set_title("Predicted and Observed Data") ax1.legend( [ "Observed (real)", "Observed (imag)", "L2-Model (real)", "L2-Model (imag)", "Sparse (real)", "Sparse (imag)", ], loc="upper left", ) plt.show() .. rst-class:: sphx-glr-horizontal * .. image-sg:: /content/tutorials/07-fdem/images/sphx_glr_plot_inv_1_em1dfm_002.png :alt: True and Recovered Models :srcset: /content/tutorials/07-fdem/images/sphx_glr_plot_inv_1_em1dfm_002.png :class: sphx-glr-multi-img * .. image-sg:: /content/tutorials/07-fdem/images/sphx_glr_plot_inv_1_em1dfm_003.png :alt: Predicted and Observed Data :srcset: /content/tutorials/07-fdem/images/sphx_glr_plot_inv_1_em1dfm_003.png :class: sphx-glr-multi-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 27.884 seconds) **Estimated memory usage:** 288 MB .. _sphx_glr_download_content_tutorials_07-fdem_plot_inv_1_em1dfm.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_1_em1dfm.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_inv_1_em1dfm.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: plot_inv_1_em1dfm.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_