Actual source code: qarnoldi.c

slepc-3.21.0 2024-03-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc quadratic eigensolver: "qarnoldi"

 13:    Method: Q-Arnoldi

 15:    Algorithm:

 17:        Quadratic Arnoldi with Krylov-Schur type restart.

 19:    References:

 21:        [1] K. Meerbergen, "The Quadratic Arnoldi method for the solution
 22:            of the quadratic eigenvalue problem", SIAM J. Matrix Anal.
 23:            Appl. 30(4):1462-1482, 2008.
 24: */

 26: #include <slepc/private/pepimpl.h>
 27: #include <petscblaslapack.h>

 29: typedef struct {
 30:   PetscReal keep;         /* restart parameter */
 31:   PetscBool lock;         /* locking/non-locking variant */
 32: } PEP_QARNOLDI;

 34: static PetscErrorCode PEPSetUp_QArnoldi(PEP pep)
 35: {
 36:   PEP_QARNOLDI   *ctx = (PEP_QARNOLDI*)pep->data;
 37:   PetscBool      flg;

 39:   PetscFunctionBegin;
 40:   PEPCheckQuadratic(pep);
 41:   PEPCheckShiftSinvert(pep);
 42:   PetscCall(PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd));
 43:   PetscCheck(ctx->lock || pep->mpd>=pep->ncv,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
 44:   if (pep->max_it==PETSC_DEFAULT) pep->max_it = PetscMax(100,4*pep->n/pep->ncv);
 45:   if (!pep->which) PetscCall(PEPSetWhichEigenpairs_Default(pep));
 46:   PetscCheck(pep->which!=PEP_ALL,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues");

 48:   PetscCall(STGetTransform(pep->st,&flg));
 49:   PetscCheck(flg,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag set, see STSetTransform()");

 51:   /* set default extraction */
 52:   if (!pep->extract) pep->extract = PEP_EXTRACT_NONE;
 53:   PEPCheckUnsupported(pep,PEP_FEATURE_NONMONOMIAL | PEP_FEATURE_EXTRACT);

 55:   if (!ctx->keep) ctx->keep = 0.5;

 57:   PetscCall(PEPAllocateSolution(pep,0));
 58:   PetscCall(PEPSetWorkVecs(pep,4));

 60:   PetscCall(DSSetType(pep->ds,DSNHEP));
 61:   PetscCall(DSSetExtraRow(pep->ds,PETSC_TRUE));
 62:   PetscCall(DSAllocate(pep->ds,pep->ncv+1));
 63:   PetscFunctionReturn(PETSC_SUCCESS);
 64: }

 66: static PetscErrorCode PEPExtractVectors_QArnoldi(PEP pep)
 67: {
 68:   PetscInt       k=pep->nconv;
 69:   Mat            X,X0;

 71:   PetscFunctionBegin;
 72:   if (pep->nconv==0) PetscFunctionReturn(PETSC_SUCCESS);
 73:   PetscCall(DSVectors(pep->ds,DS_MAT_X,NULL,NULL));

 75:   /* update vectors V = V*X */
 76:   PetscCall(DSGetMat(pep->ds,DS_MAT_X,&X));
 77:   PetscCall(MatDenseGetSubMatrix(X,0,k,0,k,&X0));
 78:   PetscCall(BVMultInPlace(pep->V,X0,0,k));
 79:   PetscCall(MatDenseRestoreSubMatrix(X,&X0));
 80:   PetscCall(DSRestoreMat(pep->ds,DS_MAT_X,&X));
 81:   PetscFunctionReturn(PETSC_SUCCESS);
 82: }

 84: /*
 85:   Compute a step of Classical Gram-Schmidt orthogonalization
 86: */
 87: 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)
 88: {
 89:   PetscBLASInt   ione = 1,j_1 = j+1;
 90:   PetscReal      x,y;
 91:   PetscScalar    dot,one = 1.0,zero = 0.0;

 93:   PetscFunctionBegin;
 94:   /* compute norm of v and w */
 95:   if (onorm) {
 96:     PetscCall(VecNorm(v,NORM_2,&x));
 97:     PetscCall(VecNorm(w,NORM_2,&y));
 98:     *onorm = SlepcAbs(x,y);
 99:   }

101:   /* orthogonalize: compute h */
102:   PetscCall(BVDotVec(V,v,h));
103:   PetscCall(BVDotVec(V,w,work));
104:   if (j>0) PetscCallBLAS("BLASgemv",BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione));
105:   PetscCall(VecDot(w,t,&dot));
106:   h[j] += dot;

108:   /* orthogonalize: update v and w */
109:   PetscCall(BVMultVec(V,-1.0,1.0,v,h));
110:   if (j>0) {
111:     PetscCallBLAS("BLASgemv",BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione));
112:     PetscCall(BVMultVec(V,-1.0,1.0,w,work));
113:   }
114:   PetscCall(VecAXPY(w,-h[j],t));

