mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
WellState FIBO: Return Segment Vectors from report()
This commit extends class WellStateFullyImplictBlackoil to report segment-related quantities as Opm::data::Segment objects (included in Opm::data::WellRates objects). All wells have at least a top segment in the context of WellState FIBO, so there is a meaningful value to report for each well. We put the extraction of segment-related quantities into a new helper function WellStateFullyImplicitBlackoil::reportSegmentResults() to avoid cluttering up the body of report() more than absolutely needed. The primary use-case for this is assigning appropriate values to items 8 through 11 of restart vector RSEG. In turn, this will enable restoring these quantities from a restart file.
This commit is contained in:
parent
519e6ad174
commit
2418df701f
@ -70,6 +70,7 @@ list (APPEND TEST_SOURCE_FILES
|
||||
tests/test_stoppedwells.cpp
|
||||
tests/test_relpermdiagnostics.cpp
|
||||
tests/test_norne_pvt.cpp
|
||||
tests/test_wellstatefullyimplicitblackoil.cpp
|
||||
)
|
||||
|
||||
if(MPI_FOUND)
|
||||
|
@ -579,6 +579,12 @@ namespace Opm
|
||||
}
|
||||
}
|
||||
assert(local_comp_index == this->wells_->well_connpos[ w + 1 ] - this->wells_->well_connpos[ w ]);
|
||||
|
||||
const auto nSeg = this->numSegments(w);
|
||||
for (auto segID = 0*nSeg; segID < nSeg; ++segID) {
|
||||
well.segments[segID + 1] =
|
||||
this->reportSegmentResults(pu, w, segID);
|
||||
}
|
||||
}
|
||||
|
||||
return res;
|
||||
@ -935,6 +941,49 @@ namespace Opm
|
||||
// Well potentials
|
||||
std::vector<double> well_potentials_;
|
||||
|
||||
::Opm::data::Segment
|
||||
reportSegmentResults(const PhaseUsage& pu,
|
||||
const int wellID,
|
||||
const int segmentID) const
|
||||
{
|
||||
auto segRes = ::Opm::data::Segment{};
|
||||
|
||||
const auto segDoF =
|
||||
this->topSegmentIndex(wellID) + segmentID;
|
||||
|
||||
const auto* rate =
|
||||
&this->segRates()[segDoF * this->numPhases()];
|
||||
|
||||
segRes.pressure = this->segPress()[segDoF];
|
||||
|
||||
if (pu.phase_used[Water]) {
|
||||
segRes.rates.set(data::Rates::opt::wat,
|
||||
rate[pu.phase_pos[Water]]);
|
||||
}
|
||||
|
||||
if (pu.phase_used[Oil]) {
|
||||
segRes.rates.set(data::Rates::opt::oil,
|
||||
rate[pu.phase_pos[Oil]]);
|
||||
}
|
||||
|
||||
if (pu.phase_used[Gas]) {
|
||||
segRes.rates.set(data::Rates::opt::gas,
|
||||
rate[pu.phase_pos[Gas]]);
|
||||
}
|
||||
|
||||
segRes.segNumber = segmentID + 1;
|
||||
|
||||
return segRes;
|
||||
}
|
||||
|
||||
int numSegments(const int wellID) const
|
||||
{
|
||||
const auto topSeg = this->topSegmentIndex(wellID);
|
||||
|
||||
return (wellID + 1 == this->numWells()) // Last well?
|
||||
? (this->numSegment() - topSeg)
|
||||
: (this->topSegmentIndex(wellID + 1) - topSeg);
|
||||
}
|
||||
};
|
||||
|
||||
} // namespace Opm
|
||||
|
285
tests/test_wellstatefullyimplicitblackoil.cpp
Normal file
285
tests/test_wellstatefullyimplicitblackoil.cpp
Normal file
@ -0,0 +1,285 @@
|
||||
/*
|
||||
Copyright 2018 Equinor ASA.
|
||||
|
||||
This file is part of the Open Porous Media Project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#include <config.h>
|
||||
|
||||
#define BOOST_TEST_MODULE WellStateFIBOTest
|
||||
|
||||
#include <opm/autodiff/WellStateFullyImplicitBlackoil.hpp>
|
||||
|
||||
#include <boost/test/unit_test.hpp>
|
||||
|
||||
#include <opm/parser/eclipse/Parser/Parser.hpp>
|
||||
#include <opm/parser/eclipse/Parser/ParseContext.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
|
||||
#include <opm/parser/eclipse/Units/Units.hpp>
|
||||
|
||||
#include <opm/core/props/BlackoilPhases.hpp>
|
||||
#include <opm/core/props/phaseUsageFromDeck.hpp>
|
||||
#include <opm/core/wells/WellsManager.hpp>
|
||||
#include <opm/core/wells.h>
|
||||
|
||||
#include <opm/grid/GridManager.hpp>
|
||||
|
||||
#include <cstddef>
|
||||
#include <string>
|
||||
|
||||
struct Setup
|
||||
{
|
||||
Setup(const std::string& filename)
|
||||
: Setup(Opm::Parser{}.parseFile(filename))
|
||||
{}
|
||||
|
||||
Setup(const Opm::Deck& deck)
|
||||
: es (deck)
|
||||
, pu (Opm::phaseUsageFromDeck(es))
|
||||
, grid (es.getInputGrid())
|
||||
, sched(deck, es)
|
||||
{}
|
||||
|
||||
Opm::EclipseState es;
|
||||
Opm::PhaseUsage pu;
|
||||
Opm::GridManager grid;
|
||||
Opm::Schedule sched;
|
||||
};
|
||||
|
||||
namespace {
|
||||
Opm::WellStateFullyImplicitBlackoil
|
||||
buildWellState(const Setup& setup, const std::size_t timeStep)
|
||||
{
|
||||
auto state = Opm::WellStateFullyImplicitBlackoil{};
|
||||
auto state0 = Opm::WellStateFullyImplicitBlackoil{}; // Empty.
|
||||
|
||||
const auto cpress =
|
||||
std::vector<double>(setup.grid.c_grid()->number_of_cells,
|
||||
100.0*Opm::unit::barsa);
|
||||
|
||||
const Opm::WellsManager wmgr{
|
||||
setup.es, setup.sched, timeStep, *setup.grid.c_grid()
|
||||
};
|
||||
|
||||
state.init(wmgr.c_wells(), cpress, setup.sched,
|
||||
setup.sched.getWells(timeStep),
|
||||
timeStep, nullptr, setup.pu);
|
||||
|
||||
state.initWellStateMSWell(wmgr.c_wells(), setup.sched.getWells(timeStep),
|
||||
timeStep, setup.pu, state0);
|
||||
|
||||
return state;
|
||||
}
|
||||
|
||||
void setSegPress(const std::vector<const Opm::Well*>& wells,
|
||||
const std::size_t tstep,
|
||||
Opm::WellStateFullyImplicitBlackoil& wstate)
|
||||
{
|
||||
const auto nWell = wells.size();
|
||||
|
||||
auto& segPress = wstate.segPress();
|
||||
|
||||
for (auto wellID = 0*nWell; wellID < nWell; ++wellID) {
|
||||
const auto* well = wells[wellID];
|
||||
const auto topSegIx = wstate.topSegmentIndex(wellID);
|
||||
const auto pressTop = 100.0 * wellID;
|
||||
|
||||
segPress[topSegIx] = pressTop;
|
||||
|
||||
if (! well->isMultiSegment(tstep)) {
|
||||
continue;
|
||||
}
|
||||
|
||||
const auto nSeg = well->getWellSegments(tstep).size();
|
||||
|
||||
for (auto segID = 0*nSeg + 1; segID < nSeg; ++segID) {
|
||||
segPress[topSegIx + segID] = pressTop + 1.0*segID;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
void setSegRates(const std::vector<const Opm::Well*>& wells,
|
||||
const std::size_t tstep,
|
||||
const Opm::PhaseUsage& pu,
|
||||
Opm::WellStateFullyImplicitBlackoil& wstate)
|
||||
{
|
||||
const auto wat = pu.phase_used[Opm::BlackoilPhases::Aqua];
|
||||
const auto iw = wat ? pu.phase_pos[Opm::BlackoilPhases::Aqua] : -1;
|
||||
|
||||
const auto oil = pu.phase_used[Opm::BlackoilPhases::Liquid];
|
||||
const auto io = oil ? pu.phase_pos[Opm::BlackoilPhases::Liquid] : -1;
|
||||
|
||||
const auto gas = pu.phase_used[Opm::BlackoilPhases::Vapour];
|
||||
const auto ig = gas ? pu.phase_pos[Opm::BlackoilPhases::Vapour] : -1;
|
||||
|
||||
const auto np = wstate.numPhases();
|
||||
|
||||
const auto nWell = wells.size();
|
||||
|
||||
auto& segRates = wstate.segRates();
|
||||
|
||||
for (auto wellID = 0*nWell; wellID < nWell; ++wellID) {
|
||||
const auto* well = wells[wellID];
|
||||
const auto topSegIx = wstate.topSegmentIndex(wellID);
|
||||
const auto rateTop = 1000.0 * wellID;
|
||||
|
||||
if (wat) { segRates[np*topSegIx + iw] = rateTop; }
|
||||
if (oil) { segRates[np*topSegIx + io] = rateTop; }
|
||||
if (gas) { segRates[np*topSegIx + ig] = rateTop; }
|
||||
|
||||
if (! well->isMultiSegment(tstep)) {
|
||||
continue;
|
||||
}
|
||||
|
||||
const auto nSeg = well->getWellSegments(tstep).size();
|
||||
|
||||
for (auto segID = 0*nSeg + 1; segID < nSeg; ++segID) {
|
||||
auto* rates = &segRates[(topSegIx + segID) * np];
|
||||
|
||||
if (wat) { rates[iw] = rateTop + 100.0*segID; }
|
||||
if (oil) { rates[io] = rateTop + 200.0*segID; }
|
||||
if (gas) { rates[ig] = rateTop + 400.0*segID; }
|
||||
}
|
||||
}
|
||||
}
|
||||
} // Anonymous
|
||||
|
||||
BOOST_AUTO_TEST_SUITE(Segment)
|
||||
|
||||
BOOST_AUTO_TEST_CASE(Linearisation)
|
||||
{
|
||||
const Setup setup{ "msw.data" };
|
||||
const auto tstep = std::size_t{0};
|
||||
|
||||
const auto wstate = buildWellState(setup, tstep);
|
||||
|
||||
BOOST_CHECK_EQUAL(wstate.numSegment(), 6 + 1);
|
||||
|
||||
const auto& wells = setup.sched.getWells(tstep);
|
||||
BOOST_CHECK_EQUAL(wells.size(), 2);
|
||||
|
||||
const auto prod01_first = wells[0]->name() == "PROD01";
|
||||
|
||||
BOOST_CHECK_EQUAL(wstate.topSegmentIndex(0), 0);
|
||||
BOOST_CHECK_EQUAL(wstate.topSegmentIndex(1),
|
||||
prod01_first ? 6 : 1);
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE(Pressure)
|
||||
{
|
||||
const Setup setup{ "msw.data" };
|
||||
const auto tstep = std::size_t{0};
|
||||
|
||||
auto wstate = buildWellState(setup, tstep);
|
||||
|
||||
const auto& wells = setup.sched.getWells(tstep);
|
||||
const auto prod01_first = wells[0]->name() == "PROD01";
|
||||
|
||||
setSegPress(wells, tstep, wstate);
|
||||
|
||||
const auto rpt = wstate.report(setup.pu, setup.grid.c_grid()->global_cell);
|
||||
|
||||
{
|
||||
const auto& xw = rpt.at("INJE01");
|
||||
|
||||
BOOST_CHECK_EQUAL(xw.segments.size(), 1); // Top Segment
|
||||
|
||||
const auto& xseg = xw.segments.at(1);
|
||||
|
||||
BOOST_CHECK_EQUAL(xseg.segNumber, 1);
|
||||
BOOST_CHECK_CLOSE(xseg.pressure, prod01_first ? 100.0 : 0.0, 1.0e-10);
|
||||
}
|
||||
|
||||
{
|
||||
const auto expect_nSeg = 6;
|
||||
const auto& xw = rpt.at("PROD01");
|
||||
|
||||
BOOST_CHECK_EQUAL(xw.segments.size(), expect_nSeg);
|
||||
|
||||
const auto pressTop = prod01_first ? 0.0 : 100.0;
|
||||
|
||||
for (auto segID = 0; segID < expect_nSeg; ++segID) {
|
||||
const auto& xseg = xw.segments.at(segID + 1);
|
||||
|
||||
BOOST_CHECK_EQUAL(xseg.segNumber, segID + 1);
|
||||
BOOST_CHECK_CLOSE(xseg.pressure, pressTop + 1.0*segID, 1.0e-10);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE(Rates)
|
||||
{
|
||||
const Setup setup{ "msw.data" };
|
||||
const auto tstep = std::size_t{0};
|
||||
|
||||
auto wstate = buildWellState(setup, tstep);
|
||||
|
||||
const auto& wells = setup.sched.getWells(tstep);
|
||||
const auto prod01_first = wells[0]->name() == "PROD01";
|
||||
|
||||
const auto& pu = setup.pu;
|
||||
|
||||
setSegRates(wells, tstep, pu, wstate);
|
||||
|
||||
const auto rpt = wstate.report(pu, setup.grid.c_grid()->global_cell);
|
||||
|
||||
const auto wat = pu.phase_used[Opm::BlackoilPhases::Aqua];
|
||||
const auto oil = pu.phase_used[Opm::BlackoilPhases::Liquid];
|
||||
const auto gas = pu.phase_used[Opm::BlackoilPhases::Vapour];
|
||||
|
||||
BOOST_CHECK(wat && oil && gas);
|
||||
|
||||
{
|
||||
const auto rateTop = prod01_first ? 1000.0 : 0.0;
|
||||
|
||||
const auto& xw = rpt.at("INJE01");
|
||||
|
||||
BOOST_CHECK_EQUAL(xw.segments.size(), 1); // Top Segment
|
||||
|
||||
const auto& xseg = xw.segments.at(1);
|
||||
|
||||
BOOST_CHECK_EQUAL(xseg.segNumber, 1);
|
||||
BOOST_CHECK_CLOSE(xseg.rates.get(Opm::data::Rates::opt::wat), rateTop, 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(xseg.rates.get(Opm::data::Rates::opt::oil), rateTop, 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(xseg.rates.get(Opm::data::Rates::opt::gas), rateTop, 1.0e-10);
|
||||
}
|
||||
|
||||
{
|
||||
const auto expect_nSeg = 6;
|
||||
const auto& xw = rpt.at("PROD01");
|
||||
|
||||
BOOST_CHECK_EQUAL(xw.segments.size(), expect_nSeg);
|
||||
|
||||
const auto rateTop = prod01_first ? 0.0 : 1000.0;
|
||||
|
||||
for (auto segID = 0; segID < expect_nSeg; ++segID) {
|
||||
const auto& xseg = xw.segments.at(segID + 1);
|
||||
|
||||
BOOST_CHECK_EQUAL(xseg.segNumber, segID + 1);
|
||||
|
||||
BOOST_CHECK_CLOSE(xseg.rates.get(Opm::data::Rates::opt::wat),
|
||||
rateTop + 100.0*segID, 1.0e-10);
|
||||
|
||||
BOOST_CHECK_CLOSE(xseg.rates.get(Opm::data::Rates::opt::oil),
|
||||
rateTop + 200.0*segID, 1.0e-10);
|
||||
|
||||
BOOST_CHECK_CLOSE(xseg.rates.get(Opm::data::Rates::opt::gas),
|
||||
rateTop + 400.0*segID, 1.0e-10);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_SUITE_END()
|
Loading…
Reference in New Issue
Block a user