Actual source code: qarnoldi.c

slepc-3.7.4 2017-05-17
Report Typos and Errors
  1: /*

  3:    SLEPc quadratic eigensolver: "qarnoldi"

  5:    Method: Q-Arnoldi

  7:    Algorithm:

  9:        Quadratic Arnoldi with Krylov-Schur type restart.

 11:    References:

 13:        [1] K. Meerbergen, "The Quadratic Arnoldi method for the solution
 14:            of the quadratic eigenvalue problem", SIAM J. Matrix Anal.
 15:            Appl. 30(4):1462-1482, 2008.

 17:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 18:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 19:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 21:    This file is part of SLEPc.

 23:    SLEPc is free software: you can redistribute it and/or modify it under  the
 24:    terms of version 3 of the GNU Lesser General Public License as published by
 25:    the Free Software Foundation.

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

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

 37: #include <slepc/private/pepimpl.h>    /*I "slepcpep.h" I*/
 38: #include <petscblaslapack.h>

 40: typedef struct {
 41:   PetscReal keep;         /* restart parameter */
 42:   PetscBool lock;         /* locking/non-locking variant */
 43: } PEP_QARNOLDI;

 47: PetscErrorCode PEPSetUp_QArnoldi(PEP pep)
 48: {
 50:   PEP_QARNOLDI   *ctx = (PEP_QARNOLDI*)pep->data;
 51:   PetscBool      shift,sinv,flg;

 54:   pep->lineariz = PETSC_TRUE;
 55:   PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
 56:   if (!ctx->lock && pep->mpd<pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
 57:   if (!pep->max_it) pep->max_it = PetscMax(100,4*pep->n/pep->ncv);
 58:   /* Set STSHIFT as the default ST */
 59:   if (!((PetscObject)pep->st)->type_name) {
 60:     STSetType(pep->st,STSHIFT);
 61:   }
 62:   PetscObjectTypeCompare((PetscObject)pep->st,STSHIFT,&shift);
 63:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
 64:   if (!shift && !sinv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Only STSHIFT and STSINVERT spectral transformations can be used");
 65:   if (!pep->which) {
 66:     if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
 67:     else pep->which = PEP_LARGEST_MAGNITUDE;
 68:   }

 70:   if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver only available for quadratic problems");
 71:   if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
 72:   STGetTransform(pep->st,&flg);
 73:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag set, see STSetTransform()");

 75:   /* set default extraction */
 76:   if (!pep->extract) {
 77:     pep->extract = PEP_EXTRACT_NONE;
 78:   }
 79:   if (pep->extract!=PEP_EXTRACT_NONE) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver does not support requested extraction");
 80:  
 81:   if (!ctx->keep) ctx->keep = 0.5;

 83:   PEPAllocateSolution(pep,0);
 84:   PEPSetWorkVecs(pep,4);

 86:   DSSetType(pep->ds,DSNHEP);
 87:   DSSetExtraRow(pep->ds,PETSC_TRUE);
 88:   DSAllocate(pep->ds,pep->ncv+1);

 90:   /* process starting vector */
 91:   if (pep->nini>-2) {
 92:     BVSetRandomColumn(pep->V,0);
 93:     BVSetRandomColumn(pep->V,1);
 94:   } else {
 95:     BVInsertVec(pep->V,0,pep->IS[0]);
 96:     BVInsertVec(pep->V,1,pep->IS[1]);
 97:   }
 98:   if (pep->nini<0) {
 99:     SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
100:   }
101:   return(0);
102: }

106: PetscErrorCode PEPExtractVectors_QArnoldi(PEP pep)
107: {
109:   PetscInt       i,k=pep->nconv,ldds;
110:   PetscScalar    *X,*pX0;
111:   Mat            X0;

114:   if (pep->nconv==0) return(0);
115:   DSGetLeadingDimension(pep->ds,&ldds);
116:   DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
117:   DSGetArray(pep->ds,DS_MAT_X,&X);

119:   /* update vectors V = V*X */ 
120:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&X0);
121:   MatDenseGetArray(X0,&pX0);
122:   for (i=0;i<k;i++) {
123:     PetscMemcpy(pX0+i*k,X+i*ldds,k*sizeof(PetscScalar));
124:   }
125:   MatDenseRestoreArray(X0,&pX0);
126:   BVMultInPlace(pep->V,X0,0,k);
127:   MatDestroy(&X0);
128:   DSRestoreArray(pep->ds,DS_MAT_X,&X);
129:   return(0);
130: }

134: /*
135:   Compute a step of Classical Gram-Schmidt orthogonalization
136: */
137: static PetscErrorCode PEPQArnoldiCGS(PEP pep,PetscScalar *H,PetscBLASInt ldh,PetscScalar *h,PetscBLASInt j,BV V,Vec t,Vec v,Vec w,PetscReal *onorm,PetscReal *norm,PetscScalar *work)
138: {
140:   PetscBLASInt   ione = 1,j_1 = j+1;
141:   PetscReal      x,y;
142:   PetscScalar    dot,one = 1.0,zero = 0.0;

145:   /* compute norm of v and w */
146:   if (onorm) {
147:     VecNorm(v,NORM_2,&x);
148:     VecNorm(w,NORM_2,&y);
149:     *onorm = PetscSqrtReal(x*x+y*y);
150:   }

152:   /* orthogonalize: compute h */
153:   BVDotVec(V,v,h);
154:   BVDotVec(V,w,work);
155:   if (j>0)
156:     PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione));
157:   VecDot(w,t,&dot);
158:   h[j] += dot;

