Actual source code: qslice.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: {
 48:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
 49:   PEP_SR         sr=ctx->sr;
 50:   PEP_shift      s;
 51:   PetscInt       i;

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

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

 81:   PEPQSliceResetSR(pep);
 82:   PetscFree(ctx->inertias);
 83:   PetscFree(ctx->shifts);
 84:   return(0);
 85: }

 87: /*
 88:   PEPQSliceAllocateSolution - Allocate memory storage for common variables such
 89:   as eigenvalues and eigenvectors.
 90: */
 91: static PetscErrorCode PEPQSliceAllocateSolution(PEP pep)
 92: {
 94:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
 95:   PetscInt       k;
 96:   PetscLogDouble cnt;
 97:   BVType         type;
 98:   Vec            t;
 99:   PEP_SR         sr = ctx->sr;

102:   /* allocate space for eigenvalues and friends */
103:   k = PetscMax(1,sr->numEigs);
104:   PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm);
105:   PetscCalloc4(k,&sr->eigr,k,&sr->eigi,k,&sr->errest,k,&sr->perm);
106:   PetscFree(sr->qinfo);
107:   PetscCalloc1(k,&sr->qinfo);
108:   cnt = 2*k*sizeof(PetscScalar) + 2*k*sizeof(PetscReal) + k*sizeof(PetscInt);
109:   PetscLogObjectMemory((PetscObject)pep,cnt);

111:   /* allocate sr->V and transfer options from pep->V */
112:   BVDestroy(&sr->V);
113:   BVCreate(PetscObjectComm((PetscObject)pep),&sr->V);
114:   PetscLogObjectParent((PetscObject)pep,(PetscObject)sr->V);
115:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
116:   if (!((PetscObject)(pep->V))->type_name) {
117:     BVSetType(sr->V,BVSVEC);
118:   } else {
119:     BVGetType(pep->V,&type);
120:     BVSetType(sr->V,type);
121:   }
122:   STMatCreateVecsEmpty(pep->st,&t,NULL);
123:   BVSetSizesFromVec(sr->V,t,k+1);
124:   VecDestroy(&t);
125:   sr->ld = k;
126:   PetscFree(sr->S);
127:   PetscMalloc1((k+1)*sr->ld*(pep->nmat-1),&sr->S);
128:   return(0);
129: }

