add mpi serialization for ColumnSchema

This commit is contained in:
Arne Morten Kvarving 2019-11-29 10:18:34 +01:00
parent e08a5c3046
commit 46c5afa490
3 changed files with 59 additions and 0 deletions

View File

@ -25,6 +25,7 @@
#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 <opm/parser/eclipse/EclipseState/Tables/ColumnSchema.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/Rock2dTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/Rock2dtrTable.hpp>
#include <dune/common/parallel/mpitraits.hh>
@ -240,6 +241,18 @@ std::size_t packSize(const Rock2dtrTable& data, Dune::MPIHelper::MPICommunicator
packSize(data.pressureValues(), comm);
}
std::size_t packSize(const ColumnSchema& data, Dune::MPIHelper::MPICommunicator comm)
{
std::size_t res = packSize(data.name(), comm) +
packSize(data.order(), comm) +
packSize(data.getDefaultMode(), comm);
if (data.getDefaultMode() == Table::DEFAULT_CONST) {
res += packSize(data.getDefaultValue(), comm);
}
return res;
}
////// pack routines
template<class T>
@ -486,6 +499,16 @@ void pack(const Rock2dtrTable& data, std::vector<char>& buffer, int& position,
pack(data.pressureValues(), buffer, position, comm);
}
void pack(const ColumnSchema& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
pack(data.name(), buffer, position, comm);
pack(data.order(), buffer, position, comm);
pack(data.getDefaultMode(), buffer, position, comm);
if (data.getDefaultMode() == Table::DEFAULT_CONST)
pack(data.getDefaultValue(), buffer, position, comm);
}
/// unpack routines
template<class T>
@ -754,6 +777,23 @@ void unpack(Rock2dtrTable& data, std::vector<char>& buffer, int& position,
data = Rock2dtrTable(transMultValues, pressureValues);
}
void unpack(ColumnSchema& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
std::string name;
Table::ColumnOrderEnum order;
Table::DefaultAction action;
unpack(name, buffer, position, comm);
unpack(order, buffer, position, comm);
unpack(action, buffer, position, comm);
if (action == Table::DEFAULT_CONST) {
double value;
unpack(value, buffer, position, comm);
data = ColumnSchema(name, order, value);
} else
data = ColumnSchema(name, order, action);
}
} // end namespace Mpi
RestartValue loadParallelRestart(const EclipseIO* eclIO, SummaryState& summaryState,
const std::vector<Opm::RestartKey>& solutionKeys,

View File

@ -36,6 +36,7 @@
namespace Opm
{
class ColumnSchema;
class EDITNNC;
class NNC;
struct NNCdata;
@ -184,6 +185,7 @@ void unpack(char* str, std::size_t length, std::vector<char>& buffer, int& posit
void unpack(T& data, std::vector<char>& buffer, int& position, \
Dune::MPIHelper::MPICommunicator comm);
ADD_PACK_PROTOTYPES(ColumnSchema)
ADD_PACK_PROTOTYPES(data::CellData)
ADD_PACK_PROTOTYPES(data::Connection)
ADD_PACK_PROTOTYPES(data::Rates)

View File

@ -28,6 +28,7 @@
#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/parser/eclipse/EclipseState/Tables/ColumnSchema.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/Rock2dTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/Rock2dtrTable.hpp>
#include <opm/output/eclipse/RestartValue.hpp>
@ -297,6 +298,22 @@ BOOST_AUTO_TEST_CASE(Rock2dtrTable)
}
BOOST_AUTO_TEST_CASE(ColumnSchema)
{
#if HAVE_MPI
Opm::ColumnSchema val1("test1", Opm::Table::INCREASING,
Opm::Table::DEFAULT_LINEAR);
auto val2 = PackUnpack(val1);
BOOST_CHECK(std::get<1>(val2) == std::get<2>(val2));
BOOST_CHECK(val1 == std::get<0>(val2));
Opm::ColumnSchema val3("test2", Opm::Table::DECREASING, 1.0);
val2 = PackUnpack(val3);
BOOST_CHECK(std::get<1>(val2) == std::get<2>(val2));
BOOST_CHECK(val3 == std::get<0>(val2));
#endif
}
bool init_unit_test_func()
{
return true;