Actual source code: stoar.c

slepc-3.18.2 2023-01-26
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: "stoar"

 13:    Method: S-TOAR

 15:    Algorithm:

 17:        Symmetric Two-Level Orthogonal Arnoldi.

 19:    References:

 21:        [1] C. Campos and J.E. Roman, "Restarted Q-Arnoldi-type methods
 22:            exploiting symmetry in quadratic eigenvalue problems", BIT
 23:            Numer. Math. 56(4):1213-1236, 2016.
 24: */

 26: #include <slepc/private/pepimpl.h>
 27: #include "../src/pep/impls/krylov/pepkrylov.h"
 28: #include <slepcblaslapack.h>

 30: static PetscBool  cited = PETSC_FALSE;
 31: static const char citation[] =
 32:   "@Article{slepc-stoar,\n"
 33:   "   author = \"C. Campos and J. E. Roman\",\n"
 34:   "   title = \"Restarted {Q-Arnoldi-type} methods exploiting symmetry in quadratic eigenvalue problems\",\n"
 35:   "   journal = \"{BIT} Numer. Math.\",\n"
 36:   "   volume = \"56\",\n"
 37:   "   number = \"4\",\n"
 38:   "   pages = \"1213--1236\",\n"
 39:   "   year = \"2016,\"\n"
 40:   "   doi = \"https://doi.org/10.1007/s10543-016-0601-5\"\n"
 41:   "}\n";

 43: typedef struct {
 44:   PetscReal   scal[2];
 45:   Mat         A[2];
 46:   Vec         t;
 47: } PEP_STOAR_MATSHELL;

 49: static PetscErrorCode MatMult_STOAR(Mat A,Vec x,Vec y)
 50: {
 51:   PEP_STOAR_MATSHELL *ctx;

 53:   MatShellGetContext(A,&ctx);
 54:   MatMult(ctx->A[0],x,y);
 55:   VecScale(y,ctx->scal[0]);
 56:   if (ctx->scal[1]) {
 57:     MatMult(ctx->A[1],x,ctx->t);
 58:     VecAXPY(y,ctx->scal[1],ctx->t);
 59:   }
 60:   return 0;
 61: }

 63: static PetscErrorCode MatDestroy_STOAR(Mat A)
 64: {
 65:   PEP_STOAR_MATSHELL *ctx;

 67:   MatShellGetContext(A,&ctx);
 68:   VecDestroy(&ctx->t);
 69:   PetscFree(ctx);
 70:   return 0;
 71: }

 73: PetscErrorCode PEPSTOARSetUpInnerMatrix(PEP pep,Mat *B)
 74: {
 75:   Mat                pB[4],Bs[3],D[3];
 76:   PetscInt           i,j,n,m;
 77:   PEP_STOAR_MATSHELL *ctxMat[3];
 78:   PEP_STOAR          *ctx=(PEP_STOAR*)pep->data;

 80:   for (i=0;i<3;i++) {
 81:     STGetMatrixTransformed(pep->st,i,&D[i]); /* D[2] = M */
 82:   }
 83:   MatGetLocalSize(D[2],&m,&n);

 85:   for (j=0;j<3;j++) {
 86:     PetscNew(ctxMat+j);
 87:     MatCreateShell(PetscObjectComm((PetscObject)pep),m,n,PETSC_DETERMINE,PETSC_DETERMINE,ctxMat[j],&Bs[j]);
 88:     MatShellSetOperation(Bs[j],MATOP_MULT,(void(*)(void))MatMult_STOAR);
 89:     MatShellSetOperation(Bs[j],MATOP_DESTROY,(void(*)(void))MatDestroy_STOAR);
 90:   }
 91:   for (i=0;i<4;i++) pB[i] = NULL;
 92:   if (ctx->alpha) {
 93:     ctxMat[0]->A[0] = D[0]; ctxMat[0]->scal[0] = ctx->alpha; ctxMat[0]->scal[1] = 0.0;
 94:     ctxMat[2]->A[0] = D[2]; ctxMat[2]->scal[0] = -ctx->alpha*pep->sfactor*pep->sfactor; ctxMat[2]->scal[1] = 0.0;
 95:     pB[0] = Bs[0]; pB[3] = Bs[2];
 96:   }
 97:   if (ctx->beta) {
 98:     i = (ctx->alpha)?1:0;
 99:     ctxMat[0]->scal[1] = 0.0;
100:     ctxMat[0]->A[i] = D[1]; ctxMat[0]->scal[i] = -ctx->beta*pep->sfactor;
101:     ctxMat[1]->A[0] = D[2]; ctxMat[1]->scal[0] = -ctx->beta*pep->sfactor*pep->sfactor; ctxMat[1]->scal[1] = 0.0;
102:     pB[0] = Bs[0]; pB[1] = pB[2] = Bs[1];
103:   }
104:   BVCreateVec(pep->V,&ctxMat[0]->t);
105:   MatCreateNest(PetscObjectComm((PetscObject)pep),2,NULL,2,NULL,pB,B);
106:   for (j=0;j<3;j++) MatDestroy(&Bs[j]);
107:   return 0;
108: }

