Reference documentation for deal.II version 8.4.2
precondition.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_h
17 #define dealii__precondition_h
18 
19 // This file contains simple preconditioners.
20 
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/smartpointer.h>
23 #include <deal.II/base/utilities.h>
24 #include <deal.II/base/parallel.h>
25 #include <deal.II/base/template_constraints.h>
26 #include <deal.II/lac/tridiagonal_matrix.h>
27 #include <deal.II/lac/solver_cg.h>
28 #include <deal.II/lac/vector_memory.h>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 // forward declarations
33 
34 template <typename number> class Vector;
35 template <typename number> class SparseMatrix;
36 namespace parallel
37 {
38  namespace distributed
39  {
40  template <typename number> class Vector;
41  }
42 }
43 
44 
45 
73 {
74 public:
79 
85  {
90  };
91 
96 
101  template <typename MatrixType>
102  void initialize (const MatrixType &matrix,
103  const AdditionalData &additional_data = AdditionalData());
104 
108  template<class VectorType>
109  void vmult (VectorType &, const VectorType &) const;
110 
115  template<class VectorType>
116  void Tvmult (VectorType &, const VectorType &) const;
117 
121  template<class VectorType>
122  void vmult_add (VectorType &, const VectorType &) const;
123 
128  template<class VectorType>
129  void Tvmult_add (VectorType &, const VectorType &) const;
130 
135  void clear () {}
136 
144  size_type m () const;
145 
153  size_type n () const;
154 
155 private:
159  size_type n_rows;
160 
164  size_type n_columns;
165 };
166 
167 
168 
183 {
184 public:
189 
194  {
195  public:
200  AdditionalData (const double relaxation = 1.);
201 
205  double relaxation;
206  };
207 
213 
217  void initialize (const AdditionalData &parameters);
218 
224  template <typename MatrixType>
225  void initialize (const MatrixType &matrix,
226  const AdditionalData &parameters);
227 
231  template<class VectorType>
232  void vmult (VectorType &, const VectorType &) const;
233 
238  template<class VectorType>
239  void Tvmult (VectorType &, const VectorType &) const;
243  template<class VectorType>
244  void vmult_add (VectorType &, const VectorType &) const;
245 
250  template<class VectorType>
251  void Tvmult_add (VectorType &, const VectorType &) const;
252 
257  void clear () {}
258 
266  size_type m () const;
267 
275  size_type n () const;
276 
277 private:
281  double relaxation;
282 
286  size_type n_rows;
287 
291  size_type n_columns;
292 };
293 
294 
295 
336 template<typename MatrixType = SparseMatrix<double>, class VectorType = Vector<double> >
338 {
339 public:
343  typedef void ( MatrixType::* function_ptr)(VectorType &, const VectorType &) const;
344 
350  PreconditionUseMatrix(const MatrixType &M,
351  const function_ptr method);
352 
357  void vmult (VectorType &dst,
358  const VectorType &src) const;
359 
360 private:
364  const MatrixType &matrix;
365 
369  const function_ptr precondition;
370 };
371 
372 
373 
382 template<typename MatrixType = SparseMatrix<double> >
384 {
385 public:
389  typedef typename MatrixType::size_type size_type;
390 
395  {
396  public:
400  AdditionalData (const double relaxation = 1.);
401 
405  double relaxation;
406  };
407 
413  void initialize (const MatrixType &A,
414  const AdditionalData &parameters = AdditionalData());
415 
419  void clear();
420 
425  size_type m () const;
426 
431  size_type n () const;
432 
433 protected:
438 
442  double relaxation;
443 };
444 
445 
446 
474 template <typename MatrixType = SparseMatrix<double> >
475 class PreconditionJacobi : public PreconditionRelaxation<MatrixType>
476 {
477 public:
481  template<class VectorType>
482  void vmult (VectorType &, const VectorType &) const;
483 
488  template<class VectorType>
489  void Tvmult (VectorType &, const VectorType &) const;
490 
494  template<class VectorType>
495  void step (VectorType &x, const VectorType &rhs) const;
496 
500  template<class VectorType>
501  void Tstep (VectorType &x, const VectorType &rhs) const;
502 };
503 
504 
551 template <typename MatrixType = SparseMatrix<double> >
552 class PreconditionSOR : public PreconditionRelaxation<MatrixType>
553 {
554 public:
558  template<class VectorType>
559  void vmult (VectorType &, const VectorType &) const;
560 
564  template<class VectorType>
565  void Tvmult (VectorType &, const VectorType &) const;
566 
570  template<class VectorType>
571  void step (VectorType &x, const VectorType &rhs) const;
572 
576  template<class VectorType>
577  void Tstep (VectorType &x, const VectorType &rhs) const;
578 };
579 
580 
581 
609 template <typename MatrixType = SparseMatrix<double> >
610 class PreconditionSSOR : public PreconditionRelaxation<MatrixType>
611 {
612 public:
616  typedef typename MatrixType::size_type size_type;
617 
622 
623 
629  void initialize (const MatrixType &A,
630  const typename BaseClass::AdditionalData &parameters = typename BaseClass::AdditionalData());
631 
635  template<class VectorType>
636  void vmult (VectorType &, const VectorType &) const;
637 
642  template<class VectorType>
643  void Tvmult (VectorType &, const VectorType &) const;
644 
645 
649  template<class VectorType>
650  void step (VectorType &x, const VectorType &rhs) const;
651 
655  template<class VectorType>
656  void Tstep (VectorType &x, const VectorType &rhs) const;
657 
658 private:
663  std::vector<std::size_t> pos_right_of_diagonal;
664 };
665 
666 
699 template <typename MatrixType = SparseMatrix<double> >
700 class PreconditionPSOR : public PreconditionRelaxation<MatrixType>
701 {
702 public:
706  typedef typename MatrixType::size_type size_type;
707 
712  {
713  public:
724  AdditionalData (const std::vector<size_type> &permutation,
725  const std::vector<size_type> &inverse_permutation,
727  &parameters = typename PreconditionRelaxation<MatrixType>::AdditionalData());
728 
732  const std::vector<size_type> &permutation;
736  const std::vector<size_type> &inverse_permutation;
741  };
742 
754  void initialize (const MatrixType &A,
755  const std::vector<size_type> &permutation,
756  const std::vector<size_type> &inverse_permutation,
759 
770  void initialize (const MatrixType &A,
771  const AdditionalData &additional_data);
772 
776  template<class VectorType>
777  void vmult (VectorType &, const VectorType &) const;
778 
782  template<class VectorType>
783  void Tvmult (VectorType &, const VectorType &) const;
784 private:
788  const std::vector<size_type> *permutation;
792  const std::vector<size_type> *inverse_permutation;
793 };
794 
795 
796 
811 template <typename MatrixType=SparseMatrix<double>, class VectorType=Vector<double> >
813 {
814 public:
819 
825  {
829  AdditionalData (const unsigned int degree = 0,
830  const double smoothing_range = 0.,
831  const bool nonzero_starting = false,
832  const unsigned int eig_cg_n_iterations = 8,
833  const double eig_cg_residual = 1e-2,
834  const double max_eigenvalue = 1);
835 
842  unsigned int degree;
843 
856 
867 
873  unsigned int eig_cg_n_iterations;
874 
880 
887 
892  };
893 
895 
907  void initialize (const MatrixType &matrix,
908  const AdditionalData &additional_data = AdditionalData());
909 
914  void vmult (VectorType &dst,
915  const VectorType &src) const;
916 
921  void Tvmult (VectorType &dst,
922  const VectorType &src) const;
923 
927  void clear ();
928 
933  size_type m () const;
934 
939  size_type n () const;
940 
941 private:
942 
947 
951  mutable VectorType update1;
952 
956  mutable VectorType update2;
957 
962 
966  double theta;
967 
972  double delta;
973 
978 };
979 
980 
981 
983 /* ---------------------------------- Inline functions ------------------- */
984 
985 #ifndef DOXYGEN
986 
987 inline
989  :
990  n_rows (0),
991  n_columns (0)
992 {}
993 
994 template <typename MatrixType>
995 inline void
996 PreconditionIdentity::initialize (const MatrixType &matrix,
998 {
999  n_rows = matrix.m();
1000  n_columns = matrix.n();
1001 }
1002 
1003 
1004 template<class VectorType>
1005 inline void
1006 PreconditionIdentity::vmult (VectorType &dst, const VectorType &src) const
1007 {
1008  dst = src;
1009 }
1010 
1011 
1012 
1013 template<class VectorType>
1014 inline void
1015 PreconditionIdentity::Tvmult (VectorType &dst, const VectorType &src) const
1016 {
1017  dst = src;
1018 }
1019 
1020 template<class VectorType>
1021 inline void
1022 PreconditionIdentity::vmult_add (VectorType &dst, const VectorType &src) const
1023 {
1024  dst.add(src);
1025 }
1026 
1027 
1028 
1029 template<class VectorType>
1030 inline void
1031 PreconditionIdentity::Tvmult_add (VectorType &dst, const VectorType &src) const
1032 {
1033  dst.add(src);
1034 }
1035 
1037 PreconditionIdentity::m () const
1038 {
1039  Assert(n_rows != 0, ExcNotInitialized());
1040  return n_rows;
1041 }
1042 
1044 PreconditionIdentity::n () const
1045 {
1046  Assert(n_columns != 0, ExcNotInitialized());
1047  return n_columns;
1048 }
1049 
1050 //---------------------------------------------------------------------------
1051 
1052 inline
1054  :
1055  relaxation(relaxation)
1056 {}
1057 
1058 
1059 inline
1061  :
1062  relaxation(0),
1063  n_rows (0),
1064  n_columns (0)
1065 {
1066  AdditionalData add_data;
1067  relaxation=add_data.relaxation;
1068 }
1069 
1070 
1071 
1072 inline void
1074 (const PreconditionRichardson::AdditionalData &parameters)
1075 {
1076  relaxation = parameters.relaxation;
1077 }
1078 
1079 
1080 
1081 template <typename MatrixType>
1082 inline void
1084 (const MatrixType &matrix,
1085  const PreconditionRichardson::AdditionalData &parameters)
1086 {
1087  relaxation = parameters.relaxation;
1088  n_rows = matrix.m();
1089  n_columns = matrix.n();
1090 }
1091 
1092 
1093 
1094 template<class VectorType>
1095 inline void
1096 PreconditionRichardson::vmult (VectorType &dst, const VectorType &src) const
1097 {
1098 #ifdef DEAL_II_WITH_CXX11
1099  static_assert(
1100  std::is_same<size_type, typename VectorType::size_type>::value,
1101  "PreconditionRichardson and VectorType must have the same size_type.");
1102 #endif // DEAL_II_WITH_CXX11
1103 
1104  dst.equ(relaxation,src);
1105 }
1106 
1107 
1108 
1109 template<class VectorType>
1110 inline void
1111 PreconditionRichardson::Tvmult (VectorType &dst, const VectorType &src) const
1112 {
1113 #ifdef DEAL_II_WITH_CXX11
1114  static_assert(
1115  std::is_same<size_type, typename VectorType::size_type>::value,
1116  "PreconditionRichardson and VectorType must have the same size_type.");
1117 #endif // DEAL_II_WITH_CXX11
1118 
1119  dst.equ(relaxation,src);
1120 }
1121 
1122 template<class VectorType>
1123 inline void
1124 PreconditionRichardson::vmult_add (VectorType &dst, const VectorType &src) const
1125 {
1126 #ifdef DEAL_II_WITH_CXX11
1127  static_assert(
1128  std::is_same<size_type, typename VectorType::size_type>::value,
1129  "PreconditionRichardson and VectorType must have the same size_type.");
1130 #endif // DEAL_II_WITH_CXX11
1131 
1132  dst.add(relaxation,src);
1133 }
1134 
1135 
1136 
1137 template<class VectorType>
1138 inline void
1139 PreconditionRichardson::Tvmult_add (VectorType &dst, const VectorType &src) const
1140 {
1141 #ifdef DEAL_II_WITH_CXX11
1142  static_assert(
1143  std::is_same<size_type, typename VectorType::size_type>::value,
1144  "PreconditionRichardson and VectorType must have the same size_type.");
1145 #endif // DEAL_II_WITH_CXX11
1146 
1147  dst.add(relaxation,src);
1148 }
1149 
1152 {
1153  Assert(n_rows != 0, ExcNotInitialized());
1154  return n_rows;
1155 }
1156 
1159 {
1160  Assert(n_columns != 0, ExcNotInitialized());
1161  return n_columns;
1162 }
1163 
1164 //---------------------------------------------------------------------------
1165 
1166 template <typename MatrixType>
1167 inline void
1169  const AdditionalData &parameters)
1170 {
1171  A = &rA;
1172  relaxation = parameters.relaxation;
1173 }
1174 
1175 
1176 template <typename MatrixType>
1177 inline void
1179 {
1180  A = 0;
1181 }
1182 
1183 template <typename MatrixType>
1186 {
1187  Assert (A!=0, ExcNotInitialized());
1188  return A->m();
1189 }
1190 
1191 template <typename MatrixType>
1194 {
1195  Assert (A!=0, ExcNotInitialized());
1196  return A->n();
1197 }
1198 
1199 //---------------------------------------------------------------------------
1200 
1201 template <typename MatrixType>
1202 template<class VectorType>
1203 inline void
1204 PreconditionJacobi<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1205 {
1206 #ifdef DEAL_II_WITH_CXX11
1207  static_assert(
1208  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1209  "PreconditionJacobi and VectorType must have the same size_type.");
1210 #endif // DEAL_II_WITH_CXX11
1211 
1212  Assert (this->A!=0, ExcNotInitialized());
1213  this->A->precondition_Jacobi (dst, src, this->relaxation);
1214 }
1215 
1216 
1217 
1218 template <typename MatrixType>
1219 template<class VectorType>
1220 inline void
1221 PreconditionJacobi<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1222 {
1223 #ifdef DEAL_II_WITH_CXX11
1224  static_assert(
1225  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1226  "PreconditionJacobi and VectorType must have the same size_type.");
1227 #endif // DEAL_II_WITH_CXX11
1228 
1229  Assert (this->A!=0, ExcNotInitialized());
1230  this->A->precondition_Jacobi (dst, src, this->relaxation);
1231 }
1232 
1233 
1234 
1235 template <typename MatrixType>
1236 template<class VectorType>
1237 inline void
1238 PreconditionJacobi<MatrixType>::step (VectorType &dst, const VectorType &src) const
1239 {
1240 #ifdef DEAL_II_WITH_CXX11
1241  static_assert(
1242  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1243  "PreconditionJacobi and VectorType must have the same size_type.");
1244 #endif // DEAL_II_WITH_CXX11
1245 
1246  Assert (this->A!=0, ExcNotInitialized());
1247  this->A->Jacobi_step (dst, src, this->relaxation);
1248 }
1249 
1250 
1251 
1252 template <typename MatrixType>
1253 template<class VectorType>
1254 inline void
1255 PreconditionJacobi<MatrixType>::Tstep (VectorType &dst, const VectorType &src) const
1256 {
1257 #ifdef DEAL_II_WITH_CXX11
1258  static_assert(
1259  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1260  "PreconditionJacobi and VectorType must have the same size_type.");
1261 #endif // DEAL_II_WITH_CXX11
1262 
1263  step (dst, src);
1264 }
1265 
1266 
1267 
1268 //---------------------------------------------------------------------------
1269 
1270 template <typename MatrixType>
1271 template<class VectorType>
1272 inline void
1273 PreconditionSOR<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1274 {
1275 #ifdef DEAL_II_WITH_CXX11
1276  static_assert(
1277  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1278  "PreconditionSOR and VectorType must have the same size_type.");
1279 #endif // DEAL_II_WITH_CXX11
1280 
1281  Assert (this->A!=0, ExcNotInitialized());
1282  this->A->precondition_SOR (dst, src, this->relaxation);
1283 }
1284 
1285 
1286 
1287 template <typename MatrixType>
1288 template<class VectorType>
1289 inline void
1290 PreconditionSOR<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1291 {
1292 #ifdef DEAL_II_WITH_CXX11
1293  static_assert(
1294  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1295  "PreconditionSOR and VectorType must have the same size_type.");
1296 #endif // DEAL_II_WITH_CXX11
1297 
1298  Assert (this->A!=0, ExcNotInitialized());
1299  this->A->precondition_TSOR (dst, src, this->relaxation);
1300 }
1301 
1302 
1303 
1304 template <typename MatrixType>
1305 template<class VectorType>
1306 inline void
1307 PreconditionSOR<MatrixType>::step (VectorType &dst, const VectorType &src) const
1308 {
1309 #ifdef DEAL_II_WITH_CXX11
1310  static_assert(
1311  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1312  "PreconditionSOR and VectorType must have the same size_type.");
1313 #endif // DEAL_II_WITH_CXX11
1314 
1315  Assert (this->A!=0, ExcNotInitialized());
1316  this->A->SOR_step (dst, src, this->relaxation);
1317 }
1318 
1319 
1320 
1321 template <typename MatrixType>
1322 template<class VectorType>
1323 inline void
1324 PreconditionSOR<MatrixType>::Tstep (VectorType &dst, const VectorType &src) const
1325 {
1326 #ifdef DEAL_II_WITH_CXX11
1327  static_assert(
1328  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1329  "PreconditionSOR and VectorType must have the same size_type.");
1330 #endif // DEAL_II_WITH_CXX11
1331 
1332  Assert (this->A!=0, ExcNotInitialized());
1333  this->A->TSOR_step (dst, src, this->relaxation);
1334 }
1335 
1336 
1337 
1338 //---------------------------------------------------------------------------
1339 
1340 template <typename MatrixType>
1341 inline void
1342 PreconditionSSOR<MatrixType>::initialize (const MatrixType &rA,
1343  const typename BaseClass::AdditionalData &parameters)
1344 {
1345  this->PreconditionRelaxation<MatrixType>::initialize (rA, parameters);
1346 
1347  // in case we have a SparseMatrix class, we can extract information about
1348  // the diagonal.
1350  dynamic_cast<const SparseMatrix<typename MatrixType::value_type> *>(&*this->A);
1351 
1352  // calculate the positions first after the diagonal.
1353  if (mat != 0)
1354  {
1355  const size_type n = this->A->n();
1356  pos_right_of_diagonal.resize(n, static_cast<std::size_t>(-1));
1357  for (size_type row=0; row<n; ++row)
1358  {
1359  // find the first element in this line which is on the right of the
1360  // diagonal. we need to precondition with the elements on the left
1361  // only. note: the first entry in each line denotes the diagonal
1362  // element, which we need not check.
1364  it = mat->begin(row)+1;
1365  for ( ; it < mat->end(row); ++it)
1366  if (it->column() > row)
1367  break;
1368  pos_right_of_diagonal[row] = it - mat->begin();
1369  }
1370  }
1371 }
1372 
1373 
1374 template <typename MatrixType>
1375 template<class VectorType>
1376 inline void
1377 PreconditionSSOR<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1378 {
1379 #ifdef DEAL_II_WITH_CXX11
1380  static_assert(
1381  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1382  "PreconditionSSOR and VectorType must have the same size_type.");
1383 #endif // DEAL_II_WITH_CXX11
1384 
1385  Assert (this->A!=0, ExcNotInitialized());
1386  this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1387 }
1388 
1389 
1390 
1391 template <typename MatrixType>
1392 template<class VectorType>
1393 inline void
1394 PreconditionSSOR<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1395 {
1396 #ifdef DEAL_II_WITH_CXX11
1397  static_assert(
1398  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1399  "PreconditionSSOR and VectorType must have the same size_type.");
1400 #endif // DEAL_II_WITH_CXX11
1401 
1402  Assert (this->A!=0, ExcNotInitialized());
1403  this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1404 }
1405 
1406 
1407 
1408 template <typename MatrixType>
1409 template<class VectorType>
1410 inline void
1411 PreconditionSSOR<MatrixType>::step (VectorType &dst, const VectorType &src) const
1412 {
1413 #ifdef DEAL_II_WITH_CXX11
1414  static_assert(
1415  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1416  "PreconditionSSOR and VectorType must have the same size_type.");
1417 #endif // DEAL_II_WITH_CXX11
1418 
1419  Assert (this->A!=0, ExcNotInitialized());
1420  this->A->SSOR_step (dst, src, this->relaxation);
1421 }
1422 
1423 
1424 
1425 template <typename MatrixType>
1426 template<class VectorType>
1427 inline void
1428 PreconditionSSOR<MatrixType>::Tstep (VectorType &dst, const VectorType &src) const
1429 {
1430 #ifdef DEAL_II_WITH_CXX11
1431  static_assert(
1432  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1433  "PreconditionSSOR and VectorType must have the same size_type.");
1434 #endif // DEAL_II_WITH_CXX11
1435 
1436  step (dst, src);
1437 }
1438 
1439 
1440 
1441 //---------------------------------------------------------------------------
1442 
1443 template <typename MatrixType>
1444 inline void
1446 (const MatrixType &rA,
1447  const std::vector<size_type> &p,
1448  const std::vector<size_type> &ip,
1449  const typename PreconditionRelaxation<MatrixType>::AdditionalData &parameters)
1450 {
1451  permutation = &p;
1452  inverse_permutation = &ip;
1454 }
1455 
1456 
1457 template <typename MatrixType>
1458 inline void
1459 PreconditionPSOR<MatrixType>::initialize (const MatrixType &A,
1460  const AdditionalData &additional_data)
1461 {
1462  initialize(A,
1463  additional_data.permutation,
1464  additional_data.inverse_permutation,
1465  additional_data.parameters);
1466 }
1467 
1468 
1469 template <typename MatrixType>
1470 template <typename VectorType>
1471 inline void
1472 PreconditionPSOR<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1473 {
1474 #ifdef DEAL_II_WITH_CXX11
1475  static_assert(
1476  std::is_same<typename PreconditionPSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1477  "PreconditionPSOR and VectorType must have the same size_type.");
1478 #endif // DEAL_II_WITH_CXX11
1479 
1480  Assert (this->A!=0, ExcNotInitialized());
1481  dst = src;
1482  this->A->PSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1483 }
1484 
1485 
1486 
1487 template <typename MatrixType>
1488 template<class VectorType>
1489 inline void
1490 PreconditionPSOR<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1491 {
1492 #ifdef DEAL_II_WITH_CXX11
1493  static_assert(
1494  std::is_same<typename PreconditionPSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1495  "PreconditionPSOR and VectorType must have the same size_type.");
1496 #endif // DEAL_II_WITH_CXX11
1497 
1498  Assert (this->A!=0, ExcNotInitialized());
1499  dst = src;
1500  this->A->TPSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1501 }
1502 
1503 template <typename MatrixType>
1505 (const std::vector<size_type> &permutation,
1506  const std::vector<size_type> &inverse_permutation,
1507  const typename PreconditionRelaxation<MatrixType>::AdditionalData &parameters)
1508  :
1509  permutation(permutation),
1510  inverse_permutation(inverse_permutation),
1511  parameters(parameters)
1512 {
1513 
1514 }
1515 
1516 
1517 //---------------------------------------------------------------------------
1518 
1519 
1520 template<typename MatrixType, class VectorType>
1522  const function_ptr method)
1523  :
1524  matrix(M), precondition(method)
1525 {}
1526 
1527 
1528 
1529 template<typename MatrixType, class VectorType>
1530 void
1532  const VectorType &src) const
1533 {
1534  (matrix.*precondition)(dst, src);
1535 }
1536 
1537 //---------------------------------------------------------------------------
1538 
1539 template<typename MatrixType>
1540 inline
1542 AdditionalData (const double relaxation)
1543  :
1544  relaxation (relaxation)
1545 {}
1546 
1547 
1548 
1549 //---------------------------------------------------------------------------
1550 
1551 namespace internal
1552 {
1553  namespace PreconditionChebyshev
1554  {
1555  // for deal.II vectors, perform updates for Chebyshev preconditioner all
1556  // at once to reduce memory transfer. Here, we select between general
1557  // vectors and deal.II vectors where we expand the loop over the (local)
1558  // size of the vector
1559 
1560  // generic part for non-deal.II vectors
1561  template <typename VectorType>
1562  inline
1563  void
1564  vector_updates (const VectorType &src,
1565  const VectorType &matrix_diagonal_inverse,
1566  const bool start_zero,
1567  const double factor1,
1568  const double factor2,
1569  VectorType &update1,
1570  VectorType &update2,
1571  VectorType &dst)
1572  {
1573  if (start_zero)
1574  {
1575  dst.equ (factor2, src);
1576  dst.scale (matrix_diagonal_inverse);
1577  update1.equ(-1.,dst);
1578  }
1579  else
1580  {
1581  update2 -= src;
1582  update2.scale (matrix_diagonal_inverse);
1583  if (factor1 == 0.)
1584  update1.equ(factor2, update2);
1585  else
1586  update1.sadd(factor1, factor2, update2);
1587  dst -= update1;
1588  }
1589  }
1590 
1591  // worker routine for deal.II vectors. Because of vectorization, we need
1592  // to put the loop into an extra structure because the virtual function of
1593  // VectorUpdatesRange prevents the compiler from applying vectorization.
1594  template <typename Number>
1595  struct VectorUpdater
1596  {
1597  VectorUpdater (const Number *src,
1598  const Number *matrix_diagonal_inverse,
1599  const bool start_zero,
1600  const Number factor1,
1601  const Number factor2,
1602  Number *update1,
1603  Number *update2,
1604  Number *dst)
1605  :
1606  src (src),
1607  matrix_diagonal_inverse (matrix_diagonal_inverse),
1608  do_startup (factor1 == Number()),
1609  start_zero (start_zero),
1610  factor1 (factor1),
1611  factor2 (factor2),
1612  update1 (update1),
1613  update2 (update2),
1614  dst (dst)
1615  {}
1616 
1617  void
1618  apply_to_subrange (const std::size_t begin,
1619  const std::size_t end) const
1620  {
1621  // To circumvent a bug in gcc
1622  // (https://gcc.gnu.org/bugzilla/show_bug.cgi?id=63945), we create copies
1623  // of the variables factor1 and factor2 and do not check based on
1624  // factor1.
1625  const Number factor1 = this->factor1;
1626  const Number factor2 = this->factor2;
1627  if (do_startup)
1628  {
1629  if (start_zero)
1630  DEAL_II_OPENMP_SIMD_PRAGMA
1631  for (std::size_t i=begin; i<end; ++i)
1632  {
1633  dst[i] = factor2 * src[i] * matrix_diagonal_inverse[i];
1634  update1[i] = -dst[i];
1635  }
1636  else
1637  DEAL_II_OPENMP_SIMD_PRAGMA
1638  for (std::size_t i=begin; i<end; ++i)
1639  {
1640  update1[i] = ((update2[i]-src[i]) *
1641  factor2*matrix_diagonal_inverse[i]);
1642  dst[i] -= update1[i];
1643  }
1644  }
1645  else
1646  DEAL_II_OPENMP_SIMD_PRAGMA
1647  for (std::size_t i=begin; i<end; ++i)
1648  {
1649  const Number update =
1650  factor1 * update1[i] + factor2 *
1651  ((update2[i] - src[i]) * matrix_diagonal_inverse[i]);
1652  update1[i] = update;
1653  dst[i] -= update;
1654  }
1655  }
1656 
1657  const Number *src;
1658  const Number *matrix_diagonal_inverse;
1659  const bool do_startup;
1660  const bool start_zero;
1661  const Number factor1;
1662  const Number factor2;
1663  mutable Number *update1;
1664  mutable Number *update2;
1665  mutable Number *dst;
1666  };
1667 
1668  template<typename Number>
1669  struct VectorUpdatesRange : public parallel::ParallelForInteger
1670  {
1671  VectorUpdatesRange(const VectorUpdater<Number> &updater,
1672  const std::size_t size)
1673  :
1674  updater (updater)
1675  {
1676  if (size < internal::Vector::minimum_parallel_grain_size)
1677  apply_to_subrange (0, size);
1678  else
1679  apply_parallel (0, size,
1680  internal::Vector::minimum_parallel_grain_size);
1681  }
1682 
1683  ~VectorUpdatesRange() {}
1684 
1685  virtual void
1686  apply_to_subrange (const std::size_t begin,
1687  const std::size_t end) const
1688  {
1689  updater.apply_to_subrange(begin, end);
1690  }
1691 
1692  const VectorUpdater<Number> &updater;
1693  };
1694 
1695  // selection for deal.II vector
1696  template <typename Number>
1697  inline
1698  void
1699  vector_updates (const ::Vector<Number> &src,
1700  const ::Vector<Number> &matrix_diagonal_inverse,
1701  const bool start_zero,
1702  const double factor1,
1703  const double factor2,
1704  ::Vector<Number> &update1,
1705  ::Vector<Number> &update2,
1706  ::Vector<Number> &dst)
1707  {
1708  VectorUpdater<Number> upd(src.begin(), matrix_diagonal_inverse.begin(),
1709  start_zero, factor1, factor2,
1710  update1.begin(), update2.begin(), dst.begin());
1711  VectorUpdatesRange<Number>(upd, src.size());
1712  }
1713 
1714  // selection for parallel deal.II vector
1715  template <typename Number>
1716  inline
1717  void
1718  vector_updates (const parallel::distributed::Vector<Number> &src,
1719  const parallel::distributed::Vector<Number> &matrix_diagonal_inverse,
1720  const bool start_zero,
1721  const double factor1,
1722  const double factor2,
1726  {
1727  VectorUpdater<Number> upd(src.begin(), matrix_diagonal_inverse.begin(),
1728  start_zero, factor1, factor2,
1729  update1.begin(), update2.begin(), dst.begin());
1730  VectorUpdatesRange<Number>(upd, src.local_size());
1731  }
1732 
1733  template <typename VectorType>
1734  struct DiagonalPreconditioner
1735  {
1736  DiagonalPreconditioner (const VectorType &vector)
1737  :
1738  diagonal_vector(vector)
1739  {}
1740 
1741  void vmult (VectorType &dst,
1742  const VectorType &src) const
1743  {
1744  dst = src;
1745  dst.scale(diagonal_vector);
1746  }
1747 
1748  const VectorType &diagonal_vector;
1749  };
1750 
1751  struct EigenvalueTracker
1752  {
1753  public:
1754  void slot(const std::vector<double> &eigenvalues)
1755  {
1756  values = eigenvalues;
1757  }
1758 
1759  std::vector<double> values;
1760  };
1761  }
1762 }
1763 
1764 
1765 
1766 template <typename MatrixType, class VectorType>
1767 inline
1769 AdditionalData (const unsigned int degree,
1770  const double smoothing_range,
1771  const bool nonzero_starting,
1772  const unsigned int eig_cg_n_iterations,
1773  const double eig_cg_residual,
1774  const double max_eigenvalue)
1775  :
1776  degree (degree),
1777  smoothing_range (smoothing_range),
1778  nonzero_starting (nonzero_starting),
1779  eig_cg_n_iterations (eig_cg_n_iterations),
1780  eig_cg_residual (eig_cg_residual),
1781  max_eigenvalue (max_eigenvalue)
1782 {}
1783 
1784 
1785 
1786 template <typename MatrixType, class VectorType>
1787 inline
1789  :
1790  is_initialized (false)
1791 {
1792 #ifdef DEAL_II_WITH_CXX11
1793  static_assert(
1794  std::is_same<size_type, typename VectorType::size_type>::value,
1795  "PreconditionChebyshev and VectorType must have the same size_type.");
1796 #endif // DEAL_II_WITH_CXX11
1797 }
1798 
1799 
1800 
1801 template <typename MatrixType, class VectorType>
1802 inline
1803 void
1805 (const MatrixType &matrix,
1806  const AdditionalData &additional_data)
1807 {
1808  matrix_ptr = &matrix;
1809  data = additional_data;
1810  if (data.matrix_diagonal_inverse.size() != matrix.m())
1811  {
1812  Assert(data.matrix_diagonal_inverse.size() == 0,
1813  ExcMessage("Matrix diagonal vector set but not sized correctly"));
1814  data.matrix_diagonal_inverse.reinit(matrix.m());
1815  for (unsigned int i=0; i<matrix.m(); ++i)
1816  data.matrix_diagonal_inverse(i) = 1./matrix.el(i,i);
1817  }
1818 
1819 
1820  // calculate largest eigenvalue using a hand-tuned CG iteration on the
1821  // matrix weighted by its diagonal. we start with a vector that consists of
1822  // ones only, weighted by the length.
1823  double max_eigenvalue, min_eigenvalue;
1824  if (data.eig_cg_n_iterations > 0)
1825  {
1826  Assert (additional_data.eig_cg_n_iterations > 2,
1827  ExcMessage ("Need to set at least two iterations to find eigenvalues."));
1828 
1829  // set a very strict tolerance to force at least two iterations
1830  ReductionControl control (data.eig_cg_n_iterations, 1e-35, 1e-10);
1832  VectorType *rhs = memory.alloc();
1833  VectorType *dummy = memory.alloc();
1834  rhs->reinit(data.matrix_diagonal_inverse);
1835  dummy->reinit(data.matrix_diagonal_inverse);
1836 
1837  // heuristically, a right hand side close to a constant has been shown
1838  // to quickly reveal the largest eigenvalue. however, avoid to use the
1839  // exact constant because that might be not in the range space of some
1840  // matrices (purely Neumann matrices with constant mode filtered out by
1841  // orthogonal projection in the matrix-vector product)
1842  *rhs = 1./std::sqrt(static_cast<double>(matrix.m()));
1843  if (rhs->locally_owned_elements().is_element(0))
1844  (*rhs)(0) = 0.;
1845  rhs->compress(VectorOperation::insert);
1846 
1847  internal::PreconditionChebyshev::EigenvalueTracker eigenvalue_tracker;
1848  SolverCG<VectorType> solver (control, memory);
1849  solver.connect_eigenvalues_slot(std_cxx11::bind(&internal::PreconditionChebyshev::EigenvalueTracker::slot,
1850  &eigenvalue_tracker,
1851  std_cxx11::_1));
1852  internal::PreconditionChebyshev::DiagonalPreconditioner<VectorType>
1853  preconditioner(data.matrix_diagonal_inverse);
1854  try
1855  {
1856  solver.solve(matrix, *dummy, *rhs, preconditioner);
1857  }
1859  {
1860  }
1861 
1862  memory.free(dummy);
1863  memory.free(rhs);
1864 
1865  // read the eigenvalues from the attached eigenvalue tracker
1866  if (eigenvalue_tracker.values.empty())
1867  min_eigenvalue = max_eigenvalue = 1;
1868  else
1869  {
1870  min_eigenvalue = eigenvalue_tracker.values.front();
1871  max_eigenvalue = eigenvalue_tracker.values.back();
1872  }
1873 
1874  // include a safety factor since the CG method will in general not be
1875  // converged
1876  max_eigenvalue *= 1.2;
1877  }
1878  else
1879  {
1880  max_eigenvalue = data.max_eigenvalue;
1881  min_eigenvalue = data.max_eigenvalue/data.smoothing_range;
1882  }
1883 
1884  const double alpha = (data.smoothing_range > 1. ?
1885  max_eigenvalue / data.smoothing_range :
1886  std::min(0.9*max_eigenvalue, min_eigenvalue));
1887  delta = (max_eigenvalue-alpha)*0.5;
1888  theta = (max_eigenvalue+alpha)*0.5;
1889 
1890  update1.reinit (data.matrix_diagonal_inverse, true);
1891  update2.reinit (data.matrix_diagonal_inverse, true);
1892 
1893  is_initialized = true;
1894 }
1895 
1896 
1897 
1898 template <typename MatrixType, class VectorType>
1899 inline
1900 void
1902  const VectorType &src) const
1903 {
1904  Assert (is_initialized, ExcMessage("Preconditioner not initialized"));
1905  double rhok = delta / theta, sigma = theta / delta;
1906  if (data.nonzero_starting && !dst.all_zero())
1907  {
1908  matrix_ptr->vmult (update2, dst);
1909  internal::PreconditionChebyshev::vector_updates
1910  (src, data.matrix_diagonal_inverse, false, 0., 1./theta, update1,
1911  update2, dst);
1912  }
1913  else
1914  internal::PreconditionChebyshev::vector_updates
1915  (src, data.matrix_diagonal_inverse, true, 0., 1./theta, update1,
1916  update2, dst);
1917 
1918  for (unsigned int k=0; k<data.degree; ++k)
1919  {
1920  matrix_ptr->vmult (update2, dst);
1921  const double rhokp = 1./(2.*sigma-rhok);
1922  const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
1923  rhok = rhokp;
1924  internal::PreconditionChebyshev::vector_updates
1925  (src, data.matrix_diagonal_inverse, false, factor1, factor2, update1,
1926  update2, dst);
1927  }
1928 }
1929 
1930 
1931 
1932 template <typename MatrixType, class VectorType>
1933 inline
1934 void
1936  const VectorType &src) const
1937 {
1938  Assert (is_initialized, ExcMessage("Preconditioner not initialized"));
1939  double rhok = delta / theta, sigma = theta / delta;
1940  if (data.nonzero_starting && !dst.all_zero())
1941  {
1942  matrix_ptr->Tvmult (update2, dst);
1943  internal::PreconditionChebyshev::vector_updates
1944  (src, data.matrix_diagonal_inverse, false, 0., 1./theta, update1,
1945  update2, dst);
1946  }
1947  else
1948  internal::PreconditionChebyshev::vector_updates
1949  (src, data.matrix_diagonal_inverse, true, 0., 1./theta, update1,
1950  update2, dst);
1951 
1952  for (unsigned int k=0; k<data.degree; ++k)
1953  {
1954  matrix_ptr->Tvmult (update2, dst);
1955  const double rhokp = 1./(2.*sigma-rhok);
1956  const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
1957  rhok = rhokp;
1958  internal::PreconditionChebyshev::vector_updates
1959  (src, data.matrix_diagonal_inverse, false, factor1, factor2, update1,
1960  update2, dst);
1961  }
1962 }
1963 
1964 
1965 
1966 template <typename MatrixType, typename VectorType>
1967 inline
1969 {
1970  is_initialized = false;
1971  matrix_ptr = 0;
1972  data.matrix_diagonal_inverse.reinit(0);
1973  update1.reinit(0);
1974  update2.reinit(0);
1975 }
1976 
1977 
1978 template <typename MatrixType, typename VectorType>
1979 inline
1982 {
1983  Assert (matrix_ptr!=0, ExcNotInitialized());
1984  return matrix_ptr->m();
1985 }
1986 
1987 
1988 template <typename MatrixType, typename VectorType>
1989 inline
1992 {
1993  Assert (matrix_ptr!=0, ExcNotInitialized());
1994  return matrix_ptr->n();
1995 }
1996 
1997 #endif // DOXYGEN
1998 
1999 DEAL_II_NAMESPACE_CLOSE
2000 
2001 #endif
void Tvmult(VectorType &, const VectorType &) const
SmartPointer< const MatrixType, PreconditionChebyshev< MatrixType, VectorType > > matrix_ptr
Definition: precondition.h:946
AdditionalData data
Definition: precondition.h:961
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
void initialize(const MatrixType &A, const typename BaseClass::AdditionalData &parameters=typename BaseClass::AdditionalData())
const_iterator end() const
void Tstep(VectorType &x, const VectorType &rhs) const
boost::signals2::connection connect_eigenvalues_slot(const std_cxx11::function< void(const std::vector< double > &)> &slot, const bool every_iteration=false)
void vmult(VectorType &dst, const VectorType &src) const
MatrixType::size_type size_type
Definition: precondition.h:616
void Tvmult_add(VectorType &, const VectorType &) const
AdditionalData(const unsigned int degree=0, const double smoothing_range=0., const bool nonzero_starting=false, const unsigned int eig_cg_n_iterations=8, const double eig_cg_residual=1e-2, const double max_eigenvalue=1)
::ExceptionBase & ExcMessage(std::string arg1)
PreconditionRelaxation< MatrixType >::AdditionalData parameters
Definition: precondition.h:740
void vmult(VectorType &, const VectorType &) const
size_type m() const
const MatrixType & matrix
Definition: precondition.h:364
size_type n() const
AdditionalData(const double relaxation=1.)
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const double relaxation=1.)
size_type m() const
types::global_dof_index size_type
Definition: precondition.h:188
types::global_dof_index size_type
Definition: precondition.h:818
const std::vector< size_type > * permutation
Definition: precondition.h:788
const function_ptr precondition
Definition: precondition.h:369
void Tvmult(VectorType &dst, const VectorType &src) const
void vmult(VectorType &, const VectorType &) const
void initialize(const AdditionalData &parameters)
void Tvmult(VectorType &, const VectorType &) const
size_type m() const
types::global_dof_index size_type
Definition: precondition.h:78
size_type n() const
void vmult(VectorType &, const VectorType &) const
unsigned int global_dof_index
Definition: types.h:88
virtual VectorType * alloc()
const std::vector< size_type > & inverse_permutation
Definition: precondition.h:736
#define Assert(cond, exc)
Definition: exceptions.h:294
void step(VectorType &x, const VectorType &rhs) const
MatrixType::size_type size_type
Definition: precondition.h:706
virtual void free(const VectorType *const)
void vmult_add(VectorType &, const VectorType &) const
void Tvmult(VectorType &, const VectorType &) const
size_type n() const
void initialize(const MatrixType &A, const AdditionalData &parameters=AdditionalData())
void Tvmult(VectorType &, const VectorType &) const
iterator begin()
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
MatrixType::size_type size_type
Definition: precondition.h:389
PreconditionRelaxation< MatrixType > BaseClass
Definition: precondition.h:621
void Tstep(VectorType &x, const VectorType &rhs) const
void vmult(VectorType &, const VectorType &) const
void Tstep(VectorType &x, const VectorType &rhs) const
AdditionalData(const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename PreconditionRelaxation< MatrixType >::AdditionalData &parameters=typename PreconditionRelaxation< MatrixType >::AdditionalData())
void Tvmult_add(VectorType &, const VectorType &) const
void Tvmult(VectorType &, const VectorType &) const
void vmult(VectorType &, const VectorType &) const
size_type n() const
SmartPointer< const MatrixType, PreconditionRelaxation< MatrixType > > A
Definition: precondition.h:437
void Tvmult(VectorType &, const VectorType &) const
void step(VectorType &x, const VectorType &rhs) const
const std::vector< size_type > * inverse_permutation
Definition: precondition.h:792
void step(VectorType &x, const VectorType &rhs) const
size_type m() const
void vmult(VectorType &dst, const VectorType &src) const
const_iterator begin() const
void vmult(VectorType &, const VectorType &) const
const std::vector< size_type > & permutation
Definition: precondition.h:732
PreconditionUseMatrix(const MatrixType &M, const function_ptr method)
void vmult_add(VectorType &, const VectorType &) const
size_type local_size() const
std::vector< std::size_t > pos_right_of_diagonal
Definition: precondition.h:663
void initialize(const MatrixType &A, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename PreconditionRelaxation< MatrixType >::AdditionalData &parameters=typename PreconditionRelaxation< MatrixType >::AdditionalData())