gf
¶
gf.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.
-
t
¶ Time vector of the GF trace.
Returns: Time vector Return type: numpy.ndarray
-
-
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_ttt()
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.Parameters:
-
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: GFTrace
-
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 None - nsamples (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: - args (tuple) –
-
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 times - weights (
numpy.ndarray
) – Trace weights - itmin (int or None) – Start time index (start time is
itmin * dt
), defaults to None - nsamples (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: - args (tuple(numpy.ndarray)) –
-
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.Parameters:
-
get_stored_phase
(phase_id)[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: - timing (str or
-
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_ttt
(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.
-
gf.builder
¶
gf.seismosizer
¶
gf.targets
¶
gf.meta
¶
-
class
InterpolationMethod
(dummy) → 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
¶ builtins.float
(pyrocko.guts.Timestamp
), default:0.0
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:{}
-
♦
-
class
ComponentScheme
(dummy) → 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
StringID
(dummy) → str[source]¶ Bases:
pyrocko.guts.StringPattern
Any
str
matching pattern'^[A-Za-z][A-Za-z0-9._]{0,64}$'
.
-
class
ScopeType
(dummy) → str[source]¶ Bases:
pyrocko.guts.StringChoice
Any
str
out of['global', 'regional', 'local']
.
-
class
NearfieldTermsType
(dummy) → str[source]¶ Bases:
pyrocko.guts.StringChoice
Any
str
out of['complete', 'incomplete', 'missing']
.
-
class
QuantityType
(dummy) → str[source]¶ Bases:
pyrocko.guts.StringChoice
Any
str
out of['displacement', 'velocity', 'acceleration', 'pressure', 'tilt', 'pore_pressure', 'darcy_velocity', 'vertical_tilt']
.
-
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
list
ofstr
objects, default:[]
-
♦
publisher
¶ str
, optional
-
♦
keywords
¶ str
, optional
-
♦
note
¶ str
, optional
-
♦
abstract
¶ str
, optional
-
♦
-
class
PhaseSelect
(dummy) → 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 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: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_slowness
¶ bool
, default:False
-
♦
select
¶ builtins.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
Receiver
(**kwargs)[source]¶ Bases:
pyrocko.model.location.Location
Undocumented.
-
♦
codes
¶ tuple
of 3str
objects, optionalnetwork, station, and location codes
-
♦
Bases:
Exception
-
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
RectangularRegion
(**kwargs)[source]¶ Bases:
pyrocko.gf.meta.Region
Undocumented.
-
♦
lat_min
¶ float
-
♦
lat_max
¶ float
-
♦
lon_min
¶ float
-
♦
lon_max
¶ float
-
♦
-
class
CircularRegion
(**kwargs)[source]¶ Bases:
pyrocko.gf.meta.Region
Undocumented.
-
♦
lat
¶ float
-
♦
lon
¶ float
-
♦
radius
¶ float
-
♦
-
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 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, receiver_distance, component)
ConfigTypeB
- cylindrical or spherical symmetry, 1D earth model, variable receiver depth- Symmetries 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 positions- Cartesian source volume around a reference point
- High level index variables:
(ireceiver, source_depth, source_east_shift, source_north_shift, component)
-
♦
id
¶ builtins.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
¶ builtins.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.
-
♦
modelling_code_id
¶ builtins.str
(StringID
), optionalIdentifier of the backend used to compute the store.
str
, optionalComma-separated list of author names.
str
, optionalAuthor’s contact email address.
-
♦
created_time
¶ builtins.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
¶ builtins.str
(ScopeType
), optionalDistance range scope of the store (choices:
'global'
,'regional'
,'local'
).
-
♦
waveform_type
¶ builtins.str
(WaveType
), optionalWave type stored (choices:
'full waveform'
,'bodywave'
,'P wave'
,'S wave'
,'surface wave'
).
-
♦
nearfield_terms
¶ builtins.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
¶ builtins.str
(ComponentScheme
), default:'elastic10'
GF component scheme (choices:
'elastic2'
,'elastic8'
,'elastic10'
,'elastic18'
,'elastic5'
,'poroelastic10'
).
-
♦
stored_quantity
¶ builtins.str
(QuantityType
), optionalPhysical quantity of stored values (choices:
'displacement'
,'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.
- lat – surface origin for coordinate system of
-
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.
- lat – surface origin for coordinate system of
-
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.
- lat – surface origin for coordinate system of
-
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.
- lat – surface origin for coordinate system of
-
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 depth - High level index variables:
(source_depth, receiver_distance, receiver_depth, 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].
- Symmetries like in
-
class
ConfigTypeC
(**kwargs)[source]¶ Bases:
pyrocko.gf.meta.Config
No symmetrical constraints but fixed receiver positions.
- Cartesian 3D source volume around a reference point
- High level index variables:
(ireceiver, 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
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
¶ builtins.str
(pyrocko.guts.StringChoice
), optional, default:'cos'
-
♦
-
class
WaveformType
(dummy) → 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
WaveformType
(dummy) → 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
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
¶ builtins.str
(WaveformType
), default:'dis'
-
♦
sample_rate
¶ float
, optional
-
♦
-
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]
-
effective_latlon
¶ Property holding the offset-corrected lat/lon pair of the location.
-
effective_lat
¶ Property holding the offset-corrected latitude of the location.
-
effective_lon
¶ Property holding the offset-corrected longitude of the location.
-
♦