Reference documentation for deal.II version 9.0.0
precondition.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2018 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  {
90  AdditionalData () = default;
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 
917  {
921  AdditionalData (const unsigned int degree = 0,
922  const double smoothing_range = 0.,
923  const bool nonzero_starting = false,
924  const unsigned int eig_cg_n_iterations = 8,
925  const double eig_cg_residual = 1e-2,
926  const double max_eigenvalue = 1);
927 
939  unsigned int degree;
940 
953 
966  bool nonzero_starting DEAL_II_DEPRECATED;
967 
973  unsigned int eig_cg_n_iterations;
974 
980 
987 
993  VectorType matrix_diagonal_inverse DEAL_II_DEPRECATED;
994 
998  std::shared_ptr<PreconditionerType> preconditioner;
999  };
1000 
1001 
1003 
1015  void initialize (const MatrixType &matrix,
1016  const AdditionalData &additional_data = AdditionalData());
1017 
1022  void vmult (VectorType &dst,
1023  const VectorType &src) const;
1024 
1029  void Tvmult (VectorType &dst,
1030  const VectorType &src) const;
1031 
1035  void step (VectorType &dst,
1036  const VectorType &src) const;
1037 
1041  void Tstep (VectorType &dst,
1042  const VectorType &src) const;
1043 
1047  void clear ();
1048 
1053  size_type m () const;
1054 
1059  size_type n () const;
1060 
1061 private:
1062 
1067 
1071  mutable VectorType update1;
1072 
1076  mutable VectorType update2;
1077 
1081  mutable VectorType update3;
1082 
1088 
1092  double theta;
1093 
1098  double delta;
1099 
1105 
1111 
1116  void do_chebyshev_loop(VectorType &dst,
1117  const VectorType &src) const;
1118 
1125  void do_transpose_chebyshev_loop(VectorType &dst,
1126  const VectorType &src) const;
1127 
1134  void estimate_eigenvalues(const VectorType &src) const;
1135 };
1136 
1137 
1138 
1140 /* ---------------------------------- Inline functions ------------------- */
1141 
1142 #ifndef DOXYGEN
1143 
1144 inline
1146  :
1147  n_rows (0),
1148  n_columns (0)
1149 {}
1150 
1151 template <typename MatrixType>
1152 inline void
1153 PreconditionIdentity::initialize (const MatrixType &matrix,
1155 {
1156  n_rows = matrix.m();
1157  n_columns = matrix.n();
1158 }
1159 
1160 
1161 template <class VectorType>
1162 inline void
1163 PreconditionIdentity::vmult (VectorType &dst, const VectorType &src) const
1164 {
1165  dst = src;
1166 }
1167 
1168 
1169 
1170 template <class VectorType>
1171 inline void
1172 PreconditionIdentity::Tvmult (VectorType &dst, const VectorType &src) const
1173 {
1174  dst = src;
1175 }
1176 
1177 template <class VectorType>
1178 inline void
1179 PreconditionIdentity::vmult_add (VectorType &dst, const VectorType &src) const
1180 {
1181  dst += src;
1182 }
1183 
1184 
1185 
1186 template <class VectorType>
1187 inline void
1188 PreconditionIdentity::Tvmult_add (VectorType &dst, const VectorType &src) const
1189 {
1190  dst += src;
1191 }
1192 
1194 PreconditionIdentity::m () const
1195 {
1196  Assert(n_rows != 0, ExcNotInitialized());
1197  return n_rows;
1198 }
1199 
1201 PreconditionIdentity::n () const
1202 {
1204  return n_columns;
1205 }
1206 
1207 //---------------------------------------------------------------------------
1208 
1209 inline
1211  :
1212  relaxation(relaxation)
1213 {}
1214 
1215 
1216 inline
1218  :
1219  relaxation(0),
1220  n_rows (0),
1221  n_columns (0)
1222 {
1223  AdditionalData add_data;
1224  relaxation=add_data.relaxation;
1225 }
1226 
1227 
1228 
1229 inline void
1231 (const PreconditionRichardson::AdditionalData &parameters)
1232 {
1233  relaxation = parameters.relaxation;
1234 }
1235 
1236 
1237 
1238 template <typename MatrixType>
1239 inline void
1241 (const MatrixType &matrix,
1242  const PreconditionRichardson::AdditionalData &parameters)
1243 {
1244  relaxation = parameters.relaxation;
1245  n_rows = matrix.m();
1246  n_columns = matrix.n();
1247 }
1248 
1249 
1250 
1251 template <class VectorType>
1252 inline void
1253 PreconditionRichardson::vmult (VectorType &dst, const VectorType &src) const
1254 {
1255  static_assert(
1256  std::is_same<size_type, typename VectorType::size_type>::value,
1257  "PreconditionRichardson and VectorType must have the same size_type.");
1258 
1259  dst.equ(relaxation,src);
1260 }
1261 
1262 
1263 
1264 template <class VectorType>
1265 inline void
1266 PreconditionRichardson::Tvmult (VectorType &dst, const VectorType &src) const
1267 {
1268  static_assert(
1269  std::is_same<size_type, typename VectorType::size_type>::value,
1270  "PreconditionRichardson and VectorType must have the same size_type.");
1271 
1272  dst.equ(relaxation,src);
1273 }
1274 
1275 template <class VectorType>
1276 inline void
1277 PreconditionRichardson::vmult_add (VectorType &dst, const VectorType &src) const
1278 {
1279  static_assert(
1280  std::is_same<size_type, typename VectorType::size_type>::value,
1281  "PreconditionRichardson and VectorType must have the same size_type.");
1282 
1283  dst.add(relaxation,src);
1284 }
1285 
1286 
1287 
1288 template <class VectorType>
1289 inline void
1290 PreconditionRichardson::Tvmult_add (VectorType &dst, const VectorType &src) const
1291 {
1292  static_assert(
1293  std::is_same<size_type, typename VectorType::size_type>::value,
1294  "PreconditionRichardson and VectorType must have the same size_type.");
1295 
1296  dst.add(relaxation,src);
1297 }
1298 
1301 {
1302  Assert(n_rows != 0, ExcNotInitialized());
1303  return n_rows;
1304 }
1305 
1308 {
1310  return n_columns;
1311 }
1312 
1313 //---------------------------------------------------------------------------
1314 
1315 template <typename MatrixType>
1316 inline void
1318  const AdditionalData &parameters)
1319 {
1320  A = &rA;
1321  relaxation = parameters.relaxation;
1322 }
1323 
1324 
1325 template <typename MatrixType>
1326 inline void
1328 {
1329  A = nullptr;
1330 }
1331 
1332 template <typename MatrixType>
1335 {
1336  Assert (A!=nullptr, ExcNotInitialized());
1337  return A->m();
1338 }
1339 
1340 template <typename MatrixType>
1343 {
1344  Assert (A!=nullptr, ExcNotInitialized());
1345  return A->n();
1346 }
1347 
1348 //---------------------------------------------------------------------------
1349 
1350 template <typename MatrixType>
1351 template <class VectorType>
1352 inline void
1353 PreconditionJacobi<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1354 {
1355  static_assert(
1356  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1357  "PreconditionJacobi and VectorType must have the same size_type.");
1358 
1359  Assert (this->A!=nullptr, ExcNotInitialized());
1360  this->A->precondition_Jacobi (dst, src, this->relaxation);
1361 }
1362 
1363 
1364 
1365 template <typename MatrixType>
1366 template <class VectorType>
1367 inline void
1368 PreconditionJacobi<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1369 {
1370  static_assert(
1371  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1372  "PreconditionJacobi and VectorType must have the same size_type.");
1373 
1374  Assert (this->A!=nullptr, ExcNotInitialized());
1375  this->A->precondition_Jacobi (dst, src, this->relaxation);
1376 }
1377 
1378 
1379 
1380 template <typename MatrixType>
1381 template <class VectorType>
1382 inline void
1383 PreconditionJacobi<MatrixType>::step (VectorType &dst, const VectorType &src) const
1384 {
1385  static_assert(
1386  std::is_same<typename PreconditionJacobi<MatrixType>::size_type, typename VectorType::size_type>::value,
1387  "PreconditionJacobi and VectorType must have the same size_type.");
1388 
1389  Assert (this->A!=nullptr, ExcNotInitialized());
1390  this->A->Jacobi_step (dst, src, this->relaxation);
1391 }
1392 
1393 
1394 
1395 template <typename MatrixType>
1396 template <class VectorType>
1397 inline void
1398 PreconditionJacobi<MatrixType>::Tstep (VectorType &dst, const VectorType &src) const
1399 {
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 
1404  step (dst, src);
1405 }
1406 
1407 
1408 
1409 //---------------------------------------------------------------------------
1410 
1411 template <typename MatrixType>
1412 template <class VectorType>
1413 inline void
1414 PreconditionSOR<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1415 {
1416  static_assert(
1417  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1418  "PreconditionSOR and VectorType must have the same size_type.");
1419 
1420  Assert (this->A!=nullptr, ExcNotInitialized());
1421  this->A->precondition_SOR (dst, src, this->relaxation);
1422 }
1423 
1424 
1425 
1426 template <typename MatrixType>
1427 template <class VectorType>
1428 inline void
1429 PreconditionSOR<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1430 {
1431  static_assert(
1432  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1433  "PreconditionSOR and VectorType must have the same size_type.");
1434 
1435  Assert (this->A!=nullptr, ExcNotInitialized());
1436  this->A->precondition_TSOR (dst, src, this->relaxation);
1437 }
1438 
1439 
1440 
1441 template <typename MatrixType>
1442 template <class VectorType>
1443 inline void
1444 PreconditionSOR<MatrixType>::step (VectorType &dst, const VectorType &src) const
1445 {
1446  static_assert(
1447  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1448  "PreconditionSOR and VectorType must have the same size_type.");
1449 
1450  Assert (this->A!=nullptr, ExcNotInitialized());
1451  this->A->SOR_step (dst, src, this->relaxation);
1452 }
1453 
1454 
1455 
1456 template <typename MatrixType>
1457 template <class VectorType>
1458 inline void
1459 PreconditionSOR<MatrixType>::Tstep (VectorType &dst, const VectorType &src) const
1460 {
1461  static_assert(
1462  std::is_same<typename PreconditionSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1463  "PreconditionSOR and VectorType must have the same size_type.");
1464 
1465  Assert (this->A!=nullptr, ExcNotInitialized());
1466  this->A->TSOR_step (dst, src, this->relaxation);
1467 }
1468 
1469 
1470 
1471 //---------------------------------------------------------------------------
1472 
1473 template <typename MatrixType>
1474 inline void
1475 PreconditionSSOR<MatrixType>::initialize (const MatrixType &rA,
1476  const typename BaseClass::AdditionalData &parameters)
1477 {
1478  this->PreconditionRelaxation<MatrixType>::initialize (rA, parameters);
1479 
1480  // in case we have a SparseMatrix class, we can extract information about
1481  // the diagonal.
1483  dynamic_cast<const SparseMatrix<typename MatrixType::value_type> *>(&*this->A);
1484 
1485  // calculate the positions first after the diagonal.
1486  if (mat != nullptr)
1487  {
1488  const size_type n = this->A->n();
1489  pos_right_of_diagonal.resize(n, static_cast<std::size_t>(-1));
1490  for (size_type row=0; row<n; ++row)
1491  {
1492  // find the first element in this line which is on the right of the
1493  // diagonal. we need to precondition with the elements on the left
1494  // only. note: the first entry in each line denotes the diagonal
1495  // element, which we need not check.
1497  it = mat->begin(row)+1;
1498  for ( ; it < mat->end(row); ++it)
1499  if (it->column() > row)
1500  break;
1501  pos_right_of_diagonal[row] = it - mat->begin();
1502  }
1503  }
1504 }
1505 
1506 
1507 template <typename MatrixType>
1508 template <class VectorType>
1509 inline void
1510 PreconditionSSOR<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1511 {
1512  static_assert(
1513  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1514  "PreconditionSSOR and VectorType must have the same size_type.");
1515 
1516  Assert (this->A!=nullptr, ExcNotInitialized());
1517  this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1518 }
1519 
1520 
1521 
1522 template <typename MatrixType>
1523 template <class VectorType>
1524 inline void
1525 PreconditionSSOR<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1526 {
1527  static_assert(
1528  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1529  "PreconditionSSOR and VectorType must have the same size_type.");
1530 
1531  Assert (this->A!=nullptr, ExcNotInitialized());
1532  this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1533 }
1534 
1535 
1536 
1537 template <typename MatrixType>
1538 template <class VectorType>
1539 inline void
1540 PreconditionSSOR<MatrixType>::step (VectorType &dst, const VectorType &src) const
1541 {
1542  static_assert(
1543  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1544  "PreconditionSSOR and VectorType must have the same size_type.");
1545 
1546  Assert (this->A!=nullptr, ExcNotInitialized());
1547  this->A->SSOR_step (dst, src, this->relaxation);
1548 }
1549 
1550 
1551 
1552 template <typename MatrixType>
1553 template <class VectorType>
1554 inline void
1555 PreconditionSSOR<MatrixType>::Tstep (VectorType &dst, const VectorType &src) const
1556 {
1557  static_assert(
1558  std::is_same<typename PreconditionSSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1559  "PreconditionSSOR and VectorType must have the same size_type.");
1560 
1561  step (dst, src);
1562 }
1563 
1564 
1565 
1566 //---------------------------------------------------------------------------
1567 
1568 template <typename MatrixType>
1569 inline void
1571 (const MatrixType &rA,
1572  const std::vector<size_type> &p,
1573  const std::vector<size_type> &ip,
1574  const typename PreconditionRelaxation<MatrixType>::AdditionalData &parameters)
1575 {
1576  permutation = &p;
1577  inverse_permutation = &ip;
1579 }
1580 
1581 
1582 template <typename MatrixType>
1583 inline void
1584 PreconditionPSOR<MatrixType>::initialize (const MatrixType &A,
1585  const AdditionalData &additional_data)
1586 {
1587  initialize(A,
1588  additional_data.permutation,
1589  additional_data.inverse_permutation,
1590  additional_data.parameters);
1591 }
1592 
1593 
1594 template <typename MatrixType>
1595 template <typename VectorType>
1596 inline void
1597 PreconditionPSOR<MatrixType>::vmult (VectorType &dst, const VectorType &src) const
1598 {
1599  static_assert(
1600  std::is_same<typename PreconditionPSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1601  "PreconditionPSOR and VectorType must have the same size_type.");
1602 
1603  Assert (this->A!=nullptr, ExcNotInitialized());
1604  dst = src;
1605  this->A->PSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1606 }
1607 
1608 
1609 
1610 template <typename MatrixType>
1611 template <class VectorType>
1612 inline void
1613 PreconditionPSOR<MatrixType>::Tvmult (VectorType &dst, const VectorType &src) const
1614 {
1615  static_assert(
1616  std::is_same<typename PreconditionPSOR<MatrixType>::size_type, typename VectorType::size_type>::value,
1617  "PreconditionPSOR and VectorType must have the same size_type.");
1618 
1619  Assert (this->A!=nullptr, ExcNotInitialized());
1620  dst = src;
1621  this->A->TPSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1622 }
1623 
1624 template <typename MatrixType>
1626 (const std::vector<size_type> &permutation,
1627  const std::vector<size_type> &inverse_permutation,
1628  const typename PreconditionRelaxation<MatrixType>::AdditionalData &parameters)
1629  :
1630  permutation(permutation),
1631  inverse_permutation(inverse_permutation),
1632  parameters(parameters)
1633 {
1634 
1635 }
1636 
1637 
1638 //---------------------------------------------------------------------------
1639 
1640 
1641 template <typename MatrixType, class VectorType>
1643  const function_ptr method)
1644  :
1645  matrix(M), precondition(method)
1646 {}
1647 
1648 
1649 
1650 template <typename MatrixType, class VectorType>
1651 void
1653  const VectorType &src) const
1654 {
1655  (matrix.*precondition)(dst, src);
1656 }
1657 
1658 //---------------------------------------------------------------------------
1659 
1660 template <typename MatrixType>
1661 inline
1663 AdditionalData (const double relaxation)
1664  :
1665  relaxation (relaxation)
1666 {}
1667 
1668 
1669 
1670 //---------------------------------------------------------------------------
1671 
1672 namespace internal
1673 {
1674  namespace PreconditionChebyshevImplementation
1675  {
1676  // for deal.II vectors, perform updates for Chebyshev preconditioner all
1677  // at once to reduce memory transfer. Here, we select between general
1678  // vectors and deal.II vectors where we expand the loop over the (local)
1679  // size of the vector
1680 
1681  // generic part for non-deal.II vectors
1682  template <typename VectorType, typename PreconditionerType>
1683  inline
1684  void
1685  vector_updates (const VectorType &src,
1686  const PreconditionerType &preconditioner,
1687  const bool start_zero,
1688  const double factor1,
1689  const double factor2,
1690  VectorType &update1,
1691  VectorType &update2,
1692  VectorType &update3,
1693  VectorType &dst)
1694  {
1695  if (start_zero)
1696  {
1697  update1.equ (factor2, src);
1698  preconditioner.vmult(dst, update1);
1699  update1.equ(-1.,dst);
1700  }
1701  else
1702  {
1703  update2 -= src;
1704  preconditioner.vmult(update3, update2);
1705  update2 = update3;
1706  if (factor1 == 0.)
1707  update1.equ(factor2, update2);
1708  else
1709  update1.sadd(factor1, factor2, update2);
1710  dst -= update1;
1711  }
1712  }
1713 
1714  // worker routine for deal.II vectors. Because of vectorization, we need
1715  // to put the loop into an extra structure because the virtual function of
1716  // VectorUpdatesRange prevents the compiler from applying vectorization.
1717  template <typename Number>
1718  struct VectorUpdater
1719  {
1720  VectorUpdater (const Number *src,
1721  const Number *matrix_diagonal_inverse,
1722  const bool start_zero,
1723  const Number factor1,
1724  const Number factor2,
1725  Number *update1,
1726  Number *update2,
1727  Number *dst)
1728  :
1729  src (src),
1730  matrix_diagonal_inverse (matrix_diagonal_inverse),
1731  do_startup (factor1 == Number()),
1732  start_zero (start_zero),
1733  factor1 (factor1),
1734  factor2 (factor2),
1735  update1 (update1),
1736  update2 (update2),
1737  dst (dst)
1738  {}
1739 
1740  void
1741  apply_to_subrange (const std::size_t begin,
1742  const std::size_t end) const
1743  {
1744  // To circumvent a bug in gcc
1745  // (https://gcc.gnu.org/bugzilla/show_bug.cgi?id=63945), we create copies
1746  // of the variables factor1 and factor2 and do not check based on
1747  // factor1.
1748  const Number factor1 = this->factor1;
1749  const Number factor2 = this->factor2;
1750  if (do_startup)
1751  {
1752  if (start_zero)
1753  DEAL_II_OPENMP_SIMD_PRAGMA
1754  for (std::size_t i=begin; i<end; ++i)
1755  {
1756  dst[i] = factor2 * src[i] * matrix_diagonal_inverse[i];
1757  update1[i] = -dst[i];
1758  }
1759  else
1760  DEAL_II_OPENMP_SIMD_PRAGMA
1761  for (std::size_t i=begin; i<end; ++i)
1762  {
1763  update1[i] = ((update2[i]-src[i]) *
1764  factor2*matrix_diagonal_inverse[i]);
1765  dst[i] -= update1[i];
1766  }
1767  }
1768  else
1769  DEAL_II_OPENMP_SIMD_PRAGMA
1770  for (std::size_t i=begin; i<end; ++i)
1771  {
1772  const Number update =
1773  factor1 * update1[i] + factor2 *
1774  ((update2[i] - src[i]) * matrix_diagonal_inverse[i]);
1775  update1[i] = update;
1776  dst[i] -= update;
1777  }
1778  }
1779 
1780  const Number *src;
1781  const Number *matrix_diagonal_inverse;
1782  const bool do_startup;
1783  const bool start_zero;
1784  const Number factor1;
1785  const Number factor2;
1786  mutable Number *update1;
1787  mutable Number *update2;
1788  mutable Number *dst;
1789  };
1790 
1791  template <typename Number>
1792  struct VectorUpdatesRange : public parallel::ParallelForInteger
1793  {
1794  VectorUpdatesRange(const VectorUpdater<Number> &updater,
1795  const std::size_t size)
1796  :
1797  updater (updater)
1798  {
1799  if (size < internal::VectorImplementation::minimum_parallel_grain_size)
1800  apply_to_subrange (0, size);
1801  else
1802  apply_parallel (0, size,
1803  internal::VectorImplementation::minimum_parallel_grain_size);
1804  }
1805 
1806  ~VectorUpdatesRange() = default;
1807 
1808  virtual void
1809  apply_to_subrange (const std::size_t begin,
1810  const std::size_t end) const
1811  {
1812  updater.apply_to_subrange(begin, end);
1813  }
1814 
1815  const VectorUpdater<Number> &updater;
1816  };
1817 
1818  // selection for diagonal matrix around deal.II vector
1819  template <typename Number>
1820  inline
1821  void
1822  vector_updates (const ::Vector<Number> &src,
1823  const DiagonalMatrix< ::Vector<Number> > &jacobi,
1824  const bool start_zero,
1825  const double factor1,
1826  const double factor2,
1827  ::Vector<Number> &update1,
1828  ::Vector<Number> &update2,
1829  ::Vector<Number> &,
1830  ::Vector<Number> &dst)
1831  {
1832  VectorUpdater<Number> upd(src.begin(), jacobi.get_vector().begin(),
1833  start_zero, factor1, factor2,
1834  update1.begin(), update2.begin(), dst.begin());
1835  VectorUpdatesRange<Number>(upd, src.size());
1836  }
1837 
1838  // selection for diagonal matrix around parallel deal.II vector
1839  template <typename Number>
1840  inline
1841  void
1842  vector_updates (const LinearAlgebra::distributed::Vector<Number> &src,
1844  const bool start_zero,
1845  const double factor1,
1846  const double factor2,
1851  {
1852  VectorUpdater<Number> upd(src.begin(), jacobi.get_vector().begin(),
1853  start_zero, factor1, factor2,
1854  update1.begin(), update2.begin(), dst.begin());
1855  VectorUpdatesRange<Number>(upd, src.local_size());
1856  }
1857 
1858  template <typename MatrixType, typename VectorType, typename PreconditionerType>
1859  inline
1860  void
1861  initialize_preconditioner(const MatrixType &matrix,
1862  std::shared_ptr<PreconditionerType> &preconditioner,
1863  VectorType &)
1864  {
1865  (void)matrix;
1866  (void)preconditioner;
1867  AssertThrow(preconditioner.get() != nullptr, ExcNotInitialized());
1868  }
1869 
1870  template <typename MatrixType, typename VectorType>
1871  inline
1872  void
1873  initialize_preconditioner(const MatrixType &matrix,
1874  std::shared_ptr<DiagonalMatrix<VectorType> > &preconditioner,
1875  VectorType &diagonal_inverse)
1876  {
1877  if (preconditioner.get() == nullptr ||
1878  preconditioner->m() != matrix.m())
1879  {
1880  if (preconditioner.get() == nullptr)
1881  preconditioner = std::make_shared<DiagonalMatrix<VectorType>>();
1882 
1883  Assert(preconditioner->m() == 0,
1884  ExcMessage("Preconditioner appears to be initialized but not sized correctly"));
1885 
1886  // Check if we can initialize from vector that then gets set to zero
1887  // as the matrix will own the memory
1888  preconditioner->reinit(diagonal_inverse);
1889  {
1890  VectorType empty_vector;
1891  diagonal_inverse.reinit(empty_vector);
1892  }
1893 
1894  // This part only works in serial
1895  if (preconditioner->m() != matrix.m())
1896  {
1897  preconditioner->get_vector().reinit(matrix.m());
1898  for (typename VectorType::size_type i=0; i<matrix.m(); ++i)
1899  preconditioner->get_vector()(i) = 1./matrix.el(i,i);
1900  }
1901  }
1902  }
1903 
1904  template <typename VectorType>
1905  void set_initial_guess(VectorType &vector)
1906  {
1907  vector = 1./std::sqrt(static_cast<double>(vector.size()));
1908  if (vector.locally_owned_elements().is_element(0))
1909  vector(0) = 0.;
1910  }
1911 
1912  template <typename Number>
1913  void set_initial_guess(::Vector<Number> &vector)
1914  {
1915  // Choose a high-frequency mode consisting of numbers between 0 and 1
1916  // that is cheap to compute (cheaper than random numbers) but avoids
1917  // obviously re-occurring numbers in multi-component systems by choosing
1918  // a period of 11
1919  for (unsigned int i=0; i<vector.size(); ++i)
1920  vector(i) = i%11;
1921 
1922  const Number mean_value = vector.mean_value();
1923  vector.add(-mean_value);
1924  }
1925 
1926  template <typename Number>
1927  void set_initial_guess(::LinearAlgebra::distributed::Vector<Number> &vector)
1928  {
1929  // Choose a high-frequency mode consisting of numbers between 0 and 1
1930  // that is cheap to compute (cheaper than random numbers) but avoids
1931  // obviously re-occurring numbers in multi-component systems by choosing
1932  // a period of 11.
1933  // Make initial guess robust with respect to number of processors
1934  // by operating on the global index.
1935  types::global_dof_index first_local_range = 0;
1936  if (!vector.locally_owned_elements().is_empty())
1937  first_local_range = vector.locally_owned_elements().nth_index_in_set(0);
1938  for (unsigned int i=0; i<vector.local_size(); ++i)
1939  vector.local_element(i) = (i+first_local_range)%11;
1940 
1941  const Number mean_value = vector.mean_value();
1942  vector.add(-mean_value);
1943  }
1944 
1945  struct EigenvalueTracker
1946  {
1947  public:
1948  void slot(const std::vector<double> &eigenvalues)
1949  {
1950  values = eigenvalues;
1951  }
1952 
1953  std::vector<double> values;
1954  };
1955  }
1956 }
1957 
1958 
1959 
1960 // avoid warning about deprecated variable nonzero_starting
1961 
1962 template <typename MatrixType, class VectorType, typename PreconditionerType>
1963 inline
1965 AdditionalData (const unsigned int degree,
1966  const double smoothing_range,
1967  const bool nonzero_starting,
1968  const unsigned int eig_cg_n_iterations,
1969  const double eig_cg_residual,
1970  const double max_eigenvalue)
1971  :
1972  degree (degree),
1973  smoothing_range (smoothing_range),
1974  nonzero_starting (nonzero_starting),
1975  eig_cg_n_iterations (eig_cg_n_iterations),
1976  eig_cg_residual (eig_cg_residual),
1977  max_eigenvalue (max_eigenvalue)
1978 {}
1979 
1980 
1981 
1982 template <typename MatrixType, typename VectorType, typename PreconditionerType>
1983 inline
1985  :
1986  theta (1.),
1987  delta (1.),
1988  eigenvalues_are_initialized (false)
1989 {
1990  static_assert(
1991  std::is_same<size_type, typename VectorType::size_type>::value,
1992  "PreconditionChebyshev and VectorType must have the same size_type.");
1993 }
1994 
1995 
1996 // avoid warning about deprecated variable AdditionalData::matrix_diagonal_inverse
1997 
1998 template <typename MatrixType, typename VectorType, typename PreconditionerType>
1999 inline
2000 void
2002 (const MatrixType &matrix,
2003  const AdditionalData &additional_data)
2004 {
2005  matrix_ptr = &matrix;
2006  data = additional_data;
2007  internal::PreconditionChebyshevImplementation::initialize_preconditioner(matrix,
2008  data.preconditioner,
2009  data.matrix_diagonal_inverse);
2010  eigenvalues_are_initialized = false;
2011 }
2012 
2013 
2014 
2015 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2016 inline
2017 void
2019 {
2020  eigenvalues_are_initialized = false;
2021  theta = delta = 1.0;
2022  matrix_ptr = nullptr;
2023  {
2024  VectorType empty_vector;
2025  data.matrix_diagonal_inverse.reinit(empty_vector);
2026  update1.reinit(empty_vector);
2027  update2.reinit(empty_vector);
2028  update3.reinit(empty_vector);
2029  }
2030  data.preconditioner.reset();
2031 }
2032 
2033 
2034 
2035 
2036 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2037 inline
2038 void
2040 (const VectorType &src) const
2041 {
2042  Assert(eigenvalues_are_initialized == false, ExcInternalError());
2043  Assert(data.preconditioner.get() != nullptr, ExcNotInitialized());
2044 
2045  update1.reinit(src);
2046  update2.reinit(src, true);
2047 
2048  // calculate largest eigenvalue using a hand-tuned CG iteration on the
2049  // matrix weighted by its diagonal. we start with a vector that consists of
2050  // ones only, weighted by the length.
2051  double max_eigenvalue, min_eigenvalue;
2052  if (data.eig_cg_n_iterations > 0)
2053  {
2054  Assert (data.eig_cg_n_iterations > 2,
2055  ExcMessage ("Need to set at least two iterations to find eigenvalues."));
2056 
2057  // set a very strict tolerance to force at least two iterations
2058  ReductionControl control (data.eig_cg_n_iterations,
2059  std::sqrt(std::numeric_limits<typename VectorType::value_type>::epsilon()),
2060  1e-10, false, false);
2061 
2062  internal::PreconditionChebyshevImplementation::EigenvalueTracker eigenvalue_tracker;
2063  SolverCG<VectorType> solver (control);
2064  solver.connect_eigenvalues_slot(std::bind(&internal::PreconditionChebyshevImplementation::EigenvalueTracker::slot,
2065  &eigenvalue_tracker,
2066  std::placeholders::_1));
2067 
2068  // set an initial guess which is close to the constant vector but where
2069  // one entry is different to trigger high frequencies
2070  internal::PreconditionChebyshevImplementation::set_initial_guess(update2);
2071 
2072  try
2073  {
2074  solver.solve(*matrix_ptr, update1, update2, *data.preconditioner);
2075  }
2077  {
2078  }
2079 
2080  // read the eigenvalues from the attached eigenvalue tracker
2081  if (eigenvalue_tracker.values.empty())
2082  min_eigenvalue = max_eigenvalue = 1;
2083  else
2084  {
2085  min_eigenvalue = eigenvalue_tracker.values.front();
2086 
2087  // include a safety factor since the CG method will in general not
2088  // be converged
2089  max_eigenvalue = 1.2*eigenvalue_tracker.values.back();
2090  }
2091  }
2092  else
2093  {
2094  max_eigenvalue = data.max_eigenvalue;
2095  min_eigenvalue = data.max_eigenvalue/data.smoothing_range;
2096  }
2097 
2098  const double alpha = (data.smoothing_range > 1. ?
2099  max_eigenvalue / data.smoothing_range :
2100  std::min(0.9*max_eigenvalue, min_eigenvalue));
2101 
2102  // in case the user set the degree to invalid unsigned int, we have to
2103  // determine the number of necessary iterations from the Chebyshev error
2104  // estimate, given the target tolerance specified by smoothing_range. This
2105  // estimate is based on the error formula given in section 5.1 of
2106  // R. S. Varga, Matrix iterative analysis, 2nd ed., Springer, 2009
2107  if (data.degree == numbers::invalid_unsigned_int)
2108  {
2109  const double actual_range = max_eigenvalue / alpha;
2110  const double sigma = (1.-std::sqrt(1./actual_range))/(1.+std::sqrt(1./actual_range));
2111  const double eps = data.smoothing_range;
2113  = 1+std::log(1./eps+std::sqrt(1./eps/eps-1))/std::log(1./sigma);
2114  }
2115 
2116  const_cast<PreconditionChebyshev<MatrixType,VectorType,PreconditionerType> *>(this)->delta = (max_eigenvalue-alpha)*0.5;
2117  const_cast<PreconditionChebyshev<MatrixType,VectorType,PreconditionerType> *>(this)->theta = (max_eigenvalue+alpha)*0.5;
2118 
2119  // We do not need the third auxiliary vector in case we have a
2120  // DiagonalMatrix as preconditioner and use deal.II's own vectors
2121  if (std::is_same<PreconditionerType,DiagonalMatrix<VectorType> >::value == false ||
2122  (std::is_same<VectorType,::Vector<typename VectorType::value_type> >::value == false
2123  &&
2124  std::is_same<VectorType,LinearAlgebra::distributed::Vector<typename VectorType::value_type> >::value == false
2125  ))
2126  update3.reinit (src, true);
2127 
2128  const_cast<PreconditionChebyshev<MatrixType,VectorType,PreconditionerType> *>(this)->eigenvalues_are_initialized = true;
2129 }
2130 
2131 
2132 
2133 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2134 inline
2135 void
2137 ::do_chebyshev_loop(VectorType &dst,
2138  const VectorType &src) const
2139 {
2140  // if delta is zero, we do not need to iterate because the updates will be
2141  // zero
2142  if (std::abs(delta) < 1e-40)
2143  return;
2144 
2145  double rhok = delta / theta, sigma = theta / delta;
2146  for (unsigned int k=0; k<data.degree; ++k)
2147  {
2148  matrix_ptr->vmult (update2, dst);
2149  const double rhokp = 1./(2.*sigma-rhok);
2150  const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
2151  rhok = rhokp;
2152  internal::PreconditionChebyshevImplementation::vector_updates
2153  (src, *data.preconditioner, false, factor1, factor2, update1, update2, update3, dst);
2154  }
2155 }
2156 
2157 
2158 
2159 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2160 inline
2161 void
2163 ::do_transpose_chebyshev_loop(VectorType &dst,
2164  const VectorType &src) const
2165 {
2166  double rhok = delta / theta, sigma = theta / delta;
2167  for (unsigned int k=0; k<data.degree; ++k)
2168  {
2169  matrix_ptr->Tvmult (update2, dst);
2170  const double rhokp = 1./(2.*sigma-rhok);
2171  const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
2172  rhok = rhokp;
2173  internal::PreconditionChebyshevImplementation::vector_updates
2174  (src, *data.preconditioner, false, factor1, factor2, update1, update2, update3, dst);
2175  }
2176 }
2177 
2178 
2179 
2180 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2181 inline
2182 void
2184 ::vmult (VectorType &dst,
2185  const VectorType &src) const
2186 {
2187  Threads::Mutex::ScopedLock lock(mutex);
2188  if (eigenvalues_are_initialized == false)
2189  estimate_eigenvalues(src);
2190 
2191  internal::PreconditionChebyshevImplementation::vector_updates
2192  (src, *data.preconditioner, true, 0., 1./theta, update1, update2, update3, dst);
2193 
2194  do_chebyshev_loop(dst, src);
2195 }
2196 
2197 
2198 
2199 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2200 inline
2201 void
2203 ::Tvmult (VectorType &dst,
2204  const VectorType &src) const
2205 {
2206  Threads::Mutex::ScopedLock lock(mutex);
2207  if (eigenvalues_are_initialized == false)
2208  estimate_eigenvalues(src);
2209 
2210  internal::PreconditionChebyshevImplementation::vector_updates
2211  (src, *data.preconditioner, true, 0., 1./theta, update1, update2, update3, dst);
2212 
2213  do_transpose_chebyshev_loop(dst, src);
2214 }
2215 
2216 
2217 
2218 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2219 inline
2220 void
2222 ::step (VectorType &dst,
2223  const VectorType &src) const
2224 {
2225  Threads::Mutex::ScopedLock lock(mutex);
2226  if (eigenvalues_are_initialized == false)
2227  estimate_eigenvalues(src);
2228 
2229  matrix_ptr->vmult (update2, dst);
2230  internal::PreconditionChebyshevImplementation::vector_updates
2231  (src, *data.preconditioner, false, 0., 1./theta, update1, update2, update3, dst);
2232 
2233  do_chebyshev_loop(dst, src);
2234 }
2235 
2236 
2237 
2238 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2239 inline
2240 void
2242 ::Tstep (VectorType &dst,
2243  const VectorType &src) const
2244 {
2245  Threads::Mutex::ScopedLock lock(mutex);
2246  if (eigenvalues_are_initialized == false)
2247  estimate_eigenvalues(src);
2248 
2249  matrix_ptr->Tvmult (update2, dst);
2250  internal::PreconditionChebyshevImplementation::vector_updates
2251  (src, *data.preconditioner, false, 0., 1./theta, update1, update2, update3, dst);
2252 
2253  do_transpose_chebyshev_loop(dst, src);
2254 }
2255 
2256 
2257 
2258 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2259 inline
2262 {
2263  Assert (matrix_ptr!=nullptr, ExcNotInitialized());
2264  return matrix_ptr->m();
2265 }
2266 
2267 
2268 
2269 template <typename MatrixType, typename VectorType, typename PreconditionerType>
2270 inline
2273 {
2274  Assert (matrix_ptr!=nullptr, ExcNotInitialized());
2275  return matrix_ptr->n();
2276 }
2277 
2278 #endif // DOXYGEN
2279 
2280 DEAL_II_NAMESPACE_CLOSE
2281 
2282 #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:173
size_type m() const
void initialize(const MatrixType &A, const typename BaseClass::AdditionalData &parameters=typename BaseClass::AdditionalData())
const_iterator end() const
Contents is actually a matrix.
size_type nth_index_in_set(const unsigned int local_index) const
Definition: index_set.h:1769
void Tstep(VectorType &x, const VectorType &rhs) const
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)
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
std::shared_ptr< PreconditionerType > preconditioner
Definition: precondition.h:998
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
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
std::array< std::pair< Number, Tensor< 1, dim, Number > >, dim > jacobi(::SymmetricTensor< 2, dim, Number > A)
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
const std::vector< size_type > & inverse_permutation
Definition: precondition.h:752
#define Assert(cond, exc)
Definition: exceptions.h:1142
virtual void add(const Number a) override
void step(VectorType &x, const VectorType &rhs) const
void(MatrixType::* function_ptr)(VectorType &, const VectorType &) const
Definition: precondition.h:344
virtual ::IndexSet locally_owned_elements() const override
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
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
virtual Number mean_value() const override
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
bool is_empty() const
Definition: index_set.h:1708
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())