File: st.html

package info (click to toggle)
slepc 3.24.1%2Bdfsg1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 122,028 kB
  • sloc: ansic: 104,353; javascript: 12,732; python: 5,958; f90: 3,312; cpp: 1,528; makefile: 761; xml: 679; sh: 347
file content (1180 lines) | stat: -rw-r--r-- 122,120 bytes parent folder | download
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

<!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>ST: Spectral Transformation &#8212; 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/st';</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="SVD: Singular Value Decomposition" href="svd.html" />
    <link rel="prev" title="EPS: Eigenvalue Problem Solver" href="eps.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 current active"><a class="current reference internal" href="#">ST: Spectral Transformation</a></li>
<li class="toctree-l2"><a class="reference internal" href="svd.html">SVD: Singular Value Decomposition</a></li>
<li class="toctree-l2"><a class="reference internal" href="pep.html">PEP: Polynomial Eigenvalue Problems</a></li>
<li class="toctree-l2"><a class="reference internal" href="nep.html">NEP: Nonlinear Eigenvalue Problems</a></li>
<li class="toctree-l2"><a class="reference internal" href="mfn.html">MFN: Matrix Function</a></li>
<li class="toctree-l2"><a class="reference internal" href="lme.html">LME: Linear Matrix Equation</a></li>
<li class="toctree-l2"><a class="reference internal" href="aux.html">Auxiliary Classes</a></li>
<li class="toctree-l2"><a class="reference internal" href="extra.html">Additional Information</a></li>
</ul>
</details></li>
<li class="toctree-l1 has-children"><a class="reference internal" href="../hands-on/index.html">Hands-on exercises</a><details><summary><span class="toctree-toggle" role="presentation"><i class="fa-solid fa-chevron-down"></i></span></summary><ul>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on0.html">Exercise 0: Hello World</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on1.html">Exercise 1: Standard Symmetric Eigenvalue Problem</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on2.html">Exercise 2: Standard Non-Symmetric Eigenvalue Problem</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on3.html">Exercise 3: Generalized Eigenvalue Problem Stored in a File</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on4.html">Exercise 4: Singular Value Decomposition</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on5.html">Exercise 5: Problem without Explicit Matrix Storage</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on6.html">Exercise 6: Parallel Execution</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on7.html">Exercise 7: Use of Deflation Subspaces</a></li>
<li class="toctree-l2"><a class="reference internal" href="../hands-on/hands-on8.html">Exercise 8: Quadratic Eigenvalue Problem</a></li>
</ul>
</details></li>
<li class="toctree-l1"><a class="reference internal" href="../faq.html">Frequently Asked Questions</a></li>
<li class="toctree-l1"><a class="reference internal" href="../presentations.html">Presentations</a></li>
<li class="toctree-l1"><a class="reference internal" href="../license.html">License</a></li>
</ul>
</div>
</nav></div>
    </div>
  
  
  <div class="sidebar-primary-items__end sidebar-primary__section">
      <div class="sidebar-primary-item">
<div id="ethical-ad-placement"
      class="flat"
      data-ea-publisher="readthedocs"
      data-ea-type="readthedocs-sidebar"
      data-ea-manual="true">
</div></div>
  </div>


      </div>
      
      <main id="main-content" class="bd-main" role="main">
        
        
          <div class="bd-content">
            <div class="bd-article-container">
              
              <div class="bd-header-article d-print-none">
<div class="header-article-items header-article__inner">
  
    <div class="header-article-items__start">
      
        <div class="header-article-item">

<nav aria-label="Breadcrumb" class="d-print-none">
  <ul class="bd-breadcrumbs">
    
    <li class="breadcrumb-item breadcrumb-home">
      <a href="../../index.html" class="nav-link" aria-label="Home">
        <i class="fa-solid fa-home"></i>
      </a>
    </li>
    
    <li class="breadcrumb-item"><a href="../index.html" class="nav-link">Documentation</a></li>
    
    
    <li class="breadcrumb-item"><a href="index.html" class="nav-link">SLEPc Users Manual</a></li>
    
    <li class="breadcrumb-item active" aria-current="page"><span class="ellipsis">ST: Spectral Transformation</span></li>
  </ul>
</nav>
</div>
      
    </div>
  
  
</div>
</div>
              
              
              
                
<div id="searchbox"></div>
                <article class="bd-article">
                  
  <section class="tex2jax_ignore mathjax_ignore" id="st-spectral-transformation">
