add mpi serialization for NNC and EDITNNC

This commit is contained in:
Arne Morten Kvarving 2019-11-28 15:12:16 +01:00
parent d9268980f3
commit fb7d4b7cb9
3 changed files with 70 additions and 0 deletions

View File

@ -22,6 +22,8 @@
#endif
#include "ParallelRestart.hpp"
#include <opm/parser/eclipse/EclipseState/Grid/NNC.hpp>
#include <opm/parser/eclipse/EclipseState/Edit/EDITNNC.hpp>
#include <opm/parser/eclipse/EclipseState/SimulationConfig/ThresholdPressure.hpp>
#include <dune/common/parallel/mpitraits.hh>
@ -214,6 +216,16 @@ std::size_t packSize(const ThresholdPressure& data, Dune::MPIHelper::MPICommunic
packSize(data.pressureTable(), comm);
}
std::size_t packSize(const NNC& data, Dune::MPIHelper::MPICommunicator comm)
{
return packSize(data.data(), comm);
}
std::size_t packSize(const EDITNNC& data, Dune::MPIHelper::MPICommunicator comm)
{
return packSize(data.data(), comm);
}
////// pack routines
template<class T>
@ -434,6 +446,18 @@ void pack(const ThresholdPressure& data, std::vector<char>& buffer, int& positio
pack(data.pressureTable(), buffer, position, comm);
}
void pack(const NNC& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
pack(data.data(), buffer, position, comm);
}
void pack(const EDITNNC& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
pack(data.data(), buffer, position, comm);
}
/// unpack routines
template<class T>
@ -666,6 +690,22 @@ void unpack(ThresholdPressure& data, std::vector<char>& buffer, int& position,
data = ThresholdPressure(active, restart, thpTable, pTable);
}
void unpack(NNC& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
std::vector<NNCdata> res;
unpack(res, buffer, position, comm);
data = NNC(res);
}
void unpack(EDITNNC& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
std::vector<NNCdata> res;
unpack(res, buffer, position, comm);
data = EDITNNC(res);
}
} // end namespace Mpi
RestartValue loadParallelRestart(const EclipseIO* eclIO, SummaryState& summaryState,
const std::vector<Opm::RestartKey>& solutionKeys,

View File

@ -36,6 +36,9 @@
namespace Opm
{
class EDITNNC;
class NNC;
struct NNCdata;
class ThresholdPressure;
namespace Mpi
@ -186,6 +189,9 @@ ADD_PACK_PROTOTYPES(data::Segment)
ADD_PACK_PROTOTYPES(data::Solution)
ADD_PACK_PROTOTYPES(data::Well)
ADD_PACK_PROTOTYPES(data::WellRates)
ADD_PACK_PROTOTYPES(EDITNNC)
ADD_PACK_PROTOTYPES(NNC)
ADD_PACK_PROTOTYPES(NNCdata)
ADD_PACK_PROTOTYPES(RestartKey)
ADD_PACK_PROTOTYPES(RestartValue)
ADD_PACK_PROTOTYPES(std::string)

View File

@ -25,6 +25,8 @@
#include <boost/test/unit_test.hpp>
#include <opm/simulators/utils/ParallelRestart.hpp>
#include <opm/parser/eclipse/EclipseState/Edit/EDITNNC.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/NNC.hpp>
#include <opm/parser/eclipse/EclipseState/SimulationConfig/ThresholdPressure.hpp>
#include <opm/output/eclipse/RestartValue.hpp>
@ -249,6 +251,28 @@ BOOST_AUTO_TEST_CASE(ThresholdPressure)
}
BOOST_AUTO_TEST_CASE(EDITNNC)
{
#if HAVE_MPI
Opm::EDITNNC val1({{1,2,1.0},{2,3,2.0}});
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(NNC)
{
#if HAVE_MPI
Opm::NNC val1({{1,2,1.0},{2,3,2.0}});
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;