160:   /* orthogonalize: update v and w */
161:   BVMultVec(V,-1.0,1.0,v,h);
162:   if (j>0) {
163:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione));
164:     BVMultVec(V,-1.0,1.0,w,work);
165:   }
166:   VecAXPY(w,-h[j],t);

168:   /* compute norm of v and w */
169:   if (norm) {
170:     VecNorm(v,NORM_2,&x);
171:     VecNorm(w,NORM_2,&y);
172:     *norm = PetscSqrtReal(x*x+y*y);
173:   }
174:   return(0);
175: }

179: /*
180:   Compute a run of Q-Arnoldi iterations
181: */
182: static PetscErrorCode PEPQArnoldi(PEP pep,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
183: {
184:   PetscErrorCode     ierr;
185:   PetscInt           i,j,l,m = *M;
186:   Vec                t = pep->work[2],u = pep->work[3];
187:   BVOrthogRefineType refinement;
188:   PetscReal          norm=0.0,onorm,eta;
189:   PetscScalar        *c = work + m;

192:   BVGetOrthogonalization(pep->V,NULL,&refinement,&eta,NULL);
193:   BVInsertVec(pep->V,k,v);
194:   for (j=k;j<m;j++) {
195:     /* apply operator */
196:     VecCopy(w,t);
197:     if (pep->Dr) {
198:       VecPointwiseMult(v,v,pep->Dr);
199:     }
200:     STMatMult(pep->st,0,v,u);
201:     VecCopy(t,v);
202:     if (pep->Dr) {
203:       VecPointwiseMult(t,t,pep->Dr);
204:     }
205:     STMatMult(pep->st,1,t,w);
206:     VecAXPY(u,pep->sfactor,w);
207:     STMatSolve(pep->st,u,w);
208:     VecScale(w,-1.0/(pep->sfactor*pep->sfactor));
209:     if (pep->Dr) {
210:       VecPointwiseDivide(w,w,pep->Dr);
211:     }
212:     VecCopy(v,t);
213:     BVSetActiveColumns(pep->V,0,j+1);

215:     /* orthogonalize */
216:     switch (refinement) {
217:       case BV_ORTHOG_REFINE_NEVER:
218:         PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,&norm,work);
219:         *breakdown = PETSC_FALSE;
220:         break;
221:       case BV_ORTHOG_REFINE_ALWAYS:
222:         PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,NULL,work);
223:         PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,&onorm,&norm,work);
224:         for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
225:         if (norm < eta * onorm) *breakdown = PETSC_TRUE;
226:         else *breakdown = PETSC_FALSE;
227:         break;
228:       case BV_ORTHOG_REFINE_IFNEEDED:
229:         PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,&onorm,&norm,work);
230:         /* ||q|| < eta ||h|| */
231:         l = 1;
232:         while (l<3 && norm < eta * onorm) {
233:           l++;
234:           onorm = norm;
235:           PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,NULL,&norm,work);
236:           for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
237:         }
238:         if (norm < eta * onorm) *breakdown = PETSC_TRUE;
239:         else *breakdown = PETSC_FALSE;
240:         break;
241:       default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of ip->orth_ref");
242:     }
243:     VecScale(v,1.0/norm);
244:     VecScale(w,1.0/norm);

