Actual source code: primme.c

  1: /*
  2:        This file implements a wrapper to the PRIMME library

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.
  9:       
 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <petscsys.h>
 25: #include <private/epsimpl.h>    /*I "slepceps.h" I*/
 26: #include <private/stimpl.h>     /*I "slepcst.h" I*/

 28: PetscErrorCode EPSSolve_PRIMME(EPS);

 31: #include <primme.h>

 34: typedef struct {
 35:   primme_params primme;           /* param struc */
 36:   primme_preset_method method;    /* primme method */
 37:   Mat A;                          /* problem matrix */
 38:   EPS eps;                        /* EPS current context */
 39:   KSP ksp;                        /* preconditioner */
 40:   Vec x,y;                        /* auxiliar vectors */
 41: } EPS_PRIMME;

 43: EPSPRIMMEMethod methodN[] = {
 44:   EPS_PRIMME_DYNAMIC,
 45:   EPS_PRIMME_DEFAULT_MIN_TIME,
 46:   EPS_PRIMME_DEFAULT_MIN_MATVECS,
 47:   EPS_PRIMME_ARNOLDI,
 48:   EPS_PRIMME_GD,
 49:   EPS_PRIMME_GD_PLUSK,
 50:   EPS_PRIMME_GD_OLSEN_PLUSK,
 51:   EPS_PRIMME_JD_OLSEN_PLUSK,
 52:   EPS_PRIMME_RQI,
 53:   EPS_PRIMME_JDQR,
 54:   EPS_PRIMME_JDQMR,
 55:   EPS_PRIMME_JDQMR_ETOL,
 56:   EPS_PRIMME_SUBSPACE_ITERATION,
 57:   EPS_PRIMME_LOBPCG_ORTHOBASIS,
 58:   EPS_PRIMME_LOBPCG_ORTHOBASISW
 59: };

 61: static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme);
 62: static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme);

 64: void par_GlobalSumDouble(void *sendBuf,void *recvBuf,int *count,primme_params *primme)
 65: {
 67:   MPI_Allreduce((double*)sendBuf,(double*)recvBuf,*count,MPI_DOUBLE,MPI_SUM,((PetscObject)(primme->commInfo))->comm);CHKERRABORT(((PetscObject)(primme->commInfo))->comm,ierr);
 68: }

 72: PetscErrorCode EPSSetUp_PRIMME(EPS eps)
 73: {
 75:   PetscMPIInt    numProcs,procID;
 76:   EPS_PRIMME     *ops = (EPS_PRIMME *)eps->data;
 77:   primme_params  *primme = &(((EPS_PRIMME *)eps->data)->primme);
 78:   PetscBool      t;

 81:   MPI_Comm_size(((PetscObject)eps)->comm,&numProcs);
 82:   MPI_Comm_rank(((PetscObject)eps)->comm,&procID);
 83: 
 84:   /* Check some constraints and set some default values */
 85:   if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
 86:   STGetOperators(eps->OP,&ops->A,PETSC_NULL);
 87:   if (!ops->A) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"The problem matrix has to be specified first");
 88:   if (!eps->ishermitian)
 89:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME is only available for Hermitian problems");
 90:   if (eps->isgeneralized)
 91:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME is not available for generalized problems");
 92:   if (!eps->which) eps->which = EPS_LARGEST_REAL;

 94:   /* Change the default sigma to inf if necessary */
 95:   if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
 96:       eps->which == EPS_LARGEST_IMAGINARY) {
 97:     STSetDefaultShift(eps->OP,3e300);
 98:   }

100:   /* Avoid setting the automatic shift when a target is set */
101:   STSetDefaultShift(eps->OP,0.0);

103:   STSetUp(eps->OP);
104:   PetscTypeCompare((PetscObject)eps->OP,STPRECOND,&t);
105:   if (!t) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME only works with STPRECOND");

107:   /* Transfer SLEPc options to PRIMME options */
108:   primme->n = eps->n;
109:   primme->nLocal = eps->nloc;
110:   primme->numEvals = eps->nev;
111:   primme->matrix = ops;
112:   primme->commInfo = eps;
113:   primme->maxMatvecs = eps->max_it;
114:   primme->eps = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
115:   primme->numProcs = numProcs;
116:   primme->procID = procID;
117:   primme->printLevel = 0;
118:   primme->correctionParams.precondition = 1;

