Actual source code: ptoar.c

slepc-3.13.0 2020-03-31
Report Typos and Errors
  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>    /*I "slepcpep.h" I*/
 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      shift,sinv,flg;
 55:   PetscInt       i;

 58:   pep->lineariz = PETSC_TRUE;
 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) pep->max_it = PetscMax(100,2*(pep->nmat-1)*pep->n/pep->ncv);

 63:   PetscObjectTypeCompare((PetscObject)pep->st,STSHIFT,&shift);
 64:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
 65:   if (!shift && !sinv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Only STSHIFT and STSINVERT spectral transformations can be used");
 66:   if (!pep->which) {
 67:     if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
 68:     else pep->which = PEP_LARGEST_MAGNITUDE;
 69:   }
 70:   if (pep->which==PEP_ALL) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of pep->which");
 71:   if (pep->problem_type!=PEP_GENERAL) {
 72:     PetscInfo(pep,"Problem type ignored, performing a non-symmetric linearization\n");
 73:   }

 75:   if (!ctx->keep) ctx->keep = 0.5;

 77:   PEPAllocateSolution(pep,pep->nmat-1);
 78:   PEPSetWorkVecs(pep,3);
 79:   DSSetType(pep->ds,DSNHEP);
 80:   DSSetExtraRow(pep->ds,PETSC_TRUE);
 81:   DSAllocate(pep->ds,pep->ncv+1);

 83:   PEPBasisCoefficients(pep,pep->pbc);
 84:   STGetTransform(pep->st,&flg);
 85:   if (!flg) {
 86:     PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
 87:     PetscLogObjectMemory((PetscObject)pep,pep->nmat*sizeof(PetscScalar));
 88:     if (sinv) {
 89:       PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL);
 90:     } else {
 91:       for (i=0;i<pep->nmat-1;i++) pep->solvematcoeffs[i] = 0.0;
 92:       pep->solvematcoeffs[pep->nmat-1] = 1.0;
 93:     }
 94:   }
 95:   BVDestroy(&ctx->V);
 96:   BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
 97:   return(0);
 98: }

100: /*
101:   Extend the TOAR basis by applying the the matrix operator
102:   over a vector which is decomposed in the TOAR way
103:   Input:
104:     - pbc: array containing the polynomial basis coefficients
105:     - S,V: define the latest Arnoldi vector (nv vectors in V)
106:   Output:
107:     - t: new vector extending the TOAR basis
108:     - r: temporary coefficients to compute the TOAR coefficients
109:          for the new Arnoldi vector
110:   Workspace: t_ (two vectors)
111: */
112: 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_)
113: {
115:   PetscInt       nmat=pep->nmat,deg=nmat-1,k,j,off=0,lss;
116:   Vec            v=t_[0],ve=t_[1],q=t_[2];
117:   PetscScalar    alpha=1.0,*ss,a;
118:   PetscReal      *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat;
119:   PetscBool      flg;

122:   BVSetActiveColumns(pep->V,0,nv);
123:   STGetTransform(pep->st,&flg);
124:   if (sinvert) {
125:     for (j=0;j<nv;j++) {
126:       if (deg>1) r[lr+j] = S[j]/ca[0];
127:       if (deg>2) r[2*lr+j] = (S[ls+j]+(sigma-cb[1])*r[lr+j])/ca[1];
128:     }
129:     for (k=2;k<deg-1;k++) {
130:       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];
131:     }
132:     k = deg-1;
133:     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];
134:     ss = r; lss = lr; off = 1; alpha = -1.0; a = pep->sfactor;
135:   } else {
136:     ss = S; lss = ls; off = 0; alpha = -ca[deg-1]; a = 1.0;
137:   }
138:   BVMultVec(V,1.0,0.0,v,ss+off*lss);
139:   if (pep->Dr) { /* balancing */
140:     VecPointwiseMult(v,v,pep->Dr);
141:   }
142:   STMatMult(pep->st,off,v,q);
143:   VecScale(q,a);
144:   for (j=1+off;j<deg+off-1;j++) {
145:     BVMultVec(V,1.0,0.0,v,ss+j*lss);
146:     if (pep->Dr) {
147:       VecPointwiseMult(v,v,pep->Dr);
148:     }
149:     STMatMult(pep->st,j,v,t);
150:     a *= pep->sfactor;
151:     VecAXPY(q,a,t);
152:   }
153:   if (sinvert) {
154:     BVMultVec(V,1.0,0.0,v,ss);
155:     if (pep->Dr) {
156:       VecPointwiseMult(v,v,pep->Dr);
157:     }
158:     STMatMult(pep->st,deg,v,t);
159:     a *= pep->sfactor;
160:     VecAXPY(q,a,t);
161:   } else {
162:     BVMultVec(V,1.0,0.0,ve,ss+(deg-1)*lss);
163:     if (pep->Dr) {
164:       VecPointwiseMult(ve,ve,pep->Dr);
165:     }
166:     a *= pep->sfactor;
167:     STMatMult(pep->st,deg-1,ve,t);
168:     VecAXPY(q,a,t);
169:     a *= pep->sfactor;
170:   }
171:   if (flg || !sinvert) alpha /= a;
172:   STMatSolve(pep->st,q,t);
173:   VecScale(t,alpha);
174:   if (!sinvert) {
175:     if (cg[deg-1]!=0) { VecAXPY(t,cg[deg-1],v); }
176:     if (cb[deg-1]!=0) { VecAXPY(t,cb[deg-1],ve); }
177:   }
178:   if (pep->Dr) {
179:     VecPointwiseDivide(t,t,pep->Dr);
180:   }
181:   return(0);
182: }

