16#ifndef dealii_cuda_sparse_matrix_h
17#define dealii_cuda_sparse_matrix_h
25#ifdef DEAL_II_WITH_CUDA
47 template <
typename Number>
86 const ::SparseMatrix<Number> &sparse_matrix_host);
123 const ::SparseMatrix<Number> &sparse_matrix_host);
162 template <
class StreamType>
164 print(StreamType &out,
165 const bool across =
false,
166 const bool diagonal_first =
true)
const;
190 const unsigned int precision = 3,
191 const bool scientific =
true,
192 const unsigned int width = 0,
193 const char * zero_string =
" ",
194 const double denominator = 1.)
const;
326 std::tuple<Number *, int *, int *, cusparseMatDescr_t, cusparseSpMatDescr_t>
354 std::unique_ptr<Number[], void (*)(Number *)>
val_dev;
379 template <
typename Number>
388 template <
typename Number>
397 template <
typename Number>
406 template <
typename Number>
407 template <
class StreamType>
411 const bool diagonal_first)
const
417 std::vector<int> rows(n_rows + 1);
418 std::vector<int> cols(nnz);
419 std::vector<double> val(nnz);
424 bool has_diagonal =
false;
425 Number diagonal = Number();
432 for (
size_type j = rows[i]; j < rows[i + 1] && cols[j] <= i; ++j)
439 out <<
' ' << i <<
',' << i <<
':' << diagonal;
441 out <<
'(' << i <<
',' << i <<
") " << diagonal
447 for (
size_type j = rows[i]; j < rows[i + 1]; ++j)
449 if (has_diagonal && i == cols[j])
452 out <<
' ' << i <<
',' << cols[j] <<
':' << val[j];
454 out <<
'(' << i <<
',' << cols[j] <<
") " << val[j] << std::endl;
463 template <
typename Number>
466 const unsigned int precision,
467 const bool scientific,
468 const unsigned int width_,
469 const char * zero_string,
470 const double denominator)
const
476 std::vector<int> rows(n_rows + 1);
477 std::vector<int> cols(nnz);
478 std::vector<Number> val(nnz);
483 unsigned int width = width_;
485 std::ios::fmtflags old_flags = out.flags();
486 unsigned int old_precision = out.precision(precision);
490 out.setf(std::ios::scientific, std::ios::floatfield);
492 width = precision + 7;
496 out.setf(std::ios::fixed, std::ios::floatfield);
498 width = precision + 2;
508 out << std::setw(width) << val[j] * Number(denominator) <<
' ';
512 out << std::setw(width) << zero_string <<
' ';
519 out.precision(old_precision);
520 out.flags(old_flags);
void vmult_add(LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &src) const
void print(StreamType &out, const bool across=false, const bool diagonal_first=true) const
std::tuple< Number *, int *, int *, cusparseMatDescr_t, cusparseSpMatDescr_t > get_cusparse_matrix() const
void Tvmult_add(LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &src) const
std::unique_ptr< int[], void(*)(int *)> column_index_dev
cusparseHandle_t cusparse_handle
void Tvmult(LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &src) const
std::unique_ptr< int[], void(*)(int *)> row_ptr_dev
void reinit(Utilities::CUDA::Handle &handle, const ::SparseMatrix< Number > &sparse_matrix_host)
Number frobenius_norm() const
Number matrix_norm_square(const LinearAlgebra::CUDAWrappers::Vector< Number > &v) const
SparseMatrix & operator*=(const Number factor)
cusparseSpMatDescr_t sp_descr
Number residual(LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &x, const LinearAlgebra::CUDAWrappers::Vector< Number > &b) const
SparseMatrix & operator/=(const Number factor)
SparseMatrix & operator=(CUDAWrappers::SparseMatrix< Number > &&)
std::unique_ptr< Number[], void(*)(Number *)> val_dev
void vmult(LinearAlgebra::CUDAWrappers::Vector< Number > &dst, const LinearAlgebra::CUDAWrappers::Vector< Number > &src) const
SparseMatrix(const CUDAWrappers::SparseMatrix< Number > &)=delete
std::size_t n_nonzero_elements() const
void print_formatted(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1.) const
Number linfty_norm() const
SparseMatrix & operator=(const CUDAWrappers::SparseMatrix< Number > &)=delete
Number matrix_scalar_product(const LinearAlgebra::CUDAWrappers::Vector< Number > &u, const LinearAlgebra::CUDAWrappers::Vector< Number > &v) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcIO()
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotInitialized()
#define AssertThrow(cond, exc)
void copy_to_host(const ArrayView< const T, MemorySpace::CUDA > &in, ArrayView< T, MemorySpace::Host > &out)