Support Distributing Guiderates in Parallel Restart

In preparation of distributing group-related information (e.g.,
the active constraint).
This commit is contained in:
Bård Skaflestad 2020-08-28 08:31:52 +02:00
parent 399ff94bce
commit a0373ed428
3 changed files with 82 additions and 0 deletions

View File

@ -209,6 +209,14 @@ HANDLE_AS_POD(data::CurrentControl)
HANDLE_AS_POD(data::Rates)
HANDLE_AS_POD(data::Segment)
std::size_t packSize(const data::GuideRateValue&, Dune::MPIHelper::MPICommunicator comm)
{
const auto nItem = static_cast<std::size_t>(data::GuideRateValue::Item::NumItems);
return packSize(std::array<int , nItem>{}, comm)
+ packSize(std::array<double, nItem>{}, comm);
}
std::size_t packSize(const data::Well& data, Dune::MPIHelper::MPICommunicator comm)
{
std::size_t size = packSize(data.rates, comm);
@ -218,6 +226,7 @@ std::size_t packSize(const data::Well& data, Dune::MPIHelper::MPICommunicator co
size += packSize(data.connections, comm);
size += packSize(data.segments, comm);
size += packSize(data.current_control, comm);
size += packSize(data.guide_rates, comm);
return size;
}
@ -436,6 +445,28 @@ void pack(const std::unordered_map<T1,T2,H,P,A>& data, std::vector<char>& buffer
}
}
void pack(const data::GuideRateValue& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
using Item = data::GuideRateValue::Item;
const auto nItem = static_cast<std::size_t>(Item::NumItems);
auto has = std::array<int , nItem>{}; has.fill(0);
auto val = std::array<double, nItem>{}; val.fill(0.0);
for (auto itemID = 0*nItem; itemID < nItem; ++itemID) {
const auto item = static_cast<Item>(itemID);
if (data.has(item)) {
has[itemID] = 1;
val[itemID] = data.get(item);
}
}
pack(has, buffer, position, comm);
pack(val, buffer, position, comm);
}
void pack(const data::Well& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
@ -447,6 +478,7 @@ void pack(const data::Well& data, std::vector<char>& buffer, int& position,
pack(data.connections, buffer, position, comm);
pack(data.segments, buffer, position, comm);
pack(data.current_control, buffer, position, comm);
pack(data.guide_rates, buffer, position, comm);
}
void pack(const RestartKey& data, std::vector<char>& buffer, int& position,
@ -707,6 +739,26 @@ void unpack(data::Well& data, std::vector<char>& buffer, int& position,
unpack(data.connections, buffer, position, comm);
unpack(data.segments, buffer, position, comm);
unpack(data.current_control, buffer, position, comm);
unpack(data.guide_rates, buffer, position, comm);
}
void unpack(data::GuideRateValue& data, std::vector<char>& buffer, int& position,
Dune::MPIHelper::MPICommunicator comm)
{
using Item = data::GuideRateValue::Item;
const auto nItem = static_cast<std::size_t>(Item::NumItems);
auto has = std::array<int , nItem>{};
auto val = std::array<double, nItem>{};
unpack(has, buffer, position, comm);
unpack(val, buffer, position, comm);
for (auto itemID = 0*nItem; itemID < nItem; ++itemID) {
if (has[itemID] != 0) {
data.set(static_cast<Item>(itemID), val[itemID]);
}
}
}
void unpack(RestartKey& data, std::vector<char>& buffer, int& position,

View File

@ -306,6 +306,7 @@ ADD_PACK_PROTOTYPES(data::CurrentControl)
ADD_PACK_PROTOTYPES(data::Rates)
ADD_PACK_PROTOTYPES(data::Segment)
ADD_PACK_PROTOTYPES(data::Solution)
ADD_PACK_PROTOTYPES(data::GuideRateValue)
ADD_PACK_PROTOTYPES(data::Well)
ADD_PACK_PROTOTYPES(data::WellRates)
ADD_PACK_PROTOTYPES(RestartKey)

View File

@ -24,6 +24,8 @@
#include <boost/test/unit_test.hpp>
#include <tuple>
#include <opm/common/OpmLog/Location.hpp>
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <opm/parser/eclipse/Deck/DeckItem.hpp>
@ -187,6 +189,15 @@ Opm::data::CurrentControl getCurrentControl()
return curr;
}
Opm::data::GuideRateValue getWellGuideRate()
{
using Item = Opm::data::GuideRateValue::Item;
return Opm::data::GuideRateValue{}.set(Item::Oil , 1.23)
.set(Item::Gas , 2.34)
.set(Item::Water, 3.45)
.set(Item::ResV , 4.56);
}
Opm::data::Well getWell()
{
@ -199,6 +210,7 @@ Opm::data::Well getWell()
well1.connections.push_back(getConnection());
well1.segments.insert({0, getSegment()});
well1.current_control = getCurrentControl();
well1.guide_rates = getWellGuideRate();
return well1;
}
@ -256,6 +268,23 @@ BOOST_AUTO_TEST_CASE(Rates)
DO_CHECKS(data::Rates)
}
BOOST_AUTO_TEST_CASE(dataGuideRateValue)
{
using Item = Opm::data::GuideRateValue::Item;
const auto val1 = Opm::data::GuideRateValue{}
.set(Item::Oil , 999.888)
.set(Item::Gas , 8888.777)
.set(Item::ResV, 12345.678);
const auto val2 = PackUnpack(val1);
BOOST_CHECK_MESSAGE(! std::get<0>(val2).has(Item::Water),
"Water Must Not Appear From "
"Serializing GuideRateValues");
DO_CHECKS(data::GuideRateValue)
}
BOOST_AUTO_TEST_CASE(dataConnection)
{