Package Core

Module z

exception mtpy.core.z.MT_Z_Error[source]
class mtpy.core.z.ResPhase(z_array=None, z_err_array=None, freq=None, **kwargs)[source]

resistivity and phase container

Attributes:
phase
phase_det
phase_det_err
phase_err
phase_err_xx
phase_err_xy
phase_err_yx
phase_err_yy
phase_xx
phase_xy
phase_yx
phase_yy
res_det
res_det_err
res_err_xx
res_err_xy
res_err_yx
res_err_yy
res_xx
res_xy
res_yx
res_yy
resistivity
resistivity_err

Methods

compute_resistivity_phase(self[, z_array, …]) compute resistivity and phase from z and z_err
set_res_phase(self, res_array, phase_array, freq) Set values for resistivity (res - in Ohm m) and phase (phase - in degrees), including error propagation.
compute_resistivity_phase(self, z_array=None, z_err_array=None, freq=None)[source]

compute resistivity and phase from z and z_err

set_res_phase(self, res_array, phase_array, freq, res_err_array=None, phase_err_array=None)[source]

Set values for resistivity (res - in Ohm m) and phase (phase - in degrees), including error propagation.

Parameters:
  • res_array (np.ndarray(num_freq, 2, 2)) – resistivity array in Ohm-m
  • phase_array (np.ndarray(num_freq, 2, 2)) – phase array in degrees
  • freq (np.ndarray(num_freq)) – frequency array in Hz
  • res_err_array (np.ndarray(num_freq, 2, 2)) – resistivity error array in Ohm-m
  • phase_err_array (np.ndarray(num_freq, 2, 2)) – phase error array in degrees
class mtpy.core.z.Tipper(tipper_array=None, tipper_err_array=None, freq=None)[source]

Tipper class –> generates a Tipper-object.

Errors are given as standard deviations (sqrt(VAR))

Parameters:
  • tipper_array (np.ndarray((nf, 1, 2), dtype='complex')) – tipper array in the shape of [Tx, Ty] default is None
  • tipper_err_array (np.ndarray((nf, 1, 2))) – array of estimated tipper errors in the shape of [Tx, Ty]. Must be the same shape as tipper_array. default is None
  • freq (np.ndarray(nf)) – array of frequencies corresponding to the tipper elements. Must be same length as tipper_array. default is None
Attributes Description
freq array of frequencies corresponding to elements of z
rotation_angle angle of which data is rotated by
tipper tipper array
tipper_err tipper error array
Methods Description
mag_direction computes magnitude and direction of real and imaginary induction arrows.
amp_phase computes amplitude and phase of Tx and Ty.
rotate rotates the data by the given angle
Attributes:
amplitude
amplitude_err
angle_err
angle_imag
angle_real
freq
mag_err
mag_imag
mag_real
phase
phase_err
tipper
tipper_err

Methods

compute_amp_phase(self) Sets attributes:
compute_mag_direction(self) computes the magnitude and direction of the real and imaginary induction vectors.
rotate(self, alpha) Rotate Tipper array.
set_amp_phase(self, r_array, phi_array) Set values for amplitude(r) and argument (phi - in degrees).
set_mag_direction(self, mag_real, ang_real, …) computes the tipper from the magnitude and direction of the real and imaginary components.
compute_amp_phase(self)[source]
Sets attributes:
  • amplitude
  • phase
  • amplitude_err
  • phase_err

values for resistivity are in in Ohm m and phase in degrees.

compute_mag_direction(self)[source]

computes the magnitude and direction of the real and imaginary induction vectors.

rotate(self, alpha)[source]

Rotate Tipper array.

Rotation angle must be given in degrees. All angles are referenced to geographic North=0, positive in clockwise direction. (Mathematically negative!)

In non-rotated state, ‘X’ refs to North and ‘Y’ to East direction.

Updates the attributes:
  • tipper
  • tipper_err
  • rotation_angle
set_amp_phase(self, r_array, phi_array)[source]

Set values for amplitude(r) and argument (phi - in degrees).

Updates the attributes:
  • tipper
  • tipper_err
set_mag_direction(self, mag_real, ang_real, mag_imag, ang_imag)[source]

computes the tipper from the magnitude and direction of the real and imaginary components.

Updates tipper

No error propagation yet

class mtpy.core.z.Z(z_array=None, z_err_array=None, freq=None)[source]

Z class - generates an impedance tensor (Z) object.

Z is a complex array of the form (n_freq, 2, 2), with indices in the following order:

  • Zxx: (0,0)
  • Zxy: (0,1)
  • Zyx: (1,0)
  • Zyy: (1,1)

All errors are given as standard deviations (sqrt(VAR))

Parameters:
  • z_array (numpy.ndarray(n_freq, 2, 2)) – array containing complex impedance values
  • z_err_array (numpy.ndarray(n_freq, 2, 2)) – array containing error values (standard deviation) of impedance tensor elements
  • freq (np.ndarray(n_freq)) – array of frequency values corresponding to impedance tensor elements.
Attributes Description
freq array of frequencies corresponding to elements of z
rotation_angle angle of which data is rotated by
z impedance tensor
z_err estimated errors of impedance tensor
resistivity apparent resisitivity estimated from z in Ohm-m
resistivity_err apparent resisitivity error
phase impedance phase (deg)
phase_err error in impedance phase
Methods Description
det calculates determinant of z with errors
invariants calculates the invariants of z
inverse calculates the inverse of z
remove_distortion removes distortion given a distortion matrix
remove_ss removes static shift by assumin Z = S * Z_0
norm calculates the norm of Z
only1d zeros diagonal components and computes the absolute valued mean of the off-diagonal components.
only2d zeros diagonal components
res_phase computes resistivity and phase
rotate rotates z positive clockwise, angle assumes North is 0.
set_res_phase recalculates z and z_err, needs attribute freq
skew calculates the invariant skew (off diagonal trace)
trace calculates the trace of z
Example:
>>> import mtpy.core.z as mtz
>>> import numpy as np
>>> z_test = np.array([[0+0j, 1+1j], [-1-1j, 0+0j]])
>>> z_object = mtz.Z(z_array=z_test, freq=[1])
>>> z_object.rotate(45)
>>> z_object.resistivity
Attributes:
det

Return the determinant of Z

det_err

Return the determinant of Z error

freq

Frequencies for each impedance tensor element

invariants

Return a dictionary of Z-invariants.

inverse

Return the inverse of Z.

norm

Return the 2-/Frobenius-norm of Z

norm_err

Return the 2-/Frobenius-norm of Z error

only_1d

Return Z in 1D form.

only_2d

Return Z in 2D form.

phase
phase_det
phase_det_err
phase_err
phase_err_xx
phase_err_xy
phase_err_yx
phase_err_yy
phase_xx
phase_xy
phase_yx
phase_yy
res_det
res_det_err
res_err_xx
res_err_xy
res_err_yx
res_err_yy
res_xx
res_xy
res_yx
res_yy
resistivity
resistivity_err
skew

Returns the skew of Z as defined by Z[0, 1] + Z[1, 0]

skew_err

Returns the skew error of Z as defined by Z_err[0, 1] + Z_err[1, 0]

trace

Return the trace of Z

trace_err

Return the trace of Z

z

Impedance tensor

z_err

Methods

compute_resistivity_phase(self[, z_array, …]) compute resistivity and phase from z and z_err
remove_distortion(self, distortion_tensor[, …]) Remove distortion D form an observed impedance tensor Z to obtain the uperturbed “correct” Z0: Z = D * Z0
remove_ss(self[, reduce_res_factor_x, …]) Remove the static shift by providing the respective correction factors for the resistivity in the x and y components.
rotate(self, alpha) Rotate Z array by angle alpha.
set_res_phase(self, res_array, phase_array, freq) Set values for resistivity (res - in Ohm m) and phase (phase - in degrees), including error propagation.
det

