Interactive diagnostics on the SAIGUP dataset

This example sets up interactive diagnostics on a corner point grid. The example uses the same geological dataset as in saigupField1phExample with a number of arbitrary wells. The example is short, as the primary purpose is to create the interactive figures produced by interactiveDiagnostics.

Contents

Set up grid and rock

mrstModule add diagnostics gridtools mrst-gui
grdecl = fullfile(ROOTDIR, 'examples', 'data', 'SAIGUP', 'SAIGUP.GRDECL');

if ~exist(grdecl, 'file'),
   error('SAIGUP model data is not available.')
end

grdecl = readGRDECL(grdecl);

actnum        = grdecl.ACTNUM;
grdecl.ACTNUM = ones(prod(grdecl.cartDims),1);
G             = processGRDECL(grdecl, 'checkgrid', false);
G             = computeGeometry(G(1));

rock = grdecl2Rock(grdecl, G.cells.indexMap);
is_pos                = rock.perm(:, 3) > 0;
rock.perm(~is_pos, 3) = min(rock.perm(is_pos, 3));
rock.perm = convertFrom(rock.perm, milli*darcy);

Set up wells

The wells are simply placed by a double for loop over the logical indices. The resulting wells give producers along the edges of the domain and injectors near the center.

pv = poreVolume(G, rock);
ijk = gridLogicalIndices(G);
W = [];
gcz = G.cells.centroids;
[pi,ii] = deal(1);
for i = 5:12:G.cartDims(1)
    for j = 5:25:G.cartDims(2)
        c = ijk{1} == i & ijk{2} == j;
        if any(c)
            c = find(c);
            x = gcz(c(1), 1); y = gcz(c(1), 2);
            if x > 500 && x < 2000 && y < 7000 && y > 1000
                % Set up rate controlled injectors for cells in the middle
                % of the domain
                val = sum(pv)/(1000*day);
                type = 'rate';
                name = ['I' num2str(ii)];
                ii = ii + 1;
            else
                % Set up BHP controlled producers at the boundary of the
                % domain
                val = 250*barsa;
                type = 'bhp';
                name = ['P' num2str(pi)];
                pi = pi + 1;
            end

            W = addWell(W, G, rock, c, 'Type', type, 'Val', val, 'Name', name);
        end
    end
end
% Plot the resulting well setup
clf
plotCellData(G, rock.poro),
plotWell(G, W)
view(-100, 25)

Set up a reservoir state

This is not strictly needed for the diagnostics application, but it allows us to see the reservoir component distribution around producers grouped by time of flight. The distribution of phases is done based on a simple hydrostatic approximation by using cell centroids. Some phase mixing is present.

state = initResSol(G, 200*barsa, [0 0 0]);

gcz = G.cells.centroids(:, 3);
height = max(gcz) - min(gcz);
pos = 1 - (gcz - min(gcz))./height;

oil = pos > 0.4 & pos < 0.8;
gas = pos > 0.7;
wat = pos < 0.45;

state.s(wat, 1) = 1;
state.s(oil, 2) = 1;
state.s(gas, 3) = 1;

state.s = bsxfun(@rdivide, state.s, sum(state.s, 2));

clf;
rgb = @(x) x(:, [2 3 1]);
plotCellData(G, rgb(state.s), 'EdgeColor', 'k', 'EdgeAlpha', 0.1)
view(-100, 25)
axis tight

Run the interactive diagnostics

We setup the interactive plot and print the help page to show functionality.

help interactiveDiagnostics
clf;
interactiveDiagnostics(G, rock, W, 'state', state);
view(-100, 25)
axis tight
 Launch an interactive diagnostics session
 
  SYNOPSIS:
    interactiveDiagnostics(G, rock, W);
    interactiveDiagnostics(G, rock, W, 'state', state)
    interactiveDiagnostics(G, rock, W, dataset, 'state', state)
 
  DESCRIPTION:
    This function launches an interactive session for doing flow
    diagnostics. The functionality differs slightly based on the input
    arguments given:
      - If a dataset is given, this set of cellwise data will be available
      for visualization.
      - If a state is given, this state will allow for visualization of
      the component ratios inside a drainage volume. The flux from this
      state can also be used to calculate time of flight if "computeFlux"
      is disabled, for instance if the user has some external means of
      computing fluxes.
 
    Once the initialization is complete, two windows will be produced:
      - A plotting window, showing the reservoir along with the wells and
      visualized quantitites. 
          In the plotting window, it is possible to click wells to get
          additional information, such as the allocation factors per
          perforation in the well, the phase distribution grouped by time
          of flight and the corresonding pore volumes swept/drained.
 
      - A controller window which is used to alter the state of the
      plotting window:
      Wells can be selected for display (if a well is selected, the
      corresponding drainage (producer) or flooding (injector) volumes will
      be visualized. A simple playback function can be used to show
      propagation of time of flight. Different quantitites can be
      visualized to get a better understanding of the system.
 
      A set of buttons provide (experimental) well editing, access to
      visualization of the well pair connections, Phi/F diagram with
      Lorenz' coefficient etc.
 
 
  REQUIRED PARAMETERS:
 
    G    - Valid grid structure.
 
    rock - Rock with valid permeability and porosity fields.
 
    W    - A set of wells which are compatible with the incompTPFA solver.
 
 
  OPTIONAL PARAMETERS (supplied in 'key'/value pairs ('pn'/pv ...)):
 
   'state' - Reservoir state containing fluid saturations and optionally
             flux and pressure (if 'computeFlux' is false)
 
   'tracerfluid' - Fluid used for tracer computation. This defaults to a
             trivial fluid.
  
   'fluid'  - Fluid used for mobility calculations if 'useMobilityArrival'
              is enabled.
 
   'LinSolve' - Linear solver for pressure systems. Defaults to mldivide
   (backslash)
 
   'computeFlux' - If set to false, fluxes are extracted from the provided
                   state keyword argument. This requires a state to be
                   provided and can be used if the fluxes are computed
                   externally (for instance from a expensive full-physics
                   simulation)
 
   'useMobilityArrival' - If the well plot showing nearby saturations
                          should plot mobility instead of saturations. This
                          may be interesting in some cases because the
                          mobile fluids are more likely to be extracted.
                          However, this plot is often dominated by very
                          mobile gas regions.
 
  RETURNS:
 
    Nothing. Creates two figures.
 
  EXAMPLE:
      G = computeGeometry(cartGrid([10, 10, 2])); 
      rock = struct('poro', ones(G.cells.num, 1), 'perm', ones(G.cells.num, 1)*darcy)
  
      W = verticalWell([], G, rock, 1,  1, [], 'Val', 1)
      W = verticalWell(W, G, rock,  1,  10, [], 'Val', 1)
      W = verticalWell(W, G, rock,  10, 1, [], 'Val', 1)
      W = verticalWell(W, G, rock,  10, 10, [], 'Val', 1)
      W = verticalWell(W, G, rock,  5,  5, [], 'Val', 0)
      interactiveDiagnostics(G, rock, W);
  
  SEE ALSO:
    Diagnostics examples