API documentation¶
Created on Tue Apr 16 20:29:19 2019
@author: jkl
-
class
emagpy.Problem.
Problem
¶ Class defining an inversion problem.
-
buildANN
(fmodel, bounds, noise=0.05, iplot=False, nsample=100, dump=None, epochs=500)¶ Build and train the artificial neural network on synthetic values derived from observed ECa values.
Parameters: fmodel : function
Function use to generated the synthetic data.
bounds : list of tuple
List of tuple with the bounds of each parameters to pass to the fmodel function.
noise : float, optional
Noise level to apply on the synthetic data generated.
iplot : bool, optional
If True, the validation graph will be plotted.
nsample : int, optional
Number of samples to be synthetically generated.
dump : function, optional
Function to dump information.
epochs : int, optional
Number of epochs.
-
calibrate
(fnameECa, fnameEC=None, fnameResMod=None, forwardModel='CS', ax=None, apply=False, meshType=None, dump=None, nbins=None, binInt=None, calib=None)¶ Calibrate ECa with given EC profile.
Parameters: fnameECa : str
Path of the .csv file with the ECa data collected on the calibration points.
fnameEC : str
Path of the .csv file with the EC profile data. One row per location corresponding to the rows of fnameECa. The header should be the corresponding depths in meters positive downards.
fnameResMod : str
File name of resmod, R2 format e.g. f001_mod.dat.
forwardModel : str, optional
Forward model to use. Either CS (default), FSlin or FSeq.
ax : matplotlib.Axes
If specified the graph will be plotted against this axis.
apply : bool, optional
If True the ECa values will be calibrated. If False, the relationship will just be plotted.
dump : function, optional
Display different parts of the calibration.
binInt : int, optional
Bin interval in metres, over which to average resistivity model.
nbins : int, optional
Number of bins to average the resistivity model over.
calib : str, optional
If specified, the corresponding GF correction will be applied prior to ERT calibration. This is needed if you apply the correction to your main dataset as well! See gfCorrection().
-
computeApparentChange
(ref=0)¶ Subtract the apparent conductivities of the reference survey to all other surveys. By default the reference survey is the first survey. The survey can then be inverted with invertGN().
Parameters: ref : int, optional
Index of the reference survey. By defaut the first survey is used as reference.
-
computeChange
(ref=0)¶ Compute the difference between inverted models compared to the ref model.
Parameters: ref : int, optional
Index of the reference model. By defaut the first model (corresponding to the first survey) is used as reference.
-
computeDOI
(conds=None, depths=None, nlayers=50)¶ Compute a depth of investigation (DOI) for each 1D EC model. Sensitivity cutoff at 0.3, i.e. 70% of signal comes from above the DOI.
Parameters: forwardModel : str,
forward model
nlayers : int,
number of layers for model. Depth index. Default is first depth.
step : int,
number of steps to use for each range
fixedParam : list of int and None,
array of whether parameters out to be fixed in grid parameter search or not
bnds : list of float, optional
If specified, will create bounds for the inversion parameters
topPer : int,
Top X percentage of models to be used for model boundaries
-
computeSens
(forwardModel='CS', coils=None, models=[], depths=[])¶ Compute normalised local sensitivity using perturbation method.
Parameters: forwardModel : str, optional
- Type of forward model:
- CS : Cumulative sensitivity (default)
- FS : Full Maxwell solution with low-induction number (LIN) approximation
- FSeq : Full Maxwell solution without LIN approximation (see Andrade 2016)
coils : list of str, optional
If None, then the default attribute of the object will be used (foward mode on inverted solution). If specified, the coil spacing, orientation and height above the ground will be set. In this case you need to assign at models and depths (full forward mode). The ECa values generated will be incorporated as a new Survey object.
models : list of numpy.array of float
List of array of shape Nsample x Nlayer with conductiivty in mS/m. If empty, self.models will be used.
depths : list of numpy.array of float
List of array of shape Nsample x (Nlayer - 1) with the depth (positive number) of the bottom of the layer in meters relative to the surface. If empty self.depths will be used.
Returns: senss : list of numpy.array of float
List of matrix of size Nsample x Ncoils x Nlayers containing the normalised local sensitivity.
-
computeStat
(timef=None)¶ Compute geometrical statistics of consective points: azimuth and bearing of walking direction, time between consective measurements and distance between consecutive measurements. Results added to the main dataframe.
Parameters: timef : str, optional
Time format of the ‘time’ column of the dataframe (if available). To be passed to pd.to_datetime(). If None, it will be inferred. e.g. ‘%Y-%m-%d %H:%M:%S’
Notes
Requires projected spatial data (using convertFromNMEA() if needed).
-
convertFromNMEA
(targetProjection='EPSG:27700')¶ Convert NMEA string to selected CRS projection.
Parameters: targetProjection : str, optional
Target CRS, in EPSG number: e.g. targetProjection=’EPSG:27700’ for the British Grid.
-
createMergedSurvey
(fnames, method='nearest', how='add', targetProjection=None, unit='ppt')¶ Create a unique survey from different files by spatially interpolating the values. This can be useful when two surveys (Hi and Lo mode, vertical and horizontal) were taken on the same site successively. The method adds the ‘x’ and ‘y’ positions of all surveys and use interpolation to compute the missing values.
-
createSurvey
(fname, freq=None, hx=None, targetProjection=None, unit='ppt')¶ Create a survey object.
Parameters: fname : str
Path to the csv file with the data.
freq : float, optional
Frequency for all the coils (can also be specified for each coil in the file).
hx : float, optional
Height of the instrument above the ground (can also be specified for each coil in the file).
targetProjection : str, optional
If specified, a conversion from NMEA string in ‘Latitude’ and ‘Longitude’ columns will be performed according to EPSG code: e.g. ‘EPSG:27700’.
unit : str, optional
Unit for the _quad and _inph columns. By default assume to be in ppt (part per thousand). ppm (part per million) can also be specified. Note that ECa columns, if present are assumed to be in mS/m.
-
createTimeLapseSurvey
(fnames, targetProjection=None, unit='ppt')¶ Import multiple surveys.
Parameters: fnames : list of str
List of file to be parsed or directory where the files are.
targetProjection : str, optional
If specified, a conversion from NMEA string in ‘Latitude’ and ‘Longitude’ columns will be performed according to EPSG code: e.g. ‘EPSG:27700’.
unit : str, optional
Unit for the _quad and _inph columns. By default assume to be in ppt (part per thousand). ppm (part per million) can also be specified. Note that ECa columns, if present are assumed to be in mS/m.
-
crossOverPointsError
(index=0, coil=None, ax=None, dump=<built-in function print>)¶ Build an error model based on the cross-over points.
Parameters: index : int, optional
Survey index to fit the model on. Default is the first.
coil : str, optional
Name of the coil.
ax : Matplotlib.Axes, optional
Matplotlib axis on which the plot is plotted against if specified.
dump : function, optional
Function to print the output.
-
filterBearing
(phiMin, phiMax)¶ Keep measurements in a certain bearing range between phiMin and phiMax.
Parameters: phiMin : float, optional
Minimum angle, in degrees.
phiMax : float, optional
Maximum angle, in degrees.
-
filterDiff
(coil=None, thresh=5)¶ Keep consecutive measurements when the difference between them is smaller than thresh.
Parameters: thresh : float, optional
Value of absolute consecutive difference above which the second data point will be discarded.
coil : str, optional
Coil on which to apply the processing.
-
filterPercentile
(coil=None, qmin=None, qmax=None)¶ Filter out measurements based on percentile.
Parameters: coil : str, optional
Coil on which apply the filtering.
qmin : float, optional
Minimum percentile value below-which measurements are discarded.
qmax : float, optional
Maximum percentila value above-which measurements are discarded.
-
filterRange
(vmin=None, vmax=None)¶ Filter out measurements that are not between vmin and vmax.
Parameters: vmin : float, optional
Minimal ECa value, default is minimum observed.
vmax : float, optional
Maximum ECa value, default is maximum observed.
-
filterRepeated
(tolerance=0.2)¶ Remove consecutive points when the distance between them is below tolerance.
Parameters: tolerance : float, optional
Minimum distance away previous point in order to be retained.
-
forward
(forwardModel='CS', coils=None, noise=0.0, models=[], depths=[])¶ Compute the forward response.
Parameters: forwardModel : str, optional
- Type of forward model:
- CS : Cumulative sensitivity (default)
- FS : Full Maxwell solution with low-induction number (LIN) approximation
- FSeq : Full Maxwell solution without LIN approximation (see Andrade 2016)
coils : list of str, optional
If None, then the default attribute of the object will be used (foward mode on inverted solution). If specified, the coil spacing, orientation and height above the ground will be set. In this case you need to assign at models and depths (full forward mode). The ECa values generated will be incorporated as a new Survey object.
noise : float, optional
Percentage of noise to add on the generated apparent conductivities.
models : list of numpy.array of float
List of array of shape Nsample x Nlayer with conductiivty in mS/m. If empty, self.models will be used.
depths : list of numpy.array of float
List of array of shape Nsample x (Nlayer - 1) with the depth (positive number) of the bottom of the layer in meters relative to the surface. If empty self.depths will be used.
Returns: df : pandas.DataFrame
With the apparent ECa in the same format as input for the Survey class.
If coils argument is specified a new Survey object will be added as well.
-
getRMSE
(forwardModel=None)¶ Returns RMSPE for all coils (columns) and all surveys (row) as percentage.
np.sqrt(np.sum(((sim-obs)/obs)**2)/len(obs))*100
Parameters: forwardModel : str, optional
- Type of forward model:
- CS : Cumulative sensitivity (default)
- FS : Full Maxwell solution with low-induction number (LIN) approximation
- FSeq : Full Maxwell solution without LIN approximation (see Andrade 2016)
If None (default), the forward model used for the inversion is used.
-
gfCorrection
(calib)¶ Apply a correction to convert the calibrated ECa taking using F-0m or F-1m on CMD Explorer and Mini-Explorer to LIN ECa.
GF instruments directly map the quadrature values measured to ECa using a linear calibration. This allows to have ECa values representative of the ground EC even when the device is operated at 1 m above the ground for instance. However, this calibration gets in the way when modelling the EM response based on physical equations for the inversion. Hence, we recommend to apply a correction and convert back the ‘calibrated ECa’ to LIN ECa. This function contains the retro-engineered coefficients of the GF calibration. The ECa values are first uncalibrated back to quadrature values and then converted back to ECa using the LIN approximation.
Parameters: calib : str
Name of the calibration used. Either ‘F-0m’ of ‘F-1m’.
-
gridData
(nx=100, ny=100, method='nearest', xmin=None, xmax=None, ymin=None, ymax=None)¶ Grid data on the same grid for all surveys.
Parameters: nx : int, optional
Number of points in x direction.
ny : int, optional
Number of points in y direction.
method : str, optional
Interpolation method (nearest, cubic or linear see scipy.interpolate.griddata). Default is nearest.
xmin : float, optional
Mininum X value.
xmax : float, optional
Maximum X value.
ymin : float, optional
Minimium Y value.
ymax : float, optional
Maximum Y value
-
gridParamSearch
(forwardModel, nlayers=2, step=25, misfitMax=0.1, regularization='l1', fixedParams=None, bnds=None)¶ Using a grid based parameter search method this returns a list of best models for a specified number of layers, the minimum and maximum parameter bounds for the the top x percentage of models is also returned. This method can be used to ‘invert’ data or provide initial model parameter and parameter bounds for McMC methods.
Parameters: forwardModel : str
Forward model name. Either ‘CS’,’FSlin’ of ‘FSeq’.
nlayers : int, optional
Number of layers for model. Depth index. Default is first depth.
step : int, optional
Number of steps to use for each range.
fixedParam : list of int, optional
Array of whether parameters out to be fixed in grid parameter search or not. TODO not explicit what this does
bnds : list of float, optional
If specified, will create bounds for the inversion parameters
maxMisfit : int, optional
Maximum allowable misfit. TODO then what? happens if it’s over?’
TODO regularization is not in the docstring
TODO the docstring parameters are not ordered the same way as the argument of the method
-
importGF
(fnameLo=None, fnameHi=None, device='CMD Mini-Explorer', hx=0, calib=None, targetProjection=None)¶ Import GF instrument data with Lo and Hi file mode. If spatial data a regridding will be performed to match the data.
Parameters: fnameLo : str
Name of the file with the Lo settings.
fnameHi : str
Name of the file with the Hi settings.
device : str, optional
Type of device. Default is Mini-Explorer.
hx : float, optional
Height of the device above the ground in meters. Note that this is different from the ‘calib’ used. Data can be collected at 1 m (hx=1) but using the ‘F-0m’ calibration.
calib : str, optional
Calibration used. Either ‘F-0m’ or ‘F-1m’. If specified, the gfCorrection() function will be called and ECa values will be converted to LIN ECa (this is recommended for inversion).
targetProjection : str, optional
If both Lo and Hi dataframe contains ‘Latitude’ with NMEA values a conversion first is done using self.convertFromNMEA() before being regrid using nearest neightbours.
-
importModel
(fnames)¶ Import a save model from previous inversion.
Parameters: fname : str or list of str
Path(s) of the .csv file.
-
interpData
(surveyIndex=-1, method='nearest')¶ All locations from all or specific surveys are merged and data are interpolated using nearest neighbours.
Parameters: surveyIndex : int, optional
If -1, all locations from all surveys are merged together. Otherwise, only the locations of the survey with the selected index is used.
method : str, optional
Either, ‘nearest’, ‘linear’ or ‘cubic’ method for interpolation.
-
invert
(forwardModel='CS', method='L-BFGS-B', regularization='l2', alpha=0.07, beta=0.0, gamma=0.0, dump=None, bnds=None, options={}, Lscaling=False, rep=100, noise=0.05, nsample=100, annplot=False, threed=False, njobs=1)¶ Invert the apparent conductivity measurements.
Parameters: forwardModel : str, optional
- Type of forward model:
- CS : Cumulative sensitivity (default)
- FSlin : Full Maxwell solution with low-induction number (LIN) approximation
- FSeq : Full Maxwell solution without LIN approximation (see Andrade et al., 2016)
method : str, optional
Name of the optimization method either L-BFGS-B, TNC, CG or Nelder-Mead to be passed to scipy.optimize.minmize() or ROPE, SCEUA, DREAM, MCMC for a MCMC-based solver based on the spotpy Python package. Alternatively ‘ANN’ can be used (requires tensorflow), it will train an artificial neural network on synthetic data and use it for inversion. Note that smoothing (alpha, beta, gamma) and regularization are not supported by ANN, MCMC, DREAM and SCUEA (but well by ROPE). Another option is the use of the Gauss-Newton algorithm which can be faster in some situation. Maximum number of GN iteration can be specified as options={maxiter:3}. Default is 1.
regularization : str, optional
Type of regularization, either l1 (blocky model) or l2 (smooth model)
alpha : float, optional
Smoothing factor for the inversion.
beta : float, optional
Smoothing factor for neightbouring profile.
gamma : float, optional
Smoothing factor between surveys (for time-lapse only).
dump : function, optional
Function to print the progression. Default is print.
bnds : list of float, optional
If specified, will create bounds for the inversion. Doesn’t work with Nelder-Mead solver.
options : dict, optional
Additional dictionary arguments will be passed to scipy.optimize.minimize().
Lscaling : bool, optional
Experimental feature If True the regularization matrix will be weighted based on centroids of layers differences.
rep : int, optional
Number of sample for the MCMC-based methods.
noise : float, optional
If ANN method is used, describe the noise applied to synthetic data in training phase. Values between 0 and 1 (100% noise).
nsample : int, optional
If ANN method is used, describe the size of the synthetic data generated in trainig phase.
annplot : bool, optional
If True, the validation plot will be plotted.
threed : bool, optional
If True, the beta parameters will serve to do a quasi3D inversion.
njobs : int, optional
If -1 all CPUs are used. If 1 is given, no parallel computing code is used at all, which is useful for debugging. For n_jobs below -1, (n_cpus + 1 + n_jobs) are used. Thus for n_jobs = -2, all CPUs but one are used.
-
invertGN
(alpha=0.07, alpha_ref=None, dump=None)¶ Fast inversion usign Gauss-Newton and cumulative sensitivity.
Parameters: alpha : float, optional
Smoothing factor.
alpha_ref : float, optional
Only used for difference inversion to contrain the bottom of the profile to not changing (see Annex in Whalley et al., 2017).
dump : function, optional
Function to output the running inversion.
-
lcurve
(isurvey=0, irow=0, alphas=None, ax=None)¶ Compute an L-curve given different values of alphas.
Parameters: isurvey : int, optional
Index of survey to be used, by default the first one.
irow : int, optional
Index of measurements to be used inside the survey. First by default.
alpha : list or array-like, optional
List or array of values of alphas to build the L-curve.
ax : matplotlib.Axes, optional
If specified, the graph will be plotted agains this axis.
-
plotCrossOverMap
(index=0, coil=None, ax=None)¶ Plot the map of the cross-over points for error model.
Parameters: index : int, optional
Survey index to fit the model on. Default is the first.
coil : str, optional
Name of the coil.
ax : Matplotlib.Axes, optional
Matplotlib axis on which the plot is plotted against if specified.
-
removeCoil
(icoil)¶ Remove coil by index.
Parameters: icoil : int
Index of the coil to be remmoved.
-
resMod2EC
(fnameECa, fnameResMod, meshType=None, binInt=None, nbins=None, calib=None)¶ Convert mesh data to dfec array to be used in calibrate.
Parameters: fnameECa : str
Path of the .csv file with the ECa data collected on the calibration points.
fnameResMod : str
Path of the .dat file with the restivity model.
binInt : int, optional
Bin interval in metres, over which to average resistivity model.
nbins : int, optional
Number of bins to average the resistivity model over.
calib str, optional
If specified, will apply a GF correction. Note that the main dataset needs to be corrected as well.
-
rollingMean
(window=3)¶ Perform a rolling mean on the data.
Parameters: window : int, optional
Size of the windows for rolling mean.
-
saveInvData
(outputdir)¶ Save inverted data as one .csv per file with columns for layer conductivity (in mS/m) and columns for depth (in meters). The filename will be the same as the survey name prefixed with ‘inv_’.
Parameters: outputdir : str
Path where the .csv files will be saved.
-
saveMap
(fname, index=0, coil=None, nx=100, ny=100, method='linear', xmin=None, xmax=None, ymin=None, ymax=None, color=True, cmap='viridis_r', vmin=None, vmax=None, nlevel=14)¶ Save a georeferenced raster TIFF file.
Parameters: fname : str
Path of where to save the .tiff file.
index : int, optional
Survey number, by default, the first survey is chosen.
coil : str, optional
Name of the coil to plot. By default, the first coil is plotted.
nx : int, optional
Number of points in x direction.
ny : int, optional
Number of points in y direction.
method : str, optional
Interpolation method (nearest, cubic or linear see scipy.interpolate.griddata) or IDW (default).
xmin : float, optional
Mininum X value.
xmax : float, optional
Maximum X value.
ymin : float, optional
Minimium Y value.
ymax : float, optional
Maximum Y value
color : bool, optional
If True a colormap will be applied.
cmap : str, optional
If color == True, name of the colormap. Default is viridis.
vmin : float, optional
Minimum value for colomap.
vmax : float, optional
Maximum value for colormap.
nlevel : int, optional
Number of level in the colormap. Default 7.
-
saveSlice
(fname, index=0, islice=0, nx=100, ny=100, method='linear', xmin=None, xmax=None, ymin=None, ymax=None, color=True, cmap='viridis', vmin=None, vmax=None, nlevel=14)¶ Save a georeferenced raster TIFF file for the specified inverted depths.
Parameters: fname : str
Path of where to save the .tiff file.
index : int, optional
Survey index. Default is first.
islice : int, optional
Depth index. Default is first depth.
nx : int, optional
Number of points in x direction.
ny : int, optional
Number of points in y direction.
method : str, optional
Interpolation method (nearest, cubic or linear see scipy.interpolate.griddata) or IDW (default).
xmin : float, optional
Mininum X value.
xmax : float, optional
Maximum X value.
ymin : float, optional
Minimium Y value.
ymax : float, optional
Maximum Y value
color : bool, optional
If True a colormap will be applied.
cmap : str, optional
If color == True, name of the colormap. Default is viridis.
vmin : float, optional
Minimum value for colomap.
vmax : float, optional
Maximum value for colormap.
-
saveVTK
(fname, index=0, maxDepth=None, elev=False)¶ Writes a vtk file.
Parameters: fname : str
Path where to save the file.
index : int, optional
Index of the survey to save.
maxDepth : float, optional
Maximum positively defined depth of the bottom infinite layer.
elev : bool, optional
If True, topography will be added.
-
setInit
(depths0, conds0=None, fixedDepths=None, fixedConds=None)¶ Set the initial depths and conductivity for the inversion. Must be set after all filtering and just before the inversion.
Parameters: depths0 : list or array
Depth as positive number of the bottom of each layer. There is N-1 depths for N layers as the last layer is infinite.
conds0 : list or array, optional
Starting conductivity in mS/m of each layer. By default a homogeneous conductivity of 20 mS/m is defined.
fixedDepths : list of type bool, optional
Boolean array of same length as depths0. True if depth is fixed. False if it’s a parameter. By default all depths are fixed.
fixedConds : list of type bool, optional
Boolean array of same length as conds0. True if conductivity if fixed. False if it’s a parameter. By default all conductivity are variable.’
-
setModels
(depths, models)¶ Set the models (depth-specific EC) and the depths. Use for forward modelling.
Parameters: depths : list of array
List of array. Each array is npos x ndepths. Each depth is a positive number representing the bottom of each layer. There is not depth 0.
models : list of array
List of array of size npos x nlayers. nlayer should be equals exactly to ndepths + 1. All specified in mS/m.
-
setProjection
(targetProjection='EPSG:27700')¶ Set surveys projection to the targetProjection.
Parameters: targetProjection : str, optional
Target CRS, in EPSG number: e.g. targetProjection=’EPSG:27700’ for the British Grid.
-
show
(index=0, coil='all', ax=None, vmin=None, vmax=None, dist=True)¶ Show the raw data of the survey.
Parameters: index : int, optional
Survey number, by default, the first survey is chosen.
coil : str, optional
Specify which coil to plot. Default is all coils available.
ax : matplotlib.Axes, optional
If supplied, the graph will be plotted against ax.
vmin : float, optional
Minimal Y value.
vmax : float, optional
Maximial Y value.
dist : bool, optional
If True the true distance between points will be computed else the sample index is used as X index.
-
show3D
(index=0, pl=None, vmin=None, vmax=None, maxDepth=None, cmap='viridis_r', elev=False, edges=False, background_color=(0.8, 0.8, 0.8), pvslices=([], [][]), pvthreshold=None, pvgrid=False, pvcontour=[])¶ Show inverted model in 3D with pyvista (pip install pyvista).
Parameters: index : int, optional
Index of the survey to plot.
pl : pyvista.Plotter, optional
If specified, the graph will be plotted against it.
vmin : float, optional
Minimum value of the colorbar.
vmax : float, optional
Maximum value of the colorbar.
maxDepth : float, optional
Maximum negative depths of the graph.
cmap : str, optional
Name of the Matplotlib colormap to use.
elev : bool, optional
If True, each inverted profile will be adjusted according to elevation.
edges : bool, optional
If True, edges will be displayed.
background_color : tuple, optional
Background color assigned to pyvista plotter object when created.
pvslices : tuple of list of float, optional
Determine the X, Y, Z slices. e.g.: ([3], [], [-3, -4]) will add a slice normal to X in 3 and two slices normal to Z in -3 and -4.
pvthreshold : list of two floats, optional
Keep values between pvthreshold[0] and pvthreshold[1].
pvgrid : bool, optional
Show grid or not.
pvcontour : list of float, optional
Values of the isosurface to be plotted.
-
showDepths
(index=0, idepth=0, contour=False, vmin=None, vmax=None, cmap='viridis_r', ax=None, pts=False, levels=[])¶ Show depth slice.
Parameters: index : int, optional
Survey index. Default is first.
idepth : int, optional
Depth index. Default is first depth.
contour : bool, optional
If True then there will be contouring.
vmin : float, optional
Minimum value for colorscale.
vmax : float, optional
Maximum value for colorscale.
cmap : str, optional
Name of colormap. Default is viridis_r.
ax : Matplotlib.Axes, optional
If specified, the graph will be plotted against it.
pts : boolean, optional
If True (default) the data points will be plotted over the contour.
levels list of float, optional
If contour == True, levels are the contour intervals.
-
showDist
(index=0, coil=None, nbins=20, vmin=None, vmax=None, ax=None)¶ Display a histogram of the recorded values.
Parameters: index : int, optional
Survey number, by default, the first survey is chosen.
coil : str or list of str, optional
By default all coils are stacked on each other. Otherwise specific coils (dataframe columns) can be specified.
nbins : int, optional
Number of bins to use for the histogram.
vmin : float, optional
Value of minimum bin limit.
vmax : float, optional
Value of maximum bin limit.
ax : matplotlib.Axes, optional
If given, the graph will be plotted against this axis.
-
showMap
(index=0, coil=None, contour=False, ax=None, vmin=None, vmax=None, pts=False, cmap='viridis_r', xlab='x', ylab='y', levels=[])¶ Display a map of the measurements.
Parameters: index : int, optional
Survey number, by default, the first survey is chosen.
coil : str, optional
Name of the coil to plot. By default, the first coil is plotted.
contour : bool, optional
If True filled contour will be plotted using tricontourf.
ax : Matplotlib.Axes, optional
If specified, the graph will be plotted against the axis.
vmin : float, optional
Minimum of the colorscale.
vmax : float, optional
Maximum of the colorscale.
pts : bool, optional
If True the measurements location will be plotted on the graph.
xlab : str, optional
X label.
ylab : str, optional
Y label.
levels list of float, optional
If contour == True, levels are the contour intervals.
-
showMisfit
(index=0, coil='all', forwardModel=None, ax=None)¶ Show Misfit after inversion.
Parameters: index : int, optional
Index of the survey to plot.
coil : str, optional
Which coil to plot. Default is all.
forwardModel : str, optional
- Type of forward model:
- CS : Cumulative sensitivity (default)
- FS : Full Maxwell solution with low-induction number (LIN) approximation
- FSeq : Full Maxwell solution without LIN approximation (see Andrade 2016)
If None (default), the forward model used for the inversion is used.
ax : matplotlib.Axes, optional
If specified the graph will be plotted on this axis.
-
showOne2one
(index=0, coil='all', forwardModel=None, ax=None, vmin=None, vmax=None)¶ Show one to one plot with inversion results.
Parameters: index : int, optional
Index of the survey to plot.
coil : str, optional
Which coil to plot. Default is all.
forwardModel : str, optional
- Type of forward model:
- CS : Cumulative sensitivity (default)
- FS : Full Maxwell solution with low-induction number (LIN) approximation
- FSeq : Full Maxwell solution without LIN approximation (see Andrade 2016)
If None (default), the forward model used for the inversion is used.
ax : matplotlib.Axes, optional
If specified the graph will be plotted on this axis.
vmin : float, optional
Minimum ECa on the graph.
vmax : float, optional
Maximum ECa on the graph.
-
showProfile
(index=0, ipos=0, ax=None, vmin=None, vmax=None, maxDepth=None, errorbar=False)¶ Show specific inverted profile.
Parameters: index : int, optional
Index of the survey to plot.
ipos : int, optional
Index of the sample in the survey.
ax : Matplotlib.Axes, optional
If specified, the graph will be plotted against this axis.
rmse : bool, optional
If True, the RMSE for each transect will be plotted on a second axis. Note that misfit can also be shown with showMisfit().
errorbar : bool, optional
If True and inversion is MCMC-based, standard deviation bar are drawn for the predicted depths and conductivities.
-
showPseudo
(index=0, coil='all', ax=None, vmin=None, vmax=None, dist=True, cmap='viridis_r')¶ Show the raw data of the survey.
Parameters: index : int, optional
Survey number, by default, the first survey is chosen.
coil : str, optional
Specify which coil to plot. Default is all coils available.
ax : matplotlib.Axes, optional
If supplied, the graph will be plotted against ax.
vmin : float, optional
Minimal Y value.
vmax : float, optional
Maximial Y value.
dist : bool, optional
If True the true distance between points will be computed else the sample index is used as X index.
cmap : str, optional
Name of the colormap.
-
showResults
(index=0, ax=None, vmin=None, vmax=None, maxDepth=None, padding=1, cmap='viridis_r', dist=True, contour=False, rmse=False, errorbar=False, overlay=False, elev=False, doi=False, levels=[])¶ Show inverted model.
Parameters: index : int, optional
Index of the survey to plot.
ax : Matplotlib.Axes, optional
If specified, the graph will be plotted against this axis.
vmin : float, optional
Minimum value of the colorbar.
vmax : float, optional
Maximum value of the colorbar.
maxDepth : float, optional
Maximum negative depths of the graph.
padding : float, optional
DONT’T KNOW
cmap : str, optional
Name of the Matplotlib colormap to use.
dist : bool, optional
If True, true distances are used for X. Otherwise measurement index is used.
contour : bool, optional
If True a contour plot will be plotted.
rmse : bool, optional
If True, the RMSE for each transect will be plotted on a second axis. Note that misfit can also be shown with showMisfit().
errorbar : bool, optional
If True and inversion is MCMC-based, standard deviation bar are drawn for the predicted depths.
overlay : bool, optional
If True, a white transparent overlay is applied depending on the conductivity standard deviation from MCMC-based inversion
elev : bool, optional
If True, each inverted profile will be adjusted according to elevation.
doi : bool, optional
If True and computeDOI() was called, the estimated DOI from above which 70% of the deeper coil configuration is coming from will be plotted on top of the graph as a red dotted line.
levels list of float, optional
If contour == True, levels are the contour intervals.
-
showSlice
(index=0, islice=0, contour=False, vmin=None, vmax=None, cmap='viridis_r', ax=None, pts=False)¶ Show depth slice of EC (if islice > 0) and depth (if islice < 0).
Parameters: index : int, optional
Survey index. Default is first.
islice : int, optional
Layer index (if islice > 0). Default is first layer. If islice < 0, the depths will be display instead (e.g. islice = -1 will display the depth of the bottom of the first layer).
contour : bool, optional
If True then there will be contouring.
vmin : float, optional
Minimum value for colorscale.
vmax : float, optional
Maximum value for colorscale.
cmap : str, optional
Name of colormap. Default is viridis_r.
ax : Matplotlib.Axes, optional
If specified, the graph will be plotted against it.
pts : boolean, optional
If True (default) the data points will be plotted over the contour.
-
tcorrEC
(tdepths, tprofile, a=0.02)¶ Temperature correction for inverted models using a 2% increase of EC per degC.
EC_t = EC * (1 - a * (t - 25))
where a == 0.02 (2% by default)
Parameters: tdepths : array-like
Depths in meters of the temperature sensors (negative downards).
tprofile : array-like
Temperature values corresponding in degree Celsius.
a : float, optional
Correction coefficient. By default a 2% (a=0.02) of EC per degC is assumed.
-
tcorrECa
(tdepths, tprofile)¶ Temperature correction based on temperature profile. An ‘apparent’ temperature is computed for each coil configuration using the CS function and the observed ECa is corrected according to a 2% increase in EC per degC.
Parameters: tdepths : list of arrays
Depths in meters of the temperature sensors (negative downards).
tprofile : list of arrays
Temperature values corresponding in degree Celsius.
-
trimSurveys
()¶ Will trim all surveys to get them ready for difference inversion where all datasets must have the same number of measurements.
-
write2vtk
()¶ Write .vtk cloud points with the inverted models.
-
Created on Tue Apr 16 20:27:12 2019
@author: jkl
-
class
emagpy.Survey.
Survey
(fname=None, freq=None, hx=None, targetProjection=None, unit='ppt')¶ Create a Survey object containing the raw EMI data.
Parameters: fname : str
Path of the .csv file to import.
freq : float, optional
The instrument frequency in Hz. Can be specified per coil in the .csv.
hx : float, optional
The height of the instrument above the ground. Can be specified per coil in the .csv.
targetProjection : str, optional
If specified, the ‘Latitude’ and ‘Longitude’ NMEA string will be converted to the targeted grid e.g. : ‘EPSG:27700’.
unit : str, optional
Unit for the _quad and _inph columns. By default assume to be in ppt (part per thousand). ppm (part per million) can also be specified. Note that ECa columns, if present are assumed to be in mS/m.
-
computeStat
(timef=None)¶ Compute geometrical statistics of consective points: azimuth and bearing of walking direction, time between consective measurements and distance between consecutive measurements. Results added to the main dataframe.
Parameters: timef : str, optional
Time format of the ‘time’ column of the dataframe (if available). To be passed to pd.to_datetime(). If None, it will be inferred. e.g. ‘%Y-%m-%d %H:%M:%S’
Notes
Requires projected spatial data (using convertFromNMEA() if needed).
-
convertFromNMEA
(targetProjection='EPSG:27700')¶ Convert coordinates string (NMEA or DMS) to selected CRS projection.
Parameters: targetProjection : str, optional
Target CRS, in EPSG number: e.g. targetProjection=’EPSG:27700’ for the British Grid.
-
crossOverPointsDrift
(coil=None, ax=None, dump=<built-in function print>, minDist=1, apply=False)¶ Build an error model based on the cross-over points.
Parameters: coil : str, optional
Name of the coil.
ax : Matplotlib.Axes, optional
Matplotlib axis on which the plot is plotted against if specified.
dump : function, optional
Output function for information.
minDist : float, optional
Point at less than minDist from each other are considered identical (cross-over). Default is 1 meter.
apply : bool, optional
If True, the drift correction will be applied.
-
crossOverPointsError
(coil=None, ax=None, dump=<built-in function print>, minDist=1)¶ Build an error model based on the cross-over points.
Parameters: coil : str, optional
Name of the coil.
ax : Matplotlib.Axes, optional
Matplotlib axis on which the plot is plotted against if specified.
dump : function, optional
Output function for information.
minDist : float, optional
Point at less than minDist from each other are considered identical (cross-over). Default is 1 meter.
-
driftCorrection
(xStation=None, yStation=None, coils='all', radius=1, fit='all', ax=None, apply=False)¶ Compute drift correction from EMI given a station point and a radius.
Parameters: xStation : float, optional
X position of the drift station. Default from first point.
yStation : float, optional
Y position of the drift station. Default from first point.
coil : str or list of str, optional
Name of coil for the analysis. Default is ‘all’.
radius : float, optional
Radius around the station point inside which data will be averaged. The default is 1.
fit : str, optional
Type of fit. Either ‘all’ if one drift correction is applied on all data (default) or ‘each’ if one fit is done between each time the user came back to the drift point.
ax : matplotlib.Axes, optional
If specified, the drift graph will be plotted against. The default is None.
apply : bool, optional
If True the drift correction will be applied. The default is False.
-
filterDiff
(coil=None, thresh=5)¶ Keep consecutive measurements when the difference between them is smaller than thresh.
Parameters: thresh : float, optional
Value of absolute consecutive difference above which the second data point will be discarded.
coil : str, optional
Coil on which to apply the processing.
-
filterPercentile
(coil=None, qmin=None, qmax=None)¶ Filter out measurements based on percentile.
Parameters: coil : str, optional
Coil on which apply the filtering.
qmin : float, optional
Minimum percentile value below-which measurements are discarded.
qmax : float, optional
Maximum percentila value above-which measurements are discarded.
-
filterRange
(vmin=None, vmax=None)¶ Filter out measurements that are not between vmin and vmax.
Parameters: vmin : float, optional
Minimal ECa value, default is minimum observed.
vmax : float, optional
Maximum ECa value, default is maximum observed.
-
filterRepeated
(tolerance=0.2)¶ Remove consecutive points when the distance between them is below tolerance.
Parameters: tolerance : float, optional
Minimum distance away previous point in order to be retained.
-
fitlerBearing
(phiMin, phiMax)¶ Keep measurements in a certain bearing range between phiMin and phiMax.
Parameters: phiMin : float, optional
Minimum angle, in degrees.
phiMax : float, optional
Maximum angle, in degrees.
-
gfCorrection
(calib)¶ Converting GF calibrated ECa to LIN ECa.
GF instruments directly map the quadrature values measured to ECa using a linear calibration. This allows to have ECa values representative of the ground EC even when the device is operated at 1 m above the ground for instance. However, this calibration gets in the way when modelling the EM response based on physical equations for the inversion. Hence, we recommend to apply a correction and convert back the ‘calibrated ECa’ to LIN ECa. This function contains the retro-engineered coefficients of the GF calibration. The ECa values are first uncalibrated back to quadrature values and then converted back to ECa using the LIN approximation.
Parameters: calib : str
Name of the calibration used. Either ‘F-0m’ of ‘F-1m’.
-
gridData
(nx=100, ny=100, method='nearest', xmin=None, xmax=None, ymin=None, ymax=None)¶ Grid data.
Parameters: nx : int, optional
Number of points in x direction.
ny : int, optional
Number of points in y direction.
xmin : float, optional
Mininum X value.
xmax : float, optional
Maximum X value.
ymin : float, optional
Minimium Y value.
ymax : float, optional
Maximum Y value
method : str, optional
Interpolation method (nearest, cubic or linear see scipy.interpolate.griddata) or IDW (default).
-
importGF
(fnameLo=None, fnameHi=None, device='CMD Mini-Explorer', hx=0, calib=None, targetProjection=None)¶ Import GF instrument data with Lo and Hi file mode. If spatial data a regridding will be performed to match the data.
Parameters: fnameLo : str
Name of the file with the Lo settings.
fnameHi : str, optional
Name of the file with the Hi settings.
device : str, optional
Type of device. Default is Mini-Explorer.
hx : float, optional
Height of the device above the ground in meters. Note that this is different from the ‘calib’ used. Data can be collected at 1 m (hx=1) but using the ‘F-0m’ calibration.
calib : str, optional
Calibration used. Either ‘F-0m’ or ‘F-1m’. If specified, the gfCorrection() function will be called and ECa values will be converted to LIN ECa (this is recommended for inversion).
targetProjection : str, optional
If both Lo and Hi dataframe contains ‘Latitude’ with NMEA values a conversion first is done using self.convertFromNMEA() before being regrid using nearest neightbours. If not supplied, regridding is done using WGS84 latitude and longitude.
-
plotCrossOverMap
(coil=None, ax=None, minDist=1)¶ Plot the map of the cross-over points for error model.
Parameters: coil : str, optional
Name of the coil.
ax : Matplotlib.Axes, optional
Matplotlib axis on which the plot is plotted against if specified.
minDist : float, optional
Point at less than minDist from each other are considered identical (cross-over). Default is 1 meter.
-
readDF
(df, name=None, sensor=None, targetProjection=None, unit='ppt')¶ Parse dataframe.
Parameters: df : pandas.DataFrame
A pandas dataframe where each column contains ECa per coil and each row a different positions.
name : str, optional
Name of the survey.
sensor : str, optional
A string describing the sensor (for metadata only).
targetProjection : str, optional
EPSG string describing the projection of a ‘Latitude’ and ‘Longitude’ column is found in the dataframe. e.g. ‘EPSG:27700’ for the British grid.
unit : str, optional
Unit for the _quad and _inph columns. By default assume to be in ppt (part per thousand). ppm (part per million) can also be specified. Note that ECa columns, if present are assumed to be in mS/m.
-
readFile
(fname, sensor=None, targetProjection=None, unit='ppt')¶ Read a .csv file.
Parameters: fname : str
Filename.
sensor : str, optional
Name of the sensor (for metadata only).
targetProjection : str, optional
EPSG string describing the projection of a ‘Latitude’ and ‘Longitude’ column is found in the dataframe. e.g. ‘EPSG:27700’ for the British grid.
unit : str, optional
Unit for the _quad and _inph columns. By default assume to be in ppt (part per thousand). ppm (part per million) can also be specified. Note that ECa columns, if present are assumed to be in mS/m.
-
removeCoil
(icoil)¶ Remove coil by index.
Parameters: icoil : int
Index of the coil to be remmoved.
-
rollingMean
(window=3)¶ Perform a rolling mean on the data.
Parameters: window : int, optional
Size of the windows for rolling mean.
-
saveMap
(fname, coil=None, nx=100, ny=100, method='linear', xmin=None, xmax=None, ymin=None, ymax=None, color=True, cmap='viridis_r', vmin=None, vmax=None, nlevel=14)¶ Save a georeferenced raster TIFF file.
Parameters: fname : str
Path of where to save the .tiff file.
coil : str, optional
Name of the coil to plot. By default, the first coil is plotted.
nx : int, optional
Number of points in x direction.
ny : int, optional
Number of points in y direction.
method : str, optional
Interpolation method (nearest, cubic or linear see scipy.interpolate.griddata) or IDW (default).
xmin : float, optional
Mininum X value.
xmax : float, optional
Maximum X value.
ymin : float, optional
Minimium Y value.
ymax : float, optional
Maximum Y value
color : bool, optional
If True a colormap will be applied.
cmap : str, optional
If color == True, name of the colormap. Default is viridis.
vmin : float, optional
Minimum value for colomap.
vmax : float, optional
Maximum value for colormap.
nlevel : int, optional
Number of level in the colormap. Default 7.
-
show
(coil='all', ax=None, vmin=None, vmax=None, dist=False)¶ Show the data.
Parameters: coil : str, optional
Specify which coil to plot. Default is all coils available.
ax : matplotlib.Axes, optional
If supplied, the graph will be plotted against ax.
vmin : float, optional
Minimal Y value.
vmax : float, optional
Maximial Y value.
dist : bool, optional
If True the true distance between points will be computed else the sample index is used as X index.
-
showDist
(coil=None, nbins=20, vmin=None, vmax=None, ax=None)¶ Display a histogram of the recorded values.
Parameters: coil : str or list of str, optional
By default all coils are stacked on each other. Otherwise specific coils (dataframe columns) can be specified.
nbins : int, optional
Number of bins to use for the histogram.
vmin : float, optional
Value of minimum bin limit.
vmax : float, optional
Value of maximum bin limit.
ax : matplotlib.Axes, optional
If given, the graph will be plotted against this axis.
-
showMap
(coil=None, contour=False, ax=None, vmin=None, vmax=None, pts=False, cmap='viridis_r', xlab='x', ylab='y', levels=[])¶ Display a map of the measurements.
Parameters: coil : str, optional
Name of the coil to plot. By default, the first coil is plotted.
contour : bool, optional
If True filled contour will be plotted using tricontourf.
ax : Matplotlib.Axes, optional
If specified, the graph will be plotted against the axis.
vmin : float, optional
Minimum of the colorscale.
vmax : float, optional
Maximum of the colorscale.
pts : bool, optional
If True the measurements location will be plotted on the graph.
xlab : str, optional
X label.
ylab : str, optional
Y label.
levels list of float, optional
If contour == True, levels are the contour intervals.
-
tcorrECa
(tdepths, tprofile)¶ Temperature correction using XXXX formula.
Parameters: tdepths : array-like
Depths in meters of the temperature sensors (negative downards).
tprofile : array-like
Temperature values corresponding in degree Celsius.
-
-
emagpy.Survey.
clipConvexHull
(xdata, ydata, x, y, z)¶ set to nan data outside the convex hull of xdata, yxdata
- Parameters:
- xdata, ydata (arrays) : x and y position of the data collected x, y (arrays) : x and y position of the computed interpolation points z (arrays) : value of the interpolation
- Return:
- znew (array) : a copy of z with nan when outside convexHull
-
emagpy.Survey.
convertFromCoord
(df, targetProjection=None)¶ Convert coordinates string (NMEA or DMS) to selected CRS projection.
Parameters: df : pandas.DataFrame
Dataframe with a ‘Latitude’ and ‘Longitude’ columns that contains NMEA string.
targetProjection : str, optional
Target CRS, in EPSG number: e.g. targetProjection=’EPSG:27700’ for the British Grid. If None, only WGS84 decimal degree will be returned.