/* Copyright 2021 Total SE 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 #include #include #include #include #include #if HAVE_DUNE_ALUGRID #include #endif // HAVE_DUNE_ALUGRID #include #include #include #include #include #include #include #include #include namespace { std::pair, int> countDomains(std::vector partition_vector) { auto maxPos = std::max_element(partition_vector.begin(), partition_vector.end()); const auto num_domains = (maxPos == partition_vector.end()) ? 0 : *maxPos + 1; return { std::move(partition_vector), num_domains }; } template struct HasZoltanPartitioning : public std::false_type {}; template struct HasZoltanPartitioning< GridType, std::void_t().zoltanPartitionWithoutScatter (std::declval*>(), std::declval(), std::declval(), std::declval()))> > : public std::true_type {}; } // anonymous namespace namespace Opm { template std::pair, int> partitionCells(const Grid& grid, const std::vector& wells, const std::string& method, const int num_local_domains, const double partition_imbalance) { if (method == "zoltan") { if constexpr (HasZoltanPartitioning::value) { return partitionCellsZoltan(grid, wells, num_local_domains, partition_imbalance); } else { OPM_THROW(std::runtime_error, "Zoltan requested for local domain partitioning, " "but is not available for the current grid type."); } } else if (method == "simple") { const int num_cells = detail::countLocalInteriorCells(grid); return partitionCellsSimple(num_cells, num_local_domains); } else if (method.size() > 10 && method.substr(method.size() - 10, 10) == ".partition") { // Method string ends with ".partition", interpret as filename for partitioning. const int num_cells = detail::countLocalInteriorCells(grid); return partitionCellsFromFile(method, num_cells); } else { OPM_THROW(std::runtime_error, "Unknown local domain partitioning method requested: " + method); } } // Read from file, containing one number per cell. std::pair, int> partitionCellsFromFile(const std::string& filename, const int num_cells) { // Read file into single vector. std::ifstream is(filename); const std::vector cellpart{std::istream_iterator(is), std::istream_iterator()}; if (cellpart.size() != size_t(num_cells)) { auto msg = fmt::format("Partition file contains {} entries, but there are {} cells.", cellpart.size(), num_cells); throw std::runtime_error(msg); } // Create and return the output domain vector. const int num_domains = (*std::max_element(cellpart.begin(), cellpart.end())) + 1; return { cellpart, num_domains }; } // Trivially simple partitioner std::pair, int> partitionCellsSimple(const int num_cells, const int num_domains) { // Build the partitions. const int dom_sz = num_cells / num_domains; std::vector bounds(num_domains + 1, dom_sz); bounds[0] = 0; for (int i = 0; i < num_cells % num_domains; ++i) { ++bounds[i + 1]; } std::partial_sum(bounds.begin(), bounds.end(), bounds.begin()); assert(bounds[num_domains] == num_cells); std::vector part(num_cells); for (int i = 0; i < num_domains; ++i) { std::fill(part.begin() + bounds[i], part.begin() + bounds[i + 1], i); } return { part, num_domains }; } template std::pair, int> partitionCellsZoltan(const Grid& grid, const std::vector& wells, const int num_domains, const double domain_imbalance) { auto partition_vector = grid.zoltanPartitionWithoutScatter (&wells, nullptr, num_domains, domain_imbalance); return countDomains(std::move(partition_vector)); } template std::pair,int> partitionCells(const Dune::CpGrid&, const std::vector&, const std::string&, const int, const double); template std::pair,int> partitionCells>(const Dune::PolyhedralGrid<3,3,double>&, const std::vector&, const std::string&, const int, const double); #if HAVE_DUNE_ALUGRID #if HAVE_MPI using ALUGrid3CN = Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridMPIComm>; #else using ALUGrid3CN = Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridNoComm>; #endif //HAVE_MPI template std::pair,int> partitionCells(const ALUGrid3CN&, const std::vector&, const std::string&, const int, const double); #endif } // namespace Opm