Actual source code: qslice.c

slepc-3.18.1 2022-11-02
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 with spectrum slicing for symmetric quadratic eigenproblems

 15:    Algorithm:

 17:        Symmetric Two-Level Orthogonal Arnoldi.

 19:    References:

 21:        [1] C. Campos and J.E. Roman, "Inertia-based spectrum slicing
 22:            for symmetric quadratic eigenvalue problems", Numer. Linear
 23:            Algebra Appl. 27(4):e2293, 2020.
 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-slice-qep,\n"
 33:   "   author = \"C. Campos and J. E. Roman\",\n"
 34:   "   title = \"Inertia-based spectrum slicing for symmetric quadratic eigenvalue problems\",\n"
 35:   "   journal = \"Numer. Linear Algebra Appl.\",\n"
 36:   "   volume = \"27\",\n"
 37:   "   number = \"4\",\n"
 38:   "   pages = \"e2293\",\n"
 39:   "   year = \"2020,\"\n"
 40:   "   doi = \"https://doi.org/10.1002/nla.2293\"\n"
 41:   "}\n";

 43: #define SLICE_PTOL PETSC_SQRT_MACHINE_EPSILON

 45: static PetscErrorCode PEPQSliceResetSR(PEP pep)
 46: {
 47:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
 48:   PEP_SR         sr=ctx->sr;
 49:   PEP_shift      s;
 50:   PetscInt       i;

 52:   if (sr) {
 53:     /* Reviewing list of shifts to free memory */
 54:     s = sr->s0;
 55:     if (s) {
 56:       while (s->neighb[1]) {
 57:         s = s->neighb[1];
 58:         PetscFree(s->neighb[0]);
 59:       }
 60:       PetscFree(s);
 61:     }
 62:     PetscFree(sr->S);
 63:     for (i=0;i<pep->nconv;i++) PetscFree(sr->qinfo[i].q);
 64:     PetscFree(sr->qinfo);
 65:     for (i=0;i<3;i++) VecDestroy(&sr->v[i]);
 66:     EPSDestroy(&sr->eps);
 67:     PetscFree(sr);
 68:   }
 69:   ctx->sr = NULL;
 70:   return 0;
 71: }

 73: PetscErrorCode PEPReset_STOAR_QSlice(PEP pep)
 74: {
 75:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;

 77:   PEPQSliceResetSR(pep);
 78:   PetscFree(ctx->inertias);
 79:   PetscFree(ctx->shifts);
 80:   return 0;
 81: }

 83: /*
 84:   PEPQSliceAllocateSolution - Allocate memory storage for common variables such
 85:   as eigenvalues and eigenvectors.
 86: */
 87: static PetscErrorCode PEPQSliceAllocateSolution(PEP pep)
 88: {
 89:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
 90:   PetscInt       k;
 91:   BVType         type;
 92:   Vec            t;
 93:   PEP_SR         sr = ctx->sr;

 95:   /* allocate space for eigenvalues and friends */
 96:   k = PetscMax(1,sr->numEigs);
 97:   PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm);
 98:   PetscCalloc4(k,&sr->eigr,k,&sr->eigi,k,&sr->errest,k,&sr->perm);
 99:   PetscFree(sr->qinfo);
100:   PetscCalloc1(k,&sr->qinfo);

102:   /* allocate sr->V and transfer options from pep->V */
103:   BVDestroy(&sr->V);
104:   BVCreate(PetscObjectComm((PetscObject)pep),&sr->V);
105:   if (!pep->V) PEPGetBV(pep,&pep->V);
106:   if (!((PetscObject)(pep->V))->type_name) BVSetType(sr->V,BVSVEC);
107:   else {
108:     BVGetType(pep->V,&type);
109:     BVSetType(sr->V,type);
110:   }
111:   STMatCreateVecsEmpty(pep->st,&t,NULL);
112:   BVSetSizesFromVec(sr->V,t,k+1);
113:   VecDestroy(&t);
114:   sr->ld = k;
115:   PetscFree(sr->S);
116:   PetscMalloc1((k+1)*sr->ld*(pep->nmat-1),&sr->S);
117:   return 0;
118: }

