Actual source code: nleigs.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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 nonlinear eigensolver: "nleigs"

 13:    Method: NLEIGS

 15:    Algorithm:

 17:        Fully rational Krylov method for nonlinear eigenvalue problems.

 19:    References:

 21:        [1] S. Guttel et al., "NLEIGS: A class of robust fully rational Krylov
 22:            method for nonlinear eigenvalue problems", SIAM J. Sci. Comput.
 23:            36(6):A2842-A2864, 2014.
 24: */

 26: #include <slepc/private/nepimpl.h>         /*I "slepcnep.h" I*/
 27: #include <slepcblaslapack.h>

 29: #define  LBPOINTS  100   /* default value of the maximum number of Leja-Bagby points */
 30: #define  NDPOINTS  1e4   /* number of discretization points */

 32: typedef struct {
 33:   PetscInt       nmat;      /* number of interpolation points */
 34:   PetscScalar    *s,*xi;    /* Leja-Bagby points */
 35:   PetscScalar    *beta;     /* scaling factors */
 36:   Mat            *D;        /* divided difference matrices */
 37:   PetscScalar    *coeffD;   /* coefficients for divided differences in split form */
 38:   PetscInt       nshifts;   /* provided number of shifts */
 39:   PetscScalar    *shifts;   /* user-provided shifts for the Rational Krylov variant */
 40:   PetscInt       nshiftsw;  /* actual number of shifts (1 if Krylov-Schur) */
 41:   PetscReal      ddtol;     /* tolerance for divided difference convergence */
 42:   PetscInt       ddmaxit;   /* maximum number of divided difference terms */
 43:   PetscReal      keep;      /* restart parameter */
 44:   PetscBool      lock;      /* locking/non-locking variant */
 45:   PetscInt       idxrk;     /* index of next shift to use */
 46:   KSP            *ksp;      /* ksp array for storing shift factorizations */
 47:   Vec            vrn;       /* random vector with normally distributed value */
 48:   void           *singularitiesctx;
 49:   PetscErrorCode (*computesingularities)(NEP,PetscInt*,PetscScalar*,void*);
 50: } NEP_NLEIGS;

 52: typedef struct {
 53:   PetscInt    nmat,maxnmat;
 54:   PetscScalar *coeff;
 55:   Mat         *A;
 56:   Vec         t;
 57: } ShellMatCtx;

 59: PETSC_STATIC_INLINE PetscErrorCode NEPNLEIGSSetShifts(NEP nep)
 60: {
 61:   NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;

 64:   if (!ctx->nshifts) {
 65:     ctx->shifts = &nep->target;
 66:     ctx->nshiftsw = 1;
 67:   } else ctx->nshiftsw = ctx->nshifts;
 68:   return(0);
 69: }

 71: static PetscErrorCode NEPNLEIGSBackTransform(PetscObject ob,PetscInt n,PetscScalar *valr,PetscScalar *vali)
 72: {
 73:   NEP         nep;
 74:   PetscInt    j;
 75: #if !defined(PETSC_USE_COMPLEX)
 76:   PetscScalar t;
 77: #endif

 80:   nep = (NEP)ob;
 81: #if !defined(PETSC_USE_COMPLEX)
 82:   for (j=0;j<n;j++) {
 83:     if (vali[j] == 0) valr[j] = 1.0 / valr[j] + nep->target;
 84:     else {
 85:       t = valr[j] * valr[j] + vali[j] * vali[j];
 86:       valr[j] = valr[j] / t + nep->target;
 87:       vali[j] = - vali[j] / t;
 88:     }
 89:   }
 90: #else
 91:   for (j=0;j<n;j++) {
 92:     valr[j] = 1.0 / valr[j] + nep->target;
 93:   }
 94: #endif
 95:   return(0);
 96: }

 98: /* Computes the roots of a polynomial */
 99: static PetscErrorCode NEPNLEIGSAuxiliarPRootFinder(PetscInt deg,PetscScalar *polcoeffs,PetscScalar *wr,PetscScalar *wi,PetscBool *avail)
100: {
102:   PetscScalar    *C,*work;
103:   PetscBLASInt   n_,info,lwork;
104:   PetscInt       i;
105: #if defined(PETSC_USE_COMPLEX)
106:   PetscReal      *rwork;
107: #endif

110: #if defined(PETSC_MISSING_LAPACK_GEEV)
111:   *avail = PETSC_FALSE;
112: #else
113:   *avail = PETSC_TRUE;
114:   if (deg>0) {
115:     PetscCalloc1(deg*deg,&C);
116:     PetscBLASIntCast(deg,&n_);
117:     for (i=0;i<deg-1;i++) {
118:       C[(deg+1)*i+1]   = 1.0;
119:       C[(deg-1)*deg+i] = -polcoeffs[deg-i]/polcoeffs[0];
120:     }
121:     C[deg*deg+-1] = -polcoeffs[1]/polcoeffs[0];
122:     PetscBLASIntCast(3*deg,&lwork);
123: #if !defined(PETSC_USE_COMPLEX)
124:     PetscMalloc1(lwork,&work);
125:     PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&n_,C,&n_,wr,wi,NULL,&n_,NULL,&n_,work,&lwork,&info));
126:     if (info) *avail = PETSC_FALSE;
127:     PetscFree(work); 
128: #else
129:     PetscMalloc2(2*deg,&rwork,lwork,&work);
130:     PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&n_,C,&n_,wr,NULL,&n_,NULL,&n_,work,&lwork,rwork,&info));
131:     if (info) *avail = PETSC_FALSE;
132:     PetscFree2(rwork,work);
133: #endif
134:     PetscFree(C);
135:   }
136: #endif
137:   return(0);
138: }

140: static PetscErrorCode NEPNLEIGSAuxiliarRmDuplicates(PetscInt nin,PetscScalar *pin,PetscInt *nout,PetscScalar *pout)
141: {
142:   PetscInt i,j;

145:   for (i=0;i<nin;i++) {
146:     pout[(*nout)++] = pin[i];
147:     for (j=0;j<*nout-1;j++)
148:       if (PetscAbsScalar(pin[i]-pout[j])<PETSC_MACHINE_EPSILON*100) {
149:         (*nout)--;
150:         break;
151:       }
152:   }
153:   return(0);
154: }

156: static PetscErrorCode NEPNLEIGSFNSingularities(FN f,PetscInt *nisol,PetscScalar **isol,PetscBool *rational)
157: {
159:   FNCombineType  ctype;
160:   FN             f1,f2;
161:   PetscInt       i,nq,nisol1,nisol2;
162:   PetscScalar    *qcoeff,*wr,*wi,*isol1,*isol2;
163:   PetscBool      flg,avail,rat1,rat2;

166:   *rational = PETSC_FALSE;
167:   PetscObjectTypeCompare((PetscObject)f,FNRATIONAL,&flg);
168:   if (flg) {
169:     *rational = PETSC_TRUE;
170:     FNRationalGetDenominator(f,&nq,&qcoeff);
171:     if (nq>1) {
172:       PetscMalloc2(nq-1,&wr,nq-1,&wi);
173:       NEPNLEIGSAuxiliarPRootFinder(nq-1,qcoeff,wr,wi,&avail);
174:       if (avail) {
175:         PetscCalloc1(nq-1,isol);
176:         *nisol = 0;
177:         for (i=0;i<nq-1;i++) 
178: #if !defined(PETSC_USE_COMPLEX)
179:           if (wi[i]==0)
180: #endif 
181:             (*isol)[(*nisol)++] = wr[i];
182:         nq = *nisol; *nisol = 0;
183:         for (i=0;i<nq;i++) wr[i] = (*isol)[i];
184:         NEPNLEIGSAuxiliarRmDuplicates(nq,wr,nisol,*isol);
185:         PetscFree2(wr,wi);
186:       } else { *nisol=0; *isol = NULL; }
187:     } else { *nisol = 0; *isol = NULL; }
188:     PetscFree(qcoeff);
189:   }
190:   PetscObjectTypeCompare((PetscObject)f,FNCOMBINE,&flg);
191:   if (flg) {
192:     FNCombineGetChildren(f,&ctype,&f1,&f2);
193:     if (ctype != FN_COMBINE_COMPOSE && ctype != FN_COMBINE_DIVIDE) {
194:       NEPNLEIGSFNSingularities(f1,&nisol1,&isol1,&rat1);
195:       NEPNLEIGSFNSingularities(f2,&nisol2,&isol2,&rat2);
196:       if (nisol1+nisol2>0) {
197:         PetscCalloc1(nisol1+nisol2,isol);
198:         *nisol = 0; 
199:         NEPNLEIGSAuxiliarRmDuplicates(nisol1,isol1,nisol,*isol);
200:         NEPNLEIGSAuxiliarRmDuplicates(nisol2,isol2,nisol,*isol);
201:       }
202:       *rational = (rat1&&rat2)?PETSC_TRUE:PETSC_FALSE;
203:       PetscFree(isol1);
204:       PetscFree(isol2);
205:     }
206:   }
207:   return(0);
208: }

210: static PetscErrorCode NEPNLEIGSRationalSingularities(NEP nep,PetscInt *ndptx,PetscScalar *dxi,PetscBool *rational)
211: {
213:   PetscInt       nt,i,nisol;
214:   FN             f;
215:   PetscScalar    *isol;
216:   PetscBool      rat;

219:   *rational = PETSC_TRUE;
220:   *ndptx = 0;
221:   NEPGetSplitOperatorInfo(nep,&nt,NULL);
222:   for (i=0;i<nt;i++) {
223:     NEPGetSplitOperatorTerm(nep,i,NULL,&f);
224:     NEPNLEIGSFNSingularities(f,&nisol,&isol,&rat);
225:     if (nisol) {
226:       NEPNLEIGSAuxiliarRmDuplicates(nisol,isol,ndptx,dxi);
227:       PetscFree(isol);
228:     }
229:     *rational = ((*rational)&&rat)?PETSC_TRUE:PETSC_FALSE;
230:   }
231:   return(0);
232: }

