Actual source code: nleigs.c

slepc-3.19.2 2023-09-05
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc 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>
 27: #include <slepcblaslapack.h>
 28: #include "nleigs.h"

 30: PetscErrorCode NEPNLEIGSBackTransform(PetscObject ob,PetscInt n,PetscScalar *valr,PetscScalar *vali)
 31: {
 32:   NEP         nep;
 33:   PetscInt    j;
 34: #if !defined(PETSC_USE_COMPLEX)
 35:   PetscScalar t;
 36: #endif

 38:   PetscFunctionBegin;
 39:   nep = (NEP)ob;
 40: #if !defined(PETSC_USE_COMPLEX)
 41:   for (j=0;j<n;j++) {
 42:     if (vali[j] == 0) valr[j] = 1.0 / valr[j] + nep->target;
 43:     else {
 44:       t = valr[j] * valr[j] + vali[j] * vali[j];
 45:       valr[j] = valr[j] / t + nep->target;
 46:       vali[j] = - vali[j] / t;
 47:     }
 48:   }
 49: #else
 50:   for (j=0;j<n;j++) {
 51:     valr[j] = 1.0 / valr[j] + nep->target;
 52:   }
 53: #endif
 54:   PetscFunctionReturn(PETSC_SUCCESS);
 55: }

 57: /* Computes the roots of a polynomial */
 58: static PetscErrorCode NEPNLEIGSAuxiliarPRootFinder(PetscInt deg,PetscScalar *polcoeffs,PetscScalar *wr,PetscScalar *wi,PetscBool *avail)
 59: {
 60:   PetscScalar    *C;
 61:   PetscBLASInt   n_,lwork;
 62:   PetscInt       i;
 63: #if defined(PETSC_USE_COMPLEX)
 64:   PetscReal      *rwork=NULL;
 65: #endif
 66:   PetscScalar    *work;
 67:   PetscBLASInt   info;

 69:   PetscFunctionBegin;
 70:   *avail = PETSC_TRUE;
 71:   if (deg>0) {
 72:     PetscCall(PetscCalloc1(deg*deg,&C));
 73:     PetscCall(PetscBLASIntCast(deg,&n_));
 74:     for (i=0;i<deg-1;i++) {
 75:       C[(deg+1)*i+1]   = 1.0;
 76:       C[(deg-1)*deg+i] = -polcoeffs[deg-i]/polcoeffs[0];
 77:     }
 78:     C[deg*deg+-1] = -polcoeffs[1]/polcoeffs[0];
 79:     PetscCall(PetscBLASIntCast(3*deg,&lwork));

 81:     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
 82: #if !defined(PETSC_USE_COMPLEX)
 83:     PetscCall(PetscMalloc1(lwork,&work));
 84:     PetscCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&n_,C,&n_,wr,wi,NULL,&n_,NULL,&n_,work,&lwork,&info));
 85:     if (info) *avail = PETSC_FALSE;
 86:     PetscCall(PetscFree(work));
 87: #else
 88:     PetscCall(PetscMalloc2(2*deg,&rwork,lwork,&work));
 89:     PetscCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&n_,C,&n_,wr,NULL,&n_,NULL,&n_,work,&lwork,rwork,&info));
 90:     if (info) *avail = PETSC_FALSE;
 91:     PetscCall(PetscFree2(rwork,work));
 92: #endif
 93:     PetscCall(PetscFPTrapPop());
 94:     PetscCall(PetscFree(C));
 95:   }
 96:   PetscFunctionReturn(PETSC_SUCCESS);
 97: }

 99: static PetscErrorCode NEPNLEIGSAuxiliarRmDuplicates(PetscInt nin,PetscScalar *pin,PetscInt *nout,PetscScalar *pout,PetscInt max)
100: {
101:   PetscInt i,j;

103:   PetscFunctionBegin;
104:   for (i=0;i<nin;i++) {
105:     if (max && *nout>=max) break;
106:     pout[(*nout)++] = pin[i];
107:     for (j=0;j<*nout-1;j++)
108:       if (PetscAbsScalar(pin[i]-pout[j])<PETSC_MACHINE_EPSILON*100) {
109:         (*nout)--;
110:         break;
111:       }
112:   }
113:   PetscFunctionReturn(PETSC_SUCCESS);
114: }

116: static PetscErrorCode NEPNLEIGSFNSingularities(FN f,PetscInt *nisol,PetscScalar **isol,PetscBool *rational)
117: {
118:   FNCombineType  ctype;
119:   FN             f1,f2;
120:   PetscInt       i,nq,nisol1,nisol2;
121:   PetscScalar    *qcoeff,*wr,*wi,*isol1,*isol2;
122:   PetscBool      flg,avail,rat1,rat2;

124:   PetscFunctionBegin;
125:   *rational = PETSC_FALSE;
126:   PetscCall(PetscObjectTypeCompare((PetscObject)f,FNRATIONAL,&flg));
127:   if (flg) {
128:     *rational = PETSC_TRUE;
129:     PetscCall(FNRationalGetDenominator(f,&nq,&qcoeff));
130:     if (nq>1) {
131:       PetscCall(PetscMalloc2(nq-1,&wr,nq-1,&wi));
132:       PetscCall(NEPNLEIGSAuxiliarPRootFinder(nq-1,qcoeff,wr,wi,&avail));
133:       if (avail) {
134:         PetscCall(PetscCalloc1(nq-1,isol));
135:         *nisol = 0;
136:         for (i=0;i<nq-1;i++)
137: #if !defined(PETSC_USE_COMPLEX)
138:           if (wi[i]==0)
139: #endif
140:             (*isol)[(*nisol)++] = wr[i];
141:         nq = *nisol; *nisol = 0;
142:         for (i=0;i<nq;i++) wr[i] = (*isol)[i];
143:         PetscCall(NEPNLEIGSAuxiliarRmDuplicates(nq,wr,nisol,*isol,0));
144:         PetscCall(PetscFree2(wr,wi));
145:       } else { *nisol=0; *isol = NULL; }
146:     } else { *nisol = 0; *isol = NULL; }
147:     PetscCall(PetscFree(qcoeff));
148:   }
149:   PetscCall(PetscObjectTypeCompare((PetscObject)f,FNCOMBINE,&flg));
150:   if (flg) {
151:     PetscCall(FNCombineGetChildren(f,&ctype,&f1,&f2));
152:     if (ctype != FN_COMBINE_COMPOSE && ctype != FN_COMBINE_DIVIDE) {
153:       PetscCall(NEPNLEIGSFNSingularities(f1,&nisol1,&isol1,&rat1));
154:       PetscCall(NEPNLEIGSFNSingularities(f2,&nisol2,&isol2,&rat2));
155:       if (nisol1+nisol2>0) {
156:         PetscCall(PetscCalloc1(nisol1+nisol2,isol));
157:         *nisol = 0;
158:         PetscCall(NEPNLEIGSAuxiliarRmDuplicates(nisol1,isol1,nisol,*isol,0));
159:         PetscCall(NEPNLEIGSAuxiliarRmDuplicates(nisol2,isol2,nisol,*isol,0));
160:       }
161:       *rational = (rat1&&rat2)?PETSC_TRUE:PETSC_FALSE;
162:       PetscCall(PetscFree(isol1));
163:       PetscCall(PetscFree(isol2));
164:     }
165:   }
166:   PetscFunctionReturn(PETSC_SUCCESS);
167: }

169: static PetscErrorCode NEPNLEIGSRationalSingularities(NEP nep,PetscInt *ndptx,PetscScalar *dxi,PetscBool *rational)
170: {
171:   PetscInt       nt,i,nisol;
172:   FN             f;
173:   PetscScalar    *isol;
174:   PetscBool      rat;

176:   PetscFunctionBegin;
177:   *rational = PETSC_TRUE;
178:   *ndptx = 0;
179:   PetscCall(NEPGetSplitOperatorInfo(nep,&nt,NULL));
180:   for (i=0;i<nt;i++) {
181:     PetscCall(NEPGetSplitOperatorTerm(nep,i,NULL,&f));
182:     PetscCall(NEPNLEIGSFNSingularities(f,&nisol,&isol,&rat));
183:     if (nisol) {
184:       PetscCall(NEPNLEIGSAuxiliarRmDuplicates(nisol,isol,ndptx,dxi,0));
185:       PetscCall(PetscFree(isol));
186:     }
187:     *rational = ((*rational)&&rat)?PETSC_TRUE:PETSC_FALSE;
188:   }
189:   PetscFunctionReturn(PETSC_SUCCESS);
190: }

192: /*  Adaptive Anderson-Antoulas algorithm */
193: static PetscErrorCode NEPNLEIGSAAAComputation(NEP nep,PetscInt ndpt,PetscScalar *ds,PetscScalar *F,PetscInt *ndptx,PetscScalar *dxi)
194: {
195:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
196:   PetscScalar    mean=0.0,*z,*f,*C,*A,*VT,*work,*ww,szero=0.0,sone=1.0;
197:   PetscScalar    *N,*D;
198:   PetscReal      *S,norm,err,*R;
199:   PetscInt       i,k,j,idx=0,cont;
200:   PetscBLASInt   n_,m_,lda_,lwork,info,one=1;
201: #if defined(PETSC_USE_COMPLEX)
202:   PetscReal      *rwork;
203: #endif

205:   PetscFunctionBegin;
206:   PetscCall(PetscBLASIntCast(8*ndpt,&lwork));
207:   PetscCall(PetscMalloc5(ndpt,&R,ndpt,&z,ndpt,&f,ndpt*ndpt,&C,ndpt,&ww));
208:   PetscCall(PetscMalloc6(ndpt*ndpt,&A,ndpt,&S,ndpt*ndpt,&VT,lwork,&work,ndpt,&D,ndpt,&N));
209: #if defined(PETSC_USE_COMPLEX)
210:   PetscCall(PetscMalloc1(8*ndpt,&rwork));
211: #endif
212:   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
213:   norm = 0.0;
214:   for (i=0;i<ndpt;i++) {
215:     mean += F[i];
216:     norm = PetscMax(PetscAbsScalar(F[i]),norm);
217:   }
218:   mean /= ndpt;
219:   PetscCall(PetscBLASIntCast(ndpt,&lda_));
220:   for (i=0;i<ndpt;i++) R[i] = PetscAbsScalar(F[i]-mean);
221:   /* next support point */
222:   err = 0.0;
223:   for (i=0;i<ndpt;i++) if (R[i]>=err) {idx = i; err = R[i];}
224:   for (k=0;k<ndpt-1;k++) {
225:     z[k] = ds[idx]; f[k] = F[idx]; R[idx] = -1.0;
226:     /* next column of Cauchy matrix */
227:     for (i=0;i<ndpt;i++) {
228:       C[i+k*ndpt] = 1.0/(ds[i]-ds[idx]);
229:     }

231:     PetscCall(PetscArrayzero(A,ndpt*ndpt));
232:     cont = 0;
233:     for (i=0;i<ndpt;i++) {
234:       if (R[i]!=-1.0) {
235:         for (j=0;j<=k;j++)A[cont+j*ndpt] = C[i+j*ndpt]*F[i]-C[i+j*ndpt]*f[j];
236:         cont++;
237:       }
238:     }
239:     PetscCall(PetscBLASIntCast(cont,&m_));
240:     PetscCall(PetscBLASIntCast(k+1,&n_));
241: #if defined(PETSC_USE_COMPLEX)
242:     PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","A",&m_,&n_,A,&lda_,S,NULL,&lda_,VT,&lda_,work,&lwork,rwork,&info));
243: #else
244:     PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","A",&m_,&n_,A,&lda_,S,NULL,&lda_,VT,&lda_,work,&lwork,&info));
245: #endif
246:     SlepcCheckLapackInfo("gesvd",info);
247:     for (i=0;i<=k;i++) {
248:       ww[i] = PetscConj(VT[i*ndpt+k]);
249:       D[i] = ww[i]*f[i];
250:     }
251:     PetscCallBLAS("BLASgemv",BLASgemv_("N",&lda_,&n_,&sone,C,&lda_,D,&one,&szero,N,&one));
252:     PetscCallBLAS("BLASgemv",BLASgemv_("N",&lda_,&n_,&sone,C,&lda_,ww,&one,&szero,D,&one));
253:     for (i=0;i<ndpt;i++) if (R[i]>=0) R[i] = PetscAbsScalar(F[i]-N[i]/D[i]);
254:     /* next support point */
255:     err = 0.0;
256:     for (i=0;i<ndpt;i++) if (R[i]>=err) {idx = i; err = R[i];}
257:     if (err <= ctx->ddtol*norm) break;
258:   }

260:   PetscCheck(k<ndpt-1,PetscObjectComm((PetscObject)nep),PETSC_ERR_CONV_FAILED,"Failed to determine singularities automatically in general problem");
261:   /* poles */
262:   PetscCall(PetscArrayzero(C,ndpt*ndpt));
263:   PetscCall(PetscArrayzero(A,ndpt*ndpt));
264:   for (i=0;i<=k;i++) {
265:     C[i+ndpt*i] = 1.0;
266:     A[(i+1)*ndpt] = ww[i];
267:     A[i+1] = 1.0;
268:     A[i+1+(i+1)*ndpt] = z[i];
269:   }
270:   C[0] = 0.0; C[k+1+(k+1)*ndpt] = 1.0;
271:   n_++;
272: #if defined(PETSC_USE_COMPLEX)
273:   PetscCallBLAS("LAPACKggev",LAPACKggev_("N","N",&n_,A,&lda_,C,&lda_,D,N,NULL,&lda_,NULL,&lda_,work,&lwork,rwork,&info));
274: #else
275:   PetscCallBLAS("LAPACKggev",LAPACKggev_("N","N",&n_,A,&lda_,C,&lda_,D,VT,N,NULL,&lda_,NULL,&lda_,work,&lwork,&info));
276: #endif
277:   SlepcCheckLapackInfo("ggev",info);
278:   cont = 0.0;
279:   for (i=0;i<n_;i++) if (N[i]!=0.0) {
280:     dxi[cont++] = D[i]/N[i];
281:   }
282:   *ndptx = cont;
283:   PetscCall(PetscFPTrapPop());
284:   PetscCall(PetscFree5(R,z,f,C,ww));
285:   PetscCall(PetscFree6(A,S,VT,work,D,N));
286: #if defined(PETSC_USE_COMPLEX)
287:   PetscCall(PetscFree(rwork));
288: #endif
289:   PetscFunctionReturn(PETSC_SUCCESS);
290: }