246:     H[j+1+ldh*j] = norm;
247:     if (j<m-1) {
248:       BVInsertVec(pep->V,j+1,v);
249:     }
250:   }
251:   *beta = norm;
252:   return(0);
253: }

257: PetscErrorCode PEPSolve_QArnoldi(PEP pep)
258: {
260:   PEP_QARNOLDI   *ctx = (PEP_QARNOLDI*)pep->data;
261:   PetscInt       j,k,l,lwork,nv,ld,newn,nconv;
262:   Vec            v=pep->work[0],w=pep->work[1];
263:   Mat            Q;
264:   PetscScalar    *S,*work;
265:   PetscReal      beta=0.0,norm,x,y;
266:   PetscBool      breakdown=PETSC_FALSE,sinv;

269:   DSGetLeadingDimension(pep->ds,&ld);
270:   lwork = 7*pep->ncv;
271:   PetscMalloc1(lwork,&work);
272:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
273:   RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
274:   STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);

276:   /* Get the starting Arnoldi vector */
277:   BVCopyVec(pep->V,0,v);
278:   BVCopyVec(pep->V,1,w);
279:   VecNorm(v,NORM_2,&x);
280:   VecNorm(w,NORM_2,&y);
281:   norm = PetscSqrtReal(x*x+y*y);
282:   VecScale(v,1.0/norm);
283:   VecScale(w,1.0/norm);

285:    /* Restart loop */
286:   l = 0;
287:   while (pep->reason == PEP_CONVERGED_ITERATING) {
288:     pep->its++;

290:     /* Compute an nv-step Arnoldi factorization */
291:     nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
292:     DSGetArray(pep->ds,DS_MAT_A,&S);
293:     PEPQArnoldi(pep,S,ld,pep->nconv+l,&nv,v,w,&beta,&breakdown,work);
294:     DSRestoreArray(pep->ds,DS_MAT_A,&S);
295:     DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
296:     if (l==0) {
297:       DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
298:     } else {
299:       DSSetState(pep->ds,DS_STATE_RAW);
300:     }
301:     BVSetActiveColumns(pep->V,pep->nconv,nv);

303:     /* Solve projected problem */
304:     DSSolve(pep->ds,pep->eigr,pep->eigi);
305:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
306:     DSUpdateExtraRow(pep->ds);

308:     /* Check convergence */
309:     PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
310:     (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
311:     nconv = k;

313:     /* Update l */
314:     if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
315:     else l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
316:     if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */

318:     if (pep->reason == PEP_CONVERGED_ITERATING) {
319:       if (breakdown) {
320:         /* Stop if breakdown */
321:         PetscInfo2(pep,"Breakdown Quadratic Arnoldi method (it=%D norm=%g)\n",pep->its,(double)beta);
322:         pep->reason = PEP_DIVERGED_BREAKDOWN;
323:       } else {
324:         /* Prepare the Rayleigh quotient for restart */
325:         DSTruncate(pep->ds,k+l);
326:         DSGetDimensions(pep->ds,&newn,NULL,NULL,NULL,NULL);
327:         l = newn-k;
328:       }
329:     }
330:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
331:     DSGetMat(pep->ds,DS_MAT_Q,&Q);
332:     BVMultInPlace(pep->V,Q,pep->nconv,k+l);
333:     MatDestroy(&Q);

335:     pep->nconv = k;
336:     PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
337:   }

339:   for (j=0;j<pep->nconv;j++) {
340:     pep->eigr[j] *= pep->sfactor;
341:     pep->eigi[j] *= pep->sfactor;
342:   }

344:   STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
345:   RGPopScale(pep->rg);

347:   /* truncate Schur decomposition and change the state to raw so that
348:      DSVectors() computes eigenvectors from scratch */
349:   DSSetDimensions(pep->ds,pep->nconv,0,0,0);
350:   DSSetState(pep->ds,DS_STATE_RAW);
351:   PetscFree(work);
352:   return(0);
353: }

