Actual source code: stoar.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  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      shift,sinv,flg,lindep;
127:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
128:   PetscInt       ld,i;
129:   PetscReal      norm,*omega;

132:   pep->lineariz = PETSC_TRUE;
133:   PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
134:   if (!ctx->lock && pep->mpd<pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
135:   if (!pep->max_it) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
136:   /* Set STSHIFT as the default ST */
137:   if (!((PetscObject)pep->st)->type_name) {
138:     STSetType(pep->st,STSHIFT);
139:   }
140:   PetscObjectTypeCompare((PetscObject)pep->st,STSHIFT,&shift);
141:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
142:   if (!shift && !sinv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Only STSHIFT and STSINVERT spectral transformations can be used");
143:   if (!pep->which) {
144:     if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
145:     else pep->which = PEP_LARGEST_MAGNITUDE;
146:   }
147:   if (pep->problem_type!=PEP_HERMITIAN) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");

149:   if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver only available for quadratic problems");
150:   if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
151:   STGetTransform(pep->st,&flg);
152:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag set, see STSetTransform()");

154:   PEPAllocateSolution(pep,2);
155:   PEPSetWorkVecs(pep,4);
156:   ld = pep->ncv+2;
157:   DSSetType(pep->ds,DSGHIEP);
158:   DSSetCompact(pep->ds,PETSC_TRUE);
159:   DSAllocate(pep->ds,ld);
160:   STGetNumMatrices(pep->st,&ctx->d);
161:   ctx->d--;
162:   ctx->ld = ld;
163:   PetscCalloc1(ctx->d*ld*ld,&ctx->S);
164:   PetscCalloc1(2*ld*ld,&ctx->qB);

166:   /* process starting vector */
167:   if (pep->nini>-2) {
168:     BVSetRandomColumn(pep->V,0);
169:     BVSetRandomColumn(pep->V,1);
170:   } else {
171:     BVInsertVec(pep->V,0,pep->IS[0]);
172:     BVInsertVec(pep->V,1,pep->IS[1]);
173:   }
174:   BVOrthogonalizeColumn(pep->V,0,NULL,&norm,&lindep);
175:   if (!lindep) {
176:     BVScaleColumn(pep->V,0,1.0/norm);
177:     ctx->S[0] = norm;
178:     PEPSTOARqKqMupdates(pep,0,pep->work);
179:   } else SETERRQ(PetscObjectComm((PetscObject)pep),1,"Problem with initial vector");
180:   BVOrthogonalizeColumn(pep->V,1,ctx->S+ld,&norm,&lindep);
181:   if (!lindep) {
182:     BVScaleColumn(pep->V,1,1.0/norm);
183:     ctx->S[1] = norm;
184:     PEPSTOARqKqMupdates(pep,1,pep->work);
185:   } else SETERRQ(PetscObjectComm((PetscObject)pep),1,"Problem with initial vector");

187:   PEPSTOARNorm(pep,0,&norm);
188:   for (i=0;i<2;i++) { ctx->S[i+ld] /= norm; ctx->S[i] /= norm; }
189:   DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
190:   omega[0] = (norm>0)?1.0:-1.0;
191:   DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
192:   if (pep->nini<0) {
193:     SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
194:   }
195:   return(0);
196: }

200: /*
201:   Computes GS orthogonalization  x = [z;x] - [Sp;Sq]*y,
202:   where y = Omega\([Sp;Sq]'*[qK zeros(size(qK,1)) ;zeros(size(qK,1)) qM]*[z;x]).
203:   n: Column from S to be orthogonalized against previous columns.
204: */
205: static PetscErrorCode PEPSTOAROrth2(PEP pep,PetscInt k,PetscReal *Omega,PetscScalar *y)
206: {
208:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
209:   PetscBLASInt   n_,lds_,k_,one=1,ld_;
210:   PetscScalar    *S=ctx->S,sonem=-1.0,sone=1.0,szero=0.0,*tp,*tq,*xp,*xq,*c,*qK,*qM;
211:   PetscInt       i,lds=ctx->d*ctx->ld,n,j;

214:   qK = ctx->qB;
215:   qM = ctx->qB+ctx->ld*ctx->ld;
216:   n = k+2;
217:   PetscMalloc3(n,&tp,n,&tq,k,&c);
218:   PetscBLASIntCast(n,&n_); /* Size of qK and qM */
219:   PetscBLASIntCast(ctx->ld,&ld_);
220:   PetscBLASIntCast(lds,&lds_);
221:   PetscBLASIntCast(k,&k_); /* Number of vectors to orthogonalize against */
222:   xp = S+k*lds;
223:   xq = S+ctx->ld+k*lds;
224:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qK,&ld_,xp,&one,&szero,tp,&one));
225:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qM,&ld_,xq,&one,&szero,tq,&one));
226:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,ctx->S,&lds_,tp,&one,&szero,y,&one));
227:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+ctx->ld,&lds_,tq,&one,&sone,y,&one));
228:   for (i=0;i<k;i++) y[i] /= Omega[i];
229:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S,&lds_,y,&one,&sone,xp,&one));
230:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+ctx->ld,&lds_,y,&one,&sone,xq,&one));
231:   /* three times */
232:   for (j=0;j<2;j++) {
233:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qK,&ld_,xp,&one,&szero,tp,&one));
234:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qM,&ld_,xq,&one,&szero,tq,&one));
235:     PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,ctx->S,&lds_,tp,&one,&szero,c,&one));
236:     PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+ctx->ld,&lds_,tq,&one,&sone,c,&one));
237:     for (i=0;i<k;i++) c[i] /= Omega[i];
238:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S,&lds_,c,&one,&sone,xp,&one));
239:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+ctx->ld,&lds_,c,&one,&sone,xq,&one));
240:     for (i=0;i<k;i++) y[i] += c[i];
241:   }
242:   PetscFree3(tp,tq,c);
243:   return(0);
244: }

