Actual source code: stoar.c

slepc-3.6.0 2015-06-12
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",
 15:            submitted, 2015.

 17:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 18:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 19:    Copyright (c) 2002-2015, 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>

 43: /*
 44:   Compute B-norm of v=[v1;v2] whith  B=diag(-pep->T[0],pep->T[2]) 
 45: */
 46: static PetscErrorCode PEPSTOARNorm(PEP pep,PetscInt j,PetscReal *norm)
 47: {
 49:   PEP_TOAR      *ctx = (PEP_TOAR*)pep->data;
 50:   PetscBLASInt   n_,one=1,ld_;
 51:   PetscScalar    sone=1.0,szero=0.0,*sp,*sq,*w1,*w2,*qK,*qM;
 52:   PetscInt       n,i,lds=ctx->d*ctx->ld;

 55:   qK = ctx->qB;
 56:   qM = ctx->qB+ctx->ld*ctx->ld;
 57:   n = j+2;
 58:   PetscMalloc2(n,&w1,n,&w2);
 59:   sp = ctx->S+lds*j;
 60:   sq = sp+ctx->ld;
 61:   PetscBLASIntCast(n,&n_);
 62:   PetscBLASIntCast(ctx->ld,&ld_);
 63:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qK,&ld_,sp,&one,&szero,w1,&one));
 64:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qM,&ld_,sq,&one,&szero,w2,&one));
 65:   *norm = 0.0;
 66:   for (i=0;i<n;i++) *norm += PetscRealPart(w1[i]*PetscConj(sp[i])+w2[i]*PetscConj(sq[i]));
 67:   *norm = (*norm>0.0)?PetscSqrtReal(*norm):-PetscSqrtReal(-*norm);
 68:   PetscFree2(w1,w2);
 69:   return(0);
 70: }

 74: static PetscErrorCode PEPSTOARqKqMupdates(PEP pep,PetscInt j,Vec *wv,PetscInt nwv)
 75: {
 77:   PEP_TOAR      *ctx = (PEP_TOAR*)pep->data;
 78:   PetscInt       i,ld=ctx->ld;
 79:   PetscScalar    *qK,*qM;
 80:   Vec            vj,v1,v2;

 83:   qK = ctx->qB;
 84:   qM = ctx->qB+ctx->ld*ctx->ld;
 85:   if (!wv||nwv<2) {
 86:     if (!wv) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",3);
 87:     else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",4);
 88:   }
 89:   v1 = wv[0];
 90:   v2 = wv[1];
 91:   BVGetColumn(pep->V,j,&vj);
 92:   STMatMult(pep->st,0,vj,v1);
 93:   STMatMult(pep->st,2,vj,v2);
 94:   BVRestoreColumn(pep->V,j,&vj);
 95:   for (i=0;i<=j;i++) {
 96:     BVGetColumn(pep->V,i,&vj);
 97:     VecDot(v1,vj,qK+j*ld+i);
 98:     VecDot(v2,vj,qM+j*ld+i);
 99:     *(qM+j*ld+i) *= pep->sfactor*pep->sfactor;
100:     BVRestoreColumn(pep->V,i,&vj);
101:   }
102:   for (i=0;i<j;i++) {
103:     qK[i+j*ld] = -qK[i+ld*j];
104:     qK[j+i*ld] = PetscConj(qK[i+j*ld]);
105:     qM[j+i*ld] = PetscConj(qM[i+j*ld]);
106:   }
107:   qK[j+j*ld] = -PetscRealPart(qK[j+ld*j]);
108:   qM[j+j*ld] = PetscRealPart(qM[j+ld*j]);
109:   return(0);
110: }