292: /*  Singularities using Adaptive Anderson-Antoulas algorithm */
293: static PetscErrorCode NEPNLEIGSAAASingularities(NEP nep,PetscInt ndpt,PetscScalar *ds,PetscInt *ndptx,PetscScalar *dxi)
294: {
295:   Vec            u,v,w;
296:   PetscRandom    rand=NULL;
297:   PetscScalar    *F,*isol;
298:   PetscInt       i,k,nisol,nt;
299:   Mat            T;
300:   FN             f;

302:   PetscFunctionBegin;
303:   PetscCall(PetscMalloc1(ndpt,&F));
304:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
305:     PetscCall(PetscMalloc1(ndpt,&isol));
306:     *ndptx = 0;
307:     PetscCall(NEPGetSplitOperatorInfo(nep,&nt,NULL));
308:     nisol = *ndptx;
309:     for (k=0;k<nt;k++) {
310:       PetscCall(NEPGetSplitOperatorTerm(nep,k,NULL,&f));
311:       for (i=0;i<ndpt;i++) PetscCall(FNEvaluateFunction(f,ds[i],&F[i]));
312:       PetscCall(NEPNLEIGSAAAComputation(nep,ndpt,ds,F,&nisol,isol));
313:       if (nisol) PetscCall(NEPNLEIGSAuxiliarRmDuplicates(nisol,isol,ndptx,dxi,ndpt));
314:     }
315:     PetscCall(PetscFree(isol));
316:   } else {
317:     PetscCall(MatCreateVecs(nep->function,&u,NULL));
318:     PetscCall(VecDuplicate(u,&v));
319:     PetscCall(VecDuplicate(u,&w));
320:     if (nep->V) PetscCall(BVGetRandomContext(nep->V,&rand));
321:     PetscCall(VecSetRandom(u,rand));
322:     PetscCall(VecNormalize(u,NULL));
323:     PetscCall(VecSetRandom(v,rand));
324:     PetscCall(VecNormalize(v,NULL));
325:     T = nep->function;
326:     for (i=0;i<ndpt;i++) {
327:       PetscCall(NEPComputeFunction(nep,ds[i],T,T));
328:       PetscCall(MatMult(T,v,w));
329:       PetscCall(VecDot(w,u,&F[i]));
330:     }
331:     PetscCall(NEPNLEIGSAAAComputation(nep,ndpt,ds,F,ndptx,dxi));
332:     PetscCall(VecDestroy(&u));
333:     PetscCall(VecDestroy(&v));
334:     PetscCall(VecDestroy(&w));
335:   }
336:   PetscCall(PetscFree(F));
337:   PetscFunctionReturn(PETSC_SUCCESS);
338: }

340: static PetscErrorCode NEPNLEIGSLejaBagbyPoints(NEP nep)
341: {
342:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
343:   PetscInt       i,k,ndpt=NDPOINTS,ndptx=NDPOINTS;
344:   PetscScalar    *ds,*dsi,*dxi,*nrs,*nrxi,*s=ctx->s,*xi=ctx->xi,*beta=ctx->beta;
345:   PetscReal      maxnrs,minnrxi;
346:   PetscBool      rational;
347: #if !defined(PETSC_USE_COMPLEX)
348:   PetscReal      a,b,h;
349: #endif

351:   PetscFunctionBegin;
352:   if (!ctx->computesingularities && nep->problem_type!=NEP_RATIONAL) ndpt = ndptx = LBPOINTS;
353:   PetscCall(PetscMalloc5(ndpt+1,&ds,ndpt+1,&dsi,ndpt,&dxi,ndpt+1,&nrs,ndpt,&nrxi));

355:   /* Discretize the target region boundary */
356:   PetscCall(RGComputeContour(nep->rg,ndpt,ds,dsi));
357: #if !defined(PETSC_USE_COMPLEX)
358:   for (i=0;i<ndpt;i++) if (dsi[i]!=0.0) break;
359:   if (i<ndpt) {
360:     PetscCheck(nep->problem_type==NEP_RATIONAL,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NLEIGS with real arithmetic requires the target set to be included in the real axis");
361:     /* Select a segment in the real axis */
362:     PetscCall(RGComputeBoundingBox(nep->rg,&a,&b,NULL,NULL));
363:     PetscCheck(a>-PETSC_MAX_REAL && b<PETSC_MAX_REAL,PetscObjectComm((PetscObject)nep),PETSC_ERR_USER_INPUT,"NLEIGS requires a bounded target set");
364:     h = (b-a)/ndpt;
365:     for (i=0;i<ndpt;i++) {ds[i] = a+h*i; dsi[i] = 0.0;}
366:   }
367: #endif
368:   /* Discretize the singularity region */
369:   if (ctx->computesingularities) PetscCall((ctx->computesingularities)(nep,&ndptx,dxi,ctx->singularitiesctx));
370:   else {
371:     if (nep->problem_type==NEP_RATIONAL) {
372:       PetscCall(NEPNLEIGSRationalSingularities(nep,&ndptx,dxi,&rational));
373:       PetscCheck(rational,PetscObjectComm((PetscObject)nep),PETSC_ERR_CONV_FAILED,"Failed to determine singularities automatically in rational problem; consider solving the problem as general");
374:     } else {
375:       /* AAA algorithm */
376:       PetscCall(NEPNLEIGSAAASingularities(nep,ndpt,ds,&ndptx,dxi));
377:     }
378:   }
379:   /* Look for Leja-Bagby points in the discretization sets */
380:   s[0]  = ds[0];
381:   xi[0] = (ndptx>0)?dxi[0]:PETSC_INFINITY;
382:   PetscCheck(PetscAbsScalar(xi[0])>=10*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)nep),PETSC_ERR_USER_INPUT,"Singularity point 0 is nearly zero: %g; consider removing the singularity or shifting the problem",(double)PetscAbsScalar(xi[0]));
383:   beta[0] = 1.0; /* scaling factors are also computed here */
384:   for (i=0;i<ndpt;i++) {
385:     nrs[i] = 1.0;
386:     nrxi[i] = 1.0;
387:   }
388:   for (k=1;k<ctx->ddmaxit;k++) {
389:     maxnrs = 0.0;
390:     minnrxi = PETSC_MAX_REAL;
391:     for (i=0;i<ndpt;i++) {
392:       nrs[i] *= ((ds[i]-s[k-1])/(1.0-ds[i]/xi[k-1]))/beta[k-1];
393:       if (PetscAbsScalar(nrs[i])>maxnrs) {maxnrs = PetscAbsScalar(nrs[i]); s[k] = ds[i];}
394:     }
395:     if (ndptx>k) {
396:       for (i=1;i<ndptx;i++) {
397:         nrxi[i] *= ((dxi[i]-s[k-1])/(1.0-dxi[i]/xi[k-1]))/beta[k-1];
398:         if (PetscAbsScalar(nrxi[i])<minnrxi) {minnrxi = PetscAbsScalar(nrxi[i]); xi[k] = dxi[i];}
399:       }
400:       PetscCheck(PetscAbsScalar(xi[k])>=10*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)nep),PETSC_ERR_USER_INPUT,"Singularity point %" PetscInt_FMT " is nearly zero: %g; consider removing the singularity or shifting the problem",k,(double)PetscAbsScalar(xi[k]));
401:     } else xi[k] = PETSC_INFINITY;
402:     beta[k] = maxnrs;
403:   }
404:   PetscCall(PetscFree5(ds,dsi,dxi,nrs,nrxi));
405:   PetscFunctionReturn(PETSC_SUCCESS);
406: }

408: PetscErrorCode NEPNLEIGSEvalNRTFunct(NEP nep,PetscInt k,PetscScalar sigma,PetscScalar *b)
409: {
410:   NEP_NLEIGS  *ctx=(NEP_NLEIGS*)nep->data;
411:   PetscInt    i;
412:   PetscScalar *beta=ctx->beta,*s=ctx->s,*xi=ctx->xi;

414:   PetscFunctionBegin;
415:   b[0] = 1.0/beta[0];
416:   for (i=0;i<k;i++) {
417:     b[i+1] = ((sigma-s[i])*b[i])/(beta[i+1]*(1.0-sigma/xi[i]));
418:   }
419:   PetscFunctionReturn(PETSC_SUCCESS);
420: }

422: static PetscErrorCode MatMult_Fun(Mat A,Vec x,Vec y)
423: {
424:   NEP_NLEIGS_MATSHELL *ctx;
425:   PetscInt            i;

427:   PetscFunctionBeginUser;
428:   PetscCall(MatShellGetContext(A,&ctx));
429:   PetscCall(MatMult(ctx->A[0],x,y));
430:   if (ctx->coeff[0]!=1.0) PetscCall(VecScale(y,ctx->coeff[0]));
431:   for (i=1;i<ctx->nmat;i++) {
432:     PetscCall(MatMult(ctx->A[i],x,ctx->t));
433:     PetscCall(VecAXPY(y,ctx->coeff[i],ctx->t));
434:   }
435:   PetscFunctionReturn(PETSC_SUCCESS);
436: }

438: static PetscErrorCode MatMultTranspose_Fun(Mat A,Vec x,Vec y)
439: {
440:   NEP_NLEIGS_MATSHELL *ctx;
441:   PetscInt            i;

443:   PetscFunctionBeginUser;
444:   PetscCall(MatShellGetContext(A,&ctx));
445:   PetscCall(MatMultTranspose(ctx->A[0],x,y));
446:   if (ctx->coeff[0]!=1.0) PetscCall(VecScale(y,ctx->coeff[0]));
447:   for (i=1;i<ctx->nmat;i++) {
448:     PetscCall(MatMultTranspose(ctx->A[i],x,ctx->t));
449:     PetscCall(VecAXPY(y,ctx->coeff[i],ctx->t));
450:   }
451:   PetscFunctionReturn(PETSC_SUCCESS);
452: }

454: static PetscErrorCode MatGetDiagonal_Fun(Mat A,Vec diag)
455: {
456:   NEP_NLEIGS_MATSHELL *ctx;
457:   PetscInt            i;

459:   PetscFunctionBeginUser;
460:   PetscCall(MatShellGetContext(A,&ctx));
461:   PetscCall(MatGetDiagonal(ctx->A[0],diag));
462:   if (ctx->coeff[0]!=1.0) PetscCall(VecScale(diag,ctx->coeff[0]));
463:   for (i=1;i<ctx->nmat;i++) {
464:     PetscCall(MatGetDiagonal(ctx->A[i],ctx->t));
465:     PetscCall(VecAXPY(diag,ctx->coeff[i],ctx->t));
466:   }
467:   PetscFunctionReturn(PETSC_SUCCESS);
468: }

470: static PetscErrorCode MatDuplicate_Fun(Mat A,MatDuplicateOption op,Mat *B)
471: {
472:   PetscInt            m,n,M,N,i;
473:   NEP_NLEIGS_MATSHELL *ctxnew,*ctx;
474:   void                (*fun)(void);

476:   PetscFunctionBeginUser;
477:   PetscCall(MatShellGetContext(A,&ctx));
478:   PetscCall(PetscNew(&ctxnew));
479:   ctxnew->nmat = ctx->nmat;
480:   ctxnew->maxnmat = ctx->maxnmat;
481:   PetscCall(PetscMalloc2(ctxnew->maxnmat,&ctxnew->A,ctxnew->maxnmat,&ctxnew->coeff));
482:   for (i=0;i<ctx->nmat;i++) {
483:     PetscCall(PetscObjectReference((PetscObject)ctx->A[i]));
484:     ctxnew->A[i] = ctx->A[i];
485:     ctxnew->coeff[i] = ctx->coeff[i];
486:   }
487:   PetscCall(MatGetSize(ctx->A[0],&M,&N));
488:   PetscCall(MatGetLocalSize(ctx->A[0],&m,&n));
489:   PetscCall(VecDuplicate(ctx->t,&ctxnew->t));
490:   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)A),m,n,M,N,(void*)ctxnew,B));
491:   PetscCall(MatShellSetManageScalingShifts(*B));
492:   PetscCall(MatShellGetOperation(A,MATOP_MULT,&fun));
493:   PetscCall(MatShellSetOperation(*B,MATOP_MULT,fun));
494:   PetscCall(MatShellGetOperation(A,MATOP_MULT_TRANSPOSE,&fun));
495:   PetscCall(MatShellSetOperation(*B,MATOP_MULT_TRANSPOSE,fun));
496:   PetscCall(MatShellGetOperation(A,MATOP_GET_DIAGONAL,&fun));
497:   PetscCall(MatShellSetOperation(*B,MATOP_GET_DIAGONAL,fun));
498:   PetscCall(MatShellGetOperation(A,MATOP_DUPLICATE,&fun));
499:   PetscCall(MatShellSetOperation(*B,MATOP_DUPLICATE,fun));
500:   PetscCall(MatShellGetOperation(A,MATOP_DESTROY,&fun));
501:   PetscCall(MatShellSetOperation(*B,MATOP_DESTROY,fun));
502:   PetscCall(MatShellGetOperation(A,MATOP_AXPY,&fun));
503:   PetscCall(MatShellSetOperation(*B,MATOP_AXPY,fun));
504:   PetscFunctionReturn(PETSC_SUCCESS);
505: }