116:   /* compute norm of v and w */
117:   if (norm) {
118:     PetscCall(VecNorm(v,NORM_2,&x));
119:     PetscCall(VecNorm(w,NORM_2,&y));
120:     *norm = SlepcAbs(x,y);
121:   }
122:   PetscFunctionReturn(PETSC_SUCCESS);
123: }

125: /*
126:   Compute a run of Q-Arnoldi iterations
127: */
128: static PetscErrorCode PEPQArnoldi(PEP pep,Mat A,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
129: {
130:   PetscInt           i,j,l,m = *M,ldh;
131:   Vec                t = pep->work[2],u = pep->work[3];
132:   BVOrthogRefineType refinement;
133:   PetscReal          norm=0.0,onorm,eta;
134:   PetscScalar        *H,*c = work + m;

136:   PetscFunctionBegin;
137:   *beta = 0.0;
138:   PetscCall(MatDenseGetArray(A,&H));
139:   PetscCall(MatDenseGetLDA(A,&ldh));
140:   PetscCall(BVGetOrthogonalization(pep->V,NULL,&refinement,&eta,NULL));
141:   PetscCall(BVInsertVec(pep->V,k,v));
142:   for (j=k;j<m;j++) {
143:     /* apply operator */
144:     PetscCall(VecCopy(w,t));
145:     if (pep->Dr) PetscCall(VecPointwiseMult(v,v,pep->Dr));
146:     PetscCall(STMatMult(pep->st,0,v,u));
147:     PetscCall(VecCopy(t,v));
148:     if (pep->Dr) PetscCall(VecPointwiseMult(t,t,pep->Dr));
149:     PetscCall(STMatMult(pep->st,1,t,w));
150:     PetscCall(VecAXPY(u,pep->sfactor,w));
151:     PetscCall(STMatSolve(pep->st,u,w));
152:     PetscCall(VecScale(w,-1.0/(pep->sfactor*pep->sfactor)));
153:     if (pep->Dr) PetscCall(VecPointwiseDivide(w,w,pep->Dr));
154:     PetscCall(VecCopy(v,t));
155:     PetscCall(BVSetActiveColumns(pep->V,0,j+1));

157:     /* orthogonalize */
158:     switch (refinement) {
159:       case BV_ORTHOG_REFINE_NEVER:
160:         PetscCall(PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,&norm,work));
161:         *breakdown = PETSC_FALSE;
162:         break;
163:       case BV_ORTHOG_REFINE_ALWAYS:
164:         PetscCall(PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,NULL,work));
165:         PetscCall(PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,&onorm,&norm,work));
166:         for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
167:         if (norm < eta * onorm) *breakdown = PETSC_TRUE;
168:         else *breakdown = PETSC_FALSE;
169:         break;
170:       case BV_ORTHOG_REFINE_IFNEEDED:
171:         PetscCall(PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,&onorm,&norm,work));
172:         /* ||q|| < eta ||h|| */
173:         l = 1;
174:         while (l<3 && norm < eta * onorm) {
175:           l++;
176:           onorm = norm;
177:           PetscCall(PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,NULL,&norm,work));
178:           for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
179:         }
180:         if (norm < eta * onorm) *breakdown = PETSC_TRUE;
181:         else *breakdown = PETSC_FALSE;
182:         break;
183:     }
184:     PetscCall(VecScale(v,1.0/norm));
185:     PetscCall(VecScale(w,1.0/norm));

187:     H[j+1+ldh*j] = norm;
188:     if (j<m-1) PetscCall(BVInsertVec(pep->V,j+1,v));
189:   }
190:   *beta = norm;
191:   PetscCall(MatDenseRestoreArray(A,&H));
192:   PetscFunctionReturn(PETSC_SUCCESS);
193: }

