spe10: Access to the SPE10 benchmark case

Functions and utilities for creating and manipulating test cases based on the 10th SPE Comparative Solution Project (SPE10). Both the model 1 and the model 2 subproblems are available for download, with routines for setting up and working the grid and petrophysical properties.

Contents

SPE10

Files
getSPE10_model_1_fluid - Construct ADI Fluid Object for Model 1 of Tenth SPE CSP getSPE10_model_1_relperm - Define Oil/Gas Relative Permeability Properties for Model 1 of tenth SPE CSP getSPE10_model_1_rock - Define rock properties for Model 1 of tenth SPE CSP getSPE10rock - Define rock properties for Model 2 of tenth SPE CSP getSPE10setup - Initialise properties for Model 2 of tenth SPE Comparative Solution Project grdeclBox - Make a GRDECL structure for simple corner-point grid, possibly faulted. make_spe10_data - Create on-disk (MAT file) representation of SPE 10 ‘rock’ data. setupSPE10_AD - Undocumented Utility Function SPE10_rock - Define rock properties for Model 2 of tenth SPE CSP SPE10_setup - Initialise properties for Model 2 of tenth SPE Comparative Solution Project
SPE10_rock(varargin)

Define rock properties for Model 2 of tenth SPE CSP

This function is DEPRECATED and will be removed in a future version of MRST. Please switch to using function getSPE10rock instead.

Synopsis:

rock = SPE10_rock
rock = SPE10_rock(layers)
rock = SPE10_rock(I, J, K)
Parameters:
  • layers

    Which of the 85 model layers to include in a specific test. OPTIONAL. If unspecified, all 85 layers (for a total of 60-by-220-by-85 == 1,122,000 cells) are included.

    Some possible values are
    layers = ( 1 : 35).’; % Tarbert formation layers = (36 : 85).’; % Upper Ness formation
  • I,J,K

    Global Cartesian indices identifying subset of rock data. Useful, e.g., in extracting lateral sections. Arrays must satisfy

    ALL((0 < I) & (I <= 60)) && … ALL((0 < J) & (J <= 220)) && … ALL((0 < K) & (K <= 85))

    OPTIONAL. If unspecified, treated as I=1:60,J=1:220,K=1:85 (i.e., the entire model; 1,122,000 cells). In other words, equivalently to an unspecified ‘layers’ input.

Returns:

rock – Rock structure having fields ‘perm’ and ‘poro’ pertaining to the specified layers or cell subset.

Note

The permeability data is returned in units of milli*darcy.

It is the caller’s responsibility to convert this data into MRST’s strict SI-only unit conventions. Function ‘convertFrom’ provides one conversion method.

Example:

rock = SPE10_rock(85)
SPE10_setup(varargin)

Initialise properties for Model 2 of tenth SPE Comparative Solution Project

This function is DEPRECATED and will be removed in a future version of MRST. Please switch to using function getSPE10setup instead.

Synopsis:

[G, W, rock] = SPE10_setup
[G, W, rock] = SPE10_setup(layers)
[G, W, rock] = SPE10_setup(layers, wloc)
Parameters:
  • layers

    Which of the 85 model layers to include in a specific test. OPTIONAL. If unspecified or empty, all 85 layers (a total of 60-by-220-by-85 == 1,122,000 cells) are included.

    Some possible values are
    layers = ( 1 : 35).’; % Tarbert formation layers = (36 : 85).’; % Upper Ness formation
  • wloc

    Location of the five wells. OPTIONAL. If unspecified, we use the default location

    wloc = [ 1, 60, 1, 60, 30;
    1, 1, 220, 220, 110];
Returns:
  • G – MRST grid structure as described in grid_structure.
  • W – Well structure. Injector at 500 Bar, producers at 200 Bar. Inner product ‘ip_tpf’.
  • rock – Rock structure having fields ‘perm’ and ‘poros’ pertaining to the specified layers.

Example:

[G, W, rock] = SPE10_setup(85)
getSPE10_model_1_fluid()

Construct ADI Fluid Object for Model 1 of Tenth SPE CSP