234: static PetscErrorCode NEPNLEIGSLejaBagbyPoints(NEP nep)
235: {
237:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
238:   PetscInt       i,k,ndpt=NDPOINTS,ndptx=NDPOINTS;
239:   PetscScalar    *ds,*dsi,*dxi,*nrs,*nrxi,*s=ctx->s,*xi=ctx->xi,*beta=ctx->beta;
240:   PetscReal      maxnrs,minnrxi;
241:   PetscBool      rational;
242: #if !defined(PETSC_USE_COMPLEX)
243:   PetscReal      a,b,h;
244: #endif

247:   PetscMalloc5(ndpt+1,&ds,ndpt+1,&dsi,ndpt,&dxi,ndpt+1,&nrs,ndpt,&nrxi);

249:   /* Discretize the target region boundary */
250:   RGComputeContour(nep->rg,ndpt,ds,dsi);
251: #if !defined(PETSC_USE_COMPLEX)
252:   for (i=0;i<ndpt;i++) if (dsi[i]!=0.0) break;
253:   if (i<ndpt) {
254:     if (nep->problem_type==NEP_RATIONAL) {
255:       /* Select a segment in the real axis */
256:       RGComputeBoundingBox(nep->rg,&a,&b,NULL,NULL);
257:       if (a<=-PETSC_MAX_REAL || b>=PETSC_MAX_REAL) SETERRQ(PetscObjectComm((PetscObject)nep),1,"NLEIGS requires a bounded target set");
258:       h = (b-a)/ndpt;
259:       for (i=0;i<ndpt;i++) {ds[i] = a+h*i; dsi[i] = 0.0;}
260:     } else SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NLEIGS with real arithmetic requires the target set to be included in the real axis");
261:   }
262: #endif
263:   /* Discretize the singularity region */
264:   if (ctx->computesingularities) {
265:     (ctx->computesingularities)(nep,&ndptx,dxi,ctx->singularitiesctx);
266:   } else {
267:     if (nep->problem_type==NEP_RATIONAL) {
268:       NEPNLEIGSRationalSingularities(nep,&ndptx,dxi,&rational);
269:       if (!rational) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Failed to determine singularities automatically in rational problem; consider solving the problem as general");
270:     } else ndptx = 0;
271:   }

273:   /* Look for Leja-Bagby points in the discretization sets */
274:   s[0]    = ds[0];
275:   xi[0]   = (ndptx>0)?dxi[0]:PETSC_INFINITY;
276:   if (PetscAbsScalar(xi[0])<10*PETSC_MACHINE_EPSILON) SETERRQ2(PetscObjectComm((PetscObject)nep),1,"Singularity point %D is nearly zero: %g; consider removing the singularity or shifting the problem",0,(double)PetscAbsScalar(xi[0]));
277:   beta[0] = 1.0; /* scaling factors are also computed here */
278:   for (i=0;i<ndpt;i++) {
279:     nrs[i] = 1.0;
280:     nrxi[i] = 1.0;
281:   }
282:   for (k=1;k<ctx->ddmaxit;k++) {
283:     maxnrs = 0.0;
284:     minnrxi = PETSC_MAX_REAL;
285:     for (i=0;i<ndpt;i++) {
286:       nrs[i] *= ((ds[i]-s[k-1])/(1.0-ds[i]/xi[k-1]))/beta[k-1];
287:       if (PetscAbsScalar(nrs[i])>maxnrs) {maxnrs = PetscAbsScalar(nrs[i]); s[k] = ds[i];}
288:     }
289:     if (ndptx>k) {
290:       for (i=1;i<ndptx;i++) {
291:         nrxi[i] *= ((dxi[i]-s[k-1])/(1.0-dxi[i]/xi[k-1]))/beta[k-1];
292:         if (PetscAbsScalar(nrxi[i])<minnrxi) {minnrxi = PetscAbsScalar(nrxi[i]); xi[k] = dxi[i];}
293:       }
294:       if (PetscAbsScalar(xi[k])<10*PETSC_MACHINE_EPSILON) SETERRQ2(PetscObjectComm((PetscObject)nep),1,"Singularity point %D is nearly zero: %g; consider removing the singularity or shifting the problem",k,(double)PetscAbsScalar(xi[k]));
295:     } else xi[k] = PETSC_INFINITY;
296:     beta[k] = maxnrs;
297:   }
298:   PetscFree5(ds,dsi,dxi,nrs,nrxi);
299:   return(0);
300: }

302: static PetscErrorCode NEPNLEIGSEvalNRTFunct(NEP nep,PetscInt k,PetscScalar sigma,PetscScalar *b)
303: {
304:   NEP_NLEIGS  *ctx=(NEP_NLEIGS*)nep->data;
305:   PetscInt    i;
306:   PetscScalar *beta=ctx->beta,*s=ctx->s,*xi=ctx->xi;

309:   b[0] = 1.0/beta[0];
310:   for (i=0;i<k;i++) {
311:     b[i+1] = ((sigma-s[i])*b[i])/(beta[i+1]*(1.0-sigma/xi[i]));
312:   }
313:   return(0);
314: }

316: static PetscErrorCode MatMult_Fun(Mat A,Vec x,Vec y)
317: {
319:   ShellMatCtx    *ctx;
320:   PetscInt       i;

323:   MatShellGetContext(A,(void**)&ctx);
324:   MatMult(ctx->A[0],x,y);
325:   if (ctx->coeff[0]!=1.0) { VecScale(y,ctx->coeff[0]); }
326:   for (i=1;i<ctx->nmat;i++) {
327:     MatMult(ctx->A[i],x,ctx->t);
328:     VecAXPY(y,ctx->coeff[i],ctx->t);
329:   }
330:   return(0);
331: }

333: static PetscErrorCode MatMultTranspose_Fun(Mat A,Vec x,Vec y)
334: {
336:   ShellMatCtx    *ctx;
337:   PetscInt       i;

340:   MatShellGetContext(A,(void**)&ctx);
341:   MatMultTranspose(ctx->A[0],x,y);
342:   if (ctx->coeff[0]!=1.0) { VecScale(y,ctx->coeff[0]); }
343:   for (i=1;i<ctx->nmat;i++) {
344:     MatMultTranspose(ctx->A[i],x,ctx->t);
345:     VecAXPY(y,ctx->coeff[i],ctx->t);
346:   }
347:   return(0);
348: }

350: static PetscErrorCode MatGetDiagonal_Fun(Mat A,Vec diag)
351: {
353:   ShellMatCtx    *ctx;
354:   PetscInt       i;

357:   MatShellGetContext(A,(void**)&ctx);
358:   MatGetDiagonal(ctx->A[0],diag);
359:   if (ctx->coeff[0]!=1.0) { VecScale(diag,ctx->coeff[0]); }
360:   for (i=1;i<ctx->nmat;i++) {
361:     MatGetDiagonal(ctx->A[i],ctx->t);
362:     VecAXPY(diag,ctx->coeff[i],ctx->t);
363:   }
364:   return(0);
365: }

367: static PetscErrorCode MatDuplicate_Fun(Mat A,MatDuplicateOption op,Mat *B)
368: {
369:   PetscInt       n,i;
370:   ShellMatCtx    *ctxnew,*ctx;
371:   void           (*fun)();

375:   MatShellGetContext(A,(void**)&ctx);
376:   PetscNew(&ctxnew);
377:   ctxnew->nmat = ctx->nmat;
378:   ctxnew->maxnmat = ctx->maxnmat;
379:   PetscMalloc2(ctxnew->maxnmat,&ctxnew->A,ctxnew->maxnmat,&ctxnew->coeff);
380:   for (i=0;i<ctx->nmat;i++) {
381:     PetscObjectReference((PetscObject)ctx->A[i]);
382:     ctxnew->A[i] = ctx->A[i];
383:     ctxnew->coeff[i] = ctx->coeff[i];
384:   }
385:   MatGetSize(ctx->A[0],&n,NULL);
386:   VecDuplicate(ctx->t,&ctxnew->t);
387:   MatCreateShell(PETSC_COMM_WORLD,n,n,n,n,(void*)ctxnew,B);
388:   MatShellGetOperation(A,MATOP_MULT,&fun);
389:   MatShellSetOperation(*B,MATOP_MULT,fun);
390:   MatShellGetOperation(A,MATOP_MULT_TRANSPOSE,&fun);
391:   MatShellSetOperation(*B,MATOP_MULT_TRANSPOSE,fun);
392:   MatShellGetOperation(A,MATOP_GET_DIAGONAL,&fun);
393:   MatShellSetOperation(*B,MATOP_GET_DIAGONAL,fun);
394:   MatShellGetOperation(A,MATOP_DUPLICATE,&fun);
395:   MatShellSetOperation(*B,MATOP_DUPLICATE,fun);
396:   MatShellGetOperation(A,MATOP_DESTROY,&fun);
397:   MatShellSetOperation(*B,MATOP_DESTROY,fun);
398:   MatShellGetOperation(A,MATOP_AXPY,&fun);
399:   MatShellSetOperation(*B,MATOP_AXPY,fun);
400:   return(0);
401: }

403: static PetscErrorCode MatDestroy_Fun(Mat A)
404: {
405:   ShellMatCtx    *ctx;
407:   PetscInt       i;

410:   if (A) {
411:     MatShellGetContext(A,(void**)&ctx);
412:     for (i=0;i<ctx->nmat;i++) {
413:       MatDestroy(&ctx->A[i]);
414:     }
415:     VecDestroy(&ctx->t);
416:     PetscFree2(ctx->A,ctx->coeff);
417:     PetscFree(ctx);
418:   }
419:   return(0);
420: }

422: static PetscErrorCode MatAXPY_Fun(Mat Y,PetscScalar a,Mat X,MatStructure str)
423: {
424:   ShellMatCtx    *ctxY,*ctxX;
426:   PetscInt       i,j;
427:   PetscBool      found;

430:   MatShellGetContext(Y,(void**)&ctxY);
431:   MatShellGetContext(X,(void**)&ctxX);
432:   for (i=0;i<ctxX->nmat;i++) {
433:     found = PETSC_FALSE;
434:     for (j=0;!found&&j<ctxY->nmat;j++) {
435:       if (ctxX->A[i]==ctxY->A[j]) {
436:         found = PETSC_TRUE;
437:         ctxY->coeff[j] += a*ctxX->coeff[i];
438:       }
439:     }
440:     if (!found) {
441:       ctxY->coeff[ctxY->nmat] = a*ctxX->coeff[i];
442:       ctxY->A[ctxY->nmat++] = ctxX->A[i];
443:       PetscObjectReference((PetscObject)ctxX->A[i]);
444:     }
445:   }
446:   return(0);
447: }

449: static PetscErrorCode MatScale_Fun(Mat M,PetscScalar a)
450: {
451:   ShellMatCtx    *ctx;
453:   PetscInt       i;

456:   MatShellGetContext(M,(void**)&ctx);
457:   for (i=0;i<ctx->nmat;i++) ctx->coeff[i] *= a;
458:   return(0);
459: }

461: static PetscErrorCode NLEIGSMatToMatShellArray(Mat M,Mat *Ms,PetscInt maxnmat)
462: {
464:   ShellMatCtx    *ctx;
465:   PetscInt       n;
466:   PetscBool      has;

469:   MatHasOperation(M,MATOP_DUPLICATE,&has);
470:   if (!has) SETERRQ(PetscObjectComm((PetscObject)M),1,"MatDuplicate operation required");
471:   PetscNew(&ctx);
472:   ctx->maxnmat = maxnmat;
473:   PetscMalloc2(ctx->maxnmat,&ctx->A,ctx->maxnmat,&ctx->coeff);
474:   MatDuplicate(M,MAT_COPY_VALUES,&ctx->A[0]);
475:   ctx->nmat = 1;
476:   ctx->coeff[0] = 1.0;
477:   MatCreateVecs(M,&ctx->t,NULL);
478:   MatGetSize(M,&n,NULL);
479:   MatCreateShell(PetscObjectComm((PetscObject)M),n,n,n,n,(void*)ctx,Ms);
480:   MatShellSetOperation(*Ms,MATOP_MULT,(void(*)())MatMult_Fun);
481:   MatShellSetOperation(*Ms,MATOP_MULT_TRANSPOSE,(void(*)())MatMultTranspose_Fun);
482:   MatShellSetOperation(*Ms,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Fun);
483:   MatShellSetOperation(*Ms,MATOP_DUPLICATE,(void(*)())MatDuplicate_Fun);
484:   MatShellSetOperation(*Ms,MATOP_DESTROY,(void(*)())MatDestroy_Fun);
485:   MatShellSetOperation(*Ms,MATOP_AXPY,(void(*)())MatAXPY_Fun);
486:   MatShellSetOperation(*Ms,MATOP_SCALE,(void(*)())MatScale_Fun);
487:   return(0);
488: }

