Quake
The Quake submodule was introduced in SeisBase v0.3.0 to isolate handling of discrete earthquake events from handling of continuous geophysical data. While the channel data are similar, fully describing an earthquake event requires many additional Types (objects) and more information (fields) in channel descriptors.
Types
SeisBase.Quake.EQMag
— TypeEQMag
Earthquake magnitude container object
Field | Type | Meaning |
---|---|---|
val | Float32 | numeric magnitude value (note: Float32!) |
scale | String | magnitude scale (freeform) |
gap | Float64 | azimuthal gap (°) |
nst | Int64 | number of stations used in magnitude calculation |
src | String | magnitude source |
SeisBase.Quake.EQLoc
— TypeEQLoc
QuakeML-compliant earthquake location
Field | Type | Meaning | SeisBase conventions/behavior |
---|---|---|---|
lat | Float64 | latitude | °N = + |
lon | Float64 | longitude | °E = + |
dep | Float64 | depth | km; down = + |
dx | Float64 | x error | uses units of data source (typically km) |
dy | Float64 | y error | uses units of data source (typically km) |
dz | Float64 | z error | uses units of data source (typically km) |
dt | Float64 | ot error | uses units of data source (typically s) |
se | Float64 | std error | uses units of data source (typically s) |
rms | Float64 | rms pick err | uses units of data source (typically s) |
gap | Float64 | azimuthal gap | uses units of data source (typically °) |
dmin | Float64 | min sta dist | uses units of data source (typically km) |
dmax | Float64 | max sta dist | uses units of data source (typically km) |
nst | Int64 | # of stations | |
flags | UInt8 | boolean flags | access flag[n] with >>(<<(flags,n-1),7) |
datum | String | geog. datum | |
typ | String | location type | freeform (e.g. "centroid", "hypocenter") |
sig | String | significance | freeform (e.g. "95%", "2σ") |
/ confidence | |||
src | String | source | freeform (e.g. "HYPOELLIPSE", "HypoDD") |
flags (0x01 = true, 0x00 = false)
- x fixed?
- y fixed?
- z fixed?
- t fixed?
SeisBase.Quake.EventChannel
— TypeEventChannel
A single channel of trace data (digital seismograms) associated with a discrete event (earthquake).
See EventTraceData
SeisBase.Quake.EventTraceData
— TypeEventTraceData
A custom structure designed to describe trace data (digital seismograms) associated with a discrete event (earthquake).
EventChannel
A single channel of trace data (digital seismograms) associated with a discrete event (earthquake).
Fields: EventTraceData, EventChannel, SeisEvent.data
Field | Description |
---|---|
:n | Number of channels [1] |
:id | Channel id. Uses NET.STA.LOC.CHA format when possible |
:name | Freeform channel name |
:loc | Location (position) vector; any subtype of InstrumentPosition |
:fs | Sampling frequency in Hz; fs=0.0 for irregularly-sampled data. |
:gain | Scalar gain |
:resp | Instrument response; any subtype of InstrumentResponse |
:units | String describing data units. UCUM standards are assumed. |
:az | Source azimuth |
:baz | Backazimuth to source |
:dist | Source-receiver distance |
:pha | Seismic phase catalog |
:src | Freeform string describing data source. |
:misc | Dictionary for non-critical information. |
:notes | Timestamped notes; includes automatically-logged information. |
:t | Matrix of time gaps in integer μs, formatted [Sample# Length] |
:x | Time-series data |
- Not present in EventChannel objects.
See also: PhaseCat
, SeisPha
, SeisData
SeisBase.Quake.SeisEvent
— TypeSeisEvent
A structure for discrete seismic events, comprising three structures:
- :hdr, a SeisHdr for the event descriptor
- :source, a SeisSrc for descrition of the seismic source process
- :data, an EventTraceData structure for channel data, including phases
See also: SeisHdr
, SeisSrc
, EventTraceData
SeisBase.Quake.SeisHdr
— TypeSeisHdr: header information for seismic events
S = SeisHdr()
Initialize an empty SeisHdr object. Fields can be initialized at creation with keywords, e.g., SeisHdr(ot=DateTime("2012-01-03T03:49:45"), int=(0x02, "MMI")).
Field | Default | Type | Meaning |
---|---|---|---|
id | "" | String | Event ID |
int | (0x00, "") | Tuple | (Intensity, Intensity Scale) |
loc | () | EQLoc | Hypocenter data |
mag | () | EQMag | Magnitude data |
misc | () | Dict{String,Any}() | Non-essential info |
ot | (unix epoch) | DateTime | Origin time |
notes | [] | Array{String,1} | Timestamped notes, logging |
src | "" | String | Data source (URL/filename) |
typ | "" | String | Event type |
See also: EQLoc
, EQMag
SeisBase.Quake.SeisPha
— TypeSeisPha()
IRIS-style seismic phase and pick container
Field | Type | Meaning | SeisBase conventions/behavior |
---|---|---|---|
amp | Float64 | amplitude | uses units of data source |
d | Float64 | distance | no unit conversion; can be m, km, or ° |
ia | Float64 | incidence angle | uses units of data source |
res | Float64 | pick residual | |
rp | Float64 | ray parameter | |
ta | Float64 | takeoff angle | |
tt | Float64 | travel time | |
unc | Float64 | uncertainty | |
pol | Char | polarity | |
qual | Char | pick quality | not (re)calculated |
SeisBase.Quake.SeisSrc
— TypeSeisSrc: container for descriptions of a seismic source process
S = SeisSrc()
Initialize an empty SeisSrc object. Fields can be initialized at creation with keywords; for example, S = SeisSrc(m0 = 1.6e22).
Field | Type | Meaning |
---|---|---|
id | String | source process ID |
eid | String | event ID (note, generally :id != :eid) |
m0 | Float64 | scalar seismic moment |
mt | Array{Float64,1} | seismic moment tensor |
dm | Array{Float64,1} | seismic moment tensor misfit |
npol | Int64 | number of polarities in focal mechanism |
gap | Float64 | max azimuthal gap in focal mechanism |
pax | Array{Float64,2} | principal axes |
planes | Array{Float64,2} | nodal planes |
src | String | data source string (filename or URL) |
st | SourceTime | source-time subfield |
misc | Dict{String,Any} | dictionary of non-essential information |
notes | Array{String,1} | notes and automated logging |
See also: EQLoc
, EQMag
, SourceTime
SeisBase.Quake.SourceTime
— TypeSourceTime()
QuakeML-compliant seismic source-time parameterization.
Field | Type | Meaning | SeisBase conventions/behavior |
---|---|---|---|
desc | String | description | |
dur | Float64 | duration | |
rise | Float64 | rise time | |
decay | Float64 | decay time |
Web Queries
Keyword descriptions for web queries appear at the end of this section.
SeisBase.Quake.FDSNevq
— Function(H,R) = FDSNevq(ot)
Multi-server query for the events with the closest origin time to ot
. Returns an Array{SeisHdr,1} in H with event headers and an Array{SeisSrc,1} in R in H with corresponding source process info.
Keywords
KW | Default | T [1] | Meaning |
---|---|---|---|
evw | [600., 600.] | Float64 | search window in seconds [2] |
mag | [6.0, 9.9] | Float64 | search magitude range |
nev | 0 | Integer | events per query [3] |
rad | [] | Float64 | radius search |
reg | [] | Float64 | geographic search region |
src [4] | "IRIS" | String | data source; ?seis_www lists |
to | 30 | Int64 | timeout (s) for web requests |
v | 0 | Integer | verbosity |
Array{T, 1}
forevw
,mag
,rad
,reg
;T
for others- search range is always
ot-|evw[1]| ≤ t ≤ ot+|evw[2]|
- if
nev=0
, all matches are returned. - In an event query, keyword
src
can be a comma-delineated list, like"IRIS, INGV, NCEDC"
.
Notes
- Specify
ot
as a string formatted YYYY-MM-DDThh:mm:ss in UTC (e.g. "2001-02-08T18:54:32"). - Incomplete string queries are read to the nearest fully-specified time constraint; thus,
FDSNevq("2001-02-08")
returns the nearest event to 2001-02-08T00:00:00. - If no event is found in the specified search window, FDSNevq exits with an error.
- For FDSNevq, keyword
src
can be a comma-delineated list of sources, provided each has a value in?seis_www
; for example,src="IRIS, INGV, NCEDC"
is valid.
See also: SeisBase.KW
, ?seis_www
SeisBase.Quake.FDSNevt
— FunctionFDSNevt(ot::String, chans::String)
Get header and trace data for the event closest to origin time ot
on channels chans
. Returns a SeisEvent structure.
Keywords
KW | Default | T [1] | Meaning |
---|---|---|---|
evw | [600., 600.] | Float64 | search window in seconds [2] |
fmt | "miniseed" | String | request data format |
len | 120.0 | Float64 | desired trace length [s] |
mag | [6.0, 9.9] | Float64 | search magitude range |
model | "iasp91" | String | velocity model for phases |
nd | 1 | Real | number of days per subrequest |
opts | "" | String | user-specified options[3] |
pha | "P" | String | phases to get [4] |
rad | [] | Float64 | radius search |
reg | [] | Float64 | geographic search region |
src | "IRIS" | String | data source; ?seis_www lists |
to | 30 | Int64 | timeout (s) for web requests |
v | 0 | Integer | verbosity |
w | false | Bool | write requests to disk? |
- KW is
Array{T, 1}
forevw
,mag
,rad
,reg
, typeT
for others - Search range is always
ot-|evw[1]| ≤ t ≤ ot+|evw[2]|
- Format like an http request string, e.g. "szsrecs=true&repo=realtime" for FDSN. String shouldn't begin with an ampersand.
- Comma-separated String, like
"P, pP"
; use"ttall"
for all phases
Notes
- Specify
ot
as a string formatted YYYY-MM-DDThh:mm:ss in UTC (e.g. "2001-02-08T18:54:32"). - Incomplete string queries are read to the nearest fully-specified time constraint; thus,
FDSNevq("2001-02-08")
returns the nearest event to 2001-02-08T00:00:00. - If no event is found in the specified search window, FDSNevt exits with an error.
- Unlike
FDSNevq
, number of events cannot be specified andsrc
must be a single source String in?seis_www
.
See also: distaz!
, FDSNevq
, FDSNsta
SeisBase.Quake.get_pha!
— Functionget_pha!(Ev::SeisEvent[, keywords])
Command-line interface to IRIS online travel time calculator, which calls TauP [1-2]. Returns a matrix of strings.
Keywords:
- pha: comma-separated String of phases ("P, S, SP")
- model: velocity model ("iasp91")
- to: timeout in seconds
- v: verbosity
References
- TauP manual: http://www.seis.sc.edu/downloads/TauP/taup.pdf
- Crotwell, H. P., Owens, T. J., & Ritsema, J. (1999). The TauP Toolkit:
Flexible seismic travel-time and ray-path utilities, SRL 70(2), 154-160.
Web Query Keywords
KW | Default | T [1] | Meaning |
---|---|---|---|
evw | [600.0, 600.0] | A{F,1} | search window in seconds [2] |
fmt | "miniseed" | S | request data format |
len | 120.0 | I | desired trace length [s] |
mag | [6.0, 9.9] | A{F,1} | magnitude range for queries |
model | "iasp91" | S | Earth velocity model for phase times |
nd | 1 | I | number of days per subrequest |
nev | 0 | I | number of events returned per query [3] |
opts | "" | S | user-specified options [4] |
pha | "P" | S | phases to get [5] |
rad | [] | A{F,1} | radial search region [6] |
reg | [] | A{F,1} | rectangular search region [7] |
src | "IRIS" | S | data source; type ?seis_www for list |
to | 30 | I | read timeout for web requests [s] |
v | 0 | I | verbosity |
w | false | B | write requests to disk? [8] |
Table Footnotes
- Types: A = Array, B = Boolean, C = Char, DT = DateTime, F = Float, I = Integer, S = String, U8 = Unsigned 8-bit integer (UInt8)
- search range is always $ot-|evw[1]| ≤ t ≤ ot+|evw[2]|$
- nev=0 returns all events in the query
- String is passed as-is, e.g. "szsrecs=true&repo=realtime" for FDSN. String should not begin with an ampersand.
- Comma-separated String, like
"P, pP"
; use"ttall"
for all phases - Specify region [centerlat, centerlon, minradius, maxradius, depmin, depmax], with lat, lon, and radius in decimal degrees (°) and depth in km with + = down. Depths are only used for earthquake searches.
- Specify region [latmin, latmax, lonmin, lonmax, depmin, depmax], with lat, lon in decimal degrees (°) and depth in km with + = down. Depths are only used for earthquake searches.
- If w=true, a file name is automatically generated from the request parameters, in addition to parsing data to a SeisData structure. Files are created from the raw download even if data processing fails, in contrast to get_data(... wsac=true).
Example
Get seismic and strainmeter records for the P-wave of the Tohoku-Oki great earthquake on two borehole stations and write to native SeisData format:
S = FDSNevt("201103110547", "PB.B004..EH?,PB.B004..BS?,PB.B001..BS?,PB.B001..EH?")
wseis("201103110547_evt.seis", S)
Utility Functions
SeisBase.Quake.distaz!
— Functiondistaz!(Ev::SeisEvent)
Compute Δ, Θ by the Haversine formula.
Updates Ev.data
with distance, azimuth, and backazimuth for each channel, written to Ev.data.dist, Ev.data.az, and Ev.data.baz, respectively.
SeisBase.Quake.gcdist
— FunctionG = gcdist(src, rec)
Compute great circle distance, azimuth, and backazimuth from single source s
with coordinates [s_lat, s_lon]
to receivers r
with coordinates [r_lat r_lon].
For a single source, pass src
as a Float64 vector of the form [s_lat, s_lon]
; gcdist will return an Array{Float64,2} of the form
[Δ₁ θ₁ β₁
Δ₂ θ₂ β₂
⋮ ⋮ ⋮
Δn θn βn]
for receivers 1:n
.
For multiple sources, pass src
as an Array{Float64,2} with each row containing one (lat, lon) pair. This returns a three-dimensional matrix where each two-dimensional slice takes the form
[Δᵢ₁ θᵢ₁ βᵢ₁
⋮ ⋮ ⋮
Δᵢn θᵢn βᵢn]
for source i
at receivers 1:n
.
SeisBase.Quake.show_phases
— Functionshow_phases(io::IO, PC::PhaseCat)
Formatted display of seismic phases in dictionary P.
SeisBase.Quake.fill_sac_evh!
— Functionfill_sac_evh!(Ev::SeisEvent, fname::String; k=i)
Fill (overwrite) values in Ev.hdr
with data from SAC file fname
. Keyword k=i
specifies the reference channel i
from which the absolute origin time Ev.hdr.ot
is set. Potentially affects header fields :id
, :loc
(subfields .lat, .lon, .dep), and :ot
.
Reading Earthquake Data Files
SeisBase.read_quake
— FunctionEv = read_quake(fmt, file)
Read data in file format fmt
from file
into SeisEvent object Ev
.
- Formats: suds, qml, uw
- Keywords: full, v
Note: because earthquake data are usually discrete, self-contained files, no "in-place" version of read_quake
exists, and read_quake
doesn't support wildcards in the file string.
Supported File Formats
File Format | String | Notes |
---|---|---|
PC-SUDS | suds | |
QuakeML | qml, quakeml | only reads first event from file |
UW | uw |
See also: read_data
, get_data
, read_meta
, UW.readuwevt