Synopsis:

fluid = getSPE10_model_1_fluid
Parameters:None.
Returns:fluid – Simplified ADI fluid object, suitable for use with simulator function simulateScheduleAD, that represents the petrochemical properties of the fluids in Model 1 of the Tenth SPE CSP. The fluids are immiscible and incompressible with constant viscosities, but there is a non-trivial relative permeability relation between the oil and gas phases.
getSPE10_model_1_relperm()

Define Oil/Gas Relative Permeability Properties for Model 1 of tenth SPE CSP

Synopsis:

kr_deck = getSPE10_model_1_relperm
Parameters:None.
Returns:kr_deck – Stripped-down simulation deck for ADI-fluid (relperm) construction. This input-deck contains just enough basic information to form relative permeability curves from the CSP benchmark data through the use of module ‘ad-props’s function initDeckADIFluid.
getSPE10_model_1_rock()

Define rock properties for Model 1 of tenth SPE CSP

Synopsis:

rock = getSPE10_model_1_rock
Parameters:None.
Returns:rock – Rock structure having fields ‘perm’ and ‘poro’ for the entire 100-by-1-by-20 cell geometry of the first benchmark dataset of the tenth SPE comparative solution project. The permeability distribution is a correlated geostatistically generated field whereas the porosity is homogenous (phi=0.2) throughout the formation.

Note

The permeability data is returned in strict SI units (metres squared).

getSPE10rock(varargin)

Define rock properties for Model 2 of tenth SPE CSP

Synopsis:

rock = getSPE10rock
rock = getSPE10rock(layers)
rock = getSPE10rock(I, J, K)
Parameters:
  • layers

    Which of the 85 model layers to include in a specific test. OPTIONAL. If unspecified, all 85 layers (for a total of 60-by-220-by-85 == 1,122,000 cells) are included.

    Some possible values are
    layers = ( 1 : 35).’; % Tarbert formation layers = (36 : 85).’; % Upper Ness formation
  • I,J,K

    Global Cartesian indices identifying subset of rock data. Useful, e.g., in extracting lateral sections. Arrays must satisfy

    ALL((0 < I) & (I <= 60)) && … ALL((0 < J) & (J <= 220)) && … ALL((0 < K) & (K <= 85))

    OPTIONAL. If unspecified, treated as I=1:60,J=1:220,K=1:85 (i.e., the entire model; 1,122,000 cells). In other words, equivalently to an unspecified ‘layers’ input.

Returns:

rock – Rock structure having fields ‘perm’ and ‘poro’ pertaining to the specified layers or cell subset.

Note

The permeability data is returned in strict SI units (metres squared).

Example:

rock = getSPE10rock(85)
getSPE10setup(varargin)

Initialise properties for Model 2 of tenth SPE Comparative Solution Project

Synopsis:

[G, W, rock] = getSPE10setup
[G, W, rock] = getSPE10setup(layers)
[G, W, rock] = getSPE10setup(layers, wloc)
Parameters:
  • layers

    Which of the 85 model layers to include in a specific test. OPTIONAL. If unspecified or empty, all 85 layers (a total of 60-by-220-by-85 == 1,122,000 cells) are included.

    Some possible values are
    layers = ( 1 : 35).’; % Tarbert formation layers = (36 : 85).’; % Upper Ness formation
  • wloc

    Location of the five wells. OPTIONAL. If unspecified, we use the default location

    wloc = [ 1, 60, 1, 60, 30;
    1, 1, 220, 220, 110];
Returns:
  • G – MRST grid structure as described in grid_structure.
  • W – Well structure. Injector at 500 Bar, producers at 200 Bar. Inner product ‘ip_tpf’.
  • rock – Rock structure having fields ‘perm’ and ‘poro’ pertaining to the specified layers. Permeability values are in strict SI (metres squared).

Example:

[G, W, rock] = getSPE10setup(85)

See also

getSPE10rock, makeSPE10DataAvailable.

grdeclBox(dims, physDims, drop, varargin)

Make a GRDECL structure for simple corner-point grid, possibly faulted.

