AIMBAT package

Submodules

algiccs.py

AIMBAT

AIMBAT (Automated and Interactive Measurement of Body wave Arrival Times) is an open-source software package for efficiently measuring teleseismic body wave arrival times for large seismic arrays (Lou et al., 2012). It is based on a widely used method called MCCC (Multi-Channel Cross-Correlation) developed by VanDecar and Crosson (1990). The package is automated in the sense of initially aligning seismograms for MCCC which is achieved by an ICCS (Iterative Cross Correlation and Stack) algorithm. Meanwhile, a graphical user interface is built to perform seismogram quality control interactively. Therefore, user processing time is reduced while valuable input from a user’s expertise is retained. As a byproduct, SAC (Goldstein et al., 2003) plotting and phase picking functionalities are replicated and enhanced.

algmccc.py

Python module for the MCCC (Multi-Channel Cross-Correlation) algorithm (VanDecar and Cross, 1990). Code transcribed from the original MCCC 3.0 fortran version using the same function names:

corread, corrmax, corrcff, correrr, corrmat, corrnow, corrwgt, corrite
copyright:Xiaoting Lou, John VanDecar
license:GNU General Public License, Version 3 (GPLv3) http://www.gnu.org/licenses/gpl.html
pysmo.aimbat.algmccc.WriteFileOriginal(mcpara, solist, solution, outvar, outcc, itmean)
pysmo.aimbat.algmccc.WriteFileWithDelay(mcpara, solist, solution, outvar, outcc, t0_times, delay_times, itmean)
pysmo.aimbat.algmccc.corrcff(ccmatrix)

Calculate mean and standard deviation of correlation coefficients.

pysmo.aimbat.algmccc.corrcff_fish(ccmatrix)

Calculate mean and standard deviation of correlation coefficients using Fisher’s transform:

z = 0.5 * ln((1+r)/(1-r))

Input: ccmatrix is the matrix of correlation coefficients. Output: ccmean is transformed back but not ccstd.

Problem: zero division if correlation coefficient is 1.

pysmo.aimbat.algmccc.corread(saclist, ipick, timewindow, taperwindow, tapertype)

Read data within timewindow+taperwindow (same length for each trace) for cross-correlation.

pysmo.aimbat.algmccc.correrr(dtmatrix, invmodel)

Calculate the rms misfit between cross correlated derived relative delay times and least-squares solution:

res_ij = dt_ij - (t_i - t_j)
pysmo.aimbat.algmccc.corrite(solist, mcpara, reftimes, solution, outvar, outcc)

Write output file, set output time picks.

pysmo.aimbat.algmccc.corrmat(windata, reftimes, mcpara)
Build matrices of cross-correlation derived relative delay times for least-squares solution.
A * t = dt

invmatrix A: sparse n*(n-1)/2+1 by n coefficient matrix including zero mean constraint invdata dt: cross-correlation derived relative delay times between every pair of stations invmodel t: optimized relative delay times for each station

pysmo.aimbat.algmccc.corrmax(datai, timei, dataj, timej, mcpara)
Calculate cross-correlation derived relative delay times by calling one of the xcorr functions.
dt_ij = t_i - t_j - tau_max

where t_i and t_j are initial time picks for the i-th and j-th traces, and tau_max is the time lag at maximum correlation

pysmo.aimbat.algmccc.corrnow(invmatrix, invdata)
Solve A * t = dt by lease-squares without weighting:
t = inv(A’A) * A’ * dt = 1/n * A’ * dt, where A’A = nI
pysmo.aimbat.algmccc.corrwgt(invmatrix, invdata, ccmatrix, resmatrix, wgtscheme='correlation', exwt=1000.0)
Solve A * t = dt by weighted least-squares:
t = inv(A’WA) * A’ * W * dt

W: n*(n-1)/2+1 by n*(n-1)/2+1 diagonal weighting matrix exwt: weight for the extra equation of zero-meaning constraint

pysmo.aimbat.algmccc.eventListName(evlist='event.list', phase='S', isol='PDE')

Read evlist (either event.list file or gsac.event list) for hypocenter and origin time. Create output filename.

