Actual source code: ptoar.c
slepc-3.16.3 2022-04-11
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: PetscFree(pep->solvematcoeffs);
80: PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
81: PetscLogObjectMemory((PetscObject)pep,pep->nmat*sizeof(PetscScalar));
82: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
83: if (sinv) {
84: PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL);
85: } else {
86: for (i=0;i<pep->nmat-1;i++) pep->solvematcoeffs[i] = 0.0;
87: pep->solvematcoeffs[pep->nmat-1] = 1.0;
88: }
89: }
90: BVDestroy(&ctx->V);
91: BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
92: return(0);
93: }
95: /*
96: Extend the TOAR basis by applying the the matrix operator
97: over a vector which is decomposed in the TOAR way
98: Input:
99: - pbc: array containing the polynomial basis coefficients
100: - S,V: define the latest Arnoldi vector (nv vectors in V)
101: Output:
102: - t: new vector extending the TOAR basis
103: - r: temporary coefficients to compute the TOAR coefficients
104: for the new Arnoldi vector
105: Workspace: t_ (two vectors)
106: */
107: 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_)
108: {
110: PetscInt nmat=pep->nmat,deg=nmat-1,k,j,off=0,lss;
111: Vec v=t_[0],ve=t_[1],q=t_[2];
112: PetscScalar alpha=1.0,*ss,a;
113: PetscReal *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat;
114: PetscBool flg;
117: BVSetActiveColumns(pep->V,0,nv);
118: STGetTransform(pep->st,&flg);
119: if (sinvert) {
120: for (j=0;j<nv;j++) {
121: if (deg>1) r[lr+j] = S[j]/ca[0];
122: if (deg>2) r[2*lr+j] = (S[ls+j]+(sigma-cb[1])*r[lr+j])/ca[1];
123: }
124: for (k=2;k<deg-1;k++) {
125: 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];
126: }
127: k = deg-1;
128: 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];
129: ss = r; lss = lr; off = 1; alpha = -1.0; a = pep->sfactor;
130: } else {
131: ss = S; lss = ls; off = 0; alpha = -ca[deg-1]; a = 1.0;
132: }
133: BVMultVec(V,1.0,0.0,v,ss+off*lss);
134: if (pep->Dr) { /* balancing */
135: VecPointwiseMult(v,v,pep->Dr);
136: }
137: STMatMult(pep->st,off,v,q);
138: VecScale(q,a);
139: for (j=1+off;j<deg+off-1;j++) {
140: BVMultVec(V,1.0,0.0,v,ss+j*lss);
141: if (pep->Dr) {
142: VecPointwiseMult(v,v,pep->Dr);
143: }
144: STMatMult(pep->st,j,v,t);
145: a *= pep->sfactor;
146: VecAXPY(q,a,t);
147: }
148: if (sinvert) {
149: BVMultVec(V,1.0,0.0,v,ss);
150: if (pep->Dr) {
151: VecPointwiseMult(v,v,pep->Dr);
152: }
153: STMatMult(pep->st,deg,v,t);
154: a *= pep->sfactor;
155: VecAXPY(q,a,t);
156: } else {
157: BVMultVec(V,1.0,0.0,ve,ss+(deg-1)*lss);
158: if (pep->Dr) {
159: VecPointwiseMult(ve,ve,pep->Dr);
160: }
161: a *= pep->sfactor;
162: STMatMult(pep->st,deg-1,ve,t);
163: VecAXPY(q,a,t);
164: a *= pep->sfactor;
165: }
166: if (flg || !sinvert) alpha /= a;
167: STMatSolve(pep->st,q,t);
168: VecScale(t,alpha);
169: if (!sinvert) {
170: if (cg[deg-1]!=0) { VecAXPY(t,cg[deg-1],v); }
171: if (cb[deg-1]!=0) { VecAXPY(t,cb[deg-1],ve); }
172: }
173: if (pep->Dr) {
174: VecPointwiseDivide(t,t,pep->Dr);
175: }
176: return(0);
177: }
179: /*
180: Compute TOAR coefficients of the blocks of the new Arnoldi vector computed
181: */
182: static PetscErrorCode PEPTOARCoefficients(PEP pep,PetscBool sinvert,PetscScalar sigma,PetscInt nv,PetscScalar *S,PetscInt ls,PetscScalar *r,PetscInt lr,PetscScalar *x)
183: {
184: PetscInt k,j,nmat=pep->nmat,d=nmat-1;
185: PetscReal *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat;
186: PetscScalar t=1.0,tp=0.0,tt;
189: if (sinvert) {
190: for (k=1;k<d;k++) {
191: tt = t;
192: t = ((sigma-cb[k-1])*t-cg[k-1]*tp)/ca[k-1]; /* k-th basis polynomial */
193: tp = tt;
194: for (j=0;j<=nv;j++) r[k*lr+j] += t*x[j];
195: }
196: } else {
197: for (j=0;j<=nv;j++) r[j] = (cb[0]-sigma)*S[j]+ca[0]*S[ls+j];
198: for (k=1;k<d-1;k++) {
199: 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];
200: }
201: if (sigma!=0.0) for (j=0;j<=nv;j++) r[(d-1)*lr+j] -= sigma*S[(d-1)*ls+j];
202: }
203: return(0);
204: }
206: /*
207: Compute a run of Arnoldi iterations dim(work)=ld
208: */
209: static PetscErrorCode PEPTOARrun(PEP pep,PetscScalar sigma,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,PetscBool *breakdown,Vec *t_)
210: {
212: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
213: PetscInt j,m=*M,deg=pep->nmat-1,ld;
214: PetscInt lds,nqt,l;
215: Vec t;
216: PetscReal norm;
217: PetscBool flg,sinvert=PETSC_FALSE,lindep;
218: PetscScalar *x,*S;
219: Mat MS;
222: BVTensorGetFactors(ctx->V,NULL,&MS);
223: MatDenseGetArray(MS,&S);
224: BVGetSizes(pep->V,NULL,NULL,&ld);
225: lds = ld*deg;
226: BVGetActiveColumns(pep->V,&l,&nqt);
227: STGetTransform(pep->st,&flg);
228: if (!flg) {
229: /* spectral transformation handled by the solver */
230: PetscObjectTypeCompareAny((PetscObject)pep->st,&flg,STSINVERT,STSHIFT,"");
231: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"ST type not supported for TOAR without transforming matrices");
232: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinvert);
233: }
234: BVSetActiveColumns(ctx->V,0,m);
235: for (j=k;j<m;j++) {
236: /* apply operator */
237: BVGetColumn(pep->V,nqt,&t);
238: PEPTOARExtendBasis(pep,sinvert,sigma,S+j*lds,ld,nqt,pep->V,t,S+(j+1)*lds,ld,t_);
239: BVRestoreColumn(pep->V,nqt,&t);
241: /* orthogonalize */
242: if (sinvert) x = S+(j+1)*lds;
243: else x = S+(deg-1)*ld+(j+1)*lds;
244: BVOrthogonalizeColumn(pep->V,nqt,x,&norm,&lindep);
245: if (!lindep) {
246: x[nqt] = norm;
247: BVScaleColumn(pep->V,nqt,1.0/norm);
248: nqt++;
249: }
251: PEPTOARCoefficients(pep,sinvert,sigma,nqt-1,S+j*lds,ld,S+(j+1)*lds,ld,x);
253: /* level-2 orthogonalization */
254: BVOrthogonalizeColumn(ctx->V,j+1,H+j*ldh,&norm,breakdown);
255: H[j+1+ldh*j] = norm;
256: if (*breakdown) {
257: *M = j+1;
258: break;
259: }
260: BVScaleColumn(ctx->V,j+1,1.0/norm);
261: BVSetActiveColumns(pep->V,l,nqt);
262: }
263: BVSetActiveColumns(ctx->V,0,*M);
264: MatDenseRestoreArray(MS,&S);
265: BVTensorRestoreFactors(ctx->V,NULL,&MS);
266: return(0);
267: }
269: /*
270: Computes T_j = phi_idx(T). In T_j and T_p are phi_{idx-1}(T)
271: and phi_{idx-2}(T) respectively or null if idx=0,1.
272: Tp and Tj are input/output arguments
273: */
274: static PetscErrorCode PEPEvaluateBasisM(PEP pep,PetscInt k,PetscScalar *T,PetscInt ldt,PetscInt idx,PetscScalar **Tp,PetscScalar **Tj)
275: {
277: PetscInt i;
278: PetscReal *ca,*cb,*cg;
279: PetscScalar *pt,g,a;
280: PetscBLASInt k_,ldt_;
283: if (idx==0) {
284: PetscArrayzero(*Tj,k*k);
285: PetscArrayzero(*Tp,k*k);
286: for (i=0;i<k;i++) (*Tj)[i+i*k] = 1.0;
287: } else {
288: PetscBLASIntCast(ldt,&ldt_);
289: PetscBLASIntCast(k,&k_);
290: ca = pep->pbc; cb = pep->pbc+pep->nmat; cg = pep->pbc+2*pep->nmat;
291: for (i=0;i<k;i++) T[i*ldt+i] -= cb[idx-1];
292: a = 1/ca[idx-1];
293: g = (idx==1)?0.0:-cg[idx-1]/ca[idx-1];
294: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&a,T,&ldt_,*Tj,&k_,&g,*Tp,&k_));
295: pt = *Tj; *Tj = *Tp; *Tp = pt;
296: for (i=0;i<k;i++) T[i*ldt+i] += cb[idx-1];
297: }
298: return(0);
299: }
301: static PetscErrorCode PEPExtractInvariantPair(PEP pep,PetscScalar sigma,PetscInt sr,PetscInt k,PetscScalar *S,PetscInt ld,PetscInt deg,PetscScalar *H,PetscInt ldh)
302: {
304: PetscInt i,j,jj,lds,ldt,d=pep->nmat-1,idxcpy=0;
305: PetscScalar *At,*Bt,*Hj,*Hp,*T,sone=1.0,g,a,*pM,*work;
306: PetscBLASInt k_,sr_,lds_,ldh_,info,*p,lwork,ldt_;
307: PetscBool transf=PETSC_FALSE,flg;
308: PetscReal norm,maxnrm,*rwork;
309: BV *R,Y;
310: Mat M,*A;
313: if (k==0) return(0);
314: lds = deg*ld;
315: PetscCalloc6(k,&p,sr*k,&At,k*k,&Bt,k*k,&Hj,k*k,&Hp,sr*k,&work);
316: PetscBLASIntCast(sr,&sr_);
317: PetscBLASIntCast(k,&k_);
318: PetscBLASIntCast(lds,&lds_);
319: PetscBLASIntCast(ldh,&ldh_);
320: STGetTransform(pep->st,&flg);
321: if (!flg) {
322: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&flg);
323: if (flg || sigma!=0.0) transf=PETSC_TRUE;
324: }
325: if (transf) {
326: PetscMalloc1(k*k,&T);
327: ldt = k;
328: for (i=0;i<k;i++) {
329: PetscArraycpy(T+k*i,H+i*ldh,k);
330: }
331: if (flg) {
332: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
333: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&k_,&k_,T,&k_,p,&info));
334: SlepcCheckLapackInfo("getrf",info);
335: PetscBLASIntCast(sr*k,&lwork);
336: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&k_,T,&k_,p,work,&lwork,&info));
337: SlepcCheckLapackInfo("getri",info);
338: PetscFPTrapPop();
339: }
340: if (sigma!=0.0) for (i=0;i<k;i++) T[i+k*i] += sigma;
341: } else {
342: T = H; ldt = ldh;
343: }
344: PetscBLASIntCast(ldt,&ldt_);
345: switch (pep->extract) {
346: case PEP_EXTRACT_NONE:
347: break;
348: case PEP_EXTRACT_NORM:
349: if (pep->basis == PEP_BASIS_MONOMIAL) {
350: PetscBLASIntCast(ldt,&ldt_);
351: PetscMalloc1(k,&rwork);
352: norm = LAPACKlange_("F",&k_,&k_,T,&ldt_,rwork);
353: PetscFree(rwork);
354: if (norm>1.0) idxcpy = d-1;
355: } else {
356: PetscBLASIntCast(ldt,&ldt_);
357: PetscMalloc1(k,&rwork);
358: maxnrm = 0.0;
359: for (i=0;i<pep->nmat-1;i++) {
360: PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
361: norm = LAPACKlange_("F",&k_,&k_,Hj,&k_,rwork);
362: if (norm > maxnrm) {
363: idxcpy = i;
364: maxnrm = norm;
365: }
366: }
367: PetscFree(rwork);
368: }
369: if (idxcpy>0) {
370: /* copy block idxcpy of S to the first one */
371: for (j=0;j<k;j++) {
372: PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr);
373: }
374: }
375: break;
376: case PEP_EXTRACT_RESIDUAL:
377: STGetTransform(pep->st,&flg);
378: if (flg) {
379: PetscMalloc1(pep->nmat,&A);
380: for (i=0;i<pep->nmat;i++) {
381: STGetMatrixTransformed(pep->st,i,A+i);
382: }
383: } else A = pep->A;
384: PetscMalloc1(pep->nmat-1,&R);
385: for (i=0;i<pep->nmat-1;i++) {
386: BVDuplicateResize(pep->V,k,R+i);
387: }
388: BVDuplicateResize(pep->V,sr,&Y);
389: MatCreateSeqDense(PETSC_COMM_SELF,sr,k,NULL,&M);
390: g = 0.0; a = 1.0;
391: BVSetActiveColumns(pep->V,0,sr);
392: for (j=0;j<pep->nmat;j++) {
393: BVMatMult(pep->V,A[j],Y);
394: PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
395: for (i=0;i<pep->nmat-1;i++) {
396: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&k_,&k_,&a,S+i*ld,&lds_,Hj,&k_,&g,At,&sr_));
397: MatDenseGetArray(M,&pM);
398: for (jj=0;jj<k;jj++) {
399: PetscArraycpy(pM+jj*sr,At+jj*sr,sr);
400: }
401: MatDenseRestoreArray(M,&pM);
402: BVMult(R[i],1.0,(i==0)?0.0:1.0,Y,M);
403: }
404: }
406: /* frobenius norm */
407: maxnrm = 0.0;
408: for (i=0;i<pep->nmat-1;i++) {
409: BVNorm(R[i],NORM_FROBENIUS,&norm);
410: if (maxnrm > norm) {
411: maxnrm = norm;
412: idxcpy = i;
413: }
414: }
415: if (idxcpy>0) {
416: /* copy block idxcpy of S to the first one */
417: for (j=0;j<k;j++) {
418: PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr);
419: }
420: }
421: if (flg) { PetscFree(A); }
422: for (i=0;i<pep->nmat-1;i++) {
423: BVDestroy(&R[i]);
424: }
425: PetscFree(R);
426: BVDestroy(&Y);
427: MatDestroy(&M);
428: break;
429: case PEP_EXTRACT_STRUCTURED:
430: for (j=0;j<k;j++) Bt[j+j*k] = 1.0;
431: for (j=0;j<sr;j++) {
432: for (i=0;i<k;i++) At[j*k+i] = PetscConj(S[i*lds+j]);
433: }
434: PEPEvaluateBasisM(pep,k,T,ldt,0,&Hp,&Hj);
435: for (i=1;i<deg;i++) {
436: PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
437: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&sr_,&k_,&sone,Hj,&k_,S+i*ld,&lds_,&sone,At,&k_));
438: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&k_,&k_,&sone,Hj,&k_,Hj,&k_,&sone,Bt,&k_));
439: }
440: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
441: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&k_,&sr_,Bt,&k_,p,At,&k_,&info));
442: PetscFPTrapPop();
443: SlepcCheckLapackInfo("gesv",info);
444: for (j=0;j<sr;j++) {
445: for (i=0;i<k;i++) S[i*lds+j] = PetscConj(At[j*k+i]);
446: }
447: break;
448: }
449: if (transf) { PetscFree(T); }
450: PetscFree6(p,At,Bt,Hj,Hp,work);
451: return(0);
452: }
454: PetscErrorCode PEPSolve_TOAR(PEP pep)
455: {
457: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
458: PetscInt i,j,k,l,nv=0,ld,lds,ldds,nq=0,nconv=0;
459: PetscInt nmat=pep->nmat,deg=nmat-1;
460: PetscScalar *S,*H,sigma;
461: PetscReal beta;
462: PetscBool breakdown=PETSC_FALSE,flg,falselock=PETSC_FALSE,sinv=PETSC_FALSE;
463: Mat MS,MQ;
466: PetscCitationsRegister(citation,&cited);
467: if (ctx->lock) {
468: /* undocumented option to use a cheaper locking instead of the true locking */
469: PetscOptionsGetBool(NULL,NULL,"-pep_toar_falselocking",&falselock,NULL);
470: }
471: DSGetLeadingDimension(pep->ds,&ldds);
472: STGetShift(pep->st,&sigma);
474: /* update polynomial basis coefficients */
475: STGetTransform(pep->st,&flg);
476: if (pep->sfactor!=1.0) {
477: for (i=0;i<nmat;i++) {
478: pep->pbc[nmat+i] /= pep->sfactor;
479: pep->pbc[2*nmat+i] /= pep->sfactor*pep->sfactor;
480: }
481: if (!flg) {
482: pep->target /= pep->sfactor;
483: RGPushScale(pep->rg,1.0/pep->sfactor);
484: STScaleShift(pep->st,1.0/pep->sfactor);
485: sigma /= pep->sfactor;
486: } else {
487: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
488: pep->target = sinv?pep->target*pep->sfactor:pep->target/pep->sfactor;
489: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
490: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
491: }
492: }
494: if (flg) sigma = 0.0;
496: /* clean projected matrix (including the extra-arrow) */
497: DSGetArray(pep->ds,DS_MAT_A,&H);
498: PetscArrayzero(H,ldds*ldds);
499: DSRestoreArray(pep->ds,DS_MAT_A,&H);
501: /* Get the starting Arnoldi vector */
502: BVTensorBuildFirstColumn(ctx->V,pep->nini);
504: /* restart loop */
505: l = 0;
506: while (pep->reason == PEP_CONVERGED_ITERATING) {
507: pep->its++;
509: /* compute an nv-step Lanczos factorization */
510: nv = PetscMax(PetscMin(nconv+pep->mpd,pep->ncv),nv);
511: DSGetArray(pep->ds,DS_MAT_A,&H);
512: PEPTOARrun(pep,sigma,H,ldds,pep->nconv+l,&nv,&breakdown,pep->work);
513: beta = PetscAbsScalar(H[(nv-1)*ldds+nv]);
514: DSRestoreArray(pep->ds,DS_MAT_A,&H);
515: DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l);
516: if (l==0) {
517: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
518: } else {
519: DSSetState(pep->ds,DS_STATE_RAW);
520: }
521: BVSetActiveColumns(ctx->V,pep->nconv,nv);
523: /* solve projected problem */
524: DSSolve(pep->ds,pep->eigr,pep->eigi);
525: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
526: DSUpdateExtraRow(pep->ds);
527: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
529: /* check convergence */
530: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
531: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
533: /* update l */
534: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
535: else {
536: l = (nv==k)?0:PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
537: DSGetTruncateSize(pep->ds,k,nv,&l);
538: if (!breakdown) {
539: /* prepare the Rayleigh quotient for restart */
540: DSTruncate(pep->ds,k+l,PETSC_FALSE);
541: }
542: }
543: nconv = k;
544: if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
545: if (l) { PetscInfo1(pep,"Preparing to restart keeping l=%D vectors\n",l); }
547: /* update S */
548: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
549: BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
550: MatDestroy(&MQ);
552: /* copy last column of S */
553: BVCopyColumn(ctx->V,nv,k+l);
555: if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
556: /* stop if breakdown */
557: PetscInfo2(pep,"Breakdown TOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
558: pep->reason = PEP_DIVERGED_BREAKDOWN;
559: }
560: if (pep->reason != PEP_CONVERGED_ITERATING) l--;
561: /* truncate S */
562: BVGetActiveColumns(pep->V,NULL,&nq);
563: if (k+l+deg<=nq) {
564: BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
565: if (!falselock && ctx->lock) {
566: BVTensorCompress(ctx->V,k-pep->nconv);
567: } else {
568: BVTensorCompress(ctx->V,0);
569: }
570: }
571: pep->nconv = k;
572: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
573: }
574: if (pep->nconv>0) {
575: /* {V*S_nconv^i}_{i=0}^{d-1} has rank nconv instead of nconv+d-1. Force zeros in each S_nconv^i block */
576: BVSetActiveColumns(ctx->V,0,pep->nconv);
577: BVGetActiveColumns(pep->V,NULL,&nq);
578: BVSetActiveColumns(pep->V,0,nq);
579: if (nq>pep->nconv) {
580: BVTensorCompress(ctx->V,pep->nconv);
581: BVSetActiveColumns(pep->V,0,pep->nconv);
582: nq = pep->nconv;
583: }
585: /* perform Newton refinement if required */
586: if (pep->refine==PEP_REFINE_MULTIPLE && pep->rits>0) {
587: /* extract invariant pair */
588: BVTensorGetFactors(ctx->V,NULL,&MS);
589: MatDenseGetArray(MS,&S);
590: DSGetArray(pep->ds,DS_MAT_A,&H);
591: BVGetSizes(pep->V,NULL,NULL,&ld);
592: lds = deg*ld;
593: PEPExtractInvariantPair(pep,sigma,nq,pep->nconv,S,ld,deg,H,ldds);
594: DSRestoreArray(pep->ds,DS_MAT_A,&H);
595: DSSetDimensions(pep->ds,pep->nconv,0,0);
596: DSSetState(pep->ds,DS_STATE_RAW);
597: PEPNewtonRefinement_TOAR(pep,sigma,&pep->rits,NULL,pep->nconv,S,lds);
598: DSSolve(pep->ds,pep->eigr,pep->eigi);
599: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
600: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
601: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
602: BVMultInPlace(ctx->V,MQ,0,pep->nconv);
603: MatDestroy(&MQ);
604: MatDenseRestoreArray(MS,&S);
605: BVTensorRestoreFactors(ctx->V,NULL,&MS);
606: }
607: }
608: STGetTransform(pep->st,&flg);
609: if (pep->refine!=PEP_REFINE_MULTIPLE || pep->rits==0) {
610: if (!flg && pep->ops->backtransform) {
611: (*pep->ops->backtransform)(pep);
612: }
613: if (pep->sfactor!=1.0) {
614: for (j=0;j<pep->nconv;j++) {
615: pep->eigr[j] *= pep->sfactor;
616: pep->eigi[j] *= pep->sfactor;
617: }
618: /* restore original values */
619: for (i=0;i<pep->nmat;i++) {
620: pep->pbc[pep->nmat+i] *= pep->sfactor;
621: pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
622: }
623: }
624: }
625: /* restore original values */
626: if (!flg) {
627: pep->target *= pep->sfactor;
628: STScaleShift(pep->st,pep->sfactor);
629: } else {
630: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
631: pep->target = (sinv)?pep->target/pep->sfactor:pep->target*pep->sfactor;
632: }
633: if (pep->sfactor!=1.0) { RGPopScale(pep->rg); }
635: /* change the state to raw so that DSVectors() computes eigenvectors from scratch */
636: DSSetDimensions(pep->ds,pep->nconv,0,0);
637: DSSetState(pep->ds,DS_STATE_RAW);
638: return(0);
639: }
641: static PetscErrorCode PEPTOARSetRestart_TOAR(PEP pep,PetscReal keep)
642: {
643: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
646: if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
647: else {
648: 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]");
649: ctx->keep = keep;
650: }
651: return(0);
652: }
654: /*@
655: PEPTOARSetRestart - Sets the restart parameter for the TOAR
656: method, in particular the proportion of basis vectors that must be kept
657: after restart.
659: Logically Collective on pep
661: Input Parameters:
662: + pep - the eigenproblem solver context
663: - keep - the number of vectors to be kept at restart
665: Options Database Key:
666: . -pep_toar_restart - Sets the restart parameter
668: Notes:
669: Allowed values are in the range [0.1,0.9]. The default is 0.5.
671: Level: advanced
673: .seealso: PEPTOARGetRestart()
674: @*/
675: PetscErrorCode PEPTOARSetRestart(PEP pep,PetscReal keep)
676: {
682: PetscTryMethod(pep,"PEPTOARSetRestart_C",(PEP,PetscReal),(pep,keep));
683: return(0);
684: }
686: static PetscErrorCode PEPTOARGetRestart_TOAR(PEP pep,PetscReal *keep)
687: {
688: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
691: *keep = ctx->keep;
692: return(0);
693: }
695: /*@
696: PEPTOARGetRestart - Gets the restart parameter used in the TOAR method.
698: Not Collective
700: Input Parameter:
701: . pep - the eigenproblem solver context
703: Output Parameter:
704: . keep - the restart parameter
706: Level: advanced
708: .seealso: PEPTOARSetRestart()
709: @*/
710: PetscErrorCode PEPTOARGetRestart(PEP pep,PetscReal *keep)
711: {
717: PetscUseMethod(pep,"PEPTOARGetRestart_C",(PEP,PetscReal*),(pep,keep));
718: return(0);
719: }
721: static PetscErrorCode PEPTOARSetLocking_TOAR(PEP pep,PetscBool lock)
722: {
723: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
726: ctx->lock = lock;
727: return(0);
728: }
730: /*@
731: PEPTOARSetLocking - Choose between locking and non-locking variants of
732: the TOAR method.
734: Logically Collective on pep
736: Input Parameters:
737: + pep - the eigenproblem solver context
738: - lock - true if the locking variant must be selected
740: Options Database Key:
741: . -pep_toar_locking - Sets the locking flag
743: Notes:
744: The default is to lock converged eigenpairs when the method restarts.
745: This behaviour can be changed so that all directions are kept in the
746: working subspace even if already converged to working accuracy (the
747: non-locking variant).
749: Level: advanced
751: .seealso: PEPTOARGetLocking()
752: @*/
753: PetscErrorCode PEPTOARSetLocking(PEP pep,PetscBool lock)
754: {
760: PetscTryMethod(pep,"PEPTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
761: return(0);
762: }
764: static PetscErrorCode PEPTOARGetLocking_TOAR(PEP pep,PetscBool *lock)
765: {
766: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
769: *lock = ctx->lock;
770: return(0);
771: }
773: /*@
774: PEPTOARGetLocking - Gets the locking flag used in the TOAR method.
776: Not Collective
778: Input Parameter:
779: . pep - the eigenproblem solver context
781: Output Parameter:
782: . lock - the locking flag
784: Level: advanced
786: .seealso: PEPTOARSetLocking()
787: @*/
788: PetscErrorCode PEPTOARGetLocking(PEP pep,PetscBool *lock)
789: {
795: PetscUseMethod(pep,"PEPTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
796: return(0);
797: }
799: PetscErrorCode PEPSetFromOptions_TOAR(PetscOptionItems *PetscOptionsObject,PEP pep)
800: {
802: PetscBool flg,lock;
803: PetscReal keep;
806: PetscOptionsHead(PetscOptionsObject,"PEP TOAR Options");
808: PetscOptionsReal("-pep_toar_restart","Proportion of vectors kept after restart","PEPTOARSetRestart",0.5,&keep,&flg);
809: if (flg) { PEPTOARSetRestart(pep,keep); }
811: PetscOptionsBool("-pep_toar_locking","Choose between locking and non-locking variants","PEPTOARSetLocking",PETSC_FALSE,&lock,&flg);
812: if (flg) { PEPTOARSetLocking(pep,lock); }
814: PetscOptionsTail();
815: return(0);
816: }
818: PetscErrorCode PEPView_TOAR(PEP pep,PetscViewer viewer)
819: {
821: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
822: PetscBool isascii;
825: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
826: if (isascii) {
827: PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
828: PetscViewerASCIIPrintf(viewer," using the %slocking variant\n",ctx->lock?"":"non-");
829: }
830: return(0);
831: }
833: PetscErrorCode PEPDestroy_TOAR(PEP pep)
834: {
836: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
839: BVDestroy(&ctx->V);
840: PetscFree(pep->data);
841: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",NULL);
842: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",NULL);
843: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",NULL);
844: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",NULL);
845: return(0);
846: }
848: SLEPC_EXTERN PetscErrorCode PEPCreate_TOAR(PEP pep)
849: {
850: PEP_TOAR *ctx;
854: PetscNewLog(pep,&ctx);
855: pep->data = (void*)ctx;
857: pep->lineariz = PETSC_TRUE;
858: ctx->lock = PETSC_TRUE;
860: pep->ops->solve = PEPSolve_TOAR;
861: pep->ops->setup = PEPSetUp_TOAR;
862: pep->ops->setfromoptions = PEPSetFromOptions_TOAR;
863: pep->ops->destroy = PEPDestroy_TOAR;
864: pep->ops->view = PEPView_TOAR;
865: pep->ops->backtransform = PEPBackTransform_Default;
866: pep->ops->computevectors = PEPComputeVectors_Default;
867: pep->ops->extractvectors = PEPExtractVectors_TOAR;
869: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",PEPTOARSetRestart_TOAR);
870: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",PEPTOARGetRestart_TOAR);
871: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",PEPTOARSetLocking_TOAR);
872: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",PEPTOARGetLocking_TOAR);
873: return(0);
874: }