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
|
<!DOCTYPE html>
<html lang="en" data-content_root="../../" >
<head>
<meta charset="utf-8" />
<meta name="viewport" content="width=device-width, initial-scale=1.0" /><meta name="viewport" content="width=device-width, initial-scale=1" />
<title>EPS: Eigenvalue Problem Solver — SLEPc 3.24.1 documentation</title>
<script data-cfasync="false">
document.documentElement.dataset.mode = localStorage.getItem("mode") || "";
document.documentElement.dataset.theme = localStorage.getItem("theme") || "";
</script>
<!--
this give us a css class that will be invisible only if js is disabled
-->
<noscript>
<style>
.pst-js-only { display: none !important; }
</style>
</noscript>
<!-- Loaded before other Sphinx assets -->
<link href="../../_static/styles/theme.css?digest=8878045cc6db502f8baf" rel="stylesheet" />
<link href="../../_static/styles/pydata-sphinx-theme.css?digest=8878045cc6db502f8baf" rel="stylesheet" />
<link rel="stylesheet" type="text/css" href="../../_static/pygments.css?v=8f2a1f02" />
<link rel="stylesheet" type="text/css" href="../../_static/copybutton.css?v=76b2166b" />
<link rel="stylesheet" type="text/css" href="../../_static/togglebutton.css?v=13237357" />
<link rel="stylesheet" type="text/css" href="../../_static/sphinx-design.min.css?v=95c83b7e" />
<link rel="stylesheet" type="text/css" href="../../_static/css/slepc.css?v=d285b177" />
<!-- So that users can add custom icons -->
<script src="../../_static/scripts/fontawesome.js?digest=8878045cc6db502f8baf"></script>
<!-- Pre-loaded scripts that we'll load fully later -->
<link rel="preload" as="script" href="../../_static/scripts/bootstrap.js?digest=8878045cc6db502f8baf" />
<link rel="preload" as="script" href="../../_static/scripts/pydata-sphinx-theme.js?digest=8878045cc6db502f8baf" />
<script src="../../_static/documentation_options.js?v=d1c46438"></script>
<script src="../../_static/doctools.js?v=9a2dae69"></script>
<script src="../../_static/sphinx_highlight.js?v=dc90522c"></script>
<script src="../../_static/clipboard.min.js?v=a7894cd8"></script>
<script src="../../_static/copybutton.js?v=a56c686a"></script>
<script>let toggleHintShow = 'Click to show';</script>
<script>let toggleHintHide = 'Click to hide';</script>
<script>let toggleOpenOnPrint = 'true';</script>
<script src="../../_static/togglebutton.js?v=4a39c7ea"></script>
<script src="../../_static/design-tabs.js?v=f930bc37"></script>
<script>var togglebuttonSelector = '.toggle, .admonition.dropdown';</script>
<script>var togglebuttonSelector = '.toggle, .admonition.dropdown';</script>
<script>window.MathJax = {"options": {"processHtmlClass": "tex2jax_process|mathjax_process|math|output_area"}}</script>
<script defer="defer" src="https://cdn.jsdelivr.net/npm/mathjax@3/es5/tex-mml-chtml.js"></script>
<script>DOCUMENTATION_OPTIONS.pagename = 'documentation/manual/eps';</script>
<link rel="icon" href="../../_static/favicon-slepc.ico"/>
<link rel="index" title="Index" href="../../genindex.html" />
<link rel="search" title="Search" href="../../search.html" />
<link rel="next" title="ST: Spectral Transformation" href="st.html" />
<link rel="prev" title="Getting Started" href="intro.html" />
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<meta name="docsearch:language" content="en"/>
<meta name="docsearch:version" content="3.24" />
</head>
<body data-bs-spy="scroll" data-bs-target=".bd-toc-nav" data-offset="180" data-bs-root-margin="0px 0px -60%" data-default-mode="">
<div id="pst-skip-link" class="skip-link d-print-none"><a href="#main-content">Skip to main content</a></div>
<div id="pst-scroll-pixel-helper"></div>
<button type="button" class="btn rounded-pill" id="pst-back-to-top">
<i class="fa-solid fa-arrow-up"></i>Back to top</button>
<dialog id="pst-search-dialog">
<form class="bd-search d-flex align-items-center"
action="../../search.html"
method="get">
<i class="fa-solid fa-magnifying-glass"></i>
<input type="search"
class="form-control"
name="q"
placeholder="Search the docs ..."
aria-label="Search the docs ..."
autocomplete="off"
autocorrect="off"
autocapitalize="off"
spellcheck="false"/>
<span class="search-button__kbd-shortcut"><kbd class="kbd-shortcut__modifier">Ctrl</kbd>+<kbd>K</kbd></span>
</form>
</dialog>
<div class="pst-async-banner-revealer d-none">
<aside id="bd-header-version-warning" class="d-none d-print-none" aria-label="Version warning"></aside>
</div>
<header class="bd-header navbar navbar-expand-lg bd-navbar d-print-none">
<div class="bd-header__inner bd-page-width">
<button class="pst-navbar-icon sidebar-toggle primary-toggle" aria-label="Site navigation">
<span class="fa-solid fa-bars"></span>
</button>
<div class="col-lg-3 navbar-header-items__start">
<div class="navbar-item">
<a class="navbar-brand logo" href="../../index.html">
<img src="../../_static/logo-slepc.gif" class="logo__image only-light" alt="SLEPc Home"/>
<img src="../../_static/logo-slepc.gif" class="logo__image only-dark pst-js-only" alt="SLEPc Home"/>
</a></div>
</div>
<div class="col-lg-9 navbar-header-items">
<div class="me-auto navbar-header-items__center">
<div class="navbar-item">
<nav>
<ul class="bd-navbar-elements navbar-nav">
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../about/index.html">
About
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../installation/index.html">
Installation
</a>
</li>
<li class="nav-item current active">
<a class="nav-link nav-internal" href="../index.html">
Documentation
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../manualpages/index.html">
C/Fortran API
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../slepc4py/index.html">
slepc4py API
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../material/index.html">
Material
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../contact/index.html">
Contact
</a>
</li>
</ul>
</nav></div>
</div>
<div class="navbar-header-items__end">
<div class="navbar-item navbar-persistent--container">
<button class="btn search-button-field search-button__button pst-js-only" title="Search" aria-label="Search" data-bs-placement="bottom" data-bs-toggle="tooltip">
<i class="fa-solid fa-magnifying-glass"></i>
<span class="search-button__default-text">Search</span>
<span class="search-button__kbd-shortcut"><kbd class="kbd-shortcut__modifier">Ctrl</kbd>+<kbd class="kbd-shortcut__modifier">K</kbd></span>
</button>
</div>
<div class="navbar-item">
<button class="btn btn-sm nav-link pst-navbar-icon theme-switch-button pst-js-only" aria-label="Color mode" data-bs-title="Color mode" data-bs-placement="bottom" data-bs-toggle="tooltip">
<i class="theme-switch fa-solid fa-sun fa-lg" data-mode="light" title="Light"></i>
<i class="theme-switch fa-solid fa-moon fa-lg" data-mode="dark" title="Dark"></i>
<i class="theme-switch fa-solid fa-circle-half-stroke fa-lg" data-mode="auto" title="System Settings"></i>
</button></div>
<div class="navbar-item"><ul class="navbar-icon-links"
aria-label="Icon Links">
<li class="nav-item">
<a href="https://gitlab.com/slepc/slepc" title="GitLab" class="nav-link pst-navbar-icon" rel="noopener" target="_blank" data-bs-toggle="tooltip" data-bs-placement="bottom"><i class="fa-brands fa-square-gitlab fa-lg" aria-hidden="true"></i>
<span class="sr-only">GitLab</span></a>
</li>
<li class="nav-item">
<a href="https://www.upv.es" title="UPV" class="nav-link pst-navbar-icon" rel="noopener" target="_blank" data-bs-toggle="tooltip" data-bs-placement="bottom"><img src="https://www.upv.es/favicon.ico" class="icon-link-image" alt="UPV"/></a>
</li>
<li class="nav-item">
<a href="https://slepc.upv.es/release/_static/rss/slepc-news.xml" title="Feed" class="nav-link pst-navbar-icon" rel="noopener" target="_blank" data-bs-toggle="tooltip" data-bs-placement="bottom"><i class="fa-solid fa-square-rss fa-lg" aria-hidden="true"></i>
<span class="sr-only">Feed</span></a>
</li>
</ul></div>
</div>
</div>
<div class="navbar-persistent--mobile">
<button class="btn search-button-field search-button__button pst-js-only" title="Search" aria-label="Search" data-bs-placement="bottom" data-bs-toggle="tooltip">
<i class="fa-solid fa-magnifying-glass"></i>
<span class="search-button__default-text">Search</span>
<span class="search-button__kbd-shortcut"><kbd class="kbd-shortcut__modifier">Ctrl</kbd>+<kbd class="kbd-shortcut__modifier">K</kbd></span>
</button>
</div>
<button class="pst-navbar-icon sidebar-toggle secondary-toggle" aria-label="On this page">
<span class="fa-solid fa-outdent"></span>
</button>
</div>
</header>
<div class="bd-container">
<div class="bd-container__inner bd-page-width">
<dialog id="pst-primary-sidebar-modal"></dialog>
<div id="pst-primary-sidebar" class="bd-sidebar-primary bd-sidebar">
<div class="sidebar-header-items sidebar-primary__section">
<div class="sidebar-header-items__center">
<div class="navbar-item">
<nav>
<ul class="bd-navbar-elements navbar-nav">
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../about/index.html">
About
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../installation/index.html">
Installation
</a>
</li>
<li class="nav-item current active">
<a class="nav-link nav-internal" href="../index.html">
Documentation
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../manualpages/index.html">
C/Fortran API
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../slepc4py/index.html">
slepc4py API
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../material/index.html">
Material
</a>
</li>
<li class="nav-item ">
<a class="nav-link nav-internal" href="../../contact/index.html">
Contact
</a>
</li>
</ul>
</nav></div>
</div>
<div class="sidebar-header-items__end">
<div class="navbar-item">
<button class="btn btn-sm nav-link pst-navbar-icon theme-switch-button pst-js-only" aria-label="Color mode" data-bs-title="Color mode" data-bs-placement="bottom" data-bs-toggle="tooltip">
<i class="theme-switch fa-solid fa-sun fa-lg" data-mode="light" title="Light"></i>
<i class="theme-switch fa-solid fa-moon fa-lg" data-mode="dark" title="Dark"></i>
<i class="theme-switch fa-solid fa-circle-half-stroke fa-lg" data-mode="auto" title="System Settings"></i>
</button></div>
<div class="navbar-item"><ul class="navbar-icon-links"
aria-label="Icon Links">
<li class="nav-item">
<a href="https://gitlab.com/slepc/slepc" title="GitLab" class="nav-link pst-navbar-icon" rel="noopener" target="_blank" data-bs-toggle="tooltip" data-bs-placement="bottom"><i class="fa-brands fa-square-gitlab fa-lg" aria-hidden="true"></i>
<span class="sr-only">GitLab</span></a>
</li>
<li class="nav-item">
<a href="https://www.upv.es" title="UPV" class="nav-link pst-navbar-icon" rel="noopener" target="_blank" data-bs-toggle="tooltip" data-bs-placement="bottom"><img src="https://www.upv.es/favicon.ico" class="icon-link-image" alt="UPV"/></a>
</li>
<li class="nav-item">
<a href="https://slepc.upv.es/release/_static/rss/slepc-news.xml" title="Feed" class="nav-link pst-navbar-icon" rel="noopener" target="_blank" data-bs-toggle="tooltip" data-bs-placement="bottom"><i class="fa-solid fa-square-rss fa-lg" aria-hidden="true"></i>
<span class="sr-only">Feed</span></a>
</li>
</ul></div>
</div>
</div>
<div class="sidebar-primary-items__start sidebar-primary__section">
<div class="sidebar-primary-item">
<nav class="bd-docs-nav bd-links"
aria-label="Section Navigation">
<p class="bd-links__title" role="heading" aria-level="1">Section Navigation</p>
<div class="bd-toc-item navbar-nav"><ul class="current nav bd-sidenav">
<li class="toctree-l1 current active has-children"><a class="reference internal" href="index.html">SLEPc Users Manual</a><details open="open"><summary><span class="toctree-toggle" role="presentation"><i class="fa-solid fa-chevron-down"></i></span></summary><ul class="current">
<li class="toctree-l2"><a class="reference internal" href="intro.html">Getting Started</a></li>
<li class="toctree-l2 current active"><a class="current reference internal" href="#">EPS: Eigenvalue Problem Solver</a></li>
<li class="toctree-l2"><a class="reference internal" href="st.html">ST: Spectral Transformation</a></li>
<li class="toctree-l2"><a class="reference internal" href="svd.html">SVD: Singular Value Decomposition</a></li>
<li class="toctree-l2"><a class="reference internal" href="pep.html">PEP: Polynomial Eigenvalue Problems</a></li>
<li class="toctree-l2"><a class="reference internal" href="nep.html">NEP: Nonlinear Eigenvalue Problems</a></li>
<li class="toctree-l2"><a class="reference internal" href="mfn.html">MFN: Matrix Function</a></li>
<li class="toctree-l2"><a class="reference internal" href="lme.html">LME: Linear Matrix Equation</a></li>
<li class="toctree-l2"><a class="reference internal" href="aux.html">Auxiliary Classes</a></li>
<li class="toctree-l2"><a class="reference internal" href="extra.html">Additional Information</a></li>
</ul>
</details></li>
<li class="toctree-l1 has-children"><a class="reference internal" href="../hands-on/index.html">Hands-on exercises</a><details><summary><span class="toctree-toggle" role="presentation"><i class="fa-solid fa-chevron-down"></i></span></summary><ul>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on0.html">Exercise 0: Hello World</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on1.html">Exercise 1: Standard Symmetric Eigenvalue Problem</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on2.html">Exercise 2: Standard Non-Symmetric Eigenvalue Problem</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on3.html">Exercise 3: Generalized Eigenvalue Problem Stored in a File</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on4.html">Exercise 4: Singular Value Decomposition</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on5.html">Exercise 5: Problem without Explicit Matrix Storage</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on6.html">Exercise 6: Parallel Execution</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on7.html">Exercise 7: Use of Deflation Subspaces</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on8.html">Exercise 8: Quadratic Eigenvalue Problem</a></li>
</ul>
</details></li>
<li class="toctree-l1"><a class="reference internal" href="../faq.html">Frequently Asked Questions</a></li>
<li class="toctree-l1"><a class="reference internal" href="../presentations.html">Presentations</a></li>
<li class="toctree-l1"><a class="reference internal" href="../license.html">License</a></li>
</ul>
</div>
</nav></div>
</div>
<div class="sidebar-primary-items__end sidebar-primary__section">
<div class="sidebar-primary-item">
<div id="ethical-ad-placement"
class="flat"
data-ea-publisher="readthedocs"
data-ea-type="readthedocs-sidebar"
data-ea-manual="true">
</div></div>
</div>
</div>
<main id="main-content" class="bd-main" role="main">
<div class="bd-content">
<div class="bd-article-container">
<div class="bd-header-article d-print-none">
<div class="header-article-items header-article__inner">
<div class="header-article-items__start">
<div class="header-article-item">
<nav aria-label="Breadcrumb" class="d-print-none">
<ul class="bd-breadcrumbs">
<li class="breadcrumb-item breadcrumb-home">
<a href="../../index.html" class="nav-link" aria-label="Home">
<i class="fa-solid fa-home"></i>
</a>
</li>
<li class="breadcrumb-item"><a href="../index.html" class="nav-link">Documentation</a></li>
<li class="breadcrumb-item"><a href="index.html" class="nav-link">SLEPc Users Manual</a></li>
<li class="breadcrumb-item active" aria-current="page"><span class="ellipsis">EPS: Eigenvalue Problem Solver</span></li>
</ul>
</nav>
</div>
</div>
</div>
</div>
<div id="searchbox"></div>
<article class="bd-article">
<section class="tex2jax_ignore mathjax_ignore" id="eps-eigenvalue-problem-solver">
<span id="ch-eps"></span><h1>EPS: Eigenvalue Problem Solver<a class="headerlink" href="#eps-eigenvalue-problem-solver" title="Link to this heading">#</a></h1>
<p>The Eigenvalue Problem Solver (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>) is the main object provided by SLEPc. It is used to specify a linear eigenvalue problem, either in standard or generalized form, and provides uniform and efficient access to all of the linear eigensolvers included in the package. Conceptually, the level of abstraction occupied by <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> is similar to other solvers in PETSc such as <a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSP/" title="(in PETSc v3.24)"><span class="xref std std-doc">KSP</span></a> for solving systems of linear equations.</p>
<section id="sec-eig">
<h2>Eigenvalue Problems<a class="headerlink" href="#sec-eig" title="Link to this heading">#</a></h2>
<p>In this section, we briefly present some basic concepts about eigenvalue problems as well as general techniques used to solve them. The description is not intended to be exhaustive. The objective is simply to define terms that will be referred to throughout the rest of the manual. Readers who are familiar with the terminology and the solution approach can skip this section. For a more comprehensive description, we refer the reader to monographs such as <span id="id1">[<a class="reference internal" href="svd.html#id15" title="Z. Bai, J. Demmel, J. Dongarra, A. Ruhe, and H. van der Vorst, editors. Templates for the Solution of Algebraic Eigenvalue Problems: A Practical Guide. Society for Industrial and Applied Mathematics, Philadelphia, PA, 2000.">Bai <em>et al.</em>, 2000</a>, <a class="reference internal" href="../../manualpages/EPS/EPS_POWER_SHIFT_WILKINSON.html#id11" title="B. N. Parlett. The Symmetric Eigenvalue Problem. Prentice-Hall, Englewood Cliffs, NJ, 1980. Reissued with revisions by SIAM, Philadelphia, 1998.">Parlett, 1980</a>, <a class="reference internal" href="#id23" title="Y. Saad. Numerical Methods for Large Eigenvalue Problems: Theory and Algorithms. John Wiley and Sons, New York, 1992.">Saad, 1992</a>, <a class="reference internal" href="#id24" title="G. W. Stewart. Matrix Algorithms. Volume II: Eigensystems. Society for Industrial and Applied Mathematics, Philadelphia, PA, USA, 2001.">Stewart, 2001</a>]</span>. A historical perspective of the topic can be found in <span id="id2">[<a class="reference internal" href="#id17" title="G. H. Golub and H. A. van der Vorst. Eigenvalue computation in the 20th century. J. Comput. Appl. Math., 123(1-2):35–65, 2000. doi:10.1016/S0377-0427(00)00413-1.">Golub and van der Vorst, 2000</a>]</span>. See also the SLEPc <a class="reference internal" href="../index.html#str"><span class="std std-ref">technical reports</span></a>.</p>
<p>In the standard formulation, the linear eigenvalue problem consists in determining <span class="math notranslate nohighlight">\(\lambda\in\mathbb{C}\)</span> for which the equation</p>
<div class="math notranslate nohighlight" id="equation-eq-eigstd">
<span class="eqno">(1)<a class="headerlink" href="#equation-eq-eigstd" title="Link to this equation">#</a></span>\[ Ax=\lambda x\]</div>
<p>has nontrivial solution, where <span class="math notranslate nohighlight">\(A\in\mathbb{C}^{n\times n}\)</span> and <span class="math notranslate nohighlight">\(x\in\mathbb{C}^n\)</span>. The scalar <span class="math notranslate nohighlight">\(\lambda\)</span> and the vector <span class="math notranslate nohighlight">\(x\neq 0\)</span> are called eigenvalue and (right) eigenvector, respectively. Note that they can be complex even when the matrix is real. If <span class="math notranslate nohighlight">\(\lambda\)</span> is an eigenvalue of <span class="math notranslate nohighlight">\(A\)</span> then <span class="math notranslate nohighlight">\(\bar{\lambda}\)</span> is an eigenvalue of its conjugate transpose, <span class="math notranslate nohighlight">\(A^*\)</span>, or equivalently</p>
<div class="math notranslate nohighlight" id="equation-eq-eigstdleft">
<span class="eqno">(2)<a class="headerlink" href="#equation-eq-eigstdleft" title="Link to this equation">#</a></span>\[y^*\!A=\lambda\, y^*,\]</div>
<p>where <span class="math notranslate nohighlight">\(y\neq 0\)</span> is called the left eigenvector.</p>
<p>In many applications, the problem is formulated as</p>
<div class="math notranslate nohighlight" id="equation-eq-eiggen">
<span class="eqno">(3)<a class="headerlink" href="#equation-eq-eiggen" title="Link to this equation">#</a></span>\[Ax=\lambda Bx,\]</div>
<p>where <span class="math notranslate nohighlight">\(B\in\mathbb{C}^{n\times n}\)</span>, which is known as the generalized eigenvalue problem. Often, this problem is solved by reformulating it in standard form, for example <span class="math notranslate nohighlight">\(B^{-1}Ax=\lambda x\)</span> if <span class="math notranslate nohighlight">\(B\)</span> is non-singular.</p>
<p>SLEPc focuses on the solution of problems in which the matrices are large and sparse. Hence, only methods that preserve sparsity are considered. These methods obtain the solution from the information generated by the application of the operator to various vectors (the operator is a simple function of matrices <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span>), that is, matrices are only used in matrix-vector products. This not only maintains sparsity but allows the solution of problems in which matrices are not available explicitly.</p>
<p>In practical analyses, from the <span class="math notranslate nohighlight">\(n\)</span> possible solutions, typically only a few eigenpairs <span class="math notranslate nohighlight">\((\lambda,x)\)</span> are considered relevant, either in the extremities of the spectrum, in an interval, or in a region of the complex plane. Depending on the application, either eigenvalues or eigenvectors (or both) are required. In some cases, left eigenvectors are also of interest.</p>
<section id="projection-methods">
<h3>Projection Methods<a class="headerlink" href="#projection-methods" title="Link to this heading">#</a></h3>
<p>Most eigensolvers provided by SLEPc perform a Rayleigh-Ritz projection for extracting the spectral approximations, that is, they project the problem onto a low-dimensional subspace that is built appropriately. Suppose that an orthogonal basis of this subspace is given by <span class="math notranslate nohighlight">\(V_j=[v_1,v_2,\ldots,v_j]\)</span>. If the solutions of the projected (reduced) problem <span class="math notranslate nohighlight">\(M_js=\theta s\)</span> (i.e., <span class="math notranslate nohighlight">\(V_j^*AV_j=M_j\)</span>) are assumed to be <span class="math notranslate nohighlight">\((\theta_i,s_i)\)</span>, <span class="math notranslate nohighlight">\(i=1,2,\ldots,j\)</span>, then the approximate eigenpairs <span class="math notranslate nohighlight">\((\tilde{\lambda}_i,\tilde{x}_i)\)</span> of the original problem (Ritz value and Ritz vector) are obtained as</p>
<div class="math notranslate nohighlight" id="equation-eq-ritz-pair">
<span class="eqno">(4)<a class="headerlink" href="#equation-eq-ritz-pair" title="Link to this equation">#</a></span>\[\tilde{\lambda}_i=\theta_i,\qquad
\tilde{x}_i=V_js_i.\]</div>
<p>Starting from this general idea, eigensolvers differ from each other in which subspace is used, how it is built and other technicalities aimed at improving convergence, reducing storage requirements, etc.</p>
<p>The subspace</p>
<div class="math notranslate nohighlight" id="equation-eq-krylov">
<span class="eqno">(5)<a class="headerlink" href="#equation-eq-krylov" title="Link to this equation">#</a></span>\[\mathcal{K}_m(A,v)\equiv\mathrm{span}\left\{v,Av,A^2v,\ldots,A^{m-1}v\right\},\]</div>
<p>is called the <span class="math notranslate nohighlight">\(m\)</span>-th Krylov subspace corresponding to <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(v\)</span>. Methods that use subspaces of this kind to carry out the projection are called Krylov methods. One example of such methods is the Arnoldi algorithm: starting with <span class="math notranslate nohighlight">\(v_1\)</span>, <span class="math notranslate nohighlight">\(\|v_1\|_2=1\)</span>, the Arnoldi basis generation process can be expressed by the recurrence</p>
<div class="math notranslate nohighlight" id="equation-eq-arnoldi">
<span class="eqno">(6)<a class="headerlink" href="#equation-eq-arnoldi" title="Link to this equation">#</a></span>\[h_{j+1,j}v_{j+1}=w_j:=Av_j-\sum_{i=1}^jh_{i,j}v_i,\]</div>
<p>where <span class="math notranslate nohighlight">\(h_{i,j}\)</span> are the scalar coefficients obtained in the Gram-Schmidt orthogonalization of <span class="math notranslate nohighlight">\(Av_j\)</span> with respect to <span class="math notranslate nohighlight">\(v_i\)</span>, <span class="math notranslate nohighlight">\(i=1,2,\ldots,j\)</span>, and <span class="math notranslate nohighlight">\(h_{j+1,j}=\|w_j\|_2\)</span>. Then, the columns of <span class="math notranslate nohighlight">\(V_j\)</span> span the Krylov subspace <span class="math notranslate nohighlight">\(\mathcal{K}_j(A,v_1)\)</span> and <span class="math notranslate nohighlight">\(Ax=\lambda x\)</span> is projected onto <span class="math notranslate nohighlight">\(H_js=\theta s\)</span>, where <span class="math notranslate nohighlight">\(H_j\)</span> is an upper Hessenberg matrix with elements <span class="math notranslate nohighlight">\(h_{i,j}\)</span>, which are 0 for <span class="math notranslate nohighlight">\(i\geq j+2\)</span>. The related Lanczos algorithms obtain a projected matrix that is tridiagonal.</p>
<p>A generalization to the above methods are the block Krylov strategies, in which the starting vector <span class="math notranslate nohighlight">\(v_1\)</span> is replaced by a full rank <span class="math notranslate nohighlight">\(n\times p\)</span> matrix <span class="math notranslate nohighlight">\(V_1\)</span>, which allows for better convergence properties when there are multiple eigenvalues and can provide better data management on modern computer architectures. Block tridiagonal and block Hessenberg matrices are then obtained as projections.</p>
<p>It is generally assumed (and observed) that the Lanczos and Arnoldi algorithms find solutions at the extremities of the spectrum. Their convergence pattern, however, is strongly related to the eigenvalue distribution. Slow convergence may be experienced in the presence of tightly clustered eigenvalues. The maximum allowable <span class="math notranslate nohighlight">\(j\)</span> may be reached without having achieved convergence for all desired solutions. Then, restarting is usually a useful technique and different strategies exist for that purpose. However, convergence can still be very slow and acceleration strategies must be applied. Usually, these techniques consist in computing eigenpairs of a transformed operator and then recovering the solution of the original problem. The aim of these transformations is twofold. On one hand, they make it possible to obtain eigenvalues other than those lying in the periphery of the spectrum. On the other hand, the separation of the eigenvalues of interest is improved in the transformed spectrum thus leading to faster convergence. The most commonly used spectral transformation is called shift-and-invert, which works with operator <span class="math notranslate nohighlight">\((A-\sigma I)^{-1}\)</span>. It allows the computation of eigenvalues closest to <span class="math notranslate nohighlight">\(\sigma\)</span> with very good separation properties. When using this approach, a linear system of equations, <span class="math notranslate nohighlight">\((A-\sigma I)y=x\)</span>, must be solved in each iteration of the eigenvalue process.</p>
</section>
<section id="preconditioned-eigensolvers">
<h3>Preconditioned Eigensolvers<a class="headerlink" href="#preconditioned-eigensolvers" title="Link to this heading">#</a></h3>
<p>In many applications, Krylov eigensolvers perform very well because Krylov subspaces are optimal in a certain theoretical sense. However, these methods may not be appropriate in some situations such as the computation of interior eigenvalues. The spectral transformation mentioned above may not be a viable solution or it may be too costly. For these reasons, other types of eigensolvers such as Davidson and Jacobi-Davidson rely on a different way of expanding the subspace. Instead of satisfying the Krylov relation, these methods compute the new basis vector by the so-called correction equation. The resulting subspace may be richer in the direction of the desired eigenvectors. These solvers may be competitive especially for computing interior eigenvalues. From a practical point of view, the correction equation may be seen as a cheap replacement for the shift-and-invert system of equations, <span class="math notranslate nohighlight">\((A-\sigma I)y=x\)</span>. By cheap we mean that it may be solved inaccurately without compromising robustness, via a preconditioned iterative linear solver. For this reason, these are known as <em>preconditioned</em> eigensolvers.</p>
</section>
<section id="related-problems">
<h3>Related Problems<a class="headerlink" href="#related-problems" title="Link to this heading">#</a></h3>
<p>In many applications such as the analysis of damped vibrating systems the problem to be solved is a <em>polynomial eigenvalue problem</em> (PEP), or more generally a <em>nonlinear eigenvalue problem</em> (NEP). For these, the reader is referred to chapters <a class="reference internal" href="pep.html#ch-pep"><span class="std std-ref">PEP: Polynomial Eigenvalue Problems</span></a> and <a class="reference internal" href="nep.html#ch-nep"><span class="std std-ref">NEP: Nonlinear Eigenvalue Problems</span></a>. Another linear algebra problem that is very closely related to the eigenvalue problem is the <em>singular value decomposition</em> (SVD), see chapter <a class="reference internal" href="svd.html#ch-svd"><span class="std std-ref">SVD: Singular Value Decomposition</span></a>.</p>
</section>
</section>
<section id="basic-usage">
<h2>Basic Usage<a class="headerlink" href="#basic-usage" title="Link to this heading">#</a></h2>
<p>The <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> module in SLEPc is used in a similar way as PETSc modules such as <a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSP/" title="(in PETSc v3.24)"><span class="xref std std-doc">KSP</span></a>. All the information related to an eigenvalue problem is handled via a context variable. The usual object management functions are available (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSCreate.html">EPSCreate</a>()</span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSDestroy.html">EPSDestroy</a>()</span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSView.html">EPSView</a>()</span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetFromOptions.html">EPSSetFromOptions</a>()</span></code>). In addition, the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object provides functions for setting several parameters such as the number of eigenvalues to compute, the dimension of the subspace, the portion of the spectrum of interest, the requested tolerance or the maximum number of iterations allowed.</p>
<p>The solution of the problem is obtained in several steps. First of all, the matrices associated with the eigenproblem are specified via <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetOperators.html">EPSSetOperators</a>()</span></code> and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetProblemType.html">EPSSetProblemType</a>()</span></code> is used to specify the type of problem. Then, a call to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a>()</span></code> is done that invokes the subroutine for the selected eigensolver. <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetConverged.html">EPSGetConverged</a>()</span></code> can be used afterwards to determine how many of the requested eigenpairs have converged to working accuracy. <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetEigenpair.html">EPSGetEigenpair</a>()</span></code> is finally used to retrieve the eigenvalues and eigenvectors.</p>
<p>As an illustration of basic usage, see listing <a class="reference internal" href="#fig-ex-eps"><span class="std std-ref">Example code for basic solution with EPS</span></a>, that implements the solution of a simple standard eigenvalue problem. Code for setting up the matrix <code class="docutils notranslate"><span class="pre">A</span></code> is not shown and error-checking code is omitted.</p>
<div class="literal-block-wrapper docutils container" id="fig-ex-eps">
<div class="code-block-caption"><span class="caption-text">Example code for basic solution with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code></span><a class="headerlink" href="#fig-ex-eps" title="Link to this code">#</a></div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="linenos"> 1</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">;</span><span class="w"> </span><span class="cm">/* eigensolver context */</span>
<span class="linenos"> 2</span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/Mat/">Mat</a></span><span class="w"> </span><span class="n">A</span><span class="p">;</span><span class="w"> </span><span class="cm">/* matrix of Ax=kx */</span>
<span class="linenos"> 3</span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">xr</span><span class="p">,</span><span class="w"> </span><span class="n">xi</span><span class="p">;</span><span class="w"> </span><span class="cm">/* eigenvector, x */</span>
<span class="linenos"> 4</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscScalar/">PetscScalar</a></span><span class="w"> </span><span class="n">kr</span><span class="p">,</span><span class="w"> </span><span class="n">ki</span><span class="p">;</span><span class="w"> </span><span class="cm">/* eigenvalue, k */</span>
<span class="linenos"> 5</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">j</span><span class="p">,</span><span class="w"> </span><span class="n">nconv</span><span class="p">;</span>
<span class="linenos"> 6</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="n">error</span><span class="p">;</span>
<span class="linenos"> 7</span>
<span class="linenos"> 8</span><span class="n"><a href="../../manualpages/EPS/EPSCreate.html">EPSCreate</a></span><span class="p">(</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PETSC_COMM_WORLD/">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="w"> </span><span class="o">&</span><span class="n">eps</span><span class="p">);</span>
<span class="linenos"> 9</span><span class="n"><a href="../../manualpages/EPS/EPSSetOperators.html">EPSSetOperators</a></span><span class="p">(</span><span class="n">eps</span><span class="p">,</span><span class="w"> </span><span class="n">A</span><span class="p">,</span><span class="w"> </span><span class="nb">NULL</span><span class="p">);</span>
<span class="linenos">10</span><span class="n"><a href="../../manualpages/EPS/EPSSetProblemType.html">EPSSetProblemType</a></span><span class="p">(</span><span class="n">eps</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../../manualpages/EPS/EPS_NHEP.html">EPS_NHEP</a></span><span class="p">);</span>
<span class="linenos">11</span><span class="n"><a href="../../manualpages/EPS/EPSSetFromOptions.html">EPSSetFromOptions</a></span><span class="p">(</span><span class="n">eps</span><span class="p">);</span>
<span class="linenos">12</span><span class="n"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a></span><span class="p">(</span><span class="n">eps</span><span class="p">);</span>
<span class="linenos">13</span><span class="n"><a href="../../manualpages/EPS/EPSGetConverged.html">EPSGetConverged</a></span><span class="p">(</span><span class="n">eps</span><span class="p">,</span><span class="w"> </span><span class="o">&</span><span class="n">nconv</span><span class="p">);</span>
<span class="linenos">14</span><span class="k">for</span><span class="w"> </span><span class="p">(</span><span class="n">j</span><span class="o">=</span><span class="mi">0</span><span class="p">;</span><span class="n">j</span><span class="o"><</span><span class="n">nconv</span><span class="p">;</span><span class="n">j</span><span class="o">++</span><span class="p">)</span><span class="w"> </span><span class="p">{</span>
<span class="linenos">15</span><span class="w"> </span><span class="n"><a href="../../manualpages/EPS/EPSGetEigenpair.html">EPSGetEigenpair</a></span><span class="p">(</span><span class="n">eps</span><span class="p">,</span><span class="w"> </span><span class="n">j</span><span class="p">,</span><span class="w"> </span><span class="o">&</span><span class="n">kr</span><span class="p">,</span><span class="w"> </span><span class="o">&</span><span class="n">ki</span><span class="p">,</span><span class="w"> </span><span class="n">xr</span><span class="p">,</span><span class="w"> </span><span class="n">xi</span><span class="p">);</span>
<span class="linenos">16</span><span class="w"> </span><span class="n"><a href="../../manualpages/EPS/EPSComputeError.html">EPSComputeError</a></span><span class="p">(</span><span class="n">eps</span><span class="p">,</span><span class="w"> </span><span class="n">j</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../../manualpages/EPS/EPSErrorType.html">EPS_ERROR_RELATIVE</a></span><span class="p">,</span><span class="w"> </span><span class="o">&</span><span class="n">error</span><span class="p">);</span>
<span class="linenos">17</span><span class="p">}</span>
<span class="linenos">18</span><span class="n"><a href="../../manualpages/EPS/EPSDestroy.html">EPSDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">eps</span><span class="p">);</span>
</pre></div>
</div>
</div>
<p>All the operations of the program are done over a single <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object. This solver context is created in line 8 with the function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSCreate.html">EPSCreate</a></span><span class="p">(</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/MPI_Comm/">MPI_Comm</a></span><span class="w"> </span><span class="n">comm</span><span class="p">,</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="o">*</span><span class="n">eps</span><span class="p">);</span>
</pre></div>
</div>
<p>Here <code class="docutils notranslate"><span class="pre">comm</span></code> is the MPI communicator, and <code class="docutils notranslate"><span class="pre">eps</span></code> is the newly formed solver context. The communicator indicates which processes are involved in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object. Most of the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> operations are collective, meaning that all the processes collaborate to perform the operation in parallel.</p>
<p>Before actually solving an eigenvalue problem with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>, the user must specify the matrices associated with the problem, as in line 9, with the following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetOperators.html">EPSSetOperators</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/Mat/">Mat</a></span><span class="w"> </span><span class="n">A</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/Mat/">Mat</a></span><span class="w"> </span><span class="n">B</span><span class="p">);</span>
</pre></div>
</div>
<p>The example specifies a standard eigenproblem. In the case of a generalized problem, it would be necessary also to provide matrix <span class="math notranslate nohighlight">\(B\)</span> as the third argument to the call. The matrices specified in this call can be in any PETSc format. In particular, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> allows the user to solve matrix-free problems by specifying matrices created via <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatCreateShell/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatCreateShell</span></a>(). A more detailed discussion of this issue is given in section <a class="reference internal" href="extra.html#sec-supported"><span class="std std-ref">Supported Matrix Types</span></a>.</p>
<p>After setting the problem matrices, the problem type is set with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetProblemType.html">EPSSetProblemType</a>()</span></code>. This is not strictly necessary since if this step is skipped then the problem type is assumed to be non-symmetric. More details are given in section <a class="reference internal" href="#sec-defprob"><span class="std std-ref">Defining the Problem</span></a>. At this point, the value of the different options could optionally be set by means of a function call such as <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetTolerances.html">EPSSetTolerances</a>()</span></code> (explained later in this chapter). After this, a call to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetFromOptions.html">EPSSetFromOptions</a>()</span></code> should be made as in line 11 of <a class="reference internal" href="#fig-ex-eps"><span class="std std-ref">Example code for basic solution with EPS</span></a>:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetFromOptions.html">EPSSetFromOptions</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">);</span>
</pre></div>
</div>
<p>The effect of this call is that options specified at runtime in the command line are passed to the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object appropriately. In this way, the user can easily experiment with different combinations of options without having to recompile. All the available options as well as the associated function calls are described later in this chapter.</p>
<p>Line 12 launches the solution algorithm, simply with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">);</span>
</pre></div>
</div>
<p>The subroutine that is actually invoked depends on which solver has been selected by the user.</p>
<p>After the call to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a>()</span></code> has finished, all the data associated with the solution of the eigenproblem are kept internally. This information can be retrieved with different function calls, as in lines 13 to 17 of <a class="reference internal" href="#fig-ex-eps"><span class="std std-ref">Example code for basic solution with EPS</span></a>. This part is described in detail in section <a class="reference internal" href="#sec-retrsol"><span class="std std-ref">Retrieving the Solution</span></a>.</p>
<p>Once the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> context is no longer needed, it should be destroyed with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSDestroy.html">EPSDestroy</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="o">*</span><span class="n">eps</span><span class="p">);</span>
</pre></div>
</div>
<p>The above procedure is sufficient for general use of the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> package. As in the case of the <a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSP/" title="(in PETSc v3.24)"><span class="xref std std-doc">KSP</span></a> solver, the user can optionally explicitly call <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetUp.html">EPSSetUp</a>()</span></code> before calling <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a>()</span></code> to perform any setup required for the eigensolver.</p>
<p>Internally, the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object works with an <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> object (spectral transformation, described in chapter <a class="reference internal" href="st.html#ch-st"><span class="std std-ref">ST: Spectral Transformation</span></a>). To allow application programmers to set any of the spectral transformation options directly within the code, the following routine is provided to extract the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> context:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSGetST.html">EPSGetST</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="o">*</span><span class="n">st</span><span class="p">);</span>
</pre></div>
</div>
<p>With the function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSView.html">EPSView</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Viewer/PetscViewer/">PetscViewer</a></span><span class="w"> </span><span class="n">viewer</span><span class="p">);</span>
</pre></div>
</div>
<p>it is possible to examine the actual values of the different settings of the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object, including also those related to the associated <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> object. This is useful for making sure that the solver is using the settings that the user wants.</p>
</section>
<section id="sec-defprob">
<h2>Defining the Problem<a class="headerlink" href="#sec-defprob" title="Link to this heading">#</a></h2>
<p>SLEPc is able to cope with different kinds of problems. Currently supported problem types are listed in table <a class="reference internal" href="#tab-ptype"><span class="std std-ref">Problem types considered in EPS</span></a>. An eigenproblem is generalized (<span class="math notranslate nohighlight">\(Ax=\lambda Bx\)</span>) if the user has specified two matrices (see <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetOperators.html">EPSSetOperators</a>()</span></code> discussed above), otherwise it is standard (<span class="math notranslate nohighlight">\(Ax=\lambda x\)</span>). A standard eigenproblem is Hermitian if matrix <span class="math notranslate nohighlight">\(A\)</span> is Hermitian (i.e., <span class="math notranslate nohighlight">\(A=A^*\)</span>) or, equivalently in the case of real matrices, if matrix <span class="math notranslate nohighlight">\(A\)</span> is symmetric (i.e., <span class="math notranslate nohighlight">\(A=A^T\)</span>). A generalized eigenproblem is Hermitian if matrix <span class="math notranslate nohighlight">\(A\)</span> is Hermitian (symmetric) and <span class="math notranslate nohighlight">\(B\)</span> is Hermitian (symmetric) and positive (semi-)definite. If <span class="math notranslate nohighlight">\(B\)</span> is not positive (semi-)definite then the problem cannot be considered Hermitian but symmetry can still be exploited to some extent in some solvers (problem type <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHIEP.html">EPS_GHIEP</a></span></code>). A special case of generalized non-Hermitian problem is when <span class="math notranslate nohighlight">\(A\)</span> is non-Hermitian but <span class="math notranslate nohighlight">\(B\)</span> is Hermitian and positive (semi-)definite, see section <a class="reference internal" href="st.html#sec-symm"><span class="std std-ref">Preserving the Symmetry in Generalized Eigenproblems</span></a> and <a class="reference internal" href="st.html#sec-purif"><span class="std std-ref">Purification of Eigenvectors</span></a> for discussion. The remaining entries in table <a class="reference internal" href="#tab-ptype"><span class="std std-ref">Problem types considered in EPS</span></a> correspond to structured eigenvalue problems, which are discussed in section <a class="reference internal" href="#sec-structured"><span class="std std-ref">Structured Eigenvalue Problems</span></a>.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-ptype">
<caption><span class="caption-text">Problem types considered in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code></span><a class="headerlink" href="#tab-ptype" title="Link to this table">#</a></caption>
<thead>
<tr class="row-odd"><th class="head"><p>Problem Type</p></th>
<th class="head"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSProblemType.html">EPSProblemType</a></span></code></p></th>
<th class="head"><p>Command line key</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Hermitian</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_hermitian</span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Generalized Hermitian</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_gen_hermitian</span></code></p></td>
</tr>
<tr class="row-even"><td><p>Non-Hermitian</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_NHEP.html">EPS_NHEP</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_non_hermitian</span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Generalized Non-Hermitian</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GNHEP.html">EPS_GNHEP</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_gen_non_hermitian</span></code></p></td>
</tr>
<tr class="row-even"><td><p>GNHEP with positive (semi-)definite <span class="math notranslate nohighlight">\(B\)</span></p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_PGNHEP.html">EPS_PGNHEP</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_pos_gen_non_hermitian</span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Generalized Hermitian indefinite</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHIEP.html">EPS_GHIEP</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_gen_indefinite</span></code></p></td>
</tr>
<tr class="row-even"><td><p>Bethe-Salpeter</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_BSE.html">EPS_BSE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_bse</span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Hamiltonian</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HAMILT.html">EPS_HAMILT</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_hamiltonian</span></code></p></td>
</tr>
</tbody>
</table>
</div>
<p>The problem type can be specified at run time with the corresponding command line key or, more usually, within the program with the function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetProblemType.html">EPSSetProblemType</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="../../manualpages/EPS/EPSProblemType.html">EPSProblemType</a></span><span class="w"> </span><span class="n">type</span><span class="p">);</span>
</pre></div>
</div>
<p>By default, SLEPc assumes that the problem is non-Hermitian. Some eigensolvers are able to exploit symmetry, that is, they compute a solution for Hermitian problems with less storage and/or computational cost than other methods that ignore this property. Also, symmetric solvers may be more accurate. On the other hand, some eigensolvers in SLEPc only have a symmetric version and will abort if the problem is non-Hermitian. In the case of generalized eigenproblems some considerations apply regarding symmetry, especially in the case of singular <span class="math notranslate nohighlight">\(B\)</span>. This topic is covered in section <a class="reference internal" href="st.html#sec-symm"><span class="std std-ref">Preserving the Symmetry in Generalized Eigenproblems</span></a> and <a class="reference internal" href="st.html#sec-purif"><span class="std std-ref">Purification of Eigenvectors</span></a>. Similarly, if your eigenproblem has a particular algebraic structure listed in table <a class="reference internal" href="#tab-ptype"><span class="std std-ref">Problem types considered in EPS</span></a>, solving it with a structured eigensolver as discussed in section <a class="reference internal" href="#sec-structured"><span class="std std-ref">Structured Eigenvalue Problems</span></a> will result in more accuracy and better efficiency. For all these reasons, the user is strongly recommended to always specify the problem type in the source code.</p>
<p>The characteristics of the problem can be determined with the functions <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSIsGeneralized.html">EPSIsGeneralized</a>()</span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSIsHermitian.html">EPSIsHermitian</a>()</span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSIsPositive.html">EPSIsPositive</a>()</span></code>, and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSIsStructured.html">EPSIsStructured</a>()</span></code>.</p>
<p>The user can specify how many eigenvalues (and eigenvectors) to compute. The default is to compute only one. The following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetDimensions.html">EPSSetDimensions</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">nev</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">ncv</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">mpd</span><span class="p">);</span>
</pre></div>
</div>
<p>allows the specification of the number of eigenvalues to compute, <code class="docutils notranslate"><span class="pre">nev</span></code>. The next argument can be set to prescribe the number of column vectors to be used by the solution algorithm, <code class="docutils notranslate"><span class="pre">ncv</span></code>, that is, the largest dimension of the working subspace. The last argument has to do with a more advanced usage, as explained in section <a class="reference internal" href="#sec-large-nev"><span class="std std-ref">Computing a Large Portion of the Spectrum</span></a>. These parameters can also be set at run time with the options <code class="docutils notranslate"><span class="pre">-eps_nev</span></code>, <code class="docutils notranslate"><span class="pre">-eps_ncv</span></code> and <code class="docutils notranslate"><span class="pre">-eps_mpd</span></code>. For example, the command line</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./program<span class="w"> </span>-eps_nev<span class="w"> </span><span class="m">10</span><span class="w"> </span>-eps_ncv<span class="w"> </span><span class="m">24</span>
</pre></div>
</div>
<p>requests 10 eigenvalues and instructs to use 24 column vectors. Note that <code class="docutils notranslate"><span class="pre">ncv</span></code> must be at least equal to <code class="docutils notranslate"><span class="pre">nev</span></code>, although in general it is recommended (depending on the method) to work with a larger subspace, for instance <span class="math notranslate nohighlight">\(\mathtt{ncv}\geq2\cdot\mathtt{nev}\)</span> or even more. The case that the user requests a relatively large number of eigenpairs is discussed in section <a class="reference internal" href="#sec-large-nev"><span class="std std-ref">Computing a Large Portion of the Spectrum</span></a>.</p>
<p>Instead of specifying the number of wanted eigenvalues <code class="docutils notranslate"><span class="pre">nev</span></code>, it is also possible to specify a threshold with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetThreshold.html">EPSSetThreshold</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="n">thres</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscBool/">PetscBool</a></span><span class="w"> </span><span class="n">rel</span><span class="p">);</span>
</pre></div>
</div>
<p>This usage is discussed for the case of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code> in section <a class="reference internal" href="svd.html#sec-thres"><span class="std std-ref">Using a threshold to specify wanted singular values</span></a>. For details about the differences in case of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>, we refer to the manual page of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetThreshold.html">EPSSetThreshold</a>()</span></code>.</p>
<section id="sec-which">
<h3>Eigenvalues of Interest<a class="headerlink" href="#sec-which" title="Link to this heading">#</a></h3>
<p>For the selection of the portion of the spectrum of interest, there are several alternatives. In real symmetric or complex Hermitian problems, one may want to compute the largest or smallest eigenvalues in magnitude, or the leftmost or rightmost ones, or even all eigenvalues in a given interval. In other problems, in which the eigenvalues can be complex, then one can select eigenvalues depending on the magnitude, or the real part or even the imaginary part. Sometimes the eigenvalues of interest are those closest to a given target value, <span class="math notranslate nohighlight">\(\tau\)</span>, measuring the distance either in the ordinary way or along the real (or imaginary) axis. In some other cases, wanted eigenvalues must be found in a given region of the complex plane. Table <a class="reference internal" href="#tab-portion"><span class="std std-ref">Available possibilities for selection of the eigenvalues of interest</span></a> summarizes all the available possibilities to be selected with the following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetWhichEigenpairs.html">EPSSetWhichEigenpairs</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="../../manualpages/EPS/EPSWhich.html">EPSWhich</a></span><span class="w"> </span><span class="n">which</span><span class="p">);</span>
</pre></div>
</div>
<p>This setting is used both for configuring how the eigensolver seeks eigenvalues (note that not all these possibilities are available for all the solvers) and also for sorting the computed values. The default is to compute the largest magnitude eigenvalues, except for those solvers in which this option is not available. There is another exception related to the use of some spectral transformations, see chapter <a class="reference internal" href="st.html#ch-st"><span class="std std-ref">ST: Spectral Transformation</span></a>.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-portion">
<caption><span class="caption-text">Available possibilities for selection of the eigenvalues of interest</span><a class="headerlink" href="#tab-portion" title="Link to this table">#</a></caption>
<thead>
<tr class="row-odd"><th class="head"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPSWhich</a></span></code></p></th>
<th class="head"><p>Command line key</p></th>
<th class="head"><p>Sorting criterion</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_LARGEST_MAGNITUDE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_largest_magnitude</span></code></p></td>
<td><p>Largest <span class="math notranslate nohighlight">\(|\lambda|\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_SMALLEST_MAGNITUDE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_smallest_magnitude</span></code></p></td>
<td><p>Smallest <span class="math notranslate nohighlight">\(|\lambda|\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_LARGEST_REAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_largest_real</span></code></p></td>
<td><p>Largest <span class="math notranslate nohighlight">\(\mathrm{Re}(\lambda)\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_SMALLEST_REAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_smallest_real</span></code></p></td>
<td><p>Smallest <span class="math notranslate nohighlight">\(\mathrm{Re}(\lambda)\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_LARGEST_IMAGINARY</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_largest_imaginary</span></code></p></td>
<td><p>Largest <span class="math notranslate nohighlight">\(\mathrm{Im}(\lambda)\)</span><a class="footnote-reference brackets" href="#eps-real" id="id3" role="doc-noteref"><span class="fn-bracket">[</span>1<span class="fn-bracket">]</span></a></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_SMALLEST_IMAGINARY</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_smallest_imaginary</span></code></p></td>
<td><p>Smallest <span class="math notranslate nohighlight">\(\mathrm{Im}(\lambda)\)</span><a class="footnote-reference brackets" href="#eps-real" id="id4" role="doc-noteref"><span class="fn-bracket">[</span>1<span class="fn-bracket">]</span></a></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_TARGET_MAGNITUDE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_target_magnitude</span></code></p></td>
<td><p>Smallest <span class="math notranslate nohighlight">\(|\lambda-\tau|\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_TARGET_REAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_target_real</span></code></p></td>
<td><p>Smallest <span class="math notranslate nohighlight">\(|\mathrm{Re}(\lambda-\tau)|\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_TARGET_IMAGINARY</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_target_imaginary</span></code></p></td>
<td><p>Smallest <span class="math notranslate nohighlight">\(|\mathrm{Im}(\lambda-\tau)|\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_ALL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_all</span></code></p></td>
<td><p>All <span class="math notranslate nohighlight">\(\lambda\in[a,b]\)</span> or <span class="math notranslate nohighlight">\(\lambda\in\Omega\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_WHICH_USER</a></span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
<td><p><em>user-defined</em></p></td>
</tr>
</tbody>
</table>
</div>
<p>For the sorting criteria relative to a target value <span class="math notranslate nohighlight">\(\tau\)</span>, the following function must be called in order to specify such value:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetTarget.html">EPSSetTarget</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscScalar/">PetscScalar</a></span><span class="w"> </span><span class="n">target</span><span class="p">);</span>
</pre></div>
</div>
<p>or, alternatively, with the command-line key <code class="docutils notranslate"><span class="pre">-eps_target</span></code>. Note that, since the target is defined as a <a class="reference external" href="https://petsc.org/release/manualpages/Sys/PetscScalar/" title="(in PETSc v3.24)"><span class="xref std std-doc">PetscScalar</span></a>, complex values of <span class="math notranslate nohighlight">\(\tau\)</span> are allowed only in the case of complex scalar builds of the SLEPc library.</p>
<p>The use of a target value makes sense especially when the eigenvalues of interest are located in the interior of the spectrum. Since these eigenvalues are usually more difficult to compute, the eigensolver by itself may not be able to obtain them, and additional tools are normally required. There are two possibilities for this:</p>
<ul class="simple">
<li><p>To use harmonic extraction (see section <a class="reference internal" href="#sec-harmonic"><span class="std std-ref">Computing Interior Eigenvalues with Harmonic Extraction</span></a>), a variant of some solvers that allows a better approximation of interior eigenvalues without changing the way the subspace is built.</p></li>
<li><p>To use a spectral transformation such as shift-and-invert (see chapter <a class="reference internal" href="st.html#ch-st"><span class="std std-ref">ST: Spectral Transformation</span></a>), where the subspace is built from a transformed problem (usually much more costly).</p></li>
</ul>
<p>The special case of computing all eigenvalues in an interval is discussed in the next chapter (sections <a class="reference internal" href="st.html#sec-filter"><span class="std std-ref">Polynomial Filtering</span></a> and <a class="reference internal" href="st.html#sec-slice"><span class="std std-ref">Spectrum Slicing</span></a>), since it is related also to spectral transformations. In this case, instead of a target value the user has to specify the computational interval with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetInterval.html">EPSSetInterval</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="n">a</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="n">b</span><span class="p">);</span>
</pre></div>
</div>
<p>which is equivalent to <code class="docutils notranslate"><span class="pre">-eps_interval</span> <span class="pre">a,b</span></code>.</p>
<p>There is also support for specifying a region of the complex plane so that the eigensolver finds eigenvalues within that region only. This possibility is described in section <a class="reference internal" href="#sec-region"><span class="std std-ref">Specifying a Region for Filtering</span></a>. If <em>all</em> eigenvalues inside the region are required, then a contour-integral method must be used, as described in <span id="id5">[<a class="reference internal" href="../../manualpages/EPS/EPSCISSSetThreshold.html#id66" title="Y. Maeda, T. Sakurai, and J. E. Roman. Contour integral spectrum slicing method in SLEPc. Technical Report STR-11, Universitat Politècnica de València, 2016. URL: https://slepc.upv.es/documentation.">Maeda <em>et al.</em>, 2016</a>]</span>.</p>
<p>Finally, we mention the possibility of defining an arbitrary sorting criterion by means of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSWhich.html">EPS_WHICH_USER</a></span></code> in combination with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetEigenvalueComparison.html">EPSSetEigenvalueComparison</a>()</span></code>.</p>
<p>The selection criteria discussed above are based solely on the eigenvalue. In some special situations, it is necessary to establish a user-defined criterion that also makes use of the eigenvector when deciding which are the most wanted eigenpairs. For these cases, use <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetArbitrarySelection.html">EPSSetArbitrarySelection</a>()</span></code>.</p>
</section>
<section id="left-eigenvectors">
<h3>Left Eigenvectors<a class="headerlink" href="#left-eigenvectors" title="Link to this heading">#</a></h3>
<p>In addition to right eigenvectors, some solvers are able to compute also left eigenvectors, as defined in equation <a class="reference internal" href="#equation-eq-eigstdleft">(2)</a>. The algorithmic variants that compute both left and right eigenvectors are usually called <em>two-sided</em>. By default, SLEPc computes right eigenvectors only. To compute also left eigenvectors, the user should set a flag by calling the following function before <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a>()</span></code>:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetTwoSided.html">EPSSetTwoSided</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscBool/">PetscBool</a></span><span class="w"> </span><span class="n">twosided</span><span class="p">);</span>
</pre></div>
</div>
<p>Note that in some problems such as Hermitian, generalized Hermitian, or some structured eigenproblems, the left eigenvector can be obtained trivially from the right eigenvector. Otherwise, the user must set the two-sided flag <em>before</em> the solver starts to compute, and this is restricted to just a few solvers, see table <a class="reference internal" href="#tab-support"><span class="std std-ref">Supported problem types for all eigensolvers available in SLEPc</span></a>.</p>
</section>
</section>
<section id="selecting-the-eigensolver">
<h2>Selecting the Eigensolver<a class="headerlink" href="#selecting-the-eigensolver" title="Link to this heading">#</a></h2>
<p>The available methods for solving the eigenvalue problems are the following:</p>
<ul class="simple">
<li><p>Basic methods (not recommended except for simple problems):</p>
<ul>
<li><p>Power Iteration with deflation. When combined with shift-and-invert (see chapter <a class="reference internal" href="st.html#ch-st"><span class="std std-ref">ST: Spectral Transformation</span></a>), it is equivalent to the inverse iteration. Also, this solver embeds the Rayleigh Quotient iteration (RQI) by allowing variable shifts. Additionally, it provides the nonlinear inverse iteration method for the case that the problem matrix is a nonlinear operator (for this advanced usage, see <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSPowerSetNonlinear.html">EPSPowerSetNonlinear</a>()</span></code>).</p></li>
<li><p>Subspace Iteration with Rayleigh-Ritz projection and locking.</p></li>
<li><p>Arnoldi method with explicit restart and deflation.</p></li>
<li><p>Lanczos with explicit restart, deflation, and different reorthogonalization strategies.</p></li>
</ul>
</li>
<li><p>Krylov-Schur, a variation of Arnoldi with a very effective restarting technique. In the case of symmetric problems, this is equivalent to the thick-restart Lanczos method.</p></li>
<li><p>Generalized Davidson, a simple iteration based on subspace expansion with the preconditioned residual.</p></li>
<li><p>Jacobi-Davidson, a preconditioned eigensolver with an effective correction equation.</p></li>
<li><p>RQCG, a basic conjugate gradient iteration for the minimization of the Rayleigh quotient.</p></li>
<li><p>LOBPCG, the locally-optimal block preconditioned conjugate gradient.</p></li>
<li><p>CISS, a contour-integral solver that allows computing all eigenvalues in a given region.</p></li>
<li><p>Lyapunov inverse iteration, to compute rightmost eigenvalues.</p></li>
</ul>
<div class="pst-scrollable-table-container"><table class="table" id="tab-solvers">
<caption><span class="caption-text">Eigenvalue solvers available in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> module</span><a class="headerlink" href="#tab-solvers" title="Link to this table">#</a></caption>
<thead>
<tr class="row-odd"><th class="head"><p>Method</p></th>
<th class="head"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSType.html">EPSType</a></span></code></p></th>
<th class="head"><p>Options Database</p></th>
<th class="head"><p>Default</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Power / Inverse / RQI</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSPOWER.html">EPSPOWER</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">power</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Subspace Iteration</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSUBSPACE.html">EPSSUBSPACE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">subspace</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Arnoldi</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSARNOLDI.html">EPSARNOLDI</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">arnoldi</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Lanczos</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSLANCZOS.html">EPSLANCZOS</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">lanczos</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Krylov-Schur</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSKRYLOVSCHUR.html">EPSKRYLOVSCHUR</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">krylovschur</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\star\)</span></p></td>
</tr>
<tr class="row-odd"><td><p>Generalized Davidson</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGD.html">EPSGD</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">gd</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Jacobi-Davidson</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSJD.html">EPSJD</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">jd</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Rayleigh quotient CG</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSRQCG.html">EPSRQCG</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">rqcg</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>LOBPCG</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSLOBPCG.html">EPSLOBPCG</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">lobpcg</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Contour integral SS</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSCISS.html">EPSCISS</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">ciss</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Lyapunov Inverse Iteration</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSLYAPII.html">EPSLYAPII</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">lyapii</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Wrapper to LAPACK</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSLAPACK.html">EPSLAPACK</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">lapack</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Wrapper to ARPACK</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSARPACK.html">EPSARPACK</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">arpack</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Wrapper to BLOPEX</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSBLOPEX.html">EPSBLOPEX</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">blopex</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Wrapper to PRIMME</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSPRIMME.html">EPSPRIMME</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">primme</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Wrapper to FEAST</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSFEAST.html">EPSFEAST</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">feast</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Wrapper to ScaLAPACK</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSCALAPACK.html">EPSSCALAPACK</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">scalapack</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Wrapper to ELPA</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSELPA.html">EPSELPA</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">elpa</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Wrapper to ELEMENTAL</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSELEMENTAL.html">EPSELEMENTAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">elemental</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Wrapper to EVSL</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSEVSL.html">EPSEVSL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">evsl</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Wrapper to CHASE</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSCHASE.html">EPSCHASE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">chase</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
</tbody>
</table>
</div>
<p>The default solver is Krylov-Schur. A detailed description of the implemented algorithms is provided in the <a class="reference internal" href="../index.html#str"><span class="std std-ref">SLEPc Technical Reports</span></a>. In addition to these methods, SLEPc also provides wrappers to external packages such as ARPACK, or PRIMME. A complete list of these interfaces can be found in section <a class="reference internal" href="extra.html#sec-wrap"><span class="std std-ref">Wrappers to External Libraries</span></a>. Note that some of these packages (LAPACK, ScaLAPACK, ELPA, ELEMENTAL) perform <em>dense</em> computations and hence return the full eigendecomposition (furthermore, take into account that LAPACK is a sequential library so the corresponding solver should be used only for debugging purposes with small problem sizes).</p>
<p>The solution method can be specified procedurally or via the command line. The application programmer can set it by means of:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetType.html">EPSSetType</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="../../manualpages/EPS/EPSType.html">EPSType</a></span><span class="w"> </span><span class="n">method</span><span class="p">);</span>
</pre></div>
</div>
<p>while the user writes the options database command <code class="docutils notranslate"><span class="pre">-eps_type</span></code> followed by the name of the method (see table <a class="reference internal" href="#tab-solvers"><span class="std std-ref">Eigenvalue solvers available in the EPS module</span></a>).</p>
<p>Not all the methods can be used for all problem types. Table <a class="reference internal" href="#tab-support"><span class="std std-ref">Supported problem types for all eigensolvers available in SLEPc</span></a> summarizes the scope of each eigensolver by listing which portion of the spectrum can be selected (as defined in table <a class="reference internal" href="#tab-portion"><span class="std std-ref">Available possibilities for selection of the eigenvalues of interest</span></a>) and which problem types are supported (as defined in table <a class="reference internal" href="#tab-ptype"><span class="std std-ref">Problem types considered in EPS</span></a>). Note that the structured problem types are not considered here, see section <a class="reference internal" href="#sec-structured"><span class="std std-ref">Structured Eigenvalue Problems</span></a>. The table also indicates whether the solvers are available or not in the complex version of SLEPc, and if they have a two-sided variant.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-support">
<caption><span class="caption-text">Supported problem types for all eigensolvers available in SLEPc</span><a class="headerlink" href="#tab-support" title="Link to this table">#</a></caption>
<thead>
<tr class="row-odd"><th class="head text-center"><p>Method</p></th>
<th class="head text-center"><p>Portion of spectrum</p></th>
<th class="head text-center"><p>Problem type</p></th>
<th class="head text-center"><p>Real/complex</p></th>
<th class="head text-center"><p>Two-sided</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">power</span></code></p></td>
<td class="text-center"><p>Largest <span class="math notranslate nohighlight">\(|\lambda|\)</span></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p>yes</p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">subspace</span></code></p></td>
<td class="text-center"><p>Largest <span class="math notranslate nohighlight">\(|\lambda|\)</span></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">arnoldi</span></code></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">lanczos</span></code></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code></p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">krylovschur</span></code></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p>yes</p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">gd</span></code></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">jd</span></code></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">rqcg</span></code></p></td>
<td class="text-center"><p>Smallest <span class="math notranslate nohighlight">\(\mathrm{Re}(\lambda)\)</span></p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code></p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">lobpcg</span></code></p></td>
<td class="text-center"><p>Largest and smallest <span class="math notranslate nohighlight">\(\mathrm{Re}(\lambda)\)</span></p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code></p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">ciss</span></code></p></td>
<td class="text-center"><p>All <span class="math notranslate nohighlight">\(\lambda\)</span> in region</p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">lyapii</span></code></p></td>
<td class="text-center"><p>Largest <span class="math notranslate nohighlight">\(\mathrm{Re}(\lambda)\)</span></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">lapack</span></code></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p>yes</p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">arpack</span></code></p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">blopex</span></code></p></td>
<td class="text-center"><p>Smallest <span class="math notranslate nohighlight">\(\mathrm{Re}(\lambda)\)</span></p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code></p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">primme</span></code></p></td>
<td class="text-center"><p>Largest and smallest <span class="math notranslate nohighlight">\(\mathrm{Re}(\lambda)\)</span></p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code></p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">feast</span></code></p></td>
<td class="text-center"><p>All <span class="math notranslate nohighlight">\(\mathrm{Re}(\lambda)\)</span> in an interval</p></td>
<td class="text-center"><p>any</p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">scalapack</span></code></p></td>
<td class="text-center"><p>All <span class="math notranslate nohighlight">\(\lambda\)</span></p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code></p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">elpa</span></code></p></td>
<td class="text-center"><p>All <span class="math notranslate nohighlight">\(\lambda\)</span></p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code></p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">elemental</span></code></p></td>
<td class="text-center"><p>All <span class="math notranslate nohighlight">\(\lambda\)</span></p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code></p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">evsl</span></code></p></td>
<td class="text-center"><p>All <span class="math notranslate nohighlight">\(\lambda\)</span> in interval</p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code></p></td>
<td class="text-center"><p>real</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td class="text-center"><p><code class="docutils notranslate"><span class="pre">chase</span></code></p></td>
<td class="text-center"><p>Smallest <span class="math notranslate nohighlight">\(\mathrm{Re}(\lambda)\)</span></p></td>
<td class="text-center"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code></p></td>
<td class="text-center"><p>both</p></td>
<td class="text-center"><p><code class="docutils notranslate"> </code></p></td>
</tr>
</tbody>
</table>
</div>
</section>
<section id="sec-retrsol">
<h2>Retrieving the Solution<a class="headerlink" href="#sec-retrsol" title="Link to this heading">#</a></h2>
<p>Once the call to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a>()</span></code> is complete, all the data associated with the solution of the eigenproblem are kept internally in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object. This information can be obtained by the calling program by means of a set of functions described in this section.</p>
<p>As explained below, the number of computed solutions depends on the convergence and, therefore, it may be different from the number of solutions requested by the user. So the first task is to find out how many solutions are available, with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSGetConverged.html">EPSGetConverged</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="o">*</span><span class="n">nconv</span><span class="p">);</span>
</pre></div>
</div>
<p>Usually, the number of converged solutions, <code class="docutils notranslate"><span class="pre">nconv</span></code>, will be equal to <code class="docutils notranslate"><span class="pre">nev</span></code>, but in general it can be a number ranging from 0 to <code class="docutils notranslate"><span class="pre">ncv</span></code> (here, <code class="docutils notranslate"><span class="pre">nev</span></code> and <code class="docutils notranslate"><span class="pre">ncv</span></code> are the arguments of function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetDimensions.html">EPSSetDimensions</a>()</span></code>).</p>
<section id="the-computed-solution">
<h3>The Computed Solution<a class="headerlink" href="#the-computed-solution" title="Link to this heading">#</a></h3>
<p>The user may be interested in the eigenvalues, or the eigenvectors, or both. The next function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSGetEigenpair.html">EPSGetEigenpair</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">j</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscScalar/">PetscScalar</a></span><span class="w"> </span><span class="o">*</span><span class="n">kr</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscScalar/">PetscScalar</a></span><span class="w"> </span><span class="o">*</span><span class="n">ki</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">xr</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">xi</span><span class="p">);</span>
</pre></div>
</div>
<p>returns the <span class="math notranslate nohighlight">\(j\)</span>-th computed eigenvalue/eigenvector pair. Typically, this function is called inside a loop for each value of <code class="docutils notranslate"><span class="pre">j</span></code> from 0 to <code class="docutils notranslate"><span class="pre">nconv</span></code>–1. Note that eigenvalues are ordered according to the same criterion specified with function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetWhichEigenpairs.html">EPSSetWhichEigenpairs</a>()</span></code> for selecting the portion of the spectrum of interest. The meaning of the last 4 arguments depends on whether SLEPc has been compiled for real or complex scalars, as detailed below. The eigenvectors are normalized so that they have a unit 2-norm, except for problem type <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code> in which case returned eigenvectors have a unit <span class="math notranslate nohighlight">\(B\)</span>-norm.</p>
<p>In case they are available, the left eigenvectors can be extracted with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSGetLeftEigenvector.html">EPSGetLeftEigenvector</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">j</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">yr</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">yi</span><span class="p">);</span>
</pre></div>
</div>
<section id="real-vs-complex-slepc">
<h4>Real vs Complex SLEPc<a class="headerlink" href="#real-vs-complex-slepc" title="Link to this heading">#</a></h4>
<p>In case of PETSc/SLEPc built with <strong>real scalars</strong>, all <a class="reference external" href="https://petsc.org/release/manualpages/Mat/Mat/" title="(in PETSc v3.24)"><span class="xref std std-doc">Mat</span></a> and <a class="reference external" href="https://petsc.org/release/manualpages/Vec/Vec/" title="(in PETSc v3.24)"><span class="xref std std-doc">Vec</span></a> objects are real. The computed approximate solution returned by the function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetEigenpair.html">EPSGetEigenpair</a></span></code> is stored in the following way: <code class="docutils notranslate"><span class="pre">kr</span></code> and <code class="docutils notranslate"><span class="pre">ki</span></code> contain the real and imaginary parts of the eigenvalue, respectively, and <code class="docutils notranslate"><span class="pre">xr</span></code> and <code class="docutils notranslate"><span class="pre">xi</span></code> contain the associated eigenvector. Two cases can be distinguished:</p>
<ul class="simple">
<li><p>When <code class="docutils notranslate"><span class="pre">ki</span></code> is zero, it means that the <span class="math notranslate nohighlight">\(j\)</span>-th eigenvalue is a real number. In this case, <code class="docutils notranslate"><span class="pre">kr</span></code> is the eigenvalue and <code class="docutils notranslate"><span class="pre">xr</span></code> is the corresponding eigenvector. The vector <code class="docutils notranslate"><span class="pre">xi</span></code> is set to all zeros.</p></li>
<li><p>If <code class="docutils notranslate"><span class="pre">ki</span></code> is different from zero, then the <span class="math notranslate nohighlight">\(j\)</span>-th eigenvalue is a complex number and, therefore, it is part of a complex conjugate pair. Thus, the <span class="math notranslate nohighlight">\(j\)</span>-th eigenvalue is <code class="docutils notranslate"><span class="pre">kr</span></code><span class="math notranslate nohighlight">\(+\,i\cdot\)</span><code class="docutils notranslate"><span class="pre">ki</span></code>. With respect to the eigenvector, <code class="docutils notranslate"><span class="pre">xr</span></code> stores the real part of the eigenvector and <code class="docutils notranslate"><span class="pre">xi</span></code> the imaginary part, that is, the <span class="math notranslate nohighlight">\(j\)</span>-th eigenvector is <code class="docutils notranslate"><span class="pre">xr</span></code><span class="math notranslate nohighlight">\(+\,i\cdot\)</span><code class="docutils notranslate"><span class="pre">xi</span></code>. The <span class="math notranslate nohighlight">\((j+1)\)</span>-th eigenvalue (and eigenvector) will be the corresponding complex conjugate and will be returned when function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetEigenpair.html">EPSGetEigenpair</a></span></code> is invoked with index <code class="docutils notranslate"><span class="pre">j</span></code>+1. Note that the sign of the imaginary part is returned correctly in all cases (users need not change signs).</p></li>
</ul>
<p>In case of PETSc/SLEPc built with <strong>complex scalars</strong>, all <a class="reference external" href="https://petsc.org/release/manualpages/Mat/Mat/" title="(in PETSc v3.24)"><span class="xref std std-doc">Mat</span></a> and <a class="reference external" href="https://petsc.org/release/manualpages/Vec/Vec/" title="(in PETSc v3.24)"><span class="xref std std-doc">Vec</span></a> objects are complex. The computed solution returned by function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetEigenpair.html">EPSGetEigenpair</a>()</span></code> is the following: <code class="docutils notranslate"><span class="pre">kr</span></code> contains the (complex) eigenvalue and <code class="docutils notranslate"><span class="pre">xr</span></code> contains the corresponding (complex) eigenvector. In this case, <code class="docutils notranslate"><span class="pre">ki</span></code> and <code class="docutils notranslate"><span class="pre">xi</span></code> are not used (set to all zeros).</p>
</section>
</section>
<section id="sec-errbnd">
<h3>Reliability of the Computed Solution<a class="headerlink" href="#sec-errbnd" title="Link to this heading">#</a></h3>
<p>In this subsection, we discuss how a-posteriori error bounds can be obtained in order to assess the accuracy of the computed solutions. These bounds are based on the so-called residual vector, defined as</p>
<div class="math notranslate nohighlight" id="equation-eq-residual-vector">
<span class="eqno">(7)<a class="headerlink" href="#equation-eq-residual-vector" title="Link to this equation">#</a></span>\[r=A\tilde{x}-\tilde{\lambda}\tilde{x},\]</div>
<p>or <span class="math notranslate nohighlight">\(r=A\tilde{x}-\tilde{\lambda}B\tilde{x}\)</span> in the case of a generalized problem, where <span class="math notranslate nohighlight">\(\tilde{\lambda}\)</span> and <span class="math notranslate nohighlight">\(\tilde{x}\)</span> represent any of the <code class="docutils notranslate"><span class="pre">nconv</span></code> computed eigenpairs delivered by <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetEigenpair.html">EPSGetEigenpair</a>()</span></code> (note that this function returns a normalized <span class="math notranslate nohighlight">\(\tilde{x}\)</span>).</p>
<p>In the case of Hermitian problems, it is possible to demonstrate the following property (see for example <span id="id6">[<a class="reference internal" href="#id23" title="Y. Saad. Numerical Methods for Large Eigenvalue Problems: Theory and Algorithms. John Wiley and Sons, New York, 1992.">Saad, 1992</a>, ch. 3]</span>):</p>
<div class="math notranslate nohighlight" id="equation-eq-reserr">
<span class="eqno">(8)<a class="headerlink" href="#equation-eq-reserr" title="Link to this equation">#</a></span>\[|\lambda-\tilde{\lambda}|\leq \|r\|_2,\]</div>
<p>where <span class="math notranslate nohighlight">\(\lambda\)</span> is an exact eigenvalue. Therefore, the 2-norm of the residual vector can be used as a bound for the absolute error in the eigenvalue.</p>
<p>In the case of non-Hermitian problems, the situation is worse because no simple relation such as equation <a class="reference internal" href="#equation-eq-reserr">(8)</a> is available. This means that in this case the residual norms may still give an indication of the actual error but the user should be aware that they may sometimes be completely wrong, especially in the case of highly non-normal matrices. A better bound would involve also the residual norm of the left eigenvector.</p>
<p>With respect to eigenvectors, we have a similar scenario in the sense that bounds for the error may be established in the Hermitian case only, for example the following one:</p>
<div class="math notranslate nohighlight" id="equation-eq-eigvector-residual">
<span class="eqno">(9)<a class="headerlink" href="#equation-eq-eigvector-residual" title="Link to this equation">#</a></span>\[\sin \theta(x,\tilde{x})\leq \frac{\|r\|_2}{\delta},\]</div>
<p>where <span class="math notranslate nohighlight">\(\theta(x,\tilde{x})\)</span> is the angle between the computed and exact eigenvectors, and <span class="math notranslate nohighlight">\(\delta\)</span> is the distance from <span class="math notranslate nohighlight">\(\tilde{\lambda}\)</span> to the rest of the spectrum. This bound is not provided by SLEPc because <span class="math notranslate nohighlight">\(\delta\)</span> is not available. The above expression is given here simply to warn the user about the fact that accuracy of eigenvectors may be deficient in the case of clustered eigenvalues.</p>
<p>In the case of non-Hermitian problems, SLEPc provides the alternative of retrieving an orthonormal basis of an invariant subspace instead of getting individual eigenvectors. This is done with the following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSGetInvariantSubspace.html">EPSGetInvariantSubspace</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">v</span><span class="p">[]);</span>
</pre></div>
</div>
<p>This is sufficient in some applications and is safer from the numerical point of view.</p>
<section id="computation-of-bounds">
<h4>Computation of Bounds<a class="headerlink" href="#computation-of-bounds" title="Link to this heading">#</a></h4>
<p>It is sometimes useful to compute error bounds based on the norm of the residual <span class="math notranslate nohighlight">\(r_j\)</span>, to assess the accuracy of the computed solution. The bound can be made in absolute terms, as in equation <a class="reference internal" href="#equation-eq-reserr">(8)</a>, or alternatively the error can be expressed relative to the eigenvalue or to the matrix norms. For this, the following function can be used:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSComputeError.html">EPSComputeError</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">j</span><span class="p">,</span><span class="n"><a href="../../manualpages/EPS/EPSErrorType.html">EPSErrorType</a></span><span class="w"> </span><span class="n">type</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="o">*</span><span class="n">error</span><span class="p">);</span>
</pre></div>
</div>
<p>The types of errors that can be computed are summarized in table <a class="reference internal" href="#tab-errors"><span class="std std-ref">Available expressions for computing error bounds</span></a>. The way in which the error is computed is unrelated to the error estimation used internally in the solver for convergence checking, as described below. Also note that in the case of two-sided eigensolvers, the error bounds are based on <span class="math notranslate nohighlight">\(\max\{\|\ell_j\|,\|r_j\|\}\)</span>, where the left residual <span class="math notranslate nohighlight">\(\ell_j\)</span> is defined as <span class="math notranslate nohighlight">\(\ell_j=A^*\tilde{y}-\bar{\tilde\lambda}\tilde{y}\)</span>.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-errors">
<caption><span class="caption-text">Available expressions for computing error bounds</span><a class="headerlink" href="#tab-errors" title="Link to this table">#</a></caption>
<thead>
<tr class="row-odd"><th class="head"><p>Error type</p></th>
<th class="head"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSErrorType.html">EPSErrorType</a></span></code></p></th>
<th class="head"><p>Command line key</p></th>
<th class="head"><p>Error bound</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Absolute error</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSErrorType.html">EPS_ERROR_ABSOLUTE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_error_absolute</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\|r\|\)</span></p></td>
</tr>
<tr class="row-odd"><td><p>Relative error</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSErrorType.html">EPS_ERROR_RELATIVE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_error_relative</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\|r\|/|\lambda|\)</span></p></td>
</tr>
<tr class="row-even"><td><p>Backward error</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSErrorType.html">EPS_ERROR_BACKWARD</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_error_backward</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\|r\|/(\|A\|+|\lambda|\|B\|)\)</span></p></td>
</tr>
</tbody>
</table>
</div>
</section>
</section>
<section id="sec-monitor">
<h3>Controlling and Monitoring Convergence<a class="headerlink" href="#sec-monitor" title="Link to this heading">#</a></h3>
<p>All the eigensolvers provided by SLEPc are iterative in nature, meaning that the solutions are (usually) improved at each iteration until they are sufficiently accurate, that is, until convergence is achieved. The number of iterations required by the process can be obtained with the function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSGetIterationNumber.html">EPSGetIterationNumber</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="o">*</span><span class="n">its</span><span class="p">);</span>
</pre></div>
</div>
<p>which returns in argument <code class="docutils notranslate"><span class="pre">its</span></code> either the iteration number at which convergence was successfully reached, or the iteration at which a problem was detected.</p>
<p>The user specifies when a solution should be considered sufficiently accurate by means of a tolerance. An approximate eigenvalue is considered to be converged if the error estimate associated with it is below the specified tolerance. The default value of the tolerance is <span class="math notranslate nohighlight">\(10^{-8}\)</span> (in case of PETSc/SLEPc built for double precision) and can be changed at run time with <code class="docutils notranslate"><span class="pre">-eps_tol</span> <span class="pre"><tol></span></code> or inside the program with the function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetTolerances.html">EPSSetTolerances</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="n">tol</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">max_it</span><span class="p">);</span>
</pre></div>
</div>
<p>The third parameter of this function allows the programmer to modify the maximum number of iterations allowed to the solution algorithm, which can also be set via <code class="docutils notranslate"><span class="pre">-eps_max_it</span> <span class="pre"><its></span></code>.</p>
<section id="convergence-check">
<h4>Convergence Check<a class="headerlink" href="#convergence-check" title="Link to this heading">#</a></h4>
<p>The error estimates used for the convergence test are based on the residual norm, as discussed in section <a class="reference internal" href="#sec-errbnd"><span class="std std-ref">Reliability of the Computed Solution</span></a>. Most eigensolvers explicitly compute the residual of the relevant eigenpairs during the iteration, but Krylov solvers use a cheap formula instead, allowing to track many eigenpairs simultaneously. When using a spectral transformation, this formula may give too optimistic bounds (corresponding to the residual of the transformed problem, not the original problem). In such cases, the users can force the computation of the residual with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetTrueResidual.html">EPSSetTrueResidual</a>()</span></code>.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-convergence">
<caption><span class="caption-text">Available possibilities for the convergence criterion</span><a class="headerlink" href="#tab-convergence" title="Link to this table">#</a></caption>
<thead>
<tr class="row-odd"><th class="head"><p>Convergence criterion</p></th>
<th class="head"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSConv.html">EPSConv</a></span></code></p></th>
<th class="head"><p>Command line key</p></th>
<th class="head"><p>Error bound</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Absolute</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSConv.html">EPS_CONV_ABS</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_conv_abs</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\|r\|\)</span></p></td>
</tr>
<tr class="row-odd"><td><p>Relative to eigenvalue</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSConv.html">EPS_CONV_REL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_conv_rel</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\|r\|/|\lambda|\)</span></p></td>
</tr>
<tr class="row-even"><td><p>Relative to matrix norms</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSConv.html">EPS_CONV_NORM</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_conv_norm</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\|r\|/(\|A\|+|\lambda|\|B\|)\)</span></p></td>
</tr>
<tr class="row-odd"><td><p>User-defined</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSConv.html">EPS_CONV_USER</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-eps_conv_user</span></code></p></td>
<td><p><em>user function</em></p></td>
</tr>
</tbody>
</table>
</div>
<p>From the residual norm, the error bound can be computed in different ways, see table <a class="reference internal" href="#tab-convergence"><span class="std std-ref">Available possibilities for the convergence criterion</span></a>. This can be set via the corresponding command-line switch or with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetConvergenceTest.html">EPSSetConvergenceTest</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="../../manualpages/EPS/EPSConv.html">EPSConv</a></span><span class="w"> </span><span class="n">conv</span><span class="p">);</span>
</pre></div>
</div>
<p>The default is to use the criterion relative to the eigenvalue (note: for computing eigenvalues close to the origin this criterion will likely give very poor accuracy, so the user is advised to use <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSConv.html">EPS_CONV_ABS</a></span></code> in that case). Finally, a custom convergence criterion may be established by specifying a user function (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetConvergenceTestFunction.html">EPSSetConvergenceTestFunction</a>()</span></code>).</p>
<p>Error estimates used internally by eigensolvers for checking convergence may be different from the error bounds provided by <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSComputeError.html">EPSComputeError</a>()</span></code>. At the end of the solution process, error estimates are available via <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetErrorEstimate.html">EPSGetErrorEstimate</a>()</span></code>.</p>
<p>By default, the eigensolver will stop iterating when the current number of eigenpairs satisfying the convergence test is equal to (or greater than) the number of requested eigenpairs (or if the maximum number of iterations has been reached). However, it is also possible to provide a user-defined stopping test that may decide to quit earlier, see <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetStoppingTest.html">EPSSetStoppingTest</a>()</span></code>.</p>
</section>
<section id="monitors">
<h4>Monitors<a class="headerlink" href="#monitors" title="Link to this heading">#</a></h4>
<p>Error estimates can be displayed during execution of the solution algorithm, as a way of monitoring convergence. There are several such monitors available. The user can activate them via the options database (see examples below), or within the code with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSMonitorSet.html">EPSMonitorSet</a>()</span></code>. By default, the solvers run silently without displaying information about the iteration. Also, application programmers can provide their own routines to perform the monitoring by using the function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSMonitorSet.html">EPSMonitorSet</a>()</span></code>.</p>
<p>The most basic monitor prints one approximate eigenvalue together with its associated error estimate in each iteration. The shown eigenvalue is the first unconverged one.</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex9<span class="w"> </span>-eps_nev<span class="w"> </span><span class="m">1</span><span class="w"> </span>-eps_tol<span class="w"> </span>1e-6<span class="w"> </span>-eps_monitor
<span class="go"> 1 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) -0.0695109+2.10989i (2.38956768e-01)</span>
<span class="go"> 2 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) -0.0231046+2.14902i (1.09212525e-01)</span>
<span class="go"> 3 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) -0.000633399+2.14178i (2.67086904e-02)</span>
<span class="go"> 4 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) 9.89074e-05+2.13924i (6.62097793e-03)</span>
<span class="go"> 5 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) -0.000149404+2.13976i (1.53444214e-02)</span>
<span class="go"> 6 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) 0.000183676+2.13939i (2.85521004e-03)</span>
<span class="go"> 7 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) 0.000192479+2.13938i (9.97563492e-04)</span>
<span class="go"> 8 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) 0.000192534+2.13938i (1.77259863e-04)</span>
<span class="go"> 9 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) 0.000192557+2.13938i (2.82539990e-05)</span>
<span class="go"> 10 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=0 first unconverged value (error) 0.000192559+2.13938i (2.51440008e-06)</span>
<span class="go"> 11 <a href="../../manualpages/EPS/EPS.html">EPS</a> nconv=2 first unconverged value (error) -0.671923+2.52712i (8.92724972e-05)</span>
</pre></div>
</div>
<p>Graphical monitoring (in an X display) is also available with <code class="docutils notranslate"><span class="pre">-eps_monitor</span> <span class="pre">draw::draw_lg</span></code>. Figure <a class="reference internal" href="#fig-plot-monitor"><span class="std std-ref">Default convergence monitor</span></a> shows the result of the following sample command line:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex9<span class="w"> </span>-n<span class="w"> </span><span class="m">200</span><span class="w"> </span>-eps_nev<span class="w"> </span><span class="m">12</span><span class="w"> </span>-eps_tol<span class="w"> </span>1e-12<span class="w"> </span>-eps_monitor<span class="w"> </span>draw::draw_lg<span class="w"> </span>-draw_pause<span class="w"> </span>.2
</pre></div>
</div>
<p>Again, only the error estimate of one eigenvalue is drawn. The spikes in the last part of the plot indicate convergence of one eigenvalue and switching to the next.</p>
<figure class="align-default" id="fig-plot-monitor">
<a class="reference internal image-reference" href="../../_images/monitor.png"><img alt="SLEPc convergence monitor" src="../../_images/monitor.png" style="width: 174.0px; height: 174.0px;" />
</a>
<figcaption>
<p><span class="caption-text">Default convergence monitor</span><a class="headerlink" href="#fig-plot-monitor" title="Link to this image">#</a></p>
</figcaption>
</figure>
<p>The two previously mentioned monitors have an alternative version (<code class="docutils notranslate"><span class="pre">*_all</span></code>) that processes all eigenvalues instead of just the first one. Figure <a class="reference internal" href="#fig-plot-monitorall"><span class="std std-ref">Simultaneous convergence monitor for all eigenvalues</span></a> corresponds to the same example but with <code class="docutils notranslate"><span class="pre">-eps_monitor_all</span> <span class="pre">draw::draw_lg</span></code>. Note that these variants have a side effect: they force the computation of all error estimates even if the method would not normally do so.</p>
<figure class="align-default" id="fig-plot-monitorall">
<a class="reference internal image-reference" href="../../_images/monitorall.png"><img alt="SLEPc convergence monitor (all eigenvalues)" src="../../_images/monitorall.png" style="width: 174.0px; height: 174.0px;" />
</a>
<figcaption>
<p><span class="caption-text">Simultaneous convergence monitor for all eigenvalues</span><a class="headerlink" href="#fig-plot-monitorall" title="Link to this image">#</a></p>
</figcaption>
</figure>
<p>A less verbose monitor is <code class="docutils notranslate"><span class="pre">-eps_monitor_conv</span></code>, which simply displays the iteration number at which convergence takes place. Note that several monitors can be used at the same time.</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex9<span class="w"> </span>-n<span class="w"> </span><span class="m">200</span><span class="w"> </span>-eps_nev<span class="w"> </span><span class="m">8</span><span class="w"> </span>-eps_tol<span class="w"> </span>1e-12<span class="w"> </span>-eps_monitor_conv
<span class="go"> 79 <a href="../../manualpages/EPS/EPS.html">EPS</a> converged value (error) #0 4.64001e-06+2.13951i (8.26091148e-13)</span>
<span class="go"> 79 <a href="../../manualpages/EPS/EPS.html">EPS</a> converged value (error) #1 4.64001e-06-2.13951i (8.26091148e-13)</span>
<span class="go"> 93 <a href="../../manualpages/EPS/EPS.html">EPS</a> converged value (error) #2 -0.674926+2.52867i (6.85260521e-13)</span>
<span class="go"> 93 <a href="../../manualpages/EPS/EPS.html">EPS</a> converged value (error) #3 -0.674926-2.52867i (6.85260521e-13)</span>
<span class="go"> 94 <a href="../../manualpages/EPS/EPS.html">EPS</a> converged value (error) #4 -1.79963+3.03259i (5.48825386e-13)</span>
<span class="go"> 94 <a href="../../manualpages/EPS/EPS.html">EPS</a> converged value (error) #5 -1.79963-3.03259i (5.48825386e-13)</span>
<span class="go"> 98 <a href="../../manualpages/EPS/EPS.html">EPS</a> converged value (error) #6 -3.37383+3.55626i (2.74909207e-13)</span>
<span class="go"> 98 <a href="../../manualpages/EPS/EPS.html">EPS</a> converged value (error) #7 -3.37383-3.55626i (2.74909207e-13)</span>
</pre></div>
</div>
</section>
</section>
<section id="sec-epsviewers">
<h3>Viewing the Solution<a class="headerlink" href="#sec-epsviewers" title="Link to this heading">#</a></h3>
<p>The computed solution (eigenvalues and eigenvectors) can be viewed in different ways, exploiting the flexibility of <a class="reference external" href="https://petsc.org/release/manualpages/Viewer/PetscViewer/" title="(in PETSc v3.24)"><span class="xref std std-doc">PetscViewer</span></a>s. The API functions for this are <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSValuesView.html">EPSValuesView</a>()</span></code> and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSVectorsView.html">EPSVectorsView</a>()</span></code>. We next illustrate their usage via the command line.</p>
<p>The command-line option <code class="docutils notranslate"><span class="pre">-eps_view_values</span></code> shows the computed eigenvalues on the standard output at the end of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a>()</span></code>. It admits an argument to specify <a class="reference external" href="https://petsc.org/release/manualpages/Viewer/PetscViewer/" title="(in PETSc v3.24)"><span class="xref std std-doc">PetscViewer</span></a> options, for instance the following will create a Matlab command file <code class="docutils notranslate"><span class="pre">myeigenvalues.m</span></code> to load the eigenvalues in Matlab:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex1<span class="w"> </span>-n<span class="w"> </span><span class="m">120</span><span class="w"> </span>-eps_nev<span class="w"> </span><span class="m">8</span><span class="w"> </span>-eps_view_values<span class="w"> </span>:myeigenvalues.m:ascii_matlab
</pre></div>
</div>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="go">Lambda_EPS_0xb430f0_0 = [</span>
<span class="go">3.9993259306070224e+00</span>
<span class="go">3.9973041767976509e+00</span>
<span class="go">3.9939361013742269e+00</span>
<span class="go">3.9892239746533216e+00</span>
<span class="go">3.9831709729353331e+00</span>
<span class="go">3.9757811763634532e+00</span>
<span class="go">3.9670595661733632e+00</span>
<span class="go">3.9570120213355646e+00</span>
<span class="go">];</span>
</pre></div>
</div>
<p>One particular instance of this option is <code class="docutils notranslate"><span class="pre">-eps_view_values</span> <span class="pre">draw</span></code>, that will plot the computed approximations of the eigenvalues on an X window. See Figure <a class="reference internal" href="#fig-plot-eigs"><span class="std std-ref">Eigenvalues plot</span></a> for an example.</p>
<figure class="align-default" id="fig-plot-eigs">
<a class="reference internal image-reference" href="../../_images/ploteigs.png"><img alt="SLEPc eigenvalues plot" src="../../_images/ploteigs.png" style="width: 174.0px; height: 174.0px;" />
</a>
<figcaption>
<p><span class="caption-text">Eigenvalues plot</span><a class="headerlink" href="#fig-plot-eigs" title="Link to this image">#</a></p>
</figcaption>
</figure>
<p>Similarly, eigenvectors may be viewed with <code class="docutils notranslate"><span class="pre">-eps_view_vectors</span></code>, either in text form, in Matlab format, in binary format, or as a draw. All eigenvectors are viewed, one after the other. As an example, the next line will dump eigenvectors to the binary file <code class="docutils notranslate"><span class="pre">evec.bin</span></code>:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex1<span class="w"> </span>-n<span class="w"> </span><span class="m">120</span><span class="w"> </span>-eps_nev<span class="w"> </span><span class="m">8</span><span class="w"> </span>-eps_view_vectors<span class="w"> </span>binary:evec.bin
</pre></div>
</div>
<p>Two more related functions are available: <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSErrorView.html">EPSErrorView</a>()</span></code> and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSConvergedReasonView.html">EPSConvergedReasonView</a>()</span></code>. These will show computed errors and the converged reason (plus number of iterations), respectively. Again, we illustrate its use via the command line. The option <code class="docutils notranslate"><span class="pre">-eps_error_relative</span></code> will show eigenvalues whose relative error are below the tolerance. The different types of errors have their corresponding options, see table <a class="reference internal" href="#tab-errors"><span class="std std-ref">Available expressions for computing error bounds</span></a>. A more detailed output can be obtained as follows:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex1<span class="w"> </span>-n<span class="w"> </span><span class="m">120</span><span class="w"> </span>-eps_nev<span class="w"> </span><span class="m">8</span><span class="w"> </span>-eps_error_relative<span class="w"> </span>::ascii_info_detail
<span class="go"> ---------------------- --------------------</span>
<span class="go"> k ||Ax-kx||/||kx||</span>
<span class="go"> ---------------------- --------------------</span>
<span class="go"> 3.999326 1.26221e-09</span>
<span class="go"> 3.997304 3.82982e-10</span>
<span class="go"> 3.993936 2.76971e-09</span>
<span class="go"> 3.989224 4.94104e-10</span>
<span class="go"> 3.983171 6.19307e-10</span>
<span class="go"> 3.975781 5.9628e-10</span>
<span class="go"> 3.967060 2.32347e-09</span>
<span class="go"> 3.957012 6.12436e-09</span>
<span class="go"> ---------------------- --------------------</span>
</pre></div>
</div>
<p>Finally, the option for showing the converged reason is:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex1<span class="w"> </span>-n<span class="w"> </span><span class="m">120</span><span class="w"> </span>-eps_nev<span class="w"> </span><span class="m">8</span><span class="w"> </span>-eps_converged_reason
<span class="go"> Linear eigensolve converged (8 eigenpairs) due to CONVERGED_TOL; iterations 14</span>
</pre></div>
</div>
</section>
</section>
<section id="advanced-usage">
<h2>Advanced Usage<a class="headerlink" href="#advanced-usage" title="Link to this heading">#</a></h2>
<p>This section includes the description of advanced features of the eigensolver object. Default settings are appropriate for most applications and modification is unnecessary for normal usage.</p>
<section id="initial-guesses">
<h3>Initial Guesses<a class="headerlink" href="#initial-guesses" title="Link to this heading">#</a></h3>
<p>In this subsection, we consider the possibility of providing initial guesses so that the eigensolver can exploit this information to get the answer faster.</p>
<p>Most of the algorithms implemented in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> iteratively build and improve a basis of a certain subspace, which will eventually become an eigenspace corresponding to the wanted eigenvalues. In some solvers such as those of Krylov type, this basis is constructed starting from an initial vector, <span class="math notranslate nohighlight">\(v_1\)</span>, whereas in other solvers such as those of Davidson type, an arbitrary subspace can be used to start the method. By default, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> initializes the starting vector or the initial subspace randomly. This default is a reasonable choice. However, it is also possible to supply an initial subspace with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetInitialSpace.html">EPSSetInitialSpace</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">n</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">is</span><span class="p">[]);</span>
</pre></div>
</div>
<p>In some cases, a suitable initial space can accelerate convergence significantly, for instance when the eigenvalue calculation is one of a sequence of closely related problems, where the eigenspace of one problem is fed as the initial guess for the next problem.</p>
<p>Note that if the eigensolver supports only a single initial vector, but several guesses are provided, then all except the first one will be discarded. One could still build a vector that is rich in the directions of all guesses, by taking a linear combination of them, but this is less effective than using a solver that considers all guesses as a subspace.</p>
</section>
<section id="dealing-with-deflation-subspaces">
<h3>Dealing with Deflation Subspaces<a class="headerlink" href="#dealing-with-deflation-subspaces" title="Link to this heading">#</a></h3>
<p>In some applications, when solving an eigenvalue problem the user wishes to use a priori knowledge about the solution. This is the case when an invariant subspace has already been computed (e.g., in a previous <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a>()</span></code> call) or when a basis of the null-space is known.</p>
<p>Consider the following example. Given a graph <span class="math notranslate nohighlight">\(G\)</span>, with vertex set <span class="math notranslate nohighlight">\(V\)</span> and edges <span class="math notranslate nohighlight">\(E\)</span>, the Laplacian matrix of <span class="math notranslate nohighlight">\(G\)</span> is a sparse symmetric positive semidefinite matrix <span class="math notranslate nohighlight">\(L\)</span> with elements</p>
<div class="math notranslate nohighlight" id="equation-eq-laplacian">
<span class="eqno">(10)<a class="headerlink" href="#equation-eq-laplacian" title="Link to this equation">#</a></span>\[\begin{split}l_{ij}=\left\{\begin{array}{cl}
d(v_i) & \mathrm{if}\;i=j\\
-1 & \mathrm{if}\;e_{ij}\in E\\
0&\mathrm{otherwise}
\end{array}\right.\end{split}\]</div>
<p>where <span class="math notranslate nohighlight">\(d(v_i)\)</span> is the degree of vertex <span class="math notranslate nohighlight">\(v_i\)</span>. This matrix is singular since all row sums are equal to zero. The constant vector is an eigenvector with zero eigenvalue, and if the graph is connected then all other eigenvalues are positive. The so-called Fiedler vector is the eigenvector associated with the smallest nonzero eigenvalue and can be used in heuristics for a number of graph manipulations such as partitioning. One possible way of computing this vector with SLEPc is to instruct the eigensolver to search for the smallest eigenvalue (with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetWhichEigenpairs.html">EPSSetWhichEigenpairs</a>()</span></code> or by using a spectral transformation as described in next chapter) but preventing it from computing the already known eigenvalue. For this, the user must provide a basis for the invariant subspace (in this case just vector <span class="math notranslate nohighlight">\([1,1,\ldots,1]^T\)</span>) so that the eigensolver can <em>deflate</em> this subspace. This process is very similar to what eigensolvers normally do with invariant subspaces associated with eigenvalues as they converge. In other words, when a deflation space has been specified, the eigensolver works with the restriction of the problem to the orthogonal complement of this subspace.</p>
<p>The following function can be used to provide the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object with some basis vectors corresponding to a subspace that should be deflated during the solution process:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetDeflationSpace.html">EPSSetDeflationSpace</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">n</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">defl</span><span class="p">[])</span>
</pre></div>
</div>
<p>The value <code class="docutils notranslate"><span class="pre">n</span></code> indicates how many vectors are passed in argument <code class="docutils notranslate"><span class="pre">defl</span></code>.</p>
<p>The deflation space can be any subspace but typically it is most useful in the case of an invariant subspace or a null-space. In any case, SLEPc internally checks to see if all (or part of) the provided subspace is a null-space of the associated linear system (see section <a class="reference internal" href="st.html#sec-lin"><span class="std std-ref">Solution of Linear Systems</span></a>). In this case, this null-space is attached to the coefficient matrix of the linear solver (see PETSc’s function <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatSetNullSpace/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatSetNullSpace</span></a>()) to enable the solution of singular systems. In practice, this allows the computation of eigenvalues of singular pencils (i.e., when <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span> share a common null-space).</p>
</section>
<section id="sec-orthog">
<h3>Orthogonalization<a class="headerlink" href="#sec-orthog" title="Link to this heading">#</a></h3>
<p>Internally, eigensolvers in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> often need to orthogonalize a vector against a set of vectors (for instance, when building an orthonormal basis of a Krylov subspace). This operation is carried out typically by a Gram-Schmidt orthogonalization procedure. The user is able to adjust several options related to this algorithm, although the default behavior is good for most cases, and we strongly suggest not to change any of these settings. This topic is covered in detail in <span id="id7">[<a class="reference internal" href="#id68" title="V. Hernandez, J. E. Roman, A. Tomas, and V. Vidal. Orthogonalization routines in SLEPc. Technical Report STR-1, Universitat Politècnica de València, 2007. URL: https://slepc.upv.es/documentation.">Hernandez <em>et al.</em>, 2007</a>]</span>.</p>
</section>
<section id="sec-region">
<h3>Specifying a Region for Filtering<a class="headerlink" href="#sec-region" title="Link to this heading">#</a></h3>
<p>Some solvers in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> (and other solver classes as well) can take into consideration a user-defined region of the complex plane, e.g., an ellipse. A region can be used for two purposes:</p>
<ul class="simple">
<li><p>To instruct the eigensolver to compute all eigenvalues lying in that region. This is available only in eigensolvers based on the contour integral technique (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSCISS.html">EPSCISS</a></span></code>).</p></li>
<li><p>To filter out eigenvalues outside the region. In this way, eigenvalues lying inside the region get higher priority during the iteration and are more likely to be returned as computed solutions.</p></li>
</ul>
<p>Regions are specified by means of an <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/RG/RG.html">RG</a></span></code> object. This object is handled internally in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> solver, as other auxiliary objects, and can be extracted with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetRG.html">EPSGetRG</a>()</span></code> to set the options that define the region. These options can also be set in the command line. The following example computes largest magnitude eigenvalues, but restricting to an ellipse of radius 0.5 centered at the origin (with vertical scale 0.1):</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex1<span class="w"> </span>-rg_type<span class="w"> </span>ellipse<span class="w"> </span>-rg_ellipse_center<span class="w"> </span><span class="m">0</span><span class="w"> </span>-rg_ellipse_radius<span class="w"> </span><span class="m">0</span>.5<span class="w"> </span>-rg_ellipse_vscale<span class="w"> </span><span class="m">0</span>.1
</pre></div>
</div>
<p>If one wants to use the region to specify where eigenvalues should <em>not</em> be computed, then the region must be the complement of the specified one. The next command line computes the smallest eigenvalues not contained in the ellipse:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex1<span class="w"> </span>-eps_smallest_magnitude<span class="w"> </span>-rg_type<span class="w"> </span>ellipse<span class="w"> </span>-rg_ellipse_center<span class="w"> </span><span class="m">0</span><span class="w"> </span>-rg_ellipse_radius<span class="w"> </span><span class="m">0</span>.5<span class="w"> </span>-rg_ellipse_vscale<span class="w"> </span><span class="m">0</span>.1<span class="w"> </span>-rg_complement
</pre></div>
</div>
<p>Additional details of the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/RG/RG.html">RG</a></span></code> class can be found in section <a class="reference internal" href="aux.html#sec-rg"><span class="std std-ref">RG: Region</span></a>.</p>
</section>
<section id="sec-large-nev">
<h3>Computing a Large Portion of the Spectrum<a class="headerlink" href="#sec-large-nev" title="Link to this heading">#</a></h3>
<p>We now consider the case when the user requests a relatively large number of eigenpairs (the related case of computing all eigenvalues in a given interval is addressed in section <a class="reference internal" href="st.html#sec-slice"><span class="std std-ref">Spectrum Slicing</span></a>). To fix ideas, suppose that the problem size (the dimension of the matrix, denoted as <code class="docutils notranslate"><span class="pre">n</span></code>), is in the order of 100,000’s, and the user wants <code class="docutils notranslate"><span class="pre">nev</span></code> to be approximately 5,000 (recall the notation of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetDimensions.html">EPSSetDimensions</a>()</span></code> in section <a class="reference internal" href="#sec-defprob"><span class="std std-ref">Defining the Problem</span></a>).</p>
<p>The first comment is that for such large values of <code class="docutils notranslate"><span class="pre">nev</span></code>, the rule of thumb suggested in section <a class="reference internal" href="#sec-defprob"><span class="std std-ref">Defining the Problem</span></a> for selecting the value of <code class="docutils notranslate"><span class="pre">ncv</span></code> (<span class="math notranslate nohighlight">\(\mathtt{ncv}\geq2\cdot\mathtt{nev}\)</span>) may be inappropriate. For small values of <code class="docutils notranslate"><span class="pre">nev</span></code>, this rule of thumb is intended to provide the solver with a sufficiently large subspace. But for large values of <code class="docutils notranslate"><span class="pre">nev</span></code>, it may be enough setting <code class="docutils notranslate"><span class="pre">ncv</span></code> to be slightly larger than <code class="docutils notranslate"><span class="pre">nev</span></code>.</p>
<p>The second thing to take into account has to do with costs, both in terms of storage and in terms of computational effort. This issue is dependent on the particular eigensolver used, but generally speaking the user can simplify to the following points:</p>
<ol class="arabic simple">
<li><p>It is necessary to store a basis of the subspace, that is, <code class="docutils notranslate"><span class="pre">ncv</span></code> vectors of length <code class="docutils notranslate"><span class="pre">n</span></code>.</p></li>
<li><p>A considerable part of the computation is devoted to orthogonalization of the basis vectors, whose cost is roughly of order <span class="math notranslate nohighlight">\(\mathtt{ncv}^2\cdot\mathtt{n}\)</span>.</p></li>
<li><p>Within the eigensolution process, a projected eigenproblem of order <code class="docutils notranslate"><span class="pre">ncv</span></code> is built. At least one dense matrix of this dimension has to be stored.</p></li>
<li><p>Solving the projected eigenproblem has a computational cost of order <span class="math notranslate nohighlight">\(\mathtt{ncv}^3\)</span>. Typically, such problems need to be solved many times within the eigensolver iteration.</p></li>
</ol>
<p>It is clear that a large value of <code class="docutils notranslate"><span class="pre">ncv</span></code> implies a high storage requirement (points 1 and 3, especially point 1), and a high computational cost (points 2 and 4, especially point 2). However, in a scenario of such big eigenproblems, it is customary to solve the problem in parallel with many processors. In that case, it turns out that the basis vectors are stored in a distributed way and the associated operations are parallelized, so that points 1 and 2 become benign as long as sufficient processors are used. Then points 3 and 4 become really critical since in the current SLEPc version the projected eigenproblem (and its associated operations) are not treated in parallel. In conclusion, the user must be aware that using a large <code class="docutils notranslate"><span class="pre">ncv</span></code> value introduces a serial step in the computation with high cost, that cannot be amortized by increasing the number of processors.</p>
<p>From SLEPc 3.0.0, another parameter <code class="docutils notranslate"><span class="pre">mpd</span></code> has been introduced to alleviate this problem. The name <code class="docutils notranslate"><span class="pre">mpd</span></code> stands for maximum projected dimension. The idea is to bound the size of the projected eigenproblem so that steps 3 and 4 work with a dimension of <code class="docutils notranslate"><span class="pre">mpd</span></code> at most, while steps 1 and 2 still work with a bigger dimension, up to <code class="docutils notranslate"><span class="pre">ncv</span></code>. Suppose we want to compute <code class="docutils notranslate"><span class="pre">nev</span></code>=5000. Setting <code class="docutils notranslate"><span class="pre">ncv</span></code>=10000 or even <code class="docutils notranslate"><span class="pre">ncv</span></code>=6000 would be prohibitively expensive, for the reasons explained above. But if we set e.g. <code class="docutils notranslate"><span class="pre">mpd</span></code>=600 then the overhead of steps 3 and 4 will be considerably diminished. Of course, this reduces the potential of approximation at each outer iteration of the algorithm, but with more iterations the same result should be obtained. The benefits will be specially noticeable in the setting of parallel computation with many processors.</p>
<p>Note that it is not necessary to set both <code class="docutils notranslate"><span class="pre">ncv</span></code> and <code class="docutils notranslate"><span class="pre">mpd</span></code>. For instance, one can do</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./program<span class="w"> </span>-eps_nev<span class="w"> </span><span class="m">5000</span><span class="w"> </span>-eps_mpd<span class="w"> </span><span class="m">600</span>
</pre></div>
</div>
</section>
<section id="sec-harmonic">
<h3>Computing Interior Eigenvalues with Harmonic Extraction<a class="headerlink" href="#sec-harmonic" title="Link to this heading">#</a></h3>
<p>The standard Rayleigh-Ritz projection procedure described in section <a class="reference internal" href="#sec-eig"><span class="std std-ref">Eigenvalue Problems</span></a> is most appropriate for approximating eigenvalues located at the periphery of the spectrum, especially those of largest magnitude. Most eigensolvers in SLEPc are restarted, meaning that the projection is carried out repeatedly with increasingly good subspaces. An effective restarting mechanism, such as that implemented in Krylov-Schur, improves the subspace by realizing a filtering effect that tries to eliminate components in the direction of unwanted eigenvectors. In that way, it is possible to compute eigenvalues located anywhere in the spectrum, even in its interior.</p>
<p>Even though in theory eigensolvers could be able to approximate interior eigenvalues with a standard extraction technique, in practice convergence difficulties may arise that prevent success. The problem comes from the property that Ritz values (the approximate eigenvalues provided by the standard projection procedure) converge from the interior to the periphery of the spectrum. That is, the Ritz values that stabilize first are those in the periphery, so convergence of interior ones requires the previous convergence of all eigenvalues between them and the periphery. Furthermore, this convergence behavior usually implies that restarting is carried out with bad approximations, so the restart is ineffective and global convergence is severely damaged.</p>
<p>Harmonic projection is a variation that uses a target value, <span class="math notranslate nohighlight">\(\tau\)</span>, around which the user wants to compute eigenvalues (see, e.g., <span id="id8">[<a class="reference internal" href="#id34" title="R. B. Morgan and M. Zeng. A harmonic restarted Arnoldi algorithm for calculating eigenvalues and determining multiplicity. Linear Algebra Appl., 415(1):96–113, 2006. doi:10.1016/j.laa.2005.07.024.">Morgan and Zeng, 2006</a>]</span>). The theory establishes that harmonic Ritz values converge in such a way that eigenvalues closest to the target stabilize first, and also that no unconverged value is ever close to the target, so restarting is safe in this case. As a conclusion, eigensolvers with harmonic extraction may be effective in computing interior eigenvalues. Whether it works or not in practical cases depends on the particular distribution of the spectrum.</p>
<p>In order to use harmonic extraction in SLEPc, it is necessary to indicate it explicitly, and provide the target value as described in section <a class="reference internal" href="#sec-defprob"><span class="std std-ref">Defining the Problem</span></a> (default is <span class="math notranslate nohighlight">\(\tau=0\)</span>). The type of extraction can be set with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetExtraction.html">EPSSetExtraction</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="../../manualpages/EPS/EPSExtraction.html">EPSExtraction</a></span><span class="w"> </span><span class="n">extr</span><span class="p">);</span>
</pre></div>
</div>
<p>Available possibilities are <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_RITZ.html">EPS_RITZ</a></span></code> for standard projection and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HARMONIC.html">EPS_HARMONIC</a></span></code> for harmonic projection (other alternatives such as refined extraction are still experimental).</p>
<p>A command line example would be:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex5<span class="w"> </span>-m<span class="w"> </span><span class="m">45</span><span class="w"> </span>-eps_harmonic<span class="w"> </span>-eps_target<span class="w"> </span><span class="m">0</span>.8<span class="w"> </span>-eps_ncv<span class="w"> </span><span class="m">60</span>
</pre></div>
</div>
<p>The example computes the eigenvalue closest to <span class="math notranslate nohighlight">\(\tau=0.8\)</span> of a real non-symmetric matrix of order 1035. Note that <code class="docutils notranslate"><span class="pre">ncv</span></code> has been set to a larger value than would be necessary for computing the largest magnitude eigenvalues. In general, users should expect a much slower convergence when computing interior eigenvalues compared to extreme eigenvalues. Increasing the value of <code class="docutils notranslate"><span class="pre">ncv</span></code> may help.</p>
<p>Currently, harmonic extraction is available in the default <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> solver, that is, Krylov-Schur, and also in Arnoldi, GD, and JD.</p>
</section>
<section id="sec-balancing">
<h3>Balancing for Non-Hermitian Problems<a class="headerlink" href="#sec-balancing" title="Link to this heading">#</a></h3>
<p>In problems where the matrix has a large norm, <span class="math notranslate nohighlight">\(\|A\|_2\)</span>, the roundoff errors introduced by the eigensolver may be large. The goal of balancing is to apply a simple similarity transformation, <span class="math notranslate nohighlight">\(DAD^{-1}\)</span>, that keeps the eigenvalues unaltered but reduces the matrix norm, thus enhancing the accuracy of the computed eigenpairs. Obviously, this makes sense only in the non-Hermitian case. The matrix <span class="math notranslate nohighlight">\(D\)</span> is chosen to be diagonal, so balancing amounts to scaling the matrix rows and columns appropriately.</p>
<p>In SLEPc, the matrix <span class="math notranslate nohighlight">\(DAD^{-1}\)</span> is not formed explicitly. Instead, the operators of table <a class="reference internal" href="st.html#tab-op"><span class="std std-ref">Operators used in each spectral transformation mode</span></a> are preceded by a multiplication by <span class="math notranslate nohighlight">\(D^{-1}\)</span> and followed by a multiplication by <span class="math notranslate nohighlight">\(D\)</span>. This allows for balancing in the case of problems with an implicit matrix.</p>
<p>A simple and effective Krylov balancing technique, described in <span id="id9">[<a class="reference internal" href="../../manualpages/EPS/EPSSetBalance.html#id26" title="T.-Y. Chen and J. W. Demmel. Balancing sparse matrices for computing eigenvalues. Linear Algebra Appl., 309(1–3):261–287, 2000. doi:10.1016/S0024-3795(00)00014-8.">Chen and Demmel, 2000</a>]</span>, is implemented in SLEPc. The user calls the following subroutine to activate it:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetBalance.html">EPSSetBalance</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="../../manualpages/EPS/EPSBalance.html">EPSBalance</a></span><span class="w"> </span><span class="n">bal</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">its</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="n">cutoff</span><span class="p">);</span>
</pre></div>
</div>
<p>Two variants are available, one-sided and two-sided, and there is also the possibility for the user to provide a pre-computed <span class="math notranslate nohighlight">\(D\)</span> matrix.</p>
</section>
<section id="sec-structured">
<h3>Structured Eigenvalue Problems<a class="headerlink" href="#sec-structured" title="Link to this heading">#</a></h3>
<p>Structured eigenvalue problems are those whose defining matrices are structured, i.e., their <span class="math notranslate nohighlight">\(n^2\)</span> entries depend on less than <span class="math notranslate nohighlight">\(n^2\)</span> parameters. Symmetry is the most obvious structure, and it is supported in SLEPc solvers via the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code> and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code> problem types, see table <a class="reference internal" href="#tab-ptype"><span class="std std-ref">Problem types considered in EPS</span></a>. The last entries listed in that table address other types of structured eigenproblems, which are discussed in this subsection. Preserving the algebraic structure can help preserve physically relevant symmetries in the eigenvalues of the matrix and may improve the accuracy and efficiency of the eigensolver. For example, in quadratic eigenvalue problems arising from gyroscopic systems (see section <a class="reference internal" href="pep.html#sec-qep"><span class="std std-ref">Quadratic Eigenvalue Problems</span></a>), eigenvalues appear in quadruples <span class="math notranslate nohighlight">\(\{\lambda,-\lambda,\bar\lambda,-\bar\lambda\}\)</span>, i.e., the spectrum is symmetric with respect to both the real and imaginary axes. This problem can be linearized to a <span class="math notranslate nohighlight">\(2n\times 2n\)</span> skew-Hamiltonian/Hamiltonian pencil with the same eigenvalues. A structure-preserving eigensolver will give a more accurate answer because it enforces the structure throughout the computation.</p>
<p>Unless otherwise stated, the structured eigenproblems discussed below are only supported in the default <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> solver, Krylov-Schur. The idea is that the user creates the structured matrix with a helper function such as <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/Sys/MatCreateBSE.html">MatCreateBSE</a>()</span></code> (see below), and then selects the appropriate problem type with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetProblemType.html">EPSSetProblemType</a>()</span></code>, in this case <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_BSE.html">EPS_BSE</a></span></code>. This will instruct the solver to exploit the problem structure. Alternatively, one can solve the problem as <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_NHEP.html">EPS_NHEP</a></span></code>, in which case the solver will neglect the structure.</p>
<div class="admonition warning">
<p class="admonition-title">Warning</p>
<p>Check below the section <a class="reference internal" href="#sec-structured-vectors"><span class="std std-ref">Extracting Eigenvectors of Structured Eigenproblems</span></a>, since the trivial approach may give vectors with disordered entries.</p>
</div>
<section id="bethe-salpeter">
<h4>Bethe-Salpeter<a class="headerlink" href="#bethe-salpeter" title="Link to this heading">#</a></h4>
<p>One structured eigenproblem that has raised interest recently is related to the Bethe–Salpeter equation, which is relevant in many state-of-the-art computational physics codes. For instance, in the Yambo software <span id="id10">[<a class="reference internal" href="#id57" title="D. Sangalli, A. Ferretti, H. Miranda, C. Attaccalite, I. Marri, E. Cannuccia, P. Melo, M. Marsili, F. Paleari, A. Marrazzo, G. Prandini, P. Bonfà, M. O. Atambo, F. Affinito, M. Palummo, A. Molina-Sánchez, C. Hogan, M. Grüning, D. Varsano, and A. Marini. Many-body perturbation theory calculations using the yambo code. J. Condens. Matter Phys., 31(32):325902, 2019. doi:10.1088/1361-648x/ab15d0.">Sangalli <em>et al.</em>, 2019</a>]</span> it is used to evaluate optical properties of materials. It is commonly formulated as an eigenvalue problem with a block-structured matrix,</p>
<div class="math notranslate nohighlight" id="equation-eq-bse">
<span class="eqno">(11)<a class="headerlink" href="#equation-eq-bse" title="Link to this equation">#</a></span>\[\begin{split}H = \begin{bmatrix}
R & C \\
-C^* & -R^T
\end{bmatrix},\end{split}\]</div>
<p>where <span class="math notranslate nohighlight">\(R\)</span> is Hermitian and <span class="math notranslate nohighlight">\(C\)</span> is complex symmetric. For a good enough approximation of the optical absorption spectrum, it is sufficient to compute a few eigenvalues with a customized version of Krylov-Schur. In this problem, eigenvalues are real and come in pairs <span class="math notranslate nohighlight">\(\{\lambda,-\lambda\}\)</span>. The eigenvalues of interest are those with the smallest magnitude, which in this case lie in the middle of the spectrum. Usually, both right and left eigenvectors are required, but the left eigenvectors can be obtained inexpensively once the corresponding right ones are known.</p>
<p>The helper function to generate the matrix <span class="math notranslate nohighlight">\(H\)</span> of equation <a class="reference internal" href="#equation-eq-bse">(11)</a> from the blocks <span class="math notranslate nohighlight">\(R\)</span> and <span class="math notranslate nohighlight">\(C\)</span> is <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/Sys/MatCreateBSE.html">MatCreateBSE</a>()</span></code>, and the associated problem type is <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_BSE.html">EPS_BSE</a></span></code> (or <code class="docutils notranslate"><span class="pre">-eps_bse</span></code> from the command line). It is possible to select a few variants of the solver with the function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSKrylovSchurSetBSEType.html">EPSKrylovSchurSetBSEType</a>()</span></code>.</p>
<p>Further details about the implementation of the SLEPc solvers for the BSE can be found in <span id="id11">[<a class="reference internal" href="../../manualpages/EPS/EPSKrylovSchurBSEType.html#id48" title="F. Alvarruiz, B. Mellado-Pinto, and J. E. Roman. Variants of thick-restart Lanczos for the Bethe-Salpeter eigenvalue problem. arXiv:2503.20920 : retrieved 27 May 2025, 2025. doi:10.48550/arXiv.2503.20920.">Alvarruiz <em>et al.</em>, 2025</a>]</span>.</p>
</section>
<section id="hamiltonian">
<h4>Hamiltonian<a class="headerlink" href="#hamiltonian" title="Link to this heading">#</a></h4>
<p>The Hamiltonian structure is relevant in many applications, particularly in control theory. A (complex) Hamiltonian matrix has the block structure</p>
<div class="math notranslate nohighlight" id="equation-eq-hamilt">
<span class="eqno">(12)<a class="headerlink" href="#equation-eq-hamilt" title="Link to this equation">#</a></span>\[\begin{split}H = \begin{bmatrix}
A & B \\
C & -A^*
\end{bmatrix},\end{split}\]</div>
<p>where <span class="math notranslate nohighlight">\(A\)</span>, <span class="math notranslate nohighlight">\(B\)</span> and <span class="math notranslate nohighlight">\(C\)</span> are either real with <span class="math notranslate nohighlight">\(B=B^T\)</span>, <span class="math notranslate nohighlight">\(C=C^T\)</span>, or complex with <span class="math notranslate nohighlight">\(B=B^*\)</span>, <span class="math notranslate nohighlight">\(C=C^*\)</span>. In the real case, eigenvalues appear in pairs <span class="math notranslate nohighlight">\(\{\lambda,-\lambda\}\)</span> and for complex eigenvalues in quadruples <span class="math notranslate nohighlight">\(\{\lambda,-\lambda,\bar\lambda,-\bar\lambda\}\)</span>. For a complex Hamiltonian matrix, if <span class="math notranslate nohighlight">\(\lambda\)</span> is an eigenvalue, then <span class="math notranslate nohighlight">\(-\bar\lambda\)</span> is also an eigenvalue. A structure-preserving eigensolver has been implemented in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> (in Krylov-Schur), which is activated by selecting the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HAMILT.html">EPS_HAMILT</a></span></code> problem type (or <code class="docutils notranslate"><span class="pre">-eps_hamiltonian</span></code> from the command line). Note that matrix <span class="math notranslate nohighlight">\(H\)</span> <a class="reference internal" href="#equation-eq-hamilt">(12)</a> must be created with the helper function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/Sys/MatCreateHamiltonian.html">MatCreateHamiltonian</a>()</span></code> in order to use this solver.</p>
<div class="admonition warning">
<p class="admonition-title">Warning</p>
<p>The structure-preserving eigensolver for Hamiltonian eigenvalue problems should be considered experimental. Depending on the problem, it may become numerically unstable after some iterations, in which case the solver will abort, returning less eigenvalues than requested.</p>
</div>
</section>
<section id="sec-structured-vectors">
<h4>Extracting Eigenvectors of Structured Eigenproblems<a class="headerlink" href="#sec-structured-vectors" title="Link to this heading">#</a></h4>
<p>In the case of structured eigenproblems, one has to create the eigenvectors in a specific way, before passing them to functions such as <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetEigenvector.html">EPSGetEigenvector</a>()</span></code> or <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetLeftEigenvector.html">EPSGetLeftEigenvector</a>()</span></code>. Otherwise, the obtained vector may have disordered entries when run with several processes, compared to sequential runs.</p>
<p>To understand the problem, it is important to note that the block matrices <a class="reference internal" href="#equation-eq-bse">(11)</a> or <a class="reference internal" href="#equation-eq-hamilt">(12)</a> are represented as PETSc matrices of type <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MATNEST/" title="(in PETSc v3.24)"><span class="xref std std-doc">MATNEST</span></a>, which implies that the individual blocks are stored independently (in parallel). The vectors associated with those matrices have the same distribution. For instance, if we are using two MPI processes, then process 0 will store the first half of the top part of the vector, and the first half of the bottom part, while process 1 will store the second half of both parts. Hence, the entries assigned to the different processes are not consecutive, but interleaved, so if we access that vector as a standard <a class="reference external" href="https://petsc.org/release/manualpages/Vec/Vec/" title="(in PETSc v3.24)"><span class="xref std std-doc">Vec</span></a> the order of entries will be different when we change the number of processes.</p>
<p>The solution is to work with the top and bottom subvectors separately. There are two possible approaches for this, as explained next.</p>
<p>The first solution is to create a vector of type <a class="reference external" href="https://petsc.org/release/manualpages/Vec/VECNEST/" title="(in PETSc v3.24)"><span class="xref std std-doc">VECNEST</span></a> compatible with the <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MATNEST/" title="(in PETSc v3.24)"><span class="xref std std-doc">MATNEST</span></a>. For instance in BSE, <a class="reference internal" href="#equation-eq-bse">(11)</a>, we can create it from the <span class="math notranslate nohighlight">\(R\)</span> and <span class="math notranslate nohighlight">\(C\)</span> blocks, as follows (assuming both matrices have the same local and global sizes):</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="w"> </span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/MatCreateVecs/">MatCreateVecs</a></span><span class="p">(</span><span class="n">R</span><span class="p">,</span><span class="o">&</span><span class="n">x1</span><span class="p">,</span><span class="o">&</span><span class="n">x2</span><span class="p">);</span>
<span class="w"> </span><span class="n">aux</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">x1</span><span class="p">;</span>
<span class="w"> </span><span class="n">aux</span><span class="p">[</span><span class="mi">1</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">x2</span><span class="p">;</span>
<span class="w"> </span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/VecCreateNest/">VecCreateNest</a></span><span class="p">(</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PETSC_COMM_WORLD/">PETSC_COMM_WORLD</a></span><span class="p">,</span><span class="mi">2</span><span class="p">,</span><span class="nb">NULL</span><span class="p">,</span><span class="n">aux</span><span class="p">,</span><span class="o">&</span><span class="n">x</span><span class="p">);</span>
</pre></div>
</div>
<p>We would pass <code class="docutils notranslate"><span class="pre">x</span></code> to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetEigenvector.html">EPSGetEigenvector</a>()</span></code>, and then the individual blocks of the eigenvector can be easily extracted with <a class="reference external" href="https://petsc.org/release/manualpages/Vec/VecNestGetSubVec/" title="(in PETSc v3.24)"><span class="xref std std-doc">VecNestGetSubVec</span></a>().</p>
<p>The second alternative is to retrieve the index sets <a class="reference external" href="https://petsc.org/release/manualpages/IS/IS/" title="(in PETSc v3.24)"><span class="xref std std-doc">IS</span></a> that define the splitting</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="w"> </span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/MatNestGetISs/">MatNestGetISs</a></span><span class="p">(</span><span class="n">H</span><span class="p">,</span><span class="n">is</span><span class="p">,</span><span class="nb">NULL</span><span class="p">);</span>
</pre></div>
</div>
<p>and then, after calling <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSGetEigenvector.html">EPSGetEigenvector</a>()</span></code> on a standard vector <code class="docutils notranslate"><span class="pre">x</span></code>, extract the subvectors using the index sets:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="w"> </span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/VecGetSubVector/">VecGetSubVector</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">is</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span><span class="o">&</span><span class="n">x1</span><span class="p">);</span>
<span class="w"> </span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/VecGetSubVector/">VecGetSubVector</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">is</span><span class="p">[</span><span class="mi">1</span><span class="p">],</span><span class="o">&</span><span class="n">x2</span><span class="p">);</span>
<span class="w"> </span><span class="cm">/* do something with x1 and x2 */</span>
<span class="w"> </span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/VecRestoreSubVector/">VecRestoreSubVector</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">is</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span><span class="o">&</span><span class="n">x1</span><span class="p">);</span>
<span class="w"> </span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/VecRestoreSubVector/">VecRestoreSubVector</a></span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">is</span><span class="p">[</span><span class="mi">1</span><span class="p">],</span><span class="o">&</span><span class="n">x2</span><span class="p">);</span>
</pre></div>
</div>
<p class="rubric">References</p>
<div class="docutils container" id="id12">
<div role="list" class="citation-list">
<div class="citation" id="id58" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id11">Alv25</a><span class="fn-bracket">]</span></span>
<p>F. Alvarruiz, B. Mellado-Pinto, and J. E. Roman. Variants of thick-restart Lanczos for the Bethe-Salpeter eigenvalue problem. <em>arXiv:2503.20920 : retrieved 27 May 2025</em>, 2025. <a class="reference external" href="https://doi.org/10.48550/arXiv.2503.20920">doi:10.48550/arXiv.2503.20920</a>.</p>
</div>
<div class="citation" id="id14" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id1">Bai00</a><span class="fn-bracket">]</span></span>
<p>Z. Bai, J. Demmel, J. Dongarra, A. Ruhe, and H. van der Vorst, editors. <em>Templates for the Solution of Algebraic Eigenvalue Problems: A Practical Guide</em>. Society for Industrial and Applied Mathematics, Philadelphia, PA, 2000.</p>
</div>
<div class="citation" id="id36" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id9">Che00</a><span class="fn-bracket">]</span></span>
<p>T.-Y. Chen and J. W. Demmel. Balancing sparse matrices for computing eigenvalues. <em>Linear Algebra Appl.</em>, 309(1–3):261–287, 2000. <a class="reference external" href="https://doi.org/10.1016/S0024-3795(00)00014-8">doi:10.1016/S0024-3795(00)00014-8</a>.</p>
</div>
<div class="citation" id="id17" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id2">Gol00</a><span class="fn-bracket">]</span></span>
<p>G. H. Golub and H. A. van der Vorst. Eigenvalue computation in the 20th century. <em>J. Comput. Appl. Math.</em>, 123(1-2):35–65, 2000. <a class="reference external" href="https://doi.org/10.1016/S0377-0427(00)00413-1">doi:10.1016/S0377-0427(00)00413-1</a>.</p>
</div>
<div class="citation" id="id68" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id7">Her07a</a><span class="fn-bracket">]</span></span>
<p>V. Hernandez, J. E. Roman, A. Tomas, and V. Vidal. Orthogonalization routines in SLEPc. Technical Report STR-1, Universitat Politècnica de València, 2007. URL: <a class="reference external" href="https://slepc.upv.es/documentation">https://slepc.upv.es/documentation</a>.</p>
</div>
<div class="citation" id="id76" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id5">Mae16</a><span class="fn-bracket">]</span></span>
<p>Y. Maeda, T. Sakurai, and J. E. Roman. Contour integral spectrum slicing method in SLEPc. Technical Report STR-11, Universitat Politècnica de València, 2016. URL: <a class="reference external" href="https://slepc.upv.es/documentation">https://slepc.upv.es/documentation</a>.</p>
</div>
<div class="citation" id="id34" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id8">Mor06</a><span class="fn-bracket">]</span></span>
<p>R. B. Morgan and M. Zeng. A harmonic restarted Arnoldi algorithm for calculating eigenvalues and determining multiplicity. <em>Linear Algebra Appl.</em>, 415(1):96–113, 2006. <a class="reference external" href="https://doi.org/10.1016/j.laa.2005.07.024">doi:10.1016/j.laa.2005.07.024</a>.</p>
</div>
<div class="citation" id="id20" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id1">Par80</a><span class="fn-bracket">]</span></span>
<p>B. N. Parlett. <em>The Symmetric Eigenvalue Problem</em>. Prentice-Hall, Englewood Cliffs, NJ, 1980. Reissued with revisions by SIAM, Philadelphia, 1998.</p>
</div>
<div class="citation" id="id23" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span>Saa92<span class="fn-bracket">]</span></span>
<span class="backrefs">(<a role="doc-backlink" href="#id1">1</a>,<a role="doc-backlink" href="#id6">2</a>)</span>
<p>Y. Saad. <em>Numerical Methods for Large Eigenvalue Problems: Theory and Algorithms</em>. John Wiley and Sons, New York, 1992.</p>
</div>
<div class="citation" id="id57" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id10">San19</a><span class="fn-bracket">]</span></span>
<p>D. Sangalli, A. Ferretti, H. Miranda, C. Attaccalite, I. Marri, E. Cannuccia, P. Melo, M. Marsili, F. Paleari, A. Marrazzo, G. Prandini, P. Bonfà, M. O. Atambo, F. Affinito, M. Palummo, A. Molina-Sánchez, C. Hogan, M. Grüning, D. Varsano, and A. Marini. Many-body perturbation theory calculations using the yambo code. <em>J. Condens. Matter Phys.</em>, 31(32):325902, 2019. <a class="reference external" href="https://doi.org/10.1088/1361-648x/ab15d0">doi:10.1088/1361-648x/ab15d0</a>.</p>
</div>
<div class="citation" id="id24" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id1">Ste01a</a><span class="fn-bracket">]</span></span>
<p>G. W. Stewart. <em>Matrix Algorithms. Volume II: Eigensystems</em>. Society for Industrial and Applied Mathematics, Philadelphia, PA, USA, 2001.</p>
</div>
</div>
</div>
<p class="rubric">Footnotes</p>
</section>
</section>
</section>
</section>
<hr class="footnotes docutils" />
<aside class="footnote-list brackets">
<aside class="footnote brackets" id="eps-real" role="doc-footnote">
<span class="label"><span class="fn-bracket">[</span>1<span class="fn-bracket">]</span></span>
<span class="backrefs">(<a role="doc-backlink" href="#id3">1</a>,<a role="doc-backlink" href="#id4">2</a>)</span>
<p>If SLEPc is compiled for real scalars, then the absolute value of the imaginary part, <span class="math notranslate nohighlight">\(|\mathrm{Im}(\lambda)|\)</span>, is used for eigenvalue selection and sorting.</p>
</aside>
</aside>
</article>
<footer class="prev-next-footer d-print-none">
<div class="prev-next-area">
<a class="left-prev"
href="intro.html"
title="previous page">
<i class="fa-solid fa-angle-left"></i>
<div class="prev-next-info">
<p class="prev-next-subtitle">previous</p>
<p class="prev-next-title">Getting Started</p>
</div>
</a>
<a class="right-next"
href="st.html"
title="next page">
<div class="prev-next-info">
<p class="prev-next-subtitle">next</p>
<p class="prev-next-title">ST: Spectral Transformation</p>
</div>
<i class="fa-solid fa-angle-right"></i>
</a>
</div>
</footer>
</div>
<dialog id="pst-secondary-sidebar-modal"></dialog>
<div id="pst-secondary-sidebar" class="bd-sidebar-secondary bd-toc"><div class="sidebar-secondary-items sidebar-secondary__inner">
<div class="sidebar-secondary-item">
<div
id="pst-page-navigation-heading-2"
class="page-toc tocsection onthispage">
<i class="fa-solid fa-list"></i> On this page
</div>
<nav class="bd-toc-nav page-toc" aria-labelledby="pst-page-navigation-heading-2">
<ul class="visible nav section-nav flex-column">
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-eig">Eigenvalue Problems</a><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#projection-methods">Projection Methods</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#preconditioned-eigensolvers">Preconditioned Eigensolvers</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#related-problems">Related Problems</a></li>
</ul>
</li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#basic-usage">Basic Usage</a></li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-defprob">Defining the Problem</a><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-which">Eigenvalues of Interest</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#left-eigenvectors">Left Eigenvectors</a></li>
</ul>
</li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#selecting-the-eigensolver">Selecting the Eigensolver</a></li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-retrsol">Retrieving the Solution</a><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#the-computed-solution">The Computed Solution</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#real-vs-complex-slepc">Real vs Complex SLEPc</a></li>
</ul>
</li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-errbnd">Reliability of the Computed Solution</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#computation-of-bounds">Computation of Bounds</a></li>
</ul>
</li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-monitor">Controlling and Monitoring Convergence</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#convergence-check">Convergence Check</a></li>
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#monitors">Monitors</a></li>
</ul>
</li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-epsviewers">Viewing the Solution</a></li>
</ul>
</li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#advanced-usage">Advanced Usage</a><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#initial-guesses">Initial Guesses</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#dealing-with-deflation-subspaces">Dealing with Deflation Subspaces</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-orthog">Orthogonalization</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-region">Specifying a Region for Filtering</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-large-nev">Computing a Large Portion of the Spectrum</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-harmonic">Computing Interior Eigenvalues with Harmonic Extraction</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-balancing">Balancing for Non-Hermitian Problems</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-structured">Structured Eigenvalue Problems</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#bethe-salpeter">Bethe-Salpeter</a></li>
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#hamiltonian">Hamiltonian</a></li>
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-structured-vectors">Extracting Eigenvectors of Structured Eigenproblems</a></li>
</ul>
</li>
</ul>
</li>
</ul>
</nav></div>
<div class="sidebar-secondary-item">
<div class="tocsection editthispage">
<a href="https://gitlab.com/slepc/slepc/-/edit/release/doc/source/documentation/manual/eps.md">
<i class="fa-solid fa-pencil"></i>
Edit on GitLab
</a>
</div>
</div>
<div class="sidebar-secondary-item">
<div role="note" aria-label="source link">
<h3>This Page</h3>
<ul class="this-page-menu">
<li><a href="../../_sources/documentation/manual/eps.md.txt"
rel="nofollow">Show Source</a></li>
</ul>
</div></div>
</div></div>
</div>
<footer class="bd-footer-content">
</footer>
</main>
</div>
</div>
<!-- Scripts loaded after <body> so the DOM is not blocked -->
<script defer src="../../_static/scripts/bootstrap.js?digest=8878045cc6db502f8baf"></script>
<script defer src="../../_static/scripts/pydata-sphinx-theme.js?digest=8878045cc6db502f8baf"></script>
<footer class="bd-footer">
<div class="bd-footer__inner bd-page-width">
<div class="footer-items__start">
<div class="footer-item">
<p class="copyright">
© Copyright 2002-2025, Universitat Politecnica de Valencia, Spain.
<br/>
</p>
</div>
<div class="footer-item">
<p class="sphinx-version">
Created using <a href="https://www.sphinx-doc.org/">Sphinx</a> 7.3.7.
<br/>
</p>
</div>
</div>
<div class="footer-items__end">
<div class="footer-item">
<p class="theme-version">
<!-- # L10n: Setting the PST URL as an argument as this does not need to be localized -->
Built with the <a href="https://pydata-sphinx-theme.readthedocs.io/en/stable/index.html">PyData Sphinx Theme</a> 0.16.1.
</p></div>
</div>
</div>
</footer>
</body>
</html>
|