def PrestackInversion(data, theta, wav, m0=None, linearization='akirich', explicit=False, simultaneous=False, epsI=None, epsR=None, dottest=False, returnres=False, epsRL1=None, kind='centered', **kwargs_solver): r"""Pre-stack linearized seismic inversion. Invert pre-stack seismic operator to retrieve a set of elastic property profiles from band-limited seismic pre-stack data (i.e., angle gathers). Depending on the choice of input parameters, inversion can be trace-by-trace with explicit operator or global with either explicit or linear operator. Parameters ---------- data : :obj:`np.ndarray` Band-limited seismic post-stack data of size :math:`[(n_{lins} \times) n_{t0} \times n_{\theta} (\times n_x \times n_y)]` theta : :obj:`np.ndarray` Incident angles in degrees wav : :obj:`np.ndarray` Wavelet in time domain (must had odd number of elements and centered to zero) m0 : :obj:`np.ndarray`, optional Background model of size :math:`[n_{t0} \times n_{m} (\times n_x \times n_y)]` linearization : :obj:`str` or :obj:`list`, optional choice of linearization, ``akirich``: Aki-Richards, ``fatti``: Fatti ``PS``: PS or a combination of them (required only when ``m0`` is ``None``). explicit : :obj:`bool`, optional Create a chained linear operator (``False``, preferred for large data) or a ``MatrixMult`` linear operator with dense matrix (``True``, preferred for small data) simultaneous : :obj:`bool`, optional Simultaneously invert entire data (``True``) or invert trace-by-trace (``False``) when using ``explicit`` operator (note that the entire data is always inverted when working with linear operator) epsI : :obj:`float` or :obj:`list`, optional Damping factor(s) for Tikhonov regularization term. If a list of :math:`n_{m}` elements is provided, the regularization term will have different strenght for each elastic property epsR : :obj:`float`, optional Damping factor for additional Laplacian regularization term dottest : :obj:`bool`, optional Apply dot-test returnres : :obj:`bool`, optional Return residuals epsRL1 : :obj:`float`, optional Damping factor for additional blockiness regularization term kind : :obj:`str`, optional Derivative kind (``forward`` or ``centered``). **kwargs_solver Arbitrary keyword arguments for :py:func:`scipy.linalg.lstsq` solver (if ``explicit=True`` and ``epsR=None``) or :py:func:`scipy.sparse.linalg.lsqr` solver (if ``explicit=False`` and/or ``epsR`` is not ``None``)) Returns ------- minv : :obj:`np.ndarray` Inverted model of size :math:`[n_{t0} \times n_{m} (\times n_x \times n_y)]` datar : :obj:`np.ndarray` Residual data (i.e., data - background data) of size :math:`[n_{t0} \times n_{\theta} (\times n_x \times n_y)]` Notes ----- The different choices of cost functions and solvers used in the seismic pre-stack inversion module follow the same convention of the seismic post-stack inversion module. Refer to :py:func:`pylops.avo.poststack.PoststackInversion` for more details. """ ncp = get_array_module(data) # find out dimensions if m0 is None and linearization is None: raise NotImplementedError('either m0 or linearization ' 'must be provided') elif m0 is None: if isinstance(linearization, str): nm = _linearizations[linearization] else: nm = _linearizations[linearization[0]] else: nm = m0.shape[1] data_shape = data.shape data_ndim = data.ndim n_lins = 1 multi = 0 if not isinstance(linearization, str): n_lins = data_shape[0] data_shape = data_shape[1:] data_ndim -= 1 multi = 1 if data_ndim == 2: dims = 1 nt0, ntheta = data_shape nspat = None nspatprod = nx = 1 elif data_ndim == 3: dims = 2 nt0, ntheta, nx = data_shape nspat = (nx, ) nspatprod = nx else: dims = 3 nt0, ntheta, nx, ny = data_shape nspat = (nx, ny) nspatprod = nx * ny data = data.reshape(nt0, ntheta, nspatprod) # check if background model and data have same shape if m0 is not None: if nt0 != m0.shape[0] or\ (dims >= 2 and nx != m0.shape[2]) or\ (dims == 3 and ny != m0.shape[3]): raise ValueError('data and m0 must have same time and space axes') # create operator if isinstance(linearization, str): # single operator PPop = PrestackLinearModelling(wav, theta, nt0=nt0, spatdims=nspat, linearization=linearization, explicit=explicit, kind=kind) else: # multiple operators if not isinstance(wav, (list, tuple)): wav = [ wav, ] * n_lins PPop = [ PrestackLinearModelling(w, theta, nt0=nt0, spatdims=nspat, linearization=lin, explicit=explicit) for w, lin in zip(wav, linearization) ] if explicit: PPop = MatrixMult(np.vstack([Op.A for Op in PPop]), dims=nspat, dtype=PPop[0].A.dtype) else: PPop = VStack(PPop) if dottest: Dottest(PPop, n_lins * nt0 * ntheta * nspatprod, nt0 * nm * nspatprod, raiseerror=True, verb=True, backend=get_module_name(ncp)) # swap axes for explicit operator if explicit: data = data.swapaxes(0 + multi, 1 + multi) if m0 is not None: m0 = m0.swapaxes(0, 1) # invert model if epsR is None: # create and remove background data from original data datar = data.flatten() if m0 is None else \ data.flatten() - PPop * m0.flatten() # inversion without spatial regularization if explicit: if epsI is None and not simultaneous: # solve unregularized equations indipendently trace-by-trace minv = \ get_lstsq(data)(PPop.A, datar.reshape(n_lins*nt0*ntheta, nspatprod).squeeze(), **kwargs_solver)[0] elif epsI is None and simultaneous: # solve unregularized equations simultaneously if ncp == np: minv = lsqr(PPop, datar, **kwargs_solver)[0] else: minv = cgls(PPop, datar, x0=ncp.zeros(int(PPop.shape[1]), PPop.dtype), **kwargs_solver)[0] elif epsI is not None: # create regularized normal equations PP = ncp.dot(PPop.A.T, PPop.A) + \ epsI * ncp.eye(nt0*nm, dtype=PPop.A.dtype) datarn = np.dot(PPop.A.T, datar.reshape(nt0 * ntheta, nspatprod)) if not simultaneous: # solve regularized normal eqs. trace-by-trace minv = get_lstsq(data)(PP, datarn, **kwargs_solver)[0] else: # solve regularized normal equations simultaneously PPop_reg = MatrixMult(PP, dims=nspatprod) if ncp == np: minv = lsqr(PPop_reg, datarn.ravel(), **kwargs_solver)[0] else: minv = cgls(PPop_reg, datarn.ravel(), x0=ncp.zeros(int(PPop_reg.shape[1]), PPop_reg.dtype), **kwargs_solver)[0] #else: # # create regularized normal eqs. and solve them simultaneously # PP = np.dot(PPop.A.T, PPop.A) + epsI * np.eye(nt0*nm) # datarn = PPop.A.T * datar.reshape(nt0*ntheta, nspatprod) # PPop_reg = MatrixMult(PP, dims=ntheta*nspatprod) # minv = lstsq(PPop_reg, datarn.flatten(), **kwargs_solver)[0] else: # solve unregularized normal equations simultaneously with lop if ncp == np: minv = lsqr(PPop, datar, **kwargs_solver)[0] else: minv = cgls(PPop, datar, x0=ncp.zeros(int(PPop.shape[1]), PPop.dtype), **kwargs_solver)[0] else: # Create Thicknov regularization if epsI is not None: if isinstance(epsI, (list, tuple)): if len(epsI) != nm: raise ValueError('epsI must be a scalar or a list of' 'size nm') RegI = Diagonal(np.array(epsI), dims=(nt0, nm, nspatprod), dir=1) else: RegI = epsI * Identity(nt0 * nm * nspatprod) if epsRL1 is None: # L2 inversion with spatial regularization if dims == 1: Regop = SecondDerivative(nt0 * nm, dtype=PPop.dtype, dims=(nt0, nm)) elif dims == 2: Regop = Laplacian((nt0, nm, nx), dirs=(0, 2), dtype=PPop.dtype) else: Regop = Laplacian((nt0, nm, nx, ny), dirs=(2, 3), dtype=PPop.dtype) if epsI is None: Regop = (Regop, ) epsR = (epsR, ) else: Regop = (Regop, RegI) epsR = (epsR, 1) minv = \ RegularizedInversion(PPop, Regop, data.ravel(), x0=m0.flatten() if m0 is not None else None, epsRs=epsR, returninfo=False, **kwargs_solver) else: # Blockiness-promoting inversion with spatial regularization if dims == 1: RegL1op = FirstDerivative(nt0 * nm, dtype=PPop.dtype) RegL2op = None elif dims == 2: RegL1op = FirstDerivative(nt0 * nx * nm, dims=(nt0, nm, nx), dir=0, dtype=PPop.dtype) RegL2op = SecondDerivative(nt0 * nx * nm, dims=(nt0, nm, nx), dir=2, dtype=PPop.dtype) else: RegL1op = FirstDerivative(nt0 * nx * ny * nm, dims=(nt0, nm, nx, ny), dir=0, dtype=PPop.dtype) RegL2op = Laplacian((nt0, nm, nx, ny), dirs=(2, 3), dtype=PPop.dtype) if dims == 1: if epsI is not None: RegL2op = (RegI, ) epsR = (1, ) else: if epsI is None: RegL2op = (RegL2op, ) epsR = (epsR, ) else: RegL2op = (RegL2op, RegI) epsR = (epsR, 1) epsRL1 = (epsRL1, ) if 'mu' in kwargs_solver.keys(): mu = kwargs_solver['mu'] kwargs_solver.pop('mu') else: mu = 1. if 'niter_outer' in kwargs_solver.keys(): niter_outer = kwargs_solver['niter_outer'] kwargs_solver.pop('niter_outer') else: niter_outer = 3 if 'niter_inner' in kwargs_solver.keys(): niter_inner = kwargs_solver['niter_inner'] kwargs_solver.pop('niter_inner') else: niter_inner = 5 minv = SplitBregman(PPop, (RegL1op, ), data.ravel(), RegsL2=RegL2op, epsRL1s=epsRL1, epsRL2s=epsR, mu=mu, niter_outer=niter_outer, niter_inner=niter_inner, x0=None if m0 is None else m0.flatten(), **kwargs_solver)[0] # compute residual if returnres: if epsR is None: datar -= PPop * minv.ravel() else: datar = data.ravel() - PPop * minv.ravel() # re-swap axes for explicit operator if explicit: if m0 is not None: m0 = m0.swapaxes(0, 1) # reshape inverted model and residual data if dims == 1: if explicit: minv = minv.reshape(nm, nt0).swapaxes(0, 1) if returnres: datar = datar.reshape(n_lins, ntheta, nt0).squeeze().swapaxes( 0 + multi, 1 + multi) else: minv = minv.reshape(nt0, nm) if returnres: datar = datar.reshape(n_lins, nt0, ntheta).squeeze() elif dims == 2: if explicit: minv = minv.reshape(nm, nt0, nx).swapaxes(0, 1) if returnres: datar = datar.reshape(n_lins, ntheta, nt0, nx).squeeze().swapaxes( 0 + multi, 1 + multi) else: minv = minv.reshape(nt0, nm, nx) if returnres: datar = datar.reshape(n_lins, nt0, ntheta, nx).squeeze() else: if explicit: minv = minv.reshape(nm, nt0, nx, ny).swapaxes(0, 1) if returnres: datar = datar.reshape(n_lins, ntheta, nt0, nx, ny).squeeze().swapaxes( 0 + multi, 1 + multi) else: minv = minv.reshape(nt0, nm, nx, ny) if returnres: datar = datar.reshape(n_lins, nt0, ntheta, nx, ny).squeeze() if m0 is not None and epsR is None: minv = minv + m0 if returnres: return minv, datar else: return minv
def PoststackInversion(data, wav, m0=None, explicit=False, simultaneous=False, epsI=None, epsR=None, dottest=False, **kwargs_solver): r"""Post-stack linearized seismic inversion. Invert post-stack seismic operator to retrieve an acoustic impedance profile from band-limited seismic post-stack data. Depending on the choice of input parameters, inversion can be trace-by-trace with explicit operator or global with either explicit or linear operator. Parameters ---------- data : :obj:`np.ndarray` Band-limited seismic post-stack data of size :math:`[n_{t0} \times n_x]` wav : :obj:`np.ndarray` Wavelet in time domain (must had odd number of elements and centered to zero) m0 : :obj:`np.ndarray`, optional Background model of size :math:`[n_{t0} \times n_x]` explicit : :obj:`bool`, optional Create a chained linear operator (``False``, preffered for large data) or a ``MatrixMult`` linear operator with dense matrix (``True``, preffered for small data) simultaneous : :obj:`bool`, optional Simultaneously invert entire data (``True``) or invert trace-by-trace (``False``) when using ``explicit`` operator (note that the entire data is always inverted when working with linear operator) epsI : :obj:`float`, optional Damping factor for Tikhonov regularization term epsR : :obj:`float`, optional Damping factor for additional Laplacian regularization term dottest : :obj:`bool`, optional Apply dot-test **kwargs_solver Arbitrary keyword arguments for :py:func:`scipy.linalg.lstsq` solver (if ``explicit=True`` and ``epsR=None``) or :py:func:`scipy.sparse.linalg.lsqr` solver (if ``explicit=False`` and/or ``epsR`` is not ``None``)) Returns ------- minv : :obj:`np.ndarray` Inverted model of size :math:`[n_{t0} \times n_x]` dr : :obj:`np.ndarray` Residual data (i.e., data - background data) of size :math:`[n_{t0} \times n_x]` Notes ----- The cost function and solver used in the seismic post-stack inversion module depends on the choice of ``explicit``, ``simultaneous``, ``epsI``, and ``epsR`` parameters: * ``explicit=True``, ``epsI=None`` and ``epsR=None``: the explicit solver :py:func:`scipy.linalg.lstsq` is used if ``simultaneous=False`` (or the iterative solver :py:func:`scipy.sparse.linalg.lsqr` is used if ``simultaneous=True``) * ``explicit=True`` with ``epsI`` and ``epsR=None``: the regularized normal equations :math:`\mathbf{W}^T\mathbf{d} = (\mathbf{W}^T \mathbf{W} + \epsilon_I^2 \mathbf{I}) \mathbf{AI}` are instead fed into the :py:func:`scipy.linalg.lstsq` solver if ``simultaneous=False`` (or the iterative solver :py:func:`scipy.sparse.linalg.lsqr` if ``simultaneous=True``) * ``explicit=False`` and ``epsR=None``: the iterative solver :py:func:`scipy.sparse.linalg.lsqr` is used * ``explicit=False`` with ``epsR``: the iterative solver :py:func:`pylops.optimization.leastsquares.RegularizedInversion` is used Note that the convergence of iterative solvers such as :py:func:`scipy.sparse.linalg.lsqr` can be very slow for this type of operator. It is suggested to take a two steps approach with first a trace-by-trace inversion using the explicit operator, followed by a regularized global inversion using the outcome of the previous inversion as initial guess. """ if data.ndim == 1: dims = 1 nt0 = data.size nspat = None nother = nx = 1 elif data.ndim == 2: dims = 2 nt0, nx = data.shape nspat = (nx, ) nother = nx else: dims = 3 nt0, nx, ny = data.shape nspat = (nx, ny) nother = nx * ny data = data.reshape(nt0, nother) # check if background model and data have same shape if m0 is not None and data.shape != m0.shape: raise ValueError('data and m0 must have same shape') # create operator PPop = PoststackLinearModelling(wav, nt0=nt0, ndims=nspat, explicit=explicit) if dottest: assert Dottest(PPop, nt0 * nother, nt0 * nother, verb=True) # create and remove background data from original data datar = data.flatten( ) if m0 is None else data.flatten() - PPop * m0.flatten() # invert model if epsR is None: # inversion without spatial regularization if explicit: if epsI is None and not simultaneous: # solve unregularized equations indipendently trace-by-trace minv = lstsq(PPop.A, datar.reshape(nt0, nother).squeeze(), **kwargs_solver)[0] elif epsI is None and simultaneous: # solve unregularized equations simultaneously minv = lsqr(PPop, datar, **kwargs_solver)[0] elif epsI is not None: # create regularized normal equations PP = np.dot(PPop.A.T, PPop.A) + epsI * np.eye(nt0) datar = np.dot(PPop.A.T, datar.reshape(nt0, nother)) if not simultaneous: # solve regularized normal equations indipendently trace-by-trace minv = lstsq(PP, datar.reshape(nt0, nother), **kwargs_solver)[0] else: # solve regularized normal equations simultaneously PPop_reg = MatrixMult(PP, dims=nother) minv = lsqr(PPop_reg, datar.flatten(), **kwargs_solver)[0] else: # create regularized normal equations and solve them simultaneously PP = np.dot(PPop.A, PPop.A) + epsI * np.eye(nt0) datar = PPop.A.T * datar.reshape(nt0, nother) PPop_reg = MatrixMult(PP, dims=nother) minv = lstsq(PPop_reg, datar.flatten(), **kwargs_solver)[0] else: # solve unregularized normal equations simultaneously with lop minv = lsqr(PPop, datar, **kwargs_solver)[0] else: # inversion with spatial regularization if dims == 1: Regop = SecondDerivative(nt0, dtype=PPop.dtype) elif dims == 2: Regop = Laplacian((nt0, nx), dtype=PPop.dtype) else: Regop = Laplacian((nt0, nx, ny), dirs=(1, 2), dtype=PPop.dtype) minv = RegularizedInversion(PPop, [Regop], data.flatten(), x0=m0.flatten(), epsRs=[epsR], returninfo=False, **kwargs_solver) if dims == 1: minv = minv.squeeze() datar = datar.squeeze() elif dims == 2: minv = minv.reshape(nt0, nx) datar = datar.reshape(nt0, nx) else: minv = minv.reshape(nt0, nx, ny) datar = datar.reshape(nt0, nx, ny) if m0 is not None and epsR is None: minv = minv + m0 return minv, datar
def SeismicInterpolation(data, nrec, iava, iava1=None, kind='fk', nffts=None, sampling=None, spataxis=None, spat1axis=None, taxis=None, paxis=None, p1axis=None, centeredh=True, nwins=None, nwin=None, nover=None, design=False, engine='numba', dottest=False, **kwargs_solver): r"""Seismic interpolation (or regularization). Interpolate seismic data from irregular to regular spatial grid. Depending on the size of the input ``data``, interpolation is either 2- or 3-dimensional. In case of 3-dimensional interpolation, data can be irregularly sampled in either one or both spatial directions. Parameters ---------- data : :obj:`np.ndarray` Irregularly sampled seismic data of size :math:`[n_{r_y} (\times n_{r_x} \times n_t)]` nrec : :obj:`int` or :obj:`tuple` Number of elements in the regularly sampled (reconstructed) spatial array, :math:`n_{R_y}` for 2-dimensional data and :math:`(n_{R_y}, n_{R_x})` for 3-dimensional data iava : :obj:`list` or :obj:`numpy.ndarray` Integer (or floating) indices of locations of available samples in first dimension of regularly sampled spatial grid of interpolated signal. The :class:`pylops.basicoperators.Restriction` operator is used in case of integer indices, while the :class:`pylops.signalprocessing.Iterp` operator is used in case of floating indices. iava1 : :obj:`list` or :obj:`numpy.ndarray`, optional Integer (or floating) indices of locations of available samples in second dimension of regularly sampled spatial grid of interpolated signal. Can be used only in case of 3-dimensional data. kind : :obj:`str`, optional Type of inversion: ``fk`` (default), ``spatial``, ``radon-linear``, ``chirpradon-linear``, ``radon-parabolic`` or , ``radon-hyperbolic`` and ``sliding`` nffts : :obj:`int` or :obj:`tuple`, optional nffts : :obj:`tuple`, optional Number of samples in Fourier Transform for each direction. Required if ``kind='fk'`` sampling : :obj:`tuple`, optional Sampling steps ``dy`` (, ``dx``) and ``dt``. Required if ``kind='fk'`` or ``kind='radon-linear'`` spataxis : :obj:`np.ndarray`, optional First spatial axis. Required for ``kind='radon-linear'``, ``kind='chirpradon-linear'``, ``kind='radon-parabolic'``, ``kind='radon-hyperbolic'``, can also be provided instead of ``sampling`` for ``kind='fk'`` spat1axis : :obj:`np.ndarray`, optional Second spatial axis. Required for ``kind='radon-linear'``, ``kind='chirpradon-linear'``, ``kind='radon-parabolic'``, ``kind='radon-hyperbolic'``, can also be provided instead of ``sampling`` for ``kind='fk'`` taxis : :obj:`np.ndarray`, optional Time axis. Required for ``kind='radon-linear'``, ``kind='chirpradon-linear'``, ``kind='radon-parabolic'``, ``kind='radon-hyperbolic'``, can also be provided instead of ``sampling`` for ``kind='fk'`` paxis : :obj:`np.ndarray`, optional First Radon axis. Required for ``kind='radon-linear'``, ``kind='chirpradon-linear'``, ``kind='radon-parabolic'``, ``kind='radon-hyperbolic'`` and ``kind='sliding'`` p1axis : :obj:`np.ndarray`, optional Second Radon axis. Required for ``kind='radon-linear'``, ``kind='chirpradon-linear'``, ``kind='radon-parabolic'``, ``kind='radon-hyperbolic'`` and ``kind='sliding'`` centeredh : :obj:`bool`, optional Assume centered spatial axis (``True``) or not (``False``). Required for ``kind='radon-linear'``, ``kind='radon-parabolic'`` and ``kind='radon-hyperbolic'`` nwins : :obj:`int` or :obj:`tuple`, optional Number of windows. Required for ``kind='sliding'`` nwin : :obj:`int` or :obj:`tuple`, optional Number of samples of window. Required for ``kind='sliding'`` nover : :obj:`int` or :obj:`tuple`, optional Number of samples of overlapping part of window. Required for ``kind='sliding'`` design : :obj:`bool`, optional Print number of sliding window (``True``) or not (``False``) when using ``kind='sliding'`` engine : :obj:`str`, optional Engine used for Radon computations (``numpy/numba`` for ``Radon2D`` and ``Radon3D`` or ``numpy/fftw`` for ``ChirpRadon2D`` and ``ChirpRadon3D`` or ) dottest : :obj:`bool`, optional Apply dot-test **kwargs_solver Arbitrary keyword arguments for :py:func:`pylops.optimization.leastsquares.RegularizedInversion` solver if ``kind='spatial'`` or :py:func:`pylops.optimization.sparsity.FISTA` solver otherwise Returns ------- recdata : :obj:`np.ndarray` Reconstructed data of size :math:`[n_{R_y} (\times n_{R_x} \times n_t)]` recprec : :obj:`np.ndarray` Reconstructed data in the sparse or preconditioned domain in case of ``kind='fk'``, ``kind='radon-linear'``, ``kind='radon-parabolic'``, ``kind='radon-hyperbolic'`` and ``kind='sliding'`` cost : :obj:`np.ndarray` Cost function norm Raises ------ KeyError If ``kind`` is neither ``spatial``, ``fl``, ``radon-linear``, ``radon-parabolic``, ``radon-hyperbolic`` nor ``sliding`` Notes ----- The problem of seismic data interpolation (or regularization) can be formally written as .. math:: \mathbf{y} = \mathbf{R} \mathbf{x} where a restriction or interpolation operator is applied along the spatial direction(s). Here :math:`\mathbf{y} = [\mathbf{y}_{R1}^T, \mathbf{y}_{R2}^T,..., \mathbf{y}_{RN^T}]^T` where each vector :math:`\mathbf{y}_{Ri}` contains all time samples recorded in the seismic data at the specific receiver :math:`R_i`. Similarly, :math:`\mathbf{x} = [\mathbf{x}_{r1}^T, \mathbf{x}_{r2}^T,..., \mathbf{x}_{rM}^T]`, contains all traces at the regularly and finely sampled receiver locations :math:`r_i`. Several alternative approaches can be taken to solve such a problem. They mostly differ in the choice of the regularization (or preconditining) used to mitigate the ill-posedness of the problem: * ``spatial``: least-squares inversion in the original time-space domain with an additional spatial smoothing regularization term, corresponding to the cost function :math:`J = ||\mathbf{y} - \mathbf{R} \mathbf{x}||_2 + \epsilon_\nabla \nabla ||\mathbf{x}||_2` where :math:`\nabla` is a second order space derivative implemented via :class:`pylops.basicoperators.SecondDerivative` in 2-dimensional case and :class:`pylops.basicoperators.Laplacian` in 3-dimensional case * ``fk``: L1 inversion in frequency-wavenumber preconditioned domain corresponding to the cost function :math:`J = ||\mathbf{y} - \mathbf{R} \mathbf{F} \mathbf{x}||_2` where :math:`\mathbf{F}` is frequency-wavenumber transform implemented via :class:`pylops.signalprocessing.FFT2D` in 2-dimensional case and :class:`pylops.signalprocessing.FFTND` in 3-dimensional case * ``radon-linear``: L1 inversion in linear Radon preconditioned domain using the same cost function as ``fk`` but with :math:`\mathbf{F}` being a Radon transform implemented via :class:`pylops.signalprocessing.Radon2D` in 2-dimensional case and :class:`pylops.signalprocessing.Radon3D` in 3-dimensional case * ``radon-parabolic``: L1 inversion in parabolic Radon preconditioned domain * ``radon-hyperbolic``: L1 inversion in hyperbolic Radon preconditioned domain * ``sliding``: L1 inversion in sliding-linear Radon preconditioned domain using the same cost function as ``fk`` but with :math:`\mathbf{F}` being a sliding Radon transform implemented via :class:`pylops.signalprocessing.Sliding2D` in 2-dimensional case and :class:`pylops.signalprocessing.Sliding3D` in 3-dimensional case """ ncp = get_array_module(data) dtype = data.dtype ndims = data.ndim if ndims == 1 or ndims > 3: raise ValueError('data must have 2 or 3 dimensions') if ndims == 2: dimsd = data.shape dims = (nrec, dimsd[1]) else: dimsd = data.shape dims = (nrec[0], nrec[1], dimsd[2]) # sampling if taxis is not None: dt = taxis[1] - taxis[0] if spataxis is not None: dspat = np.abs(spataxis[1] - spataxis[0]) if spat1axis is not None: dspat1 = np.abs(spat1axis[1] - spat1axis[0]) # create restriction/interpolation operator if iava.dtype == float: Rop = Interp(np.prod(dims), iava, dims=dims, dir=0, kind='linear', dtype=dtype) if ndims == 3 and iava1 is not None: dims1 = (len(iava), nrec[1], dimsd[2]) Rop1 = Interp(np.prod(dims1), iava1, dims=dims1, dir=1, kind='linear', dtype=dtype) Rop = Rop1 * Rop else: Rop = Restriction(np.prod(dims), iava, dims=dims, dir=0, dtype=dtype) if ndims == 3 and iava1 is not None: dims1 = (len(iava), nrec[1], dimsd[2]) Rop1 = Restriction(np.prod(dims1), iava1, dims=dims1, dir=1, dtype=dtype) Rop = Rop1 * Rop # create other operators for inversion if kind == 'spatial': prec = False dotcflag = 0 if ndims == 3 and iava1 is not None: Regop = Laplacian(dims=dims, dirs=(0, 1), dtype=dtype) else: Regop = SecondDerivative(np.prod(dims), dims=(dims), dir=0, dtype=dtype) SIop = Rop elif kind == 'fk': prec = True dimsp = nffts dotcflag = 1 if ndims == 3: if sampling is None: if spataxis is None or spat1axis is None or taxis is None: raise ValueError('Provide either sampling or spataxis, ' 'spat1axis and taxis for kind=%s' % kind) else: sampling = (np.abs(spataxis[1] - spataxis[1]), np.abs(spat1axis[1] - spat1axis[1]), np.abs(taxis[1] - taxis[1])) Pop = FFTND(dims=dims, nffts=nffts, sampling=sampling) Pop = Pop.H else: if sampling is None: if spataxis is None or taxis is None: raise ValueError('Provide either sampling or spataxis, ' 'and taxis for kind=%s' % kind) else: sampling = (np.abs(spataxis[1] - spataxis[1]), np.abs(taxis[1] - taxis[1])) Pop = FFT2D(dims=dims, nffts=nffts, sampling=sampling) Pop = Pop.H SIop = Rop * Pop elif 'chirpradon' in kind: prec = True dotcflag = 0 if ndims == 3: Pop = ChirpRadon3D( taxis, spataxis, spat1axis, (np.max(paxis) * dspat / dt, np.max(p1axis) * dspat1 / dt)).H dimsp = (spataxis.size, spat1axis.size, taxis.size) else: Pop = ChirpRadon2D(taxis, spataxis, np.max(paxis) * dspat / dt).H dimsp = (spataxis.size, taxis.size) SIop = Rop * Pop elif 'radon' in kind: prec = True dotcflag = 0 kindradon = kind.split('-')[-1] if ndims == 3: Pop = Radon3D(taxis, spataxis, spat1axis, paxis, p1axis, centeredh=centeredh, kind=kindradon, engine=engine) dimsp = (paxis.size, p1axis.size, taxis.size) else: Pop = Radon2D(taxis, spataxis, paxis, centeredh=centeredh, kind=kindradon, engine=engine) dimsp = (paxis.size, taxis.size) SIop = Rop * Pop elif kind == 'sliding': prec = True dotcflag = 0 if ndims == 3: nspat, nspat1 = spataxis.size, spat1axis.size spataxis_local = np.linspace(-dspat * nwin[0] // 2, dspat * nwin[0] // 2, nwin[0]) spat1axis_local = np.linspace(-dspat1 * nwin[1] // 2, dspat1 * nwin[1] // 2, nwin[1]) dimsslid = (nspat, nspat1, taxis.size) if ncp == np: npaxis, np1axis = paxis.size, p1axis.size Op = Radon3D(taxis, spataxis_local, spat1axis_local, paxis, p1axis, centeredh=True, kind='linear', engine=engine) else: npaxis, np1axis = nwin[0], nwin[1] Op = ChirpRadon3D(taxis, spataxis_local, spat1axis_local, (np.max(paxis) * dspat / dt, np.max(p1axis) * dspat1 / dt)).H dimsp = (nwins[0] * npaxis, nwins[1] * np1axis, dimsslid[2]) Pop = Sliding3D(Op, dimsp, dimsslid, nwin, nover, (npaxis, np1axis), tapertype='cosine') # to be able to reshape correctly the preconditioned model dimsp = (nwins[0], nwins[1], npaxis, np1axis, dimsslid[2]) else: nspat = spataxis.size spataxis_local = np.linspace(-dspat * nwin // 2, dspat * nwin // 2, nwin) dimsslid = (nspat, taxis.size) if ncp == np: npaxis = paxis.size Op = Radon2D(taxis, spataxis_local, paxis, centeredh=True, kind='linear', engine=engine) else: npaxis = nwin Op = ChirpRadon2D(taxis, spataxis_local, np.max(paxis) * dspat / dt).H dimsp = (nwins * npaxis, dimsslid[1]) Pop = Sliding2D(Op, dimsp, dimsslid, nwin, nover, tapertype='cosine', design=design) SIop = Rop * Pop else: raise KeyError('kind must be spatial, fk, radon-linear, ' 'radon-parabolic, radon-hyperbolic or sliding') # dot-test if dottest: Dottest(SIop, np.prod(dimsd), np.prod(dimsp) if prec else np.prod(dims), complexflag=dotcflag, raiseerror=True, verb=True) # inversion if kind == 'spatial': recdata = \ RegularizedInversion(SIop, [Regop], data.flatten(), **kwargs_solver) if isinstance(recdata, tuple): recdata = recdata[0] recdata = recdata.reshape(dims) recprec = None cost = None else: recprec = FISTA(SIop, data.flatten(), **kwargs_solver) if len(recprec) == 3: cost = recprec[2] else: cost = None recprec = recprec[0] recdata = np.real(Pop * recprec) recprec = recprec.reshape(dimsp) recdata = recdata.reshape(dims) return recdata, recprec, cost
def PoststackInversion(data, wav, m0=None, explicit=False, simultaneous=False, epsI=None, epsR=None, dottest=False, epsRL1=None, **kwargs_solver): r"""Post-stack linearized seismic inversion. Invert post-stack seismic operator to retrieve an elastic parameter of choice from band-limited seismic post-stack data. Depending on the choice of input parameters, inversion can be trace-by-trace with explicit operator or global with either explicit or linear operator. Parameters ---------- data : :obj:`np.ndarray` Band-limited seismic post-stack data of size :math:`[n_{t0} (\times n_x \times n_y)]` wav : :obj:`np.ndarray` Wavelet in time domain (must have odd number of elements and centered to zero). If 1d, assume stationary wavelet for the entire time axis. If 2d of size :math:`[n_{t0} \times n_h]` use as non-stationary wavelet m0 : :obj:`np.ndarray`, optional Background model of size :math:`[n_{t0} (\times n_x \times n_y)]` explicit : :obj:`bool`, optional Create a chained linear operator (``False``, preferred for large data) or a ``MatrixMult`` linear operator with dense matrix (``True``, preferred for small data) simultaneous : :obj:`bool`, optional Simultaneously invert entire data (``True``) or invert trace-by-trace (``False``) when using ``explicit`` operator (note that the entire data is always inverted when working with linear operator) epsI : :obj:`float`, optional Damping factor for Tikhonov regularization term epsR : :obj:`float`, optional Damping factor for additional Laplacian regularization term dottest : :obj:`bool`, optional Apply dot-test epsRL1 : :obj:`float`, optional Damping factor for additional blockiness regularization term **kwargs_solver Arbitrary keyword arguments for :py:func:`scipy.linalg.lstsq` solver (if ``explicit=True`` and ``epsR=None``) or :py:func:`scipy.sparse.linalg.lsqr` solver (if ``explicit=False`` and/or ``epsR`` is not ``None``) Returns ------- minv : :obj:`np.ndarray` Inverted model of size :math:`[n_{t0} (\times n_x \times n_y)]` datar : :obj:`np.ndarray` Residual data (i.e., data - background data) of size :math:`[n_{t0} (\times n_x \times n_y)]` Notes ----- The cost function and solver used in the seismic post-stack inversion module depends on the choice of ``explicit``, ``simultaneous``, ``epsI``, and ``epsR`` parameters: * ``explicit=True``, ``epsI=None`` and ``epsR=None``: the explicit solver :py:func:`scipy.linalg.lstsq` is used if ``simultaneous=False`` (or the iterative solver :py:func:`scipy.sparse.linalg.lsqr` is used if ``simultaneous=True``) * ``explicit=True`` with ``epsI`` and ``epsR=None``: the regularized normal equations :math:`\mathbf{W}^T\mathbf{d} = (\mathbf{W}^T \mathbf{W} + \epsilon_I^2 \mathbf{I}) \mathbf{AI}` are instead fed into the :py:func:`scipy.linalg.lstsq` solver if ``simultaneous=False`` (or the iterative solver :py:func:`scipy.sparse.linalg.lsqr` if ``simultaneous=True``) * ``explicit=False`` and ``epsR=None``: the iterative solver :py:func:`scipy.sparse.linalg.lsqr` is used * ``explicit=False`` with ``epsR`` and ``epsRL1=None``: the iterative solver :py:func:`pylops.optimization.leastsquares.RegularizedInversion` is used to solve the spatially regularized problem. * ``explicit=False`` with ``epsR`` and ``epsRL1``: the iterative solver :py:func:`pylops.optimization.sparsity.SplitBregman` is used to solve the blockiness-promoting (in vertical direction) and spatially regularized (in additional horizontal directions) problem. Note that the convergence of iterative solvers such as :py:func:`scipy.sparse.linalg.lsqr` can be very slow for this type of operator. It is suggested to take a two steps approach with first a trace-by-trace inversion using the explicit operator, followed by a regularized global inversion using the outcome of the previous inversion as initial guess. """ ncp = get_array_module(wav) # check if background model and data have same shape if m0 is not None and data.shape != m0.shape: raise ValueError('data and m0 must have same shape') # find out dimensions if data.ndim == 1: dims = 1 nt0 = data.size nspat = None nspatprod = nx = 1 elif data.ndim == 2: dims = 2 nt0, nx = data.shape nspat = (nx, ) nspatprod = nx else: dims = 3 nt0, nx, ny = data.shape nspat = (nx, ny) nspatprod = nx * ny data = data.reshape(nt0, nspatprod) # create operator PPop = PoststackLinearModelling(wav, nt0=nt0, spatdims=nspat, explicit=explicit) if dottest: Dottest(PPop, nt0 * nspatprod, nt0 * nspatprod, raiseerror=True, backend=get_module_name(ncp), verb=True) # create and remove background data from original data datar = data.flatten() if m0 is None else \ data.flatten() - PPop * m0.flatten() # invert model if epsR is None: # inversion without spatial regularization if explicit: if epsI is None and not simultaneous: # solve unregularized equations indipendently trace-by-trace minv = \ get_lstsq(data)(PPop.A, datar.reshape(nt0, nspatprod).squeeze(), **kwargs_solver)[0] elif epsI is None and simultaneous: # solve unregularized equations simultaneously if ncp == np: minv = lsqr(PPop, datar, **kwargs_solver)[0] else: minv = \ cgls(PPop, datar, x0=ncp.zeros(int(PPop.shape[1]), PPop.dtype), **kwargs_solver)[0] elif epsI is not None: # create regularized normal equations PP = ncp.dot(PPop.A.T, PPop.A) + \ epsI * ncp.eye(nt0, dtype=PPop.A.dtype) datarn = ncp.dot(PPop.A.T, datar.reshape(nt0, nspatprod)) if not simultaneous: # solve regularized normal eqs. trace-by-trace minv = get_lstsq(data)(PP, datarn, **kwargs_solver)[0] else: # solve regularized normal equations simultaneously PPop_reg = MatrixMult(PP, dims=nspatprod) if ncp == np: minv = lsqr(PPop_reg, datar.ravel(), **kwargs_solver)[0] else: minv = cgls(PPop_reg, datar.ravel(), x0=ncp.zeros(int(PPop_reg.shape[1]), PPop_reg.dtype), **kwargs_solver)[0] else: # create regularized normal eqs. and solve them simultaneously PP = ncp.dot(PPop.A.T, PPop.A) + \ epsI * ncp.eye(nt0, dtype=PPop.A.dtype) datarn = PPop.A.T * datar.reshape(nt0, nspatprod) PPop_reg = MatrixMult(PP, dims=nspatprod) minv = \ get_lstsq(data)(PPop_reg.A, datarn.flatten(), **kwargs_solver)[0] else: # solve unregularized normal equations simultaneously with lop if ncp == np: minv = lsqr(PPop, datar, **kwargs_solver)[0] else: minv = \ cgls(PPop, datar, x0=ncp.zeros(int(PPop.shape[1]), PPop.dtype), **kwargs_solver)[0] else: if epsRL1 is None: # L2 inversion with spatial regularization if dims == 1: Regop = SecondDerivative(nt0, dtype=PPop.dtype) elif dims == 2: Regop = Laplacian((nt0, nx), dtype=PPop.dtype) else: Regop = Laplacian((nt0, nx, ny), dirs=(1, 2), dtype=PPop.dtype) minv = RegularizedInversion( PPop, [Regop], data.flatten(), x0=None if m0 is None else m0.flatten(), epsRs=[epsR], returninfo=False, **kwargs_solver) else: # Blockiness-promoting inversion with spatial regularization if dims == 1: RegL1op = FirstDerivative(nt0, kind='forward', dtype=PPop.dtype) RegL2op = None elif dims == 2: RegL1op = FirstDerivative(nt0 * nx, dims=(nt0, nx), dir=0, kind='forward', dtype=PPop.dtype) RegL2op = SecondDerivative(nt0 * nx, dims=(nt0, nx), dir=1, dtype=PPop.dtype) else: RegL1op = FirstDerivative(nt0 * nx * ny, dims=(nt0, nx, ny), dir=0, kind='forward', dtype=PPop.dtype) RegL2op = Laplacian((nt0, nx, ny), dirs=(1, 2), dtype=PPop.dtype) if 'mu' in kwargs_solver.keys(): mu = kwargs_solver['mu'] kwargs_solver.pop('mu') else: mu = 1. if 'niter_outer' in kwargs_solver.keys(): niter_outer = kwargs_solver['niter_outer'] kwargs_solver.pop('niter_outer') else: niter_outer = 3 if 'niter_inner' in kwargs_solver.keys(): niter_inner = kwargs_solver['niter_inner'] kwargs_solver.pop('niter_inner') else: niter_inner = 5 if not isinstance(epsRL1, (list, tuple)): epsRL1 = list([epsRL1]) if not isinstance(epsR, (list, tuple)): epsR = list([epsR]) minv = SplitBregman(PPop, [RegL1op], data.ravel(), RegsL2=[RegL2op], epsRL1s=epsRL1, epsRL2s=epsR, mu=mu, niter_outer=niter_outer, niter_inner=niter_inner, x0=None if m0 is None else m0.flatten(), **kwargs_solver)[0] # compute residual if epsR is None: datar -= PPop * minv.ravel() else: datar = data.ravel() - PPop * minv.ravel() # reshape inverted model and residual data if dims == 1: minv = minv.squeeze() datar = datar.squeeze() elif dims == 2: minv = minv.reshape(nt0, nx) datar = datar.reshape(nt0, nx) else: minv = minv.reshape(nt0, nx, ny) datar = datar.reshape(nt0, nx, ny) if m0 is not None and epsR is None: minv = minv + m0 return minv, datar