120: /* Convergence test to compute positive Ritz values */
121: static PetscErrorCode ConvergedPositive(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
122: {
123:   *errest = (PetscRealPart(eigr)>0.0)?0.0:res;
124:   return 0;
125: }

127: static PetscErrorCode PEPQSliceMatGetInertia(PEP pep,PetscReal shift,PetscInt *inertia,PetscInt *zeros)
128: {
129:   KSP            ksp,kspr;
130:   PC             pc;
131:   Mat            F;
132:   PetscBool      flg;

134:   if (!pep->solvematcoeffs) PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
135:   if (shift==PETSC_MAX_REAL) { /* Inertia of matrix A[2] */
136:     pep->solvematcoeffs[0] = 0.0; pep->solvematcoeffs[1] = 0.0; pep->solvematcoeffs[2] = 1.0;
137:   } else PEPEvaluateBasis(pep,shift,0,pep->solvematcoeffs,NULL);
138:   STMatSetUp(pep->st,pep->sfactor,pep->solvematcoeffs);
139:   STGetKSP(pep->st,&ksp);
140:   KSPGetPC(ksp,&pc);
141:   PetscObjectTypeCompare((PetscObject)pc,PCREDUNDANT,&flg);
142:   if (flg) {
143:     PCRedundantGetKSP(pc,&kspr);
144:     KSPGetPC(kspr,&pc);
145:   }
146:   PCFactorGetMatrix(pc,&F);
147:   MatGetInertia(F,inertia,zeros,NULL);
148:   return 0;
149: }

151: static PetscErrorCode PEPQSliceGetInertia(PEP pep,PetscReal shift,PetscInt *inertia,PetscInt *zeros,PetscInt correction)
152: {
153:   KSP            ksp;
154:   Mat            P;
155:   PetscReal      nzshift=0.0,dot;
156:   PetscRandom    rand;
157:   PetscInt       nconv;
158:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
159:   PEP_SR         sr=ctx->sr;

161:   if (shift >= PETSC_MAX_REAL) { /* Right-open interval */
162:     *inertia = 0;
163:   } else if (shift <= PETSC_MIN_REAL) {
164:     *inertia = 0;
165:     if (zeros) *zeros = 0;
166:   } else {
167:     /* If the shift is zero, perturb it to a very small positive value.
168:        The goal is that the nonzero pattern is the same in all cases and reuse
169:        the symbolic factorizations */
170:     nzshift = (shift==0.0)? 10.0/PETSC_MAX_REAL: shift;
171:     PEPQSliceMatGetInertia(pep,nzshift,inertia,zeros);
172:     STSetShift(pep->st,nzshift);
173:   }
174:   if (!correction) {
175:     if (shift >= PETSC_MAX_REAL) *inertia = 2*pep->n;
176:     else if (shift>PETSC_MIN_REAL) {
177:       STGetKSP(pep->st,&ksp);
178:       KSPGetOperators(ksp,&P,NULL);
179:       if (*inertia!=pep->n && !sr->v[0]) {
180:         MatCreateVecs(P,&sr->v[0],NULL);
181:         VecDuplicate(sr->v[0],&sr->v[1]);
182:         VecDuplicate(sr->v[0],&sr->v[2]);
183:         BVGetRandomContext(pep->V,&rand);
184:         VecSetRandom(sr->v[0],rand);
185:       }
186:       if (*inertia<pep->n && *inertia>0) {
187:         if (!sr->eps) {
188:           EPSCreate(PetscObjectComm((PetscObject)pep),&sr->eps);
189:           EPSSetProblemType(sr->eps,EPS_HEP);
190:           EPSSetWhichEigenpairs(sr->eps,EPS_LARGEST_REAL);
191:         }
192:         EPSSetConvergenceTestFunction(sr->eps,ConvergedPositive,NULL,NULL);
193:         EPSSetOperators(sr->eps,P,NULL);
194:         EPSSolve(sr->eps);
195:         EPSGetConverged(sr->eps,&nconv);
197:         EPSGetEigenpair(sr->eps,0,NULL,NULL,sr->v[0],sr->v[1]);
198:       }
199:       if (*inertia!=pep->n) {
200:         MatMult(pep->A[1],sr->v[0],sr->v[1]);
201:         MatMult(pep->A[2],sr->v[0],sr->v[2]);
202:         VecAXPY(sr->v[1],2*nzshift,sr->v[2]);
203:         VecDotRealPart(sr->v[1],sr->v[0],&dot);
204:         if (dot>0.0) *inertia = 2*pep->n-*inertia;
205:       }
206:     }
207:   } else if (correction<0) *inertia = 2*pep->n-*inertia;
208:   return 0;
209: }

211: /*
212:    Check eigenvalue type - used only in non-hyperbolic problems.
213:    All computed eigenvalues must have the same definite type (positive or negative).
214:    If ini=TRUE the type is available in omega, otherwise we compute an eigenvalue
215:    closest to shift and determine its type.
216:  */
217: static PetscErrorCode PEPQSliceCheckEigenvalueType(PEP pep,PetscReal shift,PetscReal omega,PetscBool ini)
218: {
219:   PEP            pep2;
220:   ST             st;
221:   PetscInt       nconv;
222:   PetscScalar    lambda;
223:   PetscReal      dot;
224:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
225:   PEP_SR         sr=ctx->sr;

227:   if (!ini) {
229:   } else {
230:     PEPCreate(PetscObjectComm((PetscObject)pep),&pep2);
231:     PEPSetOptionsPrefix(pep2,((PetscObject)pep)->prefix);
232:     PEPAppendOptionsPrefix(pep2,"pep_eigenvalue_type_");
233:     PEPSetTolerances(pep2,PETSC_DEFAULT,pep->max_it/4);
234:     PEPSetType(pep2,PEPTOAR);
235:     PEPSetOperators(pep2,pep->nmat,pep->A);
236:     PEPSetWhichEigenpairs(pep2,PEP_TARGET_MAGNITUDE);
237:     PEPGetRG(pep2,&pep2->rg);
238:     RGSetType(pep2->rg,RGINTERVAL);
239: #if defined(PETSC_USE_COMPLEX)
240:     RGIntervalSetEndpoints(pep2->rg,pep->inta,pep->intb,-PETSC_SQRT_MACHINE_EPSILON,PETSC_SQRT_MACHINE_EPSILON);
241: #else
242:     RGIntervalSetEndpoints(pep2->rg,pep->inta,pep->intb,0.0,0.0);
243: #endif
244:     pep2->target = shift;
245:     st = pep2->st;
246:     pep2->st = pep->st;
247:     PEPSolve(pep2);
248:     PEPGetConverged(pep2,&nconv);
249:     if (nconv) {
250:       PEPGetEigenpair(pep2,0,&lambda,NULL,pep2->work[0],NULL);
251:       MatMult(pep->A[1],pep2->work[0],pep2->work[1]);
252:       MatMult(pep->A[2],pep2->work[0],pep2->work[2]);
253:       VecAXPY(pep2->work[1],2.0*lambda*pep->sfactor,pep2->work[2]);
254:       VecDotRealPart(pep2->work[1],pep2->work[0],&dot);
255:       PetscInfo(pep,"lambda=%g, %s type\n",(double)PetscRealPart(lambda),(dot>0.0)?"positive":"negative");
256:       if (!sr->type) sr->type = (dot>0.0)?1:-1;
258:     }
259:     pep2->st = st;
260:     PEPDestroy(&pep2);
261:   }
262:   return 0;
263: }

265: static inline PetscErrorCode PEPQSliceDiscriminant(PEP pep,Vec u,Vec w,PetscReal *d,PetscReal *smas,PetscReal *smenos)
266: {
267:   PetscReal      ap,bp,cp,dis;

269:   MatMult(pep->A[0],u,w);
270:   VecDotRealPart(w,u,&cp);
271:   MatMult(pep->A[1],u,w);
272:   VecDotRealPart(w,u,&bp);
273:   MatMult(pep->A[2],u,w);
274:   VecDotRealPart(w,u,&ap);
275:   dis = bp*bp-4*ap*cp;
276:   if (dis>=0.0 && smas) {
277:     if (ap>0) *smas = (-bp+PetscSqrtReal(dis))/(2*ap);
278:     else if (ap<0) *smas = (-bp-PetscSqrtReal(dis))/(2*ap);
279:     else {
280:       if (bp >0) *smas = -cp/bp;
281:       else *smas = PETSC_MAX_REAL;
282:     }
283:   }
284:   if (dis>=0.0 && smenos) {
285:     if (ap>0) *smenos = (-bp-PetscSqrtReal(dis))/(2*ap);
286:     else if (ap<0) *smenos = (-bp+PetscSqrtReal(dis))/(2*ap);
287:     else {
288:       if (bp<0) *smenos = -cp/bp;
289:       else *smenos = PETSC_MAX_REAL;
290:     }
291:   }
292:   if (d) *d = dis;
293:   return 0;
294: }

296: static inline PetscErrorCode PEPQSliceEvaluateQEP(PEP pep,PetscScalar x,Mat M,MatStructure str)
297: {
298:   MatCopy(pep->A[0],M,SAME_NONZERO_PATTERN);
299:   MatAXPY(M,x,pep->A[1],str);
300:   MatAXPY(M,x*x,pep->A[2],str);
301:   return 0;
302: }

304: /*@
305:    PEPCheckDefiniteQEP - Determines if a symmetric/Hermitian quadratic eigenvalue problem
306:    is definite or not.

308:    Logically Collective on pep

310:    Input Parameter:
311: .  pep  - eigensolver context

313:    Output Parameters:
314: +  xi - first computed parameter
315: .  mu - second computed parameter
316: .  definite - flag indicating that the problem is definite
317: -  hyperbolic - flag indicating that the problem is hyperbolic

319:    Notes:
320:    This function is intended for quadratic eigenvalue problems, Q(lambda)=A*lambda^2+B*lambda+C,
321:    with symmetric (or Hermitian) coefficient matrices A,B,C.

323:    On output, the flag 'definite' may have the values -1 (meaning that the QEP is not
324:    definite), 1 (if the problem is definite), or 0 if the algorithm was not able to
325:    determine whether the problem is definite or not.

327:    If definite=1, the output flag 'hyperbolic' informs in a similar way about whether the
328:    problem is hyperbolic or not.

330:    If definite=1, the computed values xi and mu satisfy Q(xi)<0 and Q(mu)>0, as
331:    obtained via the method proposed in [Niendorf and Voss, LAA 2010]. Furthermore, if
332:    hyperbolic=1 then only xi is computed.

334:    Level: advanced

336: .seealso: PEPSetProblemType()
337: @*/
338: PetscErrorCode PEPCheckDefiniteQEP(PEP pep,PetscReal *xi,PetscReal *mu,PetscInt *definite,PetscInt *hyperbolic)
339: {
340:   PetscRandom    rand;
341:   Vec            u,w;
342:   PetscReal      d=0.0,s=0.0,sp,mut=0.0,omg=0.0,omgp;
343:   PetscInt       k,its=10,hyp=0,check=0,nconv,inertia,n;
344:   Mat            M=NULL;
345:   MatStructure   str;
346:   EPS            eps;
347:   PetscBool      transform,ptypehyp;

350:   ptypehyp = (pep->problem_type==PEP_HYPERBOLIC)? PETSC_TRUE: PETSC_FALSE;
351:   if (!pep->st) PEPGetST(pep,&pep->st);
352:   PEPSetDefaultST(pep);
353:   STSetMatrices(pep->st,pep->nmat,pep->A);
354:   MatGetSize(pep->A[0],&n,NULL);
355:   STGetTransform(pep->st,&transform);
356:   STSetTransform(pep->st,PETSC_FALSE);
357:   STSetUp(pep->st);
358:   MatCreateVecs(pep->A[0],&u,&w);
359:   PEPGetBV(pep,&pep->V);
360:   BVGetRandomContext(pep->V,&rand);
361:   VecSetRandom(u,rand);
362:   VecNormalize(u,NULL);
363:   PEPQSliceDiscriminant(pep,u,w,&d,&s,NULL);
364:   if (d<0.0) check = -1;
365:   if (!check) {
366:     EPSCreate(PetscObjectComm((PetscObject)pep),&eps);
367:     EPSSetProblemType(eps,EPS_HEP);
368:     EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
369:     EPSSetTolerances(eps,PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON),PETSC_DECIDE);
370:     MatDuplicate(pep->A[0],MAT_DO_NOT_COPY_VALUES,&M);
371:     STGetMatStructure(pep->st,&str);
372:   }
373:   for (k=0;k<its&&!check;k++) {
374:     PEPQSliceEvaluateQEP(pep,s,M,str);
375:     EPSSetOperators(eps,M,NULL);
376:     EPSSolve(eps);
377:     EPSGetConverged(eps,&nconv);
378:     if (!nconv) break;
379:     EPSGetEigenpair(eps,0,NULL,NULL,u,w);
380:     sp = s;
381:     PEPQSliceDiscriminant(pep,u,w,&d,&s,&omg);
382:     if (d<0.0) {check = -1; break;}
383:     if (PetscAbsReal((s-sp)/s)<100*PETSC_MACHINE_EPSILON) break;
384:     if (s>sp) {hyp = -1;}
385:     mut = 2*s-sp;
386:     PEPQSliceMatGetInertia(pep,mut,&inertia,NULL);
387:     if (inertia == n) {check = 1; break;}
388:   }
389:   for (;k<its&&!check;k++) {
390:     mut = (s-omg)/2;
391:     PEPQSliceMatGetInertia(pep,mut,&inertia,NULL);
392:     if (inertia == n) {check = 1; break;}
393:     if (PetscAbsReal((s-omg)/omg)<100*PETSC_MACHINE_EPSILON) break;
394:     PEPQSliceEvaluateQEP(pep,omg,M,str);
395:     EPSSetOperators(eps,M,NULL);
396:     EPSSolve(eps);
397:     EPSGetConverged(eps,&nconv);
398:     if (!nconv) break;
399:     EPSGetEigenpair(eps,0,NULL,NULL,u,w);
400:     omgp = omg;
401:     PEPQSliceDiscriminant(pep,u,w,&d,NULL,&omg);
402:     if (d<0.0) {check = -1; break;}
403:     if (omg<omgp) hyp = -1;
404:   }
405:   if (check==1) *xi = mut;
407:   if (check==1 && hyp==0) {
408:     PEPQSliceMatGetInertia(pep,PETSC_MAX_REAL,&inertia,NULL);
409:     if (inertia == 0) hyp = 1;
410:     else hyp = -1;
411:   }
412:   if (check==1 && hyp!=1) {
413:     check = 0;
414:     EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
415:     for (;k<its&&!check;k++) {
416:       PEPQSliceEvaluateQEP(pep,s,M,str);
417:       EPSSetOperators(eps,M,NULL);
418:       EPSSolve(eps);
419:       EPSGetConverged(eps,&nconv);
420:       if (!nconv) break;
421:       EPSGetEigenpair(eps,0,NULL,NULL,u,w);
422:       sp = s;
423:       PEPQSliceDiscriminant(pep,u,w,&d,&s,&omg);
424:       if (d<0.0) {check = -1; break;}
425:       if (PetscAbsReal((s-sp)/s)<100*PETSC_MACHINE_EPSILON) break;
426:       mut = 2*s-sp;
427:       PEPQSliceMatGetInertia(pep,mut,&inertia,NULL);
428:       if (inertia == 0) {check = 1; break;}
429:     }
430:     for (;k<its&&!check;k++) {
431:       mut = (s-omg)/2;
432:       PEPQSliceMatGetInertia(pep,mut,&inertia,NULL);
433:       if (inertia == 0) {check = 1; break;}
434:       if (PetscAbsReal((s-omg)/omg)<100*PETSC_MACHINE_EPSILON) break;
435:       PEPQSliceEvaluateQEP(pep,omg,M,str);
436:       EPSSetOperators(eps,M,NULL);
437:       EPSSolve(eps);
438:       EPSGetConverged(eps,&nconv);
439:       if (!nconv) break;
440:       EPSGetEigenpair(eps,0,NULL,NULL,u,w);
441:       PEPQSliceDiscriminant(pep,u,w,&d,NULL,&omg);
442:       if (d<0.0) {check = -1; break;}
443:     }
444:   }
445:   if (check==1) *mu = mut;
446:   *definite = check;
447:   *hyperbolic = hyp;
448:   if (M) MatDestroy(&M);
449:   VecDestroy(&u);
450:   VecDestroy(&w);
451:   EPSDestroy(&eps);
452:   STSetTransform(pep->st,transform);
453:   return 0;
454: }