357: static PetscErrorCode PEPQArnoldiSetRestart_QArnoldi(PEP pep,PetscReal keep)
358: {
359:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

362:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
363:   else {
364:     if (keep<0.1 || keep>0.9) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
365:     ctx->keep = keep;
366:   }
367:   return(0);
368: }

372: /*@
373:    PEPQArnoldiSetRestart - Sets the restart parameter for the Q-Arnoldi
374:    method, in particular the proportion of basis vectors that must be kept
375:    after restart.

377:    Logically Collective on PEP

379:    Input Parameters:
380: +  pep  - the eigenproblem solver context
381: -  keep - the number of vectors to be kept at restart

383:    Options Database Key:
384: .  -pep_qarnoldi_restart - Sets the restart parameter

386:    Notes:
387:    Allowed values are in the range [0.1,0.9]. The default is 0.5.

389:    Level: advanced

391: .seealso: PEPQArnoldiGetRestart()
392: @*/
393: PetscErrorCode PEPQArnoldiSetRestart(PEP pep,PetscReal keep)
394: {

400:   PetscTryMethod(pep,"PEPQArnoldiSetRestart_C",(PEP,PetscReal),(pep,keep));
401:   return(0);
402: }

406: static PetscErrorCode PEPQArnoldiGetRestart_QArnoldi(PEP pep,PetscReal *keep)
407: {
408:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

411:   *keep = ctx->keep;
412:   return(0);
413: }

417: /*@
418:    PEPQArnoldiGetRestart - Gets the restart parameter used in the Q-Arnoldi method.

420:    Not Collective

422:    Input Parameter:
423: .  pep - the eigenproblem solver context

425:    Output Parameter:
426: .  keep - the restart parameter

428:    Level: advanced

430: .seealso: PEPQArnoldiSetRestart()
431: @*/
432: PetscErrorCode PEPQArnoldiGetRestart(PEP pep,PetscReal *keep)
433: {

439:   PetscUseMethod(pep,"PEPQArnoldiGetRestart_C",(PEP,PetscReal*),(pep,keep));
440:   return(0);
441: }

445: static PetscErrorCode PEPQArnoldiSetLocking_QArnoldi(PEP pep,PetscBool lock)
446: {
447:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

450:   ctx->lock = lock;
451:   return(0);
452: }

456: /*@
457:    PEPQArnoldiSetLocking - Choose between locking and non-locking variants of
458:    the Q-Arnoldi method.

460:    Logically Collective on PEP

462:    Input Parameters:
463: +  pep  - the eigenproblem solver context
464: -  lock - true if the locking variant must be selected

466:    Options Database Key:
467: .  -pep_qarnoldi_locking - Sets the locking flag

469:    Notes:
470:    The default is to keep all directions in the working subspace even if
471:    already converged to working accuracy (the non-locking variant).
472:    This behaviour can be changed so that converged eigenpairs are locked
473:    when the method restarts.

475:    Note that the default behaviour is the opposite to Krylov solvers in EPS.

477:    Level: advanced

479: .seealso: PEPQArnoldiGetLocking()
480: @*/
481: PetscErrorCode PEPQArnoldiSetLocking(PEP pep,PetscBool lock)
482: {

488:   PetscTryMethod(pep,"PEPQArnoldiSetLocking_C",(PEP,PetscBool),(pep,lock));
489:   return(0);
490: }

