API Reference

Low Level Custom Operators

FwiFlow.poisson_opFunction
poisson_op(c::Union{PyObject, Float64}, g::Union{PyObject, Float64}, 
    h::Union{PyObject, Float64}, ρ::Union{PyObject, Float64}, index::Union{Integer, PyObject}=0)

Solves the Poisson equation ($\mathbf{x}=[z\quad x]^T$)

\[\begin{aligned} -\nabla\cdot\left(c(\mathbf{x}) \nabla \left(u(\mathbf{x}) -\rho \begin{bmatrix}z \\ 0\end{bmatrix} \right)\right) &= g(\mathbf{x}) & \mathbf{x}\in \Omega\\ \frac{\partial u(x)}{\partial n} &= 0 & \mathbf{x}\in \Omega\\ \end{aligned}\]

Here $\Omega=[0,n_zh]\times [0, n_xh]$. The equation is solved using finite difference method, where the step size in each direction is $h$. Mathematically, the solution to the PDE is determined up to a constant. Numerically, we discretize the equation with the scheme

\[(A+E_{11})\mathbf{u} = \mathbf{f}\]

where $A$ is the finite difference coefficient matrix,

\[(E_{11})_{ij} = \left\{ \begin{matrix}1 & i=j=1 \\ 0 & \text{ otherwise }\end{matrix}\right.\]
  • index : Int32, when index=1, SparseLU is used to solve the linear system; otherwise the function invokes algebraic multigrid method from amgcl.
source
FwiFlow.laplacian_opFunction
laplacian_op(coef::Union{PyObject, Array{Float64}}, f::Union{PyObject, Array{Float64}}, 
        h::Union{PyObject, Float64}, ρ::Union{PyObject, Float64})

Computes the Laplacian of function $f(\mathbf{x})$; here ($\mathbf{x}=[z\quad x]^T$)

\[-\nabla\cdot\left(c(\mathbf{x}) \nabla \left(u(\mathbf{x}) -\rho \begin{bmatrix}z \\ 0\end{bmatrix} \right)\right)\]
source
FwiFlow.fwi_obs_opFunction
fwi_obs_op(lambda::Union{PyObject, Array{Float64}},mu::Union{PyObject, Array{Float64}},
    den::Union{PyObject, Array{Float64}},stf::Union{PyObject, Array{Float64}},
    gpu_id::Union{PyObject, Integer},shot_ids::Union{PyObject, Array{T}},para_fname::String) where T<:Integer

Generates the observation data and store them as files which will be used by fwi_op For the meaning of parameters, see fwi_op.

source
FwiFlow.fwi_opFunction
fwi_op(lambda::Union{PyObject, Array{Float64}},mu::Union{PyObject, Array{Float64}},
    den::Union{PyObject, Array{Float64}},stf::Union{PyObject, Array{Float64}},
    gpu_id::Union{PyObject, Integer},shot_ids::Union{PyObject, Array{T}},para_fname::String) where T<:Integer

Computes the FWI loss function.

  • lambda : Lame's first parameter (unit: MPa)
  • mu : Lame's second parameter (shear modulus, unit: MPa)
  • den : Density
  • stf : Source time functions
  • gpu_id : The ID of GPU to run this FWI operator
  • shot_ids : The source function IDs (determining the location of sources)
  • para_fname : Parameter file location
source
FwiFlow.sat_opFunction
sat_op(s0::Union{PyObject, Array{Float64}},pt::Union{PyObject, Array{Float64}},
    permi::Union{PyObject, Array{Float64}},poro::Union{PyObject, Array{Float64}},
    qw::Union{PyObject, Array{Float64}},qo::Union{PyObject, Array{Float64}},
    muw::Union{PyObject, Float64},muo::Union{PyObject, Float64},
    sref::Union{PyObject, Array{Float64}},dt::Union{PyObject, Float64},h::Union{PyObject, Float64})

Solves the following discretized equation

\[\phi (S_2^{n + 1} - S_2^n) - \nabla \cdot \left( {{m_2}(S_2^{n + 1})K\nabla \Psi _2^n} \right) \Delta t= \left(q_2^n + q_1^n \frac{m_2(S^{n+1}_2)}{m_1(S^{n+1}_2)}\right) \Delta t\]

where

\[m_2(s) = \frac{s^2}{\mu_w}\qquad m_1(s) = \frac{(1-s)^2}{\mu_o}\]

This is a nonlinear equation and is solved with the Newton-Raphson method.

  • s0 : $n_z\times n_x$, saturation of fluid, i.e., $S_2^n$
  • pt : $n_z\times n_x$, potential of fluid, i.e., $\Psi_2^n$
  • permi : $n_z\times n_x$, permeability, i.e., $K$
  • poro : $n_z\times n_x$, porosity, i.e., $\phi$
  • qw : $n_z\times n_x$, injection or production rate of fluid 1, $q_2^n$
  • qo : $n_z\times n_x$, injection or production rate of fluid 2, $q_1^n$
  • muw : viscosity of fluid 1, i.e., $\mu_w$
  • muo : viscosity of fluid 2, i.e., $\mu_o$
  • sref : $n_z\times n_x$, initial guess for $S_2^{n+1}$
  • dt : Time step size
  • h : Spatial step size
source
FwiFlow.sat_op2Function
sat_op2(s0::Union{PyObject, Array{Float64}},
dporodt::Union{PyObject, Array{Float64}},
pt::Union{PyObject, Array{Float64}},
permi::Union{PyObject, Array{Float64}},
poro::Union{PyObject, Array{Float64}},
qw::Union{PyObject, Array{Float64}},
qo::Union{PyObject, Array{Float64}},
muw::Union{PyObject, Float64},
muo::Union{PyObject, Float64},
sref::Union{PyObject, Array{Float64}},
dt::Union{PyObject, Float64},
h::Union{PyObject, Float64})

Solves the following discretized equation

\[\phi (S_2^{n + 1} - S_2^n) + \Delta t \dot \phi S_2^{n+1} - \nabla \cdot \left( {{m_2}(S_2^{n + 1})K\nabla \Psi _2^n} \right) \Delta t= \left(q_2^n + q_1^n \frac{m_2(S^{n+1}_2)}{m_1(S^{n+1}_2)}\right) \Delta t\]

where

\[m_2(s) = \frac{s^2}{\mu_w}\qquad m_1(s) = \frac{(1-s)^2}{\mu_o}\]

This is a nonlinear equation and is solved with the Newton-Raphson method.

  • s0 : $n_z\times n_x$, saturation of fluid, i.e., $S_2^n$
  • dporodt : $n_z\times n_x$, rate of porosity, $\dot \phi$
  • pt : $n_z\times n_x$, potential of fluid, i.e., $\Psi_2^n$
  • permi : $n_z\times n_x$, permeability, i.e., $K$
  • poro : $n_z\times n_x$, porosity, i.e., $\phi$
  • qw : $n_z\times n_x$, injection or production rate of fluid 1, $q_2^n$
  • qo : $n_z\times n_x$, injection or production rate of fluid 2, $q_1^n$
  • muw : viscosity of fluid 1, i.e., $\mu_w$
  • muo : viscosity of fluid 2, i.e., $\mu_o$
  • sref : $n_z\times n_x$, initial guess for $S_2^{n+1}$
  • dt : Time step size
  • h : Spatial step size
source
FwiFlow.upwlap_opFunction
upwlap_op(perm::Union{PyObject, Array{Float64}},
    mobi::Union{PyObject, Array{Float64}},
    func::Union{PyObject, Array{Float64}},
    h::Union{PyObject, Float64},
    rhograv::Union{PyObject, Float64})

Computes the Laplacian of function $f(\mathbf{x})$; here $\mathbf{x}=[z\quad x]^T$.

\[\nabla\cdot\left(m(\mathbf{x})K(\mathbf{x}) \nabla \left(f(\mathbf{x}) -\rho \begin{bmatrix}z \\ 0\end{bmatrix} \right)\right)\]

The permeability on the computational grid is computed with Harmonic mean; the mobility is computed with upwind scheme.

  • perm : $n_z\times n_x$, permeability of fluid, i.e., $K$
  • mobi : $n_z\times n_x$, mobility of fluid, i.e., $m$
  • func : $n_z\times n_x$, potential of fluid, i.e., $f$
  • h : Float64, spatial step size
  • rhograv : Float64, i.e., $\rho$
source
FwiFlow.upwps_opFunction
upwps_op(perm::Union{PyObject, Array{Float64}},mobi::Union{PyObject, Array{Float64}},
src,funcref,h::Union{PyObject, Float64},rhograv::Union{PyObject, Float64},index::Union{PyObject, Integer})

Solves the Poisson equation

\[-\nabla\cdot\left(m(\mathbf{x})K(\mathbf{x}) \nabla \left(u(\mathbf{x}) -\rho \begin{bmatrix}z \\ 0\end{bmatrix} \right)\right) = g(\mathbf{x}) \]

See upwps_op for detailed description.

  • perm : $n_z\times n_x$, permeability of fluid, i.e., $K$
  • mobi : $n_z\times n_x$, mobility of fluid, i.e., $m$
  • src : $n_z\times n_x$, source function, i.e., $g(\mathbf{x})$
  • funcref : $n_z\times n_x$, currently it is not not used
  • h : Float64, spatial step size
  • rhograv : Float64, i.e., $\rho$
  • index : Int32, when index=1, SparseLU is used to solve the linear system; otherwise the function invokes algebraic multigrid method from amgcl.
source
FwiFlow.time_fractional_opFunction
time_fractional_op(α::Union{Float64, PyObject}, f_fun::Function, 
    T::Float64, u0::Union{Float64, Array{Float64}, PyObject}, 
    NT::Int64, θ=missing)

Returns a $(NT+1)\times \texttt{size}(u_0)$ solution array. The function solves the following time fractional differential equation with explicit scheme

\[{}_0^CD_t^\alpha u(t) = f(t, u, \theta)\]
source
FwiFlow.time_fractional_t_opFunction
time_fractional_t_op(i::Union{Integer, PyObject}, ta::PyObject, 
    α::Union{Float64, PyObject}, Δt::Union{Float64, PyObject})

Returns the coefficients for the time fractional derivative

\[{}_0^CD_t^\alpha f(t) = \frac{1}{\Gamma(1-\alpha)}\int_0^t \frac{f'(\tau)d\tau}{(t-\tau)^\alpha}\]

The discretization scheme used here is

\[\begin{aligned} & {}_0^CD_\tau^\alpha u(\tau_n) \\ = &\frac{\Delta \tau^{-\alpha}}{\Gamma(2-\alpha)}\left[G_0 u_n - \sum_{k=1}^{n-1}(G_{n-k-1}-G_{n-k})u_k + G_n u_0 \right] + \mathcal{O}(\Delta \tau^{2-\alpha}) \end{aligned}\]

Here

\[G_m = (m+1)^{1-\alpha} - m^{1-\alpha}, \quad m\geq 0, \quad 0<\alpha<1\]

The function returns

  • c : $\frac{\Delta \tau^{-\alpha}}{\Gamma(2-\alpha)}$
  • cum : $- \sum_{k=1}^{n-1}(G_{n-k-1}-G_{n-k})u_k + G_n u_0$
source

FWI

FwiFlow.FWIType
FWI(nz::Int64, nx::Int64, dz::Float64, dx::Float64, nSteps::Int64, dt::Float64;
        ind_src_z::Array{Int64, 1}, ind_src_x::Array{Int64, 1}, ind_rec_z::Array{Int64, 1}, ind_rec_x::Array{Int64, 1},
        kwargs...)

Creates a FWI structure that holds geometry data and settings. Here nz = n + 1, nx = m + 1, dz and dx are mesh sizes. nSteps is the number of total iterations. dt is the time step. ind_src_z and ind_src_x are arrays of source locations. ind_rec_z and ind_rec_x are arrays of receiver locations.

source

Utility Functions