490: static PetscErrorCode NEPNLEIGSNormEstimation(NEP nep,Mat M,PetscReal *norm,Vec *w)
491: {
492:   PetscScalar    *z,*x,*y;
493:   PetscReal      tr;
494:   Vec            X=w[0],Y=w[1];
495:   PetscInt       n,i;
496:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
497:   PetscRandom    rand;

501:   if (!ctx->vrn) {
502:     /* generate a random vector with normally distributed entries with the Box-Muller transform */
503:     BVGetRandomContext(nep->V,&rand);
504:     MatCreateVecs(M,&ctx->vrn,NULL);
505:     VecSetRandom(X,rand);
506:     VecSetRandom(Y,rand);
507:     VecGetLocalSize(ctx->vrn,&n);
508:     VecGetArray(ctx->vrn,&z);
509:     VecGetArray(X,&x);
510:     VecGetArray(Y,&y);
511:     for (i=0;i<n;i++) {
512: #if defined(PETSC_USE_COMPLEX)
513:       z[i] = PetscSqrtReal(-2.0*PetscLogReal(PetscRealPart(x[i])))*PetscCosReal(2.0*PETSC_PI*PetscRealPart(y[i]));
514:       z[i] += PETSC_i*(PetscSqrtReal(-2.0*PetscLogReal(PetscImaginaryPart(x[i])))*PetscCosReal(2.0*PETSC_PI*PetscImaginaryPart(y[i])));
515: #else
516:       z[i] = PetscSqrtReal(-2.0*PetscLogReal(x[i]))*PetscCosReal(2.0*PETSC_PI*y[i]);
517: #endif
518:     }
519:     VecRestoreArray(ctx->vrn,&z);
520:     VecRestoreArray(X,&x);
521:     VecRestoreArray(Y,&y);
522:     VecNorm(ctx->vrn,NORM_2,&tr);
523:     VecScale(ctx->vrn,1/tr);
524:   }
525:   /* matrix-free norm estimator of Ipsen http://www4.ncsu.edu/~ipsen/ps/slides_ima.pdf */
526:   MatGetSize(M,&n,NULL);
527:   MatMult(M,ctx->vrn,X);
528:   VecNorm(X,NORM_2,norm);
529:   *norm *= PetscSqrtReal((PetscReal)n);
530:   return(0);
531: }

533: static PetscErrorCode NEPNLEIGSDividedDifferences_split(NEP nep)
534: {
536:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
537:   PetscInt       k,j,i,maxnmat,nmax;
538:   PetscReal      norm0,norm,*matnorm;
539:   PetscScalar    *s=ctx->s,*beta=ctx->beta,*xi=ctx->xi,*b,alpha,*coeffs,*pK,*pH,sone=1.0;
540:   Mat            T,Ts,K,H;
541:   PetscBool      shell,hasmnorm,matrix=PETSC_TRUE;
542:   PetscBLASInt   n_;

545:   nmax = ctx->ddmaxit;
546:   PetscMalloc1(nep->nt*nmax,&ctx->coeffD);
547:   PetscMalloc3(nmax+1,&b,nmax+1,&coeffs,nep->nt,&matnorm);
548:   for (j=0;j<nep->nt;j++) {
549:     MatHasOperation(nep->A[j],MATOP_NORM,&hasmnorm);
550:     if (!hasmnorm) break;
551:     MatNorm(nep->A[j],NORM_INFINITY,matnorm+j);
552:   }
553:   /* Try matrix functions scheme */
554:   PetscCalloc2(nmax*nmax,&pK,nmax*nmax,&pH);
555:   for (i=0;i<nmax-1;i++) {
556:     pK[(nmax+1)*i]   = 1.0;
557:     pK[(nmax+1)*i+1] = beta[i+1]/xi[i];
558:     pH[(nmax+1)*i]   = s[i];
559:     pH[(nmax+1)*i+1] = beta[i+1];
560:   }
561:   pH[nmax*nmax-1] = s[nmax-1];
562:   pK[nmax*nmax-1] = 1.0;
563:   PetscBLASIntCast(nmax,&n_);
564:   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","L","N","U",&n_,&n_,&sone,pK,&n_,pH,&n_));
565:   /* The matrix to be used is in H. K will be a work-space matrix */
566:   MatCreateSeqDense(PETSC_COMM_SELF,nmax,nmax,pH,&H);
567:   MatCreateSeqDense(PETSC_COMM_SELF,nmax,nmax,pK,&K);
568:   for (j=0;matrix&&j<nep->nt;j++) {
569:     PetscPushErrorHandler(PetscIgnoreErrorHandler,NULL);
570:     FNEvaluateFunctionMat(nep->f[j],H,K);
571:     PetscPopErrorHandler();
572:     if (!ierr) { 
573:       for (i=0;i<nmax;i++) { ctx->coeffD[j+i*nep->nt] = pK[i]*beta[0]; }
574:     } else {
575:       matrix = PETSC_FALSE;
576:       PetscFPTrapPop();
577:     }
578:   }
579:   MatDestroy(&H);
580:   MatDestroy(&K);
581:   if (!matrix) {
582:     for (j=0;j<nep->nt;j++) {
583:       FNEvaluateFunction(nep->f[j],s[0],ctx->coeffD+j);
584:       ctx->coeffD[j] *= beta[0];
585:     }
586:   }
587:   if (hasmnorm) {
588:     norm0 = 0.0;
589:     for (j=0;j<nep->nt;j++) norm0 += matnorm[j]*PetscAbsScalar(ctx->coeffD[j]);
590:   } else {
591:     norm0 = 0.0;
592:     for (j=0;j<nep->nt;j++) norm0 = PetscMax(PetscAbsScalar(ctx->coeffD[j]),norm0);
593:   }
594:   ctx->nmat = ctx->ddmaxit;
595:   for (k=1;k<ctx->ddmaxit;k++) {
596:     if (!matrix) {
597:       NEPNLEIGSEvalNRTFunct(nep,k,s[k],b);
598:       for (i=0;i<nep->nt;i++) {
599:         FNEvaluateFunction(nep->f[i],s[k],ctx->coeffD+k*nep->nt+i);
600:         for (j=0;j<k;j++) {
601:           ctx->coeffD[k*nep->nt+i] -= b[j]*ctx->coeffD[i+nep->nt*j];
602:         }
603:         ctx->coeffD[k*nep->nt+i] /= b[k];
604:       }
605:     }
606:     if (hasmnorm) {
607:       norm = 0.0;
608:       for (j=0;j<nep->nt;j++) norm += matnorm[j]*PetscAbsScalar(ctx->coeffD[k*nep->nt+j]);
609:     } else {
610:       norm = 0.0;
611:       for (j=0;j<nep->nt;j++) norm = PetscMax(PetscAbsScalar(ctx->coeffD[k*nep->nt+j]),norm);
612:     }
613:     if (norm/norm0 < ctx->ddtol) {
614:       ctx->nmat = k+1;
615:       break;
616:     }
617:   }
618:   if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
619:   PetscObjectTypeCompare((PetscObject)nep->A[0],MATSHELL,&shell);
620:   maxnmat = PetscMax(ctx->ddmaxit,nep->nt);
621:   for (i=0;i<ctx->nshiftsw;i++) {
622:     NEPNLEIGSEvalNRTFunct(nep,ctx->nmat-1,ctx->shifts[i],coeffs);
623:     if (!shell) {
624:       MatDuplicate(nep->A[0],MAT_COPY_VALUES,&T);
625:     } else {
626:       NLEIGSMatToMatShellArray(nep->A[0],&T,maxnmat);
627:     }
628:     alpha = 0.0;
629:     for (j=0;j<ctx->nmat;j++) alpha += coeffs[j]*ctx->coeffD[j*nep->nt];
630:     MatScale(T,alpha);
631:     for (k=1;k<nep->nt;k++) {
632:       alpha = 0.0;
633:       for (j=0;j<ctx->nmat;j++) alpha += coeffs[j]*ctx->coeffD[j*nep->nt+k];
634:       if (shell) {
635:         NLEIGSMatToMatShellArray(nep->A[k],&Ts,maxnmat);
636:       }
637:       MatAXPY(T,alpha,shell?Ts:nep->A[k],nep->mstr);
638:       if (shell) {
639:         MatDestroy(&Ts);
640:       }
641:     }
642:     KSPSetOperators(ctx->ksp[i],T,T);
643:     KSPSetUp(ctx->ksp[i]);
644:     MatDestroy(&T);
645:   }
646:   PetscFree3(b,coeffs,matnorm);
647:   PetscFree2(pK,pH);
648:   return(0);
649: }

651: static PetscErrorCode NEPNLEIGSDividedDifferences_callback(NEP nep)
652: {
654:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
655:   PetscInt       k,j,i,maxnmat;
656:   PetscReal      norm0,norm;
657:   PetscScalar    *s=ctx->s,*beta=ctx->beta,*b,*coeffs;
658:   Mat            *D=ctx->D,T;
659:   PetscBool      shell,has,vec=PETSC_FALSE;
660:   Vec            w[2];

663:   PetscMalloc2(ctx->ddmaxit+1,&b,ctx->ddmaxit+1,&coeffs);
664:   T = nep->function;
665:   NEPComputeFunction(nep,s[0],T,T);
666:   PetscObjectTypeCompare((PetscObject)T,MATSHELL,&shell);
667:   maxnmat = PetscMax(ctx->ddmaxit,nep->nt);
668:   if (!shell) {
669:     MatDuplicate(T,MAT_COPY_VALUES,&D[0]);
670:   } else {
671:     NLEIGSMatToMatShellArray(T,&D[0],maxnmat);
672:   }
673:   if (beta[0]!=1.0) {
674:     MatScale(D[0],1.0/beta[0]);
675:   }
676:   MatHasOperation(D[0],MATOP_NORM,&has);
677:   if (has) {
678:     MatNorm(D[0],NORM_FROBENIUS,&norm0);
679:   } else {
680:     MatCreateVecs(D[0],NULL,&w[0]);
681:     VecDuplicate(w[0],&w[1]);
682:     vec = PETSC_TRUE;
683:     NEPNLEIGSNormEstimation(nep,D[0],&norm0,w);
684:   }
685:   ctx->nmat = ctx->ddmaxit;
686:   for (k=1;k<ctx->ddmaxit;k++) {
687:     NEPNLEIGSEvalNRTFunct(nep,k,s[k],b);
688:     NEPComputeFunction(nep,s[k],T,T);
689:     if (!shell) {
690:       MatDuplicate(T,MAT_COPY_VALUES,&D[k]);
691:     } else {
692:       NLEIGSMatToMatShellArray(T,&D[k],maxnmat);
693:     }
694:     for (j=0;j<k;j++) {
695:       MatAXPY(D[k],-b[j],D[j],nep->mstr);
696:     }
697:     MatScale(D[k],1.0/b[k]);
698:     MatHasOperation(D[k],MATOP_NORM,&has);
699:     if (has) {
700:       MatNorm(D[k],NORM_FROBENIUS,&norm);
701:     } else {
702:       if(!vec) {
703:         MatCreateVecs(D[k],NULL,&w[0]);
704:         VecDuplicate(w[0],&w[1]);
705:         vec = PETSC_TRUE;
706:       }
707:       NEPNLEIGSNormEstimation(nep,D[k],&norm,w);
708:     }
709:     if (norm/norm0 < ctx->ddtol) {
710:       ctx->nmat = k+1;
711:       break;
712:     }
713:   }
714:   if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
715:   for (i=0;i<ctx->nshiftsw;i++) {
716:     NEPNLEIGSEvalNRTFunct(nep,ctx->nmat-1,ctx->shifts[i],coeffs);
717:     MatDuplicate(ctx->D[0],MAT_COPY_VALUES,&T);
718:     if (coeffs[0]!=1.0) { MatScale(T,coeffs[0]); }
719:     for (j=1;j<ctx->nmat;j++) {
720:       MatAXPY(T,coeffs[j],ctx->D[j],nep->mstr);
721:     }
722:     KSPSetOperators(ctx->ksp[i],T,T);
723:     KSPSetUp(ctx->ksp[i]);
724:     MatDestroy(&T);
725:   }
726:   PetscFree2(b,coeffs);
727:   if (vec) {
728:     VecDestroy(&w[0]);
729:     VecDestroy(&w[1]);
730:   }
731:   return(0);
732: }