456: /*
457:    Dummy backtransform operation
458:  */
459: static PetscErrorCode PEPBackTransform_Skip(PEP pep)
460: {
461:   return 0;
462: }

464: PetscErrorCode PEPSetUp_STOAR_QSlice(PEP pep)
465: {
466:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
467:   PEP_SR         sr;
468:   PetscInt       ld,i,zeros=0;
469:   SlepcSC        sc;
470:   PetscReal      r;

472:   PEPCheckSinvertCayley(pep);
475:   PEPCheckUnsupportedCondition(pep,PEP_FEATURE_STOPPING,PETSC_TRUE," (with spectrum slicing)");
476:   if (pep->tol==PETSC_DEFAULT) {
477: #if defined(PETSC_USE_REAL_SINGLE)
478:     pep->tol = SLEPC_DEFAULT_TOL;
479: #else
480:     /* use tighter tolerance */
481:     pep->tol = SLEPC_DEFAULT_TOL*1e-2;
482: #endif
483:   }
484:   if (ctx->nev==1) ctx->nev = PetscMin(20,pep->n);  /* nev not set, use default value */
486:   pep->ops->backtransform = PEPBackTransform_Skip;
487:   if (pep->max_it==PETSC_DEFAULT) pep->max_it = 100;

489:   /* create spectrum slicing context and initialize it */
490:   PEPQSliceResetSR(pep);
491:   PetscNew(&sr);
492:   ctx->sr   = sr;
493:   sr->itsKs = 0;
494:   sr->nleap = 0;
495:   sr->sPres = NULL;

497:   if (pep->solvematcoeffs) PetscFree(pep->solvematcoeffs);
498:   PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
499:   if (!pep->st) PEPGetST(pep,&pep->st);
500:   STSetTransform(pep->st,PETSC_FALSE);
501:   STSetUp(pep->st);

503:   ctx->hyperbolic = (pep->problem_type==PEP_HYPERBOLIC)? PETSC_TRUE: PETSC_FALSE;

505:   /* check presence of ends and finding direction */
506:   if (pep->inta > PETSC_MIN_REAL || pep->intb >= PETSC_MAX_REAL) {
507:     sr->int0 = pep->inta;
508:     sr->int1 = pep->intb;
509:     sr->dir = 1;
510:     if (pep->intb >= PETSC_MAX_REAL) { /* Right-open interval */
511:       sr->hasEnd = PETSC_FALSE;
512:     } else sr->hasEnd = PETSC_TRUE;
513:   } else {
514:     sr->int0 = pep->intb;
515:     sr->int1 = pep->inta;
516:     sr->dir = -1;
517:     sr->hasEnd = PetscNot(pep->inta <= PETSC_MIN_REAL);
518:   }

520:   /* compute inertia0 */
521:   PEPQSliceGetInertia(pep,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL,ctx->hyperbolic?0:1);
523:   if (!ctx->hyperbolic && ctx->checket) PEPQSliceCheckEigenvalueType(pep,sr->int0,0.0,PETSC_TRUE);

525:   /* compute inertia1 */
526:   PEPQSliceGetInertia(pep,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL,ctx->hyperbolic?0:1);
528:   if (!ctx->hyperbolic && ctx->checket && sr->hasEnd) {
529:     PEPQSliceCheckEigenvalueType(pep,sr->int1,0.0,PETSC_TRUE);
532:     if (sr->dir*(sr->inertia1-sr->inertia0)<0) {
533:       sr->intcorr = -1;
534:       sr->inertia0 = 2*pep->n-sr->inertia0;
535:       sr->inertia1 = 2*pep->n-sr->inertia1;
536:     } else sr->intcorr = 1;
537:   } else {
538:     if (sr->inertia0<=pep->n && sr->inertia1<=pep->n) sr->intcorr = 1;
539:     else if (sr->inertia0>=pep->n && sr->inertia1>=pep->n) sr->intcorr = -1;
540:   }

542:   if (sr->hasEnd) {
543:     sr->dir = -sr->dir; r = sr->int0; sr->int0 = sr->int1; sr->int1 = r;
544:     i = sr->inertia0; sr->inertia0 = sr->inertia1; sr->inertia1 = i;
545:   }

547:   /* number of eigenvalues in interval */
548:   sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
549:   PetscInfo(pep,"QSlice setup: allocating for %" PetscInt_FMT " eigenvalues in [%g,%g]\n",sr->numEigs,(double)pep->inta,(double)pep->intb);
550:   if (sr->numEigs) {
551:     PEPQSliceAllocateSolution(pep);
552:     PEPSetDimensions_Default(pep,ctx->nev,&ctx->ncv,&ctx->mpd);
553:     pep->nev = ctx->nev; pep->ncv = ctx->ncv; pep->mpd = ctx->mpd;
554:     ld   = ctx->ncv+2;
555:     DSSetType(pep->ds,DSGHIEP);
556:     DSSetCompact(pep->ds,PETSC_TRUE);
557:     DSSetExtraRow(pep->ds,PETSC_TRUE);
558:     DSAllocate(pep->ds,ld);
559:     DSGetSlepcSC(pep->ds,&sc);
560:     sc->rg            = NULL;
561:     sc->comparison    = SlepcCompareLargestMagnitude;
562:     sc->comparisonctx = NULL;
563:     sc->map           = NULL;
564:     sc->mapobj        = NULL;
565:   } else {pep->ncv = 0; pep->nev = 0; pep->mpd = 0;}
566:   return 0;
567: }

