1 /*******************************************************************************
2  * Copyright (c) 2018, College of William & Mary
3  * All rights reserved.
4  *
5  * Redistribution and use in source and binary forms, with or without
6  * modification, are permitted provided that the following conditions are met:
7  *     * Redistributions of source code must retain the above copyright
8  *       notice, this list of conditions and the following disclaimer.
9  *     * Redistributions in binary form must reproduce the above copyright
10  *       notice, this list of conditions and the following disclaimer in the
11  *       documentation and/or other materials provided with the distribution.
12  *     * Neither the name of the College of William & Mary nor the
13  *       names of its contributors may be used to endorse or promote products
14  *       derived from this software without specific prior written permission.
15  *
16  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
17  * ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
18  * WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
19  * DISCLAIMED. IN NO EVENT SHALL THE COLLEGE OF WILLIAM & MARY BE LIABLE FOR ANY
20  * DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
21  * (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
22  * LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
23  * ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
24  * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
25  * SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
26  *
27  * PRIMME: https://github.com/primme/primme
28  * Contact: Andreas Stathopoulos, a n d r e a s _at_ c s . w m . e d u
29  **********************************************************************
30  * File: primme_interface.c
31  *
32  * Purpose - Contains interface functions to PRIMME named primme_*
33  *           If desired, the user can call any of these functions for
34  *           initializing parameters, setting up the method,
35  *           allocating memory, or even checking a given set of parameters.
36  *
37  ******************************************************************************/
38 
39 #ifndef THIS_FILE
40 #define THIS_FILE "../eigs/primme_interface.c"
41 #endif
42 
43 #include <stdlib.h>   /* mallocs, free */
44 #include <stdio.h>
45 #include <math.h>
46 #include <string.h>  /* strcmp */
47 #include <limits.h>
48 #include "numerical.h"
49 #include "primme_interface.h"
50 
51 /* Only define these functions ones */
52 #ifdef USE_DOUBLE
53 #include "notemplate.h"
54 
55 typedef void *ptr_v;
56 typedef const char *str_v;
57 typedef union {
58    void (*matFunc_v)(void *, PRIMME_INT *, void *, PRIMME_INT *, int *,
59          struct primme_params *, int *);
60    void (*globalSumRealFunc_v)(
61          void *, void *, int *, struct primme_params *, int *);
62    void (*broadcastRealFunc_v)(void *, int *, struct primme_params *, int *);
63    void (*convTestFun_v)(
64          double *, void *, double *, int *, struct primme_params *, int *);
65    void (*monitorFun_v)(void *basisEvals, int *basisSize, int *basisFlags,
66          int *iblock, int *blockSize, void *basisNorms, int *numConverged,
67          void *lockedEvals, int *numLocked, int *lockedFlags, void *lockedNorms,
68          int *inner_its, void *LSRes, const char *msg, double *time,
69          primme_event *event, struct primme_params *primme, int *err);
70 } value_t;
71 
72 /*****************************************************************************
73  * Initialize handles also the allocation of primme structure
74  *****************************************************************************/
primme_params_create(void)75 primme_params * primme_params_create(void) {
76 
77    primme_params *primme = NULL;
78    if (MALLOC_PRIMME(1, &primme) == 0)
79       primme_initialize(primme);
80     return primme;
81 }
82 
83 /*****************************************************************************
84  *  * Free the internally allocated work arrays of the primme structure
85  *   *****************************************************************************/
primme_params_destroy(primme_params * primme)86 int primme_params_destroy(primme_params *primme) {
87     free(primme);
88     return 0;
89 }
90 
91 
92 /*******************************************************************************
93  * Subroutine primme_initialize - Set primme_params members to default values.
94  *
95  * INPUT/OUTPUT PARAMETERS
96  * ----------------------------------
97  * primme  Structure containing various solver parameters and statistics
98  *
99  ******************************************************************************/
100 
primme_initialize(primme_params * primme)101 void primme_initialize(primme_params *primme) {
102 
103    /* Essential parameters */
104    primme->n                       = 0;
105    primme->numEvals                = 1;
106    primme->target                  = primme_smallest;
107    primme->aNorm                   = 0.0L;
108    primme->BNorm                   = 0.0L;
109    primme->invBNorm                = 0.0L;
110    primme->eps                     = 0.0;
111 
112    /* Matvec and preconditioner */
113    primme->matrixMatvec            = NULL;
114    primme->matrixMatvec_type       = primme_op_default;
115    primme->applyPreconditioner     = NULL;
116    primme->applyPreconditioner_type= primme_op_default;
117    primme->massMatrixMatvec        = NULL;
118    primme->massMatrixMatvec_type   = primme_op_default;
119 
120    /* Shifts for interior eigenvalues*/
121    primme->numTargetShifts         = 0;
122    primme->targetShifts            = NULL;
123 
124    /* Parallel computing parameters */
125    primme->numProcs                = 1;
126    primme->procID                  = 0;
127    primme->nLocal                  = -1;
128    primme->commInfo                = NULL;
129    primme->globalSumReal           = NULL;
130    primme->globalSumReal_type      = primme_op_default;
131    primme->broadcastReal           = NULL;
132    primme->broadcastReal_type      = primme_op_default;
133 
134    /* Initial guesses/constraints */
135    primme->initSize                = 0;
136    primme->numOrthoConst           = 0;
137 
138    primme->projectionParams.projection = primme_proj_default;
139 
140    primme->initBasisMode                       = primme_init_default;
141 
142    /* Eigensolver parameters (outer) */
143    primme->locking                             = -1;
144    primme->dynamicMethodSwitch                 = -1;
145    primme->maxBasisSize                        = 0;
146    primme->minRestartSize                      = 0;
147    primme->maxBlockSize                        = 0;
148    primme->maxMatvecs                          = INT_MAX;
149    primme->maxOuterIterations                  = INT_MAX;
150    primme->restartingParams.maxPrevRetain      = -1;
151    primme->orth                                = primme_orth_default;
152    primme->internalPrecision                   = primme_op_default;
153 
154    /* correction parameters (inner) */
155    primme->correctionParams.precondition       = -1;
156    primme->correctionParams.robustShifts       = 0;
157    primme->correctionParams.maxInnerIterations = -INT_MAX;
158    primme->correctionParams.projectors.LeftQ   = 0;
159    primme->correctionParams.projectors.LeftX   = 0;
160    primme->correctionParams.projectors.RightQ  = 0;
161    primme->correctionParams.projectors.RightX  = 0;
162    primme->correctionParams.projectors.SkewQ   = 0;
163    primme->correctionParams.projectors.SkewX   = 0;
164    primme->correctionParams.relTolBase         = 0;
165    primme->correctionParams.convTest           = primme_adaptive_ETolerance;
166 
167    /* Printing and reporting */
168    primme->outputFile                          = stdout;
169    primme->printLevel                          = 1;
170    primme->stats.numOuterIterations            = 0;
171    primme->stats.numRestarts                   = 0;
172    primme->stats.numMatvecs                    = 0;
173    primme->stats.numPreconds                   = 0;
174    primme->stats.numGlobalSum                  = 0;
175    primme->stats.flopsDense                    = 0.0;
176    primme->stats.timeDense                     = 0;
177    primme->stats.volumeGlobalSum               = 0;
178    primme->stats.numBroadcast                  = 0;
179    primme->stats.volumeBroadcast               = 0;
180    primme->stats.numOrthoInnerProds            = 0.0;
181    primme->stats.elapsedTime                   = 0.0;
182    primme->stats.timeMatvec                    = 0.0;
183    primme->stats.timePrecond                   = 0.0;
184    primme->stats.timeOrtho                     = 0.0;
185    primme->stats.timeGlobalSum                 = 0.0;
186    primme->stats.timeBroadcast                 = 0.0;
187    primme->stats.estimateMinEVal               = -HUGE_VAL;
188    primme->stats.estimateMaxEVal               = HUGE_VAL;
189    primme->stats.estimateLargestSVal           = -HUGE_VAL;
190    primme->stats.estimateBNorm                 = -HUGE_VAL;
191    primme->stats.estimateInvBNorm              = -HUGE_VAL;
192    primme->stats.maxConvTol                    = 0.0;
193    primme->stats.estimateResidualError         = 0.0;
194    primme->stats.lockingIssue                  = 0;
195 
196    /* Optional user defined structures */
197    primme->matrix                  = NULL;
198    primme->preconditioner          = NULL;
199    primme->massMatrix              = NULL;
200 
201    /* Internally used variables */
202    primme->iseed[0] = -1;   /* To set iseed, we first need procID           */
203    primme->iseed[1] = -1;   /* Thus we set all iseeds to -1                 */
204    primme->iseed[2] = -1;   /* Unless users provide their own iseeds,       */
205    primme->iseed[3] = -1;   /* PRIMME will set thse later uniquely per proc */
206    primme->ShiftsForPreconditioner = NULL;
207    primme->convTestFun             = NULL;
208    primme->convTestFun_type        = primme_op_default;
209    primme->convtest                = NULL;
210    primme->ldevecs                 = -1;
211    primme->ldOPs                   = -1;
212    primme->monitorFun              = NULL;
213    primme->monitorFun_type         = primme_op_default;
214    primme->monitor                 = NULL;
215    primme->queue                   = NULL;
216    primme->profile                 = NULL;
217 }
218 
219 /*******************************************************************************
220  * Subroutine primme_free - Free memory resources allocated by Xprimme.
221  *
222  * INPUT/OUTPUT PARAMETERS
223  * ----------------------------------
224  * primme  Structure containing various solver parameters and statistics
225  *
226  ******************************************************************************/
227 
primme_free(primme_params * primme)228 void primme_free(primme_params *primme) {
229    (void)primme;
230    /* No function */
231 }
232 
233 /******************************************************************************
234  * Subroutine primme_set_method - Set the eigensolver parameters to implement a
235  *    method requested by the user.
236  *    A choice of 15 preset methods is provided. These implement 12 well
237  *    known methods. The two default methods are shells for easy access to
238  *    the best performing GD+k and JDQMR_ETol with expertly chosen parameters.
239  *
240  *    The DYNAMIC method makes runtime measurements and switches dynamically
241  *    between DEFAULT_MIN_TIME and DEFAULT_MIN_MATVEC keeping the one that
242  *    performs the best. Because it usually achieves best runtime over all
243  *    methods, it is recommended for the general user.
244  *
245  *    For most methods the user may specify the maxBasisSize, restart size,
246  *    block size, etc. If any (or all) of these parameters are not specified,
247  *    they will be given default values that are appropriate for the method.
248  *
249  *    primme_set_method() will override those parameters in primme that
250  *    are needed to implement the method.
251  *
252  *    Note: Spectral transformations can be applied by simply providing
253  *         (A-shift I)^{-1} as the matvec.
254  *
255  *    For additional information see the readme file
256  *
257  * INPUT
258  * -----
259  *    method    One of the following 12 enum methods:
260  *
261  *        DYNAMIC,                 : Switches dynamically to the best method
262  *        DEFAULT_MIN_TIME,        : Currently set at JDQMR_ETol
263  *        DEFAULT_MIN_MATVECS,     : Currently set at GD+block
264  *        Arnoldi,                 : obviously not an efficient choice
265  *        GD,                      : classical block Generalized Davidson
266  *        GD_plusK,                : GD+k block GD with recurrence restarting
267  *        GD_Olsen_plusK,          : GD+k with approximate Olsen precond.
268  *        JD_Olsen_plusK,          : GD+k, exact Olsen (two precond per step)
269  *        RQI,                     : Rayleigh Quotient Iteration. Also INVIT,
270  *                                 :   but for INVIT provide targetShifts
271  *        JDQR,                    : Original block, Jacobi Davidson
272  *        JDQMR,                   : Our block JDQMR method (similar to JDCG)
273  *        JDQMR_ETol,              : Slight, but efficient JDQMR modification
274  *        STEEPEST_DESCENT,      : equiv. to GD(block,2*block)
275  *        LOBPCG_OrthoBasis,       : equiv. to GD(nev,3*nev)+nev
276  *        LOBPCG_OrthoBasis_Window : equiv. to GD(block,3*block)+block nev>block
277  *
278  *
279  * INPUT/OUTPUT
280  * ------
281  *    primme    The main structure to be used by Primme, with parameters
282  *              reflecting the user's choice of method
283  *
284  *
285  * return value         Note: the return value is a LONG INT
286  *
287  *      = 0 parameters for "method" have been set.
288  *
289  *      < 0 no such method exists. If not defined by the user, defaults have
290  *          been set for maxBasisSize, minRestartSize, and maxBlockSize.
291  *
292  ******************************************************************************/
primme_set_method(primme_preset_method method,primme_params * primme)293 int primme_set_method(primme_preset_method method, primme_params *primme) {
294 
295    /* Set default method as DYNAMIC */
296    if (method == PRIMME_DEFAULT_METHOD)
297       method = PRIMME_DYNAMIC;
298 
299    /* From our experience, these two methods yield the smallest matvecs/time */
300    /* DYNAMIC will make some timings before it settles on one of the two     */
301    if (method == PRIMME_DEFAULT_MIN_MATVECS) {
302       method = PRIMME_GD_Olsen_plusK;
303    }
304    else if (method == PRIMME_DEFAULT_MIN_TIME) {
305       /* JDQMR works better than JDQMR_ETol in interior problems. */
306       if (primme->target == primme_smallest || primme->target == primme_largest) {
307          method = PRIMME_JDQMR_ETol;
308       }
309       else {
310          method = PRIMME_JDQMR;
311       }
312    }
313    if (method == PRIMME_DYNAMIC) {
314       primme->dynamicMethodSwitch = 1;
315    }
316    else {
317       primme->dynamicMethodSwitch = 0;
318    }
319 
320    if (primme->maxBlockSize == 0) {
321       primme->maxBlockSize = 1;
322    }
323    if (primme->correctionParams.precondition == -1) {
324       primme->correctionParams.precondition = primme->applyPreconditioner ? 1 : 0;
325    }
326 
327    if (method == PRIMME_Arnoldi) {
328       primme->restartingParams.maxPrevRetain      = 0;
329       primme->correctionParams.precondition       = 0;
330       primme->correctionParams.maxInnerIterations = 0;
331    }
332    else if (method == PRIMME_GD) {
333       primme->restartingParams.maxPrevRetain      = 0;
334       primme->correctionParams.robustShifts       = 1;
335       primme->correctionParams.maxInnerIterations = 0;
336       primme->correctionParams.projectors.RightX  = 0;
337       primme->correctionParams.projectors.SkewX   = 0;
338    }
339    else if (method == PRIMME_GD_plusK) {
340       if (primme->restartingParams.maxPrevRetain <= 0) {
341          if ((primme->maxBlockSize == 1 && primme->numEvals > 1) ||
342                primme->massMatrixMatvec) {
343             primme->restartingParams.maxPrevRetain = 2;
344          }
345          else {
346             primme->restartingParams.maxPrevRetain = primme->maxBlockSize;
347          }
348       }
349       primme->correctionParams.maxInnerIterations = 0;
350       primme->correctionParams.projectors.RightX  = 0;
351       primme->correctionParams.projectors.SkewX   = 0;
352    }
353    else if (method == PRIMME_GD_Olsen_plusK) {
354       if (primme->restartingParams.maxPrevRetain <= 0) {
355          if ((primme->maxBlockSize == 1 && primme->numEvals > 1) ||
356                primme->massMatrixMatvec) {
357             primme->restartingParams.maxPrevRetain = 2;
358          }
359          else {
360             primme->restartingParams.maxPrevRetain = primme->maxBlockSize;
361          }
362       }
363       primme->correctionParams.maxInnerIterations = 0;
364       primme->correctionParams.projectors.RightX  = 1;
365       primme->correctionParams.projectors.SkewX   = 0;
366    }
367    else if (method == PRIMME_JD_Olsen_plusK) {
368       if (primme->restartingParams.maxPrevRetain <= 0) {
369          if ((primme->maxBlockSize == 1 && primme->numEvals > 1) ||
370                primme->massMatrixMatvec) {
371             primme->restartingParams.maxPrevRetain = 2;
372          }
373          else {
374             primme->restartingParams.maxPrevRetain = primme->maxBlockSize;
375          }
376       }
377       primme->correctionParams.robustShifts       = 1;
378       primme->correctionParams.maxInnerIterations = 0;
379       primme->correctionParams.projectors.RightX  = 1;
380       primme->correctionParams.projectors.SkewX   = 1;
381    }
382    else if (method == PRIMME_RQI) {
383      /* This is accelerated RQI as basisSize may be > 2                 */
384      /* NOTE: If numTargetShifts > 0 and targetShifts[*] are provided,  *
385       * the interior problem solved uses these shifts in the correction *
386       * equation. Therefore RQI becomes INVIT in that case.             */
387       primme->locking                             = 1;
388       primme->restartingParams.maxPrevRetain      = 0;
389       primme->correctionParams.robustShifts       = 1;
390       primme->correctionParams.maxInnerIterations = -1;
391       primme->correctionParams.projectors.LeftQ   = 1;
392       primme->correctionParams.projectors.LeftX   = 1;
393       primme->correctionParams.projectors.RightQ  = 0;
394       primme->correctionParams.projectors.RightX  = 1;
395       primme->correctionParams.projectors.SkewQ   = 0;
396       primme->correctionParams.projectors.SkewX   = 0;
397       primme->correctionParams.convTest           = primme_full_LTolerance;
398    }
399    else if (method == PRIMME_JDQR) {
400       primme->locking                             = 1;
401       primme->restartingParams.maxPrevRetain      = 1;
402       primme->correctionParams.robustShifts       = 0;
403       if (primme->correctionParams.maxInnerIterations == -INT_MAX) {
404          primme->correctionParams.maxInnerIterations = 10;
405       }
406       primme->correctionParams.projectors.LeftQ   = 0;
407       primme->correctionParams.projectors.LeftX   = 1;
408       primme->correctionParams.projectors.RightQ  = 1;
409       primme->correctionParams.projectors.RightX  = 1;
410       primme->correctionParams.projectors.SkewQ   = 1;
411       primme->correctionParams.projectors.SkewX   = 1;
412       primme->correctionParams.relTolBase         = 1.5;
413       primme->correctionParams.convTest           = primme_full_LTolerance;
414    }
415    else if (method == PRIMME_JDQMR) {
416       if (primme->restartingParams.maxPrevRetain < 0) {
417          primme->restartingParams.maxPrevRetain   = 1;
418       }
419       primme->correctionParams.maxInnerIterations = -1;
420       if (primme->correctionParams.precondition) {
421          primme->correctionParams.projectors.LeftQ   = 1;
422       }
423       else {
424          primme->correctionParams.projectors.LeftQ   = 0;
425       }
426       primme->correctionParams.projectors.LeftX   = 1;
427       primme->correctionParams.projectors.RightQ  = 0;
428       primme->correctionParams.projectors.RightX  = 0;
429       primme->correctionParams.projectors.SkewQ   = 0;
430       primme->correctionParams.projectors.SkewX   = 1;
431       primme->correctionParams.convTest           = primme_adaptive;
432    }
433    else if (method == PRIMME_JDQMR_ETol) {
434       if (primme->restartingParams.maxPrevRetain < 0) {
435          primme->restartingParams.maxPrevRetain   = 1;
436       }
437       primme->correctionParams.maxInnerIterations = -1;
438       if (primme->correctionParams.precondition) {
439          primme->correctionParams.projectors.LeftQ   = 1;
440       }
441       else {
442          primme->correctionParams.projectors.LeftQ   = 0;
443       }
444       primme->correctionParams.projectors.LeftX   = 1;
445       primme->correctionParams.projectors.RightQ  = 0;
446       primme->correctionParams.projectors.RightX  = 0;
447       primme->correctionParams.projectors.SkewQ   = 0;
448       primme->correctionParams.projectors.SkewX   = 0;
449       primme->correctionParams.convTest           = primme_adaptive_ETolerance;
450    }
451    else if (method == PRIMME_STEEPEST_DESCENT) {
452       primme->locking                             = 1;
453       primme->maxBasisSize                        = primme->numEvals*2;
454       primme->minRestartSize                      = primme->numEvals;
455       primme->maxBlockSize                        = primme->numEvals;
456       primme->restartingParams.maxPrevRetain      = 0;
457       primme->correctionParams.robustShifts       = 0;
458       primme->correctionParams.maxInnerIterations = 0;
459       primme->correctionParams.projectors.RightX  = 1;
460       primme->correctionParams.projectors.SkewX   = 0;
461    }
462    else if (method == PRIMME_LOBPCG_OrthoBasis) {
463       primme->maxBasisSize                        = primme->numEvals*3;
464       primme->minRestartSize                      = primme->numEvals;
465       primme->maxBlockSize                        = primme->numEvals;
466       primme->restartingParams.maxPrevRetain      = primme->numEvals;
467       primme->correctionParams.robustShifts       = 0;
468       primme->correctionParams.maxInnerIterations = 0;
469       primme->correctionParams.projectors.RightX  = 1;
470       primme->correctionParams.projectors.SkewX   = 0;
471       primme->initBasisMode                       = primme_init_random;
472    }
473    else if (method == PRIMME_LOBPCG_OrthoBasis_Window) {
474       /* Observed needing to restart with two vectors at least to converge    */
475       /* in some tests like for instance testi-10-LOBPCG_OrthoBasis_Window-3- */
476       /* primme_closest_geq-primme_proj_refined.F                             */
477       if (primme->maxBlockSize == 1
478             && (primme->target == primme_closest_leq
479                || primme->target == primme_closest_geq)) {
480          primme->maxBasisSize                        = 4;
481          primme->minRestartSize                      = 2;
482          primme->restartingParams.maxPrevRetain      = 1;
483       }
484       else {
485          primme->maxBasisSize                        = primme->maxBlockSize*3;
486          primme->minRestartSize                      = primme->maxBlockSize;
487          primme->restartingParams.maxPrevRetain      = primme->maxBlockSize;
488       }
489       primme->correctionParams.robustShifts       = 0;
490       primme->correctionParams.maxInnerIterations = 0;
491       primme->correctionParams.projectors.RightX  = 1;
492       primme->correctionParams.projectors.SkewX   = 0;
493       primme->initBasisMode                       = primme_init_random;
494    }
495    else if (method == PRIMME_DYNAMIC) {
496       if (primme->restartingParams.maxPrevRetain <= 0) {
497          if ((primme->maxBlockSize == 1 && primme->numEvals > 1) ||
498                primme->massMatrixMatvec) {
499             primme->restartingParams.maxPrevRetain = 2;
500          }
501          else {
502             primme->restartingParams.maxPrevRetain = primme->maxBlockSize;
503          }
504       }
505       primme->correctionParams.maxInnerIterations = -1;
506       if (primme->correctionParams.precondition) {
507          primme->correctionParams.projectors.LeftQ   = 1;
508       }
509       else {
510          primme->correctionParams.projectors.LeftQ   = 0;
511       }
512       primme->correctionParams.projectors.LeftX   = 1;
513       primme->correctionParams.projectors.RightQ  = 0;
514       primme->correctionParams.projectors.RightX  = 0;
515       primme->correctionParams.projectors.SkewQ   = 0;
516       primme->correctionParams.projectors.SkewX   = 0;
517       if (primme->target == primme_smallest ||
518             primme->target == primme_largest) {
519          primme->correctionParams.convTest = primme_adaptive_ETolerance;
520       } else {
521          primme->correctionParams.convTest = primme_adaptive;
522       }
523    } else {
524       return -1;
525    }
526 
527    primme_set_defaults(primme);
528 
529    return 0;
530 
531 }
532 
533 /******************************************************************************
534  * Subroutine primme_set_defaults - Set valid values for options that still
535  *    has the initial invalid value set in primme_initialize.
536  *
537  * INPUT/OUTPUT PARAMETERS
538  * ----------------------------------
539  * primme  Structure containing various solver parameters and statistics
540  *
541  ******************************************************************************/
542 
primme_set_defaults(primme_params * primme)543 void primme_set_defaults(primme_params *primme) {
544    if (primme->dynamicMethodSwitch < 0) {
545       primme_set_method(PRIMME_DYNAMIC, primme);
546    }
547 
548    if (primme->ldevecs == -1 && primme->nLocal != -1)
549       primme->ldevecs = primme->nLocal;
550    if (primme->projectionParams.projection == primme_proj_default)
551       primme->projectionParams.projection = primme_proj_RR;
552    if (primme->initBasisMode == primme_init_default)
553       primme->initBasisMode = primme_init_krylov;
554 
555    /* Now that most of the parameters have been set, set defaults  */
556    /* for basisSize, restartSize (for those methods that need it)  */
557    /* For interior, larger basisSize and restartSize are advisable */
558    /* If maxBlockSize is provided, assign at least 4*blocksize     */
559    /* and consider also minRestartSize and maxPrevRetain           */
560    if (primme->maxBasisSize == 0) {
561       if (primme->target==primme_smallest || primme->target==primme_largest)
562          primme->maxBasisSize =
563                min(primme->n - primme->numOrthoConst,
564                      max(max(15, 4 * primme->maxBlockSize +
565                                        primme->restartingParams.maxPrevRetain),
566                            (int)2.5 * primme->minRestartSize +
567                                  primme->restartingParams.maxPrevRetain));
568       else
569          primme->maxBasisSize =
570                min(primme->n - primme->numOrthoConst,
571                      max(max(35, 5 * primme->maxBlockSize +
572                                        primme->restartingParams.maxPrevRetain),
573                            (int)1.7 * primme->minRestartSize +
574                                  primme->restartingParams.maxPrevRetain));
575    }
576 
577    if (primme->minRestartSize == 0) {
578       if (primme->n <= 3)
579          primme->minRestartSize = primme->n - primme->numOrthoConst;
580       else if (primme->target == primme_smallest ||
581                primme->target == primme_largest)
582          primme->minRestartSize = (int) (0.5 + 0.4*primme->maxBasisSize);
583       else
584          primme->minRestartSize = (int) (0.5 + 0.6*primme->maxBasisSize);
585 
586       /* Adjust so that an integer number of blocks are added between restarts*/
587       /* restart=basis-block*ceil((basis-restart-prevRetain)/block)-prevRetain*/
588       if (primme->maxBlockSize > 1) {
589          if (primme->restartingParams.maxPrevRetain > 0)
590             primme->minRestartSize = primme->maxBasisSize-primme->maxBlockSize*
591             (1 + (int) ((primme->maxBasisSize - primme->minRestartSize - 1
592                          -primme->restartingParams.maxPrevRetain ) / (double)
593             primme->maxBlockSize) ) - primme->restartingParams.maxPrevRetain ;
594          else
595             primme->minRestartSize = primme->maxBasisSize-primme->maxBlockSize*
596             (1 + (int) ((primme->maxBasisSize - primme->minRestartSize - 1)
597                         /(double) primme->maxBlockSize) );
598       }
599    }
600 
601    /* --------------------------------------------------------------------- */
602    /* Decide on whether to use locking (hard locking), or not (soft locking)*/
603    /* --------------------------------------------------------------------- */
604    if (primme->locking >= 0) {
605       /* Honor the user setup (do nothing) */
606    }
607    else if (primme->target != primme_smallest && primme->target != primme_largest) {
608        primme->locking = 1;
609    }
610    else if (primme->numEvals > primme->minRestartSize) {
611       /* use locking when not enough vectors to restart with */
612       primme->locking = 1;
613    }
614    else {
615       primme->locking = 0;
616    }
617 }
618 
619 /*******************************************************************************
620  * Subroutine primme_display_params - Displays the current configuration of
621  *    primme data structure.
622  *
623  * INPUT PARAMETERS
624  * ----------------------------------
625  * primme  Structure containing various solver parameters and statistics
626  *
627  ******************************************************************************/
628 
primme_display_params(primme_params primme)629 void primme_display_params(primme_params primme) {
630 
631    fprintf(primme.outputFile,
632            "// ---------------------------------------------------\n"
633            "//                 primme configuration               \n"
634            "// ---------------------------------------------------\n");
635 
636    primme_display_params_prefix("primme", primme);
637    fflush(primme.outputFile);
638 }
639 
640 /*******************************************************************************
641  * Subroutine primme_params_prefix - Displays the current configuration of
642  *    primme data structure. The options are printed as
643  *    prefix.optionName = value.
644  *
645  * INPUT PARAMETERS
646  * ----------------------------------
647  * prefix  Name of the structure
648  * primme  Structure containing various solver parameters and statistics
649  *
650  ******************************************************************************/
651 
primme_display_params_prefix(const char * prefix,primme_params primme)652 void primme_display_params_prefix(const char* prefix, primme_params primme) {
653 
654    int i;
655    FILE *outputFile = primme.outputFile;
656 
657 #define PRINT(P,L) fprintf(outputFile, "%s." #P " = " #L "\n", prefix, primme. P);
658 #define PRINTIF(P,V) if (primme. P == V) fprintf(outputFile, "%s." #P " = " #V "\n", prefix);
659 #define PRINTParams(P,S,L) fprintf(outputFile, "%s." #P "." #S " = " #L "\n", \
660                                     prefix, primme. P ## Params. S);
661 #define PRINTParamsIF(P,S,V) if (primme. P ## Params. S == V) \
662                                  fprintf(outputFile, "%s." #P "." #S " = " #V "\n", prefix);
663 #define PRINT_PRIMME_INT(P) fprintf(outputFile, "%s." #P " = %" PRIMME_INT_P "\n", prefix, primme. P);
664 
665    PRINT_PRIMME_INT(n);
666    PRINT_PRIMME_INT(nLocal);
667    PRINT(numProcs, %d);
668    PRINT(procID, %d);
669 
670    fprintf(outputFile, "\n// Output and reporting\n");
671    PRINT(printLevel, %d);
672 
673    fprintf(outputFile, "\n// Solver parameters\n");
674    PRINT(numEvals, %d);
675    PRINT(aNorm, %e);
676    PRINT(BNorm, %e);
677    PRINT(invBNorm, %e);
678    PRINT(eps, %e);
679    PRINT(maxBasisSize, %d);
680    PRINT(minRestartSize, %d);
681    PRINT(maxBlockSize, %d);
682    PRINT_PRIMME_INT(maxOuterIterations);
683    PRINT_PRIMME_INT(maxMatvecs);
684 
685    PRINTIF(target, primme_smallest);
686    PRINTIF(target, primme_largest);
687    PRINTIF(target, primme_closest_geq);
688    PRINTIF(target, primme_closest_leq);
689    PRINTIF(target, primme_closest_abs);
690    PRINTIF(target, primme_largest_abs);
691 
692    PRINTParamsIF(projection, projection, primme_proj_default);
693    PRINTParamsIF(projection, projection, primme_proj_RR);
694    PRINTParamsIF(projection, projection, primme_proj_harmonic);
695    PRINTParamsIF(projection, projection, primme_proj_refined);
696 
697    PRINTIF(initBasisMode, primme_init_default);
698    PRINTIF(initBasisMode, primme_init_krylov);
699    PRINTIF(initBasisMode, primme_init_random);
700    PRINTIF(initBasisMode, primme_init_user);
701 
702    PRINT(numTargetShifts, %d);
703    if (primme.numTargetShifts > 0 && primme.targetShifts) {
704       fprintf(outputFile, "%s.targetShifts =", prefix);
705       for (i=0; i<primme.numTargetShifts;i++) {
706          fprintf(outputFile, " %e",primme.targetShifts[i]);
707       }
708       fprintf(outputFile, "\n");
709    }
710 
711    PRINT(dynamicMethodSwitch, %d);
712    PRINT(locking, %d);
713    PRINT(initSize, %d);
714    PRINT(numOrthoConst, %d);
715    PRINT_PRIMME_INT(ldevecs);
716    PRINT_PRIMME_INT(ldOPs);
717    fprintf(outputFile, "%s.iseed =", prefix);
718    for (i=0; i<4;i++) {
719       fprintf(outputFile, " %" PRIMME_INT_P, primme.iseed[i]);
720    }
721    fprintf(outputFile, "\n");
722    PRINTIF(orth, primme_orth_implicit_I);
723    PRINTIF(orth, primme_orth_explicit_I);
724 
725    PRINTIF(internalPrecision, primme_op_half);
726    PRINTIF(internalPrecision, primme_op_float);
727    PRINTIF(internalPrecision, primme_op_double);
728    PRINTIF(internalPrecision, primme_op_quad);
729 
730    PRINTParams(restarting, maxPrevRetain, %d);
731 
732    fprintf(outputFile, "\n// Correction parameters\n");
733    PRINTParams(correction, precondition, %d);
734    PRINTParams(correction, robustShifts, %d);
735    PRINTParams(correction, maxInnerIterations, %d);
736    PRINTParams(correction, relTolBase, %g);
737 
738    PRINTParamsIF(correction, convTest, primme_full_LTolerance);
739    PRINTParamsIF(correction, convTest, primme_decreasing_LTolerance);
740    PRINTParamsIF(correction, convTest, primme_adaptive_ETolerance);
741    PRINTParamsIF(correction, convTest, primme_adaptive);
742 
743    fprintf(outputFile, "\n// projectors for JD cor.eq.\n");
744    PRINTParams(correction, projectors.LeftQ , %d);
745    PRINTParams(correction, projectors.LeftX , %d);
746    PRINTParams(correction, projectors.RightQ, %d);
747    PRINTParams(correction, projectors.SkewQ , %d);
748    PRINTParams(correction, projectors.RightX, %d);
749    PRINTParams(correction, projectors.SkewX , %d);
750    fprintf(outputFile, "// ---------------------------------------------------\n");
751 
752 #undef PRINT
753 #undef PRINTIF
754 #undef PRINTParams
755 #undef PRINTParamsIF
756 }
757 
758 /*******************************************************************************
759  * Subroutine primme_get_member - get the value of a parameter in primme_params
760  *
761  * INPUT PARAMETERS
762  * ----------------
763  * primme  Structure containing various solver parameters and statistics
764  * label   reference to the parameter
765  *
766  * OUTPUT PARAMETERS
767  * -----------------
768  * value   value of the parameter
769  *
770  * RETURN
771  * ------
772  * error code  zero if ok
773  *
774  ******************************************************************************/
775 
primme_get_member(primme_params * primme,primme_params_label label,void * value)776 int primme_get_member(primme_params *primme, primme_params_label label,
777       void *value) {
778 
779    int i;
780    value_t *v = (value_t*)value;
781 
782    switch (label) {
783       case PRIMME_n:
784               *(PRIMME_INT*)value = primme->n;
785       break;
786       case PRIMME_matrixMatvec:
787               v->matFunc_v = primme->matrixMatvec;
788       break;
789       case PRIMME_matrixMatvec_type:
790               *(PRIMME_INT*)value = primme->matrixMatvec_type;
791       break;
792       case PRIMME_massMatrixMatvec:
793               v->matFunc_v = primme->massMatrixMatvec;
794       break;
795       case PRIMME_massMatrixMatvec_type:
796               *(PRIMME_INT*)value = primme->massMatrixMatvec_type;
797       break;
798       case PRIMME_applyPreconditioner:
799               v->matFunc_v = primme->applyPreconditioner;
800       break;
801       case PRIMME_applyPreconditioner_type:
802               *(PRIMME_INT*)value = primme->applyPreconditioner_type;
803       break;
804       case PRIMME_numProcs:
805               *(PRIMME_INT*)value = primme->numProcs;
806       break;
807       case PRIMME_procID:
808               *(PRIMME_INT*)value = primme->procID;
809       break;
810       case PRIMME_commInfo:
811               *(ptr_v*)value = primme->commInfo;
812       break;
813       case PRIMME_nLocal:
814               *(PRIMME_INT*)value = primme->nLocal;
815       break;
816       case PRIMME_globalSumReal:
817               v->globalSumRealFunc_v = primme->globalSumReal;
818       break;
819       case PRIMME_broadcastReal:
820               v->broadcastRealFunc_v = primme->broadcastReal;
821       break;
822       case PRIMME_numEvals:
823               *(PRIMME_INT*)value = primme->numEvals;
824       break;
825       case PRIMME_target:
826               *(PRIMME_INT*)value = primme->target;
827       break;
828       case PRIMME_numTargetShifts:
829               *(PRIMME_INT*)value = primme->numTargetShifts;
830       break;
831       case PRIMME_targetShifts:
832               *(ptr_v*)value = primme->targetShifts;
833       break;
834       case PRIMME_ShiftsForPreconditioner:
835               *(ptr_v*)value = primme->ShiftsForPreconditioner;
836       break;
837       case PRIMME_locking:
838               *(PRIMME_INT*)value = primme->locking;
839       break;
840       case PRIMME_initSize:
841               *(PRIMME_INT*)value = primme->initSize;
842       break;
843       case PRIMME_numOrthoConst:
844               *(PRIMME_INT*)value = primme->numOrthoConst;
845       break;
846       case PRIMME_dynamicMethodSwitch:
847               *(PRIMME_INT*)value = primme->dynamicMethodSwitch;
848       break;
849       case PRIMME_maxBasisSize:
850               *(PRIMME_INT*)value = primme->maxBasisSize;
851       break;
852       case PRIMME_minRestartSize:
853               *(PRIMME_INT*)value = primme->minRestartSize;
854       break;
855       case PRIMME_maxBlockSize:
856               *(PRIMME_INT*)value = primme->maxBlockSize;
857       break;
858       case PRIMME_maxMatvecs:
859               *(PRIMME_INT*)value = primme->maxMatvecs;
860       break;
861       case PRIMME_maxOuterIterations:
862               *(PRIMME_INT*)value = primme->maxOuterIterations;
863       break;
864       case PRIMME_iseed:
865          for (i=0; i< 4; i++) {
866             ((PRIMME_INT*)value)[i] = primme->iseed[i];
867          }
868       break;
869       case PRIMME_aNorm:
870               *(double*)value = primme->aNorm;
871       break;
872       case PRIMME_BNorm:
873               *(double*)value = primme->BNorm;
874       break;
875       case PRIMME_invBNorm:
876               *(double*)value = primme->invBNorm;
877       break;
878       case PRIMME_eps:
879               *(double*)value = primme->eps;
880       break;
881       case PRIMME_orth:
882               *(PRIMME_INT*)value = primme->orth;
883       break;
884       case PRIMME_internalPrecision:
885               *(PRIMME_INT*)value = primme->internalPrecision;
886       break;
887       case PRIMME_printLevel:
888               *(PRIMME_INT*)value = primme->printLevel;
889       break;
890       case PRIMME_outputFile:
891               *(FILE**)value = primme->outputFile;
892       break;
893       case PRIMME_matrix:
894               *(ptr_v*)value = primme->matrix;
895       break;
896       case PRIMME_massMatrix:
897               *(ptr_v*)value = primme->massMatrix;
898       break;
899       case PRIMME_preconditioner:
900               *(ptr_v*)value = primme->preconditioner;
901       break;
902       case PRIMME_initBasisMode:
903               *(PRIMME_INT*)value = primme->initBasisMode;
904       break;
905       case PRIMME_projectionParams_projection:
906               *(PRIMME_INT*)value = primme->projectionParams.projection;
907       break;
908       case PRIMME_restartingParams_maxPrevRetain:
909               *(PRIMME_INT*)value = primme->restartingParams.maxPrevRetain;
910       break;
911       case PRIMME_correctionParams_precondition:
912               *(PRIMME_INT*)value = primme->correctionParams.precondition;
913       break;
914       case PRIMME_correctionParams_robustShifts:
915               *(PRIMME_INT*)value = primme->correctionParams.robustShifts;
916       break;
917       case PRIMME_correctionParams_maxInnerIterations:
918               *(PRIMME_INT*)value = primme->correctionParams.maxInnerIterations;
919       break;
920       case PRIMME_correctionParams_projectors_LeftQ:
921               *(PRIMME_INT*)value = primme->correctionParams.projectors.LeftQ;
922       break;
923       case PRIMME_correctionParams_projectors_LeftX:
924               *(PRIMME_INT*)value = primme->correctionParams.projectors.LeftX;
925       break;
926       case PRIMME_correctionParams_projectors_RightQ:
927               *(PRIMME_INT*)value = primme->correctionParams.projectors.RightQ;
928       break;
929       case PRIMME_correctionParams_projectors_RightX:
930               *(PRIMME_INT*)value = primme->correctionParams.projectors.RightX;
931       break;
932       case PRIMME_correctionParams_projectors_SkewQ:
933               *(PRIMME_INT*)value = primme->correctionParams.projectors.SkewQ;
934       break;
935       case PRIMME_correctionParams_projectors_SkewX:
936               *(PRIMME_INT*)value = primme->correctionParams.projectors.SkewX;
937       break;
938       case PRIMME_correctionParams_convTest:
939               *(PRIMME_INT*)value = primme->correctionParams.convTest;
940       break;
941       case PRIMME_correctionParams_relTolBase:
942               *(double*)value = primme->correctionParams.relTolBase;
943       break;
944       case PRIMME_stats_numOuterIterations:
945               *(PRIMME_INT*)value = primme->stats.numOuterIterations;
946       break;
947       case PRIMME_stats_numRestarts:
948               *(PRIMME_INT*)value = primme->stats.numRestarts;
949       break;
950       case PRIMME_stats_numMatvecs:
951               *(PRIMME_INT*)value = primme->stats.numMatvecs;
952       break;
953       case PRIMME_stats_numPreconds:
954               *(PRIMME_INT*)value = primme->stats.numPreconds;
955       break;
956       case PRIMME_stats_numGlobalSum:
957               *(PRIMME_INT*)value = primme->stats.numGlobalSum;
958       break;
959       case PRIMME_stats_numBroadcast:
960               *(PRIMME_INT*)value = primme->stats.numBroadcast;
961       break;
962       case PRIMME_stats_volumeGlobalSum:
963               *(PRIMME_INT*)value = primme->stats.volumeGlobalSum;
964       break;
965       case PRIMME_stats_volumeBroadcast:
966               *(PRIMME_INT*)value = primme->stats.volumeBroadcast;
967       break;
968       case PRIMME_stats_flopsDense:
969               *(double*)value = primme->stats.flopsDense;
970       break;
971       case PRIMME_stats_numOrthoInnerProds:
972               *(double*)value = primme->stats.numOrthoInnerProds;
973       break;
974       case PRIMME_stats_elapsedTime:
975               *(double*)value = primme->stats.elapsedTime;
976       break;
977       case PRIMME_stats_timeMatvec:
978               *(double*)value = primme->stats.timeMatvec;
979       break;
980       case PRIMME_stats_timePrecond:
981               *(double*)value = primme->stats.timePrecond;
982       break;
983       case PRIMME_stats_timeOrtho:
984               *(double*)value = primme->stats.timeOrtho;
985       break;
986       case PRIMME_stats_timeGlobalSum:
987               *(double*)value = primme->stats.timeGlobalSum;
988       break;
989       case PRIMME_stats_timeBroadcast:
990               *(double*)value = primme->stats.timeBroadcast;
991       break;
992       case PRIMME_stats_timeDense:
993               *(double*)value = primme->stats.timeDense;
994       break;
995       case PRIMME_stats_estimateMinEVal:
996               *(double*)value = primme->stats.estimateMinEVal;
997       break;
998       case PRIMME_stats_estimateMaxEVal:
999               *(double*)value = primme->stats.estimateMaxEVal;
1000       break;
1001       case PRIMME_stats_estimateLargestSVal:
1002               *(double*)value = primme->stats.estimateLargestSVal;
1003       break;
1004       case PRIMME_stats_estimateBNorm:
1005               *(double*)value = primme->stats.estimateBNorm;
1006       break;
1007       case PRIMME_stats_estimateInvBNorm:
1008               *(double*)value = primme->stats.estimateInvBNorm;
1009       break;
1010       case PRIMME_stats_maxConvTol:
1011               *(double*)value = primme->stats.maxConvTol;
1012       break;
1013       case PRIMME_stats_lockingIssue:
1014               *(PRIMME_INT*)value = primme->stats.lockingIssue;
1015       break;
1016       case PRIMME_ldevecs:
1017               *(PRIMME_INT*)value = primme->ldevecs;
1018       break;
1019       case PRIMME_ldOPs:
1020               *(PRIMME_INT*)value = primme->ldOPs;
1021       break;
1022       case PRIMME_convTestFun:
1023               v->convTestFun_v = primme->convTestFun;
1024       break;
1025       case PRIMME_convTestFun_type:
1026               *(PRIMME_INT*)value = primme->convTestFun_type;
1027       break;
1028       case PRIMME_convtest:
1029               *(ptr_v*)value = primme->convtest;
1030       break;
1031       case PRIMME_monitorFun:
1032               v->monitorFun_v = primme->monitorFun;
1033       break;
1034       case PRIMME_monitorFun_type:
1035               *(PRIMME_INT*)value = primme->monitorFun_type;
1036       break;
1037       case PRIMME_monitor:
1038               *(ptr_v*)value = primme->monitor;
1039       break;
1040       case PRIMME_queue:
1041               *(ptr_v*)value = primme->queue;
1042       break;
1043       case PRIMME_profile:
1044               *(str_v*)value = primme->profile;
1045       break;
1046       default :
1047       return 1;
1048    }
1049    return 0;
1050 }
1051 
1052 /*******************************************************************************
1053  * Subroutine primme_set_member - set the value to a parameter in primme_params
1054  *
1055  * INPUT PARAMETERS
1056  * ----------------
1057  * label   reference to the parameter
1058  * value   value of the parameter
1059  *
1060  * INPUT/OUTPUT PARAMETERS
1061  * -----------------
1062  * primme  Structure containing various solver parameters and statistics
1063  *
1064  * RETURN
1065  * ------
1066  * error code  zero if ok
1067  *
1068  ******************************************************************************/
1069 
primme_set_member(primme_params * primme,primme_params_label label,void * value)1070 int primme_set_member(primme_params *primme, primme_params_label label,
1071       void *value) {
1072    int i;
1073 
1074    // Workaround to avoid warnings assigning void* pointers to function pointers
1075    value_t v = *(value_t*)&value;
1076    switch (label) {
1077       case PRIMME_n:
1078               primme->n = *(PRIMME_INT*)value;
1079       break;
1080       case PRIMME_matrixMatvec:
1081               primme->matrixMatvec = v.matFunc_v;
1082       break;
1083       case PRIMME_matrixMatvec_type:
1084               primme->matrixMatvec_type = (primme_op_datatype)*(PRIMME_INT*)value;
1085       break;
1086       case PRIMME_massMatrixMatvec:
1087               primme->massMatrixMatvec = v.matFunc_v;
1088       break;
1089       case PRIMME_massMatrixMatvec_type:
1090               primme->massMatrixMatvec_type = (primme_op_datatype)*(PRIMME_INT*)value;
1091       break;
1092       case PRIMME_applyPreconditioner:
1093               primme->applyPreconditioner = v.matFunc_v;
1094       break;
1095       case PRIMME_applyPreconditioner_type:
1096               primme->applyPreconditioner_type = (primme_op_datatype)*(PRIMME_INT*)value;
1097       break;
1098       case PRIMME_numProcs:
1099               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1100               primme->numProcs = (int)*(PRIMME_INT*)value;
1101       break;
1102       case PRIMME_procID:
1103               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1104               primme->procID = (int)*(PRIMME_INT*)value;
1105       break;
1106       case PRIMME_commInfo:
1107               primme->commInfo = (ptr_v)value;
1108       break;
1109       case PRIMME_nLocal:
1110               primme->nLocal = *(PRIMME_INT*)value;
1111       break;
1112       case PRIMME_globalSumReal:
1113               primme->globalSumReal = v.globalSumRealFunc_v;
1114       break;
1115       case PRIMME_globalSumReal_type:
1116               primme->globalSumReal_type = (primme_op_datatype)*(PRIMME_INT*)value;
1117       break;
1118       case PRIMME_broadcastReal:
1119               primme->broadcastReal = v.broadcastRealFunc_v;
1120       break;
1121       case PRIMME_broadcastReal_type:
1122               primme->broadcastReal_type = (primme_op_datatype)*(PRIMME_INT*)value;
1123       break;
1124       case PRIMME_numEvals:
1125               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1126               primme->numEvals = (int)*(PRIMME_INT*)value;
1127       break;
1128       case PRIMME_target:
1129               primme->target = (primme_target)*(PRIMME_INT*)value;
1130       break;
1131       case PRIMME_numTargetShifts:
1132               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1133               primme->numTargetShifts = (int)*(PRIMME_INT*)value;
1134       break;
1135       case PRIMME_targetShifts:
1136               primme->targetShifts = (double*)value;
1137       break;
1138       case PRIMME_ShiftsForPreconditioner:
1139               primme->ShiftsForPreconditioner = (double*)value;
1140       break;
1141       case PRIMME_locking:
1142               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1143               primme->locking = (int)*(PRIMME_INT*)value;
1144       break;
1145       case PRIMME_initSize:
1146               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1147               primme->initSize = (int)*(PRIMME_INT*)value;
1148       break;
1149       case PRIMME_numOrthoConst:
1150               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1151               primme->numOrthoConst = (int)*(PRIMME_INT*)value;
1152       break;
1153       case PRIMME_dynamicMethodSwitch:
1154               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1155               primme->dynamicMethodSwitch = (int)*(PRIMME_INT*)value;
1156       break;
1157       case PRIMME_maxBasisSize:
1158               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1159               primme->maxBasisSize = (int)*(PRIMME_INT*)value;
1160       break;
1161       case PRIMME_minRestartSize:
1162               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1163               primme->minRestartSize = (int)*(PRIMME_INT*)value;
1164       break;
1165       case PRIMME_maxBlockSize:
1166               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1167               primme->maxBlockSize = (int)*(PRIMME_INT*)value;
1168       break;
1169       case PRIMME_maxMatvecs:
1170               primme->maxMatvecs = *(PRIMME_INT*)value;
1171       break;
1172       case PRIMME_maxOuterIterations:
1173               primme->maxOuterIterations = *(PRIMME_INT*)value;
1174       break;
1175       case PRIMME_iseed:
1176          for (i=0; i< 4; i++) {
1177             primme->iseed[i] = ((PRIMME_INT*)value)[i];
1178          }
1179       break;
1180       case PRIMME_aNorm:
1181               primme->aNorm = *(double*)value;
1182       break;
1183       case PRIMME_BNorm:
1184               primme->BNorm = *(double*)value;
1185       break;
1186       case PRIMME_invBNorm:
1187               primme->invBNorm = *(double*)value;
1188       break;
1189       case PRIMME_eps:
1190               primme->eps = *(double*)value;
1191       break;
1192       case PRIMME_orth:
1193               primme->orth = (primme_orth)*(PRIMME_INT*)value;
1194       break;
1195       case PRIMME_internalPrecision:
1196               primme->internalPrecision = (primme_op_datatype)*(PRIMME_INT*)value;
1197       break;
1198       case PRIMME_printLevel:
1199               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1200               primme->printLevel = (int)*(PRIMME_INT*)value;
1201       break;
1202       case PRIMME_outputFile:
1203               primme->outputFile = (FILE*)value;
1204       break;
1205       case PRIMME_matrix:
1206               primme->matrix = (ptr_v)value;
1207       break;
1208       case PRIMME_massMatrix:
1209               primme->massMatrix = (ptr_v)value;
1210       break;
1211       case PRIMME_preconditioner:
1212               primme->preconditioner = (ptr_v)value;
1213       break;
1214       case PRIMME_initBasisMode:
1215               primme->initBasisMode = (primme_init)*(PRIMME_INT*)value;
1216       break;
1217       case PRIMME_projectionParams_projection:
1218               primme->projectionParams.projection = (primme_projection)*(PRIMME_INT*)value;
1219       break;
1220       case PRIMME_restartingParams_maxPrevRetain:
1221               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1222               primme->restartingParams.maxPrevRetain = (int)*(PRIMME_INT*)value;
1223       break;
1224       case PRIMME_correctionParams_precondition:
1225               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1226               primme->correctionParams.precondition = (int)*(PRIMME_INT*)value;
1227       break;
1228       case PRIMME_correctionParams_robustShifts:
1229               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1230               primme->correctionParams.robustShifts = (int)*(PRIMME_INT*)value;
1231       break;
1232       case PRIMME_correctionParams_maxInnerIterations:
1233               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1234               primme->correctionParams.maxInnerIterations = (int)*(PRIMME_INT*)value;
1235       break;
1236       case PRIMME_correctionParams_projectors_LeftQ:
1237               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1238               primme->correctionParams.projectors.LeftQ = (int)*(PRIMME_INT*)value;
1239       break;
1240       case PRIMME_correctionParams_projectors_LeftX:
1241               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1242               primme->correctionParams.projectors.LeftX = (int)*(PRIMME_INT*)value;
1243       break;
1244       case PRIMME_correctionParams_projectors_RightQ:
1245               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1246               primme->correctionParams.projectors.RightQ = (int)*(PRIMME_INT*)value;
1247       break;
1248       case PRIMME_correctionParams_projectors_RightX:
1249               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1250               primme->correctionParams.projectors.RightX = (int)*(PRIMME_INT*)value;
1251       break;
1252       case PRIMME_correctionParams_projectors_SkewQ:
1253               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1254               primme->correctionParams.projectors.SkewQ = (int)*(PRIMME_INT*)value;
1255       break;
1256       case PRIMME_correctionParams_projectors_SkewX:
1257               if (*(PRIMME_INT*)value > INT_MAX) return 1; else
1258               primme->correctionParams.projectors.SkewX = (int)*(PRIMME_INT*)value;
1259       break;
1260       case PRIMME_correctionParams_convTest:
1261               primme->correctionParams.convTest = (primme_convergencetest)*(PRIMME_INT*)value;
1262       break;
1263       case PRIMME_correctionParams_relTolBase:
1264               primme->correctionParams.relTolBase = *(double*)value;
1265       break;
1266       case PRIMME_stats_numOuterIterations:
1267               primme->stats.numOuterIterations = *(PRIMME_INT*)value;
1268       break;
1269       case PRIMME_stats_numRestarts:
1270               primme->stats.numRestarts = *(PRIMME_INT*)value;
1271       break;
1272       case PRIMME_stats_numMatvecs:
1273               primme->stats.numMatvecs = *(PRIMME_INT*)value;
1274       break;
1275       case PRIMME_stats_numPreconds:
1276               primme->stats.numPreconds = *(PRIMME_INT*)value;
1277       break;
1278       case PRIMME_stats_volumeGlobalSum:
1279               primme->stats.volumeGlobalSum = *(PRIMME_INT*)value;
1280       break;
1281       case PRIMME_stats_volumeBroadcast:
1282               primme->stats.volumeBroadcast = *(PRIMME_INT*)value;
1283       break;
1284       case PRIMME_stats_flopsDense:
1285               primme->stats.flopsDense = *(double*)value;
1286       break;
1287       case PRIMME_stats_numOrthoInnerProds:
1288               primme->stats.numOrthoInnerProds = *(double*)value;
1289       break;
1290       case PRIMME_stats_elapsedTime:
1291               primme->stats.elapsedTime = *(double*)value;
1292       break;
1293       case PRIMME_stats_timeMatvec:
1294               primme->stats.timeMatvec = *(double*)value;
1295       break;
1296       case PRIMME_stats_timePrecond:
1297               primme->stats.timePrecond = *(double*)value;
1298       break;
1299       case PRIMME_stats_timeOrtho:
1300               primme->stats.timeOrtho = *(double*)value;
1301       break;
1302       case PRIMME_stats_timeGlobalSum:
1303               primme->stats.timeGlobalSum = *(double*)value;
1304       break;
1305       case PRIMME_stats_timeBroadcast:
1306               primme->stats.timeBroadcast = *(double*)value;
1307       break;
1308       case PRIMME_stats_timeDense:
1309               primme->stats.timeDense = *(double*)value;
1310       break;
1311       case PRIMME_stats_estimateMinEVal:
1312               primme->stats.estimateMinEVal = *(double*)value;
1313       break;
1314       case PRIMME_stats_estimateMaxEVal:
1315               primme->stats.estimateMaxEVal = *(double*)value;
1316       break;
1317       case PRIMME_stats_estimateLargestSVal:
1318               primme->stats.estimateLargestSVal = *(double*)value;
1319       break;
1320       case PRIMME_stats_estimateBNorm:
1321               primme->stats.estimateBNorm = *(double*)value;
1322       break;
1323       case PRIMME_stats_estimateInvBNorm:
1324               primme->stats.estimateInvBNorm = *(double*)value;
1325       break;
1326       case PRIMME_stats_maxConvTol:
1327               primme->stats.maxConvTol = *(double*)value;
1328       break;
1329       case PRIMME_stats_lockingIssue:
1330               primme->stats.lockingIssue = *(PRIMME_INT*)value;
1331       break;
1332       case PRIMME_convTestFun:
1333               primme->convTestFun = v.convTestFun_v;
1334       break;
1335       case PRIMME_convTestFun_type:
1336               primme->convTestFun_type = (primme_op_datatype)*(PRIMME_INT*)value;
1337       break;
1338       case PRIMME_convtest:
1339               primme->convtest = (ptr_v)value;
1340       break;
1341       case PRIMME_ldevecs:
1342               primme->ldevecs = *(PRIMME_INT*)value;
1343       break;
1344       case PRIMME_ldOPs:
1345               primme->ldOPs = *(PRIMME_INT*)value;
1346       break;
1347       case PRIMME_monitorFun:
1348               primme->monitorFun = v.monitorFun_v;
1349       break;
1350       case PRIMME_monitorFun_type:
1351               primme->monitorFun_type = (primme_op_datatype)*(PRIMME_INT*)value;
1352       break;
1353       case PRIMME_monitor:
1354               primme->monitor = (ptr_v)value;
1355       break;
1356       case PRIMME_queue:
1357               primme->queue = (ptr_v)value;
1358       break;
1359       case PRIMME_profile:
1360               primme->profile = (str_v)value;
1361       break;
1362       default :
1363       return 1;
1364    }
1365    return 0;
1366 }
1367 
1368 /*******************************************************************************
1369  * Subroutine primme_member_info - return the label value or the label name,
1370  *    the type and the arity of some member of primme_params.
1371  *
1372  * INPUT/OUTPUT PARAMETERS
1373  * ----------------
1374  * label       reference to the parameter by the value in primme_params_label
1375  * label_name  reference by the string associated to the parameter
1376  *
1377  * OUTPUT PARAMETERS
1378  * -----------------
1379  * type        kind of value
1380  * arity       number of elements
1381  *
1382  * RETURN
1383  * ------
1384  * error code  zero if ok
1385  *
1386  ******************************************************************************/
1387 
primme_member_info(primme_params_label * label_,const char ** label_name_,primme_type * type,int * arity)1388 int primme_member_info(primme_params_label *label_, const char** label_name_,
1389       primme_type *type, int *arity) {
1390    primme_params_label label = (primme_params_label)1000;
1391    const char *label_name;
1392 
1393    /* Quick exit when neither label nor label_name is given */
1394 
1395    if (label_ == NULL && (label_name_ == NULL || *label_name_ == NULL)) {
1396       return 1;
1397    }
1398 
1399    /* Get the label from label_name_ and label_name from label_ */
1400 
1401 #define IF_IS(F,O) \
1402    if ((label_name_ && *label_name_ && strcmp(#F, *label_name_) == 0) \
1403          || (label_ && *label_ == PRIMME_ ## O)) { \
1404       label = PRIMME_ ## O; \
1405       label_name = #F; \
1406    }
1407 
1408    IF_IS(n                            , n);
1409    IF_IS(matrixMatvec                 , matrixMatvec);
1410    IF_IS(matrixMatvec_type            , matrixMatvec_type);
1411    IF_IS(massMatrixMatvec             , massMatrixMatvec);
1412    IF_IS(massMatrixMatvec_type        , massMatrixMatvec_type);
1413    IF_IS(applyPreconditioner          , applyPreconditioner);
1414    IF_IS(applyPreconditioner_type     , applyPreconditioner_type);
1415    IF_IS(numProcs                     , numProcs);
1416    IF_IS(procID                       , procID);
1417    IF_IS(commInfo                     , commInfo);
1418    IF_IS(nLocal                       , nLocal);
1419    IF_IS(globalSumReal                , globalSumReal);
1420    IF_IS(broadcastReal                , broadcastReal);
1421    IF_IS(numEvals                     , numEvals);
1422    IF_IS(target                       , target);
1423    IF_IS(numTargetShifts              , numTargetShifts);
1424    IF_IS(targetShifts                 , targetShifts);
1425    IF_IS(locking                      , locking);
1426    IF_IS(initSize                     , initSize);
1427    IF_IS(numOrthoConst                , numOrthoConst);
1428    IF_IS(dynamicMethodSwitch          , dynamicMethodSwitch);
1429    IF_IS(maxBasisSize                 , maxBasisSize);
1430    IF_IS(minRestartSize               , minRestartSize);
1431    IF_IS(maxBlockSize                 , maxBlockSize);
1432    IF_IS(maxMatvecs                   , maxMatvecs);
1433    IF_IS(maxOuterIterations           , maxOuterIterations);
1434    IF_IS(iseed                        , iseed);
1435    IF_IS(aNorm                        , aNorm);
1436    IF_IS(BNorm                        , BNorm);
1437    IF_IS(invBNorm                     , invBNorm);
1438    IF_IS(eps                          , eps);
1439    IF_IS(orth                         , orth);
1440    IF_IS(internalPrecision            , internalPrecision);
1441    IF_IS(printLevel                   , printLevel);
1442    IF_IS(outputFile                   , outputFile);
1443    IF_IS(matrix                       , matrix);
1444    IF_IS(massMatrix                   , massMatrix);
1445    IF_IS(preconditioner               , preconditioner);
1446    IF_IS(ShiftsForPreconditioner      , ShiftsForPreconditioner);
1447    IF_IS(initBasisMode                , initBasisMode);
1448    IF_IS(projection_projection        , projectionParams_projection);
1449    IF_IS(restarting_maxPrevRetain     , restartingParams_maxPrevRetain);
1450    IF_IS(correction_precondition      , correctionParams_precondition);
1451    IF_IS(correction_robustShifts      , correctionParams_robustShifts);
1452    IF_IS(correction_maxInnerIterations, correctionParams_maxInnerIterations);
1453    IF_IS(correction_projectors_LeftQ  , correctionParams_projectors_LeftQ);
1454    IF_IS(correction_projectors_LeftX  , correctionParams_projectors_LeftX);
1455    IF_IS(correction_projectors_RightQ , correctionParams_projectors_RightQ);
1456    IF_IS(correction_projectors_RightX , correctionParams_projectors_RightX);
1457    IF_IS(correction_projectors_SkewQ  , correctionParams_projectors_SkewQ);
1458    IF_IS(correction_projectors_SkewX  , correctionParams_projectors_SkewX);
1459    IF_IS(correction_convTest          , correctionParams_convTest);
1460    IF_IS(correction_relTolBase        , correctionParams_relTolBase);
1461    IF_IS(stats_numOuterIterations     , stats_numOuterIterations);
1462    IF_IS(stats_numRestarts            , stats_numRestarts);
1463    IF_IS(stats_numMatvecs             , stats_numMatvecs);
1464    IF_IS(stats_numPreconds            , stats_numPreconds);
1465    IF_IS(stats_numGlobalSum           , stats_numGlobalSum);
1466    IF_IS(stats_volumeGlobalSum        , stats_volumeGlobalSum);
1467    IF_IS(stats_numBroadcast           , stats_numBroadcast);
1468    IF_IS(stats_volumeBroadcast        , stats_volumeBroadcast);
1469    IF_IS(stats_flopsDense             , stats_flopsDense);
1470    IF_IS(stats_numOrthoInnerProds     , stats_numOrthoInnerProds);
1471    IF_IS(stats_elapsedTime            , stats_elapsedTime);
1472    IF_IS(stats_timeMatvec             , stats_timeMatvec);
1473    IF_IS(stats_timePrecond            , stats_timePrecond);
1474    IF_IS(stats_timeOrtho              , stats_timeOrtho);
1475    IF_IS(stats_timeGlobalSum          , stats_timeGlobalSum);
1476    IF_IS(stats_timeBroadcast          , stats_timeBroadcast);
1477    IF_IS(stats_timeDense              , stats_timeDense);
1478    IF_IS(stats_estimateMinEVal        , stats_estimateMinEVal);
1479    IF_IS(stats_estimateMaxEVal        , stats_estimateMaxEVal);
1480    IF_IS(stats_estimateLargestSVal    , stats_estimateLargestSVal);
1481    IF_IS(stats_estimateBNorm          , stats_estimateBNorm);
1482    IF_IS(stats_estimateInvBNorm       , stats_estimateInvBNorm);
1483    IF_IS(stats_maxConvTol             , stats_maxConvTol);
1484    IF_IS(stats_lockingIssue           , stats_lockingIssue);
1485    IF_IS(convTestFun                  , convTestFun);
1486    IF_IS(convTestFun_type             , convTestFun_type);
1487    IF_IS(convtest                     , convtest);
1488    IF_IS(ldevecs                      , ldevecs);
1489    IF_IS(ldOPs                        , ldOPs);
1490    IF_IS(monitorFun                   , monitorFun);
1491    IF_IS(monitorFun_type              , monitorFun_type);
1492    IF_IS(monitor                      , monitor);
1493    IF_IS(queue                        , queue);
1494    IF_IS(profile                      , profile);
1495 #undef IF_IS
1496 
1497    /* Return label/label_name */
1498 
1499    if (label_) *label_ = label;
1500    if (label_name_) *label_name_ = label_name;
1501 
1502    /* Return type and arity */
1503 
1504    switch(label) {
1505       /* members with type int */
1506 
1507       case PRIMME_matrixMatvec_type:
1508       case PRIMME_applyPreconditioner_type:
1509       case PRIMME_globalSumReal_type:
1510       case PRIMME_broadcastReal_type:
1511       case PRIMME_massMatrixMatvec_type:
1512       case PRIMME_n:
1513       case PRIMME_numEvals:
1514       case PRIMME_target:
1515       case PRIMME_locking:
1516       case PRIMME_initSize:
1517       case PRIMME_numOrthoConst:
1518       case PRIMME_dynamicMethodSwitch:
1519       case PRIMME_maxBasisSize:
1520       case PRIMME_minRestartSize:
1521       case PRIMME_maxBlockSize:
1522       case PRIMME_maxMatvecs:
1523       case PRIMME_maxOuterIterations:
1524       case PRIMME_initBasisMode:
1525       case PRIMME_orth:
1526       case PRIMME_internalPrecision:
1527       case PRIMME_projectionParams_projection:
1528       case PRIMME_restartingParams_maxPrevRetain:
1529       case PRIMME_correctionParams_precondition:
1530       case PRIMME_correctionParams_robustShifts:
1531       case PRIMME_correctionParams_maxInnerIterations:
1532       case PRIMME_correctionParams_projectors_LeftQ:
1533       case PRIMME_correctionParams_projectors_LeftX:
1534       case PRIMME_correctionParams_projectors_RightQ:
1535       case PRIMME_correctionParams_projectors_RightX:
1536       case PRIMME_correctionParams_projectors_SkewQ:
1537       case PRIMME_correctionParams_projectors_SkewX:
1538       case PRIMME_correctionParams_convTest:
1539       case PRIMME_stats_numOuterIterations:
1540       case PRIMME_stats_numRestarts:
1541       case PRIMME_stats_numMatvecs:
1542       case PRIMME_stats_numPreconds:
1543       case PRIMME_stats_numGlobalSum:
1544       case PRIMME_stats_volumeGlobalSum:
1545       case PRIMME_stats_numBroadcast:
1546       case PRIMME_stats_volumeBroadcast:
1547       case PRIMME_stats_lockingIssue:
1548       case PRIMME_numProcs:
1549       case PRIMME_procID:
1550       case PRIMME_nLocal:
1551       case PRIMME_numTargetShifts:
1552       case PRIMME_printLevel:
1553       case PRIMME_ldevecs:
1554       case PRIMME_ldOPs:
1555       case PRIMME_monitorFun_type:
1556       case PRIMME_convTestFun_type:
1557       if (type) *type = primme_int;
1558       if (arity) *arity = 1;
1559       break;
1560 
1561       case PRIMME_iseed:
1562       if (type) *type = primme_int;
1563       if (arity) *arity = 4;
1564       break;
1565 
1566       /* members with type double */
1567 
1568       case PRIMME_aNorm:
1569       case PRIMME_BNorm:
1570       case PRIMME_invBNorm:
1571       case PRIMME_eps:
1572       case PRIMME_correctionParams_relTolBase:
1573       case PRIMME_stats_numOrthoInnerProds:
1574       case PRIMME_stats_timeMatvec:
1575       case PRIMME_stats_timePrecond:
1576       case PRIMME_stats_timeOrtho:
1577       case PRIMME_stats_timeGlobalSum:
1578       case PRIMME_stats_timeBroadcast:
1579       case PRIMME_stats_flopsDense:
1580       case PRIMME_stats_timeDense:
1581       case PRIMME_stats_elapsedTime:
1582       case PRIMME_stats_estimateMinEVal:
1583       case PRIMME_stats_estimateMaxEVal:
1584       case PRIMME_stats_estimateLargestSVal:
1585       case PRIMME_stats_estimateBNorm:
1586       case PRIMME_stats_estimateInvBNorm:
1587       case PRIMME_stats_maxConvTol:
1588       if (type) *type = primme_double;
1589       if (arity) *arity = 1;
1590       break;
1591 
1592       case PRIMME_targetShifts:
1593       case PRIMME_ShiftsForPreconditioner:
1594       if (type) *type = primme_double;
1595       if (arity) *arity = 0;
1596       break;
1597 
1598       /* members with type pointer */
1599 
1600       case PRIMME_matrixMatvec:
1601       case PRIMME_applyPreconditioner:
1602       case PRIMME_commInfo:
1603       case PRIMME_globalSumReal:
1604       case PRIMME_broadcastReal:
1605       case PRIMME_massMatrixMatvec:
1606       case PRIMME_outputFile:
1607       case PRIMME_matrix:
1608       case PRIMME_massMatrix:
1609       case PRIMME_preconditioner:
1610       case PRIMME_convTestFun:
1611       case PRIMME_convtest:
1612       case PRIMME_monitorFun:
1613       case PRIMME_monitor:
1614       case PRIMME_queue:
1615       if (type) *type = primme_pointer;
1616       if (arity) *arity = 1;
1617       break;
1618 
1619       case PRIMME_profile:
1620       if (type) *type = primme_string;
1621       if (arity) *arity = 1;
1622       break;
1623 
1624       default:
1625       return 1;
1626    }
1627 
1628    return 0;
1629 }
1630 
1631 /*******************************************************************************
1632  * Subroutine primme_constant_info - return the value of a primme enum constant.
1633  *
1634  * INPUT PARAMETERS
1635  * ----------------
1636  * label_name  name of the constant
1637  *
1638  * OUTPUT PARAMETERS
1639  * -----------------
1640  * value       numerical value of the constant
1641  *
1642  * RETURN
1643  * ------
1644  * error code  zero if ok
1645  *
1646  ******************************************************************************/
1647 
primme_constant_info(const char * label_name,int * value)1648 int primme_constant_info(const char* label_name, int *value) {
1649 
1650 #define IF_IS(F) if (strcmp(#F, label_name) == 0) {*value = (int)F; return 0;}
1651 
1652    /* primme_preset_method */
1653 
1654    IF_IS(PRIMME_DEFAULT_METHOD);
1655    IF_IS(PRIMME_DYNAMIC);
1656    IF_IS(PRIMME_DEFAULT_MIN_TIME);
1657    IF_IS(PRIMME_DEFAULT_MIN_MATVECS);
1658    IF_IS(PRIMME_Arnoldi);
1659    IF_IS(PRIMME_GD);
1660    IF_IS(PRIMME_GD_plusK);
1661    IF_IS(PRIMME_GD_Olsen_plusK);
1662    IF_IS(PRIMME_JD_Olsen_plusK);
1663    IF_IS(PRIMME_RQI);
1664    IF_IS(PRIMME_JDQR);
1665    IF_IS(PRIMME_JDQMR);
1666    IF_IS(PRIMME_JDQMR_ETol);
1667    IF_IS(PRIMME_STEEPEST_DESCENT);
1668    IF_IS(PRIMME_LOBPCG_OrthoBasis);
1669    IF_IS(PRIMME_LOBPCG_OrthoBasis_Window);
1670 
1671    /* enum members for targeting; restarting and innertest */
1672 
1673    IF_IS(primme_smallest);
1674    IF_IS(primme_largest);
1675    IF_IS(primme_closest_geq);
1676    IF_IS(primme_closest_leq);
1677    IF_IS(primme_closest_abs);
1678    IF_IS(primme_largest_abs);
1679    IF_IS(primme_proj_default);
1680    IF_IS(primme_proj_RR);
1681    IF_IS(primme_proj_harmonic);
1682    IF_IS(primme_proj_refined);
1683    IF_IS(primme_init_default);
1684    IF_IS(primme_init_krylov);
1685    IF_IS(primme_init_random);
1686    IF_IS(primme_init_user);
1687    IF_IS(primme_full_LTolerance);
1688    IF_IS(primme_decreasing_LTolerance);
1689    IF_IS(primme_adaptive_ETolerance);
1690    IF_IS(primme_adaptive);
1691 
1692    /* enum member for event */
1693 
1694    IF_IS(primme_event_outer_iteration);
1695    IF_IS(primme_event_inner_iteration);
1696    IF_IS(primme_event_restart);
1697    IF_IS(primme_event_reset);
1698    IF_IS(primme_event_converged);
1699    IF_IS(primme_event_locked);
1700    IF_IS(primme_event_message);
1701    IF_IS(primme_event_profile);
1702 
1703    /* enum member from orth */
1704 
1705    IF_IS(primme_orth_default);
1706    IF_IS(primme_orth_explicit_I);
1707    IF_IS(primme_orth_implicit_I);
1708 
1709    /* enum member from op_datatype */
1710    IF_IS(primme_op_default);
1711    IF_IS(primme_op_quad);
1712    IF_IS(primme_op_double);
1713    IF_IS(primme_op_float);
1714    IF_IS(primme_op_half);
1715    IF_IS(primme_op_int);
1716 #undef IF_IS
1717 
1718    /* return error if label not found */
1719 
1720    return 1;
1721 }
1722 
1723 /*******************************************************************************
1724  * Subroutine primme_enum_member_info - return the value of a string
1725  * representing an enum constant, or vice versa.
1726  *
1727  * INPUT/OUTPUT PARAMETERS
1728  * -----------------------
1729  * label       member to which the constant relates
1730  * value       (in) if *value >= 0, value to get the associated string,
1731  *             (in/out) if *value < 0, return the value associated to
1732  *             value_name.
1733  * value_name  (in) if *value_name > 0, string for which to seek the value
1734  *             (in/out) if *value_name == 0, return the associated to value.
1735  *
1736  * RETURN
1737  * ------
1738  * error code   0: OK
1739  *             -1: Invalid input
1740  *             -2: either value or value_name was not found
1741  *
1742  ******************************************************************************/
1743 
primme_enum_member_info(primme_params_label label,int * value,const char ** value_name)1744 int primme_enum_member_info(
1745       primme_params_label label, int *value, const char **value_name) {
1746 
1747    if (!value || !value_name || (*value >= 0 && *value_name) ||
1748          (*value < 0 && !*value_name)) {
1749       return -1;
1750    }
1751 
1752 #define IF_IS(F)                                                               \
1753    if (*value == (int)(F) || (*value_name && strcmp(#F, *value_name) == 0)) {  \
1754       *value = (int)(F);                                                       \
1755       *value_name = #F;                                                        \
1756       return 0;                                                                \
1757    }
1758 
1759    switch(label) {
1760    // Hack: Check method
1761    case PRIMME_commInfo:
1762    IF_IS(PRIMME_DEFAULT_METHOD);
1763    IF_IS(PRIMME_DYNAMIC);
1764    IF_IS(PRIMME_DEFAULT_MIN_TIME);
1765    IF_IS(PRIMME_DEFAULT_MIN_MATVECS);
1766    IF_IS(PRIMME_Arnoldi);
1767    IF_IS(PRIMME_GD);
1768    IF_IS(PRIMME_GD_plusK);
1769    IF_IS(PRIMME_GD_Olsen_plusK);
1770    IF_IS(PRIMME_JD_Olsen_plusK);
1771    IF_IS(PRIMME_RQI);
1772    IF_IS(PRIMME_JDQR);
1773    IF_IS(PRIMME_JDQMR);
1774    IF_IS(PRIMME_JDQMR_ETol);
1775    IF_IS(PRIMME_STEEPEST_DESCENT);
1776    IF_IS(PRIMME_LOBPCG_OrthoBasis);
1777    IF_IS(PRIMME_LOBPCG_OrthoBasis_Window);
1778    break;
1779 
1780    case PRIMME_target:
1781    IF_IS(primme_smallest);
1782    IF_IS(primme_largest);
1783    IF_IS(primme_closest_geq);
1784    IF_IS(primme_closest_leq);
1785    IF_IS(primme_closest_abs);
1786    IF_IS(primme_largest_abs);
1787    break;
1788 
1789    case PRIMME_projectionParams_projection:
1790    IF_IS(primme_proj_default);
1791    IF_IS(primme_proj_RR);
1792    IF_IS(primme_proj_harmonic);
1793    IF_IS(primme_proj_refined);
1794    break;
1795 
1796    case PRIMME_initBasisMode:
1797    IF_IS(primme_init_default);
1798    IF_IS(primme_init_krylov);
1799    IF_IS(primme_init_random);
1800    IF_IS(primme_init_user);
1801    break;
1802 
1803    case PRIMME_correctionParams_convTest:
1804    IF_IS(primme_full_LTolerance);
1805    IF_IS(primme_decreasing_LTolerance);
1806    IF_IS(primme_adaptive_ETolerance);
1807    IF_IS(primme_adaptive);
1808    break;
1809 
1810    case PRIMME_orth:
1811    IF_IS(primme_orth_default);
1812    IF_IS(primme_orth_explicit_I);
1813    IF_IS(primme_orth_implicit_I);
1814    break;
1815 
1816    case PRIMME_matrixMatvec_type:
1817    case PRIMME_applyPreconditioner_type:
1818    case PRIMME_globalSumReal_type:
1819    case PRIMME_broadcastReal_type:
1820    case PRIMME_massMatrixMatvec_type:
1821    IF_IS(primme_op_default);
1822    IF_IS(primme_op_quad);
1823    IF_IS(primme_op_double);
1824    IF_IS(primme_op_float);
1825    IF_IS(primme_op_half);
1826    IF_IS(primme_op_int);
1827    break;
1828 
1829    default: break;
1830    }
1831 #undef IF_IS
1832 
1833    /* return error if label not found */
1834 
1835    return -2;
1836 }
1837 
1838 #endif /* USE_DOUBLE */
1839