131: /* Convergence test to compute positive Ritz values */
132: static PetscErrorCode ConvergedPositive(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
133: {
135:   *errest = (PetscRealPart(eigr)>0.0)?0.0:res;
136:   return(0);
137: }

139: static PetscErrorCode PEPQSliceMatGetInertia(PEP pep,PetscReal shift,PetscInt *inertia,PetscInt *zeros)
140: {
141:   KSP            ksp,kspr;
142:   PC             pc;
143:   Mat            F;
144:   PetscBool      flg;

148:   if (!pep->solvematcoeffs) {
149:     PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
150:   }
151:   if (shift==PETSC_MAX_REAL) { /* Inertia of matrix A[2] */
152:     pep->solvematcoeffs[0] = 0.0; pep->solvematcoeffs[1] = 0.0; pep->solvematcoeffs[2] = 1.0;
153:   } else {
154:     PEPEvaluateBasis(pep,shift,0,pep->solvematcoeffs,NULL);
155:   }
156:   STMatSetUp(pep->st,pep->sfactor,pep->solvematcoeffs);
157:   STGetKSP(pep->st,&ksp);
158:   KSPGetPC(ksp,&pc);
159:   PetscObjectTypeCompare((PetscObject)pc,PCREDUNDANT,&flg);
160:   if (flg) {
161:     PCRedundantGetKSP(pc,&kspr);
162:     KSPGetPC(kspr,&pc);
163:   }
164:   PCFactorGetMatrix(pc,&F);
165:   MatGetInertia(F,inertia,zeros,NULL);
166:   return(0);
167: }

169: static PetscErrorCode PEPQSliceGetInertia(PEP pep,PetscReal shift,PetscInt *inertia,PetscInt *zeros,PetscInt correction)
170: {
172:   KSP            ksp;
173:   Mat            P;
174:   PetscReal      nzshift=0.0;
175:   PetscScalar    dot;
176:   PetscRandom    rand;
177:   PetscInt       nconv;
178:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
179:   PEP_SR         sr=ctx->sr;

182:   if (shift >= PETSC_MAX_REAL) { /* Right-open interval */
183:     *inertia = 0;
184:   } else if (shift <= PETSC_MIN_REAL) {
185:     *inertia = 0;
186:     if (zeros) *zeros = 0;
187:   } else {
188:     /* If the shift is zero, perturb it to a very small positive value.
189:        The goal is that the nonzero pattern is the same in all cases and reuse
190:        the symbolic factorizations */
191:     nzshift = (shift==0.0)? 10.0/PETSC_MAX_REAL: shift;
192:     PEPQSliceMatGetInertia(pep,nzshift,inertia,zeros);
193:     STSetShift(pep->st,nzshift);
194:   }
195:   if (!correction) {
196:     if (shift >= PETSC_MAX_REAL) *inertia = 2*pep->n;
197:     else if (shift>PETSC_MIN_REAL) {
198:       STGetKSP(pep->st,&ksp);
199:       KSPGetOperators(ksp,&P,NULL);
200:       if (*inertia!=pep->n && !sr->v[0]) {
201:         MatCreateVecs(P,&sr->v[0],NULL);
202:         VecDuplicate(sr->v[0],&sr->v[1]);
203:         VecDuplicate(sr->v[0],&sr->v[2]);
204:         BVGetRandomContext(pep->V,&rand);
205:         VecSetRandom(sr->v[0],rand);
206:       }
207:       if (*inertia<pep->n && *inertia>0) {
208:         if (!sr->eps) {
209:           EPSCreate(PetscObjectComm((PetscObject)pep),&sr->eps);
210:           EPSSetProblemType(sr->eps,EPS_HEP);
211:           EPSSetWhichEigenpairs(sr->eps,EPS_LARGEST_REAL);
212:         }
213:         EPSSetConvergenceTestFunction(sr->eps,ConvergedPositive,NULL,NULL);
214:         EPSSetOperators(sr->eps,P,NULL);
215:         EPSSolve(sr->eps);
216:         EPSGetConverged(sr->eps,&nconv);
217:         if (!nconv) SETERRQ1(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",nzshift);
218:         EPSGetEigenpair(sr->eps,0,NULL,NULL,sr->v[0],sr->v[1]);
219:       }
220:       if (*inertia!=pep->n) {
221:         MatMult(pep->A[1],sr->v[0],sr->v[1]);
222:         MatMult(pep->A[2],sr->v[0],sr->v[2]);
223:         VecAXPY(sr->v[1],2*nzshift,sr->v[2]);
224:         VecDot(sr->v[1],sr->v[0],&dot);
225:         if (PetscRealPart(dot)>0.0) *inertia = 2*pep->n-*inertia;
226:       }
227:     }
228:   } else if (correction<0) *inertia = 2*pep->n-*inertia;
229:   return(0);
230: }

232: /*
233:    Check eigenvalue type - used only in non-hyperbolic problems.
234:    All computed eigenvalues must have the same definite type (positive or negative).
235:    If ini=TRUE the type is available in omega, otherwise we compute an eigenvalue
236:    closest to shift and determine its type.
237:  */
238: static PetscErrorCode PEPQSliceCheckEigenvalueType(PEP pep,PetscReal shift,PetscReal omega,PetscBool ini)
239: {
241:   PEP            pep2;
242:   ST             st;
243:   PetscInt       nconv;
244:   PetscScalar    lambda,dot;
245:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
246:   PEP_SR         sr=ctx->sr;

249:   if (!ini) {
250:     if (-(omega/(shift*ctx->alpha+ctx->beta))*sr->type<0) SETERRQ1(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"Different positive/negative type detected in eigenvalue %g",(double)shift);
251:   } else {
252:     PEPCreate(PetscObjectComm((PetscObject)pep),&pep2);
253:     PEPSetOptionsPrefix(pep2,((PetscObject)pep)->prefix);
254:     PEPAppendOptionsPrefix(pep2,"pep_eigenvalue_type_");
255:     PEPSetTolerances(pep2,PETSC_DEFAULT,pep->max_it/4);
256:     PEPSetType(pep2,PEPTOAR);
257:     PEPSetOperators(pep2,pep->nmat,pep->A);
258:     PEPSetWhichEigenpairs(pep2,PEP_TARGET_MAGNITUDE);
259:     PEPGetRG(pep2,&pep2->rg);
260:     RGSetType(pep2->rg,RGINTERVAL);
261: #if defined(PETSC_USE_COMPLEX)
262:     RGIntervalSetEndpoints(pep2->rg,pep->inta,pep->intb,-PETSC_SQRT_MACHINE_EPSILON,PETSC_SQRT_MACHINE_EPSILON);
263: #else
264:     RGIntervalSetEndpoints(pep2->rg,pep->inta,pep->intb,0.0,0.0);
265: #endif
266:     pep2->target = shift;
267:     st = pep2->st;
268:     pep2->st = pep->st;
269:     PEPSolve(pep2);
270:     PEPGetConverged(pep2,&nconv);
271:     if (nconv) {
272:       PEPGetEigenpair(pep2,0,&lambda,NULL,pep2->work[0],NULL);
273:       MatMult(pep->A[1],pep2->work[0],pep2->work[1]);
274:       MatMult(pep->A[2],pep2->work[0],pep2->work[2]);
275:       VecAXPY(pep2->work[1],2.0*lambda*pep->sfactor,pep2->work[2]);
276:       VecDot(pep2->work[1],pep2->work[0],&dot);
277:       PetscInfo2(pep,"lambda=%g, %s type\n",(double)PetscRealPart(lambda),(PetscRealPart(dot)>0.0)?"positive":"negative");
278:       if (!sr->type) sr->type = (PetscRealPart(dot)>0.0)?1:-1;
279:       else {
280:         if (sr->type*PetscRealPart(dot)<0.0) SETERRQ1(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"Different positive/negative type detected in eigenvalue %g",(double)PetscRealPart(lambda));
281:       }
282:     }
283:     pep2->st = st;
284:     PEPDestroy(&pep2);
285:   }
286:   return(0);
287: }

289: PETSC_STATIC_INLINE PetscErrorCode PEPQSliceDiscriminant(PEP pep,Vec u,Vec w,PetscReal *d,PetscReal *smas,PetscReal *smenos)
290: {
291:   PetscReal      ap,bp,cp,dis;
292:   PetscScalar    ts;

296:   MatMult(pep->A[0],u,w);
297:   VecDot(w,u,&ts);
298:   cp = PetscRealPart(ts);
299:   MatMult(pep->A[1],u,w);
300:   VecDot(w,u,&ts);
301:   bp = PetscRealPart(ts);
302:   MatMult(pep->A[2],u,w);
303:   VecDot(w,u,&ts);
304:   ap = PetscRealPart(ts);
305:   dis = bp*bp-4*ap*cp;
306:   if (dis>=0.0 && smas) {
307:     if (ap>0) *smas = (-bp+PetscSqrtReal(dis))/(2*ap);
308:     else if (ap<0) *smas = (-bp-PetscSqrtReal(dis))/(2*ap);
309:     else {
310:       if (bp >0) *smas = -cp/bp;
311:       else *smas = PETSC_MAX_REAL;
312:     }
313:   }
314:   if (dis>=0.0 && smenos) {
315:     if (ap>0) *smenos = (-bp-PetscSqrtReal(dis))/(2*ap);
316:     else if (ap<0) *smenos = (-bp+PetscSqrtReal(dis))/(2*ap);
317:     else {
318:       if (bp<0) *smenos = -cp/bp;
319:       else *smenos = PETSC_MAX_REAL;
320:     }
321:   }
322:   if (d) *d = dis;
323:   return(0);
324: }

326: PETSC_STATIC_INLINE PetscErrorCode PEPQSliceEvaluateQEP(PEP pep,PetscScalar x,Mat M,MatStructure str)
327: {

331:   MatCopy(pep->A[0],M,SAME_NONZERO_PATTERN);
332:   MatAXPY(M,x,pep->A[1],str);
333:   MatAXPY(M,x*x,pep->A[2],str);
334:   return(0);
335: }

337: /*@
338:    PEPCheckDefiniteQEP - Determines if a symmetric/Hermitian quadratic eigenvalue problem
339:    is definite or not.

341:    Logically Collective on pep

343:    Input Parameter:
344: .  pep  - eigensolver context

346:    Output Parameters:
347: +  xi - first computed parameter
348: .  mu - second computed parameter
349: .  definite - flag indicating that the problem is definite
350: -  hyperbolic - flag indicating that the problem is hyperbolic

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

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

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

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

367:    Level: advanced
368: @*/
369: PetscErrorCode PEPCheckDefiniteQEP(PEP pep,PetscReal *xi,PetscReal *mu,PetscInt *definite,PetscInt *hyperbolic)
370: {
372:   PetscRandom    rand;
373:   Vec            u,w;
374:   PetscReal      d=0.0,s=0.0,sp,mut=0.0,omg=0.0,omgp;
375:   PetscInt       k,its=10,hyp=0,check=0,nconv,inertia,n;
376:   Mat            M=NULL;
377:   MatStructure   str;
378:   EPS            eps;
379:   PetscBool      transform,ptypehyp;

382:   if (pep->problem_type!=PEP_HERMITIAN && pep->problem_type!=PEP_HYPERBOLIC) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Only available for Hermitian (or hyperbolic) problems");
383:   ptypehyp = (pep->problem_type==PEP_HYPERBOLIC)? PETSC_TRUE: PETSC_FALSE;
384:   if (!pep->st) { PEPGetST(pep,&pep->st); }
385:   PEPSetDefaultST(pep);
386:   STSetMatrices(pep->st,pep->nmat,pep->A);
387:   MatGetSize(pep->A[0],&n,NULL);
388:   STGetTransform(pep->st,&transform);
389:   STSetTransform(pep->st,PETSC_FALSE);
390:   STSetUp(pep->st);
391:   MatCreateVecs(pep->A[0],&u,&w);
392:   PEPGetBV(pep,&pep->V);
393:   BVGetRandomContext(pep->V,&rand);
394:   VecSetRandom(u,rand);
395:   VecNormalize(u,NULL);
396:   PEPQSliceDiscriminant(pep,u,w,&d,&s,NULL);
397:   if (d<0.0) check = -1;
398:   if (!check) {
399:     EPSCreate(PetscObjectComm((PetscObject)pep),&eps);
400:     EPSSetProblemType(eps,EPS_HEP);
401:     EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
402:     EPSSetTolerances(eps,PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON),PETSC_DECIDE);
403:     MatDuplicate(pep->A[0],MAT_DO_NOT_COPY_VALUES,&M);
404:     STGetMatStructure(pep->st,&str);
405:   }
406:   for (k=0;k<its&&!check;k++) {
407:     PEPQSliceEvaluateQEP(pep,s,M,str);
408:     EPSSetOperators(eps,M,NULL);
409:     EPSSolve(eps);
410:     EPSGetConverged(eps,&nconv);
411:     if (!nconv) break;
412:     EPSGetEigenpair(eps,0,NULL,NULL,u,w);
413:     sp = s;
414:     PEPQSliceDiscriminant(pep,u,w,&d,&s,&omg);
415:     if (d<0.0) {check = -1; break;}
416:     if (PetscAbsReal((s-sp)/s)<100*PETSC_MACHINE_EPSILON) break;
417:     if (s>sp) {hyp = -1;}
418:     mut = 2*s-sp;
419:      PEPQSliceMatGetInertia(pep,mut,&inertia,NULL);
420:     if (inertia == n) {check = 1; break;}
421:   }
422:   for (;k<its&&!check;k++) {
423:     mut = (s-omg)/2;
424:      PEPQSliceMatGetInertia(pep,mut,&inertia,NULL);
425:     if (inertia == n) {check = 1; break;}
426:     if (PetscAbsReal((s-omg)/omg)<100*PETSC_MACHINE_EPSILON) break;
427:     PEPQSliceEvaluateQEP(pep,omg,M,str);
428:     EPSSetOperators(eps,M,NULL);
429:     EPSSolve(eps);
430:     EPSGetConverged(eps,&nconv);
431:     if (!nconv) break;
432:     EPSGetEigenpair(eps,0,NULL,NULL,u,w);
433:     omgp = omg;
434:     PEPQSliceDiscriminant(pep,u,w,&d,NULL,&omg);
435:     if (d<0.0) {check = -1; break;}
436:     if (omg<omgp) hyp = -1;
437:   }
438:   if (check==1) *xi = mut;
439:   if (hyp==-1 && ptypehyp) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Problem does not satisfy hyperbolic test; consider removing the hyperbolicity flag");
440:   if (check==1 && hyp==0) {
441:      PEPQSliceMatGetInertia(pep,PETSC_MAX_REAL,&inertia,NULL);
442:     if (inertia == 0) hyp = 1;
443:     else hyp = -1;
444:   }
445:   if (check==1 && hyp!=1) {
446:     check = 0;
447:     EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
448:     for (;k<its&&!check;k++) {
449:       PEPQSliceEvaluateQEP(pep,s,M,str);
450:       EPSSetOperators(eps,M,NULL);
451:       EPSSolve(eps);
452:       EPSGetConverged(eps,&nconv);
453:       if (!nconv) break;
454:       EPSGetEigenpair(eps,0,NULL,NULL,u,w);
455:       sp = s;
456:       PEPQSliceDiscriminant(pep,u,w,&d,&s,&omg);
457:       if (d<0.0) {check = -1; break;}
458:       if (PetscAbsReal((s-sp)/s)<100*PETSC_MACHINE_EPSILON) break;
459:       mut = 2*s-sp;
460:        PEPQSliceMatGetInertia(pep,mut,&inertia,NULL);
461:       if (inertia == 0) {check = 1; break;}
462:     }
463:     for (;k<its&&!check;k++) {
464:       mut = (s-omg)/2;
465:        PEPQSliceMatGetInertia(pep,mut,&inertia,NULL);
466:       if (inertia == 0) {check = 1; break;}
467:       if (PetscAbsReal((s-omg)/omg)<100*PETSC_MACHINE_EPSILON) break;
468:       PEPQSliceEvaluateQEP(pep,omg,M,str);
469:       EPSSetOperators(eps,M,NULL);
470:       EPSSolve(eps);
471:       EPSGetConverged(eps,&nconv);
472:       if (!nconv) break;
473:       EPSGetEigenpair(eps,0,NULL,NULL,u,w);
474:       PEPQSliceDiscriminant(pep,u,w,&d,NULL,&omg);
475:       if (d<0.0) {check = -1; break;}
476:     }
477:   }
478:   if (check==1) *mu = mut;
479:   *definite = check;
480:   *hyperbolic = hyp;
481:   if (M) { MatDestroy(&M); }
482:   VecDestroy(&u);
483:   VecDestroy(&w);
484:   EPSDestroy(&eps);
485:   STSetTransform(pep->st,transform);
486:   return(0);
487: }

489: /*
490:    Dummy backtransform operation
491:  */
492: static PetscErrorCode PEPBackTransform_Skip(PEP pep)
493: {
495:   return(0);
496: }

498: PetscErrorCode PEPSetUp_STOAR_QSlice(PEP pep)
499: {
501:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
502:   PEP_SR         sr;
503:   PetscInt       ld,i,zeros=0;
504:   SlepcSC        sc;
505:   PetscReal      r;

508:   PEPCheckSinvertCayley(pep);
509:   if (pep->inta==pep->intb) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues unless you provide a computational interval with PEPSetInterval()");
510:   if (pep->intb >= PETSC_MAX_REAL && pep->inta <= PETSC_MIN_REAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"The defined computational interval should have at least one of their sides bounded");
511:   PEPCheckUnsupportedCondition(pep,PEP_FEATURE_STOPPING,PETSC_TRUE," (with spectrum slicing)");
512:   if (pep->tol==PETSC_DEFAULT) {
513: #if defined(PETSC_USE_REAL_SINGLE)
514:     pep->tol = SLEPC_DEFAULT_TOL;
515: #else
516:     /* use tighter tolerance */
517:     pep->tol = SLEPC_DEFAULT_TOL*1e-2;
518: #endif
519:   }
520:   if (ctx->nev==1) ctx->nev = PetscMin(20,pep->n);  /* nev not set, use default value */
521:   if (pep->n>10 && ctx->nev<10) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"nev cannot be less than 10 in spectrum slicing runs");
522:   pep->ops->backtransform = PEPBackTransform_Skip;
523:   if (pep->max_it==PETSC_DEFAULT) pep->max_it = 100;

525:   /* create spectrum slicing context and initialize it */
526:   PEPQSliceResetSR(pep);
527:   PetscNewLog(pep,&sr);
528:   ctx->sr   = sr;
529:   sr->itsKs = 0;
530:   sr->nleap = 0;
531:   sr->sPres = NULL;

533:   if (pep->solvematcoeffs) { PetscFree(pep->solvematcoeffs); }
534:   PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
535:   if (!pep->st) { PEPGetST(pep,&pep->st); }
536:   STSetTransform(pep->st,PETSC_FALSE);
537:   STSetUp(pep->st);

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

541:   /* check presence of ends and finding direction */
542:   if (pep->inta > PETSC_MIN_REAL || pep->intb >= PETSC_MAX_REAL) {
543:     sr->int0 = pep->inta;
544:     sr->int1 = pep->intb;
545:     sr->dir = 1;
546:     if (pep->intb >= PETSC_MAX_REAL) { /* Right-open interval */
547:       sr->hasEnd = PETSC_FALSE;
548:     } else sr->hasEnd = PETSC_TRUE;
549:   } else {
550:     sr->int0 = pep->intb;
551:     sr->int1 = pep->inta;
552:     sr->dir = -1;
553:     sr->hasEnd = PetscNot(pep->inta <= PETSC_MIN_REAL);
554:   }

556:   /* compute inertia0 */
557:   PEPQSliceGetInertia(pep,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL,ctx->hyperbolic?0:1);
558:   if (zeros && (sr->int0==pep->inta || sr->int0==pep->intb)) SETERRQ(((PetscObject)pep)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in the interval endpoint");
559:   if (!ctx->hyperbolic && ctx->checket) {
560:     PEPQSliceCheckEigenvalueType(pep,sr->int0,0.0,PETSC_TRUE);
561:   }

563:   /* compute inertia1 */
564:   PEPQSliceGetInertia(pep,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL,ctx->hyperbolic?0:1);
565:   if (zeros) SETERRQ(((PetscObject)pep)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
566:   if (!ctx->hyperbolic && ctx->checket && sr->hasEnd) {
567:     PEPQSliceCheckEigenvalueType(pep,sr->int1,0.0,PETSC_TRUE);
568:     if (!sr->type && (sr->inertia1-sr->inertia0)) SETERRQ(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"No information of eigenvalue type in Interval");
569:     if (sr->type && !(sr->inertia1-sr->inertia0)) SETERRQ(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"Different positive/negative type detected");
570:     if (sr->dir*(sr->inertia1-sr->inertia0)<0) {
571:       sr->intcorr = -1;
572:       sr->inertia0 = 2*pep->n-sr->inertia0;
573:       sr->inertia1 = 2*pep->n-sr->inertia1;
574:     } else sr->intcorr = 1;
575:   } else {
576:     if (sr->inertia0<=pep->n && sr->inertia1<=pep->n) sr->intcorr = 1;
577:     else if (sr->inertia0>=pep->n && sr->inertia1>=pep->n) sr->intcorr = -1;
578:   }

580:   if (sr->hasEnd) {
581:     sr->dir = -sr->dir; r = sr->int0; sr->int0 = sr->int1; sr->int1 = r;
582:     i = sr->inertia0; sr->inertia0 = sr->inertia1; sr->inertia1 = i;
583:   }

585:   /* number of eigenvalues in interval */
586:   sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
587:   PetscInfo3(pep,"QSlice setup: allocating for %D eigenvalues in [%g,%g]\n",sr->numEigs,(double)pep->inta,(double)pep->intb);
588:   if (sr->numEigs) {
589:     PEPQSliceAllocateSolution(pep);
590:     PEPSetDimensions_Default(pep,ctx->nev,&ctx->ncv,&ctx->mpd);
591:     pep->nev = ctx->nev; pep->ncv = ctx->ncv; pep->mpd = ctx->mpd;
592:     ld   = ctx->ncv+2;
593:     DSSetType(pep->ds,DSGHIEP);
594:     DSSetCompact(pep->ds,PETSC_TRUE);
595:     DSSetExtraRow(pep->ds,PETSC_TRUE);
596:     DSAllocate(pep->ds,ld);
597:     DSGetSlepcSC(pep->ds,&sc);
598:     sc->rg            = NULL;
599:     sc->comparison    = SlepcCompareLargestMagnitude;
600:     sc->comparisonctx = NULL;
601:     sc->map           = NULL;
602:     sc->mapobj        = NULL;
603:   } else {pep->ncv = 0; pep->nev = 0; pep->mpd = 0;}
604:   return(0);
605: }

607: /*
608:    Fills the fields of a shift structure
609: */
610: static PetscErrorCode PEPCreateShift(PEP pep,PetscReal val,PEP_shift neighb0,PEP_shift neighb1)
611: {
613:   PEP_shift      s,*pending2;
614:   PetscInt       i;
615:   PEP_SR         sr;
616:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;

619:   sr = ctx->sr;
620:   PetscNewLog(pep,&s);
621:   s->value = val;
622:   s->neighb[0] = neighb0;
623:   if (neighb0) neighb0->neighb[1] = s;
624:   s->neighb[1] = neighb1;
625:   if (neighb1) neighb1->neighb[0] = s;
626:   s->comp[0] = PETSC_FALSE;
627:   s->comp[1] = PETSC_FALSE;
628:   s->index = -1;
629:   s->neigs = 0;
630:   s->nconv[0] = s->nconv[1] = 0;
631:   s->nsch[0] = s->nsch[1]=0;
632:   /* Inserts in the stack of pending shifts */
633:   /* If needed, the array is resized */
634:   if (sr->nPend >= sr->maxPend) {
635:     sr->maxPend *= 2;
636:     PetscMalloc1(sr->maxPend,&pending2);
637:     PetscLogObjectMemory((PetscObject)pep,sr->maxPend*sizeof(PEP_shift*));
638:     for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
639:     PetscFree(sr->pending);
640:     sr->pending = pending2;
641:   }
642:   sr->pending[sr->nPend++]=s;
643:   return(0);
644: }

646: /* Provides next shift to be computed */
647: static PetscErrorCode PEPExtractShift(PEP pep)
648: {
650:   PetscInt       iner,zeros=0;
651:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
652:   PEP_SR         sr;
653:   PetscReal      newShift,aux;
654:   PEP_shift      sPres;

657:   sr = ctx->sr;
658:   if (sr->nPend > 0) {
659:     if (sr->dirch) {
660:       aux = sr->int1; sr->int1 = sr->int0; sr->int0 = aux;
661:       iner = sr->inertia1; sr->inertia1 = sr->inertia0; sr->inertia0 = iner;
662:       sr->dir *= -1;
663:       PetscFree(sr->s0->neighb[1]);
664:       PetscFree(sr->s0);
665:       sr->nPend--;
666:       PEPCreateShift(pep,sr->int0,NULL,NULL);
667:       sr->sPrev = NULL;
668:       sr->sPres = sr->pending[--sr->nPend];
669:       pep->target = sr->sPres->value;
670:       sr->s0 = sr->sPres;
671:       pep->reason = PEP_CONVERGED_ITERATING;
672:     } else {
673:       sr->sPrev = sr->sPres;
674:       sr->sPres = sr->pending[--sr->nPend];
675:     }
676:     sPres = sr->sPres;
677:     PEPQSliceGetInertia(pep,sPres->value,&iner,ctx->detect?&zeros:NULL,sr->intcorr);
678:     if (zeros) {
679:       newShift = sPres->value*(1.0+SLICE_PTOL);
680:       if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
681:       else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
682:       PEPQSliceGetInertia(pep,newShift,&iner,&zeros,sr->intcorr);
683:       if (zeros) SETERRQ1(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",newShift);
684:       sPres->value = newShift;
685:     }
686:     sr->sPres->inertia = iner;
687:     pep->target = sr->sPres->value;
688:     pep->reason = PEP_CONVERGED_ITERATING;
689:     pep->its = 0;
690:   } else sr->sPres = NULL;
691:   return(0);
692: }

694: /*
695:   Obtains value of subsequent shift
696: */
697: static PetscErrorCode PEPGetNewShiftValue(PEP pep,PetscInt side,PetscReal *newS)
698: {
699:   PetscReal lambda,d_prev;
700:   PetscInt  i,idxP;
701:   PEP_SR    sr;
702:   PEP_shift sPres,s;
703:   PEP_STOAR *ctx=(PEP_STOAR*)pep->data;

706:   sr = ctx->sr;
707:   sPres = sr->sPres;
708:   if (sPres->neighb[side]) {
709:   /* Completing a previous interval */
710:     if (!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0) { /* One of the ends might be too far from eigenvalues */
711:       if (side) *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[sr->indexEig-1]]))/2;
712:       else *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[0]]))/2;
713:     } else *newS=(sPres->value + sPres->neighb[side]->value)/2;
714:   } else { /* (Only for side=1). Creating a new interval. */
715:     if (sPres->neigs==0) {/* No value has been accepted*/
716:       if (sPres->neighb[0]) {
717:         /* Multiplying by 10 the previous distance */
718:         *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
719:         sr->nleap++;
720:         /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
721:         if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Unable to compute the wanted eigenvalues with open interval");
722:       } else { /* First shift */
723:         if (pep->nconv != 0) {
724:           /* Unaccepted values give information for next shift */
725:           idxP=0;/* Number of values left from shift */
726:           for (i=0;i<pep->nconv;i++) {
727:             lambda = PetscRealPart(pep->eigr[i]);
728:             if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
729:             else break;
730:           }
731:           /* Avoiding subtraction of eigenvalues (might be the same).*/
732:           if (idxP>0) {
733:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(pep->eigr[0]))/(idxP+0.3);
734:           } else {
735:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(pep->eigr[pep->nconv-1]))/(pep->nconv+0.3);
736:           }
737:           *newS = sPres->value + ((sr->dir)*d_prev*pep->nev)/2;
738:           sr->dirch = PETSC_FALSE;
739:         } else { /* No values found, no information for next shift */
740:           if (!sr->dirch) {
741:             sr->dirch = PETSC_TRUE;
742:             *newS = sr->int1;
743:           } else SETERRQ(PetscObjectComm((PetscObject)pep),1,"First shift renders no information");
744:         }
745:       }
746:     } else { /* Accepted values found */
747:       sr->dirch = PETSC_FALSE;
748:       sr->nleap = 0;
749:       /* Average distance of values in previous subinterval */
750:       s = sPres->neighb[0];
751:       while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
752:         s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
753:       }
754:       if (s) {
755:         d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
756:       } else { /* First shift. Average distance obtained with values in this shift */
757:         /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
758:         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)) {
759:           d_prev =  PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
760:         } else {
761:           d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
762:         }
763:       }
764:       /* Average distance is used for next shift by adding it to value on the right or to shift */
765:       if ((sr->dir)*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
766:         *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(pep->nev))/2;
767:       } else { /* Last accepted value is on the left of shift. Adding to shift */
768:         *newS = sPres->value + ((sr->dir)*d_prev*(pep->nev))/2;
769:       }
770:     }
771:     /* End of interval can not be surpassed */
772:     if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
773:   }/* of neighb[side]==null */
774:   return(0);
775: }