569: /*
570:    Fills the fields of a shift structure
571: */
572: static PetscErrorCode PEPCreateShift(PEP pep,PetscReal val,PEP_shift neighb0,PEP_shift neighb1)
573: {
574:   PEP_shift      s,*pending2;
575:   PetscInt       i;
576:   PEP_SR         sr;
577:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;

579:   sr = ctx->sr;
580:   PetscNew(&s);
581:   s->value = val;
582:   s->neighb[0] = neighb0;
583:   if (neighb0) neighb0->neighb[1] = s;
584:   s->neighb[1] = neighb1;
585:   if (neighb1) neighb1->neighb[0] = s;
586:   s->comp[0] = PETSC_FALSE;
587:   s->comp[1] = PETSC_FALSE;
588:   s->index = -1;
589:   s->neigs = 0;
590:   s->nconv[0] = s->nconv[1] = 0;
591:   s->nsch[0] = s->nsch[1]=0;
592:   /* Inserts in the stack of pending shifts */
593:   /* If needed, the array is resized */
594:   if (sr->nPend >= sr->maxPend) {
595:     sr->maxPend *= 2;
596:     PetscMalloc1(sr->maxPend,&pending2);
597:     for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
598:     PetscFree(sr->pending);
599:     sr->pending = pending2;
600:   }
601:   sr->pending[sr->nPend++]=s;
602:   return 0;
603: }

605: /* Provides next shift to be computed */
606: static PetscErrorCode PEPExtractShift(PEP pep)
607: {
608:   PetscInt       iner,zeros=0;
609:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
610:   PEP_SR         sr;
611:   PetscReal      newShift,aux;
612:   PEP_shift      sPres;

614:   sr = ctx->sr;
615:   if (sr->nPend > 0) {
616:     if (sr->dirch) {
617:       aux = sr->int1; sr->int1 = sr->int0; sr->int0 = aux;
618:       iner = sr->inertia1; sr->inertia1 = sr->inertia0; sr->inertia0 = iner;
619:       sr->dir *= -1;
620:       PetscFree(sr->s0->neighb[1]);
621:       PetscFree(sr->s0);
622:       sr->nPend--;
623:       PEPCreateShift(pep,sr->int0,NULL,NULL);
624:       sr->sPrev = NULL;
625:       sr->sPres = sr->pending[--sr->nPend];
626:       pep->target = sr->sPres->value;
627:       sr->s0 = sr->sPres;
628:       pep->reason = PEP_CONVERGED_ITERATING;
629:     } else {
630:       sr->sPrev = sr->sPres;
631:       sr->sPres = sr->pending[--sr->nPend];
632:     }
633:     sPres = sr->sPres;
634:     PEPQSliceGetInertia(pep,sPres->value,&iner,ctx->detect?&zeros:NULL,sr->intcorr);
635:     if (zeros) {
636:       newShift = sPres->value*(1.0+SLICE_PTOL);
637:       if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
638:       else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
639:       PEPQSliceGetInertia(pep,newShift,&iner,&zeros,sr->intcorr);
641:       sPres->value = newShift;
642:     }
643:     sr->sPres->inertia = iner;
644:     pep->target = sr->sPres->value;
645:     pep->reason = PEP_CONVERGED_ITERATING;
646:     pep->its = 0;
647:   } else sr->sPres = NULL;
648:   return 0;
649: }

651: /*
652:   Obtains value of subsequent shift
653: */
654: static PetscErrorCode PEPGetNewShiftValue(PEP pep,PetscInt side,PetscReal *newS)
655: {
656:   PetscReal lambda,d_prev;
657:   PetscInt  i,idxP;
658:   PEP_SR    sr;
659:   PEP_shift sPres,s;
660:   PEP_STOAR *ctx=(PEP_STOAR*)pep->data;

662:   sr = ctx->sr;
663:   sPres = sr->sPres;
664:   if (sPres->neighb[side]) {
665:   /* Completing a previous interval */
666:     if (!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0) { /* One of the ends might be too far from eigenvalues */
667:       if (side) *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[sr->indexEig-1]]))/2;
668:       else *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[0]]))/2;
669:     } else *newS=(sPres->value + sPres->neighb[side]->value)/2;
670:   } else { /* (Only for side=1). Creating a new interval. */
671:     if (sPres->neigs==0) {/* No value has been accepted*/
672:       if (sPres->neighb[0]) {
673:         /* Multiplying by 10 the previous distance */
674:         *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
675:         sr->nleap++;
676:         /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
678:       } else { /* First shift */
679:         if (pep->nconv != 0) {
680:           /* Unaccepted values give information for next shift */
681:           idxP=0;/* Number of values left from shift */
682:           for (i=0;i<pep->nconv;i++) {
683:             lambda = PetscRealPart(pep->eigr[i]);
684:             if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
685:             else break;
686:           }
687:           /* Avoiding subtraction of eigenvalues (might be the same).*/
688:           if (idxP>0) {
689:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(pep->eigr[0]))/(idxP+0.3);
690:           } else {
691:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(pep->eigr[pep->nconv-1]))/(pep->nconv+0.3);
692:           }
693:           *newS = sPres->value + ((sr->dir)*d_prev*pep->nev)/2;
694:           sr->dirch = PETSC_FALSE;
695:         } else { /* No values found, no information for next shift */
697:           sr->dirch = PETSC_TRUE;
698:           *newS = sr->int1;
699:         }
700:       }
701:     } else { /* Accepted values found */
702:       sr->dirch = PETSC_FALSE;
703:       sr->nleap = 0;
704:       /* Average distance of values in previous subinterval */
705:       s = sPres->neighb[0];
706:       while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
707:         s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
708:       }
709:       if (s) {
710:         d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
711:       } else { /* First shift. Average distance obtained with values in this shift */
712:         /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
713:         if ((sr->dir)*(PetscRealPart(sr->eigr[0])-sPres->value)>0 && PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0]))/PetscRealPart(sr->eigr[0])) > PetscSqrtReal(pep->tol)) {
714:           d_prev =  PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
715:         } else {
716:           d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
717:         }
718:       }
719:       /* Average distance is used for next shift by adding it to value on the right or to shift */
720:       if ((sr->dir)*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
721:         *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(pep->nev))/2;
722:       } else { /* Last accepted value is on the left of shift. Adding to shift */
723:         *newS = sPres->value + ((sr->dir)*d_prev*(pep->nev))/2;
724:       }
725:     }
726:     /* End of interval can not be surpassed */
727:     if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
728:   }/* of neighb[side]==null */
729:   return 0;
730: }

732: /*
733:   Function for sorting an array of real values
734: */
735: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
736: {
737:   PetscReal re;
738:   PetscInt  i,j,tmp;

740:   if (!prev) for (i=0;i<nr;i++) perm[i] = i;
741:   /* Insertion sort */
742:   for (i=1;i<nr;i++) {
743:     re = PetscRealPart(r[perm[i]]);
744:     j = i-1;
745:     while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
746:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
747:     }
748:   }
749:   return 0;
750: }