pysmo.aimbat.algmccc.getOptions()

Parse arguments and options from command line. No default value is given here because it will override values from configuration file.

pysmo.aimbat.algmccc.getParams(gsac, mcpara, opts=None)

Get parameters for running MCCC. Hierarchy: default config < gsac < .mcccrc < command line options

pysmo.aimbat.algmccc.getWindow(stkdh, ipick, twhdrs, taperwidth=0.1)

Get timewindow and taperwindow from gsac.stkdh

pysmo.aimbat.algmccc.main()
pysmo.aimbat.algmccc.mccc(gsac, mcpara)

Run MCCC.

pysmo.aimbat.algmccc.rcdef()

Default values for the inherited .mcccrc file of MCCC 3.0

pysmo.aimbat.algmccc.rcread(rcfile='.mcccrc')

Read .mcccrc file and return ipick, time window and taper window.

pysmo.aimbat.algmccc.rcwrite(ipick, timewindow, taperwindow, rcfile='.mcccrc')

Write to .mcccrc file. Convert timewindow –> window, inset, taper.

filtering.py

pysmo.aimbat.filtering.filtering_time_freq(originalTime, originalSignalTime, delta, filterType, highFreq, lowFreq, order, runReversePass=False)
pysmo.aimbat.filtering.filtering_time_signal(originalSignalTime, delta, lowFreq, highFreq, filterType, order, MULTIPLE, runReversePass=False)
pysmo.aimbat.filtering.get_filter_params(delta, lowFreq, highFreq, filterType, order, MULTIPLE=3000)
pysmo.aimbat.filtering.time_to_freq(originalTime, originalSignalTime, delta)

pickphase.py

Python module for plot and pick phase (SAC PPK) on seismograms in one axes.

Differences from plotphase.py:
  • User interaction: set time picks and time window
  • Plot: always plot time picks
  • Plot: always use integer numbers (plot within +/-0.5) as ybases,
    but not dist/az/baz (even when sorted by d/a/b)
  • Plot: can plot seismograms in multiple pages (page navigation).
  • Normalization: can normalize within time window
Keyboard and mouse actions:
  • Click mouse to select a span to zoom in seismograms.
  • Press ‘z’ to go back to last window span.
  • Press ‘w’ to save the current xlimit as time window.
  • Press ‘t[0-9]’ to set time picks like SAC PPK.
Program structure:
PickPhaseMenu
||

PickPhase Button Front + Button Prev + Button Next + Button Last + Button Save + Button Quit

copyright:Xiaoting Lou
license:GNU General Public License, Version 3 (GPLv3) http://www.gnu.org/licenses/gpl.html
class pysmo.aimbat.pickphase.PickPhase(sacdh, opts, axpp, ybase, color='b', linew=1, alpha=1)

Bases: object

Plot one single seismogram with given attributes. See self.on_press for options on setting time picks and time window.

changeBase(newbase)

Change ybase of a seismogram.

changeColor()

Change color of a seismogram based on selection status.

connect()
disconnect()
disconnectPick()
labelStation()

label the seismogram with file name or net.sta

makeTime()

Create array x as time series and get reference time.

on_pick(event)

Click a seismogram to show file name.

on_press(event)

Key press event. Valid only if axpp contains event (within 0.5 from ybase). Options: ——– (1) t + digits 0-9: set a time pick in SAC header. (2) w: set the current xlim() as time window.

plotPicks()

Plot time picks at ybase +/- 0.5

plotWave()

Plot wiggled or filled waveform, which is normalized (if not stacking) and shifted to ybase. Fill both plus and negative side of signal but with different transparency.

If opts.fill == 0: no fill. If opts.fill > 0: alpha of negative side is a quarter of plus side. If opts.fill < 0: alpha of plus side is a quarter of negative side.
plotWindow()

Plot time window (xmin,xmax) with color fill.

resetWindow()

Reset time window when a span is selected.

updateY(xxlim)

Update ynorm for wave wiggle from given xlim.

class pysmo.aimbat.pickphase.PickPhaseMenu(gsac, opts, axs)

