import numpy as np
from simpeg.electromagnetics import natural_source as nsem
from simpeg import maps
import matplotlib.pyplot as plt
from matplotlib.colors import LogNorm
import matplotlib
import matplotlib.gridspec as gridspec
from simpeg.utils import plot_1d_layer_model
from discretize import TensorMesh
from simpeg import (
maps,
data,
data_misfit,
regularization,
optimization,
inverse_problem,
inversion,
directives,
utils,
)
Step 1: Import Yellowstone MT data in mth5
format¶
from mtpy import MTCollection
mc = MTCollection()
mc.open_collection("../../data/transfer_functions/yellowstone_mt_collection.h5")
/Users/sgkang09/Project/mtpy-v2/mtpy/modeling/simpeg/recipes/inversion_1d.py:19: FutureWarning: Importing `SimPEG` is deprecated. please import from `simpeg`.
from SimPEG.electromagnetics import natural_source as nsem
# Extract 2D line
ii = 1
kk = -7
inds_sr = ['SR' in station for station in mc.dataframe['station'].values]
station_ids = mc.dataframe['station'].values[inds_sr][ii:kk]
lat = mc.dataframe['latitude'].values[inds_sr][ii:kk]
lon = mc.dataframe['longitude'].values[inds_sr][ii:kk]
plt.plot(mc.dataframe['longitude'].values, mc.dataframe['latitude'].values, '.')
plt.plot(lon, lat, '.')
plt.gca().set_aspect(1)

