add mpi serialization for Connection

This commit is contained in:
Arne Morten Kvarving
2019-12-10 14:18:02 +01:00
parent cadd8dc414
commit 89efc489e4
3 changed files with 116 additions and 1 deletions

View File

@@ -36,6 +36,7 @@
#include <opm/parser/eclipse/EclipseState/Schedule/TimeMap.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/VFPInjTable.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/VFPProdTable.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/Connection.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/WellFoamProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/WellPolymerProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/WellTracerProperties.hpp>
@@ -966,6 +967,31 @@ std::size_t packSize(const UDAValue& data,
packSize(data.get<std::string>(), comm));
}
std::size_t packSize(const Connection& data,
Dune::MPIHelper::MPICommunicator comm)
{
return packSize(data.dir(), comm) +
packSize(data.depth(), comm) +
packSize(data.state(), comm) +
packSize(data.satTableId(), comm) +
packSize(data.complnum(), comm) +
packSize(data.CF(), comm) +
packSize(data.Kh(), comm) +
packSize(data.rw(), comm) +
packSize(data.r0(), comm) +
packSize(data.skinFactor(), comm) +
packSize(data.getI(), comm) +
packSize(data.getJ(), comm) +
packSize(data.getK(), comm) +
packSize(data.getSeqIndex(), comm) +
packSize(data.getSegDistStart(), comm) +
packSize(data.getSegDistEnd(), comm) +
packSize(data.getDefaultSatTabId(), comm) +
packSize(data.getCompSegSeqIndex(), comm) +
packSize(data.segment(), comm) +
packSize(data.wellPi(), comm);
}
////// pack routines
template<class T>
@@ -1943,6 +1969,32 @@ void pack(const UDAValue& data,
pack(data.get<std::string>(), buffer, position, comm);
}
void pack(const Connection& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
pack(data.dir(), buffer, position, comm);
pack(data.depth(), buffer, position, comm);
pack(data.state(), buffer, position, comm);
pack(data.satTableId(), buffer, position, comm);
pack(data.complnum(), buffer, position, comm);
pack(data.CF(), buffer, position, comm);
pack(data.Kh(), buffer, position, comm);
pack(data.rw(), buffer, position, comm);
pack(data.r0(), buffer, position, comm);
pack(data.skinFactor(), buffer, position, comm);
pack(data.getI(), buffer, position, comm);
pack(data.getJ(), buffer, position, comm);
pack(data.getK(), buffer, position, comm);
pack(data.getSeqIndex(), buffer, position, comm);
pack(data.getSegDistStart(), buffer, position, comm);
pack(data.getSegDistEnd(), buffer, position, comm);
pack(data.getDefaultSatTabId(), buffer, position, comm);
pack(data.getCompSegSeqIndex(), buffer, position, comm);
pack(data.segment(), buffer, position, comm);
pack(data.wellPi(), buffer, position, comm);
}
/// unpack routines
template<class T>
@@ -3273,6 +3325,51 @@ void unpack(UDAValue& data,
}
}
void unpack(Connection& data,
std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
Connection::Direction dir;
double depth;
Connection::State state;
int satTableId, complnum;
double CF, Kh, rw, r0, skinFactor;
int I, J, K;
size_t seqIndex;
double segDistStart, segDistEnd;
bool defaultSatTabId;
size_t compSegSeqIndex;
int segment;
double wellPi;
unpack(dir, buffer, position, comm);
unpack(depth, buffer, position, comm);
unpack(state, buffer, position, comm);
unpack(satTableId, buffer, position, comm);
unpack(complnum, buffer, position, comm);
unpack(CF, buffer, position, comm);
unpack(Kh, buffer, position, comm);
unpack(rw, buffer, position, comm);
unpack(r0, buffer, position, comm);
unpack(skinFactor, buffer, position, comm);
unpack(I, buffer, position, comm);
unpack(J, buffer, position, comm);
unpack(K, buffer, position, comm);
unpack(seqIndex, buffer, position, comm);
unpack(segDistStart, buffer, position, comm);
unpack(segDistEnd, buffer, position, comm);
unpack(defaultSatTabId, buffer, position, comm);
unpack(compSegSeqIndex, buffer, position, comm);
unpack(segment, buffer, position, comm);
unpack(wellPi, buffer, position, comm);
data = Connection(dir, depth, state, satTableId,
complnum, CF, Kh, rw, r0,
skinFactor, {I,J,K}, seqIndex,
segDistStart, segDistEnd,
defaultSatTabId, compSegSeqIndex,
segment, wellPi);
}
#define INSTANTIATE_PACK_VECTOR(T) \
template std::size_t packSize(const std::vector<T>& data, \
Dune::MPIHelper::MPICommunicator comm); \

View File

@@ -58,6 +58,7 @@ namespace Opm
class Actdims;
class Aqudims;
class ColumnSchema;
class Connection;
class DENSITYRecord;
class DensityTable;
class EclHysterConfig;
@@ -510,6 +511,7 @@ void unpack(char* str, std::size_t length, std::vector<char>& buffer, int& posit
ADD_PACK_PROTOTYPES(Actdims)
ADD_PACK_PROTOTYPES(Aqudims)
ADD_PACK_PROTOTYPES(ColumnSchema)
ADD_PACK_PROTOTYPES(Connection)
ADD_PACK_PROTOTYPES(data::CellData)
ADD_PACK_PROTOTYPES(data::Connection)
ADD_PACK_PROTOTYPES(data::Rates)

View File

@@ -42,6 +42,7 @@
#include <opm/parser/eclipse/EclipseState/Schedule/TimeMap.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/VFPInjTable.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/VFPProdTable.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/Connection.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/WellFoamProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/WellPolymerProperties.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well/WellTracerProperties.hpp>
@@ -316,7 +317,7 @@ BOOST_AUTO_TEST_CASE(Rates)
}
BOOST_AUTO_TEST_CASE(Connection)
BOOST_AUTO_TEST_CASE(dataConnection)
{
#if HAVE_MPI
Opm::data::Connection con1 = getConnection();
@@ -1363,6 +1364,21 @@ BOOST_AUTO_TEST_CASE(UDAValue)
}
BOOST_AUTO_TEST_CASE(Connection)
{
#ifdef HAVE_MPI
Opm::Connection val1(Opm::Connection::Direction::Y,
1.0, Opm::Connection::State::SHUT,
2, 3, 4.0, 5.0, 6.0, 7.0, 8.0,
{9, 10, 11}, 12, 13.0, 14.0, true,
15, 16, 17.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;