Bases: object

Plot a group of seismogram gathers. Set up axes attributes. Create Button Save to save SAC headers to files.

connect()
disconnect(canvas)
finish()
fron(event)
getXLimit()

Get x limit (relative to reference time)

getYLimit()

Get y limit

initIndex()

Initialize indices for page navigation.

labelSelection()

Label selection status with transform (transAxes, transData).

last(event)
next(event)
on_select(xmin, xmax)

Mouse event: select span.

on_zoom(event)

Zoom back to previous xlim when event is in event.inaxes.

plotPicks()
plotSeis()
plotSpan()

Create a SpanSelector for zoom in and zoom out.

plotWave()

Plot waveforms for this page.

prev(event)
quit(event)
replot(ipage)

Finish plotting of current page and move to prev/next.

setLabels()

Set axes labels and page label

setLimits()

Set axes limits

shdo(event)
shfp(event)
shod(event)
zoba(event)
pysmo.aimbat.pickphase.getAxes(opts)

Get axes for plotting

pysmo.aimbat.pickphase.getDataOpts()

Get SAC Data and Options

pysmo.aimbat.pickphase.getOptions()

Parse arguments and options.

pysmo.aimbat.pickphase.main()
pysmo.aimbat.pickphase.sacppk_standalone()
pysmo.aimbat.pickphase.sortSeis(gsac, opts)

Sort seismograms by file indices, quality factors, or a given header

plotphase.py

Python module for plotting multiple seismograms in one axes.
Plot only: no data/attributes of SAC files are changed.
Keyboard and mouse actions:
Click mouse to select a span to zoom in seismograms. Press the ‘z’ key to go back to last window span.
Requried options:
One of azim_on, bazim_on, dist_on, index_on, zero_on must be chosen. Default is index_on. Correspong to: paz, pbaz, prs, p1 and p2. Default: p1

These functions can be acheived by scripts created for each mode.

Optional options:
pick_on, color_on, stack_on, std_on
Program structure:
SingleSeisGather
||
SingleSeis + baseIndex + baseZero + baseDist + baseAzim + baseBAzim
| | | |

V V V V V

sacp1 sacp2 sacprs sacpaz sacpbaz

copyright:Xiaoting Lou
license:GNU General Public License, Version 3 (GPLv3) http://www.gnu.org/licenses/gpl.html
class pysmo.aimbat.plotphase.SingleSeis(sacdh, opts, axss, ybase, color='b', linew=1, alpha=1)

Bases: object

Plot a single seismogram with given attributes.

connect()
disconnect()
makeTime()

Create array x as time series and get reference time.

onpick(event)
plotPicks()

Plot time picks. Not called by default. Only works for baseIndex mode because axvline is not used to plot time picks.

plotWave()

Plot wiggled or filled waveform, which is normalized (if not stacking) and shifted to ybase. Fill both plus and negative side of signal but with different transparency.

If opts.fill == 0: no fill. If opts.fill > 0: alpha of negative side is a quarter of plus side. If opts.fill < 0: alpha of plus side is a quarter of negative side.
class pysmo.aimbat.plotphase.SingleSeisGather(saclist, opts, axss)

Bases: object

Plot a group of seismograms.

baseAzim()

Set baseline of seismograms as azimuth.

baseBAzim()

Set baseline of seismograms as back azimuth.

baseDist()

Set baseline of seismograms as epicentral distance in degree.

baseDistkm()

Set baseline of seismograms as epicentral distance in km.

baseIndex()

Set baseline of seismograms as file indices.

baseZero()

Set baseline of seismograms as zeros (stack all). Do not normalize seismogram by setting opts.ynorm<0.

connect()
disconnect()
getAzim()

Get azimuth as ybases for waveforms.

getBAzim()

Get back azimuth as ybases for waveforms.

getDist(degree=True)

Get epicentral distances in degree/km as ybases for waveforms.

getIndex()

Get file indices as ybases for waveforms.

getPlot()

Get plotting attributes

getXLimit()

Get x limit (relative to reference time)

getYLimit()

Get y limit

getZero()

Get zeros as ybases for waveforms.

