/usr/include/deal.II/lac/precondition.h is in libdeal.ii-dev 6.3.1-1.1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542 1543 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574 1575 1576 1577 1578 1579 1580 1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 1662 1663 1664 1665 1666 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717 1718 1719 1720 1721 1722 1723 1724 1725 1726 1727 1728 1729 1730 1731 1732 1733 1734 1735 1736 1737 1738 1739 1740 1741 1742 1743 1744 1745 1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756 1757 1758 1759 1760 1761 1762 1763 1764 1765 1766 1767 1768 1769 1770 1771 1772 1773 1774 1775 | //---------------------------------------------------------------------------
// $Id: precondition.h 20819 2010-03-12 23:15:56Z allmaras $
// Version: $Name$
//
// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2009, 2010 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
// to the file deal.II/doc/license.html for the text and
// further information on this license.
//
//---------------------------------------------------------------------------
#ifndef __deal2__precondition_h
#define __deal2__precondition_h
// This file contains simple preconditioners.
#include <base/config.h>
#include <base/smartpointer.h>
#include <base/template_constraints.h>
#include <lac/tridiagonal_matrix.h>
#include <lac/vector_memory.h>
DEAL_II_NAMESPACE_OPEN
template <typename number> class Vector;
template <typename number> class SparseMatrix;
/*! @addtogroup Preconditioners
*@{
*/
/**
* No preconditioning. This class helps you, if you want to use a
* linear solver without preconditioning. All solvers in LAC require a
* preconditioner. Therefore, you must use the identity provided here
* to avoid preconditioning. It can be used in the following way:
*
@verbatim
SolverControl solver_control (1000, 1e-12);
SolverCG<> cg (solver_control);
cg.solve (system_matrix, solution, system_rhs,
PreconditionIdentity());
@endverbatim
*
* See the step-3 tutorial program for an example and
* additional explanations.
*
* Alternatively, the IdentityMatrix class can be used to precondition
* in this way.
*
* @author Guido Kanschat, 1999
*/
class PreconditionIdentity : public Subscriptor
{
public:
/**
* This function is only
* present to
* provide the interface of
* a precondtioner to be
* handed to a smoother.
* This does nothing.
*/
struct AdditionalData
{
/**
* Constructor.
*/
AdditionalData (){}
};
/**
* The matrix
* argument is ignored and here
* just for compatibility with
* more complex preconditioners.
*/
template <class MATRIX>
void initialize (const MATRIX &matrix,
const AdditionalData &additional_data = AdditionalData());
/**
* Apply preconditioner.
*/
template<class VECTOR>
void vmult (VECTOR&, const VECTOR&) const;
/**
* Apply transpose
* preconditioner. Since this is
* the identity, this function is
* the same as
* vmult().
*/
template<class VECTOR>
void Tvmult (VECTOR&, const VECTOR&) const;
/**
* Apply preconditioner, adding to the previous value.
*/
template<class VECTOR>
void vmult_add (VECTOR&, const VECTOR&) const;
/**
* Apply transpose
* preconditioner, adding. Since this is
* the identity, this function is
* the same as
* vmult_add().
*/
template<class VECTOR>
void Tvmult_add (VECTOR&, const VECTOR&) const;
/**
* This function is only
* present to
* provide the interface of
* a precondtioner to be
* handed to a smoother.
* This does nothing.
*/
void clear (){}
};
/**
* Preconditioning with Richardson's method. This preconditioner just
* scales the vector with a constant relaxation factor provided by the
* AdditionalData object.
*
* In Krylov-space methods, this preconditioner should not have any
* effect. Using SolverRichardson, the two relaxation parameters will
* be just multiplied. Still, this class is useful in multigrid
* smoother objects (MGSmootherRelaxation).
*
* @author Guido Kanschat, 2005
*/
class PreconditionRichardson : public Subscriptor
{
public:
/**
* Parameters for Richardson
* preconditioner.
*/
class AdditionalData
{
public:
/**
* Constructor. Block size
* must be given since there
* is no reasonable default
* parameter.
*/
AdditionalData (const double relaxation = 1.);
/**
* Relaxation parameter.
*/
double relaxation;
};
/**
* Constructor, sets the
* relaxation parameter to its
* default.
*/
PreconditionRichardson();
/**
* Change the relaxaton parameter.
*/
void initialize (const AdditionalData ¶meters);
/**
* Change the relaxaton parameter
* in a way consistent with other
* preconditioners. The matrix
* argument is ignored and here
* just for compatibility with
* more complex preconditioners.
*/
template <class MATRIX>
void initialize (const MATRIX&,
const AdditionalData ¶meters);
/**
* Apply preconditioner.
*/
template<class VECTOR>
void vmult (VECTOR&, const VECTOR&) const;
/**
* Apply transpose
* preconditioner. Since this is
* the identity, this function is
* the same as
* vmult().
*/
template<class VECTOR>
void Tvmult (VECTOR&, const VECTOR&) const;
/**
* Apply preconditioner, adding to the previous value.
*/
template<class VECTOR>
void vmult_add (VECTOR&, const VECTOR&) const;
/**
* Apply transpose
* preconditioner, adding. Since this is
* the identity, this function is
* the same as
* vmult_add().
*/
template<class VECTOR>
void Tvmult_add (VECTOR&, const VECTOR&) const;
/**
* This function is only
* present to
* provide the interface of
* a precondtioner to be
* handed to a smoother.
* This does nothing.
*/
void clear (){}
private:
/**
* The relaxation parameter
* multiplied with the vectors.
*/
double relaxation;
};
/**
* Preconditioner using a matrix-builtin function.
* This class forms a preconditioner suitable for the LAC solver
* classes. Since many preconditioning methods are based on matrix
* entries, these have to be implemented as member functions of the
* underlying matrix implementation. This class now is intended to
* allow easy access to these member functions from LAC solver
* classes.
*
* It seems that all builtin preconditioners have a relaxation
* parameter, so please use PreconditionRelaxation for these.
*
* You will usually not want to create a named object of this type,
* although possible. The most common use is like this:
* @code
* SolverGMRES<SparseMatrix<double>,
* Vector<double> > gmres(control,memory,500);
*
* gmres.solve (matrix, solution, right_hand_side,
* PreconditionUseMatrix<SparseMatrix<double>,Vector<double> >
* (matrix,&SparseMatrix<double>::template precondition_Jacobi<double>));
* @endcode
* This creates an unnamed object to be passed as the fourth parameter to
* the solver function of the SolverGMRES class. It assumes that the
* SparseMatrix class has a function <tt>precondition_Jacobi</tt> taking two
* vectors (source and destination) as parameters (Actually, there is no
* function like that, the existing function takes a third parameter,
* denoting the relaxation parameter; this example is therefore only meant to
* illustrate the general idea).
*
* Note that due to the default template parameters, the above example
* could be written shorter as follows:
* @code
* ...
* gmres.solve (matrix, solution, right_hand_side,
* PreconditionUseMatrix<>
* (matrix,&SparseMatrix<double>::template precondition_Jacobi<double>));
* @endcode
*
* @author Guido Kanschat, Wolfgang Bangerth, 1999
*/
template<class MATRIX = SparseMatrix<double>, class VECTOR = Vector<double> >
class PreconditionUseMatrix : public Subscriptor
{
public:
/**
* Type of the preconditioning
* function of the matrix.
*/
typedef void ( MATRIX::* function_ptr)(VECTOR&, const VECTOR&) const;
/**
* Constructor.
* This constructor stores a
* reference to the matrix object
* for later use and selects a
* preconditioning method, which
* must be a member function of
* that matrix.
*/
PreconditionUseMatrix(const MATRIX &M,
const function_ptr method);
/**
* Execute preconditioning. Calls the
* function passed to the constructor
* of this object with the two
* arguments given here.
*/
void vmult (VECTOR &dst,
const VECTOR &src) const;
private:
/**
* Pointer to the matrix in use.
*/
const MATRIX &matrix;
/**
* Pointer to the preconditioning
* function.
*/
const function_ptr precondition;
};
/**
* Base class for other preconditioners.
* Here, only some common features Jacobi, SOR and SSOR preconditioners
* are implemented. For preconditioning, refer to derived classes.
*
* @author Guido Kanschat, 2000
*/
template<class MATRIX = SparseMatrix<double> >
class PreconditionRelaxation : public Subscriptor
{
public:
/**
* Class for parameters.
*/
class AdditionalData
{
public:
/**
* Constructor.
*/
AdditionalData (const double relaxation = 1.);
/**
* Relaxation parameter.
*/
double relaxation;
};
/**
* Initialize matrix and
* relaxation parameter. The
* matrix is just stored in the
* preconditioner object. The
* relaxation parameter should be
* larger than zero and smaller
* than 2 for numerical
* reasons. It defaults to 1.
*/
void initialize (const MATRIX &A,
const AdditionalData & parameters = AdditionalData());
/**
* Release the matrix and reset
* its pointer.
*/
void clear();
protected:
/**
* Pointer to the matrix object.
*/
SmartPointer<const MATRIX, PreconditionRelaxation<MATRIX> > A;
/**
* Relaxation parameter.
*/
double relaxation;
};
/**
* Jacobi preconditioner using matrix built-in function. The
* <tt>MATRIX</tt> class used is required to have a function
* <tt>precondition_Jacobi(VECTOR&, const VECTOR&, double</tt>)
*
* @code
* // Declare related objects
*
* SparseMatrix<double> A;
* Vector<double> x;
* Vector<double> b;
* SolverCG<> solver(...);
*
* //...initialize and build A
*
* // Define and initialize preconditioner
*
* PreconditionJacobi<SparseMatrix<double> > precondition;
* precondition.initialize (A, .6);
*
* solver.solve (A, x, b, precondition);
* @endcode
*
* @author Guido Kanschat, 2000
*/
template <class MATRIX = SparseMatrix<double> >
class PreconditionJacobi : public PreconditionRelaxation<MATRIX>
{
public:
/**
* Apply preconditioner.
*/
template<class VECTOR>
void vmult (VECTOR&, const VECTOR&) const;
/**
* Apply transpose
* preconditioner. Since this is
* a symmetric preconditioner,
* this function is the same as
* vmult().
*/
template<class VECTOR>
void Tvmult (VECTOR&, const VECTOR&) const;
/**
* Perform one step of the
* preconditioned Richardson
* iteration.
*/
template<class VECTOR>
void step (VECTOR& x, const VECTOR& rhs) const;
/**
* Perform one transposed step of
* the preconditioned Richardson
* iteration.
*/
template<class VECTOR>
void Tstep (VECTOR& x, const VECTOR& rhs) const;
};
/**
* SOR preconditioner using matrix built-in function.
*
* Assuming the matrix <i>A = D + L + U</i> is split into its diagonal
* <i>D</i> as well as the strict lower and upper triangles <i>L</i>
* and <i>U</i>, then the SOR preconditioner with relaxation parameter
* <i>r</i> is
* @f[
* P^{-1} = r (D+rL)^{-1}.
* @f]
* It is this operator <i>P<sup>-1</sup></i>, which is implemented by
* vmult() through forward substitution. Analogously, Tvmult()
* implements the operation of <i>r(D+rU)<sup>-1</sup></i>.
*
* The SOR iteration itself can be directly written as
* @f[
* x^{k+1} = x^k - r D^{-1} \bigl(L x^{k+1} + U x^k - b\bigr).
* @f]
* Using the right hand side <i>b</i> and the previous iterate
* <i>x</i>, this is the operation implemented by step().
*
* The MATRIX
* class used is required to have functions
* <tt>precondition_SOR(VECTOR&, const VECTOR&, double)</tt> and
* <tt>precondition_TSOR(VECTOR&, const VECTOR&, double)</tt>.
*
* @code
* // Declare related objects
*
* SparseMatrix<double> A;
* Vector<double> x;
* Vector<double> b;
* SolverCG<> solver(...);
*
* //...initialize and build A
*
* // Define and initialize preconditioner
*
* PreconditionSOR<SparseMatrix<double> > precondition;
* precondition.initialize (A, .6);
*
* solver.solve (A, x, b, precondition);
* @endcode
*
* @author Guido Kanschat, 2000
*/
template <class MATRIX = SparseMatrix<double> >
class PreconditionSOR : public PreconditionRelaxation<MATRIX>
{
public:
/**
* Apply preconditioner.
*/
template<class VECTOR>
void vmult (VECTOR&, const VECTOR&) const;
/**
* Apply transpose
* preconditioner.
*/
template<class VECTOR>
void Tvmult (VECTOR&, const VECTOR&) const;
/**
* Perform one step of the
* preconditioned Richardson
* iteration.
*/
template<class VECTOR>
void step (VECTOR& x, const VECTOR& rhs) const;
/**
* Perform one transposed step of
* the preconditioned Richardson
* iteration.
*/
template<class VECTOR>
void Tstep (VECTOR& x, const VECTOR& rhs) const;
};
/**
* SSOR preconditioner using matrix built-in function. The
* <tt>MATRIX</tt> class used is required to have a function
* <tt>precondition_SSOR(VECTOR&, const VECTOR&, double)</tt>
*
* @code
* // Declare related objects
*
* SparseMatrix<double> A;
* Vector<double> x;
* Vector<double> b;
* SolverCG<> solver(...);
*
* //...initialize and build A
*
* // Define and initialize preconditioner
*
* PreconditionSSOR<SparseMatrix<double> > precondition;
* precondition.initialize (A, .6);
*
* solver.solve (A, x, b, precondition);
* @endcode
*
* @author Guido Kanschat, 2000
*/
template <class MATRIX = SparseMatrix<double> >
class PreconditionSSOR : public PreconditionRelaxation<MATRIX>
{
public:
/**
* A typedef to the base class.
*/
typedef PreconditionRelaxation<MATRIX> BaseClass;
/**
* Initialize matrix and
* relaxation parameter. The
* matrix is just stored in the
* preconditioner object. The
* relaxation parameter should be
* larger than zero and smaller
* than 2 for numerical
* reasons. It defaults to 1.
*/
void initialize (const MATRIX &A,
const typename BaseClass::AdditionalData ¶meters = typename BaseClass::AdditionalData());
/**
* Apply preconditioner.
*/
template<class VECTOR>
void vmult (VECTOR&, const VECTOR&) const;
/**
* Apply transpose
* preconditioner. Since this is
* a symmetric preconditioner,
* this function is the same as
* vmult().
*/
template<class VECTOR>
void Tvmult (VECTOR&, const VECTOR&) const;
/**
* Perform one step of the
* preconditioned Richardson
* iteration
*/
template<class VECTOR>
void step (VECTOR& x, const VECTOR& rhs) const;
/**
* Perform one transposed step of
* the preconditioned Richardson
* iteration.
*/
template<class VECTOR>
void Tstep (VECTOR& x, const VECTOR& rhs) const;
private:
/**
* An array that stores for each matrix
* row where the first position after
* the diagonal is located.
*/
std::vector<unsigned int> pos_right_of_diagonal;
};
/**
* Permuted SOR preconditioner using matrix built-in function. The
* <tt>MATRIX</tt> class used is required to have functions
* <tt>PSOR(VECTOR&, const VECTOR&, double)</tt> and
* <tt>TPSOR(VECTOR&, const VECTOR&, double)</tt>.
*
* @code
* // Declare related objects
*
* SparseMatrix<double> A;
* Vector<double> x;
* Vector<double> b;
* SolverCG<> solver(...);
*
* //...initialize and build A
*
* std::vector<unsigned int> permutation(x.size());
* std::vector<unsigned int> inverse_permutation(x.size());
*
* //...fill permutation and its inverse with reasonable values
*
* // Define and initialize preconditioner
*
* PreconditionPSOR<SparseMatrix<double> > precondition;
* precondition.initialize (A, permutation, inverse_permutation, .6);
*
* solver.solve (A, x, b, precondition);
* @endcode
*
* @author Guido Kanschat, 2003
*/
template <class MATRIX = SparseMatrix<double> >
class PreconditionPSOR : public PreconditionRelaxation<MATRIX>
{
public:
/**
* Initialize matrix and
* relaxation parameter. The
* matrix is just stored in the
* preconditioner object.
*
* The permutation vector is
* stored as a
* pointer. Therefore, it has to
* be assured that the lifetime
* of the vector exceeds the
* lifetime of the
* preconditioner.
*
* The relaxation parameter
* should be larger than zero and
* smaller than 2 for numerical
* reasons. It defaults to 1.
*/
void initialize (const MATRIX &A,
const std::vector<unsigned int> &permutation,
const std::vector<unsigned int> &inverse_permutation,
const typename PreconditionRelaxation<MATRIX>::AdditionalData &
parameters = typename PreconditionRelaxation<MATRIX>::AdditionalData());
/**
* Apply preconditioner.
*/
template<class VECTOR>
void vmult (VECTOR&, const VECTOR&) const;
/**
* Apply transpose
* preconditioner.
*/
template<class VECTOR>
void Tvmult (VECTOR&, const VECTOR&) const;
private:
/**
* Storage for the permutation vector.
*/
const std::vector<unsigned int>* permutation;
/**
* Storage for the inverse
* permutation vector.
*/
const std::vector<unsigned int>* inverse_permutation;
};
/**
* Preconditioner using an iterative solver. This preconditioner uses
* a fully initialized LAC iterative solver for the approximate
* inverse of the matrix. Naturally, this solver needs another
* preconditionig method.
*
* Usually, the use of ReductionControl is preferred over the use of
* the basic SolverControl in defining this solver.
*
* Krylov space methods like SolverCG or SolverBicgstab
* become inefficient if soution down to machine accuracy is
* needed. This is due to the fact, that round-off errors spoil the
* orthogonality of the vector sequences. Therefore, a nested
* iteration of two methods is proposed: The outer method is
* SolverRichardson, since it is robust with respect to round-of
* errors. The inner loop is an appropriate Krylov space method, since
* it is fast.
*
* @code
* // Declare related objects
*
* SparseMatrix<double> A;
* Vector<double> x;
* Vector<double> b;
* GrowingVectorMemory<Vector<double> > mem;
* ReductionControl inner_control (10, 1.e-30, 1.e-2)
* SolverCG<Vector<double> > inner_iteration (inner_control, mem);
* PreconditionSSOR <SparseMatrix<double> > inner_precondition;
* inner_precondition.initialize (A, 1.2);
*
* PreconditionLACSolver precondition;
* precondition.initialize (inner_iteration, A, inner_precondition);
*
* SolverControl outer_control(100, 1.e-16);
* SolverRichardson<Vector<double> > outer_iteration;
*
* outer_iteration.solve (A, x, b, precondition);
* @endcode
*
* Each time we call the inner loop, reduction of the residual by a
* factor <tt>1.e-2</tt> is attempted. Since the right hand side vector of
* the inner iteration is the residual of the outer loop, the relative
* errors are far from machine accuracy, even if the errors of the
* outer loop are in the range of machine accuracy.
*
* @author Guido Kanschat, 1999
*/
template<class SOLVER, class MATRIX = SparseMatrix<double>, class PRECONDITION = PreconditionIdentity>
class PreconditionLACSolver : public Subscriptor
{
public:
/**
* Constructor. All work is done
* in initialize.
*/
PreconditionLACSolver ();
/**
* Initialization
* function. Provide a solver
* object, a matrix, and another
* preconditioner for this.
*/
void initialize (SOLVER&,
const MATRIX&,
const PRECONDITION&);
/**
* Execute preconditioning.
*/
template<class VECTOR>
void vmult (VECTOR&, const VECTOR&) const;
private:
/**
* The solver object to use.
*/
SmartPointer<SOLVER,PreconditionLACSolver<SOLVER,MATRIX,PRECONDITION> > solver;
/**
* The matrix in use.
*/
SmartPointer<const MATRIX,PreconditionLACSolver<SOLVER,MATRIX,PRECONDITION> > matrix;
/**
* The preconditioner to use.
*/
SmartPointer<const PRECONDITION,PreconditionLACSolver<SOLVER,MATRIX,PRECONDITION> > precondition;
};
/**
* @deprecated Use ProductMatrix instead.
*
* Matrix with preconditioner.
* Given a matrix $A$ and a preconditioner $P$, this class implements a new matrix
* with the matrix-vector product $PA$. It needs an auxiliary vector for that.
*
* By this time, this is considered a temporary object to be plugged
* into eigenvalue solvers. Therefore, no SmartPointer is used for
* <tt>A</tt> and <tt>P</tt>.
*
* @author Guido Kanschat, 2000
*/
template<class MATRIX, class PRECOND, class VECTOR>
class PreconditionedMatrix : public Subscriptor
{
public:
/**
* Constructor. Provide matrix,
* preconditioner and a memory
* pool to obtain the auxiliary
* vector.
*/
PreconditionedMatrix (const MATRIX & A,
const PRECOND & P,
VectorMemory<VECTOR> & mem);
/**
* Preconditioned
* matrix-vector-product.
*/
void vmult (VECTOR &dst, const VECTOR &src) const;
/**
* Transposed preconditioned
* matrix-vector-product.
*/
void Tvmult (VECTOR &dst, const VECTOR &src) const;
/**
* Residual $b-PAx$.
*/
double residual (VECTOR &dst, const VECTOR &src, const VECTOR &rhs) const;
private:
/**
* Storage for the matrix.
*/
const MATRIX &A;
/**
* Storage for preconditioner.
*/
const PRECOND &P;
/**
* Memory pool for vectors.
*/
VectorMemory<VECTOR> &mem;
};
/**
* Preconditioning with a Chebyshev polynomial for symmetric positive
* definite matrices. This preconditioner is similar to a Jacobi
* preconditioner if the degree variable is set to one, otherwise some
* higher order polynomial corrections are used. This preconditioner needs
* access to the diagonal of the matrix its acts on and needs a respective
* <tt>vmult</tt> implemention. However, it does not need to explicitly know
* the matrix entries.
*
* This class is useful e.g. in multigrid smoother objects, since it is
* trivially %parallel (assuming that matrix-vector products are %parallel).
*
* @author Martin Kronbichler, 2009
*/
template <class MATRIX=SparseMatrix<double>, class VECTOR=Vector<double> >
class PreconditionChebyshev : public Subscriptor
{
public:
/**
* Standardized data struct to
* pipe additional parameters
* to the preconditioner.
*/
struct AdditionalData
{
/**
* Constructor.
*/
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);
/**
* This determines the degree of the
* Chebyshev polynomial. The degree
* of the polynomial gives the number
* of matrix-vector products to be
* performed for one application of
* the vmult() operation. Degree zero
* corresponds to a damped Jacobi
* method.
*/
unsigned int degree;
/**
* This sets the range between the
* largest eigenvalue in the matrix
* and the smallest eigenvalue to be
* treated. If the parameter is zero,
* an estimate for the largest and
* for the smallest eigenvalue will
* be calculated
* internally. Otherwise, the
* Chebyshev polynomial will act in
* the interval
* $[\lambda_\mathrm{max}/
* \tt{smoothing\_range},
* \lambda_\mathrm{max}]$, where
* $\lambda_\mathrm{max}$ is an
* estimate of the maximum eigenvalue
* of the matrix. A choice of
* <tt>smoothing_range</tt> between 5
* and 20 is useful in case the
* preconditioner is used as a
* smoother in multigrid.
*/
double smoothing_range;
/**
* When this flag is set to
* <tt>true</tt>, it enables the
* method <tt>vmult(dst, src)</tt> to
* use non-zero data in the vector
* <tt>dst</tt>, appending to it the
* Chebyshev corrections. This can be
* useful in some situations
* (e.g. when used for high-frequency
* error smoothing in a multigrid
* algorithm), but not the way the
* solver classes expect a
* preconditioner to work (where one
* ignores the content in
* <tt>dst</tt> for the
* preconditioner application).
*/
bool nonzero_starting;
/**
* Maximum number of CG iterations
* performed for finding the maximum
* eigenvalue.
*/
unsigned int eig_cg_n_iterations;
/**
* Tolerance for CG iterations
* performed for finding the maximum
* eigenvalue.
*/
double eig_cg_residual;
/**
* Stores the inverse of the diagonal
* of the underlying matrix.
*/
VECTOR matrix_diagonal_inverse;
};
PreconditionChebyshev ();
/**
* Initialize function. Takes the
* matrix which is used to form the
* preconditioner, and additional
* flags if there are any. This
* function works only if the input
* matrix has an operator
* <tt>el(i,i)</tt> for accessing all
* the elements in the
* diagonal. Alternatively, the
* diagonal can be supplied with the
* help of the AdditionalData field.
*
* This function calculates an
* estimate of the eigenvalue range
* of the matrix weighted by its
* diagonal using a modified CG
* iteration.
*/
void initialize (const MATRIX &matrix,
const AdditionalData &additional_data = AdditionalData());
/**
* Computes the action of the
* preconditioner on <tt>src</tt>,
* storing the result in
* <tt>dst</tt>.
*/
void vmult (VECTOR &dst,
const VECTOR &src) const;
/**
* Computes the action of the
* transposed preconditioner on
* <tt>src</tt>, storing the result
* in <tt>dst</tt>.
*/
void Tvmult (VECTOR &dst,
const VECTOR &src) const;
/**
* Resets the preconditioner.
*/
void clear ();
private:
/**
* A pointer to the underlying
* matrix.
*/
SmartPointer<const MATRIX,PreconditionChebyshev<MATRIX,VECTOR> > matrix_ptr;
/**
* Internal vector used for the
* <tt>vmult</tt> operation.
*/
mutable VECTOR update1;
/**
* Internal vector used for the
* <tt>vmult</tt> operation.
*/
mutable VECTOR update2;
/**
* Stores the additional data
* provided to the initialize
* function.
*/
AdditionalData data;
/**
* Average of the largest and
* smallest eigenvalue under
* consideration.
*/
double theta;
/**
* Half the interval length between
* the largest and smallest
* eigenvalue under consideration.
*/
double delta;
/**
* Stores whether the preconditioner
* has been set up.
*/
bool is_initialized;
};
/*@}*/
/* ---------------------------------- Inline functions ------------------- */
#ifndef DOXYGEN
template <class MATRIX>
inline void
PreconditionIdentity::initialize (
const MATRIX&,
const PreconditionIdentity::AdditionalData&)
{}
template<class VECTOR>
inline void
PreconditionIdentity::vmult (VECTOR &dst, const VECTOR &src) const
{
dst = src;
}
template<class VECTOR>
inline void
PreconditionIdentity::Tvmult (VECTOR &dst, const VECTOR &src) const
{
dst = src;
}
template<class VECTOR>
inline void
PreconditionIdentity::vmult_add (VECTOR &dst, const VECTOR &src) const
{
dst.add(src);
}
template<class VECTOR>
inline void
PreconditionIdentity::Tvmult_add (VECTOR &dst, const VECTOR &src) const
{
dst.add(src);
}
//---------------------------------------------------------------------------
inline
PreconditionRichardson::AdditionalData::AdditionalData (
const double relaxation)
:
relaxation(relaxation)
{}
inline
PreconditionRichardson::PreconditionRichardson ()
:
relaxation(0)
{
AdditionalData add_data;
relaxation=add_data.relaxation;
}
inline void
PreconditionRichardson::initialize (
const PreconditionRichardson::AdditionalData ¶meters)
{
relaxation = parameters.relaxation;
}
template <class MATRIX>
inline void
PreconditionRichardson::initialize (
const MATRIX&,
const PreconditionRichardson::AdditionalData ¶meters)
{
relaxation = parameters.relaxation;
}
template<class VECTOR>
inline void
PreconditionRichardson::vmult (VECTOR &dst, const VECTOR &src) const
{
dst.equ(relaxation,src);
}
template<class VECTOR>
inline void
PreconditionRichardson::Tvmult (VECTOR &dst, const VECTOR &src) const
{
dst.equ(relaxation,src);
}
template<class VECTOR>
inline void
PreconditionRichardson::vmult_add (VECTOR &dst, const VECTOR &src) const
{
dst.add(relaxation,src);
}
template<class VECTOR>
inline void
PreconditionRichardson::Tvmult_add (VECTOR &dst, const VECTOR &src) const
{
dst.add(relaxation,src);
}
//---------------------------------------------------------------------------
template <class MATRIX>
inline void
PreconditionRelaxation<MATRIX>::initialize (const MATRIX &rA,
const AdditionalData ¶meters)
{
A = &rA;
relaxation = parameters.relaxation;
}
template <class MATRIX>
inline void
PreconditionRelaxation<MATRIX>::clear ()
{
A = 0;
}
//---------------------------------------------------------------------------
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionJacobi<MATRIX>::vmult (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->precondition_Jacobi (dst, src, this->relaxation);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionJacobi<MATRIX>::Tvmult (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->precondition_Jacobi (dst, src, this->relaxation);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionJacobi<MATRIX>::step (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->Jacobi_step (dst, src, this->relaxation);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionJacobi<MATRIX>::Tstep (VECTOR &dst, const VECTOR &src) const
{
step (dst, src);
}
//---------------------------------------------------------------------------
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionSOR<MATRIX>::vmult (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->precondition_SOR (dst, src, this->relaxation);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionSOR<MATRIX>::Tvmult (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->precondition_TSOR (dst, src, this->relaxation);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionSOR<MATRIX>::step (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->SOR_step (dst, src, this->relaxation);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionSOR<MATRIX>::Tstep (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->TSOR_step (dst, src, this->relaxation);
}
//---------------------------------------------------------------------------
template <class MATRIX>
inline void
PreconditionSSOR<MATRIX>::initialize (const MATRIX &rA,
const typename BaseClass::AdditionalData ¶meters)
{
this->PreconditionRelaxation<MATRIX>::initialize (rA, parameters);
// in case we have a SparseMatrix class,
// we can extract information about the
// diagonal.
const SparseMatrix<typename MATRIX::value_type> * mat =
dynamic_cast<const SparseMatrix<typename MATRIX::value_type> *>(&*this->A);
// calculate the positions first after
// the diagonal.
if (mat != 0)
{
const std::size_t * rowstart_ptr =
mat->get_sparsity_pattern().get_rowstart_indices();
const unsigned int * const colnums =
mat->get_sparsity_pattern().get_column_numbers();
const unsigned int n = this->A->n();
pos_right_of_diagonal.resize(n);
for (unsigned int row=0; row<n; ++row, ++rowstart_ptr)
{
// find the first element in this line
// which is on the right of the diagonal.
// we need to precondition with the
// elements on the left only.
// note: the first entry in each
// line denotes the diagonal element,
// which we need not check.
pos_right_of_diagonal[row] =
std::lower_bound (&colnums[*rowstart_ptr+1],
&colnums[*(rowstart_ptr+1)],
row)
- colnums;
}
}
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionSSOR<MATRIX>::vmult (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionSSOR<MATRIX>::Tvmult (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionSSOR<MATRIX>::step (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
this->A->SSOR_step (dst, src, this->relaxation);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionSSOR<MATRIX>::Tstep (VECTOR &dst, const VECTOR &src) const
{
step (dst, src);
}
//---------------------------------------------------------------------------
template <class MATRIX>
inline void
PreconditionPSOR<MATRIX>::initialize (
const MATRIX &rA,
const std::vector<unsigned int> &p,
const std::vector<unsigned int> &ip,
const typename PreconditionRelaxation<MATRIX>::AdditionalData ¶meters)
{
permutation = &p;
inverse_permutation = &ip;
PreconditionRelaxation<MATRIX>::initialize(rA, parameters);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionPSOR<MATRIX>::vmult (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
dst = src;
this->A->PSOR (dst, *permutation, *inverse_permutation, this->relaxation);
}
template <class MATRIX>
template<class VECTOR>
inline void
PreconditionPSOR<MATRIX>::Tvmult (VECTOR &dst, const VECTOR &src) const
{
Assert (this->A!=0, ExcNotInitialized());
dst = src;
this->A->TPSOR (dst, *permutation, *inverse_permutation, this->relaxation);
}
//---------------------------------------------------------------------------
template<class MATRIX, class VECTOR>
PreconditionUseMatrix<MATRIX,VECTOR>::PreconditionUseMatrix(const MATRIX &M,
const function_ptr method)
:
matrix(M), precondition(method)
{}
template<class MATRIX, class VECTOR>
void
PreconditionUseMatrix<MATRIX,VECTOR>::vmult (VECTOR &dst,
const VECTOR &src) const
{
(matrix.*precondition)(dst, src);
}
//---------------------------------------------------------------------------
template<class MATRIX>
inline
PreconditionRelaxation<MATRIX>::AdditionalData::
AdditionalData (const double relaxation)
:
relaxation (relaxation)
{}
//////////////////////////////////////////////////////////////////////
template<class SOLVER, class MATRIX, class PRECONDITION>
PreconditionLACSolver<SOLVER,MATRIX,PRECONDITION>
::PreconditionLACSolver ()
:
solver(0), matrix(0), precondition(0)
{}
template<class SOLVER, class MATRIX, class PRECONDITION>
void
PreconditionLACSolver<SOLVER,MATRIX,PRECONDITION>
::initialize (SOLVER &s,
const MATRIX &m,
const PRECONDITION &p)
{
solver = &s;
matrix = &m;
precondition = &p;
}
template<class SOLVER, class MATRIX, class PRECONDITION>
template<class VECTOR>
void
PreconditionLACSolver<SOLVER,MATRIX,PRECONDITION>::vmult (VECTOR &dst,
const VECTOR &src) const
{
Assert (solver !=0 && matrix != 0 && precondition != 0,
ExcNotInitialized());
solver->solve(*matrix, dst, src, *precondition);
}
//////////////////////////////////////////////////////////////////////
template<class MATRIX, class PRECOND, class VECTOR>
inline
PreconditionedMatrix<MATRIX, PRECOND, VECTOR>
::PreconditionedMatrix (const MATRIX & A,
const PRECOND &P,
VectorMemory<VECTOR> & mem):
A(A), P(P), mem(mem)
{}
template<class MATRIX, class PRECOND, class VECTOR>
inline void
PreconditionedMatrix<MATRIX, PRECOND, VECTOR>
::vmult (VECTOR &dst,
const VECTOR &src) const
{
VECTOR* h = mem.alloc();
h->reinit(src);
A.vmult(*h, src);
P.vmult(dst, *h);
mem.free(h);
}
template<class MATRIX, class PRECOND, class VECTOR>
inline void
PreconditionedMatrix<MATRIX, PRECOND, VECTOR>
::Tvmult (VECTOR &dst,
const VECTOR &src) const
{
VECTOR* h = mem.alloc();
h->reinit(src);
A.Tvmult(*h, src);
P.Tvmult(dst, *h);
mem.free(h);
}
template<class MATRIX, class PRECOND, class VECTOR>
inline double
PreconditionedMatrix<MATRIX, PRECOND, VECTOR>
::residual (VECTOR &dst,
const VECTOR &src,
const VECTOR &rhs) const
{
VECTOR* h = mem.alloc();
h->reinit(src);
A.vmult(*h, src);
P.vmult(dst, *h);
mem.free(h);
dst.sadd(-1.,1.,rhs);
return dst.l2_norm ();
}
//---------------------------------------------------------------------------
template <class MATRIX, class VECTOR>
inline
PreconditionChebyshev<MATRIX,VECTOR>::AdditionalData::
AdditionalData (const unsigned int degree,
const double smoothing_range,
const bool nonzero_starting,
const unsigned int eig_cg_n_iterations,
const double eig_cg_residual)
:
degree (degree),
smoothing_range (smoothing_range),
nonzero_starting (nonzero_starting),
eig_cg_n_iterations (eig_cg_n_iterations),
eig_cg_residual (eig_cg_residual)
{}
template <class MATRIX, class VECTOR>
inline
PreconditionChebyshev<MATRIX,VECTOR>::PreconditionChebyshev ()
:
is_initialized (false)
{}
template <class MATRIX, class VECTOR>
inline
void
PreconditionChebyshev<MATRIX,VECTOR>::initialize (const MATRIX &matrix,
const AdditionalData &additional_data)
{
Assert (matrix.m() == matrix.n(), ExcMessage("Matrix not quadratic."));
matrix_ptr = &matrix;
data = additional_data;
if (data.matrix_diagonal_inverse.size() != matrix.m())
{
data.matrix_diagonal_inverse.reinit(matrix.m());
for (unsigned int i=0; i<matrix.m(); ++i)
data.matrix_diagonal_inverse(i) = 1./matrix.el(i,i);
}
update1.reinit (data.matrix_diagonal_inverse, true);
update2.reinit (data.matrix_diagonal_inverse, true);
// calculate largest eigenvalue using a
// hand-tuned CG iteration on the matrix
// weighted by its diagonal. we start
// with a vector that consists of ones
// only, weighted by the length.
//
// TODO: can we obtain this with the
// regular CG implementation? we would need
// to read the logfile in that case, which
// does not seem feasible.
double max_eigenvalue, min_eigenvalue;
{
double eigen_beta_alpha = 0;
std::vector<double> diagonal;
std::vector<double> offdiagonal;
VECTOR rhs, g;
rhs.reinit(data.matrix_diagonal_inverse, true);
rhs = 1./std::sqrt(static_cast<double>(matrix.m()));
g.reinit(data.matrix_diagonal_inverse, true);
unsigned int it=0;
double res,gh,alpha,beta;
g.equ(-1.,rhs);
res = g.l2_norm();
update2.equ (-1., g);
gh = res*res;
while (true)
{
it++;
matrix.vmult (update1, update2);
update1.scale (data.matrix_diagonal_inverse);
alpha = update2 * update1;
alpha = gh/alpha;
g.add (alpha, update1);
res = g.l2_norm();
if (it > data.eig_cg_n_iterations || res < data.eig_cg_residual)
break;
beta = gh;
gh = res*res;
beta = gh/beta;
update2.sadd (beta, -1., g);
diagonal.push_back (1./alpha + eigen_beta_alpha);
eigen_beta_alpha = beta/alpha;
offdiagonal.push_back(std::sqrt(beta)/alpha);
}
TridiagonalMatrix<double> T(diagonal.size(), true);
for (unsigned int i=0;i<diagonal.size();++i)
{
T(i,i) = diagonal[i];
if (i< diagonal.size()-1)
T(i,i+1) = offdiagonal[i];
}
T.compute_eigenvalues();
min_eigenvalue = T.eigenvalue(0);
max_eigenvalue = T.eigenvalue(T.n()-1);
}
// include a safety factor since the CG
// method will in general not be converged
const double beta = 1.2 * max_eigenvalue;
const double alpha = (data.smoothing_range > 0 ?
max_eigenvalue / data.smoothing_range :
max_eigenvalue / min_eigenvalue);
delta = (beta-alpha)*0.5;
theta = (beta+alpha)*0.5;
is_initialized = true;
}
template <class MATRIX, class VECTOR>
inline
void
PreconditionChebyshev<MATRIX,VECTOR>::vmult (VECTOR &dst,
const VECTOR &src) const
{
Assert (is_initialized, ExcMessage("Preconditioner not initialized"));
double rhok = delta / theta, sigma = theta / delta;
if (data.nonzero_starting && !dst.all_zero())
{
matrix_ptr->vmult (update1, dst);
update1 -= src;
update1 /= theta;
update1.scale (data.matrix_diagonal_inverse);
dst -= update1;
}
else
{
dst.equ (1./theta, src);
dst.scale (data.matrix_diagonal_inverse);
update1.equ(-1.,dst);
}
for (unsigned int k=0; k<data.degree; ++k)
{
matrix_ptr->vmult (update2, dst);
update2 -= src;
update2.scale (data.matrix_diagonal_inverse);
const double rhokp = 1./(2.*sigma-rhok);
const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
rhok = rhokp;
update1.sadd (factor1, factor2, update2);
dst -= update1;
}
}
template <class MATRIX, class VECTOR>
inline
void
PreconditionChebyshev<MATRIX,VECTOR>::Tvmult (VECTOR &dst,
const VECTOR &src) const
{
Assert (is_initialized, ExcMessage("Preconditioner not initialized"));
double rhok = delta / theta, sigma = theta / delta;
if (data.nonzero_starting && !dst.all_zero())
{
matrix_ptr->Tvmult (update1, dst);
update1 -= src;
update1 /= theta;
update1.scale (data.matrix_diagonal_inverse);
dst -= update1;
}
else
{
dst.equ (1./theta, src);
dst.scale (data.matrix_diagonal_inverse);
update1.equ(-1.,dst);
}
for (unsigned int k=0; k<data.degree-1; ++k)
{
matrix_ptr->Tvmult (update2, dst);
update2 -= src;
update2.scale (data.matrix_diagonal_inverse);
const double rhokp = 1./(2.*sigma-rhok);
const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
rhok = rhokp;
update1.sadd (factor1, factor2, update2);
dst -= update1;
}
}
template <class MATRIX, class VECTOR>
inline
void PreconditionChebyshev<MATRIX,VECTOR>::clear ()
{
is_initialized = false;
matrix_ptr = 0;
data.matrix_diagonal_inverse.reinit(0);
update1.reinit(0);
update2.reinit(0);
}
#endif // DOXYGEN
DEAL_II_NAMESPACE_CLOSE
#endif
|