110: PetscErrorCode PEPSetUp_STOAR(PEP pep)
111: {
112:   PetscBool         sinv,flg;
113:   PEP_STOAR         *ctx = (PEP_STOAR*)pep->data;
114:   PetscInt          ld,i;
115:   PetscReal         eta;
116:   BVOrthogType      otype;
117:   BVOrthogBlockType obtype;

119:   PEPCheckHermitian(pep);
120:   PEPCheckQuadratic(pep);
121:   PEPCheckShiftSinvert(pep);
122:   /* spectrum slicing requires special treatment of default values */
123:   if (pep->which==PEP_ALL) {
124:     pep->ops->solve = PEPSolve_STOAR_QSlice;
125:     pep->ops->extractvectors = NULL;
126:     pep->ops->setdefaultst   = NULL;
127:     PEPSetUp_STOAR_QSlice(pep);
128:   } else {
129:     PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
131:     if (pep->max_it==PETSC_DEFAULT) pep->max_it = PetscMax(100,2*(pep->nmat-1)*pep->n/pep->ncv);
132:     pep->ops->solve = PEPSolve_STOAR;
133:     ld   = pep->ncv+2;
134:     DSSetType(pep->ds,DSGHIEP);
135:     DSSetCompact(pep->ds,PETSC_TRUE);
136:     DSSetExtraRow(pep->ds,PETSC_TRUE);
137:     DSAllocate(pep->ds,ld);
138:     PEPBasisCoefficients(pep,pep->pbc);
139:     STGetTransform(pep->st,&flg);
140:     if (!flg) {
141:       PetscFree(pep->solvematcoeffs);
142:       PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
143:       PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
144:       if (sinv) PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL);
145:       else {
146:         for (i=0;i<pep->nmat-1;i++) pep->solvematcoeffs[i] = 0.0;
147:         pep->solvematcoeffs[pep->nmat-1] = 1.0;
148:       }
149:     }
150:   }
151:   if (!pep->which) PEPSetWhichEigenpairs_Default(pep);
152:   PEPCheckUnsupported(pep,PEP_FEATURE_NONMONOMIAL | PEP_FEATURE_REGION);

154:   PEPAllocateSolution(pep,2);
155:   PEPSetWorkVecs(pep,4);
156:   BVDestroy(&ctx->V);
157:   BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
158:   BVGetOrthogonalization(pep->V,&otype,NULL,&eta,&obtype);
159:   BVSetOrthogonalization(ctx->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
160:   return 0;
161: }

163: /*
164:   Compute a run of Lanczos iterations. dim(work)=(ctx->ld)*4
165: */
166: static PetscErrorCode PEPSTOARrun(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,Vec *t_)
167: {
168:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;
169:   PetscInt       i,j,m=*M,l,lock;
170:   PetscInt       lds,d,ld,offq,nqt,ldds;
171:   Vec            v=t_[0],t=t_[1],q=t_[2];
172:   PetscReal      norm,sym=0.0,fro=0.0,*f;
173:   PetscScalar    *y,*S,*x,sigma;
174:   PetscBLASInt   j_,one=1;
175:   PetscBool      lindep,flg,sinvert=PETSC_FALSE;
176:   Mat            MS;

178:   PetscMalloc1(*M,&y);
179:   BVGetSizes(pep->V,NULL,NULL,&ld);
180:   BVTensorGetDegree(ctx->V,&d);
181:   BVGetActiveColumns(pep->V,&lock,&nqt);
182:   lds = d*ld;
183:   offq = ld;
184:   DSGetLeadingDimension(pep->ds,&ldds);
185:   *breakdown = PETSC_FALSE; /* ----- */
186:   DSGetDimensions(pep->ds,NULL,&l,NULL,NULL);
187:   BVSetActiveColumns(ctx->V,0,m);
188:   BVSetActiveColumns(pep->V,0,nqt);
189:   STGetTransform(pep->st,&flg);
190:   if (!flg) {
191:     /* spectral transformation handled by the solver */
192:     PetscObjectTypeCompareAny((PetscObject)pep->st,&flg,STSINVERT,STSHIFT,"");
194:     PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinvert);
195:     STGetShift(pep->st,&sigma);
196:   }
197:   for (j=k;j<m;j++) {
198:     /* apply operator */
199:     BVTensorGetFactors(ctx->V,NULL,&MS);
200:     MatDenseGetArray(MS,&S);
201:     BVGetColumn(pep->V,nqt,&t);
202:     BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
203:     if (!sinvert) {
204:       STMatMult(pep->st,0,v,q);
205:       BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
206:       STMatMult(pep->st,1,v,t);
207:       VecAXPY(q,pep->sfactor,t);
208:       if (ctx->beta && ctx->alpha) {
209:         STMatMult(pep->st,2,v,t);
210:         VecAXPY(q,-pep->sfactor*pep->sfactor*ctx->beta/ctx->alpha,t);
211:       }
212:       STMatSolve(pep->st,q,t);
213:       VecScale(t,-1.0/(pep->sfactor*pep->sfactor));
214:     } else {
215:       STMatMult(pep->st,1,v,q);
216:       STMatMult(pep->st,2,v,t);
217:       VecAXPY(q,sigma*pep->sfactor,t);
218:       VecScale(q,pep->sfactor);
219:       BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
220:       STMatMult(pep->st,2,v,t);
221:       VecAXPY(q,pep->sfactor*pep->sfactor,t);
222:       STMatSolve(pep->st,q,t);
223:       VecScale(t,-1.0);
224:     }
225:     BVRestoreColumn(pep->V,nqt,&t);

227:     /* orthogonalize */
228:     if (!sinvert) x = S+offq+(j+1)*lds;
229:     else x = S+(j+1)*lds;
230:     BVOrthogonalizeColumn(pep->V,nqt,x,&norm,&lindep);

232:     if (!lindep) {
233:       if (!sinvert) *(S+offq+(j+1)*lds+nqt) = norm;
234:       else *(S+(j+1)*lds+nqt) = norm;
235:       BVScaleColumn(pep->V,nqt,1.0/norm);
236:       nqt++;
237:     }
238:     if (!sinvert) {
239:       for (i=0;i<=nqt-1;i++) *(S+(j+1)*lds+i) = *(S+offq+j*lds+i);
240:       if (ctx->beta && ctx->alpha) {
241:         for (i=0;i<=nqt-1;i++) *(S+(j+1)*lds+offq+i) -= *(S+(j+1)*lds+i)*ctx->beta/ctx->alpha;
242:       }
243:     } else for (i=0;i<nqt;i++) *(S+(j+1)*lds+offq+i) = *(S+j*lds+i)+sigma*(*(S+(j+1)*lds+i));
244:     BVSetActiveColumns(pep->V,0,nqt);
245:     MatDenseRestoreArray(MS,&S);
246:     BVTensorRestoreFactors(ctx->V,NULL,&MS);

248:     /* level-2 orthogonalization */
249:     BVOrthogonalizeColumn(ctx->V,j+1,y,&norm,&lindep);
250:     a[j] = PetscRealPart(y[j]);
251:     omega[j+1] = (norm > 0)?1.0:-1.0;
252:     BVScaleColumn(ctx->V,j+1,1.0/norm);
253:     b[j] = PetscAbsReal(norm);

255:     /* check symmetry */
256:     DSGetArrayReal(pep->ds,DS_MAT_T,&f);
257:     if (j==k) {
258:       for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ldds+i]);
259:       for (i=0;i<l;i++) y[i] = 0.0;
260:     }
261:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
262:     if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
263:     PetscBLASIntCast(j,&j_);
264:     sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
265:     fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
266:     if (j>0) fro = SlepcAbs(fro,b[j-1]);
267:     if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
268:       *symmlost = PETSC_TRUE;
269:       *M=j;
270:       break;
271:     }
272:   }
273:   BVSetActiveColumns(pep->V,lock,nqt);
274:   BVSetActiveColumns(ctx->V,0,*M);
275:   PetscFree(y);
276:   return 0;
277: }

