pysep.utils.io#

Read utilities for Pysep

Module Contents#

Functions#

read_yaml(fid)

Read a YAML file and return a dictionary

read_sem(fid[, origintime, source, stations, ...])

Specfem3D outputs seismograms to ASCII (.sem? or .sem.ascii) files.

read_sem_cartesian(fid, source, stations[, location, ...])

Specfem2D and Specfem3D may have domains defined in a Cartesian coordinate

read_events_plus(fid, format, **kwargs)

Addition to the base ObsPy.read_events() function that, in addition to the

read_specfem3d_cmtsolution_cartesian(path_to_cmtsolution)

Create a barebones ObsPy Event object from a SPECFEM3D CMTSOLUTION file with

read_specfem2d_source(path_to_source[, origin_time])

Create a barebones ObsPy Event object from a SPECFEM2D Source file, which

read_forcesolution(path_to_forcesolution[, origin_time])

Create a barebones Source object from a FORCESOLUTION Source file,

_get_resource_id(name, res_type[, tag])

Helper function to create consistent resource ids, from ObsPy. Used to

read_stations(path_to_stations)

Convert a SPECFEM STATIONS file into an ObsPy Inventory object.

read_event_file(fid)

Read an event input file which is just a text file that should contain

write_sem(st, unit[, path, time_offset])

Write an ObsPy Stream in the two-column ASCII format that Specfem uses

write_pysep_station_file(inv, event[, fid, ...])

Write a list of station codes, distances, etc. useful for understanding

write_stations_file(inv[, fid, order_by, ...])

Write a SPECFEM STATIONS file given an ObsPy inventory object

write_cat_to_event_list(cat[, fid_out])

Writes out an event Catalog (ObsPy object) information to an ASCII file that

pysep.utils.io.read_yaml(fid)[source]#

Read a YAML file and return a dictionary

Parameters:

fid (str) – YAML file to read from

Return type:

dict

Returns:

YAML keys and variables in a dictionary

pysep.utils.io.read_sem(fid, origintime=None, source=None, stations=None, location='', precision=4, source_format='CMTSOLUTION')[source]#

Specfem3D outputs seismograms to ASCII (.sem? or .sem.ascii) files. Converts SPECFEM synthetics into ObsPy Stream objects with the correct header information. If source and stations files are also provided, PySEP will write appropriate SAC headers to the underlying data.

Parameters:
  • fid (str) – path of the given ascii file

  • origintime (obspy.UTCDateTime) – UTCDatetime object for the origintime of the event

  • source (str) – optional SPECFEM source file (e.g., CMTSOLUTION, SOURCE) defining the event which generated the synthetics. Used to grab event information and append as SAC headers to the ObsPy Stream

  • stations (str) – optional STATIONS file defining the station locations for the SPECFEM generated synthetics, used to generate SAC headers

  • location (str) – location value for a given station/component

  • precision (int) – dt precision determined by differencing two adjancent time steps in the underlying ascii text file.

Rtype st:

obspy.Stream.stream

Return st:

stream containing header and data info taken from ascii file

pysep.utils.io.read_sem_cartesian(fid, source, stations, location='', precision=4)[source]#

Specfem2D and Specfem3D may have domains defined in a Cartesian coordinate system. Because of this, the read_sem() function will fail because the intermediate ObsPy objects and functions expect geographic coordinates. This function bypasses these checks with some barebones objects which mimic their behavior. Only used for RecSec to plot record sections.

TODO can we combine this with cap_sac.append_sac_headers()? Currently the

code block at the bottom (manually appending header) is redundant and should use cap_sac?

Note

RecSec requires SAC header values kevnm, dist, az, baz, stlo, stla, evlo, evla

Parameters:
  • fid (str) – path of the given ascii file

  • source (str) – SOURCE or CMTSOLUTION file defining the event which generated the synthetics. Used to grab event information.

  • stations (str) – STATIONS file defining the station locations for the SPECFEM generated synthetics

  • location (str) – location value for a given station/component

Rtype st:

obspy.Stream.stream

Return st:

stream containing header and data info taken from ascii file

pysep.utils.io.read_events_plus(fid, format, **kwargs)[source]#

Addition to the base ObsPy.read_events() function that, in addition to the acceptable formats read by ObsPy, can also read the following: * SPECFEM2D SOURCE * SPECFEM3D/3D_GLOBE FORCESOLUTION * SPECFEM3D/3D_GLOBE CMTSOLUTION (both geographic and non-geographic)

See the following link for acceptable ObsPy formats: See the following link for acceptable ObsPy formats: https://docs.obspy.org/packages/autogen/obspy.core.event.read_events.html

Parameters:
  • fid (str) – full path to the event file to be read

  • format (str) – Expected format of the file (case-insensitive), available are - SOURCE - FORCESOLUTION - CMTSOLUTION - any of ObsPy’s accepted arguments for ObsPy.read_events()

Return type:

obspy.core.catalog.Catalog

Returns:

Catalog which should only contain one event, read from the fid for the given fmt (format)

pysep.utils.io.read_specfem3d_cmtsolution_cartesian(path_to_cmtsolution)[source]#

Create a barebones ObsPy Event object from a SPECFEM3D CMTSOLUTION file with coordinates defined in cartesian. Required because ObsPy read_events() will throw a ValueError when CMTSOLUTIONS do not have geographically defined coordinates.

