/* Copyright 2019 Equinor AS. This file is part of the Open Porous Media project (OPM). OPM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OPM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OPM. If not, see . */ #include #if HAVE_MPI #include #endif #include "ParallelRestart.hpp" #include #include #include #include #include #include #include #include #include #include #define HANDLE_AS_POD(T) \ std::size_t packSize(const T& data, Dune::MPIHelper::MPICommunicator comm) \ { \ return packSize(data, comm, std::integral_constant()); \ } \ void pack(const T& data, std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm) \ { \ pack(data, buffer, position, comm, std::integral_constant()); \ } \ void unpack(T& data, std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm) \ { \ unpack(data, buffer, position, comm, std::integral_constant()); \ } namespace { template std::pair, std::vector> splitDynState(const Opm::DynamicState& state) { std::vector unique; for (const auto& w : state.data()) { if (std::find(unique.begin(), unique.end(), w) == unique.end()) unique.push_back(w); } std::vector idxVec; idxVec.reserve(state.data().size()+1); for (const auto& w : state.data()) { auto uIt = std::find(unique.begin(), unique.end(), w); idxVec.push_back(uIt-unique.begin()); } idxVec.push_back(state.initialRange()); return std::make_pair(unique, idxVec); } template void reconstructDynState(const std::vector& unique, const std::vector& idxVec, Opm::DynamicState& result) { std::vector ptrData; for (size_t i = 0; i < idxVec.size()-1; ++i) { ptrData.push_back(unique[idxVec[i]]); } result = Opm::DynamicState(ptrData, idxVec.back()); } } namespace Opm { namespace Mpi { template std::size_t packSize(const T*, std::size_t, Dune::MPIHelper::MPICommunicator, std::integral_constant) { OPM_THROW(std::logic_error, "Packing not (yet) supported for this non-pod type."); } template std::size_t packSize(const T*, std::size_t l, Dune::MPIHelper::MPICommunicator comm, std::integral_constant) { #if HAVE_MPI int size; MPI_Pack_size(1, Dune::MPITraits::getType(), comm, &size); std::size_t totalSize = size; MPI_Pack_size(l, Dune::MPITraits::getType(), comm, &size); return totalSize + size; #else (void) comm; return l-l; #endif } template std::size_t packSize(const T* data, std::size_t l, Dune::MPIHelper::MPICommunicator comm) { return packSize(data, l, comm, typename std::is_pod::type()); } template std::size_t packSize(const std::pair& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.first, comm) + packSize(data.second, comm); } template std::size_t packSize(const std::vector& data, Dune::MPIHelper::MPICommunicator comm) { if (std::is_pod::value) // size written automatically return packSize(data.data(), data.size(), comm); std::size_t size = packSize(data.size(), comm); for (const auto& entry: data) size += packSize(entry, comm); return size; } template std::size_t packSize(const std::vector& data, Dune::MPIHelper::MPICommunicator comm) { bool entry; return packSize(data.size(), comm) + data.size()*packSize(entry,comm); } template typename std::enable_if::value, std::size_t>::type pack_size_tuple_entry(const Tuple&, Dune::MPIHelper::MPICommunicator) { return 0; } template typename std::enable_if::value, std::size_t>::type pack_size_tuple_entry(const Tuple& tuple, Dune::MPIHelper::MPICommunicator comm) { return packSize(std::get(tuple), comm) + pack_size_tuple_entry(tuple, comm); } template std::size_t packSize(const std::tuple& data, Dune::MPIHelper::MPICommunicator comm) { return pack_size_tuple_entry(data, comm); } template std::size_t packSize(const std::unordered_set& data, Dune::MPIHelper::MPICommunicator comm) { std::size_t totalSize = packSize(data.size(), comm); for (const auto& entry : data) { totalSize += packSize(entry, comm); } return totalSize; } template std::size_t packSize(const std::set& data, Dune::MPIHelper::MPICommunicator comm) { std::size_t totalSize = packSize(data.size(), comm); for (const auto& entry : data) { totalSize += packSize(entry, comm); } return totalSize; } template std::size_t packSize(const OrderedMap& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.getIndex(), comm) + packSize(data.getStorage(), comm); } template std::size_t packSize(const DynamicState& data, Dune::MPIHelper::MPICommunicator comm) { auto split = splitDynState(data); return packSize(split.first, comm) + packSize(split.second, comm); } std::size_t packSize(const char* str, Dune::MPIHelper::MPICommunicator comm) { #if HAVE_MPI int size; MPI_Pack_size(1, Dune::MPITraits::getType(), comm, &size); int totalSize = size; MPI_Pack_size(strlen(str)+1, MPI_CHAR, comm, &size); return totalSize + size; #else (void) str; (void) comm; return 0; #endif } std::size_t packSize(const std::string& str, Dune::MPIHelper::MPICommunicator comm) { return packSize(str.c_str(), comm); } template std::size_t packSize(const std::map& data, Dune::MPIHelper::MPICommunicator comm) { std::size_t totalSize = packSize(data.size(), comm); for (const auto& entry: data) { totalSize += packSize(entry, comm); } return totalSize; } template std::size_t packSize(const std::unordered_map& data, Dune::MPIHelper::MPICommunicator comm) { std::size_t totalSize = packSize(data.size(), comm); for (const auto& entry: data) { totalSize += packSize(entry, comm); } return totalSize; } template std::size_t packSize(const std::array& data, Dune::MPIHelper::MPICommunicator comm) { return N*packSize(data[0], comm); } HANDLE_AS_POD(data::Connection) HANDLE_AS_POD(data::CurrentControl) HANDLE_AS_POD(data::Rates) HANDLE_AS_POD(data::Segment) std::size_t packSize(const data::Well& data, Dune::MPIHelper::MPICommunicator comm) { std::size_t size = packSize(data.rates, comm); size += packSize(data.bhp, comm) + packSize(data.thp, comm); size += packSize(data.temperature, comm); size += packSize(data.control, comm); size += packSize(data.connections, comm); size += packSize(data.segments, comm); size += packSize(data.current_control, comm); return size; } std::size_t packSize(const data::CellData& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.dim, comm) + packSize(data.data, comm) + packSize(data.target, comm); } std::size_t packSize(const RestartKey& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.key, comm) + packSize(data.dim, comm) + packSize(data.required, comm); } std::size_t packSize(const data::Solution& data, Dune::MPIHelper::MPICommunicator comm) { // Needs explicit conversion to a supported base type holding the data // to prevent throwing. return packSize(static_cast&>(data), comm); } std::size_t packSize(const data::WellRates& data, Dune::MPIHelper::MPICommunicator comm) { // Needs explicit conversion to a supported base type holding the data // to prevent throwing. return packSize(static_cast&>(data), comm); } std::size_t packSize(const RestartValue& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.solution, comm) + packSize(data.wells, comm) + packSize(data.extra, comm); } std::size_t packSize(const WellType& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.producer(), comm) + packSize(data.preferred_phase(), comm); } template std::size_t packSize(const std::map& data, Dune::MPIHelper::MPICommunicator comm); std::size_t packSize(const VFPInjTable& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.getTableNum(), comm) + packSize(data.getDatumDepth(), comm) + packSize(data.getFloType(), comm) + packSize(data.getFloAxis(), comm) + packSize(data.getTHPAxis(), comm) + packSize(data.getTable(), 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) + packSize(data.getTable(), comm); } std::size_t packSize(const UDAValue& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.get_dim(), comm) + packSize(data.is(), comm) + (data.is() ? packSize(data.get(), comm) : packSize(data.get(), comm)); } template std::size_t packSize(const std::shared_ptr& data, Dune::MPIHelper::MPICommunicator comm) { std::size_t size = packSize(bool(), comm); if (data) size += packSize(*data, comm); return size; } template std::size_t packSize(const std::unique_ptr& data, Dune::MPIHelper::MPICommunicator comm) { std::size_t size = packSize(bool(), comm); if (data) size += packSize(*data, comm); return size; } std::size_t packSize(const Dimension& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.getSIScalingRaw(), comm) + packSize(data.getSIOffset(), comm); } std::size_t packSize(const UnitSystem& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.getName(), comm) + packSize(data.getType(), comm) + packSize(data.getDimensions(), comm) + packSize(data.use_count(), comm); } std::size_t packSize(const Well& data, Dune::MPIHelper::MPICommunicator comm) { std::size_t size = packSize(data.name(), comm) + packSize(data.groupName(), comm) + packSize(data.firstTimeStep(), comm) + packSize(data.seqIndex(), comm) + packSize(data.getHeadI(), comm) + packSize(data.getHeadJ(), comm) + packSize(data.getRefDepth(), comm) + packSize(data.wellType(), comm) + packSize(data.units(), comm) + packSize(data.udqUndefined(), comm) + packSize(data.getStatus(), comm) + packSize(data.getDrainageRadius(), comm) + packSize(data.getAllowCrossFlow(), comm) + packSize(data.getAutomaticShutIn(), comm) + packSize(data.wellGuideRate(), comm) + packSize(data.getEfficiencyFactor(), comm) + packSize(data.getSolventFraction(), comm) + packSize(data.predictionMode(), comm) + packSize(data.getEconLimits(), comm) + packSize(data.getFoamProperties(), comm) + packSize(data.getPolymerProperties(), comm) + packSize(data.getBrineProperties(), comm) + packSize(data.getTracerProperties(), comm) + packSize(data.getConnections(), comm) + packSize(data.getProductionProperties(), comm) + packSize(data.getInjectionProperties(), comm) + packSize(data.hasSegments(), comm); if (data.hasSegments()) size += packSize(data.getSegments(), comm); return size; } template std::size_t packSize(const IOrderSet& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.index(), comm) + packSize(data.data(), comm); } std::size_t packSize(const Group::GroupInjectionProperties& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.phase, comm) + packSize(data.cmode, comm) + packSize(data.surface_max_rate, comm) + packSize(data.resv_max_rate, comm) + packSize(data.target_reinj_fraction, comm) + packSize(data.target_void_fraction, comm) + packSize(data.reinj_group, comm) + packSize(data.voidage_group, comm) + packSize(data.injection_controls, comm); } std::size_t packSize(const Group::GroupProductionProperties& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.cmode, comm) + packSize(data.exceed_action, comm) + packSize(data.oil_target, comm) + packSize(data.water_target, comm) + packSize(data.gas_target, comm) + packSize(data.liquid_target, comm) + packSize(data.guide_rate, comm) + packSize(data.guide_rate_def, comm) + packSize(data.resv_target, comm) + packSize(data.production_controls, comm); } std::size_t packSize(const Group& data, Dune::MPIHelper::MPICommunicator comm) { return packSize(data.name(), comm) + packSize(data.insert_index(), comm) + packSize(data.initStep(), comm) + packSize(data.udqUndefined(), comm) + packSize(data.units(), comm) + packSize(data.type(), comm) + packSize(data.getGroupEfficiencyFactor(), comm) + packSize(data.getTransferGroupEfficiencyFactor(), comm) + packSize(data.isAvailableForGroupControl(), comm) + packSize(data.getGroupNetVFPTable(), comm) + packSize(data.parent(), comm) + packSize(data.iwells(), comm) + packSize(data.igroups(), comm) + packSize(data.injectionProperties(), comm) + packSize(data.productionProperties(), comm); } ////// pack routines template void pack(const T*, std::size_t, std::vector&, int&, Dune::MPIHelper::MPICommunicator, std::integral_constant) { OPM_THROW(std::logic_error, "Packing not (yet) supported for this non-pod type."); } template void pack(const T* data, std::size_t l, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm, std::integral_constant) { #if HAVE_MPI MPI_Pack(&l, 1, Dune::MPITraits::getType(), buffer.data(), buffer.size(), &position, comm); MPI_Pack(data, l, Dune::MPITraits::getType(), buffer.data(), buffer.size(), &position, comm); #else (void) data; (void) comm; (void) l; (void) buffer; (void) position; #endif } template void pack(const T* data, std::size_t l, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data, l, buffer, position, comm, typename std::is_pod::type()); } template void pack(const std::pair& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.first, buffer, position, comm); pack(data.second, buffer, position, comm); } template void pack(const std::vector& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { if (std::is_pod::value) { // size written automatically pack(data.data(), data.size(), buffer, position, comm); return; } pack(data.size(), buffer, position, comm); for (const auto& entry: data) pack(entry, buffer, position, comm); } template void pack(const std::set& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.size(), buffer, position, comm); for (const auto& entry : data) { pack(entry, buffer, position, comm); } } template void pack(const std::unordered_set& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.size(), buffer, position, comm); for (const auto& entry : data) { pack(entry, buffer, position, comm); } } template void pack(const std::array& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { for (const T& entry : data) pack(entry, buffer, position, comm); } template void pack(const std::vector& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.size(), buffer, position, comm); for (const auto& entry : data) { bool b = entry; pack(b, buffer, position, comm); } } template typename std::enable_if::value, void>::type pack_tuple_entry(const Tuple&, std::vector&, int&, Dune::MPIHelper::MPICommunicator) { } template typename std::enable_if::value, void>::type pack_tuple_entry(const Tuple& tuple, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(std::get(tuple), buffer, position, comm); pack_tuple_entry(tuple, buffer, position, comm); } template void pack(const std::tuple& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack_tuple_entry(data, buffer, position, comm); } template void pack(const OrderedMap& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.getIndex(), buffer, position, comm); pack(data.getStorage(), buffer, position, comm); } template void pack(const DynamicState& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { auto split = splitDynState(data); pack(split.first, buffer, position, comm); pack(split.second, buffer, position, comm); } void pack(const char* str, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { #if HAVE_MPI std::size_t length = strlen(str)+1; MPI_Pack(&length, 1, Dune::MPITraits::getType(), buffer.data(), buffer.size(), &position, comm); MPI_Pack(str, strlen(str)+1, MPI_CHAR, buffer.data(), buffer.size(), &position, comm); #else (void) str; (void) comm; (void) buffer; (void) position; #endif } void pack(const std::string& str, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(str.c_str(), buffer, position, comm); } template void pack(const std::map& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.size(), buffer, position, comm); for (const auto& entry: data) { pack(entry, buffer, position, comm); } } template void pack(const std::unordered_map& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.size(), buffer, position, comm); for (const auto& entry: data) { pack(entry, buffer, position, comm); } } template void pack(const std::map& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm); void pack(const data::Well& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.rates, buffer, position, comm); pack(data.bhp, buffer, position, comm); pack(data.thp, buffer, position, comm); pack(data.temperature, buffer, position, comm); pack(data.control, buffer, position, comm); pack(data.connections, buffer, position, comm); pack(data.segments, buffer, position, comm); pack(data.current_control, buffer, position, comm); } void pack(const RestartKey& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.key, buffer, position, comm); pack(data.dim, buffer, position, comm); pack(data.required, buffer, position, comm); } void pack(const data::CellData& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.dim, buffer, position, comm); pack(data.data, buffer, position, comm); pack(data.target, buffer, position, comm); } void pack(const data::Solution& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { // Needs explicit conversion to a supported base type holding the data // to prevent throwing. pack(static_cast&>(data), buffer, position, comm); } void pack(const data::WellRates& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { // Needs explicit conversion to a supported base type holding the data // to prevent throwing. pack(static_cast&>(data), buffer, position, comm); } void pack(const RestartValue& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.solution, buffer, position, comm); pack(data.wells, buffer, position, comm); pack(data.extra, buffer, position, comm); } void pack(const WellType& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.producer(), buffer, position, comm); pack(data.preferred_phase(), buffer, position, comm); } void pack(const VFPInjTable& data, std::vector& 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.getFloAxis(), buffer, position, comm); pack(data.getTHPAxis(), buffer, position, comm); pack(data.getTable(), buffer, position, comm); } void pack(const VFPProdTable& data, std::vector& 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); pack(data.getTable(), buffer, position, comm); } void pack(const UDAValue& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.get_dim(), buffer, position, comm); pack(data.is(), buffer, position, comm); if (data.is()) pack(data.get(), buffer, position, comm); else pack(data.get(), buffer, position, comm); } template void pack(const std::shared_ptr& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data != nullptr, buffer, position, comm); if (data) pack(*data, buffer, position, comm); } template void pack(const std::unique_ptr& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data != nullptr, buffer, position, comm); if (data) pack(*data, buffer, position, comm); } void pack(const Dimension& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.getSIScalingRaw(), buffer, position, comm); pack(data.getSIOffset(), buffer, position, comm); } void pack(const UnitSystem& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.getName(), buffer, position, comm); pack(data.getType(), buffer, position, comm); pack(data.getDimensions(), buffer, position, comm); pack(data.use_count(), buffer, position, comm); } void pack(const Well& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.name(), buffer, position, comm); pack(data.groupName(), buffer, position, comm); pack(data.firstTimeStep(), buffer, position, comm); pack(data.seqIndex(), buffer, position, comm); pack(data.getHeadI(), buffer, position, comm); pack(data.getHeadJ(), buffer, position, comm); pack(data.getRefDepth(), buffer, position, comm); pack(data.wellType(), buffer, position, comm); pack(data.units(), buffer, position, comm); pack(data.udqUndefined(), buffer, position, comm); pack(data.getStatus(), buffer, position, comm); pack(data.getDrainageRadius(), buffer, position, comm); pack(data.getAllowCrossFlow(), buffer, position, comm); pack(data.getAutomaticShutIn(), buffer, position, comm); pack(data.wellGuideRate(), buffer, position, comm); pack(data.getEfficiencyFactor(), buffer, position, comm); pack(data.getSolventFraction(), buffer, position, comm); pack(data.predictionMode(), buffer, position, comm); pack(data.getEconLimits(), buffer, position, comm); pack(data.getFoamProperties(), buffer, position, comm); pack(data.getPolymerProperties(), buffer, position, comm); pack(data.getBrineProperties(), buffer, position, comm); pack(data.getTracerProperties(), buffer, position, comm); pack(data.getConnections(), buffer, position, comm); pack(data.getProductionProperties(), buffer, position, comm); pack(data.getInjectionProperties(), buffer, position, comm); pack(data.hasSegments(), buffer, position, comm); if (data.hasSegments()) pack(data.getSegments(), buffer, position, comm); } template void pack(const IOrderSet& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.index(), buffer, position, comm); pack(data.data(), buffer, position, comm); } void pack(const Group::GroupInjectionProperties& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.phase, buffer, position, comm); pack(data.cmode, buffer, position, comm); pack(data.surface_max_rate, buffer, position, comm); pack(data.resv_max_rate, buffer, position, comm); pack(data.target_reinj_fraction, buffer, position, comm); pack(data.target_void_fraction, buffer, position, comm); pack(data.reinj_group, buffer, position, comm); pack(data.voidage_group, buffer, position, comm); pack(data.injection_controls, buffer, position, comm); } void pack(const Group::GroupProductionProperties& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.cmode, buffer, position, comm); pack(data.exceed_action, buffer, position, comm); pack(data.oil_target, buffer, position, comm); pack(data.water_target, buffer, position, comm); pack(data.gas_target, buffer, position, comm); pack(data.liquid_target, buffer, position, comm); pack(data.guide_rate, buffer, position, comm); pack(data.guide_rate_def, buffer, position, comm); pack(data.resv_target, buffer, position, comm); pack(data.production_controls, buffer, position, comm); } void pack(const Group& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { pack(data.name(), buffer, position, comm); pack(data.insert_index(), buffer, position, comm); pack(data.initStep(), buffer, position, comm); pack(data.udqUndefined(), buffer, position, comm); pack(data.units(), buffer, position, comm); pack(data.type(), buffer, position, comm); pack(data.getGroupEfficiencyFactor(), buffer, position, comm); pack(data.getTransferGroupEfficiencyFactor(), buffer, position, comm); pack(data.isAvailableForGroupControl(), buffer, position, comm); pack(data.getGroupNetVFPTable(), buffer, position, comm); pack(data.parent(), buffer, position, comm); pack(data.iwells(), buffer, position, comm); pack(data.igroups(), buffer, position, comm); pack(data.injectionProperties(), buffer, position, comm); pack(data.productionProperties(), buffer, position, comm); } /// unpack routines template void unpack(T*, const std::size_t&, std::vector&, int&, Dune::MPIHelper::MPICommunicator, std::integral_constant) { OPM_THROW(std::logic_error, "Packing not (yet) supported for this non-pod type."); } template void unpack(T* data, const std::size_t& l, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm, std::integral_constant) { #if HAVE_MPI MPI_Unpack(buffer.data(), buffer.size(), &position, data, l, Dune::MPITraits::getType(), comm); #else (void) data; (void) comm; (void) l; (void) buffer; (void) position; #endif } template void unpack(T* data, const std::size_t& l, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack(data, l, buffer, position, comm, typename std::is_pod::type()); } template void unpack(std::pair& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack(data.first, buffer, position, comm); unpack(data.second, buffer, position, comm); } template void unpack(std::vector& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::size_t length = 0; unpack(length, buffer, position, comm); data.resize(length); if (std::is_pod::value) { unpack(data.data(), data.size(), buffer, position, comm); return; } for (auto& entry: data) unpack(entry, buffer, position, comm); } template void unpack(std::vector& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { size_t size; unpack(size, buffer, position, comm); data.clear(); data.reserve(size); for (size_t i = 0; i < size; ++i) { bool entry; unpack(entry, buffer, position, comm); data.push_back(entry); } } template typename std::enable_if::value, void>::type unpack_tuple_entry(Tuple&, std::vector&, int&, Dune::MPIHelper::MPICommunicator) { } template typename std::enable_if::value, void>::type unpack_tuple_entry(Tuple& tuple, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack(std::get(tuple), buffer, position, comm); unpack_tuple_entry(tuple, buffer, position, comm); } template void unpack(std::tuple& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack_tuple_entry(data, buffer, position, comm); } template void unpack(std::set& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::size_t size = 0; unpack(size, buffer, position, comm); for (;size>0; size--) { K entry; unpack(entry, buffer, position, comm); data.insert(entry); } } template void unpack(std::unordered_set& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::size_t size=0; unpack(size, buffer, position, comm); for (;size>0; size--) { T entry; unpack(entry, buffer, position, comm); data.insert(entry); } } template void unpack(std::array& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { for (T& entry : data) unpack(entry, buffer, position, comm); } template void unpack(OrderedMap& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { typename OrderedMap::index_type index; typename OrderedMap::storage_type storage; unpack(index, buffer, position, comm); unpack(storage, buffer, position, comm); data = OrderedMap(index, storage); } template void unpack(DynamicState& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::vector unique; std::vector indices; Opm::Mpi::unpack(unique, buffer, position, comm); Opm::Mpi::unpack(indices, buffer, position, comm); reconstructDynState(unique, indices, data); } void unpack(char* str, std::size_t length, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { #if HAVE_MPI MPI_Unpack(buffer.data(), buffer.size(), &position, const_cast(str), length, MPI_CHAR, comm); #else (void) str; (void) comm; (void) length; (void) buffer; (void) position; #endif } void unpack(std::string& str, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::size_t length=0; unpack(length, buffer, position, comm); std::vector cStr(length, '\0'); unpack(cStr.data(), length, buffer, position, comm); str.clear(); str.append(cStr.data()); } template void unpack(std::map& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::size_t size=0; unpack(size, buffer, position, comm); for (;size>0; size--) { std::pair entry; unpack(entry, buffer, position, comm); data.insert(entry); } } template void unpack(std::unordered_map& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::size_t size=0; unpack(size, buffer, position, comm); for (;size>0; size--) { std::pair entry; unpack(entry, buffer, position, comm); data.insert(entry); } } void unpack(data::Well& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack(data.rates, buffer, position, comm); unpack(data.bhp, buffer, position, comm); unpack(data.thp, buffer, position, comm); unpack(data.temperature, buffer, position, comm); unpack(data.control, buffer, position, comm); unpack(data.connections, buffer, position, comm); unpack(data.segments, buffer, position, comm); unpack(data.current_control, buffer, position, comm); } void unpack(RestartKey& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack(data.key, buffer, position, comm); unpack(data.dim, buffer, position, comm); unpack(data.required, buffer, position, comm); } void unpack(data::CellData& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack(data.dim, buffer, position, comm); unpack(data.data, buffer, position, comm); unpack(data.target, buffer, position, comm); } void unpack(data::Solution& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { // Needs explicit conversion to a supported base type holding the data // to prevent throwing. unpack(static_cast&>(data), buffer, position, comm); } void unpack(data::WellRates& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { // Needs explicit conversion to a supported base type holding the data // to prevent throwing. unpack(static_cast&>(data), buffer, position, comm); } void unpack(RestartValue& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack(data.solution, buffer, position, comm); unpack(data.wells, buffer, position, comm); unpack(data.extra, buffer, position, comm); } void unpack(WellType& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { Phase preferred_phase; bool producer; unpack(producer, buffer, position, comm); unpack(preferred_phase, buffer, position, comm); data = WellType( producer, preferred_phase ); } void unpack(VFPInjTable& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { int tableNum; double datumDepth; VFPInjTable::FLO_TYPE floType; std::vector floAxis, thpAxis; VFPInjTable::array_type table; unpack(tableNum, buffer, position, comm); unpack(datumDepth, buffer, position, comm); unpack(floType, buffer, position, comm); unpack(floAxis, buffer, position, comm); unpack(thpAxis, buffer, position, comm); unpack(table, buffer, position, comm); data = VFPInjTable(tableNum, datumDepth, floType, floAxis, thpAxis, table); } void unpack(VFPProdTable& data, std::vector& 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 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); unpack(table, buffer, position, comm); data = VFPProdTable(tableNum, datumDepth, floType, wfrType, gfrType, alqType, floAxis, thpAxis, wfrAxis, gfrAxis, alqAxis, table); } void unpack(UDAValue& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { bool isDouble; Dimension dim; unpack(dim, buffer, position, comm); unpack(isDouble, buffer, position, comm); if (isDouble) { double val; unpack(val, buffer, position, comm); data = UDAValue(val, dim); } else { std::string val; unpack(val, buffer, position, comm); data = UDAValue(val, dim); } } template void unpack(std::shared_ptr& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { bool hasVal; unpack(hasVal, buffer, position, comm); if (hasVal) { data = std::make_shared(); unpack(*data, buffer, position, comm); } } template void unpack(std::unique_ptr& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { bool hasVal; unpack(hasVal, buffer, position, comm); if (hasVal) { data.reset(new T); unpack(*data, buffer, position, comm); } } void unpack(Dimension& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { double siScaling, siOffset; unpack(siScaling, buffer, position, comm); unpack(siOffset, buffer, position, comm); data = Dimension(siScaling, siOffset); } void unpack(UnitSystem& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::string name; UnitSystem::UnitType type; std::map dimensions; size_t use_count; unpack(name, buffer, position, comm); unpack(type, buffer, position, comm); unpack(dimensions, buffer, position, comm); unpack(use_count, buffer, position, comm); data = UnitSystem(name, type, dimensions, use_count); } void unpack(Well& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::string name, groupName; std::size_t firstTimeStep, seqIndex; int headI, headJ; double ref_depth; WellType wtype; UnitSystem units; double udq_undefined; Well::Status status; double drainageRadius; bool allowCrossFlow, automaticShutIn; Well::WellGuideRate guideRate; double efficiencyFactor; double solventFraction; bool prediction_mode; auto econLimits = std::make_shared(); auto foamProperties = std::make_shared(); auto polymerProperties = std::make_shared(); auto brineProperties = std::make_shared(); auto tracerProperties = std::make_shared(); auto connection = std::make_shared(); auto production = std::make_shared(); auto injection = std::make_shared(); std::shared_ptr segments; unpack(name, buffer, position, comm); unpack(groupName, buffer, position, comm); unpack(firstTimeStep, buffer, position, comm); unpack(seqIndex, buffer, position, comm); unpack(headI, buffer, position, comm); unpack(headJ, buffer, position, comm); unpack(ref_depth, buffer, position, comm); unpack(wtype, buffer, position, comm); unpack(units, buffer, position, comm); unpack(udq_undefined, buffer, position, comm); unpack(status, buffer, position, comm); unpack(drainageRadius, buffer, position, comm); unpack(allowCrossFlow, buffer, position, comm); unpack(automaticShutIn, buffer, position, comm); unpack(guideRate, buffer, position, comm); unpack(efficiencyFactor, buffer, position, comm); unpack(solventFraction, buffer, position, comm); unpack(prediction_mode, buffer, position, comm); unpack(*econLimits, buffer, position, comm); unpack(*foamProperties, buffer, position, comm); unpack(*polymerProperties, buffer, position, comm); unpack(*brineProperties, buffer, position, comm); unpack(*tracerProperties, buffer, position, comm); unpack(*connection, buffer, position, comm); unpack(*production, buffer, position, comm); unpack(*injection, buffer, position, comm); bool hasSegments; unpack(hasSegments, buffer, position, comm); if (hasSegments) { segments = std::make_shared(); unpack(*segments, buffer, position, comm); } data = Well(name, groupName, firstTimeStep, seqIndex, headI, headJ, ref_depth, wtype, units, udq_undefined, status, drainageRadius, allowCrossFlow, automaticShutIn, guideRate, efficiencyFactor, solventFraction, prediction_mode, econLimits, foamProperties, polymerProperties, brineProperties, tracerProperties, connection, production, injection, segments); } template void unpack(IOrderSet& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { typename IOrderSet::index_type index; typename IOrderSet::storage_type storage; unpack(index, buffer, position, comm); unpack(storage, buffer, position, comm); data = IOrderSet(index, storage); } template void unpack(std::map& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm); void unpack(Group::GroupInjectionProperties& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack(data.phase, buffer, position, comm); unpack(data.cmode, buffer, position, comm); unpack(data.surface_max_rate, buffer, position, comm); unpack(data.resv_max_rate, buffer, position, comm); unpack(data.target_reinj_fraction, buffer, position, comm); unpack(data.target_void_fraction, buffer, position, comm); unpack(data.reinj_group, buffer, position, comm); unpack(data.voidage_group, buffer, position, comm); unpack(data.injection_controls, buffer, position, comm); } void unpack(Group::GroupProductionProperties& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { unpack(data.cmode, buffer, position, comm); unpack(data.exceed_action, buffer, position, comm); unpack(data.oil_target, buffer, position, comm); unpack(data.water_target, buffer, position, comm); unpack(data.gas_target, buffer, position, comm); unpack(data.liquid_target, buffer, position, comm); unpack(data.guide_rate, buffer, position, comm); unpack(data.guide_rate_def, buffer, position, comm); unpack(data.resv_target, buffer, position, comm); unpack(data.production_controls, buffer, position, comm); } void unpack(Group& data, std::vector& buffer, int& position, Dune::MPIHelper::MPICommunicator comm) { std::string name; std::size_t insert_index, initStep; double udqUndefined; UnitSystem units; Group::GroupType type; double groupEfficiencyFactor; bool transferGroupEfficiencyFactor; bool availableForGroupControl; int groupNetVFPTable; std::string parent; IOrderSet wells, groups; std::map injection; Group::GroupProductionProperties production; unpack(name, buffer, position, comm); unpack(insert_index, buffer, position, comm); unpack(initStep, buffer, position, comm); unpack(udqUndefined, buffer, position, comm); unpack(units, buffer, position, comm); unpack(type, buffer, position, comm); unpack(groupEfficiencyFactor, buffer, position, comm); unpack(transferGroupEfficiencyFactor, buffer, position, comm); unpack(availableForGroupControl, buffer, position, comm); unpack(groupNetVFPTable, buffer, position, comm); unpack(parent, buffer, position, comm); unpack(wells, buffer, position, comm); unpack(groups, buffer, position, comm); unpack(injection, buffer, position, comm); unpack(production, buffer, position, comm); data = Group(name, insert_index, initStep, udqUndefined, units, type, groupEfficiencyFactor, transferGroupEfficiencyFactor, availableForGroupControl, groupNetVFPTable, parent, wells, groups, injection, production); } #define INSTANTIATE_PACK_VECTOR(...) \ template std::size_t packSize(const std::vector<__VA_ARGS__>& data, \ Dune::MPIHelper::MPICommunicator comm); \ template void pack(const std::vector<__VA_ARGS__>& data, \ std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm); \ template void unpack(std::vector<__VA_ARGS__>& data, \ std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm); INSTANTIATE_PACK_VECTOR(double) INSTANTIATE_PACK_VECTOR(std::vector) INSTANTIATE_PACK_VECTOR(bool) INSTANTIATE_PACK_VECTOR(char) INSTANTIATE_PACK_VECTOR(int) INSTANTIATE_PACK_VECTOR(size_t) INSTANTIATE_PACK_VECTOR(std::time_t) INSTANTIATE_PACK_VECTOR(std::array) INSTANTIATE_PACK_VECTOR(std::pair) INSTANTIATE_PACK_VECTOR(std::shared_ptr) INSTANTIATE_PACK_VECTOR(std::shared_ptr) INSTANTIATE_PACK_VECTOR(std::shared_ptr) INSTANTIATE_PACK_VECTOR(std::shared_ptr) INSTANTIATE_PACK_VECTOR(std::map) INSTANTIATE_PACK_VECTOR(std::pair>) INSTANTIATE_PACK_VECTOR(std::pair>) INSTANTIATE_PACK_VECTOR(std::pair>) INSTANTIATE_PACK_VECTOR(std::pair) INSTANTIATE_PACK_VECTOR(std::pair) INSTANTIATE_PACK_VECTOR(std::string) #undef INSTANTIATE_PACK_VECTOR #define INSTANTIATE_PACK_SET(...) \ template std::size_t packSize(const std::set<__VA_ARGS__>& data, \ Dune::MPIHelper::MPICommunicator comm); \ template void pack(const std::set<__VA_ARGS__>& data, \ std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm); \ template void unpack(std::set<__VA_ARGS__>& data, \ std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm); INSTANTIATE_PACK_SET(std::string) #undef INSTANTIATE_PACK_SET #define INSTANTIATE_PACK_SHARED_PTR(...) \ template std::size_t packSize(const std::shared_ptr<__VA_ARGS__>& data, \ Dune::MPIHelper::MPICommunicator comm); \ template void pack(const std::shared_ptr<__VA_ARGS__>& data, \ std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm); \ template void unpack(std::shared_ptr<__VA_ARGS__>& data, \ std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm); INSTANTIATE_PACK_SHARED_PTR(VFPInjTable) INSTANTIATE_PACK_SHARED_PTR(Well) #undef INSTANTIATE_PACK_SHARED_PTR #define INSTANTIATE_PACK(...) \ template std::size_t packSize(const __VA_ARGS__& data, \ Dune::MPIHelper::MPICommunicator comm); \ template void pack(const __VA_ARGS__& data, \ std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm); \ template void unpack(__VA_ARGS__& data, \ std::vector& buffer, int& position, \ Dune::MPIHelper::MPICommunicator comm); INSTANTIATE_PACK(double) INSTANTIATE_PACK(std::size_t) INSTANTIATE_PACK(bool) INSTANTIATE_PACK(int) INSTANTIATE_PACK(std::array) INSTANTIATE_PACK(std::array) INSTANTIATE_PACK(unsigned char) INSTANTIATE_PACK(std::map,std::pair>) INSTANTIATE_PACK(std::map) INSTANTIATE_PACK(std::map>) INSTANTIATE_PACK(std::map>) INSTANTIATE_PACK(std::map,int>>) INSTANTIATE_PACK(std::map) INSTANTIATE_PACK(std::map) INSTANTIATE_PACK(std::unordered_map) INSTANTIATE_PACK(std::unordered_map) INSTANTIATE_PACK(std::unordered_set) INSTANTIATE_PACK(std::pair) INSTANTIATE_PACK(std::pair) #undef INSTANTIATE_PACK } // end namespace Mpi RestartValue loadParallelRestart(const EclipseIO* eclIO, SummaryState& summaryState, const std::vector& solutionKeys, const std::vector& extraKeys, Dune::CollectiveCommunication comm) { #if HAVE_MPI data::Solution sol; data::Wells wells; RestartValue restartValues(sol, wells); if (eclIO) { assert(comm.rank() == 0); restartValues = eclIO->loadRestart(summaryState, solutionKeys, extraKeys); int packedSize = Mpi::packSize(restartValues, comm); std::vector buffer(packedSize); int position=0; Mpi::pack(restartValues, buffer, position, comm); comm.broadcast(&position, 1, 0); comm.broadcast(buffer.data(), position, 0); } else { int bufferSize{}; comm.broadcast(&bufferSize, 1, 0); std::vector buffer(bufferSize); comm.broadcast(buffer.data(), bufferSize, 0); int position{}; Mpi::unpack(restartValues, buffer, position, comm); } return restartValues; #else (void) comm; return eclIO->loadRestart(summaryState, solutionKeys, extraKeys); #endif } } // end namespace Opm