279: #if 0
280: static PetscErrorCode PEPSTOARpreKConvergence(PEP pep,PetscInt nv,PetscReal *norm,Vec *w)
281: {
282:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;
283:   PetscBLASInt   n_,one=1;
284:   PetscInt       lds=2*ctx->ld;
285:   PetscReal      t1,t2;
286:   PetscScalar    *S=ctx->S;

288:   PetscBLASIntCast(nv+2,&n_);
289:   t1 = BLASnrm2_(&n_,S+nv*2*ctx->ld,&one);
290:   t2 = BLASnrm2_(&n_,S+(nv*2+1)*ctx->ld,&one);
291:   *norm = SlepcAbs(t1,t2);
292:   BVSetActiveColumns(pep->V,0,nv+2);
293:   BVMultVec(pep->V,1.0,0.0,w[1],S+nv*lds);
294:   STMatMult(pep->st,0,w[1],w[2]);
295:   VecNorm(w[2],NORM_2,&t1);
296:   BVMultVec(pep->V,1.0,0.0,w[1],S+ctx->ld+nv*lds);
297:   STMatMult(pep->st,2,w[1],w[2]);
298:   VecNorm(w[2],NORM_2,&t2);
299:   t2 *= pep->sfactor*pep->sfactor;
300:   *norm = PetscMax(*norm,SlepcAbs(t1,t2));
301:   return 0;
302: }
303: #endif