pysep.utils.io.read_specfem2d_source(path_to_source, origin_time=None)[source]#

Create a barebones ObsPy Event object from a SPECFEM2D Source file, which only contains information required by Pyatoa.

Only requires access to: event.preferred_origin(), event.preferred_magnitude() and event.preferred_moment_tensor(). Moment tensor is wrapped in try-except so we only need origin and magnitude.

Modified from: https://docs.obspy.org/master/_modules/obspy/io/cmtsolution/

core.html#_internal_read_cmtsolution

Note

Source files do not provide origin times so we just provide an arbitrary value but allow user to set time

pysep.utils.io.read_forcesolution(path_to_forcesolution, origin_time='2000-01-01T00:00:00')[source]#

Create a barebones Source object from a FORCESOLUTION Source file, which mimics the behavior of the more complex ObsPy Event object and can be used in the same way as an Event object.

Note

Designed to read FORCESOLUTION files from SPECFEM3D/3D_GLOBE, which all have slightly different formats/keys

Parameters:
  • path_to_forcesolution (str) – path to the FORCESOLUTION file

  • origin_time (str) – FORCESOLUTION files do not natively contain any information on origin time, which is required by ObsPy Event objects. The User can provide this information if it is important, or a default value of 2000-01-01T00:00:00 will be provided as a dummy variable to keep ObsPy happy

Return type:

obspy.core.event.Event

Returns:

Barebones ObsPy Event object which contains hypocentral location of FORCE, and the origin time defined by origin_time

Raises:

KeyError – if the minimum required keys are not found in the file defined by path_to_source

pysep.utils.io._get_resource_id(name, res_type, tag=None)[source]#

Helper function to create consistent resource ids, from ObsPy. Used to create resource ID’s when generating Event objects

pysep.utils.io.read_stations(path_to_stations)[source]#

Convert a SPECFEM STATIONS file into an ObsPy Inventory object.

Specfem3D STATION files contain no channel or location information, so the inventory can only go down to the station level.

Note

This assumes a row structure for the station file is STA, NET, LAT [deg], LON [deg], ELEVATION [m], BURIAL [m]

Parameters:

path_to_stations (str) – the path to the STATIONS file that is associated with the Specfem3D DATA directory

Return type:

obspy.core.inventory.Inventory

Returns:

a station-level Inventory object

Raises:

ValueError – if latitude and longitude values are not in geographic coordinates (i.e., in cartesian coordinates). Thrown by the init of the Station class.

pysep.utils.io.read_event_file(fid)[source]#

Read an event input file which is just a text file that should contain information on different events and their hypocenters

Parameters:

fid (str) – event input file

Return type:

list of dict

Returns:

parsed in event information

pysep.utils.io.write_sem(st, unit, path='./', time_offset=0)[source]#

Write an ObsPy Stream in the two-column ASCII format that Specfem uses

Parameters:
  • st (obspy.core.stream.Stream) – stream containing synthetic data to be written

  • unit (str) – the units of the synthetic data, used for file extension, must be ‘d’, ‘v’, ‘a’ for displacement, velocity, acceleration, resp.

  • path (str) – path to save data to, defaults to cwd

  • time_offset (float) – optional argument to offset the time array. Sign matters e.g. time_offset=-20 means t0=-20

pysep.utils.io.write_pysep_station_file(inv, event, fid='./station_list.txt', order_station_list_by=None)[source]#

Write a list of station codes, distances, etc. useful for understanding characterization of all collected stations

Parameters:
  • event (obspy.core.event.Event) – optional user-provided event object which will force a skip over QuakeML/event searching

  • inv (obspy.core.inventory.Inventory) – optional user-provided inventory object which will force a skip over StationXML/inventory searching

  • fid (str) – name of the file to write to. defaults to ./station_list.txt

  • order_station_list_by (str) – how to order the stations written to file. Available are: network, station, latitude, longitude, elevation, burial. If not given, order is determined by the internal sorting of the input Inventory

pysep.utils.io.write_stations_file(inv, fid='./STATIONS', order_by=None, use_elevation=False, burials=None)[source]#

Write a SPECFEM STATIONS file given an ObsPy inventory object

Note

If topography is implemented in your mesh, elevation values should be set to 0 which means ‘surface’ in SPECFEM.

Note

burial information is not contained in the ObsPy inventory so it is always set to a constant value, which can be given as input. default 0

Parameters:
  • inv (obspy.core.inventory.Inventory) – Inventory object with station locations to write

  • fid (str) – path and file id to save the STATIONS file to.

  • order_by (str) – how to order the stations written to file. Available are: network, station, latitude, longitude, elevation, burial. If not given, order is determined by the internal sorting of the input Inventory

  • use_elevation (bool) – if True, sets the actual elevation value from the inventory. if False, sets elevation to 0.

  • burials (list of float) – If the User has burial information they want to be used as the last column. Length of burials must match the number of stations in the inventory when traversing by network and station

pysep.utils.io.write_cat_to_event_list(cat, fid_out='event_input.txt')[source]#

Writes out an event Catalog (ObsPy object) information to an ASCII file that PySEP can use for collecting data. The format of the output file is

ORIGIN_TIME LONGITUDE LATITUDE DEPTH[KM] MAGNITUDE

Parameters:
  • cat (obspy.core.catalog.Catalog) – Catalog of events to write out

  • fid_out (str) – name of the output text file to be written. Defaults to ‘event_input.txt’