Reservoir simulation in Python. It currently supports 2D rectangular grids and isotropic permeability.
-
ressim.py
: main module containing classes to define the grid, and to model and solve the pressure and saturation equations. -
utils.py
: module containing useful function definitions, e.g. linear and quadratic mobility functions, fractional flow function, etc.
See main_*.py
for examples.
A module for reservoir simulation in Python.
Grid(self, nx, ny, lx=1.0, ly=1.0)
A class to define a rectangular grid. Given nx, ny, lx, ly
, it maintains an update of attributes vol, dx, dy, ncell, shape
.
Attributes
-
nx, ny
:int
Grid resolution. -
lx, ly
:float
Grid physical dimensions. Default islx=1.0, ly=1.0
, i.e. unit square. -
vol
:float
Cell volume. -
dx, dy
:float
Cell dimensions. -
ncell
:int
Number of cells. -
shape
:int
Grid shape, i.e.(ny, nx)
.
PressureEquation(self, grid=None, q=None, k=None, diri=None, lamb_fn=None, s=None)
A class to model and solve the pressure equation,
with no-flow boundary conditions (closed reservoir).
Attributes
-
grid
:Grid
Simulation grid. -
q
:ndarray, shape (ny, nx) | (ny*nx,)
Integrated source term. -
k
:ndarray, shape (ny, nx)
Permeability -
diri
:list
of(int, float)
tuples Dirichlet boundary conditions, e.g.[(i1, val1), (i2, val2), ...]
means pressure valuesval1
at celli1
,val2
at celli2
, etc. Defaults to[(ny*nx/2, 0.0)]
, i.e. zero pressure at center of the grid. -
lamb_fn
:callable
Total mobility functionlamb_fn(s)
-
s
:ndarray, (ny, nx) | (ny*nx,)
Water saturation -
p
:ndarray, (ny, nx)
Pressure -
v
:dict
ofndarray
'x'
:ndarray, shape (ny, nx+1)
Flux in x-direction'y'
:ndarray, shape (ny+1, nx)
Flux in y-direction
Methods
-
step()
: Solve the pressure equation to obtain pressure and flux. Updateself.p
andself.v
. -
solve(mat, q)
: Method to solve the system of linear equations. Default isscipy.sparse.linalg.spsolve(mat, q)
. You can override this method to use a different solver.
SaturationEquation(self, grid=None, q=None, phi=None, s=None, f_fn=None, v=None, df_fn=None)
A class to model and solve the (water) saturation equation under water injection,
Attributes
-
grid
:Grid
Simulation grid -
q
:ndarray, (ny, nx) | (ny*nx,)
Integrated source term. -
phi
:ndarray, (ny, nx) | (ny*nx,)
Porosity -
f_fn
:callable
Water fractional flow functionf_fn(s)
-
v
:dict
ofndarray
'x'
:ndarray, (ny, nx+1)
Flux in x-direction'y'
:ndarray, (ny+1, nx)
Flux in y-direction
-
df_fn
:callable
(optional) Derivative (element-wise) of water fractional flow functiondf_fn(s)
. It is used to compute the jacobian of the residual function. IfNone
, the jacobian is approximated by the solver (which can be slow). -
s
:ndarray, (ny, nx) | (ny*nx,)
Water saturation
Methods
-
step(dt)
: Solve saturation forward in time bydt
. Updateself.s
. -
solve(residual, s0, residual_jac=None)
: Method to perform the minimization of the residual. Default isscipy.optimize.nonlin.nonlin_solve(residual, s0, jacobian=residual_jac)
. Ifresidual_jac
isNone
, defaults to'krylov'
. You can override this method to use a different solver.
Useful functions for reservoir simulation tasks.
linear_mobility(s, mu_w, mu_o, s_wir, s_oir, deriv=False)
Function to compute water and oil mobility with a linear model.
Parameters
-
s
:ndarray, (ny, nx) | (ny*nx,)
Water saturation -
mu_w
:float
Viscosity of water -
mu_o
:float
Viscosity of oil -
s_wir
:float
Irreducible water saturation -
s_oir
:float
Irreducible oil saturation -
deriv
:bool
IfTrue
, also return derivatives
Returns
if deriv=False:
lamb_w, lamb_o
:ndarray, (ny, nx) | (ny*nx,)
lamb_w
: water mobilitylamb_o
: oil mobility
if deriv=True:
lamb_w, lamb_o, dlamb_w, dlamb_o
:ndarray, (ny, nx) | (ny*nx,)
lamb_w
: water mobilitylamb_o
: oil mobilitydlamb_w
: derivative of water mobilitydlamb_o
: derivative of oil mobility
quadratic_mobility(s, mu_w, mu_o, s_wir, s_oir, deriv=False)
Function to compute water and oil mobility with a quadratic model.
Parameters
-
s
:ndarray, (ny, nx) | (ny*nx,)
Water saturation -
mu_w
:float
Viscosity of water -
mu_o
:float
Viscosity of oil -
s_wir
:float
Irreducible water saturation -
s_oir
:float
Irreducible oil saturation -
deriv
:bool
IfTrue
, also return derivatives
Returns
if deriv=False:
lamb_w, lamb_o
:ndarray, (ny, nx) | (ny*nx,)
lamb_w
: water mobilitylamb_o
: oil mobility
if deriv=True:
lamb_w, lamb_o, dlamb_w, dlamb_o
:ndarray, (ny, nx) | (ny*nx,)
lamb_w
: water mobilitylamb_o
: oil mobilitydlamb_w
: derivative of water mobilitydlamb_o
: derivative of oil mobility
f_fn(s, mobi_fn)
Water fractional flow function.
Parameters
-
s
:ndarray, (ny, nx) | (ny*nx,)
Water saturation -
mobi_fn
:callable
Mobility functionlamb_w, lamb_o = mobi_fn(s)
where:lamb_w
: water mobilitylamb_o
: oil mobility
Returns
ndarray, (ny, nx) | (ny*nx,)
Fractional flow
df_fn(s, mobi_fn)
Derivative function (element-wise) of water fractional flow.
Parameters
-
s
:ndarray, (ny, nx) | (ny*nx,)
Water saturation -
mobi_fn
:callable
Mobility functionlamb_w, lamb_o, dlamb_w, dlamb_o = mobi_fn(s, deriv=True)
where:lamb_w
: water mobilitylamb_o
: oil mobilitydlamb_w
: derivative of water mobilitydlamb_o
: derivative of oil mobility
Returns
ndarray, (ny, nx) | (ny*nx,)
Fractional flow derivative
lamb_fn(s, mobi_fn)
Total mobility function.
Parameters
-
s
:ndarray, (ny, nx) | (ny*nx,)
Water saturation -
mobi_fn
:callable
Mobility functionlamb_w, lamb_o = mobi_fn(s)
where:lamb_w
: water mobilitylamb_o
: oil mobility
Returns
ndarray, (ny, nx) | (ny*nx,)
Total mobility
Aarnes J.E., Gimse T., Lie KA. (2007) An Introduction to the Numerics of Flow in Porous Media using Matlab. In: Hasle G., Lie KA., Quak E. (eds) Geometric Modelling, Numerical Simulation, and Optimization. Springer, Berlin, Heidelberg