Reference documentation for deal.II version 8.5.1
precondition.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2017 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/base/thread_management.h>
27 #include <deal.II/lac/diagonal_matrix.h>
28 #include <deal.II/lac/solver_cg.h>
29 #include <deal.II/lac/vector_memory.h>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 // forward declarations
34 
35 template <typename number> class Vector;
36 template <typename number> class SparseMatrix;
37 namespace LinearAlgebra
38 {
39  namespace distributed
40  {
41  template <typename number> class Vector;
42  }
43 }
44 
45 
46 
74 {
75 public:
80 
86  {
91  };
92 
97 
102  template <typename MatrixType>
103  void initialize (const MatrixType &matrix,
104  const AdditionalData &additional_data = AdditionalData());
105 
109  template<class VectorType>
110  void vmult (VectorType &, const VectorType &) const;
111 
116  template<class VectorType>
117  void Tvmult (VectorType &, const VectorType &) const;
118 
122  template<class VectorType>
123  void vmult_add (VectorType &, const VectorType &) const;
124 
129  template<class VectorType>
130  void Tvmult_add (VectorType &, const VectorType &) const;
131 
136  void clear () {}
137 
145  size_type m () const;
146 
154  size_type n () const;
155 
156 private:
161 
166 };
167 
168 
169 
184 {
185 public:
190 
195  {
196  public:
201  AdditionalData (const double relaxation = 1.);
202 
206  double relaxation;
207  };
208 
214 
218  void initialize (const AdditionalData &parameters);
219 
225  template <typename MatrixType>
226  void initialize (const MatrixType &matrix,
227  const AdditionalData &parameters);
228 
232  template<class VectorType>
233  void vmult (VectorType &, const VectorType &) const;
234 
239  template<class VectorType>
240  void Tvmult (VectorType &, const VectorType &) const;
244  template<class VectorType>
245  void vmult_add (VectorType &, const VectorType &) const;
246 
251  template<class VectorType>
252  void Tvmult_add (VectorType &, const VectorType &) const;
253 
258  void clear () {}
259 
267  size_type m () const;
268 
276  size_type n () const;
277 
278 private:
282  double relaxation;
283 
288 
293 };
294 
295 
296 
337 template<typename MatrixType = SparseMatrix<double>, class VectorType = Vector<double> >
339 {
340 public:
344  typedef void ( MatrixType::* function_ptr)(VectorType &, const VectorType &) const;
345 
351  PreconditionUseMatrix(const MatrixType &M,
352  const function_ptr method);
353 
358  void vmult (VectorType &dst,
359  const VectorType &src) const;
360 
361 private:
365  const MatrixType &matrix;
366 
371 };
372 
373 
374 
383 template<typename MatrixType = SparseMatrix<double> >
385 {
386 public:
390  typedef typename MatrixType::size_type size_type;
391 
396  {
397  public:
401  AdditionalData (const double relaxation = 1.);
402 
406  double relaxation;
407  };
408 
414  void initialize (const MatrixType &A,
415  const AdditionalData &parameters = AdditionalData());
416 
420  void clear();
421 
426  size_type m () const;
427 
432  size_type n () const;
433 
434 protected:
439 
443  double relaxation;
444 };
445 
446 
447 
475 template <typename MatrixType = SparseMatrix<double> >
476 class PreconditionJacobi : public PreconditionRelaxation<MatrixType>
477 {
478 public:
483 
487  template<class VectorType>
488  void vmult (VectorType &, const VectorType &) const;
489 
494  template<class VectorType>
495  void Tvmult (VectorType &, const VectorType &) const;
496 
500  template<class VectorType>
501  void step (VectorType &x, const VectorType &rhs) const;
502 
506  template<class VectorType>
507  void Tstep (VectorType &x, const VectorType &rhs) const;
508 };
509 
510 
557 template <typename MatrixType = SparseMatrix<double> >
558 class PreconditionSOR : public PreconditionRelaxation<MatrixType>
559 {
560 public:
565 
569  template<class VectorType>
570  void vmult (VectorType &, const VectorType &) const;
571 
575  template<class VectorType>
576  void Tvmult (VectorType &, const VectorType &) const;
577 
581  template<class VectorType>
582  void step (VectorType &x, const VectorType &rhs) const;
583 
587  template<class VectorType>
588  void Tstep (VectorType &x, const VectorType &rhs) const;
589 };
590 
591 
592 
620 template <typename MatrixType = SparseMatrix<double> >
621 class PreconditionSSOR : public PreconditionRelaxation<MatrixType>
622 {
623 public:
628 
632  typedef typename MatrixType::size_type size_type;
633 
638 
639 
645  void initialize (const MatrixType &A,
646  const typename BaseClass::AdditionalData &parameters = typename BaseClass::AdditionalData());
647 
651  template<class VectorType>
652  void vmult (VectorType &, const VectorType &) const;
653 
658  template<class VectorType>
659  void Tvmult (VectorType &, const VectorType &) const;
660 
661 
665  template<class VectorType>
666  void step (VectorType &x, const VectorType &rhs) const;
667 
671  template<class VectorType>
672  void Tstep (VectorType &x, const VectorType &rhs) const;
673 
674 private:
679  std::vector<std::size_t> pos_right_of_diagonal;
680 };
681 
682 
715 template <typename MatrixType = SparseMatrix<double> >
716 class PreconditionPSOR : public PreconditionRelaxation<MatrixType>
717 {
718 public:
722  typedef typename MatrixType::size_type size_type;
723 
728  {
729  public:
740  AdditionalData (const std::vector<size_type> &permutation,
741  const std::vector<size_type> &inverse_permutation,
744 
748  const std::vector<size_type> &permutation;
752  const std::vector<size_type> &inverse_permutation;
757  };
758 
770  void initialize (const MatrixType &A,
771  const std::vector<size_type> &permutation,
772  const std::vector<size_type> &inverse_permutation,
775 
786  void initialize (const MatrixType &A,
787  const AdditionalData &additional_data);
788 
792  template<class VectorType>
793  void vmult (VectorType &, const VectorType &) const;
794 
798  template<class VectorType>
799  void Tvmult (VectorType &, const VectorType &) const;
800 private:
804  const std::vector<size_type> *permutation;
808  const std::vector<size_type> *inverse_permutation;
809 };
810 
811 
812 
901 template <typename MatrixType=SparseMatrix<double>, typename VectorType=Vector<double>, typename PreconditionerType=DiagonalMatrix<VectorType> >
903 {
904 public:
909 
910  // avoid warning about use of deprecated variables
911  DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
912 
918  {
922  AdditionalData (const unsigned int degree = 0,
923  const double smoothing_range = 0.,
924  const bool nonzero_starting = false,
925  const unsigned int eig_cg_n_iterations = 8,
926  const double eig_cg_residual = 1e-2,
927  const double max_eigenvalue = 1);
928 
940  unsigned int degree;
941 
954 
967  bool nonzero_starting DEAL_II_DEPRECATED;
968 
974  unsigned int eig_cg_n_iterations;
975 
981 
988 
994  VectorType matrix_diagonal_inverse DEAL_II_DEPRECATED;
995 
999  std_cxx11::shared_ptr<PreconditionerType> preconditioner;
1000  };
1001 
1002  DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
1003 
1005 
1017  void initialize (const MatrixType &matrix,
1018  const AdditionalData &additional_data = AdditionalData());
1019 
1024  void vmult (VectorType &dst,
1025  const VectorType &src) const;
1026 
1031  void Tvmult (VectorType &dst,
1032  const VectorType &src) const;
1033 
1037  void step (VectorType &dst,
1038  const VectorType &src) const;
1039 
1043  void Tstep (VectorType &dst,
1044  const VectorType &src) const;
1045 
1049  void clear ();
1050 
1055  size_type m () const;
1056 
1061  size_type n () const;
1062 
1063 private:
1064 
1069 
1073  mutable VectorType update1;
1074 
1078  mutable VectorType update2;
1079 
1083  mutable VectorType update3;
1084 
1090 
1094  double theta;
1095 
1100  double delta;
1101 
1107 
1113 
1118  void do_chebyshev_loop(VectorType &dst,
1119  const VectorType &src) const;
1120 
1127  void do_transpose_chebyshev_loop(VectorType &dst,
1128  const VectorType &src) const;
1129 
1136  void estimate_eigenvalues(const VectorType &src) const;
1137 };
1138 
1139 
1140 
1142 /* ---------------------------------- Inline functions ------------------- */
1143 
1144 #ifndef DOXYGEN
1145 
1146 inline
1148  :
1149  n_rows (0),
1150  n_columns (0)
1151 {}
1152 
1153 template <typename MatrixType>
1154 inline void
1155 PreconditionIdentity::initialize (const MatrixType &matrix,
1157 {
1158  n_rows = matrix.m();
1159  n_columns = matrix.n();
1160 }
1161 
1162 
1163 template<class VectorType>
1164 inline void
1165 PreconditionIdentity::vmult (VectorType &dst, const VectorType &src) const
1166 {
1167  dst = src;
1168 }
1169 
1170 
1171 
1172 template<class VectorType>
1173 inline void
1174 PreconditionIdentity::Tvmult (VectorType &dst, const VectorType &src) const
1175 {
1176  dst = src;
1177 }
1178 
1179 template<class VectorType>
1180 inline void
1181 PreconditionIdentity::vmult_add (VectorType &dst, const VectorType &src) const
1182 {
1183  dst += src;
1184 }
1185 
1186 
1187 
1188 template<class VectorType>
1189 inline void
1190 PreconditionIdentity::Tvmult_add (VectorType &dst, const VectorType &src) const
1191 {
1192  dst += src;
1193 }
1194 
1196 PreconditionIdentity::m () const
1197 {
1198  Assert(n_rows != 0, ExcNotInitialized());
1199  return n_rows;
1200 }
1201 
1203 PreconditionIdentity::n () const
1204 {
1206  return n_columns;
1207 }
1208 
1209 //---------------------------------------------------------------------------
1210 
1211 inline
1213  :
1214  relaxation(relaxation)
1215 {}
1216 
1217 
1218 inline
1220  :
1221  relaxation(0),
1222  n_rows (0),
1223  n_columns (0)
1224 {
1225  AdditionalData add_data;
1226  relaxation=add_data.relaxation;
1227 }
1228 
1229 
1230 
1231 inline void
1233 (const PreconditionRichardson::AdditionalData &parameters)
1234 {
1235  relaxation = parameters.relaxation;
1236 }
1237 
1238 
1239 
1240 template <typename MatrixType>
1241 inline void
1243 (const MatrixType &matrix,
1244  const PreconditionRichardson::AdditionalData &parameters)
1245 {
1246  relaxation = parameters.relaxation;
1247  n_rows = matrix.m();
1248  n_columns = matrix.n();
1249 }
1250 
1251 
1252 
1253 template<class VectorType>
1254 inline void
1255 PreconditionRichardson::vmult (VectorType &dst, const VectorType &src) const
1256 {
1257 #ifdef DEAL_II_WITH_CXX11
1258  static_assert(
1259  std::is_same<size_type, typename VectorType::size_type>::value,
1260  "PreconditionRichardson and VectorType must have the same size_type.");
1261 #endif // DEAL_II_WITH_CXX11
1262 
1263  dst.equ(relaxation,src);
1264 }
1265 
1266 
1267 
1268 template<class VectorType>
1269 inline void
1270 PreconditionRichardson::Tvmult (VectorType &dst, const VectorType &src) const
1271 {
1272 #ifdef DEAL_II_WITH_CXX11
1273  static_assert(
1274  std::is_same<size_type, typename VectorType::size_type>::value,
1275  "PreconditionRichardson and VectorType must have the same size_type.");
1276 #endif // DEAL_II_WITH_CXX11
1277 
1278  dst.equ(relaxation,src);
1279 }
1280 
1281 template<class VectorType>
1282 inline void
1283 PreconditionRichardson::vmult_add (VectorType &dst, const VectorType &src) const
1284 {
1285 #ifdef DEAL_II_WITH_CXX11
1286  static_assert(
1287  std::is_same<size_type, typename VectorType::size_type>::value,
1288  "PreconditionRichardson and VectorType must have the same size_type.");
1289 #endif // DEAL_II_WITH_CXX11
1290 
1291  dst.add(relaxation,src);
1292 }
1293 
1294 
1295 
1296 template<class VectorType>
1297 inline void
1298 PreconditionRichardson::Tvmult_add (VectorType &dst, const VectorType &src) const
1299 {
1300 #ifdef DEAL_II_WITH_CXX11
1301  static_assert(
1302  std::is_same<size_type, typename VectorType::size_type>::value,
1303  "PreconditionRichardson and VectorType must have the same size_type.");
1304 #endif // DEAL_II_WITH_CXX11
1305 
1306  dst.add(relaxation,src);
1307 }
1308 
1311 {
1312  Assert(n_rows != 0, ExcNotInitialized());
1313  return n_rows;
1314 }
1315 
1318 {
1320  return n_columns;
1321 }
1322 
1323 //---------------------------------------------------------------------------
1324 
1325 template <typename MatrixType>
1326 inline void
1328  const AdditionalData &parameters)
1329 {
1330  A = &rA;
1331  relaxation = parameters.relaxation;
1332 }
1333 
1334 
1335 template <typename MatrixType>
1336 inline void
1338 {
1339  A = 0;
1340 }
1341 
1342 template <typename MatrixType>
1345 {
1346  Assert (A!=0, ExcNotInitialized());
1347  return A->m();
1348 }
1349 
1350 template <typename MatrixType>
1353 {
1354  Assert (A!=0, ExcNotInitialized());
1355  return A->n();
1356 }
1357 
1358 //---------------------------------------------------------------------------
1359 
1360 template <typename MatrixType>
1361 template<class VectorType>
1362 inline void
1363 PreconditionJacobi<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1364 {
1365 #ifdef DEAL_II_WITH_CXX11
1366  static_assert(
1367  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1368  "PreconditionJacobi and VectorType must have the same size_type.");
1369 #endif // DEAL_II_WITH_CXX11
1370 
1371  Assert (this->A!=0, ExcNotInitialized());
1372  this->A->precondition_Jacobi (dst, src, this->relaxation);
1373 }
1374 
1375 
1376 
1377 template <typename MatrixType>
1378 template<class VectorType>
1379 inline void
1380 PreconditionJacobi<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1381 {
1382 #ifdef DEAL_II_WITH_CXX11
1383  static_assert(
1384  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1385  "PreconditionJacobi and VectorType must have the same size_type.");
1386 #endif // DEAL_II_WITH_CXX11
1387 
1388  Assert (this->A!=0, ExcNotInitialized());
1389  this->A->precondition_Jacobi (dst, src, this->relaxation);
1390 }
1391 
1392 
1393 
1394 template <typename MatrixType>
1395 template<class VectorType>
1396 inline void
1397 PreconditionJacobi<MatrixType>::step (VectorType &dst, const VectorType &src) const
1398 {
1399 #ifdef DEAL_II_WITH_CXX11
1400  static_assert(
1401  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1402  "PreconditionJacobi and VectorType must have the same size_type.");
1403 #endif // DEAL_II_WITH_CXX11
1404 
1405  Assert (this->A!=0, ExcNotInitialized());
1406  this->A->Jacobi_step (dst, src, this->relaxation);
1407 }
1408 
1409 
1410 
1411 template <typename MatrixType>
1412 template<class VectorType>
1413 inline void
1414 PreconditionJacobi<MatrixType>::Tstep (VectorType &dst, const VectorType &src) const
1415 {
1416 #ifdef DEAL_II_WITH_CXX11
1417  static_assert(
1418  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1419  "PreconditionJacobi and VectorType must have the same size_type.");
1420 #endif // DEAL_II_WITH_CXX11
1421 
1422  step (dst, src);
1423 }
1424 
1425 
1426 
1427 //---------------------------------------------------------------------------
1428 
1429 template <typename MatrixType>
1430 template<class VectorType>
1431 inline void
1432 PreconditionSOR<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1433 {
1434 #ifdef DEAL_II_WITH_CXX11
1435  static_assert(
1436  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1437  "PreconditionSOR and VectorType must have the same size_type.");
1438 #endif // DEAL_II_WITH_CXX11
1439 
1440  Assert (this->A!=0, ExcNotInitialized());
1441  this->A->precondition_SOR (dst, src, this->relaxation);
1442 }
1443 
1444 
1445 
1446 template <typename MatrixType>
1447 template<class VectorType>
1448 inline void
1449 PreconditionSOR<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1450 {
1451 #ifdef DEAL_II_WITH_CXX11
1452  static_assert(
1453  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1454  "PreconditionSOR and VectorType must have the same size_type.");
1455 #endif // DEAL_II_WITH_CXX11
1456 
1457  Assert (this->A!=0, ExcNotInitialized());
1458  this->A->precondition_TSOR (dst, src, this->relaxation);
1459 }
1460 
1461 
1462 
1463 template <typename MatrixType>
1464 template<class VectorType>
1465 inline void
1466 PreconditionSOR<MatrixType>::step (VectorType &dst, const VectorType &src) const
1467 {
1468 #ifdef DEAL_II_WITH_CXX11
1469  static_assert(
1470  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1471  "PreconditionSOR and VectorType must have the same size_type.");
1472 #endif // DEAL_II_WITH_CXX11
1473 
1474  Assert (this->A!=0, ExcNotInitialized());
1475  this->A->SOR_step (dst, src, this->relaxation);
1476 }
1477 
1478 
1479 
1480 template <typename MatrixType>
1481 template<class VectorType>
1482 inline void
1483 PreconditionSOR<MatrixType>::Tstep (VectorType &dst, const VectorType &src) const
1484 {
1485 #ifdef DEAL_II_WITH_CXX11
1486  static_assert(
1487  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1488  "PreconditionSOR and VectorType must have the same size_type.");
1489 #endif // DEAL_II_WITH_CXX11
1490 
1491  Assert (this->A!=0, ExcNotInitialized());
1492  this->A->TSOR_step (dst, src, this->relaxation);
1493 }
1494 
1495 
1496 
1497 //---------------------------------------------------------------------------
1498 
1499 template <typename MatrixType>
1500 inline void
1501 PreconditionSSOR<MatrixType>::initialize (const MatrixType &rA,
1502  const typename BaseClass::AdditionalData &parameters)
1503 {
1504  this->PreconditionRelaxation<MatrixType>::initialize (rA, parameters);
1505 
1506  // in case we have a SparseMatrix class, we can extract information about
1507  // the diagonal.
1509  dynamic_cast<const SparseMatrix<typename MatrixType::value_type> *>(&*this->A);
1510 
1511  // calculate the positions first after the diagonal.
1512  if (mat != 0)
1513  {
1514  const size_type n = this->A->n();
1515  pos_right_of_diagonal.resize(n, static_cast<std::size_t>(-1));
1516  for (size_type row=0; row<n; ++row)
1517  {
1518  // find the first element in this line which is on the right of the
1519  // diagonal. we need to precondition with the elements on the left
1520  // only. note: the first entry in each line denotes the diagonal
1521  // element, which we need not check.
1523  it = mat->begin(row)+1;
1524  for ( ; it < mat->end(row); ++it)
1525  if (it->column() > row)
1526  break;
1527  pos_right_of_diagonal[row] = it - mat->begin();
1528  }
1529  }
1530 }
1531 
1532 
1533 template <typename MatrixType>
1534 template<class VectorType>
1535 inline void
1536 PreconditionSSOR<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1537 {
1538 #ifdef DEAL_II_WITH_CXX11
1539  static_assert(
1540  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1541  "PreconditionSSOR and VectorType must have the same size_type.");
1542 #endif // DEAL_II_WITH_CXX11
1543 
1544  Assert (this->A!=0, ExcNotInitialized());
1545  this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1546 }
1547 
1548 
1549 
1550 template <typename MatrixType>
1551 template<class VectorType>
1552 inline void
1553 PreconditionSSOR<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1554 {
1555 #ifdef DEAL_II_WITH_CXX11
1556  static_assert(
1557  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1558  "PreconditionSSOR and VectorType must have the same size_type.");
1559 #endif // DEAL_II_WITH_CXX11
1560 
1561  Assert (this->A!=0, ExcNotInitialized());
1562  this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1563 }
1564 
1565 
1566 
1567 template <typename MatrixType>
1568 template<class VectorType>
1569 inline void
1570 PreconditionSSOR<MatrixType>::step (VectorType &dst, const VectorType &src) const
1571 {
1572 #ifdef DEAL_II_WITH_CXX11
1573  static_assert(
1574  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1575  "PreconditionSSOR and VectorType must have the same size_type.");
1576 #endif // DEAL_II_WITH_CXX11
1577 
1578  Assert (this->A!=0, ExcNotInitialized());
1579  this->A->SSOR_step (dst, src, this->relaxation);
1580 }
1581 
1582 
1583 
1584 template <typename MatrixType>
1585 template<class VectorType>
1586 inline void
1587 PreconditionSSOR<MatrixType>::Tstep (VectorType &dst, const VectorType &src) const
1588 {
1589 #ifdef DEAL_II_WITH_CXX11
1590  static_assert(
1591  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1592  "PreconditionSSOR and VectorType must have the same size_type.");
1593 #endif // DEAL_II_WITH_CXX11
1594 
1595  step (dst, src);
1596 }
1597 
1598 
1599 
1600 //---------------------------------------------------------------------------
1601 
1602 template <typename MatrixType>
1603 inline void
1605 (const MatrixType &rA,
1606  const std::vector<size_type> &p,
1607  const std::vector<size_type> &ip,
1608  const typename PreconditionRelaxation<MatrixType>::AdditionalData &parameters)
1609 {
1610  permutation = &p;
1611  inverse_permutation = &ip;
1613 }
1614 
1615 
1616 template <typename MatrixType>
1617 inline void
1618 PreconditionPSOR<MatrixType>::initialize (const MatrixType &A,
1619  const AdditionalData &additional_data)
1620 {
1621  initialize(A,
1622  additional_data.permutation,
1623  additional_data.inverse_permutation,
1624  additional_data.parameters);
1625 }
1626 
1627 
1628 template <typename MatrixType>
1629 template <typename VectorType>
1630 inline void
1631 PreconditionPSOR<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1632 {
1633 #ifdef DEAL_II_WITH_CXX11
1634  static_assert(
1635  std::is_same<typename PreconditionPSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1636  "PreconditionPSOR and VectorType must have the same size_type.");
1637 #endif // DEAL_II_WITH_CXX11
1638 
1639  Assert (this->A!=0, ExcNotInitialized());
1640  dst = src;
1641  this->A->PSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1642 }
1643 
1644 
1645 
1646 template <typename MatrixType>
1647 template<class VectorType>
1648 inline void
1649 PreconditionPSOR<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1650 {
1651 #ifdef DEAL_II_WITH_CXX11
1652  static_assert(
1653  std::is_same<typename PreconditionPSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1654  "PreconditionPSOR and VectorType must have the same size_type.");
1655 #endif // DEAL_II_WITH_CXX11
1656 
1657  Assert (this->A!=0, ExcNotInitialized());
1658  dst = src;
1659  this->A->TPSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1660 }
1661 
1662 template <typename MatrixType>
1664 (const std::vector<size_type> &permutation,
1665  const std::vector<size_type> &inverse_permutation,
1666  const typename PreconditionRelaxation<MatrixType>::AdditionalData &parameters)
1667  :
1668  permutation(permutation),
1669  inverse_permutation(inverse_permutation),
1670  parameters(parameters)
1671 {
1672 
1673 }
1674 
1675 
1676 //---------------------------------------------------------------------------
1677 
1678 
1679 template<typename MatrixType, class VectorType>
1681  const function_ptr method)
1682  :
1683  matrix(M), precondition(method)
1684 {}
1685 
1686 
1687 
1688 template<typename MatrixType, class VectorType>
1689 void
1691  const VectorType &src) const
1692 {
1693  (matrix.*precondition)(dst, src);
1694 }
1695 
1696 //---------------------------------------------------------------------------
1697 
1698 template<typename MatrixType>
1699 inline
1701 AdditionalData (const double relaxation)
1702  :
1703  relaxation (relaxation)
1704 {}
1705 
1706 
1707 
1708 //---------------------------------------------------------------------------
1709 
1710 namespace internal
1711 {
1712  namespace PreconditionChebyshev
1713  {
1714  // for deal.II vectors, perform updates for Chebyshev preconditioner all
1715  // at once to reduce memory transfer. Here, we select between general
1716  // vectors and deal.II vectors where we expand the loop over the (local)
1717  // size of the vector
1718 
1719  // generic part for non-deal.II vectors
1720  template <typename VectorType, typename PreconditionerType>
1721  inline
1722  void
1723  vector_updates (const VectorType &src,
1724  const PreconditionerType &preconditioner,
1725  const bool start_zero,
1726  const double factor1,
1727  const double factor2,
1728  VectorType &update1,
1729  VectorType &update2,
1730  VectorType &update3,
1731  VectorType &dst)
1732  {
1733  if (start_zero)
1734  {
1735  update1.equ (factor2, src);
1736  preconditioner.vmult(dst, update1);
1737  update1.equ(-1.,dst);
1738  }
1739  else
1740  {
1741  update2 -= src;
1742  preconditioner.vmult(update3, update2);
1743  update2 = update3;
1744  if (factor1 == 0.)
1745  update1.equ(factor2, update2);
1746  else
1747  update1.sadd(factor1, factor2, update2);
1748  dst -= update1;
1749  }
1750  }
1751 
1752  // worker routine for deal.II vectors. Because of vectorization, we need
1753  // to put the loop into an extra structure because the virtual function of
1754  // VectorUpdatesRange prevents the compiler from applying vectorization.
1755  template <typename Number>
1756  struct VectorUpdater
1757  {
1758  VectorUpdater (const Number *src,
1759  const Number *matrix_diagonal_inverse,
1760  const bool start_zero,
1761  const Number factor1,
1762  const Number factor2,
1763  Number *update1,
1764  Number *update2,
1765  Number *dst)
1766  :
1767  src (src),
1768  matrix_diagonal_inverse (matrix_diagonal_inverse),
1769  do_startup (factor1 == Number()),
1770  start_zero (start_zero),
1771  factor1 (factor1),
1772  factor2 (factor2),
1773  update1 (update1),
1774  update2 (update2),
1775  dst (dst)
1776  {}
1777 
1778  void
1779  apply_to_subrange (const std::size_t begin,
1780  const std::size_t end) const
1781  {
1782  // To circumvent a bug in gcc
1783  // (https://gcc.gnu.org/bugzilla/show_bug.cgi?id=63945), we create copies
1784  // of the variables factor1 and factor2 and do not check based on
1785  // factor1.
1786  const Number factor1 = this->factor1;
1787  const Number factor2 = this->factor2;
1788  if (do_startup)
1789  {
1790  if (start_zero)
1791  DEAL_II_OPENMP_SIMD_PRAGMA
1792  for (std::size_t i=begin; i<end; ++i)
1793  {
1794  dst[i] = factor2 * src[i] * matrix_diagonal_inverse[i];
1795  update1[i] = -dst[i];
1796  }
1797  else
1798  DEAL_II_OPENMP_SIMD_PRAGMA
1799  for (std::size_t i=begin; i<end; ++i)
1800  {
1801  update1[i] = ((update2[i]-src[i]) *
1802  factor2*matrix_diagonal_inverse[i]);
1803  dst[i] -= update1[i];
1804  }
1805  }
1806  else
1807  DEAL_II_OPENMP_SIMD_PRAGMA
1808  for (std::size_t i=begin; i<end; ++i)
1809  {
1810  const Number update =
1811  factor1 * update1[i] + factor2 *
1812  ((update2[i] - src[i]) * matrix_diagonal_inverse[i]);
1813  update1[i] = update;
1814  dst[i] -= update;
1815  }
1816  }
1817 
1818  const Number *src;
1819  const Number *matrix_diagonal_inverse;
1820  const bool do_startup;
1821  const bool start_zero;
1822  const Number factor1;
1823  const Number factor2;
1824  mutable Number *update1;
1825  mutable Number *update2;
1826  mutable Number *dst;
1827  };
1828 
1829  template<typename Number>
1830  struct VectorUpdatesRange : public parallel::ParallelForInteger
1831  {
1832  VectorUpdatesRange(const VectorUpdater<Number> &updater,
1833  const std::size_t size)
1834  :
1835  updater (updater)
1836  {
1837  if (size < internal::Vector::minimum_parallel_grain_size)
1838  apply_to_subrange (0, size);
1839  else
1840  apply_parallel (0, size,
1841  internal::Vector::minimum_parallel_grain_size);
1842  }
1843 
1844  ~VectorUpdatesRange() {}
1845 
1846  virtual void
1847  apply_to_subrange (const std::size_t begin,
1848  const std::size_t end) const
1849  {
1850  updater.apply_to_subrange(begin, end);
1851  }
1852 
1853  const VectorUpdater<Number> &updater;
1854  };
1855 
1856  // selection for diagonal matrix around deal.II vector
1857  template <typename Number>
1858  inline
1859  void
1860  vector_updates (const ::Vector<Number> &src,
1861  const DiagonalMatrix< ::Vector<Number> > &jacobi,
1862  const bool start_zero,
1863  const double factor1,
1864  const double factor2,
1865  ::Vector<Number> &update1,
1866  ::Vector<Number> &update2,
1867  ::Vector<Number> &,
1868  ::Vector<Number> &dst)
1869  {
1870  VectorUpdater<Number> upd(src.begin(), jacobi.get_vector().begin(),
1871  start_zero, factor1, factor2,
1872  update1.begin(), update2.begin(), dst.begin());
1873  VectorUpdatesRange<Number>(upd, src.size());
1874  }
1875 
1876  // selection for diagonal matrix around parallel deal.II vector
1877  template <typename Number>
1878  inline
1879  void
1880  vector_updates (const LinearAlgebra::distributed::Vector<Number> &src,
1882  const bool start_zero,
1883  const double factor1,
1884  const double factor2,
1889  {
1890  VectorUpdater<Number> upd(src.begin(), jacobi.get_vector().begin(),
1891  start_zero, factor1, factor2,
1892  update1.begin(), update2.begin(), dst.begin());
1893  VectorUpdatesRange<Number>(upd, src.local_size());
1894  }
1895 
1896  template <typename MatrixType, typename VectorType, typename PreconditionerType>
1897  inline
1898  void
1899  initialize_preconditioner(const MatrixType &matrix,
1900  std_cxx11::shared_ptr<PreconditionerType> &preconditioner,
1901  VectorType &)
1902  {
1903  (void)matrix;
1904  (void)preconditioner;
1905  AssertThrow(preconditioner.get() != NULL, ExcNotInitialized());
1906  }
1907 
1908  template <typename MatrixType, typename VectorType>
1909  inline
1910  void
1911  initialize_preconditioner(const MatrixType &matrix,
1912  std_cxx11::shared_ptr<DiagonalMatrix<VectorType> > &preconditioner,
1913  VectorType &diagonal_inverse)
1914  {
1915  if (preconditioner.get() == NULL ||
1916  preconditioner->m() != matrix.m())
1917  {
1918  if (preconditioner.get() == NULL)
1919  preconditioner.reset(new DiagonalMatrix<VectorType>());
1920 
1921  Assert(preconditioner->m() == 0,
1922  ExcMessage("Preconditioner appears to be initialized but not sized correctly"));
1923 
1924  // Check if we can initialize from vector that then gets set to zero
1925  // as the matrix will own the memory
1926  preconditioner->reinit(diagonal_inverse);
1927  {
1928  VectorType empty_vector;
1929  diagonal_inverse.reinit(empty_vector);
1930  }
1931 
1932  // This part only works in serial
1933  if (preconditioner->m() != matrix.m())
1934  {
1935  preconditioner->get_vector().reinit(matrix.m());
1936  for (typename VectorType::size_type i=0; i<matrix.m(); ++i)
1937  preconditioner->get_vector()(i) = 1./matrix.el(i,i);
1938  }
1939  }
1940  }
1941 
1942  template <typename VectorType>
1943  void set_initial_guess(VectorType &vector)
1944  {
1945  vector = 1./std::sqrt(static_cast<double>(vector.size()));
1946  if (vector.locally_owned_elements().is_element(0))
1947  vector(0) = 0.;
1948  }
1949 
1950  template <typename Number>
1951  void set_initial_guess(::Vector<Number> &vector)
1952  {
1953  // Choose a high-frequency mode consisting of numbers between 0 and 1
1954  // that is cheap to compute (cheaper than random numbers) but avoids
1955  // obviously re-occurring numbers in multi-component systems by choosing
1956  // a period of 11
1957  for (unsigned int i=0; i<vector.size(); ++i)
1958  vector(i) = i%11;
1959  }
1960 
1961  template <typename Number>
1962  void set_initial_guess(::LinearAlgebra::distributed::Vector<Number> &vector)
1963  {
1964  // Choose a high-frequency mode consisting of numbers between 0 and 1
1965  // that is cheap to compute (cheaper than random numbers) but avoids
1966  // obviously re-occurring numbers in multi-component systems by choosing
1967  // a period of 11
1968  for (unsigned int i=0; i<vector.local_size(); ++i)
1969  vector.local_element(i) = i%11;
1970  }
1971 
1972  struct EigenvalueTracker
1973  {
1974  public:
1975  void slot(const std::vector<double> &eigenvalues)
1976  {
1977  values = eigenvalues;
1978  }
1979 
1980  std::vector<double> values;
1981  };
1982  }
1983 }
1984 
1985 
1986 
1987 // avoid warning about deprecated variable nonzero_starting
1988 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
1989 
1990 template <typename MatrixType, class VectorType, typename PreconditionerType>
1991 inline
1993 AdditionalData (const unsigned int degree,
1994  const double smoothing_range,
1995  const bool nonzero_starting,
1996  const unsigned int eig_cg_n_iterations,
1997  const double eig_cg_residual,
1998  const double max_eigenvalue)
1999  :
2000  degree (degree),
2001  smoothing_range (smoothing_range),
2002  nonzero_starting (nonzero_starting),
2003  eig_cg_n_iterations (eig_cg_n_iterations),
2004  eig_cg_residual (eig_cg_residual),
2005  max_eigenvalue (max_eigenvalue)
2006 {}
2007 
2008 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
2009 
2010 
2011 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2012 inline
2014  :
2015  theta (1.),
2016  delta (1.),
2017  is_initialized (false)
2018 {
2019 #ifdef DEAL_II_WITH_CXX11
2020  static_assert(
2021  std::is_same<size_type, typename VectorType::size_type>::value,
2022  "PreconditionChebyshev and VectorType must have the same size_type.");
2023 #endif // DEAL_II_WITH_CXX11
2024 }
2025 
2026 
2027 // avoid warning about deprecated variable AdditionalData::matrix_diagonal_inverse
2028 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
2029 
2030 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2031 inline
2032 void
2034 (const MatrixType &matrix,
2035  const AdditionalData &additional_data)
2036 {
2037  matrix_ptr = &matrix;
2038  data = additional_data;
2039  internal::PreconditionChebyshev::initialize_preconditioner(matrix,
2040  data.preconditioner,
2041  data.matrix_diagonal_inverse);
2042  is_initialized = false;
2043 }
2044 
2045 
2046 
2047 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2048 inline
2049 void
2051 {
2052  is_initialized = false;
2053  theta = delta = 1.0;
2054  matrix_ptr = 0;
2055  {
2056  VectorType empty_vector;
2057  data.matrix_diagonal_inverse.reinit(empty_vector);
2058  update1.reinit(empty_vector);
2059  update2.reinit(empty_vector);
2060  update3.reinit(empty_vector);
2061  }
2062  data.preconditioner.reset();
2063 }
2064 
2065 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
2066 
2067 
2068 
2069 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2070 inline
2071 void
2073 (const VectorType &src) const
2074 {
2075  Assert(is_initialized == false, ExcInternalError());
2076  Assert(data.preconditioner.get() != NULL, ExcNotInitialized());
2077 
2078  update1.reinit(src);
2079  update2.reinit(src, true);
2080 
2081  // calculate largest eigenvalue using a hand-tuned CG iteration on the
2082  // matrix weighted by its diagonal. we start with a vector that consists of
2083  // ones only, weighted by the length.
2084  double max_eigenvalue, min_eigenvalue;
2085  if (data.eig_cg_n_iterations > 0)
2086  {
2087  Assert (data.eig_cg_n_iterations > 2,
2088  ExcMessage ("Need to set at least two iterations to find eigenvalues."));
2089 
2090  // set a very strict tolerance to force at least two iterations
2091  ReductionControl control (data.eig_cg_n_iterations,
2092  std::sqrt(std::numeric_limits<typename VectorType::value_type>::epsilon()),
2093  1e-10, false, false);
2094 
2095  internal::PreconditionChebyshev::EigenvalueTracker eigenvalue_tracker;
2096  SolverCG<VectorType> solver (control);
2097  solver.connect_eigenvalues_slot(std_cxx11::bind(&internal::PreconditionChebyshev::EigenvalueTracker::slot,
2098  &eigenvalue_tracker,
2099  std_cxx11::_1));
2100 
2101  // set an initial guess which is close to the constant vector but where
2102  // one entry is different to trigger high frequencies
2103  internal::PreconditionChebyshev::set_initial_guess(update2);
2104 
2105  try
2106  {
2107  solver.solve(*matrix_ptr, update1, update2, *data.preconditioner);
2108  }
2110  {
2111  }
2112 
2113  // read the eigenvalues from the attached eigenvalue tracker
2114  if (eigenvalue_tracker.values.empty())
2115  min_eigenvalue = max_eigenvalue = 1;
2116  else
2117  {
2118  min_eigenvalue = eigenvalue_tracker.values.front();
2119 
2120  // include a safety factor since the CG method will in general not
2121  // be converged
2122  max_eigenvalue = 1.2*eigenvalue_tracker.values.back();
2123  }
2124  }
2125  else
2126  {
2127  max_eigenvalue = data.max_eigenvalue;
2128  min_eigenvalue = data.max_eigenvalue/data.smoothing_range;
2129  }
2130 
2131  const double alpha = (data.smoothing_range > 1. ?
2132  max_eigenvalue / data.smoothing_range :
2133  std::min(0.9*max_eigenvalue, min_eigenvalue));
2134 
2135  // in case the user set the degree to invalid unsigned int, we have to
2136  // determine the number of necessary iterations from the Chebyshev error
2137  // estimate, given the target tolerance specified by smoothing_range. This
2138  // estimate is based on the error formula given in section 5.1 of
2139  // R. S. Varga, Matrix iterative analysis, 2nd ed., Springer, 2009
2140  if (data.degree == numbers::invalid_unsigned_int)
2141  {
2142  const double actual_range = max_eigenvalue / alpha;
2143  const double sigma = (1.-std::sqrt(1./actual_range))/(1.+std::sqrt(1./actual_range));
2144  const double eps = data.smoothing_range;
2146  = 1+std::log(1./eps+std::sqrt(1./eps/eps-1))/std::log(1./sigma);
2147  }
2148 
2149  const_cast<PreconditionChebyshev<MatrixType,VectorType,PreconditionerType> *>(this)->delta = (max_eigenvalue-alpha)*0.5;
2150  const_cast<PreconditionChebyshev<MatrixType,VectorType,PreconditionerType> *>(this)->theta = (max_eigenvalue+alpha)*0.5;
2151 
2152  // We do not need the third auxiliary vector in case we have a
2153  // DiagonalMatrix as preconditioner and use deal.II's own vectors
2154  if (types_are_equal<PreconditionerType,DiagonalMatrix<VectorType> >::value == false ||
2155  (types_are_equal<VectorType,::Vector<typename VectorType::value_type> >::value == false
2156  &&
2158  ))
2159  update3.reinit (src, true);
2160 
2161  const_cast<PreconditionChebyshev<MatrixType,VectorType,PreconditionerType> *>(this)->is_initialized = true;
2162 }
2163 
2164 
2165 
2166 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2167 inline
2168 void
2170 ::do_chebyshev_loop(VectorType &dst,
2171  const VectorType &src) const
2172 {
2173  // if delta is zero, we do not need to iterate because the updates will be
2174  // zero
2175  if (std::abs(delta) < 1e-40)
2176  return;
2177 
2178  double rhok = delta / theta, sigma = theta / delta;
2179  for (unsigned int k=0; k<data.degree; ++k)
2180  {
2181  matrix_ptr->vmult (update2, dst);
2182  const double rhokp = 1./(2.*sigma-rhok);
2183  const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
2184  rhok = rhokp;
2185  internal::PreconditionChebyshev::vector_updates
2186  (src, *data.preconditioner, false, factor1, factor2, update1, update2, update3, dst);
2187  }
2188 }
2189 
2190 
2191 
2192 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2193 inline
2194 void
2196 ::do_transpose_chebyshev_loop(VectorType &dst,
2197  const VectorType &src) const
2198 {
2199  double rhok = delta / theta, sigma = theta / delta;
2200  for (unsigned int k=0; k<data.degree; ++k)
2201  {
2202  matrix_ptr->Tvmult (update2, dst);
2203  const double rhokp = 1./(2.*sigma-rhok);
2204  const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
2205  rhok = rhokp;
2206  internal::PreconditionChebyshev::vector_updates
2207  (src, *data.preconditioner, false, factor1, factor2, update1, update2, update3, dst);
2208  }
2209 }
2210 
2211 
2212 
2213 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2214 inline
2215 void
2217 ::vmult (VectorType &dst,
2218  const VectorType &src) const
2219 {
2220  Threads::Mutex::ScopedLock lock(mutex);
2221  if (is_initialized == false)
2222  estimate_eigenvalues(src);
2223 
2224  internal::PreconditionChebyshev::vector_updates
2225  (src, *data.preconditioner, true, 0., 1./theta, update1, update2, update3, dst);
2226 
2227  do_chebyshev_loop(dst, src);
2228 }
2229 
2230 
2231 
2232 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2233 inline
2234 void
2236 ::Tvmult (VectorType &dst,
2237  const VectorType &src) const
2238 {
2239  Threads::Mutex::ScopedLock lock(mutex);
2240  if (is_initialized == false)
2241  estimate_eigenvalues(src);
2242 
2243  internal::PreconditionChebyshev::vector_updates
2244  (src, *data.preconditioner, true, 0., 1./theta, update1, update2, update3, dst);
2245 
2246  do_transpose_chebyshev_loop(dst, src);
2247 }
2248 
2249 
2250 
2251 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2252 inline
2253 void
2255 ::step (VectorType &dst,
2256  const VectorType &src) const
2257 {
2258  Threads::Mutex::ScopedLock lock(mutex);
2259  if (is_initialized == false)
2260  estimate_eigenvalues(src);
2261 
2262  if (!dst.all_zero())
2263  {
2264  matrix_ptr->vmult (update2, dst);
2265  internal::PreconditionChebyshev::vector_updates
2266  (src, *data.preconditioner, false, 0., 1./theta, update1, update2, update3, dst);
2267  }
2268  else
2269  internal::PreconditionChebyshev::vector_updates
2270  (src, *data.preconditioner, true, 0., 1./theta, update1, update2, update3, dst);
2271 
2272  do_chebyshev_loop(dst, src);
2273 }
2274 
2275 
2276 
2277 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2278 inline
2279 void
2281 ::Tstep (VectorType &dst,
2282  const VectorType &src) const
2283 {
2284  Threads::Mutex::ScopedLock lock(mutex);
2285  if (is_initialized == false)
2286  estimate_eigenvalues(src);
2287 
2288  if (!dst.all_zero())
2289  {
2290  matrix_ptr->Tvmult (update2, dst);
2291  internal::PreconditionChebyshev::vector_updates
2292  (src, *data.preconditioner, false, 0., 1./theta, update1, update2, update3, dst);
2293  }
2294  else
2295  internal::PreconditionChebyshev::vector_updates
2296  (src, *data.preconditioner, true, 0., 1./theta, update1, update2, update3, dst);
2297 
2298  do_transpose_chebyshev_loop(dst, src);
2299 }
2300 
2301 
2302 
2303 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2304 inline
2307 {
2308  Assert (matrix_ptr!=0, ExcNotInitialized());
2309  return matrix_ptr->m();
2310 }
2311 
2312 
2313 
2314 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2315 inline
2318 {
2319  Assert (matrix_ptr!=0, ExcNotInitialized());
2320  return matrix_ptr->n();
2321 }
2322 
2323 #endif // DOXYGEN
2324 
2325 DEAL_II_NAMESPACE_CLOSE
2326 
2327 #endif
void reinit(const size_type size, const bool omit_zeroing_entries=false)
void Tvmult(VectorType &, const VectorType &) const
static const unsigned int invalid_unsigned_int
Definition: types.h:170
size_type m() const
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
void vmult(VectorType &dst, const VectorType &src) const
MatrixType::size_type size_type
Definition: precondition.h:632
void Tvmult_add(VectorType &, const VectorType &) const
PreconditionRelaxation< MatrixType >::AdditionalData parameters
Definition: precondition.h:756
void do_chebyshev_loop(VectorType &dst, const VectorType &src) const
void vmult(VectorType &, const VectorType &) const
size_type m() const
PreconditionRelaxation< MatrixType >::AdditionalData AdditionalData
Definition: precondition.h:627
const MatrixType & matrix
Definition: precondition.h:365
AdditionalData data
AdditionalData(const double relaxation=1.)
static ::ExceptionBase & ExcNotInitialized()
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
SmartPointer< const MatrixType, PreconditionChebyshev< MatrixType, VectorType, PreconditionerType > > matrix_ptr
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
AdditionalData(const double relaxation=1.)
size_type m() const
types::global_dof_index size_type
Definition: precondition.h:189
const std::vector< size_type > * permutation
Definition: precondition.h:804
const function_ptr precondition
Definition: precondition.h:370
void vmult(VectorType &, const VectorType &) const
void initialize(const AdditionalData &parameters)
void Tvmult(VectorType &, const VectorType &) const
static ::ExceptionBase & ExcMessage(std::string arg1)
void Tvmult(VectorType &dst, const VectorType &src) const
types::global_dof_index size_type
Definition: precondition.h:79
size_type n() const
void vmult(VectorType &, const VectorType &) const
unsigned int global_dof_index
Definition: types.h:88
std_cxx11::shared_ptr< PreconditionerType > preconditioner
Definition: precondition.h:999
const std::vector< size_type > & inverse_permutation
Definition: precondition.h:752
#define Assert(cond, exc)
Definition: exceptions.h:313
void step(VectorType &x, const VectorType &rhs) const
std::size_t size() const
void(MatrixType::* function_ptr)(VectorType &, const VectorType &) const
Definition: precondition.h:344
MatrixType::size_type size_type
Definition: precondition.h:722
types::global_dof_index size_type
Definition: precondition.h:908
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())
Threads::Mutex mutex
void Tvmult(VectorType &, const VectorType &) const
iterator begin()
void estimate_eigenvalues(const VectorType &src) 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)
MatrixType::size_type size_type
Definition: precondition.h:390
PreconditionRelaxation< MatrixType > BaseClass
Definition: precondition.h:637
void Tstep(VectorType &x, const VectorType &rhs) const
void vmult(VectorType &, const VectorType &) const
void Tstep(VectorType &x, const VectorType &rhs) const
void vmult(VectorType &dst, const VectorType &src) 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:438
void Tvmult(VectorType &, const VectorType &) const
void step(VectorType &x, const VectorType &rhs) const
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
size_type n() const
void do_transpose_chebyshev_loop(VectorType &dst, const VectorType &src) const
const std::vector< size_type > * inverse_permutation
Definition: precondition.h:808
PreconditionRelaxation< MatrixType >::AdditionalData AdditionalData
Definition: precondition.h:482
void Tstep(VectorType &dst, const VectorType &src) const
void step(VectorType &dst, const VectorType &src) const
void step(VectorType &x, const VectorType &rhs) const
Number local_element(const size_type local_index) const
size_type m() const
const_iterator begin() const
void vmult(VectorType &, const VectorType &) const
const std::vector< size_type > & permutation
Definition: precondition.h:748
PreconditionUseMatrix(const MatrixType &M, const function_ptr method)
void vmult_add(VectorType &, const VectorType &) const
PreconditionRelaxation< MatrixType >::AdditionalData AdditionalData
Definition: precondition.h:564
std::vector< std::size_t > pos_right_of_diagonal
Definition: precondition.h:679
static ::ExceptionBase & ExcInternalError()
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())