// -*- 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. */ /*! * \file * \copydoc Opm::EclCpGridVanguard */ #ifndef EWOMS_ECL_CP_GRID_VANGUARD_HH #define EWOMS_ECL_CP_GRID_VANGUARD_HH #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #if HAVE_MPI namespace Opm { namespace details { class MPIPartitionFromFile { public: explicit MPIPartitionFromFile(const std::filesystem::path& partitionFile) : partitionFile_(partitionFile) {} std::vector operator()(const Dune::CpGrid& grid) const; private: std::filesystem::path partitionFile_{}; }; inline std::vector MPIPartitionFromFile::operator()(const Dune::CpGrid& grid) const { std::ifstream pfile { this->partitionFile_ }; auto partition = std::vector { std::istream_iterator { pfile }, std::istream_iterator {} }; const auto nc = static_cast::size_type>(grid.size(0)); if (partition.size() == nc) { // Input is one process ID for each active cell return partition; } if (partition.size() == 3 * nc) { // Partition file is of the form // // Process_ID Cartesian_Idx NLDD_Domain // // with one row for each active cell. Select first column. auto g2l = std::unordered_map{}; auto locCell = 0; for (const auto& globCell : grid.globalCell()) { g2l.insert_or_assign(globCell, locCell++); } auto filtered_partition = std::vector(nc); for (auto c = 0*nc; c < nc; ++c) { auto pos = g2l.find(partition[3*c + 1]); if (pos != g2l.end()) { filtered_partition[pos->second] = partition[3*c + 0]; } } return filtered_partition; } throw std::invalid_argument { fmt::format("Partition file '{}' with {} values does " "not match CpGrid instance with {} cells", this->partitionFile_.generic_string(), partition.size(), nc) }; } }} // namespace Opm::details #endif // HAVE_MPI namespace Opm { template class EclCpGridVanguard; } namespace Opm::Properties { namespace TTag { struct EclCpGridVanguard { using InheritsFrom = std::tuple; }; } // declare the properties template struct Vanguard { using type = EclCpGridVanguard; }; template struct Grid { using type = Dune::CpGrid; }; template struct EquilGrid { using type = GetPropType; }; } // namespace Opm::Properties namespace Opm { /*! * \ingroup EclBlackOilSimulator * * \brief Helper class for grid instantiation of ECL file-format using problems. * * This class uses Dune::CpGrid as the simulation grid. */ template class EclCpGridVanguard : public EclBaseVanguard , public EclGenericCpGridVanguard, GetPropType, GetPropType> { friend class EclBaseVanguard; using ParentType = EclBaseVanguard; using Scalar = GetPropType; using Simulator = GetPropType; using ElementMapper = GetPropType; public: using Grid = GetPropType; using CartesianIndexMapper = Dune::CartesianIndexMapper; using EquilGrid = GetPropType; using GridView = GetPropType; using TransmissibilityType = EclTransmissibility; static constexpr int dimensionworld = Grid::dimensionworld; using Indices = GetPropType; static constexpr bool waterEnabled = Indices::waterEnabled; static constexpr bool gasEnabled = Indices::gasEnabled; static constexpr bool oilEnabled = Indices::oilEnabled; private: using Element = typename GridView::template Codim<0>::Entity; public: EclCpGridVanguard(Simulator& simulator) : EclBaseVanguard(simulator) { this->checkConsistency(); this->callImplementationInit(); } /*! * Checking consistency of simulator */ void checkConsistency() { const auto& runspec = this->eclState().runspec(); const auto& config = this->eclState().getSimulationConfig(); const auto& phases = runspec.phases(); // check for correct module setup if (config.isThermal()) { if (getPropValue() == false) { throw std::runtime_error("Input specifies energy while simulator has disabled it, try xxx_energy"); } } else { if (getPropValue() == true) { throw std::runtime_error("Input specifies no energy while simulator has energy, try run without _energy"); } } if (config.isDiffusive()) { if (getPropValue() == false) { throw std::runtime_error("Input specifies diffusion while simulator has disabled it, try xxx_diffusion"); } } if (runspec.micp()) { if (getPropValue() == false) { throw std::runtime_error("Input specifies MICP while simulator has it disabled"); } } if (phases.active(Phase::BRINE)) { if (getPropValue() == false) { throw std::runtime_error("Input specifies Brine while simulator has it disabled"); } } if (phases.active(Phase::POLYMER)) { if (getPropValue() == false) { throw std::runtime_error("Input specifies Polymer while simulator has it disabled"); } } // checking for correct phases is more difficult TODO! if (phases.active(Phase::ZFRACTION)) { if (getPropValue() == false) { throw std::runtime_error("Input specifies ExBo while simulator has it disabled"); } } if (phases.active(Phase::FOAM)) { if (getPropValue() == false) { throw std::runtime_error("Input specifies Foam while simulator has it disabled"); } } if (phases.active(Phase::SOLVENT)) { if (getPropValue() == false) { throw std::runtime_error("Input specifies Solvent while simulator has it disabled"); } } if(phases.active(Phase::WATER)){ if(waterEnabled == false){ throw std::runtime_error("Input specifies water while simulator has it disabled"); } } if(phases.active(Phase::GAS)){ if(gasEnabled == false){ throw std::runtime_error("Input specifies gas while simulator has it disabled"); } } if(phases.active(Phase::OIL)){ if(oilEnabled == false){ throw std::runtime_error("Input specifies oil while simulator has it disabled"); } } } /*! * \brief Free the memory occupied by the global transmissibility object. * * After writing the initial solution, this array should not be necessary anymore. */ void releaseGlobalTransmissibilities() { globalTrans_.reset(); } const TransmissibilityType& globalTransmissibility() const { assert( globalTrans_ != nullptr ); return *globalTrans_; } void releaseGlobalTransmissibility() { globalTrans_.reset(); } /*! * \brief Distribute the simulation grid over multiple processes * * (For parallel simulation runs.) */ void loadBalance() { #if HAVE_MPI if (const auto& extPFile = this->externalPartitionFile(); !extPFile.empty() && (extPFile != "none")) { this->setExternalLoadBalancer(details::MPIPartitionFromFile { extPFile }); } this->doLoadBalance_(this->edgeWeightsMethod(), this->ownersFirst(), this->serialPartitioning(), this->enableDistributedWells(), this->zoltanImbalanceTol(), this->gridView(), this->schedule(), this->eclState(), this->parallelWells_, this->numJacobiBlocks()); #endif this->updateGridView_(); this->updateCartesianToCompressedMapping_(); this->updateCellDepths_(); this->updateCellThickness_(); #if HAVE_MPI this->distributeFieldProps_(this->eclState()); #endif } unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const { return elemIndex; } unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const { return elemIndex; } /*! * \brief Get function to query cell centroids for a distributed grid. * * Currently this only non-empty for a loadbalanced CpGrid. * It is a function return the centroid for the given element * index. */ std::function(int)> cellCentroids() const { return this->cellCentroids_(this->cartesianIndexMapper(), true); } const std::vector& globalCell() { return this->grid().globalCell(); } protected: void createGrids_() { this->doCreateGrids_(this->eclState()); } void allocTrans() override { OPM_TIMEBLOCK(allocateTrans); globalTrans_.reset(new TransmissibilityType(this->eclState(), this->gridView(), this->cartesianIndexMapper(), this->grid(), this->cellCentroids(), getPropValue(), getPropValue(), getPropValue())); globalTrans_->update(false); } double getTransmissibility(unsigned I, unsigned J) const override { return globalTrans_->transmissibility(I,J); } #if HAVE_MPI const std::string& zoltanParams() const override { return this->zoltanParams_; } #endif // removing some connection located in inactive grid cells void filterConnections_() { this->doFilterConnections_(this->schedule()); } std::unique_ptr globalTrans_; }; } // namespace Opm #endif