752: /* Stores the pairs obtained since the last shift in the global arrays */
753: static PetscErrorCode PEPStoreEigenpairs(PEP pep)
754: {
755:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
756:   PetscReal      lambda,err,*errest;
757:   PetscInt       i,*aux,count=0,ndef,ld,nconv=pep->nconv,d=pep->nmat-1,idx;
758:   PetscBool      iscayley,divide=PETSC_FALSE;
759:   PEP_SR         sr = ctx->sr;
760:   PEP_shift      sPres;
761:   Vec            w,vomega;
762:   Mat            MS;
763:   BV             tV;
764:   PetscScalar    *S,*eigr,*tS,*omega;

766:   sPres = sr->sPres;
767:   sPres->index = sr->indexEig;

769:   if (nconv>sr->ndef0+sr->ndef1) {
770:     /* Back-transform */
771:     STBackTransform(pep->st,nconv,pep->eigr,pep->eigi);
772:     for (i=0;i<nconv;i++) {
773: #if defined(PETSC_USE_COMPLEX)
774:       if (PetscImaginaryPart(pep->eigr[i])) pep->eigr[i] = sr->int0-sr->dir;
775: #else
776:       if (pep->eigi[i]) pep->eigr[i] = sr->int0-sr->dir;
777: #endif
778:     }
779:     PetscObjectTypeCompare((PetscObject)pep->st,STCAYLEY,&iscayley);
780:     /* Sort eigenvalues */
781:     sortRealEigenvalues(pep->eigr,pep->perm,nconv,PETSC_FALSE,sr->dir);
782:     VecCreateSeq(PETSC_COMM_SELF,nconv,&vomega);
783:     BVGetSignature(ctx->V,vomega);
784:     VecGetArray(vomega,&omega);
785:     BVGetSizes(pep->V,NULL,NULL,&ld);
786:     BVTensorGetFactors(ctx->V,NULL,&MS);
787:     MatDenseGetArray(MS,&S);
788:     /* Values stored in global array */
789:     PetscCalloc4(nconv,&eigr,nconv,&errest,nconv*nconv*d,&tS,nconv,&aux);
790:     ndef = sr->ndef0+sr->ndef1;
791:     for (i=0;i<nconv;i++) {
792:       lambda = PetscRealPart(pep->eigr[pep->perm[i]]);
793:       err = pep->errest[pep->perm[i]];
794:       if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
796:         PEPQSliceCheckEigenvalueType(pep,lambda,PetscRealPart(omega[pep->perm[i]]),PETSC_FALSE);
797:         eigr[count] = lambda;
798:         errest[count] = err;
799:         if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) sPres->nconv[0]++;
800:         if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) sPres->nconv[1]++;
801:         PetscArraycpy(tS+count*(d*nconv),S+pep->perm[i]*(d*ld),nconv);
802:         PetscArraycpy(tS+count*(d*nconv)+nconv,S+pep->perm[i]*(d*ld)+ld,nconv);
803:         count++;
804:       }
805:     }
806:     VecRestoreArray(vomega,&omega);
807:     VecDestroy(&vomega);
808:     for (i=0;i<count;i++) {
809:       PetscArraycpy(S+i*(d*ld),tS+i*nconv*d,nconv);
810:       PetscArraycpy(S+i*(d*ld)+ld,tS+i*nconv*d+nconv,nconv);
811:     }
812:     MatDenseRestoreArray(MS,&S);
813:     BVTensorRestoreFactors(ctx->V,NULL,&MS);
814:     BVSetActiveColumns(ctx->V,0,count);
815:     BVTensorCompress(ctx->V,count);
816:     if (sr->sPres->nconv[0] && sr->sPres->nconv[1]) {
817:       divide = PETSC_TRUE;
818:       BVTensorGetFactors(ctx->V,NULL,&MS);
819:       MatDenseGetArray(MS,&S);
820:       PetscArrayzero(tS,nconv*nconv*d);
821:       for (i=0;i<count;i++) {
822:         PetscArraycpy(tS+i*nconv*d,S+i*(d*ld),count);
823:         PetscArraycpy(tS+i*nconv*d+nconv,S+i*(d*ld)+ld,count);
824:       }
825:       MatDenseRestoreArray(MS,&S);
826:       BVTensorRestoreFactors(ctx->V,NULL,&MS);
827:       BVSetActiveColumns(pep->V,0,count);
828:       BVDuplicateResize(pep->V,count,&tV);
829:       BVCopy(pep->V,tV);
830:     }
831:     if (sr->sPres->nconv[0]) {
832:       if (divide) {
833:         BVSetActiveColumns(ctx->V,0,sr->sPres->nconv[0]);
834:         BVTensorCompress(ctx->V,sr->sPres->nconv[0]);
835:       }
836:       for (i=0;i<sr->ndef0;i++) aux[i] = sr->idxDef0[i];
837:       for (i=sr->ndef0;i<sr->sPres->nconv[0];i++) aux[i] = sr->indexEig+i-sr->ndef0;
838:       BVTensorGetFactors(ctx->V,NULL,&MS);
839:       MatDenseGetArray(MS,&S);
840:       for (i=0;i<sr->sPres->nconv[0];i++) {
841:         sr->eigr[aux[i]] = eigr[i];
842:         sr->errest[aux[i]] = errest[i];
843:         BVGetColumn(pep->V,i,&w);
844:         BVInsertVec(sr->V,aux[i],w);
845:         BVRestoreColumn(pep->V,i,&w);
846:         idx = sr->ld*d*aux[i];
847:         PetscArrayzero(sr->S+idx,sr->ld*d);
848:         PetscArraycpy(sr->S+idx,S+i*(ld*d),sr->sPres->nconv[0]);
849:         PetscArraycpy(sr->S+idx+sr->ld,S+i*(ld*d)+ld,sr->sPres->nconv[0]);
850:         PetscFree(sr->qinfo[aux[i]].q);
851:         PetscMalloc1(sr->sPres->nconv[0],&sr->qinfo[aux[i]].q);
852:         PetscArraycpy(sr->qinfo[aux[i]].q,aux,sr->sPres->nconv[0]);
853:         sr->qinfo[aux[i]].nq = sr->sPres->nconv[0];
854:       }
855:       MatDenseRestoreArray(MS,&S);
856:       BVTensorRestoreFactors(ctx->V,NULL,&MS);
857:     }

859:     if (sr->sPres->nconv[1]) {
860:       if (divide) {
861:         BVTensorGetFactors(ctx->V,NULL,&MS);
862:         MatDenseGetArray(MS,&S);
863:         for (i=0;i<sr->sPres->nconv[1];i++) {
864:           PetscArraycpy(S+i*(d*ld),tS+(sr->sPres->nconv[0]+i)*nconv*d,count);
865:           PetscArraycpy(S+i*(d*ld)+ld,tS+(sr->sPres->nconv[0]+i)*nconv*d+nconv,count);
866:         }
867:         MatDenseRestoreArray(MS,&S);
868:         BVTensorRestoreFactors(ctx->V,NULL,&MS);
869:         BVSetActiveColumns(pep->V,0,count);
870:         BVCopy(tV,pep->V);
871:         BVSetActiveColumns(ctx->V,0,sr->sPres->nconv[1]);
872:         BVTensorCompress(ctx->V,sr->sPres->nconv[1]);
873:       }
874:       for (i=0;i<sr->ndef1;i++) aux[i] = sr->idxDef1[i];
875:       for (i=sr->ndef1;i<sr->sPres->nconv[1];i++) aux[i] = sr->indexEig+sr->sPres->nconv[0]-sr->ndef0+i-sr->ndef1;
876:       BVTensorGetFactors(ctx->V,NULL,&MS);
877:       MatDenseGetArray(MS,&S);
878:       for (i=0;i<sr->sPres->nconv[1];i++) {
879:         sr->eigr[aux[i]] = eigr[sr->sPres->nconv[0]+i];
880:         sr->errest[aux[i]] = errest[sr->sPres->nconv[0]+i];
881:         BVGetColumn(pep->V,i,&w);
882:         BVInsertVec(sr->V,aux[i],w);
883:         BVRestoreColumn(pep->V,i,&w);
884:         idx = sr->ld*d*aux[i];
885:         PetscArrayzero(sr->S+idx,sr->ld*d);
886:         PetscArraycpy(sr->S+idx,S+i*(ld*d),sr->sPres->nconv[1]);
887:         PetscArraycpy(sr->S+idx+sr->ld,S+i*(ld*d)+ld,sr->sPres->nconv[1]);
888:         PetscFree(sr->qinfo[aux[i]].q);
889:         PetscMalloc1(sr->sPres->nconv[1],&sr->qinfo[aux[i]].q);
890:         PetscArraycpy(sr->qinfo[aux[i]].q,aux,sr->sPres->nconv[1]);
891:         sr->qinfo[aux[i]].nq = sr->sPres->nconv[1];
892:       }
893:       MatDenseRestoreArray(MS,&S);
894:       BVTensorRestoreFactors(ctx->V,NULL,&MS);
895:     }
896:     sPres->neigs = count-sr->ndef0-sr->ndef1;
897:     sr->indexEig += sPres->neigs;
898:     sPres->nconv[0]-= sr->ndef0;
899:     sPres->nconv[1]-= sr->ndef1;
900:     PetscFree4(eigr,errest,tS,aux);
901:   } else {
902:     sPres->neigs = 0;
903:     sPres->nconv[0]= 0;
904:     sPres->nconv[1]= 0;
905:   }
906:   /* Global ordering array updating */
907:   sortRealEigenvalues(sr->eigr,sr->perm,sr->indexEig,PETSC_FALSE,sr->dir);
908:   /* Check for completion */
909:   sPres->comp[0] = PetscNot(sPres->nconv[0] < sPres->nsch[0]);
910:   sPres->comp[1] = PetscNot(sPres->nconv[1] < sPres->nsch[1]);
912:   if (divide) BVDestroy(&tV);
913:   return 0;
914: }