734: /*
735:    NEPKrylovConvergence - This is the analogue to EPSKrylovConvergence.
736: */
737: static PetscErrorCode NEPNLEIGSKrylovConvergence(NEP nep,PetscScalar *S,PetscInt ld,PetscInt nq,PetscScalar *H,PetscBool getall,PetscInt kini,PetscInt nits,PetscScalar betak,PetscReal betah,PetscInt *kout,Vec *w)
738: {
740:   PetscInt       k,newk,marker,inside;
741:   PetscScalar    re,im;
742:   PetscReal      resnorm,tt;
743:   PetscBool      istrivial;
744:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;

747:   RGIsTrivial(nep->rg,&istrivial);
748:   marker = -1;
749:   if (nep->trackall) getall = PETSC_TRUE;
750:   for (k=kini;k<kini+nits;k++) {
751:     /* eigenvalue */
752:     re = nep->eigr[k];
753:     im = nep->eigi[k];
754:     if (!istrivial) {
755:       if (!ctx->nshifts) {
756:         NEPNLEIGSBackTransform((PetscObject)nep,1,&re,&im);
757:       }
758:       RGCheckInside(nep->rg,1,&re,&im,&inside);
759:       if (marker==-1 && inside<0) marker = k;
760:     }
761:     newk = k;
762:     DSVectors(nep->ds,DS_MAT_X,&newk,&resnorm);
763:     tt = (ctx->nshifts)?SlepcAbsEigenvalue(betak-nep->eigr[k]*betah,nep->eigi[k]*betah):betah;
764:     resnorm *=  PetscAbsReal(tt);
765:     /* error estimate */
766:     (*nep->converged)(nep,nep->eigr[k],nep->eigi[k],resnorm,&nep->errest[k],nep->convergedctx);
767:     if (marker==-1 && nep->errest[k] >= nep->tol) marker = k;
768:     if (newk==k+1) {
769:       nep->errest[k+1] = nep->errest[k];
770:       k++;
771:     }
772:     if (marker!=-1 && !getall) break;
773:   }
774:   if (marker!=-1) k = marker;
775:   *kout = k;
776:   return(0);
777: }

779: PetscErrorCode NEPSetUp_NLEIGS(NEP nep)
780: {
782:   PetscInt       k,in;
783:   PetscScalar    zero=0.0;
784:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
785:   SlepcSC        sc;
786:   PetscBool      istrivial;

789:   NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd);
790:   if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
791:   if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
792:   if (!ctx->ddmaxit) ctx->ddmaxit = LBPOINTS;
793:   RGIsTrivial(nep->rg,&istrivial);
794:   if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPNLEIGS requires a nontrivial region defining the target set");
795:   if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;

797:   /* Initialize the NLEIGS context structure */
798:   k = ctx->ddmaxit;
799:   PetscMalloc4(k,&ctx->s,k,&ctx->xi,k,&ctx->beta,k,&ctx->D);
800:   nep->data = ctx;
801:   if (nep->tol==PETSC_DEFAULT) nep->tol = SLEPC_DEFAULT_TOL;
802:   if (ctx->ddtol==PETSC_DEFAULT) ctx->ddtol = nep->tol/10.0;
803:   if (!ctx->keep) ctx->keep = 0.5;

805:   /* Compute Leja-Bagby points and scaling values */
806:   NEPNLEIGSLejaBagbyPoints(nep);
807:   if (nep->problem_type!=NEP_RATIONAL) {
808:     RGCheckInside(nep->rg,1,&nep->target,&zero,&in);
809:     if (in<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The target is not inside the target set");
810:   }

812:   /* Compute the divided difference matrices */
813:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
814:     NEPNLEIGSDividedDifferences_split(nep);
815:   } else {
816:     NEPNLEIGSDividedDifferences_callback(nep);
817:   }
818:   NEPAllocateSolution(nep,ctx->nmat);
819:   NEPSetWorkVecs(nep,4);

821:   /* set-up DS and transfer split operator functions */
822:   DSSetType(nep->ds,ctx->nshifts?DSGNHEP:DSNHEP);
823:   DSAllocate(nep->ds,nep->ncv+1);
824:   DSGetSlepcSC(nep->ds,&sc);
825:   if (!ctx->nshifts) {
826:     sc->map = NEPNLEIGSBackTransform;
827:     DSSetExtraRow(nep->ds,PETSC_TRUE);
828:   }
829:   sc->mapobj        = (PetscObject)nep;
830:   sc->rg            = nep->rg;
831:   sc->comparison    = nep->sc->comparison;
832:   sc->comparisonctx = nep->sc->comparisonctx;
833:   return(0);
834: }

836: /*
837:   Norm of [sp;sq]
838: */
839: static PetscErrorCode NEPTOARSNorm2(PetscInt n,PetscScalar *S,PetscReal *norm)
840: {
842:   PetscBLASInt   n_,one=1;

845:   PetscBLASIntCast(n,&n_);
846:   *norm = BLASnrm2_(&n_,S,&one);
847:   return(0);
848: }

850: /*
851:  Computes GS orthogonalization   [z;x] - [Sp;Sq]*y,
852:  where y = ([Sp;Sq]'*[z;x]).
853:    k: Column from S to be orthogonalized against previous columns.
854:    Sq = Sp+ld
855:    dim(work)=k;
856: */
857: static PetscErrorCode NEPTOAROrth2(NEP nep,PetscScalar *S,PetscInt ld,PetscInt deg,PetscInt k,PetscScalar *y,PetscReal *norm,PetscBool *lindep,PetscScalar *work)
858: {
860:   PetscBLASInt   n_,lds_,k_,one=1;
861:   PetscScalar    sonem=-1.0,sone=1.0,szero=0.0,*x0,*x,*c;
862:   PetscInt       i,lds=deg*ld,n;
863:   PetscReal      eta,onorm;

866:   BVGetOrthogonalization(nep->V,NULL,NULL,&eta,NULL);
867:   n = k+deg-1;
868:   PetscBLASIntCast(n,&n_);
869:   PetscBLASIntCast(deg*ld,&lds_);
870:   PetscBLASIntCast(k,&k_); /* Number of vectors to orthogonalize against them */
871:   c = work;
872:   x0 = S+k*lds;
873:   PetscStackCall("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S,&lds_,x0,&one,&szero,y,&one));
874:   for (i=1;i<deg;i++) {
875:     x = S+i*ld+k*lds;
876:     PetscStackCall("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+i*ld,&lds_,x,&one,&sone,y,&one));
877:   }
878:   for (i=0;i<deg;i++) {
879:     x= S+i*ld+k*lds;
880:     PetscStackCall("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+i*ld,&lds_,y,&one,&sone,x,&one));
881:   }
882:   NEPTOARSNorm2(lds,S+k*lds,&onorm);
883:   /* twice */
884:   PetscStackCall("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S,&lds_,x0,&one,&szero,c,&one));
885:   for (i=1;i<deg;i++) {
886:     x = S+i*ld+k*lds;
887:     PetscStackCall("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+i*ld,&lds_,x,&one,&sone,c,&one));
888:   }
889:   for (i=0;i<deg;i++) {
890:     x= S+i*ld+k*lds;
891:     PetscStackCall("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+i*ld,&lds_,c,&one,&sone,x,&one));
892:   }
893:   for (i=0;i<k;i++) y[i] += c[i];
894:   if (norm) {
895:     NEPTOARSNorm2(lds,S+k*lds,norm);
896:     if (lindep) *lindep = (*norm < eta * onorm)?PETSC_TRUE:PETSC_FALSE;
897:   }
898:   return(0);
899: }

901: /*
902:   Extend the TOAR basis by applying the the matrix operator
903:   over a vector which is decomposed on the TOAR way
904:   Input:
905:     - S,V: define the latest Arnoldi vector (nv vectors in V)
906:   Output:
907:     - t: new vector extending the TOAR basis
908:     - r: temporally coefficients to compute the TOAR coefficients
909:          for the new Arnoldi vector
910:   Workspace: t_ (two vectors)
911: */
912: static PetscErrorCode NEPTOARExtendBasis(NEP nep,PetscInt idxrktg,PetscScalar *S,PetscInt ls,PetscInt nv,BV W,BV V,Vec t,PetscScalar *r,PetscInt lr,Vec *t_)
913: {
915:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
916:   PetscInt       deg=ctx->nmat-1,k,j;
917:   Vec            v=t_[0],q=t_[1],w;
918:   PetscScalar    *beta=ctx->beta,*s=ctx->s,*xi=ctx->xi,*coeffs,sigma;

921:   if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
922:   sigma = ctx->shifts[idxrktg];
923:   BVSetActiveColumns(nep->V,0,nv);
924:   PetscMalloc1(ctx->nmat,&coeffs);
925:   if (PetscAbsScalar(s[deg-2]-sigma)<100*PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1,"Breakdown in NLEIGS");
926:   /* i-part stored in (i-1) position */
927:   for (j=0;j<nv;j++) {
928:     r[(deg-2)*lr+j] = (S[(deg-2)*ls+j]+(beta[deg-1]/xi[deg-2])*S[(deg-1)*ls+j])/(s[deg-2]-sigma);
929:   }
930:   BVSetActiveColumns(W,0,deg);
931:   BVGetColumn(W,deg-1,&w);
932:   BVMultVec(V,1.0/beta[deg],0,w,S+(deg-1)*ls);
933:   BVRestoreColumn(W,deg-1,&w);
934:   BVGetColumn(W,deg-2,&w);
935:   BVMultVec(V,1.0,0.0,w,r+(deg-2)*lr);
936:   BVRestoreColumn(W,deg-2,&w);
937:   for (k=deg-2;k>0;k--) {
938:     if (PetscAbsScalar(s[k-1]-sigma)<100*PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1,"Breakdown in NLEIGS");
939:     for (j=0;j<nv;j++) r[(k-1)*lr+j] = (S[(k-1)*ls+j]+(beta[k]/xi[k-1])*S[k*ls+j]-beta[k]*(1.0-sigma/xi[k-1])*r[(k)*lr+j])/(s[k-1]-sigma);
940:     BVGetColumn(W,k-1,&w);
941:     BVMultVec(V,1.0,0.0,w,r+(k-1)*lr);
942:     BVRestoreColumn(W,k-1,&w);
943:   }
944:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
945:     for (j=0;j<ctx->nmat-2;j++) coeffs[j] = ctx->coeffD[nep->nt*j];
946:     coeffs[ctx->nmat-2] = ctx->coeffD[nep->nt*(ctx->nmat-1)];
947:     BVMultVec(W,1.0,0.0,v,coeffs);
948:     MatMult(nep->A[0],v,q);
949:     for (k=1;k<nep->nt;k++) {
950:       for (j=0;j<ctx->nmat-2;j++) coeffs[j] = ctx->coeffD[nep->nt*j+k];
951:       coeffs[ctx->nmat-2] = ctx->coeffD[nep->nt*(ctx->nmat-1)+k];
952:       BVMultVec(W,1.0,0,v,coeffs);
953:       MatMult(nep->A[k],v,t);
954:       VecAXPY(q,1.0,t);
955:     }
956:     KSPSolve(ctx->ksp[idxrktg],q,t);
957:     VecScale(t,-1.0);
958:   } else {
959:     for (k=0;k<deg-1;k++) {
960:       BVGetColumn(W,k,&w);
961:       MatMult(ctx->D[k],w,q);
962:       BVRestoreColumn(W,k,&w);
963:       BVInsertVec(W,k,q);
964:     }
965:     BVGetColumn(W,deg-1,&w);
966:     MatMult(ctx->D[deg],w,q);
967:     BVRestoreColumn(W,k,&w);
968:     BVInsertVec(W,k,q);
969:     for (j=0;j<ctx->nmat-1;j++) coeffs[j] = 1.0;
970:     BVMultVec(W,1.0,0.0,q,coeffs);
971:     KSPSolve(ctx->ksp[idxrktg],q,t);
972:     VecScale(t,-1.0);
973:   }
974:   PetscFree(coeffs);
975:   return(0);
976: }

978: /*
979:   Compute TOAR coefficients of the blocks of the new Arnoldi vector computed
980: */
981: static PetscErrorCode NEPTOARCoefficients(NEP nep,PetscScalar sigma,PetscInt nv,PetscScalar *S,PetscInt ls,PetscScalar *r,PetscInt lr,PetscScalar *x,PetscScalar *work)
982: {
984:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
985:   PetscInt       k,j,d=ctx->nmat-1;
986:   PetscScalar    *t=work;

989:   NEPNLEIGSEvalNRTFunct(nep,d-1,sigma,t);
990:   for (k=0;k<d-1;k++) {
991:     for (j=0;j<=nv;j++) r[k*lr+j] += t[k]*x[j];
992:   }
993:   for (j=0;j<=nv;j++) r[(d-1)*lr+j] = t[d-1]*x[j];
994:   return(0);
995: }

997: /*
998:   Compute continuation vector coefficients for the Rational-Krylov run.
999:   dim(work) >= (end-ini)*(end-ini+1) + end+1 + 2*(end-ini+1), dim(t) = end.
1000: */
1001: static PetscErrorCode NEPNLEIGS_RKcontinuation(NEP nep,PetscInt ini,PetscInt end,PetscScalar *K,PetscScalar *H,PetscInt ld,PetscScalar sigma,PetscScalar *S,PetscInt lds,PetscScalar *cont,PetscScalar *t,PetscScalar *work)
1002: {
1003: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(SLEPC_MISSING_LAPACK_LARF)
1005:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/LARF - Lapack routines are unavailable");
1006: #else
1008:   PetscScalar    *x,*W,*tau,sone=1.0,szero=0.0;
1009:   PetscInt       i,j,n1,n,nwu=0;
1010:   PetscBLASInt   info,n_,n1_,one=1,dim,lds_;
1011:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;

1014:   if (!ctx->nshifts || !end) {
1015:     t[0] = 1;
1016:     PetscMemcpy(cont,S+end*lds,lds*sizeof(PetscScalar));
1017:   } else {
1018:     n   = end-ini;
1019:     n1  = n+1;
1020:     x   = work+nwu;
1021:     nwu += end+1;
1022:     tau = work+nwu;
1023:     nwu += n;
1024:     W   = work+nwu;
1025:     nwu += n1*n;
1026:     for (j=ini;j<end;j++) {
1027:       for (i=ini;i<=end;i++) W[(j-ini)*n1+i-ini] = K[j*ld+i] -H[j*ld+i]*sigma;
1028:     }
1029:     PetscBLASIntCast(n,&n_);
1030:     PetscBLASIntCast(n1,&n1_);
1031:     PetscBLASIntCast(end+1,&dim);
1032:     PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&n1_,&n_,W,&n1_,tau,work+nwu,&n1_,&info));
1033:     SlepcCheckLapackInfo("geqrf",info);
1034:     for (i=0;i<end;i++) t[i] = 0.0;
1035:     t[end] = 1.0;
1036:     for (j=n-1;j>=0;j--) {
1037:       for (i=0;i<ini+j;i++) x[i] = 0.0;
1038:       x[ini+j] = 1.0;
1039:       for (i=j+1;i<n1;i++) x[i+ini] = W[i+n1*j];
1040:       tau[j] = PetscConj(tau[j]);
1041:       PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&dim,&one,x,&one,tau+j,t,&dim,work+nwu));
1042:     }
1043:     PetscBLASIntCast(lds,&lds_);
1044:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&n1_,&sone,S,&lds_,t,&one,&szero,cont,&one));
1045:   }
1046:   return(0);
1047: #endif
1048: }

