API reference

This file is part of the ResIPy project (https://gitlab.com/hkex/resipy). @licence: GPLv3 @author: ResIPy authors and contributors

The ‘Project’ class wraps all main interactions between R* executables and other filtering or meshing part of the code. It’s the entry point for the user.

class resipy.Project.Project(dirname='', typ='R2')

Master class to handle all processing around the inversion codes.

Parameters:
dirnamestr, optional

Path of the working directory. Can also be set using R2.setwd().

typstr, optional

Either R2 or R3t for 3D. Complex equivalents are cR2 and cR3t. Automatically infered when creating the survey.

addData(index=0, **kwargs)

Adds data to the survey - used usually to add reciprocal datasets

Parameters:
index: int

Survey index to add data to.

**kwargs: Keyword arguments to be passed to Survey.addData()
addFilteredIP()

Add filtered IP to the dataframes.

addFlatError(percent=2.5)

Add a flat percentage error to resistivity data (for each survey in the class). This action is irreversable. Use case is for simulating modelling errors in unconventional surveys.

resError = res*(percent/100) + resError

Parameters:
percentfloat

Error in percent.

addRegion(xz, res0=100, phase0=1, blocky=False, fixed=False, ax=None, iplot=False)

Add region according to a polyline defined by xz and assign it the starting resistivity res0.

Parameters:
xzarray

Array with two columns for the x and y coordinates.

res0float, optional

Resistivity values of the defined area.

phase0float, optional

Read only if you choose the cR2 option. Phase value of the defined area in mrad

blockybool, optional

If True the boundary of the region will be blocky if inversion is block inversion.

fixedbool, optional

If True, the inversion will keep the starting resistivity of this region.

axmatplotlib.axes.Axes

If not None, the region will be plotted against this axes.

iplotbool, optional

If True , the updated mesh with the region will be plotted.

checkTxSign()

Checking and correcting the polarity of the transfer resistances (flat 2D surveys only !).

computeAttribute(formula, name, dump=None)

Compute a new attribute for each meshResults.

NOTE: this function present a security risk as it allows execution of python code inputed by the user. It should be used with caution.

Parameters:
formulastr

Formula as a string. All attribute already in the mesh object are available from the x dictionary and can be sed as x[‘nameOfAttribute’]. e.g. : “1/x[‘Resistivity’]” will compute the inverse of the ‘Resistivity’ attribute if this attribute exists. Operators available are addition (+), soustraction (-), division (/) multiplication (*) and exponent (**). Parenthesis are also taken into account. numpy function such as np.sqrt() or np.cos() are also available. e.g. : “(np.sqrt(x[‘Resistivity’]) + 100)*1000” Note: use a mixture of double and single quote such as double for the entire string and single for indexing the dictionnary x[‘keyName’].

namestr

Name of the new attribute computed.

dumpfunction, optional

Function to write stdout.

computeFineMeshDepth()

Compute the Fine Mesh Depth (FMD) based on electrode positions and the larger dipole spacing. Express as a positive number, it represents the relative vertical distance to extend the fine mesh region below the surface.

computeModelError(rmTree=False, dump=None)

Compute modelling error associated with the mesh for a half space problem. This is computed on a flat mesh with the same meshing parameters as the project mesh. In the case of non-conventional surveys a different modelling error scheme will be used.

Parameters:
rmTreebool

Remove the working directory used for the error modelling. Default is True.

dumpbool

Function to direct the output of the modelling error process, by default the outputs are printed to the console.

computeReciprocal(alg='Bisection Search', forceSign=False)

Compute Reciprocals and store them in self.surveys[0:n].df.

Parameters:
algstr, optional

Algorithm used to compute reciprocals. Choose between ‘Bisection Search’, ‘Pandas Merge’ or ‘Array Expansion’. The default is ‘Bisection Search’, other string are casted to ‘Array Expansion’.

forceSign: bool, optional

Force reciprocal and forward measurements to have the same polarity regarding the calculation of the reciprocal errors. Default is False.

Notes

If the relevant cython code cannot be found then the function falls back to using the pandas merge approach.

computeVol(attr='Resistivity(ohm.m)', vmin=None, vmax=None, index=0)

Given a 3D dataset, calculates volume based on mesh type. Note: The majority of this function’s code is also seen in the meshTools’ cellArea function as well. But that calculates all the mesh not the mesh results.

Parameters:
attr: str, optional

Attribute displayed in mesh results where apropriate cell vertices are assigned a location. These locations are utilized for volume calculations

vmin: float, optional

User assigned minimum threshold value of attribute

vmax: float, optional

User asigned maximum threshold value of attribute

index: int, optional

meshResults index. In case of multi-file inversion/modeling (e.g., timelapse)

Returns:
volume: float

The volume of selected mesh region in cubic meters

convertLocalGrid()

Converts UTM grid to local grid for stability in mesh generation and the finite element solution.

Parameters:
resetbool, optional

reset the grid to original

create3DSurvey(fname, lineSpacing=1, zigzag=False, ftype='Syscal', name=None, parser=None, dump=None)

Create a 3D survey based on 2D regularly spaced surveys.

Note: If one electrode has similar position between multiple lines (shared electrode), ResIPy will keep online one position of the electrode and share its label accross the lines.

Parameters:
fnamelist of str

List of 2D filenames in the right order for the grid or directory name (the files will be sorted alphabetically in this last case).

lineSpacingfloat, optional

Spacing in meter between each line.

zigzagbool, optional

If True then one survey out of two will be flipped. #TODO not implemented yet

ftypestr, optional

Type of file to be parsed. Either ‘Syscal’,’ProtocolDC’,’ResInv’, ‘BGS Prime’, ‘ProtocolIP’, ‘Sting’, ‘ABEM-Lund’, ‘Lippmann’ or ‘ARES’.

namestr, optional

Name of the merged 3D survey.

createBatchSurvey(dirname, ftype='Syscal', info={}, spacing=None, parser=None, isurveys=[], dump=None, debug=False)

Read multiples files from a folders (sorted by alphabetical order).

Parameters:
dirnamestr

Directory with files to be parsed.

ftypestr, optional

Type of file to be parsed. Either ‘Syscal’,’ProtocolDC’,’ResInv’, ‘BGS Prime’, ‘ProtocolIP’, ‘Sting’, ‘ABEM-Lund’, ‘Lippmann’ or ‘ARES’.

infodict, optional

Dictionnary of info about the survey.

spacingfloat, optional

Electrode spacing to be passed to the parser function.

parserfunction, optional

A parser function to be passed to Survey constructor.

isurveyslist, optional

List of surveys index that will be used for error modelling and so reciprocal measurements. By default all surveys are used.

dumpfunction, optional

Function to dump the information message when importing the files.

debugbool, optional

If True informations about reciprocal computation, default filtering and so on will be displayed.

createMergedSurveys(fname, ftype='Protocol DC', delimiter=',', dump=None, debug=False)

Create one or multiple surveys by merging multiple files referenced in singular human readable file with a specific format, survey details. See notes below.

Parameters:
fnamestr

File path to .csv file.

ftype: str, optional

File type, see ResIPy docs for types of file type avialable. Will be overwritten if specified in the survey details file.

delimiter: str, optional

Delimiter used in merge file. Defaults to ‘,’.

Notes

Format of survey details file should have up to 3 columns with these names. fpath, sid (optional), and ftype (optional).

fpath: str

Path to survey file, best to use full system path

ftype: str, optional

File type, should correspond to the file ‘ftype’ for each file (they might be different for some reason, though best avoided). If not passed the file type can be set be the variable at the top of this function. Type of files avialable are: Either ‘Syscal’,’ProtocolDC’, ‘ResInv’,’BGS Prime’, ‘ProtocolIP’, ‘Sting’, ‘ABEM-Lund’, ‘Lippmann’ or ‘ARES’.

sid: int, optional

Survey index, used for making timelapse surveys from multiple files

createMesh(typ='default', buried=None, surface=None, cl_factor=2, cl=-1, dump=None, res0=100, show_output=False, fmd=None, remote=None, refine=0, **kwargs)

Create a mesh.

Parameters:
typstr, optional

Type of mesh. For 2D:

  • ‘quad’: quadrilateral (fast to build)

  • ‘trian’: triangular (fast to run) - default

  • ‘circle’: close circular mesh

For 3D:
  • ‘tetra’: tetrahedral mesh for half-space - default

  • ‘cylinder’: column build with tetrahedral

  • ‘prism’: column build with prism

  • ‘tank’: closed geometry with tetrahedra

buriednumpy.array, optional

Boolean array of electrodes that are buried. Should be the same length as R2.elec

surfacenumpy.array, optional

Array with two or three columns x, y (optional) and elevation for additional surface points.

cl_factorfloat, optional

Characteristic length factor. Only used for triangular mesh to allow mesh to be refined close the electrodes and then expand.

clfloat, optional

Characteristic length that define the mesh size around the electrodes.

dumpfunction, optional
Function to which pass the output during mesh generation. print()

is the default.

res0float, optional

Starting resistivity for mesh elements.

show_outputbool, optional

If True, the output of gmsh will be shown on screen.

fmdfloat, optional

Depth of fine region specifies as a positive number relative to the mesh surface.

remotebool, optional

Boolean array of electrodes that are remote (ie not real). Should be the same length as Project.elec.

refineint, optional

Number times the mesh will be refined. Refinement split the triangles or the tetrahedra but keep the same number of parameter for the inversion. This helps having a more accurate forward response and a faster inversion (as the number of elements does not increase). Only available for triangles or tetrahedral mesh.

kwargs-, optional

Keyword arguments to be passed to mesh generation schemes Specific for ‘tank mesh’:

  • originlist of 3 floats

    Origin in X,Y,Z of one of the tank corner.

  • dimensionlist of 3 floats

    Dimension from the corner on how to extend the tank.

Specific for ‘cylinder mesh’:
  • zlim : list of 2 int

For the bottom and top of the column along the Z axis.

createModel(ax=None, dump=None, typ='poly', addAction=None)

Interactive model creation for forward modelling.

Parameters:
axmatplotlib.axes.Axes, optional

Axes to which the graph will be plotted.

dumpfunction, optional

Function that outputs messages from the interactive model creation.

typstr

Type of selection either poly for polyline or rect for rectangle.

addActionfunction

Function to be called once the selection is finished (design for GUI purpose).

Returns:
figmatplotlib.figure

If ax is None, will return a figure.

createModelErrorMesh(**kwargs)

Create an homogeneous mesh to compute modelling error.

Same arguments as Project.createMesh().

createModelMesh(**kwargs)

Create a triangular mesh given the designed geometry by R2.designModel().

Parameters:
All parameters to be passed are similar to `R2.createMesh()`.
createMultiMesh(runParallel=True, **kwargs)

Create multiple meshes from avalable Projects in self.projs.

Parameters:
runParallelbool, optional

if True, mesh generation will run in multiple threads.

kwargs-

Keyword arguments to be passed to mesh generation schemes

createPseudo3DSurvey(dirname, lineSpacing=1, ftype='Syscal', parser=None, **kwargs)
Create a pseudo 3D survey based on 2D surveys. Multiple 2D Projects to be turned into a single pseudo 3D survey.

THIS WILL NEED CORRECT ELECTRODE LAYOUT - DONE IN self._updatePseudo3DSurvey()

Parameters:
dirnamelist of str

List of 2D filenames in the right order for the grid or directory name (the files will be sorted alphabetically in this last case).

lineSpacingfloat, optional

Spacing in meter between each line.

ftypestr, optional

Type of the survey to choose which parser to use.

kwargs-

Keyword arguments to be passed to Project.createBatchSurvey()

createSequence(params=[('dpdp1', 1, 8)], seqIdx=None, *kwargs)

Creates a forward modelling sequence, see examples below for usage.

Parameters:
paramslist of tuple, optional

Each tuple is the form (<array_name>, param1, param2, …) Types of sequences available are : ‘dpdp1’,’dpdp2’,’wenner_alpha’, ‘wenner_beta’, ‘wenner_gamma’, ‘schlum1’, ‘schlum2’, ‘multigrad’, ‘custSeq’. if ‘custSeq’ is chosen, param1 should be a string of file path to a .csv file containing a custom sequence with 4 columns (a, b, m, n) containing forward model sequence.

seqIdx: list of array like, optional

Each entry in list contains electrode indices (not label and string) for a given electrode string which is to be sequenced. The advantage of a list means that sequences can be of different lengths.

Examples

>>> k = Project()
>>> k.setElec(np.c_[np.linspace(0,5.75, 24), np.zeros((24, 2))])
>>> k.createMesh(typ='trian')
>>> k.createSequence([('dpdp1', 1, 8), ('wenner_alpha', 1), ('wenner_alpha', 2)]) # dipole-dipole sequence
>>> k.createSequence([('custSeq', '<path to sequence file>/sequence.csv')]) # importing a custom sequence
>>> seqIdx = [[0,1,2,3],[4,5,6,7],[8,9,10,11,12]]
createSequenceXBH()

Custom scheme for boreholes (not yet developed)

createSurvey(fname='', ftype='Syscal', info={}, spacing=None, parser=None, debug=True, estMemory=True, **kwargs)

Read electrodes and quadrupoles data and return a survey object.

Parameters:
fnamestr

Filename to be parsed.

ftypestr, optional

Type of file to be parsed. Either ‘Syscal’,’ProtocolDC’,’ResInv’, ‘BGS Prime’, ‘ProtocolIP’, ‘Sting’, ‘ABEM-Lund’, ‘Lippmann’ or ‘ARES’.

infodict, optional

Dictionnary of info about the survey.

spacingfloat, optional

Electrode spacing to be passed to the parser function.

parserfunction, optional

A parser function to be passed to Survey constructor.

debugbool, optional

If True, information about the reciprocal measurements, default filtering, etc. will be displayed.

estMemory: bool, optional

If true, estimate the amount of RAM required to do the inversion. Default is True.

**kwargs: Keyword arguments to be passed to Survey()
createTimeLapseSurvey(dirname, ftype='Syscal', info={}, spacing=None, parser=None, isurveys=[], dump=None, debug=False)

Read electrodes and quadrupoles data and return a survey object.

Parameters:
dirnamestr or list of str

Directory with files to be parsed or list of file to be parsed.

ftypestr, optional

Type of file to be parsed. Either ‘Syscal’,’ProtocolDC’,’ResInv’, ‘BGS Prime’, ‘ProtocolIP’, ‘Sting’, ‘ABEM-Lund’, ‘Lippmann’ or ‘ARES’.

infodict, optional

Dictionnary of info about the survey. Put inverse_type = 1 to allow for a changing number of measurements between surveys.

spacingfloat, optional

Electrode spacing to be passed to the parser function.

parserfunction, optional

A parser function to be passed to Survey constructor.

isurveyslist, optional

List of surveys index that will be used for error modelling and so reciprocal measurements. By default all surveys are used.

dumpfunction, optional

Function to dump information message when importing the files.

debugbool, optional

If True informations about reciprocal computation, default filtering and so on will be displayed.

designModel(ax=None, dump=<built-in function print>, typ='poly', addAction=None, fmd=None)

Interactive model design for forward modelling (triangular only). As opposite to R2.createModel(). R2.designModel() allows to draw mesh region before meshing. This allows to have straight boundaries for triangular mesh.

Parameters:
axmatplotlib.axes.Axes, optional

Axes to which the graph will be plotted.

dumpfunction, optional

Function that outputs messages from the interactive model creation.

typstr

Type of selection either poly for polyline or rect for rectangle.

addActionfunction

Function to be called once the selection is finished (design for GUI purpose).

fmdfloat, optional

Depth of of interest specifies as a relative positive number from the surface.

Returns:
figmatplotlib.figure

If ax is None, will return a figure.

detectStrings(tolerance=15, max_itr=None)

Automatically detect electrode strings. If all electrodes lie on the same Y plane then the function assumes the problem is a borehole problem.

Parameters:
tolerancefloat, optional

Maximum (+/-) bearing (ie directional angle) tolerance each subsiquent electrode may have. The default is 15.

max_itrint, optional

Maximum number of searches that can be performed to find colinear nieghbouring electrodes. The default is number of electrodes plus 10.

Returns:
list (of list)

Each entry in the list corresponds to an electrode string, and is a list of integers which are the indices of the respective electrodes in self.elec

Raises:
ValueError

if the change in x and y direction for 2 neighbouring electrodes is the same. ie no 2 electrodes can occupy the same x y position in this code.

ValueError

if the change maxium number of searches is exceeded.

elec2distance(yDominant=False, iMoveElec=False)

Convert 3D electrode XY coordinates into just X coordinates. Use for 2D lines only! If self.elec has been set then each survey will use the electrodes set in the R2 master class. If not then the R2 master class will take on the elec values set for the first survey in a sequence.

Parameters:
yDominant: bool, optional

If electrodes are prodimently spaced in the y direction then set yDominant to True.

iMoveElec: bool, optional

If moving electrodes are present then set to True, so that the same electrode positions are not given to each survey.

elec2horidist()

Convert 2D xz data into a true horizontal distance. Assumes that survey was done with a tape measure and the X distances are not true horizontal distance but distances measured along the ground.

estimateError(a_wgt=0.01, b_wgt=0.02)

Estimate reciprocal error data for data with no reciprocals for each survey, using the same routine present in R2. This allows for the additional inclusion of modelling errors. It could be used when the user want to assign invidual errors based on a_wgt/b_wgt. This action is irreversable.

Parameters:
a_wgt: float, optional

a_wgt documented in the R2 documentation

b_wgt: float, optional

b_wgt documented in the R2 documentation

exportData(outputname=None, ftype='protocol', err=False, recip=False)

Export preconditioned data used by ResIPy into another format ( different to project.saveData).

Parameters:
outputnamestr, optional

Outputname. The default is None. If not given then the function falls back on the survey name. If set then a number will be appended to the file name in the case of time-lapse (or batch) surveys.

ftypestr, optional

Export File type, choose from either protocol, srv, csv, ResInv. The default is ‘protocol’.

errbool, optional

Flag to include errors. The default is False.

recipbool, optional

Flag to include reciprocals. The default is False.

exportElec(outputname=None, _forceLocal=False)

Export electrodes as a different file format, with coordinate conversion if Project.coordLocal set to True.

Parameters:
outputnamestr, optional
Output path with extension. Available mesh format are:
  • .vtk (Paraview)

  • .dat (R* codes)

  • .csv (general)

If not provided the electrode coordinates are saved in the working directory as electrodes.vtk.

_forceLocal: bool, optional

Forces outputs in the local grid format is self.coordLocal == True.

exportMesh(outputname=None, voxel=False, _forceLocal=False)

Export mesh as a different file format, with coordinate conversion if Project.coordLocal set to True.

Parameters:
outputnamestr, optional
Output path with extension. Available mesh format are:
  • .vtk (Paraview)

  • .node (Tetgen)

  • .dat (R* codes)

If not provided the mesh is saved in the working directory as mesh.vtk.

voxel: bool, optional

If True, mesh will be converted to voxel format.

_forceLocal: bool, optional

Forces outputs in the local grid format is self.coordLocal == True.

exportMeshResults(dirname=None, ftype='.vtk', voxel=False, _forceLocal=False)

Save mesh files of inversion results to a specified directory. If results are in a local grid, they will be converted back into thier respective utm or nationa grid coordinates.

Parameters:
dirname: str

Directory in which results will be saved. Default is the working directory.

ftype: str

Type of file extension.

voxel: bool, optional

Force export to be in a voxel format.

_forceLocal: bool, optional

Force output to be in a local grid if self.coordLocal is True. Meant for use with the UI.

filterAppResist(index=-1, vmin=None, vmax=None)

Filter measurements by apparent resistivity for surface surveys

Parameters:
vminfloat, optional

Minimum value.

vmaxfloat, optional

Maximum value.

indexint, optional

Index of the survey on which to apply the processing. If the processing is to be applied to all surveys then specifiy index=-1 (default).

filterContRes(index=-1, vmin=None, vmax=None)

Filter measurements by contact resistance if avialable.

Parameters:
vminfloat, optional

Minimum value.

vmaxfloat, optional

Maximum value.

indexint, optional

Index of the survey on which to apply the processing. If the processing is to be applied to all surveys then specifiy index=-1 (default).

filterDCA(index=-1, dump=None)

Execute DCA filtering. Decay Curve Analysis (DCA) based on. Flores Orozco, A., Gallistl, J., Bücker, M., & Williams, K. H. (2017)., Decay curve analysis for data error quantification in time-domain induced polarization imaging., Geophysics, 83(2), 1–48. https://doi.org/10.1190/geo2016-0714.1

Parameters:
indexint, optional

Index of the survey to use for processing. Default index == -1 will apply the processing to all surveys.

dumpfunction, optional

Function onto pass the progress.

filterDummy(index=-1)

Remove measurements where abs(a-b) != abs(m-n) (likely to be dummy measurements added for speed).

Parameters:
indexint, optional

Index of the survey to process. If index == -1 (default) then the processing is applied on all survey independantly.

filterElec(elec=[], index=-1)

Filter out measurements associated with specific electrodes.

Parameters:
eleclist

List of electrode number to be removed.

indexint, optional

Index of the survey on which to apply the processing. If the processing is to be applied to all surveys then specifiy index=-1 (default).

filterInvError(index=-1, vmin=None, vmax=None)

Remove measurements where inversion error is outside of the defined range.

Parameters:
vminfloat, optional

minimum value of normalized error below which data will be discarded.

vmaxfloat, optional

maximum value of normalized error above which data will be discarded.

indexint, optional

Index of the survey to process. If index == -1 (default) then the processing is applied on all survey independantly.

filterManual(index=-1, ax=None, **kwargs)

Interactive manually filters the data visually. The manually selected points index are stored in Survey.iselect or Survey.eselect``if it is an electrodes. Use `Survey.filterData() to filter them out for a single survey. Or R2._filterSimilarQuads() to filter quadrupoles amongs all R2.surveys.

filterNegative()

Remove negative apparent resistivty values

filterNested(index=-1)

Removes nested measurements: Where M or N are in between A and B.

Parameters:
indexint, optional

Index of the survey to fit. If index == -1 (default) then the fit is done on all surveys independantly. If ìndex == -2 then the fit is done on the combined surveys.

filterRangeIP(index=-1, phimin=None, phimax=None)

Filter IP data according to a specified range.

Parameters:
indexint, optional

Index of the survey to fit. If index == -1 (default) then the fit is done on all surveys independantly.

phiminfloat, optional

Minimium phase angle [mrad].

phimaxfloat, optional

Maximum phase angle [mrad].

filterRecip(percent=20, index=-1)

Filter on reciprocal errors.

Parameters:
indexint, optional

Index of the survey on which to apply the processing. If the processing is to be applied to all surveys then specifiy index=-1 (default).

percentfloat, optional

Percentage of reciprocal error above witch a measurement will be discarded. 20% by default.

filterRecipIP(index=0)

Remove reciprocal for IP data ONLY. Additional arguments to be passed to :func: ~resipy.Survey.filterRecipIP.

filterStack(percent=2, index=-1)

Filter on stacking (dev) errors.

Parameters:
indexint, optional

Index of the survey on which to apply the processing. If the processing is to be applied to all surveys then specifiy index=-1 (default).

percentfloat, optional

Percentage of stacking error above witch a measurement will be discarded. 2% by default.

filterTransferRes(index=-1, vmin=None, vmax=None)

Filter measurements by transfer resistance.

Parameters:
vminfloat, optional

Minimum value.

vmaxfloat, optional

Maximum value.

indexint, optional

Index of the survey on which to apply the processing. If the processing is to be applied to all surveys then specifiy index=-1 (default).

filterUnpaired(index=-1)

Remove quadrupoles that don’t have reciprocals. This might remove dummy measurements added for sequence optimization.

Parameters:
indexint, optional

Index of the survey on which to apply the processing. If the processing is to be applied to all surveys then specifiy index=-1 (default).

filterZeroMeasSurveys()

Filter out badly behaved surveys, where after all other QC no measurements are actually left.

fitErrorLME(index=-1, ax=None, rpath=None, iplot=True)

Fit a linear mixed effect (LME) model by having the electrodes as as grouping variables.

Parameters:
Index of the survey to fit. If `index == -1` (default) then the fit

is done on all surveys independantly. If ìndex == -2 then the fit is done on the combined surveys.

axmatplotlib.Axes, optional

If specified, the graph will be plotted against this axis, otherwise a new figure will be created.

rpathstr, optional

Path of the directory with R (for Windows only).

iplotbool, optional

If True plot it.

fitErrorLin(index=-1, ax=None)

Fit a linear relationship to the resistivity data.

Parameters:
indexint, optional

Index of the survey to fit. If index == -1 (default) then the fit is done on all surveys independantly. If ìndex == -2 then the fit is done on the combined surveys.

axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

fitErrorParabolaIP(index=-1, ax=None)

Plot the reciprocal phase errors with a parabola fit.

Parameters:
indexint, optional

Index of the survey to fit. If index == -1 (default) then the fit is done on all surveys independantly. If ìndex == -2 then the fit is done on the combined surveys.

axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

fitErrorPwl(index=-1, ax=None)

Fit an power law to the resistivity data.

Parameters:
indexint, optional

Index of the survey to fit. If index == -1 (default) then the fit is done on all surveys independantly. If ìndex == -2 then the fit is done on the combined surveys.

axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

fitErrorPwlIP(index=-1, ax=None)

Plot the reciprocal phase errors with a power-law fit.

Parameters:
indexint, optional

Index of the survey to fit. If index == -1 (default) then the fit is done on all surveys independantly. If ìndex == -2 then the fit is done on the combined surveys.

axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

fixLegendItems(ax)

Display legend items with survey names in the case of fitting individual error models to multiple data sets.

Parameters:
axmatplotlib axes

Axes handle with reciprocal error plots in it.

forward(noise=0.0, noiseIP=0.0, iplot=False, dump=None)

Operates forward modelling.

Parameters:
noisefloat, optional 0% <= noise <= 100%

Noise level in percent from a Gaussian distribution that should be applied on the forward apparent resistivities obtained.

noiseIPfloat, optional

Absolute noise level in mrad from a Gaussian distribution that should be applied on the forward phase values obtained.

iplotbool, optional

If True will plot the pseudo section after the forward modelling.

dumpfunction, optional

Function to print information messages when running the forward model.

generateElec(nb=24, dx=0.5, dz=0, nline=1, lineSpacing=2)

Generate electrodes positions for 2D and 3D surveys.

Parameters:
nbint, optional

Number of electrodes per line. For 3D survey, if multiple lines, the total number of electrodes will be nb x nline.

dxfloat, optional

Spacing in meters between electrodes in the X direction.

dzfloat, optional

Increment in the Z direction (elevation) between consecutive electrodes.

nlineint, optional

Number of lines. 1 for 2D and multiple for 3D.

lineSpacingfloat, optional

Spacing between lines (3D only).

getR2out()

Reat the .out file and parse its content.

Returns:
Dataframe with the dataset name, and the RMS decrease for each iteration.
getResults(dirname=None)

Collect inverted results after running the inversion and adding them to R2.meshResults list.

Parameters:
dirnamestr, optional

If specified, dirname will be used as the working directory (this is needed for R2.loadResults()). Default is self.dirname.

hasElecString()

Determine if a electrode strings are present in the electrode labels

Returns:
bool

True if strings present in electrode label

importElec(fname='')

Import electrodes positions.

Parameters:
fnamestr

Path of the CSV (or file containing the electrodes positions. It should contains 3 columns maximum with the X, Y, Z positions of the electrodes.

importMesh(file_path, mesh_type=None, node_pos=None, elec=None, order_nodes=True, res0=100)

Import mesh from .vtk / .msh / .dat, rather than having ResIPy create one for you.

Parameters:
file_pathstr

File path mapping to the mesh file

mesh_typestr

Not used anymore.

node_posarray like, optional

Array of ints referencing the electrode nodes. If left as none no electrodes will be added to the mesh class. Consider using mesh.moveElecNodes() to add nodes to mesh using their xyz coordinates.

elecarray, optional

N*3 numpy array of electrode x,y,z coordinates. Electrode node positions will be computed by finding the nearest nodes to the relevant coordinates.

res0float, optional

Starting resistivity for mesh elements.

importPseudo3DElec(fname='')

Import electrodes positions. The label columns should include line number separated by space (like in 3D):

label,x,y,z 1 3,0,0,0 1 4,1,1,0 1 5,1,2,1

Parameters:
fnamestr

Path of the CSV file containing the electrodes positions. It should contains 3 columns maximum with the X, Y, Z positions of the electrodes.

importSequence(fname='')

Import sequence for forward modelling.

Parameters:
fnamestr

Path of the CSV file to be imported. The file must have 4 columns with headers (a, b, m, n) containing 4 electrodes numbers.

invert(param={}, iplot=False, dump=None, err=None, modErr=False, parallel=False, iMoveElec=False, ncores=None, rmDirTree=True, modelDOI=False)

Invert the data, first generate R2.in file, then run inversion using appropriate wrapper, then return results.

Parameters:
paramdict, optional

Dictionary of parameters for inversion. Will be passed to R2.write2in().

iplotbool, optional

If True, will plot the results of the inversion using R2.showResults().

dumpfunction, optinal

Function to print the output of the inversion. To be passed to R2.runR2().

errbool, optional

If ‘True’ reciprocal error model will be used to weight the inversion. Can be used to overide self.err. Default is None (in which case self.err sets the inclusion of error weighting).

modErrbool, optional

If True, the model error will be compute and added before inversion.

parallelbool, optional

If True, batch and time-lapse survey will be inverted in //. No output will be display during inversion.

iMoveElecbool, optional

If True, then different electrode location will be used for the different surveys. Electrodes location are specified in the Survey object. Only for parallel inversion for now.

ncoresint, optional

If parallel==True then ncores is the number of cores to use (by default all the cores available are used).

rmDirTreebool, optional

Remove excess directories and files created during parallel inversion

modelDOIbool, optional

If True, the Depth of Investigation will be model by reinverting the data on with an initial res0 different of an order of magnitude. Note that this option is only available for single survey.

invertPseudo3D(invLog=None, runParallel=False, **kwargs)

Run pseudo3D inversions.

Parameters:
invLogfunction, optional

Passes project inversion outputs.

runParallelbool

if True, inversions will run in parallel based on number of CPU cores.

kwargs-

Keyword arguments to be passed to invert().

loadProject(fname)

Load data from project file.

Parameters:
fnamestr

Path where the file will be saved.

loadResults(invdir)

Given working directory, will attempt to load the results of an already run inversion.

Parameters:
invdirstr

Path to the inversion directory.

matchSurveys()

Will trim all surveys to get them ready for difference inversion so that each PAIRS of (background, surveyX) have the same number of quadrupoles. We do not take all quadrupoles in common among all surveys as this is not needed and if there is a small survey, it would reduce all other larger surveys.

mergeElec(dist=-1)

Merge electrodes that have less than a certain distance to eache other

Parameters:
distfloat, optional

maximum distance of close electrodes in meters (electrodes that have a distance less than dist will be merged) -1 flag for auto calculation where dist = average(electrode spacing)/100

Returns:
bool

True: merge successful False: merge unsuccessful

modelDOI(dump=None)

Will rerun the inversion with a background constrain (alpha_s) with the normal background and then a background 10 times more resistive. From the two different inversion a senstivity limit will be computed.

normElecIdx()

Normalise electrode indexes to start at 1 in consective and ascending order.

postProcTl()

Post processing for time-lapse surveys.

Compute the absolute and the relative difference in resistivity between inverted surveys.

resetRegions()

Just reset all regions already drawn. Shouldn’t be needed as the self.runR2() automatically use a homogenous model when starting for inversion. The only purpose of this is to use an inhomogeous starting model to invert data from forward modelling.

runParallel(dirname=None, dump=None, iMoveElec=False, ncores=None, rmDirTree=True)

Run several instances of R2 in parallel according to the number of cores available.

Parameters:
dirnamestr, optional

Path of the working directory.

dumpfunction, optional

Function to be passed to R2.runR2() for printing output during inversion.

iMoveElecbool, optional

If True will move electrodes according to their position in each Survey object.

ncoresint, optional

Number or cores to use. If None, the maximum number of cores available will be used.

rmDirTree: bool, optional

Remove excess directories and files created during parallel. Default is True.

runR2(dirname='', dump=None)

Run the executable in charge of the inversion.

Parameters:
dirnamestr, optional

Path of the directory where to run the inversion code.

dumpfunction, optional

Function to print the output of the invrsion code while running.

saveCwd(outputdir)

Save all ouputs (_res.dat, .vtk, …) from the working directory generated during inversion to the designated directory.

Parameters:
outputdirstr

Path to the directory to save the files.

saveData(outputdir)

Saves all data generated by ResIPy in the current working directory to a specified folder. Serves same function as Project.saveCwd() but retained for backwards compability.

Parameters:
outputdirstr

Path to the directory to save the files.

saveErrorData(fname)

Save quadruople, resistance, phase and their respective reciprocal errors as .csv file.

Parameters:
fnamestr

Path where to save the file.

saveFilteredData(fname, savetyp='Res2DInv (*.dat)')

Save filtered data in formats to be used outside ResIPy (e.g. Res2DInv).

Parameters:
fnamestr

Path where to save the file.

savetypstr, optional

Saving format. To be determined in GUI. Default: Res2DInv (.dat)

saveForwardModelResult(fname)

Save the result of a forward model run to a specific file name/

Parameters:
fnamestr

path to file.

saveInvPlots(outputdir=None, **kwargs)

Save all plots to output (or working directory). Parameters are passed to the showResults() method.

Parameters:
outputdirstr, optional

Path of the output directory. Default is the working directory.

saveMesh(outputname=None)

Save mesh as .vtk to be viewed in paraview.

Parameters:
outputnamestr, optional
Output path with extension. Available mesh format are:
  • .vtk (Paraview)

  • .node (Tetgen)

  • .dat (R* codes)

If not provided the mesh is saved in the working directory as mesh.vtk.

saveMeshVtk(outputname=None)

Save mesh as .vtk to be viewed in paraview.

Parameters:
outputnamestr, optional

Output path of the .vtk produced. By default the mesh is saved in the working directory self.dirname as mesh.vtk.

saveProject(fname)

Save the current project will all dataset in custom ResIPy format (.resipy) for future importation.

saveSequence(fname='')

Save sequence as .csv file.

Parameters:
fnamestr, optional

Path where to save the sequence.

saveVtks(dirname=None)

Save vtk files of inversion results to a specified directory.

Parameters:
dirname: str

Directory in which results will be saved. Default is the working directory.

setBorehole(val=False)

Set all surveys in borehole type if True is passed.

setCoordConv(flag=False, x0=None, y0=None, a=None)

Set that the project should be imported and exported according to a coordinate conversion system. Generally UTM coordinates will cause stability issues with mesh generation and the finite element solutions used in the R* codes, hence the we use a local coordinate system. Call this function before creating the mesh.

Parameters:
flagbool, optional

Flag if coordinate conversion system is required. The default is False. If flag is false but Project.coordLocal is True then the function can be used to rvert electrode coordinates back to their initial position.

x0float, optional

Reference X coordinate in UTM system. The default is None. If left as none, then the minimum X coordinate of the electrodes will be used.

y0float, optional

Reference Y coordinate in UTM system. The default is None. If left as none, then the minimum Y coordinate of the electrodes will be used.

afloat, optional

Rotation angle in degrees. The default is 0.

setElec(elec, elecList=None, _uiOverideCoordLocal=False)

Set electrodes. Automatically identified remote electrode.

Parameters:
elecnumpy array

Array of NxM dimensions. N = number of electrodes, M = 2 for x,z or M = 3 if x,y,z coordinates are supplied.

elecListlist, optional

If not None then elec is ignored in favor of elecList. This option is to be used in the advanced use case where electrodes move which each survey. Each entry of the list is a numpy array the same format of ‘elec’, which is then assigned to each survey class.

_uiOverideCoordLocal: bool, optional

UI command only, set to to True when updating the project electrodes in the UI.

static setNcores(ncores)

Set the number of cores to use for big calculations, for now this value only to mesh generation and calculations done on 3D meshes.

Parameters:
ncoresint

Number of cores to use

Raises
——
EnvironmentError

Raised if ncores is more than that avialable

setPseudo3DElec(elec)
Set pseudo 3D electrodes (with an electrode label as:

<line number> <electrode number>).

Parameters:
elecListlist of dataframes, optional

List of electrodes dataframes - each df must have 2D like XYZ (rotated to have y=0).

setRefModel(res0)

Set the reference model according to a previous inversion, avoids the need to invert reference model again for timelapse workflows. In contrast to R2.setStartingRes() which assign resistivity to group of elements, this method requires a vector of the same length as the number of elements. This enables, notably to manually perform consecutive background constrained inversion.

Parameters:
res0: array like

Array of resistivity values, ideally from a previous inversion. The length of this array should be the same as the number of elements.

setStartingRes(regionValues={}, zoneValues={}, fixedValues={}, ipValues={})

Assign starting resitivity values.

Parameters:
regionValuesdict, optional

Dictionnary with key being the region number and the value being the resistivity in [Ohm.m].

zoneValuesdict, optional

Dictionnary with key being the region number and the zone number. There would be no smoothing between the zones if ‘block inversion’ is selected (inversion_type = 4).

fixedValuesdict, optional

Dictionnary with key being the region number and a boolean value if we want to fix the resistivity of the zone to the starting one. Note that it only works for triangular mesh for now.

ipValuesdict, optional

Dictionnary with key being the region number and the values beeing the phase [mrad].

setTitle(linetitle)

Set the title of the survey name when inverting data. Input is a string.

setwd(dirname)

Set the working directory.

Parameters:
dirnamestr

Path of the working directory.

showError(index=0, ax=None)

Plot the reciprocal errors.

Parameters:
indexint, optional

Index of the survey to plot. If index == -1 then all combined data of all survey will be plotted together. Default is to plot the first survey (index==0).

axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

showErrorDist(index=0, ax=None)

Calculate and plots reciprocal error probablity histogram. Good data will have a bell shape (normal) distribution where most datapoints have near zero reciprocal error.

Parameters:
indexint, optional

Index of the survey to plot. Default is first survey index == 0. If index == -2 then the error distribution of the combined data will be plotted.

axMatplotlib.Axes

If specified, the graph will be plotted against it.

showErrorIP(index=0, ax=None)

Plot the reciprocal phase discrepancies against the reciprocal mean transfer resistance.

Parameters:
indexint, optional

Index of the survey to show. Default is the first survey index == 0. If ìndex == -2 then the combined data from all surveys are shown.

axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib.Figure, optional

If ax is not specified, the function will return a figure object.

showHeatmap(index=0, ax=None)

Plot a phase heatmap (x = M, y = A and value = -phi) based on: Orozco, A. F., K. H. Williams, and A. Kemna (2013), Time-lapse spectral induced polarization imaging of stimulated uranium bioremediation, Near Surf. Geophys., 11(5), 531–544, doi:10.3997/1873-0604.2013020)

Parameters:
indexint, optional

Index of the survey to fit. Default is the first survey index == 0.

axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

showInParaview(index=0, paraview_loc=None)

Open paraview to display the .vtk file.

Parameters:
index: int, optional

Timestep to be shown in paraview (for an individual survey this 1).

paraview_loc: str, optional

Windows ONLY maps to the excuatable paraview.exe. The program will attempt to find the location of the paraview install if not given.

showInvError(index=0, ax=None)

Display inversion error by measurment numbers.

Parameters:
indexint, optional

Index of survey (if time-lapse or batch). Default index == 0.

axmatplotlib axis

If provided, the graph will be plotted against this axis.

showIter(index=-2, ax=None, modelDOI=False, cropMaxDepth=False)

Dispay temporary inverted section after each iteration.

Parameters:
indexint, optional

Iteration number to show.

axmatplotib axis, optional

If specified, the graph will be plotted along ax.

modelDOIbool, optional

As modelDOI() is always computed using R2 (not cR2), this tells the method to look for an R2 looking iteration file.

cropMaxDepthbool, optional

if True, below max depth will be cropped

showMesh(ax=None, **kwargs)

Display the mesh and handles best default values.

showMeshInParaview(paraview_loc=None)

Show the mesh in paraview (mostly useful for 3D surveys.

Parameters:
paraview_loc: str, optional

Windows ONLY maps to the excuatable paraview.exe. The program will attempt to find the location of the paraview install if not given.

showParam()

Print parameters in R2.param dictionary.

showPseudo(index=0, vmin=None, vmax=None, ax=None, **kwargs)

Plot pseudo-section with dots.

Parameters:
indexint, optional

Index of the survey to be used for the pseudo-section (in case of timelapse or batch).

vminfloat, optional

Minimum value for colorscale.

vmaxfloat, optional

Maximum value for colorscale.

axmatplotlib.Axes, optional

If specified, axis along which to plot the graph.

**kwargsoptional

Passed to Survey.showPseudo().

showPseudo3DMesh(ax=None, color_map='Greys', meshList=None, cropMesh=True, color_bar=False, returnMesh=False, cropMaxDepth=False, clipCorners=False, pvshow=True, **kwargs)

Show 2D meshes in 3D view

Parameters:
axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

color_mapstr, optional

Name of the colormap to be used.

meshListlist of Mesh classes

If not None, pseudo 3D meshes will be plotted by default.

cropMeshbool, optional

If True, 2D mesh will be bound to electrodes and zlim.

color_barBoolean, optional

True to plot colorbar.

returnMesh: bool, optional

if True method returns a merged mesh.

cropMaxDepthbool, optional

If True, region below fine mesh depth (fmd) will be cropped.

clipCornersbool, optional

If ‘True’, triangles from bottom corners will be cropped (only if the whole mesh is not shown).

pvshowbool, optional

Set to False to not call Plotter.show() (useful for subplots).

kwargs-

Keyword arguments to be passed to Mesh.show() class.

showPseudo3DResults(cropMesh=False, **kwargs)

Show 2D Inversions in 3D view

Parameters:
cropMeshbool, optional

If True, 2D mesh will be bound to electrodes and zlim.

kwargs-

Keyword arguments to be passed to showPseudo3DMesh().

showPseudoIP(index=0, vmin=None, vmax=None, ax=None, **kwargs)

Plot pseudo-section with dots for IP data.

Parameters:
indexint, optional

Index of the survey to be used for the pseudo-section (in case of timelapse or batch).

vminfloat, optional

Minimum value for colorscale.

vmaxfloat, optional

Maximum value for colorscale.

axmatplotlib.Axes, optional

If specified, axis along which to plot the graph.

**kwargsoptional

Passed to Survey.showPseudoIP().

showPseudoInvError(index=0, ax=None, vmin=None, vmax=None, elec=True)

Plot pseudo section of errors from file f001_err.dat.

Parameters:
indexint, optional

Index of the survey (if time-lapse or batch). Default index == 0.

axmatplotlib axis

If specified, the graph will be plotted against ax.

vminfloat, optional

Min value.

vmaxfloat, optional

Max value.

elecbool, optional

If True, the electrodes are displayed and can be used for filtering.

showPseudoInvErrorIP(index=0, ax=None, vmin=None, vmax=None)

Display normalized phase error.

Parameters:
indexint, optional

Index of the survey (if time-lapse or batch). Default index == 0.

axmatplotlib axis

If specified, the graph will be plotted against ax.

vminfloat, optional

Min value.

vmaxfloat, optional

Max value.

showRMS(index=0, ax=None)

Show the RMS decrease for each iteration.

Parameters:
indexint, optional

Index of the dataset for which to plot the RMS.

axmatplotlib axis, optional

If provided, the graph will be plotted against it.

showResults(index=0, ax=None, edge_color='none', attr='', sens=True, color_map='viridis', zlim=None, clabel=None, doi=False, doiSens=False, contour=False, cropMaxDepth=True, clipContour=True, clipCorners=False, use_pyvista=True, background_color=(0.8, 0.8, 0.8), pvslices=([], [], []), pvspline=None, pvthreshold=None, pvgrid=True, pvcontour=[], pvdelaunay3d=False, volume=None, **kwargs)

Show the inverteds section.

Parameters:
indexint, optional

Index of the inverted section (mainly in the case of time-lapse inversion). If index == -1, then all 2D survey will be plotted on a 3D grid.

axmatplotlib axis, optional

If specified, the inverted graph will be plotted agains ax.

edge_colorstr, optional

Color of the edges of the mesh.

attrstr, optional

Name of the attribute to be plotted.

sensbool, optional

If True and if sensitivity is available, it will be plotted as a white transparent shade on top of the inverted section.

color_mapstr, optional

Name of the colormap to be used.

clabelstr, optional

Label of the colorbar (by default the label is the value of attr).

doibool, optional

If True, it will draw a dotted red line corresponding to 0.02 from the Oldenburg and Li method. Note that R2.modeDOI() needs to be run for that.

doiSensbool, optional

If True, it will draw a dashed line corresponding to 0.001 of the maximum of the log10 sensitivity.

contourbool, optional

If True, contours will be plotted.

cropMaxDepthbool, optional

If True, the mesh will be clipped with at a depth following the surface. If False, the mesh will be clipped at the maximum depth available. This doesn’t have any effect if clipContour is False.

clipContourbool, optional

If True, the contour of the area of interest will be clipped (default).

clipCornersbool, optional

If ‘True’, triangles from bottom corners will be cropped (only if the whole mesh is not shown).

use_pyvistabool, optional

(3D only) Use visual toolkit backend for displaying 3D mesh, note that pyvista must be installed for this to work.

background_colortuple, optional

(3D only) Background color assigned to pyvista plotter object when created. Not yet supported for matplotlib axis handles.

pvslicestuple of list of float, optional

(3D only) 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.

pvspline‘elec’ or numpy array, optional

(3D only) If ‘elec’ mesh will be sliced along the electrodes path otherwise an array of X, Y, Z of points on a path to slice the mesh along that path is needed.

pvthresholdlist of two floats, optional

(3D only) Keep values between pvthreshold[0] and pvthreshold[1].

pvgridbool, optional

(3D only) Show grid or not.

pvcontour: list of float, optional

(3D only) Values of the isosurface to be plotted.

pvdelaunay3dbool, optional

If True a “Delaunay 3D” triangulation filter will be applied on the mesh.

volumefloat, optional

If not ‘None’ then volume of selected region of the mesh will be printed onto the pyvista plot.

showSlice(index=0, ax=None, attr=None, axis='z', vmin=None, vmax=None)

Show slice of 3D mesh interactively.

Parameters:
indexint, optional

Index of the survey. Default is first survey (index == 0).

axmatplotlib.Axes, optional

Axis on which to plot the graph.

attrstr, optional

Attribute to plot. Default is ‘Resistivity(ohm.m)’.

axisstr, optional

Either ‘x’, ‘y’, or ‘z’ (default).

vminfloat, optional

Minimum value for colorbar.

vmaxfloat, optional

Maximum value for colorbar.

split3DGrid(elec=None, changeLabel=True)

Split self.elec to available lines based on ‘label’

Parameters:
elecdataframe, optional

Contains the electrodes information. “label” column must be provided and have “<line number> <electrode number>” format.

changeLablebool, optional

If True, the line number will be dropped from labels - Flase for GUI related surveys.

Returns:
elecListlist of dataframes

List of electrodes dataframes - each df can have a 3D like XYZ.

static topo2distance(x, y, z)

Convert 3d xy data in pure x lateral distance. Use for 2D data only!

write2in(param={}, err=None)

Create configuration file for inversion. Write mesh.dat and res0.dat.

Parameters:
paramdict

Dictionnary of parameters and values for the inversion settings.

err: bool

Used to overide self.err.

write2protocol(err=None, errTot=False, fm0=None, **kwargs)

Write a protocol.dat file for the inversion code.

Parameters:
errbool, optional

If True error columns will be written in protocol.dat provided an error model has been fitted or errors have been imported.

errTotbool, optional

If True, it will compute the modelling error due to the mesh and add it to the error from an error model.

fm0numpy.array of float, optional

Only for 3D time-lapse with reg_mode == 2 (difference inversion). Response of the inverted reference survey according to sequence of the reference survey as transfer resistance (Ohm).

**kwargsoptional

To be passed to Survey.write2protocol().

class resipy.Project.R2(*args, **kwargs)
class resipy.Project.cd(newPath)

Context manager for changing the current working directory

resipy.Project.cudaRm(invdir)

Compute Resolution and Covariance matrix for 2D problems using nVIDIA GPU.

Parameters:
invdirstring

Inversion directory used by R2.

Returns:
covarnd array

Values along the diagonal of the coviarance matrix.

rematnd array

Values along the diagonal of the Resolution matrix..

resipy.Project.getSysStat()

Return processor speed and usage, and free RAM and usage.

Returns:
cpu_speed: float

in Mhz.

cpu_usage: float

in percent.

ram_avail: float

avialable memory.

ram_usage: float

in percent.

ram_total: float

total memory in gb

resipy.Project.parallelRm(invdir)

Compute Resolution and Covariance matrix for 2D problems using multicore CPU. Behaves the same as cudaRm but uses numpy / openBlas.

Parameters:
invdirstring

Inversion directory used by R2.

Returns:
covarnd array

Values along the diagonal of the coviarance matrix.

rematnd array

Values along the diagonal of the Resolution matrix..

resipy.Project.systemCheck(dump=<built-in function print>)

Performs a simple diagnostic of the system, no input commands needed. System info is printed to screen, number of CPUs, memory and OS. This check is useful for parallel processing.

Parameters:
dumpfunction

stdout dump

Returns:
system_info: dict

Dictionary keys refer information about the system

Created on Fri Jun 1 11:21:54 2018

@author: ResIPy’s core developers

class resipy.Survey.Survey(fname=None, ftype='', df=None, elec=None, name='', spacing=None, parser=None, keepAll=True, debug=True, compRecip=True)

Class that handles geophysical data and some basic functions. One instance is created for each survey.

Parameters:
fnamestr

Name of the file where the data are.

ftypestr

Type of the data file. This setting is not read if a parser is not None.

dfpandas.DataFrame, optional

Pandas dataframe with ‘a’,’b’,’m’,’n’ columns as string and at least a column ‘resist’ as float. Can be provided alternatively to fname.

elecpandas.DataFrame, optional

Pandas dataframe with ‘x’,’y’,’z’,’buried’,’remote’ columns. ‘x’,’y’,’z’ as float, ‘buried’ and ‘remote’ as bool. Should be provided alternatively to fname.

namestr

A personal name for the survey.

spacingfloat, optional

This will be passed to the parser function to determine the electrode positions.

parserfunction, optional

If provided, it should return tuple containing elec a 3 columns array containing electrodes position and data a pandas.DataFrame with the a,b,m,n,resist columns at least and ip if present.

keepAll: bool, optional

If True will keep all the measurements even the ones without reciprocal. Note that if none of the quadrupoles have reciprocal they will all be kept anyway.

compRecip: bool, optional

Compute reciprocal errors, default is True.

addData(fname, ftype='Syscal', parser=None)

Add data to the actual survey (for instance the reciprocal if they are not in the same file).

addFilteredIP()

Add filtered IP data after IP filtering and pre-processing. This is because the IP filtering is done on a different dataframe and only merged when called this method.

addPerError(pnct=2.5)

Add a flat percentage error to resistivity data.

Parameters:
pnct: float

Error in percent

checkTxSign(inplace=True)

Check the sign of the measurement apparent resistivities. For example this is a necessary step for surveys using the syscal instrument as the polarity of measurements is not recorded.

Parameters:
inplace: bool

If True, Change the sign of transfer if they are negative. Defualt is True.

computeK()

Compute geometric factors, if any buried electrodes present then the generalised function for buried electrodes will be used, otherwise all electrodes are assumed to be on a flat surface.

Returns:
K: array like

Geometric factors for each measurement in Survey.df (dataframe)

computeKborehole(Gl=None)

Compute geometric factor for a borehole survey assuming a flat 2D surface. Gl = ground level. Calculated according to Eq 4.9 of Binley and Slater (2020).

computeKsurface()

Compute geomatrix factor (assuming flat surface) and store it in self.df[‘K’].

computeReciprocal(alg='Bisection Search')

Compute Reciprocals and store them in self.df (the dataframe)

Parameters:
algstr, optional

Algorithm used to compute reciprocals. Choose between ‘Bisection Search’, ‘Pandas Merge’ or ‘Array Expansion’. The default is ‘Bisection Search’, other string are casted to ‘Array Expansion’.

Notes

If the relevant cython code cannot be found then the function falls back to using the pandas merge approach.

computeReciprocalC()

Compute reciprocal measurements. using a bisection seach in cython! Parameters ———- forceSign: bool, optional

Force reciprocal and forward measurements to have the same polarity regarding the calculation of the reciprocal errors. Default is False.

Notes

The method first sorts the dipole AB and MN. Then efficiently searches for reciprocal pairs with a bisection search.

computeReciprocalN()

Compute reciprocal measurements using numpy array expansion.

Notes

The method first sorts the dipole AB and MN. Then creates a reciprocal quadrupole matrix. These matrices are then used with numpy.equal to produce a 2D matrix from which reciprocal are extracted.

computeReciprocalP()

Compute reciprocal measurements using pandas.merge().

Notes

This algorithm create two dataframes with ABMN and MNAB then use the pd.merge() function to match normal and their reciprocal. The resulting merged dataframe contains indexes of both ABMN and MNAB dataframes and is used to populate the irecip column.

convertLocalGrid()

Converts UTM grid to local grid for mesh stability

elec2distance()

Convert 3d xy data in pure x lateral distance. Use for 2D data only!

elec2horidist()

Convert 2D xz data into a true horizontal distance. Assumes that survey was done with a tape measure and the X distances are not true horizontal distance but distances measured along the ground.

estimateError(a_wgt=0.01, b_wgt=0.02)

Estimate reciprocal error data for data with no reciprocals, following the same routine present in R2. This allows for the additional inclusion of modelling errors.

Parameters:
a_wgt: float, optional

a_wgt documented in the R2 documentation

b_wgt: float, optional

b_wgt documented in the R2 documentation

exportSrv(fname=None)

Export .srv format for which is compatible with E4D. The e4d survey file includes the electrode locations, in addition to the scheduling matrix.

Parameters:
fname: string, optional

Where the output file will be written to. By default the file will take on the name of the survey and is written to the current working directory.

filterAppResist(vmin=None, vmax=None, debug=True)

Filter measurements by apparent resistivity for surface surveys

Parameters:
vminfloat, optional

Minimum value.

vmaxfloat, optional

Maximum value.

debugbool, optional

Print output to screen. Default is True.

filterContRes(vmin=None, vmax=None, debug=True)

Filter measurements by transfer resistance.

Parameters:
vminfloat, optional

Minimum value.

vmaxfloat, optional

Maximum value.

debugbool, optional

Print output to screen. Default is True.

filterDCA(dump=None)

Execute DCA filtering. Decay Curve Analysis (DCA) based on. Flores Orozco, A., Gallistl, J., Bücker, M., & Williams, K. H. (2017)., Decay curve analysis for data error quantification in time-domain induced polarization imaging., Geophysics, 83(2), 1–48. https://doi.org/10.1190/geo2016-0714.1

Parameters:
dumpfunction, optional

Callback function to print the progress in percent.

filterData(i2keep)

Filter out the data not retained in i2keep.

Parameters:
i2keepndarray of bool

Index where all measurement to be retained are True and the others False.

filterDefault(recompute_recip=True)

Remove NaN, Inf and duplicates values in the data frame.

filterDummy()

Remove measurements where abs(a-b) != abs(m-n) (likely to be dummy measurements added for speed).

filterElec(elec=[], debug=True)

Filter out specific electrodes given.

Parameters:
eleclist

List of electrode number to be removed.

debugbool, optional

Print output to screen. Default is True.

filterInvError(vmin=None, vmax=None, debug=False)

Filter measurements by inversion error.

Parameters:
vminfloat, optional

Minimum error.

vmaxfloat, optional

Maximum error.

debugbool, optional

Print output to screen. Default is False.

filterManual(attr='app', ax=None, log=False, label=None, vmin=None, vmax=None, elec=True, darkMode=False, flag3d=False)

Manually filters the data visually. The points manually selected are flagged in the Survey.iselect vector and can subsequently be removed by calling Survey.filterData(~Survey.iselect).

Parameters:
attrstr, optional

Columns of Survey.df to use for plotting. Default is app (apparent resistivity).

axmatplotlib axis, optional

If specified, the graph is plotted along the axis.

logbool, optional

If True then all data will be log transformed.

labelstr, optional

Label of the colorbar. If None, will be given from the ‘attr’ argument.

vminfloat, optional

Minimum value.

vmaxfloat, optional

Maximum value.

elecbool, optional

If True, the electrodes are shown and can be used for filtering.

darkModebool, optional

If true, electrodes wil be plotted in white, else black

flag3d: bool, optional

If true, function exits before

filterNegative()

Remove negative apparent resistivty values

filterNested()

Removes nested measurements: Where M or N are in between A and B

filterRangeIP(phimin, phimax)

Filter IP data according to a specified range.

Parameters:
phiminfloat

Minimium phase angle [mrad].

phimaxfloat

Maximum phase angle [mrad].

filterRecip(percent=20, debug=True)

Filter measurements based on the level reciprocal error.

Parameters:
percentfloat, optional

Percentage level of reciprocal error in which to filter the measurements. Percentage Errors > percentage will be removed. By default the value is 20.

debugbool, optional

Print output to screen. Default is True.

filterRecipIP()

Removing reciprocal measurements from dataset - only for visualization purposes on heatmap()

filterStack(percent=2, debug=True)

Filter measurements based on the stacking error.

Parameters:
percentfloat, optional

Percentage level of stacking error in which to filter the measurements. Percentage Errors > percentage will be removed. By default the value is 2.

debugbool, optional

Print output to screen. Default is True.

filterTransferRes(vmin=None, vmax=None, debug=True)

Filter measurements by transfer resistance.

Parameters:
vminfloat, optional

Minimum value.

vmaxfloat, optional

Maximum value.

debugbool, optional

Print output to screen. Default is True.

filterUnpaired()

Remove quadrupoles that don’t have a reciprocals. This might remove dummy measurements added for sequence optimization.

fitErrorLME(iplot=True, ax=None, rpath=None)

Fit a linear mixed effect (LME) model by having the electrodes as as grouping variables.

Parameters:
iplotbool, optional

If True, then a graph will be plotted.

axmatplotlib.Axes, optional

If specified, the graph will be plotted against this axis, otherwise a new figure will be created.

rpathstr, optional

Path of the directory with R (for Windows only).

fitErrorLin(ax=None)

Fit a linear relationship to the resistivity data.

Parameters:
axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

fitErrorParabolaIP(ax=None)

Plot the reciprocal phase errors with a parabola fit.

Parameters:
axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

fitErrorPwl(ax=None)

Fit an power law to the resistivity data.

Parameters:
axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

fitErrorPwlIP(ax=None)

Plot the reciprocal phase errors with a power-law fit.

Parameters:
axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

hasElecString()

Determine if a electrode strings are present in the electrode labels

Returns:
bool

True if strings present in electrode label

static logClasses3(datax, datay, func, class1=None)

Perform a log class of datay based on datay and applied function func to each bin.

Parameters:
dataxarray

x values to be from which the log-classes will be made.

datayarray

y values to which the function func will be applied inside each log-class

funcfunction

Function to be applied to the y value of each log class.

class1array, optional

Array of values for the log classes. If given, the limit of the classes will be computed as 10**class1

Returns:
mbinsarray

x-means of each bin.

vbinsarray

Value of each bin (output of the function func).

nbinsarray

Number of measurements in each bin.

normElecIdx(debug=True)

Normalise the electrode indexing sequencing to start at 1 and ascend consectively (ie 1 , 2 , 3 , 4 … )

Function firstly normalises all indexes so that the lowest electrode number is 1. Then removes jumps in the electrode indexing.

Parameters:
debugbool, optional

Output will be printed to console if True.

normElecIdxwSeq(expected=None)

Normalise the electrode indexing sequencing to start at 1 and ascend consectively (ie 1 , 2 , 3 , 4 … ). Also checks for any electrodes which are missing out of sequence if an expected sequence is given.

Function firstly normalises all indexes so that the lowest electrode number is 1. Then removes jumps in the electrode indexing.

Parameters:
expectedarray like

Expected sequence.

static numFormating(numList)

Formats numbers between -1 to 1 based on their decimals. e.g., 0.00001 would be 1e-5 while 0.1 would remain 0.1

Parameters:
numList: list of floats
setSeqIds()

Convert electrode labels to indexable integers, sets the class wide parameter ‘sequence’

showError(ax=None)

Plot the reciprocal errors.

Parameters:
axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

showErrorDist(ax=None)

Calculate and plots reciprocal error probablity histogram. Good data will have a bell shape (normal) distribution where most datapoints have near zero reciprocal error.

Parameters:
axMatplotlib.Axes

If specified, the graph will be plotted against it.

showErrorIP(ax=None)

Plot the reciprocal phase discrepancies.

Parameters:
axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

showHeatmap(ax=None)
Plot a phase heatmap (x = M, y = A and value = -phi) based on:

Orozco, A. F., K. H. Williams, and A. Kemna (2013), Time-lapse spectral induced polarization imaging of stimulated uranium bioremediation, Near Surf. Geophys., 11(5), 531–544, doi:10.3997/1873-0604.2013020)

Parameters:
axmatplotlib axis, optional

If specified, graph will be plotted on the given axis.

Returns:
figmatplotlib figure, optional

If ax is not specified, the function will return a figure object.

showPseudo(ax=None, bx=None, threed=False, **kwargs)

Plot pseudo section if 2D survey or just quadrupoles transfer resistance otherwise.

showPseudoIP(ax=None, bx=None, threed=False, **kwargs)

Plot pseudo section if 2D survey or just quadrupoles phase otherwise.

swapIndexes(old_indx, new_indx)

Replace the electrode number in a sequence matrix with another. Survey dataframe is updated after running.

Parameters:
old_idxint

Electrode number(or index) to be replaced

new_idxint

Replacement electrode number

write2protocol(outputname='', err=False, errTot=False, ip=False, res0=False, isubset=None, threed=False, fm0=None)

Write a protocol.dat file for R2, cR2, R3t, cR3t.

Parameters:
outputnamestr, optional

Path of the output file.

errbool, optional

If True, then the resError and phaseError (if IP present) will be used in the protocol.dat file.

errTotbool, optional

If True, the modelling error will be added to the error from the error model to form the total error.

ipbool, optional

If True and IP columns will be added to the file.

res0bool, optional

For time-lapse inversion, the background resistivity will be added if True.

isubsetarray_like of bool, optional

If specified, it will be used to take a subset of the original array. It can be used for difference inversion to take measurements in common between all the surveys.

threedbool, optional

If True, it’s for a 3D survey (and then add line numbers) or for a 2D survey (default).

fm0numpy.array of float, optional

Only for 3D timelapse, to compute d-d0+fm0 (if reg_mode == 2). If provided, will automatically set res0 to False.

Returns:
protocolpandas.DataFrame

Dataframe which contains the data for the protocol.dat.

resipy.Survey.fixSequence(sequence)

Make sequence consecutive.

Parameters:
sequencend array

N by 4 array which is the measurement sequence.

Returns:
None.
resipy.Survey.polyfit(x, y, deg=1)

Replacement function for numpy polyfit that works consistently (avoids SVD convergence error prsent on some W10 computers) Nb: is not as robust as numpy polyfit function.

Parameters:
xarray like

x values of data

yTYPE

y values of data

degTYPE, optional

DESCRIPTION. The default is 1.

Returns:
coefnp.array

coefficents of polynomail

Raises:
ValueError

If deg < 0 .

TypeError

array data is not as expected for fitting a polynomail.

Created on Wed May 30 10:19:09 2018, python 3.6.5 @author: jamyd91 Module handles mesh generation, display, discretisation and post processing. The convention for x y z coordinates is that the z coordinate is the elevation.

Dependencies:

numpy (standard scientific lib) matplotlib (standard scientific lib) gmshWrap(ResIPy resipy module) python3 standard libaries

class resipy.meshTools.Mesh(node_x, node_y, node_z, node_data, cell_type, original_file_path='N/A', order_nodes=True, compute_centre=True, check2D=True)

Mesh class.

Parameters:
node_xlist, 1d numpy array

x coordinates of nodes

node_ylist, 1d numpy array

coordinates of nodes

node_zlist, 1d numpy array

z coordinates of nodes

connection_matrix: M by N numpy array

nodes of element vertices in the form [[node1],[node2],[node3],…], each node id should be an integer type.

cell_typelist of ints

code referencing cell geometry (e.g. triangle) according to vtk format

original_file_pathstring, optional

file path to where the mesh file was originally imported

regionsoptional

element indexes for a material in the mesh (needs further explanation)

Returns:
Meshclass
addAttribute(values, key)

Add a new/update attribute to mesh.

Parameters:
values: array like

must have a length which matches the number of elements. Discrete values which map to the elements.

key: str

Name of the attribute, this will be used to reference to the values in other mesh functions.

addPtAttribute(values, key)

Associate attributes with the mesh nodes

Parameters:
values: array like

Must have a length which matches the number of nodes. Discrete values which map to the nodes.

key: str

Name of the attribute, this will be used to reference to the values in other mesh functions.

applyFunc(mesh_paras, material_no, new_key, function, *args)

Applies a function to a mesh by zone number and mesh parameter.

Parameters:
mesh_parasarray like

Mesh parameters from which new parameters are calculated.

material_noarray like of ints

Material type assigned to each element, should be numbered consectively from 1 to n. in the form 1 : 1 : 2 : n. …ie if you have 2 materials in the mesh then pass an array of ones and twos. zeros will be ignored.

new_keystring

Key assigned to the parameter in the df. DOES NOT default.

functionfunction

Function to be applied to mesh attributes, first argument must be the mesh parameter.

args[see function info]

All arguments to be passed through function after to modify the mesh parameter, … argument must be in the form of [(argA1,argB1),(argA2,argB2)], … where letters are the material, numbers refer to the argument number

Notes

Mesh object will now have the new attribute added once the function is run. Use the mesh.show() (or .draw()) function to see the result.

assignZone(poly_data)

Assign material/region assocations with certain elements in the mesh say if you have an area you’d like to forward model. *2D ONLY*

Parameters:
poly_datalist

list of 2 by N (x,y) numpy arrays describing each polygon.

Returns:
material_nonumpy.array

Element associations starting at 1. So 1 for the first region defined in the region_data variable, 2 for the second region defined and so on. If the element can’t be assigned to a region then it’ll be left at 0.

assignZoneAttribute(attr_list, new_key, zone=None)

Asssigns values to the mesh which depend on region / material only. E.G a single resistivity value.

Parameters:
attr_listlist

Values corresponding to a material number in the mesh. eg. if you had 3 regions in the mesh then you give [resistivity1,resistivity2,resistivity3].

new_keystring

Key identifier assigned to the attribute in the df.

zonearray or list

Integers starting at 0 or 1, and ascend in intervals of 1, which correspond to a material in the mesh returned from assign_attr_ID. Should have the same length as the number of elements in the mesh.

Notes

Mesh object will now have the new attribute added once the function is run. Use the mesh.show() (or .draw()) function to see the result.

cellArea()

Compute the element areas, or in the case of 3D meshes compute the cell volumes.

cellCentres()

A numpy-based approximation of cell centres for 2D and 3D elements. It’s calculated from the mean of cell x y z node coordinates

computeElmDepth()

Compute the depth of elements relative to the surface of the mesh.

computeNeigh()

Compute element neighbour matrix

copy()

Return a copy of mesh object.

crop(polyline)

Crop the mesh given a polyline in 2D. If 3D then the mesh will be cropped in the XY plane.

Parameters:
polylinearray of float

Array of size Nx2 with the XZ (or XY) coordinates forming the polyline. Notethat the first and last coordinates should be the same to close the polyline.

dat(file_path='mesh.dat')

Write a mesh.dat kind of file for mesh input for R2/R3t.

Parameters:
file_pathstr, optional

Path to the file. By default ‘mesh.dat’ is saved in the working directory.

datAdv(file_path='mesh.dat', iadvanced=True)

Write a mesh.dat kind of file for mesh input for R2/R3t. Advanced format which includes the neighbourhood and conductance matrix.

Parameters:
file_pathstr, optional

Path to the file. By default ‘mesh.dat’ is saved in the working directory.

iadvancedbool, optional

If True, use the advanced mesh format if False, use the normal mesh format.

downslopeID(attr='Resistivity')

Get the index of ‘downslope’ elements for a given attribute.

Parameters:
attrstring, optional

Key inside of mesh dataframe. The default is ‘Resistivity’.

Returns:
idxnd array

Indices of downslope elements.

draw(attr=None, edge_color='k', color_map='Spectral', color_bar=False, vmin=None, vmax=None)

Redraws a mesh over a prexisting figure canvas, this is intended for saving time when plotting each iteration of the resistivity inversion.

Parameters:
color_mapstring, optional

color map reference

color_barboolean, optional

True to plot colorbar

axmatplotlib axis handle, optional

Axis handle if preexisting (error will thrown up if not) figure is to be cast to.

edge_colorstring, optional

Color of the cell edges, set to None if you dont want an edge.

vminfloat, optional

Minimum limit for the color bar scale.

vmaxfloat, optional

Maximum limit for the color bar scale.

attrstring, optional

Which attribute in the mesh to plot, references a dictionary of attributes. attr is passed as the key for this dictionary.

Returns:
figurematplotlib figure

Figure handle for the plotted mesh object.

elemDist()

Work out the distance of each cell in the mesh to its closest electrode

exportMesh(fname, ftype=None, coordLocal=False, coordParam=None, voxel=False, meshParam=None)

Export mesh into with coordinate rotation data

Parameters:
fnamestr, optional

Path to output save file.

ftype: str, optional

File extension, if none, will be guessed from file name.

coordLocal: bool, optional

If True, convert mesh coordinates to a national grid system or utm.

coordParam: dict, optional

Coordinate conversion parameters, x0, y0 and a. Stored as a dictionary.

exportTetgenMesh(prefix='mesh', zone=None, debug=False, mixed=False)

Export a mesh like the tetgen format for input into E4D. This format is composed of several files. Currently only tested for 3D surface array like surveys. Assumes the sides of the mesh are parrallel to the x and y plane and that the surface has gentle to no topography.

Parameters:
prefix: string

Prefix assigned to exported files. Can include path to file.

zone: array like, optional

If using zones in the mesh, this attribute should include the zone attribute which is an array of integers identifying the zone associated with each element. By default each element is assigned to zone 1.

mixed: bool

If true attempt to compute mixed boundary condition at mesh edges.

Notes
———-
Please note routine is experimental and not garanteed to work 100%
externalNodes()

Find the external nodes of the mesh (triangle and tetrahedral meshes only).

Returns:
external_nodesndarray

Indices of external nodes.

surface_flagndarray

Flag if external node is on the surface of the mesh, 1 for surface and 0 for not.

extractSurface(return_idx=False, post_neigh_check=True)

Extract the surface of a triangle or tetrahedral mesh. Ouput of function will depend on mesh type.

Parameters:
return_idx: bool

Return the indices of elements on the top of the mesh (default is false)

post_neigh_check: bool

Perform a post processing step to check elements have neighbours, this is useful for meshes with uneven sides where false postives can crop up.

Returns:
mesh: class

2d faces of the top of the mesh, if the input mesh is a tetrahedral mesh

(x,z): tuple

1D faces of the top of the mesh, if the input is a triangular mesh

idx: array like

if return_idx == True an array of int

file_path()

Returns the file path from where the mesh was imported

filterIdx(in_elem)

Filter mesh down on the basis of element number / index

Parameters:
in_elem: array like

array of bool, True means to keep the element

attr: string

Name of attribute to threshold by

Returns:
mesh: Class

New instance of mesh class which is filtered

findIdirichlet()

Find the best node for the dirichlet node.

Returns:
idirchletint

A node far away as possible from each of the electrodes on the boundary of the mesh. (Add one if using inside of mesh(3d).dat)

static findParaview()

Run on windows to find paraview.exe command.

Returns:
found: bool

If True the program was able to find where paraview is installed If false then the program could not find paraview in the default install locations.

location: str

if found == True. The string maps to the paraview executable if found == False. then ‘n/a’ is returned.

flipYZ()

Make Y column Z column and vice versa, this is useful for when dealing with 2D meshes. Dimensions are modified in place.

meshLookUp(look_up_mesh)

Look up values from another mesh using nearest neighbour look up, assign attributes to the current mesh class.

Parameters:
look_up_mesh: class

Another mesh class.

moveElecNodes(new_x, new_y, new_z, debug=True)

Move the electrodes to different nodes which are close to the given coordinates. This is useful for timelapse surveys where the electrodes move through time, ideally this is implemented on a mesh which is refined near the surface. If no nodes are assigned to mesh, a mesh.e_nodes variable is created.

Parameters:
new_x: array l ike

new electrode x coordinates

new_y: array like

new electrode y coordinates, if 2D mesh let new_y=None and the array will automatically be assigned an array of zeros.

new_z: array-like

new electrode z coordinates

debug: bool, optional

Controls if any changes to electrode nodes will be output to console.

Returns:
node_in_mesh: np array

Array of mesh node numbers corresponding to the electrode postions/ coordinates.

Notes

Mesh.e_nodes parameter is updated (or added) after running function.

node2ElemAttr(node_array, name='NodeArray')

Maps node attributes onto mesh elements

Parameters:
node_arrayarray like

Array with same number of entries as there are nodes in the mesh.

namestr, optional

Name of node attribute, will be used to index array in mesh dataframe. The default is ‘NodeArray’.

Returns:
arrnd array

Node array mapped on mesh elements .

orderElem(param=None)

Order elements based on the parameter number. Ideally parameters should be concurrent, and fixed elements should go at the base of the connection matrix.

orderNodes(return_count=False)

Order mesh nodes in clockwise fashion

Parameters:
return_count:bool, optional

Function returns the number of reordered elements, default is False.

paraview(fname='ResIPy_mesh.vtk', loc=None)

Show mesh in paraview

Parameters:
fname: str,optional

A vtk file will be written to working directory and then displayed using paraview.

loc: str, optional

Path to paraview excutable, ignored if not using windows. If not provided the program will attempt to find where paraview is installed automatically.

pick3Dbox(ax=None, xlim=None, ylim=None, zlim=None, electrodes=True, darkMode=False)

Pick out regions on the mesh using the box widget (intended for 3D forward modelling)

Parameters:
axclass, optional

pyvista plotting object. The default is None.

xlimtuple, list, array like, optional

x limit of the selection area. The default is None and assigned automatically.

ylimtuple, list, array like, optional

as with xlim. The default is None.

zlimtuple, list, array like, optional

as with xlim. The default is None.

electrodesbool, optional

If True the electrodes are also plotted on the mesh. The default is True.

darkmode: bool, optional

Alters coloring of pyvista plot for a darker appearance

Returns:
pyvista clip box handle
Raises:
Exception

If mesh is 2D

quad2tri()

Make a triangle mesh from a quad mesh

Returns:
mesh: class

Triangle mesh class

quadMeshNp(topo=None)

Convert mesh nodes into x column indexes in the case of quad meshes. Does not currently support changes in electrode elevation!

Returns:
colx: list

X column indexes for quad mesh

refine()

Refine the mesh into smaller elements

resetParam()

Reorder parameters into consective ascending order

saveMesh(fname, ftype=None)

Save mesh into a file. Available formats are .dat, .vtk, .csv, .xyz, and .node. .xyz is a special case where the mesh should be transformed into a voxel mesh first.

Parameters:
fnamestr, optional

Path to output save file.

ftype: str, optional

File extension, if none, will be guessed from file name.

setElecNode(e_nodes, iremote=None)

Assign node numbers to electrodes.

Parameters:
e_nodes: array like

array of ints which index the electrode nodes in a mesh

iremote: array like, optional

Array of bool, if True then indexed electrode is remote.

show(color_map='Spectral', color_bar=True, xlim=None, zlim=None, ax=None, electrodes=True, sens=False, edge_color='k', contour=False, vmin=None, vmax=None, attr=None, clabel=None, hor_cbar=False, sensPrc=None, maxDepth=None, aspect='equal', darkMode=False, **kwargs)

Displays a 2d mesh and attribute.

Parameters:
color_mapstring, optional

color map reference

color_barBoolean, optional

True to plot colorbar

xlimtuple, optional

Axis x limits as (xmin, xmax).

zlimtuple, optional

Axis z limits as (zmin, zmax).

axmatplotlib axis handle, optional

Axis handle if preexisting (error will thrown up if not) figure is to be cast to.

electrodesboolean, optional

Enter true to add electrodes to plot.

sensboolean, optional

Enter true to plot sensitivities.

edge_colorstring, optional

Color of the cell edges, set to None if you dont want an edge.

contourboolean, optional

If True, plot filled with contours instead of the mesh.

vminfloat, optional

Minimum limit for the color bar scale.

vmaxfloat, optional

Maximum limit for the color bar scale.

attrstring, optional

Which attribute in the mesh to plot, references a dictionary of attributes. attr is passed as the key for this dictionary.

clabelstring, optional

Label of the colorbar. Default is the value of attr argument.

hor_cbarboolean, optional

‘True’ to make a horizontal color bar at the bottom of the plot, default is vertical color bar to the right of the plot.

sensPrcfloat, optional

Normalised (between 0 and 1) sensitivity value threshold. Default is None meaning the sensitivity is just overlay. Need sens=True to be used.

maxDepthfloat

Maximum absolute depth to be shown on the plotted figure.

aspectstring, optional

defines the aspect ratio of the plot. ‘equal’ locks the aspect ratio. ‘auto’, aspect ratio is define by plotting area.

darkModebool, optional

If True, electrodes will be plotted in white, else black

Returns:
figurematplotlib figure

Figure handle for the plotted mesh object.

Notes

Show a mesh object using matplotlib. The color map variable should be a string refering to the color map you want (default is “jet”). As we’re using the matplotlib package here any color map avialable within matplotlib package can be used to display the mesh here also. See: https://matplotlib.org/2.0.2/examples/color/colormaps_reference.html

show3D(color_map='Spectral', color_bar=True, xlim=None, ylim=None, zlim=None, ax=None, electrodes=True, sens=False, edge_color='k', alpha=1, vmax=None, vmin=None, attr=None, elec_color='k', elec_size=10.0, use_pyvista=True, background_color=(0.8, 0.8, 0.8), pvslices=([], [], []), pvspline=None, pvthreshold=None, pvgrid=True, pvcontour=[], pseudo3DContour=False, pvdelaunay3d=False, pvshow=True, volume=None, darkMode=False, cell_picking=False, clipping=True)

Shows a 3D tetrahedral mesh.

Parameters:
color_mapstring, optional

Matplotlib color map reference

color_barBoolean, optional

True to plot colorbar

xlimtuple, optional

Axis x limits as (xmin, xmax).

ylimtuple, optional

Axis y limits as (ymin, ymax).

zlimtuple, optional

Axis z limits as (zmin, zmax).

axmatplotlib axis handle, pvista plotter handle, optional

Axis handle if preexisting (error will thrown up if not) figure is to be cast to. If using pyvista then then ax is the plotter the object.

electrodesboolean, optional

Enter true to add electrodes to plot (if available in mesh class)

sensboolean, optional

Enter true to plot sensitivities (if available in mesh class). Note that for 3D this doesnt work so well.

edge_colorstring, optional

Color of the cell edges, set to None if you dont want an edge.

alpha: float, optional

Should be set between 0 and 1. Sets a transparancy value for the element faces.

vminfloat, optional

Minimum limit for the color bar scale.

vmaxfloat, optional

Maximum limit for the color bar scale.

attrstring, optional

Which attribute in the mesh to plot, references a dictionary of attributes. attr is passed as the key for this dictionary.

elec_colorstring, optional

Colour of the electrodes on the plot if electrodes = True. Default is ‘k’ for black. Can be a 3 by 1 tuple or string identifier.

elec_sizefloat, optional

Size of electrode points, default is 10, a number of 4-5 gives better visibility

use_pyvistabool, optional

Use visual toolkit backend for displaying 3D mesh, note that pyvista must be installed for this to work.

background_colortuple, optional

Background color assigned to pyvista plotter object when created. Not yet supported for matplotlib axis handles.

pvslicestuple 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.

pvspline‘elec’ or numpy array, optional

If ‘elec’ mesh will be sliced along the electrodes path otherwise an array of X, Y, Z of points on a path to slice the mesh along that path is needed.

pvthresholdlist of two floats, optional

Keep values between pvthreshold[0] and pvthreshold[1].

pvgridbool, optional

Show grid or not.

pvcontourlist of float, optional

Values of the isosurface to be plotted.

pseudo3DContour: bool, optional

If ‘True’, pseudo 3D plots will be contoured. Only use in case of pseudo 3D surveys.

pvdelaunay3dbool, optional

If True a “Delaunay 3D” triangulation filter will be applied on the mesh.

pvshowbool, optional

If False, that will prevent calling the pyvista.Plotter.show(). This is useful in case of subplots.

volumefloat, optional

If not ‘None’ then volume float number will be printed onto the pyvista plot.

darkmode: bool, optional

Alters coloring of pyvista plot for a darker appearance

cell_picking: bool, optional

Interactive picking flag, for the use within the UI only. Leave as False.

clipping: bool, optional

Flag to clip mesh in pyvista window, if tank type problem is expected set to False else set to True. Default is True. Note if False then the X Y and Z limit extents will be ignored.

Returns:
figurematplotlib figure, pyvista plotter

Figure handle for the plotted mesh object.

Notes

Show a mesh object using matplotlib. The color map variable should be a string refering to the color map you want (default is “Spectral”). As we’re using the matplotlib package here any color map avialable within matplotlib package can be used to display the mesh here also. See: https://matplotlib.org/2.0.2/examples/color/colormaps_reference.html

Plotting sensitivies using sens=True is not reccomended. The matplotlib renderer has trouble with it.

showAvailAttr(flag=True)

Show available attributes in mesh.df.

showPrismMesh(color_map='Spectral', color_bar=True, ax=None, electrodes=True, sens=False, edge_color='k', alpha=1, aspect='equal', vmax=None, vmin=None, attr=None, xlim=None, ylim=None, zlim=None)

Shows a 3D prism mesh.

Parameters:
color_mapstring, optional

color map reference

color_barBoolean, optional

True to plot colorbar

xlimtuple, optional

Axis x limits as (xmin, xmax).

ylimtuple, optional

Axis y limits as (ymin, ymax).

zlimtuple, optional

Axis z limits as (ymin, ymax).

axmatplotlib axis handle, optional

Axis handle if preexisting (error will thrown up if not) figure is to be cast to.

electrodesboolean, optional

Enter true to add electrodes to plot (if available in mesh class)

sensboolean, optional

Enter true to plot sensitivities (if available in mesh class). Note that for 3D this doesnt work so well.

edge_colorstring, optional

Color of the cell edges, set to None if you dont want an edge.

alpha: float, optional

Should be set between 0 and 1. Sets a transparancy value for the element faces.

vminfloat, optional

Minimum limit for the color bar scale.

vmaxfloat, optional

Maximum limit for the color bar scale.

attrstring, optional

Which attribute in the mesh to plot, references a dictionary of attributes. attr is passed as the key for this dictionary.

aspect = string, optional

defines the aspect ratio of the plot. ‘equal’ locks the aspect ratio. ‘auto’, aspect ratio is define by plotting area.

Returns:
figurematplotlib figure

Figure handle for the plotted mesh object.

showSlice(attr='Resistivity(log10)', axis='z', vmin=None, vmax=None, ax=None)

Show 3D mesh slice.

Parameters:
attrstr, optional

String the of the variable to plot. Default is Resistivity(log10).

axisstr, optional

Axis to be sliced either, x, y or z.

vminfloat, optional

Minimum value for colorbar.

vmaxfloat, optional

Maximum value for colorbar.

axmatplotlib.Axis

If provided, plot will be drawn on this axis.

splitTetra(param=None)

Refine tetrahedra by splitting them in six smaller tetrahedra.

splitTri(param=None)

Refine triangles by splitting them into 4 smaller triangles

summary(flag=True)

Prints summary information about the mesh

threshold(attr=None, vmin=None, vmax=None)

Threshold the mesh to certian attribute values.

Parameters:
attr: string

Name of attribute to threshold by

vmin: float

minimum value of attribute

vmax: float

maximum value of attribute

Returns:
mesh: Class

New instance of mesh class which is thresholded

toCsv(file_name='mesh.csv')

Write a .csv file of the mesh, the first 3 columns are the element centres at coordinates x,y,z and the rest of the columns are the attributes in the df

Parameters:
file_nameString, optional

The default is ‘mesh.csv’.

toVoxelMesh(elec=None, elec_type=None, **kwargs)

Returns a special kind of mesh optimised for visualasation inside of GIS software like Geovisionary.

Parameters:
elecarray like, pd.DataFrame, optional

Electrode coordinates. The default is None.

elec_typelist, optional

List of electrode types, i.e. ‘surface’ or ‘buried’. The default is None (all electrodes are surface electrodes).

Returns:
meshclass

Mesh object with voxel elements

Raises:
Exception

If no electrode coordinates provided.

toVoxet(file_name='mesh.vo', x0=None, y0=None, a=0)

Not yet working!

Parameters:
file_nameTYPE, optional

DESCRIPTION. The default is ‘mesh.vo’.

x0TYPE, optional

DESCRIPTION. The default is None.

y0TYPE, optional

DESCRIPTION. The default is None.

aTYPE, optional

DESCRIPTION. The default is 0.

Returns:
None.
transMesh(x, y, z)

Translate mesh by x y z coordinate

truncateMesh(xlim=None, ylim=None, zlim=None)

Crop the mesh to a box of given limits, like how R3t behaves when outputting inverted results.

Parameters:
xlimtuple, optional

Axis x limits as (xmin, xmax).

ylimtuple, optional

Axis y limits as (ymin, ymax).

zlimtuple, optional

Axis z limits as (ymin, ymax).

Returns:
nmesh: Class

New instance of mesh class which is truncated

type2FaceNo()

Converts vtk cell types into number of faces each element has

type2VertsNo()

Converts vtk cell types into number of vertices each element has

vtk(file_path='mesh.vtk', title=None, replace_nan=-9999)

Writes a vtk file for the mesh object, everything in the df will be written to file as attributes. We suggest using Paraview to display the mesh outside of ResIPy. It’s fast and open source :).

Parameters:
file_path: str, optional

Maps where python will write the file, if left as default then mesh.vtk will be written the current working directory.

titlestr, optional

Header string written at the top of the vtk file .

writeAttr(attr_key=None, file_name='_res.dat')

Writes a attribute to a _res.dat type file. file_name entered seperately because it will be needed for the R2 config file. The reason for this function is so you can write a forward model parameter input file.

Parameters:
attr_key: string

Key identifying the attr to be written in the mesh object df.

file_name: string, optional

Name of the _res.dat type file.

writeRindex(fname)

Write out the neighbourhood matrix for R3t.

Parameters:
fnamestr

Path to written file.

xyz(file_name='mesh.xyz', coordParam=None)

Write a .xyz file of the mesh, the first 3 columns are the element centres at coordinates x,y,z and the rest of the columns are the attributes in the df

Parameters:
file_nameString, optional

The default is ‘mesh.csv’.

resipy.meshTools.check4repeatNodes(X, Y, Z, flag=None)

Raise error if repeated nodes present

Parameters:
Xarray like

X coordinates of nodes (normally electrodes).

Yarray like

Y coordinates of nodes (normally electrodes).

Zarray like

Z coordinates of nodes (normally electrodes).

flaglist of string, optional

String flag assigned to each node. The default is None.

Returns:
None.
Raises:
ValueError

If repeated nodes detected

resipy.meshTools.cylinderMesh(elec_x, elec_y, elec_z, zlim=None, radius=None, file_path='cylinder_mesh.geo', cl=-1, cl_factor=2, finer=4, keep_files=True, show_output=True, path='exe', dump=<built-in function print>, handle=None)

Make a cylinder mesh.

Parameters:
eleclist of array_like or nx3 array

First column/list is the x coordinates of electrode positions and so on …

zlimlist, tuple, optional

Bottom and top z coordinate of column, in the form (min(z),max(z))

radius: float, optional

Radius of column. If not provided, will be infered from elec positions.

file_pathstring, optional

Name of the generated gmsh file (can include file path also).

clfloat, optional

Characteristic length, essentially describes how big the nodes associated elements will be. Usually no bigger than 5. If set as -1 (default) a characteristic length 1/4 the minimum electrode spacing is computed.

finerint, optional

Number of line between two consecutive electrodes to approximate the circle shape. No longer used.

handlevariable, optional

Will be assigned the output of ‘Popen’ in case the process needs to be killed in the UI for instance.

resipy.meshTools.dat_import(file_path='mesh.dat', order_nodes=True)

Import R2/cR2/R3t/cR3t .dat kind of mesh.

Parameters:
file_path: str

Maps to the mesh (.dat) file.

Returns:
mesh: class
resipy.meshTools.halfspaceControlPoints(elec_x, elec_y, r, min_r=None)

Create control points for surface electrodes in 3D (or 2d borehole electrodes)

resipy.meshTools.in_box(x, y, z, xmax, xmin, ymax, ymin, zmax, zmin)

Determine if a point lies inside a bounding volume.

Parameters:
xarray like, float

x coordinate of query point

yarray like, float

y coordinate of query point

zarray like, float

z coordinate of query point

volume_datalist

contains column of poly_data, in the form (polyx, polyy, polyz)

ray_castfloat, optional

determines how the far in the x axis a point is ray casted

Returns:
insideboolian, numpy array

true indexes where point is inside volume

resipy.meshTools.mergeMeshes(meshList)

Merge multiple meshes into one mesh object, useful for psuedo 3d surveys.

Parameters:
meshListlist

List of seperate mesh classes .

Returns:
mmesh: class

Merged mesh.

resipy.meshTools.moveMesh2D(meshObject, elecLocal, elecGrid)
Move mesh object to a certain place on the grid.

X, Y only. No Z relocation.

Parameters:
meshObjectclass object

Mesh class object.

elecLocal: dataframe

dataframe of electrodes on local 2D grid (i.e., y = 0) ‘remote’ column (bool) required.

elecGrid: dataframe

dataframe of electrodes on 3D grid (i.e., y != 0) ‘remote’ column (bool) required

Returns:
meshclass object

Copy of moved Mesh class object.

resipy.meshTools.points2vtk(x, y, z, file_name='points.vtk', title='points', data=None)

Function makes a .vtk file for some xyz coordinates. optional argument renames the name of the file (needs file path also) (default is “points.vtk”). title is the name of the vtk file.

Parameters:
xlist, tuple, np array

X coordinates of points.

ylist, tuple, np array

Y coordinates of points.

zlist, tuple, np array

Z coordinates of points.

file_namestring, optional

Path to saved file, defualts to ‘points.vtk’ in current working directory.

titlestring, optional

Title of vtk file.

Returns:
~.vtkfile
resipy.meshTools.prismMesh(elec_x, elec_y, elec_z, keep_files=True, show_output=True, path='exe', dump=<built-in function print>, handle=None, **kwargs)

Make a prism mesh.

Parameters:
elec_x: array like

electrode x coordinates

elec_y: array like

electrode y coordinates

elec_z: array like

electrode z coordinates

poly: list, tuple, optional

Describes polygon where the argument is 2 by 1 tuple/list. Each entry is the polygon x and y coordinates, ie (poly_x, poly_y)

z_lim: list, tuple, optional

top and bottom z coordinate of column, in the form (min(z),max(z))

radius: float, optional

radius of column

file_path: string, optional

name of the generated gmsh file (can include file path also) (optional)

cl: float, optional

characteristic length (optional), essentially describes how big the nodes assocaited elements will be. Usually no bigger than 5. If set as -1 (default) a characteristic length 1/4 the minimum electrode spacing is computed.

elemz: int, optional

Number of layers in between each electrode inside the column mesh.

handlevariable, optional

Will be assigned the output of ‘Popen’ in case the process needs to be killed in the UI for instance.

resipy.meshTools.quadMesh(elec_x, elec_z, elec_type=None, elemx=4, xgf=1.5, zf=1.1, zgf=1.5, fmd=None, pad=2, surface_x=None, surface_z=None, refine_x=None, refine_z=None, model_err=False)

Creates a quaderlateral mesh given the electrode x and y positions.

Parameters:
elec_xlist, np array

Electrode x coordinates

elec_zlist, np array

Electrode y coordinates

elec_type: list, optional

strings, where ‘electrode’ is a surface electrode; ‘buried’ is a buried electrode

elemxint, optional

Number of nodes between electrodes.

xgffloat, optional

X factor multiplier for fine zone.

zffloat, optional

Z factor multiplier in the fine zone (must be >1).

zgffloat, optional

Z factor multiplier in the coarse zone (must be >1).

fmdfloat (m), optional

Fine mesh region depth specifies as positive number (if None, half survey width is used).

padint, optional

X padding outside the fine area (tipicaly twice the number of elements between electrodes).

surface_x: array like, optional

Default is None. x coordinates of extra surface topography points for the generation of topography in the quad mesh

surface_z: array like, optional

Default is None. z coordinates of extra surface topography points for the generation of topography in the quad mesh. Note an error will be returned if len(surface_x) != len(surface_z)

refine_x: array like, optional

Default is None. Inserts points in the resulting quad mesh for more control over mesh refinement at depth. X coordinates. An error will be returned if len(refine_x) != len(refine_z).

refine_z: array like, optional

Default is None. Inserts points in the resulting quad mesh for more control over mesh refinement at depth. Z coordinates. An error will be returned if len(refine_x) != len(refine_z).

model_errbool, optional

If True all topography will be normalised to zero such that the returned mesh can be used to model foward modelling errors.

Returns:
Meshclass

Mesh object

meshxnumpy.array

Mesh x locations for R2in file.

meshznumpy.array

Mesh locations for R2in file (ie node depths).

toponumpy.array

Topography for R2in file.

elec_nodenumpy.array

x columns where the electrodes are.

resipy.meshTools.readMesh(file_path, node_pos=None, order_nodes=True)

Import user defined mesh, currently supports .msh, .vtk and .dat (native to R2/3t) format for quad, triangular and tetrahedral meshes. The type of file is guessed from the extension given to the code.

Parameters:
file_pathstring

Path to file.

node_posarray like, optional

Array of ints referencing the electrode nodes. If left as none no electrodes will be added to the mesh class. Consider using mesh.move_elec_nodes() to add nodes to mesh using their xyz coordinates.

order_nodesbool, optional

Order nodes when importing a mesh

Returns:
meshclass

mesh class used in ResIPy

resipy.meshTools.runGmsh(ewd, file_name, show_output=True, dump=<built-in function print>, threed=False, handle=None)
Parameters:
ewdstr

Directory where gmsh copy is stored.

file_namestr

Name of the .geo file without extension.

show_outputTYPE, optional

If True, output of gmsh is displayed to dump. The default is True.

threedbool, optional

If True, 3D mesh is done, else 2D. The default is False.

handlevariable, optional

Will be assigned the output of ‘Popen’ in case the process needs to be killed in the UI for instance.

Returns:
None.
resipy.meshTools.tankMesh(elec=None, origin=None, dimension=[10, 10, 10], file_path='tank_mesh.geo', cl=-1, keep_files=True, show_output=True, path='exe', dump=<built-in function print>, handle=None)

Make a tank mesh (3D closed box).

Parameters:
eleclist of array_like or nx3 array

First column/list is the x coordinates of electrode positions and so on …

originlist of float, optional

Origin of the corner where the mesh will be drawned from. If not provided and elec provided, the smaller elec position will be chosen.

dimensionlist of float, optional

Dimension of the mesh in X,Y,Z from the corner origin. The default is [10,10,10].

file_pathstring, optional

Name of the generated gmsh file (can include file path also).

clfloat, optional

Characteristic length, essentially describes how big the nodes associated elements will be. Usually no bigger than 5. If set as -1 (default) a characteristic length 1/4 the minimum electrode spacing is computed.

finerint, optional

Number of line between two consecutive electrodes to approximate the circle shape.

handlevariable, optional

Will be assigned the output of ‘Popen’ in case the process needs to be killed in the UI for instance.

resipy.meshTools.tetgen_import(file_path, order_nodes=True)

Import Tetgen mesh into ResIPy. This isa little different from other imports as the mesh is described by several files. From meshTools’ perspective What is needed is the node(.node) file and element (.ele) files, which describe the coordinates of mesh nodes and the connection matrix.

Parameters:
file_path: str

Maps to the mesh .node file. The program will automatically find the corresponding .ele file in the same directory.

order_nodes: bool, optional

Check ordering of the nodes on import, this ensures that the mesh will work inside of R3t / cR3t. Default is True.

Returns:
mesh: class
resipy.meshTools.tetraMesh(elec_x, elec_y, elec_z=None, elec_type=None, keep_files=True, interp_method='triangulate', surface_refinement=None, mesh_refinement=None, ball_refinement=True, path='exe', dump=<built-in function print>, whole_space=False, model_err=False, handle=None, show_output=True, **kwargs)

Generates a tetrahedral mesh for R3t (with topography) for field surveys. This function expects the current working directory has path: exe/gmsh.exe. Uses post processing after mesh generation to super impose topography on to a flat 3D tetrahedral mesh of a hemisphere. Handles, wholespace, 2d and halfspace problems.

Parameters:
elec_xarray like

electrode x coordinates

elec_yarray like

electrode y coordinates

elec_zarray like

electrode z coordinates

elec_typelist of strings, optional

Defines if electrodes are buried or not.

keep_filesboolean, optional

True if the gmsh input and output file is to be stored in the working directory.

interp_methodstring, default =’triangulate’ optional

Interpolation method to translate mesh nodes in the z direction. In other words the method in which topography is appended to the mesh. Here the topography is added to the mesh in post processing. The options are documented below in the notes.

surface_refinementnp.array, optional

Numpy array of shape (3,n), should follow the format np.array([x1,x2,x3,…],[y1,y2,y3,…],[z1,z2,z3,…]). Allows for extra refinement for the top surface of the mesh. The points are not added to the mesh, but considered in post processing in order to super impose topography on the mesh.

mesh_refinementdict, pd.DataFrame, optional

Dataframe (or dict) contianing ‘x’, ‘y’, ‘z’, ‘type’, and ‘cl’, columns which describe control points which can be used to refine the mesh, unlike surface_refinement, this argument allows the user granular control over the refinement of mesh elements. If not provided ResIPy attempts to add mesh_refinement for you (though not in the case of a wholespace). See further explanation in tetraMesh notes.

ball_refinementboolean, optional

If True, tells gmsh to add a ‘ball’ of refined mesh around electrodes. Default is True.

pathstring, optional

Path to exe folder (leave default unless you know what you are doing).

whole_spaceboolean, optional

flag for if the problem should be treated as a whole space porblem, in which case electrode type is ignored and all electrodes are buried in the middle of a large mesh.

model_errbool

If True, a flat mesh will be returned for the sake of estimating forward modelling errors.

dumpfunction, optional

Function to which pass the output during mesh generation. print() is the default.

**kwargsdict

Keyword arguments to be passed to functions in gmshWrap.py

Returns:
mesh3d: class

Notes

Possible arguments for interp_method:
‘bilinear’4 known points are used to compute the equation of a plane

in which the interpolated point lies. This method is reccomended if elevation data is organised in a regular grid.

‘nearest’Nearest neighbour interpolation. Z value at the interpolated

point takes on the same Z value as the closest known point. This method can work well for dense elevation data, say in the case using a point cloud from a digital elevation model. This method is the most computationally cheap to run.

‘spline’Like bilinear interpolation, a function is fitted to a

quadlerateral of known points encompassing the interpolated Z value. The computed function is a higher order than linear interpolation though. Works best for irregular spaced elevation data.

‘triangulate’: Triangulation method, best for complicated topographies

where an irregular grid of known topography points is not avialable. This is the default.

NoneNo interpolation method is used to interpolated topography on to

mesh, hence a flat mesh is returned.

Format of mesh_refinement:

The variable must have x,y,z and type columns ‘x’: x coordinate array like ‘y’: y coordinate array like ‘z’: z coordinate array like ‘type’: list, object array of point types, use the tag ‘surface’ for

surface points, use ‘buried’ for buried points.

‘cl’: array like of characteristic lengths for each refinement point in

the mesh.

resipy.meshTools.triMesh(elec_x, elec_z, elec_type=None, geom_input=None, keep_files=True, show_output=True, path='exe', dump=<built-in function print>, whole_space=False, model_err=False, handle=None, **kwargs)

Generates a triangular mesh for r2. Returns mesh class … this function expects the current working directory has path: exe/gmsh.exe. Uses gmsh version 3.0.6.

Parameters:
elec_x: array like

electrode x coordinates

elec_z: array like

electrode y coordinates

elec_type: list of strings, optional

List should be the same length as the electrode coordinate argument. Each entry in the list references the type of electrode: - ‘electrode’ = surface electrode coordinate, will be used to construct the topography in the mesh - ‘buried’ = buried electrode, placed the mesh surface - ‘remote’ = remote electrode (not actually placed in the mesh) - ‘borehole’ = borehole electrode, electrodes will be placed in the mesh with a line connecting them. borehole numbering starts at 1 and ascends numerically by 1.

geom_inputdict, optional

Allows for further customisation of the 2D mesh, its a dictionary contianing surface topography, polygons and boundaries

keep_filesboolean, optional

True if the gmsh input and output file is to be stored in the exe directory.

show_ouputboolean, optional

True if gmsh output is to be printed to console.

pathstring, optional

Path to exe folder (leave default unless you know what you are doing).

whole_space: boolean, optional

Flag for if the problem should be treated as a whole space porblem, in which case electrode type is ingored and all electrodes are buried in the middle of a large mesh.

model_errbool, optional

If True all topography will be normalised to zero such that the returned mesh can be used to model foward modelling errors.

dumpfunction, optional

Function to which pass the output during mesh generation. print() is the default.

handlevariable, optional

Will be assigned the output of ‘Popen’ in case the process needs to be killed in the UI for instance.

**kwargsoptional

Key word arguments to be passed to genGeoFile.

Returns:
mesh: class

<ResIPy> mesh class

Notes

geom_input format:

the code will cycle through numerically ordered keys (strings referencing objects in a dictionary”), currently the code expects a ‘surface’ and ‘electrode’ key for surface points and electrodes. the first borehole string should be given the key ‘borehole1’ and so on. The code stops searching for more keys when it cant find the next numeric key. Same concept goes for adding boundaries and polygons to the mesh. See below example:

geom_input = {‘surface’: [surf_x,surf_z],

‘boundary1’:[bound1x,bound1y], ‘polygon1’:[poly1x,poly1y]}

electrodes and electrode_type (if not None) format:

electrodes = [[x1,x2,x3,…],[y1,y2,y3,…]] electrode_type = [‘electrode’,’electrode’,’buried’,…]

like with geom_input, boreholes should be labelled borehole1, borehole2 and so on. The code will cycle through each borehole and internally sort them and add them to the mesh.

The code expects that all polygons, boundaries and electrodes fall within x values of the actaul survey area. So make sure your topography / surface electrode points cover the area you are surveying, otherwise some funky errors will occur in the mesh.

resipy.meshTools.voxelMesh(elec_x, elec_y, elec_z=None, elec_type=None, keep_files=True, interp_method='triangulate', surface_refinement=None, mesh_refinement=None, ball_refinement=True, path='exe', dump=<built-in function print>, whole_space=False, model_err=False, handle=None, show_output=True, **kwargs)

Generates a voxel mesh, not meant to be used for ERT processing with the R* codes but rather can be used for the purposes of mesh visualizuation post processing, and is more readily exportable into formats required by other geological software (e.g. Geovisionary).

Parameters:
elec_x: array like

electrode x coordinates

elec_y: array like

electrode y coordinates

elec_z: array like

electrode z coordinates

elec_type: list of strings, optional

Defines if electrodes are buried or not.

keep_filesboolean, optional

Not used. Kept as an argument for compatiblity with TetraMesh.

interp_method: string, default =’triangulate’ optional

Interpolation method to translate mesh nodes in the z direction. In other words the method in which topography is appended to the mesh. Here the topography is added to the mesh in post processing. The options are documented below in the notes.

surface_refinementnp.array, optional

Numpy array of shape (3,n), should follow the format np.array([x1,x2,x3,…],[y1,y2,y3,…],[z1,z2,z3,…]). Allows for extra refinement for the top surface of the mesh. The points are not added to the mesh, but considered in post processing in order to super impose topography on the mesh.

mesh_refinementdict, pd.DataFrame, optional

Not used. Kept as an argument for compatiblity with TetraMesh.

ball_refinement: boolean, optional

Not used. Kept as an argument for compatiblity with TetraMesh.

pathstring, optional

Not used. Kept as an argument for compatiblity with TetraMesh.

whole_space: boolean, optional

Flag for if the problem should be treated as a whole space porblem, in which case electrode type is ignored.

model_err: bool

If True, a flat mesh will be returned for the sake of estimating forward modelling errors.

dumpfunction, optional

Function to which pass the output during mesh generation. print() is the default.

**kwargs: dict

Keyword arguments to be passed to functions in gmshWrap.py. Pass cl=x to force the voxel mesh characteristic length.

force_regular: bool, optional

If passed as True, then return voxel mesh will be made of elements which are of the same size.

Returns:
mesh3d: class

Notes

Possible arguments for interp_method:
‘bilinear’4 known points are used to compute the equation of a plane

in which the interpolated point lies. This method is reccomended if elevation data is organised in a regular grid.

‘nearest’Nearest neighbour interpolation. Z value at the interpolated

point takes on the same Z value as the closest known point. This method can work well for dense elevation data, say in the case using a point cloud from a digital elevation model. This method is the most computationally cheap to run.

‘spline’Like bilinear interpolation, a function is fitted to a

quadlerateral of known points encompassing the interpolated Z value. The computed function is a higher order than linear interpolation though. Works best for irregular spaced elevation data.

‘triangulate’: Triangulation method, best for complicated topographies

where an irregular grid of known topography points is not avialable. This is the default.

NoneNo interpolation method is used to interpolated topography on to

mesh, hence a flat mesh is returned.

Format of mesh_refinement:

The variable must have x,y,z and type columns ‘x’: x coordinate array like ‘y’: y coordinate array like ‘z’: z coordinate array like ‘type’: list, object array of point types, use the tag ‘surface’ for

surface points, use ‘buried’ for buried points.

‘cl’: array like of characteristic lengths for each refinement point in

the mesh.

resipy.meshTools.vtk_import(file_path='mesh.vtk', order_nodes=True)

Imports a mesh file into the python workspace, can have triangular, quad or tetraheral shaped elements.

Parameters:
file_pathstring, optional

File path to mesh file. Note that a error will occur if the file format is not as expected.

order_nodesbool, optional

Order Nodes if true. Process can be resource intensive though. Reccomended if using the mesh for an inversion.

Returns:
meshclass

a <ResIPy> mesh class

resipy.meshTools.vtk_import_fmt4(file_path, order_nodes=True)

Vtk importer for newer format (not fully validated)