305: PetscErrorCode PEPSolve_STOAR(PEP pep)
306: {
307:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;
308:   PetscInt       j,k,l,nv=0,ld,ldds,t,nq=0;
309:   PetscInt       nconv=0,deg=pep->nmat-1;
310:   PetscScalar    sigma;
311:   PetscReal      beta,norm=1.0,*omega,*a,*b;
312:   PetscBool      breakdown,symmlost=PETSC_FALSE,sinv=PETSC_FALSE,falselock=PETSC_TRUE,flg;
313:   Mat            MQ,A,D;
314:   Vec            vomega;

316:   PetscCitationsRegister(citation,&cited);
317:   PEPSTOARSetUpInnerMatrix(pep,&A);
318:   BVSetMatrix(ctx->V,A,PETSC_TRUE);
319:   MatDestroy(&A);
320:   if (ctx->lock) {
321:     /* undocumented option to use a cheaper locking instead of the true locking */
322:     PetscOptionsGetBool(NULL,NULL,"-pep_stoar_falselocking",&falselock,NULL);
323:   }
324:   BVGetSizes(pep->V,NULL,NULL,&ld);
325:   STGetShift(pep->st,&sigma);
326:   STGetTransform(pep->st,&flg);
327:   if (pep->sfactor!=1.0) {
328:     if (!flg) {
329:       pep->target /= pep->sfactor;
330:       RGPushScale(pep->rg,1.0/pep->sfactor);
331:       STScaleShift(pep->st,1.0/pep->sfactor);
332:       sigma /= pep->sfactor;
333:     } else {
334:       PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
335:       pep->target = sinv?pep->target*pep->sfactor:pep->target/pep->sfactor;
336:       RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
337:       STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
338:     }
339:   }
340:   if (flg) sigma = 0.0;

342:   /* Get the starting Arnoldi vector */
343:   BVTensorBuildFirstColumn(ctx->V,pep->nini);
344:   DSSetDimensions(pep->ds,1,PETSC_DEFAULT,PETSC_DEFAULT);
345:   BVSetActiveColumns(ctx->V,0,1);
346:   DSGetMatAndColumn(pep->ds,DS_MAT_D,0,&D,&vomega);
347:   BVGetSignature(ctx->V,vomega);
348:   DSRestoreMatAndColumn(pep->ds,DS_MAT_D,0,&D,&vomega);

350:   /* Restart loop */
351:   l = 0;
352:   DSGetLeadingDimension(pep->ds,&ldds);
353:   while (pep->reason == PEP_CONVERGED_ITERATING) {
354:     pep->its++;
355:     DSGetArrayReal(pep->ds,DS_MAT_T,&a);
356:     b = a+ldds;
357:     DSGetArrayReal(pep->ds,DS_MAT_D,&omega);

359:     /* Compute an nv-step Lanczos factorization */
360:     nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
361:     DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l);
362:     PEPSTOARrun(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,pep->work);
363:     beta = b[nv-1];
364:     if (symmlost && nv==pep->nconv+l) {
365:       pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
366:       pep->nconv = nconv;
367:       if (falselock || !ctx->lock) {
368:        BVSetActiveColumns(ctx->V,0,pep->nconv);
369:        BVTensorCompress(ctx->V,0);
370:       }
371:       break;
372:     }
373:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
374:     DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
375:     DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l);
376:     if (l==0) DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
377:     else DSSetState(pep->ds,DS_STATE_RAW);

379:     /* Solve projected problem */
380:     DSSolve(pep->ds,pep->eigr,pep->eigi);
381:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
382:     DSUpdateExtraRow(pep->ds);
383:     DSSynchronize(pep->ds,pep->eigr,pep->eigi);

385:     /* Check convergence */
386:     /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
387:     norm = 1.0;
388:     DSGetDimensions(pep->ds,NULL,NULL,NULL,&t);
389:     PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
390:     (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);

392:     /* Update l */
393:     if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
394:     else {
395:       l = PetscMax(1,(PetscInt)((nv-k)/2));
396:       l = PetscMin(l,t);
397:       DSGetTruncateSize(pep->ds,k,t,&l);
398:       if (!breakdown) {
399:         /* Prepare the Rayleigh quotient for restart */
400:         DSTruncate(pep->ds,k+l,PETSC_FALSE);
401:       }
402:     }
403:     nconv = k;
404:     if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
405:     if (l) PetscInfo(pep,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l);

407:     /* Update S */
408:     DSGetMat(pep->ds,DS_MAT_Q,&MQ);
409:     BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
410:     DSRestoreMat(pep->ds,DS_MAT_Q,&MQ);

412:     /* Copy last column of S */
413:     BVCopyColumn(ctx->V,nv,k+l);
414:     BVSetActiveColumns(ctx->V,0,k+l);
415:     DSSetDimensions(pep->ds,k+l,PETSC_DEFAULT,PETSC_DEFAULT);
416:     DSGetMatAndColumn(pep->ds,DS_MAT_D,0,&D,&vomega);
417:     BVSetSignature(ctx->V,vomega);
418:     DSRestoreMatAndColumn(pep->ds,DS_MAT_D,0,&D,&vomega);

420:     if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
421:       /* stop if breakdown */
422:       PetscInfo(pep,"Breakdown TOAR method (it=%" PetscInt_FMT " norm=%g)\n",pep->its,(double)beta);
423:       pep->reason = PEP_DIVERGED_BREAKDOWN;
424:     }
425:     if (pep->reason != PEP_CONVERGED_ITERATING) l--;
426:     BVGetActiveColumns(pep->V,NULL,&nq);
427:     if (k+l+deg<=nq) {
428:       BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
429:       if (!falselock && ctx->lock) BVTensorCompress(ctx->V,k-pep->nconv);
430:       else BVTensorCompress(ctx->V,0);
431:     }
432:     pep->nconv = k;
433:     PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
434:   }

436:   if (pep->nconv>0) {
437:     BVSetActiveColumns(ctx->V,0,pep->nconv);
438:     BVGetActiveColumns(pep->V,NULL,&nq);
439:     BVSetActiveColumns(pep->V,0,nq);
440:     if (nq>pep->nconv) {
441:       BVTensorCompress(ctx->V,pep->nconv);
442:       BVSetActiveColumns(pep->V,0,pep->nconv);
443:     }
444:   }
445:   STGetTransform(pep->st,&flg);
446:   if (!flg) PetscTryTypeMethod(pep,backtransform);
447:   if (pep->sfactor!=1.0) {
448:     for (j=0;j<pep->nconv;j++) {
449:       pep->eigr[j] *= pep->sfactor;
450:       pep->eigi[j] *= pep->sfactor;
451:     }
452:   }
453:   /* restore original values */
454:   if (!flg) {
455:     pep->target *= pep->sfactor;
456:     STScaleShift(pep->st,pep->sfactor);
457:   } else {
458:     STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
459:     pep->target = (sinv)?pep->target/pep->sfactor:pep->target*pep->sfactor;
460:   }
461:   if (pep->sfactor!=1.0) RGPopScale(pep->rg);