Synopsis:

grdeclBox(dims, physDims)
grdecl = simpleGrdecl(n)
grdecl = simpleGrdecl(n, drop)
grdecl = simpleGrdecl(n, drop, 'pn1', pv1, ...)

Description:

This function creates a nx-by-ny-by-nz cell corner-point grid. The resulting geometry discretises the physical domain

[0,1] x [0,1] x [0,0.5].

which is deformed by adding sin(x) and sin(x*y) perturbations. In addition, the pillars are skewed. If a second input parameter is provided, a single fault of given drop size is added in the model at approximately x=0.5.

Parameters:
  • n – Three-element vector, [nx, ny, nz], specifying the number of cells in the ‘x’, ‘y’, and ‘z’ coordinate directions respectively.
  • drop – Length of fault drop (i.e., the physical length (in metres) by which the two blocks should be offset (along the fault)). The drop size can be a numeric value or a handle to a function of a spatial coordinate along the fault. OPTIONAL. Default value: drop = 0 (no fault).
  • 'pn'/pv

    List of ‘key’/value pairs defining optional paramters. The supported options are:

    • flat –
      Flag indicating if the sine perturbations should NOT be added. Default: FALSE (add perturbations).
    • undisturbed –
      Flag indicating that skew pillars should not be created. Default: FALSE (make skew pillars)
Returns:

grdecl – A GRDECL structure suitable for further processing by function ‘processGRDECL’ or for permanent storage on disk by function ‘writeGRDECL’.

Example:

% 1) Create a 20-by-20-by-5 faulted grid with a fault drop of 0.15 (m).
grdecl = simpleGrdecl([20, 20, 5], 0.15);

% Create the grid data structure.
G = processGRDECL(grdecl);

% Plot the resulting geometry.
figure, hg1 = plotGrid(G, 'FaceAlpha', 0.625);
view(3), grid on, axis tight

% 2) Create the same 20-by-20-by-5 faulted grid, but now with a
%    variable fault drop implemented as a function handle.
grdecl = simpleGrdecl([20, 20, 5], @(x) 0.05 * (sin(2*pi*x) - 0.5) );

% Create the grid data structure.
G = processGRDECL(grdecl);

% Plot the resulting geometry.
figure, hg2 = plotGrid(G, 'FaceAlpha', 0.625);
view(3), grid on, axis tight
make_spe10_data()

Create on-disk (MAT file) representation of SPE 10 ‘rock’ data.

Synopsis:

ok = make_spe10_data
Parameters:None.
Returns:ok – Status flag indicating successful creation of on-disk ‘rock’ data.

Note

The on-disk representation is a ‘rock’ structure stored in the file ‘spe10_rock.mat’ in the directory identified by

getDatasetPath(‘spe10’)

The permeability data is stored in SI units (metres squared).

See also

getDatasetPath, getSPE10rock, makeSPE10DataAvailable.

setupSPE10_AD(varargin)

Undocumented Utility Function

Contents

Files prodCurves - Calculate derived production data for wells. simulate_spe10 - Simulate the SPE10 base case spe10 - Define a model structure that represents a subset of the SPE10 base case.

simulate_spe10

Simulate the SPE10 base case The simulation uses a mimetic pressure solver and an implicit transport solver.

prodCurves(w, state, fluid)

Calculate derived production data for wells.

Synopsis:

[tot, opr, wpr, wc] = prodCurves(W, state, fluid)
Parameters:
  • W – Well data structure as defined by function addWell.
  • state – Solution state object as defined by functions ‘initResSol’, ‘initWellSol’ and, possibly, flow and/or transport solvers.
  • fluid – Fluid data structure.
Returns:
  • tot – Total well production. One value for each individual well.
  • opr – Oil (liquid) production rate. One value for each individual well.
  • wpr – Water production rate. One value for each individual well.
  • wc – Water cut. One value for each individual well.
LIMITATION:
This function assumes that the solution state object represents incompressible two-phase flow.
spe10(domain)

Define a model structure that represents a subset of the SPE10 base case.