507: static PetscErrorCode MatDestroy_Fun(Mat A)
508: {
509:   NEP_NLEIGS_MATSHELL *ctx;
510:   PetscInt            i;

512:   PetscFunctionBeginUser;
513:   if (A) {
514:     PetscCall(MatShellGetContext(A,&ctx));
515:     for (i=0;i<ctx->nmat;i++) PetscCall(MatDestroy(&ctx->A[i]));
516:     PetscCall(VecDestroy(&ctx->t));
517:     PetscCall(PetscFree2(ctx->A,ctx->coeff));
518:     PetscCall(PetscFree(ctx));
519:   }
520:   PetscFunctionReturn(PETSC_SUCCESS);
521: }

523: static PetscErrorCode MatAXPY_Fun(Mat Y,PetscScalar a,Mat X,MatStructure str)
524: {
525:   NEP_NLEIGS_MATSHELL *ctxY,*ctxX;
526:   PetscInt            i,j;
527:   PetscBool           found;

529:   PetscFunctionBeginUser;
530:   PetscCall(MatShellGetContext(Y,&ctxY));
531:   PetscCall(MatShellGetContext(X,&ctxX));
532:   for (i=0;i<ctxX->nmat;i++) {
533:     found = PETSC_FALSE;
534:     for (j=0;!found&&j<ctxY->nmat;j++) {
535:       if (ctxX->A[i]==ctxY->A[j]) {
536:         found = PETSC_TRUE;
537:         ctxY->coeff[j] += a*ctxX->coeff[i];
538:       }
539:     }
540:     if (!found) {
541:       ctxY->coeff[ctxY->nmat] = a*ctxX->coeff[i];
542:       ctxY->A[ctxY->nmat++] = ctxX->A[i];
543:       PetscCall(PetscObjectReference((PetscObject)ctxX->A[i]));
544:     }
545:   }
546:   PetscFunctionReturn(PETSC_SUCCESS);
547: }

549: static PetscErrorCode MatScale_Fun(Mat M,PetscScalar a)
550: {
551:   NEP_NLEIGS_MATSHELL *ctx;
552:   PetscInt            i;

554:   PetscFunctionBeginUser;
555:   PetscCall(MatShellGetContext(M,&ctx));
556:   for (i=0;i<ctx->nmat;i++) ctx->coeff[i] *= a;
557:   PetscFunctionReturn(PETSC_SUCCESS);
558: }

560: static PetscErrorCode NLEIGSMatToMatShellArray(Mat A,Mat *Ms,PetscInt maxnmat)
561: {
562:   NEP_NLEIGS_MATSHELL *ctx;
563:   PetscInt            m,n,M,N;
564:   PetscBool           has;

566:   PetscFunctionBegin;
567:   PetscCall(MatHasOperation(A,MATOP_DUPLICATE,&has));
568:   PetscCheck(has,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"MatDuplicate operation required");
569:   PetscCall(PetscNew(&ctx));
570:   ctx->maxnmat = maxnmat;
571:   PetscCall(PetscMalloc2(ctx->maxnmat,&ctx->A,ctx->maxnmat,&ctx->coeff));
572:   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&ctx->A[0]));
573:   ctx->nmat = 1;
574:   ctx->coeff[0] = 1.0;
575:   PetscCall(MatCreateVecs(A,&ctx->t,NULL));
576:   PetscCall(MatGetSize(A,&M,&N));
577:   PetscCall(MatGetLocalSize(A,&m,&n));
578:   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)A),m,n,M,N,(void*)ctx,Ms));
579:   PetscCall(MatShellSetManageScalingShifts(*Ms));
580:   PetscCall(MatShellSetOperation(*Ms,MATOP_MULT,(void(*)(void))MatMult_Fun));
581:   PetscCall(MatShellSetOperation(*Ms,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Fun));
582:   PetscCall(MatShellSetOperation(*Ms,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Fun));
583:   PetscCall(MatShellSetOperation(*Ms,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_Fun));
584:   PetscCall(MatShellSetOperation(*Ms,MATOP_DESTROY,(void(*)(void))MatDestroy_Fun));
585:   PetscCall(MatShellSetOperation(*Ms,MATOP_AXPY,(void(*)(void))MatAXPY_Fun));
586:   PetscCall(MatShellSetOperation(*Ms,MATOP_SCALE,(void(*)(void))MatScale_Fun));
587:   PetscFunctionReturn(PETSC_SUCCESS);
588: }

590: /*
591:    MatIsShellAny - returns true if any of the n matrices is a shell matrix
592:  */
593: static PetscErrorCode MatIsShellAny(Mat *A,PetscInt n,PetscBool *shell)
594: {
595:   PetscInt       i;
596:   PetscBool      flg;

598:   PetscFunctionBegin;
599:   *shell = PETSC_FALSE;
600:   for (i=0;i<n;i++) {
601:     PetscCall(MatIsShell(A[i],&flg));
602:     if (flg) { *shell = PETSC_TRUE; break; }
603:   }
604:   PetscFunctionReturn(PETSC_SUCCESS);
605: }

607: static PetscErrorCode NEPNLEIGSDividedDifferences_split(NEP nep)
608: {
609:   PetscErrorCode ierr;
610:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
611:   PetscInt       k,j,i,maxnmat,nmax;
612:   PetscReal      norm0,norm,*matnorm;
613:   PetscScalar    *s=ctx->s,*beta=ctx->beta,*xi=ctx->xi,*b,alpha,*coeffs,*pK,*pH,sone=1.0;
614:   Mat            T,P,Ts,K,H;
615:   PetscBool      shell,hasmnorm=PETSC_FALSE,matrix=PETSC_TRUE;
616:   PetscBLASInt   n_;

618:   PetscFunctionBegin;
619:   nmax = ctx->ddmaxit;
620:   PetscCall(PetscMalloc1(nep->nt*nmax,&ctx->coeffD));
621:   PetscCall(PetscMalloc3(nmax+1,&b,nmax+1,&coeffs,nep->nt,&matnorm));
622:   for (j=0;j<nep->nt;j++) {
623:     PetscCall(MatHasOperation(nep->A[j],MATOP_NORM,&hasmnorm));
624:     if (!hasmnorm) break;
625:     PetscCall(MatNorm(nep->A[j],NORM_INFINITY,matnorm+j));
626:   }
627:   /* Try matrix functions scheme */
628:   PetscCall(PetscCalloc2(nmax*nmax,&pK,nmax*nmax,&pH));
629:   for (i=0;i<nmax-1;i++) {
630:     pK[(nmax+1)*i]   = 1.0;
631:     pK[(nmax+1)*i+1] = beta[i+1]/xi[i];
632:     pH[(nmax+1)*i]   = s[i];
633:     pH[(nmax+1)*i+1] = beta[i+1];
634:   }
635:   pH[nmax*nmax-1] = s[nmax-1];
636:   pK[nmax*nmax-1] = 1.0;
637:   PetscCall(PetscBLASIntCast(nmax,&n_));
638:   PetscCallBLAS("BLAStrsm",BLAStrsm_("R","L","N","U",&n_,&n_,&sone,pK,&n_,pH,&n_));
639:   /* The matrix to be used is in H. K will be a work-space matrix */
640:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nmax,nmax,pH,&H));
641:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nmax,nmax,pK,&K));
642:   for (j=0;matrix&&j<nep->nt;j++) {
643:     PetscCall(PetscPushErrorHandler(PetscReturnErrorHandler,NULL));
644:     ierr = FNEvaluateFunctionMat(nep->f[j],H,K);
645:     PetscCall(PetscPopErrorHandler());
646:     if (!ierr) {
647:       for (i=0;i<nmax;i++) ctx->coeffD[j+i*nep->nt] = pK[i]*beta[0];
648:     } else {
649:       matrix = PETSC_FALSE;
650:       PetscCall(PetscFPTrapPop());
651:     }
652:   }
653:   PetscCall(MatDestroy(&H));
654:   PetscCall(MatDestroy(&K));
655:   if (!matrix) {
656:     for (j=0;j<nep->nt;j++) {
657:       PetscCall(FNEvaluateFunction(nep->f[j],s[0],ctx->coeffD+j));
658:       ctx->coeffD[j] *= beta[0];
659:     }
660:   }
661:   if (hasmnorm) {
662:     norm0 = 0.0;
663:     for (j=0;j<nep->nt;j++) norm0 += matnorm[j]*PetscAbsScalar(ctx->coeffD[j]);
664:   } else {
665:     norm0 = 0.0;
666:     for (j=0;j<nep->nt;j++) norm0 = PetscMax(PetscAbsScalar(ctx->coeffD[j]),norm0);
667:   }
668:   ctx->nmat = ctx->ddmaxit;
669:   for (k=1;k<ctx->ddmaxit;k++) {
670:     if (!matrix) {
671:       PetscCall(NEPNLEIGSEvalNRTFunct(nep,k,s[k],b));
672:       for (i=0;i<nep->nt;i++) {
673:         PetscCall(FNEvaluateFunction(nep->f[i],s[k],ctx->coeffD+k*nep->nt+i));
674:         for (j=0;j<k;j++) {
675:           ctx->coeffD[k*nep->nt+i] -= b[j]*ctx->coeffD[i+nep->nt*j];
676:         }
677:         ctx->coeffD[k*nep->nt+i] /= b[k];
678:       }
679:     }
680:     if (hasmnorm) {
681:       norm = 0.0;
682:       for (j=0;j<nep->nt;j++) norm += matnorm[j]*PetscAbsScalar(ctx->coeffD[k*nep->nt+j]);
683:     } else {
684:       norm = 0.0;
685:       for (j=0;j<nep->nt;j++) norm = PetscMax(PetscAbsScalar(ctx->coeffD[k*nep->nt+j]),norm);
686:     }
687:     if (k>1 && norm/norm0 < ctx->ddtol) {
688:       ctx->nmat = k+1;
689:       break;
690:     }
691:   }
692:   if (!ctx->ksp) PetscCall(NEPNLEIGSGetKSPs(nep,&ctx->nshiftsw,&ctx->ksp));
693:   PetscCall(MatIsShellAny(nep->A,nep->nt,&shell));
694:   maxnmat = PetscMax(ctx->ddmaxit,nep->nt);
695:   for (i=0;i<ctx->nshiftsw;i++) {
696:     PetscCall(NEPNLEIGSEvalNRTFunct(nep,ctx->nmat-1,ctx->shifts[i],coeffs));
697:     if (!shell) PetscCall(MatDuplicate(nep->A[0],MAT_COPY_VALUES,&T));
698:     else PetscCall(NLEIGSMatToMatShellArray(nep->A[0],&T,maxnmat));
699:     if (nep->P) { /* user-defined preconditioner */
700:       PetscCall(MatDuplicate(nep->P[0],MAT_COPY_VALUES,&P));
701:     } else P=T;
702:     alpha = 0.0;
703:     for (j=0;j<ctx->nmat;j++) alpha += coeffs[j]*ctx->coeffD[j*nep->nt];
704:     PetscCall(MatScale(T,alpha));
705:     if (nep->P) PetscCall(MatScale(P,alpha));
706:     for (k=1;k<nep->nt;k++) {
707:       alpha = 0.0;
708:       for (j=0;j<ctx->nmat;j++) alpha += coeffs[j]*ctx->coeffD[j*nep->nt+k];
709:       if (shell) PetscCall(NLEIGSMatToMatShellArray(nep->A[k],&Ts,maxnmat));
710:       PetscCall(MatAXPY(T,alpha,shell?Ts:nep->A[k],nep->mstr));
711:       if (nep->P) PetscCall(MatAXPY(P,alpha,nep->P[k],nep->mstrp));
712:       if (shell) PetscCall(MatDestroy(&Ts));
713:     }
714:     PetscCall(NEP_KSPSetOperators(ctx->ksp[i],T,P));
715:     PetscCall(KSPSetUp(ctx->ksp[i]));
716:     PetscCall(MatDestroy(&T));
717:     if (nep->P) PetscCall(MatDestroy(&P));
718:   }
719:   PetscCall(PetscFree3(b,coeffs,matnorm));
720:   PetscCall(PetscFree2(pK,pH));
721:   PetscFunctionReturn(PETSC_SUCCESS);
722: }

