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
30 template<
typename MatrixType,
typename inverse_type>
83 template<
typename MatrixType,
typename inverse_type =
typename MatrixType::value_type>
92 typedef typename MatrixType::value_type
number;
237 value_type
el(size_type i,
268 template <
typename number2>
272 const bool transpose_diagonal)
const;
285 template <
typename number2>
289 const bool transpose_diagonal)
const;
314 <<
"The blocksize " << arg1
315 <<
" and the size of the matrix " << arg2
316 <<
" do not match.");
370 template<
typename MatrixType,
typename inverse_type =
typename MatrixType::value_type>
378 typedef typename MatrixType::value_type
number;
403 const size_type row);
408 size_type row()
const;
413 size_type column()
const;
418 inverse_type value()
const;
457 const size_type row);
477 const Accessor *operator-> ()
const;
525 template <
typename number2>
531 template <
typename number2>
540 template <
typename number2>
546 template <
typename number2>
552 template <
typename number2>
558 template <
typename number2>
589 template <
typename number2>
635 template<
typename MatrixType,
typename inverse_type =
typename MatrixType::value_type>
653 typedef typename MatrixType::value_type
number;
682 template <
typename number2>
695 template <
typename number2>
706 template <
typename number2>
719 template <
typename number2>
725 template <
typename number2>
731 template <
typename number2>
749 template <
typename number2>
752 const bool transpose_diagonal,
753 const bool adding)
const;
764 template <
typename number2>
767 const bool transpose_diagonal,
768 const bool adding)
const;
793 template<
typename MatrixType,
typename inverse_type =
typename MatrixType::value_type>
806 typedef typename MatrixType::value_type
number;
843 template <
typename number2>
849 template <
typename number2>
855 template <
typename number2>
861 template <
typename number2>
870 template<
typename MatrixType,
typename inverse_type>
880 template<
typename MatrixType,
typename inverse_type>
887 const unsigned int nb = i/bs;
904 template<
typename MatrixType,
typename inverse_type>
911 b_iterator(&matrix->
inverse(0), 0, 0),
912 b_end(&matrix->
inverse(0), 0, 0)
919 if (a_block == matrix->
size())
928 ExcIndexRange(a_block, 0, matrix->
size()));
932 template<
typename MatrixType,
typename inverse_type>
938 ExcIteratorPastEnd());
940 return bs * a_block + b_iterator->row();
944 template<
typename MatrixType,
typename inverse_type>
950 ExcIteratorPastEnd());
952 return bs * a_block + b_iterator->column();
956 template<
typename MatrixType,
typename inverse_type>
962 ExcIteratorPastEnd());
964 return b_iterator->value();
968 template<
typename MatrixType,
typename inverse_type>
974 accessor(matrix, row)
978 template<
typename MatrixType,
typename inverse_type>
983 Assert (*
this != accessor.matrix->end(), ExcIteratorPastEnd());
985 ++accessor.b_iterator;
986 if (accessor.b_iterator == accessor.b_end)
990 if (accessor.a_block < accessor.matrix->size())
992 accessor.b_iterator = accessor.matrix->inverse(accessor.a_block).begin();
993 accessor.b_end = accessor.matrix->inverse(accessor.a_block).end();
1000 template<
typename MatrixType,
typename inverse_type>
1009 template<
typename MatrixType,
typename inverse_type>
1018 template<
typename MatrixType,
typename inverse_type>
1024 if (accessor.a_block == accessor.matrix->size() &&
1025 accessor.a_block == other.accessor.a_block)
1028 if (accessor.a_block != other.accessor.a_block)
1031 return (accessor.row() == other.accessor.row() &&
1032 accessor.column() == other.accessor.column());
1036 template<
typename MatrixType,
typename inverse_type>
1042 return ! (*
this == other);
1046 template<
typename MatrixType,
typename inverse_type>
1050 operator < (
const const_iterator &other)
const 1052 return (accessor.row() < other.accessor.row() ||
1053 (accessor.row() == other.accessor.row() &&
1054 accessor.column() < other.accessor.column()));
1058 template<
typename MatrixType,
typename inverse_type>
1063 return const_iterator(
this, 0);
1067 template<
typename MatrixType,
typename inverse_type>
1076 template<
typename MatrixType,
typename inverse_type>
1082 Assert (r < this->
A->m(), ExcIndexRange(r, 0, this->
A->m()));
1083 return const_iterator(
this, r);
1088 template<
typename MatrixType,
typename inverse_type>
1094 Assert (r < this->
A->m(), ExcIndexRange(r, 0, this->
A->m()));
1095 return const_iterator(
this, r+1);
1100 DEAL_II_NAMESPACE_CLOSE
FullMatrix< inverse_type >::const_iterator b_end
std::size_t memory_consumption() const
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++()
inverse_type value() const
FullMatrix< inverse_type >::const_iterator b_iterator
AdditionalData(const size_type block_size, const double relaxation=1., const bool invert_diagonal=true, const bool same_diagonal=false)
const Accessor * operator->() const
const_iterator begin() const
const_iterator end() const
DeclException2(ExcWrongBlockSize, int, int,<< "The blocksize "<< arg1<< " and the size of the matrix "<< arg2<< " do not match.")
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 invert_permuted_diagblocks(const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation)
void initialize(const MatrixType &A, const AdditionalData parameters)
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
FullMatrix< inverse_type > & inverse(size_type i)
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
MatrixType::value_type number
PreconditionBlockBase< inverse_type >::Inversion inversion
const_iterator begin() const
PreconditionBlock(bool store_diagonals=false)
const_iterator end() const
void forward_step(Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const
DeclException0(ExcInverseMatricesAlreadyExist)
types::global_dof_index size_type
MatrixType::value_type number
std::vector< size_type > permutation
const Accessor & operator*() const
unsigned int size() const
const PreconditionBlockJacobi< MatrixType, inverse_type > * matrix
SmartPointer< const MatrixType, PreconditionBlock< MatrixType, inverse_type > > A
bool operator<(const const_iterator &) const
MatrixType::value_type number