Function index
Public types and functions
IO
QuakeML.read
— Methodread(filename) -> ::EventParameters
Read a QuakeML file with name filename
from disk and return an EventParameters
object.
Example
julia> file = joinpath(dirname(pathof(QuakeML)), "..", "test", "data", "nepal_mw7.2.qml");
julia> events = QuakeML.read(file)
EventParameters
comment: Array{QuakeML.Comment}((0,))
event: Array{QuakeML.Event}((1,))
description: Missing missing
creation_info: QuakeML.CreationInfo
public_id: QuakeML.ResourceIdentifier
read(io) -> ::EventParameters
Read a QuakeML document from the stream io
.
Example
julia> io = IOBuffer("""
<quakeml xmlns="http://quakeml.org/xmlns/quakeml/1.2">
<eventParameters></eventParameters>
</quakeml>
""");
julia> QuakeML.read(io)
EventParameters
comment: Array{QuakeML.Comment}((0,))
event: Array{QuakeML.Event}((0,))
description: Missing missing
creation_info: Missing missing
public_id: QuakeML.ResourceIdentifier
QuakeML.readstring
— Functionreadstring(xml_string) -> ::EventParameters
Read the QuakeML contained in xml_string
and return a EventParameters
object.
Base.write
— Methodwrite(io, qml::EventParameters; kwargs...)
Write a set of EventParameters
to io
. kwargs
are passed to quakeml
to control the creation of the XML representing the catalogue.
Examples
(Note that "example_quakeml_file.xml"
may not exist.)
Write a file with the default settings
julia> qml = QuakeML.read("example_quakeml_file.xml");
julia> write("new_file.xml", qml)
Write a file with custom version
julia> write("new_file2.xml", qml, version="1.1")
QuakeML.quakeml
— Functionquakeml(qml::EventParameters; version="1.2") -> xml::EzXML.XMLDocument
Create an XML document from qml
, a set of events of type EventParameters
. xml
is an EzXML.XMLDocument
suitable for output.
The user may also set the nominal version
of QuakeML created.
The QuakeML document xml
may be written with write(io, xml)
or converted to a string with string(xml)
.
Accessors
QuakeML.preferred_focal_mechanism
— Functionpreferred_focal_mechanism(event; verbose=false) -> focal_mechanism
Return the preferred focal mechanism for an event
. This may be defined if there is more than one focal mechanism given for an event, and the preferred_focal_mechanism_id
field is set. If there is only one focal mechanism for this event
, then that is returned. If there is no focal mechanism associated with this event which matches the stated preferred_focal_mechanism_id
, then the first focal mechanism is returned, and a warning given is verbose
is true
.
QuakeML.preferred_magnitude
— Functionpreferred_magnitude(event; verbose=false) -> magnitude
Return the preferred magnitude for an event
. This may be defined if there is more than one magnitude given for an event, and the preferred_magnitude_id
field is set. If there is only one magnitude for this event
, then that is returned. If there is no magnitude associated with this event which matches the stated preferred_magnitude_id
, then the first magnitude is returned, and a warning given is verbose
is true
.
QuakeML.preferred_origin
— Functionpreferred_origin(event; verbose=false) -> origin
Return the preferred origin for an event
. This may be defined if there is more than one origin given for an event, and the preferred_origin_id
field is set. If there is only one origin for this event
, then that is returned. If there is no origin associated with this event which matches the stated preferred_origin_id
, then the first origin is returned and a warning is given when verbose=true
QuakeML.has_focal_mechanism
— Functionhas_focal_mechanism(event) -> ::Bool
Return true
if event
contains one or more focal mechanisms defined.
QuakeML.has_magnitude
— Functionhas_magnitude(event) --> ::Bool
Return true
if event
has one or more magnitudes defined.
QuakeML.has_origin
— Functionhas_origin(event) -> ::Bool
Return true
if event
has one or more origins defined.
Types
Where a constructor's arguments are given as Constructor(; kwargs...)
, this means that each listed field name can be given as a keyword argument and a value passed to the constructor that way.
By default, calling a constructor for a type which is required to have a public ID (URI) by the QuakeML specification creates a unique, random URI for that object. To specify your own ID for an object, provide a String
to the constructor's public_id
keyword argument; or you can later set the public_id
field directly. See ResourceReference
for details of the form that URIs must take.
QuakeML.Amplitude
— TypeAmplitude(; kwargs...)
Represents a quantification of the waveform anomaly, usually a single amplitude measurement or a measurement of the visible signal duration for duration magnitudes.
List of fields
public_id :: ResourceReference
: Resource identifier ofAmplitude
. (Required field.)genericAmplitude :: RealQuantity
: Measured amplitude value for the givenwaveform_id
. Note that this attribute can describe different physical quantities, depending on thetype
andcategory
of the amplitude. These can be, e.g., displacement, velocity, or a period. If the only amplitude information is a period, it has to specified here, not in theperiod
field. The latter can be used if the amplitude measurement contains information on, e.g., displacement and an additional period. Since the physical quantity described by this attributeis not fixed, the unit of measurement cannot be defined in advance. However, the quantity has to be specified in SI base units. The enumeration given in the fieldunit
provides the most likely units that could be needed here. For clarity, using the optionalunit
field is highly encouraged. (Required field.)type :: String
:String
that describes the type of amplitude using the nomenclature from Storchak et al. (2003). Possible values include unspecified amplitude reading ("A"
), amplitude reading for local magnitude ("AML"
), amplitude reading for body wave magnitude ("AMB"
), amplitude reading for surface wave magnitude ("AMS"
), and time of visible end of record for duration magnitude ("END"
). It has a maximum length of 32 characters.category :: AmplitudeCategory
: This field describes the way the waveform trace is evaluated to derive an amplitude value. This can be just reading a single value for a given point in time ("point"
), taking a mean value over a time interval ("mean"
), integrating the trace over a time interval ("integral"
), specifying just a time interval ("duration"
), or evaluating a period ("period"
). (SeeAmplitudeCategory
.)"point"
"mean"
"duration"
"period"
"integral"
"other"
unit :: AmplitudeUnit
: This field provides the most likely measurement units for the physical quantity described in thegeneric_Amplitude
field. Possible values are specified as combinations of SI base units. (SeeAmplitudeUnit
"m"
"s:
"m/s"
"m/(s*s)"
"m*s"
"dimensionless"
"other"
method_id :: ResourceReference
: Describes the method of amplitude determination.period :: RealQuantity
: Dominant period in thetime_window
in case of amplitude measurements. Not used for duration magnitude. Unit: s.snr :: Float64
: Signal-to-noise ratio of the spectrogram at the location the amplitude was measured.time_window :: TimeWindow
: Description of the time window used for amplitude measurement. Recommended for duration magnitudes.pick_id :: ResourceReference
: Refers to thepublic_id
of an associatedPick
object.waveform_id :: ResourceReference
: Identifies the waveform stream on which the amplitude was measured.filter_id :: ResourceReference
: Identifies the filter or filter setup used for filtering the waveform stream referenced bywaveform_id
.scaling_time :: TimeQuantity
: Scaling time for amplitude measurement.magnitude_hint :: String
: Type of magnitude the amplitude measurement is used for. For valid values seeMagnitude
. String value with a maximum length of 32 characters.evaluation_mode :: EvaluationMode
: Evaluation mode ofAmplitude
(seeEvaluationMode
).evaluation_status :: EvaluationStatus
: Evaluation status ofAmplitude
(seeEvaluationStatus
).comment :: Vector{Comment}
: Additional comments.creation_info :: CreationInfo
:CreationInfo
for theAmplitude
object.
QuakeML.Arrival
— TypeArrival(; kwargs...)
Successful association of a pick with an origin qualifies this pick as an arrival. An arrival thus connects a pick with an origin and provides additional attributes that describe this relationship. Usually qualification of a pick as an arrival for a given origin is a hypothesis, which is based on assumptions about the type of arrival (phase) as well as observed and (on the basis of an earth model) computed arrival times, or the residual, respectively. Additional pick attributes like the horizontal slowness and backazimuth of the observed wave—especially if derived from array data—may further constrain the nature of the arrival.
List of fields
public_id :: ResourceReference
: Resource identifier ofArrival
. (Required field.)pick_id :: ResourceReference
: Refers to apublic_id
of aPick
. (Required field.)phase :: Phase
: Phase identification. For possible values, please refer to the description of thePhase
type. (Required field.)time_correction :: Float64
: Time correction value. Usually, a value characteristic for the station at which the pick was detected, sometimes also characteristic for the phase type or the slowness. Unit: s.azimuth :: Float64
: Azimuth of station as seen from the epicenter. Unit: °.distance :: Float64
: Epicentral distance. Unit: °.takeoff_angle :: RealQuantity
: Angle of emerging ray at the source, measured against the downward normal direction. Unit: °.time_residual :: Float64
: Residual between observed and expected arrival time assuming proper phase identification and given theearth_model_id
of theOrigin
, taking into account thetimeCorrection
. Unit: s.horizontal_slowness_residual :: Float64
: Residual of horizontal slowness and the expected slowness given the current origin (refers to fieldhorizontal_slowness
ofPick
). Unit: s/°backazimuthResidual :: Float64
: Residual of backazimuth and the backazimuth computed for the current origin (refers to fieldbackazimuth
ofPick
). Unit: °.time_weight :: Float64
: Weight of the arrival time for computation of the associatedOrigin
. Note that the sum of all weights is not required to be unity.horizontal_slowness_weight :: Float64
: Weight of the horizontal slowness for computation of the associatedOrigin
. Note that the sum of all weights is not required to be unity.backazimuth_weight :: Float64
: Weight of the backazimuth for computation of the associatedOrigin
. Note that the sum of all weights is not required to be unity.earth_model_id :: ResourceReference
: Earth model which is used for the association ofArrival
toPick
and computation of the residuals.comment :: Vector{Comment}
: Additional comments.creation_info :: CreationInfo
:CreationInfo
for theArrival
object
QuakeML.Axis
— TypeAxis(; azimuth, plunge, length)
Describes an eigenvector of a moment tensor expressed in its principal-axes system. It uses the angles azimuth
, plunge
, and the eigenvalue length
.
List of fields
azimuth :: RealQuantity
: Azimuth of eigenvector of moment tensor expressed in principal-axes system. Measured clockwisefrom south-north direction at epicenter. Unit: °. (Required field.)plunge :: RealQuantity
: Plunge of eigenvector of moment tensor expressed in principal-axes system. Measured against downward vertical direction at epicenter. Unit: °. (Required field.)length :: RealQuantity
: Eigenvalue of moment tensor expressed in principal-axes system. Unit: N m. (Required field.)
QuakeML.Comment
— TypeComment(; text, creation_info, id)
Holds information on comments to a resource as well as author and creation time information.
List of fields
text :: String
: Text of comment. (Required field.)creation_info :: CreationInfo
:CreationInfo
for theComment
object.id :: ResourceReference
: Identifier of comment, in QuakeML URI format.
QuakeML.CompositeTime
— TypeCompositeTime(; year, month, day, hour, minute, second)
Focal times differ significantly in their precision. While focal times of instrumentally located earthquakes areestimated precisely down to seconds, historic events have only incomplete time descriptions. Sometimes, even contradictory information about the rupture time exist. The CompositeTime
type allows for such complex descriptions. If the specification is given with no greater accuracy than days (i.e., no time components are given), the date refers to local time. However, if time components are given, they have to refer to UTC.
As an example, consider a historic earthquake in California, e.g., on 28 February 1730, with no time information given. Expressed in UTC, this day extends from 1730-02-28T08:00:00Z until 1730-03-01T08:00:00Z. Such a specification would be against intuition. Therefore, for date-time specifications without time components, local time is used. In the example, the CompositeTime
fields are simply year=1730
, month=2
, and day=28
. In the corresponding time attribute of the origin, however, UTC has to be used. If the unknown time components are assumed to be zero, the value is DateTime("1730-02-28T08:00:00")
.
List of fields
year ::IntegerQuantity
: Year or range of years of the event’s focal time.month ::IntegerQuantity
: Month or range of months of the event’s focal time.day ::IntegerQuantity
: Day or range of days of the event’s focal time.hour ::IntegerQuantity
: Hour or range of hours of the event’s focal time.minute ::IntegerQuantity
: Minute or range of minutes of the event’s focal time.second :: RealQuantity
: Second and fraction of seconds or range of seconds with fraction of the event’s focal time.
QuakeML.ConfidenceEllipsoid
— TypeConfidenceEllipsoid(; kwargs...)
This type represents a description of the location uncertainty as a confidence ellipsoid with arbitrary orientationin space. The orientation of a rigid body in three-dimensional Euclidean space can be described by three parameters. We use the convention of Euler angles, which can be interpreted as a composition of three elemental rotations (i.e., rotations around a single axis). In the special case of Euler angles we use here, the angles are referred to as Tait-Bryan (or Cardan) angles. These angles may be familiar to the reader from their application in flight dynamics, and are referred to as heading (yaw, ψ), elevation (attitude, pitch, φ), and bank (roll, θ). For a definition of the angles, see Figure 4 of the QuakeML specification document at https://quake.ethz.ch/quakeml/docs/latest?action=AttachFile&do=get&target=QuakeML-BED.pdf. Through the three elemental rotations, a Cartesian system (x,y,z)
centered at the epicenter, with the south-north direction x
, the west-east direction y
, and the downward vertical direction z
, is transferred into a different Cartesian system (X,Y,Z)
centered on the confidence ellipsoid. Here, X
denotes the direction of the major axis, and Y
denotes the direction of the minor axis of the ellipsoid. Note that Figure 4 can be interpreted as a hypothetical view from the interior of the Earth to the inner face of a shell representing Earth's surface.
The three Tait-Bryan rotations are performed as follows: (i) a rotation about the Z
axis with angle ψ (heading, or azimuth); (ii) a rotation about the Y
axis with angle φ (elevation, or plunge); and (iii) a rotation about the X
axis with angle θ (bank). Note that in the case of Tait-Bryan angles, the rotations are performed about the ellipsoid's axes, not about the axes of the fixed (x,y,z)
Cartesian system.
List of fields
semi_major_axis_length :: Float64
: Largest uncertainty, corresponding to the semi-major axis of the confidence ellipsoid. Unit: m. (Required field.)semi_minor_axis_length :: Float64
: Smallest uncertainty, corresponding to the semi-minor axis of the confidence ellipsoid. Unit: m. (Required field.)semi_intermediate_axis_length :: Float64
: Uncertainty in direction orthogonal to major and minor axesof the confidence ellipsoid. Unit: m. (Required field.)major_axis_plunge :: Float64
: Plunge angle of major axis of confidence ellipsoid. Corresponds to Tait-Bryan angle φ. Unit: °. (Required field.)major_axis_azimuth :: Float64
: Azimuth angle of major axis of confidence ellipsoid. Corresponds to Tait-Bryan angle ψ. Unit: °. (Required field.)major_axis_rotation :: Float64
: This angle describes a rotation about the confidence ellipsoid's major axis which is required to define the direction of the ellipsoid's minor axis. Corresponds to Tait-Bryan angle θ. Unit: °. (Required field.)
QuakeML.CreationInfo
— TypeCreationInfo(; kwargs...)
Used to describe creation metadata (author, version, and creation time) of a resource.
List of fields
agency_id :: String
: Designation of agency that published a resource. The string has a maximum length of 64 characters.agency_uri :: ResourceReference
: URI of the agency that published a resource.author :: String
: Name describing the author of a resource. The string has a maximum length of 128 characters.author_uri :: ResourceReference
: URI of the author of a resource.creation_time :: Dates.DateTime
: Time of creation of a resource, in ISO 8601 format. It has to be given in UTC.version :: String
: Version string of a resource. The string has a maximum length of 64 characters.
QuakeML.DataUsed
— TypeDataUsed(; kwargs...)
Describes the type of data that has been used for a moment-tensor inversion.
List of fields
wave_type : DataUsedWaveType
: Type of waveform data. This can be one of the following values (seeDataUsedWaveType
):"P waves"
"body waves"
"surface waves"
"mantle waves"
"combined"
"unknown"
station_count :: Int
: Number of stations that have contributed data of the type given inwave_type
.component_count :: Int
: Number of data components of the type given inwave_type
.shortest_period :: Float64
: Shortest period present in data. Unit: s.longest_period :: Float64
: Longest period present in data. Unit: s.
QuakeML.Event
— TypeEvent(; kwargs...)
Describes a seismic event which does not necessarily need to be a tectonic earthquake. An event is usually associated with one or more origins, which contain information about focal time and geographic allocation of the event. Multiple origins can cover automatic and manual locations, a set of location from different agencies, locations generated with different location programs and earth models, etc. Furthermore, an eventis usually associated with one or more magnitudes, and with one or more focal mechanism determinations. In standard QuakeML-BED, Origin
, Magnitude
, StationMagnitude
, and FocalMechanism
are fields of Event
. In BED-RT (the real-time version) all these fields are on the same hierarchy level as child elements of EventParameters
. The association of origins, magnitudes, and focal mechanisms to a particular event is expressed using references inside Event
.
List of fields
public_id :: ResourceReference
: Resource identifier ofEvent
. (Required field.)preferred_origin_id :: ResourceReference
: Refers to thepublic_id
of thepreferred_origin
object.preferred_magnitude_id :: ResourceReference
: Refers to thepublic_id
of thepreferred_magnitude
object.preferred_focal_mechanism_id :: ResourceReference
: Refers to thepublic_id
of thepreferred_focal_mechanism
object.type :: EventType
: Describes the type of an event (Storchak et al. 2012). Allowed values are the following (seeEventType
):"not existing"
"not reported"
"earthquake"
"anthropogenic event"
"collapse"
"cavity collapse"
"mine collapse"
"building collapse"
"explosion"
"accidental explosion"
"chemical explosion"
"controlled explosion"
"experimental explosion"
"industrial explosion"
"mining explosion"
"quarry blast"
"road cut"
"blasting levee"
"nuclear explosion"
"induced or triggered event"
"rock burst"
"reservoir loading"
"fluid injection"
"fluid extraction"
"crash"
"plane crash"
"train crash"
"boat crash"
"other event"
"atmospheric event"
"sonic boom"
"sonic blast"
"acoustic noise"
"thunder"
"avalanche"
"snow avalanche"
"debris avalanche"
"hydroacoustic event"
"ice quake"
"slide"
"landslide"
"rockslide"
"meteorite"
"volcanic eruption"
type_certainty :: EventTypeCertainty
: Denotes how certain the information on event type is (Storchak et al. 2012). Allowed values are the following (seeEventTypeCertainty
):"known"
"suspected"
description :: Vector{EventDescription}
Additional event description, like earthquake name, Flinn-Engdahl region, etc.comment :: Vector{Comment}
: Comments.creation_info :: CreationInfo
:CreationInfo
for theEvent
object.origin :: Vector{Event}
: Set ofOrigin
s associated with thisEvent
. One of these may be the preferred origin, in which case preferredoriginid` should be set.magnitude :: Vector{Magnitude}
: Set ofMagnitude
s for thisEvent
. One of these may be the preferred magnitude, in which casepreferred_magnitude_id
should be set.station_magnitude :: Vector{StationMagnitude}
: Set ofStationMagnitude
s contributing to the magnitude of this event.focal_mechanism :: Vector{FocalMechanism}
: Set ofFocalMechanism
s for this event. One of these may be the preferred focal mechanism, in which casepreferred_focal_mechanism_id
should be set.pick :: Vector{Pick}
: Set ofPick
s made from this event.amplitude :: Vector{Amplitude}
: Set ofAmplitude
s measured at stations from this event.
(Note: The additional real-time fields origin_reference
, magnitude_reference
and focal_mechanism_reference
are not yet implemented.)
QuakeML.EventDescription
— TypeEventDescription(; text, type)
Free-form string with additional event description. This can be a well-known name, like "1906 San Francisco Earthquake"
. A number of categories can be given in type
.
List of fields
text :: String
: Free-form text with earthquake description. (Required field.)type :: EventDescriptionType
: Category of earthquake description. Values can be taken from the following:"felt report"
"Flinn-Engdahl region"
"local time"
"tectonic summary"
"nearest cities"
"earthquake name"
"region name"
QuakeML.EventParameters
— TypeEventParameters(; comment, event, description, creation_info, public_id)
Root type of QuakeML. EventParameters
objects contain a set of events and a QuakeML XML file can contain only one EventParameters
object.
In the bulletin-type (non real-time) model, this type serves as a container for Event
objects. In the real-time version, it can hold objects of type Event
, Origin
, Magnitude
, StationMagnitude
, FocalMechanism
, Reading
, Amplitude
, and Pick
.
List of fields
event :: Vector{Event}
: Set ofEvent
s making up a catalog or collection of events.description :: String
: Description string that can be assigned to the earthquake catalog, or collection of events.comment :: Vector{Comment}
: Additional comments.creation_info :: CreationInfo
:CreationInfo
for the earthquake catalog.public_id :: ResourceReference
: Resource identifier ofEventParameters
. (Required field.)
At present, QuakeML.jl only supports the non-real-time version of QuakeML.
QuakeML.FocalMechanism
— TypeFocalMechanism(; kwargs...)
Describes the focal mechanism of an event. It includes different descriptions like nodal planes, principal axes, and a moment tensor. The moment tensor description is provided by objects of the type MomentTensor
which can be specified as fields of FocalMechanism
.
List of fields
public_id :: ResourceReference
: Resource identifier ofFocalMechanism
. (Required field.)triggering_origin_id :: ResourceReference
: Refers to thepublic_id
of the triggering origin.nodal_planes :: NodalPlanes
: Nodal planes of the focal mechanism.principal_axes :: PrincipleAxes
: Principal axes of the focal mechanism.azimuthal_gap :: Float64
: Largest azimuthal gap in distribution of stations used for determination of focal mechanism. Unit: °.station_polarity_count :: Int
: Number of station polarities used for determination of focal mechanism.misfit :: Float64
: Fraction of misfit polarities in a first-motion focal mechanism determination. Decimal fraction between 0 and 1.station_distribution_ratio :: Float64
: Station distribution ratio (STDR) parameter. Indicates how the stations are distributed about the focal sphere (Reasenberg and Oppenheimer 1985). Decimal fraction between 0 and 1.method_id :: ResourceReference
: Resource identifier of the method used for determination of the focal mechanism.waveform_id :: Vector{ResourceReference}
: Refers to a set of waveform streams from which the focal mechanism was derived.evaluation_mode :: EvaluationMode
: Evaluation mode ofFocalMechanism
(seeEvaluationMode
).evaluation_status :: EvaluationStatus
: Evaluation status ofFocalMechanism
(seeEvaluationStatus
).comment :: Vector{Comment}
: Additional comments.creation_info :: CreationInfo
:CreationInfo
for theFocalMechanism
object.
QuakeML.IntegerQuantity
— TypeIntegerQuantity(; kwargs...)
Physical quantities that can be expressed numerically—either as integers or as floating point numbers—are represented by their measured or computed values and optional values for symmetric or upper and lower uncertainties. The interpretation of these uncertainties is not defined in the standard. They can contain statistically well-defined error measures, but the mechanism can also be used to simply describe a possible value range. Ifthe confidence level of the uncertainty is known, it can be listed in the optional field confidence_level
. Note that uncertainty
, upper_uncertainty
, and lower_uncertainty
are given as absolute values of the deviation from the main value
.
List of fields
value :: Int
: Value of the quantity. The unit is implicitly defined and depends on the context. (Required field.)uncertainty :: Int
: Uncertainty as the absolute value of symmetric deviation from the main value.lower_uncertainty :: Int
: Uncertainty as the absolute value of deviation from the mainvalue
towards smaller values.upper_uncertainty :: Int
: Uncertainty as the absolute value of deviation from the mainvalue
towards larger values.confidence_level :: Float64
: Confidence level of the uncertainty, given in percent.
QuakeML.Magnitude
— TypeMagnitude(; kwargs...)
Describes a magnitude which can, but does not need to be associated with an origin. Association with an origin is expressed with the optional field origin_id
. It is either a combination of different magnitude estimations, or it represents the reported magnitude for the given event.
List of fields
public_id :: ResourceReference
: Resource identifier ofMagnitude
. (Required field.)mag :: RealQuantity
: Resulting magnitude value from combining values of typeStationMagnitude
. If no estimations are available, this value can represent the reported magnitude. (Required field.)type :: String
: Describes the type of magnitude. This is a free-text field because it is impossible to cover all existing magnitude type designations with an enumeration. Possible values are unspecified magitude ("M"
), local magnitude ("ML"
), body wave magnitude ("Mb"
), surface wave magnitude ("MS"
), moment magnitude ("Mw"
), duration magnitude ("Md"
), coda magnitude ("Mc"
),"MH"
,"Mwp"
,"M50"
,"M100"
, etc.station_magnitude_contribution :: Vector{StationMagnitudeContribution}
: Set ofStationMagnitudeContribution
s describing the contributions of each station used to compute the magnitude.origin_id :: ResourceReference
: Reference to an origin’spublic_id
if the magnitude has an associatedOrigin
.method_id :: ResourceReference
: Identifies the method of magnitude estimation. Users should avoid giving contradictory information inmethod_id
andtype
.station_count
:: Int: Number of used stations for this magnitude computation.azimuthal_gap :: Float64
: Azimuthal gap for this magnitude computation. Unit: °.evaluation_mode :: EvaluationMode
: Evaluation mode ofMagnitude
(seeEvaluationMode
).evaluation_status :: EvaluationStatus
: Evaluation status ofMagnitude
(seeEvaluationStatus
).comment :: Vector{Comment}
: Additional comments.creation_info :: CreationInfo
:CreationInfo
for theMagnitude
object.
QuakeML.MomentTensor
— TypeMomentTensor(; kwargs...)
Represents a moment tensor solution for an event. It is an optional part of a FocalMechanism
description.
List of fields
public_id :: ResourceReference
Resource identifier ofMomentTensor
. (Required field.)derived_origin_id :: ResourceReference
: Refers to thepublic_id
of theOrigin
derived in the moment tensor inversion. (Required field.)moment_magnitude_id :: ResourceReference
: Refers to thepublic_id
of theMagnitude
object which represents the derived moment magnitude.scalar_moment :: RealQuantity
: Scalar moment as derived in moment tensor inversion. Unit: N m.tensor :: Tensor
:Tensor
object holding the moment tensor elements.variance :: Float64
: Variance of moment tensor inversion.variance_reduction :: Float64
: Variance reduction of moment tensor inversion, given in percent (Dreger 2003). This is a goodness-of-fit measure.double_couple :: Float64
: Double couple parameter obtained from moment tensor inversion (decimal fraction between 0 and 1).clvd :: Float64
: CLVD (compensated linear vector dipole) parameter obtained from moment tensor inversion (decimal fraction between 0 and 1).iso :: Float64
: Isotropic part obtained from moment tensor inversion (decimal fraction between 0 and 1).greens_function_id :: ResourceReference
: Resource identifier of the Green’s function used in moment tensor inversion.filter_id :: ResourceReference
: Resource identifier of the filter setup used in moment tensor inversion.source_time_function :: SourceTimeFunction
: Source time function used in moment-tensor inversion.data_used :: Vector{DataUsed}
: Describes waveform data used for moment-tensor inversion.method_id :: ResourceReference
: Resource identifier of the method used for moment-tensor inversion.category :: MomentTensorCategory
: Category of moment tensor inversion. Valid entries are given in the following list (seeMomentTensorCategory
):"teleseismic"
"regional"
inversion_type :: MTInversionType
: Type of moment tensor inversion. Users should avoid giving contradictory information ininversion_type
andmethod_id
. Valid entries are given in the following list (seeMTInversionType
):- general
- zero trace
- double couple
comment :: Vector{Comment}
: Additional comments.creation_info :: CreationInfo
:CreationInfo
for theMomentTensor
object.
QuakeML.NodalPlane
— TypeNodalPlane(; strike, dip, rake)
This class describes a nodal plane using the fields strike
, dip
, and rake
. For a definition of the angles see Aki and Richards (1980).
List of fields
strike :: RealQuantity
: Strike angle of nodal plane. Unit: °. (Required field.)dip :: RealQuantity
: Dip angle of nodal plane. Unit: °. (Required field.)rake :: RealQuantity
: Rake angle of nodal plane. Unit: °. (Required field.)
QuakeML.NodalPlanes
— TypeNodalPlanes(; nodal_plane1=::NodalPlane, nodal_plane2=::NodalPlane, preferred_plane=::Int)
This describes the nodal planes of a moment tensor. The field preferred_plane
can be used to define which plane is the preferred one, taking a value of 1
or 2
.
List of fields
nodal_plane1 :: NodalPlane
: First nodal plane of moment tensor.nodal_plane2 :: NodalPlane
: Second nodal plane of moment tensor.preferred_plane :: Int
: Indicator for preferred nodal plane of moment tensor. It can take integer values1
or2
.
QuakeML.Origin
— TypeOrigin(; kwargs...)
Represents the focal time and geographical location of an earthquake hypocenter, as well as additional meta-information. Origin
can have objects of type OriginUncertainty
and Arrival
as fields.
List of fields
public_id :: ResourceReference
: Resource identifier ofOrigin
. (Required field.)time
: Focal time. (Required field.)longitude :: RealQuantity
: Hypocenter longitude, with respect to the World Geodetic System 1984 (WGS84) reference system (National Imagery and Mapping Agency 2000). Unit: °. (Required field.)latitude :: RealQuantity
: Hypocenter latitude, with respect to the WGS84 reference system. Unit: °. (Required field.)depth :: RealQuantity
: Depth of hypocenter with respect to the nominal sea level given by the WGS84 geoid (Earth Gravitational Model, EGM96, Lemoine et al. 1998). Positive values indicate hypocenters below sea level. For shallow hypocenters, thedepth
value can be negative. Note: Other standards use different conventions for depth measurement. As an example, GSE 2.0, defines depth with respect to the local surface. If event data is converted from other formats to QuakeML, depth values may have to be modified accordingly. Unit: m.depth_type :: OriginDepthType
: Type of depth determination. Allowed values are the following (seeOriginDepthType
):"from location"
"from moment tensor inversion",
"from modeling of broad-band P waveforms"
"constrained by depth phases",
"constrained by direct phases"
"constrained by depth and direct phases",
"operator assigned"
"other"
time_fixed :: Bool
: Boolean flag.true
if focal time was kept fixed for computation of theOrigin
.epicenter_fixed :: Bool
: Boolean flag.true
if epicenter was kept fixed for computationofOrigin
.reference_system_id :: ResourceReference
: Identifies the reference system used for hypocenter determination. This is only necessary if a modified version of the standard (with local extensions) is used that provides a non-standard coordinate system.method_id :: ResourceReference
: Identifies the method used for locating the event.earth_model_id :: ResourceReference
: Identifies the earth model used inmethodID
.composite_time :: CompositeTime
: Supplementary information on time of rupture start. Complex descriptions of focal times of historic events are possible, see description of theCompositeTime
type. Note that even ifcomposite_time
is used, the mandatorytime
field has to be set too. It has to be set to the single point in time (with uncertainties allowed) that is most characteristic for the event.quality :: OriginQuality
: Additional parameters describing the quality of anOrigin
determination.type :: OriginType
: Describes theOrigin
type. Allowed values are the following (seeOriginType
):"hypocenter"
"centroid"
"amplitude"
"macroseismic"
"rupture start"
"rupture end"
region :: String
: Can be used to decribe the geographical region of the epicenter location. Useful if an event has multiple origins from different agencies, and these have different region designations. Note that an event-wide region can be defined in thedescription
field of anEvent
object. The user has to take care that this information corresponds to the region attribute of the preferredOrigin
. String with maximum length of 255 chars.evaluation_mode :: EvaluationMode
: Evaluation mode ofOrigin
(seeEvaluationMode
.evaluation_status :: EvaluationStatus
: Evaluation status ofOrigin
(seeEvaluationStatus
).comment :: Vector{Comment}
: Additional comments.creation_info :: CreationInfo
:CreationInfo
for theOrigin
object.
QuakeML.OriginQuality
— TypeOriginQuality(; kwargs...)
This type contains various attributes commonly used to describe the quality of an origin, e. g., errors, azimuthal coverage, etc. Origin
objects have an optional attribute of the type OriginQuality
.
List of fields
associated_phase_count :: Int
: Number of associated phases, regardless of their use for origin computation.used_phase_count :: Int
: Number of defining phases, i. e., phase observations that were actually used for computingthe origin. Note that there may be more than one defining phase per station.associated_station_count :: Int
: Number of stations at which the event was observed.used_station_count :: Int
: Number of stations from which data was used for origin computation.depth_phase_count :: Int
: Number of depth phases (typically pP, sometimes sP) used in depth computation.standard_error :: Float64
: RMS of the travel time residuals of the arrivals used for the origin computation. Unit: s.azimuthal_gap :: Float64
: Largest azimuthal gap in station distribution as seen from epicenter. For an illustration of azimuthal gap and secondary azimuthal gap (see below), see Fig. 5 of Bondár et al. (2004). Unit: °.secondary_azimuthal_gap :: Float64
: Secondary azimuthal gap in station distribution, i. e., the largest azimuthal gap a station closes. Unit: °.ground_truth_level :: String
:String
describing ground-truth level, e. g. GT0, GT5, etc. It has a maximum length of 32 characters.minimum_distance :: Float64
: Epicentral distance of station closest to the epicenter. Unit: °.maximum_distance :: Float64
: Epicentral distance of station farthest from the epicenter. Unit: °.median_distance :: Float64
: Median epicentral distance of used stations. Unit: °.
QuakeML.OriginUncertainty
— TypeOriginUncertainty(; kwargs...)
Describes the location uncertainties of an origin. The uncertainty can be described either as a simple circular horizontal uncertainty, an uncertainty ellipse according to IMS1.0, or a confidence ellipsoid. If multiple uncertainty models are given, the preferred variant can be specified in the field preferred_description
.
List of fields
horizontal_uncertainty :: Float64
: Circular confidence region, given by single value of horizontal uncertainty. Unit: m.min_horizontal_uncertainty :: Float64
: Semi-minor axis of confidence ellipse. Unit: m.max_horizontal_uncertainty :: Float64
: Semi-major axis of confidence ellipse. Unit: m.azimuth_max_horizontal_uncertainty :: Float64
: Azimuth of major axis of confidence ellipse. Measured clockwise from south-north direction at epicenter. Unit: °.confidence_ellipsoid :: ConfidenceEllipsoid
: Confidence ellipsoid (seeConfidenceEllipsoid
).preferred_description :: OriginUncertaintyDescription
: Preferred uncertainty description. Allowed values are the following (seeOriginUncertaintyDescription
:"horizontal uncertainty"
"uncertainty ellipse"
"confidence ellipsoid"
confidence_level :: Float64
: Confidence level of the uncertainty, given in percent.
QuakeML.Phase
— TypePhase(code)
Phase(; value=code)
Phase code as given in the IASPEI Standard Seismic Phase List (Storchak et al. 2003). String with a maximum length of 32 characters.
List of fields
value :: String
: Phase code. (Required field.)
QuakeML.Pick
— TypePick(; kwargs...)
A pick is the observation of an amplitude anomaly in a seismogram at a specific point in time. It is notnecessarily related to a seismic event.
List of fields
public_id :: ResourceReference
: Resource identifier ofPick
. (Required field.)time :: TimeQuantity
: Observed onset time of signal (“pick time”). (Required field.)waveform_id :: ResourceReference
: Identifes the waveform stream. (Required field.)filter_id :: ResourceReference
: Identifies the filter or filter setup used for filtering the waveform stream referenced bywaveform_id
.method_id :: ResourceReference
: Identifies the picker that produced the pick. This can be either a detection software program or aperson.horizontal_slowness :: RealQuantity
: Observed horizontal slowness of the signal. Most relevant in array measurements. Unit: s/°.backazimuth :: RealQuantity
: Observed backazimuth of the signal. Most relevant in array measurements. Unit: °.slowness_method_id :: ResourceReference
: Identifies the method that was used to determine the slowness.onset :: PickOnset
: Flag that roughly categorizes the sharpness of the onset. Allowed values are (seePickOnset
):"impulsive"
"emergent"
"questionable"
phase_hint :: Phase
: Tentative phase identification as specified by the picker.polarity :: PickPolarity
: Indicates the polarity of first motion, usually from impulsive onsets. Allowed values are (seePickPolarity
):"positive"
"negative"
"undecidable"
evaluation_mode :: EvaluationMode
: Evaluation mode ofPick
(seeEvaluationMode
).evaluation_status :: EvaluationStatus
: Evaluation status ofPick
(seeEvaluationStatus
).comment :: Vector{Comment}
: Additional comments.creation_info :: CreationInfo
:CreationInfo
for thePick
object.
QuakeML.PrincipleAxes
— TypePrincipleAxes(; t_axis, p_axis, n_axis)
List of fields
tAxis :: Axis
: T (tension) axis of a moment tensor. (Required field.)p_axis :: Axis
: P (pressure) axis of a moment tensor. (Required field.)n_axis :: Axis
: N (neutral) axis of a moment tensor.
QuakeML.RealQuantity
— TypeRealQuantity(; kwargs...)
Physical quantities that can be expressed numerically—either as integers or as floating point numbers—are represented by their measured or computed values and optional values for symmetric or upper and lower uncertainties. The interpretation of these uncertainties is not defined in the standard. They can contain statistically well-defined error measures, but the mechanism can also be used to simply describe a possible value range. Ifthe confidence level of the uncertainty is known, it can be listed in the optional field confidence_level
. Note that uncertainty
, upper_uncertainty
, and lower_uncertainty
are given as absolute values of the deviation from the main value
.
List of fields
value :: Float64
: Value of the quantity. The unit is implicitly defined and depends on the context. (Required field.)uncertainty :: Float64
: Uncertainty as the absolute value of symmetric deviation from the main value.lower_uncertainty :: Float64
: Uncertainty as the absolute value of deviation from the mainvalue
towards smaller values.upper_uncertainty :: Float64
: Uncertainty as the absolute value of deviation from the mainvalue
towards larger values.confidence_level :: Float64
: Confidence level of the uncertainty, given in percent.
QuakeML.ResourceReference
— TypeResourceReference(value)
ResourceReference(; value=::String)
String
that is used as a reference to a QuakeML resource. It must adhere to the format specificationsgiven in Sect. 3.1 of the QuakeML specificaiton. The string has a maximum length of 255 characters.
In this package, when creating objects which require a ResourceReference
(usually in a field called public_id
), a unique URI is created of the form "smi:local/XXXXXXXX-XXXX-XXXX-XXXX-XXXXXXXXXXXX"
, where X
represents a hexadecimal characters (matching r"[0-9a-f]"
). This is generated by calling QuakeML.random_reference
.
Further information
Identifiers take the generic form of:
[smi|quakeml]:〈authority-id〉/〈resource-id〉[#〈local-id〉]
They consist of an authority identifier, a unique resource identifier, and an optional local identifier. The URI schema name smi
stands for 'seismological meta-information', thus indicating a connection to a set of metadata associated with the resource.
The authority-id part must consist of at least three characters, of which the first character has to be alphanu-meric. The subsequent characters can be alphanumeric or from the following list: -
, .
, ~
, *
, '
, (
, )
. After the authority-id, a forward slash ("/"
) must follow which separates the authority-id from the resource-id. The resource-id must contain at least one character, which can be either alphanumeric, or from the eight special characters which are allowed for the authority-id. For the remaining characters of the resource-id, also the comma (","
) and semicolon (";"
) characters and characters from the following list can be used: +
, ?
, =
, #
, /
, &
. Note that the forward slash which separates authority-id and resource-id is always the first forwards lash in the resource identifier. The resource-id may be followed by a stop character ("#"
) and a local identifier which can be made up of alphanumeric characters, the comma (","
) and semicolon (";"
) characters, and the characters from the following list: -
, .
, ~
, *
, '
, (
, )
, /
, +
, =
, ?
. Local identifiers are thought to denote resources that have no own metadata description associated, but are part of a larger collection for which such metadata exists.
For even more information, see Section 3.1 of the QuakeML specification
ResourceReference
s are also called ResourceIdentifier
s.
QuakeML.SourceTimeFunction
— TypeSourceTimeFunction(; type, duration, rise_time, decay_time)
Source time function used in moment-tensor inversion.
List of fields
type :: SourceTimeFunctionType
: Type of source time function. Values can be taken from the following:"box car"
"triangle"
"trapezoid"
"unknown"
duration :: Float64
Source time function duration. Unit: s. (Required field.)rise_time :: Float64
: Source time function rise time. Unit: s.decay_time :: Float64
: Source time function decay time. Unit: s.
QuakeML.StationMagnitude
— TypeStationMagnitude(; kwargs...)
Describes the magnitude derived from a single waveform stream.
List of fields
public_id :: ResourceReference
: Resource identifier ofStationMagnitude
. (Required field.)origin_id :: ResourceReference
: Reference to an origin’spublic_id
if theStationMagnitude
has anassociatedOrigin
.mag :: RealQuantity
: Estimated magnitude. (Required field.)type :: String
: SeeMagnitude
.amplitude_id :: ResourceReference
: Identifies the data source of theStationMagnitude
. For magnitudes derived from amplitudes in waveforms (e. g., local magnitude ML),amplitude_id
points topublic_id
inAmplitude
.method_id :: ResourceReference
: SeeMagnitude
.waveform_id :: ResourceReference
: Identifies the waveform stream. This element can be helpful if no amplitude is referenced, or the amplitude is not available in the context. Otherwise, it would duplicate thewaveform_id
provided there and can be omitted.comment :: Vector{Comment}
: Additional comments.creationInfo :: CreationInfo
:CreationInfo
for theStationMagnitude
object.
QuakeML.StationMagnitudeContribution
— TypeStationMagnitudeContribution(; station_magnitude_id, residual, weight)
Describes the weighting of magnitude values froms everal StationMagnitude
objects for computing a network magnitude estimation.
List of fields
stationMagnitudeID :: ResourceReference
: Refers to thepublicID
of aStationMagnitude
object. (Required field.)residual :: Float64
: Residual of magnitude computation.weight :: Float64
: Weight of the magnitude value fromStationMagnitude
for computing the magnitude value inMagnitude
. Note that there is no rule for the sum of the weights of all station magnitude contributions to a specific network magnitude. In particular, the weights are not required to sum up to unity.
QuakeML.Tensor
— TypeTensor(mrr, mtt, mpp, mrt, mrp, mtp)
Tensor(; mrr, mtt, mpp, mrt, mrp, mtp)
The Tensor
type represents the six moment-tensor elements Mrr, Mtt, Mpp, Mrt, Mrp, Mtp in the spherical coordinate system defined by local upward vertical (r), North-South (t), and West-East (p) directions. See Aki and Richards(1980) for conversions to other coordinate systems.
List of fields
mrr :: RealQuantity
: Moment-tensor component Mrr. Unit: N m. (Required field.)mtt :: RealQuantity
: Moment-tensor component Mtt. Unit: N m. (Required field.)mpp :: RealQuantity
: Moment-tensor component Mpp. Unit: N m. (Required field.)mrt :: RealQuantity
: Moment-tensor component Mrt. Unit: N m. (Required field.)mrp :: RealQuantity
: Moment-tensor component Mrp. Unit: N m. (Required field.)mtp :: RealQuantity
: Moment-tensor component Mtp. Unit: N m. (Required field.)
QuakeML.TimeQuantity
— TypeTimeQuantity(; kwargs...)
Describes a point in time, given in ISO 8601 format, with optional symmetric or asymmetric uncertainties given in seconds. The time has to be specified in UTC.
List of fields
value :: Dates.DateTime
: Point in time (UTC), given in ISO 8601 format. (Required field.)uncertainty :: Float64
: Symmetric uncertainty of point in time. Unit: s.lower_uncertainty :: Float64
: Lower uncertainty of point in time. Unit: s.upper_uncertainty :: Float64
: Upper uncertainty of point in time. Unit: s.confidence_level :: Float64
: Confidence level of the uncertainty, given in percent.
QuakeML.TimeWindow
— TypeTimeWindow(; begin_, end_, reference)
Describes a time window for amplitude measurements, given by a central point in time, and points in time before and after this central point. Both points before and after may coincide with the central point.
List of fields
begin_ :: Float64
: Absolute value of duration of time interval beforereference
point in time window. The value may be zero, but not negative. Unit: s. (Required field.)end_ :: Float64
: Absolute value of duration of time interval afterreference
point in time window. The value may be zero, but not negative. Unit: s. (Required field.)reference :: Dates.DateTime
: Reference point in time (“central” point). It has to be given in UTC. (Required field.)
QuakeML.WaveformStreamID
— TypeWaveformStreamID(; kwargs...)
Reference to a stream description in an inventory. This is mostly equivalent to the combination of network_code
, station_code
, location_code
, and channel_code
. However, additional information, e. g., sampling rate, can be referenced by the resource uri
. It is recommended to use resource URI as a flexible, abstract, and unique stream ID that allows to describe different processing levels, or resampled/filtered products of the same initialstream, without violating the intrinsic meaning of the legacy identifiers (network, station, channel, and location codes). However, for operation in the context of legacy systems, the classical identifier components are upported.
List of fields
network_code :: String
: Network code. String with a maximum length of 8 characters. (Required field.)station_code :: String
: Station code. String with a maximum length of 8 characters. (Required field.)channel_code :: String
: Channel code. String with a maximum length of 8 characters.location_code :: String
: Location code. String with a maximum length of 8 characters.uri :: ResourceReference
: Resource identifier for the waveform stream.
Private types and functions
ID generation
QuakeML.random_reference
— Functionrandom_reference() -> ::ResourceReference
Create a new, random ResourceReference
.
Enumerated types
QuakeML.AmplitudeCategory
— TypeQuakeML.AmplitudeCategory(value)
QuakeML.AmplitudeCategory(; value)
Enumerated struct containing a single string which must be one of the following: "point"
, "mean"
, "duration"
, "period"
, "integral"
or "other"
.
Note that when a field of another type is a QuakeML.AmplitudeCategory
, it is not necessary to assign a field of type QuakeML.AmplitudeCategory
to the field. Instead, one can simply use a String
, from which a QuakeML.AmplitudeCategory
will be automatically constructed.
For this reason, QuakeML.AmplitudeCategory is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.AmplitudeCategory
end
julia> es = ExampleStruct("point")
ExampleStruct(QuakeML.AmplitudeCategory("point"))
julia> es.field = "mean"
"mean"
QuakeML.AmplitudeUnit
— TypeQuakeML.AmplitudeUnit(value)
QuakeML.AmplitudeUnit(; value)
Enumerated struct containing a single string which must be one of the following: "m"
, "s"
, "m/s"
, "m/(s*s)"
, "m*s"
, "dimensionless"
or "other"
.
Note that when a field of another type is a QuakeML.AmplitudeUnit
, it is not necessary to assign a field of type QuakeML.AmplitudeUnit
to the field. Instead, one can simply use a String
, from which a QuakeML.AmplitudeUnit
will be automatically constructed.
For this reason, QuakeML.AmplitudeUnit is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.AmplitudeUnit
end
julia> es = ExampleStruct("m")
ExampleStruct(QuakeML.AmplitudeUnit("m"))
julia> es.field = "s"
"s"
QuakeML.DataUsedWaveType
— TypeQuakeML.DataUsedWaveType(value)
QuakeML.DataUsedWaveType(; value)
Enumerated struct containing a single string which must be one of the following: "P waves"
, "body waves"
, "surface waves"
, "mantle waves"
, "combined"
or "unknown"
.
Note that when a field of another type is a QuakeML.DataUsedWaveType
, it is not necessary to assign a field of type QuakeML.DataUsedWaveType
to the field. Instead, one can simply use a String
, from which a QuakeML.DataUsedWaveType
will be automatically constructed.
For this reason, QuakeML.DataUsedWaveType is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.DataUsedWaveType
end
julia> es = ExampleStruct("P waves")
ExampleStruct(QuakeML.DataUsedWaveType("P waves"))
julia> es.field = "body waves"
"body waves"
QuakeML.EvaluationMode
— TypeQuakeML.EvaluationMode(value)
QuakeML.EvaluationMode(; value)
Enumerated struct containing a single string which must be one of the following: "manual"
or "automatic"
.
Note that when a field of another type is a QuakeML.EvaluationMode
, it is not necessary to assign a field of type QuakeML.EvaluationMode
to the field. Instead, one can simply use a String
, from which a QuakeML.EvaluationMode
will be automatically constructed.
For this reason, QuakeML.EvaluationMode is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.EvaluationMode
end
julia> es = ExampleStruct("manual")
ExampleStruct(QuakeML.EvaluationMode("manual"))
julia> es.field = "automatic"
"automatic"
QuakeML.EvaluationStatus
— TypeQuakeML.EvaluationStatus(value)
QuakeML.EvaluationStatus(; value)
Enumerated struct containing a single string which must be one of the following: "preliminary"
, "confirmed"
, "reviewed"
, "final"
or "rejected"
.
Note that when a field of another type is a QuakeML.EvaluationStatus
, it is not necessary to assign a field of type QuakeML.EvaluationStatus
to the field. Instead, one can simply use a String
, from which a QuakeML.EvaluationStatus
will be automatically constructed.
For this reason, QuakeML.EvaluationStatus is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.EvaluationStatus
end
julia> es = ExampleStruct("preliminary")
ExampleStruct(QuakeML.EvaluationStatus("preliminary"))
julia> es.field = "confirmed"
"confirmed"
QuakeML.EventDescriptionType
— TypeQuakeML.EventDescriptionType(value)
QuakeML.EventDescriptionType(; value)
Enumerated struct containing a single string which must be one of the following: "felt report"
, "Flinn-Engdahl region"
, "local time"
, "tectonic summary"
, "nearest cities"
, "earthquake name"
or "region name"
.
Note that when a field of another type is a QuakeML.EventDescriptionType
, it is not necessary to assign a field of type QuakeML.EventDescriptionType
to the field. Instead, one can simply use a String
, from which a QuakeML.EventDescriptionType
will be automatically constructed.
For this reason, QuakeML.EventDescriptionType is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.EventDescriptionType
end
julia> es = ExampleStruct("felt report")
ExampleStruct(QuakeML.EventDescriptionType("felt report"))
julia> es.field = "Flinn-Engdahl region"
"Flinn-Engdahl region"
QuakeML.EventType
— TypeQuakeML.EventType(value)
QuakeML.EventType(; value)
Enumerated struct containing a single string which must be one of the following: "not existing"
, "not reported"
, "earthquake"
, "anthropogenic event"
, "collapse"
, "cavity collapse"
, "mine collapse"
, "building collapse"
, "explosion"
, "accidental explosion"
, "chemical explosion"
, "controlled explosion"
, "experimental explosion"
, "industrial explosion"
, "mining explosion"
, "quarry blast"
, "road cut"
, "blasting levee"
, "nuclear explosion"
, "induced or triggered event"
, "rock burst"
, "reservoir loading"
, "fluid injection"
, "fluid extraction"
, "crash"
, "plane crash"
, "train crash"
, "boat crash"
, "other event"
, "atmospheric event"
, "sonic boom"
, "sonic blast"
, "acoustic noise"
, "thunder"
, "avalanche"
, "snow avalanche"
, "debris avalanche"
, "hydroacoustic event"
, "ice quake"
, "slide"
, "landslide"
, "rockslide"
, "meteorite"
or "volcanic eruption"
.
Note that when a field of another type is a QuakeML.EventType
, it is not necessary to assign a field of type QuakeML.EventType
to the field. Instead, one can simply use a String
, from which a QuakeML.EventType
will be automatically constructed.
For this reason, QuakeML.EventType is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.EventType
end
julia> es = ExampleStruct("not existing")
ExampleStruct(QuakeML.EventType("not existing"))
julia> es.field = "not reported"
"not reported"
QuakeML.EventTypeCertainty
— TypeQuakeML.EventTypeCertainty(value)
QuakeML.EventTypeCertainty(; value)
Enumerated struct containing a single string which must be one of the following: "known"
or "suspected"
.
Note that when a field of another type is a QuakeML.EventTypeCertainty
, it is not necessary to assign a field of type QuakeML.EventTypeCertainty
to the field. Instead, one can simply use a String
, from which a QuakeML.EventTypeCertainty
will be automatically constructed.
For this reason, QuakeML.EventTypeCertainty is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.EventTypeCertainty
end
julia> es = ExampleStruct("known")
ExampleStruct(QuakeML.EventTypeCertainty("known"))
julia> es.field = "suspected"
"suspected"
QuakeML.MomentTensorCategory
— TypeQuakeML.MomentTensorCategory(value)
QuakeML.MomentTensorCategory(; value)
Enumerated struct containing a single string which must be one of the following: "teleseismic"
or "regional"
.
Note that when a field of another type is a QuakeML.MomentTensorCategory
, it is not necessary to assign a field of type QuakeML.MomentTensorCategory
to the field. Instead, one can simply use a String
, from which a QuakeML.MomentTensorCategory
will be automatically constructed.
For this reason, QuakeML.MomentTensorCategory is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.MomentTensorCategory
end
julia> es = ExampleStruct("teleseismic")
ExampleStruct(QuakeML.MomentTensorCategory("teleseismic"))
julia> es.field = "regional"
"regional"
QuakeML.MTInversionType
— TypeQuakeML.MTInversionType(value)
QuakeML.MTInversionType(; value)
Enumerated struct containing a single string which must be one of the following: "general"
, "zero trace"
or "double couple"
.
Note that when a field of another type is a QuakeML.MTInversionType
, it is not necessary to assign a field of type QuakeML.MTInversionType
to the field. Instead, one can simply use a String
, from which a QuakeML.MTInversionType
will be automatically constructed.
For this reason, QuakeML.MTInversionType is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.MTInversionType
end
julia> es = ExampleStruct("general")
ExampleStruct(QuakeML.MTInversionType("general"))
julia> es.field = "zero trace"
"zero trace"
QuakeML.OriginDepthType
— TypeQuakeML.OriginDepthType(value)
QuakeML.OriginDepthType(; value)
Enumerated struct containing a single string which must be one of the following: "from location"
, "from moment tensor inversion"
, "from modeling of broad-band P waveforms"
, "constrained by depth phases"
, "constrained by direct phases"
, "constrained by depth and direct phases"
, "operator assigned"
or "other"
.
Note that when a field of another type is a QuakeML.OriginDepthType
, it is not necessary to assign a field of type QuakeML.OriginDepthType
to the field. Instead, one can simply use a String
, from which a QuakeML.OriginDepthType
will be automatically constructed.
For this reason, QuakeML.OriginDepthType is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.OriginDepthType
end
julia> es = ExampleStruct("from location")
ExampleStruct(QuakeML.OriginDepthType("from location"))
julia> es.field = "from moment tensor inversion"
"from moment tensor inversion"
QuakeML.OriginType
— TypeQuakeML.OriginType(value)
QuakeML.OriginType(; value)
Enumerated struct containing a single string which must be one of the following: "hypocenter"
, "centroid"
, "amplitude"
, "macroseismic"
, "rupture start"
or "rupture end"
.
Note that when a field of another type is a QuakeML.OriginType
, it is not necessary to assign a field of type QuakeML.OriginType
to the field. Instead, one can simply use a String
, from which a QuakeML.OriginType
will be automatically constructed.
For this reason, QuakeML.OriginType is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.OriginType
end
julia> es = ExampleStruct("hypocenter")
ExampleStruct(QuakeML.OriginType("hypocenter"))
julia> es.field = "centroid"
"centroid"
QuakeML.OriginUncertaintyDescription
— TypeQuakeML.OriginUncertaintyDescription(value)
QuakeML.OriginUncertaintyDescription(; value)
Enumerated struct containing a single string which must be one of the following: "horizontal uncertainty"
, "uncertainty ellipse"
or "confidence ellipsoid"
.
Note that when a field of another type is a QuakeML.OriginUncertaintyDescription
, it is not necessary to assign a field of type QuakeML.OriginUncertaintyDescription
to the field. Instead, one can simply use a String
, from which a QuakeML.OriginUncertaintyDescription
will be automatically constructed.
For this reason, QuakeML.OriginUncertaintyDescription is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.OriginUncertaintyDescription
end
julia> es = ExampleStruct("horizontal uncertainty")
ExampleStruct(QuakeML.OriginUncertaintyDescription("horizontal uncertainty"))
julia> es.field = "uncertainty ellipse"
"uncertainty ellipse"
QuakeML.PickOnset
— TypeQuakeML.PickOnset(value)
QuakeML.PickOnset(; value)
Enumerated struct containing a single string which must be one of the following: "emergent"
, "impulsive"
or "questionable"
.
Note that when a field of another type is a QuakeML.PickOnset
, it is not necessary to assign a field of type QuakeML.PickOnset
to the field. Instead, one can simply use a String
, from which a QuakeML.PickOnset
will be automatically constructed.
For this reason, QuakeML.PickOnset is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.PickOnset
end
julia> es = ExampleStruct("emergent")
ExampleStruct(QuakeML.PickOnset("emergent"))
julia> es.field = "impulsive"
"impulsive"
QuakeML.PickPolarity
— TypeQuakeML.PickPolarity(value)
QuakeML.PickPolarity(; value)
Enumerated struct containing a single string which must be one of the following: "positive"
, "negative"
or "undecidable"
.
Note that when a field of another type is a QuakeML.PickPolarity
, it is not necessary to assign a field of type QuakeML.PickPolarity
to the field. Instead, one can simply use a String
, from which a QuakeML.PickPolarity
will be automatically constructed.
For this reason, QuakeML.PickPolarity is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.PickPolarity
end
julia> es = ExampleStruct("positive")
ExampleStruct(QuakeML.PickPolarity("positive"))
julia> es.field = "negative"
"negative"
QuakeML.SourceTimeFunctionType
— TypeQuakeML.SourceTimeFunctionType(value)
QuakeML.SourceTimeFunctionType(; value)
Enumerated struct containing a single string which must be one of the following: "box car"
, "triangle"
, "trapezoid"
or "unknown"
.
Note that when a field of another type is a QuakeML.SourceTimeFunctionType
, it is not necessary to assign a field of type QuakeML.SourceTimeFunctionType
to the field. Instead, one can simply use a String
, from which a QuakeML.SourceTimeFunctionType
will be automatically constructed.
For this reason, QuakeML.SourceTimeFunctionType is not exported even when bringing QuakeML's types into scope by doing using QuakeML.Types
.
Example
julia> using QuakeML
julia> mutable struct ExampleStruct
field::QuakeML.SourceTimeFunctionType
end
julia> es = ExampleStruct("box car")
ExampleStruct(QuakeML.SourceTimeFunctionType("box car"))
julia> es.field = "triangle"
"triangle"