248: /*
249:   Compute a run of Lanczos iterations. dim(work)=(ctx->ld)*4
250: */
251: static PetscErrorCode PEPSTOARrun(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,PetscScalar *work,Vec *t_)
252: {
254:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
255:   PetscInt       i,j,m=*M,l;
256:   PetscInt       lds=ctx->d*ctx->ld,offq=ctx->ld;
257:   Vec            v=t_[0],t=t_[1],q=t_[2];
258:   PetscReal      norm,sym=0.0,fro=0.0,*f;
259:   PetscScalar    *y,*S=ctx->S;
260:   PetscBLASInt   j_,one=1;
261:   PetscBool      lindep;

264:   *breakdown = PETSC_FALSE; /* ----- */
265:   DSGetDimensions(pep->ds,NULL,NULL,&l,NULL,NULL);
266:   y = work;
267:   for (j=k;j<m;j++) {
268:     /* apply operator */
269:     BVSetActiveColumns(pep->V,0,j+2);
270:     BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
271:     STMatMult(pep->st,0,v,t);
272:     BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
273:     STMatMult(pep->st,1,v,q);
274:     VecAXPY(t,pep->sfactor,q);
275:     STMatSolve(pep->st,t,q);
276:     VecScale(q,-1.0/(pep->sfactor*pep->sfactor));

278:     /* orthogonalize */
279:     BVOrthogonalizeVec(pep->V,q,S+offq+(j+1)*lds,&norm,&lindep);
280:     if (lindep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STOAR does not support detection of linearly dependent TOAR vectors");
281:     *(S+offq+(j+1)*lds+j+2) = norm;
282:     VecScale(q,1.0/norm);
283:     BVInsertVec(pep->V,j+2,q);
284:     for (i=0;i<=j+1;i++) *(S+(j+1)*lds+i) = *(S+offq+j*lds+i);

286:     /* update qK and qM */
287:     PEPSTOARqKqMupdates(pep,j+2,t_);

289:     /* level-2 orthogonalization */
290:     PEPSTOAROrth2(pep,j+1,omega,y);
291:     a[j] = PetscRealPart(y[j])/omega[j];
292:     PEPSTOARNorm(pep,j+1,&norm);
293:     omega[j+1] = (norm > 0)?1.0:-1.0;
294:     for (i=0;i<=j+2;i++) {
295:       S[i+(j+1)*lds] /= norm;
296:       S[i+offq+(j+1)*lds] /= norm;
297:     }
298:     b[j] = PetscAbsReal(norm);

300:     /* check symmetry */
301:     DSGetArrayReal(pep->ds,DS_MAT_T,&f);
302:     if (j==k) {
303:       for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ctx->ld+i]);
304:       for (i=0;i<l;i++) y[i] = 0.0;
305:     }
306:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
307:     if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
308:     PetscBLASIntCast(j,&j_);
309:     sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
310:     fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
311:     if (j>0) fro = SlepcAbs(fro,b[j-1]);
312:     if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
313:       *symmlost = PETSC_TRUE;
314:       *M=j+1;
315:       break;
316:     }
317:   }
318:   return(0);
319: }