195: static PetscErrorCode PEPSolve_QArnoldi(PEP pep)
196: {
197:   PEP_QARNOLDI   *ctx = (PEP_QARNOLDI*)pep->data;
198:   PetscInt       j,k,l,lwork,nv,nconv;
199:   Vec            v=pep->work[0],w=pep->work[1];
200:   Mat            Q,S;
201:   PetscScalar    *work;
202:   PetscReal      beta,norm,x,y;
203:   PetscBool      breakdown=PETSC_FALSE,sinv;

205:   PetscFunctionBegin;
206:   lwork = 7*pep->ncv;
207:   PetscCall(PetscMalloc1(lwork,&work));
208:   PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv));
209:   PetscCall(RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor));
210:   PetscCall(STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor));

212:   /* Get the starting Arnoldi vector */
213:   for (j=0;j<2;j++) {
214:     if (j>=pep->nini) PetscCall(BVSetRandomColumn(pep->V,j));
215:   }
216:   PetscCall(BVCopyVec(pep->V,0,v));
217:   PetscCall(BVCopyVec(pep->V,1,w));
218:   PetscCall(VecNorm(v,NORM_2,&x));
219:   PetscCall(VecNorm(w,NORM_2,&y));
220:   norm = SlepcAbs(x,y);
221:   PetscCall(VecScale(v,1.0/norm));
222:   PetscCall(VecScale(w,1.0/norm));

224:   /* clean projected matrix (including the extra-arrow) */
225:   PetscCall(DSSetDimensions(pep->ds,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
226:   PetscCall(DSGetMat(pep->ds,DS_MAT_A,&S));
227:   PetscCall(MatZeroEntries(S));
228:   PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&S));

230:    /* Restart loop */
231:   l = 0;
232:   while (pep->reason == PEP_CONVERGED_ITERATING) {
233:     pep->its++;

235:     /* Compute an nv-step Arnoldi factorization */
236:     nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
237:     PetscCall(DSGetMat(pep->ds,DS_MAT_A,&S));
238:     PetscCall(PEPQArnoldi(pep,S,pep->nconv+l,&nv,v,w,&beta,&breakdown,work));
239:     PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&S));
240:     PetscCall(DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l));
241:     PetscCall(DSSetState(pep->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE));
242:     PetscCall(BVSetActiveColumns(pep->V,pep->nconv,nv));

244:     /* Solve projected problem */
245:     PetscCall(DSSolve(pep->ds,pep->eigr,pep->eigi));
246:     PetscCall(DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL));
247:     PetscCall(DSUpdateExtraRow(pep->ds));
248:     PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi));

250:     /* Check convergence */
251:     PetscCall(PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k));
252:     PetscCall((*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx));
253:     nconv = k;

255:     /* Update l */
256:     if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
257:     else {
258:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
259:       PetscCall(DSGetTruncateSize(pep->ds,k,nv,&l));
260:     }
261:     if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
262:     if (l) PetscCall(PetscInfo(pep,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l));

264:     if (pep->reason == PEP_CONVERGED_ITERATING) {
265:       if (PetscUnlikely(breakdown)) {
266:         /* Stop if breakdown */
267:         PetscCall(PetscInfo(pep,"Breakdown Quadratic Arnoldi method (it=%" PetscInt_FMT " norm=%g)\n",pep->its,(double)beta));
268:         pep->reason = PEP_DIVERGED_BREAKDOWN;
269:       } else {
270:         /* Prepare the Rayleigh quotient for restart */
271:         PetscCall(DSTruncate(pep->ds,k+l,PETSC_FALSE));
272:       }
273:     }
274:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
275:     PetscCall(DSGetMat(pep->ds,DS_MAT_Q,&Q));
276:     PetscCall(BVMultInPlace(pep->V,Q,pep->nconv,k+l));
277:     PetscCall(DSRestoreMat(pep->ds,DS_MAT_Q,&Q));

279:     pep->nconv = k;
280:     PetscCall(PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv));
281:   }
282:   PetscCall(BVSetActiveColumns(pep->V,0,pep->nconv));
283:   for (j=0;j<pep->nconv;j++) {
284:     pep->eigr[j] *= pep->sfactor;
285:     pep->eigi[j] *= pep->sfactor;
286:   }