777: /*
778:   Function for sorting an array of real values
779: */
780: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
781: {
782:   PetscReal re;
783:   PetscInt  i,j,tmp;

786:   if (!prev) for (i=0;i<nr;i++) perm[i] = i;
787:   /* Insertion sort */
788:   for (i=1;i<nr;i++) {
789:     re = PetscRealPart(r[perm[i]]);
790:     j = i-1;
791:     while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
792:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
793:     }
794:   }
795:   return(0);
796: }

798: /* Stores the pairs obtained since the last shift in the global arrays */
799: static PetscErrorCode PEPStoreEigenpairs(PEP pep)
800: {
802:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
803:   PetscReal      lambda,err,*errest;
804:   PetscInt       i,*aux,count=0,ndef,ld,nconv=pep->nconv,d=pep->nmat-1,idx;
805:   PetscBool      iscayley,divide=PETSC_FALSE;
806:   PEP_SR         sr = ctx->sr;
807:   PEP_shift      sPres;
808:   Vec            w,vomega;
809:   Mat            MS;
810:   BV             tV;
811:   PetscScalar    *S,*eigr,*tS,*omega;

814:   sPres = sr->sPres;
815:   sPres->index = sr->indexEig;

817:   if (nconv>sr->ndef0+sr->ndef1) {
818:     /* Back-transform */
819:     STBackTransform(pep->st,nconv,pep->eigr,pep->eigi);
820:     for (i=0;i<nconv;i++) {
821: #if defined(PETSC_USE_COMPLEX)
822:       if (PetscImaginaryPart(pep->eigr[i])) pep->eigr[i] = sr->int0-sr->dir;
823: #else
824:       if (pep->eigi[i]) pep->eigr[i] = sr->int0-sr->dir;
825: #endif
826:     }
827:     PetscObjectTypeCompare((PetscObject)pep->st,STCAYLEY,&iscayley);
828:     /* Sort eigenvalues */
829:     sortRealEigenvalues(pep->eigr,pep->perm,nconv,PETSC_FALSE,sr->dir);
830:     VecCreateSeq(PETSC_COMM_SELF,nconv,&vomega);
831:     BVGetSignature(ctx->V,vomega);
832:     VecGetArray(vomega,&omega);
833:     BVGetSizes(pep->V,NULL,NULL,&ld);
834:     BVTensorGetFactors(ctx->V,NULL,&MS);
835:     MatDenseGetArray(MS,&S);
836:     /* Values stored in global array */
837:     PetscCalloc4(nconv,&eigr,nconv,&errest,nconv*nconv*d,&tS,nconv,&aux);
838:     ndef = sr->ndef0+sr->ndef1;
839:     for (i=0;i<nconv;i++) {
840:       lambda = PetscRealPart(pep->eigr[pep->perm[i]]);
841:       err = pep->errest[pep->perm[i]];
842:       if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
843:         if (sr->indexEig+count-ndef>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Unexpected error in Spectrum Slicing");
844:         PEPQSliceCheckEigenvalueType(pep,lambda,PetscRealPart(omega[pep->perm[i]]),PETSC_FALSE);
845:         eigr[count] = lambda;
846:         errest[count] = err;
847:         if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) sPres->nconv[0]++;
848:         if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) sPres->nconv[1]++;
849:         PetscArraycpy(tS+count*(d*nconv),S+pep->perm[i]*(d*ld),nconv);
850:         PetscArraycpy(tS+count*(d*nconv)+nconv,S+pep->perm[i]*(d*ld)+ld,nconv);
851:         count++;
852:       }
853:     }
854:     VecRestoreArray(vomega,&omega);
855:     VecDestroy(&vomega);
856:     for (i=0;i<count;i++) {
857:       PetscArraycpy(S+i*(d*ld),tS+i*nconv*d,nconv);
858:       PetscArraycpy(S+i*(d*ld)+ld,tS+i*nconv*d+nconv,nconv);
859:     }
860:     MatDenseRestoreArray(MS,&S);
861:     BVTensorRestoreFactors(ctx->V,NULL,&MS);
862:     BVSetActiveColumns(ctx->V,0,count);
863:     BVTensorCompress(ctx->V,count);
864:     if (sr->sPres->nconv[0] && sr->sPres->nconv[1]) {
865:       divide = PETSC_TRUE;
866:       BVTensorGetFactors(ctx->V,NULL,&MS);
867:       MatDenseGetArray(MS,&S);
868:       PetscArrayzero(tS,nconv*nconv*d);
869:       for (i=0;i<count;i++) {
870:         PetscArraycpy(tS+i*nconv*d,S+i*(d*ld),count);
871:         PetscArraycpy(tS+i*nconv*d+nconv,S+i*(d*ld)+ld,count);
872:       }
873:       MatDenseRestoreArray(MS,&S);
874:       BVTensorRestoreFactors(ctx->V,NULL,&MS);
875:       BVSetActiveColumns(pep->V,0,count);
876:       BVDuplicateResize(pep->V,count,&tV);
877:       BVCopy(pep->V,tV);
878:     }
879:     if (sr->sPres->nconv[0]) {
880:       if (divide) {
881:         BVSetActiveColumns(ctx->V,0,sr->sPres->nconv[0]);
882:         BVTensorCompress(ctx->V,sr->sPres->nconv[0]);
883:       }
884:       for (i=0;i<sr->ndef0;i++) aux[i] = sr->idxDef0[i];
885:       for (i=sr->ndef0;i<sr->sPres->nconv[0];i++) aux[i] = sr->indexEig+i-sr->ndef0;
886:       BVTensorGetFactors(ctx->V,NULL,&MS);
887:       MatDenseGetArray(MS,&S);
888:       for (i=0;i<sr->sPres->nconv[0];i++) {
889:         sr->eigr[aux[i]] = eigr[i];
890:         sr->errest[aux[i]] = errest[i];
891:         BVGetColumn(pep->V,i,&w);
892:         BVInsertVec(sr->V,aux[i],w);
893:         BVRestoreColumn(pep->V,i,&w);
894:         idx = sr->ld*d*aux[i];
895:         PetscArrayzero(sr->S+idx,sr->ld*d);
896:         PetscArraycpy(sr->S+idx,S+i*(ld*d),sr->sPres->nconv[0]);
897:         PetscArraycpy(sr->S+idx+sr->ld,S+i*(ld*d)+ld,sr->sPres->nconv[0]);
898:         PetscFree(sr->qinfo[aux[i]].q);
899:         PetscMalloc1(sr->sPres->nconv[0],&sr->qinfo[aux[i]].q);
900:         PetscArraycpy(sr->qinfo[aux[i]].q,aux,sr->sPres->nconv[0]);
901:         sr->qinfo[aux[i]].nq = sr->sPres->nconv[0];
902:       }
903:       MatDenseRestoreArray(MS,&S);
904:       BVTensorRestoreFactors(ctx->V,NULL,&MS);
905:     }