120:   if (!eps->which) eps->which = EPS_LARGEST_REAL;
121:   switch(eps->which) {
122:     case EPS_LARGEST_REAL:
123:       primme->target = primme_largest;
124:       break;
125:     case EPS_SMALLEST_REAL:
126:       primme->target = primme_smallest;
127:       break;
128:     case EPS_TARGET_MAGNITUDE:
129:     case EPS_TARGET_REAL:
130:       primme->target = primme_closest_abs;
131:       primme->numTargetShifts = 1;
132:       primme->targetShifts = &eps->target;
133:       break;
134:     default:
135:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"'which' value does not supported by PRIMME");
136:       break;
137:   }
138: 
139:   if (primme_set_method(ops->method,primme) < 0)
140:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME method not valid");
141: 
142:   /* If user sets ncv, maxBasisSize is modified. If not, ncv is set as maxBasisSize */
143:   if (eps->ncv) primme->maxBasisSize = eps->ncv;
144:   else eps->ncv = primme->maxBasisSize;
145:   if (eps->ncv < eps->nev+primme->maxBlockSize)
146:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME needs ncv >= nev+maxBlockSize");
147:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }

149:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }

151:   /* Set workspace */
152:   EPSAllocateSolution(eps);

154:   /* Setup the preconditioner */
155:   ops->eps = eps;
156:   if (primme->correctionParams.precondition) {
157:     STGetKSP(eps->OP,&ops->ksp);
158:     PetscTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&t);
159:     if (!t) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME only works with KSPPREONLY");
160:     primme->preconditioner = PETSC_NULL;
161:     primme->applyPreconditioner = applyPreconditioner_PRIMME;
162:   }

164:   /* Prepare auxiliary vectors */
165:   VecCreateMPIWithArray(PETSC_COMM_WORLD,eps->nloc,eps->n,PETSC_NULL,&ops->x);
166:   VecCreateMPIWithArray(PETSC_COMM_WORLD,eps->nloc,eps->n,PETSC_NULL,&ops->y);
167: 
168:   /* dispatch solve method */
169:   if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
170:   eps->ops->solve = EPSSolve_PRIMME;
171:   return(0);
172: }

176: PetscErrorCode EPSSolve_PRIMME(EPS eps)
177: {
179:   EPS_PRIMME     *ops = (EPS_PRIMME *)eps->data;
180:   PetscScalar    *a;
181: #if defined(PETSC_USE_COMPLEX)
182:   PetscInt       i;
183:   PetscReal      *evals;
184: #endif

187:   /* Reset some parameters left from previous runs */
188:   ops->primme.aNorm    = 0.0;
189:   ops->primme.initSize = eps->nini;
190:   ops->primme.iseed[0] = -1;

192:   /* Call PRIMME solver */
193:   VecGetArray(eps->V[0],&a);
194: #if !defined(PETSC_USE_COMPLEX)
195:   dprimme(eps->eigr,a,eps->errest,&ops->primme);
196: #else
197:   /* PRIMME returns real eigenvalues, but SLEPc works with complex ones */
198:   PetscMalloc(eps->ncv*sizeof(PetscReal),&evals);
199:   zprimme(evals,(Complex_Z*)a,eps->errest,&ops->primme);
200:   for (i=0;i<eps->ncv;i++)
201:     eps->eigr[i] = evals[i];
202:   PetscFree(evals);
203: #endif
204:   VecRestoreArray(eps->V[0],&a);
205: 
206:   switch(ierr) {
207:     case 0: /* Successful */
208:       break;

210:     case -1:
211:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: Failed to open output file");
212:       break;

214:     case -2:
215:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: Insufficient integer or real workspace allocated");
216:       break;

218:     case -3:
219:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: main_iter encountered a problem");
220:       break;

222:     default:
223:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: some parameters wrong configured");
224:       break;
225:   }

227:   eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
228:   eps->reason = eps->ncv >= eps->nev ? EPS_CONVERGED_TOL : EPS_DIVERGED_ITS;
229:   eps->its = ops->primme.stats.numOuterIterations;
230:   eps->OP->applys = ops->primme.stats.numMatvecs;
231:   return(0);
232: }

236: static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme)
237: {
239:   PetscInt       i,N = primme->n;
240:   EPS_PRIMME     *ops = (EPS_PRIMME *)primme->matrix;
241:   Vec            x = ops->x,y = ops->y;
242:   Mat            A = ops->A;

245:   for (i=0;i<*blockSize;i++) {
246:     /* build vectors using 'in' an 'out' workspace */
247:     VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);
248:     VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);