Step 2: Visulize transfer function at each station¶
from ipywidgets import widgets, interact
station_names = mc.dataframe.station.values
def foo(name, component):
tf = mc.get_tf(name)
tf.plot_mt_response()
Q = interact(
foo,
name=widgets.Select(options=station_ids),
component=widgets.RadioButtons(options=['xy', 'yx', 'det'], value='xy')
)
Loading...
Step 3: Design a 1D vertical mesh¶
from geoana.em.fdem import skin_depth
tf = mc.get_tf(station_ids[0])
f_min = tf.frequency.min()
f_max = tf.frequency.max()
z_max = skin_depth(f_min, 1e-2) * 1
print (z_max/1e3)
z_min = skin_depth(f_max, 1e-2) / 5
print (z_min/1e3)
24:10:14T19:57:44 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR205. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
688.6936173947421
2.7223014526313616
n_layer = 21
z_factor = 1.2
layer_thicknesses_inv = z_min*z_factor**np.arange(n_layer-1)[::-1]
print (z_max/1e3, np.sum(layer_thicknesses_inv)/1e3)
688.6936173947421 508.22101256083084
Step 4: Run separate 1D inversions¶
def run_fixed_layer_inversion(
dobs,
standard_deviation,
rho_0,
rho_ref,
maxIter=10,
maxIterCG=30,
alpha_s=1e-10,
alpha_z=1,
beta0_ratio=1,
coolingFactor=2,
coolingRate=1,
chi_factor=1,
use_irls=False,
p_s=2,
p_z=2
):
mesh_inv = TensorMesh([(np.r_[layer_thicknesses_inv, layer_thicknesses_inv[-1]])], "N")
receivers_list = [
nsem.receivers.PointNaturalSource(component="app_res"),
nsem.receivers.PointNaturalSource(component="phase"),
]
source_list = []
for freq in frequencies:
source_list.append(nsem.sources.Planewave(receivers_list, freq))
survey = nsem.survey.Survey(source_list)
sigma_map = maps.ExpMap(nP=len(layer_thicknesses_inv)+1)
simulation = nsem.simulation_1d.Simulation1DRecursive(
survey=survey,
sigmaMap=sigma_map,
thicknesses=layer_thicknesses_inv,
)
# Define the data
data_object = data.Data(survey, dobs=dobs, standard_deviation=standard_deviation)
# Initial model
m0 = np.ones(len(layer_thicknesses_inv)+1) * np.log(1./rho_0)
# Reference model
mref = np.ones(len(layer_thicknesses_inv)+1) * np.log(1./rho_ref)
dmis = data_misfit.L2DataMisfit(simulation=simulation, data=data_object)
# Define the regularization (model objective function)
reg = regularization.Sparse(mesh_inv, alpha_s=alpha_s, alpha_x=alpha_z, reference_model=mref, mapping=maps.IdentityMap(mesh_inv))
# Define how the optimization problem is solved. Here we will use an inexact
# Gauss-Newton approach that employs the conjugate gradient solver.
opt = optimization.InexactGaussNewton(maxIter=maxIter, maxIterCG=maxIterCG)
# Define the inverse problem
inv_prob = inverse_problem.BaseInvProblem(dmis, reg, opt)
#######################################################################
# Define 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.
#
# Defining a starting value for the trade-off parameter (beta) between the data
# misfit and the regularization.
starting_beta = directives.BetaEstimate_ByEig(beta0_ratio=beta0_ratio)
# Set the rate of reduction in trade-off parameter (beta) each time the
# the inverse problem is solved. And set the number of Gauss-Newton iterations
# for each trade-off paramter value.
beta_schedule = directives.BetaSchedule(coolingFactor=coolingFactor, coolingRate=coolingRate)
save_dictionary = directives.SaveOutputDictEveryIteration()
save_dictionary.outDict = {}
# Setting a stopping criteria for the inversion.
target_misfit = directives.TargetMisfit(chifact=chi_factor)
if use_irls:
reg.norms = np.c_[p_s, p_z]
# 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-5)
# The directives are defined as a list.
directives_list = [
IRLS,
starting_beta,
save_dictionary,
]
else:
# The directives are defined as a list.
directives_list = [
starting_beta,
beta_schedule,
target_misfit,
save_dictionary
]
#####################################################################
# Running the Inversion
# ---------------------
#
# To define the inversion object, we need to define the inversion problem and
# the set of directives. We can then run the inversion.
#
# Here we combine the inverse problem and the set of directives
inv = inversion.BaseInversion(inv_prob, directives_list)
# Run the inversion
recovered_model = inv.run(m0)
return recovered_model, inv_prob.phi_d
resistivity_models = []
phid_s = []
# Relative error for apparent resistivity
relative_error_rho = 0.05
# Floor for phase
floor_phase = 2.
rho_0 = 100.
rho_ref = 100.
for i_station in range(len(station_ids)):
tf = mc.get_tf(station_ids[i_station])
dobs = np.c_[tf.Z.res_det, tf.Z.phase_det].flatten()
dobs_error = np.c_[tf.Z.res_error_det, tf.Z.phase_error_det].flatten()
frequencies = tf.frequency
rho_app = dobs.reshape((len(frequencies), 2))[:,0]
phase = dobs.reshape((len(frequencies), 2))[:,1]
standard_deviation = np.c_[abs(rho_app)*relative_error_rho, np.ones(len(phase))*floor_phase].flatten()
recovered_model, phi_d = run_fixed_layer_inversion(
dobs,
standard_deviation,
rho_0,
rho_ref,
maxIter=5,
maxIterCG=30,
alpha_s=1e-5,
alpha_z=1,
beta0_ratio=1,
coolingFactor=2,
coolingRate=1,
chi_factor=1
)
phid_s.append(phi_d)
resistivity_models.append(1./np.exp(recovered_model))
24:10:14T19:57:44 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR205. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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.03e+05 1.25e+06 0.00e+00 1.25e+06 6.58e+05 0
Intel MKL WARNING: Support of Intel(R) Streaming SIMD Extensions 4.2 (Intel(R) SSE4.2) enabled only processors has been deprecated. Intel oneAPI Math Kernel Library 2025.0 will require Intel(R) Advanced Vector Extensions (Intel(R) AVX) instructions.
1 5.16e+04 1.84e+05 1.12e+00 2.41e+05 9.00e+04 0
2 2.58e+04 4.33e+04 2.15e+00 9.88e+04 1.94e+04 0 Skip BFGS
3 1.29e+04 1.65e+04 2.49e+00 4.86e+04 9.16e+03 0 Skip BFGS
4 6.45e+03 5.92e+03 2.71e+00 2.34e+04 3.91e+03 0 Skip BFGS
5 3.23e+03 3.27e+03 2.53e+00 1.14e+04 9.17e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 1.1977e+04 <= tolF*(1+|f0|) = 1.2464e+05
1 : |xc-x_last| = 1.7440e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 9.1680e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 9.1680e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:57:48 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR209. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 9.54e+06 1.07e+08 0.00e+00 1.07e+08 5.33e+07 0
1 4.77e+06 1.68e+07 1.14e+00 2.22e+07 7.00e+06 0
2 2.38e+06 4.29e+06 2.10e+00 9.30e+06 1.52e+06 0 Skip BFGS
3 1.19e+06 1.51e+06 2.54e+00 4.53e+06 6.61e+05 0 Skip BFGS
4 5.96e+05 3.57e+05 2.73e+00 1.99e+06 2.41e+05 0 Skip BFGS
5 2.98e+05 6.20e+04 2.47e+00 7.99e+05 8.76e+04 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 1.1872e+06 <= tolF*(1+|f0|) = 1.0655e+07
0 : |xc-x_last| = 2.8710e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 8.7585e+04 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 8.7585e+04 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:57:50 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR212. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
/Users/sgkang09/Project/mtpy-v2/mtpy/core/transfer_function/z.py:584: RuntimeWarning: invalid value encountered in arcsin
return np.rad2deg(np.arcsin(self.det_error / abs(self.det)))
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 5.45e+06 4.72e+07 0.00e+00 4.72e+07 2.66e+07 0
1 2.73e+06 6.72e+06 7.18e-01 8.68e+06 3.49e+06 0
2 1.36e+06 1.36e+06 1.35e+00 3.21e+06 6.26e+05 0 Skip BFGS
3 6.82e+05 4.27e+05 1.56e+00 1.49e+06 2.44e+05 0 Skip BFGS
4 3.41e+05 9.33e+04 1.60e+00 6.38e+05 8.62e+04 0
5 1.70e+05 1.81e+04 1.36e+00 2.50e+05 3.84e+04 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 3.8737e+05 <= tolF*(1+|f0|) = 4.7172e+06
0 : |xc-x_last| = 2.5757e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 3.8425e+04 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 3.8425e+04 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:57:52 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR216. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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.42e+04 1.12e+05 0.00e+00 1.12e+05 6.84e+04 0
1 7.12e+03 1.61e+04 8.06e-01 2.18e+04 9.75e+03 0
2 3.56e+03 4.50e+03 1.46e+00 9.71e+03 2.34e+03 0 Skip BFGS
3 1.78e+03 2.36e+03 1.71e+00 5.39e+03 1.17e+03 0 Skip BFGS
4 8.90e+02 1.60e+03 1.90e+00 3.29e+03 5.17e+02 0
5 4.45e+02 1.21e+03 2.15e+00 2.17e+03 3.75e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 1.1173e+03 <= tolF*(1+|f0|) = 1.1184e+04
1 : |xc-x_last| = 8.8977e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 3.7529e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 3.7529e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:57:54 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR220. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
/Users/sgkang09/Project/mtpy-v2/mtpy/core/transfer_function/z.py:584: RuntimeWarning: invalid value encountered in arcsin
return np.rad2deg(np.arcsin(self.det_error / abs(self.det)))
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 2.17e+04 2.01e+05 0.00e+00 2.01e+05 1.18e+05 0
1 1.08e+04 2.94e+04 8.67e-01 3.88e+04 1.71e+04 0
2 5.42e+03 8.00e+03 1.62e+00 1.68e+04 3.95e+03 0 Skip BFGS
3 2.71e+03 4.10e+03 1.88e+00 9.19e+03 1.90e+03 0 Skip BFGS
4 1.35e+03 2.82e+03 2.04e+00 5.58e+03 8.68e+02 0 Skip BFGS
5 6.77e+02 2.52e+03 2.05e+00 3.91e+03 2.71e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 1.6701e+03 <= tolF*(1+|f0|) = 2.0077e+04
1 : |xc-x_last| = 8.5540e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 2.7137e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 2.7137e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:57:56 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR222. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 9.81e+05 1.15e+07 0.00e+00 1.15e+07 6.12e+06 0
1 4.91e+05 1.77e+06 9.66e-01 2.24e+06 8.51e+05 0
2 2.45e+05 3.96e+05 1.87e+00 8.54e+05 1.65e+05 0 Skip BFGS
3 1.23e+05 1.37e+05 2.23e+00 4.10e+05 6.76e+04 0 Skip BFGS
4 6.13e+04 3.74e+04 2.42e+00 1.86e+05 2.68e+04 0 Skip BFGS
5 3.07e+04 5.99e+03 2.30e+00 7.66e+04 1.02e+04 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 1.0935e+05 <= tolF*(1+|f0|) = 1.1548e+06
0 : |xc-x_last| = 2.4555e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.0155e+04 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.0155e+04 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:57:58 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR225. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 6.12e+03 3.48e+04 0.00e+00 3.48e+04 2.50e+04 0
1 3.06e+03 3.86e+03 6.29e-01 5.78e+03 2.98e+03 0
2 1.53e+03 9.60e+02 9.91e-01 2.48e+03 6.27e+02 0 Skip BFGS
3 7.65e+02 4.62e+02 1.11e+00 1.31e+03 2.65e+02 0
4 3.83e+02 2.93e+02 1.20e+00 7.52e+02 1.72e+02 0 Skip BFGS
5 1.91e+02 1.99e+02 1.31e+00 4.50e+02 1.64e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 3.0181e+02 <= tolF*(1+|f0|) = 3.4827e+03
1 : |xc-x_last| = 7.2090e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.6440e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.6440e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:00 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR227. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 8.15e+05 8.30e+06 0.00e+00 8.30e+06 4.50e+06 0
1 4.07e+05 1.25e+06 8.75e-01 1.60e+06 6.16e+05 0
2 2.04e+05 2.80e+05 1.66e+00 6.17e+05 1.21e+05 0 Skip BFGS
3 1.02e+05 9.77e+04 1.97e+00 2.99e+05 5.07e+04 0 Skip BFGS
4 5.09e+04 2.72e+04 2.15e+00 1.37e+05 2.05e+04 0 Skip BFGS
5 2.55e+04 4.53e+03 2.05e+00 5.68e+04 7.99e+03 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 7.9913e+04 <= tolF*(1+|f0|) = 8.3006e+05
0 : |xc-x_last| = 2.2884e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 7.9890e+03 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 7.9890e+03 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:02 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR229. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 6.23e+03 3.51e+04 0.00e+00 3.51e+04 2.53e+04 0
1 3.12e+03 3.92e+03 6.23e-01 5.86e+03 3.11e+03 0
2 1.56e+03 9.97e+02 1.00e+00 2.56e+03 7.24e+02 0 Skip BFGS
3 7.79e+02 4.84e+02 1.13e+00 1.36e+03 3.37e+02 0
4 3.89e+02 3.05e+02 1.20e+00 7.74e+02 2.14e+02 0
5 1.95e+02 2.15e+02 1.29e+00 4.65e+02 1.85e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 3.0860e+02 <= tolF*(1+|f0|) = 3.5067e+03
1 : |xc-x_last| = 8.1045e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.8545e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.8545e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:03 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR232. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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
/Users/sgkang09/Project/mtpy-v2/mtpy/core/transfer_function/z.py:584: RuntimeWarning: invalid value encountered in arcsin
return np.rad2deg(np.arcsin(self.det_error / abs(self.det)))
0 8.45e+04 7.02e+05 0.00e+00 7.02e+05 4.13e+05 0
1 4.22e+04 9.66e+04 7.25e-01 1.27e+05 5.64e+04 0
2 2.11e+04 2.14e+04 1.40e+00 5.09e+04 1.18e+04 0 Skip BFGS
3 1.06e+04 7.94e+03 1.63e+00 2.51e+04 5.66e+03 0 Skip BFGS
4 5.28e+03 2.16e+03 1.77e+00 1.15e+04 1.99e+03 0 Skip BFGS
5 2.64e+03 1.05e+03 1.65e+00 5.40e+03 4.86e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 6.1074e+03 <= tolF*(1+|f0|) = 7.0244e+04
1 : |xc-x_last| = 1.7727e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 4.8586e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 4.8586e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:05 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR236. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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.83e+04 1.41e+05 0.00e+00 1.41e+05 9.11e+04 0
1 9.15e+03 1.75e+04 6.71e-01 2.36e+04 1.25e+04 0
2 4.57e+03 3.99e+03 1.25e+00 9.68e+03 2.79e+03 0 Skip BFGS
3 2.29e+03 1.40e+03 1.44e+00 4.70e+03 1.20e+03 0 Skip BFGS
4 1.14e+03 6.54e+02 1.49e+00 2.36e+03 3.34e+02 0
5 5.72e+02 5.51e+02 1.45e+00 1.38e+03 1.88e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 9.7732e+02 <= tolF*(1+|f0|) = 1.4060e+04
1 : |xc-x_last| = 1.0410e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.8831e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.8831e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:07 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR238. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 8.30e+03 5.08e+04 0.00e+00 5.08e+04 3.66e+04 0
1 4.15e+03 5.67e+03 5.81e-01 8.08e+03 4.80e+03 0
2 2.08e+03 1.20e+03 1.00e+00 3.28e+03 1.05e+03 0 Skip BFGS
3 1.04e+03 4.30e+02 1.14e+00 1.62e+03 4.41e+02 0
4 5.19e+02 2.35e+02 1.17e+00 8.44e+02 1.49e+02 0
5 2.59e+02 2.05e+02 1.12e+00 4.94e+02 2.30e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 3.4916e+02 <= tolF*(1+|f0|) = 5.0816e+03
1 : |xc-x_last| = 8.3588e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 2.2964e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 2.2964e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:09 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR240. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 2.26e+05 2.25e+06 0.00e+00 2.25e+06 1.26e+06 0
1 1.13e+05 3.22e+05 7.98e-01 4.12e+05 1.72e+05 0
2 5.65e+04 7.14e+04 1.61e+00 1.62e+05 3.45e+04 0 Skip BFGS
3 2.82e+04 2.76e+04 1.84e+00 7.95e+04 1.58e+04 0 Skip BFGS
4 1.41e+04 8.54e+03 2.06e+00 3.76e+04 6.66e+03 0 Skip BFGS
5 7.06e+03 2.05e+03 2.01e+00 1.62e+04 1.93e+03 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 2.1417e+04 <= tolF*(1+|f0|) = 2.2502e+05
1 : |xc-x_last| = 2.0649e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.9347e+03 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.9347e+03 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:11 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR244. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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.17e+05 1.17e+06 0.00e+00 1.17e+06 6.60e+05 0
1 5.83e+04 1.64e+05 8.03e-01 2.11e+05 9.07e+04 0
2 2.91e+04 3.78e+04 1.59e+00 8.41e+04 1.88e+04 0 Skip BFGS
3 1.46e+04 1.42e+04 1.90e+00 4.18e+04 8.41e+03 0 Skip BFGS
4 7.29e+03 4.45e+03 2.08e+00 1.96e+04 3.89e+03 0 Skip BFGS
5 3.64e+03 1.26e+03 1.99e+00 8.50e+03 1.72e+03 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 1.1131e+04 <= tolF*(1+|f0|) = 1.1658e+05
1 : |xc-x_last| = 2.0014e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.7218e+03 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.7218e+03 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:13 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR247. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 2.84e+04 2.47e+05 0.00e+00 2.47e+05 1.53e+05 0
1 1.42e+04 3.15e+04 7.00e-01 4.15e+04 2.12e+04 0
2 7.11e+03 6.96e+03 1.35e+00 1.66e+04 4.49e+03 0 Skip BFGS
3 3.55e+03 2.53e+03 1.58e+00 8.15e+03 2.05e+03 0 Skip BFGS
4 1.78e+03 9.58e+02 1.72e+00 4.01e+03 8.23e+02 0
5 8.88e+02 7.40e+02 1.59e+00 2.15e+03 2.93e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 1.8559e+03 <= tolF*(1+|f0|) = 2.4670e+04
1 : |xc-x_last| = 1.1531e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 2.9255e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 2.9255e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:15 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR250. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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
/Users/sgkang09/Project/mtpy-v2/mtpy/core/transfer_function/z.py:584: RuntimeWarning: invalid value encountered in arcsin
return np.rad2deg(np.arcsin(self.det_error / abs(self.det)))
0 4.93e+04 4.29e+05 0.00e+00 4.29e+05 2.40e+05 0
1 2.47e+04 8.96e+04 5.21e-01 1.02e+05 2.89e+04 0
2 1.23e+04 4.60e+04 1.09e+00 5.95e+04 7.40e+03 0 Skip BFGS
3 6.17e+03 4.30e+04 1.11e+00 4.99e+04 6.51e+03 0
4 3.08e+03 4.01e+04 5.41e-01 4.17e+04 1.03e+04 1
5 1.54e+03 3.81e+04 6.21e-01 3.91e+04 1.35e+03 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 2.6450e+03 <= tolF*(1+|f0|) = 4.2905e+04
1 : |xc-x_last| = 3.3550e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.3533e+03 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.3533e+03 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:16 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR252. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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
/Users/sgkang09/Project/mtpy-v2/mtpy/core/transfer_function/z.py:584: RuntimeWarning: invalid value encountered in arcsin
return np.rad2deg(np.arcsin(self.det_error / abs(self.det)))
0 8.28e+03 4.63e+04 0.00e+00 4.63e+04 2.96e+04 0
1 4.14e+03 1.64e+04 2.51e-01 1.74e+04 4.39e+03 0
2 2.07e+03 1.31e+04 5.35e-01 1.42e+04 1.03e+03 0 Skip BFGS
3 1.04e+03 1.23e+04 6.94e-01 1.31e+04 5.87e+02 0
4 5.18e+02 1.20e+04 8.41e-01 1.25e+04 2.58e+02 0 Skip BFGS
5 2.59e+02 1.18e+04 1.11e+00 1.21e+04 2.54e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 3.7522e+02 <= tolF*(1+|f0|) = 4.6292e+03
1 : |xc-x_last| = 9.2758e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 2.5393e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 2.5393e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:18 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR255. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
/Users/sgkang09/Project/mtpy-v2/mtpy/core/transfer_function/z.py:584: RuntimeWarning: invalid value encountered in arcsin
return np.rad2deg(np.arcsin(self.det_error / abs(self.det)))
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 2.03e+04 1.59e+05 0.00e+00 1.59e+05 1.02e+05 0
1 1.02e+04 2.17e+04 6.54e-01 2.83e+04 1.38e+04 0
2 5.08e+03 6.83e+03 1.23e+00 1.31e+04 3.06e+03 0 Skip BFGS
3 2.54e+03 3.91e+03 1.45e+00 7.59e+03 1.37e+03 0 Skip BFGS
4 1.27e+03 3.00e+03 1.53e+00 4.94e+03 5.23e+02 0
5 6.35e+02 2.82e+03 1.50e+00 3.77e+03 2.44e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 1.1641e+03 <= tolF*(1+|f0|) = 1.5917e+04
1 : |xc-x_last| = 1.0132e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 2.4419e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 2.4419e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:20 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR258. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 4.18e+04 3.13e+05 0.00e+00 3.13e+05 1.88e+05 0
1 2.09e+04 4.39e+04 7.51e-01 5.96e+04 2.59e+04 0
2 1.04e+04 1.10e+04 1.40e+00 2.56e+04 6.15e+03 0 Skip BFGS
3 5.22e+03 4.32e+03 1.65e+00 1.29e+04 3.24e+03 0 Skip BFGS
4 2.61e+03 1.44e+03 1.83e+00 6.21e+03 1.32e+03 0 Skip BFGS
5 1.31e+03 8.99e+02 1.69e+00 3.11e+03 4.42e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 3.1040e+03 <= tolF*(1+|f0|) = 3.1320e+04
1 : |xc-x_last| = 1.3556e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 4.4215e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 4.4215e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:22 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR264. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 2.48e+02 6.22e+03 0.00e+00 6.22e+03 1.52e+03 0
1 1.24e+02 2.35e+03 1.44e+00 2.53e+03 1.20e+03 0
2 6.21e+01 1.42e+03 2.39e+00 1.57e+03 1.27e+03 0
3 3.10e+01 3.28e+02 3.17e+00 4.26e+02 6.78e+02 0
4 1.55e+01 1.24e+02 3.68e+00 1.81e+02 2.58e+02 0 Skip BFGS
5 7.76e+00 9.21e+01 5.00e+00 1.31e+02 4.44e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 5.0277e+01 <= tolF*(1+|f0|) = 6.2191e+02
1 : |xc-x_last| = 1.6739e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 4.4382e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 4.4382e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:23 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR265. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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.86e+04 1.58e+05 0.00e+00 1.58e+05 9.80e+04 0
1 9.29e+03 2.04e+04 8.40e-01 2.82e+04 1.33e+04 0
2 4.65e+03 4.25e+03 1.52e+00 1.13e+04 2.97e+03 0 Skip BFGS
3 2.32e+03 1.54e+03 1.68e+00 5.45e+03 1.43e+03 0 Skip BFGS
4 1.16e+03 6.84e+02 1.77e+00 2.74e+03 5.27e+02 0 Skip BFGS
5 5.81e+02 5.24e+02 1.70e+00 1.51e+03 1.72e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 1.2270e+03 <= tolF*(1+|f0|) = 1.5806e+04
1 : |xc-x_last| = 1.0325e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.7237e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.7237e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:25 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR270. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 4.05e+04 2.93e+05 0.00e+00 2.93e+05 1.76e+05 0
1 2.02e+04 4.47e+04 7.37e-01 5.96e+04 2.46e+04 0
2 1.01e+04 1.24e+04 1.37e+00 2.63e+04 5.85e+03 0 Skip BFGS
3 5.06e+03 5.43e+03 1.66e+00 1.39e+04 2.92e+03 0 Skip BFGS
4 2.53e+03 2.78e+03 1.79e+00 7.31e+03 1.31e+03 0
5 1.27e+03 2.15e+03 1.67e+00 4.26e+03 5.80e+02 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 3.0458e+03 <= tolF*(1+|f0|) = 2.9291e+04
1 : |xc-x_last| = 1.4357e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 5.8042e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 5.8042e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:27 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR271. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 6.15e+08 7.28e+09 0.00e+00 7.28e+09 3.81e+09 0
1 3.07e+08 1.09e+09 1.02e+00 1.41e+09 4.93e+08 0
2 1.54e+08 2.44e+08 1.96e+00 5.45e+08 9.15e+07 0 Skip BFGS
3 7.69e+07 7.95e+07 2.30e+00 2.56e+08 3.71e+07 0 Skip BFGS
4 3.84e+07 1.40e+07 2.36e+00 1.05e+08 1.40e+07 0
5 1.92e+07 8.68e+05 1.83e+00 3.61e+07 7.28e+06 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 6.8499e+07 <= tolF*(1+|f0|) = 7.2808e+08
0 : |xc-x_last| = 3.6357e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 7.2775e+06 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 7.2775e+06 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:29 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR274. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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
/Users/sgkang09/Project/mtpy-v2/mtpy/core/transfer_function/z.py:584: RuntimeWarning: invalid value encountered in arcsin
return np.rad2deg(np.arcsin(self.det_error / abs(self.det)))
0 4.67e+05 3.83e+06 0.00e+00 3.83e+06 2.25e+06 0
1 2.34e+05 5.33e+05 6.62e-01 6.87e+05 3.00e+05 0
2 1.17e+05 1.07e+05 1.24e+00 2.52e+05 5.42e+04 0 Skip BFGS
3 5.84e+04 3.74e+04 1.44e+00 1.21e+05 2.20e+04 0 Skip BFGS
4 2.92e+04 1.22e+04 1.53e+00 5.70e+04 8.84e+03 0 Skip BFGS
5 1.46e+04 4.42e+03 1.45e+00 2.55e+04 3.67e+03 0
------------------------- STOP! -------------------------
1 : |fc-fOld| = 3.1486e+04 <= tolF*(1+|f0|) = 3.8330e+05
1 : |xc-x_last| = 1.9307e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 3.6711e+03 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 3.6711e+03 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:31 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR278. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 4.71e+03 1.83e+04 0.00e+00 1.83e+04 1.39e+04 0
1 2.36e+03 2.87e+03 4.48e-01 3.92e+03 1.81e+03 0
2 1.18e+03 1.46e+03 7.50e-01 2.34e+03 5.96e+02 0 Skip BFGS
3 5.89e+02 1.10e+03 9.22e-01 1.64e+03 3.23e+02 0 Skip BFGS
4 2.94e+02 9.01e+02 1.14e+00 1.24e+03 1.74e+02 0
5 1.47e+02 7.77e+02 1.40e+00 9.84e+02 1.75e+02 0 Skip BFGS
------------------------- STOP! -------------------------
1 : |fc-fOld| = 2.5260e+02 <= tolF*(1+|f0|) = 1.8307e+03
1 : |xc-x_last| = 7.7160e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.7544e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.7544e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:33 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR280. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 2.93e+04 2.50e+05 0.00e+00 2.50e+05 1.47e+05 0
1 1.47e+04 3.69e+04 7.35e-01 4.77e+04 2.20e+04 0
2 7.33e+03 1.01e+04 1.42e+00 2.05e+04 5.21e+03 0 Skip BFGS
3 3.66e+03 4.95e+03 1.82e+00 1.16e+04 2.49e+03 0 Skip BFGS
4 1.83e+03 2.98e+03 2.14e+00 6.91e+03 1.40e+03 0 Skip BFGS
5 9.16e+02 2.21e+03 2.40e+00 4.42e+03 8.15e+02 0 Skip BFGS
------------------------- STOP! -------------------------
1 : |fc-fOld| = 2.4928e+03 <= tolF*(1+|f0|) = 2.5048e+04
1 : |xc-x_last| = 4.7019e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 8.1504e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 8.1504e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:34 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR282. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
/Users/sgkang09/Project/mtpy-v2/mtpy/core/transfer_function/z.py:584: RuntimeWarning: invalid value encountered in arcsin
return np.rad2deg(np.arcsin(self.det_error / abs(self.det)))
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 3.97e+03 2.17e+04 0.00e+00 2.17e+04 1.30e+04 0
1 1.99e+03 7.11e+03 6.61e-01 8.43e+03 2.12e+03 0
2 9.94e+02 5.30e+03 1.17e+00 6.46e+03 8.83e+02 0 Skip BFGS
3 4.97e+02 4.47e+03 1.78e+00 5.35e+03 5.59e+02 0 Skip BFGS
4 2.48e+02 3.93e+03 2.51e+00 4.56e+03 3.34e+02 0 Skip BFGS
5 1.24e+02 3.63e+03 3.25e+00 4.03e+03 1.65e+02 0 Skip BFGS
------------------------- STOP! -------------------------
1 : |fc-fOld| = 5.2179e+02 <= tolF*(1+|f0|) = 2.1710e+03
1 : |xc-x_last| = 8.5124e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 1.6512e+02 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 1.6512e+02 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:36 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR286. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 2.56e+05 1.79e+06 0.00e+00 1.79e+06 1.05e+06 0
1 1.28e+05 3.22e+05 5.30e-01 3.90e+05 1.73e+05 0
2 6.41e+04 8.88e+04 1.16e+00 1.63e+05 3.86e+04 0 Skip BFGS
3 3.21e+04 3.85e+04 1.58e+00 8.90e+04 1.64e+04 0 Skip BFGS
4 1.60e+04 1.86e+04 1.89e+00 4.90e+04 8.91e+03 0 Skip BFGS
5 8.01e+03 8.15e+03 2.08e+00 2.48e+04 4.32e+03 0 Skip BFGS
------------------------- STOP! -------------------------
1 : |fc-fOld| = 2.4188e+04 <= tolF*(1+|f0|) = 1.7916e+05
1 : |xc-x_last| = 1.4564e+00 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 4.3169e+03 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 4.3169e+03 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
24:10:14T19:58:38 | WARNING | line:300 |mtpy.core.mt_collection | get_tf | Found multiple transfer functions with ID SR290. Suggest setting survey, otherwise returning the TF from survey Yellowstone-Snake_River_Plain.
Running inversion with SimPEG v0.22.2.dev6+g67b3e9f1c
simpeg.InvProblem is setting bfgsH0 to the inverse of the eval2Deriv.
***Done using same Solver, and solver_opts as the Simulation1DRecursive problem***
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 4.10e+02 3.60e+03 0.00e+00 3.60e+03 1.67e+03 0
1 2.05e+02 9.42e+02 8.66e-01 1.12e+03 2.07e+03 0
2 1.03e+02 1.57e+02 1.13e+00 2.73e+02 1.39e+02 0
3 5.13e+01 9.18e+01 1.50e+00 1.69e+02 7.07e+01 0 Skip BFGS
4 2.56e+01 6.18e+01 1.91e+00 1.11e+02 5.06e+01 0 Skip BFGS
------------------------- STOP! -------------------------
1 : |fc-fOld| = 0.0000e+00 <= tolF*(1+|f0|) = 3.5974e+02
1 : |xc-x_last| = 3.8839e-01 <= tolX*(1+|x0|) = 2.2104e+00
0 : |proj(x-g)-x| = 5.0552e+01 <= tolG = 1.0000e-01
0 : |proj(x-g)-x| = 5.0552e+01 <= 1e3*eps = 1.0000e-02
1 : maxIter = 5 <= iter = 5
------------------------- DONE! -------------------------
Step 5: Plot recovered 1D resistivity models¶
RHO = np.vstack(resistivity_models)
top = np.cumsum(np.r_[0, layer_thicknesses_inv[::-1]])
bottom = np.cumsum(np.r_[layer_thicknesses_inv[::-1], layer_thicknesses_inv[0]])
depth = np.r_[top[0], bottom]
import matplotlib
matplotlib.rcParams['font.size'] = 12
dx = 0.02
fig = plt.figure(figsize=(10, 6))
for ii in range(len(lon)):
x_tmp = [lon[ii]-dx, lon[ii]+dx]
out = plt.pcolormesh(x_tmp, depth/1e3, RHO[ii,::-1].reshape([-1,1]), norm=LogNorm(vmin=1e0, vmax=10000), cmap='turbo_r', alpha=0.8)
out_misfit = plt.scatter(lon, np.ones(len(lon)) * -0.5e1, c=phid_s, s=10, vmin=60, vmax=100)
plt.ylim(200, -1e1)
cb = plt.colorbar(out)
cb.set_label("Resistvity (Ohm-m)")
plt.xlabel("Longitude (degree)")
plt.ylabel("Depth (km)")
cb_mis = plt.colorbar(out_misfit, orientation='horizontal', fraction=0.03)
cb_mis.set_label("Misfit")

fig = plt.figure(figsize=(10, 6))
inds_bad = []
for ii in range(len(lon)):
if phid_s[ii]<100:
x_tmp = [lon[ii]-dx, lon[ii]+dx]
out = plt.pcolormesh(x_tmp, depth/1e3, RHO[ii,::-1].reshape([-1,1]), norm=LogNorm(vmin=1e0, vmax=10000), cmap='turbo_r', alpha=0.8)
else:
inds_bad.append(ii)
out_misfit = plt.scatter(lon, np.ones(len(lon)) * -0.5e1, c=phid_s, s=10, vmin=60, vmax=100)
plt.ylim(200, -1e1)
cb = plt.colorbar(out)
cb.set_label("Resistvity (Ohm-m)")
plt.xlabel("Longitude (degree)")
plt.ylabel("Depth (km)")
cb_mis = plt.colorbar(out_misfit, orientation='horizontal', fraction=0.03)
cb_mis.set_label("Misfit")

# Plot station with a high misfit
# ii = 0
# tf = mc.get_tf(station_ids[inds_bad[ii]])
# tf.plot_mt_response()
# print (phid_s[inds_bad[ii]])