Reference documentation for deal.II version 8.5.1
precondition_block.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__precondition_block_h
17 #define dealii__precondition_block_h
18 
19 
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>
25 
26 #include <vector>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
80 template<typename MatrixType, typename inverse_type = typename MatrixType::value_type>
82  : public virtual Subscriptor,
83  protected PreconditionBlockBase<inverse_type>
84 {
85 private:
89  typedef typename MatrixType::value_type number;
90 
94  typedef inverse_type value_type;
95 
96 public:
101 
106  {
107  public:
113  const double relaxation = 1.,
114  const bool invert_diagonal = true,
115  const bool same_diagonal = false);
116 
120  double relaxation;
121 
126 
131 
140 
146  double threshold;
147  };
148 
149 
153  PreconditionBlock(bool store_diagonals = false);
154 
159 
167  void initialize (const MatrixType &A,
168  const AdditionalData parameters);
169 protected:
181  void initialize (const MatrixType &A,
182  const std::vector<size_type> &permutation,
183  const std::vector<size_type> &inverse_permutation,
184  const AdditionalData parameters);
185 
209  void set_permutation(const std::vector<size_type> &permutation,
210  const std::vector<size_type> &inverse_permutation);
211 
216 public:
222  void clear();
223 
227  bool empty () const;
228 
234  size_type j) const;
235 
251  void invert_diagblocks();
252 
264  template <typename number2>
265  void forward_step (Vector<number2> &dst,
266  const Vector<number2> &prev,
267  const Vector<number2> &src,
268  const bool transpose_diagonal) const;
269 
281  template <typename number2>
282  void backward_step (Vector<number2> &dst,
283  const Vector<number2> &prev,
284  const Vector<number2> &src,
285  const bool transpose_diagonal) const;
286 
287 
291  size_type block_size () const;
292 
297  std::size_t memory_consumption () const;
298 
309  int, int,
310  << "The blocksize " << arg1
311  << " and the size of the matrix " << arg2
312  << " do not match.");
313 
318 
320 
321 protected:
327 
338  double relaxation;
339 
343  std::vector<size_type> permutation;
344 
348  std::vector<size_type> inverse_permutation;
349 };
350 
351 
352 
366 template<typename MatrixType, typename inverse_type = typename MatrixType::value_type>
367 class PreconditionBlockJacobi : public virtual Subscriptor,
368  private PreconditionBlock<MatrixType, inverse_type>
369 {
370 private:
374  typedef typename MatrixType::value_type number;
375 
376 public:
381 
386  {
387  private:
391  class Accessor
392  {
393  public:
399  const size_type row);
400 
404  size_type row() const;
405 
409  size_type column() const;
410 
414  inverse_type value() const;
415 
416  protected:
421 
426 
431 
436 
441 
445  friend class const_iterator;
446  };
447 
448  public:
453  const size_type row);
454 
459 
464 
468  const Accessor &operator* () const;
469 
473  const Accessor *operator-> () const;
474 
478  bool operator == (const const_iterator &) const;
482  bool operator != (const const_iterator &) const;
483 
488  bool operator < (const const_iterator &) const;
489 
490  private:
495  };
496 
514 
522  template <typename number2>
523  void vmult (Vector<number2> &, const Vector<number2> &) const;
524 
528  template <typename number2>
529  void Tvmult (Vector<number2> &, const Vector<number2> &) const;
537  template <typename number2>
538  void vmult_add (Vector<number2> &, const Vector<number2> &) const;
539 
543  template <typename number2>
544  void Tvmult_add (Vector<number2> &, const Vector<number2> &) const;
545 
549  template <typename number2>
550  void step (Vector<number2> &dst, const Vector<number2> &rhs) const;
551 
555  template <typename number2>
556  void Tstep (Vector<number2> &dst, const Vector<number2> &rhs) const;
557 
561  const_iterator begin () const;
562 
566  const_iterator end () const;
567 
571  const_iterator begin (const size_type r) const;
572 
576  const_iterator end (const size_type r) const;
577 
578 
579 private:
586  template <typename number2>
587  void do_vmult (Vector<number2> &,
588  const Vector<number2> &,
589  bool adding) const;
590 
591  friend class Accessor;
592  friend class const_iterator;
593 };
594 
595 
596 
632 template<typename MatrixType, typename inverse_type = typename MatrixType::value_type>
633 class PreconditionBlockSOR : public virtual Subscriptor,
634  protected PreconditionBlock<MatrixType, inverse_type>
635 {
636 public:
641 
646 
650  typedef typename MatrixType::value_type number;
651 
668 
679  template <typename number2>
680  void vmult (Vector<number2> &, const Vector<number2> &) const;
681 
692  template <typename number2>
693  void vmult_add (Vector<number2> &, const Vector<number2> &) const;
694 
703  template <typename number2>
704  void Tvmult (Vector<number2> &, const Vector<number2> &) const;
705 
716  template <typename number2>
717  void Tvmult_add (Vector<number2> &, const Vector<number2> &) const;
718 
722  template <typename number2>
723  void step (Vector<number2> &dst, const Vector<number2> &rhs) const;
724 
728  template <typename number2>
729  void Tstep (Vector<number2> &dst, const Vector<number2> &rhs) const;
730 
731 protected:
735  PreconditionBlockSOR(bool store);
736 
746  template <typename number2>
747  void forward (Vector<number2> &,
748  const Vector<number2> &,
749  const bool transpose_diagonal,
750  const bool adding) const;
751 
761  template <typename number2>
762  void backward (Vector<number2> &,
763  const Vector<number2> &,
764  const bool transpose_diagonal,
765  const bool adding) const;
766 };
767 
768 
790 template<typename MatrixType, typename inverse_type = typename MatrixType::value_type>
791 class PreconditionBlockSSOR : public virtual Subscriptor,
792  private PreconditionBlockSOR<MatrixType, inverse_type>
793 {
794 public:
799 
803  typedef typename MatrixType::value_type number;
804 
809 
810  // Keep AdditionalData accessible
812 
813  // The following are the
814  // functions of the base classes
815  // which we want to keep
816  // accessible.
832 
840  template <typename number2>
841  void vmult (Vector<number2> &, const Vector<number2> &) const;
842 
846  template <typename number2>
847  void Tvmult (Vector<number2> &, const Vector<number2> &) const;
848 
852  template <typename number2>
853  void step (Vector<number2> &dst, const Vector<number2> &rhs) const;
854 
858  template <typename number2>
859  void Tstep (Vector<number2> &dst, const Vector<number2> &rhs) const;
860 };
861 
863 //---------------------------------------------------------------------------
864 
865 #ifndef DOXYGEN
866 
867 template<typename MatrixType, typename inverse_type>
868 inline bool
870 {
871  if (A == 0)
872  return true;
873  return A->empty();
874 }
875 
876 
877 template<typename MatrixType, typename inverse_type>
878 inline inverse_type
880  size_type i,
881  size_type j) const
882 {
883  const size_type bs = blocksize;
884  const unsigned int nb = i/bs;
885 
886  const FullMatrix<inverse_type> &B = this->inverse(nb);
887 
888  const size_type ib = i % bs;
889  const size_type jb = j % bs;
890 
891  if (jb + nb*bs != j)
892  {
893  return 0.;
894  }
895 
896  return B(ib, jb);
897 }
898 
899 //---------------------------------------------------------------------------
900 
901 template<typename MatrixType, typename inverse_type>
902 inline
905  const size_type row)
906  :
907  matrix(matrix),
908  b_iterator(&matrix->inverse(0), 0, 0),
909  b_end(&matrix->inverse(0), 0, 0)
910 {
911  bs = matrix->block_size();
912  a_block = row / bs;
913 
914  // This is the end accessor, which
915  // does not have a valid block.
916  if (a_block == matrix->size())
917  return;
918 
919  const size_type r = row % bs;
920 
921  b_iterator = matrix->inverse(a_block).begin(r);
922  b_end = matrix->inverse(a_block).end();
923 
924  Assert (a_block < matrix->size(),
925  ExcIndexRange(a_block, 0, matrix->size()));
926 }
927 
928 
929 template<typename MatrixType, typename inverse_type>
930 inline
933 {
934  Assert (a_block < matrix->size(),
936 
937  return bs * a_block + b_iterator->row();
938 }
939 
940 
941 template<typename MatrixType, typename inverse_type>
942 inline
945 {
946  Assert (a_block < matrix->size(),
948 
949  return bs * a_block + b_iterator->column();
950 }
951 
952 
953 template<typename MatrixType, typename inverse_type>
954 inline
955 inverse_type
957 {
958  Assert (a_block < matrix->size(),
960 
961  return b_iterator->value();
962 }
963 
964 
965 template<typename MatrixType, typename inverse_type>
966 inline
969  const size_type row)
970  :
971  accessor(matrix, row)
972 {}
973 
974 
975 template<typename MatrixType, typename inverse_type>
976 inline
979 {
980  Assert (*this != accessor.matrix->end(), ExcIteratorPastEnd());
981 
984  {
985  ++accessor.a_block;
986 
987  if (accessor.a_block < accessor.matrix->size())
988  {
989  accessor.b_iterator = accessor.matrix->inverse(accessor.a_block).begin();
990  accessor.b_end = accessor.matrix->inverse(accessor.a_block).end();
991  }
992  }
993  return *this;
994 }
995 
996 
997 template<typename MatrixType, typename inverse_type>
998 inline
1001 {
1002  return accessor;
1003 }
1004 
1005 
1006 template<typename MatrixType, typename inverse_type>
1007 inline
1010 {
1011  return &accessor;
1012 }
1013 
1014 
1015 template<typename MatrixType, typename inverse_type>
1016 inline
1017 bool
1019 operator == (const const_iterator &other) const
1020 {
1021  if (accessor.a_block == accessor.matrix->size() &&
1022  accessor.a_block == other.accessor.a_block)
1023  return true;
1024 
1025  if (accessor.a_block != other.accessor.a_block)
1026  return false;
1027 
1028  return (accessor.row() == other.accessor.row() &&
1029  accessor.column() == other.accessor.column());
1030 }
1031 
1032 
1033 template<typename MatrixType, typename inverse_type>
1034 inline
1035 bool
1037 operator != (const const_iterator &other) const
1038 {
1039  return ! (*this == other);
1040 }
1041 
1042 
1043 template<typename MatrixType, typename inverse_type>
1044 inline
1045 bool
1047 operator < (const const_iterator &other) const
1048 {
1049  return (accessor.row() < other.accessor.row() ||
1050  (accessor.row() == other.accessor.row() &&
1051  accessor.column() < other.accessor.column()));
1052 }
1053 
1054 
1055 template<typename MatrixType, typename inverse_type>
1056 inline
1059 {
1060  return const_iterator(this, 0);
1061 }
1062 
1063 
1064 template<typename MatrixType, typename inverse_type>
1065 inline
1068 {
1069  return const_iterator(this, this->size() * this->block_size());
1070 }
1071 
1072 
1073 template<typename MatrixType, typename inverse_type>
1074 inline
1077  const size_type r) const
1078 {
1079  Assert (r < this->A->m(), ExcIndexRange(r, 0, this->A->m()));
1080  return const_iterator(this, r);
1081 }
1082 
1083 
1084 
1085 template<typename MatrixType, typename inverse_type>
1086 inline
1089  const size_type r) const
1090 {
1091  Assert (r < this->A->m(), ExcIndexRange(r, 0, this->A->m()));
1092  return const_iterator(this, r+1);
1093 }
1094 
1095 #endif // DOXYGEN
1096 
1097 DEAL_II_NAMESPACE_CLOSE
1098 
1099 #endif
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)
void Tvmult_add(Vector< number2 > &, const Vector< number2 > &) const
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
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
Definition: types.h:88
types::global_dof_index size_type
#define Assert(cond, exc)
Definition: exceptions.h:313
value_type el(size_type i, size_type j) 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)
Definition: exceptions.h:541
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
void invert_diagblocks()
std::vector< size_type > permutation
void vmult(Vector< number2 > &, const Vector< number2 > &) const
bool empty() 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
void invert_permuted_diagblocks()
const PreconditionBlockJacobi< MatrixType, inverse_type > * matrix
inverse_type value_type
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