184: /*
185:   Compute TOAR coefficients of the blocks of the new Arnoldi vector computed
186: */
187: static PetscErrorCode PEPTOARCoefficients(PEP pep,PetscBool sinvert,PetscScalar sigma,PetscInt nv,PetscScalar *S,PetscInt ls,PetscScalar *r,PetscInt lr,PetscScalar *x)
188: {
189:   PetscInt    k,j,nmat=pep->nmat,d=nmat-1;
190:   PetscReal   *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat;
191:   PetscScalar t=1.0,tp=0.0,tt;

194:   if (sinvert) {
195:     for (k=1;k<d;k++) {
196:       tt = t;
197:       t = ((sigma-cb[k-1])*t-cg[k-1]*tp)/ca[k-1]; /* k-th basis polynomial */
198:       tp = tt;
199:       for (j=0;j<=nv;j++) r[k*lr+j] += t*x[j];
200:     }
201:   } else {
202:     for (j=0;j<=nv;j++) r[j] = (cb[0]-sigma)*S[j]+ca[0]*S[ls+j];
203:     for (k=1;k<d-1;k++) {
204:       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];
205:     }
206:     if (sigma!=0.0) for (j=0;j<=nv;j++) r[(d-1)*lr+j] -= sigma*S[(d-1)*ls+j];
207:   }
208:   return(0);
209: }