724: static PetscErrorCode NEPNLEIGSDividedDifferences_callback(NEP nep)
725: {
726:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
727:   PetscInt       k,j,i,maxnmat;
728:   PetscReal      norm0,norm;
729:   PetscScalar    *s=ctx->s,*beta=ctx->beta,*b,*coeffs;
730:   Mat            *D=ctx->D,*DP,T,P;
731:   PetscBool      shell,has,vec=PETSC_FALSE,precond=(nep->function_pre!=nep->function)?PETSC_TRUE:PETSC_FALSE;
732:   PetscRandom    rand=NULL;
733:   Vec            w[2];

735:   PetscFunctionBegin;
736:   PetscCall(PetscMalloc2(ctx->ddmaxit+1,&b,ctx->ddmaxit+1,&coeffs));
737:   if (nep->V) PetscCall(BVGetRandomContext(nep->V,&rand));
738:   T = nep->function;
739:   P = nep->function_pre;
740:   PetscCall(NEPComputeFunction(nep,s[0],T,P));
741:   PetscCall(MatIsShell(T,&shell));
742:   maxnmat = PetscMax(ctx->ddmaxit,nep->nt);
743:   if (!shell) PetscCall(MatDuplicate(T,MAT_COPY_VALUES,&D[0]));
744:   else PetscCall(NLEIGSMatToMatShellArray(T,&D[0],maxnmat));
745:   if (beta[0]!=1.0) PetscCall(MatScale(D[0],1.0/beta[0]));
746:   PetscCall(MatHasOperation(D[0],MATOP_NORM,&has));
747:   if (has) PetscCall(MatNorm(D[0],NORM_FROBENIUS,&norm0));
748:   else {
749:     PetscCall(MatCreateVecs(D[0],NULL,&w[0]));
750:     PetscCall(VecDuplicate(w[0],&w[1]));
751:     PetscCall(VecDuplicate(w[0],&ctx->vrn));
752:     PetscCall(VecSetRandomNormal(ctx->vrn,rand,w[0],w[1]));
753:     PetscCall(VecNormalize(ctx->vrn,NULL));
754:     vec = PETSC_TRUE;
755:     PetscCall(MatNormEstimate(D[0],ctx->vrn,w[0],&norm0));
756:   }
757:   if (precond) {
758:     PetscCall(PetscMalloc1(ctx->ddmaxit,&DP));
759:     PetscCall(MatDuplicate(P,MAT_COPY_VALUES,&DP[0]));
760:   }
761:   ctx->nmat = ctx->ddmaxit;
762:   for (k=1;k<ctx->ddmaxit;k++) {
763:     PetscCall(NEPNLEIGSEvalNRTFunct(nep,k,s[k],b));
764:     PetscCall(NEPComputeFunction(nep,s[k],T,P));
765:     if (!shell) PetscCall(MatDuplicate(T,MAT_COPY_VALUES,&D[k]));
766:     else PetscCall(NLEIGSMatToMatShellArray(T,&D[k],maxnmat));
767:     for (j=0;j<k;j++) PetscCall(MatAXPY(D[k],-b[j],D[j],nep->mstr));
768:     PetscCall(MatScale(D[k],1.0/b[k]));
769:     PetscCall(MatHasOperation(D[k],MATOP_NORM,&has));
770:     if (has) PetscCall(MatNorm(D[k],NORM_FROBENIUS,&norm));
771:     else {
772:       if (!vec) {
773:         PetscCall(MatCreateVecs(D[k],NULL,&w[0]));
774:         PetscCall(VecDuplicate(w[0],&w[1]));
775:         PetscCall(VecDuplicate(w[0],&ctx->vrn));
776:         PetscCall(VecSetRandomNormal(ctx->vrn,rand,w[0],w[1]));
777:         PetscCall(VecNormalize(ctx->vrn,NULL));
778:         vec = PETSC_TRUE;
779:       }
780:       PetscCall(MatNormEstimate(D[k],ctx->vrn,w[0],&norm));
781:     }
782:     if (precond) {
783:       PetscCall(MatDuplicate(P,MAT_COPY_VALUES,&DP[k]));
784:       for (j=0;j<k;j++) PetscCall(MatAXPY(DP[k],-b[j],DP[j],nep->mstrp));
785:       PetscCall(MatScale(DP[k],1.0/b[k]));
786:     }
787:     if (k>1 && norm/norm0 < ctx->ddtol && k>1) {
788:       ctx->nmat = k+1;
789:       break;
790:     }
791:   }
792:   if (!ctx->ksp) PetscCall(NEPNLEIGSGetKSPs(nep,&ctx->nshiftsw,&ctx->ksp));
793:   for (i=0;i<ctx->nshiftsw;i++) {
794:     PetscCall(NEPNLEIGSEvalNRTFunct(nep,ctx->nmat-1,ctx->shifts[i],coeffs));
795:     PetscCall(MatDuplicate(D[0],MAT_COPY_VALUES,&T));
796:     if (coeffs[0]!=1.0) PetscCall(MatScale(T,coeffs[0]));
797:     for (j=1;j<ctx->nmat;j++) PetscCall(MatAXPY(T,coeffs[j],D[j],nep->mstr));
798:     if (precond) {
799:       PetscCall(MatDuplicate(DP[0],MAT_COPY_VALUES,&P));
800:       if (coeffs[0]!=1.0) PetscCall(MatScale(P,coeffs[0]));
801:       for (j=1;j<ctx->nmat;j++) PetscCall(MatAXPY(P,coeffs[j],DP[j],nep->mstrp));
802:     } else P=T;
803:     PetscCall(NEP_KSPSetOperators(ctx->ksp[i],T,P));
804:     PetscCall(KSPSetUp(ctx->ksp[i]));
805:     PetscCall(MatDestroy(&T));
806:   }
807:   PetscCall(PetscFree2(b,coeffs));
808:   if (vec) {
809:     PetscCall(VecDestroy(&w[0]));
810:     PetscCall(VecDestroy(&w[1]));
811:   }
812:   if (precond) {
813:     PetscCall(MatDestroy(&P));
814:     PetscCall(MatDestroyMatrices(ctx->nmat,&DP));
815:   }
816:   PetscFunctionReturn(PETSC_SUCCESS);
817: }

819: /*
820:    NEPKrylovConvergence - This is the analogue to EPSKrylovConvergence.
821: */
822: PetscErrorCode NEPNLEIGSKrylovConvergence(NEP nep,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal betah,PetscScalar betak,PetscInt *kout,Vec *w)
823: {
824:   PetscInt       k,newk,marker,inside;
825:   PetscScalar    re,im;
826:   PetscReal      resnorm,tt;
827:   PetscBool      istrivial;
828:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;

830:   PetscFunctionBegin;
831:   PetscCall(RGIsTrivial(nep->rg,&istrivial));
832:   marker = -1;
833:   if (nep->trackall) getall = PETSC_TRUE;
834:   for (k=kini;k<kini+nits;k++) {
835:     /* eigenvalue */
836:     re = nep->eigr[k];
837:     im = nep->eigi[k];
838:     if (!istrivial) {
839:       if (!ctx->nshifts) PetscCall(NEPNLEIGSBackTransform((PetscObject)nep,1,&re,&im));
840:       PetscCall(RGCheckInside(nep->rg,1,&re,&im,&inside));
841:       if (marker==-1 && inside<0) marker = k;
842:     }
843:     newk = k;
844:     PetscCall(DSVectors(nep->ds,DS_MAT_X,&newk,&resnorm));
845:     tt = (ctx->nshifts)?SlepcAbsEigenvalue(betak-nep->eigr[k]*betah,nep->eigi[k]*betah):betah;
846:     resnorm *=  PetscAbsReal(tt);
847:     /* error estimate */
848:     PetscCall((*nep->converged)(nep,nep->eigr[k],nep->eigi[k],resnorm,&nep->errest[k],nep->convergedctx));
849:     if (marker==-1 && nep->errest[k] >= nep->tol) marker = k;
850:     if (newk==k+1) {
851:       nep->errest[k+1] = nep->errest[k];
852:       k++;
853:     }
854:     if (marker!=-1 && !getall) break;
855:   }
856:   if (marker!=-1) k = marker;
857:   *kout = k;
858:   PetscFunctionReturn(PETSC_SUCCESS);
859: }

861: PetscErrorCode NEPSetUp_NLEIGS(NEP nep)
862: {
863:   PetscInt       k,in;
864:   PetscScalar    zero=0.0;
865:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
866:   SlepcSC        sc;
867:   PetscBool      istrivial;

869:   PetscFunctionBegin;
870:   PetscCall(NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd));
871:   PetscCheck(nep->ncv<=nep->nev+nep->mpd,PetscObjectComm((PetscObject)nep),PETSC_ERR_USER_INPUT,"The value of ncv must not be larger than nev+mpd");
872:   if (nep->max_it==PETSC_DEFAULT) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
873:   if (!ctx->ddmaxit) ctx->ddmaxit = LBPOINTS;
874:   PetscCall(RGIsTrivial(nep->rg,&istrivial));
875:   PetscCheck(!istrivial,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPNLEIGS requires a nontrivial region defining the target set");
876:   if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
877:   PetscCheck(nep->which==NEP_TARGET_MAGNITUDE || nep->which==NEP_TARGET_REAL || nep->which==NEP_TARGET_IMAGINARY || nep->which==NEP_WHICH_USER,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver supports only target selection of eigenvalues");

879:   /* Initialize the NLEIGS context structure */
880:   k = ctx->ddmaxit;
881:   PetscCall(PetscMalloc4(k,&ctx->s,k,&ctx->xi,k,&ctx->beta,k,&ctx->D));
882:   nep->data = ctx;
883:   if (nep->tol==(PetscReal)PETSC_DEFAULT) nep->tol = SLEPC_DEFAULT_TOL;
884:   if (ctx->ddtol==(PetscReal)PETSC_DEFAULT) ctx->ddtol = nep->tol/10.0;
885:   if (!ctx->keep) ctx->keep = 0.5;

887:   /* Compute Leja-Bagby points and scaling values */
888:   PetscCall(NEPNLEIGSLejaBagbyPoints(nep));
889:   if (nep->problem_type!=NEP_RATIONAL) {
890:     PetscCall(RGCheckInside(nep->rg,1,&nep->target,&zero,&in));
891:     PetscCheck(in>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The target is not inside the target set");
892:   }

894:   /* Compute the divided difference matrices */
895:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) PetscCall(NEPNLEIGSDividedDifferences_split(nep));
896:   else PetscCall(NEPNLEIGSDividedDifferences_callback(nep));
897:   PetscCall(NEPAllocateSolution(nep,ctx->nmat-1));
898:   PetscCall(NEPSetWorkVecs(nep,4));
899:   if (!ctx->fullbasis) {
900:     PetscCheck(!nep->twosided,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Two-sided variant requires the full-basis option, rerun with -nep_nleigs_full_basis");
901:     /* set-up DS and transfer split operator functions */
902:     PetscCall(DSSetType(nep->ds,ctx->nshifts?DSGNHEP:DSNHEP));
903:     PetscCall(DSAllocate(nep->ds,nep->ncv+1));
904:     PetscCall(DSGetSlepcSC(nep->ds,&sc));
905:     if (!ctx->nshifts) sc->map = NEPNLEIGSBackTransform;
906:     PetscCall(DSSetExtraRow(nep->ds,PETSC_TRUE));
907:     sc->mapobj        = (PetscObject)nep;
908:     sc->rg            = nep->rg;
909:     sc->comparison    = nep->sc->comparison;
910:     sc->comparisonctx = nep->sc->comparisonctx;
911:     PetscCall(BVDestroy(&ctx->V));
912:     PetscCall(BVCreateTensor(nep->V,ctx->nmat-1,&ctx->V));
913:     nep->ops->solve          = NEPSolve_NLEIGS;
914:     nep->ops->computevectors = NEPComputeVectors_Schur;
915:   } else {
916:     PetscCall(NEPSetUp_NLEIGS_FullBasis(nep));
917:     nep->ops->solve          = NEPSolve_NLEIGS_FullBasis;
918:     nep->ops->computevectors = NULL;
919:   }
920:   PetscFunctionReturn(PETSC_SUCCESS);
921: }

923: /*
924:   Extend the TOAR basis by applying the matrix operator
925:   over a vector which is decomposed on the TOAR way
926:   Input:
927:     - S,V: define the latest Arnoldi vector (nv vectors in V)
928:   Output:
929:     - t: new vector extending the TOAR basis
930:     - r: temporally coefficients to compute the TOAR coefficients
931:          for the new Arnoldi vector
932:   Workspace: t_ (two vectors)
933: */
934: 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_)
935: {
936:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
937:   PetscInt       deg=ctx->nmat-1,k,j;
938:   Vec            v=t_[0],q=t_[1],w;
939:   PetscScalar    *beta=ctx->beta,*s=ctx->s,*xi=ctx->xi,*coeffs,sigma;

941:   PetscFunctionBegin;
942:   if (!ctx->ksp) PetscCall(NEPNLEIGSGetKSPs(nep,&ctx->nshiftsw,&ctx->ksp));
943:   sigma = ctx->shifts[idxrktg];
944:   PetscCall(BVSetActiveColumns(nep->V,0,nv));
945:   PetscCall(PetscMalloc1(ctx->nmat,&coeffs));
946:   PetscCheck(PetscAbsScalar(s[deg-2]-sigma)>100*PETSC_MACHINE_EPSILON,PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Breakdown in NLEIGS");
947:   /* i-part stored in (i-1) position */
948:   for (j=0;j<nv;j++) {
949:     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);
950:   }
951:   PetscCall(BVSetActiveColumns(W,0,deg));
952:   PetscCall(BVGetColumn(W,deg-1,&w));
953:   PetscCall(BVMultVec(V,1.0/beta[deg],0,w,S+(deg-1)*ls));
954:   PetscCall(BVRestoreColumn(W,deg-1,&w));
955:   PetscCall(BVGetColumn(W,deg-2,&w));
956:   PetscCall(BVMultVec(V,1.0,0.0,w,r+(deg-2)*lr));
957:   PetscCall(BVRestoreColumn(W,deg-2,&w));
958:   for (k=deg-2;k>0;k--) {
959:     PetscCheck(PetscAbsScalar(s[k-1]-sigma)>100*PETSC_MACHINE_EPSILON,PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Breakdown in NLEIGS");
960:     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);
961:     PetscCall(BVGetColumn(W,k-1,&w));
962:     PetscCall(BVMultVec(V,1.0,0.0,w,r+(k-1)*lr));
963:     PetscCall(BVRestoreColumn(W,k-1,&w));
964:   }
965:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
966:     for (j=0;j<ctx->nmat-2;j++) coeffs[j] = ctx->coeffD[nep->nt*j];
967:     coeffs[ctx->nmat-2] = ctx->coeffD[nep->nt*(ctx->nmat-1)];
968:     PetscCall(BVMultVec(W,1.0,0.0,v,coeffs));
969:     PetscCall(MatMult(nep->A[0],v,q));
970:     for (k=1;k<nep->nt;k++) {
971:       for (j=0;j<ctx->nmat-2;j++) coeffs[j] = ctx->coeffD[nep->nt*j+k];
972:       coeffs[ctx->nmat-2] = ctx->coeffD[nep->nt*(ctx->nmat-1)+k];
973:       PetscCall(BVMultVec(W,1.0,0,v,coeffs));
974:       PetscCall(MatMult(nep->A[k],v,t));
975:       PetscCall(VecAXPY(q,1.0,t));
976:     }
977:     PetscCall(KSPSolve(ctx->ksp[idxrktg],q,t));
978:     PetscCall(VecScale(t,-1.0));
979:   } else {
980:     for (k=0;k<deg-1;k++) {
981:       PetscCall(BVGetColumn(W,k,&w));
982:       PetscCall(MatMult(ctx->D[k],w,q));
983:       PetscCall(BVRestoreColumn(W,k,&w));
984:       PetscCall(BVInsertVec(W,k,q));
985:     }
986:     PetscCall(BVGetColumn(W,deg-1,&w));
987:     PetscCall(MatMult(ctx->D[deg],w,q));
988:     PetscCall(BVRestoreColumn(W,k,&w));
989:     PetscCall(BVInsertVec(W,k,q));
990:     for (j=0;j<ctx->nmat-1;j++) coeffs[j] = 1.0;
991:     PetscCall(BVMultVec(W,1.0,0.0,q,coeffs));
992:     PetscCall(KSPSolve(ctx->ksp[idxrktg],q,t));
993:     PetscCall(VecScale(t,-1.0));
994:   }
995:   PetscCall(PetscFree(coeffs));
996:   PetscFunctionReturn(PETSC_SUCCESS);
997: }