1050: /*
1051:   Compute a run of Arnoldi iterations
1052: */
1053: static PetscErrorCode NEPNLEIGSTOARrun(NEP nep,PetscInt *nq,PetscScalar *S,PetscInt ld,PetscScalar *K,PetscScalar *H,PetscInt ldh,BV W,BV V,PetscInt k,PetscInt *M,PetscBool *breakdown,Vec *t_)
1054: {
1056:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;
1057:   PetscInt       i,j,p,m=*M,lwa,deg=ctx->nmat-1,lds=ld*deg,nqt=*nq;
1058:   Vec            t=t_[0];
1059:   PetscReal      norm;
1060:   PetscScalar    *x,*work,*tt,sigma,*cont;
1061:   PetscBool      lindep;

1064:   lwa = PetscMax(ld,deg)+(m+1)*(m+1)+4*(m+1);
1065:   PetscMalloc4(ld,&x,lwa,&work,m+1,&tt,lds,&cont);
1066:   for (j=k;j<m;j++) {
1067:     sigma = ctx->shifts[(++(ctx->idxrk))%ctx->nshiftsw];

1069:     /* Continuation vector */
1070:     NEPNLEIGS_RKcontinuation(nep,0,j,K,H,ldh,sigma,S,lds,cont,tt,work);

1072:     /* apply operator */
1073:     BVGetColumn(nep->V,nqt,&t);
1074:     NEPTOARExtendBasis(nep,(ctx->idxrk)%ctx->nshiftsw,cont,ld,nqt,W,V,t,S+(j+1)*lds,ld,t_+1);
1075:     BVRestoreColumn(nep->V,nqt,&t);

1077:     /* orthogonalize */
1078:     BVOrthogonalizeColumn(nep->V,nqt,x,&norm,&lindep);
1079:     if (!lindep) {
1080:       x[nqt] = norm;
1081:       BVScaleColumn(nep->V,nqt,1.0/norm);
1082:       nqt++;
1083:     } else x[nqt] = 0.0;

1085:     NEPTOARCoefficients(nep,sigma,*nq,cont,ld,S+(j+1)*lds,ld,x,work);

1087:     /* Level-2 orthogonalization */
1088:     NEPTOAROrth2(nep,S,ld,deg,j+1,H+j*ldh,&norm,breakdown,work);
1089:     H[j+1+ldh*j] = norm;
1090:     if (ctx->nshifts) {
1091:       for (i=0;i<=j;i++) K[i+ldh*j] = sigma*H[i+ldh*j] + tt[i];
1092:       K[j+1+ldh*j] = sigma*H[j+1+ldh*j];
1093:     }
1094:     *nq = nqt;
1095:     if (*breakdown) {
1096:       *M = j+1;
1097:       break;
1098:     }
1099:     for (p=0;p<deg;p++) {
1100:       for (i=0;i<=j+deg;i++) {
1101:         S[i+p*ld+(j+1)*lds] /= norm;
1102:       }
1103:     }
1104:   }
1105:   PetscFree4(x,work,tt,cont);
1106:   return(0);
1107: }

