Actual source code: ptoar.c
slepc-3.15.0 2021-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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 polynomial eigensolver: "toar"
13: Method: TOAR
15: Algorithm:
17: Two-Level Orthogonal Arnoldi.
19: References:
21: [1] Y. Su, J. Zhang and Z. Bai, "A compact Arnoldi algorithm for
22: polynomial eigenvalue problems", talk presented at RANMEP 2008.
24: [2] C. Campos and J.E. Roman, "Parallel Krylov solvers for the
25: polynomial eigenvalue problem in SLEPc", SIAM J. Sci. Comput.
26: 38(5):S385-S411, 2016.
28: [3] D. Lu, Y. Su and Z. Bai, "Stability analysis of the two-level
29: orthogonal Arnoldi procedure", SIAM J. Matrix Anal. App.
30: 37(1):195-214, 2016.
31: */
33: #include <slepc/private/pepimpl.h>
34: #include "../src/pep/impls/krylov/pepkrylov.h"
35: #include <slepcblaslapack.h>
37: static PetscBool cited = PETSC_FALSE;
38: static const char citation[] =
39: "@Article{slepc-pep,\n"
40: " author = \"C. Campos and J. E. Roman\",\n"
41: " title = \"Parallel {Krylov} solvers for the polynomial eigenvalue problem in {SLEPc}\",\n"
42: " journal = \"{SIAM} J. Sci. Comput.\",\n"
43: " volume = \"38\",\n"
44: " number = \"5\",\n"
45: " pages = \"S385--S411\",\n"
46: " year = \"2016,\"\n"
47: " doi = \"https://doi.org/10.1137/15M1022458\"\n"
48: "}\n";
50: PetscErrorCode PEPSetUp_TOAR(PEP pep)
51: {
53: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
54: PetscBool sinv,flg;
55: PetscInt i;
58: PEPCheckShiftSinvert(pep);
59: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
60: if (!ctx->lock && pep->mpd<pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
61: if (pep->max_it==PETSC_DEFAULT) pep->max_it = PetscMax(100,2*(pep->nmat-1)*pep->n/pep->ncv);
62: if (!pep->which) { PEPSetWhichEigenpairs_Default(pep); }
63: if (pep->which==PEP_ALL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues");
64: if (pep->problem_type!=PEP_GENERAL) {
65: PetscInfo(pep,"Problem type ignored, performing a non-symmetric linearization\n");
66: }
68: if (!ctx->keep) ctx->keep = 0.5;
70: PEPAllocateSolution(pep,pep->nmat-1);
71: PEPSetWorkVecs(pep,3);
72: DSSetType(pep->ds,DSNHEP);
73: DSSetExtraRow(pep->ds,PETSC_TRUE);
74: DSAllocate(pep->ds,pep->ncv+1);
76: PEPBasisCoefficients(pep,pep->pbc);
77: STGetTransform(pep->st,&flg);
78: if (!flg) {
79: PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
80: PetscLogObjectMemory((PetscObject)pep,pep->nmat*sizeof(PetscScalar));
81: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
82: if (sinv) {
83: PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL);
84: } else {
85: for (i=0;i<pep->nmat-1;i++) pep->solvematcoeffs[i] = 0.0;
86: pep->solvematcoeffs[pep->nmat-1] = 1.0;
87: }
88: }
89: BVDestroy(&ctx->V);
90: BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
91: return(0);
92: }
94: /*
95: Extend the TOAR basis by applying the the matrix operator
96: over a vector which is decomposed in the TOAR way
97: Input:
98: - pbc: array containing the polynomial basis coefficients
99: - S,V: define the latest Arnoldi vector (nv vectors in V)
100: Output:
101: - t: new vector extending the TOAR basis
102: - r: temporary coefficients to compute the TOAR coefficients
103: for the new Arnoldi vector
104: Workspace: t_ (two vectors)
105: */
106: static PetscErrorCode PEPTOARExtendBasis(PEP pep,PetscBool sinvert,PetscScalar sigma,PetscScalar *S,PetscInt ls,PetscInt nv,BV V,Vec t,PetscScalar *r,PetscInt lr,Vec *t_)
107: {
109: PetscInt nmat=pep->nmat,deg=nmat-1,k,j,off=0,lss;
110: Vec v=t_[0],ve=t_[1],q=t_[2];
111: PetscScalar alpha=1.0,*ss,a;
112: PetscReal *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat;
113: PetscBool flg;
116: BVSetActiveColumns(pep->V,0,nv);
117: STGetTransform(pep->st,&flg);
118: if (sinvert) {
119: for (j=0;j<nv;j++) {
120: if (deg>1) r[lr+j] = S[j]/ca[0];
121: if (deg>2) r[2*lr+j] = (S[ls+j]+(sigma-cb[1])*r[lr+j])/ca[1];
122: }
123: for (k=2;k<deg-1;k++) {
124: for (j=0;j<nv;j++) r[(k+1)*lr+j] = (S[k*ls+j]+(sigma-cb[k])*r[k*lr+j]-cg[k]*r[(k-1)*lr+j])/ca[k];
125: }
126: k = deg-1;
127: for (j=0;j<nv;j++) r[j] = (S[k*ls+j]+(sigma-cb[k])*r[k*lr+j]-cg[k]*r[(k-1)*lr+j])/ca[k];
128: ss = r; lss = lr; off = 1; alpha = -1.0; a = pep->sfactor;
129: } else {
130: ss = S; lss = ls; off = 0; alpha = -ca[deg-1]; a = 1.0;
131: }
132: BVMultVec(V,1.0,0.0,v,ss+off*lss);
133: if (pep->Dr) { /* balancing */
134: VecPointwiseMult(v,v,pep->Dr);
135: }
136: STMatMult(pep->st,off,v,q);
137: VecScale(q,a);
138: for (j=1+off;j<deg+off-1;j++) {
139: BVMultVec(V,1.0,0.0,v,ss+j*lss);
140: if (pep->Dr) {
141: VecPointwiseMult(v,v,pep->Dr);
142: }
143: STMatMult(pep->st,j,v,t);
144: a *= pep->sfactor;
145: VecAXPY(q,a,t);
146: }
147: if (sinvert) {
148: BVMultVec(V,1.0,0.0,v,ss);
149: if (pep->Dr) {
150: VecPointwiseMult(v,v,pep->Dr);
151: }
152: STMatMult(pep->st,deg,v,t);
153: a *= pep->sfactor;
154: VecAXPY(q,a,t);
155: } else {
156: BVMultVec(V,1.0,0.0,ve,ss+(deg-1)*lss);
157: if (pep->Dr) {
158: VecPointwiseMult(ve,ve,pep->Dr);
159: }
160: a *= pep->sfactor;
161: STMatMult(pep->st,deg-1,ve,t);
162: VecAXPY(q,a,t);
163: a *= pep->sfactor;
164: }
165: if (flg || !sinvert) alpha /= a;
166: STMatSolve(pep->st,q,t);
167: VecScale(t,alpha);
168: if (!sinvert) {
169: if (cg[deg-1]!=0) { VecAXPY(t,cg[deg-1],v); }
170: if (cb[deg-1]!=0) { VecAXPY(t,cb[deg-1],ve); }
171: }
172: if (pep->Dr) {
173: VecPointwiseDivide(t,t,pep->Dr);
174: }
175: return(0);
176: }
178: /*
179: Compute TOAR coefficients of the blocks of the new Arnoldi vector computed
180: */
181: static PetscErrorCode PEPTOARCoefficients(PEP pep,PetscBool sinvert,PetscScalar sigma,PetscInt nv,PetscScalar *S,PetscInt ls,PetscScalar *r,PetscInt lr,PetscScalar *x)
182: {
183: PetscInt k,j,nmat=pep->nmat,d=nmat-1;
184: PetscReal *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat;
185: PetscScalar t=1.0,tp=0.0,tt;
188: if (sinvert) {
189: for (k=1;k<d;k++) {
190: tt = t;
191: t = ((sigma-cb[k-1])*t-cg[k-1]*tp)/ca[k-1]; /* k-th basis polynomial */
192: tp = tt;
193: for (j=0;j<=nv;j++) r[k*lr+j] += t*x[j];
194: }
195: } else {
196: for (j=0;j<=nv;j++) r[j] = (cb[0]-sigma)*S[j]+ca[0]*S[ls+j];
197: for (k=1;k<d-1;k++) {
198: for (j=0;j<=nv;j++) r[k*lr+j] = (cb[k]-sigma)*S[k*ls+j]+ca[k]*S[(k+1)*ls+j]+cg[k]*S[(k-1)*ls+j];
199: }
200: if (sigma!=0.0) for (j=0;j<=nv;j++) r[(d-1)*lr+j] -= sigma*S[(d-1)*ls+j];
201: }
202: return(0);
203: }
205: /*
206: Compute a run of Arnoldi iterations dim(work)=ld
207: */
208: static PetscErrorCode PEPTOARrun(PEP pep,PetscScalar sigma,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,PetscBool *breakdown,Vec *t_)
209: {
211: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
212: PetscInt j,m=*M,deg=pep->nmat-1,ld;
213: PetscInt lds,nqt,l;
214: Vec t;
215: PetscReal norm;
216: PetscBool flg,sinvert=PETSC_FALSE,lindep;
217: PetscScalar *x,*S;
218: Mat MS;
221: BVTensorGetFactors(ctx->V,NULL,&MS);
222: MatDenseGetArray(MS,&S);
223: BVGetSizes(pep->V,NULL,NULL,&ld);
224: lds = ld*deg;
225: BVGetActiveColumns(pep->V,&l,&nqt);
226: STGetTransform(pep->st,&flg);
227: if (!flg) {
228: /* spectral transformation handled by the solver */
229: PetscObjectTypeCompareAny((PetscObject)pep->st,&flg,STSINVERT,STSHIFT,"");
230: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"ST type not supported for TOAR without transforming matrices");
231: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinvert);
232: }
233: BVSetActiveColumns(ctx->V,0,m);
234: for (j=k;j<m;j++) {
235: /* apply operator */
236: BVGetColumn(pep->V,nqt,&t);
237: PEPTOARExtendBasis(pep,sinvert,sigma,S+j*lds,ld,nqt,pep->V,t,S+(j+1)*lds,ld,t_);
238: BVRestoreColumn(pep->V,nqt,&t);
240: /* orthogonalize */
241: if (sinvert) x = S+(j+1)*lds;
242: else x = S+(deg-1)*ld+(j+1)*lds;
243: BVOrthogonalizeColumn(pep->V,nqt,x,&norm,&lindep);
244: if (!lindep) {
245: x[nqt] = norm;
246: BVScaleColumn(pep->V,nqt,1.0/norm);
247: nqt++;
248: }
250: PEPTOARCoefficients(pep,sinvert,sigma,nqt-1,S+j*lds,ld,S+(j+1)*lds,ld,x);
252: /* level-2 orthogonalization */
253: BVOrthogonalizeColumn(ctx->V,j+1,H+j*ldh,&norm,breakdown);
254: H[j+1+ldh*j] = norm;
255: if (*breakdown) {
256: *M = j+1;
257: break;
258: }
259: BVScaleColumn(ctx->V,j+1,1.0/norm);
260: BVSetActiveColumns(pep->V,l,nqt);
261: }
262: BVSetActiveColumns(ctx->V,0,*M);
263: MatDenseRestoreArray(MS,&S);
264: BVTensorRestoreFactors(ctx->V,NULL,&MS);
265: return(0);
266: }
268: /*
269: Computes T_j = phi_idx(T). In T_j and T_p are phi_{idx-1}(T)
270: and phi_{idx-2}(T) respectively or null if idx=0,1.
271: Tp and Tj are input/output arguments
272: */
273: static PetscErrorCode PEPEvaluateBasisM(PEP pep,PetscInt k,PetscScalar *T,PetscInt ldt,PetscInt idx,PetscScalar **Tp,PetscScalar **Tj)
274: {
276: PetscInt i;
277: PetscReal *ca,*cb,*cg;
278: PetscScalar *pt,g,a;
279: PetscBLASInt k_,ldt_;
282: if (idx==0) {
283: PetscArrayzero(*Tj,k*k);
284: PetscArrayzero(*Tp,k*k);
285: for (i=0;i<k;i++) (*Tj)[i+i*k] = 1.0;
286: } else {
287: PetscBLASIntCast(ldt,&ldt_);
288: PetscBLASIntCast(k,&k_);
289: ca = pep->pbc; cb = pep->pbc+pep->nmat; cg = pep->pbc+2*pep->nmat;
290: for (i=0;i<k;i++) T[i*ldt+i] -= cb[idx-1];
291: a = 1/ca[idx-1];
292: g = (idx==1)?0.0:-cg[idx-1]/ca[idx-1];
293: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&a,T,&ldt_,*Tj,&k_,&g,*Tp,&k_));
294: pt = *Tj; *Tj = *Tp; *Tp = pt;
295: for (i=0;i<k;i++) T[i*ldt+i] += cb[idx-1];
296: }
297: return(0);
298: }
300: static PetscErrorCode PEPExtractInvariantPair(PEP pep,PetscScalar sigma,PetscInt sr,PetscInt k,PetscScalar *S,PetscInt ld,PetscInt deg,PetscScalar *H,PetscInt ldh)
301: {
303: PetscInt i,j,jj,lds,ldt,d=pep->nmat-1,idxcpy=0;
304: PetscScalar *At,*Bt,*Hj,*Hp,*T,sone=1.0,g,a,*pM,*work;
305: PetscBLASInt k_,sr_,lds_,ldh_,info,*p,lwork,ldt_;
306: PetscBool transf=PETSC_FALSE,flg;
307: PetscReal norm,maxnrm,*rwork;
308: BV *R,Y;
309: Mat M,*A;
312: if (k==0) return(0);
313: lds = deg*ld;
314: PetscCalloc6(k,&p,sr*k,&At,k*k,&Bt,k*k,&Hj,k*k,&Hp,sr*k,&work);
315: PetscBLASIntCast(sr,&sr_);
316: PetscBLASIntCast(k,&k_);
317: PetscBLASIntCast(lds,&lds_);
318: PetscBLASIntCast(ldh,&ldh_);
319: STGetTransform(pep->st,&flg);
320: if (!flg) {
321: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&flg);
322: if (flg || sigma!=0.0) transf=PETSC_TRUE;
323: }
324: if (transf) {
325: PetscMalloc1(k*k,&T);
326: ldt = k;
327: for (i=0;i<k;i++) {
328: PetscArraycpy(T+k*i,H+i*ldh,k);
329: }
330: if (flg) {
331: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
332: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&k_,&k_,T,&k_,p,&info));
333: SlepcCheckLapackInfo("getrf",info);
334: PetscBLASIntCast(sr*k,&lwork);
335: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&k_,T,&k_,p,work,&lwork,&info));
336: SlepcCheckLapackInfo("getri",info);
337: PetscFPTrapPop();
338: }
339: if (sigma!=0.0) for (i=0;i<k;i++) T[i+k*i] += sigma;
340: } else {
341: T = H; ldt = ldh;
342: }
343: PetscBLASIntCast(ldt,&ldt_);
344: switch (pep->extract) {
345: case PEP_EXTRACT_NONE:
346: break;
347: case PEP_EXTRACT_NORM:
348: if (pep->basis == PEP_BASIS_MONOMIAL) {
349: PetscBLASIntCast(ldt,&ldt_);
350: PetscMalloc1(k,&rwork);
351: norm = LAPACKlange_("F",&k_,&k_,T,&ldt_,rwork);
352: PetscFree(rwork);
353: if (norm>1.0) idxcpy = d-1;
354: } else {
355: PetscBLASIntCast(ldt,&ldt_);
356: PetscMalloc1(k,&rwork);
357: maxnrm = 0.0;
358: for (i=0;i<pep->nmat-1;i++) {
359: PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
360: norm = LAPACKlange_("F",&k_,&k_,Hj,&k_,rwork);
361: if (norm > maxnrm) {
362: idxcpy = i;
363: maxnrm = norm;
364: }
365: }
366: PetscFree(rwork);
367: }
368: if (idxcpy>0) {
369: /* copy block idxcpy of S to the first one */
370: for (j=0;j<k;j++) {
371: PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr);
372: }
373: }
374: break;
375: case PEP_EXTRACT_RESIDUAL:
376: STGetTransform(pep->st,&flg);
377: if (flg) {
378: PetscMalloc1(pep->nmat,&A);
379: for (i=0;i<pep->nmat;i++) {
380: STGetMatrixTransformed(pep->st,i,A+i);
381: }
382: } else A = pep->A;
383: PetscMalloc1(pep->nmat-1,&R);
384: for (i=0;i<pep->nmat-1;i++) {
385: BVDuplicateResize(pep->V,k,R+i);
386: }
387: BVDuplicateResize(pep->V,sr,&Y);
388: MatCreateSeqDense(PETSC_COMM_SELF,sr,k,NULL,&M);
389: g = 0.0; a = 1.0;
390: BVSetActiveColumns(pep->V,0,sr);
391: for (j=0;j<pep->nmat;j++) {
392: BVMatMult(pep->V,A[j],Y);
393: PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
394: for (i=0;i<pep->nmat-1;i++) {
395: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&k_,&k_,&a,S+i*ld,&lds_,Hj,&k_,&g,At,&sr_));
396: MatDenseGetArray(M,&pM);
397: for (jj=0;jj<k;jj++) {
398: PetscArraycpy(pM+jj*sr,At+jj*sr,sr);
399: }
400: MatDenseRestoreArray(M,&pM);
401: BVMult(R[i],1.0,(i==0)?0.0:1.0,Y,M);
402: }
403: }
405: /* frobenius norm */
406: maxnrm = 0.0;
407: for (i=0;i<pep->nmat-1;i++) {
408: BVNorm(R[i],NORM_FROBENIUS,&norm);
409: if (maxnrm > norm) {
410: maxnrm = norm;
411: idxcpy = i;
412: }
413: }
414: if (idxcpy>0) {
415: /* copy block idxcpy of S to the first one */
416: for (j=0;j<k;j++) {
417: PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr);
418: }
419: }
420: if (flg) PetscFree(A);
421: for (i=0;i<pep->nmat-1;i++) {
422: BVDestroy(&R[i]);
423: }
424: PetscFree(R);
425: BVDestroy(&Y);
426: MatDestroy(&M);
427: break;
428: case PEP_EXTRACT_STRUCTURED:
429: for (j=0;j<k;j++) Bt[j+j*k] = 1.0;
430: for (j=0;j<sr;j++) {
431: for (i=0;i<k;i++) At[j*k+i] = PetscConj(S[i*lds+j]);
432: }
433: PEPEvaluateBasisM(pep,k,T,ldt,0,&Hp,&Hj);
434: for (i=1;i<deg;i++) {
435: PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
436: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&sr_,&k_,&sone,Hj,&k_,S+i*ld,&lds_,&sone,At,&k_));
437: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&k_,&k_,&sone,Hj,&k_,Hj,&k_,&sone,Bt,&k_));
438: }
439: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
440: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&k_,&sr_,Bt,&k_,p,At,&k_,&info));
441: PetscFPTrapPop();
442: SlepcCheckLapackInfo("gesv",info);
443: for (j=0;j<sr;j++) {
444: for (i=0;i<k;i++) S[i*lds+j] = PetscConj(At[j*k+i]);
445: }
446: break;
447: }
448: if (transf) { PetscFree(T); }
449: PetscFree6(p,At,Bt,Hj,Hp,work);
450: return(0);
451: }
453: PetscErrorCode PEPSolve_TOAR(PEP pep)
454: {
456: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
457: PetscInt i,j,k,l,nv=0,ld,lds,ldds,nq=0,nconv=0;
458: PetscInt nmat=pep->nmat,deg=nmat-1;
459: PetscScalar *S,*H,sigma;
460: PetscReal beta;
461: PetscBool breakdown=PETSC_FALSE,flg,falselock=PETSC_FALSE,sinv=PETSC_FALSE;
462: Mat MS,MQ;
465: PetscCitationsRegister(citation,&cited);
466: if (ctx->lock) {
467: /* undocumented option to use a cheaper locking instead of the true locking */
468: PetscOptionsGetBool(NULL,NULL,"-pep_toar_falselocking",&falselock,NULL);
469: }
470: DSGetLeadingDimension(pep->ds,&ldds);
471: STGetShift(pep->st,&sigma);
473: /* update polynomial basis coefficients */
474: STGetTransform(pep->st,&flg);
475: if (pep->sfactor!=1.0) {
476: for (i=0;i<nmat;i++) {
477: pep->pbc[nmat+i] /= pep->sfactor;
478: pep->pbc[2*nmat+i] /= pep->sfactor*pep->sfactor;
479: }
480: if (!flg) {
481: pep->target /= pep->sfactor;
482: RGPushScale(pep->rg,1.0/pep->sfactor);
483: STScaleShift(pep->st,1.0/pep->sfactor);
484: sigma /= pep->sfactor;
485: } else {
486: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
487: pep->target = sinv?pep->target*pep->sfactor:pep->target/pep->sfactor;
488: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
489: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
490: }
491: }
493: if (flg) sigma = 0.0;
495: /* clean projected matrix (including the extra-arrow) */
496: DSGetArray(pep->ds,DS_MAT_A,&H);
497: PetscArrayzero(H,ldds*ldds);
498: DSRestoreArray(pep->ds,DS_MAT_A,&H);
500: /* Get the starting Arnoldi vector */
501: BVTensorBuildFirstColumn(ctx->V,pep->nini);
503: /* restart loop */
504: l = 0;
505: while (pep->reason == PEP_CONVERGED_ITERATING) {
506: pep->its++;
508: /* compute an nv-step Lanczos factorization */
509: nv = PetscMax(PetscMin(nconv+pep->mpd,pep->ncv),nv);
510: DSGetArray(pep->ds,DS_MAT_A,&H);
511: PEPTOARrun(pep,sigma,H,ldds,pep->nconv+l,&nv,&breakdown,pep->work);
512: beta = PetscAbsScalar(H[(nv-1)*ldds+nv]);
513: DSRestoreArray(pep->ds,DS_MAT_A,&H);
514: DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
515: if (l==0) {
516: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
517: } else {
518: DSSetState(pep->ds,DS_STATE_RAW);
519: }
520: BVSetActiveColumns(ctx->V,pep->nconv,nv);
522: /* solve projected problem */
523: DSSolve(pep->ds,pep->eigr,pep->eigi);
524: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
525: DSUpdateExtraRow(pep->ds);
526: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
528: /* check convergence */
529: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
530: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
532: /* update l */
533: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
534: else {
535: l = (nv==k)?0:PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
536: DSGetTruncateSize(pep->ds,k,nv,&l);
537: if (!breakdown) {
538: /* prepare the Rayleigh quotient for restart */
539: DSTruncate(pep->ds,k+l,PETSC_FALSE);
540: }
541: }
542: nconv = k;
543: if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
544: if (l) { PetscInfo1(pep,"Preparing to restart keeping l=%D vectors\n",l); }
546: /* update S */
547: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
548: BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
549: MatDestroy(&MQ);
551: /* copy last column of S */
552: BVCopyColumn(ctx->V,nv,k+l);
554: if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
555: /* stop if breakdown */
556: PetscInfo2(pep,"Breakdown TOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
557: pep->reason = PEP_DIVERGED_BREAKDOWN;
558: }
559: if (pep->reason != PEP_CONVERGED_ITERATING) l--;
560: /* truncate S */
561: BVGetActiveColumns(pep->V,NULL,&nq);
562: if (k+l+deg<=nq) {
563: BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
564: if (!falselock && ctx->lock) {
565: BVTensorCompress(ctx->V,k-pep->nconv);
566: } else {
567: BVTensorCompress(ctx->V,0);
568: }
569: }
570: pep->nconv = k;
571: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
572: }
573: if (pep->nconv>0) {
574: /* {V*S_nconv^i}_{i=0}^{d-1} has rank nconv instead of nconv+d-1. Force zeros in each S_nconv^i block */
575: BVSetActiveColumns(ctx->V,0,pep->nconv);
576: BVGetActiveColumns(pep->V,NULL,&nq);
577: BVSetActiveColumns(pep->V,0,nq);
578: if (nq>pep->nconv) {
579: BVTensorCompress(ctx->V,pep->nconv);
580: BVSetActiveColumns(pep->V,0,pep->nconv);
581: nq = pep->nconv;
582: }
584: /* perform Newton refinement if required */
585: if (pep->refine==PEP_REFINE_MULTIPLE && pep->rits>0) {
586: /* extract invariant pair */
587: BVTensorGetFactors(ctx->V,NULL,&MS);
588: MatDenseGetArray(MS,&S);
589: DSGetArray(pep->ds,DS_MAT_A,&H);
590: BVGetSizes(pep->V,NULL,NULL,&ld);
591: lds = deg*ld;
592: PEPExtractInvariantPair(pep,sigma,nq,pep->nconv,S,ld,deg,H,ldds);
593: DSRestoreArray(pep->ds,DS_MAT_A,&H);
594: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
595: DSSetState(pep->ds,DS_STATE_RAW);
596: PEPNewtonRefinement_TOAR(pep,sigma,&pep->rits,NULL,pep->nconv,S,lds);
597: DSSolve(pep->ds,pep->eigr,pep->eigi);
598: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
599: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
600: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
601: BVMultInPlace(ctx->V,MQ,0,pep->nconv);
602: MatDestroy(&MQ);
603: MatDenseRestoreArray(MS,&S);
604: BVTensorRestoreFactors(ctx->V,NULL,&MS);
605: }
606: }
607: STGetTransform(pep->st,&flg);
608: if (pep->refine!=PEP_REFINE_MULTIPLE || pep->rits==0) {
609: if (!flg && pep->ops->backtransform) {
610: (*pep->ops->backtransform)(pep);
611: }
612: if (pep->sfactor!=1.0) {
613: for (j=0;j<pep->nconv;j++) {
614: pep->eigr[j] *= pep->sfactor;
615: pep->eigi[j] *= pep->sfactor;
616: }
617: /* restore original values */
618: for (i=0;i<pep->nmat;i++){
619: pep->pbc[pep->nmat+i] *= pep->sfactor;
620: pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
621: }
622: }
623: }
624: /* restore original values */
625: if (!flg) {
626: pep->target *= pep->sfactor;
627: STScaleShift(pep->st,pep->sfactor);
628: } else {
629: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
630: pep->target = (sinv)?pep->target/pep->sfactor:pep->target*pep->sfactor;
631: }
632: if (pep->sfactor!=1.0) { RGPopScale(pep->rg); }
634: /* change the state to raw so that DSVectors() computes eigenvectors from scratch */
635: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
636: DSSetState(pep->ds,DS_STATE_RAW);
637: return(0);
638: }
640: static PetscErrorCode PEPTOARSetRestart_TOAR(PEP pep,PetscReal keep)
641: {
642: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
645: if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
646: else {
647: if (keep<0.1 || keep>0.9) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
648: ctx->keep = keep;
649: }
650: return(0);
651: }
653: /*@
654: PEPTOARSetRestart - Sets the restart parameter for the TOAR
655: method, in particular the proportion of basis vectors that must be kept
656: after restart.
658: Logically Collective on pep
660: Input Parameters:
661: + pep - the eigenproblem solver context
662: - keep - the number of vectors to be kept at restart
664: Options Database Key:
665: . -pep_toar_restart - Sets the restart parameter
667: Notes:
668: Allowed values are in the range [0.1,0.9]. The default is 0.5.
670: Level: advanced
672: .seealso: PEPTOARGetRestart()
673: @*/
674: PetscErrorCode PEPTOARSetRestart(PEP pep,PetscReal keep)
675: {
681: PetscTryMethod(pep,"PEPTOARSetRestart_C",(PEP,PetscReal),(pep,keep));
682: return(0);
683: }
685: static PetscErrorCode PEPTOARGetRestart_TOAR(PEP pep,PetscReal *keep)
686: {
687: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
690: *keep = ctx->keep;
691: return(0);
692: }
694: /*@
695: PEPTOARGetRestart - Gets the restart parameter used in the TOAR method.
697: Not Collective
699: Input Parameter:
700: . pep - the eigenproblem solver context
702: Output Parameter:
703: . keep - the restart parameter
705: Level: advanced
707: .seealso: PEPTOARSetRestart()
708: @*/
709: PetscErrorCode PEPTOARGetRestart(PEP pep,PetscReal *keep)
710: {
716: PetscUseMethod(pep,"PEPTOARGetRestart_C",(PEP,PetscReal*),(pep,keep));
717: return(0);
718: }
720: static PetscErrorCode PEPTOARSetLocking_TOAR(PEP pep,PetscBool lock)
721: {
722: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
725: ctx->lock = lock;
726: return(0);
727: }
729: /*@
730: PEPTOARSetLocking - Choose between locking and non-locking variants of
731: the TOAR method.
733: Logically Collective on pep
735: Input Parameters:
736: + pep - the eigenproblem solver context
737: - lock - true if the locking variant must be selected
739: Options Database Key:
740: . -pep_toar_locking - Sets the locking flag
742: Notes:
743: The default is to lock converged eigenpairs when the method restarts.
744: This behaviour can be changed so that all directions are kept in the
745: working subspace even if already converged to working accuracy (the
746: non-locking variant).
748: Level: advanced
750: .seealso: PEPTOARGetLocking()
751: @*/
752: PetscErrorCode PEPTOARSetLocking(PEP pep,PetscBool lock)
753: {
759: PetscTryMethod(pep,"PEPTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
760: return(0);
761: }
763: static PetscErrorCode PEPTOARGetLocking_TOAR(PEP pep,PetscBool *lock)
764: {
765: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
768: *lock = ctx->lock;
769: return(0);
770: }
772: /*@
773: PEPTOARGetLocking - Gets the locking flag used in the TOAR method.
775: Not Collective
777: Input Parameter:
778: . pep - the eigenproblem solver context
780: Output Parameter:
781: . lock - the locking flag
783: Level: advanced
785: .seealso: PEPTOARSetLocking()
786: @*/
787: PetscErrorCode PEPTOARGetLocking(PEP pep,PetscBool *lock)
788: {
794: PetscUseMethod(pep,"PEPTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
795: return(0);
796: }
798: PetscErrorCode PEPSetFromOptions_TOAR(PetscOptionItems *PetscOptionsObject,PEP pep)
799: {
801: PetscBool flg,lock;
802: PetscReal keep;
805: PetscOptionsHead(PetscOptionsObject,"PEP TOAR Options");
807: PetscOptionsReal("-pep_toar_restart","Proportion of vectors kept after restart","PEPTOARSetRestart",0.5,&keep,&flg);
808: if (flg) { PEPTOARSetRestart(pep,keep); }
810: PetscOptionsBool("-pep_toar_locking","Choose between locking and non-locking variants","PEPTOARSetLocking",PETSC_FALSE,&lock,&flg);
811: if (flg) { PEPTOARSetLocking(pep,lock); }
813: PetscOptionsTail();
814: return(0);
815: }
817: PetscErrorCode PEPView_TOAR(PEP pep,PetscViewer viewer)
818: {
820: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
821: PetscBool isascii;
824: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
825: if (isascii) {
826: PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
827: PetscViewerASCIIPrintf(viewer," using the %slocking variant\n",ctx->lock?"":"non-");
828: }
829: return(0);
830: }
832: PetscErrorCode PEPDestroy_TOAR(PEP pep)
833: {
835: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
838: BVDestroy(&ctx->V);
839: PetscFree(pep->data);
840: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",NULL);
841: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",NULL);
842: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",NULL);
843: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",NULL);
844: return(0);
845: }
847: SLEPC_EXTERN PetscErrorCode PEPCreate_TOAR(PEP pep)
848: {
849: PEP_TOAR *ctx;
853: PetscNewLog(pep,&ctx);
854: pep->data = (void*)ctx;
856: pep->lineariz = PETSC_TRUE;
857: ctx->lock = PETSC_TRUE;
859: pep->ops->solve = PEPSolve_TOAR;
860: pep->ops->setup = PEPSetUp_TOAR;
861: pep->ops->setfromoptions = PEPSetFromOptions_TOAR;
862: pep->ops->destroy = PEPDestroy_TOAR;
863: pep->ops->view = PEPView_TOAR;
864: pep->ops->backtransform = PEPBackTransform_Default;
865: pep->ops->computevectors = PEPComputeVectors_Default;
866: pep->ops->extractvectors = PEPExtractVectors_TOAR;
868: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",PEPTOARSetRestart_TOAR);
869: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",PEPTOARGetRestart_TOAR);
870: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",PEPTOARSetLocking_TOAR);
871: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",PEPTOARGetLocking_TOAR);
872: return(0);
873: }