999: /*
1000:   Compute TOAR coefficients of the blocks of the new Arnoldi vector computed
1001: */
1002: static PetscErrorCode NEPTOARCoefficients(NEP nep,PetscScalar sigma,PetscInt nv,PetscScalar *S,PetscInt ls,PetscScalar *r,PetscInt lr,PetscScalar *x,PetscScalar *work)
1003: {
1004:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
1005:   PetscInt       k,j,d=ctx->nmat-1;
1006:   PetscScalar    *t=work;

1008:   PetscFunctionBegin;
1009:   PetscCall(NEPNLEIGSEvalNRTFunct(nep,d-1,sigma,t));
1010:   for (k=0;k<d-1;k++) {
1011:     for (j=0;j<=nv;j++) r[k*lr+j] += t[k]*x[j];
1012:   }
1013:   for (j=0;j<=nv;j++) r[(d-1)*lr+j] = t[d-1]*x[j];
1014:   PetscFunctionReturn(PETSC_SUCCESS);
1015: }

1017: /*
1018:   Compute continuation vector coefficients for the Rational-Krylov run.
1019:   dim(work) >= (end-ini)*(end-ini+1) + end+1 + 2*(end-ini+1), dim(t) = end.
1020: */
1021: 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)
1022: {
1023:   PetscScalar    *x,*W,*tau,sone=1.0,szero=0.0;
1024:   PetscInt       i,j,n1,n,nwu=0;
1025:   PetscBLASInt   info,n_,n1_,one=1,dim,lds_;
1026:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;

1028:   PetscFunctionBegin;
1029:   if (!ctx->nshifts || !end) {
1030:     t[0] = 1;
1031:     PetscCall(PetscArraycpy(cont,S+end*lds,lds));
1032:   } else {
1033:     n   = end-ini;
1034:     n1  = n+1;
1035:     x   = work+nwu;
1036:     nwu += end+1;
1037:     tau = work+nwu;
1038:     nwu += n;
1039:     W   = work+nwu;
1040:     nwu += n1*n;
1041:     for (j=ini;j<end;j++) {
1042:       for (i=ini;i<=end;i++) W[(j-ini)*n1+i-ini] = K[j*ld+i] -H[j*ld+i]*sigma;
1043:     }
1044:     PetscCall(PetscBLASIntCast(n,&n_));
1045:     PetscCall(PetscBLASIntCast(n1,&n1_));
1046:     PetscCall(PetscBLASIntCast(end+1,&dim));
1047:     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1048:     PetscCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&n1_,&n_,W,&n1_,tau,work+nwu,&n1_,&info));
1049:     SlepcCheckLapackInfo("geqrf",info);
1050:     for (i=0;i<end;i++) t[i] = 0.0;
1051:     t[end] = 1.0;
1052:     for (j=n-1;j>=0;j--) {
1053:       for (i=0;i<ini+j;i++) x[i] = 0.0;
1054:       x[ini+j] = 1.0;
1055:       for (i=j+1;i<n1;i++) x[i+ini] = W[i+n1*j];
1056:       tau[j] = PetscConj(tau[j]);
1057:       PetscCallBLAS("LAPACKlarf",LAPACKlarf_("L",&dim,&one,x,&one,tau+j,t,&dim,work+nwu));
1058:     }
1059:     PetscCall(PetscBLASIntCast(lds,&lds_));
1060:     PetscCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&n1_,&sone,S,&lds_,t,&one,&szero,cont,&one));
1061:     PetscCall(PetscFPTrapPop());
1062:   }
1063:   PetscFunctionReturn(PETSC_SUCCESS);
1064: }

1066: /*
1067:   Compute a run of Arnoldi iterations
1068: */
1069: PetscErrorCode NEPNLEIGSTOARrun(NEP nep,Mat MK,Mat MH,BV W,PetscInt k,PetscInt *M,PetscReal *betah,PetscScalar *betak,PetscBool *breakdown,Vec *t_)
1070: {
1071:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;
1072:   PetscInt       i,j,m=*M,lwa,deg=ctx->nmat-1,lds,nqt,ld,l,ldh;
1073:   Vec            t;
1074:   PetscReal      norm=0.0;
1075:   PetscScalar    *x,*work,*tt,sigma=1.0,*cont,*S,*K=NULL,*H;
1076:   PetscBool      lindep;
1077:   Mat            MS;

1079:   PetscFunctionBegin;
1080:   *betah = 0.0; *betak = 0.0;
1081:   PetscCall(MatDenseGetArray(MH,&H));
1082:   if (MK) PetscCall(MatDenseGetArray(MK,&K));
1083:   PetscCall(MatDenseGetLDA(MH,&ldh));
1084:   PetscCall(BVTensorGetFactors(ctx->V,NULL,&MS));
1085:   PetscCall(MatDenseGetArray(MS,&S));
1086:   PetscCall(BVGetSizes(nep->V,NULL,NULL,&ld));
1087:   lds = ld*deg;
1088:   PetscCall(BVGetActiveColumns(nep->V,&l,&nqt));
1089:   lwa = PetscMax(ld,deg)+(m+1)*(m+1)+4*(m+1);
1090:   PetscCall(PetscMalloc4(ld,&x,lwa,&work,m+1,&tt,lds,&cont));
1091:   PetscCall(BVSetActiveColumns(ctx->V,0,m));
1092:   for (j=k;j<m;j++) {
1093:     sigma = ctx->shifts[(++(ctx->idxrk))%ctx->nshiftsw];

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

1098:     /* apply operator */
1099:     PetscCall(BVGetColumn(nep->V,nqt,&t));
1100:     PetscCall(NEPTOARExtendBasis(nep,(ctx->idxrk)%ctx->nshiftsw,cont,ld,nqt,W,nep->V,t,S+(j+1)*lds,ld,t_));
1101:     PetscCall(BVRestoreColumn(nep->V,nqt,&t));

1103:     /* orthogonalize */
1104:     PetscCall(BVOrthogonalizeColumn(nep->V,nqt,x,&norm,&lindep));
1105:     if (!lindep) {
1106:       x[nqt] = norm;
1107:       PetscCall(BVScaleColumn(nep->V,nqt,1.0/norm));
1108:       nqt++;
1109:     } else x[nqt] = 0.0;

1111:     PetscCall(NEPTOARCoefficients(nep,sigma,nqt-1,cont,ld,S+(j+1)*lds,ld,x,work));

1113:     /* Level-2 orthogonalization */
1114:     PetscCall(BVOrthogonalizeColumn(ctx->V,j+1,H+j*ldh,&norm,breakdown));
1115:     H[j+1+ldh*j] = norm;
1116:     if (ctx->nshifts && MK) {
1117:       for (i=0;i<=j;i++) K[i+ldh*j] = sigma*H[i+ldh*j] + tt[i];
1118:       K[j+1+ldh*j] = sigma*H[j+1+ldh*j];
1119:     }
1120:     if (*breakdown) {
1121:       *M = j+1;
1122:       break;
1123:     }
1124:     PetscCall(BVScaleColumn(ctx->V,j+1,1.0/norm));
1125:     PetscCall(BVSetActiveColumns(nep->V,l,nqt));
1126:   }
1127:   *betah = norm;
1128:   if (ctx->nshifts) *betak = norm*sigma;
1129:   PetscCall(PetscFree4(x,work,tt,cont));
1130:   PetscCall(MatDenseRestoreArray(MS,&S));
1131:   PetscCall(MatDenseRestoreArray(MH,&H));
1132:   if (MK) PetscCall(MatDenseRestoreArray(MK,&K));
1133:   PetscCall(BVTensorRestoreFactors(ctx->V,NULL,&MS));
1134:   PetscFunctionReturn(PETSC_SUCCESS);
1135: }