916: static PetscErrorCode PEPLookForDeflation(PEP pep)
917: {
918:   PetscReal val;
919:   PetscInt  i,count0=0,count1=0;
920:   PEP_shift sPres;
921:   PetscInt  ini,fin;
922:   PEP_SR    sr;
923:   PEP_STOAR *ctx=(PEP_STOAR*)pep->data;

925:   sr = ctx->sr;
926:   sPres = sr->sPres;

928:   if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
929:   else ini = 0;
930:   fin = sr->indexEig;
931:   /* Selection of ends for searching new values */
932:   if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
933:   else sPres->ext[0] = sPres->neighb[0]->value;
934:   if (!sPres->neighb[1]) {
935:     if (sr->hasEnd) sPres->ext[1] = sr->int1;
936:     else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
937:   } else sPres->ext[1] = sPres->neighb[1]->value;
938:   /* Selection of values between right and left ends */
939:   for (i=ini;i<fin;i++) {
940:     val=PetscRealPart(sr->eigr[sr->perm[i]]);
941:     /* Values to the right of left shift */
942:     if ((sr->dir)*(val - sPres->ext[1]) < 0) {
943:       if ((sr->dir)*(val - sPres->value) < 0) count0++;
944:       else count1++;
945:     } else break;
946:   }
947:   /* The number of values on each side are found */
948:   if (sPres->neighb[0]) {
949:     sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
951:   } else sPres->nsch[0] = 0;

953:   if (sPres->neighb[1]) {
954:     sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
956:   } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);

958:   /* Completing vector of indexes for deflation */
959:   for (i=0;i<count0;i++) sr->idxDef0[i] = sr->perm[ini+i];
960:   sr->ndef0 = count0;
961:   for (i=0;i<count1;i++) sr->idxDef1[i] = sr->perm[ini+count0+i];
962:   sr->ndef1 = count1;
963:   return 0;
964: }

966: /*
967:   Compute a run of Lanczos iterations
968: */
969: static PetscErrorCode PEPSTOARrun_QSlice(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,Vec *t_)
970: {
971:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;
972:   PetscInt       i,j,m=*M,l,lock;
973:   PetscInt       lds,d,ld,offq,nqt,ldds;
974:   Vec            v=t_[0],t=t_[1],q=t_[2];
975:   PetscReal      norm,sym=0.0,fro=0.0,*f;
976:   PetscScalar    *y,*S,sigma;
977:   PetscBLASInt   j_,one=1;
978:   PetscBool      lindep;
979:   Mat            MS;

981:   PetscMalloc1(*M,&y);
982:   BVGetSizes(pep->V,NULL,NULL,&ld);
983:   BVTensorGetDegree(ctx->V,&d);
984:   BVGetActiveColumns(pep->V,&lock,&nqt);
985:   lds = d*ld;
986:   offq = ld;
987:   DSGetLeadingDimension(pep->ds,&ldds);

989:   *breakdown = PETSC_FALSE; /* ----- */
990:   STGetShift(pep->st,&sigma);
991:   DSGetDimensions(pep->ds,NULL,&l,NULL,NULL);
992:   BVSetActiveColumns(ctx->V,0,m);
993:   BVSetActiveColumns(pep->V,0,nqt);
994:   for (j=k;j<m;j++) {
995:     /* apply operator */
996:     BVTensorGetFactors(ctx->V,NULL,&MS);
997:     MatDenseGetArray(MS,&S);
998:     BVGetColumn(pep->V,nqt,&t);
999:     BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
1000:     MatMult(pep->A[1],v,q);
1001:     MatMult(pep->A[2],v,t);
1002:     VecAXPY(q,sigma*pep->sfactor,t);
1003:     VecScale(q,pep->sfactor);
1004:     BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
1005:     MatMult(pep->A[2],v,t);
1006:     VecAXPY(q,pep->sfactor*pep->sfactor,t);
1007:     STMatSolve(pep->st,q,t);
1008:     VecScale(t,-1.0);
1009:     BVRestoreColumn(pep->V,nqt,&t);

1011:     /* orthogonalize */
1012:     BVOrthogonalizeColumn(pep->V,nqt,S+(j+1)*lds,&norm,&lindep);
1013:     if (!lindep) {
1014:       *(S+(j+1)*lds+nqt) = norm;
1015:       BVScaleColumn(pep->V,nqt,1.0/norm);
1016:       nqt++;
1017:     }
1018:     for (i=0;i<nqt;i++) *(S+(j+1)*lds+offq+i) = *(S+j*lds+i)+sigma*(*(S+(j+1)*lds+i));
1019:     BVSetActiveColumns(pep->V,0,nqt);
1020:     MatDenseRestoreArray(MS,&S);
1021:     BVTensorRestoreFactors(ctx->V,NULL,&MS);

1023:     /* level-2 orthogonalization */
1024:     BVOrthogonalizeColumn(ctx->V,j+1,y,&norm,&lindep);
1025:     a[j] = PetscRealPart(y[j]);
1026:     omega[j+1] = (norm > 0)?1.0:-1.0;
1027:     BVScaleColumn(ctx->V,j+1,1.0/norm);
1028:     b[j] = PetscAbsReal(norm);

1030:     /* check symmetry */
1031:     DSGetArrayReal(pep->ds,DS_MAT_T,&f);
1032:     if (j==k) {
1033:       for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ldds+i]);
1034:       for (i=0;i<l;i++) y[i] = 0.0;
1035:     }
1036:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
1037:     if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
1038:     PetscBLASIntCast(j,&j_);
1039:     sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
1040:     fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
1041:     if (j>0) fro = SlepcAbs(fro,b[j-1]);
1042:     if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
1043:       *symmlost = PETSC_TRUE;
1044:       *M=j;
1045:       break;
1046:     }
1047:   }
1048:   BVSetActiveColumns(pep->V,lock,nqt);
1049:   BVSetActiveColumns(ctx->V,0,*M);
1050:   PetscFree(y);
1051:   return 0;
1052: }

