API reference

Classes

Earth model

class pyprop8.LayeredStructureModel(layers=None, interface_depth_form=False)[source]

Construct an object to represent an earth model.

Parameters
  • layers (list or None) – Description of the layers that comprise the model. If interface_depth_form=False (default), list of tuples (thickness, vp, vs, rho) defining properties of each layer in turn, starting from the top; final entry must have thickness np.inf. Alternatively if interface_depth_form=True, list of tuples (ztop, vp, vs, rho) defining locations of interfaces within the model and the properties immediately below that interface; in this case one entry must have ztop=0. If None, ‘empty’ model is returned.

  • interface_depth_form (bool) – Select manner in which layers is specified; see above.

Each instance provides the following variables, which should be treated as read-only:

Variables
  • nlayers (int or None) – The number of layers in the model (including the infinite halfspace).

  • dz (numpy.ndarray or None) – Array of layer thicknesses, from top to bottom.

  • sigma (numpy.ndarray or None) – Array of P-wave moduli in each layer, from top to bottom.

  • mu (numpy.ndarray or None) – Array of S-wave moduli in each layer, from top to bottom.

  • rho (numpy.ndarray or None) – Array of densities in each layer, from top to bottom.

copy()[source]

Make a copy of the current model.

Return type

LayeredStructureModel

property vp

P-wave velocity in each layer (read-only).

Type

numpy.ndarray

property vs

S-wave velocity in each layer (read-only).

Type

numpy.ndarray

Receivers

class pyprop8.ListOfReceivers(xx, yy, depth=0, geometry='cartesian')[source]

Create a set pf receivers at known locations.

Parameters
  • xx (np.ndarray) –

  • yy (np.ndarray) – Arrays containing the x and y coordinates of each receiver. See geometry parameter for details on interpretation.

  • depth (float) – The depth of burial of receivers; zero for surface records. Note that all receivers must be at the same depth: if multiple depths are required, it will be necessary to split the receivers into groups and perform separate simulations for each.

  • geometry (str) – If geometry=’cartesian’, x and y coordinates are assumed to be given relative to a Cartesian grid with origin at some arbitrary point determined by the user. Alternatively, if geometry=’spherical’, the x and y coordinates are assumed to be given in degrees longitude/latitude. Note that the actual simulation remains in a Cartesian geometry with a flat-Earth approximation.

copy()[source]

Make a copy of the current receivers.

Return type

ListOfReceivers

property nstations

The total number of receivers present.

Type

int

class pyprop8.RegularlyDistributedReceivers(rmin=30, rmax=200, nr=18, phimin=0, phimax=360, nphi=36, depth=0, x0=0, y0=0, degrees=True)[source]

Create a set of receivers distributed regularly in cylindrical polar coordinates. Receivers lie on a set of equi-distant concentric circles centred on the origin, and at equally-distributed azimuths. This regularity enables faster computation and is useful if general characterisation of the wavefield is required.

Parameters
  • rmin (float) – Minimum radius to generate (inclusive; km)

  • rmax (float) – Maximum radius to generate (inclusive; km)

  • nr (int) – Number of radii to generate

  • phimin (float) – Minimum angle to generate.

  • phimax (float) – Maximum angles to generate. Angles are measured counter-clockwise from the x (East) axis, when viewed from above. They may be specified in degrees or radians (see degrees option).

  • nphi (int) – Number of angles to consider

  • depth (float) – Depth of burial of receivers (km)

  • x0 (float) –

  • y0 (float) – Coordinates of the origin point about which the receivers are to be distributed. These are used to verify the requirement that the source location coincide with the symmetry axis, and to generate appropriate Cartesian coordinates when .to_xy() is called.

  • degrees (bool) – True if angles are specified in degrees; False if in radians.

asListOfReceivers()[source]

Convert to ListOfReceivers object.

Return type

ListOfReceivers

copy()[source]

Make a copy of the the current receivers.

Return type

RegularlyDistributedReceivers

property nstations

The total number of receivers present.

Type

int

Seismic source

class pyprop8.PointSource(x, y, dep, Mxyz, F, time)[source]

Create an object representing a seismic source. Source mechanism represented as the sum of a moment tensor and a force, acting at a single point in space. Alternatively N distinct moment tensor/force pairs may be specified; in this case, N distinct simulations are to be performed.

Parameters
  • x (float) – The x-location of the source.

  • y (float) – The y-location of the source. Coordinate system must match that used for specifying receivers.

  • dep (float) – The depth of the source, km. Must be positive (or zero).

  • Mxyz (numpy.ndarray) – The moment tensor, expressed relative to a Cartesian system. Shape 3x3 or Nx3x3.

  • F (numpy.ndarray) – The force vector, expressed relative to a Cartesian system. Shape 1x3 or Nx1x3.

  • time (float or datetime.datetime) – The event time, expressed either as an instance of datetime.datetime, or as seconds relative to some arbitrary epoch.