211: /*
212:   Compute a run of Arnoldi iterations dim(work)=ld
213: */
214: static PetscErrorCode PEPTOARrun(PEP pep,PetscScalar sigma,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,PetscBool *breakdown,Vec *t_)
215: {
217:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
218:   PetscInt       j,m=*M,deg=pep->nmat-1,ld;
219:   PetscInt       lds,nqt,l;
220:   Vec            t;
221:   PetscReal      norm;
222:   PetscBool      flg,sinvert=PETSC_FALSE,lindep;
223:   PetscScalar    *x,*S;
224:   Mat            MS;

227:   BVTensorGetFactors(ctx->V,NULL,&MS);
228:   MatDenseGetArray(MS,&S);
229:   BVGetSizes(pep->V,NULL,NULL,&ld);
230:   lds = ld*deg;
231:   BVGetActiveColumns(pep->V,&l,&nqt);
232:   STGetTransform(pep->st,&flg);
233:   if (!flg) {
234:     /* spectral transformation handled by the solver */
235:     PetscObjectTypeCompareAny((PetscObject)pep->st,&flg,STSINVERT,STSHIFT,"");
236:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"ST type not supported for TOAR without transforming matrices");
237:     PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinvert);
238:   }
239:   BVSetActiveColumns(ctx->V,0,m);
240:   for (j=k;j<m;j++) {
241:     /* apply operator */
242:     BVGetColumn(pep->V,nqt,&t);
243:     PEPTOARExtendBasis(pep,sinvert,sigma,S+j*lds,ld,nqt,pep->V,t,S+(j+1)*lds,ld,t_);
244:     BVRestoreColumn(pep->V,nqt,&t);

246:     /* orthogonalize */
247:     if (sinvert) x = S+(j+1)*lds;
248:     else x = S+(deg-1)*ld+(j+1)*lds;
249:     BVOrthogonalizeColumn(pep->V,nqt,x,&norm,&lindep);
250:     if (!lindep) {
251:       x[nqt] = norm;
252:       BVScaleColumn(pep->V,nqt,1.0/norm);
253:       nqt++;
254:     }

256:     PEPTOARCoefficients(pep,sinvert,sigma,nqt-1,S+j*lds,ld,S+(j+1)*lds,ld,x);

258:     /* level-2 orthogonalization */
259:     BVOrthogonalizeColumn(ctx->V,j+1,H+j*ldh,&norm,breakdown);
260:     H[j+1+ldh*j] = norm;
261:     if (*breakdown) {
262:       *M = j+1;
263:       break;
264:     }
265:     BVScaleColumn(ctx->V,j+1,1.0/norm);
266:     BVSetActiveColumns(pep->V,l,nqt);
267:   }
268:   BVSetActiveColumns(ctx->V,0,*M);
269:   MatDenseRestoreArray(MS,&S);
270:   BVTensorRestoreFactors(ctx->V,NULL,&MS);
271:   return(0);
272: }

274: /*
275:   Computes T_j = phi_idx(T). In T_j and T_p are phi_{idx-1}(T)
276:    and phi_{idx-2}(T) respectively or null if idx=0,1.
277:    Tp and Tj are input/output arguments
278: */
279: static PetscErrorCode PEPEvaluateBasisM(PEP pep,PetscInt k,PetscScalar *T,PetscInt ldt,PetscInt idx,PetscScalar **Tp,PetscScalar **Tj)
280: {
282:   PetscInt       i;
283:   PetscReal      *ca,*cb,*cg;
284:   PetscScalar    *pt,g,a;
285:   PetscBLASInt   k_,ldt_;

288:   if (idx==0) {
289:     PetscArrayzero(*Tj,k*k);
290:     PetscArrayzero(*Tp,k*k);
291:     for (i=0;i<k;i++) (*Tj)[i+i*k] = 1.0;
292:   } else {
293:     PetscBLASIntCast(ldt,&ldt_);
294:     PetscBLASIntCast(k,&k_);
295:     ca = pep->pbc; cb = pep->pbc+pep->nmat; cg = pep->pbc+2*pep->nmat;
296:     for (i=0;i<k;i++) T[i*ldt+i] -= cb[idx-1];
297:     a = 1/ca[idx-1];
298:     g = (idx==1)?0.0:-cg[idx-1]/ca[idx-1];
299:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&a,T,&ldt_,*Tj,&k_,&g,*Tp,&k_));
300:     pt = *Tj; *Tj = *Tp; *Tp = pt;
301:     for (i=0;i<k;i++) T[i*ldt+i] += cb[idx-1];
302:   }
303:   return(0);
304: }

