Contents

1D flow diagnostics adjoint optimization example

This example demonstrates the diagnostics based well optimization framework on a very simple 1D problem with two injectors and a single producer.

We begin by defining a simple 10x1x1 grid discretized into 101 fine cells.

mrstModule add ad-fi diagnostics gridtools
N = 100 + 1;
% Define the middle of the domain
mid = floor(N/2);

G = cartGrid([N, 1, 1], [10, 1, 1]);
G = computeGeometry(G);

Setup porosity field, wells and other definitions

We divide the problem into two regions: One region $x < 5$ where the porosity is set to 0.5 and a region $x > 5$ where the porosity is 1. In the center, $x = 5$, the porosity is set to the average (0.75).

Two injectors are defined. One at the left boundary of the domain (x = 0) and another at the right end of the domain (x = 10). Both are allocated the same injection rate, even though they are in different porosity regions.

A single BHP controlled producer is placed in the exact middle of the domain, effectively making the reservoir into two distinct regions: The low porosity region swept by the first injector and the high porosity region swept by the other injector.

rock.poro = ones(G.cells.num, 1);
lowporo =  .5;
highporo = 1;

% Left domain
rock.poro(1:mid) = lowporo;
rock.poro(mid+1:end) = highporo;
rock.poro(mid) = (lowporo + highporo)/2 ;
% Trivial permeability
rock.perm = ones(G.cells.num, 1);

W = [];
% Left injector
W = verticalWell(W, G, rock, 1, 1, [], 'Val', 1/day, 'Type', 'rate', 'Name', 'Injector 1');
% Right injnector
W = verticalWell(W, G, rock, N, 1, [], 'Val', 1/day, 'Type', 'rate', 'Name', 'Injector 2');
% Add producer with zero bottom hole pressure for pressure support
W = verticalWell(W, G, rock, mid, 1, [], 'Val', 0, 'Type', 'bhp', 'Name', 'Producer');


T = computeTrans(G, rock);
pv = poreVolume(G, rock);

% Reservoir state

state0 = initResSol(G, 0*barsa, [1 0 0]);
% state0.wellSol = initWellSol(W, 0);

% Reservoir fluid and ad-system
mu = [1 1 1];
n = [1 1 1];

fluid_ad = initSimpleADIFluid('mu', mu, 'n', n);

s = initADISystem({'Oil', 'Water'}, G, rock, fluid_ad);

Define objective function and plot domain

The Lorenz coefficient is used here, which will have values between 0 (homogenous displacement, perfect sweep) and 1 (infinitely hetereogenous displacement). Because we want to improve the sweep of the well configuration, this is a good choice. The left injector has half the pore volume between it and the producer when compared to the other injector, and so the optimal configuration should set the ratio between the injectors' rates to 1/2.

We at the same time specify the minimum rates per well to a small value, which will not be violated at the optimum in this simple case.

Once the problem definition is complete we plot the wells and porosity on the 1D domain.

objective = getObjectiveDiagnostics(G, rock, 'minlorenz');

minRate = W(1).val./1000;
clf
x = G.cells.centroids(:,1);
stairs(x, rock.poro, '-k', 'LineWidth', 2);
grid on
ylim([0, 1.1])
hold on
colors = {'r', 'g', 'b'};
l = {'Porosity'};
% Plot each well
for i = 1:numel(W)
    c = W(i).cells(1);
    plot(x(c, 1), rock.poro(c), 'O', 'MarkerEdgeColor', 'k', ...
                    'MarkerFaceColor', colors{i}, 'MarkerSize', 8);
    l = [l; W(i).name];
end
legend(l, 'location', 'southeast')
title('Porosity and well placements')
ylabel('Porosity')
xlabel('X coordinates')

Optimize well rates

We optimize the well rates using a steepest descent-like implementation which accounts for well limits and uses adjoints to calculate the sensitivites of the objective function. Because the framework uses diagnostics for the objective evaluations and adjoints for the sensitivities, the cost per function evaluation is low.

We set the tolerance so the iteration has converged when the change between iterations is less than 5%. During the optimization, the progress will be plotted.

clf
[D_best, W_best, history] = optimizeTOF(G, W, fluid_ad, pv, T, s,...
                                     state0, minRate, objective, ...
                                     'plotProgress', true, ...
                                     'verbose', true, ...
                                     'deltatol', 0.05);
