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
|
<!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>NEP: Nonlinear Eigenvalue Problems — 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/nep';</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="MFN: Matrix Function" href="mfn.html" />
<link rel="prev" title="PEP: Polynomial Eigenvalue Problems" href="pep.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"><a class="reference internal" href="eps.html">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 current active"><a class="current reference internal" href="#">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">NEP: Nonlinear Eigenvalue Problems</span></li>
</ul>
</nav>
</div>
</div>
</div>
</div>
<div id="searchbox"></div>
<article class="bd-article">
<section class="tex2jax_ignore mathjax_ignore" id="nep-nonlinear-eigenvalue-problems">
<span id="ch-nep"></span><h1>NEP: Nonlinear Eigenvalue Problems<a class="headerlink" href="#nep-nonlinear-eigenvalue-problems" title="Link to this heading">#</a></h1>
<p>The Nonlinear Eigenvalue Problem (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code>) solver object covers the general case where the eigenproblem is nonlinear with respect to the eigenvalue, but it cannot be expressed in terms of a polynomial. We will write the problem as <span class="math notranslate nohighlight">\(T(\lambda)x=0\)</span>, where <span class="math notranslate nohighlight">\(T\)</span> is a matrix-valued function of the eigenvalue <span class="math notranslate nohighlight">\(\lambda\)</span>. Note that <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> does not cover the even more general case of having a nonlinear dependence on the eigenvector <span class="math notranslate nohighlight">\(x\)</span>.</p>
<p>In terms of the user interface, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> is quite similar to previously discussed solvers. The main difference is how to represent the function <span class="math notranslate nohighlight">\(T\)</span>. We will show different alternatives for this.</p>
<p>The <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> module of SLEPc has been explained with more detail in <span id="id1">[<a class="reference internal" href="../../manualpages/NEP/NEPNARNOLDI.html#id45" title="C. Campos and J. E. Roman. NEP: a module for the parallel solution of nonlinear eigenvalue problems in SLEPc. ACM Trans. Math. Software, 47(3):23:1–23:29, 2021. doi:10.1145/3447544.">Campos and Roman, 2021</a>]</span>, including an algorithmic description of the implemented solvers.</p>
<section id="sec-nep">
<h2>General Nonlinear Eigenproblems<a class="headerlink" href="#sec-nep" title="Link to this heading">#</a></h2>
<p>As in previous chapters, we first set up the notation and briefly review basic properties of the eigenvalue problems to be addressed. In this case, we focus on general nonlinear eigenproblems, that is, those that cannot be expressed in a simpler form such as a polynomial eigenproblem. These problems arise in many applications, such as the discretization of delay differential equations. Some of the methods designed to solve such problems are based on Newton-type iterations, so in some ways <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> has similarities to PETSc’s nonlinear solvers <a class="reference external" href="https://petsc.org/release/manualpages/SNES/SNES/" title="(in PETSc v3.24)"><span class="xref std std-doc">SNES</span></a>. For background material on the nonlinear eigenproblem, the reader is referred to <span id="id2">[<a class="reference internal" href="#id43" title="S. Güttel and F. Tisseur. The nonlinear eigenvalue problem. Acta Numerica, 26:1–94, 2017. doi:10.1017/S0962492917000034.">Güttel and Tisseur, 2017</a>, <a class="reference internal" href="#id33" title="V. Mehrmann and H. Voss. Nonlinear eigenvalue problems: a challenge for modern eigenvalue methods. GAMM Mitt., 27(2):121–152, 2004. doi:10.1002/gamm.201490007.">Mehrmann and Voss, 2004</a>]</span>.</p>
<p>We consider nonlinear eigenvalue problems of the form</p>
<div class="math notranslate nohighlight" id="equation-eq-nep">
<span class="eqno">(1)<a class="headerlink" href="#equation-eq-nep" title="Link to this equation">#</a></span>\[T(\lambda)x=0,\qquad x\neq 0,\]</div>
<p>where <span class="math notranslate nohighlight">\(T:\Omega\rightarrow\mathbb{C}^{n\times n}\)</span> is a matrix-valued function that is analytic on an open set of the complex plane <span class="math notranslate nohighlight">\(\Omega\subseteq\mathbb{C}\)</span>. Assuming that the problem is regular, that is, <span class="math notranslate nohighlight">\(\det T(\lambda)\)</span> does not vanish identically, any pair <span class="math notranslate nohighlight">\((\lambda,x)\)</span> satisfying equation <a class="reference internal" href="#equation-eq-nep">(1)</a> is an eigenpair, where <span class="math notranslate nohighlight">\(\lambda\in\Omega\)</span> is the eigenvalue and <span class="math notranslate nohighlight">\(x\in\mathbb{C}^n\)</span> is the eigenvector. Linear and polynomial eigenproblems are particular cases of equation <a class="reference internal" href="#equation-eq-nep">(1)</a>.</p>
<p>An example application is the rational eigenvalue problem</p>
<div class="math notranslate nohighlight" id="equation-eq-rep">
<span class="eqno">(2)<a class="headerlink" href="#equation-eq-rep" title="Link to this equation">#</a></span>\[-Kx+\lambda Mx+\sum_{j=1}^k\frac{\lambda}{\sigma_j-\lambda}C_jx=0,\]</div>
<p>arising in the study of free vibration of plates with elastically attached masses. Here, all matrices are symmetric, <span class="math notranslate nohighlight">\(K\)</span> and <span class="math notranslate nohighlight">\(M\)</span> are positive-definite and <span class="math notranslate nohighlight">\(C_j\)</span> have small rank. Another example comes from the discretization of parabolic partial differential equations with time delay <span class="math notranslate nohighlight">\(\tau\)</span>, resulting in</p>
<div class="math notranslate nohighlight" id="equation-eq-delay">
<span class="eqno">(3)<a class="headerlink" href="#equation-eq-delay" title="Link to this equation">#</a></span>\[(-\lambda I + A + e^{-\tau\lambda}B)x = 0.\]</div>
<p><strong>Split Form</strong>:
Equation <a class="reference internal" href="#equation-eq-nep">(1)</a> can always be rewritten as</p>
<div class="math notranslate nohighlight" id="equation-eq-split">
<span class="eqno">(4)<a class="headerlink" href="#equation-eq-split" title="Link to this equation">#</a></span>\[\big(A_0f_0(\lambda)+A_1f_1(\lambda)+\cdots+A_{\ell-1}f_{\ell-1}(\lambda)\big)x=
\left(\sum_{i=0}^{\ell-1}A_if_i(\lambda)\right)x = 0,\]</div>
<p>where <span class="math notranslate nohighlight">\(A_i\)</span> are <span class="math notranslate nohighlight">\(n\times n\)</span> matrices and <span class="math notranslate nohighlight">\(f_i:\Omega\rightarrow\mathbb{C}\)</span> are analytic functions. We will call equation <a class="reference internal" href="#equation-eq-split">(4)</a> the split form of the nonlinear eigenvalue problem. Often, the formulation arising from applications already has this form, as illustrated by the examples above. Also, a polynomial eigenvalue problem fits this form, where in this case the <span class="math notranslate nohighlight">\(f_i\)</span> functions are the polynomial bases of degree <span class="math notranslate nohighlight">\(i\)</span>, either monomial or non-monomial.</p>
</section>
<section id="defining-the-problem">
<h2>Defining the Problem<a class="headerlink" href="#defining-the-problem" title="Link to this heading">#</a></h2>
<p>The user interface of the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> package is quite similar to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/PEP/PEP.html">PEP</a></span></code>. As mentioned above, the main difference is the way in which the eigenproblem is defined. In equation <a class="reference internal" href="#sec-nepjac"><span class="std std-ref">Using Callback Functions</span></a>, we focus on the case where the problem is defined as in PETSc’s nonlinear solvers <a class="reference external" href="https://petsc.org/release/manualpages/SNES/SNES/" title="(in PETSc v3.24)"><span class="xref std std-doc">SNES</span></a>, that is, providing user-defined callback functions to compute the nonlinear function matrix, <span class="math notranslate nohighlight">\(T(\lambda)\)</span>, and its derivative, <span class="math notranslate nohighlight">\(T'(\lambda)\)</span>. We defer the discussion of using the split form of the nonlinear eigenproblem to section <a class="reference internal" href="#sec-nepsplit"><span class="std std-ref">Expressing the NEP in Split Form</span></a>.</p>
<section id="sec-nepjac">
<h3>Using Callback Functions<a class="headerlink" href="#sec-nepjac" title="Link to this heading">#</a></h3>
<p>A sample code for solving a nonlinear eigenproblem with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> is shown in listing <a class="reference internal" href="#fig-ex-nep"><span class="std std-ref">Example code for basic solution with NEP using callbacks</span></a>. The usual steps are performed, starting with the creation of the solver context with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPCreate.html">NEPCreate</a>()</span></code>. Then the problem matrices are defined, see discussion below. The call to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPSetFromOptions.html">NEPSetFromOptions</a>()</span></code> captures relevant options specified in the command line. The actual solver is invoked with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPSolve.html">NEPSolve</a>()</span></code>. Then, the solution is retrieved with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPGetConverged.html">NEPGetConverged</a>()</span></code> and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPGetEigenpair.html">NEPGetEigenpair</a>()</span></code>. Finally, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPDestroy.html">NEPDestroy</a>()</span></code> destroys the object.</p>
<div class="literal-block-wrapper docutils container" id="fig-ex-nep">
<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/NEP/NEP.html">NEP</a></span></code> using callbacks</span><a class="headerlink" href="#fig-ex-nep" title="Link to this code">#</a></div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</span><span class="p">;</span><span class="w"> </span><span class="cm">/* eigensolver context */</span>
<span class="n"><a href="https://petsc.org/release/manualpages/Mat/Mat/">Mat</a></span><span class="w"> </span><span class="n">F</span><span class="p">,</span><span class="w"> </span><span class="n">J</span><span class="p">;</span><span class="w"> </span><span class="cm">/* Function and Jacobian matrices */</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="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="n">ApplCtx</span><span class="w"> </span><span class="n">ctx</span><span class="p">;</span><span class="w"> </span><span class="cm">/* user-defined context */</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="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="n"><a href="../../manualpages/NEP/NEPCreate.html">NEPCreate</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">nep</span><span class="p">);</span>
<span class="cm">/* create and preallocate F and J matrices */</span>
<span class="n"><a href="../../manualpages/NEP/NEPSetFunction.html">NEPSetFunction</a></span><span class="p">(</span><span class="n">nep</span><span class="p">,</span><span class="w"> </span><span class="n">F</span><span class="p">,</span><span class="w"> </span><span class="n">F</span><span class="p">,</span><span class="w"> </span><span class="n">FormFunction</span><span class="p">,</span><span class="w"> </span><span class="o">&</span><span class="n">ctx</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/NEP/NEPSetJacobian.html">NEPSetJacobian</a></span><span class="p">(</span><span class="n">nep</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">FormJacobian</span><span class="p">,</span><span class="w"> </span><span class="o">&</span><span class="n">ctx</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/NEP/NEPSetFromOptions.html">NEPSetFromOptions</a></span><span class="p">(</span><span class="n">nep</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/NEP/NEPSolve.html">NEPSolve</a></span><span class="p">(</span><span class="n">nep</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/NEP/NEPGetConverged.html">NEPGetConverged</a></span><span class="p">(</span><span class="n">nep</span><span class="p">,</span><span class="w"> </span><span class="o">&</span><span class="n">nconv</span><span class="p">);</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="w"> </span><span class="n"><a href="../../manualpages/NEP/NEPGetEigenpair.html">NEPGetEigenpair</a></span><span class="p">(</span><span class="n">nep</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="w"> </span><span class="n"><a href="../../manualpages/NEP/NEPComputeError.html">NEPComputeError</a></span><span class="p">(</span><span class="n">nep</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/NEP/NEPErrorType.html">NEP_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="p">}</span>
<span class="n"><a href="../../manualpages/NEP/NEPDestroy.html">NEPDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">nep</span><span class="p">);</span>
</pre></div>
</div>
</div>
<p>In <a class="reference external" href="https://petsc.org/release/manualpages/SNES/SNES/" title="(in PETSc v3.24)"><span class="xref std std-doc">SNES</span></a>, the usual way to define a set of nonlinear equations <span class="math notranslate nohighlight">\(F(x)=0\)</span> is to provide two user-defined callback functions, one to compute the residual vector, <span class="math notranslate nohighlight">\(r=F(x)\)</span> for a given <span class="math notranslate nohighlight">\(x\)</span>, and another one to evaluate the Jacobian matrix, <span class="math notranslate nohighlight">\(J(x)=F'(x)\)</span>. In the case of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> there are some differences, since the function <span class="math notranslate nohighlight">\(T\)</span> depends on the parameter <span class="math notranslate nohighlight">\(\lambda\)</span> only. For a given value of <span class="math notranslate nohighlight">\(\lambda\)</span> and its associated vector <span class="math notranslate nohighlight">\(x\)</span>, the residual vector is defined as</p>
<div class="math notranslate nohighlight" id="equation-eq-nlres">
<span class="eqno">(5)<a class="headerlink" href="#equation-eq-nlres" title="Link to this equation">#</a></span>\[r=T(\lambda)x.\]</div>
<p>We require the user to provide a callback function to evaluate <span class="math notranslate nohighlight">\(T(\lambda)\)</span>, rather than computing the residual <span class="math notranslate nohighlight">\(r\)</span>. Once <span class="math notranslate nohighlight">\(T(\lambda)\)</span> has been built, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> solvers can compute its action on any vector <span class="math notranslate nohighlight">\(x\)</span>. Regarding the derivative, in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> we use <span class="math notranslate nohighlight">\(T'(\lambda)\)</span>, which will be referred to as the Jacobian matrix by analogy to <code class="docutils notranslate"><span class="pre"><a href="https://petsc.org/release/manualpages/SNES/SNES/">SNES</a></span></code>. This matrix must be computed with another callback function.</p>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>The derivative <span class="math notranslate nohighlight">\(T'(\lambda)\)</span> is optional in some cases, depending on the selected solver.</p>
</div>
<p>Hence, both callback functions must compute a matrix. The nonzero pattern of these matrices does not usually change, so they must be created at the beginning of the solution process. Then, these <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> objects are passed to the solver, together with the pointers to the callback functions, with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPSetFunction.html">NEPSetFunction</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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">F</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">P</span><span class="p">,</span><span class="n"><a href="../../manualpages/NEP/NEPFunctionFn.html">NEPFunctionFn</a></span><span class="w"> </span><span class="o">*</span><span class="n">fun</span><span class="p">,</span><span class="kt">void</span><span class="w"> </span><span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/NEP/NEPSetJacobian.html">NEPSetJacobian</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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">J</span><span class="p">,</span><span class="n"><a href="../../manualpages/NEP/NEPJacobianFn.html">NEPJacobianFn</a></span><span class="w"> </span><span class="o">*</span><span class="n">jac</span><span class="p">,</span><span class="kt">void</span><span class="w"> </span><span class="o">*</span><span class="n">ctx</span><span class="p">)</span>
</pre></div>
</div>
<p>The argument <code class="docutils notranslate"><span class="pre">ctx</span></code> is an optional user-defined context intended to contain application-specific parameters required to build <span class="math notranslate nohighlight">\(T(\lambda)\)</span> or <span class="math notranslate nohighlight">\(T'(\lambda)\)</span>, and it is received as the last argument in the callback functions. The callback routines also get an argument containing the value of <span class="math notranslate nohighlight">\(\lambda\)</span> at which <span class="math notranslate nohighlight">\(T\)</span> or <span class="math notranslate nohighlight">\(T'\)</span> must be evaluated. Note that the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPSetFunction.html">NEPSetFunction</a>()</span></code> callback takes two <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> arguments instead of one. The rationale for this is that some <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> solvers require to perform linear solves with <span class="math notranslate nohighlight">\(T(\lambda)\)</span> within the iteration (in <a class="reference external" href="https://petsc.org/release/manualpages/SNES/SNES/" title="(in PETSc v3.24)"><span class="xref std std-doc">SNES</span></a> this is done with the Jacobian), so <span class="math notranslate nohighlight">\(T(\lambda)\)</span> will be passed as the coefficient matrix to a <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> object. The second <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> argument <code class="docutils notranslate"><span class="pre">P</span></code> is the matrix from which the preconditioner is constructed (which is usually the same as <code class="docutils notranslate"><span class="pre">F</span></code>).</p>
<p>There is the possibility of solving the problem in a matrix-free fashion, that is, just implementing subroutines that compute the action of <span class="math notranslate nohighlight">\(T(\lambda)\)</span> or <span class="math notranslate nohighlight">\(T'(\lambda)\)</span> on a vector, instead of having to explicitly compute all nonzero entries of these two matrices. The example program <a class="reference external" href="https://slepc.upv.es/release/src/nep/tutorials/ex21.c.html" target="_blank">ex21.c</a> illustrates this, using the concept of <em>shell</em> matrices (see section <a class="reference internal" href="extra.html#sec-supported"><span class="std std-ref">Supported Matrix Types</span></a> for details).</p>
<section id="parameters-for-problem-definition">
<h4>Parameters for Problem Definition<a class="headerlink" href="#parameters-for-problem-definition" title="Link to this heading">#</a></h4>
<p>Once <span class="math notranslate nohighlight">\(T\)</span> and <span class="math notranslate nohighlight">\(T'\)</span> have been set up, the definition of the problem is completed with the number and location of the eigenvalues to compute, in a similar way as eigensolvers discussed in previous chapters.</p>
<p>The number of requested eigenvalues (and eigenvectors), <code class="docutils notranslate"><span class="pre">nev</span></code>, is established with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPSetDimensions.html">NEPSetDimensions</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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>By default, <code class="docutils notranslate"><span class="pre">nev</span></code>=1. The other two arguments control the dimension of the subspaces used internally (the number of column vectors, <code class="docutils notranslate"><span class="pre">ncv</span></code>, and the maximum projected dimension, <code class="docutils notranslate"><span class="pre">mpd</span></code>), although they are relevant only in eigensolvers based on subspace projection (basic algorithms ignore them). There are command-line keys for these parameters: <code class="docutils notranslate"><span class="pre">-nep_nev</span></code>, <code class="docutils notranslate"><span class="pre">-nep_ncv</span></code> and <code class="docutils notranslate"><span class="pre">-nep_mpd</span></code>.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-portionn">
<caption><span class="caption-text">Available possibilities for selection of the eigenvalues of interest in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code></span><a class="headerlink" href="#tab-portionn" 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/NEP/NEPWhich.html">NEPWhich</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/NEP/NEPWhich.html">NEP_LARGEST_MAGNITUDE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_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/NEP/NEPWhich.html">NEP_SMALLEST_MAGNITUDE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_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/NEP/NEPWhich.html">NEP_LARGEST_REAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_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/NEP/NEPWhich.html">NEP_SMALLEST_REAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_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/NEP/NEPWhich.html">NEP_LARGEST_IMAGINARY</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_largest_imaginary</span></code></p></td>
<td><p>Largest <span class="math notranslate nohighlight">\(\mathrm{Im}(\lambda)\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPWhich.html">NEP_SMALLEST_IMAGINARY</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_smallest_imaginary</span></code></p></td>
<td><p>Smallest <span class="math notranslate nohighlight">\(\mathrm{Im}(\lambda)\)</span></p></td>
</tr>
<tr class="row-even"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPWhich.html">NEP_TARGET_MAGNITUDE</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_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/NEP/NEPWhich.html">NEP_TARGET_REAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_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/NEP/NEPWhich.html">NEP_TARGET_IMAGINARY</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_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/NEP/NEPWhich.html">NEP_ALL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_all</span></code></p></td>
<td><p>All <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/NEP/NEPWhich.html">NEP_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 selection of the portion of the spectrum of interest, there are several alternatives listed in table <a class="reference internal" href="#tab-portionn"><span class="std std-ref">Available possibilities for selection of the eigenvalues of interest in NEP</span></a>, to be selected with the function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPSetWhichEigenpairs.html">NEPSetWhichEigenpairs</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</span><span class="p">,</span><span class="n"><a href="../../manualpages/NEP/NEPWhich.html">NEPWhich</a></span><span class="w"> </span><span class="n">which</span><span class="p">);</span>
</pre></div>
</div>
<p>The default is to compute the largest magnitude eigenvalues. For the sorting criteria relative to a target value, <span class="math notranslate nohighlight">\(\tau\)</span> must be specified with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPSetTarget.html">NEPSetTarget</a>()</span></code> or in the command-line with <code class="docutils notranslate"><span class="pre">-nep_target</span></code>.</p>
<p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> solvers can also work with a region of the complex plane (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/RG/RG.html">RG</a></span></code>), as discussed in section <a class="reference internal" href="eps.html#sec-region"><span class="std std-ref">Specifying a Region for Filtering</span></a> for linear problems. Some eigensolvers (NLEIGS) use the definition of the region to compute <code class="docutils notranslate"><span class="pre">nev</span></code> eigenvalues in its interior. If <em>all</em> eigenvalues inside the region are required, then a contour-integral method is required, see discussion in <span id="id3">[<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>
</section>
<section id="left-eigenvectors">
<h4>Left Eigenvectors<a class="headerlink" href="#left-eigenvectors" title="Link to this heading">#</a></h4>
<p>As in the case of linear eigensolvers, some <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> solvers have two-sided variants to compute also left eigenvectors. In the case of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code>, left eigenvectors are defined as</p>
<div class="math notranslate nohighlight" id="equation-eq-nepleft">
<span class="eqno">(6)<a class="headerlink" href="#equation-eq-nepleft" title="Link to this equation">#</a></span>\[y^*T(\lambda)=0^*,\qquad y\neq 0.\]</div>
<p>Two-sided variants can be selected with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPSetTwoSided.html">NEPSetTwoSided</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</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>
</section>
</section>
<section id="sec-nepsplit">
<h3>Expressing the NEP in Split Form<a class="headerlink" href="#sec-nepsplit" title="Link to this heading">#</a></h3>
<p>Instead of implementing callback functions for <span class="math notranslate nohighlight">\(T(\lambda)\)</span> and <span class="math notranslate nohighlight">\(T'(\lambda)\)</span>, a usually simpler alternative is to use the split form of the nonlinear eigenproblem, equation <a class="reference internal" href="#equation-eq-split">(4)</a>. Note that in split form, we have <span class="math notranslate nohighlight">\(T'(\lambda)=\sum_{i=0}^{\ell-1}A_if'_i(\lambda)\)</span>, so the derivatives of <span class="math notranslate nohighlight">\(f_i(\lambda)\)</span> are also required. As described below, we will represent each of the analytic functions <span class="math notranslate nohighlight">\(f_i\)</span> by means of an auxiliary object <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/FN/FN.html">FN</a></span></code> that holds both the function and its derivative.</p>
<div class="literal-block-wrapper docutils container" id="fig-ex-split">
<div class="code-block-caption"><span class="caption-text">Example code for defining the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> eigenproblem in split form</span><a class="headerlink" href="#fig-ex-split" title="Link to this code">#</a></div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/FN/FNCreate.html">FNCreate</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="o">&</span><span class="n">f1</span><span class="p">);</span><span class="w"> </span><span class="cm">/* f1 = -lambda */</span>
<span class="n"><a href="../../manualpages/FN/FNSetType.html">FNSetType</a></span><span class="p">(</span><span class="n">f1</span><span class="p">,</span><span class="n"><a href="../../manualpages/FN/FNRATIONAL.html">FNRATIONAL</a></span><span class="p">);</span>
<span class="n">coeffs</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="mf">-1.0</span><span class="p">;</span><span class="w"> </span><span class="n">coeffs</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="mf">0.0</span><span class="p">;</span>
<span class="n"><a href="../../manualpages/FN/FNRationalSetNumerator.html">FNRationalSetNumerator</a></span><span class="p">(</span><span class="n">f1</span><span class="p">,</span><span class="mi">2</span><span class="p">,</span><span class="n">coeffs</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/FN/FNCreate.html">FNCreate</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="o">&</span><span class="n">f2</span><span class="p">);</span><span class="w"> </span><span class="cm">/* f2 = 1 */</span>
<span class="n"><a href="../../manualpages/FN/FNSetType.html">FNSetType</a></span><span class="p">(</span><span class="n">f2</span><span class="p">,</span><span class="n"><a href="../../manualpages/FN/FNRATIONAL.html">FNRATIONAL</a></span><span class="p">);</span>
<span class="n">coeffs</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="mf">1.0</span><span class="p">;</span>
<span class="n"><a href="../../manualpages/FN/FNRationalSetNumerator.html">FNRationalSetNumerator</a></span><span class="p">(</span><span class="n">f2</span><span class="p">,</span><span class="mi">1</span><span class="p">,</span><span class="n">coeffs</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/FN/FNCreate.html">FNCreate</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="o">&</span><span class="n">f3</span><span class="p">);</span><span class="w"> </span><span class="cm">/* f3 = exp(-tau*lambda) */</span>
<span class="n"><a href="../../manualpages/FN/FNSetType.html">FNSetType</a></span><span class="p">(</span><span class="n">f3</span><span class="p">,</span><span class="n"><a href="../../manualpages/FN/FNEXP.html">FNEXP</a></span><span class="p">);</span>
<span class="n"><a href="../../manualpages/FN/FNSetScale.html">FNSetScale</a></span><span class="p">(</span><span class="n">f3</span><span class="p">,</span><span class="o">-</span><span class="n">tau</span><span class="p">,</span><span class="mf">1.0</span><span class="p">);</span>
<span class="n">mats</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">A</span><span class="p">;</span><span class="w"> </span><span class="n">funs</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">f2</span><span class="p">;</span>
<span class="n">mats</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">Id</span><span class="p">;</span><span class="w"> </span><span class="n">funs</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">f1</span><span class="p">;</span>
<span class="n">mats</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">B</span><span class="p">;</span><span class="w"> </span><span class="n">funs</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span><span class="w"> </span><span class="o">=</span><span class="w"> </span><span class="n">f3</span><span class="p">;</span>
<span class="n"><a href="../../manualpages/NEP/NEPSetSplitOperator.html">NEPSetSplitOperator</a></span><span class="p">(</span><span class="n">nep</span><span class="p">,</span><span class="mi">3</span><span class="p">,</span><span class="n">mats</span><span class="p">,</span><span class="n">funs</span><span class="p">,</span><span class="n">SUBSET_NONZERO_PATTERN</span><span class="p">);</span>
</pre></div>
</div>
</div>
<p>Hence, for the split form representation we must provide <span class="math notranslate nohighlight">\(\ell\)</span> matrices <span class="math notranslate nohighlight">\(A_i\)</span> and the corresponding functions <span class="math notranslate nohighlight">\(f_i(\lambda)\)</span>, by means of:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPSetSplitOperator.html">NEPSetSplitOperator</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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">l</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="../../manualpages/FN/FN.html">FN</a></span><span class="w"> </span><span class="n">f</span><span class="p">[],</span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/MatStructure/">MatStructure</a></span><span class="w"> </span><span class="n">str</span><span class="p">);</span>
</pre></div>
</div>
<p>Here, the <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatStructure/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatStructure</span></a> flag is used to indicate whether all matrices have the same (or subset) nonzero pattern with respect to the first one. Listing <a class="reference internal" href="#fig-ex-split"><span class="std std-ref">Example code for defining the NEP eigenproblem in split form</span></a> illustrates this usage with the problem of equation <a class="reference internal" href="#equation-eq-delay">(3)</a>, where <span class="math notranslate nohighlight">\(\ell=3\)</span> and the matrices are <span class="math notranslate nohighlight">\(I\)</span>, <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span> (note that in the code we have changed the order for efficiency reasons, since the nonzero pattern of <span class="math notranslate nohighlight">\(I\)</span> and <span class="math notranslate nohighlight">\(B\)</span> is a subset of <span class="math notranslate nohighlight">\(A\)</span>’s in this case). Two of the associated functions are polynomials (<span class="math notranslate nohighlight">\(-\lambda\)</span> and <span class="math notranslate nohighlight">\(1\)</span>) and the other one is the exponential <span class="math notranslate nohighlight">\(e^{-\tau\lambda}\)</span>.</p>
<p>Details of how to define the <span class="math notranslate nohighlight">\(f_i\)</span> functions by using the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/FN/FN.html">FN</a></span></code> class are provided in section <a class="reference internal" href="aux.html#sec-fn"><span class="std std-ref">FN: Mathematical Functions</span></a>.</p>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>Using the split form is required in order to be able to use some eigensolvers, in particular, those that project the nonlinear eigenproblem onto a low dimensional subspace and then use a dense nonlinear solver for the projected problem.</p>
</div>
<div class="pst-scrollable-table-container"><table class="table" id="tab-ntypeq">
<caption><span class="caption-text">Problem types considered in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code></span><a class="headerlink" href="#tab-ntypeq" 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/NEP/NEPProblemType.html">NEPProblemType</a></span></code></p></th>
<th class="head"><p>Command line key</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>General</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPProblemType.html">NEP_GENERAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_general</span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Rational</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPProblemType.html">NEP_RATIONAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-nep_rational</span></code></p></td>
</tr>
</tbody>
</table>
</div>
<p>When defining the problem in split form, it may also be useful to specify a problem type. For example, if the user knows that all <span class="math notranslate nohighlight">\(f_i\)</span> functions are rational, as in equation <a class="reference internal" href="#equation-eq-rep">(2)</a>, then setting the problem type to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPProblemType.html">NEP_RATIONAL</a></span></code> gives a hint to the solver that may simplify the solution process. The problem types currently supported for <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> are listed in table <a class="reference internal" href="#tab-ntypeq"><span class="std std-ref">Problem types considered in NEP</span></a>. When in doubt, use the default problem type (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPProblemType.html">NEP_GENERAL</a></span></code>).</p>
<p>The problem type can be specified at run time with the corresponding command line key or, more usually, within the program with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPSetProblemType.html">NEPSetProblemType</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</span><span class="p">,</span><span class="n"><a href="../../manualpages/NEP/NEPProblemType.html">NEPProblemType</a></span><span class="w"> </span><span class="n">type</span><span class="p">);</span>
</pre></div>
</div>
<p>Currently, the problem type is ignored in most solvers and it is taken into account only in NLEIGS for determining singularities automatically.</p>
</section>
</section>
<section id="selecting-the-solver">
<h2>Selecting the Solver<a class="headerlink" href="#selecting-the-solver" title="Link to this heading">#</a></h2>
<p>The solution method can be specified procedurally with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPSetType.html">NEPSetType</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</span><span class="p">,</span><span class="n"><a href="../../manualpages/NEP/NEPType.html">NEPType</a></span><span class="w"> </span><span class="n">method</span><span class="p">);</span>
</pre></div>
</div>
<p>or via the options database command <code class="docutils notranslate"><span class="pre">-nep_type</span></code> followed by the name of the method (see table <a class="reference internal" href="#tab-solversn"><span class="std std-ref">Nonlinear eigenvalue solvers available in the NEP module</span></a>). The methods currently available in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> are the following:</p>
<ul class="simple">
<li><p>Residual inverse iteration (RII), where in each iteration the eigenvector correction is computed as <span class="math notranslate nohighlight">\(T(\sigma)^{-1}\)</span> times the residual <span class="math notranslate nohighlight">\(r\)</span>.</p></li>
<li><p>Successive linear problems (SLP), where in each iteration a linear eigenvalue problem <span class="math notranslate nohighlight">\(T(\tilde\lambda)\tilde x=\mu T'(\tilde\lambda)\tilde x\)</span> is solved for the eigenvalue correction <span class="math notranslate nohighlight">\(\mu\)</span>.</p></li>
<li><p>Nonlinear Arnoldi, which builds an orthogonal basis <span class="math notranslate nohighlight">\(V_j\)</span> of a subspace expanded with the vectors generated by RII, then chooses the approximate eigenpair <span class="math notranslate nohighlight">\((\tilde\lambda,\tilde x)\)</span> such that <span class="math notranslate nohighlight">\(\tilde x=V_jy\)</span> and <span class="math notranslate nohighlight">\(V_j^*T(\tilde\lambda)V_jy=0\)</span>.</p></li>
<li><p>NLEIGS, which is based on a (rational) Krylov iteration operating on a companion-type linearization of a rational interpolant of the nonlinear function.</p></li>
<li><p>CISS, a contour-integral solver that allows computing all eigenvalues in a given region.</p></li>
<li><p>Polynomial interpolation, where a polynomial matrix <span class="math notranslate nohighlight">\(P(\lambda)\)</span> is built by evaluating <span class="math notranslate nohighlight">\(T(\cdot)\)</span> at a few points, then <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/PEP/PEP.html">PEP</a></span></code> is used to solve the resulting polynomial eigenproblem.</p></li>
</ul>
<div class="pst-scrollable-table-container"><table class="table" id="tab-solversn">
<caption><span class="caption-text">Nonlinear eigenvalue solvers available in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> module</span><a class="headerlink" href="#tab-solversn" 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/NEP/NEPType.html">NEPType</a></span></code></p></th>
<th class="head"><p>Options Database</p></th>
<th class="head"><p>Need <span class="math notranslate nohighlight">\(T'(\cdot)\)</span></p></th>
<th class="head"><p>Two-sided</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Residual inverse iteration</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPRII.html">NEPRII</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">rii</span></code></p></td>
<td><p>no</p></td>
<td><p></p></td>
</tr>
<tr class="row-odd"><td><p>Successive linear problems</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPSLP.html">NEPSLP</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">slp</span></code></p></td>
<td><p>yes</p></td>
<td><p>yes</p></td>
</tr>
<tr class="row-even"><td><p>Nonlinear Arnoldi</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPNARNOLDI.html">NEPNARNOLDI</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">narnoldi</span></code></p></td>
<td><p>no</p></td>
<td><p></p></td>
</tr>
<tr class="row-odd"><td><p>Rational Krylov (NLEIGS)</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPNLEIGS.html">NEPNLEIGS</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">nleigs</span></code></p></td>
<td><p>no</p></td>
<td><p>yes</p></td>
</tr>
<tr class="row-even"><td><p>Contour integral SS</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPCISS.html">NEPCISS</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">ciss</span></code></p></td>
<td><p>yes</p></td>
<td><p></p></td>
</tr>
<tr class="row-odd"><td><p>Polynomial interpolation</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPINTERPOL.html">NEPINTERPOL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">interpol</span></code></p></td>
<td><p>no</p></td>
<td><p></p></td>
</tr>
</tbody>
</table>
</div>
<p>The <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPSLP.html">NEPSLP</a></span></code> method performs a linearization that results in a (linear) generalized eigenvalue problem. This is handled by an <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object created internally. If required, this <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object can be extracted with the operation:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPSLPGetEPS.html">NEPSLPGetEPS</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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>This allows the application programmer to set any of the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> options directly within the code. These options can also be set through the command-line, simply by prefixing the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> options with <code class="docutils notranslate"><span class="pre">-nep_slp_</span></code>.</p>
<p>Similarly, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPINTERPOL.html">NEPINTERPOL</a></span></code> works with a <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/PEP/PEP.html">PEP</a></span></code> object internally, that can be retrieved by <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPInterpolGetPEP.html">NEPInterpolGetPEP</a>()</span></code>. Another relevant option of this solver is the degree of the interpolation polynomial, that can be set with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPInterpolSetInterpolation.html">NEPInterpolSetInterpolation</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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">deg</span><span class="p">);</span>
</pre></div>
</div>
<p>The polynomial interpolation solver currently uses Chebyshev polynomials of the 1st kind and requires the user to specify an interval of the real line where the eigenvalues must be computed, e.g.</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex22<span class="w"> </span>-nep_type<span class="w"> </span>interpol<span class="w"> </span>-rg_interval_endpoints<span class="w"> </span><span class="m">0</span>.1,14.0,-0.1,0.1<span class="w"> </span>-nep_nev<span class="w"> </span><span class="m">2</span><span class="w"> </span>-nep_interpol_interpolation_degree<span class="w"> </span><span class="m">15</span><span class="w"> </span>-nep_target<span class="w"> </span><span class="m">1</span>.0
</pre></div>
</div>
<p>Note that the target <span class="math notranslate nohighlight">\(\tau\)</span> must lie inside the region. For details about specifying a region, see section <a class="reference internal" href="aux.html#sec-rg"><span class="std std-ref">RG: Region</span></a>.</p>
<p>Some solvers such as <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPRII.html">NEPRII</a></span></code> and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPNARNOLDI.html">NEPNARNOLDI</a></span></code> need a <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> object to handle the solution of linear systems of equations. This <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> and can be retrieved, e.g., with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPRIIGetKSP.html">NEPRIIGetKSP</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/KSP/KSP/">KSP</a></span><span class="w"> </span><span class="o">*</span><span class="n">ksp</span><span class="p">);</span>
</pre></div>
</div>
<p>This <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> object is typically used to compute the action of <span class="math notranslate nohighlight">\(T(\sigma)^{-1}\)</span> on a given vector. In principle, <span class="math notranslate nohighlight">\(\sigma\)</span> is an approximation of an eigenvalue, but it is usually more efficient to keep this value constant, otherwise the factorization or preconditioner must be recomputed every time since eigensolvers update eigenvalue approximations in each iteration. This behavior can be changed with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPRIISetLagPreconditioner.html">NEPRIISetLagPreconditioner</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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">lag</span><span class="p">);</span>
</pre></div>
</div>
<p>Recomputing the preconditioner every 2 iterations, say, will introduce a considerable overhead, but may reduce the number of iterations significantly. Another related comment is that, when using an iterative linear solver, the requested accuracy is adapted as the outer iteration progresses, being the tolerance larger in the first solves. Again, the user can modify this behavior with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPRIISetConstCorrectionTol.html">NEPRIISetConstCorrectionTol</a>()</span></code>. Both options can also be changed at run time. As an example, consider the following command line:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex22<span class="w"> </span>-nep_type<span class="w"> </span>rii<span class="w"> </span>-nep_rii_lag_preconditioner<span class="w"> </span><span class="m">2</span><span class="w"> </span>-nep_rii_ksp_type<span class="w"> </span>bcgs<span class="w"> </span>-nep_rii_pc_type<span class="w"> </span>ilu<span class="w"> </span>-nep_rii_const_correction_tol<span class="w"> </span><span class="m">1</span><span class="w"> </span>-nep_rii_ksp_rtol<span class="w"> </span>1e-3
</pre></div>
</div>
<p>The example uses RII with BiCGStab plus ILU, where the preconditioner is updated every two outer iterations and linear systems are solved up to a tolerance of <span class="math notranslate nohighlight">\(10^{-3}\)</span>.</p>
<p>The NLEIGS solver is most appropriate for problems where <span class="math notranslate nohighlight">\(T(\cdot)\)</span> is singular at some known parts of the complex plane, for instance the case that <span class="math notranslate nohighlight">\(T(\cdot)\)</span> contains <span class="math notranslate nohighlight">\(\sqrt{\lambda}\)</span>. To treat this case effectively, the NLEIGS solver requires a discretization of the singularity set, which can be provided by the user in the form of a callback function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPNLEIGSSetSingularitiesFunction.html">NEPNLEIGSSetSingularitiesFunction</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</span><span class="p">,</span><span class="n"><a href="../../manualpages/NEP/NEPNLEIGSSingularitiesFn.html">NEPNLEIGSSingularitiesFn</a></span><span class="w"> </span><span class="o">*</span><span class="n">fun</span><span class="p">,</span><span class="kt">void</span><span class="w"> </span><span class="o">*</span><span class="n">ctx</span><span class="p">);</span>
</pre></div>
</div>
<p>Alternatively, if the problem is known to be a rational eigenvalue problem, the user can avoid the computation of singularities by just specifying the problem type with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPSetProblemType.html">NEPSetProblemType</a>()</span></code>, as explained at the end of the previous section. If none of the above functions is invoked by the user, then the NLEIGS solver attempts to determine the singularities automatically.</p>
</section>
<section id="retrieving-the-solution">
<h2>Retrieving the Solution<a class="headerlink" href="#retrieving-the-solution" title="Link to this heading">#</a></h2>
<p>The procedure for obtaining the computed eigenpairs is similar to previously discussed eigensolvers. After the call to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPSolve.html">NEPSolve</a>()</span></code>, the computed results are stored internally and a call to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPGetConverged.html">NEPGetConverged</a>()</span></code> must be issued to obtain the number of converged solutions. Then calling <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPGetEigenpair.html">NEPGetEigenpair</a>()</span></code> repeatedly will retrieve each eigenvalue-eigenvector pair.</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPGetEigenpair.html">NEPGetEigenpair</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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="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>In two-sided solvers (see last column of table <a class="reference internal" href="#tab-solversn"><span class="std std-ref">Nonlinear eigenvalue solvers available in the NEP module</span></a>), it is also possible to retrieve left eigenvectors with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPGetLeftEigenvector.html">NEPGetLeftEigenvector</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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>
<div class="admonition warning">
<p class="admonition-title">Warning</p>
<p><strong>Real vs complex scalar versions</strong>. The interface makes provision for returning a complex eigenvalue (or eigenvector) when doing the computation in a PETSc/SLEPc version built with real scalars, as is done in other eigensolvers such as <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>. However, in some cases this will not be possible. In particular, when callback functions are used and a complex eigenvalue approximation is hit, the solver will fail unless configured with complex scalars. The reason is that the definitions of callback functions only have a single <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> <code class="docutils notranslate"><span class="pre">lambda</span></code> argument and hence cannot handle complex arguments in real arithmetic.</p>
</div>
<p>The following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPComputeError.html">NEPComputeError</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</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/NEP/NEPErrorType.html">NEPErrorType</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>can be used to assess the accuracy of the computed solutions. The error is based on the 2-norm of the residual vector <span class="math notranslate nohighlight">\(r\)</span> defined in equation <a class="reference internal" href="#equation-eq-nlres">(5)</a>.</p>
<p>As in the case of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>, in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> the number of iterations carried out by the solver can be determined with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPGetIterationNumber.html">NEPGetIterationNumber</a>()</span></code>, and the tolerance and maximum number of iterations can be set with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEPSetTolerances.html">NEPSetTolerances</a>()</span></code>. Also, convergence can be monitored with either textual monitors <code class="docutils notranslate"><span class="pre">-nep_monitor</span></code>, <code class="docutils notranslate"><span class="pre">-nep_monitor_all</span></code>, <code class="docutils notranslate"><span class="pre">-nep_monitor_conv</span></code>, or graphical monitors <code class="docutils notranslate"><span class="pre">-nep_monitor</span> <span class="pre">draw::draw_lg</span></code>, <code class="docutils notranslate"><span class="pre">-nep_monitor_all</span> <span class="pre">draw::draw_lg</span></code>. See section <a class="reference internal" href="eps.html#sec-monitor"><span class="std std-ref">Controlling and Monitoring Convergence</span></a> for additional details. Similarly, there is support for viewing the computed solution as explained in section <a class="reference internal" href="eps.html#sec-epsviewers"><span class="std std-ref">Viewing the Solution</span></a>.</p>
<p>The <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/NEP/NEP.html">NEP</a></span></code> class also provides some kind of iterative refinement, similar to the one available in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/PEP/PEP.html">PEP</a></span></code>, see section <a class="reference internal" href="pep.html#sec-refine"><span class="std std-ref">Iterative Refinement</span></a>. The parameters can be set with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/NEP/NEPSetRefine.html">NEPSetRefine</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/NEP/NEP.html">NEP</a></span><span class="w"> </span><span class="n">nep</span><span class="p">,</span><span class="n"><a href="../../manualpages/NEP/NEPRefine.html">NEPRefine</a></span><span class="w"> </span><span class="n">refine</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">npart</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">its</span><span class="p">,</span><span class="n"><a href="../../manualpages/NEP/NEPRefineScheme.html">NEPRefineScheme</a></span><span class="w"> </span><span class="n">scheme</span><span class="p">);</span>
</pre></div>
</div>
<p class="rubric">References</p>
<div class="docutils container" id="id4">
<div role="list" class="citation-list">
<div class="citation" id="id46" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id1">Cam21</a><span class="fn-bracket">]</span></span>
<p>C. Campos and J. E. Roman. NEP: a module for the parallel solution of nonlinear eigenvalue problems in SLEPc. <em>ACM Trans. Math. Software</em>, 47(3):23:1–23:29, 2021. <a class="reference external" href="https://doi.org/10.1145/3447544">doi:10.1145/3447544</a>.</p>
</div>
<div class="citation" id="id43" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id2">Gut17</a><span class="fn-bracket">]</span></span>
<p>S. Güttel and F. Tisseur. The nonlinear eigenvalue problem. <em>Acta Numerica</em>, 26:1–94, 2017. <a class="reference external" href="https://doi.org/10.1017/S0962492917000034">doi:10.1017/S0962492917000034</a>.</p>
</div>
<div class="citation" id="id68" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id3">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="id33" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id2">Meh04</a><span class="fn-bracket">]</span></span>
<p>V. Mehrmann and H. Voss. Nonlinear eigenvalue problems: a challenge for modern eigenvalue methods. <em>GAMM Mitt.</em>, 27(2):121–152, 2004. <a class="reference external" href="https://doi.org/10.1002/gamm.201490007">doi:10.1002/gamm.201490007</a>.</p>
</div>
</div>
</div>
</section>
</section>
</article>
<footer class="prev-next-footer d-print-none">
<div class="prev-next-area">
<a class="left-prev"
href="pep.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">PEP: Polynomial Eigenvalue Problems</p>
</div>
</a>
<a class="right-next"
href="mfn.html"
title="next page">
<div class="prev-next-info">
<p class="prev-next-subtitle">next</p>
<p class="prev-next-title">MFN: Matrix Function</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-nep">General Nonlinear Eigenproblems</a></li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#defining-the-problem">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-nepjac">Using Callback Functions</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#parameters-for-problem-definition">Parameters for Problem Definition</a></li>
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#left-eigenvectors">Left Eigenvectors</a></li>
</ul>
</li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-nepsplit">Expressing the NEP in Split Form</a></li>
</ul>
</li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#selecting-the-solver">Selecting the Solver</a></li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#retrieving-the-solution">Retrieving the Solution</a></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/nep.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/nep.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>
|