forked from primme/primme
-
Notifications
You must be signed in to change notification settings - Fork 0
/
readme.txt
2208 lines (1429 loc) · 63.3 KB
/
readme.txt
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
Note: For hyperlinked html and pdf versions of this document see
directory "doc".
Welcome to PRIMME's documentation!
**********************************
Table Of Contents:
* PRIMME: PReconditioned Iterative MultiMethod Eigensolver
* Changelog
* Citing this code
* License Information
* Contact Information
* Directory Structure
* Making and Linking
* Considerations using an IDE
* Tested Systems
* C Library Interface
* Running
* Parameters Guide
* Interface Description
* dprimme
* zprimme
* primme_initialize
* primme_set_method
* primme_display_params
* primme_Free
* FORTRAN Library Interface
* primme_initialize_f77
* primme_set_method_f77
* primme_Free_f77
* dprimme_f77
* zprimme_f77
* primmetop_set_member_f77
* primmetop_get_member_f77
* primmetop_get_prec_shift_f77
* primme_set_member_f77
* primme_get_member_f77
* primme_get_prec_shift_f77
* Appendix
* primme_params
* Error Codes
* Preset Methods
PRIMME: PReconditioned Iterative MultiMethod Eigensolver
********************************************************
PRIMME, pronounced as *prime*, finds a number of eigenvalues and their
corresponding eigenvectors of a real symmetric, or complex hermitian
matrix A. Largest, smallest and interior eigenvalues are supported.
Preconditioning can be used to accelerate convergence. PRIMME is
written in C99, but complete interfaces are provided for Fortran 77
and MATLAB.
Changelog
=========
Changes in PRIMME 1.2.2 (released on October 13, 2015):
* Fixed wrong symbols in "libdprimme.a" and "libzprimme.a".
* "primme_set_method()" sets "JDQMR" instead of "JDQMR_ETol" for
preset methods "DEFAULT_MIN_TIME" and "DYNAMIC" when seeking
interior values.
* Fixed compilation of driver with a PETSc installation without
HYPRE.
* Included the content of the environment variable "INCLUDE" for
compiling the driver.
Changes in PRIMME 1.2.1 (released on September 7, 2015):
* Added MATLAB interface to full PRIMME functionality.
* Support for BLAS/LAPACK with 64bits integers
("-DPRIMME_BLASINT_SIZE=64").
* Simplified configuration of Make_flags and Make_links (removed
"TOP" variable and replaced defines "NUM_SUM" and "NUM_IBM" by
"F77UNDERSCORE").
* Replaced directories "DTEST" and "ZTEST" by "TEST", that has:
* "driver.c": read matrices in MatrixMarket format and PETSc
binary and call PRIMME with the parameters specified in a file;
support complex arithmetic and MPI and can use PETSc
preconditioners.
* "ex*.c" and "ex*.f": small, didactic examples of usage in C and
Fortran and in parallel (with PETSc).
* Fixed a few minor bugs and improved documentation (especially the
F77 interface).
* Using Sphinx to manage documentation.
Changes in PRIMME 1.2 (released on December 21, 2014):
* A Fortran compiler is no longer required for building the PRIMME
library. Fortran programs can still be linked to PRIMME's F77
interface.
* Fixed some uncommon issues with the F77 interface.
* PRIMME can be called now multiple times from the same program.
* Performance improvements in the QMR inner solver, especially for
complex arithmetic.
* Fixed a couple of bugs with the locking functionality.
* In certain extreme cases where all eigenvalues of a matrix were
needed.
* The order of selecting interior eigenvalues.
The above fixes have improved robustness and performance.
* PRIMME now assigns unique random seeds per parallel process for up
to 4096^3 (140 trillion) processes.
* For the "DYNAMIC" method, fixed issues with initialization and
synchronization decisions across multiple processes.
* Fixed uncommon library interface bugs, coordinated better setting
the method and the user setting of parameters, and improved the
interface in the sample programs and makefiles.
* Other performance and documentation improvements.
Citing this code
================
Please cite:
[r1] A. Stathopoulos and J. R. McCombs PRIMME: *PReconditioned
Iterative MultiMethod Eigensolver: Methods and software
description*, ACM Transaction on Mathematical Software Vol. 37,
No. 2, (2010), 21:1-21:30.
More information on the algorithms and research that led to this
software can be found in the rest of the papers. The work has been
supported by a number of grants from the National Science Foundation.
[r2] A. Stathopoulos, *Nearly optimal preconditioned methods for
hermitian eigenproblems under limited memory. Part I: Seeking one
eigenvalue*, SIAM J. Sci. Comput., Vol. 29, No. 2, (2007), 481--
514.
[r3] A. Stathopoulos and J. R. McCombs, *Nearly optimal
preconditioned methods for hermitian eigenproblems under limited
memory. Part II: Seeking many eigenvalues*, SIAM J. Sci. Comput.,
Vol. 29, No. 5, (2007), 2162-2188.
[r4] J. R. McCombs and A. Stathopoulos, *Iterative Validation of
Eigensolvers: A Scheme for Improving the Reliability of Hermitian
Eigenvalue Solvers*, SIAM J. Sci. Comput., Vol. 28, No. 6,
(2006), 2337-2358.
[r5] A. Stathopoulos, *Locking issues for finding a large number
of eigenvectors of hermitian matrices*, Tech Report: WM-
CS-2005-03, July, 2005.
License Information
===================
PRIMME is free software; you can redistribute it and/or modify it
under the terms of the GNU Lesser General Public License as published
by the Free Software Foundation; either version 2.1 of the License, or
(at your option) any later version.
PRIMME is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
License for more details.
You should have received a copy of the GNU Lesser General Public
License along with this library; if not, write to the Free Software
Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
USA
Contact Information
===================
For reporting bugs or questions about functionality contact Andreas
Stathopoulos by email, *andreas* at *cs.wm.edu*. See further
information in the webpage http://www.cs.wm.edu/~andreas/software .
Directory Structure
===================
The next directories and files should be available:
* "COPYING.txt", LGPL License;
* "Make_flags", flags to be used by makefiles to compile library
and tests;
* "Link_flags", flags needed in making and linking the test
programs;
* "PRIMMESRC/", directory with source code in the following
subdirectories:
* "COMMONSRC/", interface and common functions used by all
precision versions;
* "DSRC/", the source code for the double precision
"dprimme()";
* "ZSRC/", the source code for the double complex
precision "zprimme()";
* "MEX/", MATLAB interface for PRIMME;
* "TEST/", sample test programs in C and F77, both
sequential and parallel;
* "libprimme.a", the PRIMME library (to be made);
* "makefile" main make file;
* "readme.txt" text version of the documentation;
* "doc/" directory with the HTML and PDF versions of the
documentation.
Making and Linking
==================
"Make_flags" has the flags and compilers used to make "libprimme.a":
* *CC*, compiler program such as "gcc", "clang" or "icc".
* *CFLAGS*, compiler options such as "-g" or "-O3". Also include
some of the following options if required for the BLAS and LAPACK
libraries to be linked:
* "-DF77UNDERSCORE", if Fortran appends an underscore to function
names (usually they does).
* "-DPRIMME_BLASINT_SIZE=64", if the library integers are 64-bit
integer ("kind=8") type (usually they are not).
Note: When "-DPRIMME_BLASINT_SIZE=64" is set the code uses the type
"int64_t" supported by the C99 standard. In case the compiler
doesn't honor the standard, replace the next lines in
"PRIMMESRC/COMMONSRC/common_numerical.h":
#if !defined(PRIMME_BLASINT_SIZE)
# define PRIMME_BLASINT int
#else
# include <stdint.h>
# define GENERIC_INT(N) int ## N ## _t
# define XGENERIC_INT(N) GENERIC_INT(N)
# define PRIMME_BLASINT XGENERIC_INT(PRIMME_BLASINT_SIZE)
#endif
by the next macro definition with the proper type for an "int" of 64
bits:
#define PRIMME_BLASINT __int64
After customizing "Make_flags", type this to generate "libprimme.a":
make lib
Making can be also done at the command line:
make lib CC=clang CFLAGS='-O3'
"Link_flags" has the flags for linking with external libraries and
making the executables located in "TEST":
* *LDFLAGS*, linker flags such as "-framework Accelerate".
* *LIBS*, flags to link with libraries (BLAS and LAPACK are
required), such as "-lprimme -llapack -lblas -lgfortran -lm".
After that, type this to compile and execute a simple test:
$ make test
...
Test passed!
...
Test passed!
If it worked, try with other examples in "TEST" (see "README" in
"TEST" for more information about how to compile the driver and the
examples).
In case of linking problems check flags in *LDFLAGS* and *LIBS* and
consider to add/remove "-DF77UNDERSCORE" from *CFLAGS*. If the
execution fails consider to add/remove "-DPRIMME_BLASINT_SIZE=64" from
*CFLAGS*.
Full description of actions that *make* can take:
* *make lib*, builds "libprimme.a"; alternatively:
* *make libd*, if only "dprimme()" is of interest, build
"libdprimme.a":
* *make libz*, if only "zprimme()" is of interest, build
"libzprimme.a";
* *make test*, build and execute a simple example;
* *make clean*, removes all "*.o", "a.out", and core files from all
directories.
Considerations using an IDE
---------------------------
PRIMME can be built in other environments such as Anjuta, Eclipse,
KDevelop, Qt Creator, Visual Studio and XCode. To build the PRIMME
library do the following:
1. Create a new project and include the source files under the
directory "PRIMMESRC".
2. Add the directory "PRIMMESRC/COMMONSRC" as an include directory.
To build an example code using PRIMME make sure:
* to add a reference for PRIMME, BLAS and LAPACK libraries;
* to add the directory "PRIMMESRC/COMMONSRC" as an include
directory.
Tested Systems
==============
PRIMME is primary developed with GNU gcc, g++ and gfortran (versions
4.8 and later). Many users have reported builds on several other
platforms/compilers:
* SUSE 13.1 & 13.2
* CentOS 6.6
* Ubuntu 14.04
* MacOS X 10.9 & 10.10
* Cygwin & MinGW
* Cray XC30
* SunOS 5.9, quad processor Sun-Fire-280R, and several other
UltraSparcs
* AIX 5.2 IBM SP POWER 3+, 16-way SMP, 375 MHz nodes (seaborg at
nersc.gov)
C Library Interface
*******************
The PRIMME interface is composed of the following functions. To solve
real symmetric and Hermitian standard eigenproblems call respectively:
int dprimme(double *evals, double *evecs, double *resNorms,
primme_params *primme);
int zprimme(double *evals, Complex_Z *evecs, double *resNorms,
primme_params *primme);
Other useful functions:
void primme_initialize(primme_params *primme);
int primme_set_method(primme_preset_method method,
primme_params *params);
void primme_display_params(primme_params primme);
void primme_Free(primme_params primme);
PRIMME stores its data on the structure "primme_params". See
Parameters Guide for an introduction about its fields.
Running
=======
To use PRIMME, follow this basic steps.
1. Include:
#include "primme.h" /* header file is required to run primme */
2. Initialize a PRIMME parameters structure for default settings:
primme_params primme;
primme_initialize(&primme);
3. Set problem parameters (see also Parameters Guide), and,
optionally, set one of the "preset methods":
primme.matrixMatvec = LaplacianMatrixMatvec; /* MV product */
primme.n = 100; /* set problem dimension */
primme.numEvals = 10; /* Number of wanted eigenpairs */
ret = primme_set_method(method, &primme);
...
4. Then to solve a real symmetric standard eigenproblems call:
ret = dprimme(evals, evecs, resNorms, &primme);
To solve Hermitian standard eigenproblems call:
ret = zprimme(evals, evecs, resNorms, &primme);
The call arguments are:
* *evals*, array to return the found eigenvalues;
* *evecs*, array to return the found eigenvectors;
* *resNorms*, array to return the residual norms of the found
eigenpairs; and
* *ret*, returned error code.
5. Before exiting, free the work arrays in PRIMME:
primme_Free(&primme);
Parameters Guide
================
PRIMME stores the data on the structure "primme_params", which has the
next fields:
/* Basic */
int n; // matrix dimension
void (*matrixMatvec)(...); // matrix-vector product
int numEvals; // how many eigenpairs to find
primme_target target; // which eigenvalues to find
int numTargetShifts; // for targeting interior eigenpairs
double *targetShifts;
double eps; // tolerance of the converged eigenpairs
/* For parallel programs */
int numProcs;
int procID;
int nLocal;
void (*globalSumDouble)(...);
/* Accelerate the convergence */
void (*applyPreconditioner)(...); // precond-vector product
int initSize; // initial vectors as approximate solutions
int maxBasisSize;
int minRestartSize;
int maxBlockSize;
/* User data */
void *commInfo;
void *matrix;
void *preconditioner;
/* Advanced options */
int numOrthoConst; // orthogonal constrains to the eigenvectors
int dynamicMethodSwitch;
int locking;
int maxMatvecs;
int maxOuterIterations;
int intWorkSize;
long int realWorkSize;
int iseed[4];
int *intWork;
void *realWork;
double aNorm;
int printLevel;
FILE *outputFile;
double *ShiftsForPreconditioner;
struct restarting_params restartingParams;
struct correction_params correctionParams;
struct primme_stats stats;
struct stackTraceNode *stackTrace
PRIMME requires the user to set at least the dimension of the matrix
("n") and the matrix-vector product ("matrixMatvec"), as they define
the problem to be solved. For parallel programs, "nLocal", "procID"
and "globalSumDouble" are also required.
In addition, most users would want to specify how many eigenpairs to
find, and provide a preconditioner (if available).
It is useful to have set all these before calling
"primme_set_method()". Also, if users have a preference on
"maxBasisSize", "maxBlockSize", etc, they should also provide them
into "primme_params" prior to the "primme_set_method()" call. This
helps "primme_set_method()" make the right choice on other parameters.
It is sometimes useful to check the actual parameters that PRIMME is
going to use (before calling it) or used (on return) by printing them
with "primme_display_params()".
Interface Description
=====================
The next enumerations and functions are declared in "primme.h".
dprimme
-------
int dprimme(double *evals, double *evecs, double *resNorms, primme_params *primme)
Solve a real symmetric standard eigenproblem.
Parameters:
* **evals** -- array at least of size "numEvals" to store the
computed eigenvalues; all processes in a parallel run return
this local array with the same values.
* **resNorms** -- array at least of size "numEvals" to store
the residual norms of the computed eigenpairs; all processes
in parallel run return this local array with the same values.
* **evecs** -- array at least of size "nLocal" times
"numEvals" to store columnwise the (local part of the)
computed eigenvectors.
* **primme** -- parameters structure.
Returns:
error indicator; see Error Codes.
zprimme
-------
int zprimme(double *evals, Complex_Z *evecs, double *resNorms, primme_params *primme)
Solve a Hermitian standard eigenproblem; see function "dprimme()".
Note: PRIMME uses a structure called "Complex_Z" to define
complex numbers. "Complex_Z" is defined in
"PRIMMESRC/COMMONSRC/Complexz.h". In future versions of PRIMME,
"Complex_Z" will be replaced by "complex double" from the C99
standard. Because the two types are binary compatible, we
strongly recommend that calling programs use the C99 type to
maintain future compatibility. See examples in "TEST" such as
"ex_zseq.c" and "ex_zseqf77.c".
primme_initialize
-----------------
void primme_initialize(primme_params *primme)
Set PRIMME parameters structure to the default values.
Parameters:
* **primme** -- parameters structure.
primme_set_method
-----------------
int primme_set_method(primme_preset_method method, primme_params *primme)
Set PRIMME parameters to one of the preset configurations.
Parameters:
* **method** --
preset configuration; one of
"DYNAMIC"
"DEFAULT_MIN_TIME"
"DEFAULT_MIN_MATVECS"
"Arnoldi"
"GD"
"GD_plusK"
"GD_Olsen_plusK"
"JD_Olsen_plusK"
"RQI"
"JDQR"
"JDQMR"
"JDQMR_ETol"
"SUBSPACE_ITERATION"
"LOBPCG_OrthoBasis"
"LOBPCG_OrthoBasis_Window"
* **primme** -- parameters structure.
See also Preset Methods.
primme_display_params
---------------------
void primme_display_params(primme_params primme)
Display all printable settings of "primme" into the file descriptor
"outputFile".
Parameters:
* **primme** -- parameters structure.
primme_Free
-----------
void primme_Free(primme_params *primme)
Free memory allocated by PRIMME.
Parameters:
* **primme** -- parameters structure.
FORTRAN Library Interface
*************************
The next enumerations and functions are declared in "primme_f77.h".
ptr
Fortran datatype with the same size as a pointer. Use "integer*4"
when compiling in 32 bits and "integer*8" in 64 bits.
primme_initialize_f77
=====================
primme_initialize_f77(primme)
Set PRIMME parameters structure to the default values.
Parameters:
* **primme** (*ptr*) -- (output) parameters structure.
primme_set_method_f77
=====================
primme_set_method_f77(method, primme, ierr)
Set PRIMME parameters to one of the preset configurations.
Parameters:
* **method** (*integer*) --
(input) preset configuration. One of:
"PRIMMEF77_DYNAMIC"
"PRIMMEF77_DEFAULT_MIN_TIME"
"PRIMMEF77_DEFAULT_MIN_MATVECS"
"PRIMMEF77_Arnoldi"
"PRIMMEF77_GD"
"PRIMMEF77_GD_plusK"
"PRIMMEF77_GD_Olsen_plusK"
"PRIMMEF77_JD_Olsen_plusK"
"PRIMMEF77_RQI"
"PRIMMEF77_JDQR"
"PRIMMEF77_JDQMR"
"PRIMMEF77_JDQMR_ETol"
"PRIMMEF77_SUBSPACE_ITERATION"
"PRIMMEF77_LOBPCG_OrthoBasis"
"PRIMMEF77_LOBPCG_OrthoBasis_Window"
See "primme_preset_method".
* **primme** (*ptr*) -- (input) parameters structure.
* **ierr** (*integer*) -- (output) if 0, successful; if
negative, something went wrong.
primme_Free_f77
===============
primme_Free_f77(primme)
Free memory allocated by PRIMME.
Parameters:
* **primme** (*ptr*) -- parameters structure.
dprimme_f77
===========
dprimme_f77(evals, evecs, resNorms, primme, ierr)
Solve a real symmetric standard eigenproblem.
Parameters:
* **evals(*)** (*double precision*) -- (output) array at least
of size "numEvals" to store the computed eigenvalues; all
parallel calls return the same value in this array.
* **resNorms(*)** (*double precision*) -- (output) array at
least of size "numEvals" to store the residual norms of the
computed eigenpairs; all parallel calls return the same value
in this array.
* **evecs(*)** (*double precision*) -- (input/output) array at
least of size "nLocal" times "numEvals" to store columnwise
the (local part of the) computed eigenvectors.
* **primme** (*ptr*) -- parameters structure.
* **ierr** (*integer*) -- (output) error indicator; see Error
Codes.
zprimme_f77
===========
zprimme_f77(evals, evecs, resNorms, primme, ierr)
Solve a Hermitian standard eigenproblem. The arguments have the
same meaning as in function "dprimme_f77()".
Parameters:
* **evals(*)** (*double precision*) -- (output)
* **resNorms(*)** (*double precision*) -- (output)
* **evecs(*)** (*complex double precision*) -- (input/output)
* **primme** (*ptr*) -- (input) parameters structure.
* **ierr** (*integer*) -- (output) error indicator; see Error
Codes.
primmetop_set_member_f77
========================
primmetop_set_member_f77(primme, label, value)
Set a value in some field of the parameter structure.
Parameters:
* **primme** (*ptr*) -- (input) parameters structure.
* **label** (*integer*) --
field where to set value. One of:
"PRIMMEF77_n"
"PRIMMEF77_matrixMatvec"
"PRIMMEF77_applyPreconditioner"
"PRIMMEF77_numProcs"
"PRIMMEF77_procID"
"PRIMMEF77_commInfo"
"PRIMMEF77_nLocal"
"PRIMMEF77_globalSumDouble"
"PRIMMEF77_numEvals"
"PRIMMEF77_target"
"PRIMMEF77_numTargetShifts"
"PRIMMEF77_targetShifts"
"PRIMMEF77_locking"
"PRIMMEF77_initSize"
"PRIMMEF77_numOrthoConst"
"PRIMMEF77_maxBasisSize"
"PRIMMEF77_minRestartSize"
"PRIMMEF77_maxBlockSize"
"PRIMMEF77_maxMatvecs"
"PRIMMEF77_maxOuterIterations"
"PRIMMEF77_intWorkSize"
"PRIMMEF77_realWorkSize"
"PRIMMEF77_iseed"
"PRIMMEF77_intWork"
"PRIMMEF77_realWork"
"PRIMMEF77_aNorm"
"PRIMMEF77_eps"
"PRIMMEF77_printLevel"
"PRIMMEF77_outputFile"
"PRIMMEF77_matrix"
"PRIMMEF77_preconditioner"
"PRIMMEF77_restartingParams_scheme".
"PRIMMEF77_restartingParams_maxPrevRetain"
"PRIMMEF77_correctionParams_precondition"
"PRIMMEF77_correctionParams_robustShifts"
"PRIMMEF77_correctionParams_maxInnerIterations"
"PRIMMEF77_correctionParams_projectors_LeftQ"
"PRIMMEF77_correctionParams_projectors_LeftX"
"PRIMMEF77_correctionParams_projectors_RightQ"
"PRIMMEF77_correctionParams_projectors_RightX"
"PRIMMEF77_correctionParams_projectors_SkewQ"
"PRIMMEF77_correctionParams_projectors_SkewX"
"PRIMMEF77_correctionParams_convTest"
"PRIMMEF77_correctionParams_relTolBase"
"PRIMMEF77_stats_numOuterIterations"
"PRIMMEF77_stats_numRestarts"
"PRIMMEF77_stats_numMatvecs"
"PRIMMEF77_stats_numPreconds"
"PRIMMEF77_stats_elapsedTime"
"PRIMMEF77_dynamicMethodSwitch"
"PRIMMEF77_massMatrixMatvec"
* **value** -- (input) value to set.
Note: **Don't use** this function inside PRIMME's callback
functions, e.g., "matrixMatvec" or "applyPreconditioner", or in
functions called by these functions. In those cases use
"primme_set_member_f77()".
primmetop_get_member_f77
========================
primmetop_get_member_f77(primme, label, value)
Get the value in some field of the parameter structure.
Parameters:
* **primme** (*ptr*) -- (input) parameters structure.
* **label** (*integer*) -- (input) field where to get value.
One of the detailed in function "primmetop_set_member_f77()".
* **value** -- (output) value of the field.
Note: **Don't use** this function inside PRIMME's callback
functions, e.g., "matrixMatvec" or "applyPreconditioner", or in
functions called by these functions. In those cases use
"primme_get_member_f77()".
Note: When "label" is one of "PRIMMEF77_matrixMatvec",
"PRIMMEF77_applyPreconditioner", "PRIMMEF77_commInfo",
"PRIMMEF77_intWork", "PRIMMEF77_realWork", "PRIMMEF77_matrix" and
"PRIMMEF77_preconditioner", the returned "value" is a C pointer
("void*"). Use Fortran pointer or other extensions to deal with
it. For instance:
use iso_c_binding
MPI_Comm comm
comm = MPI_COMM_WORLD
call primme_set_member_f77(primme, PRIMMEF77_commInfo, comm)
...
subroutine par_GlobalSumDouble(x,y,k,primme)
use iso_c_binding
implicit none
...
MPI_Comm, pointer :: comm
type(c_ptr) :: pcomm
call primme_get_member_f77(primme, PRIMMEF77_commInfo, pcomm)
call c_f_pointer(pcomm, comm)
call MPI_Allreduce(x,y,k,MPI_DOUBLE,MPI_SUM,comm,ierr)
Most users would not need to retrieve these pointers in their
programs.
primmetop_get_prec_shift_f77
============================
primmetop_get_prec_shift_f77(primme, index, value)
Get the value in some position of the array
"ShiftsForPreconditioner".
Parameters:
* **primme** (*ptr*) -- (input) parameters structure.
* **index** (*integer*) -- (input) position of the array; the
first position is 1.
* **value** -- (output) value of the array at that position.
primme_set_member_f77
=====================
primme_set_member_f77(primme, label, value)
Set a value in some field of the parameter structure.
Parameters:
* **primme** (*ptr*) -- (input) parameters structure.
* **label** (*integer*) -- field where to set value. One of
the vales defined in "primmetop_set_member_f77()".
* **value** -- (input) value to set.
Note: Use this function exclusively inside PRIMME's callback
functions, e.g., "matrixMatvec" or "applyPreconditioner", or in
functions called by these functions. Otherwise, e.g., from the
main program, use the function "primmetop_set_member_f77()".
primme_get_member_f77
=====================
primme_get_member_f77(primme, label, value)
Get the value in some field of the parameter structure.
Parameters:
* **primme** (*ptr*) -- (input) parameters structure.
* **label** (*integer*) -- (input) field where to get value.
One of the detailed in function "primmetop_set_member_f77()".
* **value** -- (output) value of the field.
Note: Use this function exclusively inside PRIMME's callback
functions, e.g., "matrixMatvec" or "applyPreconditioner", or in
functions called by these functions. Otherwise, e.g., from the
main program, use the function "primmetop_get_member_f77()".
Note: When "label" is one of "PRIMMEF77_matrixMatvec",
"PRIMMEF77_applyPreconditioner", "PRIMMEF77_commInfo",
"PRIMMEF77_intWork", "PRIMMEF77_realWork", "PRIMMEF77_matrix" and
"PRIMMEF77_preconditioner", the returned "value" is a C pointer
("void*"). Use Fortran pointer or other extensions to deal with
it. For instance:
use iso_c_binding
MPI_Comm comm
comm = MPI_COMM_WORLD
call primme_set_member_f77(primme, PRIMMEF77_commInfo, comm)
...
subroutine par_GlobalSumDouble(x,y,k,primme)
use iso_c_binding
implicit none
...
MPI_Comm, pointer :: comm
type(c_ptr) :: pcomm
call primme_get_member_f77(primme, PRIMMEF77_commInfo, pcomm)
call c_f_pointer(pcomm, comm)
call MPI_Allreduce(x,y,k,MPI_DOUBLE,MPI_SUM,comm,ierr)
Most users would not need to retrieve these pointers in their
programs.
primme_get_prec_shift_f77
=========================
primme_get_prec_shift_f77(primme, index, value)
Get the value in some position of the array
"ShiftsForPreconditioner".
Parameters:
* **primme** (*ptr*) -- (input) parameters structure.
* **index** (*integer*) -- (input) position of the array; the
first position is 1.
* **value** -- (output) value of the array at that position.
Note: Use this function exclusively inside the function
"matrixMatvec", "massMatrixMatvec", or "applyPreconditioner".