-
Notifications
You must be signed in to change notification settings - Fork 40
/
primmef77.html
939 lines (787 loc) · 106 KB
/
primmef77.html
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
<!DOCTYPE html>
<html class="writer-html5" lang="en" >
<head>
<meta charset="utf-8">
<meta name="viewport" content="width=device-width, initial-scale=1.0">
<title>FORTRAN 77 Library Interface — PRIMME 3.2 documentation</title>
<link rel="stylesheet" href="_static/css/theme.css" type="text/css" />
<link rel="stylesheet" href="_static/pygments.css" type="text/css" />
<!--[if lt IE 9]>
<script src="_static/js/html5shiv.min.js"></script>
<![endif]-->
<script type="text/javascript" id="documentation_options" data-url_root="./" src="_static/documentation_options.js"></script>
<script src="_static/jquery.js"></script>
<script src="_static/underscore.js"></script>
<script src="_static/doctools.js"></script>
<script src="_static/language_data.js"></script>
<script async="async" src="https://cdnjs.cloudflare.com/ajax/libs/mathjax/2.7.7/latest.js?config=TeX-AMS-MML_HTMLorMML"></script>
<script type="text/javascript" src="_static/js/theme.js"></script>
<link rel="index" title="Index" href="genindex.html" />
<link rel="search" title="Search" href="search.html" />
<link rel="next" title="FORTRAN 90 Library Interface" href="primmef90.html" />
<link rel="prev" title="C Library Interface" href="primmec.html" />
</head>
<body class="wy-body-for-nav">
<div class="wy-grid-for-nav">
<nav data-toggle="wy-nav-shift" class="wy-nav-side">
<div class="wy-side-scroll">
<div class="wy-side-nav-search" style="background: white" >
<a href="readme.html">
<img src="_static/logo.png" class="logo" alt="Logo"/>
</a>
<div class="version">
3.2
</div>
<div role="search">
<form id="rtd-search-form" class="wy-form" action="search.html" method="get">
<input type="text" name="q" placeholder="Search docs" />
<input type="hidden" name="check_keywords" value="yes" />
<input type="hidden" name="area" value="default" />
</form>
</div>
</div>
<div class="wy-menu wy-menu-vertical" data-spy="affix" role="navigation" aria-label="main navigation">
<ul class="current">
<li class="toctree-l1"><a class="reference internal" href="intro.html">PRIMME: PReconditioned Iterative MultiMethod Eigensolver</a></li>
<li class="toctree-l1 current"><a class="reference internal" href="apieigs.html">Eigenvalue Problems</a><ul class="current">
<li class="toctree-l2"><a class="reference internal" href="primmec.html">C Library Interface</a></li>
<li class="toctree-l2 current"><a class="current reference internal" href="#">FORTRAN 77 Library Interface</a><ul>
<li class="toctree-l3"><a class="reference internal" href="#primme-initialize-f77">primme_initialize_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#primme-set-method-f77">primme_set_method_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#primme-free-f77">primme_free_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#sprimme-f77">sprimme_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#cprimme-f77">cprimme_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#dprimme-f77">dprimme_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#zprimme-f77">zprimme_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#magma-sprimme-f77">magma_sprimme_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#magma-cprimme-f77">magma_cprimme_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#magma-dprimme-f77">magma_dprimme_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#magma-zprimme-f77">magma_zprimme_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#cprimme-normal-f77">cprimme_normal_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#zprimme-normal-f77">zprimme_normal_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#magma-cprimme-normal-f77">magma_cprimme_normal_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#magma-zprimme-normal-f77">magma_zprimme_normal_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#primme-set-member-f77">primme_set_member_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#primmetop-get-member-f77">primmetop_get_member_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#primmetop-get-prec-shift-f77">primmetop_get_prec_shift_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#primme-get-member-f77">primme_get_member_f77</a></li>
<li class="toctree-l3"><a class="reference internal" href="#primme-get-prec-shift-f77">primme_get_prec_shift_f77</a></li>
</ul>
</li>
<li class="toctree-l2"><a class="reference internal" href="primmef90.html">FORTRAN 90 Library Interface</a></li>
<li class="toctree-l2"><a class="reference internal" href="pyeigsh.html">Python Interface</a></li>
<li class="toctree-l2"><a class="reference internal" href="mateigs.html">MATLAB Interface</a></li>
<li class="toctree-l2"><a class="reference internal" href="appendix.html">Parameters Description</a></li>
<li class="toctree-l2"><a class="reference internal" href="appendix.html#preset-methods">Preset Methods</a></li>
<li class="toctree-l2"><a class="reference internal" href="appendix.html#error-codes">Error Codes</a></li>
</ul>
</li>
<li class="toctree-l1"><a class="reference internal" href="apisvds.html">Singular Value Problems</a></li>
</ul>
</div>
</div>
</nav>
<section data-toggle="wy-nav-shift" class="wy-nav-content-wrap">
<nav class="wy-nav-top" aria-label="top navigation">
<i data-toggle="wy-nav-top" class="fa fa-bars"></i>
<a href="readme.html">PRIMME</a>
</nav>
<div class="wy-nav-content">
<div class="rst-content">
<div role="navigation" aria-label="breadcrumbs navigation">
<ul class="wy-breadcrumbs">
<li><a href="readme.html" class="icon icon-home"></a> »</li>
<li><a href="apieigs.html">Eigenvalue Problems</a> »</li>
<li>FORTRAN 77 Library Interface</li>
<li class="wy-breadcrumbs-aside">
<a href="_sources/primmef77.rst.txt" rel="nofollow"> View page source</a>
</li>
</ul>
<div class="rst-breadcrumbs-buttons" role="navigation" aria-label="breadcrumb navigation">
<a href="primmef90.html" class="btn btn-neutral float-right" title="FORTRAN 90 Library Interface" accesskey="n">Next <span class="fa fa-arrow-circle-right"></span></a>
<a href="primmec.html" class="btn btn-neutral float-left" title="C Library Interface" accesskey="p"><span class="fa fa-arrow-circle-left"></span> Previous</a>
</div>
<hr/>
</div>
<div role="main" class="document" itemscope="itemscope" itemtype="http://schema.org/Article">
<div itemprop="articleBody">
<div class="section" id="fortran-77-library-interface">
<h1>FORTRAN 77 Library Interface<a class="headerlink" href="#fortran-77-library-interface" title="Permalink to this headline">¶</a></h1>
<p>The next enumerations and functions are declared in <code class="docutils literal notranslate"><span class="pre">primme_f77.h</span></code>.</p>
<dl class="type">
<dt id="f/_/ptr">
<em class="property">type </em><code class="sig-name descname">ptr</code><a class="headerlink" href="#f/_/ptr" title="Permalink to this definition">¶</a></dt>
<dd><p>Fortran datatype with the same size as a pointer.
Use <code class="docutils literal notranslate"><span class="pre">integer*4</span></code> when compiling in 32 bits and <code class="docutils literal notranslate"><span class="pre">integer*8</span></code> in 64 bits.</p>
</dd></dl>
<div class="section" id="primme-initialize-f77">
<h2>primme_initialize_f77<a class="headerlink" href="#primme-initialize-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/primme_initialize_f77">
<em class="property">subroutine </em><code class="sig-name descname">primme_initialize_f77</code><span class="sig-paren">(</span><em class="sig-param">primme</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/primme_initialize_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Allocate and initialize a PRIMME parameters structure to the default values.</p>
<p>After calling <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a> (or a variant), call <a class="reference internal" href="#f/_/primme_free_f77" title="f/_/primme_free_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">primme_free_f77()</span></code></a> to release allocated resources by PRIMME.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (output) parameters structure.</p>
</dd>
</dl>
</dd></dl>
</div>
<div class="section" id="primme-set-method-f77">
<h2>primme_set_method_f77<a class="headerlink" href="#primme-set-method-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/primme_set_method_f77">
<em class="property">subroutine </em><code class="sig-name descname">primme_set_method_f77</code><span class="sig-paren">(</span><em class="sig-param">method</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/primme_set_method_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Set PRIMME parameters to one of the preset configurations.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>method</strong><em> [</em><em>integer</em><em>]</em> :: <p>(input) preset configuration. One of:</p>
<div class="line-block">
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_DYNAMIC</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_DEFAULT_MIN_TIME</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_DEFAULT_MIN_MATVECS</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_Arnoldi</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_GD</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_GD_plusK</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_GD_Olsen_plusK</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_JD_Olsen_plusK</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_RQI</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_JDQR</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_JDQMR</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_JDQMR_ETol</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_STEEPEST_DESCENT</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_LOBPCG_OrthoBasis</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_LOBPCG_OrthoBasis_Window</span></code></div>
</div>
<p>See <a class="reference internal" href="appendix.html#c.primme_preset_method" title="primme_preset_method"><code class="xref c c-type docutils literal notranslate"><span class="pre">primme_preset_method</span></code></a>.</p>
</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) if 0, successful; if negative, something went wrong.</p></li>
</ul>
</dd>
</dl>
</dd></dl>
</div>
<div class="section" id="primme-free-f77">
<h2>primme_free_f77<a class="headerlink" href="#primme-free-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/primme_free_f77">
<em class="property">subroutine </em><code class="sig-name descname">primme_free_f77</code><span class="sig-paren">(</span><em class="sig-param">primme</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/primme_free_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Free memory allocated by PRIMME and delete all values set.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input/output) parameters structure.</p>
</dd>
</dl>
</dd></dl>
</div>
<div class="section" id="sprimme-f77">
<h2>sprimme_f77<a class="headerlink" href="#sprimme-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/sprimme_f77">
<em class="property">subroutine </em><code class="sig-name descname">sprimme_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/sprimme_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a real symmetric standard or generalized eigenproblem.</p>
<p>All arrays should be hosted on CPU. The computations are performed on CPU (see <a class="reference internal" href="#f/_/magma_sprimme_f77" title="f/_/magma_sprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">magma_sprimme_f77()</span></code></a> for using GPUs).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>real</em><em>]</em> :: (input/output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 2.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="cprimme-f77">
<h2>cprimme_f77<a class="headerlink" href="#cprimme-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/cprimme_f77">
<em class="property">subroutine </em><code class="sig-name descname">cprimme_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/cprimme_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a Hermitian standard or generalized eigenproblem.</p>
<p>All arrays should be hosted on CPU. The computations are performed on CPU (see <a class="reference internal" href="#f/_/magma_cprimme_f77" title="f/_/magma_cprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">magma_cprimme_f77()</span></code></a> for using GPUs).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>complex real</em><em>]</em> :: (input/output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 2.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="dprimme-f77">
<h2>dprimme_f77<a class="headerlink" href="#dprimme-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/dprimme_f77">
<em class="property">subroutine </em><code class="sig-name descname">dprimme_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/dprimme_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a real symmetric standard or generalized eigenproblem.</p>
<p>All arrays should be hosted on CPU. The computations are performed on CPU (see <a class="reference internal" href="#f/_/magma_dprimme_f77" title="f/_/magma_dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">magma_dprimme_f77()</span></code></a> for using GPUs).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (input/output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>On input, <code class="docutils literal notranslate"><span class="pre">evecs</span></code> should start with the content of the <a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> vectors,
followed by the <a class="reference internal" href="appendix.html#c.primme_params.initSize" title="primme_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a> vectors.</p>
<p>On return, the i-th eigenvector starts at evecs(( <a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + i - 1)* <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a> + 1), with value <cite>evals(i)</cite> and associated residual 2-norm <cite>resNorms(i)</cite>.
The first vector has index i=1. The number of eigenpairs marked as converged (see <a class="reference internal" href="appendix.html#c.primme_params.eps" title="primme_params.eps"><code class="xref c c-member docutils literal notranslate"><span class="pre">eps</span></code></a>) is returned on <a class="reference internal" href="appendix.html#c.primme_params.initSize" title="primme_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a>. Since version 4.0, if the returned error code is <cite>PRIMME_MAIN_ITER_FAILURE</cite>, PRIMME may return also unconverged eigenpairs and its residual norms in <cite>evecs</cite>, <cite>evals</cite>, and <cite>resNorms</cite> starting at i = <a class="reference internal" href="appendix.html#c.primme_params.initSize" title="primme_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a> + 1 and going up to either <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> or the last <cite>resNorms(i)</cite> with non-negative value.</p>
<p>All internal operations are performed at the same precision than <code class="docutils literal notranslate"><span class="pre">evecs</span></code> unless the user sets <a class="reference internal" href="appendix.html#c.primme_params.internalPrecision" title="primme_params.internalPrecision"><code class="xref c c-member docutils literal notranslate"><span class="pre">internalPrecision</span></code></a> otherwise.</p>
<p>The type and precision of the callbacks is also the same as <cite>evecs</cite>. Although this can be changed. See details for <a class="reference internal" href="appendix.html#c.primme_params.matrixMatvec" title="primme_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a>, <a class="reference internal" href="appendix.html#c.primme_params.massMatrixMatvec" title="primme_params.massMatrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">massMatrixMatvec</span></code></a>, <a class="reference internal" href="appendix.html#c.primme_params.applyPreconditioner" title="primme_params.applyPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">applyPreconditioner</span></code></a>, <a class="reference internal" href="appendix.html#c.primme_params.globalSumReal" title="primme_params.globalSumReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">globalSumReal</span></code></a>, <a class="reference internal" href="appendix.html#c.primme_params.broadcastReal" title="primme_params.broadcastReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">broadcastReal</span></code></a>, and <a class="reference internal" href="appendix.html#c.primme_params.convTestFun" title="primme_params.convTestFun"><code class="xref c c-member docutils literal notranslate"><span class="pre">convTestFun</span></code></a>.</p>
</dd></dl>
</div>
<div class="section" id="zprimme-f77">
<h2>zprimme_f77<a class="headerlink" href="#zprimme-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/zprimme_f77">
<em class="property">subroutine </em><code class="sig-name descname">zprimme_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/zprimme_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a Hermitian standard or generalized eigenproblem.</p>
<p>All arrays should be hosted on CPU. The computations are performed on CPU (see <a class="reference internal" href="#f/_/magma_zprimme_f77" title="f/_/magma_zprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">magma_zprimme_f77()</span></code></a> for using GPUs).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>complex double precision</em><em>]</em> :: (input/output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
</dd></dl>
</div>
<div class="section" id="magma-sprimme-f77">
<h2>magma_sprimme_f77<a class="headerlink" href="#magma-sprimme-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/magma_sprimme_f77">
<em class="property">subroutine </em><code class="sig-name descname">magma_sprimme_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/magma_sprimme_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a real symmetric standard or generalized eigenproblem.</p>
<p>Most of the computations are performed on GPU (see <a class="reference internal" href="#f/_/sprimme_f77" title="f/_/sprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">sprimme_f77()</span></code></a> for using only the CPU).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>real</em><em>]</em> :: (input/output) GPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="magma-cprimme-f77">
<h2>magma_cprimme_f77<a class="headerlink" href="#magma-cprimme-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/magma_cprimme_f77">
<em class="property">subroutine </em><code class="sig-name descname">magma_cprimme_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/magma_cprimme_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a Hermitian standard or generalized eigenproblem.</p>
<p>Most of the computations are performed on GPU (see <a class="reference internal" href="#f/_/cprimme_f77" title="f/_/cprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">cprimme_f77()</span></code></a> for using only the CPU).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>complex real</em><em>]</em> :: (input/output) GPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="magma-dprimme-f77">
<h2>magma_dprimme_f77<a class="headerlink" href="#magma-dprimme-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/magma_dprimme_f77">
<em class="property">subroutine </em><code class="sig-name descname">magma_dprimme_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/magma_dprimme_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a real symmetric standard or generalized eigenproblem.</p>
<p>Most of the computations are performed on GPU (see <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a> for using only the CPU).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (input/output) GPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="magma-zprimme-f77">
<h2>magma_zprimme_f77<a class="headerlink" href="#magma-zprimme-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/magma_zprimme_f77">
<em class="property">subroutine </em><code class="sig-name descname">magma_zprimme_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/magma_zprimme_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a Hermitian standard or generalized eigenproblem.</p>
<p>Most of the computations are performed on GPU (see <a class="reference internal" href="#f/_/zprimme_f77" title="f/_/zprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">zprimme_f77()</span></code></a> for using only the CPU).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>complex double precision</em><em>]</em> :: (input/output) GPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="cprimme-normal-f77">
<h2>cprimme_normal_f77<a class="headerlink" href="#cprimme-normal-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/cprimme_normal_f77">
<em class="property">subroutine </em><code class="sig-name descname">cprimme_normal_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/cprimme_normal_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a normal standard eigenproblem, which may not be Hermitian.</p>
<p>All arrays should be hosted on CPU. The computations are performed on CPU (see <a class="reference internal" href="#f/_/magma_cprimme_normal_f77" title="f/_/magma_cprimme_normal_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">magma_cprimme_normal_f77()</span></code></a> for using GPUs).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>complex real</em><em>]</em> :: (input/output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="zprimme-normal-f77">
<h2>zprimme_normal_f77<a class="headerlink" href="#zprimme-normal-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/zprimme_normal_f77">
<em class="property">subroutine </em><code class="sig-name descname">zprimme_normal_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/zprimme_normal_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a normal standard eigenproblem, which may not be Hermitian.</p>
<p>All arrays should be hosted on CPU. The computations are performed on CPU (see <a class="reference internal" href="#f/_/magma_zprimme_normal_f77" title="f/_/magma_zprimme_normal_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">magma_zprimme_normal_f77()</span></code></a> for using GPUs).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>complex double precision</em><em>]</em> :: (input/output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="magma-cprimme-normal-f77">
<h2>magma_cprimme_normal_f77<a class="headerlink" href="#magma-cprimme-normal-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/magma_cprimme_normal_f77">
<em class="property">subroutine </em><code class="sig-name descname">magma_cprimme_normal_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/magma_cprimme_normal_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a normal standard eigenproblem, which may not be Hermitian.</p>
<p>Most of the arrays are stored on GPU, and also most of the computations are done on GPU (see <a class="reference internal" href="#f/_/cprimme_normal_f77" title="f/_/cprimme_normal_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">cprimme_normal_f77()</span></code></a> for using only the CPU).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>complex real</em><em>]</em> :: (input/output) GPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>real</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="magma-zprimme-normal-f77">
<h2>magma_zprimme_normal_f77<a class="headerlink" href="#magma-zprimme-normal-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/magma_zprimme_normal_f77">
<em class="property">subroutine </em><code class="sig-name descname">magma_zprimme_normal_f77</code><span class="sig-paren">(</span><em class="sig-param">evals</em>, <em class="sig-param">evecs</em>, <em class="sig-param">resNorms</em>, <em class="sig-param">primme</em>, <em class="sig-param">ierr</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/magma_zprimme_normal_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Solve a normal standard eigenproblem, which may not be Hermitian.</p>
<p>Most of the arrays are stored on GPU, and also most of the computations are done on GPU (see <a class="reference internal" href="#f/_/zprimme_normal_f77" title="f/_/zprimme_normal_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">zprimme_normal_f77()</span></code></a> for using only the CPU).</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>evals</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
computed eigenvalues; all parallel calls return the same value in this array.</p></li>
<li><p><strong>evecs</strong> (*)<em> [</em><em>complex double precision</em><em>]</em> :: (input/output) GPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> times (<a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code></a> + <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a>) with leading dimension <a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">ldevecs</span></code></a>
to store column-wise the (local part for this process of the) computed eigenvectors.</p></li>
<li><p><strong>resNorms</strong> (*)<em> [</em><em>double precision</em><em>]</em> :: (output) CPU array at least of size <a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numEvals</span></code></a> to store the
residual norms of the computed eigenpairs; all parallel calls return the same value in this array.</p></li>
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>ierr</strong><em> [</em><em>integer</em><em>]</em> :: (output) error indicator; see <a class="reference internal" href="appendix.html#error-codes"><span class="std std-ref">Error Codes</span></a>.</p></li>
</ul>
</dd>
</dl>
<p>Further descriptions of <cite>evals</cite>, <cite>evecs</cite>, and <cite>resNorms</cite> on notes in subroutine <a class="reference internal" href="#f/_/dprimme_f77" title="f/_/dprimme_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">dprimme_f77()</span></code></a>.</p>
<div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
</div>
<div class="section" id="primme-set-member-f77">
<h2>primme_set_member_f77<a class="headerlink" href="#primme-set-member-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/primme_set_member_f77">
<em class="property">subroutine </em><code class="sig-name descname">primme_set_member_f77</code><span class="sig-paren">(</span><em class="sig-param">primme</em>, <em class="sig-param">label</em>, <em class="sig-param">value</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/primme_set_member_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Set a value in some field of the parameter structure.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>label</strong><em> [</em><em>integer</em><em>]</em> :: <p>field where to set value. One of:</p>
<div class="line-block">
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.n" title="primme_params.n"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_n</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.matrixMatvec" title="primme_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_matrixMatvec</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.matrixMatvec_type" title="primme_params.matrixMatvec_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_matrixMatvec_type</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.applyPreconditioner" title="primme_params.applyPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_applyPreconditioner</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.applyPreconditioner_type" title="primme_params.applyPreconditioner_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_applyPreconditioner_type</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.massMatrixMatvec" title="primme_params.massMatrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_massMatrixMatvec</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.massMatrixMatvec_type" title="primme_params.massMatrixMatvec_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_massMatrixMatvec_type</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.numProcs" title="primme_params.numProcs"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_numProcs</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.procID" title="primme_params.procID"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_procID</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.commInfo" title="primme_params.commInfo"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_commInfo</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.nLocal" title="primme_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_nLocal</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.globalSumReal" title="primme_params.globalSumReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_globalSumReal</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.globalSumReal_type" title="primme_params.globalSumReal_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_globalSumReal_type</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.broadcastReal" title="primme_params.broadcastReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_broadcastReal</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.broadcastReal_type" title="primme_params.broadcastReal_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_broadcastReal_type</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.numEvals" title="primme_params.numEvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_numEvals</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.target" title="primme_params.target"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_target</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.numTargetShifts" title="primme_params.numTargetShifts"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_numTargetShifts</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.targetShifts" title="primme_params.targetShifts"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_targetShifts</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.locking" title="primme_params.locking"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_locking</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.initSize" title="primme_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_initSize</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.numOrthoConst" title="primme_params.numOrthoConst"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_numOrthoConst</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.maxBasisSize" title="primme_params.maxBasisSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_maxBasisSize</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.minRestartSize" title="primme_params.minRestartSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_minRestartSize</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.maxBlockSize" title="primme_params.maxBlockSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_maxBlockSize</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.maxMatvecs" title="primme_params.maxMatvecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_maxMatvecs</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.maxOuterIterations" title="primme_params.maxOuterIterations"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_maxOuterIterations</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.iseed" title="primme_params.iseed"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_iseed</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.aNorm" title="primme_params.aNorm"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_aNorm</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.BNorm" title="primme_params.BNorm"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_BNorm</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.invBNorm" title="primme_params.invBNorm"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_invBNorm</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.eps" title="primme_params.eps"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_eps</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.orth" title="primme_params.orth"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_orth</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.internalPrecision" title="primme_params.internalPrecision"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_internalPrecision</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.printLevel" title="primme_params.printLevel"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_printLevel</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.outputFile" title="primme_params.outputFile"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_outputFile</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.matrix" title="primme_params.matrix"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_matrix</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.massMatrix" title="primme_params.massMatrix"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_massMatrix</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.preconditioner" title="primme_params.preconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_preconditioner</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.initBasisMode" title="primme_params.initBasisMode"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_initBasisMode</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.projectionParams.projection" title="primme_params.projectionParams.projection"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_projectionParams_projection</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.restartingParams.maxPrevRetain" title="primme_params.restartingParams.maxPrevRetain"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_restartingParams_maxPrevRetain</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.precondition" title="primme_params.correctionParams.precondition"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_precondition</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.robustShifts" title="primme_params.correctionParams.robustShifts"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_robustShifts</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.maxInnerIterations" title="primme_params.correctionParams.maxInnerIterations"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_maxInnerIterations</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.projectors.LeftQ" title="primme_params.correctionParams.projectors.LeftQ"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_projectors_LeftQ</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.projectors.LeftX" title="primme_params.correctionParams.projectors.LeftX"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_projectors_LeftX</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.projectors.RightQ" title="primme_params.correctionParams.projectors.RightQ"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_projectors_RightQ</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.projectors.RightX" title="primme_params.correctionParams.projectors.RightX"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_projectors_RightX</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.projectors.SkewQ" title="primme_params.correctionParams.projectors.SkewQ"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_projectors_SkewQ</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.projectors.SkewX" title="primme_params.correctionParams.projectors.SkewX"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_projectors_SkewX</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.convTest" title="primme_params.correctionParams.convTest"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_convTest</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.correctionParams.relTolBase" title="primme_params.correctionParams.relTolBase"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_correctionParams_relTolBase</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.numOuterIterations" title="primme_params.stats.numOuterIterations"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_numOuterIterations</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.numRestarts" title="primme_params.stats.numRestarts"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_numRestarts</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.numMatvecs" title="primme_params.stats.numMatvecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_numMatvecs</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.numPreconds" title="primme_params.stats.numPreconds"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_numPreconds</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.numGlobalSum" title="primme_params.stats.numGlobalSum"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_numGlobalSum</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.volumeGlobalSum" title="primme_params.stats.volumeGlobalSum"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_volumeGlobalSum</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.numBroadcast" title="primme_params.stats.numBroadcast"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_numBroadcast</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.volumeBroadcast" title="primme_params.stats.volumeBroadcast"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_volumeBroadcast</span></code></a></div>
<div class="line"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_flopsDense</span></code></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.numOrthoInnerProds" title="primme_params.stats.numOrthoInnerProds"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_numOrthoInnerProds</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.elapsedTime" title="primme_params.stats.elapsedTime"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_elapsedTime</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.timeMatvec" title="primme_params.stats.timeMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_timeMatvec</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.timePrecond" title="primme_params.stats.timePrecond"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_timePrecond</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.timeOrtho" title="primme_params.stats.timeOrtho"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_timeOrtho</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.timeGlobalSum" title="primme_params.stats.timeGlobalSum"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_timeGlobalSum</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.timeBroadcast" title="primme_params.stats.timeBroadcast"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_timeBroadcast</span></code></a></div>
<div class="line"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_timeDense</span></code></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.estimateMinEVal" title="primme_params.stats.estimateMinEVal"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_estimateMinEVal</span></code></a></div>
<div class="line"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_estimateMaxEVal</span></code></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.estimateLargestSVal" title="primme_params.stats.estimateLargestSVal"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_estimateLargestSVal</span></code></a></div>
<div class="line"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_estimateBNorm</span></code></div>
<div class="line"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_estimateInvBNorm</span></code></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.maxConvTol" title="primme_params.stats.maxConvTol"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_maxConvTol</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.stats.lockingIssue" title="primme_params.stats.lockingIssue"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_stats_lockingIssue</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.dynamicMethodSwitch" title="primme_params.dynamicMethodSwitch"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_dynamicMethodSwitch</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.convTestFun" title="primme_params.convTestFun"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_convTestFun</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.convTestFun_type" title="primme_params.convTestFun_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_convTestFun_type</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.convtest" title="primme_params.convtest"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_convtest</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.ldevecs" title="primme_params.ldevecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_ldevecs</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.ldOPs" title="primme_params.ldOPs"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_ldOPs</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.monitorFun" title="primme_params.monitorFun"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_monitorFun</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.monitorFun_type" title="primme_params.monitorFun_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_monitorFun_type</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.monitor" title="primme_params.monitor"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_monitor</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params.queue" title="primme_params.queue"><code class="xref c c-member docutils literal notranslate"><span class="pre">PRIMME_queue</span></code></a></div>
</div>
</p></li>
<li><p><strong>value</strong> :: <p>(input) value to set.</p>
<p>If the type of the option is integer (<code class="docutils literal notranslate"><span class="pre">int</span></code>, <a class="reference internal" href="appendix.html#c.PRIMME_INT" title="PRIMME_INT"><code class="xref c c-type docutils literal notranslate"><span class="pre">PRIMME_INT</span></code></a>, <code class="docutils literal notranslate"><span class="pre">size_t</span></code>), the
type of <code class="docutils literal notranslate"><span class="pre">value</span></code> should be as long as <a class="reference internal" href="appendix.html#c.PRIMME_INT" title="PRIMME_INT"><code class="xref c c-type docutils literal notranslate"><span class="pre">PRIMME_INT</span></code></a>, which is <code class="docutils literal notranslate"><span class="pre">integer*8</span></code> by default.</p>
</p></li>
</ul>
</dd>
</dl>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p><strong>Don’t use</strong> this subroutine inside PRIMME’s callback functions, e.g., <a class="reference internal" href="appendix.html#c.primme_params.matrixMatvec" title="primme_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a> or
<a class="reference internal" href="appendix.html#c.primme_params.applyPreconditioner" title="primme_params.applyPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">applyPreconditioner</span></code></a>, or in functions called by these functions.</p>
</div>
</dd></dl>
</div>
<div class="section" id="primmetop-get-member-f77">
<h2>primmetop_get_member_f77<a class="headerlink" href="#primmetop-get-member-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/primmetop_get_member_f77">
<em class="property">subroutine </em><code class="sig-name descname">primmetop_get_member_f77</code><span class="sig-paren">(</span><em class="sig-param">primme</em>, <em class="sig-param">label</em>, <em class="sig-param">value</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/primmetop_get_member_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Get the value in some field of the parameter structure.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>label</strong><em> [</em><em>integer</em><em>]</em> :: (input) field where to get value. One of
the detailed in subroutine <code class="xref f f-func docutils literal notranslate"><span class="pre">primmetop_set_member_f77()</span></code>.</p></li>
<li><p><strong>value</strong> :: <p>(output) value of the field.</p>
<p>If the type of the option is integer (<code class="docutils literal notranslate"><span class="pre">int</span></code>, <a class="reference internal" href="appendix.html#c.PRIMME_INT" title="PRIMME_INT"><code class="xref c c-type docutils literal notranslate"><span class="pre">PRIMME_INT</span></code></a>, <code class="docutils literal notranslate"><span class="pre">size_t</span></code>), the
type of <code class="docutils literal notranslate"><span class="pre">value</span></code> should be as long as <a class="reference internal" href="appendix.html#c.PRIMME_INT" title="PRIMME_INT"><code class="xref c c-type docutils literal notranslate"><span class="pre">PRIMME_INT</span></code></a>, which is <code class="docutils literal notranslate"><span class="pre">integer*8</span></code> by default.</p>
</p></li>
</ul>
</dd>
</dl>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p><strong>Don’t use</strong> this subroutine inside PRIMME’s callback functions, e.g., <a class="reference internal" href="appendix.html#c.primme_params.matrixMatvec" title="primme_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a> or
<a class="reference internal" href="appendix.html#c.primme_params.applyPreconditioner" title="primme_params.applyPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">applyPreconditioner</span></code></a>, or in functions called by these functions. In those cases use
<a class="reference internal" href="#f/_/primme_get_member_f77" title="f/_/primme_get_member_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">primme_get_member_f77()</span></code></a>.</p>
</div>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>When <code class="docutils literal notranslate"><span class="pre">label</span></code> is one of <code class="docutils literal notranslate"><span class="pre">PRIMME_matrixMatvec</span></code>, <code class="docutils literal notranslate"><span class="pre">PRIMME_applyPreconditioner</span></code>, <code class="docutils literal notranslate"><span class="pre">PRIMME_commInfo</span></code>,
<code class="docutils literal notranslate"><span class="pre">PRIMME_matrix</span></code> and <code class="docutils literal notranslate"><span class="pre">PRIMME_preconditioner</span></code>,
the returned <code class="docutils literal notranslate"><span class="pre">value</span></code> is a C pointer (<code class="docutils literal notranslate"><span class="pre">void*</span></code>). Use Fortran pointer or other extensions to deal with it.
For instance:</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="k">use </span><span class="nb">iso_c_binding</span>
<span class="n">MPI_Comm</span> <span class="n">comm</span>
<span class="n">comm</span> <span class="o">=</span> <span class="n">MPI_COMM_WORLD</span>
<span class="k">call </span><span class="n">primme_set_member_f77</span><span class="p">(</span><span class="n">primme</span><span class="p">,</span> <span class="n">PRIMME_commInfo</span><span class="p">,</span> <span class="n">comm</span><span class="p">)</span>
<span class="p">...</span>
<span class="k">subroutine </span><span class="n">par_GlobalSumDouble</span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">y</span><span class="p">,</span><span class="n">k</span><span class="p">,</span><span class="n">primme</span><span class="p">)</span>
<span class="k">use </span><span class="nb">iso_c_binding</span>
<span class="k">implicit none</span>
<span class="p">...</span>
<span class="n">MPI_Comm</span><span class="p">,</span> <span class="k">pointer</span> <span class="kd">::</span> <span class="n">comm</span>
<span class="k">type</span><span class="p">(</span><span class="kt">c_ptr</span><span class="p">)</span> <span class="kd">::</span> <span class="n">pcomm</span>
<span class="k">call </span><span class="n">primme_get_member_f77</span><span class="p">(</span><span class="n">primme</span><span class="p">,</span> <span class="n">PRIMME_commInfo</span><span class="p">,</span> <span class="n">pcomm</span><span class="p">)</span>
<span class="k">call </span><span class="nb">c_f_pointer</span><span class="p">(</span><span class="n">pcomm</span><span class="p">,</span> <span class="n">comm</span><span class="p">)</span>
<span class="k">call </span><span class="n">MPI_Allreduce</span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">y</span><span class="p">,</span><span class="n">k</span><span class="p">,</span><span class="n">MPI_DOUBLE</span><span class="p">,</span><span class="n">MPI_SUM</span><span class="p">,</span><span class="n">comm</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
</pre></div>
</div>
<p>Most users would not need to retrieve these pointers in their programs.</p>
</div>
</dd></dl>
</div>
<div class="section" id="primmetop-get-prec-shift-f77">
<h2>primmetop_get_prec_shift_f77<a class="headerlink" href="#primmetop-get-prec-shift-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/primmetop_get_prec_shift_f77">
<em class="property">subroutine </em><code class="sig-name descname">primmetop_get_prec_shift_f77</code><span class="sig-paren">(</span><em class="sig-param">primme</em>, <em class="sig-param">index</em>, <em class="sig-param">value</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/primmetop_get_prec_shift_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Get the value in some position of the array <a class="reference internal" href="appendix.html#c.primme_params.ShiftsForPreconditioner" title="primme_params.ShiftsForPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">ShiftsForPreconditioner</span></code></a>.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>index</strong><em> [</em><em>integer</em><em>]</em> :: (input) position of the array; the first position is 1.</p></li>
<li><p><strong>value</strong> :: (output) value of the array at that position.</p></li>
</ul>
</dd>
</dl>
</dd></dl>
</div>
<div class="section" id="primme-get-member-f77">
<h2>primme_get_member_f77<a class="headerlink" href="#primme-get-member-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/primme_get_member_f77">
<em class="property">subroutine </em><code class="sig-name descname">primme_get_member_f77</code><span class="sig-paren">(</span><em class="sig-param">primme</em>, <em class="sig-param">label</em>, <em class="sig-param">value</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/primme_get_member_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Get the value in some field of the parameter structure.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>label</strong><em> [</em><em>integer</em><em>]</em> :: (input) field where to get value. One of
the detailed in subroutine <code class="xref f f-func docutils literal notranslate"><span class="pre">primmetop_set_member_f77()</span></code>.</p></li>
<li><p><strong>value</strong> :: <p>(output) value of the field.</p>
<p>If the type of the option is integer (<code class="docutils literal notranslate"><span class="pre">int</span></code>, <a class="reference internal" href="appendix.html#c.PRIMME_INT" title="PRIMME_INT"><code class="xref c c-type docutils literal notranslate"><span class="pre">PRIMME_INT</span></code></a>, <code class="docutils literal notranslate"><span class="pre">size_t</span></code>), the
type of <code class="docutils literal notranslate"><span class="pre">value</span></code> should be as long as <a class="reference internal" href="appendix.html#c.PRIMME_INT" title="PRIMME_INT"><code class="xref c c-type docutils literal notranslate"><span class="pre">PRIMME_INT</span></code></a>, which is <code class="docutils literal notranslate"><span class="pre">integer*8</span></code> by default.</p>
</p></li>
</ul>
</dd>
</dl>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>Use this subroutine exclusively inside PRIMME’s callback functions, e.g., <a class="reference internal" href="appendix.html#c.primme_params.matrixMatvec" title="primme_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a>
or <a class="reference internal" href="appendix.html#c.primme_params.applyPreconditioner" title="primme_params.applyPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">applyPreconditioner</span></code></a>, or in functions called by these functions. Otherwise, e.g.,
from the main program, use the subroutine <a class="reference internal" href="#f/_/primmetop_get_member_f77" title="f/_/primmetop_get_member_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">primmetop_get_member_f77()</span></code></a>.</p>
</div>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>When <code class="docutils literal notranslate"><span class="pre">label</span></code> is one of <code class="docutils literal notranslate"><span class="pre">PRIMME_matrixMatvec</span></code>, <code class="docutils literal notranslate"><span class="pre">PRIMME_applyPreconditioner</span></code>, <code class="docutils literal notranslate"><span class="pre">PRIMME_commInfo</span></code>,
<code class="docutils literal notranslate"><span class="pre">PRIMME_matrix</span></code> and <code class="docutils literal notranslate"><span class="pre">PRIMME_preconditioner</span></code>,
the returned <code class="docutils literal notranslate"><span class="pre">value</span></code> is a C pointer (<code class="docutils literal notranslate"><span class="pre">void*</span></code>). Use Fortran pointer or other extensions to deal with it.
For instance:</p>
<div class="highlight-fortran notranslate"><div class="highlight"><pre><span></span><span class="k">use </span><span class="nb">iso_c_binding</span>
<span class="n">MPI_Comm</span> <span class="n">comm</span>
<span class="n">comm</span> <span class="o">=</span> <span class="n">MPI_COMM_WORLD</span>
<span class="k">call </span><span class="n">primme_set_member_f77</span><span class="p">(</span><span class="n">primme</span><span class="p">,</span> <span class="n">PRIMME_commInfo</span><span class="p">,</span> <span class="n">comm</span><span class="p">)</span>
<span class="p">...</span>
<span class="k">subroutine </span><span class="n">par_GlobalSumDouble</span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">y</span><span class="p">,</span><span class="n">k</span><span class="p">,</span><span class="n">primme</span><span class="p">)</span>
<span class="k">use </span><span class="nb">iso_c_binding</span>
<span class="k">implicit none</span>
<span class="p">...</span>
<span class="n">MPI_Comm</span><span class="p">,</span> <span class="k">pointer</span> <span class="kd">::</span> <span class="n">comm</span>
<span class="k">type</span><span class="p">(</span><span class="kt">c_ptr</span><span class="p">)</span> <span class="kd">::</span> <span class="n">pcomm</span>
<span class="k">call </span><span class="n">primme_get_member_f77</span><span class="p">(</span><span class="n">primme</span><span class="p">,</span> <span class="n">PRIMME_commInfo</span><span class="p">,</span> <span class="n">pcomm</span><span class="p">)</span>
<span class="k">call </span><span class="nb">c_f_pointer</span><span class="p">(</span><span class="n">pcomm</span><span class="p">,</span> <span class="n">comm</span><span class="p">)</span>
<span class="k">call </span><span class="n">MPI_Allreduce</span><span class="p">(</span><span class="n">x</span><span class="p">,</span><span class="n">y</span><span class="p">,</span><span class="n">k</span><span class="p">,</span><span class="n">MPI_DOUBLE</span><span class="p">,</span><span class="n">MPI_SUM</span><span class="p">,</span><span class="n">comm</span><span class="p">,</span><span class="n">ierr</span><span class="p">)</span>
</pre></div>
</div>
<p>Most users would not need to retrieve these pointers in their programs.</p>
</div>
</dd></dl>
</div>
<div class="section" id="primme-get-prec-shift-f77">
<h2>primme_get_prec_shift_f77<a class="headerlink" href="#primme-get-prec-shift-f77" title="Permalink to this headline">¶</a></h2>
<dl class="subroutine">
<dt id="f/_/primme_get_prec_shift_f77">
<em class="property">subroutine </em><code class="sig-name descname">primme_get_prec_shift_f77</code><span class="sig-paren">(</span><em class="sig-param">primme</em>, <em class="sig-param">index</em>, <em class="sig-param">value</em><span class="sig-paren">)</span><a class="headerlink" href="#f/_/primme_get_prec_shift_f77" title="Permalink to this definition">¶</a></dt>
<dd><p>Get the value in some position of the array <a class="reference internal" href="appendix.html#c.primme_params.ShiftsForPreconditioner" title="primme_params.ShiftsForPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">ShiftsForPreconditioner</span></code></a>.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>primme</strong><em> [</em><a class="reference internal" href="#f/_/ptr" title="f/_/ptr"><em>ptr</em></a><em>]</em> :: (input) parameters structure.</p></li>
<li><p><strong>index</strong><em> [</em><em>integer</em><em>]</em> :: (input) position of the array; the first position is 1.</p></li>
<li><p><strong>value</strong> :: (output) value of the array at that position.</p></li>
</ul>
</dd>
</dl>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>Use this subroutine exclusively inside the subroutine <a class="reference internal" href="appendix.html#c.primme_params.matrixMatvec" title="primme_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a>, <a class="reference internal" href="appendix.html#c.primme_params.massMatrixMatvec" title="primme_params.massMatrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">massMatrixMatvec</span></code></a>, or <a class="reference internal" href="appendix.html#c.primme_params.applyPreconditioner" title="primme_params.applyPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">applyPreconditioner</span></code></a>.
Otherwise use the subroutine <a class="reference internal" href="#f/_/primmetop_get_prec_shift_f77" title="f/_/primmetop_get_prec_shift_f77"><code class="xref f f-func docutils literal notranslate"><span class="pre">primmetop_get_prec_shift_f77()</span></code></a>.</p>
</div>
</dd></dl>
</div>
</div>
</div>
</div>
<footer>
<div class="rst-footer-buttons" role="navigation" aria-label="footer navigation">
<a href="primmef90.html" class="btn btn-neutral float-right" title="FORTRAN 90 Library Interface" accesskey="n" rel="next">Next <span class="fa fa-arrow-circle-right"></span></a>
<a href="primmec.html" class="btn btn-neutral float-left" title="C Library Interface" accesskey="p" rel="prev"><span class="fa fa-arrow-circle-left"></span> Previous</a>
</div>
<hr/>
<div role="contentinfo">
<p>
© Copyright 2020, College of William & Mary
</p>
</div>
Built with <a href="http://sphinx-doc.org/">Sphinx</a> using a
<a href="https://github.com/rtfd/sphinx_rtd_theme">theme</a>
provided by <a href="https://readthedocs.org">Read the Docs</a>.
</footer>
</div>
</div>
</section>
</div>
<script type="text/javascript">
jQuery(function () {
SphinxRtdTheme.Navigation.enable(true);
});
</script>
</body>
</html>