16#ifndef dealii_precondition_block_h
17#define dealii_precondition_block_h
81template <
typename MatrixType,
82 typename inverse_type =
typename MatrixType::value_type>
90 using number =
typename MatrixType::value_type;
273 template <
typename number2>
278 const bool transpose_diagonal)
const;
291 template <
typename number2>
296 const bool transpose_diagonal)
const;
324 <<
"The blocksize " << arg1 <<
" and the size of the matrix "
325 << arg2 <<
" do not match.");
377template <
typename MatrixType,
378 typename inverse_type =
typename MatrixType::value_type>
387 using number =
typename MatrixType::value_type;
537 template <
typename number2>
544 template <
typename number2>
554 template <
typename number2>
561 template <
typename number2>
568 template <
typename number2>
575 template <
typename number2>
611 template <
typename number2>
654template <
typename MatrixType,
655 typename inverse_type =
typename MatrixType::value_type>
674 using number =
typename MatrixType::value_type;
700 template <
typename number2>
714 template <
typename number2>
726 template <
typename number2>
740 template <
typename number2>
747 template <
typename number2>
754 template <
typename number2>
773 template <
typename number2>
777 const bool transpose_diagonal,
778 const bool adding)
const;
789 template <
typename number2>
793 const bool transpose_diagonal,
794 const bool adding)
const;
817template <
typename MatrixType,
818 typename inverse_type =
typename MatrixType::value_type>
832 using number =
typename MatrixType::value_type;
868 template <
typename number2>
875 template <
typename number2>
882 template <
typename number2>
889 template <
typename number2>
899template <
typename MatrixType,
typename inverse_type>
909template <
typename MatrixType,
typename inverse_type>
914 const unsigned int nb = i / bs;
921 if (jb + nb * bs != j)
931template <
typename MatrixType,
typename inverse_type>
939 , b_iterator(&
matrix->inverse(0), 0, 0)
940 , b_end(&
matrix->inverse(0), 0, 0)
944 if (a_block ==
matrix->size())
949 b_iterator =
matrix->inverse(a_block).begin(r);
950 b_end =
matrix->inverse(a_block).end();
956template <
typename MatrixType,
typename inverse_type>
959 inverse_type>::const_iterator::Accessor::row()
const
963 return bs * a_block + b_iterator->row();
967template <
typename MatrixType,
typename inverse_type>
970 inverse_type>::const_iterator::Accessor::column()
const
974 return bs * a_block + b_iterator->column();
978template <
typename MatrixType,
typename inverse_type>
981 inverse_type>::const_iterator::Accessor::value()
const
985 return b_iterator->value();
989template <
typename MatrixType,
typename inverse_type>
998template <
typename MatrixType,
typename inverse_type>
1000 inverse_type>::const_iterator &
1021template <
typename MatrixType,
typename inverse_type>
1023 const_iterator::Accessor &
1031template <
typename MatrixType,
typename inverse_type>
1033 const_iterator::Accessor *
1035 inverse_type>::const_iterator::operator->()
const
1041template <
typename MatrixType,
typename inverse_type>
1044 const const_iterator &other)
const
1046 if (accessor.a_block == accessor.matrix->size() &&
1047 accessor.a_block == other.accessor.a_block)
1050 if (accessor.a_block != other.accessor.a_block)
1053 return (accessor.row() == other.accessor.row() &&
1054 accessor.column() == other.accessor.column());
1058template <
typename MatrixType,
typename inverse_type>
1061 const const_iterator &other)
const
1063 return !(*
this == other);
1067template <
typename MatrixType,
typename inverse_type>
1070 const const_iterator &other)
const
1072 return (accessor.row() < other.accessor.row() ||
1073 (accessor.row() == other.accessor.row() &&
1074 accessor.column() < other.accessor.column()));
1078template <
typename MatrixType,
typename inverse_type>
1087template <
typename MatrixType,
typename inverse_type>
1096template <
typename MatrixType,
typename inverse_type>
1108template <
typename MatrixType,
typename inverse_type>
typename Table< 2, number >::const_iterator const_iterator
iterator begin(const size_type r)
iterator end(const size_type r)
void log_statistics() const
Householder< number > & inverse_householder(size_type i)
unsigned int size() const
FullMatrix< number > & inverse(size_type i)
bool store_diagonals() const
LAPACKFullMatrix< number > & inverse_svd(size_type i)
const PreconditionBlockJacobi< MatrixType, inverse_type > * matrix
FullMatrix< inverse_type >::const_iterator b_end
FullMatrix< inverse_type >::const_iterator b_iterator
Accessor(const PreconditionBlockJacobi< MatrixType, inverse_type > *matrix, const size_type row)
inverse_type value() const
bool operator!=(const const_iterator &) const
const_iterator(const PreconditionBlockJacobi< MatrixType, inverse_type > *matrix, const size_type row)
const_iterator & operator++()
bool operator==(const const_iterator &) const
bool operator<(const const_iterator &) const
const Accessor & operator*() const
const Accessor * operator->() const
void vmult(Vector< number2 > &, const Vector< number2 > &) const
void step(Vector< number2 > &dst, const Vector< number2 > &rhs) const
void vmult_add(Vector< number2 > &, const Vector< number2 > &) const
const_iterator begin(const size_type r) const
const_iterator end() const
void Tvmult(Vector< number2 > &, const Vector< number2 > &) const
void do_vmult(Vector< number2 > &, const Vector< number2 > &, bool adding) const
void Tvmult_add(Vector< number2 > &, const Vector< number2 > &) const
typename MatrixType::value_type number
const_iterator begin() const
const_iterator end(const size_type r) const
friend class const_iterator
void Tstep(Vector< number2 > &dst, const Vector< number2 > &rhs) const
types::global_dof_index size_type
PreconditionBlockSOR(bool store)
void vmult_add(Vector< number2 > &, const Vector< number2 > &) const
void vmult(Vector< number2 > &, const Vector< number2 > &) const
typename MatrixType::value_type number
void step(Vector< number2 > &dst, const Vector< number2 > &rhs) const
void Tstep(Vector< number2 > &dst, const Vector< number2 > &rhs) const
void forward(Vector< number2 > &, const Vector< number2 > &, const bool transpose_diagonal, const bool adding) const
void backward(Vector< number2 > &, const Vector< number2 > &, const bool transpose_diagonal, const bool adding) const
void Tvmult(Vector< number2 > &, const Vector< number2 > &) const
void Tvmult_add(Vector< number2 > &, const Vector< number2 > &) const
void vmult(Vector< number2 > &, const Vector< number2 > &) const
void step(Vector< number2 > &dst, const Vector< number2 > &rhs) const
typename MatrixType::value_type number
void Tstep(Vector< number2 > &dst, const Vector< number2 > &rhs) const
types::global_dof_index size_type
void Tvmult(Vector< number2 > &, const Vector< number2 > &) const
AdditionalData(const size_type block_size, const double relaxation=1., const bool invert_diagonal=true, const bool same_diagonal=false)
PreconditionBlockBase< inverse_type >::Inversion inversion
PreconditionBlock(bool store_diagonals=false)
~PreconditionBlock() override=default
value_type el(size_type i, size_type j) const
SmartPointer< const MatrixType, PreconditionBlock< MatrixType, inverse_type > > A
void forward_step(Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const
void initialize(const MatrixType &A, const AdditionalData parameters)
void invert_permuted_diagblocks()
void backward_step(Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const
std::size_t memory_consumption() const
std::vector< size_type > inverse_permutation
void set_permutation(const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation)
std::vector< size_type > permutation
size_type block_size() const
typename MatrixType::value_type number
void initialize(const MatrixType &A, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const AdditionalData parameters)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DeclException0(Exception0)
static ::ExceptionBase & ExcInverseMatricesAlreadyExist()
#define Assert(cond, exc)
static ::ExceptionBase & ExcIteratorPastEnd()
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & ExcWrongBlockSize(int arg1, int arg2)
#define AssertIndexRange(index, range)
@ matrix
Contents is actually a matrix.
types::global_dof_index size_type
unsigned int global_dof_index