Actual source code: nleigs.c
slepc-3.7.2 2016-07-19
1: /*
3: SLEPc nonlinear eigensolver: "nleigs"
5: Method: NLEIGS
7: Algorithm:
9: Fully rational Krylov method for nonlinear eigenvalue problems.
11: References:
13: [1] S. Guttel et al., "NLEIGS: A class of robust fully rational Krylov
14: method for nonlinear eigenvalue problems", SIAM J. Sci. Comput.
15: 36(6):A2842-A2864, 2014.
17: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18: SLEPc - Scalable Library for Eigenvalue Problem Computations
19: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
21: This file is part of SLEPc.
23: SLEPc is free software: you can redistribute it and/or modify it under the
24: terms of version 3 of the GNU Lesser General Public License as published by
25: the Free Software Foundation.
27: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
28: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
29: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
30: more details.
32: You should have received a copy of the GNU Lesser General Public License
33: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
34: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: */
37: #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/
38: #include <slepcblaslapack.h>
40: #define MAX_LBPOINTS 100
41: #define NDPOINTS 1e4
42: #define MAX_NSHIFTS 100
44: typedef struct {
45: PetscInt nmat; /* number of interpolation points */
46: PetscScalar *s,*xi; /* Leja-Bagby points */
47: PetscScalar *beta; /* scaling factors */
48: Mat *D; /* divided difference matrices */
49: PetscScalar *coeffD; /* coefficients for divided differences in split form */
50: PetscInt nshifts; /* provided number of shifts */
51: PetscScalar *shifts; /* user-provided shifts for the Rational Krylov variant */
52: PetscInt nshiftsw; /* actual number of shifts (1 if Krylov-Schur) */
53: PetscReal ddtol; /* tolerance for divided difference convergence */
54: PetscInt ddmaxit; /* maximum number of divided difference terms */
55: BV W; /* auxiliary BV object */
56: PetscReal keep; /* restart parameter */
57: PetscBool lock; /* locking/non-locking variant */
58: PetscBool trueres; /* whether the true residual norm must be computed */
59: PetscInt idxrk; /* index of next shift to use */
60: KSP *ksp; /* ksp array for storing shift factorizations */
61: Vec vrn; /* random vector with normally distributed value */
62: void *singularitiesctx;
63: PetscErrorCode (*computesingularities)(NEP,PetscInt*,PetscScalar*,void*);
64: } NEP_NLEIGS;
66: typedef struct {
67: PetscInt nmat;
68: PetscScalar coeff[MAX_NSHIFTS];
69: Mat A[MAX_NSHIFTS];
70: Vec t;
71: } ShellMatCtx;
75: PETSC_STATIC_INLINE PetscErrorCode NEPNLEIGSSetShifts(NEP nep)
76: {
77: NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;
80: if (!ctx->nshifts) {
81: ctx->shifts = &nep->target;
82: ctx->nshiftsw = 1;
83: } else ctx->nshiftsw = ctx->nshifts;
84: return(0);
85: }
89: static PetscErrorCode NEPNLEIGSBackTransform(PetscObject ob,PetscInt n,PetscScalar *valr,PetscScalar *vali)
90: {
91: NEP nep;
92: PetscInt j;
93: #if !defined(PETSC_USE_COMPLEX)
94: PetscScalar t;
95: #endif
98: nep = (NEP)ob;
99: #if !defined(PETSC_USE_COMPLEX)
100: for (j=0;j<n;j++) {
101: if (vali[j] == 0) valr[j] = 1.0 / valr[j] + nep->target;
102: else {
103: t = valr[j] * valr[j] + vali[j] * vali[j];
104: valr[j] = valr[j] / t + nep->target;
105: vali[j] = - vali[j] / t;
106: }
107: }
108: #else
109: for (j=0;j<n;j++) {
110: valr[j] = 1.0 / valr[j] + nep->target;
111: }
112: #endif
113: return(0);
114: }
118: static PetscErrorCode NEPNLEIGSLejaBagbyPoints(NEP nep)
119: {
121: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
122: PetscInt i,k,ndpt=NDPOINTS,ndptx=NDPOINTS;
123: PetscScalar *ds,*dsi,*dxi,*nrs,*nrxi,*s=ctx->s,*xi=ctx->xi,*beta=ctx->beta;
124: PetscReal maxnrs,minnrxi;
127: PetscMalloc5(ndpt+1,&ds,ndpt+1,&dsi,ndpt,&dxi,ndpt+1,&nrs,ndpt,&nrxi);
129: /* Discretize the target region boundary */
130: RGComputeContour(nep->rg,ndpt,ds,dsi);
131: #if !defined(PETSC_USE_COMPLEX)
132: for (i=0;i<ndpt;i++) if (dsi[i]!=0.0) break;
133: if (i<ndpt) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NLEIGS with real arithmetic requires the target set to be included in the real axis");
134: #endif
135: /* Discretize the singularity region */
136: if (ctx->computesingularities) {
137: (ctx->computesingularities)(nep,&ndptx,dxi,ctx->singularitiesctx);
138: } else ndptx = 0;
140: /* Look for Leja-Bagby points in the discretization sets */
141: s[0] = ds[0];
142: xi[0] = (ndptx>0)?dxi[0]:PETSC_INFINITY;
143: beta[0] = 1.0; /* scaling factors are also computed here */
144: maxnrs = 0.0;
145: minnrxi = PETSC_MAX_REAL;
146: for (i=0;i<ndpt;i++) {
147: nrs[i] = (ds[i]-s[0])/(1.0-ds[i]/xi[0]);
148: if (PetscAbsScalar(nrs[i])>=maxnrs) {maxnrs = PetscAbsScalar(nrs[i]); s[1] = ds[i];}
149: }
150: for (i=1;i<ndptx;i++) {
151: nrxi[i] = (dxi[i]-s[0])/(1.0-dxi[i]/xi[0]);
152: if (PetscAbsScalar(nrxi[i])<=minnrxi) {minnrxi = PetscAbsScalar(nrxi[i]); xi[1] = dxi[i];}
153: }
154: if (ndptx<2) xi[1] = PETSC_INFINITY;
156: beta[1] = maxnrs;
157: for (k=2;k<ctx->ddmaxit;k++) {
158: maxnrs = 0.0;
159: minnrxi = PETSC_MAX_REAL;
160: for (i=0;i<ndpt;i++) {
161: nrs[i] *= ((ds[i]-s[k-1])/(1.0-ds[i]/xi[k-1]))/beta[k-1];
162: if (PetscAbsScalar(nrs[i])>maxnrs) {maxnrs = PetscAbsScalar(nrs[i]); s[k] = ds[i];}
163: }
164: if (ndptx>=k) {
165: for (i=1;i<ndptx;i++) {
166: nrxi[i] *= ((dxi[i]-s[k-1])/(1.0-dxi[i]/xi[k-1]))/beta[k-1];
167: if (PetscAbsScalar(nrxi[i])<minnrxi) {minnrxi = PetscAbsScalar(nrxi[i]); xi[k] = dxi[i];}
168: }
169: } else xi[k] = PETSC_INFINITY;
170: beta[k] = maxnrs;
171: }
172: PetscFree5(ds,dsi,dxi,nrs,nrxi);
173: return(0);
174: }
178: static PetscErrorCode NEPNLEIGSEvalNRTFunct(NEP nep,PetscInt k,PetscScalar sigma,PetscScalar *b)
179: {
180: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
181: PetscInt i;
182: PetscScalar *beta=ctx->beta,*s=ctx->s,*xi=ctx->xi;
185: b[0] = 1.0/beta[0];
186: for (i=0;i<k;i++) {
187: b[i+1] = ((sigma-s[i])*b[i])/(beta[i+1]*(1.0-sigma/xi[i]));
188: }
189: return(0);
190: }
194: static PetscErrorCode MatMult_Fun(Mat A,Vec x,Vec y)
195: {
197: ShellMatCtx *ctx;
198: PetscInt i;
201: MatShellGetContext(A,(void**)&ctx);
202: MatMult(ctx->A[0],x,y);
203: if (ctx->coeff[0]!=1.0) { VecScale(y,ctx->coeff[0]); }
204: for (i=1;i<ctx->nmat;i++) {
205: MatMult(ctx->A[i],x,ctx->t);
206: VecAXPY(y,ctx->coeff[i],ctx->t);
207: }
208: return(0);
209: }
213: static PetscErrorCode MatMultTranspose_Fun(Mat A,Vec x,Vec y)
214: {
216: ShellMatCtx *ctx;
217: PetscInt i;
220: MatShellGetContext(A,(void**)&ctx);
221: MatMultTranspose(ctx->A[0],x,y);
222: if (ctx->coeff[0]!=1.0) { VecScale(y,ctx->coeff[0]); }
223: for (i=1;i<ctx->nmat;i++) {
224: MatMultTranspose(ctx->A[i],x,ctx->t);
225: VecAXPY(y,ctx->coeff[i],ctx->t);
226: }
227: return(0);
228: }
232: static PetscErrorCode MatGetDiagonal_Fun(Mat A,Vec diag)
233: {
235: ShellMatCtx *ctx;
236: PetscInt i;
239: MatShellGetContext(A,(void**)&ctx);
240: MatGetDiagonal(ctx->A[0],diag);
241: if (ctx->coeff[0]!=1.0) { VecScale(diag,ctx->coeff[0]); }
242: for (i=1;i<ctx->nmat;i++) {
243: MatGetDiagonal(ctx->A[i],ctx->t);
244: VecAXPY(diag,ctx->coeff[i],ctx->t);
245: }
246: return(0);
247: }
251: static PetscErrorCode MatDuplicate_Fun(Mat A,MatDuplicateOption op,Mat *B)
252: {
253: PetscInt n,i;
254: ShellMatCtx *ctxnew,*ctx;
255: void (*fun)();
259: MatShellGetContext(A,(void**)&ctx);
260: PetscNew(&ctxnew);
261: ctxnew->nmat = ctx->nmat;
262: for (i=0;i<ctx->nmat;i++) {
263: PetscObjectReference((PetscObject)ctx->A[i]);
264: ctxnew->A[i] = ctx->A[i];
265: ctxnew->coeff[i] = ctx->coeff[i];
266: }
267: MatGetSize(ctx->A[0],&n,NULL);
268: VecDuplicate(ctx->t,&ctxnew->t);
269: MatCreateShell(PETSC_COMM_WORLD,n,n,n,n,(void*)ctxnew,B);
270: MatShellGetOperation(A,MATOP_MULT,&fun);
271: MatShellSetOperation(*B,MATOP_MULT,fun);
272: MatShellGetOperation(A,MATOP_MULT_TRANSPOSE,&fun);
273: MatShellSetOperation(*B,MATOP_MULT_TRANSPOSE,fun);
274: MatShellGetOperation(A,MATOP_GET_DIAGONAL,&fun);
275: MatShellSetOperation(*B,MATOP_GET_DIAGONAL,fun);
276: MatShellGetOperation(A,MATOP_DUPLICATE,&fun);
277: MatShellSetOperation(*B,MATOP_DUPLICATE,fun);
278: MatShellGetOperation(A,MATOP_DESTROY,&fun);
279: MatShellSetOperation(*B,MATOP_DESTROY,fun);
280: MatShellGetOperation(A,MATOP_AXPY,&fun);
281: MatShellSetOperation(*B,MATOP_AXPY,fun);
282: return(0);
283: }
287: static PetscErrorCode MatDestroy_Fun(Mat A)
288: {
289: ShellMatCtx *ctx;
291: PetscInt i;
294: if (A) {
295: MatShellGetContext(A,(void**)&ctx);
296: for (i=0;i<ctx->nmat;i++) {
297: MatDestroy(&ctx->A[i]);
298: }
299: VecDestroy(&ctx->t);
300: PetscFree(ctx);
301: }
302: return(0);
303: }
307: static PetscErrorCode MatAXPY_Fun(Mat Y,PetscScalar a,Mat X,MatStructure str)
308: {
309: ShellMatCtx *ctxY,*ctxX;
311: PetscInt i,j;
312: PetscBool found;
315: MatShellGetContext(Y,(void**)&ctxY);
316: MatShellGetContext(X,(void**)&ctxX);
317: for (i=0;i<ctxX->nmat;i++) {
318: found = PETSC_FALSE;
319: for (j=0;!found&&j<ctxY->nmat;j++) {
320: if (ctxX->A[i]==ctxY->A[j]) {
321: found = PETSC_TRUE;
322: ctxY->coeff[j] += a*ctxX->coeff[i];
323: }
324: }
325: if (!found) {
326: ctxY->coeff[ctxY->nmat] = a*ctxX->coeff[i];
327: ctxY->A[ctxY->nmat++] = ctxX->A[i];
328: PetscObjectReference((PetscObject)ctxX->A[i]);
329: }
330: }
331: return(0);
332: }
336: static PetscErrorCode MatScale_Fun(Mat M,PetscScalar a)
337: {
338: ShellMatCtx *ctx;
340: PetscInt i;
343: MatShellGetContext(M,(void**)&ctx);
344: for (i=0;i<ctx->nmat;i++) ctx->coeff[i] *= a;
345: return(0);
346: }
350: static PetscErrorCode NLEIGSMatToMatShellArray(Mat M,Mat *Ms)
351: {
353: ShellMatCtx *ctx;
354: PetscInt n;
355: PetscBool has;
356:
358: MatHasOperation(M,MATOP_DUPLICATE,&has);
359: if (!has) SETERRQ(PetscObjectComm((PetscObject)M),1,"MatDuplicate operation required");
360: PetscNew(&ctx);
361: MatDuplicate(M,MAT_COPY_VALUES,&ctx->A[0]);
362: ctx->nmat = 1;
363: ctx->coeff[0] = 1.0;
364: MatCreateVecs(M,&ctx->t,NULL);
365: MatGetSize(M,&n,NULL);
366: MatCreateShell(PetscObjectComm((PetscObject)M),n,n,n,n,(void*)ctx,Ms);
367: MatShellSetOperation(*Ms,MATOP_MULT,(void(*)())MatMult_Fun);
368: MatShellSetOperation(*Ms,MATOP_MULT_TRANSPOSE,(void(*)())MatMultTranspose_Fun);
369: MatShellSetOperation(*Ms,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Fun);
370: MatShellSetOperation(*Ms,MATOP_DUPLICATE,(void(*)())MatDuplicate_Fun);
371: MatShellSetOperation(*Ms,MATOP_DESTROY,(void(*)())MatDestroy_Fun);
372: MatShellSetOperation(*Ms,MATOP_AXPY,(void(*)())MatAXPY_Fun);
373: MatShellSetOperation(*Ms,MATOP_SCALE,(void(*)())MatScale_Fun);
374: return(0);
375: }
379: static PetscErrorCode NEPNLEIGSNormEstimation(NEP nep,Mat M,PetscReal *norm,Vec *w)
380: {
381: PetscScalar *z,*x,*y;
382: PetscReal tr;
383: Vec X=w[0],Y=w[1];
384: PetscInt n,i;
385: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
386: PetscRandom rand;
388:
390: if (!ctx->vrn) {
391: /* generate a random vector with normally distributed entries with the Box-Muller transform */
392: BVGetRandomContext(nep->V,&rand);
393: MatCreateVecs(M,&ctx->vrn,NULL);
394: VecSetRandom(X,rand);
395: VecSetRandom(Y,rand);
396: VecGetLocalSize(ctx->vrn,&n);
397: VecGetArray(ctx->vrn,&z);
398: VecGetArray(X,&x);
399: VecGetArray(Y,&y);
400: for (i=0;i<n;i++) {
401: #if defined(PETSC_USE_COMPLEX)
402: z[i] = PetscSqrtReal(-2.0*PetscLogReal(PetscRealPart(x[i])))*PetscCosReal(2.0*PETSC_PI*PetscRealPart(y[i]));
403: z[i] += PETSC_i*(PetscSqrtReal(-2.0*PetscLogReal(PetscImaginaryPart(x[i])))*PetscCosReal(2.0*PETSC_PI*PetscImaginaryPart(y[i])));
404: #else
405: z[i] = PetscSqrtReal(-2.0*PetscLogReal(x[i]))*PetscCosReal(2.0*PETSC_PI*y[i]);
406: #endif
407: }
408: VecRestoreArray(ctx->vrn,&z);
409: VecRestoreArray(X,&x);
410: VecRestoreArray(Y,&y);
411: VecNorm(ctx->vrn,NORM_2,&tr);
412: VecScale(ctx->vrn,1/tr);
413: }
414: /* matrix-free norm estimator of Ipsen http://www4.ncsu.edu/~ipsen/ps/slides_ima.pdf */
415: MatGetSize(M,&n,NULL);
416: MatMult(M,ctx->vrn,X);
417: VecNorm(X,NORM_2,norm);
418: *norm *= PetscSqrtReal((PetscReal)n);
419: return(0);
420: }
424: static PetscErrorCode NEPNLEIGSDividedDifferences_split(NEP nep)
425: {
427: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
428: PetscInt k,j,i;
429: PetscReal norm0,norm,max;
430: PetscScalar *s=ctx->s,*beta=ctx->beta,*b,alpha,*coeffs;
431: Mat T,Ts;
432: PetscBool shell;
435: PetscMalloc1(nep->nt*ctx->ddmaxit,&ctx->coeffD);
436: PetscMalloc2(ctx->ddmaxit+1,&b,ctx->ddmaxit+1,&coeffs);
437: max = 0.0;
438: for (j=0;j<nep->nt;j++) {
439: FNEvaluateFunction(nep->f[j],s[0],ctx->coeffD+j);
440: ctx->coeffD[j] /= beta[0];
441: max = PetscMax(PetscAbsScalar(ctx->coeffD[j]),max);
442: }
443: norm0 = max;
444: ctx->nmat = ctx->ddmaxit;
445: for (k=1;k<ctx->ddmaxit;k++) {
446: NEPNLEIGSEvalNRTFunct(nep,k,s[k],b);
447: max = 0.0;
448: for (i=0;i<nep->nt;i++) {
449: FNEvaluateFunction(nep->f[i],s[k],ctx->coeffD+k*nep->nt+i);
450: for (j=0;j<k;j++) {
451: ctx->coeffD[k*nep->nt+i] -= b[j]*ctx->coeffD[i+nep->nt*j];
452: }
453: ctx->coeffD[k*nep->nt+i] /= b[k];
454: max = PetscMax(PetscAbsScalar(ctx->coeffD[k*nep->nt+i]),max);
455: }
456: norm = max;
457: if (norm/norm0 < ctx->ddtol) {
458: ctx->nmat = k+1;
459: break;
460: }
461: }
462: if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
463: PetscObjectTypeCompare((PetscObject)nep->A[0],MATSHELL,&shell);
464: for (i=0;i<ctx->nshiftsw;i++) {
465: NEPNLEIGSEvalNRTFunct(nep,ctx->nmat,ctx->shifts[i],coeffs);
466: if (!shell) {
467: MatDuplicate(nep->A[0],MAT_COPY_VALUES,&T);
468: } else {
469: NLEIGSMatToMatShellArray(nep->A[0],&T);
470: }
471: alpha = 0.0;
472: for (j=0;j<ctx->nmat-1;j++) alpha += coeffs[j]*ctx->coeffD[j*nep->nt];
473: MatScale(T,alpha);
474: for (k=1;k<nep->nt;k++) {
475: alpha = 0.0;
476: for (j=0;j<ctx->nmat-1;j++) alpha += coeffs[j]*ctx->coeffD[j*nep->nt+k];
477: if (shell) {
478: NLEIGSMatToMatShellArray(nep->A[k],&Ts);
479: }
480: MatAXPY(T,alpha,shell?Ts:nep->A[k],nep->mstr);
481: if (shell) {
482: MatDestroy(&Ts);
483: }
484: }
485: KSPSetOperators(ctx->ksp[i],T,T);
486: KSPSetUp(ctx->ksp[i]);
487: MatDestroy(&T);
488: }
489: PetscFree2(b,coeffs);
490: return(0);
491: }
495: static PetscErrorCode NEPNLEIGSDividedDifferences_callback(NEP nep)
496: {
498: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
499: PetscInt k,j,i;
500: PetscReal norm0,norm;
501: PetscScalar *s=ctx->s,*beta=ctx->beta,*b,*coeffs;
502: Mat *D=ctx->D,T;
503: PetscBool shell,has,vec=PETSC_FALSE;
504: Vec w[2];
507: PetscMalloc2(ctx->ddmaxit+1,&b,ctx->ddmaxit+1,&coeffs);
508: T = nep->function;
509: NEPComputeFunction(nep,s[0],T,T);
510: PetscObjectTypeCompare((PetscObject)T,MATSHELL,&shell);
511: if (!shell) {
512: MatDuplicate(T,MAT_COPY_VALUES,&D[0]);
513: } else {
514: NLEIGSMatToMatShellArray(T,&D[0]);
515: }
516: if (beta[0]!=1.0) {
517: MatScale(D[0],1.0/beta[0]);
518: }
519: MatHasOperation(D[0],MATOP_NORM,&has);
520: if (has) {
521: MatNorm(D[0],NORM_FROBENIUS,&norm0);
522: } else {
523: MatCreateVecs(D[0],NULL,&w[0]);
524: VecDuplicate(w[0],&w[1]);
525: vec = PETSC_TRUE;
526: NEPNLEIGSNormEstimation(nep,D[0],&norm0,w);
527: }
528: ctx->nmat = ctx->ddmaxit;
529: for (k=1;k<ctx->ddmaxit;k++) {
530: NEPNLEIGSEvalNRTFunct(nep,k,s[k],b);
531: NEPComputeFunction(nep,s[k],T,T);
532: if (!shell) {
533: MatDuplicate(T,MAT_COPY_VALUES,&D[k]);
534: } else {
535: NLEIGSMatToMatShellArray(T,&D[k]);
536: }
537: for (j=0;j<k;j++) {
538: MatAXPY(D[k],-b[j],D[j],nep->mstr);
539: }
540: MatScale(D[k],1.0/b[k]);
541: MatHasOperation(D[k],MATOP_NORM,&has);
542: if (has) {
543: MatNorm(D[k],NORM_FROBENIUS,&norm);
544: } else {
545: if(!vec) {
546: MatCreateVecs(D[k],NULL,&w[0]);
547: VecDuplicate(w[0],&w[1]);
548: vec = PETSC_TRUE;
549: }
550: NEPNLEIGSNormEstimation(nep,D[k],&norm,w);
551: }
552: if (norm/norm0 < ctx->ddtol) {
553: ctx->nmat = k+1; /* increment (the last matrix is not used) */
554: MatDestroy(&D[k]);
555: break;
556: }
557: }
558: if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
559: for (i=0;i<ctx->nshiftsw;i++) {
560: NEPNLEIGSEvalNRTFunct(nep,ctx->nmat,ctx->shifts[i],coeffs);
561: MatDuplicate(ctx->D[0],MAT_COPY_VALUES,&T);
562: if (coeffs[0]!=1.0) { MatScale(T,coeffs[0]); }
563: for (j=1;j<ctx->nmat-1;j++) {
564: MatAXPY(T,coeffs[j],ctx->D[j],nep->mstr);
565: }
566: KSPSetOperators(ctx->ksp[i],T,T);
567: KSPSetUp(ctx->ksp[i]);
568: MatDestroy(&T);
569: }
570: PetscFree2(b,coeffs);
571: if (vec) {
572: VecDestroy(&w[0]);
573: VecDestroy(&w[1]);
574: }
575: return(0);
576: }
580: static PetscErrorCode NEPNLEIGSRitzVector(NEP nep,PetscScalar *S,PetscInt ld,PetscInt nq,PetscScalar *H,PetscInt k,Vec t)
581: {
583: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
584: PetscInt deg=ctx->nmat-1,ldds,n;
585: PetscBLASInt nq_,lds_,n_,one=1,ldds_;
586: PetscScalar sone=1.0,zero=0.0,*x,*y,*X;
589: DSGetDimensions(nep->ds,&n,NULL,NULL,NULL,NULL);
590: PetscMalloc1(n+nq,&y);
591: x = y+nq;
592: DSGetLeadingDimension(nep->ds,&ldds);
593: PetscBLASIntCast(n,&n_);
594: PetscBLASIntCast(nq,&nq_);
595: PetscBLASIntCast(ldds,&ldds_);
596: PetscBLASIntCast(deg*ld,&lds_);
597: DSGetArray(nep->ds,DS_MAT_X,&X);
598: if (ctx->nshifts) PetscStackCall("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,H,&ldds_,X+k*ldds,&one,&zero,x,&one));
599: else x = X+k*ldds;
600: PetscStackCall("BLASgemv",BLASgemv_("N",&nq_,&n_,&sone,S,&lds_,x,&one,&zero,y,&one));
601: DSRestoreArray(nep->ds,DS_MAT_X,&X);
602: BVSetActiveColumns(nep->V,0,nq);
603: BVMultVec(nep->V,1.0,0.0,t,y);
604: VecNormalize(t,NULL);
605: PetscFree(y);
606: return(0);
607: }
611: /*
612: NEPKrylovConvergence - This is the analogue to EPSKrylovConvergence.
613: */
614: 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)
615: {
617: PetscInt k,newk,marker,inside;
618: PetscScalar re,im;
619: PetscReal resnorm,tt;
620: PetscBool istrivial;
621: Vec t;
622: NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;
625: t = w[0];
626: RGIsTrivial(nep->rg,&istrivial);
627: marker = -1;
628: if (nep->trackall) getall = PETSC_TRUE;
629: for (k=kini;k<kini+nits;k++) {
630: /* eigenvalue */
631: re = nep->eigr[k];
632: im = nep->eigi[k];
633: if (!istrivial) {
634: if (!ctx->nshifts) {
635: NEPNLEIGSBackTransform((PetscObject)nep,1,&re,&im);
636: }
637: RGCheckInside(nep->rg,1,&re,&im,&inside);
638: if (marker==-1 && inside<0) marker = k;
639: }
640: newk = k;
641: DSVectors(nep->ds,DS_MAT_X,&newk,&resnorm);
642: tt = (ctx->nshifts)?SlepcAbsEigenvalue(betak-nep->eigr[k]*betah,nep->eigi[k]*betah):betah;
643: resnorm *= PetscAbsReal(tt);
644: /* error estimate */
645: (*nep->converged)(nep,nep->eigr[k],nep->eigi[k],resnorm,&nep->errest[k],nep->convergedctx);
646: if (ctx->trueres && (nep->errest[k] < nep->tol) ) {
647: /* check explicit residual */
648: NEPNLEIGSRitzVector(nep,S,ld,nq,H,k,t);
649: NEPComputeResidualNorm_Private(nep,re,t,w+1,&resnorm);
650: (*nep->converged)(nep,re,im,resnorm,&nep->errest[k],nep->convergedctx);
651: }
652: if (marker==-1 && nep->errest[k] >= nep->tol) marker = k;
653: if (newk==k+1) {
654: nep->errest[k+1] = nep->errest[k];
655: k++;
656: }
657: if (marker!=-1 && !getall) break;
658: }
659: if (marker!=-1) k = marker;
660: *kout = k;
661: return(0);
662: }
666: PetscErrorCode NEPSetUp_NLEIGS(NEP nep)
667: {
669: PetscInt k,in;
670: PetscScalar zero=0.0;
671: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
672: SlepcSC sc;
673: PetscBool istrivial;
676: NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd);
677: if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
678: if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
679: if (!ctx->ddmaxit) ctx->ddmaxit = MAX_LBPOINTS;
680: RGIsTrivial(nep->rg,&istrivial);
681: if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPNLEIGS requires a nontrivial region defining the target set");
682: RGCheckInside(nep->rg,1,&nep->target,&zero,&in);
683: if (in<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The target is not inside the target set");
684: if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
686: /* Initialize the NLEIGS context structure */
687: k = ctx->ddmaxit;
688: PetscMalloc4(k,&ctx->s,k,&ctx->xi,k,&ctx->beta,k,&ctx->D);
689: nep->data = ctx;
690: if (nep->tol==PETSC_DEFAULT) nep->tol = SLEPC_DEFAULT_TOL;
691: if (ctx->ddtol==PETSC_DEFAULT) ctx->ddtol = nep->tol/10.0;
692: if (!ctx->keep) ctx->keep = 0.5;
694: /* Compute Leja-Bagby points and scaling values */
695: NEPNLEIGSLejaBagbyPoints(nep);
697: /* Compute the divided difference matrices */
698: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
699: NEPNLEIGSDividedDifferences_split(nep);
700: } else {
701: NEPNLEIGSDividedDifferences_callback(nep);
702: }
703: NEPAllocateSolution(nep,ctx->nmat);
704: NEPSetWorkVecs(nep,4);
706: /* set-up DS and transfer split operator functions */
707: DSSetType(nep->ds,ctx->nshifts?DSGNHEP:DSNHEP);
708: DSAllocate(nep->ds,nep->ncv+1);
709: DSGetSlepcSC(nep->ds,&sc);
710: if (!ctx->nshifts) {
711: sc->map = NEPNLEIGSBackTransform;
712: DSSetExtraRow(nep->ds,PETSC_TRUE);
713: }
714: sc->mapobj = (PetscObject)nep;
715: sc->rg = nep->rg;
716: sc->comparison = nep->sc->comparison;
717: sc->comparisonctx = nep->sc->comparisonctx;
718: return(0);
719: }
723: /*
724: Norm of [sp;sq]
725: */
726: static PetscErrorCode NEPTOARSNorm2(PetscInt n,PetscScalar *S,PetscReal *norm)
727: {
729: PetscBLASInt n_,one=1;
732: PetscBLASIntCast(n,&n_);
733: *norm = BLASnrm2_(&n_,S,&one);
734: return(0);
735: }
739: /*
740: Computes GS orthogonalization [z;x] - [Sp;Sq]*y,
741: where y = ([Sp;Sq]'*[z;x]).
742: k: Column from S to be orthogonalized against previous columns.
743: Sq = Sp+ld
744: dim(work)=k;
745: */
746: static PetscErrorCode NEPTOAROrth2(NEP nep,PetscScalar *S,PetscInt ld,PetscInt deg,PetscInt k,PetscScalar *y,PetscReal *norm,PetscBool *lindep,PetscScalar *work)
747: {
749: PetscBLASInt n_,lds_,k_,one=1;
750: PetscScalar sonem=-1.0,sone=1.0,szero=0.0,*x0,*x,*c;
751: PetscInt i,lds=deg*ld,n;
752: PetscReal eta,onorm;
753:
755: BVGetOrthogonalization(nep->V,NULL,NULL,&eta,NULL);
756: n = k+deg-1;
757: PetscBLASIntCast(n,&n_);
758: PetscBLASIntCast(deg*ld,&lds_);
759: PetscBLASIntCast(k,&k_); /* Number of vectors to orthogonalize against them */
760: c = work;
761: x0 = S+k*lds;
762: PetscStackCall("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S,&lds_,x0,&one,&szero,y,&one));
763: for (i=1;i<deg;i++) {
764: x = S+i*ld+k*lds;
765: PetscStackCall("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+i*ld,&lds_,x,&one,&sone,y,&one));
766: }
767: for (i=0;i<deg;i++) {
768: x= S+i*ld+k*lds;
769: PetscStackCall("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+i*ld,&lds_,y,&one,&sone,x,&one));
770: }
771: NEPTOARSNorm2(lds,S+k*lds,&onorm);
772: /* twice */
773: PetscStackCall("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S,&lds_,x0,&one,&szero,c,&one));
774: for (i=1;i<deg;i++) {
775: x = S+i*ld+k*lds;
776: PetscStackCall("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+i*ld,&lds_,x,&one,&sone,c,&one));
777: }
778: for (i=0;i<deg;i++) {
779: x= S+i*ld+k*lds;
780: PetscStackCall("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+i*ld,&lds_,c,&one,&sone,x,&one));
781: }
782: for (i=0;i<k;i++) y[i] += c[i];
783: if (norm) {
784: NEPTOARSNorm2(lds,S+k*lds,norm);
785: if (lindep) *lindep = (*norm < eta * onorm)?PETSC_TRUE:PETSC_FALSE;
786: }
787: return(0);
788: }
792: /*
793: Extend the TOAR basis by applying the the matrix operator
794: over a vector which is decomposed on the TOAR way
795: Input:
796: - S,V: define the latest Arnoldi vector (nv vectors in V)
797: Output:
798: - t: new vector extending the TOAR basis
799: - r: temporally coefficients to compute the TOAR coefficients
800: for the new Arnoldi vector
801: Workspace: t_ (two vectors)
802: */
803: static PetscErrorCode NEPTOARExtendBasis(NEP nep,PetscInt idxrktg,PetscScalar *S,PetscInt ls,PetscInt nv,BV V,Vec t,PetscScalar *r,PetscInt lr,Vec *t_)
804: {
806: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
807: PetscInt deg=ctx->nmat-1,k,j;
808: Vec v=t_[0],q=t_[1],w;
809: PetscScalar *beta=ctx->beta,*s=ctx->s,*xi=ctx->xi,*coeffs,sigma;
812: if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
813: sigma = ctx->shifts[idxrktg];
814: BVSetActiveColumns(nep->V,0,nv);
815: PetscMalloc1(ctx->nmat-1,&coeffs);
816: if (PetscAbsScalar(s[deg-2]-sigma)<100*PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1,"Breakdown in NLEIGS");
817: /* i-part stored in (i-1) position */
818: for (j=0;j<nv;j++) {
819: 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);
820: }
821: BVSetActiveColumns(ctx->W,0,deg-1);
822: BVGetColumn(ctx->W,deg-2,&w);
823: BVMultVec(V,1.0,0.0,w,r+(deg-2)*lr);
824: BVRestoreColumn(ctx->W,deg-2,&w);
825: for (k=deg-2;k>0;k--) {
826: if (PetscAbsScalar(s[k-1]-sigma)<100*PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1,"Breakdown in NLEIGS");
827: 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);
828: BVGetColumn(ctx->W,k-1,&w);
829: BVMultVec(V,1.0,0.0,w,r+(k-1)*lr);
830: BVRestoreColumn(ctx->W,k-1,&w);
831: }
832: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
833: for (j=0;j<ctx->nmat-1;j++) coeffs[j] = ctx->coeffD[nep->nt*j];
834: BVMultVec(ctx->W,1.0,0.0,v,coeffs);
835: MatMult(nep->A[0],v,q);
836: for (k=1;k<nep->nt;k++) {
837: for (j=0;j<ctx->nmat-1;j++) coeffs[j] = ctx->coeffD[nep->nt*j+k];
838: BVMultVec(ctx->W,1.0,0,v,coeffs);
839: MatMult(nep->A[k],v,t);
840: VecAXPY(q,1.0,t);
841: }
842: KSPSolve(ctx->ksp[idxrktg],q,t);
843: VecScale(t,-1.0);
844: } else {
845: for (k=0;k<deg-1;k++) {
846: BVGetColumn(ctx->W,k,&w);
847: MatMult(ctx->D[k],w,q);
848: BVRestoreColumn(ctx->W,k,&w);
849: BVInsertVec(ctx->W,k,q);
850: }
851: for (j=0;j<ctx->nmat-1;j++) coeffs[j] = 1.0;
852: BVMultVec(ctx->W,1.0,0.0,q,coeffs);
853: KSPSolve(ctx->ksp[idxrktg],q,t);
854: VecScale(t,-1.0);
855: }
856: PetscFree(coeffs);
857: return(0);
858: }
862: /*
863: Compute TOAR coefficients of the blocks of the new Arnoldi vector computed
864: */
865: static PetscErrorCode NEPTOARCoefficients(NEP nep,PetscScalar sigma,PetscInt nv,PetscScalar *S,PetscInt ls,PetscScalar *r,PetscInt lr,PetscScalar *x,PetscScalar *work)
866: {
868: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
869: PetscInt k,j,d=ctx->nmat-1;
870: PetscScalar *t=work;
873: NEPNLEIGSEvalNRTFunct(nep,d-1,sigma,t);
874: for (k=0;k<d-1;k++) {
875: for (j=0;j<=nv;j++) r[k*lr+j] += t[k]*x[j];
876: }
877: for (j=0;j<=nv;j++) r[(d-1)*lr+j] = t[d-1]*x[j];
878: return(0);
879: }
883: /*
884: Compute continuation vector coefficients for the Rational-Krylov run.
885: dim(work) >= (end-ini)*(end-ini+1) + end+1 + 2*(end-ini+1), dim(t) = end.
886: */
887: 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)
888: {
889: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(SLEPC_MISSING_LAPACK_LARF)
891: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/LARF - Lapack routines are unavailable");
892: #else
894: PetscScalar *x,*W,*tau,sone=1.0,szero=0.0;
895: PetscInt i,j,n1,n,nwu=0;
896: PetscBLASInt info,n_,n1_,one=1,dim,lds_;
897: NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;
900: if (!ctx->nshifts || !end) {
901: t[0] = 1;
902: PetscMemcpy(cont,S+end*lds,lds*sizeof(PetscScalar));
903: } else {
904: n = end-ini;
905: n1 = n+1;
906: x = work+nwu;
907: nwu += end+1;
908: tau = work+nwu;
909: nwu += n;
910: W = work+nwu;
911: nwu += n1*n;
912: for (j=ini;j<end;j++) {
913: for (i=ini;i<=end;i++) W[(j-ini)*n1+i-ini] = K[j*ld+i] -H[j*ld+i]*sigma;
914: }
915: PetscBLASIntCast(n,&n_);
916: PetscBLASIntCast(n1,&n1_);
917: PetscBLASIntCast(end+1,&dim);
918: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&n1_,&n_,W,&n1_,tau,work+nwu,&n1_,&info));
919: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
920: for (i=0;i<end;i++) t[i] = 0.0;
921: t[end] = 1.0;
922: for (j=n-1;j>=0;j--) {
923: for (i=0;i<ini+j;i++) x[i] = 0.0;
924: x[ini+j] = 1.0;
925: for (i=j+1;i<n1;i++) x[i+ini] = W[i+n1*j];
926: tau[j] = PetscConj(tau[j]);
927: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&dim,&one,x,&one,tau+j,t,&dim,work+nwu));
928: }
929: PetscBLASIntCast(lds,&lds_);
930: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&n1_,&sone,S,&lds_,t,&one,&szero,cont,&one));
931: }
932: return(0);
933: #endif
934: }
938: /*
939: Compute a run of Arnoldi iterations
940: */
941: static PetscErrorCode NEPNLEIGSTOARrun(NEP nep,PetscInt *nq,PetscScalar *S,PetscInt ld,PetscScalar *K,PetscScalar *H,PetscInt ldh,BV V,PetscInt k,PetscInt *M,PetscBool *breakdown,Vec *t_)
942: {
944: NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;
945: PetscInt i,j,p,m=*M,lwa,deg=ctx->nmat-1,lds=ld*deg,nqt=*nq;
946: Vec t=t_[0];
947: PetscReal norm;
948: PetscScalar *x,*work,*tt,sigma,*cont;
949: PetscBool lindep;
952: lwa = PetscMax(ld,deg)+(m+1)*(m+1)+4*(m+1);
953: PetscMalloc4(ld,&x,lwa,&work,m+1,&tt,lds,&cont);
954: for (j=k;j<m;j++) {
955: sigma = ctx->shifts[(++(ctx->idxrk))%ctx->nshiftsw];
957: /* Continuation vector */
958: NEPNLEIGS_RKcontinuation(nep,0,j,K,H,ldh,sigma,S,lds,cont,tt,work);
959:
960: /* apply operator */
961: BVGetColumn(nep->V,nqt,&t);
962: NEPTOARExtendBasis(nep,(ctx->idxrk)%ctx->nshiftsw,cont,ld,nqt,V,t,S+(j+1)*lds,ld,t_+1);
963: BVRestoreColumn(nep->V,nqt,&t);
965: /* orthogonalize */
966: BVOrthogonalizeColumn(nep->V,nqt,x,&norm,&lindep);
967: if (!lindep) {
968: x[nqt] = norm;
969: BVScaleColumn(nep->V,nqt,1.0/norm);
970: nqt++;
971: } else x[nqt] = 0.0;
973: NEPTOARCoefficients(nep,sigma,*nq,cont,ld,S+(j+1)*lds,ld,x,work);
975: /* Level-2 orthogonalization */
976: NEPTOAROrth2(nep,S,ld,deg,j+1,H+j*ldh,&norm,breakdown,work);
977: H[j+1+ldh*j] = norm;
978: if (ctx->nshifts) {
979: for (i=0;i<=j;i++) K[i+ldh*j] = sigma*H[i+ldh*j] + tt[i];
980: K[j+1+ldh*j] = sigma*H[j+1+ldh*j];
981: }
982: *nq = nqt;
983: if (*breakdown) {
984: *M = j+1;
985: break;
986: }
987: for (p=0;p<deg;p++) {
988: for (i=0;i<=j+deg;i++) {
989: S[i+p*ld+(j+1)*lds] /= norm;
990: }
991: }
992: }
993: PetscFree4(x,work,tt,cont);
994: return(0);
995: }
999: /* dim(work)=5*ld*lds dim(rwork)=6*n */
1000: static PetscErrorCode NEPTOARTrunc(NEP nep,PetscScalar *S,PetscInt ld,PetscInt deg,PetscInt *nq,PetscInt cs1,PetscScalar *work,PetscReal *rwork)
1001: {
1003: PetscInt lwa,nwu=0,nrwu=0;
1004: PetscInt j,i,n,lds=deg*ld,rk=0,rs1;
1005: PetscScalar *M,*V,*pU,t;
1006: PetscReal *sg,tol;
1007: PetscBLASInt cs1_,rs1_,cs1tdeg,n_,info,lw_;
1008: Mat U;
1011: rs1 = *nq;
1012: n = (rs1>deg*cs1)?deg*cs1:rs1;
1013: lwa = 5*ld*lds;
1014: M = work+nwu;
1015: nwu += rs1*cs1*deg;
1016: sg = rwork+nrwu;
1017: nrwu += n;
1018: pU = work+nwu;
1019: nwu += rs1*n;
1020: V = work+nwu;
1021: nwu += deg*cs1*n;
1022: for (i=0;i<cs1;i++) {
1023: for (j=0;j<deg;j++) {
1024: PetscMemcpy(M+(i+j*cs1)*rs1,S+i*lds+j*ld,rs1*sizeof(PetscScalar));
1025: }
1026: }
1027: PetscBLASIntCast(n,&n_);
1028: PetscBLASIntCast(cs1,&cs1_);
1029: PetscBLASIntCast(rs1,&rs1_);
1030: PetscBLASIntCast(cs1*deg,&cs1tdeg);
1031: PetscBLASIntCast(lwa-nwu,&lw_);
1032: #if !defined (PETSC_USE_COMPLEX)
1033: PetscStackCall("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1tdeg,M,&rs1_,sg,pU,&rs1_,V,&n_,work+nwu,&lw_,&info));
1034: #else
1035: PetscStackCall("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1tdeg,M,&rs1_,sg,pU,&rs1_,V,&n_,work+nwu,&lw_,rwork+nrwu,&info));
1036: #endif
1037: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESVD %d",info);
1038:
1039: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1040: MatCreateSeqDense(PETSC_COMM_SELF,rs1,cs1+deg-1,pU,&U);
1041: BVSetActiveColumns(nep->V,0,rs1);
1042: BVMultInPlace(nep->V,U,0,cs1+deg-1);
1043: BVSetActiveColumns(nep->V,0,cs1+deg-1);
1044: MatDestroy(&U);
1045: tol = PetscMax(rs1,deg*cs1)*PETSC_MACHINE_EPSILON*sg[0];
1046: for (i=0;i<PetscMin(n_,cs1tdeg);i++) if (sg[i]>tol) rk++;
1047: rk = PetscMin(cs1+deg-1,rk);
1048:
1049: /* Update S */
1050: PetscMemzero(S,lds*ld*sizeof(PetscScalar));
1051: for (i=0;i<rk;i++) {
1052: t = sg[i];
1053: PetscStackCall("BLASscal",BLASscal_(&cs1tdeg,&t,V+i,&n_));
1054: }
1055: for (j=0;j<cs1;j++) {
1056: for (i=0;i<deg;i++) {
1057: PetscMemcpy(S+j*lds+i*ld,V+(cs1*i+j)*n,rk*sizeof(PetscScalar));
1058: }
1059: }
1060: *nq = rk;
1061: return(0);
1062: }
1066: /*
1067: S <- S*Q
1068: columns s-s+ncu of S
1069: rows 0-sr of S
1070: size(Q) qr x ncu
1071: dim(work)=sr*ncu
1072: */
1073: static PetscErrorCode NEPTOARSupdate(PetscScalar *S,PetscInt ld,PetscInt deg,PetscInt sr,PetscInt s,PetscInt ncu,PetscInt qr,PetscScalar *Q,PetscInt ldq,PetscScalar *work)
1074: {
1076: PetscScalar a=1.0,b=0.0;
1077: PetscBLASInt sr_,ncu_,ldq_,lds_,qr_;
1078: PetscInt j,lds=deg*ld,i;
1081: PetscBLASIntCast(sr,&sr_);
1082: PetscBLASIntCast(qr,&qr_);
1083: PetscBLASIntCast(ncu,&ncu_);
1084: PetscBLASIntCast(lds,&lds_);
1085: PetscBLASIntCast(ldq,&ldq_);
1086: for (i=0;i<deg;i++) {
1087: PetscStackCall("BLASgemm",BLASgemm_("N","N",&sr_,&ncu_,&qr_,&a,S+i*ld,&lds_,Q,&ldq_,&b,work,&sr_));
1088: for (j=0;j<ncu;j++) {
1089: PetscMemcpy(S+lds*(s+j)+i*ld,work+j*sr,sr*sizeof(PetscScalar));
1090: }
1091: }
1092: return(0);
1093: }
1097: PetscErrorCode NEPSolve_NLEIGS(NEP nep)
1098: {
1100: NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;
1101: PetscInt i,j,k=0,l,nv=0,ld,lds,off,ldds,rs1,nq=0,newn;
1102: PetscInt lwa,lrwa,nwu=0,nrwu=0,deg=ctx->nmat-1,nconv=0;
1103: PetscScalar *S,*Q,*work,*H,*pU,*K,betak=0,*Hc,*eigr,*eigi;
1104: PetscReal betah,norm,*rwork;
1105: PetscBool breakdown=PETSC_FALSE,lindep;
1106: Mat U;
1109: ld = nep->ncv+deg;
1110: lds = deg*ld;
1111: lwa = (deg+6)*ld*lds;
1112: lrwa = 7*lds;
1113: DSGetLeadingDimension(nep->ds,&ldds);
1114: PetscMalloc4(lwa,&work,lrwa,&rwork,lds*ld,&S,ldds*ldds,&Hc);
1115: PetscMemzero(S,lds*ld*sizeof(PetscScalar));
1116: if (!ctx->nshifts) {
1117: PetscMalloc2(nep->ncv,&eigr,nep->ncv,&eigi);
1118: } else { eigr = nep->eigr; eigi = nep->eigi; }
1119: BVDuplicateResize(nep->V,PetscMax(nep->nt-1,ctx->nmat-1),&ctx->W);
1121: /* Get the starting vector */
1122: for (i=0;i<deg;i++) {
1123: BVSetRandomColumn(nep->V,i);
1124: BVOrthogonalizeColumn(nep->V,i,S+i*ld,&norm,&lindep);
1125: if (!lindep) {
1126: BVScaleColumn(nep->V,i,1/norm);
1127: S[i+i*ld] = norm;
1128: nq++;
1129: }
1130: }
1131: if (!nq) SETERRQ(PetscObjectComm((PetscObject)nep),1,"NEP: Problem with initial vector");
1132: NEPTOARSNorm2(lds,S,&norm);
1133: for (j=0;j<deg;j++) {
1134: for (i=0;i<=j;i++) S[i+j*ld] /= norm;
1135: }
1137: /* Restart loop */
1138: l = 0;
1139: while (nep->reason == NEP_CONVERGED_ITERATING) {
1140: nep->its++;
1141:
1142: /* Compute an nv-step Krylov relation */
1143: nv = PetscMin(nep->nconv+nep->mpd,nep->ncv);
1144: if (ctx->nshifts) { DSGetArray(nep->ds,DS_MAT_A,&K); }
1145: DSGetArray(nep->ds,ctx->nshifts?DS_MAT_B:DS_MAT_A,&H);
1146: NEPNLEIGSTOARrun(nep,&nq,S,ld,K,H,ldds,nep->V,nep->nconv+l,&nv,&breakdown,nep->work);
1147: betah = PetscAbsScalar(H[(nv-1)*ldds+nv]);
1148: DSRestoreArray(nep->ds,ctx->nshifts?DS_MAT_B:DS_MAT_A,&H);
1149: if (ctx->nshifts) {
1150: betak = K[(nv-1)*ldds+nv];
1151: DSRestoreArray(nep->ds,DS_MAT_A,&K);
1152: }
1153: DSSetDimensions(nep->ds,nv,0,nep->nconv,nep->nconv+l);
1154: if (l==0) {
1155: DSSetState(nep->ds,DS_STATE_INTERMEDIATE);
1156: } else {
1157: DSSetState(nep->ds,DS_STATE_RAW);
1158: }
1160: /* Solve projected problem */
1161: if (ctx->nshifts) {
1162: DSGetArray(nep->ds,DS_MAT_B,&H);
1163: PetscMemcpy(Hc,H,ldds*ldds*sizeof(PetscScalar));
1164: DSRestoreArray(nep->ds,DS_MAT_B,&H);
1165: }
1166: DSSolve(nep->ds,nep->eigr,nep->eigi);
1167: DSSort(nep->ds,nep->eigr,nep->eigi,NULL,NULL,NULL);
1168: if (!ctx->nshifts) {
1169: DSUpdateExtraRow(nep->ds);
1170: }
1172: /* Check convergence */
1173: NEPNLEIGSKrylovConvergence(nep,S,ld,nq,Hc,PETSC_FALSE,nep->nconv,nv-nep->nconv,betak,betah,&k,nep->work);
1174: (*nep->stopping)(nep,nep->its,nep->max_it,k,nep->nev,&nep->reason,nep->stoppingctx);
1175: nconv = k;
1177: /* Update l */
1178: if (nep->reason != NEP_CONVERGED_ITERATING || breakdown) l = 0;
1179: else {
1180: l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
1181: if (!breakdown) {
1182: /* Prepare the Rayleigh quotient for restart */
1183: if (!ctx->nshifts) {
1184: DSTruncate(nep->ds,k+l);
1185: DSGetDimensions(nep->ds,&newn,NULL,NULL,NULL,NULL);
1186: l = newn-k;
1187: } else {
1188: DSGetArray(nep->ds,DS_MAT_Q,&Q);
1189: DSGetArray(nep->ds,DS_MAT_B,&H);
1190: DSGetArray(nep->ds,DS_MAT_A,&K);
1191: for (i=ctx->lock?k:0;i<k+l;i++) {
1192: H[k+l+i*ldds] = betah*Q[nv-1+i*ldds];
1193: K[k+l+i*ldds] = betak*Q[nv-1+i*ldds];
1194: }
1195: DSRestoreArray(nep->ds,DS_MAT_B,&H);
1196: DSRestoreArray(nep->ds,DS_MAT_A,&K);
1197: DSRestoreArray(nep->ds,DS_MAT_Q,&Q);
1198: DSSetDimensions(nep->ds,k+l,0,nep->nconv,0);
1199: }
1200: }
1201: }
1202: if (!ctx->lock && l>0) { l += k; k = 0; }
1204: /* Update S */
1205: off = nep->nconv*ldds;
1206: DSGetArray(nep->ds,ctx->nshifts?DS_MAT_Z:DS_MAT_Q,&Q);
1207: NEPTOARSupdate(S,ld,deg,nq,nep->nconv,k+l-nep->nconv,nv,Q+off,ldds,work+nwu);
1208: DSRestoreArray(nep->ds,ctx->nshifts?DS_MAT_Z:DS_MAT_Q,&Q);
1210: /* Copy last column of S */
1211: PetscMemcpy(S+lds*(k+l),S+lds*nv,lds*sizeof(PetscScalar));
1212: if (nep->reason == NEP_CONVERGED_ITERATING) {
1213: if (breakdown) {
1215: /* Stop if breakdown */
1216: PetscInfo2(nep,"Breakdown (it=%D norm=%g)\n",nep->its,(double)betah);
1217: nep->reason = NEP_DIVERGED_BREAKDOWN;
1218: } else {
1219: /* Truncate S */
1220: NEPTOARTrunc(nep,S,ld,deg,&nq,k+l+1,work+nwu,rwork+nrwu);
1221: }
1222: }
1223: nep->nconv = k;
1224: if (!ctx->nshifts) {
1225: for (i=0;i<nv;i++) { eigr[i] = nep->eigr[i]; eigi[i] = nep->eigi[i]; }
1226: NEPNLEIGSBackTransform((PetscObject)nep,nv,eigr,eigi);
1227: }
1228: NEPMonitor(nep,nep->its,nconv,eigr,eigi,nep->errest,nv);
1229: }
1230: nep->nconv = nconv;
1231: if (nep->nconv>0) {
1232: /* Extract invariant pair */
1233: NEPTOARTrunc(nep,S,ld,deg,&nq,nep->nconv,work+nwu,rwork+nrwu);
1234: /* Update vectors V = V*S or V=V*S*H */
1235: rs1 = nep->nconv;
1236: if (ctx->nshifts) {
1237: DSGetArray(nep->ds,DS_MAT_B,&H);
1238: NEPTOARSupdate(S,ld,deg,rs1,0,nep->nconv,nep->nconv,H,ldds,work+nwu);
1239: DSRestoreArray(nep->ds,DS_MAT_B,&H);
1240: }
1241: PetscMalloc1(rs1*nep->nconv,&pU);
1242: for (i=0;i<nep->nconv;i++) {
1243: PetscMemcpy(pU+i*rs1,S+i*lds,rs1*sizeof(PetscScalar));
1244: }
1245: MatCreateSeqDense(PETSC_COMM_SELF,rs1,nep->nconv,pU,&U);
1246: BVSetActiveColumns(nep->V,0,rs1);
1247: BVMultInPlace(nep->V,U,0,nep->nconv);
1248: BVSetActiveColumns(nep->V,0,nep->nconv);
1249: MatDestroy(&U);
1250: PetscFree(pU);
1251: }
1252: /* truncate Schur decomposition and change the state to raw so that
1253: DSVectors() computes eigenvectors from scratch */
1254: DSSetDimensions(nep->ds,nep->nconv,0,0,0);
1255: DSSetState(nep->ds,DS_STATE_RAW);
1257: PetscFree4(work,rwork,S,Hc);
1258: /* Map eigenvalues back to the original problem */
1259: if (!ctx->nshifts) {
1260: NEPNLEIGSBackTransform((PetscObject)nep,nep->nconv,nep->eigr,nep->eigi);
1261: PetscFree2(eigr,eigi);
1262: }
1263: BVDestroy(&ctx->W);
1264: return(0);
1265: }
1269: static PetscErrorCode NEPNLEIGSSetSingularitiesFunction_NLEIGS(NEP nep,PetscErrorCode (*fun)(NEP,PetscInt*,PetscScalar*,void*),void *ctx)
1270: {
1271: NEP_NLEIGS *nepctx=(NEP_NLEIGS*)nep->data;
1274: if (fun) nepctx->computesingularities = fun;
1275: if (ctx) nepctx->singularitiesctx = ctx;
1276: return(0);
1277: }
1281: /*@C
1282: NEPNLEIGSSetSingularitiesFunction - Sets a user function to compute a discretization
1283: of the singularity set (where T(.) is not analytic).
1285: Logically Collective on NEP
1287: Input Parameters:
1288: + nep - the NEP context
1289: . fun - user function (if NULL then NEP retains any previously set value)
1290: - ctx - [optional] user-defined context for private data for the function
1291: (may be NULL, in which case NEP retains any previously set value)
1293: Calling Sequence of fun:
1294: $ fun(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *ctx)
1296: + nep - the NEP context
1297: . maxnp - on input number of requested points in the discretization (can be set)
1298: . xi - computed values of the discretization
1299: - ctx - optional context, as set by NEPNLEIGSSetSingularitiesFunction()
1301: Note:
1302: The user-defined function can set a smaller value of maxnp if necessary.
1304: Level: intermediate
1306: .seealso: NEPNLEIGSGetSingularitiesFunction()
1307: @*/
1308: PetscErrorCode NEPNLEIGSSetSingularitiesFunction(NEP nep,PetscErrorCode (*fun)(NEP,PetscInt*,PetscScalar*,void*),void *ctx)
1309: {
1314: PetscTryMethod(nep,"NEPNLEIGSSetSingularitiesFunction_C",(NEP,PetscErrorCode(*)(NEP,PetscInt*,PetscScalar*,void*),void*),(nep,fun,ctx));
1315: return(0);
1316: }
1320: static PetscErrorCode NEPNLEIGSGetSingularitiesFunction_NLEIGS(NEP nep,PetscErrorCode (**fun)(NEP,PetscInt*,PetscScalar*,void*),void **ctx)
1321: {
1322: NEP_NLEIGS *nepctx=(NEP_NLEIGS*)nep->data;
1325: if (fun) *fun = nepctx->computesingularities;
1326: if (ctx) *ctx = nepctx->singularitiesctx;
1327: return(0);
1328: }
1332: /*@C
1333: NEPNLEIGSGetSingularitiesFunction - Returns the Function and optionally the user
1334: provided context for computing a discretization of the singularity set.
1336: Not Collective
1338: Input Parameter:
1339: . nep - the nonlinear eigensolver context
1341: Output Parameters:
1342: + fun - location to put the function (or NULL)
1343: - ctx - location to stash the function context (or NULL)
1345: Level: advanced
1347: .seealso: NEPNLEIGSSetSingularitiesFunction()
1348: @*/
1349: PetscErrorCode NEPNLEIGSGetSingularitiesFunction(NEP nep,PetscErrorCode (**fun)(NEP,PetscInt*,PetscScalar*,void*),void **ctx)
1350: {
1355: PetscUseMethod(nep,"NEPNLEIGSGetSingularitiesFunction_C",(NEP,PetscErrorCode(**)(NEP,PetscInt*,PetscScalar*,void*),void**),(nep,fun,ctx));
1356: return(0);
1357: }
1361: static PetscErrorCode NEPNLEIGSSetRestart_NLEIGS(NEP nep,PetscReal keep)
1362: {
1363: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1366: if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
1367: else {
1368: 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]");
1369: ctx->keep = keep;
1370: }
1371: return(0);
1372: }
1376: /*@
1377: NEPNLEIGSSetRestart - Sets the restart parameter for the NLEIGS
1378: method, in particular the proportion of basis vectors that must be kept
1379: after restart.
1381: Logically Collective on NEP
1383: Input Parameters:
1384: + nep - the nonlinear eigensolver context
1385: - keep - the number of vectors to be kept at restart
1387: Options Database Key:
1388: . -nep_nleigs_restart - Sets the restart parameter
1390: Notes:
1391: Allowed values are in the range [0.1,0.9]. The default is 0.5.
1393: Level: advanced
1395: .seealso: NEPNLEIGSGetRestart()
1396: @*/
1397: PetscErrorCode NEPNLEIGSSetRestart(NEP nep,PetscReal keep)
1398: {
1404: PetscTryMethod(nep,"NEPNLEIGSSetRestart_C",(NEP,PetscReal),(nep,keep));
1405: return(0);
1406: }
1410: static PetscErrorCode NEPNLEIGSGetRestart_NLEIGS(NEP nep,PetscReal *keep)
1411: {
1412: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1415: *keep = ctx->keep;
1416: return(0);
1417: }
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: {
1443: PetscUseMethod(nep,"NEPNLEIGSGetRestart_C",(NEP,PetscReal*),(nep,keep));
1444: return(0);
1445: }
1449: static PetscErrorCode NEPNLEIGSSetLocking_NLEIGS(NEP nep,PetscBool lock)
1450: {
1451: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1454: ctx->lock = lock;
1455: return(0);
1456: }
1460: /*@
1461: NEPNLEIGSSetLocking - Choose between locking and non-locking variants of
1462: the NLEIGS method.
1464: Logically Collective on NEP
1466: Input Parameters:
1467: + nep - the nonlinear eigensolver context
1468: - lock - true if the locking variant must be selected
1470: Options Database Key:
1471: . -nep_nleigs_locking - Sets the locking flag
1473: Notes:
1474: The default is to lock converged eigenpairs when the method restarts.
1475: This behaviour can be changed so that all directions are kept in the
1476: working subspace even if already converged to working accuracy (the
1477: non-locking variant).
1479: Level: advanced
1481: .seealso: NEPNLEIGSGetLocking()
1482: @*/
1483: PetscErrorCode NEPNLEIGSSetLocking(NEP nep,PetscBool lock)
1484: {
1490: PetscTryMethod(nep,"NEPNLEIGSSetLocking_C",(NEP,PetscBool),(nep,lock));
1491: return(0);
1492: }
1496: static PetscErrorCode NEPNLEIGSGetLocking_NLEIGS(NEP nep,PetscBool *lock)
1497: {
1498: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1501: *lock = ctx->lock;
1502: return(0);
1503: }
1507: /*@
1508: NEPNLEIGSGetLocking - Gets the locking flag used in the NLEIGS method.
1510: Not Collective
1512: Input Parameter:
1513: . nep - the nonlinear eigensolver context
1515: Output Parameter:
1516: . lock - the locking flag
1518: Level: advanced
1520: .seealso: NEPNLEIGSSetLocking()
1521: @*/
1522: PetscErrorCode NEPNLEIGSGetLocking(NEP nep,PetscBool *lock)
1523: {
1529: PetscUseMethod(nep,"NEPNLEIGSGetLocking_C",(NEP,PetscBool*),(nep,lock));
1530: return(0);
1531: }
1535: static PetscErrorCode NEPNLEIGSSetInterpolation_NLEIGS(NEP nep,PetscReal tol,PetscInt maxits)
1536: {
1537: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1540: if (tol == PETSC_DEFAULT) {
1541: ctx->ddtol = PETSC_DEFAULT;
1542: nep->state = NEP_STATE_INITIAL;
1543: } else {
1544: if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1545: ctx->ddtol = tol;
1546: }
1547: if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
1548: ctx->ddmaxit = 0;
1549: nep->state = NEP_STATE_INITIAL;
1550: } else {
1551: if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
1552: ctx->ddmaxit = maxits;
1553: }
1554: return(0);
1555: }
1559: /*@
1560: NEPNLEIGSSetInterpolation - Sets the tolerance and maximum iteration count used
1561: by the NLEIGS method when building the interpolation via divided differences.
1563: Logically Collective on NEP
1565: Input Parameters:
1566: + nep - the nonlinear eigensolver context
1567: . tol - the convergence tolerance
1568: - maxits - maximum number of iterations to use
1570: Options Database Key:
1571: + -nep_nleigs_interpolation_tol <tol> - Sets the convergence tolerance
1572: - -nep_nleigs_interpolation_max_it <maxits> - Sets the maximum number of iterations
1574: Notes:
1575: Use PETSC_DEFAULT for either argument to assign a reasonably good value.
1577: Level: advanced
1579: .seealso: NEPNLEIGSGetInterpolation()
1580: @*/
1581: PetscErrorCode NEPNLEIGSSetInterpolation(NEP nep,PetscReal tol,PetscInt maxits)
1582: {
1589: PetscTryMethod(nep,"NEPNLEIGSSetInterpolation_C",(NEP,PetscReal,PetscInt),(nep,tol,maxits));
1590: return(0);
1591: }
1595: static PetscErrorCode NEPNLEIGSGetInterpolation_NLEIGS(NEP nep,PetscReal *tol,PetscInt *maxits)
1596: {
1597: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1600: if (tol) *tol = ctx->ddtol;
1601: if (maxits) *maxits = ctx->ddmaxit;
1602: return(0);
1603: }
1607: /*@
1608: NEPNLEIGSGetInterpolation - Gets the tolerance and maximum iteration count used
1609: by the NLEIGS method when building the interpolation via divided differences.
1611: Not Collective
1613: Input Parameter:
1614: . nep - the nonlinear eigensolver context
1616: Output Parameter:
1617: + tol - the convergence tolerance
1618: - maxits - maximum number of iterations
1620: Level: advanced
1622: .seealso: NEPNLEIGSSetInterpolation()
1623: @*/
1624: PetscErrorCode NEPNLEIGSGetInterpolation(NEP nep,PetscReal *tol,PetscInt *maxits)
1625: {
1630: PetscTryMethod(nep,"NEPNLEIGSGetInterpolation_C",(NEP,PetscReal*,PetscInt*),(nep,tol,maxits));
1631: return(0);
1632: }
1636: static PetscErrorCode NEPNLEIGSSetTrueResidual_NLEIGS(NEP nep,PetscBool trueres)
1637: {
1638: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1641: ctx->trueres = trueres;
1642: return(0);
1643: }
1647: /*@
1648: NEPNLEIGSSetTrueResidual - Specifies if the solver must compute the true residual
1649: explicitly or not.
1651: Logically Collective on NEP
1653: Input Parameters:
1654: + nep - the nonlinear eigensolver context
1655: - trueres - whether true residuals are required or not
1657: Options Database Key:
1658: . -nep_nleigs_true_residual <boolean> - Sets/resets the boolean flag 'trueres'
1660: Notes:
1661: If the user sets trueres=PETSC_TRUE then the solver explicitly computes
1662: the true residual norm for each eigenpair approximation, and uses it for
1663: convergence testing. The default is to use the cheaper approximation
1664: available from the (rational) Krylov iteration.
1666: Level: advanced
1668: .seealso: NEPNLEIGSGetTrueResidual()
1669: @*/
1670: PetscErrorCode NEPNLEIGSSetTrueResidual(NEP nep,PetscBool trueres)
1671: {
1677: PetscTryMethod(nep,"NEPNLEIGSSetTrueResidual_C",(NEP,PetscBool),(nep,trueres));
1678: return(0);
1679: }
1683: static PetscErrorCode NEPNLEIGSGetTrueResidual_NLEIGS(NEP nep,PetscBool *trueres)
1684: {
1685: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1688: *trueres = ctx->trueres;
1689: return(0);
1690: }
1694: /*@
1695: NEPNLEIGSGetTrueResidual - Returns the flag indicating whether true
1696: residuals must be computed explicitly or not.
1698: Not Collective
1700: Input Parameter:
1701: . nep - the nonlinear eigensolver context
1703: Output Parameter:
1704: . trueres - the returned flag
1706: Level: advanced
1708: .seealso: NEPNLEIGSSetTrueResidual()
1709: @*/
1710: PetscErrorCode NEPNLEIGSGetTrueResidual(NEP nep,PetscBool *trueres)
1711: {
1717: PetscTryMethod(nep,"NEPNLEIGSGetTrueResidual_C",(NEP,PetscBool*),(nep,trueres));
1718: return(0);
1719: }
1723: static PetscErrorCode NEPNLEIGSSetRKShifts_NLEIGS(NEP nep,PetscInt ns,PetscScalar *shifts)
1724: {
1726: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1727: PetscInt i;
1730: if (ns<=0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"Number of shifts must be positive");
1731: if (ctx->nshifts) { PetscFree(ctx->shifts); }
1732: for (i=0;i<ctx->nshiftsw;i++) { KSPDestroy(&ctx->ksp[i]); }
1733: PetscFree(ctx->ksp);
1734: ctx->ksp = NULL;
1735: PetscMalloc(ns,&ctx->shifts);
1736: for (i=0;i<ns;i++) ctx->shifts[i] = shifts[i];
1737: ctx->nshifts = ns;
1738: nep->state = NEP_STATE_INITIAL;
1739: return(0);
1740: }
1744: /*@C
1745: NEPNLEIGSSetRKShifts - Sets a list of shifts to be used in the Rational
1746: Krylov method.
1748: Logically Collective on NEP
1750: Input Parameters:
1751: + nep - the nonlinear eigensolver context
1752: . ns - number of shifts
1753: - shifts - array of scalar values specifying the shifts
1755: Options Database Key:
1756: . -nep_nleigs_rk_shifts - Sets the list of shifts
1758: Notes:
1759: If only one shift is provided, the subspace is built with the simpler
1760: shift-and-invert Krylov-Schur.
1762: In the case of real scalars, complex shifts are not allowed. In the
1763: command line, a comma-separated list of complex values can be provided with
1764: the format [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
1765: -nep_nleigs_rk_shifts 1.0+2.0i,1.5+2.0i,1.0+1.5i
1767: Level: advanced
1769: .seealso: NEPNLEIGSGetRKShifts()
1770: @*/
1771: PetscErrorCode NEPNLEIGSSetRKShifts(NEP nep,PetscInt ns,PetscScalar *shifts)
1772: {
1779: PetscTryMethod(nep,"NEPNLEIGSSetRKShifts_C",(NEP,PetscInt,PetscScalar*),(nep,ns,shifts));
1780: return(0);
1781: }
1785: static PetscErrorCode NEPNLEIGSGetRKShifts_NLEIGS(NEP nep,PetscInt *ns,PetscScalar **shifts)
1786: {
1788: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1789: PetscInt i;
1792: *ns = ctx->nshifts;
1793: if (ctx->nshifts) {
1794: PetscMalloc1(ctx->nshifts,shifts);
1795: for (i=0;i<ctx->nshifts;i++) (*shifts)[i] = ctx->shifts[i];
1796: }
1797: return(0);
1798: }
1802: /*@C
1803: NEPNLEIGSGetRKShifts - Gets the list of shifts used in the Rational
1804: Krylov method.
1806: Not Collective
1808: Input Parameter:
1809: . nep - the nonlinear eigensolver context
1811: Output Parameter:
1812: + ns - number of shifts
1813: - shifts - array of shifts
1815: Level: advanced
1817: .seealso: NEPNLEIGSSetRKShifts()
1818: @*/
1819: PetscErrorCode NEPNLEIGSGetRKShifts(NEP nep,PetscInt *ns,PetscScalar **shifts)
1820: {
1827: PetscTryMethod(nep,"NEPNLEIGSGetRKShifts_C",(NEP,PetscInt*,PetscScalar**),(nep,ns,shifts));
1828: return(0);
1829: }
1831: #define SHIFTMAX 30
1835: PetscErrorCode NEPSetFromOptions_NLEIGS(PetscOptionItems *PetscOptionsObject,NEP nep)
1836: {
1838: NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;
1839: PetscInt i,k;
1840: PetscBool flg1,flg2,b;
1841: PetscReal r;
1842: PetscScalar array[SHIFTMAX];
1843: PC pc;
1844: PCType pctype;
1845: KSPType ksptype;
1848: PetscOptionsHead(PetscOptionsObject,"NEP NLEIGS Options");
1849: PetscOptionsReal("-nep_nleigs_restart","Proportion of vectors kept after restart","NEPNLEIGSSetRestart",0.5,&r,&flg1);
1850: if (flg1) {
1851: NEPNLEIGSSetRestart(nep,r);
1852: }
1853: PetscOptionsBool("-nep_nleigs_locking","Choose between locking and non-locking variants","NEPNLEIGSSetLocking",PETSC_FALSE,&b,&flg1);
1854: if (flg1) {
1855: NEPNLEIGSSetLocking(nep,b);
1856: }
1857: PetscOptionsBool("-nep_nleigs_true_residual","Compute true residuals explicitly","NEPNLEIGSSetTrueResidual",PETSC_FALSE,&b,&flg1);
1858: if (flg1) {
1859: NEPNLEIGSSetTrueResidual(nep,b);
1860: }
1861: NEPNLEIGSGetInterpolation(nep,&r,&i);
1862: if (!i) i = PETSC_DEFAULT;
1863: PetscOptionsInt("-nep_nleigs_interpolation_max_it","Maximum number of terms for interpolation via divided differences","NEPNLEIGSSetInterpolation",i,&i,&flg1);
1864: PetscOptionsReal("-nep_nleigs_interpolation_tol","Tolerance for interpolation via divided differences","NEPNLEIGSSetInterpolation",r,&r,&flg2);
1865: if (flg1 || flg2) {
1866: NEPNLEIGSSetInterpolation(nep,r,i);
1867: }
1868: k = SHIFTMAX;
1869: for (i=0;i<k;i++) array[i] = 0;
1870: PetscOptionsScalarArray("-nep_nleigs_rk_shifts","Shifts for Rational Krylov","NEPNLEIGSSetRKShifts",array,&k,&flg1);
1871: if (flg1) {
1872: NEPNLEIGSSetRKShifts(nep,k,array);
1873: }
1875: if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
1876: for (i=0;i<ctx->nshiftsw;i++) {
1877: KSPGetPC(ctx->ksp[i],&pc);
1878: KSPGetType(ctx->ksp[i],&ksptype);
1879: PCGetType(pc,&pctype);
1880: if (!pctype && !ksptype) {
1881: KSPSetType(ctx->ksp[i],KSPPREONLY);
1882: PCSetType(pc,PCLU);
1883: }
1884: KSPSetOperators(ctx->ksp[i],nep->function,nep->function_pre);
1885: KSPSetFromOptions(ctx->ksp[i]);
1886: }
1887: PetscOptionsTail();
1888: return(0);
1889: }
1893: static PetscErrorCode NEPNLEIGSGetKSPs_NLEIGS(NEP nep,KSP **ksp)
1894: {
1896: NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;
1897: PetscInt i;
1900: if (!ctx->ksp) {
1901: NEPNLEIGSSetShifts(nep);
1902: PetscMalloc1(ctx->nshiftsw,&ctx->ksp);
1903: for (i=0;i<ctx->nshiftsw;i++) {
1904: KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp[i]);
1905: KSPSetOptionsPrefix(ctx->ksp[i],((PetscObject)nep)->prefix);
1906: KSPAppendOptionsPrefix(ctx->ksp[i],"nep_nleigs_");
1907: PetscObjectIncrementTabLevel((PetscObject)ctx->ksp[i],(PetscObject)nep,1);
1908: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp[i]);
1909: KSPSetErrorIfNotConverged(ctx->ksp[i],PETSC_TRUE);
1910: }
1911: }
1912: *ksp = ctx->ksp;
1913: return(0);
1914: }
1918: /*@C
1919: NEPNLEIGSGetKSPs - Retrieve the array of linear solver objects associated with
1920: the nonlinear eigenvalue solver.
1922: Not Collective
1924: Input Parameter:
1925: . nep - nonlinear eigenvalue solver
1927: Output Parameter:
1928: . ksp - array of linear solver object
1930: Level: advanced
1931: @*/
1932: PetscErrorCode NEPNLEIGSGetKSPs(NEP nep,KSP **ksp)
1933: {
1939: PetscUseMethod(nep,"NEPNLEIGSGetKSPs_C",(NEP,KSP**),(nep,ksp));
1940: return(0);
1941: }
1945: PetscErrorCode NEPView_NLEIGS(NEP nep,PetscViewer viewer)
1946: {
1948: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1949: PetscBool isascii;
1950: PetscInt i;
1951: char str[50];
1954: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1955: if (isascii) {
1956: PetscViewerASCIIPrintf(viewer," NLEIGS: %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
1957: PetscViewerASCIIPrintf(viewer," NLEIGS: using the %slocking variant\n",ctx->lock?"":"non-");
1958: PetscViewerASCIIPrintf(viewer," NLEIGS: divided difference terms: used=%D, max=%D\n",ctx->nmat-1,ctx->ddmaxit);
1959: PetscViewerASCIIPrintf(viewer," NLEIGS: tolerance for divided difference convergence: %g\n",(double)ctx->ddtol);
1960: if (ctx->nshifts) {
1961: PetscViewerASCIIPrintf(viewer," NLEIGS: RK shifts: ");
1962: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1963: for (i=0;i<ctx->nshifts;i++) {
1964: SlepcSNPrintfScalar(str,50,ctx->shifts[i],PETSC_FALSE);
1965: PetscViewerASCIIPrintf(viewer,"%s%s",str,(i<ctx->nshifts-1)?",":"");
1966: }
1967: PetscViewerASCIIPrintf(viewer,"\n");
1968: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1969: }
1970: if (ctx->trueres) { PetscViewerASCIIPrintf(viewer," NLEIGS: computing true residuals for convergence check\n"); }
1971: if (!ctx->ksp) { NEPNLEIGSGetKSPs(nep,&ctx->ksp); }
1972: PetscViewerASCIIPushTab(viewer);
1973: KSPView(ctx->ksp[0],viewer);
1974: PetscViewerASCIIPopTab(viewer);
1975: }
1976: return(0);
1977: }
1981: PetscErrorCode NEPReset_NLEIGS(NEP nep)
1982: {
1984: PetscInt k;
1985: NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data;
1988: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
1989: PetscFree(ctx->coeffD);
1990: } else {
1991: for (k=0;k<ctx->nmat;k++) { MatDestroy(&ctx->D[k]); }
1992: }
1993: if (ctx->vrn) {
1994: VecDestroy(&ctx->vrn);
1995: }
1996: return(0);
1997: }
2001: PetscErrorCode NEPDestroy_NLEIGS(NEP nep)
2002: {
2004: PetscInt k;
2005: NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;
2008: for (k=0;k<ctx->nshiftsw;k++) { KSPDestroy(&ctx->ksp[k]); }
2009: PetscFree(ctx->ksp);
2010: if (ctx->nshifts) { PetscFree(ctx->shifts); }
2011: PetscFree4(ctx->s,ctx->xi,ctx->beta,ctx->D);
2012: PetscFree(nep->data);
2013: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetSingularitiesFunction_C",NULL);
2014: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetSingularitiesFunction_C",NULL);
2015: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRestart_C",NULL);
2016: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRestart_C",NULL);
2017: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetLocking_C",NULL);
2018: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetLocking_C",NULL);
2019: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetInterpolation_C",NULL);
2020: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetInterpolation_C",NULL);
2021: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetTrueResidual_C",NULL);
2022: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetTrueResidual_C",NULL);
2023: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRKShifts_C",NULL);
2024: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRKShifts_C",NULL);
2025: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetKSPs_C",NULL);
2026: return(0);
2027: }
2031: PETSC_EXTERN PetscErrorCode NEPCreate_NLEIGS(NEP nep)
2032: {
2034: NEP_NLEIGS *ctx;
2037: PetscNewLog(nep,&ctx);
2038: nep->data = (void*)ctx;
2039: ctx->lock = PETSC_TRUE;
2040: ctx->ddtol = PETSC_DEFAULT;
2041: ctx->ddmaxit = 0;
2042: ctx->trueres = PETSC_FALSE;
2043: ctx->nshifts = 0;
2045: nep->ops->solve = NEPSolve_NLEIGS;
2046: nep->ops->setup = NEPSetUp_NLEIGS;
2047: nep->ops->setfromoptions = NEPSetFromOptions_NLEIGS;
2048: nep->ops->view = NEPView_NLEIGS;
2049: nep->ops->destroy = NEPDestroy_NLEIGS;
2050: nep->ops->reset = NEPReset_NLEIGS;
2051: nep->ops->computevectors = NEPComputeVectors_Schur;
2052: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetSingularitiesFunction_C",NEPNLEIGSSetSingularitiesFunction_NLEIGS);
2053: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetSingularitiesFunction_C",NEPNLEIGSGetSingularitiesFunction_NLEIGS);
2054: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRestart_C",NEPNLEIGSSetRestart_NLEIGS);
2055: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRestart_C",NEPNLEIGSGetRestart_NLEIGS);
2056: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetLocking_C",NEPNLEIGSSetLocking_NLEIGS);
2057: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetLocking_C",NEPNLEIGSGetLocking_NLEIGS);
2058: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetInterpolation_C",NEPNLEIGSSetInterpolation_NLEIGS);
2059: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetInterpolation_C",NEPNLEIGSGetInterpolation_NLEIGS);
2060: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetTrueResidual_C",NEPNLEIGSSetTrueResidual_NLEIGS);
2061: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetTrueResidual_C",NEPNLEIGSGetTrueResidual_NLEIGS);
2062: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSSetRKShifts_C",NEPNLEIGSSetRKShifts_NLEIGS);
2063: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetRKShifts_C",NEPNLEIGSGetRKShifts_NLEIGS);
2064: PetscObjectComposeFunction((PetscObject)nep,"NEPNLEIGSGetKSPs_C",NEPNLEIGSGetKSPs_NLEIGS);
2065: return(0);
2066: }