Basic Flow-Solver Tutorial

The purpose of this example is to give an overview of how to set up and use a standard two-point pressure solver to solve the single-phase pressure equation

$$\nabla\cdot v = q, \qquad v=\textbf{--}\frac{K}{\mu}\nabla p,$$

for a flow driven by Dirichlet and Neumann boundary conditions. Our geological model will be simple a Cartesian grid with anisotropic, homogeneous permeability.

In this tutorial example, you will learn about:

  1. the grid structure,
  2. how to specify rock and fluid data,
  3. the structure of the data-objects used to hold solution,
  4. how to assemble and solve linear systems,
  5. useful routines for visualizing and interacting with the grids and simulation results.

Contents

Define geometry

Construct a Cartesian grid of size 10-by-10-by-4 cells, where each cell has dimension 1-by-1-by-1. Because our flow solvers are applicable for general unstructured grids, the Cartesian grid is here represented using an unstructured format, in which cells, faces, nodes, etc. are given explicitly.

nx = 10; ny = 10; nz = 4;
G = cartGrid([nx, ny, nz]);
display(G);
G = 

       cells: [1x1 struct]
       faces: [1x1 struct]
       nodes: [1x1 struct]
    cartDims: [10 10 4]
        type: {'tensorGrid'  'cartGrid'}
     griddim: 3

After the grid structure is generated, we plot the geometry.

clf, plotGrid(G); view(3)

Process geometry

Having set up the basic structure, we continue to compute centroids and volumes of the cells and centroids, normals, and areas for the faces. For a Cartesian grid, this information can trivially be computed, but is given explicitly so that the flow solver is compatible with fully unstructured grids.

G = computeGeometry(G);

Set rock and fluid data

The only parameters in the single-phase pressure equation are the permeability $K$ and the fluid viscosity $\mu$. We set the permeability to be homogeneous and anisotropic

$$ K = \left(\begin{array}{ccc}
     1000 & 0  & 0 \\ 0 & 100 & 0 \\ 0 & 0 & 10 \end{array}\right) $$

The viscosity is specified by saying that the reservoir is filled with a single fluid, for which de default viscosity value equals unity. Our flow solver is written for a general incompressible flow and requires the evaluation of a total mobility, which is provided by the fluid object.

gravity off
rock.perm = repmat([1000, 100, 10].* milli*darcy(), [G.cells.num, 1]);
fluid     = initSingleFluid('mu' ,    1*centi*poise     , ...
                            'rho', 1014*kilogram/meter^3);

Initialize reservoir simulator

To simplify communication among different flow and transport solvers, all unknowns are collected in a structure. Here this structure is initialized with uniform initial reservoir pressure equal 0 and (single-phase) saturation equal 0.0 (using the default behavior of initResSol)

resSol = initResSol(G, 0.0);
display(resSol);
resSol = 

    pressure: [400x1 double]
        flux: [1380x1 double]
           s: [400x1 double]

Impose Dirichlet boundary conditions

Our flow solvers automatically assume no-flow conditions on all outer (and inner) boundaries; other type of boundary conditions need to be specified explicitly.

Here, we impose Neumann conditions (flux of 1 m^3/day) on the global left-hand side. The fluxes must be given in units of m^3/s, and thus we need to divide by the number of seconds in a day (day()). Similarly, we set Dirichlet boundary conditions p = 0 on the global right-hand side of the grid, respectively. For a single-phase flow, we need not specify the saturation at inflow boundaries. Similarly, fluid composition over outflow faces (here, right) is ignored by pside.

bc = fluxside([], G, 'LEFT',  1*meter^3/day());
bc = pside   (bc, G, 'RIGHT', 0);
display(bc);
bc = 

     face: [80x1 int32]
     type: {1x80 cell}
    value: [80x1 double]
      sat: []

Construct linear system

Compute transmissibilities based on input grid and rock properties and use this to set up the linear system.

T = computeTrans(G, rock, 'Verbose', true);

% Solve linear system construced from T and bc to obtain solution for flow
% and pressure in the reservoir. The <matlab:help('incompTPFA')
% option> 'MatrixOutput=true' adds the system matrix A to resSol to enable
% inspection of the matrix.
resSol = incompTPFA(resSol, G, T, fluid, ...
                   'bc', bc, 'MatrixOutput', true);
display(resSol);
Computing one-sided transmissibilities...	Elapsed time is 0.008676 seconds.

resSol = 

        pressure: [400x1 double]
            flux: [1380x1 double]
               s: [400x1 double]
               A: [400x400 double]
             rhs: [400x1 double]
    facePressure: [1380x1 double]

Inspect results

The resSol object contains the matrix used to solve the system.

spy(resSol.A);

We then plot convert the computed pressure to unit 'bar' before plotting result.

clf
plotCellData(G, convertTo(resSol.pressure(1:G.cells.num), barsa()), ...
             'EdgeColor', 'k');
title('Cell Pressure [bar]')
xlabel('x'), ylabel('y'), zlabel('Depth');
view(3); shading faceted; camproj perspective; axis tight;
colorbar