907:     if (sr->sPres->nconv[1]) {
908:       if (divide) {
909:         BVTensorGetFactors(ctx->V,NULL,&MS);
910:         MatDenseGetArray(MS,&S);
911:         for (i=0;i<sr->sPres->nconv[1];i++) {
912:           PetscArraycpy(S+i*(d*ld),tS+(sr->sPres->nconv[0]+i)*nconv*d,count);
913:           PetscArraycpy(S+i*(d*ld)+ld,tS+(sr->sPres->nconv[0]+i)*nconv*d+nconv,count);
914:         }
915:         MatDenseRestoreArray(MS,&S);
916:         BVTensorRestoreFactors(ctx->V,NULL,&MS);
917:         BVSetActiveColumns(pep->V,0,count);
918:         BVCopy(tV,pep->V);
919:         BVSetActiveColumns(ctx->V,0,sr->sPres->nconv[1]);
920:         BVTensorCompress(ctx->V,sr->sPres->nconv[1]);
921:       }
922:       for (i=0;i<sr->ndef1;i++) aux[i] = sr->idxDef1[i];
923:       for (i=sr->ndef1;i<sr->sPres->nconv[1];i++) aux[i] = sr->indexEig+sr->sPres->nconv[0]-sr->ndef0+i-sr->ndef1;
924:       BVTensorGetFactors(ctx->V,NULL,&MS);
925:       MatDenseGetArray(MS,&S);
926:       for (i=0;i<sr->sPres->nconv[1];i++) {
927:         sr->eigr[aux[i]] = eigr[sr->sPres->nconv[0]+i];
928:         sr->errest[aux[i]] = errest[sr->sPres->nconv[0]+i];
929:         BVGetColumn(pep->V,i,&w);
930:         BVInsertVec(sr->V,aux[i],w);
931:         BVRestoreColumn(pep->V,i,&w);
932:         idx = sr->ld*d*aux[i];
933:         PetscArrayzero(sr->S+idx,sr->ld*d);
934:         PetscArraycpy(sr->S+idx,S+i*(ld*d),sr->sPres->nconv[1]);
935:         PetscArraycpy(sr->S+idx+sr->ld,S+i*(ld*d)+ld,sr->sPres->nconv[1]);
936:         PetscFree(sr->qinfo[aux[i]].q);
937:         PetscMalloc1(sr->sPres->nconv[1],&sr->qinfo[aux[i]].q);
938:         PetscArraycpy(sr->qinfo[aux[i]].q,aux,sr->sPres->nconv[1]);
939:         sr->qinfo[aux[i]].nq = sr->sPres->nconv[1];
940:       }
941:       MatDenseRestoreArray(MS,&S);
942:       BVTensorRestoreFactors(ctx->V,NULL,&MS);
943:     }
944:     sPres->neigs = count-sr->ndef0-sr->ndef1;
945:     sr->indexEig += sPres->neigs;
946:     sPres->nconv[0]-= sr->ndef0;
947:     sPres->nconv[1]-= sr->ndef1;
948:     PetscFree4(eigr,errest,tS,aux);
949:   } else {
950:     sPres->neigs = 0;
951:     sPres->nconv[0]= 0;
952:     sPres->nconv[1]= 0;
953:   }
954:   /* Global ordering array updating */
955:   sortRealEigenvalues(sr->eigr,sr->perm,sr->indexEig,PETSC_FALSE,sr->dir);
956:   /* Check for completion */
957:   sPres->comp[0] = PetscNot(sPres->nconv[0] < sPres->nsch[0]);
958:   sPres->comp[1] = PetscNot(sPres->nconv[1] < sPres->nsch[1]);
959:   if (sPres->nconv[0] > sPres->nsch[0] || sPres->nconv[1] > sPres->nsch[1]) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Mismatch between number of values found and information from inertia");
960:   if (divide) { BVDestroy(&tV); }
961:   return(0);
962: }