114: PetscErrorCode PEPSetUp_STOAR(PEP pep)
115: {
117:   PetscBool      sinv,flg,lindep;
118:   PEP_TOAR      *ctx = (PEP_TOAR*)pep->data;
119:   PetscInt       ld,i;
120:   PetscReal      norm;

123:   PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
124:   if (!ctx->lock && pep->mpd<pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
125:   if (!pep->max_it) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
126:   if (!pep->which) {
127:     PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
128:     if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
129:     else pep->which = PEP_LARGEST_MAGNITUDE;
130:   }
131:   if (pep->problem_type!=PEP_HERMITIAN) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");

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

138:   PEPAllocateSolution(pep,2);
139:   PEPSetWorkVecs(pep,4);
140:   ld = pep->ncv+2;
141:   DSSetType(pep->ds,DSGHIEP);
142:   DSSetCompact(pep->ds,PETSC_TRUE);
143:   DSAllocate(pep->ds,ld);
144:   STGetNumMatrices(pep->st,&ctx->d);
145:   ctx->d--;
146:   ctx->ld = ld;
147:   PetscCalloc1(ctx->d*ld*ld,&ctx->S);
148:   PetscCalloc1(2*ld*ld,&ctx->qB);

150:   /* process starting vector */
151:   if (pep->nini>-2) {  
152:     BVSetRandomColumn(pep->V,0,pep->rand);
153:     BVSetRandomColumn(pep->V,1,pep->rand);
154:   } else {
155:     BVInsertVec(pep->V,0,pep->IS[0]);
156:     BVInsertVec(pep->V,1,pep->IS[1]);
157:   }
158:   BVOrthogonalizeColumn(pep->V,0,NULL,&norm,&lindep);
159:   if (!lindep) {
160:     BVScaleColumn(pep->V,0,1.0/norm);
161:     ctx->S[0] = norm;
162:     PEPSTOARqKqMupdates(pep,0,pep->work,2);
163:   } else {
164:     SETERRQ(PetscObjectComm((PetscObject)pep),1,"Problem with initial vector");
165:   }
166:   BVOrthogonalizeColumn(pep->V,1,ctx->S+ld,&norm,&lindep);
167:   if (!lindep) {
168:     BVScaleColumn(pep->V,1,1.0/norm);
169:     ctx->S[1] = norm;
170:     PEPSTOARqKqMupdates(pep,1,pep->work,2);
171:   } else {
172:     SETERRQ(PetscObjectComm((PetscObject)pep),1,"Problem with initial vector");
173:   }

175:   PEPSTOARNorm(pep,0,&norm);
176:   for (i=0;i<2;i++) { ctx->S[i+ld] /= norm; ctx->S[i] /= norm; }
177:   if (pep->nini<0) {
178:     SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
179:   }
180:   return(0);
181: }

185: /*
186:   Computes GS orthogonalization  x = [z;x] - [Sp;Sq]*y,
187:   where y = Omega\([Sp;Sq]'*[qK zeros(size(qK,1)) ;zeros(size(qK,1)) qM]*[z;x]).
188:   n: Column from S to be orthogonalized against previous columns.
189: */
190: static PetscErrorCode PEPSTOAROrth2(PEP pep,PetscInt k,PetscReal *Omega,PetscScalar *y)
191: {
193:   PEP_TOAR      *ctx = (PEP_TOAR*)pep->data;
194:   PetscBLASInt   n_,lds_,k_,one=1,ld_;
195:   PetscScalar    *S=ctx->S,sonem=-1.0,sone=1.0,szero=0.0,*tp,*tq,*xp,*xq,*c,*qK,*qM;
196:   PetscInt       i,lds=ctx->d*ctx->ld,n,j;
197:   
199:   qK = ctx->qB;
200:   qM = ctx->qB+ctx->ld*ctx->ld;
201:   n = k+2;
202:   PetscMalloc3(n,&tp,n,&tq,k,&c);
203:   PetscBLASIntCast(n,&n_); /* Size of qK and qM */
204:   PetscBLASIntCast(ctx->ld,&ld_);
205:   PetscBLASIntCast(lds,&lds_);
206:   PetscBLASIntCast(k,&k_); /* Number of vectors to orthogonalize against */
207:   xp = S+k*lds;
208:   xq = S+ctx->ld+k*lds;
209:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qK,&ld_,xp,&one,&szero,tp,&one));
210:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qM,&ld_,xq,&one,&szero,tq,&one));
211:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,ctx->S,&lds_,tp,&one,&szero,y,&one));
212:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+ctx->ld,&lds_,tq,&one,&sone,y,&one));
213:   for (i=0;i<k;i++) y[i] /= Omega[i];
214:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S,&lds_,y,&one,&sone,xp,&one));
215:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+ctx->ld,&lds_,y,&one,&sone,xq,&one));
216:   /* three times */
217:   for (j=0;j<2;j++) {
218:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qK,&ld_,xp,&one,&szero,tp,&one));
219:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,qM,&ld_,xq,&one,&szero,tq,&one));
220:     PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,ctx->S,&lds_,tp,&one,&szero,c,&one));
221:     PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&k_,&sone,S+ctx->ld,&lds_,tq,&one,&sone,c,&one));
222:     for (i=0;i<k;i++) c[i] /= Omega[i];
223:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S,&lds_,c,&one,&sone,xp,&one));
224:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&k_,&sonem,S+ctx->ld,&lds_,c,&one,&sone,xq,&one));
225:     for (i=0;i<k;i++) y[i] += c[i];
226:   }
227:   PetscFree3(tp,tq,c);
228:   return(0);
229: }