1137: PetscErrorCode NEPSolve_NLEIGS(NEP nep)
1138: {
1139:   NEP_NLEIGS        *ctx = (NEP_NLEIGS*)nep->data;
1140:   PetscInt          i,k=0,l,nv=0,ld,lds,nq;
1141:   PetscInt          deg=ctx->nmat-1,nconv=0,dsn,dsk;
1142:   PetscScalar       *pU,betak=0,*eigr,*eigi;
1143:   const PetscScalar *S;
1144:   PetscReal         betah;
1145:   PetscBool         falselock=PETSC_FALSE,breakdown=PETSC_FALSE;
1146:   BV                W;
1147:   Mat               H,K=NULL,MS,MQ,U;

1149:   PetscFunctionBegin;
1150:   if (ctx->lock) {
1151:     /* undocumented option to use a cheaper locking instead of the true locking */
1152:     PetscCall(PetscOptionsGetBool(NULL,NULL,"-nep_nleigs_falselocking",&falselock,NULL));
1153:   }

1155:   PetscCall(BVGetSizes(nep->V,NULL,NULL,&ld));
1156:   lds = deg*ld;
1157:   if (!ctx->nshifts) PetscCall(PetscMalloc2(nep->ncv,&eigr,nep->ncv,&eigi));
1158:   else { eigr = nep->eigr; eigi = nep->eigi; }
1159:   PetscCall(BVDuplicateResize(nep->V,PetscMax(nep->nt-1,ctx->nmat-1),&W));

1161:   /* clean projected matrix (including the extra-arrow) */
1162:   PetscCall(DSSetDimensions(nep->ds,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
1163:   PetscCall(DSGetMat(nep->ds,DS_MAT_A,&H));
1164:   PetscCall(MatZeroEntries(H));
1165:   PetscCall(DSRestoreMat(nep->ds,DS_MAT_A,&H));
1166:   if (ctx->nshifts) {
1167:     PetscCall(DSGetMat(nep->ds,DS_MAT_B,&H));
1168:     PetscCall(MatZeroEntries(H));
1169:     PetscCall(DSRestoreMat(nep->ds,DS_MAT_B,&H));
1170:   }

1172:   /* Get the starting Arnoldi vector */
1173:   PetscCall(BVTensorBuildFirstColumn(ctx->V,nep->nini));

1175:   /* Restart loop */
1176:   l = 0;
1177:   while (nep->reason == NEP_CONVERGED_ITERATING) {
1178:     nep->its++;

1180:     /* Compute an nv-step Krylov relation */
1181:     nv = PetscMin(nep->nconv+nep->mpd,nep->ncv);
1182:     if (ctx->nshifts) PetscCall(DSGetMat(nep->ds,DS_MAT_A,&K));
1183:     PetscCall(DSGetMat(nep->ds,ctx->nshifts?DS_MAT_B:DS_MAT_A,&H));
1184:     PetscCall(NEPNLEIGSTOARrun(nep,K,H,W,nep->nconv+l,&nv,&betah,&betak,&breakdown,nep->work));
1185:     PetscCall(DSRestoreMat(nep->ds,ctx->nshifts?DS_MAT_B:DS_MAT_A,&H));
1186:     if (ctx->nshifts) PetscCall(DSRestoreMat(nep->ds,DS_MAT_A,&K));
1187:     PetscCall(DSSetDimensions(nep->ds,nv,nep->nconv,nep->nconv+l));
1188:     if (l==0) PetscCall(DSSetState(nep->ds,DS_STATE_INTERMEDIATE));
1189:     else PetscCall(DSSetState(nep->ds,DS_STATE_RAW));

1191:     /* Solve projected problem */
1192:     PetscCall(DSSolve(nep->ds,nep->eigr,nep->eigi));
1193:     PetscCall(DSSort(nep->ds,nep->eigr,nep->eigi,NULL,NULL,NULL));
1194:     PetscCall(DSUpdateExtraRow(nep->ds));
1195:     PetscCall(DSSynchronize(nep->ds,nep->eigr,nep->eigi));

1197:     /* Check convergence */
1198:     PetscCall(NEPNLEIGSKrylovConvergence(nep,PETSC_FALSE,nep->nconv,nv-nep->nconv,betah,betak,&k,nep->work));
1199:     PetscCall((*nep->stopping)(nep,nep->its,nep->max_it,k,nep->nev,&nep->reason,nep->stoppingctx));

1201:     /* Update l */
1202:     if (nep->reason != NEP_CONVERGED_ITERATING || breakdown) l = 0;
1203:     else {
1204:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
1205:       PetscCall(DSGetTruncateSize(nep->ds,k,nv,&l));
1206:       if (!breakdown) {
1207:         /* Prepare the Rayleigh quotient for restart */
1208:         PetscCall(DSGetDimensions(nep->ds,&dsn,NULL,&dsk,NULL));
1209:         PetscCall(DSSetDimensions(nep->ds,dsn,k,dsk));
1210:         PetscCall(DSTruncate(nep->ds,k+l,PETSC_FALSE));
1211:       }
1212:     }
1213:     nconv = k;
1214:     if (!ctx->lock && nep->reason == NEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; }
1215:     if (l) PetscCall(PetscInfo(nep,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l));

1217:     /* Update S */
1218:     PetscCall(DSGetMat(nep->ds,ctx->nshifts?DS_MAT_Z:DS_MAT_Q,&MQ));
1219:     PetscCall(BVMultInPlace(ctx->V,MQ,nep->nconv,k+l));
1220:     PetscCall(DSRestoreMat(nep->ds,ctx->nshifts?DS_MAT_Z:DS_MAT_Q,&MQ));

1222:     /* Copy last column of S */
1223:     PetscCall(BVCopyColumn(ctx->V,nv,k+l));

1225:     if (breakdown && nep->reason == NEP_CONVERGED_ITERATING) {
1226:       /* Stop if breakdown */
1227:       PetscCall(PetscInfo(nep,"Breakdown (it=%" PetscInt_FMT " norm=%g)\n",nep->its,(double)betah));
1228:       nep->reason = NEP_DIVERGED_BREAKDOWN;
1229:     }
1230:     if (nep->reason != NEP_CONVERGED_ITERATING) l--;
1231:     /* truncate S */
1232:     PetscCall(BVGetActiveColumns(nep->V,NULL,&nq));
1233:     if (k+l+deg<=nq) {
1234:       PetscCall(BVSetActiveColumns(ctx->V,nep->nconv,k+l+1));
1235:       if (!falselock && ctx->lock) PetscCall(BVTensorCompress(ctx->V,k-nep->nconv));
1236:       else PetscCall(BVTensorCompress(ctx->V,0));
1237:     }
1238:     nep->nconv = k;
1239:     if (!ctx->nshifts) {
1240:       for (i=0;i<nv;i++) { eigr[i] = nep->eigr[i]; eigi[i] = nep->eigi[i]; }
1241:       PetscCall(NEPNLEIGSBackTransform((PetscObject)nep,nv,eigr,eigi));
1242:     }
1243:     PetscCall(NEPMonitor(nep,nep->its,nconv,eigr,eigi,nep->errest,nv));
1244:   }
1245:   nep->nconv = nconv;
1246:   if (nep->nconv>0) {
1247:     PetscCall(BVSetActiveColumns(ctx->V,0,nep->nconv));
1248:     PetscCall(BVGetActiveColumns(nep->V,NULL,&nq));
1249:     PetscCall(BVSetActiveColumns(nep->V,0,nq));
1250:     if (nq>nep->nconv) {
1251:       PetscCall(BVTensorCompress(ctx->V,nep->nconv));
1252:       PetscCall(BVSetActiveColumns(nep->V,0,nep->nconv));
1253:       nq = nep->nconv;
1254:     }
1255:     if (ctx->nshifts) {
1256:       PetscCall(DSGetMat(nep->ds,DS_MAT_B,&MQ));
1257:       PetscCall(BVMultInPlace(ctx->V,MQ,0,nep->nconv));
1258:       PetscCall(DSRestoreMat(nep->ds,DS_MAT_B,&MQ));
1259:     }
1260:     PetscCall(BVTensorGetFactors(ctx->V,NULL,&MS));
1261:     PetscCall(MatDenseGetArrayRead(MS,&S));
1262:     PetscCall(PetscMalloc1(nq*nep->nconv,&pU));
1263:     for (i=0;i<nep->nconv;i++) PetscCall(PetscArraycpy(pU+i*nq,S+i*lds,nq));
1264:     PetscCall(MatDenseRestoreArrayRead(MS,&S));
1265:     PetscCall(BVTensorRestoreFactors(ctx->V,NULL,&MS));
1266:     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nq,nep->nconv,pU,&U));
1267:     PetscCall(BVSetActiveColumns(nep->V,0,nq));
1268:     PetscCall(BVMultInPlace(nep->V,U,0,nep->nconv));
1269:     PetscCall(BVSetActiveColumns(nep->V,0,nep->nconv));
1270:     PetscCall(MatDestroy(&U));
1271:     PetscCall(PetscFree(pU));
1272:     PetscCall(DSTruncate(nep->ds,nep->nconv,PETSC_TRUE));
1273:   }

1275:   /* Map eigenvalues back to the original problem */
1276:   if (!ctx->nshifts) {
1277:     PetscCall(NEPNLEIGSBackTransform((PetscObject)nep,nep->nconv,nep->eigr,nep->eigi));
1278:     PetscCall(PetscFree2(eigr,eigi));
1279:   }
1280:   PetscCall(BVDestroy(&W));
1281:   PetscFunctionReturn(PETSC_SUCCESS);
1282: }

1284: static PetscErrorCode NEPNLEIGSSetSingularitiesFunction_NLEIGS(NEP nep,PetscErrorCode (*fun)(NEP,PetscInt*,PetscScalar*,void*),void *ctx)
1285: {
1286:   NEP_NLEIGS *nepctx=(NEP_NLEIGS*)nep->data;

1288:   PetscFunctionBegin;
1289:   if (fun) nepctx->computesingularities = fun;
1290:   if (ctx) nepctx->singularitiesctx     = ctx;
1291:   nep->state = NEP_STATE_INITIAL;
1292:   PetscFunctionReturn(PETSC_SUCCESS);
1293: }

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

1299:    Logically Collective

1301:    Input Parameters:
1302: +  nep - the NEP context
1303: .  fun - user function (if NULL then NEP retains any previously set value)
1304: -  ctx - [optional] user-defined context for private data for the function
1305:          (may be NULL, in which case NEP retains any previously set value)

1307:    Calling sequence of fun:
1308: $  PetscErrorCode fun(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *ctx)
1309: +   nep   - the NEP context
1310: .   maxnp - on input number of requested points in the discretization (can be set)
1311: .   xi    - computed values of the discretization
1312: -   ctx   - optional context, as set by NEPNLEIGSSetSingularitiesFunction()

1314:    Notes:
1315:    The user-defined function can set a smaller value of maxnp if necessary.
1316:    It is wrong to return a larger value.

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

1322:    Level: intermediate

1324: .seealso: NEPNLEIGSGetSingularitiesFunction(), NEPSetProblemType()
1325: @*/
1326: PetscErrorCode NEPNLEIGSSetSingularitiesFunction(NEP nep,PetscErrorCode (*fun)(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *ctx),void *ctx)
1327: {
1328:   PetscFunctionBegin;
1330:   PetscTryMethod(nep,"NEPNLEIGSSetSingularitiesFunction_C",(NEP,PetscErrorCode(*)(NEP,PetscInt*,PetscScalar*,void*),void*),(nep,fun,ctx));
1331:   PetscFunctionReturn(PETSC_SUCCESS);
1332: }

1334: static PetscErrorCode NEPNLEIGSGetSingularitiesFunction_NLEIGS(NEP nep,PetscErrorCode (**fun)(NEP,PetscInt*,PetscScalar*,void*),void **ctx)
1335: {
1336:   NEP_NLEIGS *nepctx=(NEP_NLEIGS*)nep->data;

1338:   PetscFunctionBegin;
1339:   if (fun) *fun = nepctx->computesingularities;
1340:   if (ctx) *ctx = nepctx->singularitiesctx;
1341:   PetscFunctionReturn(PETSC_SUCCESS);
1342: }

1344: /*@C
1345:    NEPNLEIGSGetSingularitiesFunction - Returns the Function and optionally the user
1346:    provided context for computing a discretization of the singularity set.

1348:    Not Collective

1350:    Input Parameter:
1351: .  nep - the nonlinear eigensolver context

1353:    Output Parameters:
1354: +  fun - location to put the function (or NULL)
1355: -  ctx - location to stash the function context (or NULL)

1357:    Level: advanced

1359: .seealso: NEPNLEIGSSetSingularitiesFunction()
1360: @*/
1361: PetscErrorCode NEPNLEIGSGetSingularitiesFunction(NEP nep,PetscErrorCode (**fun)(NEP,PetscInt*,PetscScalar*,void*),void **ctx)
1362: {
1363:   PetscFunctionBegin;
1365:   PetscUseMethod(nep,"NEPNLEIGSGetSingularitiesFunction_C",(NEP,PetscErrorCode(**)(NEP,PetscInt*,PetscScalar*,void*),void**),(nep,fun,ctx));
1366:   PetscFunctionReturn(PETSC_SUCCESS);
1367: }