Obj: 0.173377 at initial controls 
Obj: 0.323985 (delta: 0.9, alpha: 8e-06) at iteration 0
Obj: 0.324053 (delta: 0.9, alpha: 8e-07) at iteration 1
Obj: 0.324015 (delta: 0.9, alpha: 8e-08) at iteration 2
Obj: 0.324088 (delta: 0.9, alpha: 8e-09) at iteration 3
Obj: 0.323995 (delta: 0.9, alpha: 8e-10) at iteration 4
Obj: 0.099343 (delta: -0.4, alpha: 8e-11) at iteration 5
================ VALUE IMPROVED ==================
Obj: 0.324032 (delta: 2, alpha: 8e-10) at iteration 0
Obj: 0.0253323 (delta: -0.7, alpha: 8e-11) at iteration 1
================ VALUE IMPROVED ==================
Obj: 0.324069 (delta: 1e+01, alpha: 8e-10) at iteration 0
Obj: 0.0488718 (delta: 0.9, alpha: 8e-11) at iteration 1
Obj: 0.0179193 (delta: -0.3, alpha: 8e-12) at iteration 2
================ VALUE IMPROVED ==================
Obj: 0.0563162 (delta: 2, alpha: 8e-11) at iteration 0
Obj: 0.0105036 (delta: -0.4, alpha: 8e-12) at iteration 1
================ VALUE IMPROVED ==================
Obj: 0.0637655 (delta: 5, alpha: 8e-11) at iteration 0
Obj: 0.00308722 (delta: -0.7, alpha: 8e-12) at iteration 1
================ VALUE IMPROVED ==================
Obj: 0.0697707 (delta: 2e+01, alpha: 8e-11) at iteration 0
Obj: 0.00423779 (delta: 0.4, alpha: 8e-12) at iteration 1
Obj: 0.00237375 (delta: -0.2, alpha: 8e-13) at iteration 2
================ VALUE IMPROVED ==================
Obj: 0.00495155 (delta: 1, alpha: 8e-12) at iteration 0
Obj: 0.00166025 (delta: -0.3, alpha: 8e-13) at iteration 1
================ VALUE IMPROVED ==================
Obj: 0.00566534 (delta: 2, alpha: 8e-12) at iteration 0
Obj: 0.000946724 (delta: -0.4, alpha: 8e-13) at iteration 1
================ VALUE IMPROVED ==================
Obj: 0.00637917 (delta: 6, alpha: 8e-12) at iteration 0
Obj: 0.000233168 (delta: -0.8, alpha: 8e-13) at iteration 1
================ VALUE IMPROVED ==================
Obj: 0.00710396 (delta: 3e+01, alpha: 8e-12) at iteration 0
Obj: 0.000670118 (delta: 2, alpha: 8e-13) at iteration 1
Obj: 0.00016181 (delta: -0.3, alpha: 8e-14) at iteration 2
================ VALUE IMPROVED ==================
Obj: 0.000741478 (delta: 4, alpha: 8e-13) at iteration 0
Obj: 9.92573e-05 (delta: -0.4, alpha: 8e-14) at iteration 1
================ VALUE IMPROVED ==================
Obj: 0.000804017 (delta: 7, alpha: 8e-13) at iteration 0
Obj: 0.000161809 (delta: 0.6, alpha: 8e-14) at iteration 1
Obj: 9.75872e-05 (delta: -0.02, alpha: 8e-15) at iteration 2
Optimization loop DONE, improvement 0.17 (relative improvement:  1)

Plot the time of flight of the initial and the best well rates

We plot the time of flight in the domain for each configuration. The optimized well configuration shouws equal arrival times from each injector, indicating an optimal solution.

clf
pw = @() plotWell(G, W, 'height', 0.05, 'radius', 0, 'Color', 'red');

subplot(2,2,1:2)
ind = 1;
x = G.cells.centroids(:,1);
initial = history.D(1).tof(:,1);
optimized = history.D(end).tof(:,1);

hold on
plot(x, [initial, optimized]./max(initial), 'LineWidth', 2)

grid on
title('Time of flight')
legend({'Initial TOF', 'Optimized TOF'});

subplot(2,2,3)
plotCellData(G, history.D(1).tof(:,1))
title('Initial')
[htop, htext] = pw();
set(htext, 'Interpreter', 'None')
axis tight off equal
view(20, 40)

subplot(2,2,4)
colormap winter
plotCellData(G, history.D(end).tof(:,1))
title('Optimized')
[htop, htext] = pw();
set(htext, 'Interpreter', 'None')
axis tight off equal
view(20, 40)

Plot sweep and flow capacity diagrams

The Lorenz coefficient is really a derived measure from the flow-capacity diagram, equal to the integral of the deviation from a perfect, linear flow curve where every unit of fluid injected corresponds to a unit of recovery. We compute the F and Phi quantities for the base case and the optimum and plot both.

At the same time, we plot the sweep diagrams for the same configurations.

[F_0, Phi_0] = computeFandPhi(pv, history.D(1).tof);
[F_end, Phi_end] = computeFandPhi(pv, history.D(end).tof);

clf;

subplot(2,1,1)
plot([Phi_0, Phi_end], [F_0, F_end], 'linewidth', 2)
legend({'Equal rates', 'Optimized rates'}, 'location', 'SouthEast')
xlabel('\Phi')
ylabel('F')
title('Flow-capacity diagram')
grid on

subplot(2,1,2)
[Ev_0, tD_0] = computeSweep(F_0, Phi_0);
[Ev_end, tD_end] = computeSweep(F_end, Phi_end);
plot([tD_0, tD_end], [Ev_0, Ev_end], 'linewidth', 2)
legend({'Equal rates', 'Optimized rates'}, 'location', 'SouthEast')
xlabel('\Phi')
ylabel('F')
title('Sweep')
grid on
axis tight