mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
WellContainer: add serialization support
This commit is contained in:
parent
e1942d145f
commit
7f36bac579
@ -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) {
|
||||
|
@ -113,6 +113,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();
|
||||
|
Loading…
Reference in New Issue
Block a user