<span id="ch-st"></span><h1>ST: Spectral Transformation<a class="headerlink" href="#st-spectral-transformation" title="Link to this heading">#</a></h1>
<p>The Spectral Transformation (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code>) is the object that encapsulates the functionality required for acceleration techniques based on the transformation of the spectrum. Most eigensolvers in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> work by applying an operator to a set of vectors and this operator can adopt different forms. The <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> object handles all the different possibilities in a uniform way, so that the solver can proceed without knowing which transformation has been selected. The spectral transformation can be specified at run time, together with related options such as which linear solver to use.</p>
<p>Despite being a rather unrelated concept, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> is also used to handle the preconditioners and correction-equation solvers used in preconditioned eigensolvers such as GD and JD.</p>
<p>The description in this chapter focuses on the use of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> in the context of <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code>. For usage within other solver classes, we will provide further details, e.g., shift-and-invert for polynomial eigenproblems in section <a class="reference internal" href="pep.html#sec-qst"><span class="std std-ref">Spectral Transformation for PEP</span></a>.</p>
<section id="general-description">
<h2>General Description<a class="headerlink" href="#general-description" title="Link to this heading">#</a></h2>
<p>Spectral transformations are powerful tools for adjusting the way in which eigensolvers behave when coping with a problem. The general strategy consists in transforming the original problem into a new one in which eigenvalues are mapped to a new position while eigenvectors remain unchanged. These transformations can be used with several goals in mind:</p>
<ul class="simple">
<li><p>Compute internal eigenvalues. In some applications, the eigenpairs of interest are not the extreme ones (largest magnitude, rightmost, leftmost), but those contained in a certain interval or those closest to a certain value of the complex plane.</p></li>
<li><p>Accelerate convergence. Convergence properties typically depend on how close the eigenvalues are from each other. With some spectral transformations, difficult eigenvalue distributions can be remapped in a more favorable way in terms of convergence.</p></li>
<li><p>Handle some special situations. For instance, in generalized problems when matrix <span class="math notranslate nohighlight">\(B\)</span> is singular, it may be necessary to use a spectral transformation.</p></li>
</ul>
<p>SLEPc separates spectral transformations from solution methods so that any combination of them can be specified by the user. To achieve this, most eigensolvers contained in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> are implemented in such a way that they are independent of which transformation has been selected by the user (the exception are preconditioned solvers, see below). That is, the solver algorithm has to work with a generic operator, whose actual form depends on the transformation used. After convergence, eigenvalues are transformed back appropriately.</p>
<p>For technical details of the transformations described in this chapter, the interested user is referred to <span id="id1">[<a class="reference internal" href="#id28" title="T. Ericsson and A. Ruhe. The spectral transformation Lanczos method for the numerical solution of large sparse generalized symmetric eigenvalue problems. Math. Comp., 35(152):1251–1268, 1980. doi:10.1090/S0025-5718-1980-0583502-2.">Ericsson and Ruhe, 1980</a>, <a class="reference internal" href="#id29" title="K. Meerbergen, A. Spence, and D. Roose. Shift-invert and Cayley transforms for detection of rightmost eigenvalues of nonsymmetric matrices. BIT, 34(3):409–423, 1994. doi:10.1007/BF01935650.">Meerbergen <em>et al.</em>, 1994</a>, <a class="reference internal" href="#id66" title="B. Nour-Omid, B. N. Parlett, T. Ericsson, and P. S. Jensen. How to implement the spectral transformation. Math. Comp., 48(178):663–673, 1987. doi:10.1090/S0025-5718-1987-0878698-5.">Nour-Omid <em>et al.</em>, 1987</a>, <a class="reference internal" href="#id31" title="D. S. Scott. The advantages of inverted operators in Rayleigh-Ritz approximations. SIAM J. Sci. Statist. Comput., 3(1):68–75, 1982. doi:10.1137/0903006.">Scott, 1982</a>]</span>.</p>
<p><strong>Preconditioners</strong>:
As explained in the previous chapter, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> contains preconditioned eigensolvers such as GD or JD. These solvers either apply a preconditioner at a certain step of the computation, or need to solve a correction equation with a preconditioned linear solver. One of the main goals of these solvers is to achieve a similar effect as an inverse-based spectral transformation such as shift-and-invert, but with less computational cost. For this reason, a “preconditioner” spectral transformation has been included in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> object. However, this is just a convenient way of organizing the functionality, since this fake spectral transform cannot be used with non-preconditioned eigensolvers, and conversely preconditioned eigensolvers cannot be used with conventional spectral transformations.</p>
</section>
<section id="basic-usage">
<h2>Basic Usage<a class="headerlink" href="#basic-usage" title="Link to this heading">#</a></h2>
<p>The <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> module is the analog of some PETSc modules such as <a class="reference external" href="https://petsc.org/release/manualpages/PC/PC/" title="(in PETSc v3.24)"><span class="xref std std-doc">PC</span></a>. The user does not usually need to create a stand-alone <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> object explicitly. Instead, every <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object internally sets up an associated <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code>. Therefore, the usual object management methods such as <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STCreate.html">STCreate</a>()</span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STDestroy.html">STDestroy</a>()</span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STView.html">STView</a>()</span></code>, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSetFromOptions.html">STSetFromOptions</a>()</span></code>, are not usually called by the user.</p>
<p>Although the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> context is hidden inside the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object, the user still has control over all the options, by means of the command line, or also inside the program. To allow application programmers to set any of the spectral transformation options directly within the code, the following function is provided to extract the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> context from the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> object:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSGetST.html">EPSGetST</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="o">*</span><span class="n">st</span><span class="p">);</span>
</pre></div>
</div>
<p>After this, one is able to set any options associated with the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> object. For example, to set the value of the shift, <span class="math notranslate nohighlight">\(\sigma\)</span>, the following function is available:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/ST/STSetShift.html">STSetShift</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="n">st</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscScalar/">PetscScalar</a></span><span class="w"> </span><span class="n">shift</span><span class="p">);</span>
</pre></div>
</div>
<p>This can also be done with the command line option <code class="docutils notranslate"><span class="pre">-st_shift</span> <span class="pre">&lt;shift&gt;</span></code>. Note that the argument <code class="docutils notranslate"><span class="pre">shift</span></code> is defined as a <a class="reference external" href="https://petsc.org/release/manualpages/Sys/PetscScalar/" title="(in PETSc v3.24)"><span class="xref std std-doc">PetscScalar</span></a>, and this means that complex shifts are not allowed unless the complex version of SLEPc is used.</p>
<div class="admonition warning">
<p class="admonition-title">Warning</p>
<p>Usually, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSetShift.html">STSetShift</a>()</span></code> is never called from application code, as the value of the shift is taken from the target. In all transformations except <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSHIFT.html">STSHIFT</a></span></code>, there is a direct connection between the target <span class="math notranslate nohighlight">\(\tau\)</span> (described in section <a class="reference internal" href="eps.html#sec-which"><span class="std std-ref">Eigenvalues of Interest</span></a>) and the shift <span class="math notranslate nohighlight">\(\sigma\)</span>, as will be discussed below. The normal usage is that the user sets the target and then <span class="math notranslate nohighlight">\(\sigma\)</span> is set to <span class="math notranslate nohighlight">\(\tau\)</span> automatically (though it is still possible for the user to set a different value of the shift).</p>
</div>
<p>Other object operations are available, which are not usually called by the user. The most important ones are <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STApply.html">STApply</a>()</span></code>, that applies the operator to a vector, and <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSetUp.html">STSetUp</a>()</span></code>, that prepares all the necessary data structures before the solution process starts. The term “operator” refers to one of <span class="math notranslate nohighlight">\(A\)</span>, <span class="math notranslate nohighlight">\(B^{-1}\!A\)</span>, <span class="math notranslate nohighlight">\(A-\sigma I\)</span>, … depending on which kind of spectral transformation is being used.</p>
</section>
<section id="available-transformations">
<h2>Available Transformations<a class="headerlink" href="#available-transformations" title="Link to this heading">#</a></h2>
<p>This section describes the spectral transformations that are provided in SLEPc. As in the case of eigensolvers, the spectral transformation to be used 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/ST/STSetType.html">STSetType</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="n">st</span><span class="p">,</span><span class="n"><a href="../../manualpages/ST/STType.html">STType</a></span><span class="w"> </span><span class="n">type</span><span class="p">);</span>
</pre></div>
</div>
<p>The <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> type can also be set with the command-line option <code class="docutils notranslate"><span class="pre">-st_type</span></code> followed by the name of the method (see table <a class="reference internal" href="#tab-transforms"><span class="std std-ref">Spectral transformations available in the ST package</span></a>). The first five spectral transformations are described in detail in the rest of this section. The last possibility, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSHELL.html">STSHELL</a></span></code>, uses a specific, application-provided spectral transformation. Section <a class="reference internal" href="extra.html#sec-shell"><span class="std std-ref">Extending SLEPc</span></a> in the final chapter describes how to implement one of these transformations.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-transforms">
<caption><span class="caption-text">Spectral transformations available in the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> package</span><a class="headerlink" href="#tab-transforms" title="Link to this table">#</a></caption>
<thead>
<tr class="row-odd"><th class="head text-left"><p>Spectral Transformation</p></th>
<th class="head text-left"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STType.html">STType</a></span></code></p></th>
<th class="head text-left"><p>Options Name</p></th>
<th class="head text-center"><p>Operator</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td class="text-left"><p>Shift of Origin</p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSHIFT.html">STSHIFT</a></span></code></p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre">shift</span></code></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(B^{-1}A-\sigma I\)</span></p></td>
</tr>
<tr class="row-odd"><td class="text-left"><p>Shift-and-invert</p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSINVERT.html">STSINVERT</a></span></code></p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre">sinvert</span></code></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\((A-\sigma B)^{-1}B\)</span></p></td>
</tr>
<tr class="row-even"><td class="text-left"><p>Generalized Cayley</p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STCAYLEY.html">STCAYLEY</a></span></code></p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre">cayley</span></code></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\((A-\sigma B)^{-1}(A+\nu B)\)</span></p></td>
</tr>
<tr class="row-odd"><td class="text-left"><p>Preconditioner</p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STPRECOND.html">STPRECOND</a></span></code></p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre">precond</span></code></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(K^{-1}\approx(A-\sigma B)^{-1}\)</span></p></td>
</tr>
<tr class="row-even"><td class="text-left"><p>Polynomial Filter</p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STFILTER.html">STFILTER</a></span></code></p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre">filter</span></code></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(p(A)\)</span></p></td>
</tr>
<tr class="row-odd"><td class="text-left"><p>Shell Transformation</p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSHELL.html">STSHELL</a></span></code></p></td>
<td class="text-left"><p><code class="docutils notranslate"><span class="pre">shell</span></code></p></td>
<td class="text-center"><p><em>user-defined</em></p></td>
</tr>
</tbody>
</table>
</div>
<p>The last column of table <a class="reference internal" href="#tab-transforms"><span class="std std-ref">Spectral transformations available in the ST package</span></a> shows a general form of the operator used in each case. This generic operator can adopt different particular forms depending on whether the eigenproblem is standard or generalized, or whether the value of the shift (<span class="math notranslate nohighlight">\(\sigma\)</span>) and anti-shift (<span class="math notranslate nohighlight">\(\nu\)</span>) is zero or not. All the possible combinations are listed in table <a class="reference internal" href="#tab-op"><span class="std std-ref">Operators used in each spectral transformation mode</span></a>.</p>
<div class="pst-scrollable-table-container"><table class="table" id="tab-op">
<caption><span class="caption-text">Operators used in each spectral transformation mode</span><a class="headerlink" href="#tab-op" title="Link to this table">#</a></caption>
<thead>
<tr class="row-odd"><th class="head text-left"><p><code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code></p></th>
<th class="head text-left"><p>Choice of <span class="math notranslate nohighlight">\(\sigma,\nu\)</span></p></th>
<th class="head text-center"><p>Standard problem</p></th>
<th class="head text-center"><p>Generalized problem</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td class="text-left"><p><code class="docutils notranslate"><span class="pre">shift</span></code></p></td>
<td class="text-left"><p><span class="math notranslate nohighlight">\(\sigma=0\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(A\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(B^{-1}A\)</span></p></td>
</tr>
<tr class="row-odd"><td class="text-left"><p><code class="docutils notranslate"> </code></p></td>
<td class="text-left"><p><span class="math notranslate nohighlight">\(\sigma\not=0\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(A-\sigma I\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(B^{-1}A-\sigma I\)</span></p></td>
</tr>
<tr class="row-even"><td class="text-left"><p><code class="docutils notranslate"><span class="pre">sinvert</span></code></p></td>
<td class="text-left"><p><span class="math notranslate nohighlight">\(\sigma=0\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(A^{-1}\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(A^{-1}B\)</span></p></td>
</tr>
<tr class="row-odd"><td class="text-left"><p><code class="docutils notranslate"> </code></p></td>
<td class="text-left"><p><span class="math notranslate nohighlight">\(\sigma\not=0\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\((A-\sigma I)^{-1}\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\((A-\sigma B)^{-1}B\)</span></p></td>
</tr>
<tr class="row-even"><td class="text-left"><p><code class="docutils notranslate"><span class="pre">cayley</span></code></p></td>
<td class="text-left"><p><span class="math notranslate nohighlight">\(\sigma\not=0,\nu=0\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\((A-\sigma I)^{-1}A\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\((A-\sigma B)^{-1}A\)</span></p></td>
</tr>
<tr class="row-odd"><td class="text-left"><p><code class="docutils notranslate"> </code></p></td>
<td class="text-left"><p><span class="math notranslate nohighlight">\(\sigma=0,\nu\not=0\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(I+\nu A^{-1}\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(I+\nu A^{-1}B\)</span></p></td>
</tr>
<tr class="row-even"><td class="text-left"><p><code class="docutils notranslate"> </code></p></td>
<td class="text-left"><p><span class="math notranslate nohighlight">\(\sigma\not=0,\nu\not=0\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\((A-\sigma I)^{-1}(A+\nu I)\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\((A-\sigma B)^{-1}(A+\nu B)\)</span></p></td>
</tr>
<tr class="row-odd"><td class="text-left"><p><code class="docutils notranslate"><span class="pre">precond</span></code></p></td>
<td class="text-left"><p><span class="math notranslate nohighlight">\(\sigma=0\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(K^{-1}\approx A^{-1}\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(K^{-1}\approx A^{-1}\)</span></p></td>
</tr>
<tr class="row-even"><td class="text-left"><p><code class="docutils notranslate"> </code></p></td>
<td class="text-left"><p><span class="math notranslate nohighlight">\(\sigma\not=0\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(K^{-1}\approx(A-\sigma I)^{-1}\)</span></p></td>
<td class="text-center"><p><span class="math notranslate nohighlight">\(K^{-1}\approx(A-\sigma B)^{-1}\)</span></p></td>
</tr>
</tbody>
</table>
</div>
<p>The expressions shown in table <a class="reference internal" href="#tab-op"><span class="std std-ref">Operators used in each spectral transformation mode</span></a> are not built explicitly. Instead, the appropriate operations are carried out when applying the operator to a certain vector. The inverses imply the solution of a system of linear equations that is managed by setting up an associated <a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSP/" title="(in PETSc v3.24)"><span class="xref std std-doc">KSP</span></a> object. The user can control the behavior of this object by adjusting the appropriate options, as will be illustrated with examples in section <a class="reference internal" href="#sec-lin"><span class="std std-ref">Solution of Linear Systems</span></a>.</p>
<section id="shift-of-origin">
<h3>Shift of Origin<a class="headerlink" href="#shift-of-origin" title="Link to this heading">#</a></h3>
<p>By default, no spectral transformation is performed. This is equivalent to a shift of origin (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSHIFT.html">STSHIFT</a></span></code>) with <span class="math notranslate nohighlight">\(\sigma=0\)</span>, that is, the first line of table <a class="reference internal" href="#tab-op"><span class="std std-ref">Operators used in each spectral transformation mode</span></a>. The solver works with the original expressions of the eigenvalue problems,</p>
<div class="math notranslate nohighlight" id="equation-eq-std-eig-problem">
<span class="eqno">(1)<a class="headerlink" href="#equation-eq-std-eig-problem" title="Link to this equation">#</a></span>\[Ax=\lambda x,\]</div>
<p>for standard problems, and <span class="math notranslate nohighlight">\(Ax=\lambda Bx\)</span> for generalized ones. Note that this last equation is actually treated internally as</p>
<div class="math notranslate nohighlight" id="equation-eq-gen-eig-problem">
<span class="eqno">(2)<a class="headerlink" href="#equation-eq-gen-eig-problem" title="Link to this equation">#</a></span>\[B^{-1}Ax=\lambda x.\]</div>
<p>When the eigensolver in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> requests the application of the operator to a vector, a matrix-vector multiplication by matrix <span class="math notranslate nohighlight">\(A\)</span> is carried out (in the standard case) or a matrix-vector multiplication by matrix <span class="math notranslate nohighlight">\(A\)</span> followed by a linear system solve with coefficient matrix <span class="math notranslate nohighlight">\(B\)</span> (in the generalized case). Note that in the last case, the operation will fail if matrix <span class="math notranslate nohighlight">\(B\)</span> is singular.</p>
<p>When the shift, <span class="math notranslate nohighlight">\(\sigma\)</span>, is given a value different from the default, 0, the effect is to move the whole spectrum by that exact quantity, <span class="math notranslate nohighlight">\(\sigma\)</span>, which is called <em>shift of origin</em>. To achieve this, the solver works with the shifted matrix, that is, the expressions it has to cope with are</p>
<div class="math notranslate nohighlight" id="equation-eq-std-eig-problem-shift">
<span class="eqno">(3)<a class="headerlink" href="#equation-eq-std-eig-problem-shift" title="Link to this equation">#</a></span>\[(A-\sigma I)x=\theta x,\]</div>
<p>for standard problems, and</p>
<div class="math notranslate nohighlight" id="equation-eq-gen-eig-problem-shift">
<span class="eqno">(4)<a class="headerlink" href="#equation-eq-gen-eig-problem-shift" title="Link to this equation">#</a></span>\[(B^{-1}A-\sigma I) x=\theta x,\]</div>
<p>for generalized ones. The important property that is used is that shifting does not alter the eigenvectors and that it does change the eigenvalues in a simple known way, it shifts them by <span class="math notranslate nohighlight">\(\sigma\)</span>. In both the standard and the generalized problems, the following relation holds</p>
<div class="math notranslate nohighlight" id="equation-eq-eig-problem-shift">
<span class="eqno">(5)<a class="headerlink" href="#equation-eq-eig-problem-shift" title="Link to this equation">#</a></span>\[\theta=\lambda-\sigma.\]</div>
<p>This means that after the solution process, the value <span class="math notranslate nohighlight">\(\sigma\)</span> has to be added<a class="footnote-reference brackets" href="#slepc-v3-5" id="id2" role="doc-noteref"><span class="fn-bracket">[</span>1<span class="fn-bracket">]</span></a> to the computed eigenvalues, <span class="math notranslate nohighlight">\(\theta\)</span>, in order to retrieve the solution of the original problem, <span class="math notranslate nohighlight">\(\lambda\)</span>. This is done by means of the function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STBackTransform.html">STBackTransform</a>()</span></code>, which does not need to be called directly by the user.</p>
</section>
<section id="shift-and-invert">
<h3>Shift-and-invert<a class="headerlink" href="#shift-and-invert" title="Link to this heading">#</a></h3>
<figure class="align-default" id="fig-sinvert">
<img alt="The shift-and-invert spectral transformation" src="../../_images/fig-sinvert.svg" /><figcaption>
<p><span class="caption-text">The shift-and-invert spectral transformation</span><a class="headerlink" href="#fig-sinvert" title="Link to this image">#</a></p>
</figcaption>
</figure>
<p>The shift-and-invert spectral transformation (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSINVERT.html">STSINVERT</a></span></code>) is used to enhance convergence of eigenvalues in the neighborhood of a given value. In this case, the solver deals with the expressions</p>
<div class="math notranslate nohighlight" id="equation-eq-std-eig-problem-shift-and-invert">
<span class="eqno">(6)<a class="headerlink" href="#equation-eq-std-eig-problem-shift-and-invert" title="Link to this equation">#</a></span>\[(A-\sigma I)^{-1}x=\theta x,\]</div>
<div class="math notranslate nohighlight" id="equation-eq-gen-eig-problem-shift-and-invert">
<span class="eqno">(7)<a class="headerlink" href="#equation-eq-gen-eig-problem-shift-and-invert" title="Link to this equation">#</a></span>\[(A-\sigma B)^{-1}B x=\theta x,\]</div>
<p>for standard and generalized problems, respectively. This transformation is effective for finding eigenvalues near <span class="math notranslate nohighlight">\(\sigma\)</span> since the eigenvalues <span class="math notranslate nohighlight">\(\theta\)</span> of the operator that are largest in magnitude correspond to the eigenvalues <span class="math notranslate nohighlight">\(\lambda\)</span> of the original problem that are closest to the shift <span class="math notranslate nohighlight">\(\sigma\)</span> in absolute value, as illustrated in figure <a class="reference internal" href="#fig-sinvert"><span class="std std-ref">The shift-and-invert spectral transformation</span></a> for an example with real eigenvalues. Once the wanted eigenvalues have been found, they may be transformed back to eigenvalues of the original problem. Again, the eigenvectors remain unchanged. In this case, the relation between the eigenvalues of both problems is</p>
<div class="math notranslate nohighlight" id="equation-eq-eig-problem-shift-and-invert">
<span class="eqno">(8)<a class="headerlink" href="#equation-eq-eig-problem-shift-and-invert" title="Link to this equation">#</a></span>\[\theta=1/(\lambda-\sigma).\]</div>
<p>Therefore, after the solution process, the operation to be performed in function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STBackTransform.html">STBackTransform</a>()</span></code> is <span class="math notranslate nohighlight">\(\lambda=\sigma+1/\theta\)</span> for each of the computed eigenvalues.</p>
<p>This spectral transformation is used in the spectrum slicing technique, see section <a class="reference internal" href="#sec-slice"><span class="std std-ref">Spectrum Slicing</span></a>.</p>
</section>
<section id="sec-cayley">
<h3>Cayley<a class="headerlink" href="#sec-cayley" title="Link to this heading">#</a></h3>
<p>The generalized Cayley transform (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STCAYLEY.html">STCAYLEY</a></span></code>) is defined from the expressions</p>
<div class="math notranslate nohighlight" id="equation-eq-std-eig-problem-cayley">
<span class="eqno">(9)<a class="headerlink" href="#equation-eq-std-eig-problem-cayley" title="Link to this equation">#</a></span>\[(A-\sigma I)^{-1}(A+\nu I)x=\theta x,\]</div>
<div class="math notranslate nohighlight" id="equation-eq-gen-eig-problem-cayley">
<span class="eqno">(10)<a class="headerlink" href="#equation-eq-gen-eig-problem-cayley" title="Link to this equation">#</a></span>\[(A-\sigma B)^{-1}(A+\nu B)x=\theta x,\]</div>
<p>for standard and generalized problems, respectively. Sometimes, the term Cayley transform is applied for the particular case in which <span class="math notranslate nohighlight">\(\nu=\sigma\)</span>. This is the default if <span class="math notranslate nohighlight">\(\nu\)</span> is not given a value explicitly. The value of <span class="math notranslate nohighlight">\(\nu\)</span> (the anti-shift) can be set with the following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/ST/STCayleySetAntishift.html">STCayleySetAntishift</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="n">st</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscScalar/">PetscScalar</a></span><span class="w"> </span><span class="n">nu</span><span class="p">);</span>
</pre></div>
</div>
<p>or in the command line with <code class="docutils notranslate"><span class="pre">-st_cayley_antishift</span></code>.</p>
<p>This transformation is mathematically equivalent to shift-and-invert and, therefore, it is effective for finding eigenvalues near <span class="math notranslate nohighlight">\(\sigma\)</span> as well. However, in some situations it is numerically advantageous with respect to shift-and-invert (see <span id="id3">[<a class="reference internal" href="svd.html#id15" title="Z. Bai, J. Demmel, J. Dongarra, A. Ruhe, and H. van der Vorst, editors. Templates for the Solution of Algebraic Eigenvalue Problems: A Practical Guide. Society for Industrial and Applied Mathematics, Philadelphia, PA, 2000.">Bai <em>et al.</em>, 2000</a>, 11.2, <a class="reference internal" href="#id32" title="R. B. Lehoucq and A. G. Salinger. Large-scale eigenvalue calculations for stability analysis of steady flows on massively parallel computers. Int. J. Numer. Methods Fluids, 36:309–327, 2001. doi:10.1002/fld.135.">Lehoucq and Salinger, 2001</a>]</span>).</p>
<p>In this case, the relation between the eigenvalues of both problems is</p>
<div class="math notranslate nohighlight" id="equation-eq-eig-problem-cayley">
<span class="eqno">(11)<a class="headerlink" href="#equation-eq-eig-problem-cayley" title="Link to this equation">#</a></span>\[\theta=(\lambda+\nu)/(\lambda-\sigma).\]</div>
<p>Therefore, after the solution process, the operation to be performed in function <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STBackTransform.html">STBackTransform</a>()</span></code> is <span class="math notranslate nohighlight">\(\lambda=(\theta\sigma+\nu)/(\theta-1)\)</span> for each of the computed eigenvalues.</p>
</section>
<section id="sec-precond">
<h3>Preconditioner<a class="headerlink" href="#sec-precond" title="Link to this heading">#</a></h3>
<p>As mentioned in the introduction of this chapter, the special type <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STPRECOND.html">STPRECOND</a></span></code> is used for handling preconditioners or preconditioned iterative linear solvers, which are used in the context of preconditioned eigensolvers for expanding the subspace. For instance, in the GD solver the so-called correction vector <span class="math notranslate nohighlight">\(d_i\)</span> to be added to the subspace in each iteration is computed as</p>
<div class="math notranslate nohighlight" id="equation-eq-correction-vector">
<span class="eqno">(12)<a class="headerlink" href="#equation-eq-correction-vector" title="Link to this equation">#</a></span>\[d_i=K^{-1}P_i(A-\theta_i B)x_i,\]</div>
<p>where <span class="math notranslate nohighlight">\((\theta_i,x_i)\)</span> is the current approximation of the sought-after eigenpair, and <span class="math notranslate nohighlight">\(P_i\)</span> is a projector involving <span class="math notranslate nohighlight">\(x_i\)</span> and <span class="math notranslate nohighlight">\(K^{-1}x_i\)</span>. In the above expressions, <span class="math notranslate nohighlight">\(K\)</span> is a preconditioner matrix that is built from <span class="math notranslate nohighlight">\(A-\theta_i B\)</span>. However, since <span class="math notranslate nohighlight">\(\theta_i\)</span> changes at each iteration, which would force recomputation of the preconditioner, we opt for using</p>
<div class="math notranslate nohighlight" id="equation-eq-precon">
<span class="eqno">(13)<a class="headerlink" href="#equation-eq-precon" title="Link to this equation">#</a></span>\[K^{-1}\approx (A-\sigma B)^{-1}.\]</div>
<p>Similarly, in the JD eigensolver the expansion of the subspace is carried out by solving a correction equation similar to</p>
<div class="math notranslate nohighlight" id="equation-eq-jd-correction">
<span class="eqno">(14)<a class="headerlink" href="#equation-eq-jd-correction" title="Link to this equation">#</a></span>\[(I-x_ix_i^*)(A-\theta_i B)(I-x_ix_i^*)d_i=-(A-\theta_i B)x_i,\]</div>
<p>where the system is solved approximately with a preconditioned iterative linear solver. For building the preconditioner of this linear system, the projectors <span class="math notranslate nohighlight">\(I-x_ix_i^*\)</span> are ignored, and again it is not recomputed in each iteration. Therefore, the preconditioner is built as in equation <a class="reference internal" href="#equation-eq-precon">(13)</a> as well.</p>
<p>It should be clear from the previous discussion, that <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STPRECOND.html">STPRECOND</a></span></code> does not work in the same way as the rest of spectral transformations. In particular, it does not rely on <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STBackTransform.html">STBackTransform</a>()</span></code>. It is rather a convenient mechanism for handling the preconditioner and linear solver (see examples in section <a class="reference internal" href="#sec-lin"><span class="std std-ref">Solution of Linear Systems</span></a>). The expressions shown in tables <a class="reference internal" href="#tab-transforms"><span class="std std-ref">Spectral transformations available in the ST package</span></a> and <a class="reference internal" href="#tab-op"><span class="std std-ref">Operators used in each spectral transformation mode</span></a> are just a reference to indicate from which matrix the preconditioner is built by default.</p>
<p>There is the possibility that the user overrides the default behavior, that is, to explicitly supply a matrix from which the preconditioner is to be built, with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/ST/STSetPreconditionerMat.html">STSetPreconditionerMat</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="n">st</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">mat</span><span class="p">);</span>
</pre></div>
</div>
<p>The above function can also be used in other spectral transformations such as shift-and-invert in case the user has a cheap approximation <span class="math notranslate nohighlight">\(K\)</span> of the coefficient matrix <span class="math notranslate nohighlight">\(A-\sigma B\)</span>. An alternative is to pass approximations of both <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span> so that <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> builds the preconditioner matrix internally, with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/ST/STSetSplitPreconditioner.html">STSetSplitPreconditioner</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="n">st</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">n</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/Mat/">Mat</a></span><span class="w"> </span><span class="n">Psplit</span><span class="p">[],</span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/MatStructure/">MatStructure</a></span><span class="w"> </span><span class="n">strp</span><span class="p">);</span>
</pre></div>
</div>
<p>Note that preconditioned eigensolvers in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> select <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STPRECOND.html">STPRECOND</a></span></code> by default, so the user does not need to specify it explicitly.</p>
</section>
<section id="sec-filter">
<h3>Polynomial Filtering<a class="headerlink" href="#sec-filter" title="Link to this heading">#</a></h3>
<p>The type <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STFILTER.html">STFILTER</a></span></code> is also special. It is used in the case of standard symmetric (or Hermitian) eigenvalue problems when the eigenvalues of interest are interior to the spectrum and we want to avoid the high cost associated with the matrix factorization of the shift-and-invert spectral transformation. The techniques generically known as <em>polynomial filtering</em> aim at this goal.</p>
<p>The polynomial filtering methods address the eigenvalue problem</p>
<div class="math notranslate nohighlight" id="equation-eq-polyfilt">
<span class="eqno">(15)<a class="headerlink" href="#equation-eq-polyfilt" title="Link to this equation">#</a></span>\[p(A)x=\theta x,\]</div>
<p>where <span class="math notranslate nohighlight">\(p(\cdot)\)</span> is a suitable high-degree polynomial. Once the polynomial is built, the eigensolver relies on <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STApply.html">STApply</a>()</span></code> to compute approximations of the eigenvalues <span class="math notranslate nohighlight">\(\theta\)</span> of the transformed problem. These approximations must be processed in some way in order to recover the <span class="math notranslate nohighlight">\(\lambda\)</span> eigenvalues. Note that in this case there is no <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STBackTransform.html">STBackTransform</a>()</span></code> operation. Details of the method can be found in <span id="id4">[<a class="reference internal" href="../../manualpages/ST/ST_FILTER_FILTLAN.html#id42" title="H. Fang and Y. Saad. A filtered Lanczos procedure for extreme and interior eigenvalue problems. SIAM J. Sci. Comput., 34(4):A2220–A2246, 2012. doi:10.1137/110836535.">Fang and Saad, 2012</a>]</span>.</p>
<p>Currently, SLEPc provides several types of polynomial filtering techniques, which can be selected via <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STFilterSetType.html">STFilterSetType</a>()</span></code>. Note that the external package EVSL also implements polynomial filters to compute all eigenvalues in an interval.</p>
</section>
</section>
<section id="advanced-usage">
<h2>Advanced Usage<a class="headerlink" href="#advanced-usage" title="Link to this heading">#</a></h2>
<p>Using the <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> object is very straightforward. However, when using spectral transformations many things are happening behind the scenes, mainly the solution of systems of linear equations. The user must be aware of what is going on in each case, so that it is possible to guide the solution process in the most beneficial way. This section describes several advanced aspects that can have a considerable impact on efficiency.</p>
<section id="sec-lin">
<h3>Solution of Linear Systems<a class="headerlink" href="#sec-lin" title="Link to this heading">#</a></h3>
<p>In many of the cases shown in table <a class="reference internal" href="#tab-op"><span class="std std-ref">Operators used in each spectral transformation mode</span></a>, the operator contains an inverted matrix, which means that a system of linear equations must be solved whenever the application of the operator to a vector is required. These cases are handled internally by means of a <a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSP/" title="(in PETSc v3.24)"><span class="xref std std-doc">KSP</span></a> object.</p>
<p>In the simplest case, a generalized problem is to be solved with a zero shift. Suppose you run a program that solves a generalized eigenproblem, with default options:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./program
</pre></div>
</div>
<p>In this case, 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> solver creates a <a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSP/" title="(in PETSc v3.24)"><span class="xref std std-doc">KSP</span></a> object whose coefficient matrix is <span class="math notranslate nohighlight">\(B\)</span>. By default, 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 set to use a direct solver, in particular an LU factorization. However, default settings can be changed, as illustrated below.</p>
<p>The following command-line is equivalent to the previous one:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./program<span class="w"> </span>-st_ksp_type<span class="w"> </span>preonly<span class="w"> </span>-st_pc_type<span class="w"> </span>lu
</pre></div>
</div>
<p>The two options specify the type of the linear solver and preconditioner to be used. The <code class="docutils notranslate"><span class="pre">-st_</span></code> prefix indicates that the option corresponds to the linear solver within <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code>. The combination <code class="docutils notranslate"><span class="pre">preonly</span></code><span class="math notranslate nohighlight">\(+\)</span><code class="docutils notranslate"><span class="pre">lu</span></code> instructs to use a direct solver (LU factorization, see <a class="reference external" href="https://petsc.org/release/manual/ksp/" target="_blank">PETSc documentation</a> for details), so this is the same as the default. Adding a new option changes the default behavior, for instance</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./program<span class="w"> </span>-st_ksp_type<span class="w"> </span>preonly<span class="w"> </span>-st_pc_type<span class="w"> </span>lu<span class="w"> </span>-st_pc_factor_mat_solver_type<span class="w"> </span>mumps
</pre></div>
</div>
<p>In this case, an external linear solver package is used (MUMPS, see <a class="reference external" href="https://petsc.org/release/overview/linear_solve_table/" target="_blank">PETSc documentation</a> for other available packages). Note that an external package is required for computing a matrix factorization in parallel, since PETSc itself only provides sequential direct linear solvers.</p>
<p>Instead of a direct linear solver, it is possible to use an iterative solver. This may be necessary in some cases, specially for very large problems. However, the user is warned that using an iterative linear solver makes the overall solution process less robust (see also the discussion of preconditioned eigensolvers below). As an 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>-st_ksp_type<span class="w"> </span>gmres<span class="w"> </span>-st_pc_type<span class="w"> </span>bjacobi<span class="w"> </span>-st_ksp_rtol<span class="w"> </span>1e-9
</pre></div>
</div>
<p>selects the GMRES solver with block Jacobi preconditioning. In the case of iterative solvers, it is important to use an appropriate tolerance, usually slightly more stringent for the linear solves relative to the desired accuracy of the eigenvalue calculation (<span class="math notranslate nohighlight">\(10^{-9}\)</span> in the example, compared to <span class="math notranslate nohighlight">\(10^{-8}\)</span> for the eigensolver).</p>
<p>Although the direct solver approach may seem too costly, note that the factorization is only carried out at the beginning of the eigenvalue calculation and this cost is amortized in each subsequent application of the operator. This is also the case for iterative methods with preconditioners with high-cost set-up such as ILU.</p>
<p>The application programmer is able to set the desired linear systems solver options also from within the code. In order to do this, first the context of the <a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSP/" title="(in PETSc v3.24)"><span class="xref std std-doc">KSP</span></a> object must be retrieved with the following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/ST/STGetKSP.html">STGetKSP</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="n">st</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>The above functionality is also applicable to the other spectral transformations. For instance, for the shift-and-invert technique with <span class="math notranslate nohighlight">\(\tau=10\)</span> using BiCGStab+Jacobi:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./program<span class="w"> </span>-st_type<span class="w"> </span>sinvert<span class="w"> </span>-eps_target<span class="w"> </span><span class="m">10</span><span class="w"> </span>-st_ksp_type<span class="w"> </span>bcgs<span class="w"> </span>-st_pc_type<span class="w"> </span>jacobi
</pre></div>
</div>
<p>In shift-and-invert and Cayley, unless <span class="math notranslate nohighlight">\(\sigma=0\)</span>, the coefficient matrix is not a simple matrix but an expression that can be explicitly constructed or not, depending on the user’s choice. This issue is examined in detail in section <a class="reference internal" href="#sec-explicit"><span class="std std-ref">Explicit Computation of the Coefficient Matrix</span></a> below.</p>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>By default the preconditioner is built from the matrix <span class="math notranslate nohighlight">\(A-\sigma I\)</span> (or <span class="math notranslate nohighlight">\(A-\sigma B\)</span> in generalized problems). In some special situations, the user may want to build the preconditioner from a different matrix, e.g., stemming from simplified versions of <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span>. To do this, one should use the functions mentioned in section <a class="reference internal" href="#sec-precond"><span class="std std-ref">Preconditioner</span></a>.</p>
</div>
<p>In many cases, especially if a shift-and-invert or Cayley transformation is being used, iterative methods may not be well suited for solving linear systems (because of the properties of the coefficient matrix that can be indefinite and ill-conditioned). When using an iterative linear solver, it may be helpful to run with the option <code class="docutils notranslate"><span class="pre">-st_ksp_converged_reason</span></code>, which will display the number of iterations required in each operator application. In extreme cases, the iterative solver fails, so <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSolve.html">EPSSolve</a>()</span></code> aborts with an error</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="n">PETSC</span><span class="w"> </span><span class="n">ERROR</span><span class="o">:</span><span class="w"> </span><span class="n"><a href="https://petsc.org/release/manualpages/KSP/KSP/">KSP</a></span><span class="w"> </span><span class="n">did</span><span class="w"> </span><span class="n">not</span><span class="w"> </span><span class="n">converge</span><span class="w"> </span><span class="p">(</span><span class="n">reason</span><span class="o">=</span><span class="n">DIVERGED_ITS</span><span class="p">)</span><span class="o">!</span>
</pre></div>
</div>
<p>If this happens, it is necessary to use a direct method for solving the linear systems, as explained above.</p>
<section id="the-case-of-preconditioned-eigensolvers">
<h4>The Case of Preconditioned Eigensolvers<a class="headerlink" href="#the-case-of-preconditioned-eigensolvers" title="Link to this heading">#</a></h4>
<p>The <a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSP/" title="(in PETSc v3.24)"><span class="xref std std-doc">KSP</span></a> object contained internally in <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> is also used for applying the preconditioner or solving the correction equation in preconditioned eigensolvers.</p>
<p>The GD eigensolver employs just a preconditioner. Therefore, by default it sets the <a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSP/" title="(in PETSc v3.24)"><span class="xref std std-doc">KSP</span></a> type to <code class="docutils notranslate"><span class="pre">preonly</span></code> (no other <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> is allowed) and the <a class="reference external" href="https://petsc.org/release/manualpages/PC/PC/" title="(in PETSc v3.24)"><span class="xref std std-doc">PC</span></a> type to <code class="docutils notranslate"><span class="pre">jacobi</span></code>. The user may change the preconditioner, for example as</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex5<span class="w"> </span>-eps_type<span class="w"> </span>gd<span class="w"> </span>-st_pc_type<span class="w"> </span>asm
</pre></div>
</div>
<p>The JD eigensolver uses both an iterative linear solver and a preconditioner, so both <code class="docutils notranslate"><span class="pre"><a href="https://petsc.org/release/manualpages/KSP/KSP/">KSP</a></span></code> and <code class="docutils notranslate"><span class="pre"><a href="https://petsc.org/release/manualpages/PC/PC/">PC</a></span></code> are meaningful in this case. It is important to note that, contrary to the ordinary spectral transformations where a direct linear solver is recommended, in JD using an iterative linear solver is usually better than a direct solver. Indeed, the best performance may be achieved with a few iterations of the linear solver (or a large tolerance). For instance, the next example uses JD with GMRES+Jacobi limiting to 10 the number of allowed iterations for the linear solver:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex5<span class="w"> </span>-eps_type<span class="w"> </span>jd<span class="w"> </span>-st_ksp_type<span class="w"> </span>gmres<span class="w"> </span>-st_pc_type<span class="w"> </span>jacobi<span class="w"> </span>-st_ksp_max_it<span class="w"> </span><span class="m">10</span>
</pre></div>
</div>
<p>A discussion on the different options available for the Davidson solvers can be found in <span id="id5">[<a class="reference internal" href="../../manualpages/EPS/EPSGD.html#id40" title="E. Romero and J. E. Roman. A parallel implementation of Davidson methods for large-scale eigenvalue problems in SLEPc. ACM Trans. Math. Software, 40(2):13:1–13:29, 2014. doi:10.1145/2543696.">Romero and Roman, 2014</a>]</span>.</p>
</section>
</section>
<section id="sec-explicit">
<h3>Explicit Computation of the Coefficient Matrix<a class="headerlink" href="#sec-explicit" title="Link to this heading">#</a></h3>
<p>Three possibilities can be distinguished regarding the form of the coefficient matrix of the systems of linear equations associated with the different spectral transformations. The possible coefficient matrices are:</p>
<ul class="simple">
<li><p>Simple: <span class="math notranslate nohighlight">\(B\)</span>.</p></li>
<li><p>Shifted: <span class="math notranslate nohighlight">\(A-\sigma I\)</span>.</p></li>
<li><p>Axpy: <span class="math notranslate nohighlight">\(A-\sigma B\)</span>.</p></li>
</ul>
<p>The first case has already been described and presents no difficulty. In the other two cases, there are three possible approaches:</p>
<dl class="simple myst">
<dt>“<code class="docutils notranslate"><span class="pre">shell</span></code>”</dt><dd><p>To work with the corresponding expression without forming the matrix explicitly. This is achieved by internally setting a matrix-free matrix with <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatCreateShell/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatCreateShell</span></a>().</p>
</dd>
<dt>“<code class="docutils notranslate"><span class="pre">inplace</span></code>”</dt><dd><p>To build the coefficient matrix explicitly. This is done by means of a <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatShift/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatShift</span></a>() or a <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatAXPY/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatAXPY</span></a>() operation, which overwrites matrix <span class="math notranslate nohighlight">\(A\)</span> with the corresponding expression. This alteration of matrix <span class="math notranslate nohighlight">\(A\)</span> is reversed after the eigensolution process has finished.</p>
</dd>
<dt>“<code class="docutils notranslate"><span class="pre">copy</span></code>”</dt><dd><p>To build the matrix explicitly, as in the previous option, but using a working copy of the matrix, that is, without modifying the original matrix <span class="math notranslate nohighlight">\(A\)</span>.</p>
</dd>
</dl>
<p>The default behavior is to build the coefficient matrix explicitly in a copy of <span class="math notranslate nohighlight">\(A\)</span> (option “<code class="docutils notranslate"><span class="pre">copy</span></code>”). The user can change this 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>-st_type<span class="w"> </span>sinvert<span class="w"> </span>-eps_target<span class="w"> </span><span class="m">10</span><span class="w"> </span>-st_ksp_type<span class="w"> </span>cg<span class="w"> </span>-st_pc_type<span class="w"> </span>jacobi<span class="w"> </span>-st_matmode<span class="w"> </span>shell
</pre></div>
</div>
<p>As always, the procedural equivalent is also available for specifying this option in the code of the program:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/ST/STSetMatMode.html">STSetMatMode</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="n">st</span><span class="p">,</span><span class="n"><a href="../../manualpages/ST/STMatMode.html">STMatMode</a></span><span class="w"> </span><span class="n">mode</span><span class="p">);</span>
</pre></div>
</div>
<p>The user must consider which approach is the most appropriate for the particular application. The different options have advantages and drawbacks. The “<code class="docutils notranslate"><span class="pre">shell</span></code>” approach is the simplest one but severely restricts the number of possibilities available for solving the system, in particular most of the PETSc preconditioners would not be available, including direct methods. The only preconditioners that can be used in this case are Jacobi (only if matrices <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span> have the operation <code class="docutils notranslate"><span class="pre">MATOP_GET_DIAGONAL</span></code>) or a user-defined one.</p>
<p>The second approach (“<code class="docutils notranslate"><span class="pre">inplace</span></code>”) can be much faster, specially in the generalized case. A more important advantage of this approach is that, in this case, the linear system solver can be combined with any of the preconditioners available in PETSc, including those which need to access internal matrix data-structures such as ILU. The main drawback is that, in the generalized problem, this approach probably makes sense only in the case that <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span> have the same sparse pattern, because otherwise the function <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatAXPY/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatAXPY</span></a>() might be inefficient. If the user knows that the pattern is the same (or a subset), then this can be specified with the function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/ST/STSetMatStructure.html">STSetMatStructure</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/ST/ST.html">ST</a></span><span class="w"> </span><span class="n">st</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Mat/MatStructure/">MatStructure</a></span><span class="w"> </span><span class="n">str</span><span class="p">);</span>
</pre></div>
</div>
<p>Note that when the value of the shift <span class="math notranslate nohighlight">\(\sigma\)</span> is very close to an eigenvalue, then the linear system will be ill-conditioned and using iterative methods may be problematic. On the other hand, in symmetric definite problems, the coefficient matrix will be indefinite whenever <span class="math notranslate nohighlight">\(\sigma\)</span> is a point in the interior of the spectrum.</p>
<p>The third approach (“<code class="docutils notranslate"><span class="pre">copy</span></code>”) uses more memory but avoids a potential problem that could appear in the “<code class="docutils notranslate"><span class="pre">inplace</span></code>” approach: the recovered matrix might be slightly different from the original one (due to roundoff).</p>
</section>
<section id="sec-symm">
<h3>Preserving the Symmetry in Generalized Eigenproblems<a class="headerlink" href="#sec-symm" title="Link to this heading">#</a></h3>
<p>As mentioned in section <a class="reference internal" href="eps.html#sec-defprob"><span class="std std-ref">Defining the Problem</span></a>, some eigensolvers can exploit symmetry and compute a solution for Hermitian problems with less storage and/or computational cost than other methods. Also, symmetric solvers can be more accurate in some cases. However, in the case of generalized eigenvalue problems in which both <span class="math notranslate nohighlight">\(A\)</span> and <span class="math notranslate nohighlight">\(B\)</span> are symmetric, it happens that, due to the spectral transformation, symmetry is lost since none of the transformed operators <span class="math notranslate nohighlight">\(B^{-1}\!A-\sigma I\)</span>, <span class="math notranslate nohighlight">\((A-\sigma B)^{-1}B\)</span>, etc. is symmetric (the same applies in the Hermitian case for complex matrices).</p>
<p>The solution proposed in SLEPc is based on selecting different kinds of inner products. Currently, we have the following choice of inner products:</p>
<ul class="simple">
<li><p>Standard Hermitian inner product: <span class="math notranslate nohighlight">\(\langle x,y\rangle=x^*y\)</span>.</p></li>
<li><p><span class="math notranslate nohighlight">\(B\)</span>-inner product: <span class="math notranslate nohighlight">\(\langle x,y\rangle_B=x^*\!B\,y\)</span>.</p></li>
</ul>
<p>The second one can be used for preserving the symmetry in symmetric definite generalized problems, as described below. Note that <span class="math notranslate nohighlight">\(\langle x,y\rangle_B\)</span> is a genuine inner product only if <span class="math notranslate nohighlight">\(B\)</span> is symmetric positive definite (for the case of symmetric positive semi-definite <span class="math notranslate nohighlight">\(B\)</span> see section <a class="reference internal" href="#sec-purif"><span class="std std-ref">Purification of Eigenvectors</span></a>).</p>
<p>It can be shown that <span class="math notranslate nohighlight">\(\mathbb{R}^n\)</span> with the <span class="math notranslate nohighlight">\(\langle x,y\rangle_B\)</span> inner product is isomorphic to the Euclidean <span class="math notranslate nohighlight">\(n\)</span>-space <span class="math notranslate nohighlight">\(\mathbb{R}^n\)</span> with the standard Hermitian inner product. This means that if we use <span class="math notranslate nohighlight">\(\langle x,y\rangle_B\)</span> instead of the standard inner product, we are just changing the way lengths and angles are measured, but otherwise all the algebraic properties are maintained and therefore algorithms remain correct. What is interesting to observe is that the transformed operators such as <span class="math notranslate nohighlight">\(B^{-1}\!A\)</span> or <span class="math notranslate nohighlight">\((A-\sigma B)^{-1}B\)</span> are self-adjoint with respect to <span class="math notranslate nohighlight">\(\langle x,y\rangle_B\)</span>.</p>
<figure class="align-default" id="fig-abstr">
<img alt="Abstraction used by SLEPc solvers" src="../../_images/fig-abstr.svg" /><figcaption>
<p><span class="caption-text">Abstraction used by SLEPc solvers</span><a class="headerlink" href="#fig-abstr" title="Link to this image">#</a></p>
</figcaption>
</figure>
<p>Internally, SLEPc operates with the scheme illustrated in figure <a class="reference internal" href="#fig-abstr"><span class="std std-ref">Abstraction used by SLEPc solvers</span></a>. The operations indicated by dashed arrows are implemented as virtual functions. From the user point of view, all the above explanation is transparent. The only thing he/she has to care about is to set the problem type appropriately with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetProblemType.html">EPSSetProblemType</a>()</span></code> (see section <a class="reference internal" href="eps.html#sec-defprob"><span class="std std-ref">Defining the Problem</span></a>). In the case of the Cayley transform, SLEPc is using <span class="math notranslate nohighlight">\(\langle x,y\rangle_{A+\nu B}\)</span> as the inner product for preserving symmetry.</p>
<p>Using the <span class="math notranslate nohighlight">\(B\)</span>-inner product may be attractive also in the non-symmetric case (<span class="math notranslate nohighlight">\(A\)</span> non-symmetric) as described in the next subsection.</p>
<p>Note that the above discussion is not directly applicable to <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STPRECOND.html">STPRECOND</a></span></code> and the preconditioned eigensolvers, in the sense that the goal is not to recover the symmetry of the operator. Still, the <span class="math notranslate nohighlight">\(B\)</span>-inner product is also used in generalized symmetric-definite problems.</p>
</section>
<section id="sec-purif">
<h3>Purification of Eigenvectors<a class="headerlink" href="#sec-purif" title="Link to this heading">#</a></h3>
<p>In generalized eigenproblems, the case of singular <span class="math notranslate nohighlight">\(B\)</span> deserves especial consideration. In this case the default spectral transformation (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSHIFT.html">STSHIFT</a></span></code>) cannot be used since <span class="math notranslate nohighlight">\(B^{-1}\)</span> does not exist.</p>
<p>In shift-and-invert with operator matrix <span class="math notranslate nohighlight">\(T=(A-\sigma B)^{-1}B\)</span>, when <span class="math notranslate nohighlight">\(B\)</span> is singular all the eigenvectors that belong to finite eigenvalues are also eigenvectors of <span class="math notranslate nohighlight">\(T\)</span> and belong to the range of <span class="math notranslate nohighlight">\(T\)</span>, <span class="math notranslate nohighlight">\(\mathcal{R}(T)\)</span>. In this case, the bilinear form <span class="math notranslate nohighlight">\(\langle x,y\rangle_B\)</span> introduced in section <a class="reference internal" href="#sec-symm"><span class="std std-ref">Preserving the Symmetry in Generalized Eigenproblems</span></a> is a semi-inner product and <span class="math notranslate nohighlight">\(\|x\|_B=\sqrt{\langle x,x\rangle_B}\)</span> is a semi-norm. As before, <span class="math notranslate nohighlight">\(T\)</span> is self-adjoint with respect to this inner product since <span class="math notranslate nohighlight">\(B\,T=T^*B\)</span>. Also, <span class="math notranslate nohighlight">\(\langle x,y\rangle_B\)</span> is a true inner product on <span class="math notranslate nohighlight">\(\mathcal{R}(T)\)</span>.</p>
<p>The implication of all this is that, for singular <span class="math notranslate nohighlight">\(B\)</span>, if the <span class="math notranslate nohighlight">\(B\)</span>-inner product is used throughout the eigensolver then, assuming that the initial vector has been forced to lie in <span class="math notranslate nohighlight">\(\mathcal{R}(T)\)</span>, the computed eigenvectors should be correct, i.e., they should belong to <span class="math notranslate nohighlight">\(\mathcal{R}(T)\)</span> as well. Nevertheless, finite precision arithmetic spoils this nice picture, and computed eigenvectors are easily corrupted by components of vectors in the null-space of <span class="math notranslate nohighlight">\(B\)</span>. Additional computation is required for achieving the desired property. This is usually referred to as <em>eigenvector purification</em>.</p>
<p>Although more elaborate purification strategies have been proposed (usually trying to reduce the computational effort, see <span id="id6">[<a class="reference internal" href="#id25" title="K. Meerbergen and A. Spence. Implicitly restarted Arnoldi with purification for the shift-invert transformation. Math. Comp., 66(218):667–689, 1997. doi:10.1090/S0025-5718-97-00844-2.">Meerbergen and Spence, 1997</a>, <a class="reference internal" href="#id66" title="B. Nour-Omid, B. N. Parlett, T. Ericsson, and P. S. Jensen. How to implement the spectral transformation. Math. Comp., 48(178):663–673, 1987. doi:10.1090/S0025-5718-1987-0878698-5.">Nour-Omid <em>et al.</em>, 1987</a>]</span>), the approach in SLEPc is simply to explicitly force the initial vector in the range of <span class="math notranslate nohighlight">\(T\)</span>, with <span class="math notranslate nohighlight">\(v_0\leftarrow Tv_0\)</span>, as well as the computed eigenvectors at the end, <span class="math notranslate nohighlight">\(x_i\leftarrow Tx_i\)</span>. Since this computation can be costly, it can be deactivated if the user knows that <span class="math notranslate nohighlight">\(B\)</span> is non-singular, with:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSSetPurify.html">EPSSetPurify</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscBool/">PetscBool</a></span><span class="w"> </span><span class="n">purify</span><span class="p">);</span>
</pre></div>
</div>
<p>A final comment is that eigenvector corruption happens also in the non-symmetric case. If <span class="math notranslate nohighlight">\(A\)</span> is non-symmetric but <span class="math notranslate nohighlight">\(B\)</span> is symmetric positive semi-definite, then the scheme presented above (<span class="math notranslate nohighlight">\(B\)</span>-inner product together with purification) can still be applied and is generally more successful than the straightforward approach with the standard inner product. For using this scheme in SLEPc, the user has to specify the special problem type <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_PGNHEP.html">EPS_PGNHEP</a></span></code>, see table <a class="reference internal" href="eps.html#tab-ptype"><span class="std std-ref">Problem types considered in EPS</span></a>.</p>
</section>
<section id="sec-slice">
<h3>Spectrum Slicing<a class="headerlink" href="#sec-slice" title="Link to this heading">#</a></h3>
<p>In the context of symmetric-definite generalized eigenvalue problems (<code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS_GHEP.html">EPS_GHEP</a></span></code>) it is often required to compute all eigenvalues contained in a given interval <span class="math notranslate nohighlight">\([a,b]\)</span>. This poses some difficulties, such as:</p>
<ul class="simple">
<li><p>The number of eigenvalues in the interval is not known a priori.</p></li>
<li><p>There might be many eigenvalues, in some applications a significant percentage of the spectrum (20%, say).</p></li>
<li><p>We must make certain that no eigenvalues are missed, and in particular all eigenvalues must be computed with their correct multiplicity.</p></li>
<li><p>In some applications, the interval is open in one end, i.e., either <span class="math notranslate nohighlight">\(a\)</span> or <span class="math notranslate nohighlight">\(b\)</span> can be infinite.</p></li>
</ul>
<p>One possible strategy to solve this problem is to sweep the interval from one end to the other, computing chunks of eigenvalues with a spectral transformation that updates the shift dynamically. This is generally referred to as <em>spectrum slicing</em>. The method implemented in SLEPc is similar to that proposed by <span id="id7">Grimes <em>et al.</em> [<a class="reference internal" href="#id27" title="R. G. Grimes, J. G. Lewis, and H. D. Simon. A shifted block Lanczos algorithm for solving sparse symmetric generalized eigenproblems. SIAM J. Matrix Anal. Appl., 15(1):228–272, 1994. doi:10.1137/S0895479888151111.">1994</a>]</span>, where inertia information is used to validate sub-intervals. Given a symmetric-indefinite triangular factorization</p>
<div class="math notranslate nohighlight" id="equation-eq-symindtri">
<span class="eqno">(16)<a class="headerlink" href="#equation-eq-symindtri" title="Link to this equation">#</a></span>\[A-\sigma B=LDL^T,\]</div>
<p>by Sylvester’s law of inertia we know that the number of eigenvalues on the left of <span class="math notranslate nohighlight">\(\sigma\)</span> is equal to the number of negative eigenvalues of <span class="math notranslate nohighlight">\(D\)</span>,</p>
<div class="math notranslate nohighlight" id="equation-eq-inertia">
<span class="eqno">(17)<a class="headerlink" href="#equation-eq-inertia" title="Link to this equation">#</a></span>\[\nu(A-\sigma B)=\nu(D).\]</div>
<p>A detailed description of the method available in SLEPc can be found in <span id="id8">[<a class="reference internal" href="../../manualpages/EPS/EPSKRYLOVSCHUR.html#id34" title="C. Campos and J. E. Roman. Strategies for spectrum slicing based on restarted Lanczos methods. Numer. Algorithms, 60(2):279–295, 2012. doi:10.1007/s11075-012-9564-z.">Campos and Roman, 2012</a>]</span>. The SLEPc interface hides all the complications of the algorithm. However, the user must be aware of all the restrictions for this technique to be employed:</p>
<ul class="simple">
<li><p>This is currently implemented only in Krylov-Schur.</p></li>
<li><p>The method is based on shift-and-invert, so <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/STSINVERT.html">STSINVERT</a></span></code> must be used. Furthermore, direct linear solvers are required.</p></li>
<li><p>The direct linear solver must provide the matrix inertia (see PETSc’s <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatGetInertia/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatGetInertia</span></a>()).</p></li>
</ul>
<p>An example command-line that sets up all the required options is:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex2<span class="w"> </span>-n<span class="w"> </span><span class="m">50</span><span class="w"> </span>-eps_interval<span class="w"> </span><span class="m">0</span>.4,0.8<span class="w"> </span>-st_type<span class="w"> </span>sinvert<span class="w"> </span>-st_ksp_type<span class="w"> </span>preonly<span class="w"> </span>-st_pc_type<span class="w"> </span>cholesky
</pre></div>
</div>
<p>Note that PETSc’s Cholesky factorization is not parallel, so for doing spectrum slicing in parallel it is required to use an external solver that supports inertia. For example, with MUMPS (see section <a class="reference internal" href="#sec-lin"><span class="std std-ref">Solution of Linear Systems</span></a> on how to use external linear solvers) we would do:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex2<span class="w"> </span>-n<span class="w"> </span><span class="m">50</span><span class="w"> </span>-eps_interval<span class="w"> </span><span class="m">0</span>.4,0.8<span class="w"> </span>-st_type<span class="w"> </span>sinvert<span class="w"> </span>-st_ksp_type<span class="w"> </span>preonly<span class="w"> </span>-st_pc_type<span class="w"> </span>cholesky<span class="w"> </span>-st_pc_factor_mat_solver_type<span class="w"> </span>mumps<span class="w"> </span>-st_mat_mumps_icntl_13<span class="w"> </span><span class="m">1</span>
</pre></div>
</div>
<p>The last option is required by MUMPS to compute the inertia. An alternative is to use SuperLU_DIST, in which case the required options would be:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>./ex2<span class="w"> </span>-n<span class="w"> </span><span class="m">50</span><span class="w"> </span>-eps_interval<span class="w"> </span><span class="m">0</span>.4,0.8<span class="w"> </span>-st_type<span class="w"> </span>sinvert<span class="w"> </span>-st_ksp_type<span class="w"> </span>preonly<span class="w"> </span>-st_pc_type<span class="w"> </span>cholesky<span class="w"> </span>-st_pc_factor_mat_solver_type<span class="w"> </span>superlu_dist<span class="w"> </span>-st_mat_superlu_dist_rowperm<span class="w"> </span>NOROWPERM
</pre></div>
</div>
<p>In the latter example, <a class="reference external" href="https://petsc.org/release/manualpages/Mat/MatSetOption/" title="(in PETSc v3.24)"><span class="xref std std-doc">MatSetOption</span></a>() must be used in both matrices to explicitly state that they are symmetric (or Hermitian in the complex case).</p>
<p>Apart from the above recommendations, the following must be taken into account:</p>
<ul class="simple">
<li><p>The parameters <code class="docutils notranslate"><span class="pre">nev</span></code> and <code class="docutils notranslate"><span class="pre">ncv</span></code> from <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSSetDimensions.html">EPSSetDimensions</a>()</span></code> are determined internally (user-provided values are ignored, and set to the number of eigenvalues in the interval). It is possible for the user to specify the “local” <code class="docutils notranslate"><span class="pre">nev</span></code> and <code class="docutils notranslate"><span class="pre">ncv</span></code> parameters that will be used when computing eigenvalues around each shift, with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPSKrylovSchurSetDimensions.html">EPSKrylovSchurSetDimensions</a>()</span></code>.</p></li>
<li><p>The user can also tune the computation by setting the value of <code class="docutils notranslate"><span class="pre">max_it</span></code>.</p></li>
</ul>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p><strong>Usage with Complex Scalars</strong>:
Some external packages that provide inertia information (MUMPS, Pardiso) do so only in real scalars, but not in the case of complex scalars. Hence, with complex scalars spectrum slicing is available only sequentially (with PETSc’s Cholesky factorization) or via SuperLU_DIST (as in the last example above). An alternative to spectrum slicing is to use the CISS solver with a region enclosing an interval on the real axis, see <span id="id9">[<a class="reference internal" href="../../manualpages/EPS/EPSCISSSetThreshold.html#id66" title="Y. Maeda, T. Sakurai, and J. E. Roman. Contour integral spectrum slicing method in SLEPc. Technical Report STR-11, Universitat Politècnica de València, 2016. URL: https://slepc.upv.es/documentation.">Maeda <em>et al.</em>, 2016</a>]</span> for details.</p>
</div>
<section id="use-of-multiple-communicators">
<h4>Use of Multiple Communicators<a class="headerlink" href="#use-of-multiple-communicators" title="Link to this heading">#</a></h4>
<p>Since spectrum slicing requires direct linear solves, parallel runs may suffer from bad scalability in the sense that increasing the number of MPI processes does not imply a performance gain. For this reason, SLEPc provides the option of using multiple communicators, that is, splitting the initial MPI communicator in several groups, each of them in charge of processing part of the interval.</p>
<p>The multi-communicator setting is activated with a value of <code class="docutils notranslate"><span class="pre">npart</span></code>&gt;1 in the following function:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSKrylovSchurSetPartitions.html">EPSKrylovSchurSetPartitions</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscInt/">PetscInt</a></span><span class="w"> </span><span class="n">npart</span><span class="p">);</span>
</pre></div>
</div>
<p>The interval <span class="math notranslate nohighlight">\([a,b]\)</span> is then divided in <code class="docutils notranslate"><span class="pre">npart</span></code> subintervals of equal size, and the problem of computing all eigenvalues in <span class="math notranslate nohighlight">\([a,b]\)</span> is divided in <code class="docutils notranslate"><span class="pre">npart</span></code> independent subproblems. Each subproblem is solved using only a subset of the initial <span class="math notranslate nohighlight">\(p\)</span> processes, with <span class="math notranslate nohighlight">\(\lceil p/\texttt{npart}\rceil\)</span> processes at most. A final step will gather all computed solutions so that they are available in the whole <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/EPS/EPS.html">EPS</a></span></code> communicator.</p>
<p>The division of the interval in subintervals is done blindly, and this may result in load imbalance if some subintervals contain much more eigenvalues than others. This can be prevented by passing a list of subinterval boundaries, provided that the user has a priori information to roughly determine the eigenvalue distribution:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n"><a href="../../manualpages/EPS/EPSKrylovSchurSetSubintervals.html">EPSKrylovSchurSetSubintervals</a></span><span class="p">(</span><span class="n"><a href="../../manualpages/EPS/EPS.html">EPS</a></span><span class="w"> </span><span class="n">eps</span><span class="p">,</span><span class="n"><a href="https://petsc.org/release/manualpages/Sys/PetscReal/">PetscReal</a></span><span class="w"> </span><span class="n">subint</span><span class="p">[]);</span>
</pre></div>
</div>
<p>An additional benefit of multi-communicator support is that it enables parallel spectrum slicing runs without the need to install a parallel direct solver (MUMPS), by setting the number of partitions equal to the number of MPI processes. The following command-line example uses sequential linear solves in 4 partitions, one process each:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>mpiexec<span class="w"> </span>-n<span class="w"> </span><span class="m">4</span><span class="w"> </span>./ex25<span class="w"> </span>-eps_interval<span class="w"> </span><span class="m">0</span>.4,0.8<span class="w"> </span>-eps_krylovschur_partitions<span class="w"> </span><span class="m">4</span><span class="w"> </span>-st_type<span class="w"> </span>sinvert<span class="w"> </span>-st_ksp_type<span class="w"> </span>preonly<span class="w"> </span>-st_pc_type<span class="w"> </span>cholesky
</pre></div>
</div>
<p>The analog example using MUMPS with 5 processes in each partition:</p>
<div class="highlight-console notranslate"><div class="highlight"><pre><span></span><span class="gp">$ </span>mpiexec<span class="w"> </span>-n<span class="w"> </span><span class="m">20</span><span class="w"> </span>./ex25<span class="w"> </span>-eps_interval<span class="w"> </span><span class="m">0</span>.4,0.8<span class="w"> </span>-eps_krylovschur_partitions<span class="w"> </span><span class="m">4</span><span class="w"> </span>-st_type<span class="w"> </span>sinvert<span class="w"> </span>-st_ksp_type<span class="w"> </span>preonly<span class="w"> </span>-st_pc_type<span class="w"> </span>cholesky<span class="w"> </span>-st_pc_factor_mat_solver_type<span class="w"> </span>mumps<span class="w"> </span>-st_mat_mumps_icntl_13<span class="w"> </span><span class="m">1</span>
</pre></div>
</div>
</section>
</section>
<section id="spectrum-folding">
<h3>Spectrum Folding<a class="headerlink" href="#spectrum-folding" title="Link to this heading">#</a></h3>
<p>In SLEPc versions prior to 3.5, <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code> had another type intended to perform the spectrum folding technique described below. It is no longer available with <code class="docutils notranslate"><span class="pre"><a href="../../manualpages/ST/ST.html">ST</a></span></code>, but it can be implemented directly in application code as illustrated in example <code class="docutils notranslate"><span class="pre">ex24.c</span></code>.</p>
<p>Spectrum folding involves squaring in addition to shifting. This makes sense for standard Hermitian eigenvalue problems, where the transformed problem to be addressed is</p>
<div class="math notranslate nohighlight" id="equation-eq-spectrum-folding">
<span class="eqno">(18)<a class="headerlink" href="#equation-eq-spectrum-folding" title="Link to this equation">#</a></span>\[(A-\sigma I)^2x=\theta x.\]</div>
<p>The following relation holds</p>
<div class="math notranslate nohighlight" id="equation-eq-spectrum-folding-relation">
<span class="eqno">(19)<a class="headerlink" href="#equation-eq-spectrum-folding-relation" title="Link to this equation">#</a></span>\[\theta=(\lambda-\sigma)^2.\]</div>
<p>Note that the mapping between <span class="math notranslate nohighlight">\(\lambda\)</span> and <span class="math notranslate nohighlight">\(\theta\)</span> is not injective, and hence this cannot be considered a true spectral transformation.</p>
<figure class="align-default" id="fig-fold">
<img alt="Illustration of the effect of spectrum folding" src="../../_images/fig-fold.svg" /><figcaption>
<p><span class="caption-text">Illustration of the effect of spectrum folding</span><a class="headerlink" href="#fig-fold" title="Link to this image">#</a></p>
</figcaption>
</figure>
<p>The effect is that the spectrum is folded around the value of <span class="math notranslate nohighlight">\(\sigma\)</span>. Thus, eigenvalues that are closest to the shift become the smallest eigenvalues in the folded spectrum, see figure <a class="reference internal" href="#fig-fold"><span class="std std-ref">Illustration of the effect of spectrum folding</span></a>. For this reason, spectrum folding is commonly used in combination with eigensolvers that compute the smallest eigenvalues, for instance in the context of electronic structure calculations <span id="id10">[<a class="reference internal" href="#id24" title="A. Canning, L. W. Wang, A. Williamson, and A. Zunger. Parallel empirical pseudopotential electronic structure calculations for million atom systems. J. Comput. Phys., 160(1):29–41, 2000. doi:10.1006/jcph.2000.6440.">Canning <em>et al.</em>, 2000</a>]</span>. This transformation can be an effective, low-cost alternative to shift-and-invert.</p>
<p class="rubric">References</p>
<div class="docutils container" id="id11">
<div role="list" class="citation-list">
<div class="citation" id="id13" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id3">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="id41" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id8">Cam12</a><span class="fn-bracket">]</span></span>
<p>C. Campos and J. E. Roman. Strategies for spectrum slicing based on restarted Lanczos methods. <em>Numer. Algorithms</em>, 60(2):279–295, 2012. <a class="reference external" href="https://doi.org/10.1007/s11075-012-9564-z">doi:10.1007/s11075-012-9564-z</a>.</p>
</div>
<div class="citation" id="id24" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id10">Can00</a><span class="fn-bracket">]</span></span>
<p>A. Canning, L. W. Wang, A. Williamson, and A. Zunger. Parallel empirical pseudopotential electronic structure calculations for million atom systems. <em>J. Comput. Phys.</em>, 160(1):29–41, 2000. <a class="reference external" href="https://doi.org/10.1006/jcph.2000.6440">doi:10.1006/jcph.2000.6440</a>.</p>
</div>
<div class="citation" id="id28" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id1">Eri80</a><span class="fn-bracket">]</span></span>
<p>T. Ericsson and A. Ruhe. The spectral transformation Lanczos method for the numerical solution of large sparse generalized symmetric eigenvalue problems. <em>Math. Comp.</em>, 35(152):1251–1268, 1980. <a class="reference external" href="https://doi.org/10.1090/S0025-5718-1980-0583502-2">doi:10.1090/S0025-5718-1980-0583502-2</a>.</p>
</div>
<div class="citation" id="id51" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id4">Fan12</a><span class="fn-bracket">]</span></span>
<p>H. Fang and Y. Saad. A filtered Lanczos procedure for extreme and interior eigenvalue problems. <em>SIAM J. Sci. Comput.</em>, 34(4):A2220–A2246, 2012. <a class="reference external" href="https://doi.org/10.1137/110836535">doi:10.1137/110836535</a>.</p>
</div>
<div class="citation" id="id27" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id7">Gri94</a><span class="fn-bracket">]</span></span>
<p>R. G. Grimes, J. G. Lewis, and H. D. Simon. A shifted block Lanczos algorithm for solving sparse symmetric generalized eigenproblems. <em>SIAM J. Matrix Anal. Appl.</em>, 15(1):228–272, 1994. <a class="reference external" href="https://doi.org/10.1137/S0895479888151111">doi:10.1137/S0895479888151111</a>.</p>
</div>
<div class="citation" id="id32" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id3">Leh01</a><span class="fn-bracket">]</span></span>
<p>R. B. Lehoucq and A. G. Salinger. Large-scale eigenvalue calculations for stability analysis of steady flows on massively parallel computers. <em>Int. J. Numer. Methods Fluids</em>, 36:309–327, 2001. <a class="reference external" href="https://doi.org/10.1002/fld.135">doi:10.1002/fld.135</a>.</p>
</div>
<div class="citation" id="id75" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id9">Mae16</a><span class="fn-bracket">]</span></span>
<p>Y. Maeda, T. Sakurai, and J. E. Roman. Contour integral spectrum slicing method in SLEPc. Technical Report STR-11, Universitat Politècnica de València, 2016. URL: <a class="reference external" href="https://slepc.upv.es/documentation">https://slepc.upv.es/documentation</a>.</p>
</div>
<div class="citation" id="id25" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id6">Mee97</a><span class="fn-bracket">]</span></span>
<p>K. Meerbergen and A. Spence. Implicitly restarted Arnoldi with purification for the shift-invert transformation. <em>Math. Comp.</em>, 66(218):667–689, 1997. <a class="reference external" href="https://doi.org/10.1090/S0025-5718-97-00844-2">doi:10.1090/S0025-5718-97-00844-2</a>.</p>
</div>
<div class="citation" id="id29" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id1">Mee94</a><span class="fn-bracket">]</span></span>
<p>K. Meerbergen, A. Spence, and D. Roose. Shift-invert and Cayley transforms for detection of rightmost eigenvalues of nonsymmetric matrices. <em>BIT</em>, 34(3):409–423, 1994. <a class="reference external" href="https://doi.org/10.1007/BF01935650">doi:10.1007/BF01935650</a>.</p>
</div>
<div class="citation" id="id66" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span>Nou87<span class="fn-bracket">]</span></span>
<span class="backrefs">(<a role="doc-backlink" href="#id1">1</a>,<a role="doc-backlink" href="#id6">2</a>)</span>
<p>B. Nour-Omid, B. N. Parlett, T. Ericsson, and P. S. Jensen. How to implement the spectral transformation. <em>Math. Comp.</em>, 48(178):663–673, 1987. <a class="reference external" href="https://doi.org/10.1090/S0025-5718-1987-0878698-5">doi:10.1090/S0025-5718-1987-0878698-5</a>.</p>
</div>
<div class="citation" id="id49" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id5">Rom14</a><span class="fn-bracket">]</span></span>
<p>E. Romero and J. E. Roman. A parallel implementation of Davidson methods for large-scale eigenvalue problems in SLEPc. <em>ACM Trans. Math. Software</em>, 40(2):13:1–13:29, 2014. <a class="reference external" href="https://doi.org/10.1145/2543696">doi:10.1145/2543696</a>.</p>
</div>
<div class="citation" id="id31" role="doc-biblioentry">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id1">Sco82</a><span class="fn-bracket">]</span></span>
<p>D. S. Scott. The advantages of inverted operators in Rayleigh-Ritz approximations. <em>SIAM J. Sci. Statist. Comput.</em>, 3(1):68–75, 1982. <a class="reference external" href="https://doi.org/10.1137/0903006">doi:10.1137/0903006</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="slepc-v3-5" role="doc-footnote">
<span class="label"><span class="fn-bracket">[</span><a role="doc-backlink" href="#id2">1</a><span class="fn-bracket">]</span></span>
<p>Note that the sign changed in SLEPc 3.5 with respect to previous versions.</p>
</aside>
</aside>


                </article>
              
              
              
              
              
                <footer class="prev-next-footer d-print-none">
                  
<div class="prev-next-area">
    <a class="left-prev"
       href="eps.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">EPS: Eigenvalue Problem Solver</p>
      </div>
    </a>
    <a class="right-next"
       href="svd.html"
       title="next page">
      <div class="prev-next-info">
        <p class="prev-next-subtitle">next</p>
        <p class="prev-next-title">SVD: Singular Value Decomposition</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="#general-description">General Description</a></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="#available-transformations">Available Transformations</a><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#shift-of-origin">Shift of Origin</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#shift-and-invert">Shift-and-invert</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-cayley">Cayley</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-precond">Preconditioner</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-filter">Polynomial Filtering</a></li>
</ul>
</li>
<li class="toc-h2 nav-item toc-entry"><a class="reference internal nav-link" href="#advanced-usage">Advanced Usage</a><ul class="nav section-nav flex-column">
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-lin">Solution of Linear Systems</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#the-case-of-preconditioned-eigensolvers">The Case of Preconditioned Eigensolvers</a></li>
</ul>
</li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-explicit">Explicit Computation of the Coefficient Matrix</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-symm">Preserving the Symmetry in Generalized Eigenproblems</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-purif">Purification of Eigenvectors</a></li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#sec-slice">Spectrum Slicing</a><ul class="nav section-nav flex-column">
<li class="toc-h4 nav-item toc-entry"><a class="reference internal nav-link" href="#use-of-multiple-communicators">Use of Multiple Communicators</a></li>
</ul>
</li>
<li class="toc-h3 nav-item toc-entry"><a class="reference internal nav-link" href="#spectrum-folding">Spectrum Folding</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/st.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/st.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>