964: static PetscErrorCode PEPLookForDeflation(PEP pep)
965: {
966:   PetscReal val;
967:   PetscInt  i,count0=0,count1=0;
968:   PEP_shift sPres;
969:   PetscInt  ini,fin;
970:   PEP_SR    sr;
971:   PEP_STOAR *ctx=(PEP_STOAR*)pep->data;

974:   sr = ctx->sr;
975:   sPres = sr->sPres;

977:   if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
978:   else ini = 0;
979:   fin = sr->indexEig;
980:   /* Selection of ends for searching new values */
981:   if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
982:   else sPres->ext[0] = sPres->neighb[0]->value;
983:   if (!sPres->neighb[1]) {
984:     if (sr->hasEnd) sPres->ext[1] = sr->int1;
985:     else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
986:   } else sPres->ext[1] = sPres->neighb[1]->value;
987:   /* Selection of values between right and left ends */
988:   for (i=ini;i<fin;i++) {
989:     val=PetscRealPart(sr->eigr[sr->perm[i]]);
990:     /* Values to the right of left shift */
991:     if ((sr->dir)*(val - sPres->ext[1]) < 0) {
992:       if ((sr->dir)*(val - sPres->value) < 0) count0++;
993:       else count1++;
994:     } else break;
995:   }
996:   /* The number of values on each side are found */
997:   if (sPres->neighb[0]) {
998:     sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
999:     if (sPres->nsch[0]<0) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Mismatch between number of values found and information from inertia");
1000:   } else sPres->nsch[0] = 0;

1002:   if (sPres->neighb[1]) {
1003:     sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
1004:     if (sPres->nsch[1]<0) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Mismatch between number of values found and information from inertia");
1005:   } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);

