Actual source code: stoar.c
slepc-3.7.0 2016-05-16
1: /*
3: SLEPc polynomial eigensolver: "stoar"
5: Method: S-TOAR
7: Algorithm:
9: Symmetric Two-Level Orthogonal Arnoldi.
11: References:
13: [1] C. Campos and J.E. Roman, "Restarted Q-Arnoldi-type methods
14: exploiting symmetry in quadratic eigenvalue problems", BIT
15: Numer. Math. (in press), 2016.
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/pepimpl.h> /*I "slepcpep.h" I*/
38: #include ../src/pep/impls/krylov/pepkrylov.h
39: #include <slepcblaslapack.h>
41: static PetscBool cited = PETSC_FALSE;
42: static const char citation[] =
43: "@Article{slepc-stoar,\n"
44: " author = \"C. Campos and J. E. Roman\",\n"
45: " title = \"Restarted {Q-Arnoldi-type} methods exploiting symmetry in quadratic eigenvalue problems\",\n"
46: " journal = \"{BIT} Numer. Math.\",\n"
47: " volume = \"to appear\",\n"
48: " number = \"\",\n"
49: " pages = \"\",\n"
50: " year = \"2016,\"\n"
51: " doi = \"http://dx.doi.org/10.1007/s10543-016-0601-5\"\n"
52: "}\n";
56: /*
57: Compute B-norm of v=[v1;v2] whith B=diag(-pep->T[0],pep->T[2])
58: */
59: static PetscErrorCode PEPSTOARNorm(PEP pep,PetscInt j,PetscReal *norm)
60: {
62: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
63: PetscBLASInt n_,one=1,ld_;
64: PetscScalar sone=1.0,szero=0.0,*sp,*sq,*w1,*w2,*qK,*qM;
65: PetscInt n,i,lds=ctx->d*ctx->ld;
68: qK = ctx->qB;
69: qM = ctx->qB+ctx->ld*ctx->ld;
70: n = j+2;
71: PetscMalloc2(n,&w1,n,&w2);
72: sp = ctx->S+lds*j;
73: sq = sp+ctx->ld;
74: PetscBLASIntCast(n,&n_);
75: PetscBLASIntCast(ctx->ld,&ld_);
76: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qK,&ld_,sp,&one,&szero,w1,&one));
77: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qM,&ld_,sq,&one,&szero,w2,&one));
78: *norm = 0.0;
79: for (i=0;i<n;i++) *norm += PetscRealPart(w1[i]*PetscConj(sp[i])+w2[i]*PetscConj(sq[i]));
80: *norm = (*norm>0.0)?PetscSqrtReal(*norm):-PetscSqrtReal(-*norm);
81: PetscFree2(w1,w2);
82: return(0);
83: }
87: static PetscErrorCode PEPSTOARqKqMupdates(PEP pep,PetscInt j,Vec *wv)
88: {
90: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
91: PetscInt i,ld=ctx->ld;
92: PetscScalar *qK,*qM;
93: Vec vj,v1,v2;
96: qK = ctx->qB;
97: qM = ctx->qB+ctx->ld*ctx->ld;
98: v1 = wv[0];
99: v2 = wv[1];
100: BVGetColumn(pep->V,j,&vj);
101: STMatMult(pep->st,0,vj,v1);
102: STMatMult(pep->st,2,vj,v2);
103: BVRestoreColumn(pep->V,j,&vj);
104: for (i=0;i<=j;i++) {
105: BVGetColumn(pep->V,i,&vj);
106: VecDot(v1,vj,qK+j*ld+i);
107: VecDot(v2,vj,qM+j*ld+i);
108: *(qM+j*ld+i) *= pep->sfactor*pep->sfactor;
109: BVRestoreColumn(pep->V,i,&vj);
110: }
111: for (i=0;i<j;i++) {
112: qK[i+j*ld] = -qK[i+ld*j];
113: qK[j+i*ld] = PetscConj(qK[i+j*ld]);
114: qM[j+i*ld] = PetscConj(qM[i+j*ld]);
115: }
116: qK[j+j*ld] = -PetscRealPart(qK[j+ld*j]);
117: qM[j+j*ld] = PetscRealPart(qM[j+ld*j]);
118: return(0);
119: }
123: PetscErrorCode PEPSetUp_STOAR(PEP pep)
124: {
126: PetscBool sinv,flg,lindep;
127: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
128: PetscInt ld,i;
129: PetscReal norm,*omega;
132: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
133: if (!ctx->lock && pep->mpd<pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
134: if (!pep->max_it) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
135: if (!pep->which) {
136: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
137: if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
138: else pep->which = PEP_LARGEST_MAGNITUDE;
139: }
140: if (pep->problem_type!=PEP_HERMITIAN) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
142: if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver only available for quadratic problems");
143: if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
144: STGetTransform(pep->st,&flg);
145: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag set, see STSetTransform()");
147: PEPAllocateSolution(pep,2);
148: PEPSetWorkVecs(pep,4);
149: ld = pep->ncv+2;
150: DSSetType(pep->ds,DSGHIEP);
151: DSSetCompact(pep->ds,PETSC_TRUE);
152: DSAllocate(pep->ds,ld);
153: STGetNumMatrices(pep->st,&ctx->d);
154: ctx->d--;
155: ctx->ld = ld;
156: PetscCalloc1(ctx->d*ld*ld,&ctx->S);
157: PetscCalloc1(2*ld*ld,&ctx->qB);
159: /* process starting vector */
160: if (pep->nini>-2) {
161: BVSetRandomColumn(pep->V,0);
162: BVSetRandomColumn(pep->V,1);
163: } else {
164: BVInsertVec(pep->V,0,pep->IS[0]);
165: BVInsertVec(pep->V,1,pep->IS[1]);
166: }
167: BVOrthogonalizeColumn(pep->V,0,NULL,&norm,&lindep);
168: if (!lindep) {
169: BVScaleColumn(pep->V,0,1.0/norm);
170: ctx->S[0] = norm;
171: PEPSTOARqKqMupdates(pep,0,pep->work);
172: } else SETERRQ(PetscObjectComm((PetscObject)pep),1,"Problem with initial vector");
173: BVOrthogonalizeColumn(pep->V,1,ctx->S+ld,&norm,&lindep);
174: if (!lindep) {
175: BVScaleColumn(pep->V,1,1.0/norm);
176: ctx->S[1] = norm;
177: PEPSTOARqKqMupdates(pep,1,pep->work);
178: } else SETERRQ(PetscObjectComm((PetscObject)pep),1,"Problem with initial vector");
180: PEPSTOARNorm(pep,0,&norm);
181: for (i=0;i<2;i++) { ctx->S[i+ld] /= norm; ctx->S[i] /= norm; }
182: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
183: omega[0] = (norm>0)?1.0:-1.0;
184: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
185: if (pep->nini<0) {
186: SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
187: }
188: return(0);
189: }
193: /*
194: Computes GS orthogonalization x = [z;x] - [Sp;Sq]*y,
195: where y = Omega\([Sp;Sq]'*[qK zeros(size(qK,1)) ;zeros(size(qK,1)) qM]*[z;x]).
196: n: Column from S to be orthogonalized against previous columns.
197: */
198: static PetscErrorCode PEPSTOAROrth2(PEP pep,PetscInt k,PetscReal *Omega,PetscScalar *y)
199: {
201: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
202: PetscBLASInt n_,lds_,k_,one=1,ld_;
203: PetscScalar *S=ctx->S,sonem=-1.0,sone=1.0,szero=0.0,*tp,*tq,*xp,*xq,*c,*qK,*qM;
204: PetscInt i,lds=ctx->d*ctx->ld,n,j;
207: qK = ctx->qB;
208: qM = ctx->qB+ctx->ld*ctx->ld;
209: n = k+2;
210: PetscMalloc3(n,&tp,n,&tq,k,&c);
211: PetscBLASIntCast(n,&n_); /* Size of qK and qM */
212: PetscBLASIntCast(ctx->ld,&ld_);
213: PetscBLASIntCast(lds,&lds_);
214: PetscBLASIntCast(k,&k_); /* Number of vectors to orthogonalize against */
215: xp = S+k*lds;
216: xq = S+ctx->ld+k*lds;
217: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qK,&ld_,xp,&one,&szero,tp,&one));
218: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qM,&ld_,xq,&one,&szero,tq,&one));
219: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,ctx->S,&lds_,tp,&one,&szero,y,&one));
220: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+ctx->ld,&lds_,tq,&one,&sone,y,&one));
221: for (i=0;i<k;i++) y[i] /= Omega[i];
222: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S,&lds_,y,&one,&sone,xp,&one));
223: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+ctx->ld,&lds_,y,&one,&sone,xq,&one));
224: /* three times */
225: for (j=0;j<2;j++) {
226: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qK,&ld_,xp,&one,&szero,tp,&one));
227: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qM,&ld_,xq,&one,&szero,tq,&one));
228: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,ctx->S,&lds_,tp,&one,&szero,c,&one));
229: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+ctx->ld,&lds_,tq,&one,&sone,c,&one));
230: for (i=0;i<k;i++) c[i] /= Omega[i];
231: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S,&lds_,c,&one,&sone,xp,&one));
232: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+ctx->ld,&lds_,c,&one,&sone,xq,&one));
233: for (i=0;i<k;i++) y[i] += c[i];
234: }
235: PetscFree3(tp,tq,c);
236: return(0);
237: }
241: /*
242: Compute a run of Lanczos iterations. dim(work)=(ctx->ld)*4
243: */
244: static PetscErrorCode PEPSTOARrun(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,PetscScalar *work,Vec *t_)
245: {
247: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
248: PetscInt i,j,m=*M,l;
249: PetscInt lds=ctx->d*ctx->ld,offq=ctx->ld;
250: Vec v=t_[0],t=t_[1],q=t_[2];
251: PetscReal norm,sym=0.0,fro=0.0,*f;
252: PetscScalar *y,*S=ctx->S;
253: PetscBLASInt j_,one=1;
254: PetscBool lindep;
257: *breakdown = PETSC_FALSE; /* ----- */
258: DSGetDimensions(pep->ds,NULL,NULL,&l,NULL,NULL);
259: y = work;
260: for (j=k;j<m;j++) {
261: /* apply operator */
262: BVSetActiveColumns(pep->V,0,j+2);
263: BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
264: STMatMult(pep->st,0,v,t);
265: BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
266: STMatMult(pep->st,1,v,q);
267: VecAXPY(t,pep->sfactor,q);
268: STMatSolve(pep->st,t,q);
269: VecScale(q,-1.0/(pep->sfactor*pep->sfactor));
271: /* orthogonalize */
272: BVOrthogonalizeVec(pep->V,q,S+offq+(j+1)*lds,&norm,&lindep);
273: if (lindep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STOAR does not support detection of linearly dependent TOAR vectors");
274: *(S+offq+(j+1)*lds+j+2) = norm;
275: VecScale(q,1.0/norm);
276: BVInsertVec(pep->V,j+2,q);
277: for (i=0;i<=j+1;i++) *(S+(j+1)*lds+i) = *(S+offq+j*lds+i);
279: /* update qK and qM */
280: PEPSTOARqKqMupdates(pep,j+2,t_);
282: /* level-2 orthogonalization */
283: PEPSTOAROrth2(pep,j+1,omega,y);
284: a[j] = PetscRealPart(y[j])/omega[j];
285: PEPSTOARNorm(pep,j+1,&norm);
286: omega[j+1] = (norm > 0)?1.0:-1.0;
287: for (i=0;i<=j+2;i++) {
288: S[i+(j+1)*lds] /= norm;
289: S[i+offq+(j+1)*lds] /= norm;
290: }
291: b[j] = PetscAbsReal(norm);
293: /* check symmetry */
294: DSGetArrayReal(pep->ds,DS_MAT_T,&f);
295: if (j==k) {
296: for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ctx->ld+i]);
297: for (i=0;i<l;i++) y[i] = 0.0;
298: }
299: DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
300: if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
301: PetscBLASIntCast(j,&j_);
302: sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
303: fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
304: if (j>0) fro = SlepcAbs(fro,b[j-1]);
305: if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
306: *symmlost = PETSC_TRUE;
307: *M=j+1;
308: break;
309: }
310: }
311: return(0);
312: }
316: static PetscErrorCode PEPSTOARTrunc(PEP pep,PetscInt rs1,PetscInt cs1,PetscScalar *work,PetscReal *rwork)
317: {
318: #if defined(PETSC_MISSING_LAPACK_GESVD)
320: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable");
321: #else
323: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
324: Mat G;
325: PetscInt lwa,nwu=0,nrwu=0;
326: PetscInt i,n,lds=2*ctx->ld;
327: PetscScalar *M,*V,*U,*S=ctx->S,sone=1.0,zero=0.0,t,*qK,*qM;
328: PetscReal *sg;
329: PetscBLASInt cs1_,rs1_,cs1t2,cs1p1,n_,info,lw_,lds_,ld_;
332: qK = ctx->qB;
333: qM = ctx->qB+ctx->ld*ctx->ld;
334: n = (rs1>2*cs1)?2*cs1:rs1;
335: lwa = cs1*rs1*4+n*(rs1+2*cs1)+(cs1+1)*(cs1+2);
336: M = work+nwu;
337: nwu += rs1*cs1*2;
338: U = work+nwu;
339: nwu += rs1*n;
340: V = work+nwu;
341: nwu += 2*cs1*n;
342: sg = rwork+nrwu;
343: nrwu += n;
344: for (i=0;i<cs1;i++) {
345: PetscMemcpy(M+i*rs1,S+i*lds,rs1*sizeof(PetscScalar));
346: PetscMemcpy(M+(i+cs1)*rs1,S+i*lds+ctx->ld,rs1*sizeof(PetscScalar));
347: }
348: PetscBLASIntCast(n,&n_);
349: PetscBLASIntCast(cs1,&cs1_);
350: PetscBLASIntCast(rs1,&rs1_);
351: PetscBLASIntCast(cs1*2,&cs1t2);
352: PetscBLASIntCast(cs1+1,&cs1p1);
353: PetscBLASIntCast(lds,&lds_);
354: PetscBLASIntCast(ctx->ld,&ld_);
355: PetscBLASIntCast(lwa-nwu,&lw_);
356: #if !defined(PETSC_USE_COMPLEX)
357: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1t2,M,&rs1_,sg,U,&rs1_,V,&n_,work+nwu,&lw_,&info));
358: #else
359: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1t2,M,&rs1_,sg,U,&rs1_,V,&n_,work+nwu,&lw_,rwork+nrwu,&info));
360: #endif
361: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESVD %d",info);
363: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
364: MatCreateSeqDense(PETSC_COMM_SELF,rs1,2*cs1,U,&G);
365: BVSetActiveColumns(pep->V,0,rs1);
366: BVMultInPlace(pep->V,G,0,cs1+1);
367: MatDestroy(&G);
369: /* Update S */
370: PetscMemzero(S,lds*ctx->ld*sizeof(PetscScalar));
372: for (i=0;i<cs1+1;i++) {
373: t = sg[i];
374: PetscStackCallBLAS("BLASscal",BLASscal_(&cs1t2,&t,V+i,&n_));
375: }
376: for (i=0;i<cs1;i++) {
377: PetscMemcpy(S+i*lds,V+i*n,(cs1+1)*sizeof(PetscScalar));
378: PetscMemcpy(S+ctx->ld+i*lds,V+(cs1+i)*n,(cs1+1)*sizeof(PetscScalar));
379: }
381: /* Update qM and qK */
382: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&rs1_,&cs1p1,&rs1_,&sone,qK,&ld_,U,&rs1_,&zero,work+nwu,&rs1_));
383: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&cs1p1,&cs1p1,&rs1_,&sone,U,&rs1_,work+nwu,&rs1_,&zero,qK,&ld_));
384: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&rs1_,&cs1p1,&rs1_,&sone,qM,&ld_,U,&rs1_,&zero,work+nwu,&rs1_));
385: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&cs1p1,&cs1p1,&rs1_,&sone,U,&rs1_,work+nwu,&rs1_,&zero,qM,&ld_));
386: return(0);
387: #endif
388: }
392: /*
393: S <- S*Q
394: columns s-s+ncu of S
395: rows 0-sr of S
396: size(Q) qr x ncu
397: dim(work)=sr*ncu;
398: */
399: static PetscErrorCode PEPSTOARSupdate(PetscScalar *S,PetscInt ld,PetscInt sr,PetscInt s,PetscInt ncu,PetscInt qr,PetscScalar *Q,PetscInt ldq,PetscScalar *work)
400: {
402: PetscScalar a=1.0,b=0.0;
403: PetscBLASInt sr_,ncu_,ldq_,lds_,qr_;
404: PetscInt j,lds=2*ld;
407: PetscBLASIntCast(sr,&sr_);
408: PetscBLASIntCast(qr,&qr_);
409: PetscBLASIntCast(ncu,&ncu_);
410: PetscBLASIntCast(lds,&lds_);
411: PetscBLASIntCast(ldq,&ldq_);
412: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&ncu_,&qr_,&a,S,&lds_,Q,&ldq_,&b,work,&sr_));
413: for (j=0;j<ncu;j++) {
414: PetscMemcpy(S+lds*(s+j),work+j*sr,sr*sizeof(PetscScalar));
415: }
416: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&ncu_,&qr_,&a,S+ld,&lds_,Q,&ldq_,&b,work,&sr_));
417: for (j=0;j<ncu;j++) {
418: PetscMemcpy(S+lds*(s+j)+ld,work+j*sr,sr*sizeof(PetscScalar));
419: }
420: return(0);
421: }
423: #if 0
426: static PetscErrorCode PEPSTOARpreKConvergence(PEP pep,PetscInt nv,PetscReal *norm,Vec *w)
427: {
429: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
430: PetscBLASInt n_,one=1;
431: PetscInt lds=2*ctx->ld;
432: PetscReal t1,t2;
433: PetscScalar *S=ctx->S;
436: PetscBLASIntCast(nv+2,&n_);
437: t1 = BLASnrm2_(&n_,S+nv*2*ctx->ld,&one);
438: t2 = BLASnrm2_(&n_,S+(nv*2+1)*ctx->ld,&one);
439: *norm = SlepcAbs(t1,t2);
440: BVSetActiveColumns(pep->V,0,nv+2);
441: BVMultVec(pep->V,1.0,0.0,w[1],S+nv*lds);
442: STMatMult(pep->st,0,w[1],w[2]);
443: VecNorm(w[2],NORM_2,&t1);
444: BVMultVec(pep->V,1.0,0.0,w[1],S+ctx->ld+nv*lds);
445: STMatMult(pep->st,2,w[1],w[2]);
446: VecNorm(w[2],NORM_2,&t2);
447: t2 *= pep->sfactor*pep->sfactor;
448: *norm = PetscMax(*norm,SlepcAbs(t1,t2));
449: return(0);
450: }
451: #endif
455: PetscErrorCode PEPSolve_STOAR(PEP pep)
456: {
458: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
459: PetscInt j,k,l,nv=0,ld=ctx->ld,lds=ctx->d*ctx->ld,off,ldds,t;
460: PetscInt lwa,lrwa,nwu=0,nrwu=0,nconv=0;
461: PetscScalar *S=ctx->S,*Q,*work;
462: PetscReal beta,norm=1.0,*omega,*a,*b,*r,*rwork;
463: PetscBool breakdown,symmlost=PETSC_FALSE,sinv;
466: PetscCitationsRegister(citation,&cited);
467: BVSetMatrix(pep->V,NULL,PETSC_FALSE);
468: lwa = 9*ld*ld+5*ld;
469: lrwa = 8*ld;
470: PetscMalloc2(lwa,&work,lrwa,&rwork); /* REVIEW */
471: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
472: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
473: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
475: /* Restart loop */
476: l = 0;
477: DSGetLeadingDimension(pep->ds,&ldds);
478: while (pep->reason == PEP_CONVERGED_ITERATING) {
479: pep->its++;
480: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
481: b = a+ldds;
482: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
484: /* Compute an nv-step Lanczos factorization */
485: nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
486: PEPSTOARrun(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,work+nwu,pep->work);
487: beta = b[nv-1];
488: if (symmlost) {
489: pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
490: if (nv==pep->nconv+l+1) { pep->nconv = nconv; break; }
491: }
492: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
493: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
494: DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
495: if (l==0) {
496: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
497: } else {
498: DSSetState(pep->ds,DS_STATE_RAW);
499: }
501: /* Solve projected problem */
502: DSSolve(pep->ds,pep->eigr,pep->eigi);
503: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
505: /* Check convergence */
506: /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
507: norm = 1.0;
508: DSGetDimensions(pep->ds,NULL,NULL,NULL,NULL,&t);
509: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
510: nconv = k;
511: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
513: /* Update l */
514: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
515: else {
516: l = PetscMax(1,(PetscInt)((nv-k)/2));
517: l = PetscMin(l,t);
518: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
519: if (*(a+ldds+k+l-1)!=0) {
520: if (k+l<nv-1) l = l+1;
521: else l = l-1;
522: }
523: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
524: }
525: if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
527: /* Update S */
528: off = pep->nconv*ldds;
529: DSGetArray(pep->ds,DS_MAT_Q,&Q);
530: PEPSTOARSupdate(S,ld,nv+2,pep->nconv,k+l-pep->nconv,nv,Q+off,ldds,work+nwu);
531: DSRestoreArray(pep->ds,DS_MAT_Q,&Q);
533: /* Copy last column of S */
534: PetscMemcpy(S+lds*(k+l),S+lds*nv,lds*sizeof(PetscScalar));
536: if (pep->reason == PEP_CONVERGED_ITERATING) {
537: if (breakdown) {
538: /* Stop if breakdown */
539: PetscInfo2(pep,"Breakdown STOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
540: pep->reason = PEP_DIVERGED_BREAKDOWN;
541: } else {
542: /* Prepare the Rayleigh quotient for restart */
543: DSGetArray(pep->ds,DS_MAT_Q,&Q);
544: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
545: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
546: r = a + 2*ldds;
547: for (j=k;j<k+l;j++) {
548: r[j] = PetscRealPart(Q[nv-1+j*ldds]*beta);
549: }
550: b = a+ldds;
551: b[k+l-1] = r[k+l-1];
552: omega[k+l] = omega[nv];
553: DSRestoreArray(pep->ds,DS_MAT_Q,&Q);
554: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
555: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
556: /* Truncate S */
557: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
558: PEPSTOARTrunc(pep,nv+2,k+l+1,work+nwu,rwork+nrwu);
559: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
560: }
561: }
564: pep->nconv = k;
565: PEPMonitor(pep,pep->its,pep->nconv,pep->eigr,pep->eigi,pep->errest,nv);
566: }
568: if (pep->nconv>0) {
569: /* Truncate S */
570: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
571: PEPSTOARTrunc(pep,nv+2,pep->nconv,work+nwu,rwork+nrwu);
572: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
574: /* Extraction */
575: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
576: DSSetState(pep->ds,DS_STATE_RAW);
578: for (j=0;j<pep->nconv;j++) {
579: pep->eigr[j] *= pep->sfactor;
580: pep->eigi[j] *= pep->sfactor;
581: }
582: }
583: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
584: RGPopScale(pep->rg);
586: /* truncate Schur decomposition and change the state to raw so that
587: DSVectors() computes eigenvectors from scratch */
588: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
589: DSSetState(pep->ds,DS_STATE_RAW);
590: PetscFree2(work,rwork);
591: return(0);
592: }
596: PetscErrorCode PEPSetFromOptions_STOAR(PetscOptionItems *PetscOptionsObject,PEP pep)
597: {
599: PetscBool flg,lock;
602: PetscOptionsHead(PetscOptionsObject,"PEP STOAR Options");
603: PetscOptionsBool("-pep_stoar_locking","Choose between locking and non-locking variants","PEPSTOARSetLocking",PETSC_FALSE,&lock,&flg);
604: if (flg) {
605: PEPSTOARSetLocking(pep,lock);
606: }
607: PetscOptionsTail();
608: return(0);
609: }
613: static PetscErrorCode PEPSTOARSetLocking_STOAR(PEP pep,PetscBool lock)
614: {
615: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
618: ctx->lock = lock;
619: return(0);
620: }
624: /*@
625: PEPSTOARSetLocking - Choose between locking and non-locking variants of
626: the STOAR method.
628: Logically Collective on PEP
630: Input Parameters:
631: + pep - the eigenproblem solver context
632: - lock - true if the locking variant must be selected
634: Options Database Key:
635: . -pep_stoar_locking - Sets the locking flag
637: Notes:
638: The default is to lock converged eigenpairs when the method restarts.
639: This behaviour can be changed so that all directions are kept in the
640: working subspace even if already converged to working accuracy (the
641: non-locking variant).
643: Level: advanced
645: .seealso: PEPSTOARGetLocking()
646: @*/
647: PetscErrorCode PEPSTOARSetLocking(PEP pep,PetscBool lock)
648: {
654: PetscTryMethod(pep,"PEPSTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
655: return(0);
656: }
660: static PetscErrorCode PEPSTOARGetLocking_STOAR(PEP pep,PetscBool *lock)
661: {
662: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
665: *lock = ctx->lock;
666: return(0);
667: }
671: /*@
672: PEPSTOARGetLocking - Gets the locking flag used in the STOAR method.
674: Not Collective
676: Input Parameter:
677: . pep - the eigenproblem solver context
679: Output Parameter:
680: . lock - the locking flag
682: Level: advanced
684: .seealso: PEPSTOARSetLocking()
685: @*/
686: PetscErrorCode PEPSTOARGetLocking(PEP pep,PetscBool *lock)
687: {
693: PetscUseMethod(pep,"PEPSTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
694: return(0);
695: }
699: PetscErrorCode PEPView_STOAR(PEP pep,PetscViewer viewer)
700: {
702: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
703: PetscBool isascii;
706: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
707: if (isascii) {
708: PetscViewerASCIIPrintf(viewer," STOAR: using the %slocking variant\n",ctx->lock?"":"non-");
709: }
710: return(0);
711: }
715: PetscErrorCode PEPDestroy_STOAR(PEP pep)
716: {
720: PetscFree(pep->data);
721: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",NULL);
722: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",NULL);
723: return(0);
724: }
728: PETSC_EXTERN PetscErrorCode PEPCreate_STOAR(PEP pep)
729: {
731: PEP_TOAR *ctx;
734: PetscNewLog(pep,&ctx);
735: pep->data = (void*)ctx;
736: ctx->lock = PETSC_TRUE;
738: pep->ops->solve = PEPSolve_STOAR;
739: pep->ops->setup = PEPSetUp_STOAR;
740: pep->ops->setfromoptions = PEPSetFromOptions_STOAR;
741: pep->ops->view = PEPView_STOAR;
742: pep->ops->destroy = PEPDestroy_STOAR;
743: pep->ops->backtransform = PEPBackTransform_Default;
744: pep->ops->computevectors = PEPComputeVectors_Default;
745: pep->ops->extractvectors = PEPExtractVectors_TOAR;
746: pep->ops->reset = PEPReset_TOAR;
747: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",PEPSTOARSetLocking_STOAR);
748: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",PEPSTOARGetLocking_STOAR);
749: return(0);
750: }