463:   DSTruncate(pep->ds,pep->nconv,PETSC_TRUE);
464:   return 0;
465: }

467: PetscErrorCode PEPSetFromOptions_STOAR(PEP pep,PetscOptionItems *PetscOptionsObject)
468: {
469:   PetscBool      flg,lock,b,f1,f2,f3;
470:   PetscInt       i,j,k;
471:   PetscReal      array[2]={0,0};
472:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;

474:   PetscOptionsHeadBegin(PetscOptionsObject,"PEP STOAR Options");

476:     PetscOptionsBool("-pep_stoar_locking","Choose between locking and non-locking variants","PEPSTOARSetLocking",PETSC_FALSE,&lock,&flg);
477:     if (flg) PEPSTOARSetLocking(pep,lock);

479:     b = ctx->detect;
480:     PetscOptionsBool("-pep_stoar_detect_zeros","Check zeros during factorizations at interval boundaries","PEPSTOARSetDetectZeros",ctx->detect,&b,&flg);
481:     if (flg) PEPSTOARSetDetectZeros(pep,b);

483:     i = 1;
484:     j = k = PETSC_DECIDE;
485:     PetscOptionsInt("-pep_stoar_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","PEPSTOARSetDimensions",20,&i,&f1);
486:     PetscOptionsInt("-pep_stoar_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","PEPSTOARSetDimensions",40,&j,&f2);
487:     PetscOptionsInt("-pep_stoar_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","PEPSTOARSetDimensions",40,&k,&f3);
488:     if (f1 || f2 || f3) PEPSTOARSetDimensions(pep,i,j,k);

490:     k = 2;
491:     PetscOptionsRealArray("-pep_stoar_linearization","Parameters of the linearization","PEPSTOARSetLinearization",array,&k,&flg);
492:     if (flg) PEPSTOARSetLinearization(pep,array[0],array[1]);

494:     b = ctx->checket;
495:     PetscOptionsBool("-pep_stoar_check_eigenvalue_type","Check eigenvalue type during spectrum slicing","PEPSTOARSetCheckEigenvalueType",ctx->checket,&b,&flg);
496:     if (flg) PEPSTOARSetCheckEigenvalueType(pep,b);

498:   PetscOptionsHeadEnd();
499:   return 0;
500: }

502: static PetscErrorCode PEPSTOARSetLocking_STOAR(PEP pep,PetscBool lock)
503: {
504:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

506:   ctx->lock = lock;
507:   return 0;
508: }

