Merge pull request #4440 from akva2/groupstate_serialize

GroupState: add serialization support
This commit is contained in:
Bård Skaflestad 2023-02-15 23:33:50 +01:00 committed by GitHub
commit 498ca650d5
No known key found for this signature in database
GPG Key ID: 4AEE18F83AFDEB23
4 changed files with 81 additions and 4 deletions

View File

@ -33,6 +33,25 @@ GroupState::GroupState(std::size_t np) :
num_phases(np)
{}
GroupState GroupState::serializationTestObject()
{
GroupState result(3);
result.m_production_rates = {{"test1", {1.0, 2.0}}};
result.production_controls = {{"test2", Group::ProductionCMode::LRAT}};
result.prod_red_rates = {{"test3", {3.0, 4.0, 5.0}}};
result.inj_red_rates = {{"test4", {6.0, 7.0}}};
result.inj_surface_rates = {{"test5", {8.0}}};
result.inj_resv_rates = {{"test6", {9.0, 10.0}}};
result.inj_rein_rates = {{"test7", {11.0}}};
result.inj_vrep_rate = {{"test8", 12.0}, {"test9", 13.0}};
result.m_grat_sales_target = {{"test10", 14.0}};
result.m_gpmaint_target = {{"test11", 15.0}};
result.injection_controls = {{{Phase::FOAM, "test12"}, Group::InjectionCMode::REIN}};
result.gpmaint_state.add("foo", GPMaint::State::serializationTestObject());
return result;
}
bool GroupState::operator==(const GroupState& other) const {
return this->m_production_rates == other.m_production_rates &&
this->production_controls == other.production_controls &&
@ -43,7 +62,8 @@ bool GroupState::operator==(const GroupState& other) const {
this->inj_vrep_rate == other.inj_vrep_rate &&
this->inj_surface_rates == other.inj_surface_rates &&
this->m_grat_sales_target == other.m_grat_sales_target &&
this->injection_controls == other.injection_controls;
this->injection_controls == other.injection_controls &&
this->gpmaint_state == other.gpmaint_state;
}
//-------------------------------------------------------------------------

View File

@ -33,7 +33,11 @@ namespace Opm {
class GroupState {
public:
GroupState() = default;
explicit GroupState(std::size_t num_phases);
static GroupState serializationTestObject();
bool operator==(const GroupState& other) const;
bool has_production_rates(const std::string& gname) const;
@ -159,9 +163,26 @@ public:
throw std::logic_error("Internal size mismatch when distributing groupData");
}
template<class Serializer>
void serializeOp(Serializer& serializer)
{
serializer(num_phases);
serializer(m_production_rates);
serializer(production_controls);
serializer(prod_red_rates);
serializer(inj_red_rates);
serializer(inj_surface_rates);
serializer(inj_resv_rates);
serializer(inj_rein_rates);
serializer(inj_vrep_rate);
serializer(m_grat_sales_target);
serializer(m_gpmaint_target);
serializer(injection_controls);
serializer(gpmaint_state);
}
private:
std::size_t num_phases;
std::size_t num_phases{};
std::map<std::string, std::vector<double>> m_production_rates;
std::map<std::string, Group::ProductionCMode> production_controls;
std::map<std::string, std::vector<double>> prod_red_rates;

View File

@ -45,15 +45,23 @@ namespace Opm {
template <class T>
class WellContainer {
public:
WellContainer() = default;
WellContainer(std::initializer_list<std::pair<std::string,T>> init_list) {
for (const auto& [name, value] : init_list)
this->add(name, value);
}
static WellContainer serializationTestObject(const T& data)
{
WellContainer<T> result;
result.m_data = {data};
result.index_map = {{"test1", 1}, {"test2", 4}};
return result;
}
bool empty() const {
return this->index_map.empty();
}
@ -167,6 +175,18 @@ public:
return wlist;
}
template<class Serializer>
void serializeOp(Serializer& serializer)
{
serializer(m_data);
serializer(index_map);
}
bool operator==(const WellContainer<T>& rhs) const
{
return this->m_data == rhs.m_data &&
this->index_map == rhs.index_map;
}
private:
void update_if(std::size_t index, const std::string& name, const WellContainer<T>& other) {

View File

@ -36,6 +36,7 @@
#include <opm/simulators/timestepping/SimulatorTimer.hpp>
#include <opm/simulators/timestepping/TimeStepControl.hpp>
#include <opm/simulators/utils/SerializationPackers.hpp>
#include <opm/simulators/wells/GroupState.hpp>
#include <opm/simulators/wells/SegmentState.hpp>
#include <opm/simulators/wells/SingleWellState.hpp>
@ -76,6 +77,7 @@ BOOST_AUTO_TEST_CASE(NAME) \
#define TEST_FOR_TYPE(TYPE) \
TEST_FOR_TYPE_NAMED(TYPE, TYPE)
TEST_FOR_TYPE(GroupState)
TEST_FOR_TYPE(HardcodedTimeStepControl)
TEST_FOR_TYPE(PIDAndIterationCountTimeStepControl)
TEST_FOR_TYPE(PIDTimeStepControl)
@ -113,6 +115,20 @@ BOOST_AUTO_TEST_CASE(SingleWellState)
BOOST_CHECK_MESSAGE(data_out == data_in, "Deserialized SingleWellState differ");
}
BOOST_AUTO_TEST_CASE(WellContainer)
{
auto data_out = Opm::WellContainer<double>::serializationTestObject(1.0);
Opm::Serialization::MemPacker packer;
Opm::Serializer ser(packer);
ser.pack(data_out);
const size_t pos1 = ser.position();
decltype(data_out) data_in;
ser.unpack(data_in);
const size_t pos2 = ser.position();
BOOST_CHECK_MESSAGE(pos1 == pos2, "Packed size differ from unpack size for WellContainer");
BOOST_CHECK_MESSAGE(data_out == data_in, "Deserialized WellContainer differ");
}
BOOST_AUTO_TEST_CASE(EclGenericVanguard)
{
auto in_params = Opm::EclGenericVanguard::serializationTestParams();