1054: static PetscErrorCode PEPSTOAR_QSlice(PEP pep,Mat B)
1055: {
1056:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;
1057:   PetscInt       j,k,l,nv=0,ld,ldds,t,nq=0,idx;
1058:   PetscInt       nconv=0,deg=pep->nmat-1,count0=0,count1=0;
1059:   PetscScalar    *om,sigma,*back,*S,*pQ;
1060:   PetscReal      beta,norm=1.0,*omega,*a,*b,eta,lambda;
1061:   PetscBool      breakdown,symmlost=PETSC_FALSE,sinv,falselock=PETSC_TRUE;
1062:   Mat            MS,MQ,D;
1063:   Vec            v,vomega;
1064:   PEP_SR         sr;
1065:   BVOrthogType   otype;
1066:   BVOrthogBlockType obtype;

1068:   /* Resize if needed for deflating vectors  */
1069:   sr = ctx->sr;
1070:   sigma = sr->sPres->value;
1071:   k = sr->ndef0+sr->ndef1;
1072:   pep->ncv = ctx->ncv+k;
1073:   pep->nev = ctx->nev+k;
1074:   PEPAllocateSolution(pep,3);
1075:   BVDestroy(&ctx->V);
1076:   BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
1077:   BVGetOrthogonalization(pep->V,&otype,NULL,&eta,&obtype);
1078:   BVSetOrthogonalization(ctx->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
1079:   DSAllocate(pep->ds,pep->ncv+2);
1080:   PetscMalloc1(pep->ncv,&back);
1081:   DSGetLeadingDimension(pep->ds,&ldds);
1082:   BVSetMatrix(ctx->V,B,PETSC_TRUE);
1084:   /* undocumented option to use a cheaper locking instead of the true locking */
1085:   PetscOptionsGetBool(NULL,NULL,"-pep_stoar_falselocking",&falselock,NULL);
1086:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
1087:   RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
1088:   STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);

1090:   /* Get the starting Arnoldi vector */
1091:   BVSetActiveColumns(pep->V,0,1);
1092:   BVTensorBuildFirstColumn(ctx->V,pep->nini);
1093:   BVSetActiveColumns(ctx->V,0,1);
1094:   if (k) {
1095:     /* Insert deflated vectors */
1096:     BVSetActiveColumns(pep->V,0,0);
1097:     idx = sr->ndef0?sr->idxDef0[0]:sr->idxDef1[0];
1098:     for (j=0;j<k;j++) {
1099:       BVGetColumn(pep->V,j,&v);
1100:       BVCopyVec(sr->V,sr->qinfo[idx].q[j],v);
1101:       BVRestoreColumn(pep->V,j,&v);
1102:     }
1103:     /* Update innerproduct matrix */
1104:     BVSetActiveColumns(ctx->V,0,0);
1105:     BVTensorGetFactors(ctx->V,NULL,&MS);
1106:     BVSetActiveColumns(pep->V,0,k);
1107:     BVTensorRestoreFactors(ctx->V,NULL,&MS);

1109:     BVGetSizes(pep->V,NULL,NULL,&ld);
1110:     BVTensorGetFactors(ctx->V,NULL,&MS);
1111:     MatDenseGetArray(MS,&S);
1112:     for (j=0;j<sr->ndef0;j++) {
1113:       PetscArrayzero(S+j*ld*deg,ld*deg);
1114:       PetscArraycpy(S+j*ld*deg,sr->S+sr->idxDef0[j]*sr->ld*deg,k);
1115:       PetscArraycpy(S+j*ld*deg+ld,sr->S+sr->idxDef0[j]*sr->ld*deg+sr->ld,k);
1116:       pep->eigr[j] = sr->eigr[sr->idxDef0[j]];
1117:       pep->errest[j] = sr->errest[sr->idxDef0[j]];
1118:     }
1119:     for (j=0;j<sr->ndef1;j++) {
1120:       PetscArrayzero(S+(j+sr->ndef0)*ld*deg,ld*deg);
1121:       PetscArraycpy(S+(j+sr->ndef0)*ld*deg,sr->S+sr->idxDef1[j]*sr->ld*deg,k);
1122:       PetscArraycpy(S+(j+sr->ndef0)*ld*deg+ld,sr->S+sr->idxDef1[j]*sr->ld*deg+sr->ld,k);
1123:       pep->eigr[j+sr->ndef0] = sr->eigr[sr->idxDef1[j]];
1124:       pep->errest[j+sr->ndef0] = sr->errest[sr->idxDef1[j]];
1125:     }
1126:     MatDenseRestoreArray(MS,&S);
1127:     BVTensorRestoreFactors(ctx->V,NULL,&MS);
1128:     BVSetActiveColumns(ctx->V,0,k+1);
1129:     VecCreateSeq(PETSC_COMM_SELF,k+1,&vomega);
1130:     VecGetArray(vomega,&om);
1131:     for (j=0;j<k;j++) {
1132:       BVOrthogonalizeColumn(ctx->V,j,NULL,&norm,NULL);
1133:       BVScaleColumn(ctx->V,j,1/norm);
1134:       om[j] = (norm>=0.0)?1.0:-1.0;
1135:     }
1136:     BVTensorGetFactors(ctx->V,NULL,&MS);
1137:     MatDenseGetArray(MS,&S);
1138:     for (j=0;j<deg;j++) {
1139:       BVSetRandomColumn(pep->V,k+j);
1140:       BVOrthogonalizeColumn(pep->V,k+j,S+k*ld*deg+j*ld,&norm,NULL);
1141:       BVScaleColumn(pep->V,k+j,1.0/norm);
1142:       S[k*ld*deg+j*ld+k+j] = norm;
1143:     }
1144:     MatDenseRestoreArray(MS,&S);
1145:     BVSetActiveColumns(pep->V,0,k+deg);
1146:     BVTensorRestoreFactors(ctx->V,NULL,&MS);
1147:     BVOrthogonalizeColumn(ctx->V,k,NULL,&norm,NULL);
1148:     BVScaleColumn(ctx->V,k,1.0/norm);
1149:     om[k] = (norm>=0.0)?1.0:-1.0;
1150:     VecRestoreArray(vomega,&om);
1151:     BVSetSignature(ctx->V,vomega);
1152:     DSGetArrayReal(pep->ds,DS_MAT_T,&a);
1153:     VecGetArray(vomega,&om);
1154:     for (j=0;j<k;j++) a[j] = PetscRealPart(om[j]/(pep->eigr[j]-sigma));
1155:     VecRestoreArray(vomega,&om);
1156:     VecDestroy(&vomega);
1157:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
1158:     DSGetArray(pep->ds,DS_MAT_Q,&pQ);
1159:     PetscArrayzero(pQ,ldds*k);
1160:     for (j=0;j<k;j++) pQ[j+j*ldds] = 1.0;
1161:     DSRestoreArray(pep->ds,DS_MAT_Q,&pQ);
1162:   }
1163:   BVSetActiveColumns(ctx->V,0,k+1);
1164:   DSSetDimensions(pep->ds,k+1,PETSC_DEFAULT,PETSC_DEFAULT);
1165:   DSGetMatAndColumn(pep->ds,DS_MAT_D,0,&D,&vomega);
1166:   BVGetSignature(ctx->V,vomega);
1167:   DSRestoreMatAndColumn(pep->ds,DS_MAT_D,0,&D,&vomega);

1169:   PetscInfo(pep,"Start STOAR: sigma=%g in [%g,%g], for deflation: left=%" PetscInt_FMT " right=%" PetscInt_FMT ", searching: left=%" PetscInt_FMT " right=%" PetscInt_FMT "\n",(double)sr->sPres->value,(double)(sr->sPres->neighb[0]?sr->sPres->neighb[0]->value:sr->int0),(double)(sr->sPres->neighb[1]?sr->sPres->neighb[1]->value:sr->int1),sr->ndef0,sr->ndef1,sr->sPres->nsch[0],sr->sPres->nsch[1]);

1171:   /* Restart loop */
1172:   l = 0;
1173:   pep->nconv = k;
1174:   while (pep->reason == PEP_CONVERGED_ITERATING) {
1175:     pep->its++;
1176:     DSGetArrayReal(pep->ds,DS_MAT_T,&a);
1177:     b = a+ldds;
1178:     DSGetArrayReal(pep->ds,DS_MAT_D,&omega);

1180:     /* Compute an nv-step Lanczos factorization */
1181:     nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
1182:     PEPSTOARrun_QSlice(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,pep->work);
1183:     beta = b[nv-1];
1184:     if (symmlost && nv==pep->nconv+l) {
1185:       pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
1186:       pep->nconv = nconv;
1187:       PetscInfo(pep,"Symmetry lost in STOAR sigma=%g nconv=%" PetscInt_FMT "\n",(double)sr->sPres->value,nconv);
1188:       if (falselock || !ctx->lock) {
1189:         BVSetActiveColumns(ctx->V,0,pep->nconv);
1190:         BVTensorCompress(ctx->V,0);
1191:       }
1192:       break;
1193:     }
1194:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
1195:     DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
1196:     DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l);
1197:     if (l==0) DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
1198:     else DSSetState(pep->ds,DS_STATE_RAW);

1200:     /* Solve projected problem */
1201:     DSSolve(pep->ds,pep->eigr,pep->eigi);
1202:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
1203:     DSUpdateExtraRow(pep->ds);
1204:     DSSynchronize(pep->ds,pep->eigr,pep->eigi);

1206:     /* Check convergence */
1207:     /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
1208:     norm = 1.0;
1209:     DSGetDimensions(pep->ds,NULL,NULL,NULL,&t);
1210:     PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
1211:     (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
1212:     for (j=0;j<k;j++) back[j] = pep->eigr[j];
1213:     STBackTransform(pep->st,k,back,pep->eigi);
1214:     count0=count1=0;
1215:     for (j=0;j<k;j++) {
1216:       lambda = PetscRealPart(back[j]);
1217:       if (((sr->dir)*(sr->sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sr->sPres->ext[0]) > 0)) count0++;
1218:       if (((sr->dir)*(lambda - sr->sPres->value) > 0) && ((sr->dir)*(sr->sPres->ext[1] - lambda) > 0)) count1++;
1219:     }
1220:     if ((count0-sr->ndef0 >= sr->sPres->nsch[0]) && (count1-sr->ndef1 >= sr->sPres->nsch[1])) pep->reason = PEP_CONVERGED_TOL;
1221:     /* Update l */
1222:     if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
1223:     else {
1224:       l = PetscMax(1,(PetscInt)((nv-k)/2));
1225:       l = PetscMin(l,t);
1226:       DSGetTruncateSize(pep->ds,k,t,&l);
1227:       if (!breakdown) {
1228:         /* Prepare the Rayleigh quotient for restart */
1229:         DSTruncate(pep->ds,k+l,PETSC_FALSE);
1230:       }
1231:     }
1232:     nconv = k;
1233:     if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
1234:     if (l) PetscInfo(pep,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l);

1236:     /* Update S */
1237:     DSGetMat(pep->ds,DS_MAT_Q,&MQ);
1238:     BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
1239:     DSRestoreMat(pep->ds,DS_MAT_Q,&MQ);

1241:     /* Copy last column of S */
1242:     BVCopyColumn(ctx->V,nv,k+l);
1243:     BVSetActiveColumns(ctx->V,0,k+l);
1244:     if (k+l) {
1245:       DSSetDimensions(pep->ds,k+l,PETSC_DEFAULT,PETSC_DEFAULT);
1246:       DSGetMatAndColumn(pep->ds,DS_MAT_D,0,&D,&vomega);
1247:       BVSetSignature(ctx->V,vomega);
1248:       DSRestoreMatAndColumn(pep->ds,DS_MAT_D,0,&D,&vomega);
1249:     }

1251:     if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
1252:       /* stop if breakdown */
1253:       PetscInfo(pep,"Breakdown TOAR method (it=%" PetscInt_FMT " norm=%g)\n",pep->its,(double)beta);
1254:       pep->reason = PEP_DIVERGED_BREAKDOWN;
1255:     }
1256:     if (pep->reason != PEP_CONVERGED_ITERATING) l--;
1257:     BVGetActiveColumns(pep->V,NULL,&nq);
1258:     if (k+l+deg<=nq) {
1259:       BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
1260:       if (!falselock && ctx->lock) BVTensorCompress(ctx->V,k-pep->nconv);
1261:       else BVTensorCompress(ctx->V,0);
1262:     }
1263:     pep->nconv = k;
1264:     PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
1265:   }
1266:   sr->itsKs += pep->its;
1267:   if (pep->nconv>0) {
1268:     BVSetActiveColumns(ctx->V,0,pep->nconv);
1269:     BVGetActiveColumns(pep->V,NULL,&nq);
1270:     BVSetActiveColumns(pep->V,0,nq);
1271:     if (nq>pep->nconv) {
1272:       BVTensorCompress(ctx->V,pep->nconv);
1273:       BVSetActiveColumns(pep->V,0,pep->nconv);
1274:     }
1275:     for (j=0;j<pep->nconv;j++) {
1276:       pep->eigr[j] *= pep->sfactor;
1277:       pep->eigi[j] *= pep->sfactor;
1278:     }
1279:   }
1280:   PetscInfo(pep,"Finished STOAR: nconv=%" PetscInt_FMT " (deflated=%" PetscInt_FMT ", left=%" PetscInt_FMT ", right=%" PetscInt_FMT ")\n",pep->nconv,sr->ndef0+sr->ndef1,count0-sr->ndef0,count1-sr->ndef1);
1281:   STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
1282:   RGPopScale(pep->rg);

