/* Copyright 2021 Equinor ASA 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 . */ #ifndef OPM_OPENCLCPR_HPP #define OPM_OPENCLCPR_HPP #include #include #include #include #include #include #include #include #include #include namespace Opm::Accelerator { template class BlockedMatrix; /// This class implements a Constrained Pressure Residual (CPR) preconditioner template class openclCPR : public openclPreconditioner, public CprCreation { using Base = openclPreconditioner; using Base::N; using Base::Nb; using Base::nnz; using Base::nnzb; using Base::verbosity; using Base::context; using Base::queue; using Base::events; using Base::err; private: std::vector> d_Amatrices, d_Rmatrices; // scalar matrices that represent the AMG hierarchy std::vector d_PcolIndices; std::vector d_invDiags; std::vector d_t, d_f, d_u; // intermediate vectors used during amg cycle std::unique_ptr d_rs; // use before extracting the pressure std::unique_ptr d_weights; // the quasiimpes weights, used to extract pressure std::unique_ptr> d_mat; // stores blocked matrix std::unique_ptr d_coarse_y, d_coarse_x; // stores the scalar vectors std::once_flag opencl_buffers_allocated; // only allocate OpenCL Buffers once std::unique_ptr> bilu0; // Blocked ILU0 preconditioner std::unique_ptr> coarse_solver; // coarse solver is scalar bool opencl_ilu_parallel; // whether ILU0 operation should be parallelized // Initialize and allocate matrices and vectors void init_opencl_buffers(); // Copy matrices and vectors to GPU void opencl_upload(); // apply pressure correction to vector void apply_amg(const cl::Buffer& y, cl::Buffer& x); void amg_cycle_gpu(const int level, cl::Buffer &y, cl::Buffer &x); public: openclCPR(bool opencl_ilu_parallel, int verbosity); bool analyze_matrix(BlockedMatrix* mat) override; bool analyze_matrix(BlockedMatrix* mat, BlockedMatrix* jacMat) override; // set own Opencl variables, but also that of the bilu0 preconditioner void setOpencl(std::shared_ptr& context, std::shared_ptr& queue) override; // applies blocked ilu0 // also applies amg for pressure component void apply(const cl::Buffer& y, cl::Buffer& x) override; void apply(Scalar&, Scalar&) override {} bool create_preconditioner(BlockedMatrix* mat) override; bool create_preconditioner(BlockedMatrix* mat, BlockedMatrix* jacMat) override; }; } // namespace Opm::Accelerator #endif