16 #ifndef dealii__precondition_block_h 17 #define dealii__precondition_block_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/subscriptor.h> 23 #include <deal.II/base/smartpointer.h> 24 #include <deal.II/lac/precondition_block_base.h> 28 DEAL_II_NAMESPACE_OPEN
80 template<
typename MatrixType,
typename inverse_type =
typename MatrixType::value_type>
89 typedef typename MatrixType::value_type
number;
264 template <
typename number2>
268 const bool transpose_diagonal)
const;
281 template <
typename number2>
285 const bool transpose_diagonal)
const;
310 <<
"The blocksize " << arg1
311 <<
" and the size of the matrix " << arg2
312 <<
" do not match.");
366 template<
typename MatrixType,
typename inverse_type =
typename MatrixType::value_type>
374 typedef typename MatrixType::value_type
number;
414 inverse_type
value()
const;
522 template <
typename number2>
528 template <
typename number2>
537 template <
typename number2>
543 template <
typename number2>
549 template <
typename number2>
555 template <
typename number2>
586 template <
typename number2>
591 friend class Accessor;
632 template<
typename MatrixType,
typename inverse_type =
typename MatrixType::value_type>
650 typedef typename MatrixType::value_type
number;
679 template <
typename number2>
692 template <
typename number2>
703 template <
typename number2>
716 template <
typename number2>
722 template <
typename number2>
728 template <
typename number2>
746 template <
typename number2>
749 const bool transpose_diagonal,
750 const bool adding)
const;
761 template <
typename number2>
764 const bool transpose_diagonal,
765 const bool adding)
const;
790 template<
typename MatrixType,
typename inverse_type =
typename MatrixType::value_type>
803 typedef typename MatrixType::value_type
number;
840 template <
typename number2>
846 template <
typename number2>
852 template <
typename number2>
858 template <
typename number2>
867 template<
typename MatrixType,
typename inverse_type>
877 template<
typename MatrixType,
typename inverse_type>
884 const unsigned int nb = i/bs;
901 template<
typename MatrixType,
typename inverse_type>
908 b_iterator(&matrix->inverse(0), 0, 0),
909 b_end(&matrix->inverse(0), 0, 0)
911 bs = matrix->block_size();
916 if (a_block == matrix->size())
921 b_iterator = matrix->inverse(a_block).begin(r);
922 b_end = matrix->inverse(a_block).end();
929 template<
typename MatrixType,
typename inverse_type>
941 template<
typename MatrixType,
typename inverse_type>
949 return bs * a_block + b_iterator->column();
953 template<
typename MatrixType,
typename inverse_type>
961 return b_iterator->value();
965 template<
typename MatrixType,
typename inverse_type>
975 template<
typename MatrixType,
typename inverse_type>
997 template<
typename MatrixType,
typename inverse_type>
1006 template<
typename MatrixType,
typename inverse_type>
1015 template<
typename MatrixType,
typename inverse_type>
1021 if (accessor.a_block == accessor.matrix->size() &&
1022 accessor.a_block == other.accessor.a_block)
1025 if (accessor.a_block != other.accessor.a_block)
1028 return (accessor.row() == other.accessor.row() &&
1029 accessor.column() == other.accessor.column());
1033 template<
typename MatrixType,
typename inverse_type>
1039 return ! (*
this == other);
1043 template<
typename MatrixType,
typename inverse_type>
1047 operator < (
const const_iterator &other)
const 1049 return (accessor.row() < other.accessor.row() ||
1050 (accessor.row() == other.accessor.row() &&
1051 accessor.column() < other.accessor.column()));
1055 template<
typename MatrixType,
typename inverse_type>
1060 return const_iterator(
this, 0);
1064 template<
typename MatrixType,
typename inverse_type>
1073 template<
typename MatrixType,
typename inverse_type>
1080 return const_iterator(
this, r);
1085 template<
typename MatrixType,
typename inverse_type>
1092 return const_iterator(
this, r+1);
1097 DEAL_II_NAMESPACE_CLOSE
FullMatrix< inverse_type >::const_iterator b_end
std::size_t memory_consumption() const
types::global_dof_index size_type
types::global_dof_index size_type
void set_permutation(const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation)
const_iterator & operator++()
void Tvmult_add(Vector< number2 > &, const Vector< number2 > &) const
inverse_type value() const
#define DeclException2(Exception2, type1, type2, outsequence)
FullMatrix< inverse_type >::const_iterator b_iterator
void Tstep(Vector< number2 > &dst, const Vector< number2 > &rhs) const
AdditionalData(const size_type block_size, const double relaxation=1., const bool invert_diagonal=true, const bool same_diagonal=false)
static ::ExceptionBase & ExcWrongBlockSize(int arg1, int arg2)
void step(Vector< number2 > &dst, const Vector< number2 > &rhs) const
const Accessor * operator->() const
void Tvmult(Vector< number2 > &, const Vector< number2 > &) const
const_iterator end() const
Accessor(const PreconditionBlockJacobi< MatrixType, inverse_type > *matrix, const size_type row)
void backward_step(Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const
types::global_dof_index size_type
void initialize(const MatrixType &A, const AdditionalData parameters)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void Tvmult_add(Vector< number2 > &, const Vector< number2 > &) const
std::vector< size_type > inverse_permutation
const_iterator(const PreconditionBlockJacobi< MatrixType, inverse_type > *matrix, const size_type row)
size_type block_size() const
MatrixType::value_type number
static ::ExceptionBase & ExcInverseMatricesAlreadyExist()
void step(Vector< number2 > &dst, const Vector< number2 > &rhs) const
bool operator!=(const const_iterator &) const
bool operator==(const const_iterator &) const
unsigned int global_dof_index
types::global_dof_index size_type
#define Assert(cond, exc)
value_type el(size_type i, size_type j) const
bool store_diagonals() const
void vmult(Vector< number2 > &, const Vector< number2 > &) const
void Tvmult(Vector< number2 > &, const Vector< number2 > &) const
void vmult_add(Vector< number2 > &, const Vector< number2 > &) const
#define DeclException0(Exception0)
MatrixType::value_type number
PreconditionBlockBase< inverse_type >::Inversion inversion
void vmult_add(Vector< number2 > &, const Vector< number2 > &) const
const_iterator begin() const
PreconditionBlock(bool store_diagonals=false)
void forward_step(Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const
void Tvmult(Vector< number2 > &, const Vector< number2 > &) const
static ::ExceptionBase & ExcIteratorPastEnd()
void do_vmult(Vector< number2 > &, const Vector< number2 > &, bool adding) const
void backward(Vector< number2 > &, const Vector< number2 > &, const bool transpose_diagonal, const bool adding) const
void vmult(Vector< number2 > &, const Vector< number2 > &) const
types::global_dof_index size_type
MatrixType::value_type number
void Tstep(Vector< number2 > &dst, const Vector< number2 > &rhs) const
std::vector< size_type > permutation
void vmult(Vector< number2 > &, const Vector< number2 > &) const
void Tstep(Vector< number2 > &dst, const Vector< number2 > &rhs) const
void step(Vector< number2 > &dst, const Vector< number2 > &rhs) const
const Accessor & operator*() const
unsigned int size() const
void invert_permuted_diagblocks()
const PreconditionBlockJacobi< MatrixType, inverse_type > * matrix
void forward(Vector< number2 > &, const Vector< number2 > &, const bool transpose_diagonal, const bool adding) const
SmartPointer< const MatrixType, PreconditionBlock< MatrixType, inverse_type > > A
bool operator<(const const_iterator &) const
MatrixType::value_type number