// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- // vi: set et ts=4 sw=4 sts=4: /* 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 2 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 . Consult the COPYING file in the top-level source directory of this module for the precise wording of the license and the list of copyright holders. */ #ifndef EWOMS_COLLECT_TO_IO_RANK_HH #define EWOMS_COLLECT_TO_IO_RANK_HH #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include namespace Dune { template class CartesianIndexMapper; } namespace Opm { template class CollectDataToIORank { public: using CollectiveCommunication = typename Grid::CollectiveCommunication; using P2PCommunicatorType = Dune::Point2PointCommunicator; using IndexMapType = std::vector; using IndexMapStorageType = std::vector; static constexpr int dimension = Grid::dimension; enum { ioRank = 0 }; static const bool needsReordering = !std::is_same::value; CollectDataToIORank(const Grid& grid, const EquilGrid* equilGrid, const GridView& gridView, const Dune::CartesianIndexMapper& cartMapper, const Dune::CartesianIndexMapper* equilCartMapper, const std::set& fipRegionsInterregFlow = {}); // gather solution to rank 0 for EclipseWriter void collect(const data::Solution& localCellData, const std::map, double>& localBlockData, const data::Wells& localWellData, const data::WellBlockAveragePressures& localWBPData, const data::GroupAndNetworkValues& localGroupAndNetworkData, const data::Aquifers& localAquiferData, const WellTestState& localWellTestState, const EclInterRegFlowMap& interRegFlows, const std::array, std::vector>>, 3>& localFlowsn, const std::array, std::vector>>, 3>& localFloresn); const std::map, double>& globalBlockData() const { return globalBlockData_; } const data::Solution& globalCellData() const { return globalCellData_; } const data::Wells& globalWellData() const { return globalWellData_; } const data::WellBlockAveragePressures& globalWBPData() const { return this->globalWBPData_; } const data::GroupAndNetworkValues& globalGroupAndNetworkData() const { return globalGroupAndNetworkData_; } const data::Aquifers& globalAquiferData() const { return globalAquiferData_; } const WellTestState& globalWellTestState() const { return this->globalWellTestState_; } EclInterRegFlowMap& globalInterRegFlows() { return this->globalInterRegFlows_; } const EclInterRegFlowMap& globalInterRegFlows() const { return this->globalInterRegFlows_; } const std::array, std::vector>>, 3>& globalFlowsn() const { return globalFlowsn_; } const std::array, std::vector>>, 3>& globalFloresn() const { return globalFloresn_; } bool isIORank() const { return toIORankComm_.rank() == ioRank; } bool isParallel() const { return toIORankComm_.size() > 1; } int localIdxToGlobalIdx(unsigned localIdx) const; const std::vector& localIdxToGlobalIdxMapping() const { return localIdxToGlobalIdx_; } bool doesNeedReordering() const { return needsReordering;} size_t numCells () const { return globalCartesianIndex_.size(); } const std::vector& globalRanks() const { return globalRanks_; } bool isCartIdxOnThisRank(int cartIdx) const; protected: P2PCommunicatorType toIORankComm_; EclInterRegFlowMap globalInterRegFlows_; IndexMapType globalCartesianIndex_; IndexMapType localIndexMap_; IndexMapStorageType indexMaps_; std::vector globalRanks_; data::Solution globalCellData_; std::map, double> globalBlockData_; data::Wells globalWellData_; data::WellBlockAveragePressures globalWBPData_; data::GroupAndNetworkValues globalGroupAndNetworkData_; data::Aquifers globalAquiferData_; WellTestState globalWellTestState_; std::vector localIdxToGlobalIdx_; std::array, std::vector>>, 3> globalFlowsn_; std::array, std::vector>>, 3> globalFloresn_; /// \brief sorted list of cartesian indices present- /// /// non-empty only when running in parallel std::vector sortedCartesianIdx_; }; } // end namespace Opm #endif