306: static PetscErrorCode PEPExtractInvariantPair(PEP pep,PetscScalar sigma,PetscInt sr,PetscInt k,PetscScalar *S,PetscInt ld,PetscInt deg,PetscScalar *H,PetscInt ldh)
307: {
309:   PetscInt       i,j,jj,lds,ldt,d=pep->nmat-1,idxcpy=0;
310:   PetscScalar    *At,*Bt,*Hj,*Hp,*T,sone=1.0,g,a,*pM,*work;
311:   PetscBLASInt   k_,sr_,lds_,ldh_,info,*p,lwork,ldt_;
312:   PetscBool      transf=PETSC_FALSE,flg;
313:   PetscReal      nrm,norm,maxnrm,*rwork;
314:   BV             *R,Y;
315:   Mat            M,*A;
316:   Vec            v;

319:   if (k==0) return(0);
320:   lds = deg*ld;
321:   PetscCalloc6(k,&p,sr*k,&At,k*k,&Bt,k*k,&Hj,k*k,&Hp,sr*k,&work);
322:   PetscBLASIntCast(sr,&sr_);
323:   PetscBLASIntCast(k,&k_);
324:   PetscBLASIntCast(lds,&lds_);
325:   PetscBLASIntCast(ldh,&ldh_);
326:   STGetTransform(pep->st,&flg);
327:   if (!flg) {
328:      PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&flg);
329:     if (flg || sigma!=0.0) transf=PETSC_TRUE;
330:   }
331:   if (transf) {
332:     PetscMalloc1(k*k,&T);
333:     ldt = k;
334:     for (i=0;i<k;i++) {
335:       PetscArraycpy(T+k*i,H+i*ldh,k);
336:     }
337:     if (flg) {
338:       PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&k_,&k_,T,&k_,p,&info));
339:       SlepcCheckLapackInfo("getrf",info);
340:       PetscBLASIntCast(sr*k,&lwork);
341:       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&k_,T,&k_,p,work,&lwork,&info));
342:       SlepcCheckLapackInfo("getri",info);
343:     }
344:     if (sigma!=0.0) for (i=0;i<k;i++) T[i+k*i] += sigma;
345:   } else {
346:     T = H; ldt = ldh;
347:   }
348:   PetscBLASIntCast(ldt,&ldt_);
349:   switch (pep->extract) {
350:   case PEP_EXTRACT_NONE:
351:     break;
352:   case PEP_EXTRACT_NORM:
353:     if (pep->basis == PEP_BASIS_MONOMIAL) {
354:       PetscBLASIntCast(ldt,&ldt_);
355:       PetscMalloc1(k,&rwork);
356:       norm = LAPACKlange_("F",&k_,&k_,T,&ldt_,rwork);
357:       PetscFree(rwork);
358:       if (norm>1.0) idxcpy = d-1;
359:     } else {
360:       PetscBLASIntCast(ldt,&ldt_);
361:       PetscMalloc1(k,&rwork);
362:       maxnrm = 0.0;
363:       for (i=0;i<pep->nmat-1;i++) {
364:         PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
365:         norm = LAPACKlange_("F",&k_,&k_,Hj,&k_,rwork);
366:         if (norm > maxnrm) {
367:           idxcpy = i;
368:           maxnrm = norm;
369:         }
370:       }
371:       PetscFree(rwork);
372:     }
373:     if (idxcpy>0) {
374:       /* copy block idxcpy of S to the first one */
375:       for (j=0;j<k;j++) {
376:         PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr);
377:       }
378:     }
379:     break;
380:   case PEP_EXTRACT_RESIDUAL:
381:     STGetTransform(pep->st,&flg);
382:     if (flg) {
383:       PetscMalloc1(pep->nmat,&A);
384:       for (i=0;i<pep->nmat;i++) {
385:         STGetMatrixTransformed(pep->st,i,A+i);
386:       }
387:     } else A = pep->A;
388:     PetscMalloc1(pep->nmat-1,&R);
389:     for (i=0;i<pep->nmat-1;i++) {
390:       BVDuplicateResize(pep->V,k,R+i);
391:     }
392:     BVDuplicateResize(pep->V,sr,&Y);
393:     MatCreateSeqDense(PETSC_COMM_SELF,sr,k,NULL,&M);
394:     g = 0.0; a = 1.0;
395:     BVSetActiveColumns(pep->V,0,sr);
396:     for (j=0;j<pep->nmat;j++) {
397:       BVMatMult(pep->V,A[j],Y);
398:       PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
399:       for (i=0;i<pep->nmat-1;i++) {
400:         PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&k_,&k_,&a,S+i*ld,&lds_,Hj,&k_,&g,At,&sr_));
401:         MatDenseGetArray(M,&pM);
402:         for (jj=0;jj<k;jj++) {
403:           PetscArraycpy(pM+jj*sr,At+jj*sr,sr);
404:         }
405:         MatDenseRestoreArray(M,&pM);
406:         BVMult(R[i],1.0,(i==0)?0.0:1.0,Y,M);
407:       }
408:     }

