From 20ecd37709bd5cc11817b71a57f7bd517211d8be Mon Sep 17 00:00:00 2001 From: =?UTF-8?q?Atgeirr=20Fl=C3=B8=20Rasmussen?= Date: Tue, 25 Mar 2014 11:14:11 +0100 Subject: [PATCH] Use currentControls() field of well state class. Now the well state is consulted for the controls to use when assembling the well control equations, instead of the (immutable) wells_ struct. Also, added dummy implementation of updateWellControls(). --- opm/autodiff/FullyImplicitBlackoilSolver.cpp | 31 +++++++++++++++----- opm/autodiff/FullyImplicitBlackoilSolver.hpp | 7 +++-- 2 files changed, 28 insertions(+), 10 deletions(-) diff --git a/opm/autodiff/FullyImplicitBlackoilSolver.cpp b/opm/autodiff/FullyImplicitBlackoilSolver.cpp index 1f52070c4..55d460733 100644 --- a/opm/autodiff/FullyImplicitBlackoilSolver.cpp +++ b/opm/autodiff/FullyImplicitBlackoilSolver.cpp @@ -687,7 +687,7 @@ namespace { FullyImplicitBlackoilSolver:: assemble(const V& pvdt, const BlackoilState& x , - const WellStateFullyImplicitBlackoil& xw ) + WellStateFullyImplicitBlackoil& xw ) { // Create the primary variables. const SolutionState state = variableState(x, xw); @@ -748,7 +748,8 @@ namespace { } addWellEq(state); - addWellControlEq(state); + updateWellControls(xw); + addWellControlEq(state, xw); } @@ -933,7 +934,17 @@ namespace { - void FullyImplicitBlackoilSolver::addWellControlEq(const SolutionState& state) + void FullyImplicitBlackoilSolver::updateWellControls(WellStateFullyImplicitBlackoil& /*xw*/) const + { + // Do stuff. + } + + + + + + void FullyImplicitBlackoilSolver::addWellControlEq(const SolutionState& state, + const WellStateFullyImplicitBlackoil& xw) { // Handling BHP and SURFACE_RATE wells. @@ -945,14 +956,18 @@ namespace { M rate_distr(nw, np*nw); for (int w = 0; w < nw; ++w) { const WellControls* wc = wells_.ctrls[w]; - if (well_controls_get_current_type(wc) == BHP) { - bhp_targets[w] = well_controls_get_current_target(wc); + // The current control in the well state overrides + // the current control set in the Wells struct, which + // is instead treated as a default. + const int current = xw.currentControls()[w]; + if (well_controls_iget_type(wc, current) == BHP) { + bhp_targets[w] = well_controls_iget_target(wc, current); rate_targets[w] = -1e100; - } else if (well_controls_get_current_type( wc ) == SURFACE_RATE) { + } else if (well_controls_iget_type(wc, current) == SURFACE_RATE) { bhp_targets[w] = -1e100; - rate_targets[w] = well_controls_get_current_target(wc); + rate_targets[w] = well_controls_iget_target(wc, current); { - const double * distr = well_controls_get_current_distr( wc ); + const double * distr = well_controls_iget_distr(wc, current); for (int phase = 0; phase < np; ++phase) { rate_distr.insert(w, phase*nw + w) = distr[phase]; } diff --git a/opm/autodiff/FullyImplicitBlackoilSolver.hpp b/opm/autodiff/FullyImplicitBlackoilSolver.hpp index c9c8ac487..4cd7f1528 100644 --- a/opm/autodiff/FullyImplicitBlackoilSolver.hpp +++ b/opm/autodiff/FullyImplicitBlackoilSolver.hpp @@ -166,15 +166,18 @@ namespace Opm { addOldWellEq(const SolutionState& state); void - addWellControlEq(const SolutionState& state); + addWellControlEq(const SolutionState& state, + const WellStateFullyImplicitBlackoil& xw); void addWellEq(const SolutionState& state); + void updateWellControls(WellStateFullyImplicitBlackoil& xw) const; + void assemble(const V& dtpv, const BlackoilState& x, - const WellStateFullyImplicitBlackoil& xw); + WellStateFullyImplicitBlackoil& xw); V solveJacobianSystem() const;