288:   PetscCall(STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor));
289:   PetscCall(RGPopScale(pep->rg));

291:   PetscCall(DSTruncate(pep->ds,pep->nconv,PETSC_TRUE));
292:   PetscCall(PetscFree(work));
293:   PetscFunctionReturn(PETSC_SUCCESS);
294: }

296: static PetscErrorCode PEPQArnoldiSetRestart_QArnoldi(PEP pep,PetscReal keep)
297: {
298:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

300:   PetscFunctionBegin;
301:   if (keep==(PetscReal)PETSC_DEFAULT) ctx->keep = 0.5;
302:   else {
303:     PetscCheck(keep>=0.1 && keep<=0.9,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
304:     ctx->keep = keep;
305:   }
306:   PetscFunctionReturn(PETSC_SUCCESS);
307: }

309: /*@
310:    PEPQArnoldiSetRestart - Sets the restart parameter for the Q-Arnoldi
311:    method, in particular the proportion of basis vectors that must be kept
312:    after restart.

314:    Logically Collective

316:    Input Parameters:
317: +  pep  - the eigenproblem solver context
318: -  keep - the number of vectors to be kept at restart

320:    Options Database Key:
321: .  -pep_qarnoldi_restart - Sets the restart parameter

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

326:    Level: advanced

328: .seealso: PEPQArnoldiGetRestart()
329: @*/
330: PetscErrorCode PEPQArnoldiSetRestart(PEP pep,PetscReal keep)
331: {
332:   PetscFunctionBegin;
335:   PetscTryMethod(pep,"PEPQArnoldiSetRestart_C",(PEP,PetscReal),(pep,keep));
336:   PetscFunctionReturn(PETSC_SUCCESS);
337: }

339: static PetscErrorCode PEPQArnoldiGetRestart_QArnoldi(PEP pep,PetscReal *keep)
340: {
341:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

343:   PetscFunctionBegin;
344:   *keep = ctx->keep;
345:   PetscFunctionReturn(PETSC_SUCCESS);
346: }

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

351:    Not Collective

353:    Input Parameter:
354: .  pep - the eigenproblem solver context

356:    Output Parameter:
357: .  keep - the restart parameter

359:    Level: advanced

361: .seealso: PEPQArnoldiSetRestart()
362: @*/
363: PetscErrorCode PEPQArnoldiGetRestart(PEP pep,PetscReal *keep)
364: {
365:   PetscFunctionBegin;
367:   PetscAssertPointer(keep,2);
368:   PetscUseMethod(pep,"PEPQArnoldiGetRestart_C",(PEP,PetscReal*),(pep,keep));
369:   PetscFunctionReturn(PETSC_SUCCESS);
370: }

372: static PetscErrorCode PEPQArnoldiSetLocking_QArnoldi(PEP pep,PetscBool lock)
373: {
374:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

376:   PetscFunctionBegin;
377:   ctx->lock = lock;
378:   PetscFunctionReturn(PETSC_SUCCESS);
379: }

381: /*@
382:    PEPQArnoldiSetLocking - Choose between locking and non-locking variants of
383:    the Q-Arnoldi method.

385:    Logically Collective

387:    Input Parameters:
388: +  pep  - the eigenproblem solver context
389: -  lock - true if the locking variant must be selected

391:    Options Database Key:
392: .  -pep_qarnoldi_locking - Sets the locking flag

394:    Notes:
395:    The default is to lock converged eigenpairs when the method restarts.
396:    This behaviour can be changed so that all directions are kept in the
397:    working subspace even if already converged to working accuracy (the
398:    non-locking variant).

400:    Level: advanced

402: .seealso: PEPQArnoldiGetLocking()
403: @*/
404: PetscErrorCode PEPQArnoldiSetLocking(PEP pep,PetscBool lock)
405: {
406:   PetscFunctionBegin;
409:   PetscTryMethod(pep,"PEPQArnoldiSetLocking_C",(PEP,PetscBool),(pep,lock));
410:   PetscFunctionReturn(PETSC_SUCCESS);
411: }

413: static PetscErrorCode PEPQArnoldiGetLocking_QArnoldi(PEP pep,PetscBool *lock)
414: {
415:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

417:   PetscFunctionBegin;
418:   *lock = ctx->lock;
419:   PetscFunctionReturn(PETSC_SUCCESS);
420: }

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

425:    Not Collective

427:    Input Parameter:
428: .  pep - the eigenproblem solver context

430:    Output Parameter:
431: .  lock - the locking flag

433:    Level: advanced

435: .seealso: PEPQArnoldiSetLocking()
436: @*/
437: PetscErrorCode PEPQArnoldiGetLocking(PEP pep,PetscBool *lock)
438: {
439:   PetscFunctionBegin;
441:   PetscAssertPointer(lock,2);
442:   PetscUseMethod(pep,"PEPQArnoldiGetLocking_C",(PEP,PetscBool*),(pep,lock));
443:   PetscFunctionReturn(PETSC_SUCCESS);
444: }

446: static PetscErrorCode PEPSetFromOptions_QArnoldi(PEP pep,PetscOptionItems *PetscOptionsObject)
447: {
448:   PetscBool      flg,lock;
449:   PetscReal      keep;

451:   PetscFunctionBegin;
452:   PetscOptionsHeadBegin(PetscOptionsObject,"PEP Q-Arnoldi Options");

454:     PetscCall(PetscOptionsReal("-pep_qarnoldi_restart","Proportion of vectors kept after restart","PEPQArnoldiSetRestart",0.5,&keep,&flg));
455:     if (flg) PetscCall(PEPQArnoldiSetRestart(pep,keep));

457:     PetscCall(PetscOptionsBool("-pep_qarnoldi_locking","Choose between locking and non-locking variants","PEPQArnoldiSetLocking",PETSC_FALSE,&lock,&flg));
458:     if (flg) PetscCall(PEPQArnoldiSetLocking(pep,lock));

460:   PetscOptionsHeadEnd();
461:   PetscFunctionReturn(PETSC_SUCCESS);
462: }

464: static PetscErrorCode PEPView_QArnoldi(PEP pep,PetscViewer viewer)
465: {
466:   PEP_QARNOLDI   *ctx = (PEP_QARNOLDI*)pep->data;
467:   PetscBool      isascii;

469:   PetscFunctionBegin;
470:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
471:   if (isascii) {
472:     PetscCall(PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep)));
473:     PetscCall(PetscViewerASCIIPrintf(viewer,"  using the %slocking variant\n",ctx->lock?"":"non-"));
474:   }
475:   PetscFunctionReturn(PETSC_SUCCESS);
476: }