410:     /* frobenius norm */
411:     maxnrm = 0.0;
412:     for (i=0;i<pep->nmat-1;i++) {
413:       norm = 0.0;
414:       for (j=0;j<k;j++) {
415:         BVGetColumn(R[i],j,&v);
416:         VecNorm(v,NORM_2,&nrm);
417:         BVRestoreColumn(R[i],j,&v);
418:         norm += nrm*nrm;
419:       }
420:       norm = PetscSqrtReal(norm);
421:       if (maxnrm > norm) {
422:         maxnrm = norm;
423:         idxcpy = i;
424:       }
425:     }
426:     if (idxcpy>0) {
427:       /* copy block idxcpy of S to the first one */
428:       for (j=0;j<k;j++) {
429:         PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr);
430:       }
431:     }
432:     if (flg) PetscFree(A);
433:     for (i=0;i<pep->nmat-1;i++) {
434:       BVDestroy(&R[i]);
435:     }
436:     PetscFree(R);
437:     BVDestroy(&Y);
438:     MatDestroy(&M);
439:     break;
440:   case PEP_EXTRACT_STRUCTURED:
441:     for (j=0;j<k;j++) Bt[j+j*k] = 1.0;
442:     for (j=0;j<sr;j++) {
443:       for (i=0;i<k;i++) At[j*k+i] = PetscConj(S[i*lds+j]);
444:     }
445:     PEPEvaluateBasisM(pep,k,T,ldt,0,&Hp,&Hj);
446:     for (i=1;i<deg;i++) {
447:       PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj);
448:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&sr_,&k_,&sone,Hj,&k_,S+i*ld,&lds_,&sone,At,&k_));
449:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&k_,&k_,&sone,Hj,&k_,Hj,&k_,&sone,Bt,&k_));
450:     }
451:     PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&k_,&sr_,Bt,&k_,p,At,&k_,&info));
452:     SlepcCheckLapackInfo("gesv",info);
453:     for (j=0;j<sr;j++) {
454:       for (i=0;i<k;i++) S[i*lds+j] = PetscConj(At[j*k+i]);
455:     }
456:     break;
457:   }
458:   if (transf) { PetscFree(T); }
459:   PetscFree6(p,At,Bt,Hj,Hp,work);
460:   return(0);
461: }