Return the determinant of Z

Returns:det_Z
Return type:np.ndarray(nfreq)
det_err

Return the determinant of Z error

Returns:det_Z_err
Return type:np.ndarray(nfreq)
freq

Frequencies for each impedance tensor element

Units are Hz.

invariants

Return a dictionary of Z-invariants.

inverse

Return the inverse of Z.

(no error propagtaion included yet)

norm

Return the 2-/Frobenius-norm of Z

Returns:norm
Return type:np.ndarray(nfreq)
norm_err

Return the 2-/Frobenius-norm of Z error

Returns:norm_err
Return type:np.ndarray(nfreq)
only_1d

Return Z in 1D form.

If Z is not 1D per se, the diagonal elements are set to zero, the off-diagonal elements keep their signs, but their absolute is set to the mean of the original Z off-diagonal absolutes.

only_2d

Return Z in 2D form.

If Z is not 2D per se, the diagonal elements are set to zero.

remove_distortion(self, distortion_tensor, distortion_err_tensor=None)[source]

Remove distortion D form an observed impedance tensor Z to obtain the uperturbed “correct” Z0: Z = D * Z0

Propagation of errors/uncertainties included

Parameters:
  • distortion_tensor (np.ndarray(2, 2, dtype=real)) – real distortion tensor as a 2x2
  • distortion_err_tensor – default is None
Return type:

np.ndarray(2, 2, dtype=’real’)

returns:impedance tensor with distorion removed

Return type:

np.ndarray(num_freq, 2, 2, dtype=’complex’)

returns:impedance tensor error after distortion is removed

Return type:

np.ndarray(num_freq, 2, 2, dtype=’complex’)

Example:
>>> import mtpy.core.z as mtz
>>> distortion = np.array([[1.2, .5],[.35, 2.1]])
>>> d, new_z, new_z_err = z_obj.remove_distortion(distortion)

remove_ss(self, reduce_res_factor_x=1.0, reduce_res_factor_y=1.0)[source]

Remove the static shift by providing the respective correction factors for the resistivity in the x and y components. (Factors can be determined by using the “Analysis” module for the impedance tensor)

Assume the original observed tensor Z is built by a static shift S and an unperturbated “correct” Z0 :

  • Z = S * Z0
therefore the correct Z will be :
  • Z0 = S^(-1) * Z
Parameters:
  • reduce_res_factor_x (float or iterable list or array) – static shift factor to be applied to x components (ie z[:, 0, :]). This is assumed to be in resistivity scale
  • reduce_res_factor_y (float or iterable list or array) – static shift factor to be applied to y components (ie z[:, 1, :]). This is assumed to be in resistivity scale
Returns:

static shift matrix,

Return type:

np.ndarray ((2, 2))

Returns:

corrected Z

Return type:

mtpy.core.z.Z

Note

The factors are in resistivity scale, so the entries of the matrix “S” need to be given by their square-roots!

rotate(self, alpha)[source]

Rotate Z array by angle alpha.

Rotation angle must be given in degrees. All angles are referenced to geographic North, positive in clockwise direction. (Mathematically negative!)

In non-rotated state, X refs to North and Y to East direction.

Updates the attributes
  • z
  • z_err
  • zrot
  • resistivity
  • phase
  • resistivity_err
  • phase_err
skew

Returns the skew of Z as defined by Z[0, 1] + Z[1, 0]

Note

This is not the MT skew, but simply the linear algebra skew

Returns:skew
Return type:np.ndarray(nfreq, 2, 2)
skew_err

Returns the skew error of Z as defined by Z_err[0, 1] + Z_err[1, 0]

Note

This is not the MT skew, but simply the linear algebra skew

Returns:skew_err
Return type:np.ndarray(nfreq, 2, 2)
trace

Return the trace of Z

Returns:Trace(z)
Return type:np.ndarray(nfreq, 2, 2)
trace_err

Return the trace of Z

Returns:Trace(z)
Return type:np.ndarray(nfreq, 2, 2)
z

Impedance tensor

np.ndarray(nfreq, 2, 2)

mtpy.core.z.correct4sensor_orientation(Z_prime, Bx=0, By=90, Ex=0, Ey=90, Z_prime_error=None)[source]

Correct a Z-array for wrong orientation of the sensors.

Assume, E’ is measured by sensors orientated with the angles
E’x: a E’y: b
Assume, B’ is measured by sensors orientated with the angles
B’x: c B’y: d
With those data, one obtained the impedance tensor Z’:
E’ = Z’ * B’
Now we define change-of-basis matrices T,U so that
E = T * E’ B = U * B’

=> T contains the expression of the E’-basis in terms of E (the standard basis) and U contains the expression of the B’-basis in terms of B (the standard basis) The respective expressions for E’x-basis vector and E’y-basis vector are the columns of T. The respective expressions for B’x-basis vector and B’y-basis vector are the columns of U.

We obtain the impedance tensor in default coordinates as:

E’ = Z’ * B’ => T^(-1) * E = Z’ * U^(-1) * B
=> E = T * Z’ * U^(-1) * B => Z = T * Z’ * U^(-1)
Parameters:
  • Z_prime – impedance tensor to be adjusted
  • Bx (float (angle in degrees)) – orientation of Bx relative to geographic north (0) default is 0
  • By
  • Ex (float (angle in degrees)) – orientation of Ex relative to geographic north (0) default is 0
  • Ey (float (angle in degrees)) – orientation of Ey relative to geographic north (0) default is 90
  • Z_prime_error (np.ndarray(Z_prime.shape)) – impedance tensor error (std) default is None
Dtype Z_prime:

np.ndarray(num_freq, 2, 2, dtype=’complex’)

Returns:

adjusted impedance tensor

Return type:

np.ndarray(Z_prime.shape, dtype=’complex’)

Returns:

impedance tensor standard deviation in default orientation

Return type:

np.ndarray(Z_prime.shape, dtype=’real’)

Module TS

class mtpy.core.ts.MT_TS(**kwargs)[source]

MT time series object that will read/write data in different formats including hdf5, txt, miniseed.

The foundations are based on Pandas Python package.

The data are store in the variable ts, which is a pandas dataframe with the data in the column ‘data’. This way the data can be indexed as a numpy array:

>>> MT_TS.ts['data'][0:256]

or

>>> MT_TS.ts.data[0:256]

Also, the data can be indexed by time (note needs to be exact time):

