19#ifndef CU_MATRIX_DESCRIPTION_HPP
20#define CU_MATRIX_DESCRIPTION_HPP
21#include <opm/simulators/linalg/gpuistl/detail/CuSparseResource.hpp>
22#include <opm/simulators/linalg/gpuistl/detail/cusparse_safe_call.hpp>
44 auto description = std::make_shared<GpuSparseMatrixDescription>();
47 OPM_CUSPARSE_SAFE_CALL(cusparseSetMatType(description->get(), CUSPARSE_MATRIX_TYPE_GENERAL));
48 OPM_CUSPARSE_SAFE_CALL(cusparseSetMatIndexBase(description->get(), CUSPARSE_INDEX_BASE_ZERO));
63 OPM_CUSPARSE_SAFE_CALL(cusparseSetMatFillMode(description->get(), CUSPARSE_FILL_MODE_LOWER));
64 OPM_CUSPARSE_SAFE_CALL(cusparseSetMatDiagType(description->get(), CUSPARSE_DIAG_TYPE_UNIT));
78 OPM_CUSPARSE_SAFE_CALL(cusparseSetMatFillMode(description->get(), CUSPARSE_FILL_MODE_UPPER));
79 OPM_CUSPARSE_SAFE_CALL(cusparseSetMatDiagType(description->get(), CUSPARSE_DIAG_TYPE_NON_UNIT));
The CuSparseResource class wraps a CuSparse resource in a proper RAII pattern.
Definition CuSparseResource.hpp:55
Contains wrappers to make the CuBLAS library behave as a modern C++ library with function overlading.
Definition autotuner.hpp:29
GpuSparseMatrixDescriptionPtr createLowerDiagonalDescription()
createLowerDiagonalDescription creates a lower diagonal matrix description
Definition CuMatrixDescription.hpp:60
GpuSparseMatrixDescriptionPtr createUpperDiagonalDescription()
createUpperDiagonalDescription creates an upper diagonal matrix description
Definition CuMatrixDescription.hpp:75
GpuSparseMatrixDescriptionPtr createMatrixDescription()
createMatrixDescription creates a default matrix description
Definition CuMatrixDescription.hpp:42
std::shared_ptr< CuSparseResource< cusparseMatDescr_t > > GpuSparseMatrixDescriptionPtr
Pointer to GpuSparseMatrixDescription holder.
Definition CuMatrixDescription.hpp:35