1369: static PetscErrorCode NEPNLEIGSSetRestart_NLEIGS(NEP nep,PetscReal keep)
1370: {
1371:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1373:   PetscFunctionBegin;
1374:   if (keep==(PetscReal)PETSC_DEFAULT) ctx->keep = 0.5;
1375:   else {
1376:     PetscCheck(keep>=0.1 && keep<=0.9,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
1377:     ctx->keep = keep;
1378:   }
1379:   PetscFunctionReturn(PETSC_SUCCESS);
1380: }

1382: /*@
1383:    NEPNLEIGSSetRestart - Sets the restart parameter for the NLEIGS
1384:    method, in particular the proportion of basis vectors that must be kept
1385:    after restart.

1387:    Logically Collective

1389:    Input Parameters:
1390: +  nep  - the nonlinear eigensolver context
1391: -  keep - the number of vectors to be kept at restart

1393:    Options Database Key:
1394: .  -nep_nleigs_restart - Sets the restart parameter

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

1399:    Level: advanced

1401: .seealso: NEPNLEIGSGetRestart()
1402: @*/
1403: PetscErrorCode NEPNLEIGSSetRestart(NEP nep,PetscReal keep)
1404: {
1405:   PetscFunctionBegin;
1408:   PetscTryMethod(nep,"NEPNLEIGSSetRestart_C",(NEP,PetscReal),(nep,keep));
1409:   PetscFunctionReturn(PETSC_SUCCESS);
1410: }

1412: static PetscErrorCode NEPNLEIGSGetRestart_NLEIGS(NEP nep,PetscReal *keep)
1413: {
1414:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1416:   PetscFunctionBegin;
1417:   *keep = ctx->keep;
1418:   PetscFunctionReturn(PETSC_SUCCESS);
1419: }

1421: /*@
1422:    NEPNLEIGSGetRestart - Gets the restart parameter used in the NLEIGS method.

1424:    Not Collective

1426:    Input Parameter:
1427: .  nep - the nonlinear eigensolver context

1429:    Output Parameter:
1430: .  keep - the restart parameter

1432:    Level: advanced

1434: .seealso: NEPNLEIGSSetRestart()
1435: @*/
1436: PetscErrorCode NEPNLEIGSGetRestart(NEP nep,PetscReal *keep)
1437: {
1438:   PetscFunctionBegin;
1441:   PetscUseMethod(nep,"NEPNLEIGSGetRestart_C",(NEP,PetscReal*),(nep,keep));
1442:   PetscFunctionReturn(PETSC_SUCCESS);
1443: }

1445: static PetscErrorCode NEPNLEIGSSetLocking_NLEIGS(NEP nep,PetscBool lock)
1446: {
1447:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1449:   PetscFunctionBegin;
1450:   ctx->lock = lock;
1451:   PetscFunctionReturn(PETSC_SUCCESS);
1452: }

1454: /*@
1455:    NEPNLEIGSSetLocking - Choose between locking and non-locking variants of
1456:    the NLEIGS method.

1458:    Logically Collective

1460:    Input Parameters:
1461: +  nep  - the nonlinear eigensolver context
1462: -  lock - true if the locking variant must be selected

1464:    Options Database Key:
1465: .  -nep_nleigs_locking - Sets the locking flag

1467:    Notes:
1468:    The default is to lock converged eigenpairs when the method restarts.
1469:    This behaviour can be changed so that all directions are kept in the
1470:    working subspace even if already converged to working accuracy (the
1471:    non-locking variant).

1473:    Level: advanced

1475: .seealso: NEPNLEIGSGetLocking()
1476: @*/
1477: PetscErrorCode NEPNLEIGSSetLocking(NEP nep,PetscBool lock)
1478: {
1479:   PetscFunctionBegin;
1482:   PetscTryMethod(nep,"NEPNLEIGSSetLocking_C",(NEP,PetscBool),(nep,lock));
1483:   PetscFunctionReturn(PETSC_SUCCESS);
1484: }

1486: static PetscErrorCode NEPNLEIGSGetLocking_NLEIGS(NEP nep,PetscBool *lock)
1487: {
1488:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1490:   PetscFunctionBegin;
1491:   *lock = ctx->lock;
1492:   PetscFunctionReturn(PETSC_SUCCESS);
1493: }

1495: /*@
1496:    NEPNLEIGSGetLocking - Gets the locking flag used in the NLEIGS method.

1498:    Not Collective

1500:    Input Parameter:
1501: .  nep - the nonlinear eigensolver context

1503:    Output Parameter:
1504: .  lock - the locking flag

1506:    Level: advanced

1508: .seealso: NEPNLEIGSSetLocking()
1509: @*/
1510: PetscErrorCode NEPNLEIGSGetLocking(NEP nep,PetscBool *lock)
1511: {
1512:   PetscFunctionBegin;
1515:   PetscUseMethod(nep,"NEPNLEIGSGetLocking_C",(NEP,PetscBool*),(nep,lock));
1516:   PetscFunctionReturn(PETSC_SUCCESS);
1517: }

1519: static PetscErrorCode NEPNLEIGSSetInterpolation_NLEIGS(NEP nep,PetscReal tol,PetscInt degree)
1520: {
1521:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;

1523:   PetscFunctionBegin;
1524:   if (tol == (PetscReal)PETSC_DEFAULT) {
1525:     ctx->ddtol = PETSC_DEFAULT;
1526:     nep->state = NEP_STATE_INITIAL;
1527:   } else {
1528:     PetscCheck(tol>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1529:     ctx->ddtol = tol;
1530:   }
1531:   if (degree == PETSC_DEFAULT || degree == PETSC_DECIDE) {
1532:     ctx->ddmaxit = 0;
1533:     if (nep->state) PetscCall(NEPReset(nep));
1534:     nep->state = NEP_STATE_INITIAL;
1535:   } else {
1536:     PetscCheck(degree>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of degree. Must be > 0");
1537:     if (ctx->ddmaxit != degree) {
1538:       ctx->ddmaxit = degree;
1539:       if (nep->state) PetscCall(NEPReset(nep));
1540:       nep->state = NEP_STATE_INITIAL;
1541:     }
1542:   }
1543:   PetscFunctionReturn(PETSC_SUCCESS);
1544: }

1546: /*@
1547:    NEPNLEIGSSetInterpolation - Sets the tolerance and maximum degree
1548:    when building the interpolation via divided differences.

1550:    Collective

1552:    Input Parameters:
1553: +  nep    - the nonlinear eigensolver context
1554: .  tol    - tolerance to stop computing divided differences
1555: -  degree - maximum degree of interpolation

1557:    Options Database Key:
1558: +  -nep_nleigs_interpolation_tol <tol> - Sets the tolerance to stop computing divided differences
1559: -  -nep_nleigs_interpolation_degree <degree> - Sets the maximum degree of interpolation

1561:    Notes:
1562:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

1564:    Level: advanced

1566: .seealso: NEPNLEIGSGetInterpolation()
1567: @*/
1568: PetscErrorCode NEPNLEIGSSetInterpolation(NEP nep,PetscReal tol,PetscInt degree)
1569: {
1570:   PetscFunctionBegin;
1574:   PetscTryMethod(nep,"NEPNLEIGSSetInterpolation_C",(NEP,PetscReal,PetscInt),(nep,tol,degree));
1575:   PetscFunctionReturn(PETSC_SUCCESS);
1576: }

1578: static PetscErrorCode NEPNLEIGSGetInterpolation_NLEIGS(NEP nep,PetscReal *tol,PetscInt *degree)
1579: {
1580:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1582:   PetscFunctionBegin;
1583:   if (tol)    *tol    = ctx->ddtol;
1584:   if (degree) *degree = ctx->ddmaxit;
1585:   PetscFunctionReturn(PETSC_SUCCESS);
1586: }

1588: /*@
1589:    NEPNLEIGSGetInterpolation - Gets the tolerance and maximum degree
1590:    when building the interpolation via divided differences.

1592:    Not Collective

1594:    Input Parameter:
1595: .  nep - the nonlinear eigensolver context

1597:    Output Parameters:
1598: +  tol    - tolerance to stop computing divided differences
1599: -  degree - maximum degree of interpolation

1601:    Level: advanced

1603: .seealso: NEPNLEIGSSetInterpolation()
1604: @*/
1605: PetscErrorCode NEPNLEIGSGetInterpolation(NEP nep,PetscReal *tol,PetscInt *degree)
1606: {
1607:   PetscFunctionBegin;
1609:   PetscTryMethod(nep,"NEPNLEIGSGetInterpolation_C",(NEP,PetscReal*,PetscInt*),(nep,tol,degree));
1610:   PetscFunctionReturn(PETSC_SUCCESS);
1611: }

1613: static PetscErrorCode NEPNLEIGSSetRKShifts_NLEIGS(NEP nep,PetscInt ns,PetscScalar *shifts)
1614: {
1615:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
1616:   PetscInt       i;

1618:   PetscFunctionBegin;
1619:   PetscCheck(ns>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"Number of shifts must be non-negative");
1620:   if (ctx->nshifts) PetscCall(PetscFree(ctx->shifts));
1621:   for (i=0;i<ctx->nshiftsw;i++) PetscCall(KSPDestroy(&ctx->ksp[i]));
1622:   PetscCall(PetscFree(ctx->ksp));
1623:   ctx->ksp = NULL;
1624:   if (ns) {
1625:     PetscCall(PetscMalloc1(ns,&ctx->shifts));
1626:     for (i=0;i<ns;i++) ctx->shifts[i] = shifts[i];
1627:   }
1628:   ctx->nshifts = ns;
1629:   nep->state   = NEP_STATE_INITIAL;
1630:   PetscFunctionReturn(PETSC_SUCCESS);
1631: }

1633: /*@
1634:    NEPNLEIGSSetRKShifts - Sets a list of shifts to be used in the Rational
1635:    Krylov method.

1637:    Collective

1639:    Input Parameters:
1640: +  nep    - the nonlinear eigensolver context
1641: .  ns     - number of shifts
1642: -  shifts - array of scalar values specifying the shifts

1644:    Options Database Key:
1645: .  -nep_nleigs_rk_shifts - Sets the list of shifts

1647:    Notes:
1648:    If only one shift is provided, the built subspace built is equivalent to
1649:    shift-and-invert Krylov-Schur (provided that the absolute convergence
1650:    criterion is used).

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

1657:    Use ns=0 to remove previously set shifts.

1659:    Level: advanced

1661: .seealso: NEPNLEIGSGetRKShifts()
1662: @*/
1663: PetscErrorCode NEPNLEIGSSetRKShifts(NEP nep,PetscInt ns,PetscScalar shifts[])
1664: {
1665:   PetscFunctionBegin;
1669:   PetscTryMethod(nep,"NEPNLEIGSSetRKShifts_C",(NEP,PetscInt,PetscScalar*),(nep,ns,shifts));
1670:   PetscFunctionReturn(PETSC_SUCCESS);
1671: }

1673: static PetscErrorCode NEPNLEIGSGetRKShifts_NLEIGS(NEP nep,PetscInt *ns,PetscScalar **shifts)
1674: {
1675:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
1676:   PetscInt       i;

1678:   PetscFunctionBegin;
1679:   *ns = ctx->nshifts;
1680:   if (ctx->nshifts) {
1681:     PetscCall(PetscMalloc1(ctx->nshifts,shifts));
1682:     for (i=0;i<ctx->nshifts;i++) (*shifts)[i] = ctx->shifts[i];
1683:   }
1684:   PetscFunctionReturn(PETSC_SUCCESS);
1685: }

1687: /*@C
1688:    NEPNLEIGSGetRKShifts - Gets the list of shifts used in the Rational
1689:    Krylov method.

1691:    Not Collective

1693:    Input Parameter:
1694: .  nep - the nonlinear eigensolver context

1696:    Output Parameters:
1697: +  ns     - number of shifts
1698: -  shifts - array of shifts

1700:    Note:
1701:    The user is responsible for deallocating the returned array.

1703:    Level: advanced

1705: .seealso: NEPNLEIGSSetRKShifts()
1706: @*/
1707: PetscErrorCode NEPNLEIGSGetRKShifts(NEP nep,PetscInt *ns,PetscScalar *shifts[])
1708: {
1709:   PetscFunctionBegin;
1713:   PetscTryMethod(nep,"NEPNLEIGSGetRKShifts_C",(NEP,PetscInt*,PetscScalar**),(nep,ns,shifts));
1714:   PetscFunctionReturn(PETSC_SUCCESS);
1715: }

1717: static PetscErrorCode NEPNLEIGSGetKSPs_NLEIGS(NEP nep,PetscInt *nsolve,KSP **ksp)
1718: {
1719:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;
1720:   PetscInt       i;
1721:   PC             pc;

1723:   PetscFunctionBegin;
1724:   if (!ctx->ksp) {
1725:     PetscCall(NEPNLEIGSSetShifts(nep,&ctx->nshiftsw));
1726:     PetscCall(PetscMalloc1(ctx->nshiftsw,&ctx->ksp));
1727:     for (i=0;i<ctx->nshiftsw;i++) {
1728:       PetscCall(KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp[i]));
1729:       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->ksp[i],(PetscObject)nep,1));
1730:       PetscCall(KSPSetOptionsPrefix(ctx->ksp[i],((PetscObject)nep)->prefix));
1731:       PetscCall(KSPAppendOptionsPrefix(ctx->ksp[i],"nep_nleigs_"));
1732:       PetscCall(PetscObjectSetOptions((PetscObject)ctx->ksp[i],((PetscObject)nep)->options));
1733:       PetscCall(KSPSetErrorIfNotConverged(ctx->ksp[i],PETSC_TRUE));
1734:       PetscCall(KSPSetTolerances(ctx->ksp[i],1e-3*SlepcDefaultTol(nep->tol),PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
1735:       PetscCall(KSPGetPC(ctx->ksp[i],&pc));
1736:       if ((nep->fui==NEP_USER_INTERFACE_SPLIT && nep->P) || (nep->fui==NEP_USER_INTERFACE_CALLBACK && nep->function_pre!=nep->function)) {
1737:         PetscCall(KSPSetType(ctx->ksp[i],KSPBCGS));
1738:         PetscCall(PCSetType(pc,PCBJACOBI));
1739:       } else {
1740:         PetscCall(KSPSetType(ctx->ksp[i],KSPPREONLY));
1741:         PetscCall(PCSetType(pc,PCLU));
1742:       }
1743:     }
1744:   }
1745:   if (nsolve) *nsolve = ctx->nshiftsw;
1746:   if (ksp)    *ksp    = ctx->ksp;
1747:   PetscFunctionReturn(PETSC_SUCCESS);
1748: }

1750: /*@C
1751:    NEPNLEIGSGetKSPs - Retrieve the array of linear solver objects associated with
1752:    the nonlinear eigenvalue solver.

1754:    Collective

1756:    Input Parameter:
1757: .  nep - nonlinear eigenvalue solver

1759:    Output Parameters:
1760: +  nsolve - number of returned KSP objects
1761: -  ksp - array of linear solver object

1763:    Notes:
1764:    The number of KSP objects is equal to the number of shifts provided by the user,
1765:    or 1 if the user did not provide shifts.

1767:    Level: advanced

1769: .seealso: NEPNLEIGSSetRKShifts()
1770: @*/
1771: PetscErrorCode NEPNLEIGSGetKSPs(NEP nep,PetscInt *nsolve,KSP **ksp)
1772: {
1773:   PetscFunctionBegin;
1775:   PetscUseMethod(nep,"NEPNLEIGSGetKSPs_C",(NEP,PetscInt*,KSP**),(nep,nsolve,ksp));
1776:   PetscFunctionReturn(PETSC_SUCCESS);
1777: }

1779: static PetscErrorCode NEPNLEIGSSetFullBasis_NLEIGS(NEP nep,PetscBool fullbasis)
1780: {
1781:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1783:   PetscFunctionBegin;
1784:   if (fullbasis!=ctx->fullbasis) {
1785:     ctx->fullbasis = fullbasis;
1786:     nep->state     = NEP_STATE_INITIAL;
1787:     nep->useds     = PetscNot(fullbasis);
1788:   }
1789:   PetscFunctionReturn(PETSC_SUCCESS);
1790: }

1792: /*@
1793:    NEPNLEIGSSetFullBasis - Choose between TOAR-basis (default) and full-basis
1794:    variants of the NLEIGS method.

1796:    Logically Collective

1798:    Input Parameters:
1799: +  nep  - the nonlinear eigensolver context
1800: -  fullbasis - true if the full-basis variant must be selected

1802:    Options Database Key:
1803: .  -nep_nleigs_full_basis - Sets the full-basis flag

1805:    Notes:
1806:    The default is to use a compact representation of the Krylov basis, that is,
1807:    V = (I otimes U) S, with a tensor BV. This behaviour can be changed so that
1808:    the full basis V is explicitly stored and operated with. This variant is more
1809:    expensive in terms of memory and computation, but is necessary in some cases,
1810:    particularly for two-sided computations, see NEPSetTwoSided().

1812:    In the full-basis variant, the NLEIGS solver uses an EPS object to explicitly
1813:    solve the linearized eigenproblem, see NEPNLEIGSGetEPS().

1815:    Level: advanced

1817: .seealso: NEPNLEIGSGetFullBasis(), NEPNLEIGSGetEPS(), NEPSetTwoSided(), BVCreateTensor()
1818: @*/
1819: PetscErrorCode NEPNLEIGSSetFullBasis(NEP nep,PetscBool fullbasis)
1820: {
1821:   PetscFunctionBegin;
1824:   PetscTryMethod(nep,"NEPNLEIGSSetFullBasis_C",(NEP,PetscBool),(nep,fullbasis));
1825:   PetscFunctionReturn(PETSC_SUCCESS);
1826: }

1828: static PetscErrorCode NEPNLEIGSGetFullBasis_NLEIGS(NEP nep,PetscBool *fullbasis)
1829: {
1830:   NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;

1832:   PetscFunctionBegin;
1833:   *fullbasis = ctx->fullbasis;
1834:   PetscFunctionReturn(PETSC_SUCCESS);
1835: }

1837: /*@
1838:    NEPNLEIGSGetFullBasis - Gets the flag that indicates if NLEIGS is using the
1839:    full-basis variant.

1841:    Not Collective

1843:    Input Parameter:
1844: .  nep - the nonlinear eigensolver context

1846:    Output Parameter:
1847: .  fullbasis - the flag

1849:    Level: advanced

1851: .seealso: NEPNLEIGSSetFullBasis()
1852: @*/
1853: PetscErrorCode NEPNLEIGSGetFullBasis(NEP nep,PetscBool *fullbasis)
1854: {
1855:   PetscFunctionBegin;
1858:   PetscUseMethod(nep,"NEPNLEIGSGetFullBasis_C",(NEP,PetscBool*),(nep,fullbasis));
1859:   PetscFunctionReturn(PETSC_SUCCESS);
1860: }

1862: #define SHIFTMAX 30

1864: PetscErrorCode NEPSetFromOptions_NLEIGS(NEP nep,PetscOptionItems *PetscOptionsObject)
1865: {
1866:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;
1867:   PetscInt       i=0,k;
1868:   PetscBool      flg1,flg2,b;
1869:   PetscReal      r;
1870:   PetscScalar    array[SHIFTMAX];

1872:   PetscFunctionBegin;
1873:   PetscOptionsHeadBegin(PetscOptionsObject,"NEP NLEIGS Options");

1875:     PetscCall(PetscOptionsReal("-nep_nleigs_restart","Proportion of vectors kept after restart","NEPNLEIGSSetRestart",0.5,&r,&flg1));
1876:     if (flg1) PetscCall(NEPNLEIGSSetRestart(nep,r));

1878:     PetscCall(PetscOptionsBool("-nep_nleigs_locking","Choose between locking and non-locking variants","NEPNLEIGSSetLocking",PETSC_FALSE,&b,&flg1));
1879:     if (flg1) PetscCall(NEPNLEIGSSetLocking(nep,b));

1881:     PetscCall(PetscOptionsBool("-nep_nleigs_full_basis","Choose between TOAR and full-basis variants","NEPNLEIGSSetFullBasis",PETSC_FALSE,&b,&flg1));
1882:     if (flg1) PetscCall(NEPNLEIGSSetFullBasis(nep,b));

1884:     PetscCall(NEPNLEIGSGetInterpolation(nep,&r,&i));
1885:     if (!i) i = PETSC_DEFAULT;
1886:     PetscCall(PetscOptionsInt("-nep_nleigs_interpolation_degree","Maximum number of terms for interpolation via divided differences","NEPNLEIGSSetInterpolation",i,&i,&flg1));
1887:     PetscCall(PetscOptionsReal("-nep_nleigs_interpolation_tol","Tolerance for interpolation via divided differences","NEPNLEIGSSetInterpolation",r,&r,&flg2));
1888:     if (flg1 || flg2) PetscCall(NEPNLEIGSSetInterpolation(nep,r,i));

1890:     k = SHIFTMAX;
1891:     for (i=0;i<k;i++) array[i] = 0;
1892:     PetscCall(PetscOptionsScalarArray("-nep_nleigs_rk_shifts","Shifts for Rational Krylov","NEPNLEIGSSetRKShifts",array,&k,&flg1));
1893:     if (flg1) PetscCall(NEPNLEIGSSetRKShifts(nep,k,array));

1895:   PetscOptionsHeadEnd();

1897:   if (!ctx->ksp) PetscCall(NEPNLEIGSGetKSPs(nep,&ctx->nshiftsw,&ctx->ksp));
1898:   for (i=0;i<ctx->nshiftsw;i++) PetscCall(KSPSetFromOptions(ctx->ksp[i]));

1900:   if (ctx->fullbasis) {
1901:     if (!ctx->eps) PetscCall(NEPNLEIGSGetEPS(nep,&ctx->eps));
1902:     PetscCall(EPSSetFromOptions(ctx->eps));
1903:   }
1904:   PetscFunctionReturn(PETSC_SUCCESS);
1905: }

1907: PetscErrorCode NEPView_NLEIGS(NEP nep,PetscViewer viewer)
1908: {
1909:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;
1910:   PetscBool      isascii;
1911:   PetscInt       i;
1912:   char           str[50];

1914:   PetscFunctionBegin;
1915:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1916:   if (isascii) {
1917:     PetscCall(PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep)));
1918:     if (ctx->fullbasis) PetscCall(PetscViewerASCIIPrintf(viewer,"  using the full-basis variant\n"));
1919:     else PetscCall(PetscViewerASCIIPrintf(viewer,"  using the %slocking variant\n",ctx->lock?"":"non-"));
1920:     PetscCall(PetscViewerASCIIPrintf(viewer,"  divided difference terms: used=%" PetscInt_FMT ", max=%" PetscInt_FMT "\n",ctx->nmat,ctx->ddmaxit));
1921:     PetscCall(PetscViewerASCIIPrintf(viewer,"  tolerance for divided difference convergence: %g\n",(double)ctx->ddtol));
1922:     if (ctx->nshifts) {
1923:       PetscCall(PetscViewerASCIIPrintf(viewer,"  RK shifts: "));
1924:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
1925:       for (i=0;i<ctx->nshifts;i++) {
1926:         PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->shifts[i],PETSC_FALSE));
1927:         PetscCall(PetscViewerASCIIPrintf(viewer,"%s%s",str,(i<ctx->nshifts-1)?",":""));
1928:       }
1929:       PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
1930:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
1931:     }
1932:     if (!ctx->ksp) PetscCall(NEPNLEIGSGetKSPs(nep,&ctx->nshiftsw,&ctx->ksp));
1933:     PetscCall(PetscViewerASCIIPushTab(viewer));
1934:     PetscCall(KSPView(ctx->ksp[0],viewer));
1935:     PetscCall(PetscViewerASCIIPopTab(viewer));
1936:     if (ctx->fullbasis) {
1937:       if (!ctx->eps) PetscCall(NEPNLEIGSGetEPS(nep,&ctx->eps));
1938:       PetscCall(PetscViewerASCIIPushTab(viewer));
1939:       PetscCall(EPSView(ctx->eps,viewer));
1940:       PetscCall(PetscViewerASCIIPopTab(viewer));
1941:     }
1942:   }
1943:   PetscFunctionReturn(PETSC_SUCCESS);
1944: }

