Actual source code: qarnoldi.c
slepc-3.20.1 2023-11-27
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));
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
67: static PetscErrorCode PEPExtractVectors_QArnoldi(PEP pep)
68: {
69: PetscInt k=pep->nconv;
70: Mat X,X0;
72: PetscFunctionBegin;
73: if (pep->nconv==0) PetscFunctionReturn(PETSC_SUCCESS);
74: PetscCall(DSVectors(pep->ds,DS_MAT_X,NULL,NULL));
76: /* update vectors V = V*X */
77: PetscCall(DSGetMat(pep->ds,DS_MAT_X,&X));
78: PetscCall(MatDenseGetSubMatrix(X,0,k,0,k,&X0));
79: PetscCall(BVMultInPlace(pep->V,X0,0,k));
80: PetscCall(MatDenseRestoreSubMatrix(X,&X0));
81: PetscCall(DSRestoreMat(pep->ds,DS_MAT_X,&X));
82: PetscFunctionReturn(PETSC_SUCCESS);
83: }
85: /*
86: Compute a step of Classical Gram-Schmidt orthogonalization
87: */
88: 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)
89: {
90: PetscBLASInt ione = 1,j_1 = j+1;
91: PetscReal x,y;
92: PetscScalar dot,one = 1.0,zero = 0.0;
94: PetscFunctionBegin;
95: /* compute norm of v and w */
96: if (onorm) {
97: PetscCall(VecNorm(v,NORM_2,&x));
98: PetscCall(VecNorm(w,NORM_2,&y));
99: *onorm = SlepcAbs(x,y);
100: }
102: /* orthogonalize: compute h */
103: PetscCall(BVDotVec(V,v,h));
104: PetscCall(BVDotVec(V,w,work));
105: if (j>0) PetscCallBLAS("BLASgemv",BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione));
106: PetscCall(VecDot(w,t,&dot));
107: h[j] += dot;
109: /* orthogonalize: update v and w */
110: PetscCall(BVMultVec(V,-1.0,1.0,v,h));
111: if (j>0) {
112: PetscCallBLAS("BLASgemv",BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione));
113: PetscCall(BVMultVec(V,-1.0,1.0,w,work));
114: }
115: PetscCall(VecAXPY(w,-h[j],t));
117: /* compute norm of v and w */
118: if (norm) {
119: PetscCall(VecNorm(v,NORM_2,&x));
120: PetscCall(VecNorm(w,NORM_2,&y));
121: *norm = SlepcAbs(x,y);
122: }
123: PetscFunctionReturn(PETSC_SUCCESS);
124: }
126: /*
127: Compute a run of Q-Arnoldi iterations
128: */
129: static PetscErrorCode PEPQArnoldi(PEP pep,Mat A,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
130: {
131: PetscInt i,j,l,m = *M,ldh;
132: Vec t = pep->work[2],u = pep->work[3];
133: BVOrthogRefineType refinement;
134: PetscReal norm=0.0,onorm,eta;
135: PetscScalar *H,*c = work + m;
137: PetscFunctionBegin;
138: *beta = 0.0;
139: PetscCall(MatDenseGetArray(A,&H));
140: PetscCall(MatDenseGetLDA(A,&ldh));
141: PetscCall(BVGetOrthogonalization(pep->V,NULL,&refinement,&eta,NULL));
142: PetscCall(BVInsertVec(pep->V,k,v));
143: for (j=k;j<m;j++) {
144: /* apply operator */
145: PetscCall(VecCopy(w,t));
146: if (pep->Dr) PetscCall(VecPointwiseMult(v,v,pep->Dr));
147: PetscCall(STMatMult(pep->st,0,v,u));
148: PetscCall(VecCopy(t,v));
149: if (pep->Dr) PetscCall(VecPointwiseMult(t,t,pep->Dr));
150: PetscCall(STMatMult(pep->st,1,t,w));
151: PetscCall(VecAXPY(u,pep->sfactor,w));
152: PetscCall(STMatSolve(pep->st,u,w));
153: PetscCall(VecScale(w,-1.0/(pep->sfactor*pep->sfactor)));
154: if (pep->Dr) PetscCall(VecPointwiseDivide(w,w,pep->Dr));
155: PetscCall(VecCopy(v,t));
156: PetscCall(BVSetActiveColumns(pep->V,0,j+1));
158: /* orthogonalize */
159: switch (refinement) {
160: case BV_ORTHOG_REFINE_NEVER:
161: PetscCall(PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,&norm,work));
162: *breakdown = PETSC_FALSE;
163: break;
164: case BV_ORTHOG_REFINE_ALWAYS:
165: PetscCall(PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,NULL,work));
166: PetscCall(PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,&onorm,&norm,work));
167: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
168: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
169: else *breakdown = PETSC_FALSE;
170: break;
171: case BV_ORTHOG_REFINE_IFNEEDED:
172: PetscCall(PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,&onorm,&norm,work));
173: /* ||q|| < eta ||h|| */
174: l = 1;
175: while (l<3 && norm < eta * onorm) {
176: l++;
177: onorm = norm;
178: PetscCall(PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,NULL,&norm,work));
179: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
180: }
181: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
182: else *breakdown = PETSC_FALSE;
183: break;
184: }
185: PetscCall(VecScale(v,1.0/norm));
186: PetscCall(VecScale(w,1.0/norm));
188: H[j+1+ldh*j] = norm;
189: if (j<m-1) PetscCall(BVInsertVec(pep->V,j+1,v));
190: }
191: *beta = norm;
192: PetscCall(MatDenseRestoreArray(A,&H));
193: PetscFunctionReturn(PETSC_SUCCESS);
194: }
196: static PetscErrorCode PEPSolve_QArnoldi(PEP pep)
197: {
198: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
199: PetscInt j,k,l,lwork,nv,nconv;
200: Vec v=pep->work[0],w=pep->work[1];
201: Mat Q,S;
202: PetscScalar *work;
203: PetscReal beta,norm,x,y;
204: PetscBool breakdown=PETSC_FALSE,sinv;
206: PetscFunctionBegin;
207: lwork = 7*pep->ncv;
208: PetscCall(PetscMalloc1(lwork,&work));
209: PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv));
210: PetscCall(RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor));
211: PetscCall(STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor));
213: /* Get the starting Arnoldi vector */
214: for (j=0;j<2;j++) {
215: if (j>=pep->nini) PetscCall(BVSetRandomColumn(pep->V,j));
216: }
217: PetscCall(BVCopyVec(pep->V,0,v));
218: PetscCall(BVCopyVec(pep->V,1,w));
219: PetscCall(VecNorm(v,NORM_2,&x));
220: PetscCall(VecNorm(w,NORM_2,&y));
221: norm = SlepcAbs(x,y);
222: PetscCall(VecScale(v,1.0/norm));
223: PetscCall(VecScale(w,1.0/norm));
225: /* clean projected matrix (including the extra-arrow) */
226: PetscCall(DSSetDimensions(pep->ds,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
227: PetscCall(DSGetMat(pep->ds,DS_MAT_A,&S));
228: PetscCall(MatZeroEntries(S));
229: PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&S));
231: /* Restart loop */
232: l = 0;
233: while (pep->reason == PEP_CONVERGED_ITERATING) {
234: pep->its++;
236: /* Compute an nv-step Arnoldi factorization */
237: nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
238: PetscCall(DSGetMat(pep->ds,DS_MAT_A,&S));
239: PetscCall(PEPQArnoldi(pep,S,pep->nconv+l,&nv,v,w,&beta,&breakdown,work));
240: PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&S));
241: PetscCall(DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l));
242: PetscCall(DSSetState(pep->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE));
243: PetscCall(BVSetActiveColumns(pep->V,pep->nconv,nv));
245: /* Solve projected problem */
246: PetscCall(DSSolve(pep->ds,pep->eigr,pep->eigi));
247: PetscCall(DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL));
248: PetscCall(DSUpdateExtraRow(pep->ds));
249: PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi));
251: /* Check convergence */
252: PetscCall(PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k));
253: PetscCall((*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx));
254: nconv = k;
256: /* Update l */
257: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
258: else {
259: l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
260: PetscCall(DSGetTruncateSize(pep->ds,k,nv,&l));
261: }
262: if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
263: if (l) PetscCall(PetscInfo(pep,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l));
265: if (pep->reason == PEP_CONVERGED_ITERATING) {
266: if (PetscUnlikely(breakdown)) {
267: /* Stop if breakdown */
268: PetscCall(PetscInfo(pep,"Breakdown Quadratic Arnoldi method (it=%" PetscInt_FMT " norm=%g)\n",pep->its,(double)beta));
269: pep->reason = PEP_DIVERGED_BREAKDOWN;
270: } else {
271: /* Prepare the Rayleigh quotient for restart */
272: PetscCall(DSTruncate(pep->ds,k+l,PETSC_FALSE));
273: }
274: }
275: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
276: PetscCall(DSGetMat(pep->ds,DS_MAT_Q,&Q));
277: PetscCall(BVMultInPlace(pep->V,Q,pep->nconv,k+l));
278: PetscCall(DSRestoreMat(pep->ds,DS_MAT_Q,&Q));
280: pep->nconv = k;
281: PetscCall(PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv));
282: }
283: PetscCall(BVSetActiveColumns(pep->V,0,pep->nconv));
284: for (j=0;j<pep->nconv;j++) {
285: pep->eigr[j] *= pep->sfactor;
286: pep->eigi[j] *= pep->sfactor;
287: }
289: PetscCall(STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor));
290: PetscCall(RGPopScale(pep->rg));
292: PetscCall(DSTruncate(pep->ds,pep->nconv,PETSC_TRUE));
293: PetscCall(PetscFree(work));
294: PetscFunctionReturn(PETSC_SUCCESS);
295: }
297: static PetscErrorCode PEPQArnoldiSetRestart_QArnoldi(PEP pep,PetscReal keep)
298: {
299: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
301: PetscFunctionBegin;
302: if (keep==(PetscReal)PETSC_DEFAULT) ctx->keep = 0.5;
303: else {
304: 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]");
305: ctx->keep = keep;
306: }
307: PetscFunctionReturn(PETSC_SUCCESS);
308: }
310: /*@
311: PEPQArnoldiSetRestart - Sets the restart parameter for the Q-Arnoldi
312: method, in particular the proportion of basis vectors that must be kept
313: after restart.
315: Logically Collective
317: Input Parameters:
318: + pep - the eigenproblem solver context
319: - keep - the number of vectors to be kept at restart
321: Options Database Key:
322: . -pep_qarnoldi_restart - Sets the restart parameter
324: Notes:
325: Allowed values are in the range [0.1,0.9]. The default is 0.5.
327: Level: advanced
329: .seealso: PEPQArnoldiGetRestart()
330: @*/
331: PetscErrorCode PEPQArnoldiSetRestart(PEP pep,PetscReal keep)
332: {
333: PetscFunctionBegin;
336: PetscTryMethod(pep,"PEPQArnoldiSetRestart_C",(PEP,PetscReal),(pep,keep));
337: PetscFunctionReturn(PETSC_SUCCESS);
338: }
340: static PetscErrorCode PEPQArnoldiGetRestart_QArnoldi(PEP pep,PetscReal *keep)
341: {
342: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
344: PetscFunctionBegin;
345: *keep = ctx->keep;
346: PetscFunctionReturn(PETSC_SUCCESS);
347: }
349: /*@
350: PEPQArnoldiGetRestart - Gets the restart parameter used in the Q-Arnoldi method.
352: Not Collective
354: Input Parameter:
355: . pep - the eigenproblem solver context
357: Output Parameter:
358: . keep - the restart parameter
360: Level: advanced
362: .seealso: PEPQArnoldiSetRestart()
363: @*/
364: PetscErrorCode PEPQArnoldiGetRestart(PEP pep,PetscReal *keep)
365: {
366: PetscFunctionBegin;
368: PetscAssertPointer(keep,2);
369: PetscUseMethod(pep,"PEPQArnoldiGetRestart_C",(PEP,PetscReal*),(pep,keep));
370: PetscFunctionReturn(PETSC_SUCCESS);
371: }
373: static PetscErrorCode PEPQArnoldiSetLocking_QArnoldi(PEP pep,PetscBool lock)
374: {
375: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
377: PetscFunctionBegin;
378: ctx->lock = lock;
379: PetscFunctionReturn(PETSC_SUCCESS);
380: }
382: /*@
383: PEPQArnoldiSetLocking - Choose between locking and non-locking variants of
384: the Q-Arnoldi method.
386: Logically Collective
388: Input Parameters:
389: + pep - the eigenproblem solver context
390: - lock - true if the locking variant must be selected
392: Options Database Key:
393: . -pep_qarnoldi_locking - Sets the locking flag
395: Notes:
396: The default is to lock converged eigenpairs when the method restarts.
397: This behaviour can be changed so that all directions are kept in the
398: working subspace even if already converged to working accuracy (the
399: non-locking variant).
401: Level: advanced
403: .seealso: PEPQArnoldiGetLocking()
404: @*/
405: PetscErrorCode PEPQArnoldiSetLocking(PEP pep,PetscBool lock)
406: {
407: PetscFunctionBegin;
410: PetscTryMethod(pep,"PEPQArnoldiSetLocking_C",(PEP,PetscBool),(pep,lock));
411: PetscFunctionReturn(PETSC_SUCCESS);
412: }
414: static PetscErrorCode PEPQArnoldiGetLocking_QArnoldi(PEP pep,PetscBool *lock)
415: {
416: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
418: PetscFunctionBegin;
419: *lock = ctx->lock;
420: PetscFunctionReturn(PETSC_SUCCESS);
421: }
423: /*@
424: PEPQArnoldiGetLocking - Gets the locking flag used in the Q-Arnoldi method.
426: Not Collective
428: Input Parameter:
429: . pep - the eigenproblem solver context
431: Output Parameter:
432: . lock - the locking flag
434: Level: advanced
436: .seealso: PEPQArnoldiSetLocking()
437: @*/
438: PetscErrorCode PEPQArnoldiGetLocking(PEP pep,PetscBool *lock)
439: {
440: PetscFunctionBegin;
442: PetscAssertPointer(lock,2);
443: PetscUseMethod(pep,"PEPQArnoldiGetLocking_C",(PEP,PetscBool*),(pep,lock));
444: PetscFunctionReturn(PETSC_SUCCESS);
445: }
447: static PetscErrorCode PEPSetFromOptions_QArnoldi(PEP pep,PetscOptionItems *PetscOptionsObject)
448: {
449: PetscBool flg,lock;
450: PetscReal keep;
452: PetscFunctionBegin;
453: PetscOptionsHeadBegin(PetscOptionsObject,"PEP Q-Arnoldi Options");
455: PetscCall(PetscOptionsReal("-pep_qarnoldi_restart","Proportion of vectors kept after restart","PEPQArnoldiSetRestart",0.5,&keep,&flg));
456: if (flg) PetscCall(PEPQArnoldiSetRestart(pep,keep));
458: PetscCall(PetscOptionsBool("-pep_qarnoldi_locking","Choose between locking and non-locking variants","PEPQArnoldiSetLocking",PETSC_FALSE,&lock,&flg));
459: if (flg) PetscCall(PEPQArnoldiSetLocking(pep,lock));
461: PetscOptionsHeadEnd();
462: PetscFunctionReturn(PETSC_SUCCESS);
463: }
465: static PetscErrorCode PEPView_QArnoldi(PEP pep,PetscViewer viewer)
466: {
467: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
468: PetscBool isascii;
470: PetscFunctionBegin;
471: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
472: if (isascii) {
473: PetscCall(PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep)));
474: PetscCall(PetscViewerASCIIPrintf(viewer," using the %slocking variant\n",ctx->lock?"":"non-"));
475: }
476: PetscFunctionReturn(PETSC_SUCCESS);
477: }
479: static PetscErrorCode PEPDestroy_QArnoldi(PEP pep)
480: {
481: PetscFunctionBegin;
482: PetscCall(PetscFree(pep->data));
483: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",NULL));
484: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",NULL));
485: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",NULL));
486: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",NULL));
487: PetscFunctionReturn(PETSC_SUCCESS);
488: }
490: SLEPC_EXTERN PetscErrorCode PEPCreate_QArnoldi(PEP pep)
491: {
492: PEP_QARNOLDI *ctx;
494: PetscFunctionBegin;
495: PetscCall(PetscNew(&ctx));
496: pep->data = (void*)ctx;
498: pep->lineariz = PETSC_TRUE;
499: ctx->lock = PETSC_TRUE;
501: pep->ops->solve = PEPSolve_QArnoldi;
502: pep->ops->setup = PEPSetUp_QArnoldi;
503: pep->ops->setfromoptions = PEPSetFromOptions_QArnoldi;
504: pep->ops->destroy = PEPDestroy_QArnoldi;
505: pep->ops->view = PEPView_QArnoldi;
506: pep->ops->backtransform = PEPBackTransform_Default;
507: pep->ops->computevectors = PEPComputeVectors_Default;
508: pep->ops->extractvectors = PEPExtractVectors_QArnoldi;
509: pep->ops->setdefaultst = PEPSetDefaultST_Transform;
511: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",PEPQArnoldiSetRestart_QArnoldi));
512: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",PEPQArnoldiGetRestart_QArnoldi));
513: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",PEPQArnoldiSetLocking_QArnoldi));
514: PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",PEPQArnoldiGetLocking_QArnoldi));
515: PetscFunctionReturn(PETSC_SUCCESS);
516: }