1007:   /* Completing vector of indexes for deflation */
1008:   for (i=0;i<count0;i++) sr->idxDef0[i] = sr->perm[ini+i];
1009:   sr->ndef0 = count0;
1010:   for (i=0;i<count1;i++) sr->idxDef1[i] = sr->perm[ini+count0+i];
1011:   sr->ndef1 = count1;
1012:   return(0);
1013: }

1015: /*
1016:   Compute a run of Lanczos iterations
1017: */
1018: static PetscErrorCode PEPSTOARrun_QSlice(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,Vec *t_)
1019: {
1021:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;
1022:   PetscInt       i,j,m=*M,l,lock;
1023:   PetscInt       lds,d,ld,offq,nqt,ldds;
1024:   Vec            v=t_[0],t=t_[1],q=t_[2];
1025:   PetscReal      norm,sym=0.0,fro=0.0,*f;
1026:   PetscScalar    *y,*S,sigma;
1027:   PetscBLASInt   j_,one=1;
1028:   PetscBool      lindep;
1029:   Mat            MS;

1032:   PetscMalloc1(*M,&y);
1033:   BVGetSizes(pep->V,NULL,NULL,&ld);
1034:   BVTensorGetDegree(ctx->V,&d);
1035:   BVGetActiveColumns(pep->V,&lock,&nqt);
1036:   lds = d*ld;
1037:   offq = ld;
1038:   DSGetLeadingDimension(pep->ds,&ldds);

1040:   *breakdown = PETSC_FALSE; /* ----- */
1041:   STGetShift(pep->st,&sigma);
1042:   DSGetDimensions(pep->ds,NULL,NULL,&l,NULL,NULL);
1043:   BVSetActiveColumns(ctx->V,0,m);
1044:   BVSetActiveColumns(pep->V,0,nqt);
1045:   for (j=k;j<m;j++) {
1046:     /* apply operator */
1047:     BVTensorGetFactors(ctx->V,NULL,&MS);
1048:     MatDenseGetArray(MS,&S);
1049:     BVGetColumn(pep->V,nqt,&t);
1050:     BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
1051:     MatMult(pep->A[1],v,q);
1052:     MatMult(pep->A[2],v,t);
1053:     VecAXPY(q,sigma*pep->sfactor,t);
1054:     VecScale(q,pep->sfactor);
1055:     BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
1056:     MatMult(pep->A[2],v,t);
1057:     VecAXPY(q,pep->sfactor*pep->sfactor,t);
1058:     STMatSolve(pep->st,q,t);
1059:     VecScale(t,-1.0);
1060:     BVRestoreColumn(pep->V,nqt,&t);

1062:     /* orthogonalize */
1063:     BVOrthogonalizeColumn(pep->V,nqt,S+(j+1)*lds,&norm,&lindep);
1064:     if (!lindep) {
1065:       *(S+(j+1)*lds+nqt) = norm;
1066:       BVScaleColumn(pep->V,nqt,1.0/norm);
1067:       nqt++;
1068:     }
1069:     for (i=0;i<nqt;i++) *(S+(j+1)*lds+offq+i) = *(S+j*lds+i)+sigma*(*(S+(j+1)*lds+i));
1070:     BVSetActiveColumns(pep->V,0,nqt);
1071:     MatDenseRestoreArray(MS,&S);
1072:     BVTensorRestoreFactors(ctx->V,NULL,&MS);

1074:     /* level-2 orthogonalization */
1075:     BVOrthogonalizeColumn(ctx->V,j+1,y,&norm,&lindep);
1076:     a[j] = PetscRealPart(y[j]);
1077:     omega[j+1] = (norm > 0)?1.0:-1.0;
1078:     BVScaleColumn(ctx->V,j+1,1.0/norm);
1079:     b[j] = PetscAbsReal(norm);

1081:     /* check symmetry */
1082:     DSGetArrayReal(pep->ds,DS_MAT_T,&f);
1083:     if (j==k) {
1084:       for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ldds+i]);
1085:       for (i=0;i<l;i++) y[i] = 0.0;
1086:     }
1087:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
1088:     if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
1089:     PetscBLASIntCast(j,&j_);
1090:     sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
1091:     fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
1092:     if (j>0) fro = SlepcAbs(fro,b[j-1]);
1093:     if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
1094:       *symmlost = PETSC_TRUE;
1095:       *M=j;
1096:       break;
1097:     }
1098:   }
1099:   BVSetActiveColumns(pep->V,lock,nqt);
1100:   BVSetActiveColumns(ctx->V,0,*M);
1101:   PetscFree(y);
1102:   return(0);
1103: }

1105: static PetscErrorCode PEPSTOAR_QSlice(PEP pep,Mat B)
1106: {
1108:   PEP_STOAR      *ctx = (PEP_STOAR*)pep->data;
1109:   PetscInt       j,k,l,nv=0,ld,ldds,t,nq=0,idx;
1110:   PetscInt       nconv=0,deg=pep->nmat-1,count0=0,count1=0;
1111:   PetscScalar    *om,sigma,*back,*S,*pQ;
1112:   PetscReal      beta,norm=1.0,*omega,*a,*b,eta,lambda;
1113:   PetscBool      breakdown,symmlost=PETSC_FALSE,sinv,falselock=PETSC_TRUE;
1114:   Mat            MS,MQ;
1115:   Vec            v,vomega;
1116:   PEP_SR         sr;
1117:   BVOrthogType   otype;
1118:   BVOrthogBlockType obtype;

1121:   /* Resize if needed for deflating vectors  */
1122:   sr = ctx->sr;
1123:   sigma = sr->sPres->value;
1124:   k = sr->ndef0+sr->ndef1;
1125:   pep->ncv = ctx->ncv+k;
1126:   pep->nev = ctx->nev+k;
1127:   PEPAllocateSolution(pep,3);
1128:   BVDestroy(&ctx->V);
1129:   BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
1130:   BVGetOrthogonalization(pep->V,&otype,NULL,&eta,&obtype);
1131:   BVSetOrthogonalization(ctx->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
1132:   DSAllocate(pep->ds,pep->ncv+2);
1133:   PetscMalloc1(pep->ncv,&back);
1134:   DSGetLeadingDimension(pep->ds,&ldds);
1135:   BVSetMatrix(ctx->V,B,PETSC_TRUE);
1136:   if (ctx->lock) {
1137:     /* undocumented option to use a cheaper locking instead of the true locking */
1138:     PetscOptionsGetBool(NULL,NULL,"-pep_stoar_falselocking",&falselock,NULL);
1139:   } else SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"A locking variant is needed for spectrum slicing");
1140:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
1141:   RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
1142:   STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);

