gf
¶
Pyrocko-GF: storage and calculation of synthetic seismograms
The pyrocko.gf
subpackage splits functionality into several submodules:
The
pyrocko.gf.store
module deals with the storage, retrieval and summation of Green’s functions.The
pyrocko.gf.meta
module provides data structures for the meta information associated with the Green’s function stores and implements various the Green’s function lookup indexing schemes.The
pyrocko.gf.builder
module defines a common base for Green’s function store builders.The
pyrocko.gf.seismosizer
module provides high level synthetic seismogram synthesis.
All classes defined in the pyrocko.gf.*
submodules are imported into the
pyrocko.gf
namespace, so user scripts may simply use from pyrocko
import gf
or from pyrocko.gf import *
for convenience.
store
¶
- class GFTrace(data=None, itmin=None, deltat=1.0, is_zero=False, begin_value=None, end_value=None, tmin=None)[source]¶
Bases:
object
Green’s Function trace class for handling traces from the GF store.
- property t¶
Time vector of the GF trace.
- Returns
Time vector
- Return type
- class Store(store_dir, mode='r', use_memmap=True)[source]¶
Bases:
pyrocko.gf.store.BaseStore
Green’s function disk storage and summation machine.
The Store can be used to efficiently store, retrieve, and sum Green’s function traces. A Store contains many 1D time traces sampled at even multiples of a global sampling rate, where each time trace has an individual start and end time. The traces are treated as having repeating end points, so the functions they represent can be non-constant only between begin and end time. It provides capabilities to retrieve decimated traces and to extract parts of the traces. The main purpose of this class is to provide a fast, easy to use, and flexible machanism to compute weighted delay-and-sum stacks with many Green’s function traces involved.
Individual Green’s functions are accessed through a single integer index at low level. On top of that, various indexing schemes can be implemented by providing a mapping from physical coordinates to the low level index i. E.g. for a problem with cylindrical symmetry, one might define a mapping from the coordinates (receiver_depth, source_depth, distance) to the low level index. Index translation is done in the
pyrocko.gf.meta.Config
subclass object associated with the Store. It is accessible through the store’sconfig
attribute, and contains all meta information about the store, including physical extent, geometry, sampling rate, and information about the type of the stored Green’s functions. Information about the underlying earth model can also be found there.A GF store can also contain tabulated phase arrivals. In basic cases, these can be created with the
make_travel_time_tables()
and evaluated with thet()
methods.- config¶
The
pyrocko.gf.meta.Config
derived object associated with the store which contains all meta information about the store and provides the high-level to low-level index mapping.
- store_dir¶
Path to the store’s data directory.
- mode¶
The mode in which the store is opened:
'r'
: read-only,'w'
: writeable.
- classmethod create(store_dir, config, force=False, extra=None)[source]¶
Create new GF store.
Creates a new GF store at path
store_dir
. The layout of the GF is defined with the parameters given inconfig
, which should be an object of a subclass ofConfig
. This function will refuse to overwrite an existing GF store, unlessforce
is set toTrue
. If more information, e.g. parameters used for the modelling code, earth models or other, should be saved along with the GF store, these may be provided though a dict given toextra
. The keys of this dict must be names and the values must be guts type objects.
- put(args, trace)[source]¶
Insert trace into GF store.
Store a single GF trace at (high-level) index
args
.- Parameters
args (tuple) –
Config
index tuple, e.g.(source_depth, distance, component)
as inConfigTypeA
.- Returns
GF trace at
args
- Return type
- get(args, itmin=None, nsamples=None, decimate=1, interpolation='nearest_neighbor', implementation='c')[source]¶
Retrieve GF trace from store.
Retrieve a single GF trace from the store at (high-level) index
args
. By default, the full trace is retrieved. Givenitmin
andnsamples
, only the selected portion of the trace is extracted. Ifdecimate
is an integer in the range [2,8], the trace is decimated on the fly or, if available, the trace is read from a decimated version of the GF store.- Parameters
args (tuple) –
Config
index tuple, e.g.(source_depth, distance, component)
as inConfigTypeA
.itmin (int or None) – Start time index (start time is
itmin * dt
), defaults to Nonensamples (int or None) – Number of samples, defaults to None
decimate (int) – Decimatation factor, defaults to 1
interpolation (str) – Interpolation method
['nearest_neighbor', 'multilinear', 'off']
, defaults to'nearest_neighbor'
implementation (str) – Implementation to use
['c', 'reference']
, defaults to'c'
.
- Returns
GF trace at
args
- Return type
- sum(args, delays, weights, itmin=None, nsamples=None, decimate=1, interpolation='nearest_neighbor', implementation='c', optimization='enable')[source]¶
Sum delayed and weighted GF traces.
Calculate sum of delayed and weighted GF traces.
args
is a tuple of arrays forming the (high-level) indices of the GF traces to be selected. Delays and weights for the summation are given in the arraysdelays
andweights
.itmin
andnsamples
can be given to restrict to computation to a given time interval. Ifdecimate
is an integer in the range [2,8], decimated traces are used in the summation.- Parameters
args (tuple(numpy.ndarray)) –
Config
index tuple, e.g.(source_depth, distance, component)
as inConfigTypeA
.delays (
numpy.ndarray
) – Delay timesweights (
numpy.ndarray
) – Trace weightsitmin (int or None) – Start time index (start time is
itmin * dt
), defaults to Nonensamples (int or None) – Number of samples, defaults to None
decimate (int) – Decimatation factor, defaults to 1
interpolation (str) – Interpolation method
['nearest_neighbor', 'multilinear', 'off']
, defaults to'nearest_neighbor'
implementation (str) – Implementation to use,
['c', 'alternative', 'reference']
, where'alternative'
and'reference'
use a Python implementation, defaults to ‘c’optimization (str) – Optimization mode
['enable', 'disable']
, defaults to'enable'
- Returns
Stacked GF trace.
- Return type
- make_decimated(decimate, config=None, force=False, show_progress=False)[source]¶
Create decimated version of GF store.
Create a downsampled version of the GF store. Downsampling is done for the integer factor
decimate
which should be in the range [2,8]. Ifconfig
isNone
, all traces of the GF store are decimated and held available (i.e. the index mapping of the original store is used), otherwise, a different spacial stepping can be specified by giving a modified GF store configuration inconfig
(seecreate()
). Decimated GF sub-stores are created under thedecimated
subdirectory within the GF store directory. Holding available decimated versions of the GF store can save computation time, IO bandwidth, or decrease memory footprint at the cost of increased disk space usage, when computation are done for lower frequency signals.
- get_stored_phase(phase_def, attribute='phase')[source]¶
Get stored phase from GF store.
- Returns
Phase information
- Return type
pyrocko.spit.SPTree
- t(timing, *args)[source]¶
Compute interpolated phase arrivals.
Examples:
If
test_store
is ofConfigTypeA
:test_store.t('p', (1000, 10000)) test_store.t('last{P|Pdiff}', (1000, 10000)) # The latter arrival # of P or diffracted # P phase
If
test_store
is ofConfigTypeB
:test_store.t('S', (1000, 1000, 10000)) test_store.t('first{P|p|Pdiff|sP}', (1000, 1000, 10000)) # The ` # first arrival of # the given phases is # selected
- Parameters
timing (str or
Timing
) – Timing string as described above*args (tuple) –
Config
index tuple, e.g.(source_depth, distance, component)
as inConfigTypeA
.
- Returns
Phase arrival according to
timing
- Return type
float or None
- get_stored_attribute(phase_def, attribute, *args)[source]¶
Return interpolated store attribute
- Parameters
attribute (str) – takeoff_angle / incidence_angle [deg]
*args (tuple) –
Config
index tuple, e.g.(source_depth, distance, component)
as inConfigTypeA
.
- get_many_stored_attributes(phase_def, attribute, coords)[source]¶
Return interpolated store attribute
- Parameters
attribute (str) – takeoff_angle / incidence_angle [deg]
coords (
num.array.Array
) –num.array.Array
, with columns being(source_depth, distance, component)
as inConfigTypeA
.
- make_stored_table(attribute, force=False)[source]¶
Compute tables for selected ray attributes.
- Parameters
attribute (str) – phase / takeoff_angle [deg]/ incidence_angle [deg]
Tables are computed using the 1D earth model defined in
earthmodel_1d
for each defined phase intabulated_phases
.
- make_timing_params(begin, end, snap_vred=True, force=False)[source]¶
Compute tight parameterized time ranges to include given timings.
Calculates appropriate time ranges to cover given begin and end timings over all GF points in the store. A dict with the following keys is returned:
'tmin'
: time [s], minimum of begin timing over all GF points'tmax'
: time [s], maximum of end timing over all GF points'vred'
,'tmin_vred'
: slope [m/s] and offset [s] of reduction velocity [m/s] appropriate to catch begin timing over all GF points'tlenmax_vred'
: maximum time length needed to cover all end timings, when using linear slope given with (vred
,tmin_vred
) as start
- make_travel_time_tables(force=False)[source]¶
Compute travel time tables.
Travel time tables are computed using the 1D earth model defined in
earthmodel_1d
for each defined phase intabulated_phases
. The accuracy of the tablulated times is adjusted to the sampling rate of the store.
- make_takeoff_angle_tables(force=False)[source]¶
Compute takeoff-angle tables.
Takeoff-angle tables [deg] are computed using the 1D earth model defined in
earthmodel_1d
for each defined phase intabulated_phases
. The accuracy of the tablulated times is adjusted to 0.01 times the sampling rate of the store.
- make_incidence_angle_tables(force=False)[source]¶
Compute incidence-angle tables.
Incidence-angle tables [deg] are computed using the 1D earth model defined in
earthmodel_1d
for each defined phase intabulated_phases
. The accuracy of the tablulated times is adjusted to 0.01 times the sampling rate of the store.
builder
¶
seismosizer
¶
Coordinate systems¶
Coordinate system names commonly used in source models.
Name |
Description |
---|---|
|
northing, easting, depth in [m] |
|
northing, easting in [m] |
|
latitude, longitude in [deg] |
|
longitude, latitude in [deg] |
|
latitude, longitude in [deg], depth in [m] |
- class STFMode(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringChoice
Any
str
out of['pre', 'post']
.
- class Source(**kwargs)[source]¶
Bases:
pyrocko.model.location.Location
,pyrocko.gf.seismosizer.Cloneable
Base class for all source models.
- ♦ name¶
str
, optional, default:''
- ♦ time¶
time_float (
pyrocko.guts.Timestamp
), default:str_to_time('1970-01-01 00:00:00')
source origin time.
- ♦ stf_mode¶
str
(STFMode
), default:'post'
whether to apply source time function in pre or post-processing.
- update(**kwargs)[source]¶
Change some of the source models parameters.
Example:
>>> from pyrocko import gf >>> s = gf.DCSource() >>> s.update(strike=66., dip=33.) >>> print(s) --- !pf.DCSource depth: 0.0 time: 1970-01-01 00:00:00 magnitude: 6.0 strike: 66.0 dip: 33.0 rake: 0.0
- grid(**variables)[source]¶
Create grid of source model variations.
- Returns
SourceGrid
instance.
Example:
>>> from pyrocko import gf >>> base = DCSource() >>> R = gf.Range >>> for s in base.grid(R('
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- effective_stf_pre()[source]¶
Return the STF applied before stacking of the Green’s functions.
This STF is used during discretization of the parameterized source models, i.e. to produce a temporal distribution of point sources.
Handling of the STF before stacking of the GFs is less efficient but allows to use different source time functions for different parts of the source.
- effective_stf_post()[source]¶
Return the STF applied after stacking of the Green’s fuctions.
This STF is used in the post-processing of the synthetic seismograms.
Handling of the STF after stacking of the GFs is usually more efficient but is only possible when a common STF is used for all subsources.
- class SourceWithMagnitude(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.Source
Base class for sources containing a moment magnitude.
- ♦ magnitude¶
float
, default:6.0
Moment magnitude Mw as in [Hanks and Kanamori, 1979]
- class SourceWithDerivedMagnitude(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.Source
Undocumented.
- check_conflicts()[source]¶
Check for parameter conflicts.
To be overloaded in subclasses. Raises
DerivedMagnitudeError
on conflicts.
- class ExplosionSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.SourceWithDerivedMagnitude
An isotropic explosion point source.
- ♦ magnitude¶
float
, optionalmoment magnitude Mw as in [Hanks and Kanamori, 1979]
- ♦ volume_change¶
float
, optionalvolume change of the explosion/implosion or the contracting/extending magmatic source. [m^3]
- discretized_source_class¶
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- check_conflicts()[source]¶
Check for parameter conflicts.
To be overloaded in subclasses. Raises
DerivedMagnitudeError
on conflicts.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- class RectangularExplosionSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.ExplosionSource
Rectangular or line explosion source.
- ♦ strike¶
float
, default:0.0
strike direction in [deg], measured clockwise from north
- ♦ dip¶
float
, default:90.0
dip angle in [deg], measured downward from horizontal
- ♦ length¶
float
, default:0.0
length of rectangular source area [m]
- ♦ width¶
float
, default:0.0
width of rectangular source area [m]
- ♦ anchor¶
str
(pyrocko.guts.StringChoice
), optional, default:'center'
Anchor point for positioning the plane, can be: top, center orbottom and also top_left, top_right,bottom_left,bottom_right, center_left and center right
- ♦ nucleation_x¶
float
, optionalhorizontal position of rupture nucleation in normalized fault plane coordinates (-1 = left edge, +1 = right edge)
- ♦ nucleation_y¶
float
, optionaldown-dip position of rupture nucleation in normalized fault plane coordinates (-1 = upper edge, +1 = lower edge)
- ♦ velocity¶
float
, default:3500.0
speed of explosion front [m/s]
- ♦ aggressive_oversampling¶
bool
, default:False
Aggressive oversampling for basesource discretization. When using ‘multilinear’ interpolation oversampling has practically no effect.
- discretized_source_class¶
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- class DCSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.SourceWithMagnitude
A double-couple point source.
- ♦ strike¶
float
, default:0.0
strike direction in [deg], measured clockwise from north
- ♦ dip¶
float
, default:90.0
dip angle in [deg], measured downward from horizontal
- ♦ rake¶
float
, default:0.0
rake angle in [deg], measured counter-clockwise from right-horizontal in on-plane view
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedMTSource
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- class CLVDSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.SourceWithMagnitude
A pure CLVD point source.
- ♦ azimuth¶
float
, default:0.0
azimuth direction of largest dipole, clockwise from north [deg]
- ♦ dip¶
float
, default:90.0
dip direction of largest dipole, downward from horizontal [deg]
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedMTSource
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- class VLVDSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.SourceWithDerivedMagnitude
Volumetric linear vector dipole source.
This source is a parameterization for a restricted moment tensor point source, useful to represent dyke or sill like inflation or deflation sources. The restriction is such that the moment tensor is rotational symmetric. It can be represented by a superposition of a linear vector dipole (here we use a CLVD for convenience) and an isotropic component. The restricted moment tensor has 4 degrees of freedom: 2 independent eigenvalues and 2 rotation angles orienting the the symmetry axis.
In this parameterization, the isotropic component is controlled by
volume_change
. To define the moment tensor, it must be converted to the scalar moment of the the MT’s isotropic component. For the conversion, the shear modulus at the source’s position must be known. This value is extracted from the earth model defined in the GF store in use.The CLVD part by controlled by its scalar moment :
clvd_moment
. The sign ofclvd_moment
is used to switch between a positiv or negativ CLVD (the sign of the largest eigenvalue).- ♦ azimuth¶
float
, default:0.0
azimuth direction of symmetry axis, clockwise from north [deg].
- ♦ dip¶
float
, default:90.0
dip direction of symmetry axis, downward from horizontal [deg].
- ♦ volume_change¶
float
, default:0.0
volume change of the inflation/deflation [m^3].
- ♦ clvd_moment¶
float
, default:0.0
scalar moment of the CLVD component [Nm]. The sign controls the sign of the CLVD (the sign of its largest eigenvalue).
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedMTSource
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- class MTSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.Source
A moment tensor point source.
- ♦ mnn¶
float
, default:1.0
north-north component of moment tensor in [Nm]
- ♦ mee¶
float
, default:1.0
east-east component of moment tensor in [Nm]
- ♦ mdd¶
float
, default:1.0
down-down component of moment tensor in [Nm]
- ♦ mne¶
float
, default:0.0
north-east component of moment tensor in [Nm]
- ♦ mnd¶
float
, default:0.0
north-down component of moment tensor in [Nm]
- ♦ med¶
float
, default:0.0
east-down component of moment tensor in [Nm]
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedMTSource
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- class RectangularSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.SourceWithDerivedMagnitude
Classical Haskell source model modified for bilateral rupture.
- ♦ magnitude¶
float
, optionalmoment magnitude Mw as in [Hanks and Kanamori, 1979]
- ♦ strike¶
float
, default:0.0
strike direction in [deg], measured clockwise from north
- ♦ dip¶
float
, default:90.0
dip angle in [deg], measured downward from horizontal
- ♦ rake¶
float
, default:0.0
rake angle in [deg], measured counter-clockwise from right-horizontal in on-plane view
- ♦ length¶
float
, default:0.0
length of rectangular source area [m]
- ♦ width¶
float
, default:0.0
width of rectangular source area [m]
- ♦ anchor¶
str
(pyrocko.guts.StringChoice
), optional, default:'center'
Anchor point for positioning the plane, can be:
top, center bottom, top_left, top_right,bottom_left,bottom_right, center_left, center right
.
- ♦ nucleation_x¶
float
, optionalhorizontal position of rupture nucleation in normalized fault plane coordinates (
-1.
= left edge,+1.
= right edge)
- ♦ nucleation_y¶
float
, optionaldown-dip position of rupture nucleation in normalized fault plane coordinates (
-1.
= upper edge,+1.
= lower edge)
- ♦ velocity¶
float
, default:3500.0
speed of rupture front [m/s]
- ♦ slip¶
float
, optionalSlip on the rectangular source area [m]
- ♦ opening_fraction¶
float
, default:0.0
Determines fraction of slip related to opening. (
-1
: pure tensile closing,0
: pure shear,1
: pure tensile opening)
- ♦ decimation_factor¶
int
, optional, default:1
Sub-source decimation factor, a larger decimation will make the result inaccurate but shorten the necessary computation time (use for testing puposes only).
- ♦ aggressive_oversampling¶
bool
, default:False
Aggressive oversampling for basesource discretization. When using ‘multilinear’ interpolation oversampling has practically no effect.
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedMTSource
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- check_conflicts()[source]¶
Check for parameter conflicts.
To be overloaded in subclasses. Raises
DerivedMagnitudeError
on conflicts.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- class PseudoDynamicRupture(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.SourceWithDerivedMagnitude
Combined Eikonal and Okada quasi-dynamic rupture model.
Details are described in Pseudo Dynamic Rupture - A stress-based self-similar finite source model. Note: attribute stf is not used so far, but kept for future applications.
- ♦ strike¶
float
, default:0.0
Strike direction in [deg], measured clockwise from north.
- ♦ dip¶
float
, default:0.0
Dip angle in [deg], measured downward from horizontal.
- ♦ length¶
float
, default:10000.0
Length of rectangular source area in [m].
- ♦ width¶
float
, default:5000.0
Width of rectangular source area in [m].
- ♦ anchor¶
str
(pyrocko.guts.StringChoice
), optional, default:'center'
Anchor point for positioning the plane, can be:
top, center, bottom, top_left, top_right, bottom_left, bottom_right, center_left, center_right
.
- ♦ nucleation_x¶
numpy.ndarray
(pyrocko.guts_array.Array
), default:array([0.])
Horizontal position of rupture nucleation in normalized fault plane coordinates (
-1.
= left edge,+1.
= right edge).
- ♦ nucleation_y¶
numpy.ndarray
(pyrocko.guts_array.Array
), default:array([0.])
Down-dip position of rupture nucleation in normalized fault plane coordinates (
-1.
= upper edge,+1.
= lower edge).
- ♦ nucleation_time¶
numpy.ndarray
(pyrocko.guts_array.Array
), optionalTime in [s] after origin, when nucleation points defined by
nucleation_x
andnucleation_y
rupture.
- ♦ gamma¶
float
, default:0.8
Scaling factor between rupture velocity and S-wave velocity: .
- ♦ nx¶
int
, default:2
Number of discrete source patches in x direction (along strike).
- ♦ ny¶
int
, default:2
Number of discrete source patches in y direction (down dip).
- ♦ slip¶
float
, optionalMaximum slip of the rectangular source [m]. Setting the slip the tractions/stress field will be normalized to accomodate the desired maximum slip.
- ♦ rake¶
float
, optionalRake angle in [deg], measured counter-clockwise from right-horizontal in on-plane view. Rake is translated into homogenous tractions in strike and up-dip direction.
rake
is mutually exclusive with tractions parameter.
- ♦ patches¶
list
ofpyrocko.modelling.okada.OkadaSource
objects, optionalList of all boundary elements/sub faults/fault patches.
- ♦ patch_mask¶
numpy.ndarray
(pyrocko.guts_array.Array
), optionalMask for all boundary elements/sub faults/fault patches. True leaves the patch in the calculation, False excludes the patch.
- ♦ tractions¶
pyrocko.gf.tractions.TractionField
, optionalTraction field the rupture plane is exposed to. See the:py:mod:pyrocko.gf.tractions module for more details. If
tractions=None
andrake
is givenDirectedTractions
will be used.
- ♦ coef_mat¶
numpy.ndarray
(pyrocko.guts_array.Array
), optionalCoefficient matrix linking traction and dislocation field.
- ♦ eikonal_decimation¶
int
, optional, default:1
Sub-source eikonal factor, a smaller eikonal factor will increase the accuracy of rupture front calculation but increases also the computation time.
- ♦ decimation_factor¶
int
, optional, default:1
Sub-source decimation factor, a larger decimation will make the result inaccurate but shorten the necessary computation time (use for testing puposes only).
- ♦ nthreads¶
int
, optional, default:1
Number of threads for Okada forward modelling, matrix inversion and calculation of point subsources. Note: for small/medium matrices 1 thread is most efficient.
- ♦ pure_shear¶
bool
, optional, default:False
Calculate only shear tractions and omit tensile tractions.
- ♦ smooth_rupture¶
bool
, default:True
Smooth the tractions by weighting partially ruptured fault patches.
- ♦ aggressive_oversampling¶
bool
, default:False
Aggressive oversampling for basesource discretization. When using ‘multilinear’ interpolation oversampling has practically no effect.
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedMTSource
- get_tractions()[source]¶
Get source traction vectors.
If
rake
is given, unit length directed traction vectors (DirectedTractions
) are returned, else the giventractions
are used.- Returns
Traction vectors per patch.
- Return type
ndarray
:(n_patches, 3)
.
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- check_conflicts()[source]¶
Check for parameter conflicts.
To be overloaded in subclasses. Raises
DerivedMagnitudeError
on conflicts.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- outline(cs='xyz')[source]¶
Get source outline corner coordinates.
- Parameters
cs (optional, str) – Output coordinate system.
- Returns
Corner points in desired coordinate system.
- Return type
ndarray
:(5, [2, 3])
.
- points_on_source(cs='xyz', **kwargs)[source]¶
Convert relative plane coordinates to geographical coordinates.
Given x and y coordinates (relative source coordinates between -1. and 1.) are converted to desired geographical coordinates. Coordinates need to be given as
ndarray
argumentspoints_x
andpoints_y
.- Parameters
cs (optional, str) – Output coordinate system.
- Returns
Point coordinates in desired coordinate system.
- Return type
ndarray
:(n_points, [2, 3])
.
- discretize_time(store, interpolation='nearest_neighbor', vr=None, times=None, *args, **kwargs)[source]¶
Get rupture start time for discrete points on source plane.
- Parameters
store (
Store
) – Green’s function database (needs to cover whole region of the source)interpolation (optional, str) – Interpolation method to use (choose between
'nearest_neighbor'
and'multilinear'
).vr (optional,
ndarray
) – Array, containing rupture user defined rupture velocity values.times (optional,
ndarray
) – Array, containing zeros, where rupture is starting, real positive numbers at later secondary nucleation points and -1, where time will be calculated. If not given, rupture starts at nucleation_x, nucleation_y. Times are given for discrete points with equal horizontal and vertical spacing.
- Returns
Coordinates (latlondepth), coordinates (xy), rupture velocity, rupture propagation time of discrete points.
- Return type
ndarray
:(n_points, 3)
,ndarray
:(n_points, 2)
,ndarray
:(n_points_dip, n_points_strike)
,ndarray
:(n_points_dip, n_points_strike)
.
- get_vr_time_interpolators(store, interpolation='nearest_neighbor', force=False, **kwargs)[source]¶
Get interpolators for rupture velocity and rupture time.
Additional
**kwargs
are passed todiscretize_time()
.- Parameters
store (
Store
) – Green’s function database (needs to cover whole region of the source).interpolation (optional, str) – Interpolation method to use (choose between
'nearest_neighbor'
and'multilinear'
).force (optional, bool) – Force recalculation of the interpolators (e.g. after change of nucleation point locations/times). Default is
False
.
- discretize_patches(store, interpolation='nearest_neighbor', force=False, grid_shape=(), **kwargs)[source]¶
Get rupture start time and OkadaSource elements for points on rupture.
All source elements and their corresponding center points are calculated and stored in the
patches
attribute.Additional
**kwargs
are passed todiscretize_time()
.- Parameters
store (
Store
) – Green’s function database (needs to cover whole region of the source).interpolation (optional, str) – Interpolation method to use (choose between
'nearest_neighbor'
and'multilinear'
).force (optional, bool) – Force recalculation of the vr and time interpolators ( e.g. after change of nucleation point locations/times). Default is
False
.grid_shape (optional, tuple of int) – Desired sub fault patch grid size (nlength, nwidth). Either factor or grid_shape should be set.
- discretize_basesource(store, target=None)[source]¶
Prepare source for synthetic waveform calculation.
- Parameters
- Returns
Source discretized by a set of moment tensors and times.
- Return type
- get_slip(time=None, scale_slip=True, interpolation='nearest_neighbor', **kwargs)[source]¶
Get slip per subfault patch for given time after rupture start.
- Parameters
time (optional, float > 0.) – Time after origin [s], for which slip is computed. If not given, final static slip is returned.
scale_slip (optional, bool) – If
True
andslip
given, all slip values are scaled to fit the given maximum slip.interpolation (optional, str) – Interpolation method to use (choose between
'nearest_neighbor'
and'multilinear'
).
- Returns
Inverted dislocations () for each source patch.
- Return type
ndarray
:(n_sources, 3)
- get_delta_slip(store=None, deltat=None, delta=True, interpolation='nearest_neighbor', **kwargs)[source]¶
Get slip change snapshots.
The time interval, within which the slip changes are computed is determined by the sampling rate of the Green’s function database or
deltat
. Additional**kwargs
are passed toget_slip()
.- Parameters
store (optional,
Store
) – Green’s function database (needs to cover whole region of of the source). Its sampling interval is used as time increment for slip difference calculation. Eitherdeltat
orstore
should be given.deltat (optional, float) – Time interval for slip difference calculation [s]. Either
deltat
orstore
should be given.delta (optional, bool) – If
True
, slip differences between two time steps are given. IfFalse
, cumulative slip for all time steps.interpolation (optional, str) – Interpolation method to use (choose between
'nearest_neighbor'
and'multilinear'
).
- Returns
Displacement changes() for each source patch and time; corner times, for which delta slip is computed. The order of displacement changes array is:
- Return type
- get_slip_rate(*args, **kwargs)[source]¶
Get slip rate inverted from patches.
The time interval, within which the slip rates are computed is determined by the sampling rate of the Green’s function database or
deltat
. Additional*args
and**kwargs
are passed toget_delta_slip()
.
- get_moment_rate_patches(*args, **kwargs)[source]¶
Get scalar seismic moment rate for each patch individually.
Additional
*args
and**kwargs
are passed toget_slip_rate()
.
- get_moment_rate(store, target=None, deltat=None)[source]¶
Get seismic source moment rate for the total source (STF).
- Parameters
store (
Store
) – Green’s function database (needs to cover whole region of of the source). Itsdeltat
[s] is used as time increment for slip difference calculation. Eitherdeltat
orstore
should be given.target (optional,
Target
) – Target information, needed for interpolation method.deltat (optional, float) – Time increment for slip difference calculation [s]. If not given
store.deltat
is used.
- Returns
Seismic moment rate [Nm/s] for each time; corner times, for which moment rate is computed. The order of the moment rate array is:
- Return type
- get_moment(*args, **kwargs)[source]¶
Get seismic cumulative moment.
Additional
*args
and**kwargs
are passed toget_magnitude()
.- Returns
Cumulative seismic moment in [Nm].
- Return type
- rescale_slip(magnitude=None, moment=None, **kwargs)[source]¶
Rescale source slip based on given target magnitude or seismic moment.
Rescale the maximum source slip to fit the source moment magnitude or seismic moment to the given target values. Either
magnitude
ormoment
need to be given. Additional**kwargs
are passed toget_moment()
.
- get_centroid(store, *args, **kwargs)[source]¶
Centroid of the pseudo dynamic rupture model.
The centroid location and time are derived from the locations and times of the individual patches weighted with their moment contribution. Additional
**kwargs
are passed topyrocko_moment_tensor()
.- Parameters
store (
Store
) – Green’s function database (needs to cover whole region of of the source). Itsdeltat
[s] is used as time increment for slip difference calculation. Eitherdeltat
orstore
should be given.- Returns
The centroid location and associated moment tensor.
- Return type
pyrocko.model.Event
- class DoubleDCSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.SourceWithMagnitude
Two double-couple point sources separated in space and time. Moment share between the sub-sources is controlled by the parameter mix. The position of the subsources is dependent on the moment distribution between the two sources. Depth, east and north shift are given for the centroid between the two double-couples. The subsources will positioned according to their moment shares around this centroid position. This is done according to their delta parameters, which are therefore in relation to that centroid. Note that depth of the subsources therefore can be depth+/-delta_depth. For shallow earthquakes therefore the depth has to be chosen deeper to avoid sampling above surface.
- ♦ strike1¶
float
, default:0.0
strike direction in [deg], measured clockwise from north
- ♦ dip1¶
float
, default:90.0
dip angle in [deg], measured downward from horizontal
- ♦ azimuth¶
float
, default:0.0
azimuth to second double-couple [deg], measured at first, clockwise from north
- ♦ rake1¶
float
, default:0.0
rake angle in [deg], measured counter-clockwise from right-horizontal in on-plane view
- ♦ strike2¶
float
, default:0.0
strike direction in [deg], measured clockwise from north
- ♦ dip2¶
float
, default:90.0
dip angle in [deg], measured downward from horizontal
- ♦ rake2¶
float
, default:0.0
rake angle in [deg], measured counter-clockwise from right-horizontal in on-plane view
- ♦ delta_time¶
float
, default:0.0
separation of double-couples in time (t2-t1) [s]
- ♦ delta_depth¶
float
, default:0.0
difference in depth (z2-z1) [m]
- ♦ distance¶
float
, default:0.0
distance between the two double-couples [m]
- ♦ mix¶
float
, default:0.5
how to distribute the moment to the two doublecouples mix=0 -> m1=1 and m2=0; mix=1 -> m1=0, m2=1
- ♦ stf1¶
STF
, optionalSource time function of subsource 1 (if given, overrides STF from attribute
Source.stf
)
- ♦ stf2¶
STF
, optionalSource time function of subsource 2 (if given, overrides STF from attribute
Source.stf
)
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedMTSource
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- effective_stf_post()[source]¶
Return the STF applied after stacking of the Green’s fuctions.
This STF is used in the post-processing of the synthetic seismograms.
Handling of the STF after stacking of the GFs is usually more efficient but is only possible when a common STF is used for all subsources.
- class RingfaultSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.SourceWithMagnitude
A ring fault with vertical doublecouples.
- ♦ diameter¶
float
, default:1.0
diameter of the ring in [m]
- ♦ sign¶
float
, default:1.0
inside of the ring moves up (+1) or down (-1)
- ♦ strike¶
float
, default:0.0
strike direction of the ring plane, clockwise from north, in [deg]
- ♦ dip¶
float
, default:0.0
dip angle of the ring plane from horizontal in [deg]
- ♦ npointsources¶
int
, default:360
number of point sources to use
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedMTSource
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- class CombiSource(subsources=[], **kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.Source
Composite source model.
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedMTSource
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- class SFSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.Source
A single force point source.
Supported GF schemes: ‘elastic5’.
- ♦ fn¶
float
, default:0.0
northward component of single force [N]
- ♦ fe¶
float
, default:0.0
eastward component of single force [N]
- ♦ fd¶
float
, default:0.0
downward component of single force [N]
- discretized_source_class¶
alias of
pyrocko.gf.meta.DiscretizedSFSource
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- class PorePressurePointSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.Source
Excess pore pressure point source.
For poro-elastic initial value problem where an excess pore pressure is brought into a small source volume.
- ♦ pp¶
float
, default:1.0
initial excess pore pressure in [Pa]
- discretized_source_class¶
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- class PorePressureLineSource(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.Source
Excess pore pressure line source.
The line source is centered at (north_shift, east_shift, depth).
- ♦ pp¶
float
, default:1.0
initial excess pore pressure in [Pa]
- ♦ length¶
float
, default:0.0
length of the line source [m]
- ♦ azimuth¶
float
, default:0.0
azimuth direction, clockwise from north [deg]
- ♦ dip¶
float
, default:90.0
dip direction, downward from horizontal [deg]
- discretized_source_class¶
- base_key()[source]¶
Get key to decide about source discretization / GF stack sharing.
When two source models differ only in amplitude and origin time, the discretization and the GF stacking can be done only once for a unit amplitude and a zero origin time and the amplitude and origin times of the seismograms can be applied during post-processing of the synthetic seismogram.
For any derived parameterized source model, this method is called to decide if discretization and stacking of the source should be shared. When two source models return an equal vector of values discretization is shared.
- get_factor()[source]¶
Get the scaling factor to be applied during post-processing.
Discretization of the base seismogram is usually done for a unit amplitude, because a common factor can be efficiently multiplied to final seismograms. This eliminates to do repeat the stacking when creating seismograms for a series of source models only differing in amplitude.
This method should return the scaling factor to apply in the post-processing (often this is simply the scalar moment of the source).
- class STF(effective_duration=None, **kwargs)[source]¶
Bases:
pyrocko.guts.Object
,pyrocko.gf.seismosizer.Cloneable
Base class for source time functions.
- class BoxcarSTF(effective_duration=None, **kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.STF
Boxcar type source time function.
- ♦ duration¶
float
, default:0.0
duration of the boxcar
- ♦ anchor¶
float
, default:0.0
anchor point with respect to source.time: (-1.0: left -> source duration [0, T] ~ hypocenter time, 0.0: center -> source duration [-T/2, T/2] ~ centroid time, +1.0: right -> source duration [-T, 0] ~ rupture end time)
- class TriangularSTF(effective_duration=None, **kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.STF
Triangular type source time function.
- ♦ duration¶
float
, default:0.0
baseline of the triangle
- ♦ peak_ratio¶
float
, default:0.5
fraction of time compared to duration, when the maximum amplitude is reached
- ♦ anchor¶
float
, default:0.0
anchor point with respect to source.time: (-1.0: left -> source duration [0, T] ~ hypocenter time, 0.0: center -> source duration [-T/2, T/2] ~ centroid time, +1.0: right -> source duration [-T, 0] ~ rupture end time)
- class HalfSinusoidSTF(effective_duration=None, **kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.STF
Half sinusoid type source time function.
- ♦ duration¶
float
, default:0.0
duration of the half-sinusoid (baseline)
- ♦ anchor¶
float
, default:0.0
anchor point with respect to source.time: (-1.0: left -> source duration [0, T] ~ hypocenter time, 0.0: center -> source duration [-T/2, T/2] ~ centroid time, +1.0: right -> source duration [-T, 0] ~ rupture end time)
- ♦ exponent¶
int
, default:1
set to 2 to use square of the half-period sinusoidal function.
- class ResonatorSTF(effective_duration=None, **kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.STF
Simple resonator like source time function.
- ♦ duration¶
float
, default:0.0
decay time
- ♦ frequency¶
float
, default:1.0
resonance frequency
- class Request(*args, **kwargs)[source]¶
Bases:
pyrocko.guts.Object
Synthetic seismogram computation request.
Request(**kwargs) Request(sources, targets, **kwargs)
- ♦ targets¶
list
ofpyrocko.gf.targets.Target
objects, default:[]
list of targets for which to produce synthetics.
- class ProcessingStats(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Undocumented.
- ♦ t_perc_get_store_and_receiver¶
float
, default:0.0
- ♦ t_perc_discretize_source¶
float
, default:0.0
- ♦ t_perc_make_base_seismogram¶
float
, default:0.0
- ♦ t_perc_make_same_span¶
float
, default:0.0
- ♦ t_perc_post_process¶
float
, default:0.0
- ♦ t_perc_optimize¶
float
, default:0.0
- ♦ t_perc_stack¶
float
, default:0.0
- ♦ t_perc_static_get_store¶
float
, default:0.0
- ♦ t_perc_static_discretize_basesource¶
float
, default:0.0
- ♦ t_perc_static_sum_statics¶
float
, default:0.0
- ♦ t_perc_static_post_process¶
float
, default:0.0
- ♦ t_wallclock¶
float
, default:0.0
- ♦ t_cpu¶
float
, default:0.0
- ♦ n_read_blocks¶
int
, default:0
- ♦ n_results¶
int
, default:0
- ♦ n_subrequests¶
int
, default:0
- ♦ n_stores¶
int
, default:0
- ♦ n_records_stacked¶
int
, default:0
- class Response(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Resonse object to a synthetic seismogram computation request.
- ♦ results_list¶
list
oflist
ofpyrocko.gf.meta.SeismosizerResult
objects objects, default:[]
- ♦ stats¶
- static_results()[source]¶
Return a list of requested
StaticResult
instances.
- class Engine(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Base class for synthetic seismogram calculators.
- class LocalEngine(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.Engine
Offline synthetic seismogram calculator.
- Parameters
use_env – if
True
, fillstore_superdirs
andstore_dirs
with paths set in environment variables GF_STORE_SUPERDIRS and GF_STORE_DIRS.use_config –
if
True
, fillstore_superdirs
andstore_dirs
with paths set in the user’s config file.The config file can be found at
~/.pyrocko/config.pf
gf_store_dirs: ['/home/pyrocko/gf_stores/ak135/'] gf_store_superdirs: ['/home/pyrocko/gf_stores/']
- ♦ store_superdirs¶
list
ofstr
objects, default:[]
directories which are searched for Green’s function stores
- ♦ store_dirs¶
list
ofstr
objects, default:[]
additional individual Green’s function store directories
- ♦ default_store_id¶
str
, optionaldefault store ID to be used when a request does not provide one
- get_store(store_id=None)[source]¶
Get a store from the engine.
- Parameters
store_id – identifier of the store (optional)
- Returns
Store
object
If no
store_id
is provided the store associated with thedefault_store_id
is returned. RaisesNoDefaultStoreSet
ifdefault_store_id
is undefined.
- class RemoteEngine(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.Engine
Client for remote synthetic seismogram calculator.
- ♦ site¶
str
, optional, default:'localhost'
- ♦ url¶
str
, optional, default:'%(site)s/gfws/%(service)s/%(majorversion)i/%(method)s'
- class Range(*args, **kwargs)[source]¶
Bases:
pyrocko.guts.SObject
Convenient range specification.
Equivalent ways to sepecify the range [ 0., 1000., … 10000. ]:
Range('0 .. 10k : 1k') Range(start=0., stop=10e3, step=1e3) Range(0, 10e3, 1e3) Range('0 .. 10k @ 11') Range(start=0., stop=10*km, n=11) Range(0, 10e3, n=11) Range(values=[x*1e3 for x in range(11)])
Depending on the use context, it can be possible to omit any part of the specification. E.g. in the context of extracting a subset of an already existing range, the existing range’s specification values would be filled in where missing.
The values are distributed with equal spacing, unless the
spacing
argument is modified. The values can be created offset or relative to an external base value with therelative
argument if the use context supports this.The range specification can be expressed with a short string representation:
'start .. stop @ num | spacing, relative' 'start .. stop : step | spacing, relative'
most parts of the expression can be omitted if not needed. Whitespace is allowed for readability but can also be omitted.
- ♦ start¶
float
, optional
- ♦ stop¶
float
, optional
- ♦ step¶
float
, optional
- ♦ n¶
int
, optional
- ♦ values¶
numpy.ndarray
(pyrocko.guts_array.Array
), optional
- ♦ spacing¶
str
(pyrocko.guts.StringChoice
), optional, default:'lin'
- ♦ relative¶
str
(pyrocko.guts.StringChoice
), optional, default:''
- class SourceList(**kwargs)[source]¶
Bases:
pyrocko.gf.seismosizer.SourceGroup
Undocumented.
targets
¶
- class OptimizationMethod(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringChoice
Any
str
out of['enable', 'disable']
.
- component_orientation(source, target, component)[source]¶
Get component and azimuth for standard components R, T, Z, N, and E.
- Parameters
source –
pyrocko.gf.Location
objecttarget –
pyrocko.gf.Location
objectcomponent – string
'R'
,'T'
,'Z'
,'N'
or'E'
- class Target(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.Receiver
A seismogram computation request for a single component, including its post-processing parmeters.
- ♦ quantity¶
str
(pyrocko.gf.meta.QuantityType
), optionalMeasurement quantity type. If not given, it is guessed from the channel code. For some common cases, derivatives of the stored quantities are supported by using finite difference approximations (e.g. displacement to velocity or acceleration). 4th order central FD schemes are used.
- ♦ codes¶
tuple
of 4str
objects, default:('', 'STA', '', 'Z')
network, station, location and channel codes to be set on the response trace.
- ♦ elevation¶
float
, default:0.0
station surface elevation in [m]
- ♦ store_id¶
str
(pyrocko.gf.meta.StringID
), optionalID of Green’s function store to use for the computation. If not given, the processor may use a system default.
- ♦ sample_rate¶
float
, optionalsample rate to produce. If not given the GF store’s default sample rate is used. GF store specific restrictions may apply.
- ♦ interpolation¶
str
(pyrocko.gf.meta.InterpolationMethod
), default:'nearest_neighbor'
Interpolation method between Green’s functions. Supported are
nearest_neighbor
andmultilinear
- ♦ optimization¶
str
(OptimizationMethod
), optional, default:'enable'
disable/enable optimizations in weight-delay-and-sum operation
- ♦ tmin¶
time_float (
pyrocko.guts.Timestamp
), optionaltime of first sample to request in [s]. If not given, it is determined from the Green’s functions.
- ♦ tmax¶
time_float (
pyrocko.guts.Timestamp
), optionaltime of last sample to request in [s]. If not given, it is determined from the Green’s functions.
- ♦ azimuth¶
float
, optionalazimuth of sensor component in [deg], clockwise from north. If not given, it is guessed from the channel code.
- ♦ dip¶
float
, optionaldip of sensor component in [deg], measured downward from horizontal. If not given, it is guessed from the channel code.
- class StaticTarget(*args, **kwargs)[source]¶
Bases:
pyrocko.gf.meta.MultiLocation
A computation request for a spatial multi-location target of static/geodetic quantities.
- ♦ quantity¶
str
(pyrocko.gf.meta.QuantityType
), optional, default:'displacement'
Measurement quantity type, for now only displacement issupported.
- ♦ interpolation¶
str
(pyrocko.gf.meta.InterpolationMethod
), default:'nearest_neighbor'
Interpolation method between Green’s functions. Supported are
nearest_neighbor
andmultilinear
- ♦ tsnapshot¶
time_float (
pyrocko.guts.Timestamp
), optionaltime of the desired snapshot in [s], If not given, the first sample is taken. If the desired sample exceeds the length of the Green’s function store, the last sample is taken.
- ♦ store_id¶
str
(pyrocko.gf.meta.StringID
), optionalID of Green’s function store to use for the computation. If not given, the processor may use a system default.
- property ntargets¶
Number of targets held by instance.
- class SatelliteTarget(*args, **kwargs)[source]¶
Bases:
pyrocko.gf.targets.StaticTarget
A computation request for a spatial multi-location target of static/geodetic quantities measured from a satellite instrument. The line of sight angles are provided and projecting post-processing is applied.
- ♦ theta¶
numpy.ndarray
(pyrocko.guts_array.Array
)Horizontal angle towards satellite’s line of sight in radians.
Important
is east and is north.
- ♦ phi¶
numpy.ndarray
(pyrocko.guts_array.Array
)Theta is look vector elevation angle towards satellite from horizon in radians. Matrix of theta towards satellite’s line of sight.
Important
is down and is up.
- class KiteSceneTarget(scene, **kwargs)[source]¶
Bases:
pyrocko.gf.targets.SatelliteTarget
Undocumented.
- ♦ shape¶
tuple
of 2int
objects, default:(None, None)
Shape of the displacement vectors.
- class GNSSCampaignTarget(*args, **kwargs)[source]¶
Bases:
pyrocko.gf.targets.StaticTarget
Undocumented.
tractions
¶
- class AbstractTractionField(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Base class for multiplicative traction fields (tapers).
Fields of this type a re multiplied in the
TractionComposition
- class TractionField(**kwargs)[source]¶
Bases:
pyrocko.gf.tractions.AbstractTractionField
Base class for additive traction fields.
Fields of this type are added in the
TractionComposition
- class TractionComposition(**kwargs)[source]¶
Bases:
pyrocko.gf.tractions.TractionField
Composition of traction fields.
TractionField
andAbstractTractionField
can be combined to realize a combination of different fields.- ♦ components¶
list
ofAbstractTractionField
objects, default:[]
Ordered list of tractions.
- class HomogeneousTractions(**kwargs)[source]¶
Bases:
pyrocko.gf.tractions.TractionField
Homogeneous traction field.
The traction vectors in strike, dip and normal direction are acting homogeneously on the rupture plane.
- ♦ strike¶
float
, default:1.0
Tractions in strike direction [Pa].
- ♦ dip¶
float
, default:1.0
Traction in dip direction (up) [Pa].
- ♦ normal¶
float
, default:1.0
Traction in normal direction [Pa].
- class DirectedTractions(**kwargs)[source]¶
Bases:
pyrocko.gf.tractions.TractionField
Directed traction field.
The traction vectors are following a uniform
rake
.- ♦ rake¶
float
, default:0.0
Rake angle in [deg], measured counter-clockwise from right-horizontal in on-plane view. Rake is translated into homogenous tractions in strike and up-dip direction.
- ♦ traction¶
float
, default:1.0
Traction in rake direction [Pa].
- class FractalTractions(*args, **kwargs)[source]¶
Bases:
pyrocko.gf.tractions.TractionField
Fractal traction field.
- ♦ rseed¶
int
, optionalSeed for
RandomState
.IfNone
, an random seed will be initialized.
- ♦ rake¶
float
, default:0.0
Rake angle in [deg], measured counter-clockwise from right-horizontal in on-plane view. Rake is translated into homogenous tractions in strike and up-dip direction.
- ♦ traction_max¶
float
, default:1.0
Maximum traction vector length [Pa].
- class SelfSimilarTractions(**kwargs)[source]¶
Bases:
pyrocko.gf.tractions.TractionField
Traction model following Power & Tullis (1991).
The traction vectors are calculated as a sum of 2D-cosines with a constant amplitude / wavelength ratio. The wavenumber kx and ky are constant for each cosine function. The rank defines the maximum wavenumber used for summation. So, e.g. a rank of 3 will lead to a summation of cosines with
kx = ky
in (1, 2, 3). Each cosine has an associated phases, which defines both the phase shift and also the shift from the rupture plane centre. Finally the summed cosines are translated into shear tractions based on the rake and normalized withtraction_max
.- ♦ rank¶
int
, default:1
Maximum summed cosine wavenumber/spatial frequency.
- ♦ rake¶
float
, default:0.0
Rake angle in [deg], measured counter-clockwise from right-horizontal in on-plane view. Rake is translated into homogenous tractions in strike and up-dip direction.
- ♦ traction_max¶
float
, default:1.0
Maximum traction vector length [Pa].
- ♦ phases¶
numpy.ndarray
(pyrocko.guts_array.Array
), optionalPhase shift of the cosines in [rad].
- class RectangularTaper(**kwargs)[source]¶
Bases:
pyrocko.gf.tractions.AbstractTractionField
Undocumented.
- ♦ width¶
float
, default:0.2
Width of the taper as a fraction of the plane.
- ♦ type¶
str
(pyrocko.guts.StringChoice
), default:'tukey'
Type of the taper, default: “tukey”.
- class DepthTaper(**kwargs)[source]¶
Bases:
pyrocko.gf.tractions.AbstractTractionField
Undocumented.
- ♦ depth_start¶
float
Depth where the taper begins [m].
- ♦ depth_stop¶
float
Depth where taper ends and drops to zero [m].
- ♦ type¶
str
(pyrocko.guts.StringChoice
), default:'linear'
Type of the taper, default: “linear”.
- plot_tractions(tractions, nx=15, ny=12, depth=10000.0, component='strike')[source]¶
Plot traction model for quick inspection.
- Parameters
tractions (
pyrocko.gf.tractions.TractionField
) – Traction field or traction composition to be displayed.nx (optional, int) – Number of patches along strike.
ny (optional, int) – Number of patches down dip.
depth (optional, float) – Depth of the rupture plane center in [m].
component (optional, str) – Choice of traction component to be shown. Available:
'tx'
(along strike),'ty'
(up dip),'tz'
(normal),'absolute'
(vector length).
meta
¶
- class StringID(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringPattern
Any
str
matching pattern'^[A-Za-z][A-Za-z0-9._]{0,64}$'
.
- class ScopeType(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringChoice
Any
str
out of['global', 'regional', 'local']
.
- class WaveformType(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringChoice
Any
str
out of['dis', 'vel', 'acc', 'amp_spec_dis', 'amp_spec_vel', 'amp_spec_acc', 'envelope_dis', 'envelope_vel', 'envelope_acc']
.
- class QuantityType(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringChoice
Any
str
out of['displacement', 'rotation', 'velocity', 'acceleration', 'pressure', 'tilt', 'pore_pressure', 'darcy_velocity', 'vertical_tilt']
.
- class NearfieldTermsType(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringChoice
Any
str
out of['complete', 'incomplete', 'missing']
.
- class Reference(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Undocumented.
- ♦ type¶
str
- ♦ title¶
str
- ♦ journal¶
str
, optional
- ♦ volume¶
str
, optional
- ♦ number¶
str
, optional
- ♦ pages¶
str
, optional
- ♦ year¶
str
- ♦ issn¶
str
, optional
- ♦ doi¶
str
, optional
- ♦ url¶
str
, optional
- ♦ eprint¶
str
, optional
- ♦ authors¶
list
ofstr
objects, default:[]
- ♦ publisher¶
str
, optional
- ♦ keywords¶
str
, optional
- ♦ note¶
str
, optional
- ♦ abstract¶
str
, optional
- class CircularRegion(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.Region
Undocumented.
- ♦ lat¶
float
- ♦ lon¶
float
- ♦ radius¶
float
- class RectangularRegion(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.Region
Undocumented.
- ♦ lat_min¶
float
- ♦ lat_max¶
float
- ♦ lon_min¶
float
- ♦ lon_max¶
float
- class PhaseSelect(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringChoice
Any
str
out of['', 'first', 'last']
.
- class Timing(s=None, **kwargs)[source]¶
Bases:
pyrocko.guts.SObject
Definition of a time instant relative to one or more named phase arrivals.
Instances of this class can be used e.g. in cutting and tapering operations. They can hold an absolute time or an offset to a named phase arrival or group of such arrivals.
Timings can be instantiated from a simple string defintion i.e. with
Timing(str)
wherestr
is something like'SELECT{PHASE_DEFS}[+-]OFFSET[S|%]'
where'SELECT'
is'first'
,'last'
or empty,'PHASE_DEFS'
is a'|'
-separated list of phase definitions, and'OFFSET'
is the time offset in seconds. If a'%'
is appended, it is interpreted as percent. If the an'S'
is appended to'OFFSET'
, it is interpreted as a surface slowness in [s/km].Phase definitions can be specified in either of the following ways:
'stored:PHASE_ID'
- retrieves value from stored travel time table'cake:CAKE_PHASE_DEF'
- evaluates first arrival of phase with cake (seepyrocko.cake.PhaseDef
)'vel_surface:VELOCITY'
- arrival according to surface distance / velocity [km/s]'vel:VELOCITY'
- arrival according to 3D-distance / velocity [km/s]
Examples:
'100'
: absolute time; 100 s'{stored:P}-100'
: 100 s before arrival of P phase according to stored travel time table named'P'
'{stored:P}-5.1S'
: 10% before arrival of P phase according to stored travel time table named'P'
'{stored:P}-10%'
: 10% before arrival of P phase according to stored travel time table named'P'
'{stored:A|stored:B}'
: time instant of phase arrival A, or B if A is undefined for a given geometry'first{stored:A|stored:B}'
: as above, but the earlier arrival of A and B is chosen, if both phases are defined for a given geometry'last{stored:A|stored:B}'
: as above but the later arrival is chosen'first{stored:A|stored:B|stored:C}-100'
: 100 s before first out of arrivals A, B, and C
- ♦ phase_defs¶
list
ofstr
objects, default:[]
- ♦ offset¶
float
, default:0.0
- ♦ offset_is¶
str
, optional
- ♦ select¶
str
(PhaseSelect
), default:''
Can be either
''
,'first'
, or'last'
.
- class TPDef(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Maps an arrival phase identifier to an arrival phase definition.
- ♦ definition¶
str
definition of the phase in either cake syntax as defined in
pyrocko.cake.PhaseDef
, or, if prepended with an!
, as a classic phase name, or, if it is a simple number, as a constant horizontal velocity.
- class Location(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Geographical location.
The location is given by a reference point at the earth’s surface (
lat
,lon
,elevation
) and a cartesian offset from this point (north_shift
,east_shift
,depth
). The offset corrected lat/lon coordinates of the location can be accessed though theeffective_latlon
,effective_lat
, andeffective_lon
properties.- ♦ lat¶
float
, optional, default:0.0
latitude of reference point [deg]
- ♦ lon¶
float
, optional, default:0.0
longitude of reference point [deg]
- ♦ north_shift¶
float
, optional, default:0.0
northward cartesian offset from reference point [m]
- ♦ east_shift¶
float
, optional, default:0.0
eastward cartesian offset from reference point [m]
- ♦ elevation¶
float
, optional, default:0.0
surface elevation, above sea level [m]
- ♦ depth¶
float
, default:0.0
depth, below surface [m]
- property effective_latlon¶
Property holding the offset-corrected lat/lon pair of the location.
- property effective_lat¶
Property holding the offset-corrected latitude of the location.
- property effective_lon¶
Property holding the offset-corrected longitude of the location.
- class Receiver(**kwargs)[source]¶
Bases:
pyrocko.model.location.Location
Undocumented.
- ♦ codes¶
tuple
of 3str
objects, optionalnetwork, station, and location codes
- class DiscretizedExplosionSource(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.DiscretizedSource
Undocumented.
- ♦ m0s¶
numpy.ndarray
(pyrocko.guts_array.Array
)
- classmethod combine(sources, **kwargs)[source]¶
Combine several discretized source models.
Concatenenates all point sources in the given discretized
sources
. Care must be taken when using this function that the external amplitude factors and reference times of the parameterized (undiscretized) sources match or are accounted for.
- class DiscretizedSFSource(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.DiscretizedSource
Undocumented.
- ♦ forces¶
numpy.ndarray
(pyrocko.guts_array.Array
)
- classmethod combine(sources, **kwargs)[source]¶
Combine several discretized source models.
Concatenenates all point sources in the given discretized
sources
. Care must be taken when using this function that the external amplitude factors and reference times of the parameterized (undiscretized) sources match or are accounted for.
- class DiscretizedMTSource(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.DiscretizedSource
Undocumented.
- ♦ m6s¶
numpy.ndarray
(pyrocko.guts_array.Array
)rows with (m_nn, m_ee, m_dd, m_ne, m_nd, m_ed)
- classmethod combine(sources, **kwargs)[source]¶
Combine several discretized source models.
Concatenenates all point sources in the given discretized
sources
. Care must be taken when using this function that the external amplitude factors and reference times of the parameterized (undiscretized) sources match or are accounted for.
- class DiscretizedPorePressureSource(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.DiscretizedSource
Undocumented.
- ♦ pp¶
numpy.ndarray
(pyrocko.guts_array.Array
)
- classmethod combine(sources, **kwargs)[source]¶
Combine several discretized source models.
Concatenenates all point sources in the given discretized
sources
. Care must be taken when using this function that the external amplitude factors and reference times of the parameterized (undiscretized) sources match or are accounted for.
- class ConfigTypeA(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.Config
Cylindrical symmetry, 1D earth model, single receiver depth
Problem is invariant to horizontal translations and rotations around vertical axis.
All receivers must be at the same depth (e.g. at the surface) High level index variables:
(source_depth, distance, component)
The
distance
is the surface distance between source and receiver points.
- ♦ receiver_depth¶
float
, default:0.0
Fixed receiver depth [m].
- ♦ source_depth_min¶
float
Minimum source depth [m].
- ♦ source_depth_max¶
float
Maximum source depth [m].
- ♦ source_depth_delta¶
float
Grid spacing of source depths [m]
- ♦ distance_min¶
float
Minimum source-receiver surface distance [m].
- ♦ distance_max¶
float
Maximum source-receiver surface distance [m].
- ♦ distance_delta¶
float
Grid spacing of source-receiver surface distance [m].
- class ConfigTypeB(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.Config
Cylindrical symmetry, 1D earth model, variable receiver depth
Symmetries like in
ConfigTypeA
but has additional index for receiver depthHigh level index variables:
(receiver_depth, source_depth, receiver_distance, component)
- ♦ receiver_depth_min¶
float
Minimum receiver depth [m].
- ♦ receiver_depth_max¶
float
Maximum receiver depth [m].
- ♦ receiver_depth_delta¶
float
Grid spacing of receiver depths [m]
- ♦ source_depth_min¶
float
Minimum source depth [m].
- ♦ source_depth_max¶
float
Maximum source depth [m].
- ♦ source_depth_delta¶
float
Grid spacing of source depths [m]
- ♦ distance_min¶
float
Minimum source-receiver surface distance [m].
- ♦ distance_max¶
float
Maximum source-receiver surface distance [m].
- ♦ distance_delta¶
float
Grid spacing of source-receiver surface distances [m].
- class ConfigTypeC(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.Config
No symmetrical constraints, one fixed receiver position.
Cartesian 3D source volume around a reference point
High level index variables:
(source_depth, source_east_shift, source_north_shift, component)
- ♦ source_origin¶
pyrocko.model.location.Location
Origin of the source volume grid.
- ♦ source_depth_min¶
float
Minimum source depth [m].
- ♦ source_depth_max¶
float
Maximum source depth [m].
- ♦ source_depth_delta¶
float
Source depth grid spacing [m].
- ♦ source_east_shift_min¶
float
Minimum easting of source grid [m].
- ♦ source_east_shift_max¶
float
Maximum easting of source grid [m].
- ♦ source_east_shift_delta¶
float
Source volume grid spacing in east direction [m].
- ♦ source_north_shift_min¶
float
Minimum northing of source grid [m].
- ♦ source_north_shift_max¶
float
Maximum northing of source grid [m].
- ♦ source_north_shift_delta¶
float
Source volume grid spacing in north direction [m].
- class ComponentScheme(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringChoice
Different Green’s Function component schemes are available:
Name
Description
elastic10
Elastodynamic for
ConfigTypeA
andConfigTypeB
stores, MT sources onlyelastic8
Elastodynamic for far-field only
ConfigTypeA
andConfigTypeB
stores, MT sources onlyelastic2
Elastodynamic for
ConfigTypeA
andConfigTypeB
stores, purely isotropic sources onlyelastic5
Elastodynamic for
ConfigTypeA
andConfigTypeB
stores, SF sources onlyelastic18
Elastodynamic for
ConfigTypeC
stores, MT sources onlyporoelastic10
Poroelastic for
ConfigTypeA
andConfigTypeB
stores
- class Config(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Green’s function store meta information.
Currently implemented
Store
configuration types are:ConfigTypeA
- cylindrical or spherical symmetry, 1D earth model, single receiver depthProblem is invariant to horizontal translations and rotations around vertical axis.
All receivers must be at the same depth (e.g. at the surface)
High level index variables:
(source_depth, receiver_distance, component)
ConfigTypeB
- cylindrical or spherical symmetry, 1D earth model, variable receiver depthSymmetries like in Type A but has additional index for receiver depth
High level index variables:
(source_depth, receiver_distance, receiver_depth, component)
ConfigTypeC
- no symmetrical constraints but fixed receiver positionsCartesian source volume around a reference point
High level index variables:
(ireceiver, source_depth, source_east_shift, source_north_shift, component)
- ♦ id¶
str
(StringID
)Name of the store. May consist of upper and lower-case letters, digits, dots and underscores. The name must start with a letter.
- ♦ derived_from_id¶
str
(StringID
), optionalName of the original store, if this store has been derived from another one (e.g. extracted subset).
- ♦ version¶
str
, optional, default:'1.0'
User-defined version string. Use <major>.<minor> format.
- ♦ author¶
str
, optionalComma-separated list of author names.
- ♦ author_email¶
str
, optionalAuthor’s contact email address.
- ♦ created_time¶
time_float (
pyrocko.guts.Timestamp
), optionalTime of creation of the store.
- ♦ regions¶
list
ofRegion
objects, default:[]
Geographical regions for which the store is representative.
- ♦ scope_type¶
str
(ScopeType
), optionalDistance range scope of the store (choices:
'global'
,'regional'
,'local'
).
- ♦ waveform_type¶
str
(WaveType
), optionalWave type stored (choices:
'full waveform'
,'bodywave'
,'P wave'
,'S wave'
,'surface wave'
).
- ♦ nearfield_terms¶
str
(NearfieldTermsType
), optionalInformation about the inclusion of near-field terms in the modelling (choices:
'complete'
,'incomplete'
,'missing'
).
- ♦ description¶
str
, optionalFree form textual description of the GF store.
- ♦ references¶
list
ofReference
objects, default:[]
Reference list to cite the modelling code, earth model or related work.
- ♦ earthmodel_1d¶
pyrocko.cake.LayeredModel
(Earthmodel1D
), optionalLayered earth model in ND (named discontinuity) format.
- ♦ earthmodel_receiver_1d¶
pyrocko.cake.LayeredModel
(Earthmodel1D
), optionalReceiver-side layered earth model in ND format.
- ♦ can_interpolate_source¶
bool
, optionalHint to indicate if the spatial sampling of the store is dense enough for multi-linear interpolation at the source.
- ♦ can_interpolate_receiver¶
bool
, optionalHint to indicate if the spatial sampling of the store is dense enough for multi-linear interpolation at the receiver.
- ♦ frequency_min¶
float
, optionalHint to indicate the lower bound of valid frequencies [Hz].
- ♦ frequency_max¶
float
, optionalHint to indicate the upper bound of valid frequencies [Hz].
- ♦ sample_rate¶
float
, optionalSample rate of the GF store [Hz].
- ♦ factor¶
float
, optional, default:1.0
Gain value, factored out of the stored GF samples. (may not work properly, keep at 1.0).
- ♦ component_scheme¶
str
(ComponentScheme
), default:'elastic10'
GF component scheme (choices:
'elastic2'
,'elastic8'
,'elastic10'
,'elastic18'
,'elastic5'
,'poroelastic10'
).
- ♦ stored_quantity¶
str
(QuantityType
), optionalPhysical quantity of stored values (choices:
'displacement'
,'rotation'
,'velocity'
,'acceleration'
,'pressure'
,'tilt'
,'pore_pressure'
,'darcy_velocity'
,'vertical_tilt'
). If not given, a default is used based on the GF component scheme. The default for the"elastic*"
family of component schemes is"displacement"
.
- ♦ tabulated_phases¶
list
ofTPDef
objects, default:[]
Mapping of phase names to phase definitions, for which travel time tables are available in the GF store.
- ♦ ncomponents¶
int
, optionalNumber of GF components. Use
component_scheme
instead.
- ♦ uuid¶
str
, optionalHeuristic hash value which can be used to uniquely identify the GF store for practical purposes.
- ♦ reference¶
str
, optionalStore reference name composed of the store’s
id
and the first six letters of itsuuid
.
- get_shear_moduli(lat, lon, points, interpolation=None)[source]¶
Get shear moduli at given points from contained velocity model.
- Parameters
lat – surface origin for coordinate system of
points
points – NumPy array of shape
(N, 3)
, where each row is a point(north, east, depth)
, relative to origin at(lat, lon)
interpolation – interpolation method. Choose from
('nearest_neighbor', 'multilinear')
- Returns
NumPy array of length N with extracted shear moduli at each point
The default implementation retrieves and interpolates the shear moduli from the contained 1D velocity profile.
- get_lambda_moduli(lat, lon, points, interpolation=None)[source]¶
Get lambda moduli at given points from contained velocity model.
- Parameters
lat – surface origin for coordinate system of
points
points – NumPy array of shape
(N, 3)
, where each row is a point(north, east, depth)
, relative to origin at(lat, lon)
interpolation – interpolation method. Choose from
('nearest_neighbor', 'multilinear')
- Returns
NumPy array of length N with extracted shear moduli at each point
The default implementation retrieves and interpolates the lambda moduli from the contained 1D velocity profile.
- get_bulk_moduli(lat, lon, points, interpolation=None)[source]¶
Get bulk moduli at given points from contained velocity model.
- Parameters
lat – surface origin for coordinate system of
points
points – NumPy array of shape
(N, 3)
, where each row is a point(north, east, depth)
, relative to origin at(lat, lon)
interpolation – interpolation method. Choose from
('nearest_neighbor', 'multilinear')
- Returns
NumPy array of length N with extracted shear moduli at each point
The default implementation retrieves and interpolates the lambda moduli from the contained 1D velocity profile.
- get_vs(lat, lon, points, interpolation=None)[source]¶
Get Vs at given points from contained velocity model.
- Parameters
lat – surface origin for coordinate system of
points
points – NumPy array of shape
(N, 3)
, where each row is a point(north, east, depth)
, relative to origin at(lat, lon)
interpolation – interpolation method. Choose from
('nearest_neighbor', 'multilinear')
- Returns
NumPy array of length N with extracted shear moduli at each point
The default implementation retrieves and interpolates Vs from the contained 1D velocity profile.
- get_vp(lat, lon, points, interpolation=None)[source]¶
Get Vp at given points from contained velocity model.
- Parameters
lat – surface origin for coordinate system of
points
points – NumPy array of shape
(N, 3)
, where each row is a point(north, east, depth)
, relative to origin at(lat, lon)
interpolation – interpolation method. Choose from
('nearest_neighbor', 'multilinear')
- Returns
NumPy array of length N with extracted shear moduli at each point
The default implementation retrieves and interpolates Vp from the contained 1D velocity profile.
- get_rho(lat, lon, points, interpolation=None)[source]¶
Get rho at given points from contained velocity model.
- Parameters
lat – surface origin for coordinate system of
points
points – NumPy array of shape
(N, 3)
, where each row is a point(north, east, depth)
, relative to origin at(lat, lon)
interpolation – interpolation method. Choose from
('nearest_neighbor', 'multilinear')
- Returns
NumPy array of length N with extracted shear moduli at each point
The default implementation retrieves and interpolates rho from the contained 1D velocity profile.
- class Weighting(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Undocumented.
- ♦ factor¶
float
, default:1.0
- class Taper(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Undocumented.
- ♦ tfade¶
float
, default:0.0
- ♦ shape¶
str
(pyrocko.guts.StringChoice
), optional, default:'cos'
- class ChannelSelection(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Undocumented.
- ♦ pattern¶
SimplePattern
, optional
- ♦ min_sample_rate¶
float
, optional
- ♦ max_sample_rate¶
float
, optional
- class StationSelection(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Undocumented.
- ♦ includes¶
- ♦ excludes¶
- ♦ distance_min¶
float
, optional
- ♦ distance_max¶
float
, optional
- ♦ azimuth_min¶
float
, optional
- ♦ azimuth_max¶
float
, optional
- class WaveformSelection(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Undocumented.
- ♦ channel_selection¶
ChannelSelection
, optional
- ♦ station_selection¶
StationSelection
, optional
- ♦ waveform_type¶
str
(WaveformType
), default:'dis'
- ♦ sample_rate¶
float
, optional
Bases:
Exception
- class InterpolationMethod(...) dummy for str [source]¶
Bases:
pyrocko.guts.StringChoice
Any
str
out of['nearest_neighbor', 'multilinear']
.
- class SeismosizerTrace(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Undocumented.
- ♦ codes¶
tuple
of 4str
objects, default:('', 'STA', '', 'Z')
network, station, location and channel codes
- ♦ data¶
numpy.ndarray
(pyrocko.guts_array.Array
)numpy array with data samples
- ♦ deltat¶
float
, default:1.0
sampling interval [s]
- ♦ tmin¶
time_float (
pyrocko.guts.Timestamp
), default:str_to_time('1970-01-01 00:00:00')
time of first sample as a system timestamp [s]
- class SeismosizerResult(**kwargs)[source]¶
Bases:
pyrocko.guts.Object
Undocumented.
- ♦ n_records_stacked¶
int
, optional, default:1
- ♦ t_stack¶
float
, optional, default:0.0
- class Result(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.SeismosizerResult
Undocumented.
- ♦ trace¶
SeismosizerTrace
, optional
int
, optional, default:1
- ♦ t_optimize¶
float
, optional, default:0.0
- class StaticResult(**kwargs)[source]¶
Bases:
pyrocko.gf.meta.SeismosizerResult
Undocumented.
- ♦ result¶
dict
ofnumpy.ndarray
(pyrocko.guts_array.Array
) objects, default:{}