16#ifndef dealii_block_sparse_matrix_h
17#define dealii_block_sparse_matrix_h
49template <
typename number>
220 template <
typename block_number>
229 template <
typename block_number,
typename nonblock_number>
238 template <
typename block_number,
typename nonblock_number>
247 template <
typename nonblock_number>
256 template <
typename block_number>
265 template <
typename block_number,
typename nonblock_number>
274 template <
typename block_number,
typename nonblock_number>
283 template <
typename nonblock_number>
300 template <
class BlockVectorType>
303 const BlockVectorType &src,
304 const number omega = 1.)
const;
311 template <
typename number2>
315 const number omega = 1.)
const;
344 const unsigned int precision = 3,
345 const bool scientific =
true,
346 const unsigned int width = 0,
347 const char * zero_string =
" ",
348 const double denominator = 1.)
const;
378template <
typename number>
384 for (
size_type r = 0; r < this->n_block_rows(); ++r)
385 for (
size_type c = 0; c < this->n_block_cols(); ++c)
386 this->block(r, c) = d;
393template <
typename number>
394template <
typename block_number>
399 BaseClass::vmult_block_block(dst, src);
404template <
typename number>
405template <
typename block_number,
typename nonblock_number>
410 BaseClass::vmult_block_nonblock(dst, src);
415template <
typename number>
416template <
typename block_number,
typename nonblock_number>
421 BaseClass::vmult_nonblock_block(dst, src);
426template <
typename number>
427template <
typename nonblock_number>
432 BaseClass::vmult_nonblock_nonblock(dst, src);
437template <
typename number>
438template <
typename block_number>
443 BaseClass::Tvmult_block_block(dst, src);
448template <
typename number>
449template <
typename block_number,
typename nonblock_number>
454 BaseClass::Tvmult_block_nonblock(dst, src);
459template <
typename number>
460template <
typename block_number,
typename nonblock_number>
465 BaseClass::Tvmult_nonblock_block(dst, src);
470template <
typename number>
471template <
typename nonblock_number>
476 BaseClass::Tvmult_nonblock_nonblock(dst, src);
481template <
typename number>
482template <
class BlockVectorType>
485 const BlockVectorType &src,
486 const number omega)
const
489 Assert(dst.n_blocks() == this->n_block_rows(),
491 Assert(src.n_blocks() == this->n_block_cols(),
496 for (
size_type i = 0; i < this->n_block_rows(); ++i)
497 this->block(i, i).precondition_Jacobi(dst.block(i), src.block(i), omega);
502template <
typename number>
503template <
typename number2>
507 const number omega)
const
512 Assert(this->n_block_cols() == 1,
513 ExcMessage(
"This function only works if the matrix has "
515 Assert(this->n_block_rows() == 1,
516 ExcMessage(
"This function only works if the matrix has "
521 this->block(0, 0).precondition_Jacobi(dst, src, omega);
SparseMatrix< number > BlockType
const value_type & const_reference
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, true > > const_iterator
types::global_dof_index size_type
typename BlockType::value_type value_type
const value_type * const_pointer
MatrixIterator< BlockMatrixIterators::Accessor< BlockMatrixBase, false > > iterator
BlockSparseMatrix(const BlockSparsityPattern &sparsity)
typename BaseClass::value_type value_type
typename BaseClass::iterator iterator
typename BaseClass::const_pointer const_pointer
BlockSparseMatrix()=default
virtual void reinit(const BlockSparsityPattern &sparsity)
size_type n_nonzero_elements() const
void precondition_Jacobi(Vector< number2 > &dst, const Vector< number2 > &src, const number omega=1.) const
std::size_t memory_consumption() const
void vmult(BlockVector< block_number > &dst, const Vector< nonblock_number > &src) const
void Tvmult(Vector< nonblock_number > &dst, const Vector< nonblock_number > &src) const
virtual ~BlockSparseMatrix() override
typename BaseClass::pointer pointer
BlockSparseMatrix & operator=(const BlockSparseMatrix &)
typename BaseClass::const_reference const_reference
size_type n_actually_nonzero_elements(const double threshold=0.0) const
void Tvmult(Vector< nonblock_number > &dst, const BlockVector< block_number > &src) const
typename BaseClass::size_type size_type
void vmult(Vector< nonblock_number > &dst, const BlockVector< block_number > &src) const
void vmult(Vector< nonblock_number > &dst, const Vector< nonblock_number > &src) const
SmartPointer< const BlockSparsityPattern, BlockSparseMatrix< number > > sparsity_pattern
typename BaseClass::const_iterator const_iterator
void Tvmult(BlockVector< block_number > &dst, const BlockVector< block_number > &src) const
void Tvmult(BlockVector< block_number > &dst, const Vector< nonblock_number > &src) const
void vmult(BlockVector< block_number > &dst, const BlockVector< block_number > &src) const
typename BaseClass::BlockType BlockType
typename BaseClass::reference reference
size_type get_row_length(const size_type row) 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
void precondition_Jacobi(BlockVectorType &dst, const BlockVectorType &src, const number omega=1.) const
BlockSparseMatrix & operator=(const double d)
const BlockSparsityPattern & get_sparsity_pattern() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DeclException0(Exception0)
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcBlockDimensionMismatch()
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcMessage(std::string arg1)