250:     MatMult(A,x,y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
251: 
252:     VecResetArray(x);CHKERRABORT(PETSC_COMM_WORLD,ierr);
253:     VecResetArray(y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
254:   }
255:   PetscFunctionReturnVoid();
256: }

260: static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme)
261: {
263:   PetscInt       i,N = primme->n,lits;
264:   EPS_PRIMME     *ops = (EPS_PRIMME *)primme->matrix;
265:   Vec            x = ops->x,y = ops->y;
266: 
268:   for (i=0;i<*blockSize;i++) {
269:     /* build vectors using 'in' an 'out' workspace */
270:     VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);
271:     VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);

273:     KSPSolve(ops->ksp,x,y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
274:     KSPGetIterationNumber(ops->ksp,&lits);CHKERRABORT(PETSC_COMM_WORLD,ierr);
275:     ops->eps->OP->lineariterations+= lits;
276: 
277:     VecResetArray(x);CHKERRABORT(PETSC_COMM_WORLD,ierr);
278:     VecResetArray(y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
279:   }
280:   PetscFunctionReturnVoid();
281: }

285: PetscErrorCode EPSReset_PRIMME(EPS eps)
286: {
288:   EPS_PRIMME     *ops = (EPS_PRIMME *)eps->data;

291:   primme_Free(&ops->primme);
292:   VecDestroy(&ops->x);
293:   VecDestroy(&ops->y);
294:   EPSFreeSolution(eps);
295:   return(0);
296: }

300: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
301: {

305:   PetscFree(eps->data);
306:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetBlockSize_C","",PETSC_NULL);
307:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetMethod_C","",PETSC_NULL);
308:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetBlockSize_C","",PETSC_NULL);
309:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetMethod_C","",PETSC_NULL);
310:   return(0);
311: }

315: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
316: {
317:   PetscErrorCode  ierr;
318:   PetscBool       isascii;
319:   primme_params   *primme = &((EPS_PRIMME *)eps->data)->primme;
320:   EPSPRIMMEMethod methodn;
321:   PetscMPIInt     rank;

324:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
325:   if (!isascii) {
326:     SETERRQ1(((PetscObject)eps)->comm,1,"Viewer type %s not supported for EPSPRIMME",((PetscObject)viewer)->type_name);
327:   }
328: 
329:   PetscViewerASCIIPrintf(viewer,"  PRIMME: block size=%d\n",primme->maxBlockSize);
330:   EPSPRIMMEGetMethod(eps,&methodn);
331:   PetscViewerASCIIPrintf(viewer,"  PRIMME: solver method: %s\n",EPSPRIMMEMethods[methodn]);

333:   /* Display PRIMME params */
334:   MPI_Comm_rank(((PetscObject)eps)->comm,&rank);
335:   if (!rank) primme_display_params(*primme);
336:   return(0);
337: }

341: PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps)
342: {
343:   PetscErrorCode  ierr;
344:   EPS_PRIMME      *ctx = (EPS_PRIMME *)eps->data;
345:   PetscInt        bs;
346:   EPSPRIMMEMethod meth;
347:   PetscBool       flg;

350:   PetscOptionsHead("EPS PRIMME Options");
351:   PetscOptionsInt("-eps_primme_block_size","Maximum block size","EPSPRIMMESetBlockSize",ctx->primme.maxBlockSize,&bs,&flg);
352:   if (flg) { EPSPRIMMESetBlockSize(eps,bs); }
353:   PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
354:   if (flg) { EPSPRIMMESetMethod(eps,meth); }
355:   PetscOptionsTail();
356:   return(0);
357: }

362: PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
363: {
364:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

367:   if (bs == PETSC_DEFAULT) ops->primme.maxBlockSize = 1;
368:   else if (bs <= 0) {
369:     SETERRQ(((PetscObject)eps)->comm,1,"PRIMME: wrong block size");
370:   } else ops->primme.maxBlockSize = bs;
371:   return(0);
372: }