1144:   /* Get the starting Arnoldi vector */
1145:   BVSetActiveColumns(pep->V,0,1);
1146:   BVTensorBuildFirstColumn(ctx->V,pep->nini);
1147:   BVSetActiveColumns(ctx->V,0,1);
1148:   if (k) {
1149:     /* Insert deflated vectors */
1150:     BVSetActiveColumns(pep->V,0,0);
1151:     idx = sr->ndef0?sr->idxDef0[0]:sr->idxDef1[0];
1152:     for (j=0;j<k;j++) {
1153:       BVGetColumn(pep->V,j,&v);
1154:       BVCopyVec(sr->V,sr->qinfo[idx].q[j],v);
1155:       BVRestoreColumn(pep->V,j,&v);
1156:     }
1157:     /* Update innerproduct matrix */
1158:     BVSetActiveColumns(ctx->V,0,0);
1159:     BVTensorGetFactors(ctx->V,NULL,&MS);
1160:     BVSetActiveColumns(pep->V,0,k);
1161:     BVTensorRestoreFactors(ctx->V,NULL,&MS);

1163:     BVGetSizes(pep->V,NULL,NULL,&ld);
1164:     BVTensorGetFactors(ctx->V,NULL,&MS);
1165:     MatDenseGetArray(MS,&S);
1166:     for (j=0;j<sr->ndef0;j++) {
1167:       PetscArrayzero(S+j*ld*deg,ld*deg);
1168:       PetscArraycpy(S+j*ld*deg,sr->S+sr->idxDef0[j]*sr->ld*deg,k);
1169:       PetscArraycpy(S+j*ld*deg+ld,sr->S+sr->idxDef0[j]*sr->ld*deg+sr->ld,k);
1170:       pep->eigr[j] = sr->eigr[sr->idxDef0[j]];
1171:       pep->errest[j] = sr->errest[sr->idxDef0[j]];
1172:     }
1173:     for (j=0;j<sr->ndef1;j++) {
1174:       PetscArrayzero(S+(j+sr->ndef0)*ld*deg,ld*deg);
1175:       PetscArraycpy(S+(j+sr->ndef0)*ld*deg,sr->S+sr->idxDef1[j]*sr->ld*deg,k);
1176:       PetscArraycpy(S+(j+sr->ndef0)*ld*deg+ld,sr->S+sr->idxDef1[j]*sr->ld*deg+sr->ld,k);
1177:       pep->eigr[j+sr->ndef0] = sr->eigr[sr->idxDef1[j]];
1178:       pep->errest[j+sr->ndef0] = sr->errest[sr->idxDef1[j]];
1179:     }
1180:     MatDenseRestoreArray(MS,&S);
1181:     BVTensorRestoreFactors(ctx->V,NULL,&MS);
1182:     BVSetActiveColumns(ctx->V,0,k+1);
1183:     VecCreateSeq(PETSC_COMM_SELF,k+1,&vomega);
1184:     VecGetArray(vomega,&om);
1185:     for (j=0;j<k;j++) {
1186:       BVOrthogonalizeColumn(ctx->V,j,NULL,&norm,NULL);
1187:       BVScaleColumn(ctx->V,j,1/norm);
1188:       om[j] = (norm>=0.0)?1.0:-1.0;
1189:     }
1190:     BVTensorGetFactors(ctx->V,NULL,&MS);
1191:     MatDenseGetArray(MS,&S);
1192:     for (j=0;j<deg;j++) {
1193:       BVSetRandomColumn(pep->V,k+j);
1194:       BVOrthogonalizeColumn(pep->V,k+j,S+k*ld*deg+j*ld,&norm,NULL);
1195:       BVScaleColumn(pep->V,k+j,1.0/norm);
1196:       S[k*ld*deg+j*ld+k+j] = norm;
1197:     }
1198:     MatDenseRestoreArray(MS,&S);
1199:     BVSetActiveColumns(pep->V,0,k+deg);
1200:     BVTensorRestoreFactors(ctx->V,NULL,&MS);
1201:     BVOrthogonalizeColumn(ctx->V,k,NULL,&norm,NULL);
1202:     BVScaleColumn(ctx->V,k,1.0/norm);
1203:     om[k] = (norm>=0.0)?1.0:-1.0;
1204:     VecRestoreArray(vomega,&om);
1205:     BVSetSignature(ctx->V,vomega);
1206:     DSGetArrayReal(pep->ds,DS_MAT_T,&a);
1207:     VecGetArray(vomega,&om);
1208:     for (j=0;j<k;j++) a[j] = PetscRealPart(om[j]/(pep->eigr[j]-sigma));
1209:     VecRestoreArray(vomega,&om);
1210:     VecDestroy(&vomega);
1211:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
1212:     DSGetArray(pep->ds,DS_MAT_Q,&pQ);
1213:     PetscArrayzero(pQ,ldds*k);
1214:     for (j=0;j<k;j++) pQ[j+j*ldds] = 1.0;
1215:     DSRestoreArray(pep->ds,DS_MAT_Q,&pQ);
1216:   }
1217:   BVSetActiveColumns(ctx->V,0,k+1);
1218:   DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
1219:   VecCreateSeq(PETSC_COMM_SELF,k+1,&vomega);
1220:   BVGetSignature(ctx->V,vomega);
1221:   VecGetArray(vomega,&om);
1222:   for (j=0;j<k+1;j++) omega[j] = PetscRealPart(om[j]);
1223:   VecRestoreArray(vomega,&om);
1224:   DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
1225:   VecDestroy(&vomega);

1227:   PetscInfo7(pep,"Start STOAR: sigma=%g in [%g,%g], for deflation: left=%D right=%D, searching: left=%D right=%D\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]);

1229:   /* Restart loop */
1230:   l = 0;
1231:   pep->nconv = k;
1232:   while (pep->reason == PEP_CONVERGED_ITERATING) {
1233:     pep->its++;
1234:     DSGetArrayReal(pep->ds,DS_MAT_T,&a);
1235:     b = a+ldds;
1236:     DSGetArrayReal(pep->ds,DS_MAT_D,&omega);

1238:     /* Compute an nv-step Lanczos factorization */
1239:     nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
1240:     PEPSTOARrun_QSlice(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,pep->work);
1241:     beta = b[nv-1];
1242:     if (symmlost && nv==pep->nconv+l) {
1243:       pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
1244:       pep->nconv = nconv;
1245:       PetscInfo2(pep,"Symmetry lost in STOAR sigma=%g nconv=%D\n",(double)sr->sPres->value,nconv);
1246:       if (falselock || !ctx->lock) {
1247:         BVSetActiveColumns(ctx->V,0,pep->nconv);
1248:         BVTensorCompress(ctx->V,0);
1249:       }
1250:       break;
1251:     }
1252:     DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
1253:     DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
1254:     DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
1255:     if (l==0) {
1256:       DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
1257:     } else {
1258:       DSSetState(pep->ds,DS_STATE_RAW);
1259:     }

1261:     /* Solve projected problem */
1262:     DSSolve(pep->ds,pep->eigr,pep->eigi);
1263:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
1264:     DSUpdateExtraRow(pep->ds);
1265:     DSSynchronize(pep->ds,pep->eigr,pep->eigi);

1267:     /* Check convergence */
1268:     /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
1269:     norm = 1.0;
1270:     DSGetDimensions(pep->ds,NULL,NULL,NULL,NULL,&t);
1271:     PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
1272:     (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
1273:     for (j=0;j<k;j++) back[j] = pep->eigr[j];
1274:     STBackTransform(pep->st,k,back,pep->eigi);
1275:     count0=count1=0;
1276:     for (j=0;j<k;j++) {
1277:       lambda = PetscRealPart(back[j]);
1278:       if (((sr->dir)*(sr->sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sr->sPres->ext[0]) > 0)) count0++;
1279:       if (((sr->dir)*(lambda - sr->sPres->value) > 0) && ((sr->dir)*(sr->sPres->ext[1] - lambda) > 0)) count1++;
1280:     }
1281:     if ((count0-sr->ndef0 >= sr->sPres->nsch[0]) && (count1-sr->ndef1 >= sr->sPres->nsch[1])) pep->reason = PEP_CONVERGED_TOL;
1282:     /* Update l */
1283:     if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
1284:     else {
1285:       l = PetscMax(1,(PetscInt)((nv-k)/2));
1286:       l = PetscMin(l,t);
1287:       DSGetTruncateSize(pep->ds,k,t,&l);
1288:       if (!breakdown) {
1289:         /* Prepare the Rayleigh quotient for restart */
1290:         DSTruncate(pep->ds,k+l,PETSC_FALSE);
1291:       }
1292:     }
1293:     nconv = k;
1294:     if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
1295:     if (l) { PetscInfo1(pep,"Preparing to restart keeping l=%D vectors\n",l); }

1297:     /* Update S */
1298:     DSGetMat(pep->ds,DS_MAT_Q,&MQ);
1299:     BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
1300:     MatDestroy(&MQ);

1302:     /* Copy last column of S */
1303:     BVCopyColumn(ctx->V,nv,k+l);
1304:     DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
1305:     VecCreateSeq(PETSC_COMM_SELF,k+l,&vomega);
1306:     VecGetArray(vomega,&om);
1307:     for (j=0;j<k+l;j++) om[j] = omega[j];
1308:     VecRestoreArray(vomega,&om);
1309:     BVSetActiveColumns(ctx->V,0,k+l);
1310:     BVSetSignature(ctx->V,vomega);
1311:     VecDestroy(&vomega);
1312:     DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);

1314:     if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
1315:       /* stop if breakdown */
1316:       PetscInfo2(pep,"Breakdown TOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
1317:       pep->reason = PEP_DIVERGED_BREAKDOWN;
1318:     }
1319:     if (pep->reason != PEP_CONVERGED_ITERATING) l--;
1320:     BVGetActiveColumns(pep->V,NULL,&nq);
1321:     if (k+l+deg<=nq) {
1322:       BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
1323:       if (!falselock && ctx->lock) {
1324:         BVTensorCompress(ctx->V,k-pep->nconv);
1325:       } else {
1326:         BVTensorCompress(ctx->V,0);
1327:       }
1328:     }
1329:     pep->nconv = k;
1330:     PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
1331:   }
1332:   sr->itsKs += pep->its;
1333:   if (pep->nconv>0) {
1334:     BVSetActiveColumns(ctx->V,0,pep->nconv);
1335:     BVGetActiveColumns(pep->V,NULL,&nq);
1336:     BVSetActiveColumns(pep->V,0,nq);
1337:     if (nq>pep->nconv) {
1338:       BVTensorCompress(ctx->V,pep->nconv);
1339:       BVSetActiveColumns(pep->V,0,pep->nconv);
1340:     }
1341:     for (j=0;j<pep->nconv;j++) {
1342:       pep->eigr[j] *= pep->sfactor;
1343:       pep->eigi[j] *= pep->sfactor;
1344:     }
1345:   }
1346:   PetscInfo4(pep,"Finished STOAR: nconv=%D (deflated=%D, left=%D, right=%D)\n",pep->nconv,sr->ndef0+sr->ndef1,count0-sr->ndef0,count1-sr->ndef1);
1347:   STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
1348:   RGPopScale(pep->rg);

1350:   if (pep->reason == PEP_DIVERGED_SYMMETRY_LOST && nconv<sr->ndef0+sr->ndef1) SETERRQ1(PetscObjectComm((PetscObject)pep),1,"Symmetry lost at sigma=%g",(double)sr->sPres->value);
1351:   if (pep->reason == PEP_DIVERGED_SYMMETRY_LOST && nconv==sr->ndef0+sr->ndef1) {
1352:     if (++sr->symmlost>10) SETERRQ1(PetscObjectComm((PetscObject)pep),1,"Symmetry lost at sigma=%g",(double)sr->sPres->value);
1353:   } else sr->symmlost = 0;

1355:   DSTruncate(pep->ds,pep->nconv,PETSC_TRUE);
1356:   PetscFree(back);
1357:   return(0);
1358: }

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