323: static PetscErrorCode PEPSTOARTrunc(PEP pep,PetscInt rs1,PetscInt cs1,PetscScalar *work,PetscReal *rwork)
324: {
325: #if defined(PETSC_MISSING_LAPACK_GESVD)
327:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable");
328: #else
330:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
331:   Mat            G;
332:   PetscInt       lwa,nwu=0,nrwu=0;
333:   PetscInt       i,n,lds=2*ctx->ld;
334:   PetscScalar    *M,*V,*U,*S=ctx->S,sone=1.0,zero=0.0,t,*qK,*qM;
335:   PetscReal      *sg;
336:   PetscBLASInt   cs1_,rs1_,cs1t2,cs1p1,n_,info,lw_,lds_,ld_;

339:   qK = ctx->qB;
340:   qM = ctx->qB+ctx->ld*ctx->ld;
341:   n = (rs1>2*cs1)?2*cs1:rs1;
342:   lwa = cs1*rs1*4+n*(rs1+2*cs1)+(cs1+1)*(cs1+2);
343:   M = work+nwu;
344:   nwu += rs1*cs1*2;
345:   U = work+nwu;
346:   nwu += rs1*n;
347:   V = work+nwu;
348:   nwu += 2*cs1*n;
349:   sg = rwork+nrwu;
350:   nrwu += n;
351:   for (i=0;i<cs1;i++) {
352:     PetscMemcpy(M+i*rs1,S+i*lds,rs1*sizeof(PetscScalar));
353:     PetscMemcpy(M+(i+cs1)*rs1,S+i*lds+ctx->ld,rs1*sizeof(PetscScalar));
354:   }
355:   PetscBLASIntCast(n,&n_);
356:   PetscBLASIntCast(cs1,&cs1_);
357:   PetscBLASIntCast(rs1,&rs1_);
358:   PetscBLASIntCast(cs1*2,&cs1t2);
359:   PetscBLASIntCast(cs1+1,&cs1p1);
360:   PetscBLASIntCast(lds,&lds_);
361:   PetscBLASIntCast(ctx->ld,&ld_);
362:   PetscBLASIntCast(lwa-nwu,&lw_);
363: #if !defined(PETSC_USE_COMPLEX)
364:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1t2,M,&rs1_,sg,U,&rs1_,V,&n_,work+nwu,&lw_,&info));
365: #else
366:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1t2,M,&rs1_,sg,U,&rs1_,V,&n_,work+nwu,&lw_,rwork+nrwu,&info));
367: #endif
368:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESVD %d",info);

370:   /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
371:   MatCreateSeqDense(PETSC_COMM_SELF,rs1,2*cs1,U,&G);
372:   BVSetActiveColumns(pep->V,0,rs1);
373:   BVMultInPlace(pep->V,G,0,cs1+1);
374:   MatDestroy(&G);

376:   /* Update S */
377:   PetscMemzero(S,lds*ctx->ld*sizeof(PetscScalar));

379:   for (i=0;i<cs1+1;i++) {
380:     t = sg[i];
381:     PetscStackCallBLAS("BLASscal",BLASscal_(&cs1t2,&t,V+i,&n_));
382:   }
383:   for (i=0;i<cs1;i++) {
384:     PetscMemcpy(S+i*lds,V+i*n,(cs1+1)*sizeof(PetscScalar));
385:     PetscMemcpy(S+ctx->ld+i*lds,V+(cs1+i)*n,(cs1+1)*sizeof(PetscScalar));
386:   }

388:   /* Update qM and qK */
389:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&rs1_,&cs1p1,&rs1_,&sone,qK,&ld_,U,&rs1_,&zero,work+nwu,&rs1_));
390:   PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&cs1p1,&cs1p1,&rs1_,&sone,U,&rs1_,work+nwu,&rs1_,&zero,qK,&ld_));
391:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&rs1_,&cs1p1,&rs1_,&sone,qM,&ld_,U,&rs1_,&zero,work+nwu,&rs1_));
392:   PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&cs1p1,&cs1p1,&rs1_,&sone,U,&rs1_,work+nwu,&rs1_,&zero,qM,&ld_));
393:   return(0);
394: #endif
395: }