1109: /* dim(work)=5*ld*lds dim(rwork)=6*n */
1110: static PetscErrorCode NEPTOARTrunc(NEP nep,PetscScalar *S,PetscInt ld,PetscInt deg,PetscInt *nq,PetscInt cs1,PetscScalar *work,PetscReal *rwork)
1111: {
1112: #if defined(PETSC_MISSING_LAPACK_GESVD)
1114:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable");
1115: #else
1117:   PetscInt       lwa,nwu=0,nrwu=0;
1118:   PetscInt       j,i,n,lds=deg*ld,rk=0,rs1;
1119:   PetscScalar    *M,*V,*pU,t;
1120:   PetscReal      *sg,tol;
1121:   PetscBLASInt   cs1_,rs1_,cs1tdeg,n_,info,lw_;
1122:   Mat            U;

1125:   rs1 = *nq;
1126:   n = (rs1>deg*cs1)?deg*cs1:rs1;
1127:   lwa = 5*ld*lds;
1128:   M = work+nwu;
1129:   nwu += rs1*cs1*deg;
1130:   sg = rwork+nrwu;
1131:   nrwu += n;
1132:   pU = work+nwu;
1133:   nwu += rs1*n;
1134:   V = work+nwu;
1135:   nwu += deg*cs1*n;
1136:   for (i=0;i<cs1;i++) {
1137:     for (j=0;j<deg;j++) {
1138:       PetscMemcpy(M+(i+j*cs1)*rs1,S+i*lds+j*ld,rs1*sizeof(PetscScalar));
1139:     }
1140:   }
1141:   PetscBLASIntCast(n,&n_);
1142:   PetscBLASIntCast(cs1,&cs1_);
1143:   PetscBLASIntCast(rs1,&rs1_);
1144:   PetscBLASIntCast(cs1*deg,&cs1tdeg);
1145:   PetscBLASIntCast(lwa-nwu,&lw_);
1146: #if !defined (PETSC_USE_COMPLEX)
1147:   PetscStackCall("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1tdeg,M,&rs1_,sg,pU,&rs1_,V,&n_,work+nwu,&lw_,&info));
1148: #else
1149:   PetscStackCall("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1tdeg,M,&rs1_,sg,pU,&rs1_,V,&n_,work+nwu,&lw_,rwork+nrwu,&info));
1150: #endif
1151:   SlepcCheckLapackInfo("gesvd",info);

1153:   /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1154:   MatCreateSeqDense(PETSC_COMM_SELF,rs1,cs1+deg-1,pU,&U);
1155:   BVSetActiveColumns(nep->V,0,rs1);
1156:   BVMultInPlace(nep->V,U,0,cs1+deg-1);
1157:   BVSetActiveColumns(nep->V,0,cs1+deg-1);
1158:   MatDestroy(&U);
1159:   tol = PetscMax(rs1,deg*cs1)*PETSC_MACHINE_EPSILON*sg[0];
1160:   for (i=0;i<PetscMin(n_,cs1tdeg);i++) if (sg[i]>tol) rk++;
1161:   rk = PetscMin(cs1+deg-1,rk);

1163:   /* Update S */
1164:   PetscMemzero(S,lds*ld*sizeof(PetscScalar));
1165:   for (i=0;i<rk;i++) {
1166:     t = sg[i];
1167:     PetscStackCall("BLASscal",BLASscal_(&cs1tdeg,&t,V+i,&n_));
1168:   }
1169:   for (j=0;j<cs1;j++) {
1170:     for (i=0;i<deg;i++) {
1171:       PetscMemcpy(S+j*lds+i*ld,V+(cs1*i+j)*n,rk*sizeof(PetscScalar));
1172:     }
1173:   }
1174:   *nq = rk;
1175:   return(0);
1176: #endif
1177: }

1179: /*
1180:   S <- S*Q
1181:   columns s-s+ncu of S
1182:   rows 0-sr of S
1183:   size(Q) qr x ncu
1184:   dim(work)=sr*ncu
1185: */
1186: static PetscErrorCode NEPTOARSupdate(PetscScalar *S,PetscInt ld,PetscInt deg,PetscInt sr,PetscInt s,PetscInt ncu,PetscInt qr,PetscScalar *Q,PetscInt ldq,PetscScalar *work)
1187: {
1189:   PetscScalar    a=1.0,b=0.0;
1190:   PetscBLASInt   sr_,ncu_,ldq_,lds_,qr_;
1191:   PetscInt       j,lds=deg*ld,i;

1194:   PetscBLASIntCast(sr,&sr_);
1195:   PetscBLASIntCast(qr,&qr_);
1196:   PetscBLASIntCast(ncu,&ncu_);
1197:   PetscBLASIntCast(lds,&lds_);
1198:   PetscBLASIntCast(ldq,&ldq_);
1199:   for (i=0;i<deg;i++) {
1200:     PetscStackCall("BLASgemm",BLASgemm_("N","N",&sr_,&ncu_,&qr_,&a,S+i*ld,&lds_,Q,&ldq_,&b,work,&sr_));
1201:     for (j=0;j<ncu;j++) {
1202:       PetscMemcpy(S+lds*(s+j)+i*ld,work+j*sr,sr*sizeof(PetscScalar));
1203:     }
1204:   }
1205:   return(0);
1206: }

1208: PetscErrorCode NEPSolve_NLEIGS(NEP nep)
1209: {
1211:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;
1212:   PetscInt       i,j,k=0,l,nv=0,ld,lds,off,ldds,rs1,nq=0,newn;
1213:   PetscInt       lwa,lrwa,nwu=0,nrwu=0,deg=ctx->nmat-1,nconv=0;
1214:   PetscScalar    *S,*Q,*work,*H,*pU,*K,betak=0,*Hc,*eigr,*eigi;
1215:   PetscReal      betah,norm,*rwork;
1216:   PetscBool      breakdown=PETSC_FALSE,lindep;
1217:   Mat            U;
1218:   BV             W;

1221:   ld = nep->ncv+deg;
1222:   lds = deg*ld;
1223:   lwa = (deg+6)*ld*lds;
1224:   lrwa = 7*lds;
1225:   DSGetLeadingDimension(nep->ds,&ldds);
1226:   PetscMalloc4(lwa,&work,lrwa,&rwork,lds*ld,&S,ldds*ldds,&Hc);
1227:   PetscMemzero(S,lds*ld*sizeof(PetscScalar));
1228:   if (!ctx->nshifts) {
1229:     PetscMalloc2(nep->ncv,&eigr,nep->ncv,&eigi);
1230:   } else { eigr = nep->eigr; eigi = nep->eigi; }
1231:   BVDuplicateResize(nep->V,PetscMax(nep->nt-1,ctx->nmat-1),&W);

1233:   /* Get the starting vector */
1234:   for (i=0;i<deg;i++) {
1235:     BVSetRandomColumn(nep->V,i);
1236:     BVOrthogonalizeColumn(nep->V,i,S+i*ld,&norm,&lindep);
1237:     if (!lindep) {
1238:       BVScaleColumn(nep->V,i,1/norm);
1239:       S[i+i*ld] = norm;
1240:       nq++;
1241:     }
1242:   }
1243:   if (!nq) SETERRQ(PetscObjectComm((PetscObject)nep),1,"NEP: Problem with initial vector");
1244:   NEPTOARSNorm2(lds,S,&norm);
1245:   for (j=0;j<deg;j++) {
1246:     for (i=0;i<=j;i++) S[i+j*ld] /= norm;
1247:   }

1249:   /* Restart loop */
1250:   l = 0;
1251:   while (nep->reason == NEP_CONVERGED_ITERATING) {
1252:     nep->its++;

1254:     /* Compute an nv-step Krylov relation */
1255:     nv = PetscMin(nep->nconv+nep->mpd,nep->ncv);
1256:     if (ctx->nshifts) { DSGetArray(nep->ds,DS_MAT_A,&K); }
1257:     DSGetArray(nep->ds,ctx->nshifts?DS_MAT_B:DS_MAT_A,&H);
1258:     NEPNLEIGSTOARrun(nep,&nq,S,ld,K,H,ldds,W,nep->V,nep->nconv+l,&nv,&breakdown,nep->work);
1259:     betah = PetscAbsScalar(H[(nv-1)*ldds+nv]);
1260:     DSRestoreArray(nep->ds,ctx->nshifts?DS_MAT_B:DS_MAT_A,&H);
1261:     if (ctx->nshifts) {
1262:       betak = K[(nv-1)*ldds+nv];
1263:       DSRestoreArray(nep->ds,DS_MAT_A,&K);
1264:     }
1265:     DSSetDimensions(nep->ds,nv,0,nep->nconv,nep->nconv+l);
1266:     if (l==0) {
1267:       DSSetState(nep->ds,DS_STATE_INTERMEDIATE);
1268:     } else {
1269:       DSSetState(nep->ds,DS_STATE_RAW);
1270:     }

1272:     /* Solve projected problem */
1273:     if (ctx->nshifts) {
1274:       DSGetArray(nep->ds,DS_MAT_B,&H);
1275:       PetscMemcpy(Hc,H,ldds*ldds*sizeof(PetscScalar));
1276:       DSRestoreArray(nep->ds,DS_MAT_B,&H);
1277:     }
1278:     DSSolve(nep->ds,nep->eigr,nep->eigi);
1279:     DSSort(nep->ds,nep->eigr,nep->eigi,NULL,NULL,NULL);
1280:     if (!ctx->nshifts) {
1281:       DSUpdateExtraRow(nep->ds);
1282:     }
1283:     DSSynchronize(nep->ds,nep->eigr,nep->eigi);

1285:     /* Check convergence */
1286:     NEPNLEIGSKrylovConvergence(nep,S,ld,nq,Hc,PETSC_FALSE,nep->nconv,nv-nep->nconv,betak,betah,&k,nep->work);
1287:     (*nep->stopping)(nep,nep->its,nep->max_it,k,nep->nev,&nep->reason,nep->stoppingctx);
1288:     nconv = k;

1290:     /* Update l */
1291:     if (nep->reason != NEP_CONVERGED_ITERATING || breakdown) l = 0;
1292:     else {
1293:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
1294:       if (!breakdown) {
1295:         /* Prepare the Rayleigh quotient for restart */
1296:         if (!ctx->nshifts) {
1297:           DSTruncate(nep->ds,k+l);
1298:           DSGetDimensions(nep->ds,&newn,NULL,NULL,NULL,NULL);
1299:           l = newn-k;
1300:         } else {
1301:           DSGetArray(nep->ds,DS_MAT_Q,&Q);
1302:           DSGetArray(nep->ds,DS_MAT_B,&H);
1303:           DSGetArray(nep->ds,DS_MAT_A,&K);
1304:           for (i=ctx->lock?k:0;i<k+l;i++) {
1305:             H[k+l+i*ldds] = betah*Q[nv-1+i*ldds];
1306:             K[k+l+i*ldds] = betak*Q[nv-1+i*ldds];
1307:           }
1308:           DSRestoreArray(nep->ds,DS_MAT_B,&H);
1309:           DSRestoreArray(nep->ds,DS_MAT_A,&K);
1310:           DSRestoreArray(nep->ds,DS_MAT_Q,&Q);
1311:           DSSetDimensions(nep->ds,k+l,0,nep->nconv,0);
1312:         }
1313:       }
1314:     }
1315:     if (!ctx->lock && l>0) { l += k; k = 0; }

1317:     /* Update S */
1318:     off = nep->nconv*ldds;
1319:     DSGetArray(nep->ds,ctx->nshifts?DS_MAT_Z:DS_MAT_Q,&Q);
1320:     NEPTOARSupdate(S,ld,deg,nq,nep->nconv,k+l-nep->nconv,nv,Q+off,ldds,work+nwu);
1321:     DSRestoreArray(nep->ds,ctx->nshifts?DS_MAT_Z:DS_MAT_Q,&Q);

1323:     /* Copy last column of S */
1324:     PetscMemcpy(S+lds*(k+l),S+lds*nv,lds*sizeof(PetscScalar));
1325:     if (nep->reason == NEP_CONVERGED_ITERATING) {
1326:       if (breakdown) {

1328:         /* Stop if breakdown */
1329:         PetscInfo2(nep,"Breakdown (it=%D norm=%g)\n",nep->its,(double)betah);
1330:         nep->reason = NEP_DIVERGED_BREAKDOWN;
1331:       } else {
1332:         /* Truncate S */
1333:         NEPTOARTrunc(nep,S,ld,deg,&nq,k+l+1,work+nwu,rwork+nrwu);
1334:       }
1335:     }
1336:     nep->nconv = k;
1337:     if (!ctx->nshifts) {
1338:       for (i=0;i<nv;i++) { eigr[i] = nep->eigr[i]; eigi[i] = nep->eigi[i]; }
1339:       NEPNLEIGSBackTransform((PetscObject)nep,nv,eigr,eigi);
1340:     }
1341:     NEPMonitor(nep,nep->its,nconv,eigr,eigi,nep->errest,nv);
1342:   }
1343:   nep->nconv = nconv;
1344:   if (nep->nconv>0) {
1345:     /* Extract invariant pair */
1346:     NEPTOARTrunc(nep,S,ld,deg,&nq,nep->nconv,work+nwu,rwork+nrwu);
1347:     /* Update vectors V = V*S or V=V*S*H */
1348:     rs1 = nep->nconv;
1349:     if (ctx->nshifts) {
1350:       DSGetArray(nep->ds,DS_MAT_B,&H);
1351:       NEPTOARSupdate(S,ld,deg,rs1,0,nep->nconv,nep->nconv,H,ldds,work+nwu);
1352:       DSRestoreArray(nep->ds,DS_MAT_B,&H);
1353:     }
1354:     PetscMalloc1(rs1*nep->nconv,&pU);
1355:     for (i=0;i<nep->nconv;i++) {
1356:       PetscMemcpy(pU+i*rs1,S+i*lds,rs1*sizeof(PetscScalar));
1357:     }
1358:     MatCreateSeqDense(PETSC_COMM_SELF,rs1,nep->nconv,pU,&U);
1359:     BVSetActiveColumns(nep->V,0,rs1);
1360:     BVMultInPlace(nep->V,U,0,nep->nconv);
1361:     BVSetActiveColumns(nep->V,0,nep->nconv);
1362:     MatDestroy(&U);
1363:     PetscFree(pU);
1364:   }
1365:   /* truncate Schur decomposition and change the state to raw so that
1366:      DSVectors() computes eigenvectors from scratch */
1367:   DSSetDimensions(nep->ds,nep->nconv,0,0,0);
1368:   DSSetState(nep->ds,DS_STATE_RAW);

1370:   PetscFree4(work,rwork,S,Hc);
1371:   /* Map eigenvalues back to the original problem */
1372:   if (!ctx->nshifts) {
1373:     NEPNLEIGSBackTransform((PetscObject)nep,nep->nconv,nep->eigr,nep->eigi);
1374:     PetscFree2(eigr,eigi);
1375:   }
1376:   BVDestroy(&W);
1377:   return(0);
1378: }

1380: static PetscErrorCode NEPNLEIGSSetSingularitiesFunction_NLEIGS(NEP nep,PetscErrorCode (*fun)(NEP,PetscInt*,PetscScalar*,void*),void *ctx)
1381: {
1382:   NEP_NLEIGS *nepctx=(NEP_NLEIGS*)nep->data;

1385:   if (fun) nepctx->computesingularities = fun;
1386:   if (ctx) nepctx->singularitiesctx     = ctx;
1387:   nep->state = NEP_STATE_INITIAL;
1388:   return(0);
1389: }

1391: /*@C
1392:    NEPNLEIGSSetSingularitiesFunction - Sets a user function to compute a discretization
1393:    of the singularity set (where T(.) is not analytic).

1395:    Logically Collective on NEP

1397:    Input Parameters:
1398: +  nep - the NEP context
1399: .  fun - user function (if NULL then NEP retains any previously set value)
1400: -  ctx - [optional] user-defined context for private data for the function
1401:          (may be NULL, in which case NEP retains any previously set value)

1403:    Calling Sequence of fun:
1404: $   fun(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *ctx)

1406: +   nep   - the NEP context
1407: .   maxnp - on input number of requested points in the discretization (can be set)
1408: .   xi    - computed values of the discretization
1409: -   ctx   - optional context, as set by NEPNLEIGSSetSingularitiesFunction()

1411:    Notes:
1412:    The user-defined function can set a smaller value of maxnp if necessary.
1413:    It is wrong to return a larger value.

1415:    If the problem type has been set to rational with NEPSetProblemType(),
1416:    then it is not necessary to set the singularities explicitly since the
1417:    solver will try to determine them automatically.

1419:    Level: intermediate

1421: .seealso: NEPNLEIGSGetSingularitiesFunction(), NEPSetProblemType()
1422: @*/
1423: PetscErrorCode NEPNLEIGSSetSingularitiesFunction(NEP nep,PetscErrorCode (*fun)(NEP,PetscInt*,PetscScalar*,void*),void *ctx)
1424: {

1429:   PetscTryMethod(nep,"NEPNLEIGSSetSingularitiesFunction_C",(NEP,PetscErrorCode(*)(NEP,PetscInt*,PetscScalar*,void*),void*),(nep,fun,ctx));
1430:   return(0);
1431: }

1433: static PetscErrorCode NEPNLEIGSGetSingularitiesFunction_NLEIGS(NEP nep,PetscErrorCode (**fun)(NEP,PetscInt*,PetscScalar*,void*),void **ctx)
1434: {
1435:   NEP_NLEIGS *nepctx=(NEP_NLEIGS*)nep->data;

1438:   if (fun) *fun = nepctx->computesingularities;
1439:   if (ctx) *ctx = nepctx->singularitiesctx;
1440:   return(0);
1441: }

1443: /*@C
1444:    NEPNLEIGSGetSingularitiesFunction - Returns the Function and optionally the user
1445:    provided context for computing a discretization of the singularity set.

1447:    Not Collective

1449:    Input Parameter:
1450: .  nep - the nonlinear eigensolver context

1452:    Output Parameters:
1453: +  fun - location to put the function (or NULL)
1454: -  ctx - location to stash the function context (or NULL)

1456:    Level: advanced

1458: .seealso: NEPNLEIGSSetSingularitiesFunction()
1459: @*/
1460: PetscErrorCode NEPNLEIGSGetSingularitiesFunction(NEP nep,PetscErrorCode (**fun)(NEP,PetscInt*,PetscScalar*,void*),void **ctx)
1461: {

1466:   PetscUseMethod(nep,"NEPNLEIGSGetSingularitiesFunction_C",(NEP,PetscErrorCode(**)(NEP,PetscInt*,PetscScalar*,void*),void**),(nep,fun,ctx));
1467:   return(0);
1468: }

1470: static PetscErrorCode NEPNLEIGSSetRestart_NLEIGS(NEP nep,PetscReal keep)
1471: {
1472:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1475:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
1476:   else {
1477:     if (keep<0.1 || keep>0.9) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
1478:     ctx->keep = keep;
1479:   }
1480:   return(0);
1481: }

1483: /*@
1484:    NEPNLEIGSSetRestart - Sets the restart parameter for the NLEIGS
1485:    method, in particular the proportion of basis vectors that must be kept
1486:    after restart.

1488:    Logically Collective on NEP

1490:    Input Parameters:
1491: +  nep  - the nonlinear eigensolver context
1492: -  keep - the number of vectors to be kept at restart

1494:    Options Database Key:
1495: .  -nep_nleigs_restart - Sets the restart parameter

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

1500:    Level: advanced

1502: .seealso: NEPNLEIGSGetRestart()
1503: @*/
1504: PetscErrorCode NEPNLEIGSSetRestart(NEP nep,PetscReal keep)
1505: {

1511:   PetscTryMethod(nep,"NEPNLEIGSSetRestart_C",(NEP,PetscReal),(nep,keep));
1512:   return(0);
1513: }

1515: static PetscErrorCode NEPNLEIGSGetRestart_NLEIGS(NEP nep,PetscReal *keep)
1516: {
1517:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1520:   *keep = ctx->keep;
1521:   return(0);
1522: }

1524: /*@
1525:    NEPNLEIGSGetRestart - Gets the restart parameter used in the NLEIGS method.

1527:    Not Collective

1529:    Input Parameter:
1530: .  nep - the nonlinear eigensolver context

1532:    Output Parameter:
1533: .  keep - the restart parameter

1535:    Level: advanced

1537: .seealso: NEPNLEIGSSetRestart()
1538: @*/
1539: PetscErrorCode NEPNLEIGSGetRestart(NEP nep,PetscReal *keep)
1540: {

1546:   PetscUseMethod(nep,"NEPNLEIGSGetRestart_C",(NEP,PetscReal*),(nep,keep));
1547:   return(0);
1548: }

1550: static PetscErrorCode NEPNLEIGSSetLocking_NLEIGS(NEP nep,PetscBool lock)
1551: {
1552:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1555:   ctx->lock = lock;
1556:   return(0);
1557: }

1559: /*@
1560:    NEPNLEIGSSetLocking - Choose between locking and non-locking variants of
1561:    the NLEIGS method.

1563:    Logically Collective on NEP

1565:    Input Parameters:
1566: +  nep  - the nonlinear eigensolver context
1567: -  lock - true if the locking variant must be selected

1569:    Options Database Key:
1570: .  -nep_nleigs_locking - Sets the locking flag

1572:    Notes:
1573:    The default is to lock converged eigenpairs when the method restarts.
1574:    This behaviour can be changed so that all directions are kept in the
1575:    working subspace even if already converged to working accuracy (the
1576:    non-locking variant).

1578:    Level: advanced

1580: .seealso: NEPNLEIGSGetLocking()
1581: @*/
1582: PetscErrorCode NEPNLEIGSSetLocking(NEP nep,PetscBool lock)
1583: {

1589:   PetscTryMethod(nep,"NEPNLEIGSSetLocking_C",(NEP,PetscBool),(nep,lock));
1590:   return(0);
1591: }

1593: static PetscErrorCode NEPNLEIGSGetLocking_NLEIGS(NEP nep,PetscBool *lock)
1594: {
1595:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1598:   *lock = ctx->lock;
1599:   return(0);
1600: }

1602: /*@
1603:    NEPNLEIGSGetLocking - Gets the locking flag used in the NLEIGS method.

1605:    Not Collective

1607:    Input Parameter:
1608: .  nep - the nonlinear eigensolver context

1610:    Output Parameter:
1611: .  lock - the locking flag

1613:    Level: advanced

1615: .seealso: NEPNLEIGSSetLocking()
1616: @*/
1617: PetscErrorCode NEPNLEIGSGetLocking(NEP nep,PetscBool *lock)
1618: {

1624:   PetscUseMethod(nep,"NEPNLEIGSGetLocking_C",(NEP,PetscBool*),(nep,lock));
1625:   return(0);
1626: }

1628: static PetscErrorCode NEPNLEIGSSetInterpolation_NLEIGS(NEP nep,PetscReal tol,PetscInt degree)
1629: {
1631:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;

1634:   if (tol == PETSC_DEFAULT) {
1635:     ctx->ddtol = PETSC_DEFAULT;
1636:     nep->state = NEP_STATE_INITIAL;
1637:   } else {
1638:     if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1639:     ctx->ddtol = tol;
1640:   }
1641:   if (degree == PETSC_DEFAULT || degree == PETSC_DECIDE) {
1642:     ctx->ddmaxit = 0;
1643:     if (nep->state) { NEPReset(nep); }
1644:     nep->state = NEP_STATE_INITIAL;
1645:   } else {
1646:     if (degree <= 0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of degree. Must be > 0");
1647:     if (ctx->ddmaxit != degree) {
1648:       ctx->ddmaxit = degree;
1649:       if (nep->state) { NEPReset(nep); }
1650:       nep->state = NEP_STATE_INITIAL;
1651:     }
1652:   }
1653:   return(0);
1654: }

1656: /*@
1657:    NEPNLEIGSSetInterpolation - Sets the tolerance and maximum degree
1658:    when building the interpolation via divided differences.

1660:    Logically Collective on NEP

1662:    Input Parameters:
1663: +  nep    - the nonlinear eigensolver context
1664: .  tol    - tolerance to stop computing divided differences
1665: -  degree - maximum degree of interpolation

1667:    Options Database Key:
1668: +  -nep_nleigs_interpolation_tol <tol> - Sets the tolerance to stop computing divided differences
1669: -  -nep_nleigs_interpolation_degree <degree> - Sets the maximum degree of interpolation

1671:    Notes:
1672:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

1674:    Level: advanced

1676: .seealso: NEPNLEIGSGetInterpolation()
1677: @*/
1678: PetscErrorCode NEPNLEIGSSetInterpolation(NEP nep,PetscReal tol,PetscInt degree)
1679: {

1686:   PetscTryMethod(nep,"NEPNLEIGSSetInterpolation_C",(NEP,PetscReal,PetscInt),(nep,tol,degree));
1687:   return(0);
1688: }

1690: static PetscErrorCode NEPNLEIGSGetInterpolation_NLEIGS(NEP nep,PetscReal *tol,PetscInt *degree)
1691: {
1692:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1695:   if (tol)    *tol    = ctx->ddtol;
1696:   if (degree) *degree = ctx->ddmaxit;
1697:   return(0);
1698: }

1700: /*@
1701:    NEPNLEIGSGetInterpolation - Gets the tolerance and maximum degree
1702:    when building the interpolation via divided differences.

1704:    Not Collective

1706:    Input Parameter:
1707: .  nep - the nonlinear eigensolver context

1709:    Output Parameter:
1710: +  tol    - tolerance to stop computing divided differences
1711: -  degree - maximum degree of interpolation

1713:    Level: advanced

1715: .seealso: NEPNLEIGSSetInterpolation()
1716: @*/
1717: PetscErrorCode NEPNLEIGSGetInterpolation(NEP nep,PetscReal *tol,PetscInt *degree)
1718: {

1723:   PetscTryMethod(nep,"NEPNLEIGSGetInterpolation_C",(NEP,PetscReal*,PetscInt*),(nep,tol,degree));
1724:   return(0);
1725: }

1727: static PetscErrorCode NEPNLEIGSSetRKShifts_NLEIGS(NEP nep,PetscInt ns,PetscScalar *shifts)
1728: {
1730:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
1731:   PetscInt       i;

1734:   if (ns<=0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"Number of shifts must be positive");
1735:   if (ctx->nshifts) { PetscFree(ctx->shifts); }
1736:   for (i=0;i<ctx->nshiftsw;i++) { KSPDestroy(&ctx->ksp[i]); }
1737:   PetscFree(ctx->ksp);
1738:   ctx->ksp = NULL;
1739:   PetscMalloc1(ns,&ctx->shifts);
1740:   for (i=0;i<ns;i++) ctx->shifts[i] = shifts[i];
1741:   ctx->nshifts = ns;
1742:   nep->state   = NEP_STATE_INITIAL;
1743:   return(0);
1744: }

1746: /*@C
1747:    NEPNLEIGSSetRKShifts - Sets a list of shifts to be used in the Rational
1748:    Krylov method.

1750:    Logically Collective on NEP

1752:    Input Parameters:
1753: +  nep    - the nonlinear eigensolver context
1754: .  ns     - number of shifts
1755: -  shifts - array of scalar values specifying the shifts

1757:    Options Database Key:
1758: .  -nep_nleigs_rk_shifts - Sets the list of shifts

1760:    Notes:
1761:    If only one shift is provided, the built subspace built is equivalent to
1762:    shift-and-invert Krylov-Schur (provided that the absolute convergence
1763:    criterion is used).

1765:    In the case of real scalars, complex shifts are not allowed. In the
1766:    command line, a comma-separated list of complex values can be provided with
1767:    the format [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
1768:    -nep_nleigs_rk_shifts 1.0+2.0i,1.5+2.0i,1.0+1.5i

1770:    Level: advanced

1772: .seealso: NEPNLEIGSGetRKShifts()
1773: @*/
1774: PetscErrorCode NEPNLEIGSSetRKShifts(NEP nep,PetscInt ns,PetscScalar *shifts)
1775: {

1782:   PetscTryMethod(nep,"NEPNLEIGSSetRKShifts_C",(NEP,PetscInt,PetscScalar*),(nep,ns,shifts));
1783:   return(0);
1784: }

1786: static PetscErrorCode NEPNLEIGSGetRKShifts_NLEIGS(NEP nep,PetscInt *ns,PetscScalar **shifts)
1787: {
1789:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
1790:   PetscInt       i;

1793:   *ns = ctx->nshifts;
1794:   if (ctx->nshifts) {
1795:     PetscMalloc1(ctx->nshifts,shifts);
1796:     for (i=0;i<ctx->nshifts;i++) (*shifts)[i] = ctx->shifts[i];
1797:   }
1798:   return(0);
1799: }

1801: /*@C
1802:    NEPNLEIGSGetRKShifts - Gets the list of shifts used in the Rational
1803:    Krylov method.

1805:    Not Collective

1807:    Input Parameter:
1808: .  nep - the nonlinear eigensolver context

1810:    Output Parameter:
1811: +  ns     - number of shifts
1812: -  shifts - array of shifts

1814:    Level: advanced

1816: .seealso: NEPNLEIGSSetRKShifts()
1817: @*/
1818: PetscErrorCode NEPNLEIGSGetRKShifts(NEP nep,PetscInt *ns,PetscScalar **shifts)
1819: {

1826:   PetscTryMethod(nep,"NEPNLEIGSGetRKShifts_C",(NEP,PetscInt*,PetscScalar**),(nep,ns,shifts));
1827:   return(0);
1828: }

1830: #define SHIFTMAX 30

1832: PetscErrorCode NEPSetFromOptions_NLEIGS(PetscOptionItems *PetscOptionsObject,NEP nep)
1833: {
1835:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;
1836:   PetscInt       i,k;
1837:   PetscBool      flg1,flg2,b;
1838:   PetscReal      r;
1839:   PetscScalar    array[SHIFTMAX];
1840:   PC             pc;
1841:   PCType         pctype;
1842:   KSPType        ksptype;

1845:   PetscOptionsHead(PetscOptionsObject,"NEP NLEIGS Options");

1847:     PetscOptionsReal("-nep_nleigs_restart","Proportion of vectors kept after restart","NEPNLEIGSSetRestart",0.5,&r,&flg1);
1848:     if (flg1) { NEPNLEIGSSetRestart(nep,r); }

1850:     PetscOptionsBool("-nep_nleigs_locking","Choose between locking and non-locking variants","NEPNLEIGSSetLocking",PETSC_FALSE,&b,&flg1);
1851:     if (flg1) { NEPNLEIGSSetLocking(nep,b); }

1853:     NEPNLEIGSGetInterpolation(nep,&r,&i);
1854:     if (!i) i = PETSC_DEFAULT;
1855:     PetscOptionsInt("-nep_nleigs_interpolation_degree","Maximum number of terms for interpolation via divided differences","NEPNLEIGSSetInterpolation",i,&i,&flg1);
1856:     PetscOptionsReal("-nep_nleigs_interpolation_tol","Tolerance for interpolation via divided differences","NEPNLEIGSSetInterpolation",r,&r,&flg2);
1857:     if (flg1 || flg2) { NEPNLEIGSSetInterpolation(nep,r,i); }

1859:     k = SHIFTMAX;
1860:     for (i=0;i<k;i++) array[i] = 0;
1861:     PetscOptionsScalarArray("-nep_nleigs_rk_shifts","Shifts for Rational Krylov","NEPNLEIGSSetRKShifts",array,&k,&flg1);
1862:     if (flg1) { NEPNLEIGSSetRKShifts(nep,k,array); }

1864:   PetscOptionsTail();

1866:   if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
1867:   for (i=0;i<ctx->nshiftsw;i++) {
1868:     KSPGetPC(ctx->ksp[i],&pc);
1869:     KSPGetType(ctx->ksp[i],&ksptype);
1870:     PCGetType(pc,&pctype);
1871:     if (!pctype && !ksptype) {
1872:       KSPSetType(ctx->ksp[i],KSPPREONLY);
1873:       PCSetType(pc,PCLU);
1874:     }
1875:     KSPSetFromOptions(ctx->ksp[i]);
1876:   }
1877:   return(0);
1878: }

1880: static PetscErrorCode NEPNLEIGSGetKSPs_NLEIGS(NEP nep,KSP **ksp)
1881: {
1883:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;
1884:   PetscInt       i;

1887:   if (!ctx->ksp) {
1888:     NEPNLEIGSSetShifts(nep);
1889:     PetscMalloc1(ctx->nshiftsw,&ctx->ksp);
1890:     for (i=0;i<ctx->nshiftsw;i++) {
1891:       KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp[i]);
1892:       KSPSetOptionsPrefix(ctx->ksp[i],((PetscObject)nep)->prefix);
1893:       KSPAppendOptionsPrefix(ctx->ksp[i],"nep_nleigs_");
1894:       PetscObjectIncrementTabLevel((PetscObject)ctx->ksp[i],(PetscObject)nep,1);
1895:       PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp[i]);
1896:       KSPSetErrorIfNotConverged(ctx->ksp[i],PETSC_TRUE);
1897:       KSPSetTolerances(ctx->ksp[i],SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
1898:     }
1899:   }
1900:   *ksp = ctx->ksp;
1901:   return(0);
1902: }

1904: /*@C
1905:    NEPNLEIGSGetKSPs - Retrieve the array of linear solver objects associated with
1906:    the nonlinear eigenvalue solver.

1908:    Not Collective

1910:    Input Parameter:
1911: .  nep - nonlinear eigenvalue solver

1913:    Output Parameter:
1914: .  ksp - array of linear solver object

1916:    Level: advanced
1917: @*/
1918: PetscErrorCode NEPNLEIGSGetKSPs(NEP nep,KSP **ksp)
1919: {

1925:   PetscUseMethod(nep,"NEPNLEIGSGetKSPs_C",(NEP,KSP**),(nep,ksp));
1926:   return(0);
1927: }

1929: PetscErrorCode NEPView_NLEIGS(NEP nep,PetscViewer viewer)
1930: {
1932:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
1933:   PetscBool      isascii;
1934:   PetscInt       i;
1935:   char           str[50];

1938:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1939:   if (isascii) {
1940:     PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
1941:     PetscViewerASCIIPrintf(viewer,"  using the %slocking variant\n",ctx->lock?"":"non-");
1942:     PetscViewerASCIIPrintf(viewer,"  divided difference terms: used=%D, max=%D\n",ctx->nmat,ctx->ddmaxit);
1943:     PetscViewerASCIIPrintf(viewer,"  tolerance for divided difference convergence: %g\n",(double)ctx->ddtol);
1944:     if (ctx->nshifts) {
1945:       PetscViewerASCIIPrintf(viewer,"  RK shifts: ");
1946:       PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1947:       for (i=0;i<ctx->nshifts;i++) {
1948:         SlepcSNPrintfScalar(str,50,ctx->shifts[i],PETSC_FALSE);
1949:         PetscViewerASCIIPrintf(viewer,"%s%s",str,(i<ctx->nshifts-1)?",":"");
1950:       }
1951:       PetscViewerASCIIPrintf(viewer,"\n");
1952:       PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1953:     }
1954:     if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
1955:     PetscViewerASCIIPushTab(viewer);
1956:     KSPView(ctx->ksp[0],viewer);
1957:     PetscViewerASCIIPopTab(viewer);
1958:   }
1959:   return(0);
1960: }

1962: PetscErrorCode NEPReset_NLEIGS(NEP nep)
1963: {
1965:   PetscInt       k;
1966:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;

1969:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
1970:     PetscFree(ctx->coeffD);
1971:   } else {
1972:     for (k=0;k<ctx->nmat;k++) { MatDestroy(&ctx->D[k]); }
1973:   }
1974:   PetscFree4(ctx->s,ctx->xi,ctx->beta,ctx->D);
1975:   for (k=0;k<ctx->nshiftsw;k++) { KSPReset(ctx->ksp[k]); }
1976:   if (ctx->vrn) {
1977:     VecDestroy(&ctx->vrn);
1978:   }
1979:   return(0);
1980: }

1982: PetscErrorCode NEPDestroy_NLEIGS(NEP nep)
1983: {
1985:   PetscInt       k;
1986:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;

1989:   for (k=0;k<ctx->nshiftsw;k++) { KSPDestroy(&ctx->ksp[k]); }
1990:   PetscFree(ctx->ksp);
1991:   if (ctx->nshifts) { PetscFree(ctx->shifts); }
1992:   PetscFree(nep->data);
1993:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetSingularitiesFunction_C",NULL);
1994:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetSingularitiesFunction_C",NULL);
1995:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRestart_C",NULL);
1996:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRestart_C",NULL);
1997:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetLocking_C",NULL);
1998:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetLocking_C",NULL);
1999:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetInterpolation_C",NULL);
2000:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetInterpolation_C",NULL);
2001:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRKShifts_C",NULL);
2002:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRKShifts_C",NULL);
2003:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetKSPs_C",NULL);
2004:   return(0);
2005: }

2007: PETSC_EXTERN PetscErrorCode NEPCreate_NLEIGS(NEP nep)
2008: {
2010:   NEP_NLEIGS     *ctx;

2013:   PetscNewLog(nep,&ctx);
2014:   nep->data  = (void*)ctx;
2015:   ctx->lock  = PETSC_TRUE;
2016:   ctx->ddtol = PETSC_DEFAULT;

2018:   nep->ops->solve          = NEPSolve_NLEIGS;
2019:   nep->ops->setup          = NEPSetUp_NLEIGS;
2020:   nep->ops->setfromoptions = NEPSetFromOptions_NLEIGS;
2021:   nep->ops->view           = NEPView_NLEIGS;
2022:   nep->ops->destroy        = NEPDestroy_NLEIGS;
2023:   nep->ops->reset          = NEPReset_NLEIGS;
2024:   nep->ops->computevectors = NEPComputeVectors_Schur;

2026:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetSingularitiesFunction_C",NEPNLEIGSSetSingularitiesFunction_NLEIGS);
2027:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetSingularitiesFunction_C",NEPNLEIGSGetSingularitiesFunction_NLEIGS);
2028:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRestart_C",NEPNLEIGSSetRestart_NLEIGS);
2029:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRestart_C",NEPNLEIGSGetRestart_NLEIGS);
2030:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetLocking_C",NEPNLEIGSSetLocking_NLEIGS);
2031:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetLocking_C",NEPNLEIGSGetLocking_NLEIGS);
2032:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetInterpolation_C",NEPNLEIGSSetInterpolation_NLEIGS);
2033:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetInterpolation_C",NEPNLEIGSGetInterpolation_NLEIGS);
2034:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRKShifts_C",NEPNLEIGSSetRKShifts_NLEIGS);
2035:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRKShifts_C",NEPNLEIGSGetRKShifts_NLEIGS);
2036:   PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetKSPs_C",NEPNLEIGSGetKSPs_NLEIGS);
2037:   return(0);
2038: }