>>> MT_TS.ts['2017-05-04 12:32:00.0078125':'2017-05-05 12:35:00]

Input ts as a numpy.ndarray or Pandas DataFrame

Metadata Description
azimuth clockwise angle from coordinate system N (deg)
calibration_fn file name for calibration data
component component name [ ‘ex’ | ‘ey’ | ‘hx’ | ‘hy’ | ‘hz’]
coordinate_system [ geographic | geomagnetic ]
datum datum of geographic location ex. WGS84
declination geomagnetic declination (deg)
dipole_length length of dipole (m)
data_logger data logger type
instrument_id ID number of instrument for calibration
lat latitude of station in decimal degrees
lon longitude of station in decimal degrees
n_samples number of samples in time series
sampling_rate sampling rate in samples/second
start_time_epoch_sec start time in epoch seconds
start_time_utc start time in UTC
station station name
units units of time series

Note

Currently only supports hdf5 and text files

Method Description
read_hdf5 read an hdf5 file
write_hdf5 write an hdf5 file
write_ascii_file write an ascii file
read_ascii_file read an ascii file
Example:
>>> import mtpy.core.ts as ts
>>> import numpy as np
>>> mt_ts = ts.MT_TS()
>>> mt_ts.ts = np.random.randn(1024)
>>> mt_ts.station = 'test'
>>> mt_ts.lon = 30.00
>>> mt_ts.lat = -122.00
>>> mt_ts.component = 'HX'
>>> mt_ts.units = 'counts'
>>> mt_ts.write_hdf5(r"/home/test.h5")
Attributes:
elev

elevation in elevation units

lat

Latitude in decimal degrees

lon

Longitude in decimal degrees

n_samples

number of samples

sampling_rate

sampling rate in samples/second

start_time_epoch_sec

start time in epoch seconds

start_time_utc

start time in UTC given in time format

ts

Methods

apply_addaptive_notch_filter(self[, …]) apply notch filter to the data that finds the peak around each frequency.
decimate(self[, dec_factor]) decimate the data by using scipy.signal.decimate
low_pass_filter(self[, low_pass_freq, …]) low pass the data
plot_spectra(self[, spectra_type]) Plot spectra using the spectral type
read_ascii(self, fn_ascii) Read an ascii format file with metadata
read_ascii_header(self, fn_ascii) Read an ascii metadata
read_hdf5(self, fn_hdf5[, …]) Read an hdf5 file with metadata using Pandas.
write_ascii_file(self, fn_ascii[, chunk_size]) Write an ascii format file with metadata
write_hdf5(self, fn_hdf5[, …]) Write an hdf5 file with metadata using pandas to write the file.
apply_addaptive_notch_filter(self, notches=None, notch_radius=0.5, freq_rad=0.5, rp=0.1)[source]

apply notch filter to the data that finds the peak around each frequency.

see mtpy.processing.filter.adaptive_notch_filter

Parameters:notch_dict (dictionary) – dictionary of filter parameters. if an empty dictionary is input the filter looks for 60 Hz and harmonics to filter out.
decimate(self, dec_factor=1)[source]

decimate the data by using scipy.signal.decimate

Parameters:dec_factor (int) – decimation factor
  • refills ts.data with decimated data and replaces sampling_rate
elev

elevation in elevation units

lat

Latitude in decimal degrees

lon

Longitude in decimal degrees

low_pass_filter(self, low_pass_freq=15, cutoff_freq=55)[source]

low pass the data

Parameters:
  • low_pass_freq (float) – low pass corner in Hz
  • cutoff_freq (float) – cut off frequency in Hz
  • filters ts.data
n_samples

number of samples

plot_spectra(self, spectra_type='welch', **kwargs)[source]

Plot spectra using the spectral type

Note

Only spectral type supported is welch

Parameters:

spectra_type – [ ‘welch’ ]

Example:
>>> ts_obj = mtts.MT_TS()
>>> ts_obj.read_hdf5(r"/home/MT/mt01.h5")
>>> ts_obj.plot_spectra()
read_ascii(self, fn_ascii)[source]

Read an ascii format file with metadata

Parameters:

fn_ascii (string) – full path to ascii file

Example:
>>> ts_obj.read_ascii(r"/home/ts/mt01.EX")
read_ascii_header(self, fn_ascii)[source]

Read an ascii metadata

Parameters:

fn_ascii (string) – full path to ascii file

Example:
>>> ts_obj.read_ascii_header(r"/home/ts/mt01.EX")
read_hdf5(self, fn_hdf5, compression_level=0, compression_lib='blosc')[source]

Read an hdf5 file with metadata using Pandas.

Parameters:
  • fn_hdf5 (string) – full path to hdf5 file, has .h5 extension
  • compression_level (int) – compression level of file [ 0-9 ]
  • compression_lib (string) – compression library default is blosc
Returns:

fn_hdf5

See also

Pandas.HDf5Store

sampling_rate

sampling rate in samples/second

start_time_epoch_sec

start time in epoch seconds

start_time_utc

start time in UTC given in time format

write_ascii_file(self, fn_ascii, chunk_size=4096)[source]

Write an ascii format file with metadata

Parameters:
  • fn_ascii (string) – full path to ascii file
  • chunk_size (int) – read in file by chunks for efficiency
Example:
>>> ts_obj.write_ascii_file(r"/home/ts/mt01.EX")
write_hdf5(self, fn_hdf5, compression_level=0, compression_lib='blosc')[source]

Write an hdf5 file with metadata using pandas to write the file.

Parameters:
  • fn_hdf5 (string) – full path to hdf5 file, has .h5 extension
  • compression_level (int) – compression level of file [ 0-9 ]
  • compression_lib (string) – compression library default is blosc
Returns:

fn_hdf5

See also

Pandas.HDf5Store

exception mtpy.core.ts.MT_TS_Error[source]
class mtpy.core.ts.Spectra(**kwargs)[source]

compute spectra of time series

Methods

compute_spectra(self, data, spectra_type, …) compute spectra according to input type
welch_method(self, data[, plot]) Compute the spectra using the Welch method, which is an average spectra of the data.
compute_spectra(self, data, spectra_type, **kwargs)[source]

compute spectra according to input type

welch_method(self, data, plot=True, **kwargs)[source]

Compute the spectra using the Welch method, which is an average spectra of the data. Computes short time window of length nperseg and averages them to reduce noise.

Module MT

class mtpy.core.mt.Citation(**kwargs)[source]

Information for a citation.

Holds the following information:

Attributes Type Explanation
author string Author names
title string Title of article, or publication
journal string Name of journal
doi string DOI number (doi:10.110/sf454)
year int year published

More attributes can be added by inputing a key word dictionary

>>> Citation(**{'volume':56, 'pages':'234--214'})
class mtpy.core.mt.Copyright(**kwargs)[source]

Information of copyright, mainly about how someone else can use these data. Be sure to read over the conditions_of_use.

Holds the following information:

Attributes Type Explanation
citation Citation citation of published work using these data
conditions_of_use string conditions of use of these data
release_status string release status [ open | public | proprietary]

More attributes can be added by inputing a key word dictionary

>>> Copyright(**{'owner':'University of MT', 'contact':'Cagniard'})
class mtpy.core.mt.DataQuality(**kwargs)[source]

Information on data quality.

Holds the following information:

Attributes Type Explanation
comments string comments on data quality
good_from_period float minimum period data are good
good_to_period float maximum period data are good
rating int [1-5]; 1 = poor, 5 = excellent
warrning_comments string any comments on warnings in the data
warnings_flag int [0-#of warnings]

More attributes can be added by inputing a key word dictionary

>>>DataQuality(**{‘time_series_comments’:’Periodic Noise’})

class mtpy.core.mt.FieldNotes(**kwargs)[source]

Field note information.

Holds the following information:

Attributes Type Explanation
data_quality DataQuality notes on data quality
electrode Instrument type of electrode used
data_logger Instrument type of data logger
magnetometer Instrument type of magnetotmeter

More attributes can be added by inputing a key word dictionary

>>> FieldNotes(**{'electrode_ex':'Ag-AgCl 213', 'magnetometer_hx':'102'})
class mtpy.core.mt.Instrument(**kwargs)[source]

Information on an instrument that was used.

Holds the following information:

Attributes Type Explanation
id string serial number or id number of data logger
manufacturer string company whom makes the instrument
type string Broadband, long period, something else

More attributes can be added by inputing a key word dictionary

>>> Instrument(**{'ports':'5', 'gps':'time_stamped'})
class mtpy.core.mt.Location(**kwargs)[source]

location details

Attributes:
easting
elevation
latitude
longitude
northing

Methods

project_location2ll(self) project location coordinates into meters given the reference ellipsoid, for now that is constrained to WGS84
project_location2utm(self) project location coordinates into meters given the reference ellipsoid, for now that is constrained to WGS84
project_location2ll(self)[source]

project location coordinates into meters given the reference ellipsoid, for now that is constrained to WGS84

Returns East, North, Zone

project_location2utm(self)[source]

project location coordinates into meters given the reference ellipsoid, for now that is constrained to WGS84

Returns East, North, Zone

class mtpy.core.mt.MT(fn=None, **kwargs)[source]

Basic MT container to hold all information necessary for a MT station including the following parameters.

  • Site –> information on site details (lat, lon, name, etc)
  • FieldNotes –> information on instruments, setup, etc.
  • Copyright –> information on how the data can be used and citations
  • Provenance –> where the data come from and how they are stored
  • Processing –> how the data were processed.

The most used attributes are made available from MT, namely the following.

Attribute Description
station station name
lat station latitude in decimal degrees
lon station longitude in decimal degrees
elev station elevation in meters
Z mtpy.core.z.Z object for impedance tensor
Tipper mtpy.core.z.Tipper object for tipper
pt mtpy.analysis.pt.PhaseTensor for phase tensor
east station location in UTM coordinates assuming WGS-84
north station location in UTM coordinates assuming WGS-84
utm_zone zone of UTM coordinates assuming WGS-84
rotation_angle rotation angle of the data
fn absolute path to the data file

Other information is contained with in the different class attributes. For instance survey name is in MT.Site.survey

Note

  • The best way to see what all the information is and where it is contained would be to write out a configuration file

    >>> import mtpy.core.mt as mt
    >>> mt_obj = mt.MT()
    >>> mt_obj.write_cfg_file(r"/home/mt/generic.cfg")
    
  • Currently EDI, XML, and j file are supported to read in information, and can write out EDI and XML formats. Will be extending to j and Egberts Z format.

Methods Description
read_mt_file read in a MT file [ EDI | XML | j ]
write_mt_file write a MT file [ EDI | XML ]
read_cfg_file read a configuration file
write_cfg_file write a configuration file
remove_distortion remove distortion following Bibby et al. [2005]
remove_static_shift Shifts apparent resistivity curves up or down
interpolate interpolates Z and T onto specified frequency array.

Examples

Read from an .edi File:
 
>>> import mtpy.core.mt as mt
>>> mt_obj = mt.MT(r"/home/edi_files/s01.edi")
Remove Distortion:
 
>>> import mtpy.core.mt as mt
>>> mt1 = mt.MT(fn=r"/home/mt/edi_files/mt01.edi")
>>> D, new_z = mt1.remove_distortion()
>>> mt1.write_mt_file(new_fn=r"/home/mt/edi_files/mt01_dr.edi",        >>>                    new_Z=new_z)
Remove Static Shift:
 
>>> new_z_obj = mt_obj.remove_static_shift(ss_x=.78, ss_y=1.1)
>>> # write a new edi file
>>> mt_obj.write_mt_file(new_fn=r"/home/edi_files/s01_ss.edi",
>>>                       new_Z=new_z)
>>> wrote file to: /home/edi_files/s01_ss.edi
Interpolate:
>>> new_freq = np.logspace(-3, 3, num=24)
>>> new_z_obj, new_tipper_obj = mt_obj.interpolate(new_freq)
>>> mt_obj.write_mt_file(new_Z=new_z_obj, new_Tipper=new_tipper_obj)
>>> wrote file to: /home/edi_files/s01_RW.edi
Attributes:
Tipper

mtpy.core.z.Tipper object to hold tipper information

Z

mtpy.core.z.Z object to hole impedance tensor

east

easting (m)

elev

Elevation

fn

reference to original data file

lat

Latitude

lon

Longitude

north

northing (m)

pt

mtpy.analysis.pt.PhaseTensor object to hold phase tensor

rotation_angle

rotation angle in degrees from north

station

station name

utm_zone

utm zone

Methods

interpolate(self, new_freq_array[, …]) Interpolate the impedance tensor onto different frequencies
plot_mt_response(self, \*\*kwargs) Returns a mtpy.imaging.plotresponse.PlotResponse object
read_cfg_file(self, cfg_fn) Read in a configuration file and populate attributes accordingly.
read_mt_file(self, fn[, file_type]) Read an MT response file.
remove_distortion(self[, num_freq]) remove distortion following Bibby et al.
remove_static_shift(self[, ss_x, ss_y]) Remove static shift from the apparent resistivity
write_cfg_file(self, cfg_fn) Write a configuration file for the MT sections
write_mt_file(self[, save_dir, fn_basename, …]) Write an mt file, the supported file types are EDI and XML.
Tipper

mtpy.core.z.Tipper object to hold tipper information

Z

mtpy.core.z.Z object to hole impedance tensor

east

easting (m)

elev

Elevation

fn

reference to original data file

interpolate(self, new_freq_array, interp_type='slinear', bounds_error=True, period_buffer=None)[source]

Interpolate the impedance tensor onto different frequencies

Parameters:

new_freq_array (np.ndarray) – a 1-d array of frequencies to interpolate on to. Must be with in the bounds of the existing frequency range, anything outside and an error will occur.

Returns:

a new impedance object with the corresponding frequencies and components.

Return type:

mtpy.core.z.Z

Returns:

a new tipper object with the corresponding frequencies and components.

Return type:

mtpy.core.z.Tipper

Interpolate:
>>> import mtpy.core.mt as mt
>>> edi_fn = r"/home/edi_files/mt_01.edi"
>>> mt_obj = mt.MT(edi_fn)
>>> # create a new frequency range to interpolate onto
>>> new_freq = np.logspace(-3, 3, 24)
>>> new_z_object, new_tipper_obj = mt_obj.interpolate(new_freq)
>>> mt_obj.write_mt_file(new_fn=r"/home/edi_files/mt_01_interp.edi",
>>> ...                   new_Z_obj=new_z_object,
>>> ...                   new_Tipper_obj=new_tipper_object)
lat

Latitude

lon

Longitude

north

northing (m)

plot_mt_response(self, **kwargs)[source]

Returns a mtpy.imaging.plotresponse.PlotResponse object

Plot Response:
>>> mt_obj = mt.MT(edi_file)
>>> pr = mt.plot_mt_response()
>>> # if you need more info on plot_mt_response
>>> help(pr)
pt

mtpy.analysis.pt.PhaseTensor object to hold phase tensor

read_cfg_file(self, cfg_fn)[source]

Read in a configuration file and populate attributes accordingly.

The configuration file should be in the form:
Site.Location.latitude = 46.5
Site.Location.longitude = 122.7
Site.Location.datum = ‘WGS84’

Processing.Software.name = BIRRP
Processing.Software.version = 5.2.1

Provenance.Creator.name = L. Cagniard
Provenance.Submitter.name = I. Larionov
Parameters:cfg_fn (string) – full path to configuration file

Note

The best way to make a configuration file would be to save a configuration file first from MT, then filling in the fields.

Make configuration file:
 
>>> import mtpy.core.mt as mt
>>> mt_obj = mt.MT()
>>> mt_obj.write_cfg_file(r"/mt/generic_config.cfg")
Read in configuration file:
 
>>> import mtpy.core.mt as mt
>>> mt_obj = mt.MT()
>>> mt_obj.read_cfg_file(r"/home/mt/survey_config.cfg")
read_mt_file(self, fn, file_type=None)[source]

Read an MT response file.

Note

Currently only .edi, .xml, and .j files are supported

Parameters:
  • fn (string) – full path to input file
  • file_type (string) – [‘edi’ | ‘j’ | ‘xml’ | … ] if None, automatically detects file type by the extension.
Example:
>>> import mtpy.core.mt as mt
>>> mt_obj = mt.MT()
>>> mt_obj.read_mt_file(r"/home/mt/mt01.xml")
remove_distortion(self, num_freq=None)[source]

remove distortion following Bibby et al. [2005].

Parameters:

num_freq (int) – number of frequencies to look for distortion from the highest frequency

Returns:

Distortion matrix

Return type:

np.ndarray(2, 2, dtype=real)

Returns:

Z with distortion removed

Return type:

mtpy.core.z.Z

Remove distortion and write new .edi file:
 
>>> import mtpy.core.mt as mt
>>> mt1 = mt.MT(fn=r"/home/mt/edi_files/mt01.edi")
>>> D, new_z = mt1.remove_distortion()
>>> mt1.write_mt_file(new_fn=r"/home/mt/edi_files/mt01_dr.edi",            >>>                    new_Z=new_z)
remove_static_shift(self, ss_x=1.0, ss_y=1.0)[source]

Remove static shift from the apparent resistivity

Assume the original observed tensor Z is built by a static shift S and an unperturbated “correct” Z0 :

  • Z = S * Z0
therefore the correct Z will be :
  • Z0 = S^(-1) * Z
Parameters:
  • ss_x (float) – correction factor for x component
  • ss_y (float) – correction factor for y component
Returns:

new Z object with static shift removed

Return type:

mtpy.core.z.Z

Note

The factors are in resistivity scale, so the entries of the matrix “S” need to be given by their square-roots!

Remove Static Shift:
 
>>> import mtpy.core.mt as mt
>>> mt_obj = mt.MT(r"/home/mt/mt01.edi")
>>> new_z_obj = mt.remove_static_shift(ss_x=.5, ss_y=1.2)
>>> mt_obj.write_mt_file(new_fn=r"/home/mt/mt01_ss.edi",
>>> ...                   new_Z_obj=new_z_obj)
rotation_angle

rotation angle in degrees from north

station

station name

utm_zone

utm zone

write_cfg_file(self, cfg_fn)[source]

Write a configuration file for the MT sections

Parameters:

cfg_fn (string) – full path to configuration file to write to

Return cfg_fn:

full path to configuration file

Rtype cfg_fn:

string

Write configuration file:
 
>>> import mtpy.core.mt as mt
>>> mt_obj = mt.MT()
>>> mt_obj.read_mt_file(r"/home/mt/edi_files/mt01.edi")
>>> mt_obj.write_cfg_file(r"/home/mt/survey_config.cfg")
write_mt_file(self, save_dir=None, fn_basename=None, file_type='edi', new_Z_obj=None, new_Tipper_obj=None, longitude_format='LON', latlon_format='dms')[source]

Write an mt file, the supported file types are EDI and XML.

Parameters:
  • save_dir (string) – full path save directory
  • fn_basename (string) – name of file with or without extension
  • file_type (string) – [ ‘edi’ | ‘xml’ ]
  • new_Z_obj (mtpy.core.z.Z) – new Z object
  • new_Tipper_obj (mtpy.core.z.Tipper) – new Tipper object
  • longitude_format (string) – whether to write longitude as LON or LONG. options are ‘LON’ or ‘LONG’, default ‘LON’
  • latlon_format (string) – format of latitude and longitude in output edi, degrees minutes seconds (‘dms’) or decimal degrees (‘dd’)
Returns:

full path to file

Return type:

string

Example:
>>> mt_obj.write_mt_file(file_type='xml')
exception mtpy.core.mt.MT_Error[source]
class mtpy.core.mt.Person(**kwargs)[source]

Information for a person

Holds the following information:

Attributes Type Explanation
email string email of person
name string name of person
organization string name of person’s organization
organization_url string organizations web address

More attributes can be added by inputing a key word dictionary

>>> Person(**{'phone':'650-888-6666'})
class mtpy.core.mt.Processing(**kwargs)[source]

Information for a processing

Holds the following information:

Attributes Type Explanation
email string email of person
name string name of person
organization string name of person’s organization
organization_url string organizations web address

More attributes can be added by inputing a key word dictionary

>>> Person(**{'phone':'888-867-5309'})
class mtpy.core.mt.Provenance(**kwargs)[source]

Information of the file history, how it was made

Holds the following information:

Attributes Type Explanation
creation_time string creation time of file YYYY-MM-DD,hh:mm:ss
creating_application string name of program creating the file
creator Person person whom created the file
submitter Person person whom is submitting file for archiving

More attributes can be added by inputing a key word dictionary

>>> Provenance(**{'archive':'IRIS', 'reprocessed_by':'grad_student'})
class mtpy.core.mt.Site(**kwargs)[source]

Information on the site, including location, id, etc.

Holds the following information:

Attributes Type Explanation
aqcuired_by string name of company or person whom aqcuired the data.
id string station name
Location object Location Holds location information, lat, lon, elev datum, easting, northing see Location class
start_date string YYYY-MM-DD start date of measurement
end_date string YYYY-MM-DD end date of measurement
year_collected string year data collected
survey string survey name
project string project name
run_list string list of measurment runs ex. [mt01a, mt01b]

More attributes can be added by inputing a key word dictionary

>>> Site(**{'state':'Nevada', 'Operator':'MTExperts'})
Attributes:
year_collected
class mtpy.core.mt.Software(**kwargs)[source]

software

Module EDI

class mtpy.core.edi.DataSection(edi_fn=None, edi_lines=None)[source]

DataSection contains the small metadata block that describes which channel is which. A typical block looks like:

>=MTSECT

    ex=1004.001
    ey=1005.001
    hx=1001.001
    hy=1002.001
    hz=1003.001
    nfreq=14
    sectid=par28ew
    nchan=None
    maxblks=None
Parameters:edi_fn (string) – full path to .edi file to read in.
[Rf5ecdd4b8de8-1]Changes these values to change what is written to edi file

Methods

get_data_sect(self) read in the data of the file, will detect if reading spectra or impedance.
read_data_sect(self[, data_sect_list]) read data section
write_data_sect(self[, data_sect_list, …]) write a data section
get_data_sect(self)[source]

read in the data of the file, will detect if reading spectra or impedance.

read_data_sect(self, data_sect_list=None)[source]

read data section

write_data_sect(self, data_sect_list=None, over_dict=None)[source]

write a data section

class mtpy.core.edi.DefineMeasurement(edi_fn=None, edi_lines=None)[source]

DefineMeasurement class holds information about the measurement. This includes how each channel was setup. The main block contains information on the reference location for the station. This is a bit of an archaic part and was meant for a multiple station .edi file. This section is also important if you did any forward modeling with Winglink cause it only gives the station location in this section. The other parts are how each channel was collected. An example define measurement section looks like:

>=DEFINEMEAS

    MAXCHAN=7
    MAXRUN=999
    MAXMEAS=9999
    UNITS=M
    REFTYPE=CART
    REFLAT=-30:12:49.4693
    REFLONG=139:47:50.87
    REFELEV=0

>HMEAS ID=1001.001 CHTYPE=HX X=0.0 Y=0.0 Z=0.0 AZM=0.0
>HMEAS ID=1002.001 CHTYPE=HY X=0.0 Y=0.0 Z=0.0 AZM=90.0
>HMEAS ID=1003.001 CHTYPE=HZ X=0.0 Y=0.0 Z=0.0 AZM=0.0
>EMEAS ID=1004.001 CHTYPE=EX X=0.0 Y=0.0 Z=0.0 X2=0.0 Y2=0.0
>EMEAS ID=1005.001 CHTYPE=EY X=0.0 Y=0.0 Z=0.0 X2=0.0 Y2=0.0
>HMEAS ID=1006.001 CHTYPE=HX X=0.0 Y=0.0 Z=0.0 AZM=0.0
>HMEAS ID=1007.001 CHTYPE=HY X=0.0 Y=0.0 Z=0.0 AZM=90.0
Parameters:edi_fn (string) – full path to .edi file to read in.
[R5ea4377773dd-1]Each channel with have its own define measurement and depending on whether it is an E or H channel the metadata will be different. the #### correspond to the channel number.
[R5ea4377773dd-2]Internally everything is converted to decimal degrees. Output is written as HH:MM:SS.ss so Winglink can read them in.
[R5ea4377773dd-3]

If you want to change what metadata is written into the .edi file change the items in _header_keys. Default attributes are:

  • maxchan
  • maxrun
  • maxmeas
  • reflat
  • reflon
  • refelev
  • reftype
  • units

Methods

get_measurement_dict(self) get a dictionary for the xmeas parts
get_measurement_lists(self) get measurement list including measurement setup
read_define_measurement(self[, measurement_list]) read the define measurment section of the edi file
write_define_measurement(self[, …]) write the define measurement block as a list of strings
get_measurement_dict(self)[source]

get a dictionary for the xmeas parts

get_measurement_lists(self)[source]

get measurement list including measurement setup

read_define_measurement(self, measurement_list=None)[source]

read the define measurment section of the edi file

should be a list with lines for:
  • maxchan

  • maxmeas

  • maxrun

  • refelev

  • reflat

  • reflon

  • reftype

  • units

  • dictionaries for >XMEAS with keys:
    • id
    • chtype
    • x
    • y
    • axm

    -acqchn

write_define_measurement(self, measurement_list=None, longitude_format='LON', latlon_format='dd')[source]

write the define measurement block as a list of strings

class mtpy.core.edi.EMeasurement(**kwargs)[source]

EMeasurement contains metadata for an electric field measurement

Attributes Description
id Channel number
chtype [ EX | EY ]
x x (m) north from reference point (station) of one electrode of the dipole
y y (m) east from reference point (station) of one electrode of the dipole
x2 x (m) north from reference point (station) of the other electrode of the dipole
y2 y (m) north from reference point (station) of the other electrode of the dipole
acqchan name of the channel acquired usually same as chtype
Fill Metadata:
>>> import mtpy.core.edi as mtedi
>>> e_dict = {'id': '1', 'chtype':'ex', 'x':0, 'y':0, 'x2':50, 'y2':50}
>>> e_dict['acqchn'] = 'ex'
>>> emeas = mtedi.EMeasurement(**e_dict)
class mtpy.core.edi.Edi(edi_fn=None)[source]

This class is for .edi files, mainly reading and writing. Has been tested on Winglink and Phoenix output .edi’s, which are meant to follow the archaic EDI format put forward by SEG. Can read impedance, Tipper and/or spectra data.

The Edi class contains a class for each major section of the .edi file.

Frequency and components are ordered from highest to lowest frequency.

Parameters:edi_fn (string) – full path to .edi file to be read in. default is None. If an .edi file is input, it is automatically read in and attributes of Edi are filled
Methods Description
read_edi_file Reads in an edi file and populates the associated classes and attributes.
write_edi_file Writes an .edi file following the EDI format given the apporpriate attributes are filled. Writes out in impedance and Tipper format.
_read_data Reads in the impedance and Tipper blocks, if the .edi file is in ‘spectra’ format, read_data converts the data to impedance and Tipper.
_read_mt Reads impedance and tipper data from the appropriate blocks of the .edi file.
_read_spectra Reads in spectra data and converts it to impedance and Tipper data.
Attributes Description default
Data_sect DataSection class, contains basic information on the data collected and in whether the data is in impedance or spectra.  
Define_measurement DefineMeasurement class, contains information on how the data was collected.  
edi_fn full path to edi file read in None
Header Header class, contains metadata on where, when, and who collected the data  
Info Information class, contains information on how the data was processed and how the transfer functions where estimated.  
Tipper mtpy.core.z.Tipper class, contains the tipper data  
Z mtpy.core.z.Z class, contains the impedance data  
_block_len number of data in one line. 6
_data_header_str header string for each of the data section ‘>!****{0}****!’
_num_format string format of data. ‘ 15.6e’
_t_labels labels for tipper blocks  
_z_labels labels for impedance blocks  
Change Latitude:
 
>>> import mtpy.core.edi as mtedi
>>> edi_obj = mtedi.Edi(edi_fn=r"/home/mt/mt01.edi")
>>> # change the latitude
>>> edi_obj.header.lat = 45.7869
>>> new_edi_fn = edi_obj.write_edi_file()
Attributes:
elev

Elevation in elevation units

lat

latitude in decimal degrees

lon

longitude in decimal degrees

station

station name

Methods

read_edi_file(self[, edi_fn]) Read in an edi file and fill attributes of each section’s classes.
write_edi_file(self[, new_edi_fn, …]) Write a new edi file from either an existing .edi file or from data input by the user into the attributes of Edi.
elev

Elevation in elevation units

lat

latitude in decimal degrees

lon

longitude in decimal degrees

read_edi_file(self, edi_fn=None)[source]

Read in an edi file and fill attributes of each section’s classes. Including:

  • Header
  • Info
  • Define_measurement
  • Data_sect
  • Z
  • Tipper

Note

Automatically detects if data is in spectra format. All data read in is converted to impedance and Tipper.

Parameters:

edi_fn (string) – full path to .edi file to be read in default is None

Example:
>>> import mtpy.core.Edi as mtedi
>>> edi_obj = mtedi.Edi()
>>> edi_obj.read_edi_file(edi_fn=r"/home/mt/mt01.edi")
station

station name

write_edi_file(self, new_edi_fn=None, longitude_format='LON', latlon_format='dms')[source]

Write a new edi file from either an existing .edi file or from data input by the user into the attributes of Edi.

Parameters:
  • new_edi_fn (string) – full path to new edi file. default is None, which will write to the same file as the input .edi with as: r”/home/mt/mt01_1.edi”
  • longitude_format (string) – whether to write longitude as LON or LONG. options are ‘LON’ or ‘LONG’, default ‘LON’
  • latlon_format (string) – format of latitude and longitude in output edi, degrees minutes seconds (‘dms’) or decimal degrees (‘dd’)
Returns:

full path to new edi file

Return type:

string

Example:
>>> import mtpy.core.edi as mtedi
>>> edi_obj = mtedi.Edi(edi_fn=r"/home/mt/mt01/edi")
>>> edi_obj.Header.dataid = 'mt01_rr'
>>> n_edi_fn = edi_obj.write_edi_file()
class mtpy.core.edi.HMeasurement(**kwargs)[source]

HMeasurement contains metadata for a magnetic field measurement

Attributes Description
id Channel number
chtype [ HX | HY | HZ | RHX | RHY ]
x x (m) north from reference point (station)
y y (m) east from reference point (station)
azm angle of sensor relative to north = 0
acqchan name of the channel acquired usually same as chtype
Fill Metadata:
>>> import mtpy.core.edi as mtedi
>>> h_dict = {'id': '1', 'chtype':'hx', 'x':0, 'y':0, 'azm':0}
>>> h_dict['acqchn'] = 'hx'
>>> hmeas = mtedi.HMeasurement(**h_dict)
class mtpy.core.edi.Header(edi_fn=None, **kwargs)[source]

Header class contains all the information in the header section of the .edi file. A typical header block looks like:

>HEAD

    ACQBY=None
    ACQDATE=None
    DATAID=par28ew
    ELEV=0.000
    EMPTY=1e+32
    FILEBY=WG3DForward
    FILEDATE=2016/04/11 19:37:37 UTC
    LAT=-30:12:49
    LOC=None
    LON=139:47:50
    PROGDATE=2002-04-22
    PROGVERS=WINGLINK EDI 1.0.22
    COORDINATE SYSTEM = GEOGRAPHIC NORTH
    DECLINATION = 10.0
Parameters:edi_fn (string) – full path to .edi file to be read in. default is None. If an .edi file is input attributes of Header are filled.

Many of the attributes are needed in the .edi file. They are marked with a yes for ‘In .edi’

[R60960842fb28-1]Internally everything is converted to decimal degrees. Output is written as HH:MM:SS.ss so Winglink can read them in.
[R60960842fb28-2]

If you want to change what metadata is written into the .edi file change the items in _header_keys. Default attributes are:

  • acqby
  • acqdate
  • coordinate_system
  • dataid
  • declination
  • elev
  • fileby
  • lat
  • loc
  • lon
  • filedate
  • empty
  • progdate
  • progvers
Methods Description
get_header_list get header lines from edi file
read_header read in header information from header_lines
write_header write header lines, returns a list of lines to write
Read Header:
>>> import mtpy.core.edi as mtedi
>>> header_obj = mtedi.Header(edi_fn=r"/home/mt/mt01.edi")

Methods

get_header_list(self) Get the header information from the .edi file in the form of a list, where each item is a line in the header section.
read_header(self[, header_list]) read a header information from either edi file or a list of lines containing header information.
write_header(self[, header_list, …]) Write header information to a list of lines.
get_header_list(self)[source]

Get the header information from the .edi file in the form of a list, where each item is a line in the header section.

read_header(self, header_list=None)[source]

read a header information from either edi file or a list of lines containing header information.

Parameters:

header_list (list) – should be read from an .edi file or input as [‘key_01=value_01’, ‘key_02=value_02’]

Input header_list:
 
>>> h_list = ['lat=36.7898', 'lon=120.73532', 'elev=120.0', ...
>>>           'dataid=mt01']
>>> import mtpy.core.edi as mtedi
>>> header = mtedi.Header()
>>> header.read_header(h_list)
write_header(self, header_list=None, longitude_format='LON', latlon_format='dms')[source]

Write header information to a list of lines.

param header_list:
 

should be read from an .edi file or input as [‘key_01=value_01’, ‘key_02=value_02’]

type header_list:
 

list

param longitude_format:
 

whether to write longitude as LON or LONG. options are ‘LON’ or ‘LONG’, default ‘LON’

type longitude_format:
 

string

param latlon_format:
 

format of latitude and longitude in output edi, degrees minutes seconds (‘dms’) or decimal degrees (‘dd’)

type latlon_format:
 

string

returns header_lines:
 

list of lines containing header information will be of the form:

['>HEAD
‘,
‘ key_01=value_01
‘]
if None is input then reads from input .edi file or uses attribute information to write metadata.
class mtpy.core.edi.Information(edi_fn=None, edi_lines=None)[source]

Contain, read, and write info section of .edi file

not much to really do here, but just keep it in the same format that it is read in as, except if it is in phoenix format then split the two paragraphs up so they are sequential.

Methods

get_info_list(self) get a list of lines from the info section
read_info(self[, info_list]) read information section of the .edi file
write_info(self[, info_list]) write out information
get_info_list(self)[source]

get a list of lines from the info section

read_info(self, info_list=None)[source]

read information section of the .edi file

write_info(self, info_list=None)[source]

write out information

Module EDI_Collection

Description: To compute and encapsulate the properties of a set of EDI files

Author: fei.zhang@ga.gov.au

CreateDate: 2017-04-20

class mtpy.core.edi_collection.EdiCollection(edilist=None, mt_objs=None, outdir=None, ptol=0.05)[source]

A super class to encapsulate the properties pertinent to a set of EDI files

Parameters:
  • edilist – a list of edifiles with full path, for read-only
  • outdir – computed result to be stored in outdir
  • ptol – period tolerance considered as equal, default 0.05 means 5 percent

The ptol parameter controls what freqs/periods are grouped together: 10 percent may result more double counting of freq/period data than 5 pct. (eg: MT_Datasets/WPJ_EDI)

Methods

create_measurement_csv(self, dest_dir[, …]) create csv file from the data of EDI files: IMPEDANCE, APPARENT RESISTIVITIES AND PHASES see also utils/shapefiles_creator.py
create_mt_station_gdf(self[, outshpfile]) create station location geopandas dataframe, and output to shape file
create_phase_tensor_csv(self, dest_dir[, …]) create phase tensor ellipse and tipper properties.
create_phase_tensor_csv_with_image(\*args, …) Using PlotPhaseTensorMaps class to generate csv file of phase tensor attributes, etc.
display_on_basemap(self) display MT stations which are in stored in geopandas dataframe in a base map.
display_on_image(self) display/overlay the MT properties on a background geo-referenced map image
export_edi_files(self, dest_dir[, …]) export edi files.
get_bounding_box(self[, epsgcode]) compute bounding box
get_min_max_distance(self) get the min and max distance between all possible pairs of stations.
get_period_occurance(self, aper) For a given aperiod, compute its occurance frequencies among the stations/edi :param aper: a float value of the period :return:
get_periods_by_stats(self[, percentage]) check the presence of each period in all edi files, keep a list of periods which are at least percentage present :return: a list of periods which are present in at least percentage edi files
get_phase_tensor_tippers(self, period[, …]) For a given MT period (s) value, compute the phase tensor and tippers etc.
get_station_utmzones_stats(self) A simple method to find what UTM zones these (edi files) MT stations belong to are they in a single UTM zone, which corresponds to a unique EPSG code? or do they belong to multiple UTM zones?
get_stations_distances_stats(self) get the min max statistics of the distances between stations.
plot_stations(self[, savefile, showfig]) Visualise the geopandas df of MT stations
select_periods(self[, show, period_list, …]) Use edi_collection to analyse the whole set of EDI files
show_obj(self[, dest_dir]) test call object’s methods and show it’s properties
create_measurement_csv(self, dest_dir, period_list=None, interpolate=True)[source]

create csv file from the data of EDI files: IMPEDANCE, APPARENT RESISTIVITIES AND PHASES see also utils/shapefiles_creator.py

Parameters:
  • dest_dir – output directory
  • period_list – list of periods; default=None, in which data for all available frequencies are output
  • interpolate – Boolean to indicate whether to interpolate data onto given period_list
Returns:

csvfname

create_mt_station_gdf(self, outshpfile=None)[source]

create station location geopandas dataframe, and output to shape file

Parameters:outshpfile – output file
Returns:gdf
create_phase_tensor_csv(self, dest_dir, period_list=None, interpolate=True, file_name='phase_tensor.csv')[source]

create phase tensor ellipse and tipper properties. Implementation based on mtpy.utils.shapefiles_creator.ShapeFilesCreator.create_csv_files

Parameters:
  • dest_dir – output directory
  • period_list – list of periods; default=None, in which data for all available frequencies are output
  • interpolate – Boolean to indicate whether to interpolate data onto given period_list
  • file_name – output file name
Returns:

pt_dict

create_phase_tensor_csv_with_image(*args, **kwargs)[source]

Using PlotPhaseTensorMaps class to generate csv file of phase tensor attributes, etc. Only for comparison. This method is more expensive because it will create plot object first.

Returns:
display_on_basemap(self)[source]

display MT stations which are in stored in geopandas dataframe in a base map.

Returns:plot object
display_on_image(self)[source]

display/overlay the MT properties on a background geo-referenced map image

Returns:plot object
export_edi_files(self, dest_dir, period_list=None, interpolate=True, period_buffer=None, longitude_format='LON')[source]

export edi files. :param dest_dir: output directory :param period_list: list of periods; default=None, in which data for all available

frequencies are output
Parameters:
  • interpolate – Boolean to indicate whether to interpolate data onto given period_list; otherwise a period_list is obtained from get_periods_by_stats()
  • file_name – output file name
  • period_buffer – buffer so that interpolation doesn’t stretch too far over periods. Provide a float or integer factor, greater than which interpolation will not stretch. e.g. 1.5 means only interpolate to a maximum of 1.5 times each side of each frequency value
Returns:

get_bounding_box(self, epsgcode=None)[source]

compute bounding box

Returns:bounding box in given proj coord system
get_min_max_distance(self)[source]

get the min and max distance between all possible pairs of stations.

Returns:min_dist, max_dist
get_period_occurance(self, aper)[source]

For a given aperiod, compute its occurance frequencies among the stations/edi :param aper: a float value of the period :return:

get_periods_by_stats(self, percentage=10.0)[source]

check the presence of each period in all edi files, keep a list of periods which are at least percentage present :return: a list of periods which are present in at least percentage edi files

get_phase_tensor_tippers(self, period, interpolate=True)[source]

For a given MT period (s) value, compute the phase tensor and tippers etc.

Parameters:
  • period – MT_period (s)
  • interpolate – Boolean to indicate whether to interpolate on to the given period
Returns:

dictionary pt_dict_list

pt_dict keys [‘station’, ‘freq’, ‘lon’, ‘lat’, ‘phi_min’, ‘phi_max’, ‘azimuth’, ‘skew’, ‘n_skew’, ‘elliptic’,
‘tip_mag_re’, ‘tip_mag_im’, ‘tip_ang_re’, ‘tip_ang_im’]
get_station_utmzones_stats(self)[source]

A simple method to find what UTM zones these (edi files) MT stations belong to are they in a single UTM zone, which corresponds to a unique EPSG code? or do they belong to multiple UTM zones?

Returns:a_dict like {UTMZone:Number_of_MT_sites}
get_stations_distances_stats(self)[source]

get the min max statistics of the distances between stations. useful for determining the ellipses tipper sizes etc

Returns:dict={}
plot_stations(self, savefile=None, showfig=True)[source]

Visualise the geopandas df of MT stations

Parameters:
  • savefile
  • showfig
Returns:

select_periods(self, show=True, period_list=None, percentage=10.0)[source]

Use edi_collection to analyse the whole set of EDI files

Parameters:
  • show – True or false
  • period_list
  • percentage
Returns:

select_period_list

show_obj(self, dest_dir=None)[source]

test call object’s methods and show it’s properties

Returns:
mtpy.core.edi_collection.is_num_in_seq(anum, aseq, atol=0.0001)[source]

check if anum is in a sequence by a small tolerance

Parameters:
  • anum – a number to be checked
  • aseq – a sequence or a list of values
  • atol – absolute tolerance
Returns:

True | False

Module XML

Note

This module is written to align with the tools written by Anna Kelbert <akelbert@usgs.gov>

class mtpy.core.mt_xml.MT_XML(**kwargs)[source]

Class to read and write MT information from XML format. This tries to follow the format put forward by Anna Kelbert for archiving MT response data.

A configuration file can be read in that might make it easier to write multiple files for the same survey.

See also

mtpy.core.mt_xml.XML_Config

Attributes Description
Z object of type mtpy.core.z.Z
Tipper object of type mtpy.core.z.Tipper

Note

All other attributes are of the same name and of type XML_element, where attributes are name, value and attr. Attr contains any tag information. This is left this way so that mtpy.core.mt.MT can read in the information. Use mtpy.core.mt.MT for conversion between data formats.

Methods Description
read_cfg_file Read a configuration file in the format of XML_Config
read_xml_file Read an xml file
write_xml_file Write an xml file
Example:

:: >>> import mtpy.core.mt_xml as mtxml >>> x = mtxml.read_xml_file(r”/home/mt_data/mt01.xml”) >>> x.read_cfg_file(r”/home/mt_data/survey_xml.cfg”) >>> x.write_xml_file(r”/home/mt_data/xml/mt01.xml”)

Attributes:
Tipper

get Tipper information

Z

get z information

Methods

read_cfg_file(self[, cfg_fn]) Read in a cfg file making all key = value pairs attribures of XML_Config.
read_xml_file(self, xml_fn) read in an xml file and set attributes appropriately.
write_cfg_file(self[, cfg_fn]) Write out configuration file in the style of: parent.attribute = value
write_xml_file(self, xml_fn[, cfg_fn]) write xml from edi
Tipper

get Tipper information

Z

get z information

read_xml_file(self, xml_fn)[source]

read in an xml file and set attributes appropriately.

write_xml_file(self, xml_fn, cfg_fn=None)[source]

write xml from edi

exception mtpy.core.mt_xml.MT_XML_Error[source]
class mtpy.core.mt_xml.XML_Config(**kwargs)[source]

Class to deal with configuration files for xml.

Includes all the important information for the station and how data was processed.

Key Information includes:

Name Purpose
ProductID Station name
ExternalUrl External URL to link to data
Notes Any important information on station, data collection.
TimeSeriesArchived Information on Archiving time series including URL.
Image A location to an image of the station or the MT response.
  • ProductID –> station name
    • ExternalUrl –> external url to link to data
    • Notes –> any

Methods

read_cfg_file(self[, cfg_fn]) Read in a cfg file making all key = value pairs attribures of XML_Config.
write_cfg_file(self[, cfg_fn]) Write out configuration file in the style of: parent.attribute = value
read_cfg_file(self, cfg_fn=None)[source]

Read in a cfg file making all key = value pairs attribures of XML_Config. Being sure all new attributes are XML_element objects.

The assumed structure of the xml.cfg file is similar to:

``# XML Configuration File MTpy

Attachement.Description = Original file use to produce XML Attachment.Filename = None

Copyright.Citation.Authors = None Copyright.Citation.DOI = None Copyright.Citation.Journal = None Copyright.Citation.Title = None Copyright.Citation.Volume = None Copyright.Citation.Year = None

PeriodRange(max=0)(min=0) = None``

where the heirarchy of information is separated by a . and if the information has attribures they are in the name with (key=value) syntax.

write_cfg_file(self, cfg_fn=None)[source]

Write out configuration file in the style of: parent.attribute = value

class mtpy.core.mt_xml.XML_element(name, attr, value, **kwargs)[source]
Basically an ET element. The key components are
  • ‘name’ –> name of the element
  • ‘attr’ –> attribute information of the element
  • ‘value’ –> value of the element

Used the property function here to be sure that these 3 cannot be set through the common k.value = 10, just in case there are similar names in the xml file. This seemed to be the safest to avoid those cases.

Attributes:
attr
name
value

Module JFile

class mtpy.core.jfile.JFile(j_fn=None)[source]

be able to read and write a j-file

Methods

read_header(self[, j_lines]) Parsing the header lines of a j-file to extract processing information.
read_j_file(self[, j_fn]) read_j_file will read in a *.j file output by BIRRP (better than reading lots of *.<k>r<l>.rf files)
read_metadata(self[, j_lines, j_fn]) read in the metadata of the station, or information of station logistics like: lat, lon, elevation
read_header(self, j_lines=None)[source]

Parsing the header lines of a j-file to extract processing information.

Input: - j-file as list of lines (output of readlines())

Output: - Dictionary with all parameters found

read_j_file(self, j_fn=None)[source]

read_j_file will read in a *.j file output by BIRRP (better than reading lots of *.<k>r<l>.rf files)

Input: j-filename

Output: 4-tuple - periods : N-array - Z_array : 2-tuple - values and errors - tipper_array : 2-tuple - values and errors - processing_dict : parsed processing parameters from j-file header

read_metadata(self, j_lines=None, j_fn=None)[source]

read in the metadata of the station, or information of station logistics like: lat, lon, elevation

Not really needed for a birrp output since all values are nan’s