// -*- 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::Linear::ParallelBaseBackend */ #ifndef EWOMS_PARALLEL_BASE_BACKEND_HH #define EWOMS_PARALLEL_BASE_BACKEND_HH #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include namespace Opm::Properties { namespace TTag { struct ParallelBaseLinearSolver {}; } //! Set the type of a global jacobian matrix for linear solvers that are based on //! dune-istl. template struct SparseMatrixAdapter { private: using Scalar = GetPropType; enum { numEq = getPropValue() }; using Block = Opm::MatrixBlock; public: using type = typename Opm::Linear::IstlSparseMatrixAdapter; }; } // namespace Opm::Properties namespace Opm { namespace Linear { /*! * \ingroup Linear * * \brief Provides the common code which is required by most linear solvers. * * This class provides access to all preconditioners offered by dune-istl using the * PreconditionerWrapper property: * \code * template * struct PreconditionerWrapper * { using type = Opm::Linear::PreconditionerWrapper$PRECONDITIONER; }; * \endcode * * Where the choices possible for '\c $PRECONDITIONER' are: * - \c Jacobi: A Jacobi preconditioner * - \c GaussSeidel: A Gauss-Seidel preconditioner * - \c SSOR: A symmetric successive overrelaxation (SSOR) preconditioner * - \c SOR: A successive overrelaxation (SOR) preconditioner * - \c ILUn: An ILU(n) preconditioner * - \c ILU0: An ILU(0) preconditioner. The results of this * preconditioner are the same as setting the * PreconditionerOrder property to 0 and using the ILU(n) * preconditioner. The reason for the existence of ILU0 is * that it is computationally cheaper because it does not * need to consider things which are only required for * higher orders */ template class ParallelBaseBackend { protected: using Implementation = GetPropType; using Simulator = GetPropType; using Scalar = GetPropType; using LinearSolverScalar = GetPropType; using SparseMatrixAdapter = GetPropType; using Vector = GetPropType; using BorderListCreator = GetPropType; using GridView = GetPropType; using Overlap = GetPropType; using OverlappingVector = GetPropType; using OverlappingMatrix = GetPropType; using PreconditionerWrapper = GetPropType; using SequentialPreconditioner = typename PreconditionerWrapper::SequentialPreconditioner; using ParallelPreconditioner = Opm::Linear::OverlappingPreconditioner; using ParallelScalarProduct = Opm::Linear::OverlappingScalarProduct; using ParallelOperator = Opm::Linear::OverlappingOperator; enum { dimWorld = GridView::dimensionworld }; public: ParallelBaseBackend(const Simulator& simulator) : simulator_(simulator) , gridSequenceNumber_( -1 ) , lastIterations_( -1 ) { overlappingMatrix_ = nullptr; overlappingb_ = nullptr; overlappingx_ = nullptr; } ~ParallelBaseBackend() { cleanup_(); } /*! * \brief Register all run-time parameters for the linear solver. */ static void registerParameters() { EWOMS_REGISTER_PARAM(TypeTag, Scalar, LinearSolverTolerance, "The maximum allowed error between of the linear solver"); EWOMS_REGISTER_PARAM(TypeTag, Scalar, LinearSolverAbsTolerance, "The maximum accepted error of the norm of the residual"); EWOMS_REGISTER_PARAM(TypeTag, unsigned, LinearSolverOverlapSize, "The size of the algebraic overlap for the linear solver"); EWOMS_REGISTER_PARAM(TypeTag, int, LinearSolverMaxIterations, "The maximum number of iterations of the linear solver"); EWOMS_REGISTER_PARAM(TypeTag, int, LinearSolverVerbosity, "The verbosity level of the linear solver"); PreconditionerWrapper::registerParameters(); } /*! * \brief Causes the solve() method to discared the structure of the linear system of * equations the next time it is called. */ void eraseMatrix() { cleanup_(); } /*! * \brief Set up the internal data structures required for the linear solver. * * This only specified the topology of the linear system of equations; it does does * *not* assign the values of the residual vector and its Jacobian matrix. */ void prepare(const SparseMatrixAdapter& M, const Vector& ) { // if grid has changed the sequence number has changed too int curSeqNum = simulator_.vanguard().gridSequenceNumber(); if (gridSequenceNumber_ == curSeqNum && overlappingMatrix_) // the grid has not changed since the overlappingMatrix_has been created, so // there's noting to do return; asImp_().cleanup_(); gridSequenceNumber_ = curSeqNum; BorderListCreator borderListCreator(simulator_.gridView(), simulator_.model().dofMapper()); // create the overlapping Jacobian matrix unsigned overlapSize = EWOMS_GET_PARAM(TypeTag, unsigned, LinearSolverOverlapSize); overlappingMatrix_ = new OverlappingMatrix(M.istlMatrix(), borderListCreator.borderList(), borderListCreator.blackList(), overlapSize); // create the overlapping vectors for the residual and the // solution overlappingb_ = new OverlappingVector(overlappingMatrix_->overlap()); overlappingx_ = new OverlappingVector(*overlappingb_); // writeOverlapToVTK_(); } /*! * \brief Assign values to the internal data structure for the residual vector. * * This method also cares about synchronizing that vector with the peer processes. */ void setResidual(const Vector& b) { // copy the interior values of the non-overlapping residual vector to the // overlapping one overlappingb_->assignAddBorder(b); } /*! * \brief Retrieve the synchronized internal residual vector. * * This only deals with entries which are local to the current process. */ void getResidual(Vector& b) const { // update the non-overlapping vector with the overlapping one overlappingb_->assignTo(b); } /*! * \brief Sets the values of the residual's Jacobian matrix. * * This method also synchronizes the data structure across the processes which are * involved in the simulation run. */ void setMatrix(const SparseMatrixAdapter& M) { overlappingMatrix_->assignFromNative(M.istlMatrix()); overlappingMatrix_->syncAdd(); } /*! * \brief Actually solve the linear system of equations. * * \return true if the residual reduction could be achieved, else false. */ bool solve(Vector& x) { #if ! DUNE_VERSION_NEWER(DUNE_COMMON, 2,7) Dune::FMatrixPrecision::set_singular_limit(1.e-30); Dune::FMatrixPrecision::set_absolute_limit(1.e-30); #endif (*overlappingx_) = 0.0; auto parPreCond = asImp_().preparePreconditioner_(); auto precondCleanupFn = [this]() -> void { this->asImp_().cleanupPreconditioner_(); }; auto precondCleanupGuard = Opm::make_guard(precondCleanupFn); // create the parallel scalar product and the parallel operator ParallelScalarProduct parScalarProduct(overlappingMatrix_->overlap()); ParallelOperator parOperator(*overlappingMatrix_); // retrieve the linear solver auto solver = asImp_().prepareSolver_(parOperator, parScalarProduct, *parPreCond); auto cleanupSolverFn = [this]() -> void { this->asImp_().cleanupSolver_(); }; GenericGuard solverGuard(cleanupSolverFn); // run the linear solver and have some fun auto result = asImp_().runSolver_(solver); // store number of iterations used lastIterations_ = result.second; // copy the result back to the non-overlapping vector overlappingx_->assignTo(x); // return the result of the solver return result.first; } /*! * \brief Return number of iterations used during last solve. */ size_t iterations () const { return lastIterations_; } protected: Implementation& asImp_() { return *static_cast(this); } const Implementation& asImp_() const { return *static_cast(this); } void cleanup_() { // create the overlapping Jacobian matrix and vectors delete overlappingMatrix_; delete overlappingb_; delete overlappingx_; overlappingMatrix_ = 0; overlappingb_ = 0; overlappingx_ = 0; } std::shared_ptr preparePreconditioner_() { int preconditionerIsReady = 1; try { // update sequential preconditioner precWrapper_.prepare(*overlappingMatrix_); } catch (const Dune::Exception& e) { std::cout << "Preconditioner threw exception \"" << e.what() << " on rank " << overlappingMatrix_->overlap().myRank() << "\n" << std::flush; preconditionerIsReady = 0; } // make sure that the preconditioner is also ready on all peer // ranks. preconditionerIsReady = simulator_.gridView().comm().min(preconditionerIsReady); if (!preconditionerIsReady) throw Opm::NumericalIssue("Creating the preconditioner failed"); // create the parallel preconditioner return std::make_shared(precWrapper_.get(), overlappingMatrix_->overlap()); } void cleanupPreconditioner_() { precWrapper_.cleanup(); } void writeOverlapToVTK_() { for (int lookedAtRank = 0; lookedAtRank < simulator_.gridView().comm().size(); ++lookedAtRank) { std::cout << "writing overlap for rank " << lookedAtRank << "\n" << std::flush; using VtkField = Dune::BlockVector >; int n = simulator_.gridView().size(/*codim=*/dimWorld); VtkField isInOverlap(n); VtkField rankField(n); isInOverlap = 0.0; rankField = 0.0; assert(rankField.two_norm() == 0.0); assert(isInOverlap.two_norm() == 0.0); auto vIt = simulator_.gridView().template begin(); const auto& vEndIt = simulator_.gridView().template end(); const auto& overlap = overlappingMatrix_->overlap(); for (; vIt != vEndIt; ++vIt) { int nativeIdx = simulator_.model().vertexMapper().map(*vIt); int localIdx = overlap.foreignOverlap().nativeToLocal(nativeIdx); if (localIdx < 0) continue; rankField[nativeIdx] = simulator_.gridView().comm().rank(); if (overlap.peerHasIndex(lookedAtRank, localIdx)) isInOverlap[nativeIdx] = 1.0; } using VtkWriter = Dune::VTKWriter; VtkWriter writer(simulator_.gridView(), Dune::VTK::conforming); writer.addVertexData(isInOverlap, "overlap"); writer.addVertexData(rankField, "rank"); std::ostringstream oss; oss << "overlap_rank=" << lookedAtRank; writer.write(oss.str().c_str(), Dune::VTK::ascii); } } const Simulator& simulator_; int gridSequenceNumber_; size_t lastIterations_; OverlappingMatrix *overlappingMatrix_; OverlappingVector *overlappingb_; OverlappingVector *overlappingx_; PreconditionerWrapper precWrapper_; }; }} // namespace Linear, Opm namespace Opm::Properties { //! make the linear solver shut up by default template struct LinearSolverVerbosity { static constexpr int value = 0; }; //! set the preconditioner relaxation parameter to 1.0 by default template struct PreconditionerRelaxation { using type = GetPropType; static constexpr type value = 1.0; }; //! set the preconditioner order to 0 by default template struct PreconditionerOrder { static constexpr int value = 0; }; //! by default use the same kind of floating point values for the linearization and for //! the linear solve template struct LinearSolverScalar { using type = GetPropType; }; template struct OverlappingMatrix { private: static constexpr int numEq = getPropValue(); using LinearSolverScalar = GetPropType; using MatrixBlock = Opm::MatrixBlock; using NonOverlappingMatrix = Dune::BCRSMatrix; public: using type = Opm::Linear::OverlappingBCRSMatrix; }; template struct Overlap { using type = typename GetPropType::Overlap; }; template struct OverlappingVector { static constexpr int numEq = getPropValue(); using LinearSolverScalar = GetPropType; using VectorBlock = Dune::FieldVector; using Overlap = GetPropType; using type = Opm::Linear::OverlappingBlockVector; }; template struct OverlappingScalarProduct { using OverlappingVector = GetPropType; using Overlap = GetPropType; using type = Opm::Linear::OverlappingScalarProduct; }; template struct OverlappingLinearOperator { using OverlappingMatrix = GetPropType; using OverlappingVector = GetPropType; using type = Opm::Linear::OverlappingOperator; }; #if DUNE_VERSION_NEWER(DUNE_ISTL, 2,7) template struct PreconditionerWrapper { using type = Opm::Linear::PreconditionerWrapperILU; }; #else template struct PreconditionerWrapper { using type = Opm::Linear::PreconditionerWrapperILU0; }; #endif //! set the default overlap size to 2 template struct LinearSolverOverlapSize { static constexpr int value = 2; }; //! set the default number of maximum iterations for the linear solver template struct LinearSolverMaxIterations { static constexpr int value = 1000; }; } // namespace Opm::Properties #endif