Actual source code: ptoar.c
slepc-3.14.2 2021-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&k_,&k_,T,&k_,p,&info));
332: SlepcCheckLapackInfo("getrf",info);
333: PetscBLASIntCast(sr*k,&lwork);
334: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&k_,T,&k_,p,work,&lwork,&info));
335: SlepcCheckLapackInfo("getri",info);
336: }
337: if (sigma!=0.0) for (i=0;i<k;i++) T[i+k*i] += sigma;
338: } else {
339: T = H; ldt = ldh;
340: }
341: PetscBLASIntCast(ldt,&ldt_);
342: switch (pep->extract) {
343: case PEP_EXTRACT_NONE:
344: break;
345: case PEP_EXTRACT_NORM:
346: if (pep->basis == PEP_BASIS_MONOMIAL) {
347: PetscBLASIntCast(ldt,&ldt_);
348: PetscMalloc1(k,&rwork);
349: norm = LAPACKlange_("F",&k_,&k_,T,&ldt_,rwork);
350: PetscFree(rwork);
351: if (norm>1.0) idxcpy = d-1;
352: } else {
353: PetscBLASIntCast(ldt,&ldt_);
354: PetscMalloc1(k,&rwork);
355: maxnrm = 0.0;
356: for (i=0;i<pep->nmat-1;i++) {
357: PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
358: norm = LAPACKlange_("F",&k_,&k_,Hj,&k_,rwork);
359: if (norm > maxnrm) {
360: idxcpy = i;
361: maxnrm = norm;
362: }
363: }
364: PetscFree(rwork);
365: }
366: if (idxcpy>0) {
367: /* copy block idxcpy of S to the first one */
368: for (j=0;j<k;j++) {
369: PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr);
370: }
371: }
372: break;
373: case PEP_EXTRACT_RESIDUAL:
374: STGetTransform(pep->st,&flg);
375: if (flg) {
376: PetscMalloc1(pep->nmat,&A);
377: for (i=0;i<pep->nmat;i++) {
378: STGetMatrixTransformed(pep->st,i,A+i);
379: }
380: } else A = pep->A;
381: PetscMalloc1(pep->nmat-1,&R);
382: for (i=0;i<pep->nmat-1;i++) {
383: BVDuplicateResize(pep->V,k,R+i);
384: }
385: BVDuplicateResize(pep->V,sr,&Y);
386: MatCreateSeqDense(PETSC_COMM_SELF,sr,k,NULL,&M);
387: g = 0.0; a = 1.0;
388: BVSetActiveColumns(pep->V,0,sr);
389: for (j=0;j<pep->nmat;j++) {
390: BVMatMult(pep->V,A[j],Y);
391: PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
392: for (i=0;i<pep->nmat-1;i++) {
393: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&k_,&k_,&a,S+i*ld,&lds_,Hj,&k_,&g,At,&sr_));
394: MatDenseGetArray(M,&pM);
395: for (jj=0;jj<k;jj++) {
396: PetscArraycpy(pM+jj*sr,At+jj*sr,sr);
397: }
398: MatDenseRestoreArray(M,&pM);
399: BVMult(R[i],1.0,(i==0)?0.0:1.0,Y,M);
400: }
401: }
403: /* frobenius norm */
404: maxnrm = 0.0;
405: for (i=0;i<pep->nmat-1;i++) {
406: BVNorm(R[i],NORM_FROBENIUS,&norm);
407: if (maxnrm > norm) {
408: maxnrm = norm;
409: idxcpy = i;
410: }
411: }
412: if (idxcpy>0) {
413: /* copy block idxcpy of S to the first one */
414: for (j=0;j<k;j++) {
415: PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr);
416: }
417: }
418: if (flg) PetscFree(A);
419: for (i=0;i<pep->nmat-1;i++) {
420: BVDestroy(&R[i]);
421: }
422: PetscFree(R);
423: BVDestroy(&Y);
424: MatDestroy(&M);
425: break;
426: case PEP_EXTRACT_STRUCTURED:
427: for (j=0;j<k;j++) Bt[j+j*k] = 1.0;
428: for (j=0;j<sr;j++) {
429: for (i=0;i<k;i++) At[j*k+i] = PetscConj(S[i*lds+j]);
430: }
431: PEPEvaluateBasisM(pep,k,T,ldt,0,&Hp,&Hj);
432: for (i=1;i<deg;i++) {
433: PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
434: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&sr_,&k_,&sone,Hj,&k_,S+i*ld,&lds_,&sone,At,&k_));
435: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&k_,&k_,&sone,Hj,&k_,Hj,&k_,&sone,Bt,&k_));
436: }
437: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&k_,&sr_,Bt,&k_,p,At,&k_,&info));
438: SlepcCheckLapackInfo("gesv",info);
439: for (j=0;j<sr;j++) {
440: for (i=0;i<k;i++) S[i*lds+j] = PetscConj(At[j*k+i]);
441: }
442: break;
443: }
444: if (transf) { PetscFree(T); }
445: PetscFree6(p,At,Bt,Hj,Hp,work);
446: return(0);
447: }
449: PetscErrorCode PEPSolve_TOAR(PEP pep)
450: {
452: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
453: PetscInt i,j,k,l,nv=0,ld,lds,ldds,newn,nq=0,nconv=0;
454: PetscInt nmat=pep->nmat,deg=nmat-1;
455: PetscScalar *S,*H,sigma;
456: PetscReal beta;
457: PetscBool breakdown=PETSC_FALSE,flg,falselock=PETSC_FALSE,sinv=PETSC_FALSE;
458: Mat MS,MQ;
461: PetscCitationsRegister(citation,&cited);
462: if (ctx->lock) {
463: /* undocumented option to use a cheaper locking instead of the true locking */
464: PetscOptionsGetBool(NULL,NULL,"-pep_toar_falselocking",&falselock,NULL);
465: }
466: DSGetLeadingDimension(pep->ds,&ldds);
467: STGetShift(pep->st,&sigma);
469: /* update polynomial basis coefficients */
470: STGetTransform(pep->st,&flg);
471: if (pep->sfactor!=1.0) {
472: for (i=0;i<nmat;i++) {
473: pep->pbc[nmat+i] /= pep->sfactor;
474: pep->pbc[2*nmat+i] /= pep->sfactor*pep->sfactor;
475: }
476: if (!flg) {
477: pep->target /= pep->sfactor;
478: RGPushScale(pep->rg,1.0/pep->sfactor);
479: STScaleShift(pep->st,1.0/pep->sfactor);
480: sigma /= pep->sfactor;
481: } else {
482: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
483: pep->target = sinv?pep->target*pep->sfactor:pep->target/pep->sfactor;
484: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
485: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
486: }
487: }
489: if (flg) sigma = 0.0;
491: /* clean projected matrix (including the extra-arrow) */
492: DSGetArray(pep->ds,DS_MAT_A,&H);
493: PetscArrayzero(H,ldds*ldds);
494: DSRestoreArray(pep->ds,DS_MAT_A,&H);
496: /* Get the starting Arnoldi vector */
497: BVTensorBuildFirstColumn(ctx->V,pep->nini);
499: /* restart loop */
500: l = 0;
501: while (pep->reason == PEP_CONVERGED_ITERATING) {
502: pep->its++;
504: /* compute an nv-step Lanczos factorization */
505: nv = PetscMax(PetscMin(nconv+pep->mpd,pep->ncv),nv);
506: DSGetArray(pep->ds,DS_MAT_A,&H);
507: PEPTOARrun(pep,sigma,H,ldds,pep->nconv+l,&nv,&breakdown,pep->work);
508: beta = PetscAbsScalar(H[(nv-1)*ldds+nv]);
509: DSRestoreArray(pep->ds,DS_MAT_A,&H);
510: DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
511: if (l==0) {
512: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
513: } else {
514: DSSetState(pep->ds,DS_STATE_RAW);
515: }
516: BVSetActiveColumns(ctx->V,pep->nconv,nv);
518: /* solve projected problem */
519: DSSolve(pep->ds,pep->eigr,pep->eigi);
520: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
521: DSUpdateExtraRow(pep->ds);
522: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
524: /* check convergence */
525: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
526: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
528: /* update l */
529: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
530: else {
531: l = (nv==k)?0:PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
532: if (!breakdown) {
533: /* prepare the Rayleigh quotient for restart */
534: DSTruncate(pep->ds,k+l);
535: DSGetDimensions(pep->ds,&newn,NULL,NULL,NULL,NULL);
536: l = newn-k;
537: }
538: }
539: nconv = k;
540: if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
542: /* update S */
543: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
544: BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
545: MatDestroy(&MQ);
547: /* copy last column of S */
548: BVCopyColumn(ctx->V,nv,k+l);
550: if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
551: /* stop if breakdown */
552: PetscInfo2(pep,"Breakdown TOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
553: pep->reason = PEP_DIVERGED_BREAKDOWN;
554: }
555: if (pep->reason != PEP_CONVERGED_ITERATING) l--;
556: /* truncate S */
557: BVGetActiveColumns(pep->V,NULL,&nq);
558: if (k+l+deg<=nq) {
559: BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
560: if (!falselock && ctx->lock) {
561: BVTensorCompress(ctx->V,k-pep->nconv);
562: } else {
563: BVTensorCompress(ctx->V,0);
564: }
565: }
566: pep->nconv = k;
567: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
568: }
569: if (pep->nconv>0) {
570: /* {V*S_nconv^i}_{i=0}^{d-1} has rank nconv instead of nconv+d-1. Force zeros in each S_nconv^i block */
571: BVSetActiveColumns(ctx->V,0,pep->nconv);
572: BVGetActiveColumns(pep->V,NULL,&nq);
573: BVSetActiveColumns(pep->V,0,nq);
574: if (nq>pep->nconv) {
575: BVTensorCompress(ctx->V,pep->nconv);
576: BVSetActiveColumns(pep->V,0,pep->nconv);
577: nq = pep->nconv;
578: }
580: /* perform Newton refinement if required */
581: if (pep->refine==PEP_REFINE_MULTIPLE && pep->rits>0) {
582: /* extract invariant pair */
583: BVTensorGetFactors(ctx->V,NULL,&MS);
584: MatDenseGetArray(MS,&S);
585: DSGetArray(pep->ds,DS_MAT_A,&H);
586: BVGetSizes(pep->V,NULL,NULL,&ld);
587: lds = deg*ld;
588: PEPExtractInvariantPair(pep,sigma,nq,pep->nconv,S,ld,deg,H,ldds);
589: DSRestoreArray(pep->ds,DS_MAT_A,&H);
590: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
591: DSSetState(pep->ds,DS_STATE_RAW);
592: PEPNewtonRefinement_TOAR(pep,sigma,&pep->rits,NULL,pep->nconv,S,lds);
593: DSSolve(pep->ds,pep->eigr,pep->eigi);
594: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
595: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
596: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
597: BVMultInPlace(ctx->V,MQ,0,pep->nconv);
598: MatDestroy(&MQ);
599: MatDenseRestoreArray(MS,&S);
600: BVTensorRestoreFactors(ctx->V,NULL,&MS);
601: }
602: }
603: STGetTransform(pep->st,&flg);
604: if (pep->refine!=PEP_REFINE_MULTIPLE || pep->rits==0) {
605: if (!flg && pep->ops->backtransform) {
606: (*pep->ops->backtransform)(pep);
607: }
608: if (pep->sfactor!=1.0) {
609: for (j=0;j<pep->nconv;j++) {
610: pep->eigr[j] *= pep->sfactor;
611: pep->eigi[j] *= pep->sfactor;
612: }
613: /* restore original values */
614: for (i=0;i<pep->nmat;i++){
615: pep->pbc[pep->nmat+i] *= pep->sfactor;
616: pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
617: }
618: }
619: }
620: /* restore original values */
621: if (!flg) {
622: pep->target *= pep->sfactor;
623: STScaleShift(pep->st,pep->sfactor);
624: } else {
625: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
626: pep->target = (sinv)?pep->target/pep->sfactor:pep->target*pep->sfactor;
627: }
628: if (pep->sfactor!=1.0) { RGPopScale(pep->rg); }
630: /* change the state to raw so that DSVectors() computes eigenvectors from scratch */
631: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
632: DSSetState(pep->ds,DS_STATE_RAW);
633: return(0);
634: }
636: static PetscErrorCode PEPTOARSetRestart_TOAR(PEP pep,PetscReal keep)
637: {
638: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
641: if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
642: else {
643: 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]");
644: ctx->keep = keep;
645: }
646: return(0);
647: }
649: /*@
650: PEPTOARSetRestart - Sets the restart parameter for the TOAR
651: method, in particular the proportion of basis vectors that must be kept
652: after restart.
654: Logically Collective on pep
656: Input Parameters:
657: + pep - the eigenproblem solver context
658: - keep - the number of vectors to be kept at restart
660: Options Database Key:
661: . -pep_toar_restart - Sets the restart parameter
663: Notes:
664: Allowed values are in the range [0.1,0.9]. The default is 0.5.
666: Level: advanced
668: .seealso: PEPTOARGetRestart()
669: @*/
670: PetscErrorCode PEPTOARSetRestart(PEP pep,PetscReal keep)
671: {
677: PetscTryMethod(pep,"PEPTOARSetRestart_C",(PEP,PetscReal),(pep,keep));
678: return(0);
679: }
681: static PetscErrorCode PEPTOARGetRestart_TOAR(PEP pep,PetscReal *keep)
682: {
683: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
686: *keep = ctx->keep;
687: return(0);
688: }
690: /*@
691: PEPTOARGetRestart - Gets the restart parameter used in the TOAR method.
693: Not Collective
695: Input Parameter:
696: . pep - the eigenproblem solver context
698: Output Parameter:
699: . keep - the restart parameter
701: Level: advanced
703: .seealso: PEPTOARSetRestart()
704: @*/
705: PetscErrorCode PEPTOARGetRestart(PEP pep,PetscReal *keep)
706: {
712: PetscUseMethod(pep,"PEPTOARGetRestart_C",(PEP,PetscReal*),(pep,keep));
713: return(0);
714: }
716: static PetscErrorCode PEPTOARSetLocking_TOAR(PEP pep,PetscBool lock)
717: {
718: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
721: ctx->lock = lock;
722: return(0);
723: }
725: /*@
726: PEPTOARSetLocking - Choose between locking and non-locking variants of
727: the TOAR method.
729: Logically Collective on pep
731: Input Parameters:
732: + pep - the eigenproblem solver context
733: - lock - true if the locking variant must be selected
735: Options Database Key:
736: . -pep_toar_locking - Sets the locking flag
738: Notes:
739: The default is to lock converged eigenpairs when the method restarts.
740: This behaviour can be changed so that all directions are kept in the
741: working subspace even if already converged to working accuracy (the
742: non-locking variant).
744: Level: advanced
746: .seealso: PEPTOARGetLocking()
747: @*/
748: PetscErrorCode PEPTOARSetLocking(PEP pep,PetscBool lock)
749: {
755: PetscTryMethod(pep,"PEPTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
756: return(0);
757: }
759: static PetscErrorCode PEPTOARGetLocking_TOAR(PEP pep,PetscBool *lock)
760: {
761: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
764: *lock = ctx->lock;
765: return(0);
766: }
768: /*@
769: PEPTOARGetLocking - Gets the locking flag used in the TOAR method.
771: Not Collective
773: Input Parameter:
774: . pep - the eigenproblem solver context
776: Output Parameter:
777: . lock - the locking flag
779: Level: advanced
781: .seealso: PEPTOARSetLocking()
782: @*/
783: PetscErrorCode PEPTOARGetLocking(PEP pep,PetscBool *lock)
784: {
790: PetscUseMethod(pep,"PEPTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
791: return(0);
792: }
794: PetscErrorCode PEPSetFromOptions_TOAR(PetscOptionItems *PetscOptionsObject,PEP pep)
795: {
797: PetscBool flg,lock;
798: PetscReal keep;
801: PetscOptionsHead(PetscOptionsObject,"PEP TOAR Options");
803: PetscOptionsReal("-pep_toar_restart","Proportion of vectors kept after restart","PEPTOARSetRestart",0.5,&keep,&flg);
804: if (flg) { PEPTOARSetRestart(pep,keep); }
806: PetscOptionsBool("-pep_toar_locking","Choose between locking and non-locking variants","PEPTOARSetLocking",PETSC_FALSE,&lock,&flg);
807: if (flg) { PEPTOARSetLocking(pep,lock); }
809: PetscOptionsTail();
810: return(0);
811: }
813: PetscErrorCode PEPView_TOAR(PEP pep,PetscViewer viewer)
814: {
816: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
817: PetscBool isascii;
820: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
821: if (isascii) {
822: PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
823: PetscViewerASCIIPrintf(viewer," using the %slocking variant\n",ctx->lock?"":"non-");
824: }
825: return(0);
826: }
828: PetscErrorCode PEPDestroy_TOAR(PEP pep)
829: {
831: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
834: BVDestroy(&ctx->V);
835: PetscFree(pep->data);
836: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",NULL);
837: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",NULL);
838: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",NULL);
839: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",NULL);
840: return(0);
841: }
843: SLEPC_EXTERN PetscErrorCode PEPCreate_TOAR(PEP pep)
844: {
845: PEP_TOAR *ctx;
849: PetscNewLog(pep,&ctx);
850: pep->data = (void*)ctx;
852: pep->lineariz = PETSC_TRUE;
853: ctx->lock = PETSC_TRUE;
855: pep->ops->solve = PEPSolve_TOAR;
856: pep->ops->setup = PEPSetUp_TOAR;
857: pep->ops->setfromoptions = PEPSetFromOptions_TOAR;
858: pep->ops->destroy = PEPDestroy_TOAR;
859: pep->ops->view = PEPView_TOAR;
860: pep->ops->backtransform = PEPBackTransform_Default;
861: pep->ops->computevectors = PEPComputeVectors_Default;
862: pep->ops->extractvectors = PEPExtractVectors_TOAR;
864: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",PEPTOARSetRestart_TOAR);
865: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",PEPTOARGetRestart_TOAR);
866: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",PEPTOARSetLocking_TOAR);
867: PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",PEPTOARGetLocking_TOAR);
868: return(0);
869: }