mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
changed: move all partitionCells code to partitionCells.cpp
This commit is contained in:
@@ -17,8 +17,20 @@
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#include <config.h>
|
||||
#include <opm/simulators/flow/partitionCells.hpp>
|
||||
|
||||
#include <opm/grid/CpGrid.hpp>
|
||||
#include <opm/grid/polyhedralgrid.hh>
|
||||
|
||||
#include <opm/input/eclipse/Schedule/Well/Well.hpp>
|
||||
|
||||
#include <opm/simulators/flow/countGlobalCells.hpp>
|
||||
|
||||
#if HAVE_DUNE_ALUGRID
|
||||
#include <dune/alugrid/grid.hh>
|
||||
#endif // HAVE_DUNE_ALUGRID
|
||||
|
||||
#include <fmt/format.h>
|
||||
|
||||
#include <algorithm>
|
||||
@@ -28,9 +40,64 @@
|
||||
#include <numeric>
|
||||
#include <stdexcept>
|
||||
#include <string>
|
||||
#include <type_traits>
|
||||
|
||||
namespace Opm
|
||||
{
|
||||
namespace {
|
||||
|
||||
std::pair<std::vector<int>, int>
|
||||
countDomains(std::vector<int> 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 <typename, class = void>
|
||||
struct HasZoltanPartitioning : public std::false_type {};
|
||||
|
||||
template <typename GridType>
|
||||
struct HasZoltanPartitioning<
|
||||
GridType,
|
||||
std::void_t<decltype(std::declval<const GridType&>().zoltanPartitionWithoutScatter
|
||||
(std::declval<const std::vector<Opm::Well>*>(),
|
||||
std::declval<const double*>(),
|
||||
std::declval<const int>(),
|
||||
std::declval<const double>()))>
|
||||
> : public std::true_type {};
|
||||
|
||||
} // anonymous namespace
|
||||
|
||||
namespace Opm {
|
||||
|
||||
template<class Grid>
|
||||
std::pair<std::vector<int>, int> partitionCells(const Grid& grid,
|
||||
const std::vector<Well>& wells,
|
||||
const std::string& method,
|
||||
const int num_local_domains,
|
||||
const double partition_imbalance)
|
||||
{
|
||||
if (method == "zoltan") {
|
||||
if constexpr (HasZoltanPartitioning<Grid>::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.
|
||||
@@ -71,4 +138,46 @@ namespace Opm
|
||||
}
|
||||
|
||||
|
||||
template<class Grid>
|
||||
std::pair<std::vector<int>, int> partitionCellsZoltan(const Grid& grid,
|
||||
const std::vector<Well>& 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<std::vector<int>,int>
|
||||
partitionCells<Dune::CpGrid>(const Dune::CpGrid&,
|
||||
const std::vector<Well>&,
|
||||
const std::string&,
|
||||
const int,
|
||||
const double);
|
||||
|
||||
template std::pair<std::vector<int>,int>
|
||||
partitionCells<Dune::PolyhedralGrid<3,3,double>>(const Dune::PolyhedralGrid<3,3,double>&,
|
||||
const std::vector<Well>&,
|
||||
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<std::vector<int>,int>
|
||||
partitionCells<ALUGrid3CN>(const ALUGrid3CN&,
|
||||
const std::vector<Well>&,
|
||||
const std::string&,
|
||||
const int,
|
||||
const double);
|
||||
#endif
|
||||
|
||||
|
||||
} // namespace Opm
|
||||
|
||||
Reference in New Issue
Block a user