mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
add mpi serialization for WellConnections
This commit is contained in:
@@ -37,6 +37,7 @@
|
||||
#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/WellConnections.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>
|
||||
@@ -1031,6 +1032,15 @@ std::size_t packSize(const WellEconProductionLimits& data,
|
||||
packSize(data.minReservoirFluidRate(), comm);
|
||||
}
|
||||
|
||||
std::size_t packSize(const WellConnections& data,
|
||||
Dune::MPIHelper::MPICommunicator comm)
|
||||
{
|
||||
return packSize(data.getHeadI(), comm) +
|
||||
packSize(data.getHeadJ(), comm) +
|
||||
packSize(data.getNumRemoved(), comm) +
|
||||
packSize(data.getConnections(), comm);
|
||||
}
|
||||
|
||||
////// pack routines
|
||||
|
||||
template<class T>
|
||||
@@ -2074,6 +2084,16 @@ void pack(const WellEconProductionLimits& data,
|
||||
pack(data.minReservoirFluidRate(), buffer, position, comm);
|
||||
}
|
||||
|
||||
void pack(const WellConnections& data,
|
||||
std::vector<char>& buffer, int& position,
|
||||
Dune::MPIHelper::MPICommunicator comm)
|
||||
{
|
||||
pack(data.getHeadI(), buffer, position, comm);
|
||||
pack(data.getHeadJ(), buffer, position, comm);
|
||||
pack(data.getNumRemoved(), buffer, position, comm);
|
||||
pack(data.getConnections(), buffer, position, comm);
|
||||
}
|
||||
|
||||
/// unpack routines
|
||||
|
||||
template<class T>
|
||||
@@ -3503,6 +3523,22 @@ void unpack(WellEconProductionLimits& data,
|
||||
minReservoirFluidRate);
|
||||
}
|
||||
|
||||
void unpack(WellConnections& data,
|
||||
std::vector<char>& buffer, int& position,
|
||||
Dune::MPIHelper::MPICommunicator comm)
|
||||
{
|
||||
int headI, headJ;
|
||||
size_t numRemoved;
|
||||
std::vector<Connection> connections;
|
||||
|
||||
unpack(headI, buffer, position, comm),
|
||||
unpack(headJ, buffer, position, comm),
|
||||
unpack(numRemoved, buffer, position, comm),
|
||||
unpack(connections, buffer, position, comm),
|
||||
|
||||
data = WellConnections(headI, headJ, numRemoved, connections);
|
||||
}
|
||||
|
||||
#define INSTANTIATE_PACK_VECTOR(T) \
|
||||
template std::size_t packSize(const std::vector<T>& data, \
|
||||
Dune::MPIHelper::MPICommunicator comm); \
|
||||
|
||||
@@ -114,6 +114,7 @@ class VISCREFRecord;
|
||||
class ViscrefTable;
|
||||
class WATDENTRecord;
|
||||
class WatdentTable;
|
||||
class WellConnections;
|
||||
class Welldims;
|
||||
class WellEconProductionLimits;
|
||||
class WellFoamProperties;
|
||||
@@ -582,6 +583,7 @@ ADD_PACK_PROTOTYPES(WATDENTRecord)
|
||||
ADD_PACK_PROTOTYPES(WatdentTable)
|
||||
ADD_PACK_PROTOTYPES(Well::WellInjectionProperties)
|
||||
ADD_PACK_PROTOTYPES(Well::WellGuideRate)
|
||||
ADD_PACK_PROTOTYPES(WellConnections)
|
||||
ADD_PACK_PROTOTYPES(Welldims)
|
||||
ADD_PACK_PROTOTYPES(WellEconProductionLimits)
|
||||
ADD_PACK_PROTOTYPES(WellFoamProperties)
|
||||
|
||||
@@ -1428,6 +1428,22 @@ BOOST_AUTO_TEST_CASE(WellGuideRate)
|
||||
}
|
||||
|
||||
|
||||
BOOST_AUTO_TEST_CASE(WellConnections)
|
||||
{
|
||||
#ifdef HAVE_MPI
|
||||
Opm::Connection conn(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);
|
||||
Opm::WellConnections val1(1, 2, 3, {conn, conn});
|
||||
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;
|
||||
|
||||
Reference in New Issue
Block a user