Actual source code: qarnoldi.c

slepc-3.6.1 2015-09-03
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-2015, 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      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:   if (!pep->which) {
 59:     PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
 60:     if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
 61:     else pep->which = PEP_LARGEST_MAGNITUDE;
 62:   }

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

 69:   /* set default extraction */
 70:   if (!pep->extract) {
 71:     pep->extract = PEP_EXTRACT_NONE;
 72:   }
 73:   if (pep->extract!=PEP_EXTRACT_NONE) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver does not support requested extraction");
 74:  
 75:   if (!ctx->keep) ctx->keep = 0.5;

 77:   PEPAllocateSolution(pep,0);
 78:   PEPSetWorkVecs(pep,4);

 80:   DSSetType(pep->ds,DSNHEP);
 81:   DSSetExtraRow(pep->ds,PETSC_TRUE);
 82:   DSAllocate(pep->ds,pep->ncv+1);

 84:   /* process starting vector */
 85:   if (pep->nini>-2) {
 86:     BVSetRandomColumn(pep->V,0,pep->rand);
 87:     BVSetRandomColumn(pep->V,1,pep->rand);
 88:   } else {
 89:     BVInsertVec(pep->V,0,pep->IS[0]);
 90:     BVInsertVec(pep->V,1,pep->IS[1]);
 91:   }
 92:   if (pep->nini<0) {
 93:     SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
 94:   }
 95:   return(0);
 96: }

100: PetscErrorCode PEPExtractVectors_QArnoldi(PEP pep)
101: {
103:   PetscInt       i,k=pep->nconv,ldds;
104:   PetscScalar    *X,*pX0;
105:   Mat            X0;

108:   if (pep->nconv==0) return(0);
109:   DSGetLeadingDimension(pep->ds,&ldds);
110:   DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
111:   DSGetArray(pep->ds,DS_MAT_X,&X);

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

130: /*
131:   Compute a step of Classical Gram-Schmidt orthogonalization
132: */
133: 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)
134: {
136:   PetscBLASInt   ione = 1,j_1 = j+1;
137:   PetscReal      x,y;
138:   PetscScalar    dot,one = 1.0,zero = 0.0;

141:   /* compute norm of v and w */
142:   if (onorm) {
143:     VecNorm(v,NORM_2,&x);
144:     VecNorm(w,NORM_2,&y);
145:     *onorm = PetscSqrtReal(x*x+y*y);
146:   }

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

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

164:   /* compute norm of v and w */
165:   if (norm) {
166:     VecNorm(v,NORM_2,&x);
167:     VecNorm(w,NORM_2,&y);
168:     *norm = PetscSqrtReal(x*x+y*y);
169:   }
170:   return(0);
171: }

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

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

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

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

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

265:   DSGetLeadingDimension(pep->ds,&ld);
266:   lwork = 7*pep->ncv;
267:   PetscMalloc1(lwork,&work);
268:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
269:   RGSetScale(pep->rg,sinv?1.0/pep->sfactor:pep->sfactor);
270:   STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);

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

281:    /* Restart loop */
282:   l = 0;
283:   while (pep->reason == PEP_CONVERGED_ITERATING) {
284:     pep->its++;

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

299:     /* Solve projected problem */
300:     DSSolve(pep->ds,pep->eigr,pep->eigi);
301:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
302:     DSUpdateExtraRow(pep->ds);

304:     /* Check convergence */
305:     PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
306:     if (pep->its >= pep->max_it) pep->reason = PEP_DIVERGED_ITS;
307:     if (k >= pep->nev) pep->reason = PEP_CONVERGED_TOL;
308:     nconv = k;

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

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

332:     pep->nconv = k;
333:     PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
334:   }

336:   for (j=0;j<pep->nconv;j++) {
337:     pep->eigr[j] *= pep->sfactor;
338:     pep->eigi[j] *= pep->sfactor;
339:   }

341:   STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
342:   RGSetScale(pep->rg,1.0);

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