377: /*@
378:    EPSPRIMMESetBlockSize - The maximum block size the code will try to use. 
379:    The user should set
380:    this based on the architecture specifics of the target computer, 
381:    as well as any a priori knowledge of multiplicities. The code does 
382:    NOT require BlockSize > 1 to find multiple eigenvalues.  For some 
383:    methods, keeping BlockSize = 1 yields the best overall performance.

385:    Collective on EPS

387:    Input Parameters:
388: +  eps - the eigenproblem solver context
389: -  bs - block size

391:    Options Database Key:
392: .  -eps_primme_block_size - Sets the max allowed block size value

394:    Notes:
395:    If the block size is not set, the value established by primme_initialize
396:    is used.

398:    Level: advanced
399: .seealso: EPSPRIMMEGetBlockSize()
400: @*/
401: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
402: {

408:   PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
409:   return(0);
410: }

415: PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
416: {
417:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

420:   if (bs) *bs = ops->primme.maxBlockSize;
421:   return(0);
422: }

427: /*@
428:    EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use. 

430:    Collective on EPS

432:    Input Parameters:
433: .  eps - the eigenproblem solver context
434:     
435:    Output Parameters:  
436: .  bs - returned block size 

438:    Level: advanced
439: .seealso: EPSPRIMMESetBlockSize()
440: @*/
441: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
442: {

447:   PetscTryMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
448:   return(0);
449: }

454: PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
455: {
456:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

459:   if (method == PETSC_DEFAULT) ops->method = DEFAULT_MIN_TIME;
460:   else ops->method = (primme_preset_method)method;
461:   return(0);
462: }

467: /*@
468:    EPSPRIMMESetMethod - Sets the method for the PRIMME library.

470:    Collective on EPS

472:    Input Parameters:
473: +  eps - the eigenproblem solver context
474: -  method - method that will be used by PRIMME. It must be one of:
475:     EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
476:     EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
477:     EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK, 
478:     EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR, 
479:     EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
480:     EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW

482:    Options Database Key:
483: .  -eps_primme_set_method - Sets the method for the PRIMME library.

485:    Note:
486:    If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.

488:    Level: advanced

490: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
491: @*/
492: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
493: {

499:   PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
500:   return(0);
501: }

506: PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
507: {
508:   EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;

511:   if (method) *method = (EPSPRIMMEMethod)ops->method;
512:   return(0);
513: }

518: /*@C
519:     EPSPRIMMEGetMethod - Gets the method for the PRIMME library.

521:     Mon Collective on EPS

523:    Input Parameters:
524: .  eps - the eigenproblem solver context
525:     
526:    Output Parameters: 
527: .  method - method that will be used by PRIMME. It must be one of:
528:     EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
529:     EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
530:     EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK, 
531:     EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR, 
532:     EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
533:     EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW

535:     Level: advanced

537: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
538: @*/
539: PetscErrorCode EPSPRIMMEGetMethod(EPS eps, EPSPRIMMEMethod *method)
540: {

545:   PetscTryMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
546:   return(0);
547: }

552: PetscErrorCode EPSCreate_PRIMME(EPS eps)
553: {
555:   EPS_PRIMME     *primme;

558:   STSetType(eps->OP,STPRECOND);
559:   STPrecondSetKSPHasMat(eps->OP,PETSC_TRUE);

561:   PetscNewLog(eps,EPS_PRIMME,&primme);
562:   eps->data                      = (void*)primme;
563:   eps->ops->setup                = EPSSetUp_PRIMME;
564:   eps->ops->setfromoptions       = EPSSetFromOptions_PRIMME;
565:   eps->ops->destroy              = EPSDestroy_PRIMME;
566:   eps->ops->reset                = EPSReset_PRIMME;
567:   eps->ops->view                 = EPSView_PRIMME;
568:   eps->ops->backtransform        = EPSBackTransform_Default;
569:   eps->ops->computevectors       = EPSComputeVectors_Default;

571:   primme_initialize(&primme->primme);
572:   primme->primme.matrixMatvec = multMatvec_PRIMME;
573:   primme->primme.globalSumDouble = par_GlobalSumDouble;
574:   primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;

576:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetBlockSize_C","EPSPRIMMESetBlockSize_PRIMME",EPSPRIMMESetBlockSize_PRIMME);
577:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetMethod_C","EPSPRIMMESetMethod_PRIMME",EPSPRIMMESetMethod_PRIMME);
578:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetBlockSize_C","EPSPRIMMEGetBlockSize_PRIMME",EPSPRIMMEGetBlockSize_PRIMME);
579:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetMethod_C","EPSPRIMMEGetMethod_PRIMME",EPSPRIMMEGetMethod_PRIMME);
580:   return(0);
581: }