-
Notifications
You must be signed in to change notification settings - Fork 40
/
svdsc.html
721 lines (568 loc) · 75.3 KB
/
svdsc.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
<!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>C 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 Library Interface" href="svdsf77.html" />
<link rel="prev" title="Singular Value Problems" href="apisvds.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"><a class="reference internal" href="apieigs.html">Eigenvalue Problems</a></li>
<li class="toctree-l1 current"><a class="reference internal" href="apisvds.html">Singular Value Problems</a><ul class="current">
<li class="toctree-l2 current"><a class="current reference internal" href="#">C Library Interface</a><ul>
<li class="toctree-l3"><a class="reference internal" href="#running">Running</a></li>
<li class="toctree-l3"><a class="reference internal" href="#parameters-guide">Parameters Guide</a></li>
<li class="toctree-l3"><a class="reference internal" href="#interface-description">Interface Description</a><ul>
<li class="toctree-l4"><a class="reference internal" href="#primme-svds">?primme_svds</a></li>
<li class="toctree-l4"><a class="reference internal" href="#magma-primme-svds"><cite>magma_?primme_svds</cite></a></li>
<li class="toctree-l4"><a class="reference internal" href="#primme-svds-initialize">primme_svds_initialize</a></li>
<li class="toctree-l4"><a class="reference internal" href="#primme-svds-create">primme_svds_create</a></li>
<li class="toctree-l4"><a class="reference internal" href="#primme-svds-set-method">primme_svds_set_method</a></li>
<li class="toctree-l4"><a class="reference internal" href="#primme-svds-display-params">primme_svds_display_params</a></li>
<li class="toctree-l4"><a class="reference internal" href="#primme-svds-free">primme_svds_free</a></li>
<li class="toctree-l4"><a class="reference internal" href="#primme-svds-params-destroy">primme_svds_params_destroy</a></li>
</ul>
</li>
</ul>
</li>
<li class="toctree-l2"><a class="reference internal" href="svdsf77.html">FORTRAN Library Interface</a></li>
<li class="toctree-l2"><a class="reference internal" href="svdsf90.html">FORTRAN 90 Library Interface</a></li>
<li class="toctree-l2"><a class="reference internal" href="pysvds.html">Python Interface</a></li>
<li class="toctree-l2"><a class="reference internal" href="matsvds.html">MATLAB Interface</a></li>
<li class="toctree-l2"><a class="reference internal" href="appendixsvds.html">Parameter Description</a></li>
<li class="toctree-l2"><a class="reference internal" href="appendixsvds.html#preset-methods">Preset Methods</a></li>
<li class="toctree-l2"><a class="reference internal" href="appendixsvds.html#error-codes">Error Codes</a></li>
</ul>
</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="apisvds.html">Singular Value Problems</a> »</li>
<li>C Library Interface</li>
<li class="wy-breadcrumbs-aside">
<a href="_sources/svdsc.rst.txt" rel="nofollow"> View page source</a>
</li>
</ul>
<div class="rst-breadcrumbs-buttons" role="navigation" aria-label="breadcrumb navigation">
<a href="svdsf77.html" class="btn btn-neutral float-right" title="FORTRAN Library Interface" accesskey="n">Next <span class="fa fa-arrow-circle-right"></span></a>
<a href="apisvds.html" class="btn btn-neutral float-left" title="Singular Value Problems" 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="c-library-interface">
<h1>C Library Interface<a class="headerlink" href="#c-library-interface" title="Permalink to this headline">¶</a></h1>
<div class="versionadded">
<p><span class="versionmodified added">New in version 2.0.</span></p>
</div>
<p>The PRIMME SVDS interface is composed of the following functions.
To solve real and complex singular value problems call respectively:</p>
<pre class="literal-block">int <a class="reference internal" href="#c.dprimme_svds" title="dprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">dprimme_svds</span></code></a> (double *svals, double *svecs, double *resNorms,
primme_svds_params *primme_svds)
int <a class="reference internal" href="#c.zprimme_svds" title="zprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">zprimme_svds</span></code></a> (double *svals, <a class="reference internal" href="appendix.html#c.PRIMME_COMPLEX_DOUBLE" title="PRIMME_COMPLEX_DOUBLE"><code class="xref c c-type docutils literal notranslate"><span class="pre">PRIMME_COMPLEX_DOUBLE</span></code></a> *svecs,
double *resNorms, primme_svds_params *primme_svds)</pre>
<p>There are more versions for matrix problems working in other precisions:</p>
<table class="docutils align-default">
<colgroup>
<col style="width: 18%" />
<col style="width: 41%" />
<col style="width: 41%" />
</colgroup>
<thead>
<tr class="row-odd"><th class="head"><p>Precision</p></th>
<th class="head"><p>Real</p></th>
<th class="head"><p>Complex</p></th>
</tr>
</thead>
<tbody>
<tr class="row-even"><td><p>half</p></td>
<td><p><a class="reference internal" href="#c.hprimme_svds" title="hprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">hprimme_svds()</span></code></a>
<a class="reference internal" href="#c.hsprimme_svds" title="hsprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">hsprimme_svds()</span></code></a></p></td>
<td><p><a class="reference internal" href="#c.kprimme_svds" title="kprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">kprimme_svds()</span></code></a>
<a class="reference internal" href="#c.ksprimme_svds" title="ksprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">ksprimme_svds()</span></code></a></p></td>
</tr>
<tr class="row-odd"><td><p>single</p></td>
<td><p><a class="reference internal" href="#c.sprimme_svds" title="sprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">sprimme_svds()</span></code></a></p></td>
<td><p><a class="reference internal" href="#c.cprimme_svds" title="cprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">cprimme_svds()</span></code></a></p></td>
</tr>
<tr class="row-even"><td><p>double</p></td>
<td><p><a class="reference internal" href="#c.dprimme_svds" title="dprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">dprimme_svds()</span></code></a></p></td>
<td><p><a class="reference internal" href="#c.zprimme_svds" title="zprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">zprimme_svds()</span></code></a></p></td>
</tr>
</tbody>
</table>
<p>Other useful functions:</p>
<pre class="literal-block">void <a class="reference internal" href="#c.primme_svds_initialize" title="primme_svds_initialize"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_initialize</span></code></a> (primme_svds_params *primme_svds)
int <a class="reference internal" href="#c.primme_svds_set_method" title="primme_svds_set_method"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_set_method</span></code></a> (primme_svds_preset_method method,
primme_preset_method methodStage1,
primme_preset_method methodStage2, primme_svds_params *primme_svds)
void <a class="reference internal" href="#c.primme_svds_display_params" title="primme_svds_display_params"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_display_params</span></code></a> (primme_svds_params primme_svds)
void <a class="reference internal" href="#c.primme_svds_free" title="primme_svds_free"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_free</span></code></a> (primme_svds_params *primme_svds)</pre>
<p>PRIMME SVDS stores its data on the structure <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params"><code class="xref c c-type docutils literal notranslate"><span class="pre">primme_svds_params</span></code></a>.
See <a class="reference internal" href="#svds-guide-params"><span class="std std-ref">Parameters Guide</span></a> for an introduction about its fields.</p>
<div class="section" id="running">
<h2>Running<a class="headerlink" href="#running" title="Permalink to this headline">¶</a></h2>
<p>To use PRIMME SVDS, follow these basic steps.</p>
<ol class="arabic">
<li><p>Include:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="cp">#include</span> <span class="cpf">"primme.h"</span><span class="c1"> /* header file is required to run primme */</span><span class="cp"></span>
</pre></div>
</div>
</li>
<li><p>Initialize a PRIMME SVDS parameters structure for default settings:</p>
<pre class="literal-block"><a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params"><code class="xref c c-type docutils literal notranslate"><span class="pre">primme_svds_params</span></code></a> primme_svds;
<a class="reference internal" href="#c.primme_svds_initialize" title="primme_svds_initialize"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_initialize</span></code></a> (&primme_svds);</pre>
</li>
<li><p>Set problem parameters (see also <a class="reference internal" href="#svds-guide-params"><span class="std std-ref">Parameters Guide</span></a>), and,
optionally, set one of the <a class="reference internal" href="appendixsvds.html#c.primme_svds_preset_method" title="primme_svds_preset_method"><code class="xref c c-type docutils literal notranslate"><span class="pre">preset</span> <span class="pre">methods</span></code></a>:</p>
<pre class="literal-block">primme_svds.<a class="reference internal" href="appendixsvds.html#c.primme_svds_params.matrixMatvec" title="primme_svds_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a> = matrixMatvec; /* MV product */
primme_svds.<a class="reference internal" href="appendixsvds.html#c.primme_svds_params.m" title="primme_svds_params.m"><code class="xref c c-member docutils literal notranslate"><span class="pre">m</span></code></a> = 1000; /* set the matrix dimensions */
primme_svds.<a class="reference internal" href="appendixsvds.html#c.primme_svds_params.n" title="primme_svds_params.n"><code class="xref c c-member docutils literal notranslate"><span class="pre">n</span></code></a> = 100;
primme_svds.<a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numSvals" title="primme_svds_params.numSvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numSvals</span></code></a> = 10; /* Number of singular values */
<a class="reference internal" href="#c.primme_svds_set_method" title="primme_svds_set_method"><code class="xref c c-func docutils literal notranslate"><span class="pre">primmesvds_set_method</span></code></a> (primme_svds_hybrid, PRIMME_DEFAULT_METHOD,
PRIMME_DEFAULT_METHOD, &primme_svds);
...</pre>
</li>
<li><p>Then to solve a real singular value problem call:</p>
<pre class="literal-block">ret = <a class="reference internal" href="#c.dprimme_svds" title="dprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">dprimme_svds</span></code></a> (svals, svecs, resNorms, &primme_svds);</pre>
<p>The previous is the double precision call. There is available calls for complex
double, single and complex single; check <a class="reference internal" href="#c.zprimme_svds" title="zprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">zprimme_svds()</span></code></a>, <a class="reference internal" href="#c.sprimme_svds" title="sprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">sprimme_svds()</span></code></a>
and <a class="reference internal" href="#c.cprimme_svds" title="cprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">cprimme_svds()</span></code></a>.</p>
<p>To solve complex singular value problems call:</p>
<pre class="literal-block">ret = <a class="reference internal" href="#c.zprimme_svds" title="zprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">zprimme_svds</span></code></a> (svals, svecs, resNorms, &primme_svds);</pre>
<p>The call arguments are:</p>
<ul class="simple">
<li><p><cite>svals</cite>, array to return the found singular values;</p></li>
<li><p><cite>svecs</cite>, array to return the found left and right singular vectors;</p></li>
<li><p><cite>resNorms</cite>, array to return the residual norms of the found triplets; and</p></li>
<li><p><cite>ret</cite>, returned error code.</p></li>
</ul>
</li>
<li><p>To free the work arrays in PRIMME SVDS:</p>
<pre class="literal-block"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_free</span></code> (&primme_svds);</pre>
</li>
</ol>
</div>
<div class="section" id="parameters-guide">
<span id="svds-guide-params"></span><h2>Parameters Guide<a class="headerlink" href="#parameters-guide" title="Permalink to this headline">¶</a></h2>
<p>PRIMME SVDS stores the data on the structure <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params"><code class="xref c c-type docutils literal notranslate"><span class="pre">primme_svds_params</span></code></a>, which has the next fields:</p>
<div class="line-block">
<div class="line"><em>Basic</em></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_INT</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.m" title="primme_svds_params.m"><code class="xref c c-member docutils literal notranslate"><span class="pre">m</span></code></a>, number of rows of the matrix.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_INT</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.n" title="primme_svds_params.n"><code class="xref c c-member docutils literal notranslate"><span class="pre">n</span></code></a>, number of columns of the matrix.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">(*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.matrixMatvec" title="primme_svds_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a> <code class="docutils literal notranslate"><span class="pre">)(...)</span></code>, matrix-vector product.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numSvals" title="primme_svds_params.numSvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numSvals</span></code></a>, how many singular triplets to find.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">primme_svds_target</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.target" title="primme_svds_params.target"><code class="xref c c-member docutils literal notranslate"><span class="pre">target</span></code></a>, which singular values to find.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">double</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.eps" title="primme_svds_params.eps"><code class="xref c c-member docutils literal notranslate"><span class="pre">eps</span></code></a>, tolerance of the residual norm of converged triplets.</div>
<div class="line"><br /></div>
<div class="line"><em>For parallel programs</em></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numProcs" title="primme_svds_params.numProcs"><code class="xref c c-member docutils literal notranslate"><span class="pre">numProcs</span></code></a>, number of processes</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.procID" title="primme_svds_params.procID"><code class="xref c c-member docutils literal notranslate"><span class="pre">procID</span></code></a>, rank of this process</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_INT</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.mLocal" title="primme_svds_params.mLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">mLocal</span></code></a>, number of rows stored in this process</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_INT</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.nLocal" title="primme_svds_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a>, number of columns stored in this process</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">(*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.globalSumReal" title="primme_svds_params.globalSumReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">globalSumReal</span></code></a> <code class="docutils literal notranslate"><span class="pre">)(...)</span></code>, sum reduction among processes</div>
<div class="line"><br /></div>
<div class="line"><em>Accelerate the convergence</em></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">(*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.applyPreconditioner" title="primme_svds_params.applyPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">applyPreconditioner</span></code></a> <code class="docutils literal notranslate"><span class="pre">)(...)</span></code>, preconditioner-vector product.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.initSize" title="primme_svds_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a>, initial vectors as approximate solutions.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.maxBasisSize" title="primme_svds_params.maxBasisSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">maxBasisSize</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <code class="xref c c-member docutils literal notranslate"><span class="pre">minRestartSize</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.maxBlockSize" title="primme_svds_params.maxBlockSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">maxBlockSize</span></code></a></div>
<div class="line"><br /></div>
<div class="line"><em>User data</em></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.commInfo" title="primme_svds_params.commInfo"><code class="xref c c-member docutils literal notranslate"><span class="pre">commInfo</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.matrix" title="primme_svds_params.matrix"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrix</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.preconditioner" title="primme_svds_params.preconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">preconditioner</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.convtest" title="primme_svds_params.convtest"><code class="xref c c-member docutils literal notranslate"><span class="pre">convtest</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.monitor" title="primme_svds_params.monitor"><code class="xref c c-member docutils literal notranslate"><span class="pre">monitor</span></code></a></div>
<div class="line"><br /></div>
<div class="line"><em>Advanced options</em></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numTargetShifts" title="primme_svds_params.numTargetShifts"><code class="xref c c-member docutils literal notranslate"><span class="pre">numTargetShifts</span></code></a>, for targeting interior singular values.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">double</span> <span class="pre">*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.targetShifts" title="primme_svds_params.targetShifts"><code class="xref c c-member docutils literal notranslate"><span class="pre">targetShifts</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code>, orthogonal constrains to the singular vectors.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.locking" title="primme_svds_params.locking"><code class="xref c c-member docutils literal notranslate"><span class="pre">locking</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_INT</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.maxMatvecs" title="primme_svds_params.maxMatvecs"><code class="xref c c-member docutils literal notranslate"><span class="pre">maxMatvecs</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">PRIMME_INT</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.iseed" title="primme_svds_params.iseed"><code class="xref c c-member docutils literal notranslate"><span class="pre">iseed</span></code></a> <code class="docutils literal notranslate"><span class="pre">[4]</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">double</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.aNorm" title="primme_svds_params.aNorm"><code class="xref c c-member docutils literal notranslate"><span class="pre">aNorm</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">int</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.printLevel" title="primme_svds_params.printLevel"><code class="xref c c-member docutils literal notranslate"><span class="pre">printLevel</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">FILE</span> <span class="pre">*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.outputFile" title="primme_svds_params.outputFile"><code class="xref c c-member docutils literal notranslate"><span class="pre">outputFile</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">primme_svds_operator</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.method" title="primme_svds_params.method"><code class="xref c c-member docutils literal notranslate"><span class="pre">method</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">primme_svds_operator</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.methodStage2" title="primme_svds_params.methodStage2"><code class="xref c c-member docutils literal notranslate"><span class="pre">methodStage2</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params" title="primme_params"><code class="xref c c-type docutils literal notranslate"><span class="pre">primme_params</span></code></a> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.primme" title="primme_svds_params.primme"><code class="xref c c-member docutils literal notranslate"><span class="pre">primme</span></code></a></div>
<div class="line"><a class="reference internal" href="appendix.html#c.primme_params" title="primme_params"><code class="xref c c-type docutils literal notranslate"><span class="pre">primme_params</span></code></a> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.primmeStage2" title="primme_svds_params.primmeStage2"><code class="xref c c-member docutils literal notranslate"><span class="pre">primmeStage2</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">(*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.convTestFun" title="primme_svds_params.convTestFun"><code class="xref c c-member docutils literal notranslate"><span class="pre">convTestFun</span></code></a> <code class="docutils literal notranslate"><span class="pre">)(...)</span></code>, custom convergence criterion.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">void</span> <span class="pre">(*</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.monitorFun" title="primme_svds_params.monitorFun"><code class="xref c c-member docutils literal notranslate"><span class="pre">monitorFun</span></code></a> <code class="docutils literal notranslate"><span class="pre">)(...)</span></code>, custom convergence history.</div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">primme_op_datatype</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.matrixMatvec_type" title="primme_svds_params.matrixMatvec_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec_type</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">primme_op_datatype</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.applyPreconditioner_type" title="primme_svds_params.applyPreconditioner_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">applyPreconditioner_type</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">primme_op_datatype</span></code> <code class="xref c c-member docutils literal notranslate"><span class="pre">globalSumReal_type</span></code></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">primme_op_datatype</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.broadcastReal_type" title="primme_svds_params.broadcastReal_type"><code class="xref c c-member docutils literal notranslate"><span class="pre">broadcastReal_type</span></code></a></div>
<div class="line"><code class="docutils literal notranslate"><span class="pre">primme_op_datatype</span></code> <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.internalPrecision" title="primme_svds_params.internalPrecision"><code class="xref c c-member docutils literal notranslate"><span class="pre">internalPrecision</span></code></a></div>
</div>
<p>PRIMME SVDS requires the user to set at least the matrix dimensions (<a class="reference internal" href="appendixsvds.html#c.primme_svds_params.m" title="primme_svds_params.m"><code class="xref c c-member docutils literal notranslate"><span class="pre">m</span></code></a> x <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.n" title="primme_svds_params.n"><code class="xref c c-member docutils literal notranslate"><span class="pre">n</span></code></a>) and
the matrix-vector product (<a class="reference internal" href="appendixsvds.html#c.primme_svds_params.matrixMatvec" title="primme_svds_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a>), as they define the problem to be solved.
For parallel programs, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.mLocal" title="primme_svds_params.mLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">mLocal</span></code></a>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.nLocal" title="primme_svds_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.procID" title="primme_svds_params.procID"><code class="xref c c-member docutils literal notranslate"><span class="pre">procID</span></code></a> and <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.globalSumReal" title="primme_svds_params.globalSumReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">globalSumReal</span></code></a> are also required.</p>
<p>In addition, most users would want to specify how many singular triplets to find,
and provide a preconditioner (if available).</p>
<p>It is useful to have set all these before calling <a class="reference internal" href="#c.primme_svds_set_method" title="primme_svds_set_method"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_set_method()</span></code></a>.
Also, if users have a preference on <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.maxBasisSize" title="primme_svds_params.maxBasisSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">maxBasisSize</span></code></a>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.maxBlockSize" title="primme_svds_params.maxBlockSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">maxBlockSize</span></code></a>, etc, they
should also provide them into <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params"><code class="xref c c-type docutils literal notranslate"><span class="pre">primme_svds_params</span></code></a> prior to the
<a class="reference internal" href="#c.primme_svds_set_method" title="primme_svds_set_method"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_set_method()</span></code></a> call. This helps <a class="reference internal" href="#c.primme_svds_set_method" title="primme_svds_set_method"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_set_method()</span></code></a> make
the right choice on other parameters. It is sometimes useful to check the actual
parameters that PRIMME SVDS is going to use (before calling it) or used (on return)
by printing them with <a class="reference internal" href="#c.primme_svds_display_params" title="primme_svds_display_params"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_display_params()</span></code></a>.</p>
</div>
<div class="section" id="interface-description">
<h2>Interface Description<a class="headerlink" href="#interface-description" title="Permalink to this headline">¶</a></h2>
<p>The next enumerations and functions are declared in <code class="docutils literal notranslate"><span class="pre">primme.h</span></code>.</p>
<div class="section" id="primme-svds">
<h3>?primme_svds<a class="headerlink" href="#primme-svds" title="Permalink to this headline">¶</a></h3>
<dl class="c function">
<dt id="c.hprimme_svds">
int <code class="sig-name descname">hprimme_svds</code><span class="sig-paren">(</span><a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>svecs</em>, <a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.hprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.hsprimme_svds">
int <code class="sig-name descname">hsprimme_svds</code><span class="sig-paren">(</span>float *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>svecs</em>, float *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.hsprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.kprimme_svds">
int <code class="sig-name descname">kprimme_svds</code><span class="sig-paren">(</span><a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_COMPLEX_HALF" title="PRIMME_COMPLEX_HALF">PRIMME_COMPLEX_HALF</a> *<em>svecs</em>, <a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.kprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.ksprimme_svds">
int <code class="sig-name descname">ksprimme_svds</code><span class="sig-paren">(</span>float *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_COMPLEX_HALF" title="PRIMME_COMPLEX_HALF">PRIMME_COMPLEX_HALF</a> *<em>svecs</em>, float *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.ksprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd><div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
<dl class="c function">
<dt id="c.sprimme_svds">
int <code class="sig-name descname">sprimme_svds</code><span class="sig-paren">(</span>float *<em>svals</em>, float *<em>svecs</em>, float *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.sprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.cprimme_svds">
int <code class="sig-name descname">cprimme_svds</code><span class="sig-paren">(</span>float *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_COMPLEX_FLOAT" title="PRIMME_COMPLEX_FLOAT">PRIMME_COMPLEX_FLOAT</a> *<em>svecs</em>, float *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.cprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd><div class="versionadded">
<p><span class="versionmodified added">New in version 2.0.</span></p>
</div>
</dd></dl>
<dl class="c function">
<dt id="c.dprimme_svds">
int <code class="sig-name descname">dprimme_svds</code><span class="sig-paren">(</span>double *<em>svals</em>, double *<em>svecs</em>, double *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.dprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.zprimme_svds">
int <code class="sig-name descname">zprimme_svds</code><span class="sig-paren">(</span>double *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_COMPLEX_DOUBLE" title="PRIMME_COMPLEX_DOUBLE">PRIMME_COMPLEX_DOUBLE</a> *<em>svecs</em>, double *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.zprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd><p>Solve a real singular value problem.</p>
<p>All arrays should be hosted on CPU. The computations are performed on CPU (see <a class="reference internal" href="#c.magma_dprimme_svds" title="magma_dprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">magma_dprimme_svds()</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>svals</strong> – array at least of size <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numSvals" title="primme_svds_params.numSvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numSvals</span></code></a> to store the
computed singular values; all processes in a parallel run return this local array with the same values.</p></li>
<li><p><strong>svecs</strong> – array at least of size (<a class="reference internal" href="appendixsvds.html#c.primme_svds_params.mLocal" title="primme_svds_params.mLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">mLocal</span></code></a> + <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.nLocal" title="primme_svds_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a>) times (<code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> + <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numSvals" title="primme_svds_params.numSvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numSvals</span></code></a>)
to store column-wise the (local part for this process of the) computed left singular vectors
and the right singular vectors.</p></li>
<li><p><strong>resNorms</strong> – array at least of size <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numSvals" title="primme_svds_params.numSvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numSvals</span></code></a> to store the
residual norms of the computed triplets; all processes in parallel run return this local array with
the same values.</p></li>
<li><p><strong>primme_svds</strong> – parameters structure.</p></li>
</ul>
</dd>
<dt class="field-even">Returns</dt>
<dd class="field-even"><p>error indicator; see <a class="reference internal" href="appendixsvds.html#error-codes-svds"><span class="std std-ref">Error Codes</span></a>.</p>
</dd>
</dl>
<p>On input, <code class="docutils literal notranslate"><span class="pre">svecs</span></code> should start with the content of the <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> left vectors,
followed by the <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.initSize" title="primme_svds_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a> left vectors, followed by the <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> right vectors, and
followed by the <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.initSize" title="primme_svds_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a> right vectors.</p>
<p>On return, the i-th left singular vector starts at svecs[( <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> +i)* <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.mLocal" title="primme_svds_params.mLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">mLocal</span></code></a> ].
The i-th right singular vector starts at svecs[( <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> + <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.initSize" title="primme_svds_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a> )* <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.mLocal" title="primme_svds_params.mLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">mLocal</span></code></a> + ( <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> +i)* <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.nLocal" title="primme_svds_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> ].
The first vector has i=0.</p>
<p>All internal operations are performed at the same precision than <code class="docutils literal notranslate"><span class="pre">svecs</span></code> unless the user sets <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.internalPrecision" title="primme_svds_params.internalPrecision"><code class="xref c c-member docutils literal notranslate"><span class="pre">internalPrecision</span></code></a> otherwise.
The functions <a class="reference internal" href="#c.hsprimme_svds" title="hsprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">hsprimme_svds()</span></code></a> and <a class="reference internal" href="#c.ksprimme_svds" title="ksprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">ksprimme_svds()</span></code></a> perform all computations in half precision by default and report the eigenvalues and the residual norms in single precision. These functions may help in applications that may be not built with a compiler supporting half precision.</p>
<p>The type and precision of the callbacks depends on the type and precision of <code class="docutils literal notranslate"><span class="pre">svecs</span></code>. Although this can be changed. See details for <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.matrixMatvec" title="primme_svds_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.applyPreconditioner" title="primme_svds_params.applyPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">applyPreconditioner</span></code></a>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.globalSumReal" title="primme_svds_params.globalSumReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">globalSumReal</span></code></a>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.broadcastReal" title="primme_svds_params.broadcastReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">broadcastReal</span></code></a>, and <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.convTestFun" title="primme_svds_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="magma-primme-svds">
<h3><cite>magma_?primme_svds</cite><a class="headerlink" href="#magma-primme-svds" title="Permalink to this headline">¶</a></h3>
<dl class="c function">
<dt id="c.magma_hprimme_svds">
int <code class="sig-name descname">magma_hprimme_svds</code><span class="sig-paren">(</span><a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>svecs</em>, <a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.magma_hprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.magma_hsprimme_svds">
int <code class="sig-name descname">magma_hsprimme_svds</code><span class="sig-paren">(</span>float *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>svecs</em>, float *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.magma_hsprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.magma_kprimme_svds">
int <code class="sig-name descname">magma_kprimme_svds</code><span class="sig-paren">(</span><a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_COMPLEX_HALF" title="PRIMME_COMPLEX_HALF">PRIMME_COMPLEX_HALF</a> *<em>svecs</em>, <a class="reference internal" href="appendix.html#c.PRIMME_HALF" title="PRIMME_HALF">PRIMME_HALF</a> *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.magma_kprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.magma_ksprimme_svds">
int <code class="sig-name descname">magma_ksprimme_svds</code><span class="sig-paren">(</span>float *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_COMPLEX_HALF" title="PRIMME_COMPLEX_HALF">PRIMME_COMPLEX_HALF</a> *<em>svecs</em>, float *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.magma_ksprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.magma_sprimme_svds">
int <code class="sig-name descname">magma_sprimme_svds</code><span class="sig-paren">(</span>float *<em>svals</em>, float *<em>svecs</em>, float *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.magma_sprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.magma_cprimme_svds">
int <code class="sig-name descname">magma_cprimme_svds</code><span class="sig-paren">(</span>float *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_COMPLEX_FLOAT" title="PRIMME_COMPLEX_FLOAT">PRIMME_COMPLEX_FLOAT</a> *<em>svecs</em>, float *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.magma_cprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.magma_dprimme_svds">
int <code class="sig-name descname">magma_dprimme_svds</code><span class="sig-paren">(</span>double *<em>svals</em>, double *<em>svecs</em>, double *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.magma_dprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd></dd></dl>
<dl class="c function">
<dt id="c.magma_zprimme_svds">
int <code class="sig-name descname">magma_zprimme_svds</code><span class="sig-paren">(</span>double *<em>svals</em>, <a class="reference internal" href="appendix.html#c.PRIMME_COMPLEX_DOUBLE" title="PRIMME_COMPLEX_DOUBLE">PRIMME_COMPLEX_DOUBLE</a> *<em>svecs</em>, double *<em>resNorms</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.magma_zprimme_svds" title="Permalink to this definition">¶</a><br /></dt>
<dd><p>Solve a real singular value problem.</p>
<p>Most of the computations are performed on GPU (see <a class="reference internal" href="#c.dprimme_svds" title="dprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">dprimme_svds()</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>svals</strong> – CPU array at least of size <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numSvals" title="primme_svds_params.numSvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numSvals</span></code></a> to store the
computed singular values; all processes in a parallel run return this local array with the same values.</p></li>
<li><p><strong>svecs</strong> – GPU array at least of size (<a class="reference internal" href="appendixsvds.html#c.primme_svds_params.mLocal" title="primme_svds_params.mLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">mLocal</span></code></a> + <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.nLocal" title="primme_svds_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a>) times (<code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> + <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numSvals" title="primme_svds_params.numSvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numSvals</span></code></a>)
to store column-wise the (local part for this process of the) computed left singular vectors
and the right singular vectors.</p></li>
<li><p><strong>resNorms</strong> – CPU array at least of size <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.numSvals" title="primme_svds_params.numSvals"><code class="xref c c-member docutils literal notranslate"><span class="pre">numSvals</span></code></a> to store the
residual norms of the computed triplets; all processes in parallel run return this local array with
the same values.</p></li>
<li><p><strong>primme_svds</strong> – parameters structure.</p></li>
</ul>
</dd>
<dt class="field-even">Returns</dt>
<dd class="field-even"><p>error indicator; see <a class="reference internal" href="appendixsvds.html#error-codes-svds"><span class="std std-ref">Error Codes</span></a>.</p>
</dd>
</dl>
<p>On input, <code class="docutils literal notranslate"><span class="pre">svecs</span></code> should start with the content of the <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> left vectors,
followed by the <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.initSize" title="primme_svds_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a> left vectors, followed by the <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> right vectors, and
followed by the <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.initSize" title="primme_svds_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a> right vectors.</p>
<p>On return, the i-th left singular vector starts at svecs[( <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> +i)* <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.mLocal" title="primme_svds_params.mLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">mLocal</span></code></a> ].
The i-th right singular vector starts at svecs[( <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> + <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.initSize" title="primme_svds_params.initSize"><code class="xref c c-member docutils literal notranslate"><span class="pre">initSize</span></code></a> )* <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.mLocal" title="primme_svds_params.mLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">mLocal</span></code></a> + ( <code class="xref c c-member docutils literal notranslate"><span class="pre">numOrthoConst</span></code> +i)* <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.nLocal" title="primme_svds_params.nLocal"><code class="xref c c-member docutils literal notranslate"><span class="pre">nLocal</span></code></a> ].
The first vector has i=0.</p>
<p>All internal operations are performed at the same precision than <code class="docutils literal notranslate"><span class="pre">svecs</span></code> unless the user sets <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.internalPrecision" title="primme_svds_params.internalPrecision"><code class="xref c c-member docutils literal notranslate"><span class="pre">internalPrecision</span></code></a> otherwise.
The functions <a class="reference internal" href="#c.magma_hsprimme_svds" title="magma_hsprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">magma_hsprimme_svds()</span></code></a> and <a class="reference internal" href="#c.magma_ksprimme_svds" title="magma_ksprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">magma_ksprimme_svds()</span></code></a> perform all computations in half precision by default and report the eigenvalues and the residual norms in single precision. These functions may help in applications that may be not built with a compiler supporting half precision.</p>
<p>The type and precision of the callbacks depends on the type and precision of <code class="docutils literal notranslate"><span class="pre">svecs</span></code>. Although this can be changed. See details for <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.matrixMatvec" title="primme_svds_params.matrixMatvec"><code class="xref c c-member docutils literal notranslate"><span class="pre">matrixMatvec</span></code></a>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.applyPreconditioner" title="primme_svds_params.applyPreconditioner"><code class="xref c c-member docutils literal notranslate"><span class="pre">applyPreconditioner</span></code></a>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.globalSumReal" title="primme_svds_params.globalSumReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">globalSumReal</span></code></a>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.broadcastReal" title="primme_svds_params.broadcastReal"><code class="xref c c-member docutils literal notranslate"><span class="pre">broadcastReal</span></code></a>, and <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.convTestFun" title="primme_svds_params.convTestFun"><code class="xref c c-member docutils literal notranslate"><span class="pre">convTestFun</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-svds-initialize">
<h3>primme_svds_initialize<a class="headerlink" href="#primme-svds-initialize" title="Permalink to this headline">¶</a></h3>
<dl class="c function">
<dt id="c.primme_svds_initialize">
void <code class="sig-name descname">primme_svds_initialize</code><span class="sig-paren">(</span><a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.primme_svds_initialize" title="Permalink to this definition">¶</a><br /></dt>
<dd><p>Initialize PRIMME SVDS parameters structure to the default values.</p>
<p>After calling <a class="reference internal" href="#c.dprimme_svds" title="dprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">dprimme_svds()</span></code></a> (or a variant), call <a class="reference internal" href="#c.primme_svds_free" title="primme_svds_free"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_free()</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"><ul class="simple">
<li><p><strong>primme_svds</strong> – parameters structure.</p></li>
</ul>
</dd>
</dl>
<p>Example:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n">primme_svds_params</span> <span class="n">primme_svds</span><span class="p">;</span>
<span class="n">primme_svds_initialize</span><span class="p">(</span><span class="o">&</span><span class="n">primme_svds</span><span class="p">);</span>
<span class="n">primme_svds</span><span class="p">.</span><span class="n">n</span> <span class="o">=</span> <span class="mi">100</span><span class="p">;</span>
<span class="p">...</span>
<span class="n">dprimme_svds</span><span class="p">(</span><span class="n">svals</span><span class="p">,</span> <span class="n">svecs</span><span class="p">,</span> <span class="n">rnorms</span><span class="p">,</span> <span class="o">&</span><span class="n">primme_svds</span><span class="p">);</span>
<span class="p">...</span>
<span class="n">primme_svds_free</span><span class="p">(</span><span class="o">&</span><span class="n">primme_svds</span><span class="p">);</span>
</pre></div>
</div>
<p>See the alternative function <a class="reference internal" href="#c.primme_svds_params_create" title="primme_svds_params_create"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_params_create()</span></code></a> that also allocates the structure.</p>
</dd></dl>
</div>
<div class="section" id="primme-svds-create">
<h3>primme_svds_create<a class="headerlink" href="#primme-svds-create" title="Permalink to this headline">¶</a></h3>
<dl class="c function">
<dt id="c.primme_svds_params_create">
<a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<code class="sig-name descname">primme_svds_params_create</code><span class="sig-paren">(</span>void<span class="sig-paren">)</span><a class="headerlink" href="#c.primme_svds_params_create" title="Permalink to this definition">¶</a><br /></dt>
<dd><p>Allocate and initialize a parameters structure to the default values.</p>
<p>After calling <a class="reference internal" href="#c.dprimme_svds" title="dprimme_svds"><code class="xref c c-func docutils literal notranslate"><span class="pre">dprimme_svds()</span></code></a> (or a variant), call <a class="reference internal" href="#c.primme_svds_params_destroy" title="primme_svds_params_destroy"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_params_destroy()</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"><ul class="simple">
<li><p><strong>primme_sv</strong> – parameters structure.</p></li>
</ul>
</dd>
</dl>
<p>Example:</p>
<div class="highlight-c notranslate"><div class="highlight"><pre><span></span><span class="n">primme_svds_params</span> <span class="o">*</span><span class="n">primme_svds</span> <span class="o">=</span> <span class="n">primme_svds_params_create</span><span class="p">();</span>
<span class="n">primme_svds</span><span class="o">-></span><span class="n">n</span> <span class="o">=</span> <span class="mi">100</span><span class="p">;</span>
<span class="p">...</span>
<span class="n">dprimme_svds</span><span class="p">(</span><span class="n">svals</span><span class="p">,</span> <span class="n">svecs</span><span class="p">,</span> <span class="n">rnorms</span><span class="p">,</span> <span class="n">primme_svds</span><span class="p">);</span>
<span class="p">...</span>
<span class="n">primme_svds_params_destroy</span><span class="p">(</span><span class="n">primme_svds</span><span class="p">);</span>
</pre></div>
</div>
<p>See the alternative function <a class="reference internal" href="#c.primme_svds_initialize" title="primme_svds_initialize"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_initialize()</span></code></a> that only initializes the structure.</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-svds-set-method">
<h3>primme_svds_set_method<a class="headerlink" href="#primme-svds-set-method" title="Permalink to this headline">¶</a></h3>
<dl class="c function">
<dt id="c.primme_svds_set_method">
int <code class="sig-name descname">primme_svds_set_method</code><span class="sig-paren">(</span><a class="reference internal" href="appendixsvds.html#c.primme_svds_preset_method" title="primme_svds_preset_method">primme_svds_preset_method</a> <em>method</em>, <a class="reference internal" href="appendix.html#c.primme_preset_method" title="primme_preset_method">primme_preset_method</a> <em>methodStage1</em>, <a class="reference internal" href="appendix.html#c.primme_preset_method" title="primme_preset_method">primme_preset_method</a> <em>methodStage2</em>, <a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.primme_svds_set_method" title="Permalink to this definition">¶</a><br /></dt>
<dd><p>Set PRIMME SVDS 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> – <p>preset method to compute the singular triplets; one of</p>
<ul>
<li><p><a class="reference internal" href="appendixsvds.html#c.primme_svds_preset_method.primme_svds_default" title="primme_svds_preset_method.primme_svds_default"><code class="xref c c-member docutils literal notranslate"><span class="pre">primme_svds_default</span></code></a>, currently set as <a class="reference internal" href="appendixsvds.html#c.primme_svds_preset_method.primme_svds_hybrid" title="primme_svds_preset_method.primme_svds_hybrid"><code class="xref c c-member docutils literal notranslate"><span class="pre">primme_svds_hybrid</span></code></a>.</p></li>
<li><p><a class="reference internal" href="appendixsvds.html#c.primme_svds_preset_method.primme_svds_normalequations" title="primme_svds_preset_method.primme_svds_normalequations"><code class="xref c c-member docutils literal notranslate"><span class="pre">primme_svds_normalequations</span></code></a>, compute the eigenvectors of <span class="math notranslate nohighlight">\(A^*A\)</span> or <span class="math notranslate nohighlight">\(A A^*\)</span>.</p></li>
<li><p><a class="reference internal" href="appendixsvds.html#c.primme_svds_preset_method.primme_svds_augmented" title="primme_svds_preset_method.primme_svds_augmented"><code class="xref c c-member docutils literal notranslate"><span class="pre">primme_svds_augmented</span></code></a>, compute the eigenvectors of the augmented matrix, <span class="math notranslate nohighlight">\(\left(\begin{array}{cc} 0 & A^* \\ A & 0 \end{array}\right)\)</span>.</p></li>
<li><p><a class="reference internal" href="appendixsvds.html#c.primme_svds_preset_method.primme_svds_hybrid" title="primme_svds_preset_method.primme_svds_hybrid"><code class="xref c c-member docutils literal notranslate"><span class="pre">primme_svds_hybrid</span></code></a>, start with <a class="reference internal" href="appendixsvds.html#c.primme_svds_preset_method.primme_svds_normalequations" title="primme_svds_preset_method.primme_svds_normalequations"><code class="xref c c-member docutils literal notranslate"><span class="pre">primme_svds_normalequations</span></code></a>; use the
resulting approximate singular vectors as initial vectors for
<a class="reference internal" href="appendixsvds.html#c.primme_svds_preset_method.primme_svds_augmented" title="primme_svds_preset_method.primme_svds_augmented"><code class="xref c c-member docutils literal notranslate"><span class="pre">primme_svds_augmented</span></code></a> if the required accuracy was not achieved.</p></li>
</ul>
</p></li>
<li><p><strong>methodStage1</strong> – preset method to compute the eigenpairs at the first stage; see available values at <a class="reference internal" href="primmec.html#c.primme_set_method" title="primme_set_method"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_set_method()</span></code></a>.</p></li>
<li><p><strong>methodStage2</strong> – preset method to compute the eigenpairs with
the second stage of <code class="docutils literal notranslate"><span class="pre">primme_svds_hybrid</span></code>; see available values at <a class="reference internal" href="primmec.html#c.primme_set_method" title="primme_set_method"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_set_method()</span></code></a>.</p></li>
<li><p><strong>primme_svds</strong> – parameters structure.</p></li>
</ul>
</dd>
</dl>
<p>See also <a class="reference internal" href="appendixsvds.html#methods-svds"><span class="std std-ref">Preset Methods</span></a>.</p>
</dd></dl>
</div>
<div class="section" id="primme-svds-display-params">
<h3>primme_svds_display_params<a class="headerlink" href="#primme-svds-display-params" title="Permalink to this headline">¶</a></h3>
<dl class="c function">
<dt id="c.primme_svds_display_params">
void <code class="sig-name descname">primme_svds_display_params</code><span class="sig-paren">(</span><a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> <em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.primme_svds_display_params" title="Permalink to this definition">¶</a><br /></dt>
<dd><p>Display all printable settings of <code class="docutils literal notranslate"><span class="pre">primme_svds</span></code> into the file descriptor <a class="reference internal" href="appendixsvds.html#c.primme_svds_params.outputFile" title="primme_svds_params.outputFile"><code class="xref c c-member docutils literal notranslate"><span class="pre">outputFile</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_svds</strong> – parameters structure.</p></li>
</ul>
</dd>
</dl>
</dd></dl>
</div>
<div class="section" id="primme-svds-free">
<h3>primme_svds_free<a class="headerlink" href="#primme-svds-free" title="Permalink to this headline">¶</a></h3>
<dl class="c function">
<dt id="c.primme_svds_free">
void <code class="sig-name descname">primme_svds_free</code><span class="sig-paren">(</span><a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme_svds</em><span class="sig-paren">)</span><a class="headerlink" href="#c.primme_svds_free" title="Permalink to this definition">¶</a><br /></dt>
<dd><p>Free memory allocated by PRIMME SVDS.</p>
<dl class="field-list simple">
<dt class="field-odd">Parameters</dt>
<dd class="field-odd"><ul class="simple">
<li><p><strong>primme_svds</strong> – parameters structure.</p></li>
</ul>
</dd>
</dl>
</dd></dl>
</div>
<div class="section" id="primme-svds-params-destroy">
<h3>primme_svds_params_destroy<a class="headerlink" href="#primme-svds-params-destroy" title="Permalink to this headline">¶</a></h3>
<dl class="c function">
<dt id="c.primme_svds_params_destroy">
int <code class="sig-name descname">primme_svds_params_destroy</code><span class="sig-paren">(</span><a class="reference internal" href="appendixsvds.html#c.primme_svds_params" title="primme_svds_params">primme_svds_params</a> *<em>primme</em><span class="sig-paren">)</span><a class="headerlink" href="#c.primme_svds_params_destroy" title="Permalink to this definition">¶</a><br /></dt>
<dd><p>Free memory allocated by PRIMME associated to a parameters structure created
with <a class="reference internal" href="#c.primme_svds_params_create" title="primme_svds_params_create"><code class="xref c c-func docutils literal notranslate"><span class="pre">primme_svds_params_create()</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_svds</strong> – parameters structure.</p></li>
</ul>
</dd>
<dt class="field-even">Returns</dt>
<dd class="field-even"><p>nonzero value if the call is not successful.</p>
</dd>
</dl>
<div class="versionadded">
<p><span class="versionmodified added">New in version 3.0.</span></p>
</div>
</dd></dl>
</div>
</div>
</div>
</div>
</div>
<footer>
<div class="rst-footer-buttons" role="navigation" aria-label="footer navigation">
<a href="svdsf77.html" class="btn btn-neutral float-right" title="FORTRAN Library Interface" accesskey="n" rel="next">Next <span class="fa fa-arrow-circle-right"></span></a>
<a href="apisvds.html" class="btn btn-neutral float-left" title="Singular Value Problems" 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>