233: /*
234:   Compute a run of Lanczos iterations
235: */
236: static PetscErrorCode PEPSTOARrun(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,PetscScalar *work,PetscInt nw,Vec *t_,PetscInt nwv)
237: {
239:   PEP_TOAR      *ctx = (PEP_TOAR*)pep->data;
240:   PetscInt       i,j,m=*M,nwu=0,lwa,l;
241:   PetscInt       lds=ctx->d*ctx->ld,offq=ctx->ld;
242:   Vec            v=t_[0],t=t_[1],q=t_[2];
243:   PetscReal      norm,sym=0.0,fro=0.0,*f;
244:   PetscScalar    *y,*S=ctx->S;
245:   PetscBLASInt   j_,one=1;
246:   PetscBool      lindep;

249:   *breakdown = PETSC_FALSE; /* ----- */
250:   if (!t_||nwv<3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",12);
251:   lwa = (ctx->ld)*4;
252:   if (!work||nw<lwa) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",10);
253:   DSGetDimensions(pep->ds,NULL,NULL,&l,NULL,NULL);
254:   y = work;
255:   nwu += ctx->ld;
256:   for (j=k;j<m;j++) {
257:     /* apply operator */
258:     BVSetActiveColumns(pep->V,0,j+2);
259:     BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
260:     STMatMult(pep->st,0,v,t);
261:     BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
262:     STMatMult(pep->st,1,v,q);
263:     VecAXPY(t,pep->sfactor,q);
264:     STMatSolve(pep->st,t,q);
265:     VecScale(q,-1.0/(pep->sfactor*pep->sfactor));

267:     /* orthogonalize */
268:     BVOrthogonalizeVec(pep->V,q,S+offq+(j+1)*lds,&norm,&lindep);
269:     if (lindep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STOAR does not support detection of linearly dependent TOAR vectors");
270:     *(S+offq+(j+1)*lds+j+2) = norm;
271:     VecScale(q,1.0/norm);
272:     BVInsertVec(pep->V,j+2,q);
273:     for (i=0;i<=j+1;i++) *(S+(j+1)*lds+i) = *(S+offq+j*lds+i);
274:    
275:     /* update qK and qM */
276:     PEPSTOARqKqMupdates(pep,j+2,t_,2);

278:     /* level-2 orthogonalization */
279:     PEPSTOAROrth2(pep,j+1,omega,y);
280:     a[j] = PetscRealPart(y[j])/omega[j];
281:     PEPSTOARNorm(pep,j+1,&norm);
282:     omega[j+1] = (norm > 0)?1.0:-1.0;
283:     for (i=0;i<=j+2;i++) {
284:       S[i+(j+1)*lds] /= norm;
285:       S[i+offq+(j+1)*lds] /= norm;
286:     }
287:     b[j] = PetscAbsReal(norm);

289:     /* check symmetry */
290:     DSGetArrayReal(pep->ds,DS_MAT_T,&f);
291:     if (j==k) {
292:       for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ctx->ld+i]);
293:       for (i=0;i<l;i++) y[i] = 0.0;
294:     }
295:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
296:     if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
297:     PetscBLASIntCast(j,&j_);
298:     sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
299:     fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
300:     if (j>0) fro = SlepcAbs(fro,b[j-1]);
301:     if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
302:       *symmlost = PETSC_TRUE;
303:       *M=j+1;
304:       break;
305:     }
306:   }
307:   return(0);
308: }