510: /*@
511:    PEPSTOARSetLocking - Choose between locking and non-locking variants of
512:    the STOAR method.

514:    Logically Collective on pep

516:    Input Parameters:
517: +  pep  - the eigenproblem solver context
518: -  lock - true if the locking variant must be selected

520:    Options Database Key:
521: .  -pep_stoar_locking - Sets the locking flag

523:    Notes:
524:    The default is to lock converged eigenpairs when the method restarts.
525:    This behaviour can be changed so that all directions are kept in the
526:    working subspace even if already converged to working accuracy (the
527:    non-locking variant).

529:    Level: advanced

531: .seealso: PEPSTOARGetLocking()
532: @*/
533: PetscErrorCode PEPSTOARSetLocking(PEP pep,PetscBool lock)
534: {
537:   PetscTryMethod(pep,"PEPSTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
538:   return 0;
539: }

541: static PetscErrorCode PEPSTOARGetLocking_STOAR(PEP pep,PetscBool *lock)
542: {
543:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

545:   *lock = ctx->lock;
546:   return 0;
547: }

549: /*@
550:    PEPSTOARGetLocking - Gets the locking flag used in the STOAR method.

552:    Not Collective

554:    Input Parameter:
555: .  pep - the eigenproblem solver context

557:    Output Parameter:
558: .  lock - the locking flag

560:    Level: advanced

562: .seealso: PEPSTOARSetLocking()
563: @*/
564: PetscErrorCode PEPSTOARGetLocking(PEP pep,PetscBool *lock)
565: {
568:   PetscUseMethod(pep,"PEPSTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
569:   return 0;
570: }

572: static PetscErrorCode PEPSTOARGetInertias_STOAR(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
573: {
574:   PetscInt       i,numsh;
575:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;
576:   PEP_SR         sr = ctx->sr;

580:   switch (pep->state) {
581:   case PEP_STATE_INITIAL:
582:     break;
583:   case PEP_STATE_SETUP:
584:     if (n) *n = 2;
585:     if (shifts) {
586:       PetscMalloc1(2,shifts);
587:       (*shifts)[0] = pep->inta;
588:       (*shifts)[1] = pep->intb;
589:     }
590:     if (inertias) {
591:       PetscMalloc1(2,inertias);
592:       (*inertias)[0] = (sr->dir==1)?sr->inertia0:sr->inertia1;
593:       (*inertias)[1] = (sr->dir==1)?sr->inertia1:sr->inertia0;
594:     }
595:     break;
596:   case PEP_STATE_SOLVED:
597:   case PEP_STATE_EIGENVECTORS:
598:     numsh = ctx->nshifts;
599:     if (n) *n = numsh;
600:     if (shifts) {
601:       PetscMalloc1(numsh,shifts);
602:       for (i=0;i<numsh;i++) (*shifts)[i] = ctx->shifts[i];
603:     }
604:     if (inertias) {
605:       PetscMalloc1(numsh,inertias);
606:       for (i=0;i<numsh;i++) (*inertias)[i] = ctx->inertias[i];
607:     }
608:     break;
609:   }
610:   return 0;
611: }

613: /*@C
614:    PEPSTOARGetInertias - Gets the values of the shifts and their
615:    corresponding inertias in case of doing spectrum slicing for a
616:    computational interval.

618:    Not Collective

620:    Input Parameter:
621: .  pep - the eigenproblem solver context

623:    Output Parameters:
624: +  n        - number of shifts, including the endpoints of the interval
625: .  shifts   - the values of the shifts used internally in the solver
626: -  inertias - the values of the inertia in each shift

628:    Notes:
629:    If called after PEPSolve(), all shifts used internally by the solver are
630:    returned (including both endpoints and any intermediate ones). If called
631:    before PEPSolve() and after PEPSetUp() then only the information of the
632:    endpoints of subintervals is available.

634:    This function is only available for spectrum slicing runs.

636:    The returned arrays should be freed by the user. Can pass NULL in any of
637:    the two arrays if not required.

639:    Fortran Notes:
640:    The calling sequence from Fortran is
641: .vb
642:    PEPSTOARGetInertias(pep,n,shifts,inertias,ierr)
643:    integer n
644:    double precision shifts(*)
645:    integer inertias(*)
646: .ve
647:    The arrays should be at least of length n. The value of n can be determined
648:    by an initial call
649: .vb
650:    PEPSTOARGetInertias(pep,n,PETSC_NULL_REAL,PETSC_NULL_INTEGER,ierr)
651: .ve

653:    Level: advanced

655: .seealso: PEPSetInterval()
656: @*/
657: PetscErrorCode PEPSTOARGetInertias(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
658: {
661:   PetscUseMethod(pep,"PEPSTOARGetInertias_C",(PEP,PetscInt*,PetscReal**,PetscInt**),(pep,n,shifts,inertias));
662:   return 0;
663: }

665: static PetscErrorCode PEPSTOARSetDetectZeros_STOAR(PEP pep,PetscBool detect)
666: {
667:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

669:   ctx->detect = detect;
670:   pep->state  = PEP_STATE_INITIAL;
671:   return 0;
672: }

674: /*@
675:    PEPSTOARSetDetectZeros - Sets a flag to enforce detection of
676:    zeros during the factorizations throughout the spectrum slicing computation.

678:    Logically Collective on pep

680:    Input Parameters:
681: +  pep    - the eigenproblem solver context
682: -  detect - check for zeros

684:    Options Database Key:
685: .  -pep_stoar_detect_zeros - Check for zeros; this takes an optional
686:    bool value (0/1/no/yes/true/false)

688:    Notes:
689:    A zero in the factorization indicates that a shift coincides with an eigenvalue.

691:    This flag is turned off by default, and may be necessary in some cases.
692:    This feature currently requires an external package for factorizations
693:    with support for zero detection, e.g. MUMPS.

695:    Level: advanced

697: .seealso: PEPSetInterval()
698: @*/
699: PetscErrorCode PEPSTOARSetDetectZeros(PEP pep,PetscBool detect)
700: {
703:   PetscTryMethod(pep,"PEPSTOARSetDetectZeros_C",(PEP,PetscBool),(pep,detect));
704:   return 0;
705: }

707: static PetscErrorCode PEPSTOARGetDetectZeros_STOAR(PEP pep,PetscBool *detect)
708: {
709:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

711:   *detect = ctx->detect;
712:   return 0;
713: }

715: /*@
716:    PEPSTOARGetDetectZeros - Gets the flag that enforces zero detection
717:    in spectrum slicing.

719:    Not Collective

721:    Input Parameter:
722: .  pep - the eigenproblem solver context

724:    Output Parameter:
725: .  detect - whether zeros detection is enforced during factorizations

727:    Level: advanced

729: .seealso: PEPSTOARSetDetectZeros()
730: @*/
731: PetscErrorCode PEPSTOARGetDetectZeros(PEP pep,PetscBool *detect)
732: {
735:   PetscUseMethod(pep,"PEPSTOARGetDetectZeros_C",(PEP,PetscBool*),(pep,detect));
736:   return 0;
737: }

739: static PetscErrorCode PEPSTOARSetLinearization_STOAR(PEP pep,PetscReal alpha,PetscReal beta)
740: {
741:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

744:   ctx->alpha = alpha;
745:   ctx->beta  = beta;
746:   return 0;
747: }

749: /*@
750:    PEPSTOARSetLinearization - Set the coefficients that define
751:    the linearization of a quadratic eigenproblem.

753:    Logically Collective on pep

755:    Input Parameters:
756: +  pep   - polynomial eigenvalue solver
757: .  alpha - first parameter of the linearization
758: -  beta  - second parameter of the linearization

760:    Options Database Key:
761: .  -pep_stoar_linearization <alpha,beta> - Sets the coefficients

763:    Notes:
764:    Cannot pass zero for both alpha and beta. The default values are
765:    alpha=1 and beta=0.

767:    Level: advanced

769: .seealso: PEPSTOARGetLinearization()
770: @*/
771: PetscErrorCode PEPSTOARSetLinearization(PEP pep,PetscReal alpha,PetscReal beta)
772: {
776:   PetscTryMethod(pep,"PEPSTOARSetLinearization_C",(PEP,PetscReal,PetscReal),(pep,alpha,beta));
777:   return 0;
778: }

780: static PetscErrorCode PEPSTOARGetLinearization_STOAR(PEP pep,PetscReal *alpha,PetscReal *beta)
781: {
782:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

784:   if (alpha) *alpha = ctx->alpha;
785:   if (beta)  *beta  = ctx->beta;
786:   return 0;
787: }

789: /*@
790:    PEPSTOARGetLinearization - Returns the coefficients that define
791:    the linearization of a quadratic eigenproblem.

793:    Not Collective

795:    Input Parameter:
796: .  pep  - polynomial eigenvalue solver

798:    Output Parameters:
799: +  alpha - the first parameter of the linearization
800: -  beta  - the second parameter of the linearization

802:    Level: advanced

804: .seealso: PEPSTOARSetLinearization()
805: @*/
806: PetscErrorCode PEPSTOARGetLinearization(PEP pep,PetscReal *alpha,PetscReal *beta)
807: {
809:   PetscUseMethod(pep,"PEPSTOARGetLinearization_C",(PEP,PetscReal*,PetscReal*),(pep,alpha,beta));
810:   return 0;
811: }

813: static PetscErrorCode PEPSTOARSetDimensions_STOAR(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
814: {
815:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

818:   ctx->nev = nev;
819:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
820:     ctx->ncv = PETSC_DEFAULT;
821:   } else {
823:     ctx->ncv = ncv;
824:   }
825:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
826:     ctx->mpd = PETSC_DEFAULT;
827:   } else {
829:     ctx->mpd = mpd;
830:   }
831:   pep->state = PEP_STATE_INITIAL;
832:   return 0;
833: }

835: /*@
836:    PEPSTOARSetDimensions - Sets the dimensions used for each subsolve
837:    step in case of doing spectrum slicing for a computational interval.
838:    The meaning of the parameters is the same as in PEPSetDimensions().

840:    Logically Collective on pep

842:    Input Parameters:
843: +  pep - the eigenproblem solver context
844: .  nev - number of eigenvalues to compute
845: .  ncv - the maximum dimension of the subspace to be used by the subsolve
846: -  mpd - the maximum dimension allowed for the projected problem

848:    Options Database Key:
849: +  -eps_stoar_nev <nev> - Sets the number of eigenvalues
850: .  -eps_stoar_ncv <ncv> - Sets the dimension of the subspace
851: -  -eps_stoar_mpd <mpd> - Sets the maximum projected dimension

853:    Level: advanced

855: .seealso: PEPSTOARGetDimensions(), PEPSetDimensions(), PEPSetInterval()
856: @*/
857: PetscErrorCode PEPSTOARSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
858: {
863:   PetscTryMethod(pep,"PEPSTOARSetDimensions_C",(PEP,PetscInt,PetscInt,PetscInt),(pep,nev,ncv,mpd));
864:   return 0;
865: }

867: static PetscErrorCode PEPSTOARGetDimensions_STOAR(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
868: {
869:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

871:   if (nev) *nev = ctx->nev;
872:   if (ncv) *ncv = ctx->ncv;
873:   if (mpd) *mpd = ctx->mpd;
874:   return 0;
875: }

877: /*@
878:    PEPSTOARGetDimensions - Gets the dimensions used for each subsolve
879:    step in case of doing spectrum slicing for a computational interval.

881:    Not Collective

883:    Input Parameter:
884: .  pep - the eigenproblem solver context

886:    Output Parameters:
887: +  nev - number of eigenvalues to compute
888: .  ncv - the maximum dimension of the subspace to be used by the subsolve
889: -  mpd - the maximum dimension allowed for the projected problem

891:    Level: advanced

893: .seealso: PEPSTOARSetDimensions()
894: @*/
895: PetscErrorCode PEPSTOARGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
896: {
898:   PetscUseMethod(pep,"PEPSTOARGetDimensions_C",(PEP,PetscInt*,PetscInt*,PetscInt*),(pep,nev,ncv,mpd));
899:   return 0;
900: }

902: static PetscErrorCode PEPSTOARSetCheckEigenvalueType_STOAR(PEP pep,PetscBool checket)
903: {
904:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

906:   ctx->checket = checket;
907:   pep->state   = PEP_STATE_INITIAL;
908:   return 0;
909: }

911: /*@
912:    PEPSTOARSetCheckEigenvalueType - Sets a flag to check that all the eigenvalues
913:    obtained throughout the spectrum slicing computation have the same definite type.

915:    Logically Collective on pep

917:    Input Parameters:
918: +  pep     - the eigenproblem solver context
919: -  checket - check eigenvalue type

921:    Options Database Key:
922: .  -pep_stoar_check_eigenvalue_type - Check eigenvalue type; this takes an optional
923:    bool value (0/1/no/yes/true/false)

925:    Notes:
926:    This option is relevant only for spectrum slicing computations, but it is
927:    ignored if the problem type is PEP_HYPERBOLIC.

929:    This flag is turned on by default, to guarantee that the computed eigenvalues
930:    have the same type (otherwise the computed solution might be wrong). But since
931:    the check is computationally quite expensive, the check may be turned off if
932:    the user knows for sure that all eigenvalues in the requested interval have
933:    the same type.

935:    Level: advanced

937: .seealso: PEPSetProblemType(), PEPSetInterval()
938: @*/
939: PetscErrorCode PEPSTOARSetCheckEigenvalueType(PEP pep,PetscBool checket)
940: {
943:   PetscTryMethod(pep,"PEPSTOARSetCheckEigenvalueType_C",(PEP,PetscBool),(pep,checket));
944:   return 0;
945: }

947: static PetscErrorCode PEPSTOARGetCheckEigenvalueType_STOAR(PEP pep,PetscBool *checket)
948: {
949:   PEP_STOAR *ctx = (PEP_STOAR*)pep->data;

951:   *checket = ctx->checket;
952:   return 0;
953: }

955: /*@
956:    PEPSTOARGetCheckEigenvalueType - Gets the flag for the eigenvalue type
957:    check in spectrum slicing.

959:    Not Collective

961:    Input Parameter:
962: .  pep - the eigenproblem solver context

964:    Output Parameter:
965: .  checket - whether eigenvalue type must be checked during spectrum slcing

967:    Level: advanced

969: .seealso: PEPSTOARSetCheckEigenvalueType()
970: @*/
971: PetscErrorCode PEPSTOARGetCheckEigenvalueType(PEP pep,PetscBool *checket)
972: {
975:   PetscUseMethod(pep,"PEPSTOARGetCheckEigenvalueType_C",(PEP,PetscBool*),(pep,checket));
976:   return 0;
977: }

979: PetscErrorCode PEPView_STOAR(PEP pep,PetscViewer viewer)
980: {
981:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;
982:   PetscBool      isascii;

984:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
985:   if (isascii) {
986:     PetscViewerASCIIPrintf(viewer,"  using the %slocking variant\n",ctx->lock?"":"non-");
987:     PetscViewerASCIIPrintf(viewer,"  linearization parameters: alpha=%g beta=%g\n",(double)ctx->alpha,(double)ctx->beta);
988:     if (pep->which==PEP_ALL && !ctx->hyperbolic) PetscViewerASCIIPrintf(viewer,"  checking eigenvalue type: %s\n",ctx->checket?"enabled":"disabled");
989:   }
990:   return 0;
991: }

993: PetscErrorCode PEPReset_STOAR(PEP pep)
994: {
995:   if (pep->which==PEP_ALL) PEPReset_STOAR_QSlice(pep);
996:   return 0;
997: }

999: PetscErrorCode PEPDestroy_STOAR(PEP pep)
1000: {
1001:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;

1003:   BVDestroy(&ctx->V);
1004:   PetscFree(pep->data);
1005:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",NULL);
1006:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",NULL);
1007:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDetectZeros_C",NULL);
1008:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDetectZeros_C",NULL);
1009:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetInertias_C",NULL);
1010:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDimensions_C",NULL);
1011:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDimensions_C",NULL);
1012:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLinearization_C",NULL);
1013:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLinearization_C",NULL);
1014:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetCheckEigenvalueType_C",NULL);
1015:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetCheckEigenvalueType_C",NULL);
1016:   return 0;
1017: }

1019: SLEPC_EXTERN PetscErrorCode PEPCreate_STOAR(PEP pep)
1020: {
1021:   PEP_STOAR      *ctx;

1023:   PetscNew(&ctx);
1024:   pep->data = (void*)ctx;

1026:   pep->lineariz = PETSC_TRUE;
1027:   ctx->lock     = PETSC_TRUE;
1028:   ctx->nev      = 1;
1029:   ctx->ncv      = PETSC_DEFAULT;
1030:   ctx->mpd      = PETSC_DEFAULT;
1031:   ctx->alpha    = 1.0;
1032:   ctx->beta     = 0.0;
1033:   ctx->checket  = PETSC_TRUE;

1035:   pep->ops->setup          = PEPSetUp_STOAR;
1036:   pep->ops->setfromoptions = PEPSetFromOptions_STOAR;
1037:   pep->ops->destroy        = PEPDestroy_STOAR;
1038:   pep->ops->view           = PEPView_STOAR;
1039:   pep->ops->backtransform  = PEPBackTransform_Default;
1040:   pep->ops->computevectors = PEPComputeVectors_Default;
1041:   pep->ops->extractvectors = PEPExtractVectors_TOAR;
1042:   pep->ops->reset          = PEPReset_STOAR;

1044:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",PEPSTOARSetLocking_STOAR);
1045:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",PEPSTOARGetLocking_STOAR);
1046:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDetectZeros_C",PEPSTOARSetDetectZeros_STOAR);
1047:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDetectZeros_C",PEPSTOARGetDetectZeros_STOAR);
1048:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetInertias_C",PEPSTOARGetInertias_STOAR);
1049:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDimensions_C",PEPSTOARGetDimensions_STOAR);
1050:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDimensions_C",PEPSTOARSetDimensions_STOAR);
1051:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLinearization_C",PEPSTOARSetLinearization_STOAR);
1052:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLinearization_C",PEPSTOARGetLinearization_STOAR);
1053:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetCheckEigenvalueType_C",PEPSTOARSetCheckEigenvalueType_STOAR);
1054:   PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetCheckEigenvalueType_C",PEPSTOARGetCheckEigenvalueType_STOAR);
1055:   return 0;
1056: }