copy()[source]

Make a copy of the current source.

Derivatives

class pyprop8.DerivativeSwitches(moment_tensor=False, force=False, r=False, phi=False, x=False, y=False, z=False, time=False, thickness=False, structure=None)[source]

Object to specify the set of derivatives that are sought.

Parameters
  • moment_tensor (bool) – Compute derivatives wrt six independent components of the moment tensor.

  • force (bool) – Compute derivatives wrt three components of the force vector.

  • r (bool) – Compute derivatives wrt source-receiver distance.

  • phi (bool) – Compute derivatives wrt source-receiver azimuth.

  • x (bool) – Compute derivatives wrt x (east-west) location of source.

  • y (bool) – Compute derivatives wrt y (north-south) location of source.

  • z (bool) – Compute derivatives wrt z (vertical) location of source.

  • time (bool) – Compute derivatives wrt event time.

  • thickness (bool) – Compute derivatives wrt layer thicknesses

  • structure (LayeredStructureModel or None) – Current model (only required when layer-thickness derivatives are sought).

Main computational routines

The following routines perform simulations of seismic response.

pyprop8.compute_spectra(structure, source, stations, omegas, derivatives=None, show_progress=True, stencil=<function kIntegrationStencil>, stencil_kwargs={'kmax': 2.04, 'kmin': 0, 'nk': 1200}, squeeze_outputs=True)[source]

Calculate and return velocity spectra for a given source and earth model at specified locations.

Parameters
  • structure (LayeredStructureModel) – The earth model within which calculations should be performed.

  • source (PointSource) – The source(s) for which calculations should be performed

  • stations (ListOfReceivers or RegularlyDistributedReceivers) – The locations for which spectra should be generated.

  • omegas (numpy.ndarray) – Freqencies at which spectrum is to be computed. Array may be complex and should have shape (n, ).

  • derivatives (DerivativeSwitches or None) – Determines which derivatives are computed and returned. See also discussion of return value, below.

  • show_progress (bool) – Display progress bars if available.

  • stencil (callable) –

    Function to generate quadrature points/weights for performing integral over wavenumber, k. Pass a callable satisfying:

    kk, wts = stencil(**stencil_kwargs)
    

    where kk represents sampling points and wts the associated quadrature weights.

  • stencil_kwargs (dict) – Arguments that will be passed to stencil()

  • squeeze_outputs (bool) – If true, apply numpy.squeeze() to all output arrays to eliminate dimensions of size ‘1’.

The output of compute_spectra depends on the value of the derivatives parameter.

Returns

If derivatives = None then compute_spectra returns a single array, spectra.Otherwise it returns a tuple of two arrays, (spectra, deriv).

The shapes of spectra and deriv depend on the nature of the object passed to stations.

  • If stations is an instance of ListOfReceivers, spectra will have shape (source.nsources, receivers.nstations, 3, nfreq), where source.nsources is the number of moment tensor/source vector pairs specified within the source object, receivers.nstations is the total number of receivers, and nfreq is the number of frequency points at which evaluation was requested (i.e., size(omegas)). The third dimension indexes the three components of motion: radial, transverse, and vertical. deriv will have shape (source.nsources, receivers.nstations,derivatives.nderivs,3,nfreq) where derivatives.nderivs is the total number of derivative components requested within the DerivativeSwitches object.

  • If stations is an instance of RegularlyDistributedReceivers, spectra will have shape (source.nsources, receivers.nr, receivers.nphi, 3, nfreq), where receivers.nr and receivers.nphi are the number of grid points in the radial and azimuthal directions, respectively. deriv will have shape (source.nsources, receivers.nr, receivers.nphi, derivatives.nderiv, 3, nfreq).

In all cases, if squeeze_outputs=True (default), then numpy.squeeze() will be called, discarding dimensions with only one entry.

pyprop8.compute_seismograms(structure, source, stations, nt, dt, alpha=None, source_time_function=None, pad_frac=0.5, xyz=True, derivatives=None, show_progress=True, squeeze_outputs=True, number_of_processes=1, **kwargs)[source]

Calculate and return displacement seismograms for a given source and earth model at specified locations.