labelStation()

Label stations at y axis on the right. The xcoords of the transform are axes, and the yscoords are data.

on_zoom(event)

Zoom back to previous xlim when event is in event.inaxes.

plotPicks()
plotSeis()

Plot wiggles or filled waveforms.

plotSpan()

Create a SpanSelector on axss.

plotStack()

Calculate mean stack from all traces and plot it. No taper window is added in.

pysmo.aimbat.plotphase.dopts(opts)

Default options

pysmo.aimbat.plotphase.getAxes(opts)

Get axes for plotting

pysmo.aimbat.plotphase.getDataOpts()

Get SAC Data and Options

pysmo.aimbat.plotphase.getOptions()

Parse arguments and options.

pysmo.aimbat.plotphase.main()
pysmo.aimbat.plotphase.sacp1(saclist, opts, axss)

SAC P1 style of plotting.

pysmo.aimbat.plotphase.sacp1_standalone()
pysmo.aimbat.plotphase.sacp2(saclist, opts, axss)

SAC P2 style of plotting.

pysmo.aimbat.plotphase.sacp2_standalone()
pysmo.aimbat.plotphase.sacpaz(saclist, opts, axss)

SAC plotting along azimuth.

pysmo.aimbat.plotphase.sacpaz_standalone()
pysmo.aimbat.plotphase.sacpbaz(saclist, opts, axss)

SAC plotting along backazimuth.

pysmo.aimbat.plotphase.sacpbaz_standalone()
pysmo.aimbat.plotphase.sacplot_standalone()
pysmo.aimbat.plotphase.sacprs(saclist, opts, axss)

SAC PRS style of plotting: record section.

pysmo.aimbat.plotphase.sacprs_standalone()
pysmo.aimbat.plotphase.splitAxesH(fig, rect=[0.1, 0.1, 0.6, 0.6], n=2, hspace=0, axshare=False)

Split an axes to multiple (n,1,i) horizontal axes. Share x-axis if axshare is True.

plotstations.py

class pysmo.aimbat.plotstations.PlotStations(plotname, gsac)

Bases: object

bounding_rectangle()
plot_deleted_stations(axes_handle)
plot_selected_stations(axes_handle)
plot_selected_stations_color_delay_times(axes_handle, delay_times)
plot_stations(gsac)
show_station_name(event)

plotutils.py

Python module for plotting seismograms:
functions for axes and legend control, SpanSelector, and multiple-page navigation.
copyright:Xiaoting Lou
license:GNU General Public License, Version 3 (GPLv3) http://www.gnu.org/licenses/gpl.html
class pysmo.aimbat.plotutils.TimeSelector(ax, onselect, direction, minspan=None, useblit=False, rectprops=None, onmove_callback=None, span_stays=False, button=None)

Bases: matplotlib.widgets.SpanSelector

To disable SpanSelector when pan, zoom or other interactive/navigation modes are active. Also disable it when event is out of axes, which is needed to avoid error interfering with pick_event.

ignore(event)

return True if event should be ignored

pysmo.aimbat.plotutils.axLimit(minmax, w=0.05)

Calculate axis limit with white space (default 5%) from given min/max values.

pysmo.aimbat.plotutils.dataNorm(d, w=0.05)

Calculate normalization factor for d, which can be multi-dimensional arrays. Extra white space is added.

pysmo.aimbat.plotutils.getAxes(opts)

Get axes for pickphase

pysmo.aimbat.plotutils.indexBaseTick(na, nb, pagesize, pna)

Indexing for page navigation with two lists of length na and nb.

Example

list b (nb=5) list a (na=11) [ 0, 1, 2, 3, 4] [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10] <– yindex [ 5, 4, 3, 2, 1] [-1,-2,-3,-4,-5,-6,-7,-8,-9,-10,-11] <– ybases [-5,-4,-3,-2,-1] [1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11] <– yticks

