mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
add mpi serialization for VFPProdTable
This commit is contained in:
@@ -35,6 +35,7 @@
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/OilVaporizationProperties.hpp>
|
||||
#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/SimulationConfig/SimulationConfig.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/SimulationConfig/ThresholdPressure.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Tables/Aqudims.hpp>
|
||||
@@ -911,6 +912,24 @@ std::size_t packSize(const VFPInjTable& data,
|
||||
packSize(double(), comm);
|
||||
}
|
||||
|
||||
std::size_t packSize(const VFPProdTable& data,
|
||||
Dune::MPIHelper::MPICommunicator comm)
|
||||
{
|
||||
return packSize(data.getTableNum(), comm) +
|
||||
packSize(data.getDatumDepth(), comm) +
|
||||
packSize(data.getFloType(), comm) +
|
||||
packSize(data.getWFRType(), comm) +
|
||||
packSize(data.getGFRType(), comm) +
|
||||
packSize(data.getALQType(), comm) +
|
||||
packSize(data.getFloAxis(), comm) +
|
||||
packSize(data.getTHPAxis(), comm) +
|
||||
packSize(data.getWFRAxis(), comm) +
|
||||
packSize(data.getGFRAxis(), comm) +
|
||||
packSize(data.getALQAxis(), comm) +
|
||||
data.getTable().num_elements() *
|
||||
packSize(double(), comm);
|
||||
}
|
||||
|
||||
////// pack routines
|
||||
|
||||
template<class T>
|
||||
@@ -1831,6 +1850,26 @@ void pack(const VFPInjTable& data,
|
||||
for (size_t i = 0; i < data.getTable().num_elements(); ++i)
|
||||
pack(*(data.getTable().data() + i), buffer, position, comm);
|
||||
}
|
||||
|
||||
void pack(const VFPProdTable& data,
|
||||
std::vector<char>& buffer, int& position,
|
||||
Dune::MPIHelper::MPICommunicator comm)
|
||||
{
|
||||
pack(data.getTableNum(), buffer, position, comm);
|
||||
pack(data.getDatumDepth(), buffer, position, comm);
|
||||
pack(data.getFloType(), buffer, position, comm);
|
||||
pack(data.getWFRType(), buffer, position, comm);
|
||||
pack(data.getGFRType(), buffer, position, comm);
|
||||
pack(data.getALQType(), buffer, position, comm);
|
||||
pack(data.getFloAxis(), buffer, position, comm);
|
||||
pack(data.getTHPAxis(), buffer, position, comm);
|
||||
pack(data.getWFRAxis(), buffer, position, comm);
|
||||
pack(data.getGFRAxis(), buffer, position, comm);
|
||||
pack(data.getALQAxis(), buffer, position, comm);
|
||||
for (size_t i = 0; i < data.getTable().num_elements(); ++i)
|
||||
pack(*(data.getTable().data() + i), buffer, position, comm);
|
||||
}
|
||||
|
||||
/// unpack routines
|
||||
|
||||
template<class T>
|
||||
@@ -3075,6 +3114,45 @@ void unpack(VFPInjTable& data,
|
||||
floAxis, thpAxis, table);
|
||||
}
|
||||
|
||||
void unpack(VFPProdTable& data,
|
||||
std::vector<char>& buffer, int& position,
|
||||
Dune::MPIHelper::MPICommunicator comm)
|
||||
{
|
||||
int tableNum;
|
||||
double datumDepth;
|
||||
VFPProdTable::FLO_TYPE floType;
|
||||
VFPProdTable::WFR_TYPE wfrType;
|
||||
VFPProdTable::GFR_TYPE gfrType;
|
||||
VFPProdTable::ALQ_TYPE alqType;
|
||||
std::vector<double> floAxis, thpAxis, wfrAxis, gfrAxis, alqAxis;
|
||||
VFPProdTable::array_type table;
|
||||
|
||||
unpack(tableNum, buffer, position, comm);
|
||||
unpack(datumDepth, buffer, position, comm);
|
||||
unpack(floType, buffer, position, comm);
|
||||
unpack(wfrType, buffer, position, comm);
|
||||
unpack(gfrType, buffer, position, comm);
|
||||
unpack(alqType, buffer, position, comm);
|
||||
unpack(floAxis, buffer, position, comm);
|
||||
unpack(thpAxis, buffer, position, comm);
|
||||
unpack(wfrAxis, buffer, position, comm);
|
||||
unpack(gfrAxis, buffer, position, comm);
|
||||
unpack(alqAxis, buffer, position, comm);
|
||||
VFPProdTable::extents extents;
|
||||
extents[0] = thpAxis.size();
|
||||
extents[1] = wfrAxis.size();
|
||||
extents[2] = gfrAxis.size();
|
||||
extents[3] = alqAxis.size();
|
||||
extents[4] = floAxis.size();
|
||||
table.resize(extents);
|
||||
for (size_t i = 0; i < table.num_elements(); ++i)
|
||||
unpack(*(table.data() + i), buffer, position, comm);
|
||||
|
||||
data = VFPProdTable(tableNum, datumDepth, floType, wfrType,
|
||||
gfrType, alqType, floAxis, thpAxis,
|
||||
wfrAxis, gfrAxis, alqAxis, table);
|
||||
}
|
||||
|
||||
} // end namespace Mpi
|
||||
RestartValue loadParallelRestart(const EclipseIO* eclIO, SummaryState& summaryState,
|
||||
const std::vector<Opm::RestartKey>& solutionKeys,
|
||||
|
||||
@@ -105,6 +105,7 @@ class TableSchema;
|
||||
class ThresholdPressure;
|
||||
class UDQParams;
|
||||
class VFPInjTable;
|
||||
class VFPProdTable;
|
||||
class VISCREFRecord;
|
||||
class ViscrefTable;
|
||||
class WATDENTRecord;
|
||||
@@ -564,6 +565,7 @@ ADD_PACK_PROTOTYPES(TimeMap)
|
||||
ADD_PACK_PROTOTYPES(TimeMap::StepData)
|
||||
ADD_PACK_PROTOTYPES(UDQParams)
|
||||
ADD_PACK_PROTOTYPES(VFPInjTable)
|
||||
ADD_PACK_PROTOTYPES(VFPProdTable)
|
||||
ADD_PACK_PROTOTYPES(VISCREFRecord)
|
||||
ADD_PACK_PROTOTYPES(ViscrefTable)
|
||||
ADD_PACK_PROTOTYPES(WATDENTRecord)
|
||||
|
||||
@@ -41,6 +41,7 @@
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/OilVaporizationProperties.hpp>
|
||||
#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/SimulationConfig/SimulationConfig.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/SimulationConfig/ThresholdPressure.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Tables/Aqudims.hpp>
|
||||
@@ -243,6 +244,31 @@ Opm::VFPInjTable getVFPInjTable()
|
||||
return Opm::VFPInjTable(1, 2.0, Opm::VFPInjTable::FLO_WAT, {1.0, 2.0},
|
||||
{3.0, 4.0, 5.0}, table);
|
||||
}
|
||||
|
||||
|
||||
Opm::VFPProdTable getVFPProdTable()
|
||||
{
|
||||
Opm::VFPProdTable::array_type table;
|
||||
Opm::VFPProdTable::extents shape;
|
||||
shape[0] = 1;
|
||||
shape[1] = 2;
|
||||
shape[2] = 3;
|
||||
shape[3] = 4;
|
||||
shape[4] = 5;
|
||||
table.resize(shape);
|
||||
double foo = 1.0;
|
||||
for (size_t i = 0; i < table.num_elements(); ++i)
|
||||
*(table.data() + i) = foo++;
|
||||
return Opm::VFPProdTable(1, 2.0, Opm::VFPProdTable::FLO_OIL,
|
||||
Opm::VFPProdTable::WFR_WOR,
|
||||
Opm::VFPProdTable::GFR_GLR,
|
||||
Opm::VFPProdTable::ALQ_TGLR,
|
||||
{1.0, 2.0, 3.0, 4.0, 5.0},
|
||||
{1.0},
|
||||
{1.0, 2.0},
|
||||
{1.0, 2.0, 3.0},
|
||||
{1.0, 2.0, 3.0, 4.0}, table);
|
||||
}
|
||||
#endif
|
||||
|
||||
|
||||
@@ -1250,6 +1276,17 @@ BOOST_AUTO_TEST_CASE(VFPInjTable)
|
||||
}
|
||||
|
||||
|
||||
BOOST_AUTO_TEST_CASE(VFPProdTable)
|
||||
{
|
||||
#ifdef HAVE_MPI
|
||||
Opm::VFPProdTable val1 = getVFPProdTable();
|
||||
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