pyatoa.utils.srcrcv

Functions used in determining information related to sources and receivers, or their corresponding representations in ObsPy

Module Contents

Functions

lonlat_utm(lon_or_x, lat_or_y[, utm_zone, inverse])

Convert latitude and longitude coordinates to UTM projection using PyProj

utm_zone_from_lat_lon(lat, lon)

Calculate the UTM zone longitude value using quick maffs.

gcd_and_baz(event, sta)

Calculate great circle distance and backazimuth values for a given

merge_inventories(inv_a, inv_b)

Adding inventories together duplicates network and station codes, which is

seismogram_length(distance_km[, slow_wavespeed_km_s, ...])

Dynamically determine the length of the seismogram based on source-receiver

sort_by_backazimuth(ds[, clockwise])

Its illustrative to show misfit information for an event, sorted by

event_by_distance(cat[, filter_type, filter_bounds, ...])

Sort through an obspy catalog by interevent distance. If we have a lot of

pyatoa.utils.srcrcv.lonlat_utm(lon_or_x, lat_or_y, utm_zone=None, inverse=False)[source]

Convert latitude and longitude coordinates to UTM projection using PyProj

Parameters:
  • lon_or_x (float or int) – longitude value in WGS84 or X in UTM-‘zone’ projection

  • lat_or_y (float or int) – latude value in WGS84 or Y in UTM-‘zone’ projection

  • utm_zone (int) – UTM zone for conversion from WGS84

  • inverse (bool) – if inverse == False, latlon => UTM, vice versa.

Return type:

tuple (float, float)

Returns:

(x in UTM or longitude in WGS84, y in UTM or latitude in WGS84)

pyatoa.utils.srcrcv.utm_zone_from_lat_lon(lat, lon)[source]

Calculate the UTM zone longitude value using quick maffs. Get the sign of the UTM zone based on the latitude value.

Parameters:
  • lat (float) – latitude coordinate in degrees

  • lon (float) – longitude coordinate in degrees

Return type:

int

Returns:

UTM zone number

pyatoa.utils.srcrcv.gcd_and_baz(event, sta)[source]

Calculate great circle distance and backazimuth values for a given station and event configuration

Parameters:
  • event (obspy.core.event.Event) – event object

  • sta (obspy.core.inventory.station.Station) – station object

Return type:

tuple (float, float)

Returns:

(great circle distance in km, backazimuth in degrees)

pyatoa.utils.srcrcv.merge_inventories(inv_a, inv_b)[source]

Adding inventories together duplicates network and station codes, which is kind of annoying for looping. This function will add two inventories together while minimizing the amount of redundant networks, stations, channels inside the merged inventory.

Parameters:
  • inv_a (obspy.core.inventory.Inventory) – inventory to merge into, will be returned

  • inv_b (obspy.core.inventory.Inventory) – inventory to merge into inv_a

Return type:

obspy.core.inventory.Inventory

Returns:

merged inventories

pyatoa.utils.srcrcv.seismogram_length(distance_km, slow_wavespeed_km_s=2, binsize=50, minimum_length=100)[source]

Dynamically determine the length of the seismogram based on source-receiver distance. Bin into lengths to keep some uniformity in window lengths

Parameters:
  • distance_km (float) – source-receiver distance in km

  • slow_wavespeed_km_s (int) – slowest wavespeed in model, in km/s

  • binsize (int) – bin size for rounding the length to the nearest value

  • minimum_length (int) – the shortest a seismogram can be

Return type:

int

Returns:

expected seismogram length

pyatoa.utils.srcrcv.sort_by_backazimuth(ds, clockwise=True)[source]

Its illustrative to show misfit information for an event, sorted by backazimuth. Stations with misfit information are generally sorted alphabetically, so this function just calcualtes backazimuth and returns a sorted list of station names. Can go clockwise or counter clockwise, starting from 0 degrees.

Parameters:
  • ds (pyasdf.ASDFDataSet()) – dataset containing event and station information

  • clockwise (bool) – False = counter clockwise

Return type:

list

Returns:

list of stations in order from 0deg to 360deg in direction

pyatoa.utils.srcrcv.event_by_distance(cat, filter_type=False, filter_bounds=None, random=False)[source]

Sort through an obspy catalog by interevent distance. If we have a lot of events in a catalog, it’s best to spatially vary them such that we don’t redundantly oversample a spatial region. Returns an index list for events that are most distant from one another, without repeating any used events.

Catalog filter parameters can be found here: https://docs.obspy.org/packages/autogen/obspy.core.event.Catalog.filter.html

Parameters:
  • cat (obspy.event.Catalog) – catalog to sort through

  • filter_type (str) – filter to be passed to the Catalog filter

  • filter_bounds (list of floats) – (min filter bound, max filter bound)

  • random (bool) – randomly determined starting point

Return type:

obspy.event.Catalog

Returns:

filtered catalog object