Pressure Solver: Example of a Real Field Model
In the example, we will solve the single-phase pressure equation
using the corner-point geometry from a real North Sea field.
The purpose of this example is to demonstrate how the mimetic flow solver can be applied to compute flow on a real grid model that has degenerate cell geometries and non-neighbouring connections arising from a large number of faults, eroded layers, pinch-outs, and inactive cells. A more thorough examination of the model is given in a separate example.
For confidentiallity reasons, we will not use real data for rock parameters and well paths. Instead, synthetic permeability data are generated by (make-believe) geostatistics and vertical wells are placed quite haphazardly throughout the domain.
Experimentation with this model is continued in another example, in which we solve a simplified two-phase model.
Contents
Check for existence of input model data
grdecl = fullfile(ROOTDIR, 'examples', 'grids', 'GSmodel.grdecl'); if ~exist(grdecl, 'file'), error('Model data is not available.') end
Read and process the model
We start by reading the model from a file in the Eclipse formate (GRDECL). As shown in the more thorough discussion, the model has two components, of which we will only use the first one.
grdecl = readGRDECL(grdecl);
G = processGRDECL(grdecl); clear grdecl;
G = computeGeometry(G(1));
Set rock and fluid data
The permeability is lognormal and isotropic within nine distinct layers and is generated using our simplified 'geostatistics' function and then transformed to lay in the interval 200-2000 mD. The single fluid has density 1000 kg/m^3 and viscosity 1 cP.
gravity off K = logNormLayers(G.cartDims, rand(9,1), 'sigma', 1.5); K = 200 + (K-min(K))/(max(K)-min(K))*1800; rock.perm = K(G.cells.indexMap)*milli*darcy; clear K; fluid = initSingleFluid(); clf, plotCellData(G,log10(rock.perm),'EdgeColor','k'); axis off, view(15,60), h=colorbar('horiz'); zoom(2.5), cs = [200 400 700 1000 1500 2000]; caxis(log10([min(cs) max(cs)]*milli*darcy)); set(h, 'XTick', log10(cs*milli*darcy), 'XTickLabel', num2str(round(cs)'));

Introduce wells
The reservoir is produced using a set production wells controlled by bottom-hole pressure and rate-controlled injectors. Wells are described using a Peacemann model, giving an extra set of equations that need to be assembled. For simplicity, all wells are assumed to be vertical and are assigned using the logical (i,j) subindex.
% Plot grid outline clf subplot('position',[0.02 0.02 0.96 0.96]); plotGrid(G,'FaceColor','none','EdgeAlpha',0.1); axis tight off, zoom(1.1), view(-5,58) % Set six vertical injectors, completed in each layer. nz = G.cartDims(3); I = [ 9, 26, 8, 25, 35, 10]; J = [14, 14, 35, 35, 68, 75]; R = [ 2, 2, 4, 4, 0.5,0.5]*1000*meter^3/day; radius = .1; W = []; for i = 1 : numel(I), W = verticalWell(W, G, rock, I(i), J(i), 1:nz, 'Type', 'rate', ... 'Val', R(i), 'Radius', radius, 'Comp_i', [1,0,0], ... 'name', ['I$_{', int2str(i), '}$']); end plotGrid(G, vertcat(W.cells), 'FaceColor', 'b'); prod_off = numel(W); % Set five vertical producers, completed in each layer. I = [17, 12, 25, 35, 15]; J = [23, 51, 51, 88, 94]; for i = 1 : numel(I), W = verticalWell(W, G, rock, I(i), J(i), 1:nz, 'Type', 'bhp', ... 'Val', 300*barsa(), 'Radius', radius, ... 'name', ['P$_{', int2str(i), '}$']); end plotGrid(G, vertcat(W(prod_off + 1 : end).cells), 'FaceColor', 'r'); plotWell(G,W,'height',200);

Initialize and construct linear system
Initialize solution structures and assemble linear hybrid system from input grid, rock properties, and well structure.
S = computeMimeticIP(G, rock, 'Verbose', true);
W = assembleWellSystem(G, W);
rSol = initResSol(G, 350*barsa);
wSol = initWellSol(W, 300*barsa());
Using inner product: 'ip_simple'. Computing component matrices C and D ... Elapsed time is 0.019604 seconds. Computing cell inner products ... Elapsed time is 6.877634 seconds. Assembling global inner product matrix ... Elapsed time is 0.119223 seconds.
Solve the linear system
Solve linear system construced from S and W to obtain solution for flow and pressure in the reservoir and the wells.
[rSol, wSol] = solveIncompFlow(rSol, wSol, G, S, fluid, 'wells', W);
Plot the fluxes
clf cellNo = rldecode(1:G.cells.num, double(G.cells.numFaces), 2) .'; plotCellData(G, sqrt(accumarray(cellNo, ... abs(convertTo(rSol.cellFlux, meter^3/day)))),'EdgeColor', 'k'); plotWell(G,W,'height',200,'color','c'); colorbar('horiz'), axis tight off; view(20,80) zoom(2)

Plot the pressure distribution
clf plotCellData(G, convertTo(rSol.cellPressure, barsa),'EdgeColor','k'); plotWell(G,W,'height',200,'color','c'); colorbar('horiz'), axis tight off; view(20,70) zoom(2)

#COPYRIGHT_EXAMPLE#
Last time modified: $Id: realField1phExample.m 1945 2009-03-31 09:22:39Z bska $