API Reference
Low Level Custom Operators
FwiFlow.poisson_op
— Functionpoisson_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$)
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
where $A$ is the finite difference coefficient matrix,
index
:Int32
, whenindex=1
,SparseLU
is used to solve the linear system; otherwise the function invokes algebraic multigrid method fromamgcl
.
FwiFlow.laplacian_op
— Functionlaplacian_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$)
FwiFlow.fwi_obs_op
— Functionfwi_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
.
FwiFlow.fwi_op
— Functionfwi_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
: Densitystf
: Source time functionsgpu_id
: The ID of GPU to run this FWI operatorshot_ids
: The source function IDs (determining the location of sources)para_fname
: Parameter file location
FwiFlow.sat_op
— Functionsat_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
where
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 sizeh
: Spatial step size
FwiFlow.sat_op2
— Functionsat_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
where
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 sizeh
: Spatial step size
FwiFlow.upwlap_op
— Functionupwlap_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$.
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 sizerhograv
:Float64
, i.e., $\rho$
FwiFlow.upwps_op
— Functionupwps_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
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 usedh
:Float64
, spatial step sizerhograv
:Float64
, i.e., $\rho$index
:Int32
, whenindex=1
,SparseLU
is used to solve the linear system; otherwise the function invokes algebraic multigrid method fromamgcl
.
FwiFlow.time_fractional_op
— Functiontime_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
FwiFlow.time_fractional_t_op
— Functiontime_fractional_t_op(i::Union{Integer, PyObject}, ta::PyObject,
α::Union{Float64, PyObject}, Δt::Union{Float64, PyObject})
Returns the coefficients for the time fractional derivative
The discretization scheme used here is
Here
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$
FWI
FwiFlow.FWI
— TypeFWI(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.