LCOV - code coverage report
Current view: top level - eps/impls/external/primme - primme.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 291 318 91.5 %
Date: 2024-12-18 00:42:09 Functions: 21 22 95.5 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       6             :    This file is part of SLEPc.
       7             :    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
       8             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       9             : */
      10             : /*
      11             :    This file implements a wrapper to the PRIMME package
      12             : */
      13             : 
      14             : #include <slepc/private/epsimpl.h>    /*I "slepceps.h" I*/
      15             : 
      16             : #include <primme.h>
      17             : 
      18             : #if defined(PETSC_USE_COMPLEX)
      19             : #if defined(PETSC_USE_REAL_SINGLE)
      20             : #define PRIMME_DRIVER cprimme
      21             : #else
      22             : #define PRIMME_DRIVER zprimme
      23             : #endif
      24             : #else
      25             : #if defined(PETSC_USE_REAL_SINGLE)
      26             : #define PRIMME_DRIVER sprimme
      27             : #else
      28             : #define PRIMME_DRIVER dprimme
      29             : #endif
      30             : #endif
      31             : 
      32             : #if defined(PRIMME_VERSION_MAJOR) && PRIMME_VERSION_MAJOR*100+PRIMME_VERSION_MINOR >= 202
      33             : #define SLEPC_HAVE_PRIMME2p2
      34             : #endif
      35             : 
      36             : typedef struct {
      37             :   primme_params        primme;    /* param struct */
      38             :   PetscInt             bs;        /* block size */
      39             :   primme_preset_method method;    /* primme method */
      40             :   Mat                  A,B;       /* problem matrices */
      41             :   KSP                  ksp;       /* linear solver and preconditioner */
      42             :   Vec                  x,y;       /* auxiliary vectors */
      43             :   double               target;    /* a copy of eps's target */
      44             : } EPS_PRIMME;
      45             : 
      46         852 : static void par_GlobalSumReal(void *sendBuf,void *recvBuf,int *count,primme_params *primme,int *ierr)
      47             : {
      48         852 :   if (sendBuf == recvBuf) {
      49        1704 :     *ierr = MPI_Allreduce(MPI_IN_PLACE,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
      50             :   } else {
      51           0 :     *ierr = MPI_Allreduce(sendBuf,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
      52             :   }
      53         852 : }
      54             : 
      55             : #if defined(SLEPC_HAVE_PRIMME3)
      56        2252 : static void par_broadcastReal(void *buf,int *count,primme_params *primme,int *ierr)
      57             : {
      58        2252 :   *ierr = MPI_Bcast(buf,*count,MPIU_REAL,0/*root*/,PetscObjectComm((PetscObject)primme->commInfo));
      59        2252 : }
      60             : #endif
      61             : 
      62             : #if defined(SLEPC_HAVE_PRIMME2p2)
      63        7413 : static void convTestFun(double *eval,void *evec,double *resNorm,int *isconv,primme_params *primme,int *err)
      64             : {
      65        7413 :   PetscErrorCode ierr;
      66        7413 :   EPS            eps = (EPS)primme->commInfo;
      67        7413 :   PetscScalar    eigvr = eval?*eval:0.0;
      68        7413 :   PetscReal      r = resNorm?*resNorm:HUGE_VAL,errest;
      69             : 
      70        7413 :   ierr = (*eps->converged)(eps,eigvr,0.0,r,&errest,eps->convergedctx);
      71        7413 :   if (ierr) *err = 1;
      72             :   else {
      73        7413 :     *isconv = (errest<=eps->tol?1:0);
      74        7413 :     *err = 0;
      75             :   }
      76        7413 : }
      77             : 
      78        6591 : static void monitorFun(void *basisEvals,int *basisSize,int *basisFlags,int *iblock,int *blockSize,void *basisNorms,int *numConverged,void *lockedEvals,int *numLocked,int *lockedFlags,void *lockedNorms,int *inner_its,void *LSRes,
      79             : #if defined(SLEPC_HAVE_PRIMME3)
      80             :                        const char *msg,double *time,
      81             : #endif
      82             :                        primme_event *event,struct primme_params *primme,int *err)
      83             : {
      84        6591 :   PetscErrorCode ierr = PETSC_SUCCESS;
      85        6591 :   EPS            eps = (EPS)primme->commInfo;
      86        6591 :   PetscInt       i,k,nerrest;
      87             : 
      88        6591 :   switch (*event) {
      89         795 :     case primme_event_outer_iteration:
      90             :       /* Update EPS */
      91         795 :       PetscCallVoid(PetscIntCast(primme->stats.numOuterIterations,&eps->its));
      92         795 :       eps->nconv = primme->initSize;
      93         795 :       k=0;
      94        1411 :       if (lockedEvals && numLocked) for (i=0; i<*numLocked && k<eps->ncv; i++) eps->eigr[k++] = ((PetscReal*)lockedEvals)[i];
      95         795 :       nerrest = k;
      96         795 :       if (iblock && blockSize) {
      97        2112 :         for (i=0; i<*blockSize && k+iblock[i]<eps->ncv; i++) eps->errest[k+iblock[i]] = ((PetscReal*)basisNorms)[i];
      98         795 :         nerrest = k+(*blockSize>0?1+iblock[*blockSize-1]:0);
      99             :       }
     100       16867 :       if (basisEvals && basisSize) for (i=0; i<*basisSize && k<eps->ncv; i++) eps->eigr[k++] = ((PetscReal*)basisEvals)[i];
     101             :       /* Show progress */
     102         795 :       ierr = EPSMonitor(eps,eps->its,numConverged?*numConverged:0,eps->eigr,eps->eigi,eps->errest,nerrest);
     103         795 :       break;
     104             : #if defined(SLEPC_HAVE_PRIMME3)
     105           0 :     case primme_event_message:
     106             :       /* Print PRIMME information messages */
     107           0 :       ierr = PetscInfo(eps,"%s\n",msg);
     108           0 :       break;
     109             : #endif
     110             :     default:
     111             :       break;
     112             :   }
     113        6591 :   *err = (ierr!=0)? 1: 0;
     114             : }
     115             : #endif /* SLEPC_HAVE_PRIMME2p2 */
     116             : 
     117        5940 : static void matrixMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
     118             : {
     119        5940 :   PetscInt   i;
     120        5940 :   EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
     121        5940 :   Vec        x = ops->x,y = ops->y;
     122        5940 :   Mat        A = ops->A;
     123             : 
     124        5940 :   PetscFunctionBegin;
     125       13752 :   for (i=0;i<*blockSize;i++) {
     126        7812 :     PetscCallAbort(PetscObjectComm((PetscObject)A),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
     127        7812 :     PetscCallAbort(PetscObjectComm((PetscObject)A),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
     128        7812 :     PetscCallAbort(PetscObjectComm((PetscObject)A),MatMult(A,x,y));
     129        7812 :     PetscCallAbort(PetscObjectComm((PetscObject)A),VecResetArray(x));
     130        7812 :     PetscCallAbort(PetscObjectComm((PetscObject)A),VecResetArray(y));
     131             :   }
     132        5940 :   PetscFunctionReturnVoid();
     133             : }
     134             : 
     135             : #if defined(SLEPC_HAVE_PRIMME3)
     136        5506 : static void massMatrixMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
     137             : {
     138        5506 :   PetscInt   i;
     139        5506 :   EPS_PRIMME *ops = (EPS_PRIMME*)primme->massMatrix;
     140        5506 :   Vec        x = ops->x,y = ops->y;
     141        5506 :   Mat        B = ops->B;
     142             : 
     143        5506 :   PetscFunctionBegin;
     144       11012 :   for (i=0;i<*blockSize;i++) {
     145        5506 :     PetscCallAbort(PetscObjectComm((PetscObject)B),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
     146        5506 :     PetscCallAbort(PetscObjectComm((PetscObject)B),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
     147        5506 :     PetscCallAbort(PetscObjectComm((PetscObject)B),MatMult(B,x,y));
     148        5506 :     PetscCallAbort(PetscObjectComm((PetscObject)B),VecResetArray(x));
     149        5506 :     PetscCallAbort(PetscObjectComm((PetscObject)B),VecResetArray(y));
     150             :   }
     151        5506 :   PetscFunctionReturnVoid();
     152             : }
     153             : #endif
     154             : 
     155        5336 : static void applyPreconditioner_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,struct primme_params *primme,int *ierr)
     156             : {
     157        5336 :   PetscInt   i;
     158        5336 :   EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
     159        5336 :   Vec        x = ops->x,y = ops->y;
     160             : 
     161        5336 :   PetscFunctionBegin;
     162       12375 :   for (i=0;i<*blockSize;i++) {
     163        7039 :     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
     164        7039 :     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
     165        7039 :     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),KSPSolve(ops->ksp,x,y));
     166        7039 :     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecResetArray(x));
     167        7039 :     PetscCallAbort(PetscObjectComm((PetscObject)ops->ksp),VecResetArray(y));
     168             :   }
     169        5336 :   PetscFunctionReturnVoid();
     170             : }
     171             : 
     172          10 : static PetscErrorCode EPSSetUp_PRIMME(EPS eps)
     173             : {
     174          10 :   PetscMPIInt    numProcs,procID;
     175          10 :   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
     176          10 :   primme_params  *primme = &ops->primme;
     177          10 :   PetscBool      flg;
     178             : 
     179          10 :   PetscFunctionBegin;
     180          10 :   EPSCheckHermitianDefinite(eps);
     181          10 :   EPSCheckNotStructured(eps);
     182          10 :   if (eps->nev==0) eps->nev = 1;
     183          10 :   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs));
     184          10 :   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&procID));
     185             : 
     186             :   /* Check some constraints and set some default values */
     187          10 :   if (eps->max_it==PETSC_DETERMINE) eps->max_it = PETSC_INT_MAX;
     188          10 :   PetscCall(STGetMatrix(eps->st,0,&ops->A));
     189          10 :   if (eps->isgeneralized) {
     190             : #if defined(SLEPC_HAVE_PRIMME3)
     191           1 :     PetscCall(STGetMatrix(eps->st,1,&ops->B));
     192             : #else
     193             :     SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This version of PRIMME is not available for generalized problems");
     194             : #endif
     195             :   }
     196          10 :   EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
     197          10 :   EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);
     198          10 :   if (!eps->which) eps->which = EPS_LARGEST_REAL;
     199             : #if !defined(SLEPC_HAVE_PRIMME2p2)
     200             :   if (eps->converged != EPSConvergedAbsolute) PetscCall(PetscInfo(eps,"Warning: using absolute convergence test\n"));
     201             : #else
     202          10 :   EPSCheckIgnored(eps,EPS_FEATURE_CONVERGENCE);
     203             : #endif
     204             : 
     205             :   /* Transfer SLEPc options to PRIMME options */
     206          10 :   primme_free(primme);
     207          10 :   primme_initialize(primme);
     208          10 :   primme->n                             = (PRIMME_INT)eps->n;
     209          10 :   primme->nLocal                        = (PRIMME_INT)eps->nloc;
     210          10 :   primme->numEvals                      = (int)eps->nev;
     211          10 :   primme->matrix                        = ops;
     212          10 :   primme->matrixMatvec                  = matrixMatvec_PRIMME;
     213             : #if defined(SLEPC_HAVE_PRIMME3)
     214          10 :   if (eps->isgeneralized) {
     215           1 :     primme->massMatrix                  = ops;
     216           1 :     primme->massMatrixMatvec            = massMatrixMatvec_PRIMME;
     217             :   }
     218             : #endif
     219          10 :   primme->commInfo                      = eps;
     220          10 :   primme->maxOuterIterations            = (PRIMME_INT)eps->max_it;
     221             : #if !defined(SLEPC_HAVE_PRIMME2p2)
     222             :   primme->eps                           = SlepcDefaultTol(eps->tol);
     223             : #endif
     224          10 :   primme->numProcs                      = numProcs;
     225          10 :   primme->procID                        = procID;
     226          10 :   primme->printLevel                    = 1;
     227          10 :   primme->correctionParams.precondition = 1;
     228          10 :   primme->globalSumReal                 = par_GlobalSumReal;
     229             : #if defined(SLEPC_HAVE_PRIMME3)
     230          10 :   primme->broadcastReal                 = par_broadcastReal;
     231             : #endif
     232             : #if defined(SLEPC_HAVE_PRIMME2p2)
     233          10 :   primme->convTestFun                   = convTestFun;
     234          10 :   primme->monitorFun                    = monitorFun;
     235             : #endif
     236          10 :   if (ops->bs > 0) primme->maxBlockSize = (int)ops->bs;
     237             : 
     238          10 :   switch (eps->which) {
     239           3 :     case EPS_LARGEST_REAL:
     240           3 :       primme->target = primme_largest;
     241           3 :       break;
     242           2 :     case EPS_SMALLEST_REAL:
     243           2 :       primme->target = primme_smallest;
     244           2 :       break;
     245           1 :     case EPS_LARGEST_MAGNITUDE:
     246           1 :       primme->target = primme_largest_abs;
     247           1 :       ops->target = 0.0;
     248           1 :       primme->numTargetShifts = 1;
     249           1 :       primme->targetShifts = &ops->target;
     250           1 :       break;
     251           1 :     case EPS_SMALLEST_MAGNITUDE:
     252           1 :       primme->target = primme_closest_abs;
     253           1 :       ops->target = 0.0;
     254           1 :       primme->numTargetShifts = 1;
     255           1 :       primme->targetShifts = &ops->target;
     256           1 :       break;
     257           3 :     case EPS_TARGET_MAGNITUDE:
     258             :     case EPS_TARGET_REAL:
     259           3 :       primme->target = primme_closest_abs;
     260           3 :       primme->numTargetShifts = 1;
     261           3 :       ops->target = PetscRealPart(eps->target);
     262           3 :       primme->targetShifts = &ops->target;
     263           3 :       break;
     264           0 :     default:
     265           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'which' value not supported by PRIMME");
     266             :   }
     267             : 
     268          10 :   switch (eps->extraction) {
     269           7 :     case EPS_RITZ:
     270           7 :       primme->projectionParams.projection = primme_proj_RR;
     271           7 :       break;
     272           1 :     case EPS_HARMONIC:
     273           1 :       primme->projectionParams.projection = primme_proj_harmonic;
     274           1 :       break;
     275           2 :     case EPS_REFINED:
     276           2 :       primme->projectionParams.projection = primme_proj_refined;
     277           2 :       break;
     278           0 :     default:
     279           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'extraction' value not supported by PRIMME");
     280             :   }
     281             : 
     282             :   /* If user sets mpd or ncv, maxBasisSize is modified */
     283          10 :   if (eps->mpd!=PETSC_DETERMINE) {
     284           1 :     primme->maxBasisSize = (int)eps->mpd;
     285           1 :     if (eps->ncv!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: 'ncv' is ignored by PRIMME\n"));
     286           9 :   } else if (eps->ncv!=PETSC_DETERMINE) primme->maxBasisSize = (int)eps->ncv;
     287             : 
     288          10 :   PetscCheck(primme_set_method(ops->method,primme)>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME method not valid");
     289             : 
     290          10 :   eps->mpd = (PetscInt)primme->maxBasisSize;
     291          10 :   eps->ncv = (PetscInt)(primme->locking?eps->nev:0)+primme->maxBasisSize;
     292          10 :   ops->bs  = (PetscInt)primme->maxBlockSize;
     293             : 
     294             :   /* Set workspace */
     295          10 :   PetscCall(EPSAllocateSolution(eps,0));
     296             : 
     297             :   /* Setup the preconditioner */
     298          10 :   if (primme->correctionParams.precondition) {
     299          10 :     PetscCall(STGetKSP(eps->st,&ops->ksp));
     300          10 :     PetscCall(PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg));
     301          10 :     if (!flg) PetscCall(PetscInfo(eps,"Warning: ignoring KSP, should use KSPPREONLY\n"));
     302          10 :     primme->preconditioner = NULL;
     303          10 :     primme->applyPreconditioner = applyPreconditioner_PRIMME;
     304             :   }
     305             : 
     306             :   /* Prepare auxiliary vectors */
     307          10 :   if (!ops->x) PetscCall(MatCreateVecsEmpty(ops->A,&ops->x,&ops->y));
     308          10 :   PetscFunctionReturn(PETSC_SUCCESS);
     309             : }
     310             : 
     311          10 : static PetscErrorCode EPSSolve_PRIMME(EPS eps)
     312             : {
     313          10 :   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
     314          10 :   PetscScalar    *a;
     315          10 :   PetscInt       i,ld,ierrprimme;
     316          10 :   PetscReal      *evals,*rnorms;
     317             : 
     318          10 :   PetscFunctionBegin;
     319             :   /* Reset some parameters left from previous runs */
     320             : #if defined(SLEPC_HAVE_PRIMME2p2)
     321          10 :   ops->primme.aNorm    = 0.0;
     322             : #else
     323             :   /* Force PRIMME to stop by absolute error */
     324             :   ops->primme.aNorm    = 1.0;
     325             : #endif
     326          10 :   ops->primme.initSize = (int)eps->nini;
     327          10 :   ops->primme.iseed[0] = -1;
     328          10 :   ops->primme.iseed[1] = -1;
     329          10 :   ops->primme.iseed[2] = -1;
     330          10 :   ops->primme.iseed[3] = -1;
     331          10 :   PetscCall(BVGetLeadingDimension(eps->V,&ld));
     332          10 :   ops->primme.ldevecs  = (PRIMME_INT)ld;
     333             : 
     334             :   /* Call PRIMME solver */
     335          10 :   PetscCall(BVGetArray(eps->V,&a));
     336          10 :   PetscCall(PetscMalloc2(eps->ncv,&evals,eps->ncv,&rnorms));
     337          10 :   ierrprimme = PRIMME_DRIVER(evals,a,rnorms,&ops->primme);
     338         267 :   for (i=0;i<eps->ncv;i++) eps->eigr[i] = evals[i];
     339         257 :   for (i=0;i<eps->ncv;i++) eps->errest[i] = rnorms[i];
     340          10 :   PetscCall(PetscFree2(evals,rnorms));
     341          10 :   PetscCall(BVRestoreArray(eps->V,&a));
     342             : 
     343          10 :   eps->nconv  = ops->primme.initSize >= 0 ? (PetscInt)ops->primme.initSize : 0;
     344          10 :   eps->reason = eps->nconv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
     345          10 :   PetscCall(PetscIntCast(ops->primme.stats.numOuterIterations,&eps->its));
     346             : 
     347             :   /* Process PRIMME error code */
     348          10 :   switch (ierrprimme) {
     349             :     case 0: /* no error */
     350             :       break;
     351           0 :     case -1:
     352           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": unexpected error",ierrprimme);
     353           0 :     case -2:
     354           0 :       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": allocation error",ierrprimme);
     355             :     case -3: /* stop due to maximum number of iterations or matvecs */
     356             :       break;
     357           0 :     default:
     358           0 :       PetscCheck(ierrprimme<-39,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": configuration error; check PRIMME's manual",ierrprimme);
     359           0 :       PetscCheck(ierrprimme>=-39,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": runtime error; check PRIMME's manual",ierrprimme);
     360             :   }
     361          10 :   PetscFunctionReturn(PETSC_SUCCESS);
     362             : }
     363             : 
     364           9 : static PetscErrorCode EPSReset_PRIMME(EPS eps)
     365             : {
     366           9 :   EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
     367             : 
     368           9 :   PetscFunctionBegin;
     369           9 :   primme_free(&ops->primme);
     370           9 :   PetscCall(VecDestroy(&ops->x));
     371           9 :   PetscCall(VecDestroy(&ops->y));
     372           9 :   PetscFunctionReturn(PETSC_SUCCESS);
     373             : }
     374             : 
     375           9 : static PetscErrorCode EPSDestroy_PRIMME(EPS eps)
     376             : {
     377           9 :   PetscFunctionBegin;
     378           9 :   PetscCall(PetscFree(eps->data));
     379           9 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL));
     380           9 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL));
     381           9 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL));
     382           9 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL));
     383           9 :   PetscFunctionReturn(PETSC_SUCCESS);
     384             : }
     385             : 
     386           0 : static PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
     387             : {
     388           0 :   PetscBool      isascii;
     389           0 :   EPS_PRIMME     *ctx = (EPS_PRIMME*)eps->data;
     390           0 :   PetscMPIInt    rank;
     391             : 
     392           0 :   PetscFunctionBegin;
     393           0 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
     394           0 :   if (isascii) {
     395           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  block size=%" PetscInt_FMT "\n",ctx->bs));
     396           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  solver method: %s\n",EPSPRIMMEMethods[(EPSPRIMMEMethod)ctx->method]));
     397             : 
     398             :     /* Display PRIMME params */
     399           0 :     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank));
     400           0 :     if (!rank) primme_display_params(ctx->primme);
     401             :   }
     402           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     403             : }
     404             : 
     405           8 : static PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps,PetscOptionItems *PetscOptionsObject)
     406             : {
     407           8 :   EPS_PRIMME      *ctx = (EPS_PRIMME*)eps->data;
     408           8 :   PetscInt        bs;
     409           8 :   EPSPRIMMEMethod meth;
     410           8 :   PetscBool       flg;
     411             : 
     412           8 :   PetscFunctionBegin;
     413           8 :   PetscOptionsHeadBegin(PetscOptionsObject,"EPS PRIMME Options");
     414             : 
     415           8 :     PetscCall(PetscOptionsInt("-eps_primme_blocksize","Maximum block size","EPSPRIMMESetBlockSize",ctx->bs,&bs,&flg));
     416           8 :     if (flg) PetscCall(EPSPRIMMESetBlockSize(eps,bs));
     417             : 
     418           8 :     PetscCall(PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg));
     419           8 :     if (flg) PetscCall(EPSPRIMMESetMethod(eps,meth));
     420             : 
     421           8 :   PetscOptionsHeadEnd();
     422           8 :   PetscFunctionReturn(PETSC_SUCCESS);
     423             : }
     424             : 
     425           5 : static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
     426             : {
     427           5 :   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
     428             : 
     429           5 :   PetscFunctionBegin;
     430           5 :   if (bs == PETSC_DEFAULT || bs == PETSC_DECIDE) ops->bs = 0;
     431             :   else {
     432           5 :     PetscCheck(bs>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
     433           5 :     ops->bs = bs;
     434             :   }
     435           5 :   PetscFunctionReturn(PETSC_SUCCESS);
     436             : }
     437             : 
     438             : /*@
     439             :    EPSPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.
     440             : 
     441             :    Logically Collective
     442             : 
     443             :    Input Parameters:
     444             : +  eps - the eigenproblem solver context
     445             : -  bs - block size
     446             : 
     447             :    Options Database Key:
     448             : .  -eps_primme_blocksize - Sets the max allowed block size value
     449             : 
     450             :    Notes:
     451             :    If the block size is not set, the value established by primme_initialize
     452             :    is used.
     453             : 
     454             :    The user should set the block size based on the architecture specifics
     455             :    of the target computer, as well as any a priori knowledge of multiplicities.
     456             :    The code does NOT require bs > 1 to find multiple eigenvalues. For some
     457             :    methods, keeping bs = 1 yields the best overall performance.
     458             : 
     459             :    Level: advanced
     460             : 
     461             : .seealso: EPSPRIMMEGetBlockSize()
     462             : @*/
     463           5 : PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
     464             : {
     465           5 :   PetscFunctionBegin;
     466           5 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     467          15 :   PetscValidLogicalCollectiveInt(eps,bs,2);
     468           5 :   PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
     469           5 :   PetscFunctionReturn(PETSC_SUCCESS);
     470             : }
     471             : 
     472           4 : static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
     473             : {
     474           4 :   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
     475             : 
     476           4 :   PetscFunctionBegin;
     477           4 :   *bs = ops->bs;
     478           4 :   PetscFunctionReturn(PETSC_SUCCESS);
     479             : }
     480             : 
     481             : /*@
     482             :    EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
     483             : 
     484             :    Not Collective
     485             : 
     486             :    Input Parameter:
     487             : .  eps - the eigenproblem solver context
     488             : 
     489             :    Output Parameter:
     490             : .  bs - returned block size
     491             : 
     492             :    Level: advanced
     493             : 
     494             : .seealso: EPSPRIMMESetBlockSize()
     495             : @*/
     496           4 : PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
     497             : {
     498           4 :   PetscFunctionBegin;
     499           4 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     500           4 :   PetscAssertPointer(bs,2);
     501           4 :   PetscUseMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
     502           4 :   PetscFunctionReturn(PETSC_SUCCESS);
     503             : }
     504             : 
     505           5 : static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
     506             : {
     507           5 :   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
     508             : 
     509           5 :   PetscFunctionBegin;
     510           5 :   ops->method = (primme_preset_method)method;
     511           5 :   PetscFunctionReturn(PETSC_SUCCESS);
     512             : }
     513             : 
     514             : /*@
     515             :    EPSPRIMMESetMethod - Sets the method for the PRIMME library.
     516             : 
     517             :    Logically Collective
     518             : 
     519             :    Input Parameters:
     520             : +  eps - the eigenproblem solver context
     521             : -  method - method that will be used by PRIMME
     522             : 
     523             :    Options Database Key:
     524             : .  -eps_primme_method - Sets the method for the PRIMME library
     525             : 
     526             :    Note:
     527             :    If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
     528             : 
     529             :    Level: advanced
     530             : 
     531             : .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
     532             : @*/
     533           5 : PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
     534             : {
     535           5 :   PetscFunctionBegin;
     536           5 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     537          15 :   PetscValidLogicalCollectiveEnum(eps,method,2);
     538           5 :   PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
     539           5 :   PetscFunctionReturn(PETSC_SUCCESS);
     540             : }
     541             : 
     542           4 : static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
     543             : {
     544           4 :   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
     545             : 
     546           4 :   PetscFunctionBegin;
     547           4 :   *method = (EPSPRIMMEMethod)ops->method;
     548           4 :   PetscFunctionReturn(PETSC_SUCCESS);
     549             : }
     550             : 
     551             : /*@
     552             :    EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
     553             : 
     554             :    Not Collective
     555             : 
     556             :    Input Parameter:
     557             : .  eps - the eigenproblem solver context
     558             : 
     559             :    Output Parameter:
     560             : .  method - method that will be used by PRIMME
     561             : 
     562             :    Level: advanced
     563             : 
     564             : .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
     565             : @*/
     566           4 : PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
     567             : {
     568           4 :   PetscFunctionBegin;
     569           4 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
     570           4 :   PetscAssertPointer(method,2);
     571           4 :   PetscUseMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
     572           4 :   PetscFunctionReturn(PETSC_SUCCESS);
     573             : }
     574             : 
     575           9 : SLEPC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
     576             : {
     577           9 :   EPS_PRIMME     *primme;
     578             : 
     579           9 :   PetscFunctionBegin;
     580           9 :   PetscCall(PetscNew(&primme));
     581           9 :   eps->data = (void*)primme;
     582             : 
     583           9 :   primme_initialize(&primme->primme);
     584           9 :   primme->primme.globalSumReal = par_GlobalSumReal;
     585             : #if defined(SLEPC_HAVE_PRIMME3)
     586           9 :   primme->primme.broadcastReal = par_broadcastReal;
     587             : #endif
     588             : #if defined(SLEPC_HAVE_PRIMME2p2)
     589           9 :   primme->primme.convTestFun = convTestFun;
     590           9 :   primme->primme.monitorFun = monitorFun;
     591             : #endif
     592           9 :   primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
     593             : 
     594           9 :   eps->categ = EPS_CATEGORY_PRECOND;
     595             : 
     596           9 :   eps->ops->solve          = EPSSolve_PRIMME;
     597           9 :   eps->ops->setup          = EPSSetUp_PRIMME;
     598           9 :   eps->ops->setupsort      = EPSSetUpSort_Basic;
     599           9 :   eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
     600           9 :   eps->ops->destroy        = EPSDestroy_PRIMME;
     601           9 :   eps->ops->reset          = EPSReset_PRIMME;
     602           9 :   eps->ops->view           = EPSView_PRIMME;
     603           9 :   eps->ops->backtransform  = EPSBackTransform_Default;
     604           9 :   eps->ops->setdefaultst   = EPSSetDefaultST_GMRES;
     605             : 
     606           9 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME));
     607           9 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME));
     608           9 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME));
     609           9 :   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME));
     610           9 :   PetscFunctionReturn(PETSC_SUCCESS);
     611             : }

Generated by: LCOV version 1.14