1285:   if (pep->reason == PEP_DIVERGED_SYMMETRY_LOST && nconv==sr->ndef0+sr->ndef1) {
1287:   } else sr->symmlost = 0;

1289:   DSTruncate(pep->ds,pep->nconv,PETSC_TRUE);
1290:   PetscFree(back);
1291:   return 0;
1292: }

1294: #define SWAP(a,b,t) {t=a;a=b;b=t;}

1296: static PetscErrorCode PEPQSliceGetInertias(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
1297: {
1298:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
1299:   PEP_SR          sr=ctx->sr;
1300:   PetscInt        i=0,j,tmpi;
1301:   PetscReal       v,tmpr;
1302:   PEP_shift       s;

1306:   if (!sr->s0) {  /* PEPSolve not called yet */
1307:     *n = 2;
1308:   } else {
1309:     *n = 1;
1310:     s = sr->s0;
1311:     while (s) {
1312:       (*n)++;
1313:       s = s->neighb[1];
1314:     }
1315:   }
1316:   PetscMalloc1(*n,shifts);
1317:   PetscMalloc1(*n,inertias);
1318:   if (!sr->s0) {  /* PEPSolve not called yet */
1319:     (*shifts)[0]   = sr->int0;
1320:     (*shifts)[1]   = sr->int1;
1321:     (*inertias)[0] = sr->inertia0;
1322:     (*inertias)[1] = sr->inertia1;
1323:   } else {
1324:     s = sr->s0;
1325:     while (s) {
1326:       (*shifts)[i]     = s->value;
1327:       (*inertias)[i++] = s->inertia;
1328:       s = s->neighb[1];
1329:     }
1330:     (*shifts)[i]   = sr->int1;
1331:     (*inertias)[i] = sr->inertia1;
1332:   }
1333:   /* remove possible duplicate in last position */
1334:   if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
1335:   /* sort result */
1336:   for (i=0;i<*n;i++) {
1337:     v = (*shifts)[i];
1338:     for (j=i+1;j<*n;j++) {
1339:       if (v > (*shifts)[j]) {
1340:         SWAP((*shifts)[i],(*shifts)[j],tmpr);
1341:         SWAP((*inertias)[i],(*inertias)[j],tmpi);
1342:         v = (*shifts)[i];
1343:       }
1344:     }
1345:   }
1346:   return 0;
1347: }

1349: PetscErrorCode PEPSolve_STOAR_QSlice(PEP pep)
1350: {
1351:   PetscInt       i,j,ti,deg=pep->nmat-1;
1352:   PetscReal      newS;
1353:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
1354:   PEP_SR         sr=ctx->sr;
1355:   Mat            S,B;
1356:   PetscScalar    *pS;

1358:   PetscCitationsRegister(citation,&cited);

1360:   /* Only with eigenvalues present in the interval ...*/
1361:   if (sr->numEigs==0) {
1362:     pep->reason = PEP_CONVERGED_TOL;
1363:     return 0;
1364:   }

1366:   /* Inner product matrix */
1367:   PEPSTOARSetUpInnerMatrix(pep,&B);

1369:   /* Array of pending shifts */
1370:   sr->maxPend = 100; /* Initial size */
1371:   sr->nPend = 0;
1372:   PetscMalloc1(sr->maxPend,&sr->pending);
1373:   PEPCreateShift(pep,sr->int0,NULL,NULL);
1374:   /* extract first shift */
1375:   sr->sPrev = NULL;
1376:   sr->sPres = sr->pending[--sr->nPend];
1377:   sr->sPres->inertia = sr->inertia0;
1378:   pep->target = sr->sPres->value;
1379:   sr->s0 = sr->sPres;
1380:   sr->indexEig = 0;

1382:   for (i=0;i<sr->numEigs;i++) {
1383:     sr->eigr[i]   = 0.0;
1384:     sr->eigi[i]   = 0.0;
1385:     sr->errest[i] = 0.0;
1386:     sr->perm[i]   = i;
1387:   }
1388:   /* Vectors for deflation */
1389:   PetscMalloc2(sr->numEigs,&sr->idxDef0,sr->numEigs,&sr->idxDef1);
1390:   sr->indexEig = 0;
1391:   while (sr->sPres) {
1392:     /* Search for deflation */
1393:     PEPLookForDeflation(pep);
1394:     /* KrylovSchur */
1395:     PEPSTOAR_QSlice(pep,B);

1397:     PEPStoreEigenpairs(pep);
1398:     /* Select new shift */
1399:     if (!sr->sPres->comp[1]) {
1400:       PEPGetNewShiftValue(pep,1,&newS);
1401:       PEPCreateShift(pep,newS,sr->sPres,sr->sPres->neighb[1]);
1402:     }
1403:     if (!sr->sPres->comp[0]) {
1404:       /* Completing earlier interval */
1405:       PEPGetNewShiftValue(pep,0,&newS);
1406:       PEPCreateShift(pep,newS,sr->sPres->neighb[0],sr->sPres);
1407:     }
1408:     /* Preparing for a new search of values */
1409:     PEPExtractShift(pep);
1410:   }

1412:   /* Updating pep values prior to exit */
1413:   PetscFree2(sr->idxDef0,sr->idxDef1);
1414:   PetscFree(sr->pending);
1415:   pep->nconv  = sr->indexEig;
1416:   pep->reason = PEP_CONVERGED_TOL;
1417:   pep->its    = sr->itsKs;
1418:   pep->nev    = sr->indexEig;
1419:   MatCreateSeqDense(PETSC_COMM_SELF,pep->nconv,pep->nconv,NULL,&S);
1420:   MatDenseGetArray(S,&pS);
1421:   for (i=0;i<pep->nconv;i++) {
1422:     for (j=0;j<sr->qinfo[i].nq;j++) pS[i*pep->nconv+sr->qinfo[i].q[j]] = *(sr->S+i*sr->ld*deg+j);
1423:   }
1424:   MatDenseRestoreArray(S,&pS);
1425:   BVSetActiveColumns(sr->V,0,pep->nconv);
1426:   BVMultInPlace(sr->V,S,0,pep->nconv);
1427:   MatDestroy(&S);
1428:   BVDestroy(&pep->V);
1429:   pep->V = sr->V;
1430:   PetscFree4(pep->eigr,pep->eigi,pep->errest,pep->perm);
1431:   pep->eigr   = sr->eigr;
1432:   pep->eigi   = sr->eigi;
1433:   pep->perm   = sr->perm;
1434:   pep->errest = sr->errest;
1435:   if (sr->dir<0) {
1436:     for (i=0;i<pep->nconv/2;i++) {
1437:       ti = sr->perm[i]; sr->perm[i] = sr->perm[pep->nconv-1-i]; sr->perm[pep->nconv-1-i] = ti;
1438:     }
1439:   }
1440:   PetscFree(ctx->inertias);
1441:   PetscFree(ctx->shifts);
1442:   MatDestroy(&B);
1443:   PEPQSliceGetInertias(pep,&ctx->nshifts,&ctx->shifts,&ctx->inertias);
1444:   return 0;
1445: }