399: /*
400:   S <- S*Q
401:   columns s-s+ncu of S
402:   rows 0-sr of S
403:   size(Q) qr x ncu
404:   dim(work)=sr*ncu;
405: */
406: static PetscErrorCode PEPSTOARSupdate(PetscScalar *S,PetscInt ld,PetscInt sr,PetscInt s,PetscInt ncu,PetscInt qr,PetscScalar *Q,PetscInt ldq,PetscScalar *work)
407: {
409:   PetscScalar    a=1.0,b=0.0;
410:   PetscBLASInt   sr_,ncu_,ldq_,lds_,qr_;
411:   PetscInt       j,lds=2*ld;

414:   PetscBLASIntCast(sr,&sr_);
415:   PetscBLASIntCast(qr,&qr_);
416:   PetscBLASIntCast(ncu,&ncu_);
417:   PetscBLASIntCast(lds,&lds_);
418:   PetscBLASIntCast(ldq,&ldq_);
419:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&ncu_,&qr_,&a,S,&lds_,Q,&ldq_,&b,work,&sr_));
420:   for (j=0;j<ncu;j++) {
421:     PetscMemcpy(S+lds*(s+j),work+j*sr,sr*sizeof(PetscScalar));
422:   }
423:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&ncu_,&qr_,&a,S+ld,&lds_,Q,&ldq_,&b,work,&sr_));
424:   for (j=0;j<ncu;j++) {
425:     PetscMemcpy(S+lds*(s+j)+ld,work+j*sr,sr*sizeof(PetscScalar));
426:   }
427:   return(0);
428: }

430: #if 0
433: static PetscErrorCode PEPSTOARpreKConvergence(PEP pep,PetscInt nv,PetscReal *norm,Vec *w)
434: {
436:   PEP_TOAR      *ctx = (PEP_TOAR*)pep->data;
437:   PetscBLASInt   n_,one=1;
438:   PetscInt       lds=2*ctx->ld;
439:   PetscReal      t1,t2;
440:   PetscScalar    *S=ctx->S;

443:   PetscBLASIntCast(nv+2,&n_);
444:   t1 = BLASnrm2_(&n_,S+nv*2*ctx->ld,&one);
445:   t2 = BLASnrm2_(&n_,S+(nv*2+1)*ctx->ld,&one);
446:   *norm = SlepcAbs(t1,t2);
447:   BVSetActiveColumns(pep->V,0,nv+2);
448:   BVMultVec(pep->V,1.0,0.0,w[1],S+nv*lds);
449:   STMatMult(pep->st,0,w[1],w[2]);
450:   VecNorm(w[2],NORM_2,&t1);
451:   BVMultVec(pep->V,1.0,0.0,w[1],S+ctx->ld+nv*lds);
452:   STMatMult(pep->st,2,w[1],w[2]);
453:   VecNorm(w[2],NORM_2,&t2);
454:   t2 *= pep->sfactor*pep->sfactor;
455:   *norm = PetscMax(*norm,SlepcAbs(t1,t2));
456:   return(0);
457: }
458: #endif

462: PetscErrorCode PEPSolve_STOAR(PEP pep)
463: {
465:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
466:   PetscInt       j,k,l,nv=0,ld=ctx->ld,lds=ctx->d*ctx->ld,off,ldds,t;
467:   PetscInt       lwa,lrwa,nwu=0,nrwu=0,nconv=0;
468:   PetscScalar    *S=ctx->S,*Q,*work;
469:   PetscReal      beta,norm=1.0,*omega,*a,*b,*r,*rwork;
470:   PetscBool      breakdown,symmlost=PETSC_FALSE,sinv;

473:   PetscCitationsRegister(citation,&cited);
474:   BVSetMatrix(pep->V,NULL,PETSC_FALSE);
475:   lwa = 9*ld*ld+5*ld;
476:   lrwa = 8*ld;
477:   PetscMalloc2(lwa,&work,lrwa,&rwork); /* REVIEW */
478:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
479:   RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
480:   STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);