–page -1] [—-page 0—-] [—page 1—] [—page 2—
[pnb=2] [pna=3 ]

yindex for na and nb: {-1: [[], [0, 1, 2]],

0: [[0, 1, 2], [3, 4]], 1: [[3, 4, 5, 6, 7], []], 2: [[8, 9, 10, 11, 12], []]}

yybase: {-1: [[], [5, 4, 3]],

0: [[-1, -2, -3], [2, 1]], 1: [[-4, -5, -6, -7, -8], []], 2: [[-9, -10, -11, -12, -13], []]}

yticks: {-1: [[], [-5, -4, -3]],

0: [[1, 2, 3], [-2, -1]], 1: [[4, 5, 6, 7, 8], []], 2: [[9, 10, 11, 12, 13], []]}
pysmo.aimbat.plotutils.pickLegend(ax, npick, pickcolors, pickstyles, left=True)

Plot only legend box for time picks.

pysmo.aimbat.plotutils.plotDelay(stalos, stalas, dtimes, opts)

prepdata.py

Python module for preparing time and data arrays (original and filtered in memory) for plotting.

copyright:Xiaoting Lou
license:GNU General Public License, Version 3 (GPLv3) http://www.gnu.org/licenses/gpl.html
pysmo.aimbat.prepdata.axLimit(minmax, w=0.05)

Calculate axis limit with white space (default 5%) from given min/max values.

pysmo.aimbat.prepdata.dataNorm(d, w=0.05)

Calculate normalization factor for d, which can be multi-dimensional arrays. Extra white space is added.

pysmo.aimbat.prepdata.dataNormWindow(d, t, twindow)

Calculate normalization factor in a time window.

pysmo.aimbat.prepdata.findPhase(filename)

Find phase (P or S) from component info (BH?) in file name

pysmo.aimbat.prepdata.getFilterPara(sacdh, pppara)

Get default filter parameters from ttdefaults.conf. Override defaults if already set in SAC file

pysmo.aimbat.prepdata.paraDataOpts(opts, ifiles)

Common parameters, data and options

pysmo.aimbat.prepdata.prepData(gsac, opts)

Prepare data for plotting

pysmo.aimbat.prepdata.prepStack(opts)

Prep stack from existing sacfile opts.fstack

pysmo.aimbat.prepdata.sacDataNorm(sacdh, opts)

get data normalization factor one seismogram

pysmo.aimbat.prepdata.seisApplyFilter(saclist, filtParas)

Filter seismograms by butterworth filter

pysmo.aimbat.prepdata.seisDataBaseline(gsac)

Create plotting baselines for each seismogram in selected and deselected sac lists Example

delist (n=5) selist a(na=11)

[ -5, -4, -3, -2, -1] [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10] <– yindex [ 5, 4, 3, 2, 1] [0, -1,-2,-3,-4,-5,-6,-7,-8, -9,-10] <– ybases

pysmo.aimbat.prepdata.seisDataNorm(saclist, opts)

get data normalization factor each seismogram

pysmo.aimbat.prepdata.seisSort(gsac, opts)

Sort seismograms by file indices, quality factors, time difference, or a given header.

pysmo.aimbat.prepdata.seisTimeData(saclist)

Create time and data (original and in memory) arrays

pysmo.aimbat.prepdata.seisTimeRefr(saclist, opts)

get reference time for each seismogram

pysmo.aimbat.prepdata.seisTimeWindow(saclist, twhdrs)

get time window for each seismogram

pysmo.aimbat.prepdata.seisUnApplyFilter(saclist)

Filter seismograms by butterworth filter

pysmo.aimbat.prepdata.seisWave(sacdh)

Calculate waveform X and Y for plot

pysmo.aimbat.prepdata.setFilterPara(sacdh, pppara, filterParameters)

Set filter parameters dict to sacdh.

prepplot.py

qualctrl.py

Python module for interactively measuring body wave travel times and quality control. Used by ttpick.py

PickPhaseMenuMore
||

PickPhaseMenu + Buttons: ICCS-A Sync ICCS-B MCCC SACP2

copyright:Xiaoting Lou
license:GNU General Public License, Version 3 (GPLv3) http://www.gnu.org/licenses/gpl.html
class pysmo.aimbat.qualctrl.PickPhaseMenuMore(gsac, opts, axs)

Bases: object

Pick phase for multiple seismograms and array stack Button: Sync

addEarthquakeInfo()

Set Earthquake Info * Magnitude * Location (Lat and Long) * Depth

applyFilter(event)
ccStack()

Call iccs.ccWeightStack. Change reference time pick to the input pick before ICCS and to the output pick afterwards.

ccff(event)

Run iccs with time window from final stack. Time picks: hdrfin, hdrfin.

ccim(event)
connect()

Connect button events.

disconnect()

Disconnect button events.

dismiss_sort(event)

Dismiss the sorting selection popup Window

filter_connect()
filtering(event)
getBandpassFreq(event)
getBandtype(event)
getButterOrder(event)
getFilterAxes()
getHighFreq(event)
getLowFreq(event)
getPicks()

Get time picks of stack

getReversePassOption(event)
getSortAxes()
getWindow(hdr)

Get time window twcorr (relative to hdr) from array stack, which is from last run.

initPlot()

Plot waveforms

mccc(event)

Run mccc.py

modifyFilterTextLabels()
on_zoom(event)

Zoom back to previous xlim when event is in event.inaxes.

plot2(event)

Plot P2 stack of seismograms for defined time picks (ichdrs + wpick).

plotSpan()

Span for array stack

plotStack()

Plot array stack and span

plot_stations(event)
replot()

Replot seismograms and array stack after running iccs.

replot_seismograms()

Replot seismograms and array stack after running iccs.

setLabels()

Set plot attributes

sort_connect()

write the position for the buttons into self

sort_disconnect()
sort_file(event)
sort_haz(event)
sort_hb(event)
sort_hbaz(event)
sort_hdelta(event)
sort_hdist(event)
sort_he(event)
sort_hgcarc(event)
sort_hkstnm(event)
sort_hnpts(event)
sort_hstla(event)
sort_hstlo(event)
sort_qall(event)
sort_qccc(event)
sort_qcoh(event)
sort_qsnr(event)
sorting(event)

Sort the seismograms in particular order

spreadButter()
summarize_sort()
sync(event)

Sync final time pick and time window from array stack to each trace and update current page.

syncPick()

Sync final time pick hdrfin from array stack to all traces.

syncWind()

Sync time window relative to hdrfin from array stack to all traces. Times saved to twhdrs are alway absolute.

unapplyFilter(event)
pysmo.aimbat.qualctrl.getAxes(opts)

Get axes for plotting

pysmo.aimbat.qualctrl.getDataOpts()

Get SAC Data and Options

pysmo.aimbat.qualctrl.getOptions()

Parse arguments and options.

pysmo.aimbat.qualctrl.main()
pysmo.aimbat.qualctrl.sortSeis(gsac, opts)

Sort seismograms by file indices, quality factors, time difference, or a given header.

qualsort.py

Python module for selecting and sorting seismograms by quality factors and other header variables.

copyright:Xiaoting Lou
license:GNU General Public License, Version 3 (GPLv3) http://www.gnu.org/licenses/gpl.html
pysmo.aimbat.qualsort.getOptions()

Parse arguments and options.

pysmo.aimbat.qualsort.hdrtype(sacdh, hdr)

Indentify type of a header

pysmo.aimbat.qualsort.initQual(saclist, hdrsel, qheaders)

Set initial values for selection status and quality factors.

pysmo.aimbat.qualsort.seleSeis(saclist)

Select seismograms. Return sacdh lists of selected and deleted seismograms.

selelist: selected seismograms delelist: deleted seismograms, user doe snot want them

pysmo.aimbat.qualsort.sortQual(saclist, qheaders, qweights, increase=True)

Sort quality factors by weighted averaging. Return sorted sacdh list and means of quality factors.

pysmo.aimbat.qualsort.sortSeisHeader(saclist, hdr, increase=True)

Sort saclist by a header value.

pysmo.aimbat.qualsort.sortSeisHeaderDiff(saclist, hdr0, hdr1, increase=True)

Sort saclist by header value difference (hdr1-hdr0) in increase/decrease order. Limited to t_n, user_n and kuser_n headers.

pysmo.aimbat.qualsort.sortSeisQual(saclist, qheaders, qweights, qfactors, increase=True)

Select and sort seismograms by quality factors.

sacpickle.py

Module sacpickle.py

class pysmo.aimbat.sacpickle.SacDataHdrs(ifile, delta=-1)

Bases: object

Class for individual SAC file’s data and headers.

gethdr(hdr)

Read a header variable (t_n, user_n, or kuser_n).

resampleData(delta)
savesac()

Save all data and header variables to an existing or new sacfile.

sethdr(hdr, val)

Write a header variable (t_n, user_n, or kuser_n).

sethdrs(sacobj)

Write SAC headers (t_n, user_n, and kuser_n) in python obj to SAC obj.

writeHdrs()

Write SAC headers (t_n, user_n, and kuser_n) in python obj to existing SAC file.

class pysmo.aimbat.sacpickle.SacGroup(ifiles, delta=-1)

Bases: object

Read a group of SAC files’ headers and data to python objects in memory. Get event information.

resampleData(delta)

resample data of all sacdh

pysmo.aimbat.sacpickle.date2jul(year, mon, day)

date –> julian day

pysmo.aimbat.sacpickle.fileZipMode(ifilename)

Determine if an input file is SAC, pickle or compressed pickle file

pysmo.aimbat.sacpickle.getOptions()

Parse arguments and options.

pysmo.aimbat.sacpickle.jul2date(year, jday)

julian day –> date

pysmo.aimbat.sacpickle.loadData(ifiles, opts, para)

Load data either from SAC files or (gz/bz2 compressed) pickle file. Get sampling rate from command line option or default config file. Resample data if a positive sample rate is given. Output file type is the same as input file (filemode and zipmode do not change). If filemode == ‘sac’: zipmode = None If filemode == ‘pkl’: zipmode = None/bz2/gz

pysmo.aimbat.sacpickle.main()
pysmo.aimbat.sacpickle.obj2sac(gsac)

Save headers in python objects to SAC files.

pysmo.aimbat.sacpickle.pkl2sac(pkfile, zipmode)

Save headers in python pickle to SAC files.

pysmo.aimbat.sacpickle.readPickle(picklefile, zipmode=None)

Read compressed pickle file to python objects.

pysmo.aimbat.sacpickle.resampleSeis(data, deltaold, delta)

Resample data of a seismogram if given a different positive delta

pysmo.aimbat.sacpickle.sac2obj(ifiles, delta=-1)

Convert SAC files to python objects.

pysmo.aimbat.sacpickle.sac2pkl(ifiles, pkfile='sac.pkl', delta=-1, zipmode='gz')

Convert SAC files to python pickle files.

pysmo.aimbat.sacpickle.saveData(gsac, opts)

Save pickle or sac files.

pysmo.aimbat.sacpickle.taper(data, taperwidth=0.1, tapertype='hanning')

Apply a symmetric taper to each end of data. http://www.iris.edu/software/sac/commands/taper.html Default width: 0.1/2=0.05 on each end.

pysmo.aimbat.sacpickle.taperWindow(timewindow, taperwidth=0.1)

Calculate length of taper window from time window so that: taperwidth = (taperwindow)/(taperwindow+timewindow)

pysmo.aimbat.sacpickle.windowData(saclist, nstart, ntotal, taperwidth, tapertype='hanning', datatype='data')

Cut data within a time window using given indices. Pad dat with zero if not enough sample. Use sacdh.data or sacdh.datamem based on data type

pysmo.aimbat.sacpickle.windowIndex(saclist, reftimes, timewindow=(-5.0, 5.0), taperwindow=1.0)

Calculate indices for cutting data at a time window and taper window. Indices nstart and notal bound the entire window, which is sum of time window and taper window. Only the taper window part of data is tapered: __———–__ The original MCCC code defines taper window differently: ____——-____ :param saclist: :type saclist: list of sacdh :param reftimes: :type reftimes: list of reference times for the time window :param timewindow: :type timewindow: relative time window to cut data :param taperwindow: :type taperwindow: length of taper window :param nstart: :type nstart: index of first sample point of datacut in each sacdh.data :param ntotal: :type ntotal: length of data within the time window

pysmo.aimbat.sacpickle.windowTime(saclist, nstart, ntotal, taperwidth, tapertype='hanning')

Cut time within a time window using given indices.

pysmo.aimbat.sacpickle.windowTimeData(saclist, nstart, ntotal, taperwidth, tapertype='hanning')

Cut part of the time and data based on given indices.

pysmo.aimbat.sacpickle.writePickle(d, picklefile, zipmode=None)

Write python objects to pickle file and compress with highest protocal (binary) if zipmode is not None.

pysmo.aimbat.sacpickle.zipFile(zipmode='gz')

Return file compress method: bz2 or gz.

stationmapping.py

Python module for mapping station with the help of GMT scripts.

copyright:Arnav Sankaran
license:GNU General Public License, Version 3 (GPLv3) http://www.gnu.org/licenses/gpl.html
class pysmo.aimbat.stationmapping.StationMapper(sacgroup)

Bases: object

extractData()
plotData(selectedPath, deselectedPath)
start()
writeToFile()

ttconfig.py

Module: ttconfig.py

class pysmo.aimbat.ttconfig.CCConfig

Bases: object

Class for ICCS configuration.

class pysmo.aimbat.ttconfig.MCConfig

Bases: object

Class for MCCC configuration.

class pysmo.aimbat.ttconfig.PPConfig

Bases: object

Class for SAC PlotPhase and PickPhase configurations.

class pysmo.aimbat.ttconfig.QCConfig

Bases: object

Class for QCTRL configuration.

pysmo.aimbat.ttconfig.getDefaults()

Get default parameters from a configuration file: ttdefaults.conf. The file is searched in this order:

  1. the current working directory
  2. home directory
  3. environment variable TTCONFIG
  4. the same directory as this program.
pysmo.aimbat.ttconfig.getParser()

Parse command line arguments and options.

ttguiqt.py

xcorr.py

Python module to calculate cross-correlation function of two time series with the same length.
c(k) = sum_i x(i) * y(i+k)
where k is the time shift of y relative to x.

Delay time, correlation coefficient, and polarity at maximum correlation are returned.

Added suport to correlation polarity to allow negative correlation maximum. Output ccpol=1 if positive or ccpol=-1 if negative. xlou 03/2011

copyright:Xiaoting Lou
license:GNU General Public License, Version 3 (GPLv3) http://www.gnu.org/licenses/gpl.html
pysmo.aimbat.xcorr.xcorr_fast(x, y, shift=10)

Fast cross-correlation of two time series of the same length. One level of coarse shift by downsampling the signal.

pysmo.aimbat.xcorr.xcorr_fast_polarity(x, y, shift=10)

Fast cross-correlation of two time series of the same length. One level of coarse shift by downsampling the signal. Do not correct polarity

pysmo.aimbat.xcorr.xcorr_faster(x, y, shift=10)

Faster cross-correlation only for time lags around zero.

pysmo.aimbat.xcorr.xcorr_faster_polarity(x, y, shift=10)

Faster cross-correlation only for time lags around zero. Do not correct polarity

pysmo.aimbat.xcorr.xcorr_full(x, y, shift=1)

Cross-correlation of two 1-D arrays using ‘full’ mode. Argument shift=1 is here only in order to make the same number of arguments for all xcorr functions.

pysmo.aimbat.xcorr.xcorr_full_polarity(x, y, shift=1)

Cross-correlation of two 1-D arrays using ‘full’ mode. Argument shift=1 is here only in order to make the same number of arguments for all xcorr functions. Not correct the polarity

pysmo.aimbat.xcorr.xcorr_same(x, y)

Cross-correlation of two 1-D arrays using ‘same’ mode.

pysmo.aimbat.xcorr.xcorr_select(x, y, lags)

Cross-correlation of two time series of the same length for selected lag(shift) times

pysmo.aimbat.xcorr.xcorr_select_polarity(x, y, lags)

Cross-correlation of two time series of the same length for selected lag(shift) times Do not correct polarity