/* Copyright 2016 IRIS AS Copyright 2019, 2020 Equinor ASA Copyright 2020 SINTEF Digital, Mathematics and Cybernetics 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 #include #include #include #if COMPILE_BDA_BRIDGE #include #include #endif #if HAVE_DUNE_ALUGRID #include #include #endif // HAVE_DUNE_ALUGRID #include namespace Opm { namespace detail { #ifdef HAVE_MPI void copyParValues(std::any& parallelInformation, size_t size, Dune::OwnerOverlapCopyCommunication& comm) { if (parallelInformation.type() == typeid(ParallelISTLInformation)) { const ParallelISTLInformation* parinfo = std::any_cast(¶llelInformation); assert(parinfo); parinfo->copyValuesTo(comm.indexSet(), comm.remoteIndices(), size, 1); } } #endif template void makeOverlapRowsInvalid(Matrix& matrix, const std::vector& overlapRows) { //value to set on diagonal const int numEq = Matrix::block_type::rows; typename Matrix::block_type diag_block(0.0); for (int eq = 0; eq < numEq; ++eq) diag_block[eq][eq] = 1.0; //loop over precalculated overlap rows and columns for (const auto row : overlapRows) { // Zero out row. matrix[row] = 0.0; //diagonal block set to diag(1.0). matrix[row][row] = diag_block; } } /// Return an appropriate weight function if a cpr preconditioner is asked for. template std::function getWeightsCalculator(const PropertyTree& prm, const Matrix& matrix, size_t pressureIndex, std::function trueFunc) { std::function weightsCalculator; using namespace std::string_literals; auto preconditionerType = prm.get("preconditioner.type"s, "cpr"s); if (preconditionerType == "cpr" || preconditionerType == "cprt" || preconditionerType == "cprw" || preconditionerType == "cprwt") { const bool transpose = preconditionerType == "cprt" || preconditionerType == "cprwt"; const auto weightsType = prm.get("preconditioner.weight_type"s, "quasiimpes"s); if (weightsType == "quasiimpes") { // weights will be created as default in the solver // assignment p = pressureIndex prevent compiler warning about // capturing variable with non-automatic storage duration weightsCalculator = [matrix, transpose, pressureIndex]() { return Amg::getQuasiImpesWeights(matrix, pressureIndex, transpose); }; } else if (weightsType == "trueimpes") { weightsCalculator = trueFunc; } else { OPM_THROW(std::invalid_argument, "Weights type " + weightsType + "not implemented for cpr." " Please use quasiimpes or trueimpes."); } } return weightsCalculator; } template void FlexibleSolverInfo::create(const Matrix& matrix, bool parallel, const PropertyTree& prm, size_t pressureIndex, std::function trueFunc, [[maybe_unused]] Comm& comm) { // Write sizes of linear systems on all ranks to debug log. { #if HAVE_MPI auto basic_comm = comm.communicator(); #else auto basic_comm = Dune::Communication{}; #endif // HAVE_MPI std::ostringstream os; os << "Linear system "; if (basic_comm.size() > 1) { os << fmt::format("on MPI rank: {} ", basic_comm.rank()); } os << fmt::format("blocksize: {} size: {:7d} block nonzeroes: {:9d}", Matrix::block_type::rows, matrix.N(), matrix.nonzeroes()); DeferredLogger local_logger; local_logger.debug(os.str()); auto global_logger = gatherDeferredLogger(local_logger, basic_comm); if (basic_comm.rank() == 0) { global_logger.logMessages(); } } std::function weightsCalculator = getWeightsCalculator(prm, matrix, pressureIndex, trueFunc); if (parallel) { #if HAVE_MPI if (!wellOperator_) { using ParOperatorType = Dune::OverlappingSchwarzOperator; auto pop = std::make_unique(matrix, comm); using FlexibleSolverType = Dune::FlexibleSolver; auto sol = std::make_unique(*pop, comm, prm, weightsCalculator, pressureIndex); this->pre_ = &sol->preconditioner(); this->op_ = std::move(pop); this->solver_ = std::move(sol); } else { using ParOperatorType = WellModelGhostLastMatrixAdapter; auto pop = std::make_unique(matrix, *wellOperator_, interiorCellNum_); using FlexibleSolverType = Dune::FlexibleSolver; auto sol = std::make_unique(*pop, comm, prm, weightsCalculator, pressureIndex); this->pre_ = &sol->preconditioner(); this->op_ = std::move(pop); this->solver_ = std::move(sol); } #endif } else { if (!wellOperator_) { using SeqOperatorType = Dune::MatrixAdapter; auto sop = std::make_unique(matrix); using FlexibleSolverType = Dune::FlexibleSolver; auto sol = std::make_unique(*sop, prm, weightsCalculator, pressureIndex); this->pre_ = &sol->preconditioner(); this->op_ = std::move(sop); this->solver_ = std::move(sol); } else { using SeqOperatorType = WellModelMatrixAdapter; auto sop = std::make_unique(matrix, *wellOperator_); using FlexibleSolverType = Dune::FlexibleSolver; auto sol = std::make_unique(*sop, prm, weightsCalculator, pressureIndex); this->pre_ = &sol->preconditioner(); this->op_ = std::move(sop); this->solver_ = std::move(sol); } } } template using BM = Dune::BCRSMatrix>; template using BV = Dune::BlockVector>; #if HAVE_MPI using CommunicationType = Dune::OwnerOverlapCopyCommunication; #else using CommunicationType = Dune::CollectiveCommunication; #endif #define INSTANCE_FLEX(Dim) \ template void makeOverlapRowsInvalid>(BM&, const std::vector&); \ template struct FlexibleSolverInfo,BV,CommunicationType>; INSTANCE_FLEX(1) INSTANCE_FLEX(2) INSTANCE_FLEX(3) INSTANCE_FLEX(4) INSTANCE_FLEX(5) INSTANCE_FLEX(6) } }