pymarchenko.mme.MME#
- class pymarchenko.mme.MME(R, wav, wav_c=None, dt=0.004, nt=None, dr=1.0, nfmax=None, toff=0.0, nsmooth=10, dtype='float64', saveRt=True, prescaled=False)[source]#
Marchenko Multiple Elimination
Solve multi-dimensional Marchenko Multiple Elimination problem using Neumann iterative substitution.
- Parameters:
- R
numpy.ndarray
Multi-dimensional reflection response in time or frequency domain of size \([n_s \times n_r \times n_t/n_{fmax}]\). If provided in time, R should not be of complex type. If provided in frequency, R should contain the positive time axis followed by the negative one. Note that the reflection response should have already been multiplied by 2.
- wav
numpy.ndarray
Wavelet to apply to the reflection response shot gather used as initial guess
- wav_c
int
, optional Index of center of wavelet. If
None
the middle sample is used.- dt
float
, optional Sampling of time integration axis
- nt
float
, optional Number of samples in time (not required if
R
is in time)- dr
float
, optional Sampling of receiver integration axis
- nfmax
int
, optional Index of max frequency to include in deconvolution process
- toff
float
, optional Time-offset to apply to traveltime
- nsmooth
int
, optional Number of samples of smoothing operator to apply to window
- dtype
bool
, optional Type of elements in input array.
- saveRt
bool
, optional Save
R
andR^H
to speed up the computation of adjoint ofpylops.signalprocessing.Fredholm1
(True
) or createR^H
on-the-fly (False
) Note thatsaveRt=True
will be faster but double the amount of required memory- prescaled
bool
, optional Apply scaling to
R
(False
) or not (False
) when performing spatial and temporal summations within thepylops.waveeqprocessing.MDC
operator. In caseprescaled=True
, theR
is assumed to have been pre-scaled by the user.
- R
- Raises:
- TypeError
If
t
is notnumpy.ndarray
.
Notes
Marchenko Multiple Elimination is a method that allows to produce a primary-only reflection data by repeated filtering of the data with itself and windowing [1].
The projected Marchenko equations can be written in compact matrix-vector notation as:
\[\begin{split}\mathbf{v^-} = \Theta \mathbf{R} + \Theta \mathbf{R} \mathbf{v^+} \\ \mathbf{v^+} = \Theta \mathbf{R^*} \mathbf{v^-}\end{split}\]and solved for \(\mathbf{v^+}\) via Neumann iterative substitution:
\[\mathbf{v^+} = \sum_{k=0}^\inf (\Theta \mathbf{R^*}\Theta \mathbf{R})^k \Theta \mathbf{R^*} \Theta R = \sum_{k=1}^\inf (\Theta \mathbf{R^*} \Theta \mathbf{R})^k \delta\]where \(R\) is the reflection response (or \(\delta\) is a spatio-temporal delta) and \(\Theta\) is the window. The MME algorithm requires a small time offset
toff
(on the order of the wavelet lenght) to be subtracted to the window time, whilst the TMME requires such time offset to be added.At this point the projected demultipled data \(\mathbf{U^-}\) is computed and its values at time sample \(t=2t_d\) is extracted:
\[\mathbf{U^-} = \mathbf{R} + \mathbf{R} \mathbf{v^+} = \mathbf{R} + \mathbf{R} \sum_{k=1}^\inf (\Theta \mathbf{R^*}\Theta \mathbf{R})^k \delta\]If we repeat the same procedure for all possible $t=2t_d$, the retrived dataset is deprived of all internal multiples.
[1]Zhang, L., Thorbecke, J., Wapenaar, K., and Slob, E., “Data-driven internal multiple elimination and its consequences for imaging: A comparison of strategies”, Geophysics, vol. 84, pp. S365–S372. 2019.
[2]Zhang, L., Thorbecke, J., Wapenaar, K., and Slob, E., “Transmission compensated primary reflection retrieval in the data domain and consequences for imaging”, Geophysics, vol. 84, pp. Q27–Q36. 2019.
- Attributes:
Methods
__init__
(R, wav[, wav_c, dt, nt, dr, nfmax, ...])apply_multisrc
(Rsrcs[, usematmul, trcomp, ...])Marchenko Multiple elimination for multiple shot gathers
apply_onesrc
(Rsrc[, usematmul, trcomp, ...])Marchenko Multiple elimination for one shot gather
Examples using pymarchenko.mme.MME
#
5. Marchenko Multiple Elimination