1362: static PetscErrorCode PEPQSliceGetInertias(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
1363: {
1364:   PetscErrorCode  ierr;
1365:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
1366:   PEP_SR          sr=ctx->sr;
1367:   PetscInt        i=0,j,tmpi;
1368:   PetscReal       v,tmpr;
1369:   PEP_shift       s;

1372:   if (!pep->state) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSetUp() first");
1373:   if (!sr) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see PEPSetInterval()");
1374:   if (!sr->s0) {  /* PEPSolve not called yet */
1375:     *n = 2;
1376:   } else {
1377:     *n = 1;
1378:     s = sr->s0;
1379:     while (s) {
1380:       (*n)++;
1381:       s = s->neighb[1];
1382:     }
1383:   }
1384:   PetscMalloc1(*n,shifts);
1385:   PetscMalloc1(*n,inertias);
1386:   if (!sr->s0) {  /* PEPSolve not called yet */
1387:     (*shifts)[0]   = sr->int0;
1388:     (*shifts)[1]   = sr->int1;
1389:     (*inertias)[0] = sr->inertia0;
1390:     (*inertias)[1] = sr->inertia1;
1391:   } else {
1392:     s = sr->s0;
1393:     while (s) {
1394:       (*shifts)[i]     = s->value;
1395:       (*inertias)[i++] = s->inertia;
1396:       s = s->neighb[1];
1397:     }
1398:     (*shifts)[i]   = sr->int1;
1399:     (*inertias)[i] = sr->inertia1;
1400:   }
1401:   /* remove possible duplicate in last position */
1402:   if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
1403:   /* sort result */
1404:   for (i=0;i<*n;i++) {
1405:     v = (*shifts)[i];
1406:     for (j=i+1;j<*n;j++) {
1407:       if (v > (*shifts)[j]) {
1408:         SWAP((*shifts)[i],(*shifts)[j],tmpr);
1409:         SWAP((*inertias)[i],(*inertias)[j],tmpi);
1410:         v = (*shifts)[i];
1411:       }
1412:     }
1413:   }
1414:   return(0);
1415: }

1417: PetscErrorCode PEPSolve_STOAR_QSlice(PEP pep)
1418: {
1420:   PetscInt       i,j,ti,deg=pep->nmat-1;
1421:   PetscReal      newS;
1422:   PEP_STOAR      *ctx=(PEP_STOAR*)pep->data;
1423:   PEP_SR         sr=ctx->sr;
1424:   Mat            S,B;
1425:   PetscScalar    *pS;

1428:   PetscCitationsRegister(citation,&cited);

1430:   /* Only with eigenvalues present in the interval ...*/
1431:   if (sr->numEigs==0) {
1432:     pep->reason = PEP_CONVERGED_TOL;
1433:     return(0);
1434:   }

1436:   /* Inner product matrix */
1437:   PEPSTOARSetUpInnerMatrix(pep,&B);

1439:   /* Array of pending shifts */
1440:   sr->maxPend = 100; /* Initial size */
1441:   sr->nPend = 0;
1442:   PetscMalloc1(sr->maxPend,&sr->pending);
1443:   PetscLogObjectMemory((PetscObject)pep,sr->maxPend*sizeof(PEP_shift*));
1444:   PEPCreateShift(pep,sr->int0,NULL,NULL);
1445:   /* extract first shift */
1446:   sr->sPrev = NULL;
1447:   sr->sPres = sr->pending[--sr->nPend];
1448:   sr->sPres->inertia = sr->inertia0;
1449:   pep->target = sr->sPres->value;
1450:   sr->s0 = sr->sPres;
1451:   sr->indexEig = 0;

1453:   for (i=0;i<sr->numEigs;i++) {
1454:     sr->eigr[i]   = 0.0;
1455:     sr->eigi[i]   = 0.0;
1456:     sr->errest[i] = 0.0;
1457:     sr->perm[i]   = i;
1458:   }
1459:   /* Vectors for deflation */
1460:   PetscMalloc2(sr->numEigs,&sr->idxDef0,sr->numEigs,&sr->idxDef1);
1461:   PetscLogObjectMemory((PetscObject)pep,2*sr->numEigs*sizeof(PetscInt));
1462:   sr->indexEig = 0;
1463:   while (sr->sPres) {
1464:     /* Search for deflation */
1465:     PEPLookForDeflation(pep);
1466:     /* KrylovSchur */
1467:     PEPSTOAR_QSlice(pep,B);

1469:     PEPStoreEigenpairs(pep);
1470:     /* Select new shift */
1471:     if (!sr->sPres->comp[1]) {
1472:       PEPGetNewShiftValue(pep,1,&newS);
1473:       PEPCreateShift(pep,newS,sr->sPres,sr->sPres->neighb[1]);
1474:     }
1475:     if (!sr->sPres->comp[0]) {
1476:       /* Completing earlier interval */
1477:       PEPGetNewShiftValue(pep,0,&newS);
1478:       PEPCreateShift(pep,newS,sr->sPres->neighb[0],sr->sPres);
1479:     }
1480:     /* Preparing for a new search of values */
1481:     PEPExtractShift(pep);
1482:   }

1484:   /* Updating pep values prior to exit */
1485:   PetscFree2(sr->idxDef0,sr->idxDef1);
1486:   PetscFree(sr->pending);
1487:   pep->nconv  = sr->indexEig;
1488:   pep->reason = PEP_CONVERGED_TOL;
1489:   pep->its    = sr->itsKs;
1490:   pep->nev    = sr->indexEig;
1491:   MatCreateSeqDense(PETSC_COMM_SELF,pep->nconv,pep->nconv,NULL,&S);
1492:   MatDenseGetArray(S,&pS);
1493:   for (i=0;i<pep->nconv;i++) {
1494:     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);
1495:   }
1496:   MatDenseRestoreArray(S,&pS);
1497:   BVSetActiveColumns(sr->V,0,pep->nconv);
1498:   BVMultInPlace(sr->V,S,0,pep->nconv);
1499:   MatDestroy(&S);
1500:   BVDestroy(&pep->V);
1501:   pep->V = sr->V;
1502:   PetscFree4(pep->eigr,pep->eigi,pep->errest,pep->perm);
1503:   pep->eigr   = sr->eigr;
1504:   pep->eigi   = sr->eigi;
1505:   pep->perm   = sr->perm;
1506:   pep->errest = sr->errest;
1507:   if (sr->dir<0) {
1508:     for (i=0;i<pep->nconv/2;i++) {
1509:       ti = sr->perm[i]; sr->perm[i] = sr->perm[pep->nconv-1-i]; sr->perm[pep->nconv-1-i] = ti;
1510:     }
1511:   }
1512:   PetscFree(ctx->inertias);
1513:   PetscFree(ctx->shifts);
1514:   MatDestroy(&B);
1515:   PEPQSliceGetInertias(pep,&ctx->nshifts,&ctx->shifts,&ctx->inertias);
1516:   return(0);
1517: }