354: static PetscErrorCode PEPQArnoldiSetRestart_QArnoldi(PEP pep,PetscReal keep)
355: {
356:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

359:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
360:   else {
361:     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]");
362:     ctx->keep = keep;
363:   }
364:   return(0);
365: }

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

374:    Logically Collective on PEP

376:    Input Parameters:
377: +  pep  - the eigenproblem solver context
378: -  keep - the number of vectors to be kept at restart

380:    Options Database Key:
381: .  -pep_qarnoldi_restart - Sets the restart parameter

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

386:    Level: advanced

388: .seealso: PEPQArnoldiGetRestart()
389: @*/
390: PetscErrorCode PEPQArnoldiSetRestart(PEP pep,PetscReal keep)
391: {

397:   PetscTryMethod(pep,"PEPQArnoldiSetRestart_C",(PEP,PetscReal),(pep,keep));
398:   return(0);
399: }

403: static PetscErrorCode PEPQArnoldiGetRestart_QArnoldi(PEP pep,PetscReal *keep)
404: {
405:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

408:   *keep = ctx->keep;
409:   return(0);
410: }

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

417:    Not Collective

419:    Input Parameter:
420: .  pep - the eigenproblem solver context

422:    Output Parameter:
423: .  keep - the restart parameter

425:    Level: advanced

427: .seealso: PEPQArnoldiSetRestart()
428: @*/
429: PetscErrorCode PEPQArnoldiGetRestart(PEP pep,PetscReal *keep)
430: {

436:   PetscTryMethod(pep,"PEPQArnoldiGetRestart_C",(PEP,PetscReal*),(pep,keep));
437:   return(0);
438: }

442: static PetscErrorCode PEPQArnoldiSetLocking_QArnoldi(PEP pep,PetscBool lock)
443: {
444:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

447:   ctx->lock = lock;
448:   return(0);
449: }

453: /*@
454:    PEPQArnoldiSetLocking - Choose between locking and non-locking variants of
455:    the Q-Arnoldi method.

457:    Logically Collective on PEP

459:    Input Parameters:
460: +  pep  - the eigenproblem solver context
461: -  lock - true if the locking variant must be selected

463:    Options Database Key:
464: .  -pep_qarnoldi_locking - Sets the locking flag

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

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

474:    Level: advanced

476: .seealso: PEPQArnoldiGetLocking()
477: @*/
478: PetscErrorCode PEPQArnoldiSetLocking(PEP pep,PetscBool lock)
479: {

485:   PetscTryMethod(pep,"PEPQArnoldiSetLocking_C",(PEP,PetscBool),(pep,lock));
486:   return(0);
487: }

491: static PetscErrorCode PEPQArnoldiGetLocking_QArnoldi(PEP pep,PetscBool *lock)
492: {
493:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

496:   *lock = ctx->lock;
497:   return(0);
498: }

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

505:    Not Collective

507:    Input Parameter:
508: .  pep - the eigenproblem solver context

510:    Output Parameter:
511: .  lock - the locking flag

513:    Level: advanced

515: .seealso: PEPQArnoldiSetLocking()
516: @*/
517: PetscErrorCode PEPQArnoldiGetLocking(PEP pep,PetscBool *lock)
518: {

524:   PetscTryMethod(pep,"PEPQArnoldiGetLocking_C",(PEP,PetscBool*),(pep,lock));
525:   return(0);
526: }

530: PetscErrorCode PEPSetFromOptions_QArnoldi(PetscOptions *PetscOptionsObject,PEP pep)
531: {
533:   PetscBool      flg,lock;
534:   PetscReal      keep;

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

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

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

569: PetscErrorCode PEPDestroy_QArnoldi(PEP pep)
570: {

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

584: PETSC_EXTERN PetscErrorCode PEPCreate_QArnoldi(PEP pep)
585: {
586:   PEP_QARNOLDI   *ctx;

590:   PetscNewLog(pep,&ctx);
591:   pep->data = (void*)ctx;
592:   ctx->lock = PETSC_TRUE;

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