482:   /* Restart loop */
483:   l = 0;
484:   DSGetLeadingDimension(pep->ds,&ldds);
485:   while (pep->reason == PEP_CONVERGED_ITERATING) {
486:     pep->its++;
487:     DSGetArrayReal(pep->ds,DS_MAT_T,&a);
488:     b = a+ldds;
489:     DSGetArrayReal(pep->ds,DS_MAT_D,&omega);

491:     /* Compute an nv-step Lanczos factorization */
492:     nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
493:     PEPSTOARrun(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,work+nwu,pep->work);
494:     beta = b[nv-1];
495:     if (symmlost) {
496:       pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
497:       if (nv==pep->nconv+l+1) { pep->nconv = nconv; break; }
498:     }
499:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
500:     DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
501:     DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
502:     if (l==0) {
503:       DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
504:     } else {
505:       DSSetState(pep->ds,DS_STATE_RAW);
506:     }

508:     /* Solve projected problem */
509:     DSSolve(pep->ds,pep->eigr,pep->eigi);
510:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);

512:     /* Check convergence */
513:     /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
514:     norm = 1.0;
515:     DSGetDimensions(pep->ds,NULL,NULL,NULL,NULL,&t);
516:     PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
517:     nconv = k;
518:     (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);

520:     /* Update l */
521:     if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
522:     else {
523:       l = PetscMax(1,(PetscInt)((nv-k)/2));
524:       l = PetscMin(l,t);
525:       DSGetArrayReal(pep->ds,DS_MAT_T,&a);
526:       if (*(a+ldds+k+l-1)!=0) {
527:         if (k+l<nv-1) l = l+1;
528:         else l = l-1;
529:       }
530:       DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
531:     }
532:     if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */

534:     /* Update S */
535:     off = pep->nconv*ldds;
536:     DSGetArray(pep->ds,DS_MAT_Q,&Q);
537:     PEPSTOARSupdate(S,ld,nv+2,pep->nconv,k+l-pep->nconv,nv,Q+off,ldds,work+nwu);
538:     DSRestoreArray(pep->ds,DS_MAT_Q,&Q);

540:     /* Copy last column of S */
541:     PetscMemcpy(S+lds*(k+l),S+lds*nv,lds*sizeof(PetscScalar));

543:     if (pep->reason == PEP_CONVERGED_ITERATING) {
544:       if (breakdown) {
545:         /* Stop if breakdown */
546:         PetscInfo2(pep,"Breakdown STOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
547:         pep->reason = PEP_DIVERGED_BREAKDOWN;
548:       } else {
549:         /* Prepare the Rayleigh quotient for restart */
550:         DSGetArray(pep->ds,DS_MAT_Q,&Q);
551:         DSGetArrayReal(pep->ds,DS_MAT_T,&a);
552:         DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
553:         r = a + 2*ldds;
554:         for (j=k;j<k+l;j++) {
555:           r[j] = PetscRealPart(Q[nv-1+j*ldds]*beta);
556:         }
557:         b = a+ldds;
558:         b[k+l-1] = r[k+l-1];
559:         omega[k+l] = omega[nv];
560:         DSRestoreArray(pep->ds,DS_MAT_Q,&Q);
561:         DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
562:         DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
563:         /* Truncate S */
564:         DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
565:         PEPSTOARTrunc(pep,nv+2,k+l+1,work+nwu,rwork+nrwu);
566:         DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
567:       }
568:     }


571:     pep->nconv = k;
572:     PEPMonitor(pep,pep->its,pep->nconv,pep->eigr,pep->eigi,pep->errest,nv);
573:   }

575:   if (pep->nconv>0) {
576:     /* Truncate S */
577:     DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
578:     PEPSTOARTrunc(pep,nv+2,pep->nconv,work+nwu,rwork+nrwu);
579:     DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);

581:     /* Extraction */
582:     DSSetDimensions(pep->ds,pep->nconv,0,0,0);
583:     DSSetState(pep->ds,DS_STATE_RAW);

585:     for (j=0;j<pep->nconv;j++) {
586:       pep->eigr[j] *= pep->sfactor;
587:       pep->eigi[j] *= pep->sfactor;
588:     }
589:   }
590:   STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
591:   RGPopScale(pep->rg);