Parameters
  • structure (LayeredStructureModel) – The earth model within which calculations should be performed.

  • source (PointSource) – The source(s) for which calculations should be performed

  • stations (ListOfReceivers or RegularlyDistributedReceivers) – The locations for which seismograms should be generated.

  • nt (int) – Number of time-series points to be returned.

  • dt (float) – Desired interval (seconds) between successive time series points (i.e. the inverse of sampling frequency). Total time series length is thus (nt-1)*dt seconds.

  • alpha (float or None) – Shift applied to integration contour to avoid singularities on real axis. If None, this is determined according to the rule-of-thumb given in O’Toole & Woodhouse (2011).

  • source_time_function (callable or None) – Function representing the spectrum of a source time-function to be applied to the seismograms. Function should take a single argument, representing angular frequency, and return a single complex number representing the amplitude of the spectrum at that point. If None, no source time-function is applied.

  • pad_frac (float) – Simulations are performed for a longer time series, which is then truncated to the desired length. This helps improve accuracy and stability. pad_frac determines the additional padding used, specified as a proportion of the desired time series length.

  • xyz (bool) – If True, the three components of each seismogram will correspond to motion relative to Cartesian x/y/z axes. If Fales, the three components will be expressed in polar coordinates: radial/transverse/vertical.

  • derivatives (DerivativeSwitches or None) – Determines which derivatives are computed and returned. See also discussion of return value, below.

  • show_progress (bool) – Display progress bars if available.

  • squeeze_outputs (bool) – If true, apply numpy.squeeze() to all output arrays to eliminate dimensions of size ‘1’.

  • number_of_processes (int) – The number of processes to use while performing computations. If >1, the multiprocessing module will be used, splitting up the frequency band across processs. Due to the cost of creating and managing separate processes, more is not always faster.

  • **kwargs (bool) – Any additional keyword options will be passed to compute_spectra().

The output of compute_seismograms depends on value of the derivatives parameter.

Returns

If derivatives = None then compute_spectra returns a tuple, (tt, seis). Otherwise compute_spectra returns a tuple, (tt, seis, deriv).

Here, tt is a numpy.ndarray containing the sequence of time points for which the seismograms have been evaluated. It will have shape (nt,).

The shapes of seis and deriv depend on the nature of the object passed to stations:

  • If stations is an instance of ListOfReceivers, seis will have shape (source.nsources, receivers.nstations, 3, nt), where source.nsources is the number of moment tensor/source vector pairs specified within the source object, receivers.nstations is the total number of receivers, and nt is the number of time points at which evaluation was requested. The third dimension indexes the three components of motion, governed by the xyz option. deriv will have shape (source.nsources, receivers.nstations,derivatives.nderivs,3,nt) where derivatives.nderivs is the total number of derivative components requested within the DerivativeSwitches object.

  • If stations is an instance of RegularlyDistributedReceivers, seis will have shape (source.nsources, receivers.nr, receivers.nphi, 3, nt), where receivers.nr and receivers.nphi are the number of grid points in the radial and azimuthal directions, respectively. deriv will have shape (source.nsources, receivers.nr, receivers.nphi, derivatives.nderiv, 3, nt).

In all cases, if squeeze_outputs=True (default), then numpy.squeeze() will be called, discarding dimensions with only one entry.

pyprop8.compute_static(structure, source, stations, los_vector=array([[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0]]), derivatives=None, squeeze_outputs=True, **kwargs)[source]

Calculate and return static offset measurements for a given source and earth model at specified locations.

Parameters
  • structure (LayeredStructureModel) – The earth model within which calculations should be performed.

  • source (PointSource) – The source(s) for which calculations should be performed

  • stations (ListOfReceivers or RegularlyDistributedReceivers) – The locations for which seismograms should be generated.

  • los_vector (numpy.ndarray) – Vector(s) defining ‘line(s) of sight’ along which static displacement should be measured. Should be expressed relative to a Cartesian basis and have shape (3) or (3, nlos). Default will return displacements relative to Cartesian basis.

  • derivatives (DerivativeSwitches or None) – Determines which derivatives are computed and returned. See also discussion of return value, below.

  • show_progress (bool) – Display progress bars if available.

  • squeeze_outputs (bool) – If true, apply numpy.squeeze() to all output arrays to eliminate dimensions of size ‘1’.

  • **kwargs (bool) – Any additional keyword options will be passed to compute_spectra().

The output of compute_static depends on value of the derivatives parameter.

Returns

If derivatives = None then compute_spectra returns a single array, static. Otherwise compute_spectra returns a tuple, (static, deriv).

The shapes of static and deriv depend on the nature of the object passed to stations:

  • If stations is an instance of ListOfReceivers, static will have shape (source.nsources, receivers.nstations, nlos), where source.nsources is the number of moment tensor/source vector pairs specified within the source object, receivers.nstations is the total number of receivers, and nlos is the number of lines-of-sight for which evaluation was requested. deriv will have shape (source.nsources, receivers.nstations,derivatives.nderivs,nlos) where derivatives.nderivs is the total number of derivative components requested within the DerivativeSwitches object.

  • If stations is an instance of RegularlyDistributedReceivers, seis will have shape (source.nsources, receivers.nr, receivers.nphi, nlos), where receivers.nr and receivers.nphi are the number of grid points in the radial and azimuthal directions, respectively. deriv will have shape (source.nsources, receivers.nr, receivers.nphi, derivatives.nderiv, nlos).

In all cases, if squeeze_outputs=True (default), then numpy.squeeze() will be called, discarding dimensions with only one entry.