add mpi serialization for ASTNode

This commit is contained in:
Arne Morten Kvarving 2019-12-17 14:58:39 +01:00
parent 5be1955e47
commit 2d682a0838
3 changed files with 62 additions and 0 deletions

View File

@ -31,6 +31,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/Action/ASTNode.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Events.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Group/GuideRateConfig.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/MessageLimits.hpp>
@ -1564,6 +1565,17 @@ std::size_t packSize(const Tuning& data,
packSize(data.getResetValues(), comm);
}
std::size_t packSize(const Action::ASTNode& data,
Dune::MPIHelper::MPICommunicator comm)
{
return packSize(data.type, comm) +
packSize(data.func_type, comm) +
packSize(data.func, comm) +
packSize(data.argList(), comm) +
packSize(data.getNumber(), comm) +
packSize(data.childrens(), comm);
}
////// pack routines
template<class T>
@ -3155,6 +3167,18 @@ void pack(const Tuning& data,
pack(data.getResetValues(), buffer, position, comm);
}
void pack(const Action::ASTNode& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
pack(data.type, buffer, position, comm);
pack(data.func_type, buffer, position, comm);
pack(data.func, buffer, position, comm);
pack(data.argList(), buffer, position, comm);
pack(data.getNumber(), buffer, position, comm);
pack(data.childrens(), buffer, position, comm);
}
/// unpack routines
template<class T>
@ -5421,6 +5445,25 @@ void unpack(Tuning& data, std::vector<char>& buffer, int& position,
DDSLIM, TRGDPR, XXXDPR, XXXDPR_has_value, ResetValue);
}
void unpack(Action::ASTNode& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
TokenType token;
FuncType func_type;
std::string func;
std::vector<std::string> argList;
double number;
std::vector<Action::ASTNode> children;
unpack(token, buffer, position, comm);
unpack(func_type, buffer, position, comm);
unpack(func, buffer, position, comm);
unpack(argList, buffer, position, comm);
unpack(number, buffer, position, comm);
unpack(children, buffer, position, comm);
data = Action::ASTNode(token, func_type, func, argList, number, children);
}
#define INSTANTIATE_PACK_VECTOR(T) \
template std::size_t packSize(const std::vector<T>& data, \
Dune::MPIHelper::MPICommunicator comm); \

View File

@ -63,6 +63,11 @@ namespace Opm
{
class Actdims;
namespace Action {
class ASTNode;
}
class Aqudims;
class ColumnSchema;
class Connection;
@ -598,6 +603,7 @@ void unpack(char* str, std::size_t length, std::vector<char>& buffer, int& posit
Dune::MPIHelper::MPICommunicator comm);
ADD_PACK_PROTOTYPES(Actdims)
ADD_PACK_PROTOTYPES(Action::ASTNode)
ADD_PACK_PROTOTYPES(Aqudims)
ADD_PACK_PROTOTYPES(ColumnSchema)
ADD_PACK_PROTOTYPES(Connection)

View File

@ -38,6 +38,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/Action/ASTNode.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Events.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Group/GConSale.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Group/Group.hpp>
@ -2123,6 +2124,18 @@ BOOST_AUTO_TEST_CASE(Tuning)
}
BOOST_AUTO_TEST_CASE(ASTNode)
{
#ifdef HAVE_MPI
Opm::Action::ASTNode child(number, FuncType::field, "test3", {"test2"}, 2.0, {});
Opm::Action::ASTNode val1(number, FuncType::field, "test1", {"test2"}, 1.0, {child});
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;