463: PetscErrorCode PEPSolve_TOAR(PEP pep)
464: {
466:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
467:   PetscInt       i,j,k,l,nv=0,ld,lds,ldds,newn,nq=0,nconv=0;
468:   PetscInt       nmat=pep->nmat,deg=nmat-1;
469:   PetscScalar    *S,*H,sigma;
470:   PetscReal      beta;
471:   PetscBool      breakdown=PETSC_FALSE,flg,falselock=PETSC_FALSE,sinv=PETSC_FALSE;
472:   Mat            MS,MQ;

475:   PetscCitationsRegister(citation,&cited);
476:   if (ctx->lock) {
477:     /* undocumented option to use a cheaper locking instead of the true locking */
478:     PetscOptionsGetBool(NULL,NULL,"-pep_toar_falselocking",&falselock,NULL);
479:   }
480:   DSGetLeadingDimension(pep->ds,&ldds);
481:   STGetShift(pep->st,&sigma);

483:   /* update polynomial basis coefficients */
484:   STGetTransform(pep->st,&flg);
485:   if (pep->sfactor!=1.0) {
486:     for (i=0;i<nmat;i++) {
487:       pep->pbc[nmat+i] /= pep->sfactor;
488:       pep->pbc[2*nmat+i] /= pep->sfactor*pep->sfactor;
489:     }
490:     if (!flg) {
491:       pep->target /= pep->sfactor;
492:       RGPushScale(pep->rg,1.0/pep->sfactor);
493:       STScaleShift(pep->st,1.0/pep->sfactor);
494:       sigma /= pep->sfactor;
495:     } else {
496:       PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
497:       pep->target = sinv?pep->target*pep->sfactor:pep->target/pep->sfactor;
498:       RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
499:       STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
500:     }
501:   }

503:   if (flg) sigma = 0.0;

505:   /* clean projected matrix (including the extra-arrow) */
506:   DSGetArray(pep->ds,DS_MAT_A,&H);
507:   PetscArrayzero(H,ldds*ldds);
508:   DSRestoreArray(pep->ds,DS_MAT_A,&H);

510:   /* Get the starting Arnoldi vector */
511:   BVTensorBuildFirstColumn(ctx->V,pep->nini);

513:   /* restart loop */
514:   l = 0;
515:   while (pep->reason == PEP_CONVERGED_ITERATING) {
516:     pep->its++;

518:     /* compute an nv-step Lanczos factorization */
519:     nv = PetscMax(PetscMin(nconv+pep->mpd,pep->ncv),nv);
520:     DSGetArray(pep->ds,DS_MAT_A,&H);
521:     PEPTOARrun(pep,sigma,H,ldds,pep->nconv+l,&nv,&breakdown,pep->work);
522:     beta = PetscAbsScalar(H[(nv-1)*ldds+nv]);
523:     DSRestoreArray(pep->ds,DS_MAT_A,&H);
524:     DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
525:     if (l==0) {
526:       DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
527:     } else {
528:       DSSetState(pep->ds,DS_STATE_RAW);
529:     }
530:     BVSetActiveColumns(ctx->V,pep->nconv,nv);

532:     /* solve projected problem */
533:     DSSolve(pep->ds,pep->eigr,pep->eigi);
534:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
535:     DSUpdateExtraRow(pep->ds);
536:     DSSynchronize(pep->ds,pep->eigr,pep->eigi);

538:     /* check convergence */
539:     PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
540:     (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);

542:     /* update l */
543:     if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
544:     else {
545:       l = (nv==k)?0:PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
546:       if (!breakdown) {
547:         /* prepare the Rayleigh quotient for restart */
548:         DSTruncate(pep->ds,k+l);
549:         DSGetDimensions(pep->ds,&newn,NULL,NULL,NULL,NULL);
550:         l = newn-k;
551:       }
552:     }
553:     nconv = k;
554:     if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */

556:     /* update S */
557:     DSGetMat(pep->ds,DS_MAT_Q,&MQ);
558:     BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
559:     MatDestroy(&MQ);

561:     /* copy last column of S */
562:     BVCopyColumn(ctx->V,nv,k+l);

564:     if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
565:       /* stop if breakdown */
566:       PetscInfo2(pep,"Breakdown TOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
567:       pep->reason = PEP_DIVERGED_BREAKDOWN;
568:     }
569:     if (pep->reason != PEP_CONVERGED_ITERATING) l--;
570:     /* truncate S */
571:     BVGetActiveColumns(pep->V,NULL,&nq);
572:     if (k+l+deg<=nq) {
573:       BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
574:       if (!falselock && ctx->lock) {
575:         BVTensorCompress(ctx->V,k-pep->nconv);
576:       } else {
577:         BVTensorCompress(ctx->V,0);
578:       }
579:     }
580:     pep->nconv = k;
581:     PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
582:   }
583:   if (pep->nconv>0) {
584:     /* {V*S_nconv^i}_{i=0}^{d-1} has rank nconv instead of nconv+d-1. Force zeros in each S_nconv^i block */
585:     BVSetActiveColumns(ctx->V,0,pep->nconv);
586:     BVGetActiveColumns(pep->V,NULL,&nq);
587:     BVSetActiveColumns(pep->V,0,nq);
588:     if (nq>pep->nconv) {
589:       BVTensorCompress(ctx->V,pep->nconv);
590:       BVSetActiveColumns(pep->V,0,pep->nconv);
591:       nq = pep->nconv;
592:     }

594:     /* perform Newton refinement if required */
595:     if (pep->refine==PEP_REFINE_MULTIPLE && pep->rits>0) {
596:       /* extract invariant pair */
597:       BVTensorGetFactors(ctx->V,NULL,&MS);
598:       MatDenseGetArray(MS,&S);
599:       DSGetArray(pep->ds,DS_MAT_A,&H);
600:       BVGetSizes(pep->V,NULL,NULL,&ld);
601:       lds = deg*ld;
602:       PEPExtractInvariantPair(pep,sigma,nq,pep->nconv,S,ld,deg,H,ldds);
603:       DSRestoreArray(pep->ds,DS_MAT_A,&H);
604:       DSSetDimensions(pep->ds,pep->nconv,0,0,0);
605:       DSSetState(pep->ds,DS_STATE_RAW);
606:       PEPNewtonRefinement_TOAR(pep,sigma,&pep->rits,NULL,pep->nconv,S,lds);
607:       DSSolve(pep->ds,pep->eigr,pep->eigi);
608:       DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
609:       DSSynchronize(pep->ds,pep->eigr,pep->eigi);
610:       DSGetMat(pep->ds,DS_MAT_Q,&MQ);
611:       BVMultInPlace(ctx->V,MQ,0,pep->nconv);
612:       MatDestroy(&MQ);
613:       MatDenseRestoreArray(MS,&S);
614:       BVTensorRestoreFactors(ctx->V,NULL,&MS);
615:     }
616:   }
617:   STGetTransform(pep->st,&flg);
618:   if (pep->refine!=PEP_REFINE_MULTIPLE || pep->rits==0) {
619:     if (!flg && pep->ops->backtransform) {
620:         (*pep->ops->backtransform)(pep);
621:     }
622:     if (pep->sfactor!=1.0) {
623:       for (j=0;j<pep->nconv;j++) {
624:         pep->eigr[j] *= pep->sfactor;
625:         pep->eigi[j] *= pep->sfactor;
626:       }
627:       /* restore original values */
628:       for (i=0;i<pep->nmat;i++){
629:         pep->pbc[pep->nmat+i] *= pep->sfactor;
630:         pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
631:       }
632:     }
633:   }
634:   /* restore original values */
635:   if (!flg) {
636:     pep->target *= pep->sfactor;
637:     STScaleShift(pep->st,pep->sfactor);
638:   } else {
639:     STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
640:     pep->target = (sinv)?pep->target/pep->sfactor:pep->target*pep->sfactor;
641:   }
642:   if (pep->sfactor!=1.0) { RGPopScale(pep->rg); }

644:   /* change the state to raw so that DSVectors() computes eigenvectors from scratch */
645:   DSSetDimensions(pep->ds,pep->nconv,0,0,0);
646:   DSSetState(pep->ds,DS_STATE_RAW);
647:   return(0);
648: }