593:   /* truncate Schur decomposition and change the state to raw so that
594:      DSVectors() computes eigenvectors from scratch */
595:   DSSetDimensions(pep->ds,pep->nconv,0,0,0);
596:   DSSetState(pep->ds,DS_STATE_RAW);
597:   PetscFree2(work,rwork);
598:   return(0);
599: }

603: PetscErrorCode PEPSetFromOptions_STOAR(PetscOptionItems *PetscOptionsObject,PEP pep)
604: {
606:   PetscBool      flg,lock;

609:   PetscOptionsHead(PetscOptionsObject,"PEP STOAR Options");
610:   PetscOptionsBool("-pep_stoar_locking","Choose between locking and non-locking variants","PEPSTOARSetLocking",PETSC_FALSE,&lock,&flg);
611:   if (flg) {
612:     PEPSTOARSetLocking(pep,lock);
613:   }
614:   PetscOptionsTail();
615:   return(0);
616: }

620: static PetscErrorCode PEPSTOARSetLocking_STOAR(PEP pep,PetscBool lock)
621: {
622:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;

625:   ctx->lock = lock;
626:   return(0);
627: }

631: /*@
632:    PEPSTOARSetLocking - Choose between locking and non-locking variants of
633:    the STOAR method.

635:    Logically Collective on PEP

637:    Input Parameters:
638: +  pep  - the eigenproblem solver context
639: -  lock - true if the locking variant must be selected

641:    Options Database Key:
642: .  -pep_stoar_locking - Sets the locking flag

644:    Notes:
645:    The default is to lock converged eigenpairs when the method restarts.
646:    This behaviour can be changed so that all directions are kept in the
647:    working subspace even if already converged to working accuracy (the
648:    non-locking variant).

650:    Level: advanced

652: .seealso: PEPSTOARGetLocking()
653: @*/
654: PetscErrorCode PEPSTOARSetLocking(PEP pep,PetscBool lock)
655: {

661:   PetscTryMethod(pep,"PEPSTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
662:   return(0);
663: }

667: static PetscErrorCode PEPSTOARGetLocking_STOAR(PEP pep,PetscBool *lock)
668: {
669:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;

672:   *lock = ctx->lock;
673:   return(0);
674: }

678: /*@
679:    PEPSTOARGetLocking - Gets the locking flag used in the STOAR method.

681:    Not Collective

683:    Input Parameter:
684: .  pep - the eigenproblem solver context

686:    Output Parameter:
687: .  lock - the locking flag

689:    Level: advanced

691: .seealso: PEPSTOARSetLocking()
692: @*/
693: PetscErrorCode PEPSTOARGetLocking(PEP pep,PetscBool *lock)
694: {

700:   PetscUseMethod(pep,"PEPSTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
701:   return(0);
702: }

706: PetscErrorCode PEPView_STOAR(PEP pep,PetscViewer viewer)
707: {
709:   PEP_TOAR      *ctx = (PEP_TOAR*)pep->data;
710:   PetscBool      isascii;

713:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
714:   if (isascii) {
715:     PetscViewerASCIIPrintf(viewer,"  STOAR: using the %slocking variant\n",ctx->lock?"":"non-");
716:   }
717:   return(0);
718: }

722: PetscErrorCode PEPDestroy_STOAR(PEP pep)
723: {

727:   PetscFree(pep->data);
728:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",NULL);
729:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",NULL);
730:   return(0);
731: }

735: PETSC_EXTERN PetscErrorCode PEPCreate_STOAR(PEP pep)
736: {
738:   PEP_TOAR      *ctx;

741:   PetscNewLog(pep,&ctx);
742:   pep->data = (void*)ctx;
743:   ctx->lock = PETSC_TRUE;

745:   pep->ops->solve          = PEPSolve_STOAR;
746:   pep->ops->setup          = PEPSetUp_STOAR;
747:   pep->ops->setfromoptions = PEPSetFromOptions_STOAR;
748:   pep->ops->view           = PEPView_STOAR;
749:   pep->ops->destroy        = PEPDestroy_STOAR;
750:   pep->ops->backtransform  = PEPBackTransform_Default;
751:   pep->ops->computevectors = PEPComputeVectors_Default;
752:   pep->ops->extractvectors = PEPExtractVectors_TOAR;
753:   pep->ops->reset          = PEPReset_TOAR;
754:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",PEPSTOARSetLocking_STOAR);
755:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",PEPSTOARGetLocking_STOAR);
756:   return(0);
757: }