312: static PetscErrorCode PEPSTOARTrunc(PEP pep,PetscInt rs1,PetscInt cs1,PetscScalar *work,PetscInt nw,PetscReal *rwork,PetscInt nrw)
313: {
315:   PEP_TOAR      *ctx = (PEP_TOAR*)pep->data;
316:   Mat            G;
317:   PetscInt       lwa,nwu=0,lrwa,nrwu=0;
318:   PetscInt       i,n,lds=2*ctx->ld;
319:   PetscScalar    *M,*V,*U,*S=ctx->S,sone=1.0,zero=0.0,t,*qK,*qM;
320:   PetscReal      *sg;
321:   PetscBLASInt   cs1_,rs1_,cs1t2,cs1p1,n_,info,lw_,lds_,ld_;

324:   qK = ctx->qB;
325:   qM = ctx->qB+ctx->ld*ctx->ld;
326:   n = (rs1>2*cs1)?2*cs1:rs1;
327:   lwa = cs1*rs1*4+n*(rs1+2*cs1)+(cs1+1)*(cs1+2);
328:   lrwa = n+cs1+1+5*n;
329:   if (!work||nw<lwa) {
330:     if (nw<lwa) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",6);
331:     if (!work) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",5);
332:   }
333:   if (!rwork||nrw<lrwa) {
334:     if (nrw<lrwa) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",8);
335:     if (!work) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",7);
336:   }
337:   M = work+nwu;
338:   nwu += rs1*cs1*2;
339:   U = work+nwu;
340:   nwu += rs1*n;
341:   V = work+nwu;
342:   nwu += 2*cs1*n;
343:   sg = rwork+nrwu;
344:   nrwu += n;
345:   for (i=0;i<cs1;i++) {
346:     PetscMemcpy(M+i*rs1,S+i*lds,rs1*sizeof(PetscScalar));  
347:     PetscMemcpy(M+(i+cs1)*rs1,S+i*lds+ctx->ld,rs1*sizeof(PetscScalar));  
348:   }
349:   PetscBLASIntCast(n,&n_);
350:   PetscBLASIntCast(cs1,&cs1_);
351:   PetscBLASIntCast(rs1,&rs1_);
352:   PetscBLASIntCast(cs1*2,&cs1t2);
353:   PetscBLASIntCast(cs1+1,&cs1p1);
354:   PetscBLASIntCast(lds,&lds_);
355:   PetscBLASIntCast(ctx->ld,&ld_);
356:   PetscBLASIntCast(lwa-nwu,&lw_);
357: #if !defined(PETSC_USE_COMPLEX)
358:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1t2,M,&rs1_,sg,U,&rs1_,V,&n_,work+nwu,&lw_,&info));
359: #else
360:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&rs1_,&cs1t2,M,&rs1_,sg,U,&rs1_,V,&n_,work+nwu,&lw_,rwork+nrwu,&info));  
361: #endif
362:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESVD %d",info);

364:   /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
365:   MatCreateSeqDense(PETSC_COMM_SELF,rs1,2*cs1,U,&G);
366:   BVSetActiveColumns(pep->V,0,rs1);
367:   BVMultInPlace(pep->V,G,0,cs1+1);
368:   MatDestroy(&G);
369:   
370:   /* Update S */
371:   PetscMemzero(S,lds*ctx->ld*sizeof(PetscScalar));

373:   for (i=0;i<cs1+1;i++) {
374:     t = sg[i];
375:     PetscStackCallBLAS("BLASscal",BLASscal_(&cs1t2,&t,V+i,&n_));
376:   }
377:   for (i=0;i<cs1;i++) {
378:     PetscMemcpy(S+i*lds,V+i*n,(cs1+1)*sizeof(PetscScalar));
379:     PetscMemcpy(S+ctx->ld+i*lds,V+(cs1+i)*n,(cs1+1)*sizeof(PetscScalar));
380:   }

382:   /* Update qM and qK */
383:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&rs1_,&cs1p1,&rs1_,&sone,qK,&ld_,U,&rs1_,&zero,work+nwu,&rs1_));
384:   PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&cs1p1,&cs1p1,&rs1_,&sone,U,&rs1_,work+nwu,&rs1_,&zero,qK,&ld_));
385:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&rs1_,&cs1p1,&rs1_,&sone,qM,&ld_,U,&rs1_,&zero,work+nwu,&rs1_));
386:   PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&cs1p1,&cs1p1,&rs1_,&sone,U,&rs1_,work+nwu,&rs1_,&zero,qM,&ld_));
387:   return(0);
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: */
398: static PetscErrorCode PEPSTOARSupdate(PetscScalar *S,PetscInt ld,PetscInt sr,PetscInt s,PetscInt ncu,PetscInt qr,PetscScalar *Q,PetscInt ldq,PetscScalar *work,PetscInt nw)
399: {
401:   PetscScalar    a=1.0,b=0.0;
402:   PetscBLASInt   sr_,ncu_,ldq_,lds_,qr_;
403:   PetscInt       lwa,j,lds=2*ld;

406:   lwa = sr*ncu;
407:   if (!work||nw<lwa) {
408:     if (nw<lwa) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",10);
409:     if (!work) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid argument %d",9);
410:   }
411:   PetscBLASIntCast(sr,&sr_);
412:   PetscBLASIntCast(qr,&qr_);
413:   PetscBLASIntCast(ncu,&ncu_);
414:   PetscBLASIntCast(lds,&lds_);
415:   PetscBLASIntCast(ldq,&ldq_);
416:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&ncu_,&qr_,&a,S,&lds_,Q,&ldq_,&b,work,&sr_));
417:   for (j=0;j<ncu;j++) {
418:     PetscMemcpy(S+lds*(s+j),work+j*sr,sr*sizeof(PetscScalar));
419:   }
420:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&ncu_,&qr_,&a,S+ld,&lds_,Q,&ldq_,&b,work,&sr_));
421:   for (j=0;j<ncu;j++) {
422:     PetscMemcpy(S+lds*(s+j)+ld,work+j*sr,sr*sizeof(PetscScalar));
423:   }
424:   return(0);
425: }

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

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

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

