WGState: add serialization support

This commit is contained in:
Arne Morten Kvarving 2023-02-02 11:52:08 +01:00
parent 5b0c87b946
commit 90fbdcee0a
3 changed files with 48 additions and 1 deletions

View File

@ -32,10 +32,27 @@ WGState::WGState(const PhaseUsage& pu) :
well_test_state{}
{}
WGState WGState::serializationTestObject(const ParallelWellInfo& pinfo)
{
WGState result(PhaseUsage{});
result.well_state = WellState::serializationTestObject(pinfo);
result.group_state = GroupState::serializationTestObject();
result.well_test_state = WellTestState::serializationTestObject();
return result;
}
void WGState::wtest_state(WellTestState wtest_state)
{
wtest_state.filter_wells( this->well_state.wells() );
this->well_test_state = std::move(wtest_state);
}
bool WGState::operator==(const WGState& rhs) const
{
return this->well_state == rhs.well_state &&
this->group_state == rhs.group_state &&
this->well_test_state == rhs.well_test_state;
}
}

View File

@ -26,20 +26,34 @@
namespace Opm {
class ParallelWellInfo;
/*
Microscopic class to handle well , group and well test state.
Microscopic class to handle well, group and well test state.
*/
struct PhaseUsage;
struct WGState {
WGState(const PhaseUsage& pu);
static WGState serializationTestObject(const ParallelWellInfo& pinfo);
void wtest_state(WellTestState wtest_state);
WellState well_state;
GroupState group_state;
WellTestState well_test_state;
bool operator==(const WGState&) const;
template<class Serializer>
void serializeOp(Serializer& serializer)
{
serializer(well_state);
serializer(group_state);
serializer(well_test_state);
}
};
}

View File

@ -146,6 +146,22 @@ BOOST_AUTO_TEST_CASE(WellState)
BOOST_CHECK_MESSAGE(data_out == data_in, "Deserialized WellState differ");
}
BOOST_AUTO_TEST_CASE(WGState)
{
Opm::ParallelWellInfo dummy;
auto data_out = Opm::WGState::serializationTestObject(dummy);
Opm::Serialization::MemPacker packer;
Opm::Serializer ser(packer);
ser.pack(data_out);
const size_t pos1 = ser.position();
decltype(data_out) data_in(Opm::PhaseUsage{});
data_in.well_state = Opm::WellState(dummy);
ser.unpack(data_in);
const size_t pos2 = ser.position();
BOOST_CHECK_MESSAGE(pos1 == pos2, "Packed size differ from unpack size for WGState");
BOOST_CHECK_MESSAGE(data_out == data_in, "Deserialized WGState differ");
}
BOOST_AUTO_TEST_CASE(EclGenericVanguard)
{
auto in_params = Opm::EclGenericVanguard::serializationTestParams();