Merge pull request #400 from atgeirr/serialize_SCON

Added implementation of serialize_SCON
This commit is contained in:
Joakim Hove
2018-05-14 13:26:21 +02:00
committed by GitHub
5 changed files with 145 additions and 11 deletions

View File

@@ -243,6 +243,7 @@ if(ENABLE_ECL_OUTPUT)
tests/test_Wells.cpp
tests/test_writenumwells.cpp
tests/test_serialize_ICON.cpp
tests/test_serialize_SCON.cpp
)
endif()

View File

@@ -22,6 +22,11 @@
#include <vector>
// Missing definitions (really belong in ert/ecl_well/well_const.h, but not
// defined there)
#define SCON_KH_INDEX 3
// Forward declarations
namespace Opm {
@@ -30,11 +35,14 @@ namespace Opm {
class EclipseState;
class Schedule;
class Well;
class UnitSystem;
} // Opm
namespace Opm { namespace RestartIO { namespace Helpers {
const double UNIMPLEMENTED_VALUE = 1e-100; // placeholder for values not yet available
std::vector<double>
createDoubHead(const EclipseState& es,
const Schedule& sched,
@@ -60,6 +68,12 @@ namespace Opm { namespace RestartIO { namespace Helpers {
int niconz, // Number of elements per completion in ICON, should be entry 32 from createInteHead.
const std::vector<const Well*>& sched_wells);
std::vector<double> serialize_SCON(int lookup_step, // The integer index used to look up dynamic properties, e.g. the number of well.
int ncwmax, // Max number of completions per well, should be entry 17 from createInteHead.
int nsconz, // Number of elements per completion in SCON, should be entry 33 from createInteHead.
const std::vector<const Well*>& sched_wells,
const UnitSystem& units);
}}} // Opm::RestartIO::Helpers
#endif // OPM_WRITE_RESTART_HELPERS_HPP

View File

@@ -17,11 +17,65 @@
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#include <opm/output/eclipse/WriteRestartHelpers.hpp>
#include <ert/ecl_well/well_const.h> // containts ICON_XXX_INDEX
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
#include <opm/output/eclipse/WriteRestartHelpers.hpp>
#include <opm/parser/eclipse/Units/UnitSystem.hpp>
#include <vector>
// ----------------------------------------------------------------------------
std::vector<double>
Opm::RestartIO::Helpers::
serialize_SCON(int lookup_step,
int ncwmax,
int nsconz,
const std::vector<const Well*>& sched_wells,
const UnitSystem& units)
// ----------------------------------------------------------------------------
{
const size_t well_field_size = ncwmax * nsconz;
std::vector<double> data(sched_wells.size() * well_field_size, 0);
size_t well_offset = 0;
for (const Opm::Well* well : sched_wells) {
const auto& completions = well->getCompletions( lookup_step );
size_t completion_offset = 0;
bool explicit_ctf_not_found = false;
for (const auto& completion : completions) {
const size_t offset = well_offset + completion_offset;
const auto& ctf = completion.getConnectionTransmissibilityFactorAsValueObject();
if (ctf.hasValue()) {
// CTF explicitly set in deck, overrides calculation
// from Peaceman model. We should also give the Kh
// factor, we output an explicitly invalid value
// instead. This is acceptable since it will not be
// used (the explicit CTF factor is used instead).
const double ctf_SI = ctf.getValue();
const double ctf_output = units.from_si(UnitSystem::measure::transmissibility, ctf_SI);
data[ offset + SCON_CF_INDEX ] = ctf_output;
data[ offset + SCON_KH_INDEX ] = UNIMPLEMENTED_VALUE;
} else {
// CTF not set in deck, Peaceman formula used to
// compute it. Here we should store the data for the
// connection required to recalculate the CTF (the Kh
// factor), as well as the actual CTF used by the
// simulator, but that requires access to more data
// from the simulator. As an interim measure we write
// invalid values and give a warning.
data[ offset + SCON_CF_INDEX ] = UNIMPLEMENTED_VALUE;
data[ offset + SCON_KH_INDEX ] = UNIMPLEMENTED_VALUE;
explicit_ctf_not_found = true;
}
completion_offset += nsconz;
}
if (explicit_ctf_not_found) {
OpmLog::warning("restart output completion data missing",
"Explicit connection transmissibility factors for well " + well->name() + " missing, writing dummy values to restart file.");
}
well_offset += well_field_size;
}
return data;
}
// ----------------------------------------------------------------------------
std::vector<int>
Opm::RestartIO::Helpers::

View File

@@ -19,10 +19,6 @@
#include <config.h>
#include <iostream> // @@
#include <algorithm> // @@
#include <iterator> // @@
#define BOOST_TEST_MODULE serialize_ICON_TEST
#include <boost/test/unit_test.hpp>
@@ -88,12 +84,6 @@ BOOST_AUTO_TEST_CASE( serialize_icon_test )
}
w_offset += (ICONZ * ncwmax);
}
std::copy(icondata.begin(),
icondata.end(),
std::ostream_iterator<int>(std::cout, " "));
std::cout << std::endl;
BOOST_CHECK_EQUAL(1, 1);// @@@
}
};

View File

@@ -0,0 +1,75 @@
/*
Copyright 2016 Statoil 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 serialize_SCON_TEST
#include <boost/test/unit_test.hpp>
#include <opm/parser/eclipse/Parser/Parser.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
#include <ert/ecl_well/well_const.h> // containts SCON_CF_INDEX
#include <opm/output/eclipse/WriteRestartHelpers.hpp>
BOOST_AUTO_TEST_CASE( serialize_scon_test )
{
const Opm::Deck deck(Opm::Parser{}.parseFile("FIRST_SIM.DATA"));
Opm::UnitSystem units(Opm::UnitSystem::UnitType::UNIT_TYPE_METRIC); // Unit system used in deck FIRST_SIM.DATA.
const Opm::EclipseState state(deck);
const Opm::Schedule schedule(deck, state);
const Opm::TimeMap timemap(deck);
for (size_t tstep = 0; tstep != timemap.numTimesteps(); ++tstep) {
const size_t ncwmax = schedule.getMaxNumCompletionsForWells(tstep);
const int SCONZ = 40; // normally obtained from InteHead
const auto wells = schedule.getWells(tstep);
const std::vector<double> scondata =
Opm::RestartIO::Helpers::serialize_SCON(tstep,
ncwmax,
SCONZ,
wells,
units);
size_t w_offset = 0;
for (const auto w : wells) {
size_t c_offset = 0;
for (const auto c : w->getCompletions(tstep)) {
const size_t offset = w_offset + c_offset;
const auto ctf = c.getConnectionTransmissibilityFactorAsValueObject();
const double expected =
ctf.hasValue()
? units.from_si(Opm::UnitSystem::measure::transmissibility, ctf.getValue())
: Opm::RestartIO::Helpers::UNIMPLEMENTED_VALUE;
BOOST_CHECK_EQUAL(scondata[offset + SCON_CF_INDEX],
expected);
BOOST_CHECK_EQUAL(scondata[offset + SCON_KH_INDEX],
Opm::RestartIO::Helpers::UNIMPLEMENTED_VALUE);
c_offset += SCONZ;
}
w_offset += (SCONZ * ncwmax);
}
}
};