470:   BVSetMatrix(pep->V,NULL,PETSC_FALSE);
471:   lwa = 9*ld*ld+5*ld;
472:   lrwa = 8*ld;
473:   PetscMalloc2(lwa,&work,lrwa,&rwork); /* REVIEW */
474:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
475:   RGSetScale(pep->rg,sinv?1.0/pep->sfactor:pep->sfactor);
476:   STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);

478:   /* Restart loop */
479:   l = 0;
480:   DSGetLeadingDimension(pep->ds,&ldds);
481:   while (pep->reason == PEP_CONVERGED_ITERATING) {
482:     pep->its++;
483:     DSGetArrayReal(pep->ds,DS_MAT_T,&a);
484:     b = a+ldds;
485:     r = b+ldds;
486:     DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
487:     
488:     /* Compute an nv-step Lanczos factorization */
489:     nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
490:     PEPSTOARrun(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,work+nwu,lwa-nwu,pep->work,3);
491:     beta = b[nv-1];
492:     if (symmlost) {
493:       pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
494:       if (nv==pep->nconv+l+1) { pep->nconv = nconv; break; }
495:     }
496:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
497:     DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
498:     DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
499:     if (l==0) {
500:       DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
501:     } else {
502:       DSSetState(pep->ds,DS_STATE_RAW);
503:     }

505:     /* Solve projected problem */
506:     DSSolve(pep->ds,pep->eigr,pep->eigi);
507:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);

509:     /* Check convergence */
510:     /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
511:     norm = 1.0;
512:     DSGetDimensions(pep->ds,NULL,NULL,NULL,NULL,&t);    
513:     PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
514:     nconv = k;
515:     if (pep->its >= pep->max_it) pep->reason = PEP_DIVERGED_ITS;
516:     if (k >= pep->nev) pep->reason = PEP_CONVERGED_TOL;

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

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

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

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


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

573:   if (pep->nconv>0) {
574:     /* Truncate S */
575:     DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
576:     PEPSTOARTrunc(pep,nv+2,pep->nconv,work+nwu,lwa-nwu,rwork+nrwu,lrwa-nrwu);
577:     DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
578:   
579:     /* Extraction */
580:     DSSetDimensions(pep->ds,pep->nconv,0,0,0);
581:     DSSetState(pep->ds,DS_STATE_RAW);
582:   
583:     for (j=0;j<pep->nconv;j++) {
584:       pep->eigr[j] *= pep->sfactor;
585:       pep->eigi[j] *= pep->sfactor;
586:     }
587:   }
588:   STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
589:   RGSetScale(pep->rg,1.0);

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

601: PetscErrorCode PEPSetFromOptions_STOAR(PetscOptions *PetscOptionsObject,PEP pep)
602: {
604:   PetscBool      flg,lock;

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

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

623:   ctx->lock = lock;
624:   return(0);
625: }

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

633:    Logically Collective on PEP

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

639:    Options Database Key:
640: .  -pep_stoar_locking - Sets the locking flag

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

648:    Level: advanced

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

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

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

670:   *lock = ctx->lock;
671:   return(0);
672: }

676: /*@
677:    PEPSTOARGetLocking - Gets the locking flag used in the STOAR method.

679:    Not Collective

681:    Input Parameter:
682: .  pep - the eigenproblem solver context

684:    Output Parameter:
685: .  lock - the locking flag

687:    Level: advanced

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

698:   PetscTryMethod(pep,"PEPSTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
699:   return(0);
700: }

704: PetscErrorCode PEPView_STOAR(PEP pep,PetscViewer viewer)
705: {
707:   PEP_TOAR      *ctx = (PEP_TOAR*)pep->data;
708:   PetscBool      isascii;

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

720: PetscErrorCode PEPDestroy_STOAR(PEP pep)
721: {

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

733: PETSC_EXTERN PetscErrorCode PEPCreate_STOAR(PEP pep)
734: {
736:   PEP_TOAR      *ctx;

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

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