1946: PetscErrorCode NEPReset_NLEIGS(NEP nep)
1947: {
1948:   PetscInt       k;
1949:   NEP_NLEIGS     *ctx=(NEP_NLEIGS*)nep->data;

1951:   PetscFunctionBegin;
1952:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) PetscCall(PetscFree(ctx->coeffD));
1953:   else {
1954:     for (k=0;k<ctx->nmat;k++) PetscCall(MatDestroy(&ctx->D[k]));
1955:   }
1956:   PetscCall(PetscFree4(ctx->s,ctx->xi,ctx->beta,ctx->D));
1957:   for (k=0;k<ctx->nshiftsw;k++) PetscCall(KSPReset(ctx->ksp[k]));
1958:   PetscCall(VecDestroy(&ctx->vrn));
1959:   if (ctx->fullbasis) {
1960:     PetscCall(MatDestroy(&ctx->A));
1961:     PetscCall(EPSReset(ctx->eps));
1962:     for (k=0;k<4;k++) PetscCall(VecDestroy(&ctx->w[k]));
1963:   }
1964:   PetscFunctionReturn(PETSC_SUCCESS);
1965: }

1967: PetscErrorCode NEPDestroy_NLEIGS(NEP nep)
1968: {
1969:   PetscInt       k;
1970:   NEP_NLEIGS     *ctx = (NEP_NLEIGS*)nep->data;

1972:   PetscFunctionBegin;
1973:   PetscCall(BVDestroy(&ctx->V));
1974:   for (k=0;k<ctx->nshiftsw;k++) PetscCall(KSPDestroy(&ctx->ksp[k]));
1975:   PetscCall(PetscFree(ctx->ksp));
1976:   if (ctx->nshifts) PetscCall(PetscFree(ctx->shifts));
1977:   if (ctx->fullbasis) PetscCall(EPSDestroy(&ctx->eps));
1978:   PetscCall(PetscFree(nep->data));
1979:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetSingularitiesFunction_C",NULL));
1980:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetSingularitiesFunction_C",NULL));
1981:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRestart_C",NULL));
1982:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRestart_C",NULL));
1983:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetLocking_C",NULL));
1984:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetLocking_C",NULL));
1985:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetInterpolation_C",NULL));
1986:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetInterpolation_C",NULL));
1987:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRKShifts_C",NULL));
1988:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRKShifts_C",NULL));
1989:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetKSPs_C",NULL));
1990:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetFullBasis_C",NULL));
1991:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetFullBasis_C",NULL));
1992:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetEPS_C",NULL));
1993:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetEPS_C",NULL));
1994:   PetscFunctionReturn(PETSC_SUCCESS);
1995: }

1997: SLEPC_EXTERN PetscErrorCode NEPCreate_NLEIGS(NEP nep)
1998: {
1999:   NEP_NLEIGS     *ctx;

2001:   PetscFunctionBegin;
2002:   PetscCall(PetscNew(&ctx));
2003:   nep->data  = (void*)ctx;
2004:   ctx->lock  = PETSC_TRUE;
2005:   ctx->ddtol = PETSC_DEFAULT;

2007:   nep->useds = PETSC_TRUE;

2009:   nep->ops->setup          = NEPSetUp_NLEIGS;
2010:   nep->ops->setfromoptions = NEPSetFromOptions_NLEIGS;
2011:   nep->ops->view           = NEPView_NLEIGS;
2012:   nep->ops->destroy        = NEPDestroy_NLEIGS;
2013:   nep->ops->reset          = NEPReset_NLEIGS;

2015:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetSingularitiesFunction_C",NEPNLEIGSSetSingularitiesFunction_NLEIGS));
2016:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetSingularitiesFunction_C",NEPNLEIGSGetSingularitiesFunction_NLEIGS));
2017:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRestart_C",NEPNLEIGSSetRestart_NLEIGS));
2018:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRestart_C",NEPNLEIGSGetRestart_NLEIGS));
2019:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetLocking_C",NEPNLEIGSSetLocking_NLEIGS));
2020:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetLocking_C",NEPNLEIGSGetLocking_NLEIGS));
2021:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetInterpolation_C",NEPNLEIGSSetInterpolation_NLEIGS));
2022:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetInterpolation_C",NEPNLEIGSGetInterpolation_NLEIGS));
2023:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRKShifts_C",NEPNLEIGSSetRKShifts_NLEIGS));
2024:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRKShifts_C",NEPNLEIGSGetRKShifts_NLEIGS));
2025:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetKSPs_C",NEPNLEIGSGetKSPs_NLEIGS));
2026:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetFullBasis_C",NEPNLEIGSSetFullBasis_NLEIGS));
2027:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetFullBasis_C",NEPNLEIGSGetFullBasis_NLEIGS));
2028:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetEPS_C",NEPNLEIGSSetEPS_NLEIGS));
2029:   PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetEPS_C",NEPNLEIGSGetEPS_NLEIGS));
2030:   PetscFunctionReturn(PETSC_SUCCESS);
2031: }