478: static PetscErrorCode PEPDestroy_QArnoldi(PEP pep)
479: {
480:   PetscFunctionBegin;
481:   PetscCall(PetscFree(pep->data));
482:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",NULL));
483:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",NULL));
484:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",NULL));
485:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",NULL));
486:   PetscFunctionReturn(PETSC_SUCCESS);
487: }

489: SLEPC_EXTERN PetscErrorCode PEPCreate_QArnoldi(PEP pep)
490: {
491:   PEP_QARNOLDI   *ctx;

493:   PetscFunctionBegin;
494:   PetscCall(PetscNew(&ctx));
495:   pep->data = (void*)ctx;

497:   pep->lineariz = PETSC_TRUE;
498:   ctx->lock     = PETSC_TRUE;

500:   pep->ops->solve          = PEPSolve_QArnoldi;
501:   pep->ops->setup          = PEPSetUp_QArnoldi;
502:   pep->ops->setfromoptions = PEPSetFromOptions_QArnoldi;
503:   pep->ops->destroy        = PEPDestroy_QArnoldi;
504:   pep->ops->view           = PEPView_QArnoldi;
505:   pep->ops->backtransform  = PEPBackTransform_Default;
506:   pep->ops->computevectors = PEPComputeVectors_Default;
507:   pep->ops->extractvectors = PEPExtractVectors_QArnoldi;
508:   pep->ops->setdefaultst   = PEPSetDefaultST_Transform;

510:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",PEPQArnoldiSetRestart_QArnoldi));
511:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",PEPQArnoldiGetRestart_QArnoldi));
512:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",PEPQArnoldiSetLocking_QArnoldi));
513:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",PEPQArnoldiGetLocking_QArnoldi));
514:   PetscFunctionReturn(PETSC_SUCCESS);
515: }