650: static PetscErrorCode PEPTOARSetRestart_TOAR(PEP pep,PetscReal keep)
651: {
652:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;

655:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
656:   else {
657:     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]");
658:     ctx->keep = keep;
659:   }
660:   return(0);
661: }

663: /*@
664:    PEPTOARSetRestart - Sets the restart parameter for the TOAR
665:    method, in particular the proportion of basis vectors that must be kept
666:    after restart.

668:    Logically Collective on pep

670:    Input Parameters:
671: +  pep  - the eigenproblem solver context
672: -  keep - the number of vectors to be kept at restart

674:    Options Database Key:
675: .  -pep_toar_restart - Sets the restart parameter

677:    Notes:
678:    Allowed values are in the range [0.1,0.9]. The default is 0.5.

680:    Level: advanced

682: .seealso: PEPTOARGetRestart()
683: @*/
684: PetscErrorCode PEPTOARSetRestart(PEP pep,PetscReal keep)
685: {

691:   PetscTryMethod(pep,"PEPTOARSetRestart_C",(PEP,PetscReal),(pep,keep));
692:   return(0);
693: }

695: static PetscErrorCode PEPTOARGetRestart_TOAR(PEP pep,PetscReal *keep)
696: {
697:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;

700:   *keep = ctx->keep;
701:   return(0);
702: }

704: /*@
705:    PEPTOARGetRestart - Gets the restart parameter used in the TOAR method.

707:    Not Collective

709:    Input Parameter:
710: .  pep - the eigenproblem solver context

712:    Output Parameter:
713: .  keep - the restart parameter

715:    Level: advanced

717: .seealso: PEPTOARSetRestart()
718: @*/
719: PetscErrorCode PEPTOARGetRestart(PEP pep,PetscReal *keep)
720: {

726:   PetscUseMethod(pep,"PEPTOARGetRestart_C",(PEP,PetscReal*),(pep,keep));
727:   return(0);
728: }

730: static PetscErrorCode PEPTOARSetLocking_TOAR(PEP pep,PetscBool lock)
731: {
732:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;

735:   ctx->lock = lock;
736:   return(0);
737: }

