pysource package
Submodules
FD_utils module
FDutils.Rmat(model)
Rotation matrix according to tilt and asymut.
Parameters
model (Model) – Model structure
FDutils.divs(func, sofact=1, side=- 1)
GrDivergenceadient shifted by half a grid point, only to be used in combination with grads.
FDutils.grads(func, sofact=1, side=1)
Gradient shifted by half a grid point, only to be used in combination with divs.
FD_utils.laplacian(v, irho)
Laplacian with density div( 1/rho grad) (u)
FDutils.satti(u, v, model)
Tensor factorized SSA TTI wave equation spatial derivatives.
Parameters
- u (TimeFunction) – first TTI field
* **v** (*TimeFunction*) – second TTI field
* **model** (*Model*) – Model structure
FDutils.thomsenmat(model)
Diagonal Matrices with Thomsen parameters for vectorial temporaries computation.
Parameters
model (Model) – Model structure
checkpoint module
class checkpoint.CheckpointOperator(op, **kwargs)
Devito’s concrete implementation of the ABC pyrevolve.Operator. This class wraps devito.Operator so it conforms to the pyRevolve API. pyRevolve will call apply with arguments tstart and tend. Devito calls these arguments ts and te so the following dict is used to perform the translations between different names.
Parameters
- op (Operator) – devito.Operator object that this object will wrap.
* **args** (*dict*) – If devito.Operator.apply() expects any arguments, they can be provided
here to be cached. Any calls to CheckpointOperator.apply() will
automatically include these cached arguments in the call to the
underlying devito.Operator.apply().
apply(tstart, tend)
If the devito operator requires some extra arguments in the call to apply they can be stored in the args property of this object so pyRevolve calls pyRevolve.Operator.apply() without caring about these extra arguments while this method passes them on correctly to devito.Operator
targnames(_ = {'tend': 'timeM', 'tstart': 'timem'_ )
class checkpoint.DevitoCheckpoint(objects)
Devito’s concrete implementation of the Checkpoint abstract base class provided by pyRevolve. Holds a list of symbol objects that hold data.
property dtype()
data type
get_data(timestep)
returns the data (wavefield) for the time-step timestep
getdatalocation(timestep)
returns the data (wavefield) for the time-step timestep
load()
NotImplementedError
save()
NotImplementedError
property size()
The memory consumption of the data contained in a checkpoint.
checkpoint.getsymboldata(symbol, timestep)
Return the symbol corresponding to the data at time-step timestep
geom_utils module
geomutils.geomexpr(model, u, srccoords=None, reccoords=None, wavelet=None, fw=True, nt=None)
Generates the source injection and receiver interpolation. This function is fully abstracted and does not care whether this is a forward or adjoint wave-equation. The source is the source term of the equation The receiver is the measurment term
Therefore, for the adjoint, this function has to be called as: srcrec(model, v, srccoords=rec_coords, …) because the data is the sources
Parameters
- model (Model) – Physical model
* **u** (*TimeFunction** or **tuple*) – Wavefield to inject into and read from
* **src_coords** (*Array*) – Physical coordinates of the sources
* **rec_coords** (*Array*) – Physical coordinates of the receivers
* **wavelet** (*Array*) – Data for the source
* **fw=True** – Whether the direction is forward or backward in time
* **nt** (*int*) – Number of time steps
geomutils.srcrec(model, u, srccoords=None, reccoords=None, wavelet=None, nt=None)
interface module
interface.Jadjoint(model, srccoords, wavelet, reccoords, recin, spaceorder=8, isresidual=False, checkpointing=False, ncheckpoints=None, tsub=1, returnobj=False, freqlist=[], dftsub=None, ic='as', ws=None, f0=0.015, born_fwd=False, nlind=False, misfit=None)
Jacobian (adjoint fo born modeling operator) operator on a shot record as a source (i.e data residual). Supports three modes: * Checkpinting * Frequency compression (on-the-fly DFT) * Standard zero lag cross correlation over time
Parameters
- model (Model) – Physical model
* **src_coords** (*Array*) – Coordiantes of the source(s)
* **wavelet** (*Array*) – Source signature
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **recin** (*Array*) – Receiver data
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **checkpointing** (*Bool*) – Whether or not to use checkpointing
* **n_checkpoints** (*Int*) – Number of checkpoints for checkpointing
* **maxmem** (*Float*) – Maximum memory to use for checkpointing
* **freq_list** (*List*) – List of frequencies for on-the-fly DFT
* **dft_sub** (*Int*) – Subsampling factor for on-the-fly DFT
* **ic** (*String*) – Imaging conditions (“as”, “isic” or “fwi”), defaults to “as”
* **ws** (*Array*) – Extended source spatial distribution
* **f0** (*peak frequency*) –
Returns
Adjoint jacobian on the input data (gradient)
Return type
Array
interface.Jadjointcheckpointing(model, srccoords, wavelet, reccoords, recin, spaceorder=8, isresidual=False, ncheckpoints=None, bornfwd=False, return_obj=False, ic='as', ws=None, nlind=False, f0=0.015, misfit=None)
Jacobian (adjoint fo born modeling operator) operator on a shot record as a source (i.e data residual). Outputs the gradient with Checkpointing.
Parameters
- model (Model) – Physical model
* **src_coords** (*Array*) – Coordiantes of the source(s)
* **wavelet** (*Array*) – Source signature
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **recin** (*Array*) – Receiver data
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **checkpointing** (*Bool*) – Whether or not to use checkpointing
* **n_checkpoints** (*Int*) – Number of checkpoints for checkpointing
* **maxmem** (*Float*) – Maximum memory to use for checkpointing
* **ic** (*String*) – Imaging conditions (“as”, “isic” or “fwi”), defaults to “as”
* **ws** (*Array*) – Extended source spatial distribution
* **is_residual** (*Bool*) – Whether to treat the input as the residual or as the observed data
* **born_fwd** (*Bool*) – Whether to use the forward or linearized forward modeling operator
* **nlind** (*Bool*) – Whether to remove the non linear data from the input data. This option is
only available in combination with born_fwd
* **f0** (*peak frequency*) –
Returns
Adjoint jacobian on the input data (gradient)
Return type
Array
interface.Jadjointfreq(model, srccoords, wavelet, reccoords, recin, spaceorder=8, freqlist=[], isresidual=False, returnobj=False, nlind=False, dftsub=None, ic='as', ws=None, bornfwd=False, f0=0.015, misfit=None)
Jacobian (adjoint fo born modeling operator) operator on a shot record as a source (i.e data residual). Outputs the gradient with Frequency compression (on-the-fly DFT).
Parameters
- model (Model) – Physical model
* **src_coords** (*Array*) – Coordiantes of the source(s)
* **wavelet** (*Array*) – Source signature
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **recin** (*Array*) – Receiver data
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **freq_list** (*List*) – List of frequencies for on-the-fly DFT
* **dft_sub** (*Int*) – Subsampling factor for on-the-fly DFT
* **ic** (*String*) – Imaging conditions (“as”, “isic” or “fwi”), defaults to “as”
* **ws** (*Array*) – Extended source spatial distribution
* **is_residual** (*Bool*) – Whether to treat the input as the residual or as the observed data
* **born_fwd** (*Bool*) – Whether to use the forward or linearized forward modeling operator
* **nlind** (*Bool*) – Whether to remove the non linear data from the input data. This option is
only available in combination with born_fwd
* **f0** (*peak frequency*) –
Returns
Adjoint jacobian on the input data (gradient)
Return type
Array
interface.Jadjointstandard(model, srccoords, wavelet, reccoords, recin, spaceorder=8, isresidual=False, returnobj=False, bornfwd=False, ic='as', ws=None, t_sub=1, nlind=False, f0=0.015, misfit=None)
Adjoint Jacobian (adjoint fo born modeling operator) operator on a shot record as a source (i.e data residual). Outputs the gradient with standard zero lag cross correlation over time.
Parameters
- model (Model) – Physical model
* **src_coords** (*Array*) – Coordiantes of the source(s)
* **wavelet** (*Array*) – Source signature
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **recin** (*Array*) – Receiver data
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **ic** (*String*) – Imaging conditions (“as”, “isic” or “fwi”), defaults to “as”
* **ws** (*Array*) – Extended source spatial distribution
* **is_residual** (*Bool*) – Whether to treat the input as the residual or as the observed data
* **born_fwd** (*Bool*) – Whether to use the forward or linearized forward modeling operator
* **nlind** (*Bool*) – Whether to remove the non linear data from the input data. This option is
only available in combination with born_fwd
* **f0** (*peak frequency*) –
Returns
Adjoint jacobian on the input data (gradient)
Return type
Array
interface.adjointnorec(model, reccoords, data, spaceorder=8, f0=0.015)
Adjoint/backward modeling of a shot record (receivers as source) without source sampling F^T*Pr^T*d_obs.
Parameters
- model (Model) – Physical model
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **data** (*Array*) – Shot gather
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
Adjoint wavefield
Return type
Array
interface.adjointrec(model, srccoords, reccoords, data, spaceorder=8, f0=0.015)
Adjoint/backward modeling of a shot record (receivers as source) Ps*F^T*Pr^T*d.
Parameters
- model (Model) – Physical model
* **src_coords** (*Array*) – Coordiantes of the source(s)
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **data** (*Array*) – Shot gather
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
Shot record (adjoint wavefield at source position(s))
Return type
Array
interface.adjointw(model, reccoords, data, wavelet, space_order=8, f0=0.015)
Adjoint/backward modeling of a shot record (receivers as source) for an extended source setup Pw*F^T*Pr^T*d_obs.
Parameters
- model (Model) – Physical model
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **data** (*Array*) – Shot gather
* **wavelet** (*Array*) – Time signature of the forward source for stacking along time
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
spatial distribution
Return type
Array
interface.adjointwfsrc(model, u, srccoords, spaceorder=8, f0=0.015)
Adjoint/backward modeling of a full wavefield (full wavefield as adjoint source) Ps*F^T*u.
Parameters
- model (Model) – Physical model
* **u** (*Array** or **TimeFunction*) – Time-space dependent source
* **src_coords** (*Array*) – Source coordinates
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
Shot record (sampled at source position(s))
Return type
Array
interface.adjointwfsrcnorec(model, u, spaceorder=8, f0=0.015)
Adjoint/backward modeling of a full wavefield (full wavefield as adjoint source) F^T*u.
Parameters
- model (Model) – Physical model
* **u** (*Array** or **TimeFunction*) – Time-space dependent source
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
Adjoint wavefield
Return type
Array
interface.bornrec(model, srccoords, wavelet, reccoords, spaceorder=8, ic='as', f0=0.015)
Linearized (Born) modeling of a point source for a model perturbation (square slowness) dm.
Parameters
- model (Model) – Physical model
* **src_coords** (*Array*) – Coordiantes of the source(s)
* **wavelet** (*Array*) – Source signature
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **ic** (*String*) – Imaging conditions (“as”, “isic” or “fwi”), defaults to “as”
* **f0** (*peak frequency*) –
Returns
Shot record
Return type
Array
interface.bornrecw(model, weight, wavelet, reccoords, spaceorder=8, ic='as', f0=0.015)
Linearized (Born) modeling of an extended source for a model perturbation (square slowness) dm with an extended source
Parameters
- model (Model) – Physical model
* **weight** (*Array*) – Spatial distriubtion of the extended source
* **wavelet** (*Array*) – Source signature
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **ic** (*String*) – Imaging conditions (“as”, “isic” or “fwi”), defaults to “as”
* **f0** (*peak frequency*) –
Returns
Shot record
Return type
Array
interface.forwardnorec(model, srccoords, wavelet, spaceorder=8, f0=0.015)
Forward modeling of a point source without receiver.
Parameters
- model (Model) – Physical model
* **src_coords** (*Array*) – Coordiantes of the source(s)
* **wavelet** (*Array*) – Source signature
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
Wavefield
Return type
Array
interface.forwardrec(model, srccoords, wavelet, reccoords, spaceorder=8, f0=0.015)
Forward modeling of a point source with receivers Pr*F*Ps^T*q.
Parameters
- model (Model) – Physical model
* **src_coords** (*Array*) – Coordiantes of the source(s)
* **wavelet** (*Array*) – Source signature
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
Shot record
Return type
Array
interface.forwardrecw(model, weight, wavelet, reccoords, spaceorder=8, f0=0.015)
Forward modeling of an extended source with receivers Pr*F*Pw^T*w
Parameters
- model (Model) – Physical model
* **weights** (*Array*) – Spatial distribution of the extended source.
* **wavelet** (*Array*) – Source signature
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
Shot record
Return type
Array
interface.forwardwfsrc(model, u, reccoords, spaceorder=8, f0=0.015)
Forward modeling of a full wavefield source Pr*F*u.
Parameters
- model (Model) – Physical model
* **u** (*TimeFunction** or **Array*) – Time-space dependent wavefield
* **rec_coords** (*Array*) – Coordiantes of the receiver(s)
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
Shot record
Return type
Array
interface.forwardwfsrcnorec(model, u, spaceorder=8, f0=0.015)
Forward modeling of a full wavefield source without receiver F*u.
Parameters
- model (Model) – Physical model
* **u** (*TimeFunction** or **Array*) – Time-space dependent wavefield
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
Wavefield
Return type
Array
interface.gradfwi(model, recin, reccoords, u, space_order=8, f0=0.015)
FWI gradient, i.e adjoint Jacobian on a data residual.
Parameters
- model (Model) – Physical model
* **recin** (*Array*) – Data residual
* **rec_coords** (*Array*) – Receivers coordinates
* **u** (*TimeFunction*) – Forward wavefield
* **space_order** (*Int** (**optional**)*) – Spatial discretization order, defaults to 8
* **f0** (*peak frequency*) –
Returns
FWI gradient
Return type
Array
interface.wrifunc(model, srccoords, wavelet, reccoords, recin, yin, spaceorder=8, ic='as', ws=None, tsub=1, grad='m', gradcorr=False, alphaop=False, wfun=None, eps=0, freq_list=[], wfilt=None, f0=0.015)
Time domain wavefield reconstruction inversion wrapper
kernels module
kernels.SLS2ndorder(model, p, fw=True, q=None, f0=0.015)
Viscoacoustic 2nd SLS wave equation. https://library.seg.org/doi/10.1190/geo2013-0030.1
Bulk modulus moved to rhs. The adjoint equation is directly derived as the discrete adjoint of the forward PDE which leads to a slightly different formulation than in the paper.
Parameters
- model (Model) – Physical model
* **u1** (*TimeFunction*) – Pressure field
* **u2** (*TimeFunction*) – Attenuation Memory variable
* **fw** (*Bool*) – Whether forward or backward in time propagation
* **q** (*TimeFunction** or **Expr*) – Full time-space source as a tuple (one value for each component)
* **f0** (*Peak frequency*) –
kernels.acoustic_kernel(model, u, fw=True, q=None)
Acoustic wave equation time stepper
Parameters
- model (Model) – Physical model
* **u** (*TimeFunction** or **tuple*) – wavefield (tuple if TTI)
* **fw** (*Bool*) – Whether forward or backward in time propagation
* **q** (*TimeFunction** or **Expr*) – Full time-space source
kernels.tti_kernel(model, u1, u2, fw=True, q=None)
TTI wave equation (one from my paper) time stepper
Parameters
- model (Model) – Physical model
* **u1** (*TimeFunction*) – First component (pseudo-P) of the wavefield
* **u2** (*TimeFunction*) – First component (pseudo-P) of the wavefield
* **fw** (*Bool*) – Whether forward or backward in time propagation
* **q** (*TimeFunction** or **Expr*) – Full time-space source as a tuple (one value for each component)
kernels.wave_kernel(model, u, fw=True, q=None, f0=0.015)
Pde kernel corresponding the the model for the input wavefield
Parameters
- model (Model) – Physical model
* **u** (*TimeFunction** or **tuple*) – wavefield (tuple if TTI or Viscoacoustic)
* **fw** (*Bool*) – Whether forward or backward in time propagation
* **q** (*TimeFunction** or **Expr*) – Full time-space source
* **f0** (*Peak frequency*) –
models module
class models.Model(origin, spacing, shape, m, space_order=2, nbl=40, dtype=<class 'numpy.float32'>, epsilon=None, delta=None, theta=None, phi=None, rho=None, b=None, qp=None, dm=None, fs=False, **kwargs)
The physical model used in seismic inversion
shape_pml = np.array(shape) + 2 \* self.nbl processes.
Parameters
- origin (tuple of floats) – Origin of the model in m as a tuple in (x,y,z) order.
* **spacing** (*tuple of floats*) – Grid size in m as a Tuple in (x,y,z) order.
* **shape** (*tuple of int*) – Number of grid points size in (x,y,z) order.
* **space_order** (*int*) – Order of the spatial stencil discretisation.
* **m** (*array_like** or **float*) – Squared slownes in s^2/km^2
* **nbl** (*int**, **optional*) – The number of absorbin layers for boundary damping.
* **dtype** (*np.float32** or **np.float64*) – Defaults to 32.
* **epsilon** (*array_like** or **float**, **optional*) – Thomsen epsilon parameter (0<epsilon<1).
* **delta** (*array_like** or **float*) – Thomsen delta parameter (0<delta<1), delta<epsilon.
* **theta** (*array_like** or **float*) – Tilt angle in radian.
* **phi** (*array_like** or **float*) – Asymuth angle in radian.
* **dt** (*Float*) – User provided computational time-step
property critical_dt()
Critical computational time step value from the CFL condition.
property dim()
Spatial dimension of the problem and model domain.
property dm()
Model perturbation for linearized modeling
property domain_size()
Physical size of the domain as determined by shape and spacing
property dt()
User provided dt
property dtype()
Data type for all assocaited data objects.
property is_tti()
Whether the model is TTI or isotopic
property is_viscoacoustic()
Whether the model is TTI or isotopic
property m()
Function holding the squared slowness in s^2/km^2.
property padsizes()
property physical_parameters()
List of physical parameteres
physical_params(**kwargs)
Return all set physical parameters and update to input values if provided
property space_dimensions()
Spatial dimensions of the grid
property space_order()
Spatial discretization order
property spacing()
Grid spacing for all fields in the physical model.
property spacing_map()
Map between spacing symbols and their values for each SpaceDimension.
property vp()
Symbolic representation of the velocity vp = sqrt(1 / m)
property zero_thomsen()
propagators module
propagators.adjoint(model, y, srccoords, rcvcoords, spaceorder=8, qwf=None, dftsub=None, save=False, ws=None, normv=False, wfun=None, freq_list=None, f0=0.015)
Low level propagator, to be used through interface.py Compute adjoint wavefield v = adjoint(F(m))*y and related quantities (||v||_w, v(xsrc))
propagators.born(model, srccoords, rcvcoords, wavelet, spaceorder=8, save=False, qwf=None, returnop=False, ic='as', freqlist=None, dftsub=None, ws=None, t_sub=1, nlind=False, f0=0.015)
Low level propagator, to be used through interface.py Compute linearized wavefield U = J(m)* δ m and related quantities.
propagators.forward(model, srccoords, rcvcoords, wavelet, spaceorder=8, save=False, qwf=None, returnop=False, freqlist=None, dftsub=None, ws=None, t_sub=1, f0=0.015, **kwargs)
Low level propagator, to be used through interface.py Compute forward wavefield u = A(m)^{-1}*f and related quantities (u(xrcv))
propagators.forwardgrad(model, srccoords, rcvcoords, wavelet, v, spaceorder=8, q=None, ws=None, ic='as', w=None, freq=None, f0=0.015, **kwargs)
Low level propagator, to be used through interface.py Compute forward wavefield u = A(m)^{-1}*f and related quantities (u(xrcv))
propagators.gradient(model, residual, rcvcoords, u, returnop=False, spaceorder=8, w=None, freq=None, dftsub=None, ic='as', f0=0.015, save=True)
Low level propagator, to be used through interface.py Compute the action of the adjoint Jacobian onto a residual J’* δ d.
sensitivity module
sensitivity.Loss(dsyn, dobs, dt, is_residual=False, misfit=None)
L2 loss and residual between the synthetic data dsyn and observed data dobs
Parameters
- dsyn (SparseTimeFunction** or **tuple) – Synthetic data or tuple (background, linearized) data
* **dobs** (*SparseTimeFunction*) – Observed data
* **dt** (*float*) – Time sampling rate
* **is_residual** (*bool*) – Whether input dobs is already the data residual
* **misfit** (*function*) – User provided function of the form:
misifit(dsyn, dob) -> obj, adjoint_source
sensitivity.basic_src(model, u, **kwargs)
Basic source for linearized modeling
Parameters
- model (Model) – Model containing the perturbation dm
* **u** (*TimeFunction** or **Tuple*) – Forward wavefield (tuple of fields for TTI or dft)
sensitivity.crosscorrfreq(u, v, model, freq=None, dftsub=None, **kwargs)
Standard cross-correlation imaging condition with on-th-fly-dft
Parameters
- u (TimeFunction** or **Tuple) – Forward wavefield (tuple of fields for TTI or dft)
* **v** (*TimeFunction** or **Tuple*) – Adjoint wavefield (tuple of fields for TTI)
* **model** (*Model*) – Model structure
* **freq** (*Array*) – Array of frequencies for on-the-fly DFT
* **factor** (*int*) – Subsampling factor for DFT
sensitivity.crosscorr_time(u, v, model, **kwargs)
Cross correlation of forward and adjoint wavefield
Parameters
- u (TimeFunction** or **Tuple) – Forward wavefield (tuple of fields for TTI or dft)
* **v** (*TimeFunction** or **Tuple*) – Adjoint wavefield (tuple of fields for TTI)
* **model** (*Model*) – Model structure
sensitivity.func_name(freq=None, ic='as')
Get key for imaging condition/linearized source function
sensitivity.fwi_freq(*ar, **kw)
sensitivity.fwi_src(*ar, **kw)
sensitivity.fwi_time(*ar, **kw)
sensitivity.gradexpr(gradm, u, v, model, w=None, freq=None, dftsub=None, ic='as')
Gradient update stencil
Parameters
- u (TimeFunction** or **Tuple) – Forward wavefield (tuple of fields for TTI or dft)
* **v** (*TimeFunction** or **Tuple*) – Adjoint wavefield (tuple of fields for TTI)
* **model** (*Model*) – Model structure
* **w** (*Float** or **Expr** (**optional**)*) – Weight for the gradient expression (default=1)
* **freq** (*Array*) – Array of frequencies for on-the-fly DFT
* **factor** (*int*) – Subsampling factor for DFT
* **isic** (*Bool*) – Whether or not to use inverse scattering imaging condition (not supported yet)
sensitivity.inner_grad(u, v)
Inner product of the gradient of two Function.
Parameters
- u (TimeFunction** or **Function) – First wavefield
* **v** (*TimeFunction** or **Function*) – Second wavefield
sensitivity.isic_freq(u, v, model, **kwargs)
Inverse scattering imaging condition
Parameters
- u (TimeFunction** or **Tuple) – Forward wavefield (tuple of fields for TTI or dft)
* **v** (*TimeFunction** or **Tuple*) – Adjoint wavefield (tuple of fields for TTI)
* **model** (*Model*) – Model structure
sensitivity.isic_src(model, u, **kwargs)
ISIC source for linearized modeling
Parameters
- model (Model) – Model containing the perturbation dm
* **u** (*TimeFunction** or **Tuple*) – Forward wavefield (tuple of fields for TTI or dft)
sensitivity.isic_time(u, v, model, **kwargs)
Inverse scattering imaging condition
Parameters
- u (TimeFunction** or **Tuple) – Forward wavefield (tuple of fields for TTI or dft)
* **v** (*TimeFunction** or **Tuple*) – Adjoint wavefield (tuple of fields for TTI)
* **model** (*Model*) – Model structure
sensitivity.lin_src(model, u, ic='as')
Source for linearized modeling
Parameters
- model (Model) – Model containing the perturbation dm
* **u** (*TimeFunction** or **Tuple*) – Forward wavefield (tuple of fields for TTI or dft)
* **ic** (*String*) – Imaging condition of which we compute the linearized source
sources module
class sources.PointSource(*args, **kwargs)
Symbolic data object for a set of sparse point sources
Parameters
- name (String) – Name of the symbol representing this source
* **grid** (*Grid*) – Grid object defining the computational domain.
* **coordinates** (*Array*) – Point coordinates for this source
* **data** (*(**Optional**) **Data*) – values to initialise point data
* **ntime** (*Int** (**Optional**)*) – Number of timesteps for which to allocate data
* **npoint** (*Int** (**Optional**)*) –
* **of sparse points represented by this source** (*Number*) –
* **dimension** (*Dimension** (**Optional**)*) – object for representing the number of points in this source
* **either the dimensions ntime and npoint**** or ****the fully** (*Note**,*) –
* **data array need to be provided.** (*initialised*) –
defaultassumptions( = {'commutative': True, 'complex': True, 'extendedreal': True, 'finite': True, 'hermitian': True, 'imaginary': False, 'infinite': False, 'real': True )
iscommutative( = Tru_ )
iscomplex( = Tru_ )
isextendedreal(_ = Tru_ )
isfinite( = Tru_ )
ishermitian( = Tru_ )
isimaginary( = Fals_ )
isinfinite( = Fals_ )
isreal( = Tru_ )
sources.Receiver()
alias of sources.PointSource
class sources.RickerSource(*args, **kwargs)
Symbolic object that encapsulate a set of sources with a pre-defined Ricker wavelet: http://subsurfwiki.org/wiki/Ricker_wavelet name: Name for the resulting symbol grid: Grid
object defining the computational domain. f0: Peak frequency for Ricker wavelet in kHz time: Discretized values of time in ms
defaultassumptions( = {'commutative': True, 'complex': True, 'extendedreal': True, 'finite': True, 'hermitian': True, 'imaginary': False, 'infinite': False, 'real': True )
iscommutative( = Tru_ )
iscomplex( = Tru_ )
isextendedreal(_ = Tru_ )
isfinite( = Tru_ )
ishermitian( = Tru_ )
isimaginary( = Fals_ )
isinfinite( = Fals_ )
isreal( = Tru_ )
wavelet(timev)
class sources.TimeAxis(start=None, step=None, num=None, stop=None)
Data object to store the TimeAxis. Exactly three of the four key arguments must be prescribed. Because of remainder values it is not possible to create a TimeAxis that exactly adhears to the inputs therefore start, stop, step and num values should be taken from the TimeAxis object rather than relying upon the input values. The four possible cases are: * start is None: start = step*(1 - num) + stop * step is None: step = (stop - start)/(num - 1) * num is None: num = ceil((stop - start + step)/step) and because of remainder stop = step*(num - 1) + start * stop is None: stop = step*(num - 1) + start
Parameters
- start (float**, **optional) – Start of time axis.
* **step** (*float**, **optional*) – Time interval.
* **num** (*int**, **optional*) – Number of values (Note: this is the number of intervals + 1).
Stop value is reset to correct for remainder.
* **stop** (*float**, **optional*) – End time.