1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206
|
<!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>SVD: Singular Value Decomposition — 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/svd';</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="PEP: Polynomial Eigenvalue Problems" href="pep.html" />
<link rel="prev" title="ST: Spectral Transformation" href="st.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 current active"><a class="current reference internal" href="#">SVD: Singular Value Decomposition</a></li>
<li class="toctree-l2"><a class="reference internal" href="pep.html">PEP: Polynomial Eigenvalue Problems</a></li>
<li class="toctree-l2"><a class="reference internal" href="nep.html">NEP: Nonlinear Eigenvalue Problems</a></li>
<li class="toctree-l2"><a class="reference internal" href="mfn.html">MFN: Matrix Function</a></li>
<li class="toctree-l2"><a class="reference internal" href="lme.html">LME: Linear Matrix Equation</a></li>
<li class="toctree-l2"><a class="reference internal" href="aux.html">Auxiliary Classes</a></li>
<li class="toctree-l2"><a class="reference internal" href="extra.html">Additional Information</a></li>
</ul>
</details></li>
<li class="toctree-l1 has-children"><a class="reference internal" href="../hands-on/index.html">Hands-on exercises</a><details><summary><span class="toctree-toggle" role="presentation"><i class="fa-solid fa-chevron-down"></i></span></summary><ul>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on0.html">Exercise 0: Hello World</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on1.html">Exercise 1: Standard Symmetric Eigenvalue Problem</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on2.html">Exercise 2: Standard Non-Symmetric Eigenvalue Problem</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on3.html">Exercise 3: Generalized Eigenvalue Problem Stored in a File</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on4.html">Exercise 4: Singular Value Decomposition</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on5.html">Exercise 5: Problem without Explicit Matrix Storage</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on6.html">Exercise 6: Parallel Execution</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on7.html">Exercise 7: Use of Deflation Subspaces</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on8.html">Exercise 8: Quadratic Eigenvalue Problem</a></li>
</ul>
</details></li>
<li class="toctree-l1"><a class="reference internal" href="../faq.html">Frequently Asked Questions</a></li>
<li class="toctree-l1"><a class="reference internal" href="../presentations.html">Presentations</a></li>
<li class="toctree-l1"><a class="reference internal" href="../license.html">License</a></li>
</ul>
</div>
</nav></div>
</div>
<div class="sidebar-primary-items__end sidebar-primary__section">
<div class="sidebar-primary-item">
<div id="ethical-ad-placement"
class="flat"
data-ea-publisher="readthedocs"
data-ea-type="readthedocs-sidebar"
data-ea-manual="true">
</div></div>
</div>
</div>
<main id="main-content" class="bd-main" role="main">
<div class="bd-content">
<div class="bd-article-container">
<div class="bd-header-article d-print-none">
<div class="header-article-items header-article__inner">
<div class="header-article-items__start">
<div class="header-article-item">
<nav aria-label="Breadcrumb" class="d-print-none">
<ul class="bd-breadcrumbs">
<li class="breadcrumb-item breadcrumb-home">
<a href="../../index.html" class="nav-link" aria-label="Home">
<i class="fa-solid fa-home"></i>
</a>
</li>
<li class="breadcrumb-item"><a href="../index.html" class="nav-link">Documentation</a></li>
<li class="breadcrumb-item"><a href="index.html" class="nav-link">SLEPc Users Manual</a></li>
<li class="breadcrumb-item active" aria-current="page"><span class="ellipsis">SVD: Singular Value Decomposition</span></li>
</ul>
</nav>
</div>
</div>
</div>
</div>
<div id="searchbox"></div>
<article class="bd-article">
<section class="tex2jax_ignore mathjax_ignore" id="svd-singular-value-decomposition">
<span id="ch-svd"></span><h1>SVD: Singular Value Decomposition<a class="headerlink" href="#svd-singular-value-decomposition" title="Link to this heading">#</a></h1>
<p>The Singular Value Decomposition (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code>) solver object can be used for computing a partial SVD of a rectangular matrix, and other related problems. It provides uniform and efficient access to several specific SVD solvers included in SLEPc, and also gives the possibility to compute the decomposition via the eigensolvers provided in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> package.</p>
<section id="sec-svdback">
<h2>Mathematical Background<a class="headerlink" href="#sec-svdback" title="Link to this heading">#</a></h2>
<p>This section provides some basic concepts about the singular value decomposition and other related problems. The objective is to set up the notation and also to justify some of the solution approaches, particularly those based on the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object. As in the case of eigensolvers, some of the implemented methods are described in detail in the SLEPc <a class="reference internal" href="../index.html#str"><span class="std std-ref">technical reports</span></a>.</p>
<p>For background material about the SVD, see for instance <span id="id1">[<a class="reference internal" href="#id15" title="Z. Bai, J. Demmel, J. Dongarra, A. Ruhe, and H. van der Vorst, editors. Templates for the Solution of Algebraic Eigenvalue Problems: A Practical Guide. Society for Industrial and Applied Mathematics, Philadelphia, PA, 2000.">Bai <em>et al.</em>, 2000</a>, ch. 6]</span>. Many other books such as <span id="id2">[<a class="reference internal" href="#id23" title="Å. Björck. Numerical Methods for Least Squares Problems. Society for Industrial and Applied Mathematics, Philadelphia, 1996.">Björck, 1996</a>, <a class="reference internal" href="#id32" title="P. C. Hansen. Rank-Deficient and Discrete Ill-Posed Problems: Numerical Aspects of Linear Inversion. Society for Industrial and Applied Mathematics, Philadelphia, PA, 1998.">Hansen, 1998</a>]</span> present the SVD from the perspective of its application to the solution of least squares problems and other related linear algebra problems.</p>
<section id="sec-svd">
<h3>The (Standard) Singular Value Decomposition (SVD)<a class="headerlink" href="#sec-svd" title="Link to this heading">#</a></h3>
<p>The singular value decomposition (SVD) of an <span class="math notranslate nohighlight">\(m\times n\)</span> matrix <span class="math notranslate nohighlight">\(A\)</span> can be written as</p>
<div class="math notranslate nohighlight" id="equation-eq-svd">
<span class="eqno">(1)<a class="headerlink" href="#equation-eq-svd" title="Link to this equation">#</a></span>\[A=U\Sigma V^*,\]</div>
<p>where <span class="math notranslate nohighlight">\(U=[u_1,\ldots,u_m]\)</span> is an <span class="math notranslate nohighlight">\(m\times m\)</span> unitary matrix (<span class="math notranslate nohighlight">\(U^*U=I\)</span>), <span class="math notranslate nohighlight">\(V=[v_1,\ldots,v_n]\)</span> is an <span class="math notranslate nohighlight">\(n\times n\)</span> unitary matrix (<span class="math notranslate nohighlight">\(V^*V=I\)</span>), and <span class="math notranslate nohighlight">\(\Sigma\)</span> is an <span class="math notranslate nohighlight">\(m\times n\)</span> diagonal matrix with real diagonal entries <span class="math notranslate nohighlight">\(\Sigma_{ii}=\sigma_i\)</span> for <span class="math notranslate nohighlight">\(i=1,\ldots,\min\{m,n\}\)</span>. If <span class="math notranslate nohighlight">\(A\)</span> is real, <span class="math notranslate nohighlight">\(U\)</span> and <span class="math notranslate nohighlight">\(V\)</span> are real and orthogonal. The vectors <span class="math notranslate nohighlight">\(u_i\)</span> are called the left singular vectors, the <span class="math notranslate nohighlight">\(v_i\)</span> are the right singular vectors, and the <span class="math notranslate nohighlight">\(\sigma_i\)</span> are the singular values.</p>
<figure class="align-default" id="fig-svd">
<img alt="Scheme of the thin SVD of a rectangular matrix." src="../../_images/fig-svd.svg" /><figcaption>
<p><span class="caption-text">Scheme of the thin SVD of a rectangular matrix</span><a class="headerlink" href="#fig-svd" title="Link to this image">#</a></p>
</figcaption>
</figure>
<p>In the following, we will assume that <span class="math notranslate nohighlight">\(m\geq n\)</span>. If <span class="math notranslate nohighlight">\(m<n\)</span> then <span class="math notranslate nohighlight">\(A\)</span> should be replaced by <span class="math notranslate nohighlight">\(A^*\)</span> (note that in SLEPc this is done transparently as described later in this chapter and the user need not worry about this). In the case that <span class="math notranslate nohighlight">\(m\geq n\)</span>, the top <span class="math notranslate nohighlight">\(n\)</span> rows of <span class="math notranslate nohighlight">\(\Sigma\)</span> contain <span class="math notranslate nohighlight">\(\mathrm{diag}(\sigma_1,\ldots,\sigma_n)\)</span> and its bottom <span class="math notranslate nohighlight">\(m-n\)</span> rows are zero. The relation equation <a class="reference internal" href="#equation-eq-svd">(1)</a> may also be written as <span class="math notranslate nohighlight">\(AV=U\Sigma\)</span>, or</p>
<div class="math notranslate nohighlight" id="equation-eq-svdleft">
<span class="eqno">(2)<a class="headerlink" href="#equation-eq-svdleft" title="Link to this equation">#</a></span>\[Av_i=u_i\sigma_i\;,\quad i=1,\ldots,n,\]</div>
<p>and also as <span class="math notranslate nohighlight">\(A^*U=V\Sigma^*\)</span>, or</p>
<div class="math notranslate nohighlight" id="equation-eq-svdright">
<span class="eqno">(3)<a class="headerlink" href="#equation-eq-svdright" title="Link to this equation">#</a></span>\[\begin{split}\begin{aligned}
A^*u_i&=v_i\sigma_i\;,\quad i=1,\ldots,n,\\
\end{aligned}\end{split}\]</div>
<div class="math notranslate nohighlight" id="equation-eq-svdright2">
<span class="eqno">(4)<a class="headerlink" href="#equation-eq-svdright2" title="Link to this equation">#</a></span>\[\begin{aligned}
A^*u_i&=0\;,\quad i=n+1,\ldots,m.
\end{aligned}\]</div>
<p>The last left singular vectors corresponding to equation <a class="reference internal" href="#equation-eq-svdright2">(4)</a> are often not computed, especially if <span class="math notranslate nohighlight">\(m\gg n\)</span>. In that case, the resulting factorization is sometimes called the <em>thin</em> SVD, <span class="math notranslate nohighlight">\(A=U_n\Sigma_n V_n^*\)</span>, and is depicted in figure <a class="reference internal" href="#fig-svd"><span class="std std-ref">Scheme of the thin SVD of a rectangular matrix</span></a>. This factorization can also be written as</p>
<div class="math notranslate nohighlight" id="equation-eq-svdouter">
<span class="eqno">(5)<a class="headerlink" href="#equation-eq-svdouter" title="Link to this equation">#</a></span>\[A=\sum_{i=1}^{n}\sigma_iu_iv_i^*.\]</div>
<p>Each <span class="math notranslate nohighlight">\((\sigma_i,u_i,v_i)\)</span> is called a singular triplet.</p>
<p>The singular values are real and nonnegative, <span class="math notranslate nohighlight">\(\sigma_1\geq\sigma_2\geq\ldots\geq\sigma_r>\sigma_{r+1}=\ldots=\sigma_n=0\)</span>, where <span class="math notranslate nohighlight">\(r=\mathrm{rank}(A)\)</span>. It can be shown that <span class="math notranslate nohighlight">\(\{u_1,\ldots,u_r\}\)</span> span the range of <span class="math notranslate nohighlight">\(A\)</span>, <span class="math notranslate nohighlight">\(\mathcal{R}(A)\)</span>, whereas <span class="math notranslate nohighlight">\(\{v_{r+1},\ldots,v_n\}\)</span> span the null space of <span class="math notranslate nohighlight">\(A\)</span>, <span class="math notranslate nohighlight">\(\mathcal{N}(A)\)</span>.</p>
<p>If the zero singular values are dropped from the sum in equation <a class="reference internal" href="#equation-eq-svdouter">(5)</a>, the resulting factorization, <span class="math notranslate nohighlight">\(A=\sum_{i=1}^{r}\sigma_iu_iv_i^*\)</span>, is called the <em>compact</em> SVD, <span class="math notranslate nohighlight">\(A=U_r\Sigma_r V_r^*\)</span>.</p>
<p>In the case of a very large and sparse <span class="math notranslate nohighlight">\(A\)</span>, it is usual to compute only a subset of <span class="math notranslate nohighlight">\(k\leq r\)</span> singular triplets. We will refer to this decomposition as the <em>truncated</em> SVD of <span class="math notranslate nohighlight">\(A\)</span>. It can be shown that the matrix <span class="math notranslate nohighlight">\(A_k=U_k\Sigma_k V_k^*\)</span> is the best rank-<span class="math notranslate nohighlight">\(k\)</span> approximation to matrix <span class="math notranslate nohighlight">\(A\)</span>, in the least squares sense.</p>
<p>In general, one can take an arbitrary subset of the summands in equation <a class="reference internal" href="#equation-eq-svdouter">(5)</a>, and the resulting factorization is called the <em>partial</em> SVD of <span class="math notranslate nohighlight">\(A\)</span>. As described later in this chapter, SLEPc allows the computation of a partial SVD corresponding to either the <span class="math notranslate nohighlight">\(k\)</span> largest or smallest singular triplets.</p>
<section id="equivalent-eigenvalue-problems">
<h4>Equivalent Eigenvalue Problems<a class="headerlink" href="#equivalent-eigenvalue-problems" title="Link to this heading">#</a></h4>
<p>It is possible to formulate the problem of computing the singular triplets of a matrix <span class="math notranslate nohighlight">\(A\)</span> as an eigenvalue problem involving a Hermitian matrix related to <span class="math notranslate nohighlight">\(A\)</span>. There are two possible ways of achieving this:</p>
<ol class="arabic simple">
<li><p>With the <em>cross product</em> matrix, either <span class="math notranslate nohighlight">\(A^*A\)</span> or <span class="math notranslate nohighlight">\(AA^*\)</span>.</p></li>
<li><p>With the <em>cyclic</em> matrix, <span class="math notranslate nohighlight">\(H(A)=\left[\begin{smallmatrix}0&A\\A^*&0\end{smallmatrix}\right]\)</span>.</p></li>
</ol>
<p>In SLEPc, the computation of the SVD is usually based on one of these two alternatives, either by passing one of these matrices to an <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object or by performing the computation implicitly.</p>
<p>By pre-multiplying equation <a class="reference internal" href="#equation-eq-svdleft">(2)</a> by <span class="math notranslate nohighlight">\(A^*\)</span> and then using equation <a class="reference internal" href="#equation-eq-svdright">(3)</a>, the following relation results</p>
<div class="math notranslate nohighlight" id="equation-eq-eigleft">
<span class="eqno">(6)<a class="headerlink" href="#equation-eq-eigleft" title="Link to this equation">#</a></span>\[A^*Av_i=\sigma_i^2v_i,\]</div>
<p>that is, the <span class="math notranslate nohighlight">\(v_i\)</span> are the eigenvectors of matrix <span class="math notranslate nohighlight">\(A^*A\)</span> with corresponding eigenvalues equal to <span class="math notranslate nohighlight">\(\sigma_i^2\)</span>. Note that after computing <span class="math notranslate nohighlight">\(v_i\)</span> the corresponding left singular vector, <span class="math notranslate nohighlight">\(u_i\)</span>, is readily available through equation <a class="reference internal" href="#equation-eq-svdleft">(2)</a> with just a matrix-vector product, <span class="math notranslate nohighlight">\(u_i=\frac{1}{\sigma_i}Av_i\)</span>.</p>
<p>Alternatively, one could first compute the left vectors and then the right ones. For this, pre-multiply equation <a class="reference internal" href="#equation-eq-svdright">(3)</a> by <span class="math notranslate nohighlight">\(A\)</span> and then use equation <a class="reference internal" href="#equation-eq-svdleft">(2)</a> to get</p>
<div class="math notranslate nohighlight" id="equation-eq-eigright">
<span class="eqno">(7)<a class="headerlink" href="#equation-eq-eigright" title="Link to this equation">#</a></span>\[AA^*u_i=\sigma_i^2u_i.\]</div>
<p>In this case, the right singular vectors are obtained as <span class="math notranslate nohighlight">\(v_i=\frac{1}{\sigma_i}A^*u_i\)</span>.</p>
<p>The two approaches represented in equations <a class="reference internal" href="#equation-eq-eigleft">(6)</a> and <a class="reference internal" href="#equation-eq-eigright">(7)</a> are very similar. Note however that <span class="math notranslate nohighlight">\(A^*A\)</span> is a square matrix of order <span class="math notranslate nohighlight">\(n\)</span> whereas <span class="math notranslate nohighlight">\(AA^*\)</span> is of order <span class="math notranslate nohighlight">\(m\)</span>. In cases where <span class="math notranslate nohighlight">\(m\gg n\)</span>, the computational effort will favor the <span class="math notranslate nohighlight">\(A^*A\)</span> approach. On the other hand, the eigenproblem <a class="reference internal" href="#equation-eq-eigleft">(6)</a> has <span class="math notranslate nohighlight">\(n-r\)</span> zero eigenvalues and the eigenproblem <a class="reference internal" href="#equation-eq-eigright">(7)</a> has <span class="math notranslate nohighlight">\(m-r\)</span> zero eigenvalues. Therefore, continuing with the assumption that <span class="math notranslate nohighlight">\(m\geq n\)</span>, even in the full rank case the <span class="math notranslate nohighlight">\(AA^*\)</span> approach may have a large null space resulting in difficulties if the smallest singular values are sought. In SLEPc, this will be referred to as the cross product approach and will use whichever matrix is smaller, either <span class="math notranslate nohighlight">\(A^*A\)</span> or <span class="math notranslate nohighlight">\(AA^*\)</span>.</p>
<p>Computing the SVD via the cross product approach may be adequate for determining the largest singular triplets of <span class="math notranslate nohighlight">\(A\)</span>, but the loss of accuracy can be severe for the smallest singular triplets. The cyclic matrix approach is an alternative that avoids this problem, but at the expense of significantly increasing the cost of the computation. Consider the eigendecomposition of</p>
<div class="math notranslate nohighlight" id="equation-eq-cyclic">
<span class="eqno">(8)<a class="headerlink" href="#equation-eq-cyclic" title="Link to this equation">#</a></span>\[\begin{split}H(A)=\begin{bmatrix}0&A\\A^*&0\end{bmatrix},\end{split}\]</div>
<p>which is a Hermitian matrix of order <span class="math notranslate nohighlight">\((m+n)\)</span>. It can be shown that <span class="math notranslate nohighlight">\(\pm\sigma_i\)</span> is a pair of eigenvalues of <span class="math notranslate nohighlight">\(H(A)\)</span> for <span class="math notranslate nohighlight">\(i=1,\ldots,r\)</span> and the other <span class="math notranslate nohighlight">\(m+n-2r\)</span> eigenvalues are zero. The unit eigenvectors associated with <span class="math notranslate nohighlight">\(\pm\sigma_i\)</span> are <span class="math notranslate nohighlight">\(\frac{1}{\sqrt{2}}\left[\begin{smallmatrix}\pm u_i\\v_i\end{smallmatrix}\right]\)</span>. Thus it is possible to extract the singular values and the left and right singular vectors of <span class="math notranslate nohighlight">\(A\)</span> directly from the eigenvalues and eigenvectors of <span class="math notranslate nohighlight">\(H(A)\)</span>. Note that in this case the singular values are not squared, and therefore the computed values will be more accurate (especially the small ones). The drawback in this case is that small eigenvalues are located in the interior of the spectrum.</p>
</section>
</section>
<section id="sec-gsvd">
<h3>The Generalized Singular Value Decomposition (GSVD)<a class="headerlink" href="#sec-gsvd" title="Link to this heading">#</a></h3>
<p>An extension of the SVD to the case of two matrices is the generalized singular value decomposition (GSVD), which can be applied in constrained least squares problems, among others. An overview of the problem can be found in <span id="id3">[<a class="reference internal" href="#id22" title="G. H. Golub and C. F. van Loan. Matrix Computations. The Johns Hopkins University Press, Baltimore, MD, third edition, 1996.">Golub and van Loan, 1996</a>, 8.7.3]</span>.</p>
<p>Consider two matrices, <span class="math notranslate nohighlight">\(A\in\mathbb{C}^{m\times n}\)</span> with <span class="math notranslate nohighlight">\(m\geq n\)</span> and <span class="math notranslate nohighlight">\(B\in\mathbb{C}^{p\times n}\)</span>. Note that both matrices must have the same column dimension. Then there exist two unitary matrices <span class="math notranslate nohighlight">\(U\in\mathbb{C}^{m\times m}\)</span> and <span class="math notranslate nohighlight">\(V\in\mathbb{C}^{p\times p}\)</span> and an invertible matrix <span class="math notranslate nohighlight">\(X\in\mathbb{C}^{n\times n}\)</span> such that</p>
<div class="math notranslate nohighlight" id="equation-eq-gsvd">
<span class="eqno">(9)<a class="headerlink" href="#equation-eq-gsvd" title="Link to this equation">#</a></span>\[U^*AX=C,\qquad V^*BX=S,\]</div>
<p>where <span class="math notranslate nohighlight">\(C=\mathrm{diag}(c_1,\dots,c_n)\)</span> and <span class="math notranslate nohighlight">\(S=\mathrm{diag}(s_{n-q+1},\dots,s_n)\)</span> with <span class="math notranslate nohighlight">\(q=\min(p,n)\)</span>. The values <span class="math notranslate nohighlight">\(c_i\)</span> and <span class="math notranslate nohighlight">\(s_i\)</span> are real and nonnegative, and the ratios define the generalized singular values,</p>
<div class="math notranslate nohighlight" id="equation-eq-gsvd-values">
<span class="eqno">(10)<a class="headerlink" href="#equation-eq-gsvd-values" title="Link to this equation">#</a></span>\[\sigma(A,B)\equiv\{c_1/s_1,\dots,c_q/s_q\},\]</div>
<p>and if <span class="math notranslate nohighlight">\(p<n\)</span> we can consider that the first <span class="math notranslate nohighlight">\(n-p\)</span> generalized singular values are infinite, as if <span class="math notranslate nohighlight">\(s_1=\dots=s_{n-p}=0\)</span>. Note that if <span class="math notranslate nohighlight">\(B\)</span> is the identity matrix, <span class="math notranslate nohighlight">\(X\)</span> can be taken to be unitary and then we recover the standard SVD, <span class="math notranslate nohighlight">\(\sigma(A,I)=\sigma(A)\)</span>, that is why equation <a class="reference internal" href="#equation-eq-gsvd">(9)</a> is considered a generalization of the SVD.</p>
<p>The diagonal matrices <span class="math notranslate nohighlight">\(C\)</span> and <span class="math notranslate nohighlight">\(S\)</span> satisfy <span class="math notranslate nohighlight">\(C^*C+S^*S=I\)</span>, and are related to the CS decomposition <span id="id4">[<a class="reference internal" href="#id22" title="G. H. Golub and C. F. van Loan. Matrix Computations. The Johns Hopkins University Press, Baltimore, MD, third edition, 1996.">Golub and van Loan, 1996</a>, 2.6.4]</span> associated with the orthogonal factor of the QR factorization of matrices <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span> stacked, that is, if</p>
<div class="math notranslate nohighlight" id="equation-eq-qr">
<span class="eqno">(11)<a class="headerlink" href="#equation-eq-qr" title="Link to this equation">#</a></span>\[\begin{split}Z:=\begin{bmatrix}A\\B\end{bmatrix}=\begin{bmatrix}Q_1\\Q_2\end{bmatrix}R,\end{split}\]</div>
<p>then <span class="math notranslate nohighlight">\(C\)</span> and <span class="math notranslate nohighlight">\(S\)</span> can be obtained from the singular values of <span class="math notranslate nohighlight">\(Q_1\)</span> and <span class="math notranslate nohighlight">\(Q_2\)</span>, respectively. The matrix <span class="math notranslate nohighlight">\(Z\)</span> is relevant for algorithms and is often built explicitly.</p>
<figure class="align-default" id="fig-gsvd">
<img alt="Scheme of the thin GSVD of two matrices $A$ and $B$, for the case $p < n$." src="../../_images/fig-gsvd.svg" /><figcaption>
<p><span class="caption-text">Scheme of the thin GSVD of two matrices <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span>, for the case <span class="math notranslate nohighlight">\(p < n\)</span></span><a class="headerlink" href="#fig-gsvd" title="Link to this image">#</a></p>
</figcaption>
</figure>
<p>The above description assumes that matrix <span class="math notranslate nohighlight">\(Z\)</span> has full column rank, and the rank of <span class="math notranslate nohighlight">\(B\)</span> is also <span class="math notranslate nohighlight">\(q\)</span>. In the general case where these assumptions do not hold, the structure of matrices <span class="math notranslate nohighlight">\(C\)</span> and <span class="math notranslate nohighlight">\(S\)</span> is a bit more complicated. This includes also the case where <span class="math notranslate nohighlight">\(m<n\)</span>. A detailed description of those cases can be found in <span id="id5">[<a class="reference internal" href="../../manualpages/EPS/EPSLAPACK.html#id3" title="E. Anderson, Z. Bai, C. Bischof, L. S. Blackford, J. Demmel, J. Dongarra, J. Du Croz, A. Greenbaum, S. Hammarling, A. McKenney, and D. Sorensen. LAPACK Users' Guide. Society for Industrial and Applied Mathematics, Philadelphia, PA, third edition, 1999.">Anderson <em>et al.</em>, 1999</a>, 2.3.5.3]</span>.</p>
<p>As in the case of the SVD, one can consider thin, compact, truncated, and partial versions of the GSVD. A pictorial representation of the thin GSVD is shown in figure <a class="reference internal" href="#fig-gsvd"><span class="std std-ref">Scheme of the thin GSVD of two matrices A and B, for the case p < n</span></a>.</p>
<p>The columns of <span class="math notranslate nohighlight">\(X\)</span>, <span class="math notranslate nohighlight">\(x_i\)</span>, are called the (right) generalized singular vectors. The left vectors in this case would correspond to the columns of <span class="math notranslate nohighlight">\(U\)</span> and <span class="math notranslate nohighlight">\(V\)</span>. In SLEPc, the user interface will provide these vectors stacked on top of each other, as a single <span class="math notranslate nohighlight">\((m+p)\)</span>-vector <span class="math notranslate nohighlight">\(\begin{bmatrix}u_i\\v_i\end{bmatrix}\)</span>.</p>
<section id="id6">
<h4>Equivalent Eigenvalue Problems<a class="headerlink" href="#id6" title="Link to this heading">#</a></h4>
<p>In the GSVD it is also possible to formulate the problem as an eigenvalue problem, which opens the door to approach its solution via <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>. The columns of <span class="math notranslate nohighlight">\(X\)</span> satisfy</p>
<div class="math notranslate nohighlight" id="equation-eq-gsvdeigcross">
<span class="eqno">(12)<a class="headerlink" href="#equation-eq-gsvdeigcross" title="Link to this equation">#</a></span>\[s_i^2A^*Ax_i=c_i^2B^*Bx_i,\]</div>
<p>and so if <span class="math notranslate nohighlight">\(s_i\neq 0\)</span> then <span class="math notranslate nohighlight">\(A^*Ax_i=\sigma_i^2B^*Bx_i\)</span>, a generalized eigenvalue problem for the matrix pencil <span class="math notranslate nohighlight">\((A^*A,B^*B)\)</span>. This is the analog of the cross product matrix eigenproblem of equation <a class="reference internal" href="#equation-eq-eigleft">(6)</a>.</p>
<p>The formulation that is analog to the eigenproblem associated with the cyclic matrix equation <a class="reference internal" href="#equation-eq-cyclic">(8)</a> is to solve the generalized eigenvalue problem defined by any of the matrix pairs</p>
<div class="math notranslate nohighlight" id="equation-eq-gsvdeigcyclic">
<span class="eqno">(13)<a class="headerlink" href="#equation-eq-gsvdeigcyclic" title="Link to this equation">#</a></span>\[\begin{split}\left(
\begin{bmatrix}0&A\\A^*&0\end{bmatrix},
\begin{bmatrix}I&0\\0&B^*B\end{bmatrix}
\right),
\qquad\text{or}\qquad
\left(
\begin{bmatrix}0&B\\B^*&0\end{bmatrix},
\begin{bmatrix}I&0\\0&A^*A\end{bmatrix}
\right).\end{split}\]</div>
</section>
</section>
<section id="sec-hsvd">
<h3>The Hyperbolic Singular Value Decomposition (HSVD)<a class="headerlink" href="#sec-hsvd" title="Link to this heading">#</a></h3>
<p>The hyperbolic singular value decomposition (HSVD) was introduced in <span id="id7">[<a class="reference internal" href="#id57" title="R. Onn, A. O. Steinhardt, and A. Bojanczyk. The hyperbolic singular value decomposition and applications. IEEE Trans. Signal Proces., 39(7):1575–1588, 1991. doi:10.1109/78.134396.">Onn <em>et al.</em>, 1991</a>]</span>, motivated by some signal processing applications such as the so-called covariance differencing problem. The formulation of the HSVD is similar to that of the SVD, except that <span class="math notranslate nohighlight">\(U\)</span> is orthogonal with respect to a signature matrix,</p>
<div class="math notranslate nohighlight" id="equation-eq-hsvd">
<span class="eqno">(14)<a class="headerlink" href="#equation-eq-hsvd" title="Link to this equation">#</a></span>\[A=U\Sigma V^*,\qquad U^*\Omega U=\tilde\Omega,\]</div>
<p>where <span class="math notranslate nohighlight">\(\Omega=\mathrm{diag}(\pm 1)\)</span> is an <span class="math notranslate nohighlight">\(m\times m\)</span> signature matrix provided by the user, while <span class="math notranslate nohighlight">\(\tilde\Omega\)</span> is another signature matrix obtained as part of the solution. Sometimes <span class="math notranslate nohighlight">\(U\)</span> is said to be a hyperexchange matrix, or also an <span class="math notranslate nohighlight">\((\Omega,\tilde\Omega)\)</span>-orthogonal matrix. Note that in the problem definition normally found in the literature it is <span class="math notranslate nohighlight">\(V\)</span> that is <span class="math notranslate nohighlight">\((\Omega,\tilde\Omega)\)</span>-orthogonal and not <span class="math notranslate nohighlight">\(U\)</span>. We choose this definition for consistency with respect to the generalized HSVD of two matrices. If the user wants to compute the HSVD according to the alternative definition, then it suffices to (conjugate) transpose the input matrix <span class="math notranslate nohighlight">\(A\)</span>, using for instance <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatHermitianTranspose/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatHermitianTranspose</span></a>().</p>
<p>As in the case of the SVD, the solution of the problem consists in singular triplets <span class="math notranslate nohighlight">\((\sigma_i,u_i,v_i)\)</span>, with <span class="math notranslate nohighlight">\(\sigma_i\)</span> real and nonnegative and sorted in nonincreasing order. Note that these quantities are different from those of section <a class="reference internal" href="#sec-svd"><span class="std std-ref">The (Standard) Singular Value Decomposition (SVD)</span></a>, even though we use the same notation here. With each singular triplet, there is an associated sign <span class="math notranslate nohighlight">\(\tilde\omega_i\)</span> (either 1 or <span class="math notranslate nohighlight">\(-1\)</span>), the corresponding diagonal element of <span class="math notranslate nohighlight">\(\tilde\Omega\)</span>. In SLEPc, this value is not returned by the user interface, but if required it can be easily computed as <span class="math notranslate nohighlight">\(u_i^*\Omega u_i\)</span>.</p>
<p>The relations between left and right singular vectors are slightly different from those of the standard SVD. We have <span class="math notranslate nohighlight">\(AV=U\Sigma\)</span> and <span class="math notranslate nohighlight">\(A^*\Omega U=V\Sigma^*\tilde\Omega\)</span>, so for <span class="math notranslate nohighlight">\(m\geq n\)</span> the following relations hold, together with equation <a class="reference internal" href="#equation-eq-svdleft">(2)</a>:</p>
<div class="math notranslate nohighlight" id="equation-eq-hsvdright">
<span class="eqno">(15)<a class="headerlink" href="#equation-eq-hsvdright" title="Link to this equation">#</a></span>\[\begin{split}\begin{aligned}
A^*\Omega u_i&=v_i\sigma_i\tilde\omega_i\;,\quad i=1,\ldots,n,\\
\end{aligned}\end{split}\]</div>
<div class="math notranslate nohighlight" id="equation-eq-hsvdright2">
<span class="eqno">(16)<a class="headerlink" href="#equation-eq-hsvdright2" title="Link to this equation">#</a></span>\[\begin{aligned}
A^*\Omega u_i&=0\;,\quad i=n+1,\ldots,m.
\end{aligned}\]</div>
<p>In SLEPc we will compute a partial HSVD consisting of either the largest or smallest hyperbolic singular triplets. Note that the sign <span class="math notranslate nohighlight">\(\tilde\omega_i\)</span> is not used when sorting for largest or smallest <span class="math notranslate nohighlight">\(\sigma_i\)</span>.</p>
<section id="id8">
<h4>Equivalent Eigenvalue Problems<a class="headerlink" href="#id8" title="Link to this heading">#</a></h4>
<p>Once again, we can derive cross and cyclic schemes to compute the decomposition by solving an eigenvalue problem. The cross product matrix approach has two forms, to be selected depending on whether <span class="math notranslate nohighlight">\(m\geq n\)</span> or not, as discussed in section <a class="reference internal" href="#sec-svd"><span class="std std-ref">The (Standard) Singular Value Decomposition (SVD)</span></a>. The first form,</p>
<div class="math notranslate nohighlight" id="equation-eq-heigleft">
<span class="eqno">(17)<a class="headerlink" href="#equation-eq-heigleft" title="Link to this equation">#</a></span>\[A^*\Omega Av_i=\sigma_i^2\tilde\omega_iv_i,\]</div>
<p>is derived by pre-multiplying equation <a class="reference internal" href="#equation-eq-svdleft">(2)</a> by <span class="math notranslate nohighlight">\(A^*\Omega\)</span> and then using equation <a class="reference internal" href="#equation-eq-hsvdright">(15)</a>. This eigenproblem can be solved as a HEP (cf. section <a class="reference internal" href="eps.html#sec-defprob"><span class="std std-ref">Defining the Problem</span></a>) and may have both positive and negative eigenvalues, corresponding to <span class="math notranslate nohighlight">\(\tilde\omega_i=1\)</span> and <span class="math notranslate nohighlight">\(\tilde\omega_i=-1\)</span>, respectively. Once the right vector <span class="math notranslate nohighlight">\(v_i\)</span> has been computed, the corresponding left vector can be obtained using equation <a class="reference internal" href="#equation-eq-svdleft">(2)</a> with just a matrix-vector product, <span class="math notranslate nohighlight">\(u_i=\sigma_i^{-1}Av_i\)</span>.</p>
<p>The second form of cross computes the left vectors first, by pre-multiplying equation <a class="reference internal" href="#equation-eq-hsvdright">(15)</a> by <span class="math notranslate nohighlight">\(A\)</span> and then using equation <a class="reference internal" href="#equation-eq-svdleft">(2)</a>,</p>
<div class="math notranslate nohighlight" id="equation-eq-heigright">
<span class="eqno">(18)<a class="headerlink" href="#equation-eq-heigright" title="Link to this equation">#</a></span>\[AA^*\Omega u_i=\sigma_i^2\tilde\omega_iu_i.\]</div>
<p>In this case, the right singular vectors are obtained as <span class="math notranslate nohighlight">\(v_i=(\sigma_i\tilde\omega_i)^{-1}A^*\Omega u_i\)</span>. The coefficient matrix of <a class="reference internal" href="#equation-eq-heigright">(18)</a> is non-Hermitian, so the eigenproblem has to be solved as non-Hermitian, or alternatively it can be formulated as a generalized eigenvalue problem of type <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHIEP.html">EPS_GHIEP</a></span></code> (cf. section <a class="reference internal" href="eps.html#sec-defprob"><span class="std std-ref">Defining the Problem</span></a>) for the indefinite pencil <span class="math notranslate nohighlight">\((AA^*,\Omega)\)</span>,</p>
<div class="math notranslate nohighlight" id="equation-eq-heigright2">
<span class="eqno">(19)<a class="headerlink" href="#equation-eq-heigright2" title="Link to this equation">#</a></span>\[AA^*\hat{u}_i=\sigma_i^2\tilde\omega_i\Omega \hat{u}_i,\]</div>
<p>with <span class="math notranslate nohighlight">\(\hat{u}_i=\Omega u_i\)</span>. The eigenvectors obtained from equation <a class="reference internal" href="#equation-eq-heigright">(18)</a> or equation <a class="reference internal" href="#equation-eq-heigright2">(19)</a> must be normalized so that <span class="math notranslate nohighlight">\(U^*\Omega U=\tilde\Omega\)</span> holds.</p>
<p>Finally, in the cyclic matrix approach for the HSVD we must solve a generalized eigenvalue problem defined by the matrices of order <span class="math notranslate nohighlight">\((m+n)\)</span></p>
<div class="math notranslate nohighlight" id="equation-eq-hcyclic">
<span class="eqno">(20)<a class="headerlink" href="#equation-eq-hcyclic" title="Link to this equation">#</a></span>\[\begin{split}H(A)=\begin{bmatrix}0&A\\A^*&0\end{bmatrix},\qquad
\hat\Omega=\begin{bmatrix}\Omega&0\\0&I\end{bmatrix}.\end{split}\]</div>
<p>As in the case of equation <a class="reference internal" href="#equation-eq-heigright2">(19)</a>, this problem is Hermitian-indefinite and hence it may have complex eigenvalues. However, it can be shown that nonzero eigenvalues of the matrix pair <span class="math notranslate nohighlight">\((H(A),\hat\Omega)\)</span> are either real (equal to <span class="math notranslate nohighlight">\(\pm\sigma_i\)</span> for a certain hyperbolic singular value <span class="math notranslate nohighlight">\(\sigma_i\)</span>) or purely imaginary (equal to <span class="math notranslate nohighlight">\(\pm\sigma_ij\)</span> for a certain <span class="math notranslate nohighlight">\(\sigma_i\)</span> with <span class="math notranslate nohighlight">\(j=\sqrt{-1}\)</span>). The associated eigenvectors are <span class="math notranslate nohighlight">\(\left[\begin{smallmatrix}\varsigma_i u_i\\v_i\end{smallmatrix}\right]\)</span>, where <span class="math notranslate nohighlight">\(\varsigma_i\)</span> is either <span class="math notranslate nohighlight">\(\pm 1\)</span> or <span class="math notranslate nohighlight">\(\pm j\)</span>.</p>
</section>
</section>
</section>
<section id="basic-usage">
<h2>Basic Usage<a class="headerlink" href="#basic-usage" title="Link to this heading">#</a></h2>
<p>From the perspective of the user interface, the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code> package is very similar to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>, with some differences that will be highlighted shortly.</p>
<div class="literal-block-wrapper docutils container" id="fig-ex-svd">
<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/SVD/SVD.html">SVD</a></span></code>.</span><a class="headerlink" href="#fig-ex-svd" title="Link to this code">#</a></div>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</span><span class="p">;</span><span class="w"> </span><span class="cm">/* <a href="../../manualpages/SVD/SVD.html">SVD</a> solver context */</span>
<span class="n"><a href="https://petsc.org/release/manualpages/Mat/Mat/">Mat</a></span><span class="w"> </span><span class="n">A</span><span class="p">;</span><span class="w"> </span><span class="cm">/* problem matrix */</span>
<span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">u</span><span class="p">,</span><span class="w"> </span><span class="n">v</span><span class="p">;</span><span class="w"> </span><span class="cm">/* singular vectors */</span>
<span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="n">sigma</span><span class="p">;</span><span class="w"> </span><span class="cm">/* singular value */</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/SVD/SVDCreate.html">SVDCreate</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">svd</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/SVD/SVDSetOperators.html">SVDSetOperators</a></span><span class="p">(</span><span class="n">svd</span><span class="p">,</span><span class="w"> </span><span class="n">A</span><span class="p">,</span><span class="w"> </span><span class="nb">NULL</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/SVD/SVDSetProblemType.html">SVDSetProblemType</a></span><span class="p">(</span><span class="n">svd</span><span class="p">,</span><span class="w"> </span><span class="n"><a href="../../manualpages/SVD/SVD_STANDARD.html">SVD_STANDARD</a></span><span class="p">);</span>
<span class="n"><a href="../../manualpages/SVD/SVDSetFromOptions.html">SVDSetFromOptions</a></span><span class="p">(</span><span class="n">svd</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/SVD/SVDSolve.html">SVDSolve</a></span><span class="p">(</span><span class="n">svd</span><span class="p">);</span>
<span class="n"><a href="../../manualpages/SVD/SVDGetConverged.html">SVDGetConverged</a></span><span class="p">(</span><span class="n">svd</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/SVD/SVDGetSingularTriplet.html">SVDGetSingularTriplet</a></span><span class="p">(</span><span class="n">svd</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">sigma</span><span class="p">,</span><span class="w"> </span><span class="n">u</span><span class="p">,</span><span class="w"> </span><span class="n">v</span><span class="p">);</span>
<span class="w"> </span><span class="n"><a href="../../manualpages/SVD/SVDComputeError.html">SVDComputeError</a></span><span class="p">(</span><span class="n">svd</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/SVD/SVDErrorType.html">SVD_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/SVD/SVDDestroy.html">SVDDestroy</a></span><span class="p">(</span><span class="o">&</span><span class="n">svd</span><span class="p">);</span>
</pre></div>
</div>
</div>
<p>The basic steps for computing a partial SVD with SLEPc are illustrated in listing <a class="reference internal" href="#fig-ex-svd"><span class="std std-ref">Example code for basic solution with SVD.</span></a>. The steps are more or less the same as those described in chapter <a class="reference internal" href="eps.html#ch-eps"><span class="std std-ref">EPS: Eigenvalue Problem Solver</span></a> for the eigenvalue problem. First, the solver context is created with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDCreate.html">SVDCreate</a>()</span></code>. Then the problem matrices have to be specified with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSetOperators.html">SVDSetOperators</a>()</span></code> and the type of problem can be selected via <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSetProblemType.html">SVDSetProblemType</a>()</span></code>. Then, a call to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSolve.html">SVDSolve</a>()</span></code> invokes the actual solver. After that, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDGetConverged.html">SVDGetConverged</a>()</span></code> is used to determine how many solutions have been computed, which are retrieved with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDGetSingularTriplet.html">SVDGetSingularTriplet</a>()</span></code>. Finally, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDDestroy.html">SVDDestroy</a>()</span></code> gets rid of the object.</p>
<p>If one compares this example code with the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> example in listing <a class="reference internal" href="eps.html#fig-ex-eps"><span class="std std-ref">Example code for basic solution with EPS</span></a>, the most outstanding differences are the following:</p>
<ul class="simple">
<li><p>The singular value is a <a class="reference external" href="https://petsc.org/release/manualpages/Sys/PetscReal/" title="(in PETSc v3.24)"><span class="xref std std-doc">PetscReal</span></a>, not a <a class="reference external" href="https://petsc.org/release/manualpages/Sys/PetscScalar/" title="(in PETSc v3.24)"><span class="xref std std-doc">PetscScalar</span></a>.</p></li>
<li><p>Each singular vector is defined with a single <a class="reference external" href="https://petsc.org/release/manualpages/Vec/Vec/" title="(in PETSc v3.24)"><span class="xref std std-doc">Vec</span></a> object, not two as was the case for eigenvectors.</p></li>
</ul>
</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>Defining the problem consists in specifying the problem matrix <span class="math notranslate nohighlight">\(A\)</span> (and also a second matrix <span class="math notranslate nohighlight">\(B\)</span> in case of the GSVD), the problem type, and the portion of the spectrum to be computed.</p>
<p>The problem matrices are provided with the following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDSetOperators.html">SVDSetOperators</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/Mat/">Mat</a></span><span class="w"> </span><span class="n">A</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/Mat/">Mat</a></span><span class="w"> </span><span class="n">B</span><span class="p">);</span>
</pre></div>
</div>
<p>where <code class="docutils notranslate"><span class="pre">A</span></code> can be any matrix, not necessarily square, stored in any allowed PETSc format including the matrix-free mechanism (see section <a class="reference internal" href="extra.html#sec-supported"><span class="std std-ref">Supported Matrix Types</span></a> for a detailed discussion), and <code class="docutils notranslate"><span class="pre">B</span></code> is generally set to <code class="docutils notranslate"><span class="pre">NULL</span></code> except in the case of computing the GSVD. If the problem is hyperbolic (cf. <a class="reference internal" href="#sec-hsvd"><span class="std std-ref">The Hyperbolic Singular Value Decomposition (HSVD)</span></a>) then in addition a signature matrix <span class="math notranslate nohighlight">\(\Omega\)</span> must be provided with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDSetSignature.html">SVDSetSignature</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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">omega</span><span class="p">);</span>
</pre></div>
</div>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>The signature matrix <span class="math notranslate nohighlight">\(\Omega\)</span> is represented as a vector containing values equal to <span class="math notranslate nohighlight">\(1\)</span> or <span class="math notranslate nohighlight">\(-1\)</span>.</p>
</div>
<p>The problem type can be specified with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDSetProblemType.html">SVDSetProblemType</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</span><span class="p">,</span><span class="n"><a href="../../manualpages/SVD/SVDProblemType.html">SVDProblemType</a></span><span class="w"> </span><span class="n">type</span><span class="p">);</span>
</pre></div>
</div>
<p>Note that in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code> it is currently not strictly required to call this function, since the problem type can be deduced from the information passed by the user: if only one matrix is passed to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSetOperators.html">SVDSetOperators</a>()</span></code> the problem type will default to a standard SVD (or hyperbolic SVD if a signature has been provided), while if two matrices are passed then it defaults to GSVD. Still, it is recommended to always call <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSetProblemType.html">SVDSetProblemType</a>()</span></code>. Table <a class="reference internal" href="#tab-ptypesvd"><span class="std std-ref">Problem types considered in SVD</span></a> lists the supported problem types.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-ptypesvd">
<caption><span class="caption-text">Problem types considered in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code></span><a class="headerlink" href="#tab-ptypesvd" 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/SVD/SVDProblemType.html">SVDProblemType</a></span></code></p></th>
<th class="head"><p>Command line key</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Standard SVD</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD_STANDARD.html">SVD_STANDARD</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_standard</span></code></p></td>
</tr>
<tr class="row-odd"><td><p>Generalized SVD (GSVD)</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD_GENERALIZED.html">SVD_GENERALIZED</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_generalized</span></code></p></td>
</tr>
<tr class="row-even"><td><p>Hyperbolic SVD (HSVD)</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD_HYPERBOLIC.html">SVD_HYPERBOLIC</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_hyperbolic</span></code></p></td>
</tr>
</tbody>
</table>
</div>
<p>It is important to note that all SVD solvers in SLEPc make use of both <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(A^*\)</span>, as suggested by the description in section <a class="reference internal" href="#sec-svd"><span class="std std-ref">The (Standard) Singular Value Decomposition (SVD)</span></a>. <span class="math notranslate nohighlight">\(A^*\)</span> is not explicitly passed as an argument to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSetOperators.html">SVDSetOperators</a>()</span></code>, therefore it will have to stem from <span class="math notranslate nohighlight">\(A\)</span>. There are two possibilities for this: either <span class="math notranslate nohighlight">\(A\)</span> is transposed explicitly and <span class="math notranslate nohighlight">\(A^*\)</span> is created as a distinct matrix, or <span class="math notranslate nohighlight">\(A^*\)</span> is handled implicitly via <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatMultTranspose/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatMultTranspose</span></a>() (or <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatMultHermitianTranspose/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatMultHermitianTranspose</span></a>() in the complex case) operations whenever a matrix-vector product is required in the algorithm. The default is to build <span class="math notranslate nohighlight">\(A^*\)</span> explicitly, but this behavior can be changed with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDSetImplicitTranspose.html">SVDSetImplicitTranspose</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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">impl</span><span class="p">);</span>
</pre></div>
</div>
<p>In section <a class="reference internal" href="#sec-svd"><span class="std std-ref">The (Standard) Singular Value Decomposition (SVD)</span></a>, it was mentioned that in SLEPc the cross product approach chooses the smallest of the two possible cases <span class="math notranslate nohighlight">\(A^*A\)</span> or <span class="math notranslate nohighlight">\(AA^*\)</span>, that is, <span class="math notranslate nohighlight">\(A^*A\)</span> is used if <span class="math notranslate nohighlight">\(A\)</span> is a tall, thin matrix (<span class="math notranslate nohighlight">\(m\geq n\)</span>), and <span class="math notranslate nohighlight">\(AA^*\)</span> is used if <span class="math notranslate nohighlight">\(A\)</span> is a fat, short matrix (<span class="math notranslate nohighlight">\(m<n\)</span>). In fact, what SLEPc does internally is that if <span class="math notranslate nohighlight">\(m<n\)</span> the roles of <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(A^*\)</span> are reversed. This is equivalent to transposing all the SVD factorization, so left singular vectors become right singular vectors and vice versa. This is actually done in all singular value solvers, not only the cross product approach. The objective is to simplify the number of cases to be treated internally by SLEPc, as well as to reduce the computational cost in some situations. Note that this is done transparently and the user need not worry about transposing the matrix, only to indicate how the transpose has to be handled, as explained above.</p>
<p>The user can specify how many singular values and vectors to compute. The default is to compute only one singular triplet. The next function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDSetDimensions.html">SVDSetDimensions</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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">nsv</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">ncv</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">mpd</span><span class="p">);</span>
</pre></div>
</div>
<p>allows the specification of the number of singular values to compute, <code class="docutils notranslate"><span class="pre">nsv</span></code>, and the number of column vectors to be used by the solution algorithm, <code class="docutils notranslate"><span class="pre">ncv</span></code>, that is, the largest dimension of the working subspace. These two parameters can also be set at run time with the options <code class="docutils notranslate"><span class="pre">-svd_nsv</span></code> and <code class="docutils notranslate"><span class="pre">-svd_ncv</span></code>. For example, the command line</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./program<span class="w"> </span>-svd_nsv<span class="w"> </span><span class="m">10</span><span class="w"> </span>-svd_ncv<span class="w"> </span><span class="m">24</span>
</pre></div>
</div>
<p>requests 10 singular values and instructs to use 24 column vectors. Note that <code class="docutils notranslate"><span class="pre">ncv</span></code> must be at least equal to <code class="docutils notranslate"><span class="pre">nsv</span></code>, although in general it is recommended (depending on the method) to work with a larger subspace, for instance <span class="math notranslate nohighlight">\(\mathtt{ncv}\geq2\cdot\mathtt{nsv}\)</span> or even more. As in the case of the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object, the last argument, <code class="docutils notranslate"><span class="pre">mpd</span></code>, can be used to limit the maximum dimension of the projected problem, as discussed in section <a class="reference internal" href="eps.html#sec-large-nev"><span class="std std-ref">Computing a Large Portion of the Spectrum</span></a>. Using this parameter is especially important in the case that a large number of singular values are requested.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-whichsvd">
<caption><span class="caption-text">Available possibilities for selection of the singular values of interest</span><a class="headerlink" href="#tab-whichsvd" 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/SVD/SVDWhich.html">SVDWhich</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/SVD/SVD_LARGEST.html">SVD_LARGEST</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_largest</span></code></p></td>
<td><p>Largest <span class="math notranslate nohighlight">\(\sigma\)</span></p></td>
</tr>
<tr class="row-odd"><td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD_SMALLEST.html">SVD_SMALLEST</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_smallest</span></code></p></td>
<td><p>Smallest <span class="math notranslate nohighlight">\(\sigma\)</span></p></td>
</tr>
</tbody>
</table>
</div>
<p>For the selection of the portion of the spectrum of interest, there are only two possibilities in the case of SVD: largest and smallest singular values, see table <a class="reference internal" href="#tab-whichsvd"><span class="std std-ref">Available possibilities for selection of the singular values of interest</span></a>. The default is to compute the largest ones, but this can be changed with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDSetWhichSingularTriplets.html">SVDSetWhichSingularTriplets</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</span><span class="p">,</span><span class="n"><a href="../../manualpages/SVD/SVDWhich.html">SVDWhich</a></span><span class="w"> </span><span class="n">which</span><span class="p">);</span>
</pre></div>
</div>
<p>This setting is used both for configuring how the algorithm seeks singular values and also for sorting the computed values. In contrast to the case of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>, computing singular values located in the interior part of the spectrum is difficult, the only possibility is to use an <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object combined with a spectral transformation (this possibility is explained in detail in the next section). Note that in this case, the value of <code class="docutils notranslate"><span class="pre">which</span></code> applies to the transformed spectrum.</p>
<section id="sec-thres">
<h3>Using a threshold to specify wanted singular values<a class="headerlink" href="#sec-thres" title="Link to this heading">#</a></h3>
<p>In some applications, the number of wanted singular values is not known a priori. For instance, suppose that matrix <span class="math notranslate nohighlight">\(A\)</span> is known to have low rank, and we want to approximate it by a truncated SVD containing the leading singular values. The goal is to have a good approximation, that is, discard only small singular values. The numerical rank <span class="math notranslate nohighlight">\(k\)</span> might be difficult to estimate, due to discretization error or other reasons. But since we assume that singular values decay abruptly around some unknown value <span class="math notranslate nohighlight">\(k\)</span>, we can configure SLEPc to detect this decay.</p>
<p>The threshold <span class="math notranslate nohighlight">\(\delta\)</span> can be set with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDSetThreshold.html">SVDSetThreshold</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="n">thres</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscBool/">PetscBool</a></span><span class="w"> </span><span class="n">rel</span><span class="p">);</span>
</pre></div>
</div>
<p>When <code class="docutils notranslate"><span class="pre">rel</span></code> is <a class="reference external" href="https://petsc.org/release/manualpages/Sys/PETSC_TRUE/" title="(in PETSc v3.24)"><span class="xref std std-doc">PETSC_TRUE</span></a>, the solver will compute <span class="math notranslate nohighlight">\(k\)</span> singular values, with <span class="math notranslate nohighlight">\(\sigma_j\geq\delta\cdot\sigma_1\)</span>, for <span class="math notranslate nohighlight">\(j=1,\dots,k-1\)</span>, where <span class="math notranslate nohighlight">\(k\)</span> is computed internally. In this case, the threshold is relative to the largest singular value <span class="math notranslate nohighlight">\(\sigma_1\)</span>, i.e., the matrix norm. It can be interpreted as a percentage, for instance <span class="math notranslate nohighlight">\(\delta=0.8\)</span> means compute singular values that are at least 80% of the matrix norm. Alternatively, an absolute threshold can be used by setting <code class="docutils notranslate"><span class="pre">rel</span></code> to <a class="reference external" href="https://petsc.org/release/manualpages/Sys/PETSC_FALSE/" title="(in PETSc v3.24)"><span class="xref std std-doc">PETSC_FALSE</span></a>.</p>
<figure class="align-default" id="fig-thres">
<img alt="Illustration of threshold usage with the singular values of the `rdb200` matrix. The red line represents a threshold of {math}`\delta=0.8`." src="../../_images/fig-thres.svg" /><figcaption>
<p><span class="caption-text">Illustration of threshold usage with the singular values of the <code class="docutils notranslate"><span class="pre">rdb200</span></code> matrix; the red line represents a relative threshold of <span class="math notranslate nohighlight">\(\delta=0.8\)</span></span><a class="headerlink" href="#fig-thres" title="Link to this image">#</a></p>
</figcaption>
</figure>
<p>In the low-rank example suggested above, the singular values will decay fast for a relatively small <span class="math notranslate nohighlight">\(k\)</span>, so a relative threshold <span class="math notranslate nohighlight">\(\delta=0.5\)</span> or less could do the job. But care must be taken in matrices whose singular values decay progressively, as in the example of figure <a class="reference internal" href="#fig-thres">Illustration of threshold usage with the singular values of the <code class="docutils notranslate"><span class="pre">rdb200</span></code> matrix</a>, since a small <span class="math notranslate nohighlight">\(\delta\)</span> would imply computing many singular triplets and hence a very high cost, both computationally and in memory. Since the number of computed singular values is not known a priori, the solver will need to reallocate the basis of vectors internally, to have enough room to accommodate all the singular vectors, so this option must be used with caution to avoid out-of-memory problems. The recommendation is to set the value of <code class="docutils notranslate"><span class="pre">ncv</span></code> to be larger than the estimated number of singular values, to minimize the number of reallocations. You can also use the <code class="docutils notranslate"><span class="pre">nsv</span></code> parameter in combination with the threshold to stop in case the number of computed singular triplets exceeds that value.</p>
<p>An absolute threshold is also available when computing smallest singular values, in which case the solver will compute the <span class="math notranslate nohighlight">\(k\)</span> smallest singular values, where <span class="math notranslate nohighlight">\(\sigma_j<\delta\)</span>, <span class="math notranslate nohighlight">\(j=n-k+1,\dots,n\)</span>.</p>
</section>
</section>
<section id="selecting-the-svd-solver">
<h2>Selecting the SVD Solver<a class="headerlink" href="#selecting-the-svd-solver" title="Link to this heading">#</a></h2>
<p>The available methods for computing the partial SVD are shown in table <a class="reference internal" href="#tab-svdsolvers">List of solvers available in the SVD module</a>. These methods can be classified in the following categories:</p>
<ul class="simple">
<li><p>Solvers based on <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>. These solvers set up an <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object internally, thus using the available eigensolvers for solving the SVD problem. The two possible approaches in this case are the cross product matrix and the cyclic matrix, as described in section <a class="reference internal" href="#sec-svdback"><span class="std std-ref">Mathematical Background</span></a>.</p></li>
<li><p>Specific SVD solvers. These are typically eigensolvers that have been adapted algorithmically to exploit the structure of the SVD or GSVD problems. There are two Lanczos-type solvers in this category: Lanczos and thick-restart Lanczos, see <span id="id9">[<a class="reference internal" href="#id74" title="V. Hernandez, J. E. Roman, and A. Tomas. Restarted Lanczos bidiagonalization for the SVD in SLEPc. Technical Report STR-8, Universitat Politècnica de València, 2007. URL: https://slepc.upv.es/documentation.">Hernandez <em>et al.</em>, 2007</a>]</span> for a detailed description of these methods. In this category, we could also add the randomized SVD (RSVD), a special solver that does not compute individual singular vectors accurately, but rather a low-rank approximation of <span class="math notranslate nohighlight">\(A\)</span> by means of randomization techniques.</p></li>
<li><p>The LAPACK solver. This is an interface to some LAPACK routines, analog of those in the case of eigenproblems. These routines operate in dense mode with only one processor and therefore are suitable only for moderate size problems. This solver should be used only for debugging purposes.</p></li>
<li><p>External packages: ScaLAPACK, KSVD and ELEMENTAL are dense packages and compute the complete SVD, while PRIMME offers Davidson-type methods to compute only a few singular triplets.</p></li>
</ul>
<p>The default solver is the one that uses the cross product matrix (<code class="docutils notranslate"><span class="pre">cross</span></code>), usually the fastest and most memory-efficient approach, at least for the standard SVD. See a more detailed explanation below.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-svdsolvers">
<caption><span class="caption-text">List of solvers available in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code> module; in the column of supported problems, ‘G’ means GSVD and ‘H’ means HSVD (the standard SVD is supported by all solvers)</span><a class="headerlink" href="#tab-svdsolvers" 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/SVD/SVDType.html">SVDType</a></span></code></p></th>
<th class="head"><p>Options Database</p></th>
<th class="head"><p>Supported</p></th>
<th class="head"><p>Default</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Cross Product</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDCROSS.html">SVDCROSS</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">cross</span></code></p></td>
<td><p>G,H</p></td>
<td><p><span class="math notranslate nohighlight">\(\star\)</span></p></td>
</tr>
<tr class="row-odd"><td><p>Cyclic Matrix</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDCYCLIC.html">SVDCYCLIC</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">cyclic</span></code></p></td>
<td><p>G,H</p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Lanczos</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDLANCZOS.html">SVDLANCZOS</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">lanczos</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Thick-restart Lanczos</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDTRLANCZOS.html">SVDTRLANCZOS</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">trlanczos</span></code></p></td>
<td><p>G,H</p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Randomized SVD (RSVD)</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDRANDOMIZED.html">SVDRANDOMIZED</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">randomized</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Wrapper to LAPACK</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDLAPACK.html">SVDLAPACK</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">lapack</span></code></p></td>
<td><p>G</p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Wrapper to ScaLAPACK</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSCALAPACK.html">SVDSCALAPACK</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">scalapack</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Wrapper to KSVD</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDKSVD.html">SVDKSVD</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">ksvd</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-even"><td><p>Wrapper to ELEMENTAL</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDELEMENTAL.html">SVDELEMENTAL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">elemental</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
<tr class="row-odd"><td><p>Wrapper to PRIMME</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDPRIMME.html">SVDPRIMME</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">primme</span></code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
<td><p><code class="docutils notranslate"> </code></p></td>
</tr>
</tbody>
</table>
</div>
<p>The solution method can be specified procedurally or via the command line. The application programmer can set it by means of:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDSetType.html">SVDSetType</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</span><span class="p">,</span><span class="n"><a href="../../manualpages/SVD/SVDType.html">SVDType</a></span><span class="w"> </span><span class="n">method</span><span class="p">);</span>
</pre></div>
</div>
<p>while the user writes the options database command <code class="docutils notranslate"><span class="pre">-svd_type</span></code> followed by the name of the method (see table <a class="reference internal" href="#tab-svdsolvers">List of solvers available in the SVD module</a>).</p>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>Not all solvers in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code> module support non-standard problems (GSVD and HSVD). <a class="reference internal" href="#tab-svdsolvers">List of solvers available in the SVD module</a> indicates which solvers can be used for these.</p>
</div>
<p>The <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>-based solvers deserve some additional comments. These SVD solvers work by creating an <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object internally and setting up an eigenproblem of type <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_HEP.html">EPS_HEP</a></span></code> (or <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code> in the case of the GSVD, or <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHIEP.html">EPS_GHIEP</a></span></code> in the case of the HSVD). These solvers implement the cross product matrix and the cyclic matrix approaches as described in section <a class="reference internal" href="#sec-svdback"><span class="std std-ref">Mathematical Background</span></a>. Therefore, the operator matrix associated with the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object will be <span class="math notranslate nohighlight">\(A^*A\)</span> in the case of the <code class="docutils notranslate"><span class="pre">cross</span></code> solver and <span class="math notranslate nohighlight">\(H(A)\)</span> in the case of the <code class="docutils notranslate"><span class="pre">cyclic</span></code> solver (with variations in the case of non-standard problems).</p>
<p>In the case of the <code class="docutils notranslate"><span class="pre">cross</span></code> solver, the matrix <span class="math notranslate nohighlight">\(A^*A\)</span> is not built explicitly by default, since the resulting matrix may be much less sparse than the original matrix. By default, a shell matrix is created internally in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code> object and passed to the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object. Still, the user may choose to force the computation of <span class="math notranslate nohighlight">\(A^*A\)</span> explicitly, by means of PETSc’s sparse matrix-matrix product subroutines. This is set with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDCrossSetExplicitMatrix.html">SVDCrossSetExplicitMatrix</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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">explicit</span><span class="p">);</span>
</pre></div>
</div>
<p>In the <code class="docutils notranslate"><span class="pre">cyclic</span></code> solver the user can also choose between handling <span class="math notranslate nohighlight">\(H(A)\)</span> implicitly as a shell matrix (the default), or forming it explicitly, that is, storing its elements in a distinct matrix. The function for setting this option is:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDCyclicSetExplicitMatrix.html">SVDCyclicSetExplicitMatrix</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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">explicit</span><span class="p">);</span>
</pre></div>
</div>
<p>The <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object associated with the <code class="docutils notranslate"><span class="pre">cross</span></code> and <code class="docutils notranslate"><span class="pre">cyclic</span></code> SVD solvers is created with a set of reasonable default parameters. However, it may sometimes be necessary to change some of the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> options such as the eigensolver. To allow application programmers 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, the following routines are provided to extract the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> context from the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code> object:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDCrossGetEPS.html">SVDCrossGetEPS</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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>
<span class="n"><a href="../../manualpages/SVD/SVDCyclicGetEPS.html">SVDCyclicGetEPS</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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>A more convenient way of changing <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> options is through the command-line. This is achieved 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">-svd_cross_</span></code> or <code class="docutils notranslate"><span class="pre">-svd_cyclic_</span></code> as in the following example:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./program<span class="w"> </span>-svd_type<span class="w"> </span>cross<span class="w"> </span>-svd_cross_eps_type<span class="w"> </span>gd
</pre></div>
</div>
<p>At this point, one may consider changing also the options of the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> object associated with the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object in <code class="docutils notranslate"><span class="pre">cross</span></code> and <code class="docutils notranslate"><span class="pre">cyclic</span></code> SVD solvers, for example to compute singular values located at the interior of the spectrum via a shift-and-invert transformation. This is indeed possible, but some considerations have to be taken into account. When <span class="math notranslate nohighlight">\(A^*A\)</span> or <span class="math notranslate nohighlight">\(H(A)\)</span> are managed as shell matrices, then the potential of the spectral transformation is limited seriously, because some of the required operations will not be defined if working with implicit matrices (this is discussed briefly in sections <a class="reference internal" href="extra.html#sec-supported"><span class="std std-ref">Supported Matrix Types</span></a> and <a class="reference internal" href="st.html#sec-explicit"><span class="std std-ref">Explicit Computation of the Coefficient Matrix</span></a>). The following example illustrates the computation of interior singular values with the <code class="docutils notranslate"><span class="pre">cyclic</span></code> solver with explicit <span class="math notranslate nohighlight">\(H(A)\)</span> matrix:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./program<span class="w"> </span>-svd_type<span class="w"> </span>cyclic<span class="w"> </span>-svd_cyclic_explicitmatrix<span class="w"> </span>-svd_cyclic_st_type<span class="w"> </span>sinvert<span class="w"> </span>-svd_cyclic_eps_target<span class="w"> </span><span class="m">12</span>.0<span class="w"> </span>-svd_cyclic_st_ksp_type<span class="w"> </span>preonly<span class="w"> </span>-svd_cyclic_st_pc_type<span class="w"> </span>lu
</pre></div>
</div>
<p>Similarly, in the case of GSVD the thick-restart Lanczos solver uses 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> solver internally, that can be configured by accessing it with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDTRLanczosGetKSP.html">SVDTRLanczosGetKSP</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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>or with the corresponding command-line options prefixed with <code class="docutils notranslate"><span class="pre">-svd_trlanczos_</span></code>. 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 used to solve a linear least squares problem at each Lanczos step with coefficient matrix <span class="math notranslate nohighlight">\(Z\)</span> equation <a class="reference internal" href="#equation-eq-qr">(11)</a>, which by default is a shell matrix but the user can choose to create it explicitly with the function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDTRLanczosSetExplicitMatrix.html">SVDTRLanczosSetExplicitMatrix</a>()</span></code>.</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>Once the call to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSolve.html">SVDSolve</a></span></code> is complete, all the data associated with the computed partial SVD is kept internally in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</a></span></code> object. This information can be obtained by the calling program by means of a set of functions described below.</p>
<p>As in the case of eigenproblems, the number of computed singular triplets depends on the convergence and, therefore, it may be different from the number of solutions requested by the user. So the first task is to find out how many solutions are available, with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDGetConverged.html">SVDGetConverged</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="o">*</span><span class="n">nconv</span><span class="p">);</span>
</pre></div>
</div>
<p>Usually, the number of converged solutions, <code class="docutils notranslate"><span class="pre">nconv</span></code>, will be equal to <code class="docutils notranslate"><span class="pre">nsv</span></code>, but in general it can be a number ranging from 0 to <code class="docutils notranslate"><span class="pre">ncv</span></code> (here, <code class="docutils notranslate"><span class="pre">nsv</span></code> and <code class="docutils notranslate"><span class="pre">ncv</span></code> are the arguments of function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSetDimensions.html">SVDSetDimensions</a>()</span></code>).</p>
<p>Normally, the user is interested in the singular values only, or the complete singular triplets. The function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDGetSingularTriplet.html">SVDGetSingularTriplet</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="o">*</span><span class="n">sigma</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">u</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Vec/Vec/">Vec</a></span><span class="w"> </span><span class="n">v</span><span class="p">);</span>
</pre></div>
</div>
<p>returns the <span class="math notranslate nohighlight">\(j\)</span>-th computed singular triplet, <span class="math notranslate nohighlight">\((\sigma_j,u_j,v_j)\)</span>, where both <span class="math notranslate nohighlight">\(u_j\)</span> and <span class="math notranslate nohighlight">\(v_j\)</span> are normalized to have unit norm<a class="footnote-reference brackets" href="#hsvd" id="id10" role="doc-noteref"><span class="fn-bracket">[</span>1<span class="fn-bracket">]</span></a>. Typically, this function is called inside a loop for each value of <code class="docutils notranslate"><span class="pre">j</span></code> from 0 to <code class="docutils notranslate"><span class="pre">nconv</span></code>–1. Note that singular values are ordered according to the same criterion specified with function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDSetWhichSingularTriplets.html">SVDSetWhichSingularTriplets</a>()</span></code> for selecting the portion of the spectrum of interest.</p>
<p>In some applications, it may be enough to compute only the right singular vectors. This is especially important in cases in which memory requirements are critical (remember that both <span class="math notranslate nohighlight">\(U_k\)</span> and <span class="math notranslate nohighlight">\(V_k\)</span> are dense matrices, and <span class="math notranslate nohighlight">\(U_k\)</span> may require much more storage than <span class="math notranslate nohighlight">\(V_k\)</span>, see figure <a class="reference internal" href="#fig-svd"><span class="std std-ref">Scheme of the thin SVD of a rectangular matrix</span></a>). In SLEPc, there is no general option for specifying this, but the default behavior of some solvers is to compute only right vectors and allocate/compute left vectors only in the case that the user requests them. This is done in the <code class="docutils notranslate"><span class="pre">cross</span></code> solver and in some special variants of other solvers such as one-sided Lanczos (consult <span id="id11">[<a class="reference internal" href="#id74" title="V. Hernandez, J. E. Roman, and A. Tomas. Restarted Lanczos bidiagonalization for the SVD in SLEPc. Technical Report STR-8, Universitat Politècnica de València, 2007. URL: https://slepc.upv.es/documentation.">Hernandez <em>et al.</em>, 2007</a>]</span> for specific solver options).</p>
<p>In the case of the GSVD, the <code class="docutils notranslate"><span class="pre">sigma</span></code> argument of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDGetSingularTriplet.html">SVDGetSingularTriplet</a>()</span></code> contains <span class="math notranslate nohighlight">\(\sigma_i=c_i/s_i\)</span> and the second <a class="reference external" href="https://petsc.org/release/manualpages/Vec/Vec/" title="(in PETSc v3.24)"><span class="xref std std-doc">Vec</span></a> argument (<code class="docutils notranslate"><span class="pre">v</span></code>) contains the right singular vectors (<span class="math notranslate nohighlight">\(x_i\)</span>), while the first <a class="reference external" href="https://petsc.org/release/manualpages/Vec/Vec/" title="(in PETSc v3.24)"><span class="xref std std-doc">Vec</span></a> argument (<code class="docutils notranslate"><span class="pre">u</span></code>) contains the other vectors of the decomposition stacked on top of each other, as a single <span class="math notranslate nohighlight">\((m+p)\)</span>-vector: <span class="math notranslate nohighlight">\(\left[\begin{smallmatrix}u_i\\v_i\end{smallmatrix}\right]\)</span>.</p>
<div class="admonition warning">
<p class="admonition-title">Warning</p>
<p>In the GSVD, one has to be careful when dealing with the <span class="math notranslate nohighlight">\((m+p)\)</span>-vector <code class="docutils notranslate"><span class="pre">u</span></code> in parallel runs, since the entries of the upper and lower part will get mixed unless we manage it as a <a class="reference external" href="https://petsc.org/release/manualpages/Vec/VECNEST/" title="(in PETSc v3.24)"><span class="xref std std-doc">VECNEST</span></a>. One should follow a procedure similar to the one discussed in section <a class="reference internal" href="eps.html#sec-structured-vectors"><span class="std std-ref">Extracting Eigenvectors of Structured Eigenproblems</span></a>.</p>
</div>
<section id="reliability-of-the-computed-solution">
<h3>Reliability of the Computed Solution<a class="headerlink" href="#reliability-of-the-computed-solution" title="Link to this heading">#</a></h3>
<p>In SVD computations, a-posteriori error bounds are much the same as in the case of Hermitian eigenproblems, due to the equivalence discussed in section <a class="reference internal" href="#sec-svd"><span class="std std-ref">The (Standard) Singular Value Decomposition (SVD)</span></a>. The residual vector is defined in terms of the cyclic matrix, <span class="math notranslate nohighlight">\(H(A)\)</span>, so its norm is</p>
<div class="math notranslate nohighlight" id="equation-eq-svd-residual-norm">
<span class="eqno">(21)<a class="headerlink" href="#equation-eq-svd-residual-norm" title="Link to this equation">#</a></span>\[\|r_\mathrm{SVD}\|_2=\left(\|A\tilde{v}-\tilde{\sigma}\tilde{u}\|_2^2+\|A^*\tilde{u}-\tilde{\sigma}\tilde{v}\|_2^2\right)^{\frac{1}{2}},\]</div>
<p>where <span class="math notranslate nohighlight">\(\tilde{\sigma}\)</span>, <span class="math notranslate nohighlight">\(\tilde{u}\)</span> and <span class="math notranslate nohighlight">\(\tilde{v}\)</span> represent any of the <code class="docutils notranslate"><span class="pre">nconv</span></code> computed singular triplets delivered by <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDGetSingularTriplet.html">SVDGetSingularTriplet</a>()</span></code>.</p>
<p>Given the above definition, the following relation holds</p>
<div class="math notranslate nohighlight" id="equation-eq-svd-residual-norm-relation">
<span class="eqno">(22)<a class="headerlink" href="#equation-eq-svd-residual-norm-relation" title="Link to this equation">#</a></span>\[|\sigma-\tilde{\sigma}|\leq \|r_\mathrm{SVD}\|_2,\]</div>
<p>where <span class="math notranslate nohighlight">\(\sigma\)</span> is an exact singular value. The associated error can be obtained in terms of <span class="math notranslate nohighlight">\(\|r_\mathrm{SVD}\|_2\)</span> with the following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/SVD/SVDComputeError.html">SVDComputeError</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/SVD/SVD.html">SVD</a></span><span class="w"> </span><span class="n">svd</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/SVD/SVDErrorType.html">SVDErrorType</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>In the case of the GSVD, the function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDComputeError.html">SVDComputeError</a>()</span></code> will compute a residual norm based on the two relations <a class="reference internal" href="#equation-eq-gsvd">(9)</a>,</p>
<div class="math notranslate nohighlight" id="equation-eq-gsvd-residual-norm">
<span class="eqno">(23)<a class="headerlink" href="#equation-eq-gsvd-residual-norm" title="Link to this equation">#</a></span>\[\|r_\mathrm{GSVD}\|_2=\left(\|\tilde{s}^2A^*\tilde{u}-\tilde{c}B^*B\tilde{x}\|_2^2+\|\tilde{c}^2B^*\tilde{v}-\tilde{s}A^*A\tilde{x}\|_2^2\right)^{\frac{1}{2}},\]</div>
<p>where <span class="math notranslate nohighlight">\(\tilde{x}\)</span>, <span class="math notranslate nohighlight">\(\tilde{u}\)</span>, <span class="math notranslate nohighlight">\(\tilde{v}\)</span> are the computed singular vectors corresponding to <span class="math notranslate nohighlight">\(\tilde{\sigma}\)</span>, and <span class="math notranslate nohighlight">\(\tilde{c}\)</span>, <span class="math notranslate nohighlight">\(\tilde{s}\)</span> are obtained from <span class="math notranslate nohighlight">\(\tilde{\sigma}\)</span> as <span class="math notranslate nohighlight">\(\tilde{s}=1/\sqrt{1+\tilde{\sigma}^2}\)</span> and <span class="math notranslate nohighlight">\(\tilde{c}=\tilde{\sigma}\tilde{s}\)</span>. See <span id="id12">[<a class="reference internal" href="#id56" title="F. Alvarruiz, C. Campos, and J. E. Roman. Thick-restarted joint Lanczos bidiagonalization for the GSVD. J. Comput. Appl. Math., 440:115506, 2024. doi:10.1016/j.cam.2023.115506.">Alvarruiz <em>et al.</em>, 2024</a>]</span> for details.</p>
<p>Similarly, in the HSVD we employ a modified residual</p>
<div class="math notranslate nohighlight" id="equation-eq-hsvd-residual-norm">
<span class="eqno">(24)<a class="headerlink" href="#equation-eq-hsvd-residual-norm" title="Link to this equation">#</a></span>\[\|r_\mathrm{HSVD}\|_2=\left(\|A\tilde{v}-\tilde{\sigma}\tilde{u}\|_2^2+\|A^*\Omega\tilde{u}-\tilde{\sigma}\tilde{\omega}\tilde{v}\|_2^2\right)^{\frac{1}{2}},\]</div>
<p>where <span class="math notranslate nohighlight">\(\tilde\omega\)</span> is the corresponding element of the signature <span class="math notranslate nohighlight">\(\tilde\Omega\)</span> of the definition <a class="reference internal" href="#equation-eq-hsvd">(14)</a>.</p>
</section>
<section id="controlling-and-monitoring-convergence">
<h3>Controlling and Monitoring Convergence<a class="headerlink" href="#controlling-and-monitoring-convergence" title="Link to this heading">#</a></h3>
<p>Similarly to the case of eigensolvers, in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVD.html">SVD</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/SVD/SVDGetIterationNumber.html">SVDGetIterationNumber</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/SVD/SVDSetTolerances.html">SVDSetTolerances</a>()</span></code>. Also, convergence can be monitored with command-line keys <code class="docutils notranslate"><span class="pre">-svd_monitor</span></code>, <code class="docutils notranslate"><span class="pre">-svd_monitor_all</span></code>, <code class="docutils notranslate"><span class="pre">-svd_monitor_conv</span></code>, or graphically with <code class="docutils notranslate"><span class="pre">-svd_monitor</span> <span class="pre">draw::draw_lg</span></code>, or alternatively with <code class="docutils notranslate"><span class="pre">-svd_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.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-svdconvergence">
<caption><span class="caption-text">Available possibilities for the convergence criterion, with <span class="math notranslate nohighlight">\(Z=A\)</span> in the standard SVD or as defined in equation <a class="reference internal" href="#equation-eq-qr">(11)</a> for the GSVD</span><a class="headerlink" href="#tab-svdconvergence" title="Link to this table">#</a></caption>
<thead>
<tr class="row-odd"><th class="head"><p>Convergence criterion</p></th>
<th class="head"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDConv.html">SVDConv</a></span></code></p></th>
<th class="head"><p>Command line key</p></th>
<th class="head"><p>Error bound</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>Absolute</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDConv.html">SVD_CONV_ABS</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_conv_abs</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\|r\|\)</span></p></td>
</tr>
<tr class="row-odd"><td><p>Relative to eigenvalue</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDConv.html">SVD_CONV_REL</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_conv_rel</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\|r\|/|\lambda|\)</span></p></td>
</tr>
<tr class="row-even"><td><p>Relative to matrix norms</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDConv.html">SVD_CONV_NORM</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_conv_norm</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\|r\|/\|Z\|\)</span></p></td>
</tr>
<tr class="row-odd"><td><p>Fixed number of iterations</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDConv.html">SVD_CONV_MAXIT</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_conv_maxit</span></code></p></td>
<td><p><span class="math notranslate nohighlight">\(\infty\)</span></p></td>
</tr>
<tr class="row-even"><td><p>User-defined</p></td>
<td><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDConv.html">SVD_CONV_USER</a></span></code></p></td>
<td><p><code class="docutils notranslate"><span class="pre">-svd_conv_user</span></code></p></td>
<td><p><em>user function</em></p></td>
</tr>
</tbody>
</table>
</div>
<p>As in the case of eigensolvers, the user can choose different convergence tests, based on an error bound obtained from the computed residual norm, <span class="math notranslate nohighlight">\(\|r\|\)</span>. Table <a class="reference internal" href="#tab-svdconvergence">Available possibilities for the convergence criterion</a> lists the available options. It is worth mentioning the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDConv.html">SVD_CONV_MAXIT</a></span></code> convergence criterion, which is a bit special. With this criterion, the solver will not compute any residual norms and will stop with a successful status when the maximum number of iterations is reached. This is intended for the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/SVD/SVDRANDOMIZED.html">SVDRANDOMIZED</a></span></code> solver in cases when a low-rank approximation of a matrix needs to be computed instead of accurate singular vectors.</p>
</section>
<section id="viewing-the-solution">
<h3>Viewing the Solution<a class="headerlink" href="#viewing-the-solution" title="Link to this heading">#</a></h3>
<p>There is support for different kinds of viewers for the solution, as in the case of eigensolvers. One can for instance use <code class="docutils notranslate"><span class="pre">-svd_view_values</span></code>, <code class="docutils notranslate"><span class="pre">-svd_view_vectors</span></code>, <code class="docutils notranslate"><span class="pre">-svd_error_relative</span></code>, or <code class="docutils notranslate"><span class="pre">-svd_converged_reason</span></code>. See the description in section <a class="reference internal" href="eps.html#sec-epsviewers"><span class="std std-ref">Viewing the Solution</span></a>.</p>
<p class="rubric">References</p>
<div class="docutils container" id="id13">
<div role="list" class="citation-list">
<div class="citation" id="id56" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id12">Alv24</a><span class="fn-bracket">]</span></span>
<p>F. Alvarruiz, C. Campos, and J. E. Roman. Thick-restarted joint Lanczos bidiagonalization for the GSVD. <em>J. Comput. Appl. Math.</em>, 440:115506, 2024. <a class="reference external" href="https://doi.org/10.1016/j.cam.2023.115506">doi:10.1016/j.cam.2023.115506</a>.</p>
</div>
<div class="citation" id="id14" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id5">And99</a><span class="fn-bracket">]</span></span>
<p>E. Anderson, Z. Bai, C. Bischof, L. S. Blackford, J. Demmel, J. Dongarra, J. Du Croz, A. Greenbaum, S. Hammarling, A. McKenney, and D. Sorensen. <em>LAPACK Users' Guide</em>. Society for Industrial and Applied Mathematics, Philadelphia, PA, third edition, 1999.</p>
</div>
<div class="citation" id="id15" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id1">Bai00</a><span class="fn-bracket">]</span></span>
<p>Z. Bai, J. Demmel, J. Dongarra, A. Ruhe, and H. van der Vorst, editors. <em>Templates for the Solution of Algebraic Eigenvalue Problems: A Practical Guide</em>. Society for Industrial and Applied Mathematics, Philadelphia, PA, 2000.</p>
</div>
<div class="citation" id="id23" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id2">Bjo96</a><span class="fn-bracket">]</span></span>
<p>Å. Björck. <em>Numerical Methods for Least Squares Problems</em>. Society for Industrial and Applied Mathematics, Philadelphia, 1996.</p>
</div>
<div class="citation" id="id22" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span>Gol96<span class="fn-bracket">]</span></span>
<span class="backrefs">(<a role="doc-backlink" href="#id3">1</a>,<a role="doc-backlink" href="#id4">2</a>)</span>
<p>G. H. Golub and C. F. van Loan. <em>Matrix Computations</em>. The Johns Hopkins University Press, Baltimore, MD, third edition, 1996.</p>
</div>
<div class="citation" id="id32" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id2">Han98</a><span class="fn-bracket">]</span></span>
<p>P. C. Hansen. <em>Rank-Deficient and Discrete Ill-Posed Problems: Numerical Aspects of Linear Inversion</em>. Society for Industrial and Applied Mathematics, Philadelphia, PA, 1998.</p>
</div>
<div class="citation" id="id74" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span>Her07c<span class="fn-bracket">]</span></span>
<span class="backrefs">(<a role="doc-backlink" href="#id9">1</a>,<a role="doc-backlink" href="#id11">2</a>)</span>
<p>V. Hernandez, J. E. Roman, and A. Tomas. Restarted Lanczos bidiagonalization for the SVD in SLEPc. Technical Report STR-8, Universitat Politècnica de València, 2007. URL: <a class="reference external" href="https://slepc.upv.es/documentation">https://slepc.upv.es/documentation</a>.</p>
</div>
<div class="citation" id="id57" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id7">Onn91</a><span class="fn-bracket">]</span></span>
<p>R. Onn, A. O. Steinhardt, and A. Bojanczyk. The hyperbolic singular value decomposition and applications. <em>IEEE Trans. Signal Proces.</em>, 39(7):1575–1588, 1991. <a class="reference external" href="https://doi.org/10.1109/78.134396">doi:10.1109/78.134396</a>.</p>
</div>
</div>
</div>
<p class="rubric">Footnotes</p>
</section>
</section>
</section>
<hr class="footnotes docutils" />
<aside class="footnote-list brackets">
<aside class="footnote brackets" id="hsvd" role="doc-footnote">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id10">1</a><span class="fn-bracket">]</span></span>
<p>The exception is in the HSVD, where <span class="math notranslate nohighlight">\(u_j\)</span> is normalized so that <span class="math notranslate nohighlight">\(U^*\Omega U=\tilde\Omega\)</span>.</p>
</aside>
</aside>
</article>
<footer class="prev-next-footer d-print-none">
<div class="prev-next-area">
<a class="left-prev"
href="st.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">ST: Spectral Transformation</p>
</div>
</a>
<a class="right-next"
href="pep.html"
title="next page">
<div class="prev-next-info">
<p class="prev-next-subtitle">next</p>
<p class="prev-next-title">PEP: Polynomial Eigenvalue Problems</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-svdback">Mathematical Background</a><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-svd">The (Standard) Singular Value Decomposition (SVD)</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#equivalent-eigenvalue-problems">Equivalent Eigenvalue Problems</a></li>
</ul>
</li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-gsvd">The Generalized Singular Value Decomposition (GSVD)</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#id6">Equivalent Eigenvalue Problems</a></li>
</ul>
</li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-hsvd">The Hyperbolic Singular Value Decomposition (HSVD)</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#id8">Equivalent Eigenvalue Problems</a></li>
</ul>
</li>
</ul>
</li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#basic-usage">Basic Usage</a></li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#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-thres">Using a threshold to specify wanted singular values</a></li>
</ul>
</li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#selecting-the-svd-solver">Selecting the SVD 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><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#reliability-of-the-computed-solution">Reliability of the Computed Solution</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#controlling-and-monitoring-convergence">Controlling and Monitoring Convergence</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#viewing-the-solution">Viewing the Solution</a></li>
</ul>
</li>
</ul>
</nav></div>
<div class="sidebar-secondary-item">
<div class="tocsection editthispage">
<a href="https://gitlab.com/slepc/slepc/-/edit/release/doc/source/documentation/manual/svd.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/svd.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>
|