pyadjoint
- copyright
adjTomo Dev Team (adjtomo@gmail.com), 2022 Lion Krischer (krischer@geophysik.uni-muenchen.de), 2015
- license
BSD 3-Clause (“BSD New” or “BSD Simplified”)
Subpackages
Submodules
Package Contents
Classes
Adjoint Source class to hold calculated adjoint sources |
Functions
|
Central function of Pyadjoint used to calculate adjoint sources and misfit. |
Helper function returning example data for Pyadjoint. |
|
|
Generic plotting function for adjoint sources and data. |
|
Defines two common parameters for all configuration objects and then |
|
Wrapper for getting the correct adjoint source function based on the |
Attributes
- exception pyadjoint.PyadjointError[source]
Bases:
Exception
Base class for all Pyadjoint exceptions. Will probably be used for all exceptions to not overcomplicate things as the whole package is pretty small.
- exception pyadjoint.PyadjointWarning[source]
Bases:
UserWarning
Base class for all Pyadjoint warnings.
- class pyadjoint.AdjointSource(adjsrc_type, misfit, dt, min_period, max_period, component, adjoint_source=None, windows=None, network=None, station=None, location=None, starttime=None, window_stats=None)[source]
Adjoint Source class to hold calculated adjoint sources
- __str__()
Return str(self).
- write(filename, format, **kwargs)
Write the adjoint source to a file.
- Parameters
SPECFEM
SPECFEM requires one additional parameter: the temporal offset of the first sample in seconds. The following example sets the time of the first sample in the adjoint source to
-10
.>>> adj_src.write("NET.STA.CHAN.adj", format="SPECFEM", ... time_offset=-10)
Adjoint sources can be written directly to an ASDFDataSet if provided. Optional
coordinates
parameter specifies the location of the station that generated the adjoint source>>> adj_src.write(ds, format="ASDF", time_offset=-10, ... coordinates={'latitude':19.2, ... 'longitude':13.4, ... 'elevation_in_m':2.0})
- _write_specfem(filename, time_offset)
Write the adjoint source for SPECFEM.
- _write_asdf(ds, time_offset, coordinates=None, **kwargs)
Writes the adjoint source to an ASDF file.
- Parameters
ds (str) – The ASDF data structure read in using pyasdf.
time_offset (float) – The temporal offset of the first sample in seconds. This is required if using the adjoint source as input to SPECFEM.
coordinates (list) – If given, the coordinates of the adjoint source. The ‘latitude’, ‘longitude’, and ‘elevation_in_m’ of the adjoint source must be defined.
- pyadjoint.calculate_adjoint_source(observed, synthetic, config, windows, plot=False, plot_filename=None, observed_2=None, synthetic_2=None, windows_2=None, **kwargs)[source]
Central function of Pyadjoint used to calculate adjoint sources and misfit.
This function uses the notion of observed and synthetic data to offer a nomenclature most users are familiar with. Please note that it is nonetheless independent of what the two data arrays actually represent.
The function tapers the data from
left_window_border
toright_window_border
, both in seconds since the first sample in the data arrays.- Parameters
observed (
obspy.core.trace.Trace
) – The observed data.synthetic (
obspy.core.trace.Trace
) – The synthetic data.config (configuration parameters that control measurement) –
pyadjoint.config.Config
windows (list of tuples) – [(left, right),…] representing left and right window borders to be tapered in units of seconds since first sample in data array
plot (bool or empty
matplotlib.figure.Figure
instance) – Also produce a plot of the adjoint source. This will force the adjoint source to be calculated regardless of the value ofadjoint_src
.plot_filename (str) – If given, the plot of the adjoint source will be saved there. Only used if
plot
isTrue
.observed_2 (obspy.core.trace.Trace) – second observed waveform to calculate adjoint sources from station pairs
synthetic_2 (obspy.core.trace.Trace) – second synthetic waveform to calculate adjoint sources from station pairs
windows_2 (list of tuples) – [(left, right),…] representing left and right window borders to be tapered in units of seconds since first sample in data array. Used to window observed_2 and synthetic_2
- pyadjoint.get_example_data()[source]
Helper function returning example data for Pyadjoint.
The returned data is fully preprocessed and ready to be used with Pyadjoint.
- Returns
Tuple of observed and synthetic streams
- Return type
tuple of
obspy.core.stream.Stream
objects
Example
>>> from pyadjoint import get_example_data >>> observed, synthetic = get_example_data() >>> print(observed) 3 Trace(s) in Stream: SY.DBO.S3.MXR | 2014-11-15T02:31:50.259999Z - ... | 1.0 Hz, 3600 samples SY.DBO.S3.MXT | 2014-11-15T02:31:50.259999Z - ... | 1.0 Hz, 3600 samples SY.DBO.S3.MXZ | 2014-11-15T02:31:50.259999Z - ... | 1.0 Hz, 3600 samples >>> print(synthetic) 3 Trace(s) in Stream: SY.DBO..LXR | 2014-11-15T02:31:50.259999Z - ... | 1.0 Hz, 3600 samples SY.DBO..LXT | 2014-11-15T02:31:50.259999Z - ... | 1.0 Hz, 3600 samples SY.DBO..LXZ | 2014-11-15T02:31:50.259999Z - ... | 1.0 Hz, 3600 samples
- pyadjoint.plot_adjoint_source(observed, synthetic, adjoint_source, misfit, windows, adjoint_source_name)[source]
Generic plotting function for adjoint sources and data.
Many types of adjoint sources can be represented in the same manner. This is a convenience function that can be called by different the implementations for different adjoint sources.
- Parameters
observed (
obspy.core.trace.Trace
) – The observed data.synthetic (
obspy.core.trace.Trace
) – The synthetic data.adjoint_source (numpy.ndarray) – The adjoint source.
misfit (float) – The associated misfit value.
windows (list of tuples) – [(left, right),…] representing left and right window borders to be tapered in units of seconds since first sample in data array
adjoint_source_name (str) – The name of the adjoint source.
- pyadjoint.get_config(adjsrc_type, min_period, max_period, **kwargs)[source]
Defines two common parameters for all configuration objects and then reassigns self to a sub Config class which dictates its own required parameters
- pyadjoint.get_function(adjsrc_type)[source]
Wrapper for getting the correct adjoint source function based on the adjsrc_type. Many adjoint sources share functions with different flags so this function takes care of the logic of choosing which.
- Parameters
adjsrc_type (str) – choice of adjoint source
- Return type
function
- Returns
calculate_adjoint_source function for the correct adjoint source type