494: static PetscErrorCode PEPQArnoldiGetLocking_QArnoldi(PEP pep,PetscBool *lock)
495: {
496:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

499:   *lock = ctx->lock;
500:   return(0);
501: }

505: /*@
506:    PEPQArnoldiGetLocking - Gets the locking flag used in the Q-Arnoldi method.

508:    Not Collective

510:    Input Parameter:
511: .  pep - the eigenproblem solver context

513:    Output Parameter:
514: .  lock - the locking flag

516:    Level: advanced

518: .seealso: PEPQArnoldiSetLocking()
519: @*/
520: PetscErrorCode PEPQArnoldiGetLocking(PEP pep,PetscBool *lock)
521: {

527:   PetscUseMethod(pep,"PEPQArnoldiGetLocking_C",(PEP,PetscBool*),(pep,lock));
528:   return(0);
529: }

533: PetscErrorCode PEPSetFromOptions_QArnoldi(PetscOptionItems *PetscOptionsObject,PEP pep)
534: {
536:   PetscBool      flg,lock;
537:   PetscReal      keep;

540:   PetscOptionsHead(PetscOptionsObject,"PEP Q-Arnoldi Options");
541:   PetscOptionsReal("-pep_qarnoldi_restart","Proportion of vectors kept after restart","PEPQArnoldiSetRestart",0.5,&keep,&flg);
542:   if (flg) {
543:     PEPQArnoldiSetRestart(pep,keep);
544:   }
545:   PetscOptionsBool("-pep_qarnoldi_locking","Choose between locking and non-locking variants","PEPQArnoldiSetLocking",PETSC_FALSE,&lock,&flg);
546:   if (flg) {
547:     PEPQArnoldiSetLocking(pep,lock);
548:   }
549:   PetscOptionsTail();
550:   return(0);
551: }

555: PetscErrorCode PEPView_QArnoldi(PEP pep,PetscViewer viewer)
556: {
558:   PEP_QARNOLDI   *ctx = (PEP_QARNOLDI*)pep->data;
559:   PetscBool      isascii;

562:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
563:   if (isascii) {
564:     PetscViewerASCIIPrintf(viewer,"  Q-Arnoldi: %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
565:     PetscViewerASCIIPrintf(viewer,"  Q-Arnoldi: using the %slocking variant\n",ctx->lock?"":"non-");
566:   }
567:   return(0);
568: }

572: PetscErrorCode PEPDestroy_QArnoldi(PEP pep)
573: {

577:   PetscFree(pep->data);
578:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",NULL);
579:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",NULL);
580:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",NULL);
581:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",NULL);
582:   return(0);
583: }

587: PETSC_EXTERN PetscErrorCode PEPCreate_QArnoldi(PEP pep)
588: {
589:   PEP_QARNOLDI   *ctx;

593:   PetscNewLog(pep,&ctx);
594:   pep->data = (void*)ctx;
595:   ctx->lock = PETSC_TRUE;

597:   pep->ops->solve          = PEPSolve_QArnoldi;
598:   pep->ops->setup          = PEPSetUp_QArnoldi;
599:   pep->ops->setfromoptions = PEPSetFromOptions_QArnoldi;
600:   pep->ops->destroy        = PEPDestroy_QArnoldi;
601:   pep->ops->view           = PEPView_QArnoldi;
602:   pep->ops->backtransform  = PEPBackTransform_Default;
603:   pep->ops->computevectors = PEPComputeVectors_Default;
604:   pep->ops->extractvectors = PEPExtractVectors_QArnoldi;
605:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",PEPQArnoldiSetRestart_QArnoldi);
606:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",PEPQArnoldiGetRestart_QArnoldi);
607:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",PEPQArnoldiSetLocking_QArnoldi);
608:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",PEPQArnoldiGetLocking_QArnoldi);
609:   return(0);
610: }