add mpi serialization for UDQActive

This commit is contained in:
Arne Morten Kvarving 2019-12-12 11:54:46 +01:00
parent 01d05c9585
commit ec1a18ca70
3 changed files with 145 additions and 0 deletions

View File

@ -1354,6 +1354,36 @@ std::size_t packSize(const UDQConfig& data,
packSize(data.typeCount(), comm);
}
std::size_t packSize(const UDQActive::InputRecord& data,
Dune::MPIHelper::MPICommunicator comm)
{
return packSize(data.input_index, comm) +
packSize(data.udq, comm) +
packSize(data.wgname, comm) +
packSize(data.control, comm);
}
std::size_t packSize(const UDQActive::Record& data,
Dune::MPIHelper::MPICommunicator comm)
{
return packSize(data.udq, comm) +
packSize(data.input_index, comm) +
packSize(data.use_index, comm) +
packSize(data.wgname, comm) +
packSize(data.control, comm) +
packSize(data.uad_code, comm) +
packSize(data.use_count, comm);
}
std::size_t packSize(const UDQActive& data,
Dune::MPIHelper::MPICommunicator comm)
{
return packSize(data.getInputRecords(), comm) +
packSize(data.getOutputRecords(), comm) +
packSize(data.getUdqKeys(), comm) +
packSize(data.getWgKeys(), comm);
}
////// pack routines
template<class T>
@ -2726,6 +2756,39 @@ void pack(const UDQConfig& data,
pack(data.typeCount(), buffer, position, comm);
}
void pack(const UDQActive::InputRecord& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
pack(data.input_index, buffer, position, comm);
pack(data.udq, buffer, position, comm);
pack(data.wgname, buffer, position, comm);
pack(data.control, buffer, position, comm);
}
void pack(const UDQActive::Record& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
pack(data.udq, buffer, position, comm);
pack(data.input_index, buffer, position, comm);
pack(data.use_index, buffer, position, comm);
pack(data.wgname, buffer, position, comm);
pack(data.control, buffer, position, comm);
pack(data.uad_code, buffer, position, comm);
pack(data.use_count, buffer, position, comm);
}
void pack(const UDQActive& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
pack(data.getInputRecords(), buffer, position, comm);
pack(data.getOutputRecords(), buffer, position, comm);
pack(data.getUdqKeys(), buffer, position, comm);
pack(data.getWgKeys(), buffer, position, comm);
}
/// unpack routines
template<class T>
@ -4664,6 +4727,44 @@ void unpack(UDQConfig& data,
assignmentsMap, units, inputIndex, typeCount);
}
void unpack(UDQActive::InputRecord& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
unpack(data.input_index, buffer, position, comm);
unpack(data.udq, buffer, position, comm);
unpack(data.wgname, buffer, position, comm);
unpack(data.control, buffer, position, comm);
}
void unpack(UDQActive::Record& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
unpack(data.udq, buffer, position, comm);
unpack(data.input_index, buffer, position, comm);
unpack(data.use_index, buffer, position, comm);
unpack(data.wgname, buffer, position, comm);
unpack(data.control, buffer, position, comm);
unpack(data.uad_code, buffer, position, comm);
unpack(data.use_count, buffer, position, comm);
}
void unpack(UDQActive& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
std::vector<UDQActive::InputRecord> inputRecords;
std::vector<UDQActive::Record> outputRecords;
std::unordered_map<std::string,std::size_t> udqKeys, wgKeys;
unpack(inputRecords, buffer, position, comm);
unpack(outputRecords, buffer, position, comm);
unpack(udqKeys, buffer, position, comm);
unpack(wgKeys, buffer, position, comm);
data = UDQActive(inputRecords, outputRecords, udqKeys, wgKeys);
}
#define INSTANTIATE_PACK_VECTOR(T) \
template std::size_t packSize(const std::vector<T>& data, \
Dune::MPIHelper::MPICommunicator comm); \

View File

@ -44,6 +44,7 @@
#include <opm/parser/eclipse/EclipseState/Schedule/Group/Group.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/TimeMap.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/UDQ/UDQAssign.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/UDQ/UDQActive.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/Well.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/WellTestConfig.hpp>
#include <opm/parser/eclipse/EclipseState/Util/OrderedMap.hpp>
@ -634,6 +635,9 @@ ADD_PACK_PROTOTYPES(ThresholdPressure)
ADD_PACK_PROTOTYPES(TimeMap)
ADD_PACK_PROTOTYPES(TimeMap::StepData)
ADD_PACK_PROTOTYPES(UDAValue)
ADD_PACK_PROTOTYPES(UDQActive)
ADD_PACK_PROTOTYPES(UDQActive::InputRecord)
ADD_PACK_PROTOTYPES(UDQActive::Record)
ADD_PACK_PROTOTYPES(UDQAssign)
ADD_PACK_PROTOTYPES(UDQAssign::AssignRecord)
ADD_PACK_PROTOTYPES(UDQASTNode)

View File

@ -43,6 +43,7 @@
#include <opm/parser/eclipse/EclipseState/Schedule/MSW/Valve.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/OilVaporizationProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/TimeMap.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/UDQ/UDQActive.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/UDQ/UDQAssign.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/UDQ/UDQASTNode.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/UDQ/UDQConfig.hpp>
@ -1791,6 +1792,45 @@ BOOST_AUTO_TEST_CASE(UDQConfig)
}
BOOST_AUTO_TEST_CASE(UDQActiveInputRecord)
{
#ifdef HAVE_MPI
Opm::UDQActive::InputRecord val1(1, "test1", "test2",
Opm::UDAControl::WCONPROD_ORAT);
auto val2 = PackUnpack(val1);
BOOST_CHECK(std::get<1>(val2) == std::get<2>(val2));
BOOST_CHECK(val1 == std::get<0>(val2));
#endif
}
BOOST_AUTO_TEST_CASE(UDQActiveRecord)
{
#ifdef HAVE_MPI
Opm::UDQActive::Record val1("test1", 1, 2, "test2",
Opm::UDAControl::WCONPROD_ORAT);
auto val2 = PackUnpack(val1);
BOOST_CHECK(std::get<1>(val2) == std::get<2>(val2));
BOOST_CHECK(val1 == std::get<0>(val2));
#endif
}
BOOST_AUTO_TEST_CASE(UDQActive)
{
#ifdef HAVE_MPI
Opm::UDQActive val1({Opm::UDQActive::InputRecord(1, "test1", "test2",
Opm::UDAControl::WCONPROD_ORAT)},
{Opm::UDQActive::Record("test1", 1, 2, "test2",
Opm::UDAControl::WCONPROD_ORAT)},
{{"test1", 1}}, {{"test2", 2}});
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;