739: /*@
740:    PEPTOARSetLocking - Choose between locking and non-locking variants of
741:    the TOAR method.

743:    Logically Collective on pep

745:    Input Parameters:
746: +  pep  - the eigenproblem solver context
747: -  lock - true if the locking variant must be selected

749:    Options Database Key:
750: .  -pep_toar_locking - Sets the locking flag

752:    Notes:
753:    The default is to lock converged eigenpairs when the method restarts.
754:    This behaviour can be changed so that all directions are kept in the
755:    working subspace even if already converged to working accuracy (the
756:    non-locking variant).

758:    Level: advanced

760: .seealso: PEPTOARGetLocking()
761: @*/
762: PetscErrorCode PEPTOARSetLocking(PEP pep,PetscBool lock)
763: {

769:   PetscTryMethod(pep,"PEPTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
770:   return(0);
771: }

773: static PetscErrorCode PEPTOARGetLocking_TOAR(PEP pep,PetscBool *lock)
774: {
775:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;

778:   *lock = ctx->lock;
779:   return(0);
780: }

782: /*@
783:    PEPTOARGetLocking - Gets the locking flag used in the TOAR method.

785:    Not Collective

787:    Input Parameter:
788: .  pep - the eigenproblem solver context

790:    Output Parameter:
791: .  lock - the locking flag

793:    Level: advanced

795: .seealso: PEPTOARSetLocking()
796: @*/
797: PetscErrorCode PEPTOARGetLocking(PEP pep,PetscBool *lock)
798: {

804:   PetscUseMethod(pep,"PEPTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
805:   return(0);
806: }

808: PetscErrorCode PEPSetFromOptions_TOAR(PetscOptionItems *PetscOptionsObject,PEP pep)
809: {
811:   PetscBool      flg,lock;
812:   PetscReal      keep;

815:   PetscOptionsHead(PetscOptionsObject,"PEP TOAR Options");

817:     PetscOptionsReal("-pep_toar_restart","Proportion of vectors kept after restart","PEPTOARSetRestart",0.5,&keep,&flg);
818:     if (flg) { PEPTOARSetRestart(pep,keep); }

820:     PetscOptionsBool("-pep_toar_locking","Choose between locking and non-locking variants","PEPTOARSetLocking",PETSC_FALSE,&lock,&flg);
821:     if (flg) { PEPTOARSetLocking(pep,lock); }

823:   PetscOptionsTail();
824:   return(0);
825: }

827: PetscErrorCode PEPView_TOAR(PEP pep,PetscViewer viewer)
828: {
830:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
831:   PetscBool      isascii;

834:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
835:   if (isascii) {
836:     PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
837:     PetscViewerASCIIPrintf(viewer,"  using the %slocking variant\n",ctx->lock?"":"non-");
838:   }
839:   return(0);
840: }

842: PetscErrorCode PEPDestroy_TOAR(PEP pep)
843: {
845:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;

848:   BVDestroy(&ctx->V);
849:   PetscFree(pep->data);
850:   PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",NULL);
851:   PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",NULL);
852:   PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",NULL);
853:   PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",NULL);
854:   return(0);
855: }

857: SLEPC_EXTERN PetscErrorCode PEPCreate_TOAR(PEP pep)
858: {
859:   PEP_TOAR       *ctx;

863:   PetscNewLog(pep,&ctx);
864:   pep->data = (void*)ctx;
865:   ctx->lock = PETSC_TRUE;

867:   pep->ops->solve          = PEPSolve_TOAR;
868:   pep->ops->setup          = PEPSetUp_TOAR;
869:   pep->ops->setfromoptions = PEPSetFromOptions_TOAR;
870:   pep->ops->destroy        = PEPDestroy_TOAR;
871:   pep->ops->view           = PEPView_TOAR;
872:   pep->ops->backtransform  = PEPBackTransform_Default;
873:   pep->ops->computevectors = PEPComputeVectors_Default;
874:   pep->ops->extractvectors = PEPExtractVectors_TOAR;

876:   PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",PEPTOARSetRestart_TOAR);
877:   PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",PEPTOARGetRestart_TOAR);
878:   PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",PEPTOARSetLocking_TOAR);
879:   PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",PEPTOARGetLocking_TOAR);
880:   return(0);
881: }