add mpi serialization for Events

This commit is contained in:
Arne Morten Kvarving 2019-12-06 23:26:38 +01:00
parent 74aca3fefb
commit 4be5a5207d
3 changed files with 38 additions and 0 deletions

View File

@ -30,6 +30,7 @@
#include <opm/parser/eclipse/EclipseState/IOConfig/IOConfig.hpp>
#include <opm/parser/eclipse/EclipseState/IOConfig/RestartConfig.hpp>
#include <opm/parser/eclipse/EclipseState/Edit/EDITNNC.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Events.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/OilVaporizationProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/TimeMap.hpp>
#include <opm/parser/eclipse/EclipseState/SimulationConfig/SimulationConfig.hpp>
@ -883,6 +884,12 @@ std::size_t packSize(const OilVaporizationProperties& data,
packSize(data.maxDRVDT(), comm);
}
std::size_t packSize(const Events& data,
Dune::MPIHelper::MPICommunicator comm)
{
return packSize(data.events(), comm);
}
////// pack routines
template<class T>
@ -1775,6 +1782,14 @@ void pack(const OilVaporizationProperties& data,
pack(data.maxDRSDT(), buffer, position, comm);
pack(data.maxDRSDT_allCells(), buffer, position, comm);
pack(data.maxDRVDT(), buffer, position, comm);
}
void pack(const Events& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
pack(data.events(), buffer, position, comm);
}
/// unpack routines
@ -2977,6 +2992,15 @@ void unpack(OilVaporizationProperties& data,
maxDRSDT_allCells, maxDRVDT);
}
void unpack(Events& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
DynamicVector<uint64_t> events;
unpack(events, buffer, position, comm);
data = Events(events);
}
} // end namespace Mpi
RestartValue loadParallelRestart(const EclipseIO* eclIO, SummaryState& summaryState,
const std::vector<Opm::RestartKey>& solutionKeys,

View File

@ -65,6 +65,7 @@ class EDITNNC;
class EndpointScaling;
class Equil;
class EquilRecord;
class Events;
class FoamConfig;
class FoamData;
class InitConfig;
@ -513,6 +514,7 @@ ADD_PACK_PROTOTYPES(EDITNNC)
ADD_PACK_PROTOTYPES(EndpointScaling)
ADD_PACK_PROTOTYPES(Equil)
ADD_PACK_PROTOTYPES(EquilRecord)
ADD_PACK_PROTOTYPES(Events)
ADD_PACK_PROTOTYPES(FoamConfig)
ADD_PACK_PROTOTYPES(FoamData)
ADD_PACK_PROTOTYPES(EclHysterConfig)

View File

@ -36,6 +36,7 @@
#include <opm/parser/eclipse/EclipseState/InitConfig/InitConfig.hpp>
#include <opm/parser/eclipse/EclipseState/IOConfig/IOConfig.hpp>
#include <opm/parser/eclipse/EclipseState/IOConfig/RestartConfig.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Events.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/OilVaporizationProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/TimeMap.hpp>
#include <opm/parser/eclipse/EclipseState/SimulationConfig/SimulationConfig.hpp>
@ -1186,6 +1187,17 @@ BOOST_AUTO_TEST_CASE(OilVaporizationProperties)
}
BOOST_AUTO_TEST_CASE(Events)
{
#ifdef HAVE_MPI
Opm::Events val1(Opm::DynamicVector<uint64_t>({1,2,3,4,5}));
auto val2 = PackUnpack(val1);
BOOST_CHECK(std::get<1>(val2) == std::get<2>(val2));
BOOST_CHECK(val1 == std::get<0>(val2));
#endif
}
bool init_unit_test_func()
{
return true;