Actual source code: qslice.c
slepc-3.13.4 2020-09-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: SLEPc polynomial eigensolver: "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. (in press), 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 = \"IP\",\n"
37: " number = \"x\",\n"
38: " pages = \"xx--xx\",\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: PetscBool issinv;
506: PetscReal r;
509: 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");
510: PetscObjectTypeCompareAny((PetscObject)pep->st,&issinv,STSINVERT,STCAYLEY,"");
511: if (!issinv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
512: if (pep->tol==PETSC_DEFAULT) pep->tol = SLEPC_DEFAULT_TOL*1e-2; /* use tighter tolerance */
513: if (ctx->nev==1) ctx->nev = PetscMin(20,pep->n); /* nev not set, use default value */
514: 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");
515: pep->ops->backtransform = PEPBackTransform_Skip;
516: if (pep->max_it==PETSC_DEFAULT) pep->max_it = 100;
518: /* create spectrum slicing context and initialize it */
519: PEPQSliceResetSR(pep);
520: PetscNewLog(pep,&sr);
521: ctx->sr = sr;
522: sr->itsKs = 0;
523: sr->nleap = 0;
524: sr->sPres = NULL;
526: if (pep->solvematcoeffs) { PetscFree(pep->solvematcoeffs); }
527: PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
528: if (!pep->st) { PEPGetST(pep,&pep->st); }
529: STSetTransform(pep->st,PETSC_FALSE);
530: STSetUp(pep->st);
532: ctx->hyperbolic = (pep->problem_type==PEP_HYPERBOLIC)? PETSC_TRUE: PETSC_FALSE;
534: /* check presence of ends and finding direction */
535: if (pep->inta > PETSC_MIN_REAL || pep->intb >= PETSC_MAX_REAL) {
536: sr->int0 = pep->inta;
537: sr->int1 = pep->intb;
538: sr->dir = 1;
539: if (pep->intb >= PETSC_MAX_REAL) { /* Right-open interval */
540: sr->hasEnd = PETSC_FALSE;
541: } else sr->hasEnd = PETSC_TRUE;
542: } else {
543: sr->int0 = pep->intb;
544: sr->int1 = pep->inta;
545: sr->dir = -1;
546: sr->hasEnd = PetscNot(pep->inta <= PETSC_MIN_REAL);
547: }
549: /* compute inertia0 */
550: PEPQSliceGetInertia(pep,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL,ctx->hyperbolic?0:1);
551: 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");
552: if (!ctx->hyperbolic && ctx->checket) {
553: PEPQSliceCheckEigenvalueType(pep,sr->int0,0.0,PETSC_TRUE);
554: }
556: /* compute inertia1 */
557: PEPQSliceGetInertia(pep,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL,ctx->hyperbolic?0:1);
558: if (zeros) SETERRQ(((PetscObject)pep)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
559: if (!ctx->hyperbolic && ctx->checket && sr->hasEnd) {
560: PEPQSliceCheckEigenvalueType(pep,sr->int1,0.0,PETSC_TRUE);
561: if (!sr->type && (sr->inertia1-sr->inertia0)) SETERRQ(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"No information of eigenvalue type in Interval");
562: if (sr->type && !(sr->inertia1-sr->inertia0)) SETERRQ(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"Different positive/negative type detected");
563: if (sr->dir*(sr->inertia1-sr->inertia0)<0) {
564: sr->intcorr = -1;
565: sr->inertia0 = 2*pep->n-sr->inertia0;
566: sr->inertia1 = 2*pep->n-sr->inertia1;
567: } else sr->intcorr = 1;
568: } else {
569: if (sr->inertia0<=pep->n && sr->inertia1<=pep->n) sr->intcorr = 1;
570: else if (sr->inertia0>=pep->n && sr->inertia1>=pep->n) sr->intcorr = -1;
571: }
573: if (sr->hasEnd) {
574: sr->dir = -sr->dir; r = sr->int0; sr->int0 = sr->int1; sr->int1 = r;
575: i = sr->inertia0; sr->inertia0 = sr->inertia1; sr->inertia1 = i;
576: }
578: /* number of eigenvalues in interval */
579: sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
580: PetscInfo3(pep,"QSlice setup: allocating for %D eigenvalues in [%g,%g]\n",sr->numEigs,(double)pep->inta,(double)pep->intb);
581: if (sr->numEigs) {
582: PEPQSliceAllocateSolution(pep);
583: PEPSetDimensions_Default(pep,ctx->nev,&ctx->ncv,&ctx->mpd);
584: pep->nev = ctx->nev; pep->ncv = ctx->ncv; pep->mpd = ctx->mpd;
585: ld = ctx->ncv+2;
586: DSSetType(pep->ds,DSGHIEP);
587: DSSetCompact(pep->ds,PETSC_TRUE);
588: DSAllocate(pep->ds,ld);
589: DSGetSlepcSC(pep->ds,&sc);
590: sc->rg = NULL;
591: sc->comparison = SlepcCompareLargestMagnitude;
592: sc->comparisonctx = NULL;
593: sc->map = NULL;
594: sc->mapobj = NULL;
595: }
596: return(0);
597: }
599: /*
600: Fills the fields of a shift structure
601: */
602: static PetscErrorCode PEPCreateShift(PEP pep,PetscReal val,PEP_shift neighb0,PEP_shift neighb1)
603: {
605: PEP_shift s,*pending2;
606: PetscInt i;
607: PEP_SR sr;
608: PEP_STOAR *ctx=(PEP_STOAR*)pep->data;
611: sr = ctx->sr;
612: PetscNewLog(pep,&s);
613: s->value = val;
614: s->neighb[0] = neighb0;
615: if (neighb0) neighb0->neighb[1] = s;
616: s->neighb[1] = neighb1;
617: if (neighb1) neighb1->neighb[0] = s;
618: s->comp[0] = PETSC_FALSE;
619: s->comp[1] = PETSC_FALSE;
620: s->index = -1;
621: s->neigs = 0;
622: s->nconv[0] = s->nconv[1] = 0;
623: s->nsch[0] = s->nsch[1]=0;
624: /* Inserts in the stack of pending shifts */
625: /* If needed, the array is resized */
626: if (sr->nPend >= sr->maxPend) {
627: sr->maxPend *= 2;
628: PetscMalloc1(sr->maxPend,&pending2);
629: PetscLogObjectMemory((PetscObject)pep,sizeof(PEP_shift));
630: for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
631: PetscFree(sr->pending);
632: sr->pending = pending2;
633: }
634: sr->pending[sr->nPend++]=s;
635: return(0);
636: }
638: /* Provides next shift to be computed */
639: static PetscErrorCode PEPExtractShift(PEP pep)
640: {
642: PetscInt iner,zeros=0;
643: PEP_STOAR *ctx=(PEP_STOAR*)pep->data;
644: PEP_SR sr;
645: PetscReal newShift,aux;
646: PEP_shift sPres;
649: sr = ctx->sr;
650: if (sr->nPend > 0) {
651: if (sr->dirch) {
652: aux = sr->int1; sr->int1 = sr->int0; sr->int0 = aux;
653: iner = sr->inertia1; sr->inertia1 = sr->inertia0; sr->inertia0 = iner;
654: sr->dir *= -1;
655: PetscFree(sr->s0->neighb[1]);
656: PetscFree(sr->s0);
657: sr->nPend--;
658: PEPCreateShift(pep,sr->int0,NULL,NULL);
659: sr->sPrev = NULL;
660: sr->sPres = sr->pending[--sr->nPend];
661: pep->target = sr->sPres->value;
662: sr->s0 = sr->sPres;
663: pep->reason = PEP_CONVERGED_ITERATING;
664: } else {
665: sr->sPrev = sr->sPres;
666: sr->sPres = sr->pending[--sr->nPend];
667: }
668: sPres = sr->sPres;
669: PEPQSliceGetInertia(pep,sPres->value,&iner,ctx->detect?&zeros:NULL,sr->intcorr);
670: if (zeros) {
671: newShift = sPres->value*(1.0+SLICE_PTOL);
672: if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
673: else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
674: PEPQSliceGetInertia(pep,newShift,&iner,&zeros,sr->intcorr);
675: if (zeros) SETERRQ1(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",newShift);
676: sPres->value = newShift;
677: }
678: sr->sPres->inertia = iner;
679: pep->target = sr->sPres->value;
680: pep->reason = PEP_CONVERGED_ITERATING;
681: pep->its = 0;
682: } else sr->sPres = NULL;
683: return(0);
684: }
686: /*
687: Obtains value of subsequent shift
688: */
689: static PetscErrorCode PEPGetNewShiftValue(PEP pep,PetscInt side,PetscReal *newS)
690: {
691: PetscReal lambda,d_prev;
692: PetscInt i,idxP;
693: PEP_SR sr;
694: PEP_shift sPres,s;
695: PEP_STOAR *ctx=(PEP_STOAR*)pep->data;
698: sr = ctx->sr;
699: sPres = sr->sPres;
700: if (sPres->neighb[side]) {
701: /* Completing a previous interval */
702: if (!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0) { /* One of the ends might be too far from eigenvalues */
703: if (side) *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[sr->indexEig-1]]))/2;
704: else *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[0]]))/2;
705: } else *newS=(sPres->value + sPres->neighb[side]->value)/2;
706: } else { /* (Only for side=1). Creating a new interval. */
707: if (sPres->neigs==0) {/* No value has been accepted*/
708: if (sPres->neighb[0]) {
709: /* Multiplying by 10 the previous distance */
710: *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
711: sr->nleap++;
712: /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
713: if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Unable to compute the wanted eigenvalues with open interval");
714: } else { /* First shift */
715: if (pep->nconv != 0) {
716: /* Unaccepted values give information for next shift */
717: idxP=0;/* Number of values left from shift */
718: for (i=0;i<pep->nconv;i++) {
719: lambda = PetscRealPart(pep->eigr[i]);
720: if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
721: else break;
722: }
723: /* Avoiding subtraction of eigenvalues (might be the same).*/
724: if (idxP>0) {
725: d_prev = PetscAbsReal(sPres->value - PetscRealPart(pep->eigr[0]))/(idxP+0.3);
726: } else {
727: d_prev = PetscAbsReal(sPres->value - PetscRealPart(pep->eigr[pep->nconv-1]))/(pep->nconv+0.3);
728: }
729: *newS = sPres->value + ((sr->dir)*d_prev*pep->nev)/2;
730: sr->dirch = PETSC_FALSE;
731: } else { /* No values found, no information for next shift */
732: if (!sr->dirch) {
733: sr->dirch = PETSC_TRUE;
734: *newS = sr->int1;
735: } else SETERRQ(PetscObjectComm((PetscObject)pep),1,"First shift renders no information");
736: }
737: }
738: } else { /* Accepted values found */
739: sr->dirch = PETSC_FALSE;
740: sr->nleap = 0;
741: /* Average distance of values in previous subinterval */
742: s = sPres->neighb[0];
743: while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
744: s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
745: }
746: if (s) {
747: d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
748: } else { /* First shift. Average distance obtained with values in this shift */
749: /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
750: 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)) {
751: d_prev = PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
752: } else {
753: d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
754: }
755: }
756: /* Average distance is used for next shift by adding it to value on the right or to shift */
757: if ((sr->dir)*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
758: *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(pep->nev))/2;
759: } else { /* Last accepted value is on the left of shift. Adding to shift */
760: *newS = sPres->value + ((sr->dir)*d_prev*(pep->nev))/2;
761: }
762: }
763: /* End of interval can not be surpassed */
764: if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
765: }/* of neighb[side]==null */
766: return(0);
767: }
769: /*
770: Function for sorting an array of real values
771: */
772: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
773: {
774: PetscReal re;
775: PetscInt i,j,tmp;
778: if (!prev) for (i=0;i<nr;i++) perm[i] = i;
779: /* Insertion sort */
780: for (i=1;i<nr;i++) {
781: re = PetscRealPart(r[perm[i]]);
782: j = i-1;
783: while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
784: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
785: }
786: }
787: return(0);
788: }
790: /* Stores the pairs obtained since the last shift in the global arrays */
791: static PetscErrorCode PEPStoreEigenpairs(PEP pep)
792: {
794: PEP_STOAR *ctx=(PEP_STOAR*)pep->data;
795: PetscReal lambda,err,*errest;
796: PetscInt i,*aux,count=0,ndef,ld,nconv=pep->nconv,d=pep->nmat-1,idx;
797: PetscBool iscayley,divide=PETSC_FALSE;
798: PEP_SR sr = ctx->sr;
799: PEP_shift sPres;
800: Vec w,vomega;
801: Mat MS;
802: BV tV;
803: PetscScalar *S,*eigr,*tS,*omega;
806: sPres = sr->sPres;
807: sPres->index = sr->indexEig;
809: if (nconv>sr->ndef0+sr->ndef1) {
810: /* Back-transform */
811: STBackTransform(pep->st,nconv,pep->eigr,pep->eigi);
812: for (i=0;i<nconv;i++) {
813: #if defined(PETSC_USE_COMPLEX)
814: if (PetscImaginaryPart(pep->eigr[i])) pep->eigr[i] = sr->int0-sr->dir;
815: #else
816: if (pep->eigi[i]) pep->eigr[i] = sr->int0-sr->dir;
817: #endif
818: }
819: PetscObjectTypeCompare((PetscObject)pep->st,STCAYLEY,&iscayley);
820: /* Sort eigenvalues */
821: sortRealEigenvalues(pep->eigr,pep->perm,nconv,PETSC_FALSE,sr->dir);
822: VecCreateSeq(PETSC_COMM_SELF,nconv,&vomega);
823: BVGetSignature(ctx->V,vomega);
824: VecGetArray(vomega,&omega);
825: BVGetSizes(pep->V,NULL,NULL,&ld);
826: BVTensorGetFactors(ctx->V,NULL,&MS);
827: MatDenseGetArray(MS,&S);
828: /* Values stored in global array */
829: PetscCalloc4(nconv,&eigr,nconv,&errest,nconv*nconv*d,&tS,nconv,&aux);
830: ndef = sr->ndef0+sr->ndef1;
831: for (i=0;i<nconv;i++) {
832: lambda = PetscRealPart(pep->eigr[pep->perm[i]]);
833: err = pep->errest[pep->perm[i]];
834: if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
835: if (sr->indexEig+count-ndef>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Unexpected error in Spectrum Slicing");
836: PEPQSliceCheckEigenvalueType(pep,lambda,PetscRealPart(omega[pep->perm[i]]),PETSC_FALSE);
837: eigr[count] = lambda;
838: errest[count] = err;
839: if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) sPres->nconv[0]++;
840: if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) sPres->nconv[1]++;
841: PetscArraycpy(tS+count*(d*nconv),S+pep->perm[i]*(d*ld),nconv);
842: PetscArraycpy(tS+count*(d*nconv)+nconv,S+pep->perm[i]*(d*ld)+ld,nconv);
843: count++;
844: }
845: }
846: VecRestoreArray(vomega,&omega);
847: VecDestroy(&vomega);
848: for (i=0;i<count;i++) {
849: PetscArraycpy(S+i*(d*ld),tS+i*nconv*d,nconv);
850: PetscArraycpy(S+i*(d*ld)+ld,tS+i*nconv*d+nconv,nconv);
851: }
852: MatDenseRestoreArray(MS,&S);
853: BVTensorRestoreFactors(ctx->V,NULL,&MS);
854: BVSetActiveColumns(ctx->V,0,count);
855: BVTensorCompress(ctx->V,count);
856: if (sr->sPres->nconv[0] && sr->sPres->nconv[1]) {
857: divide = PETSC_TRUE;
858: BVTensorGetFactors(ctx->V,NULL,&MS);
859: MatDenseGetArray(MS,&S);
860: PetscArrayzero(tS,nconv*nconv*d);
861: for (i=0;i<count;i++) {
862: PetscArraycpy(tS+i*nconv*d,S+i*(d*ld),count);
863: PetscArraycpy(tS+i*nconv*d+nconv,S+i*(d*ld)+ld,count);
864: }
865: MatDenseRestoreArray(MS,&S);
866: BVTensorRestoreFactors(ctx->V,NULL,&MS);
867: BVSetActiveColumns(pep->V,0,count);
868: BVDuplicateResize(pep->V,count,&tV);
869: BVCopy(pep->V,tV);
870: }
871: if (sr->sPres->nconv[0]) {
872: if (divide) {
873: BVSetActiveColumns(ctx->V,0,sr->sPres->nconv[0]);
874: BVTensorCompress(ctx->V,sr->sPres->nconv[0]);
875: }
876: for (i=0;i<sr->ndef0;i++) aux[i] = sr->idxDef0[i];
877: for (i=sr->ndef0;i<sr->sPres->nconv[0];i++) aux[i] = sr->indexEig+i-sr->ndef0;
878: BVTensorGetFactors(ctx->V,NULL,&MS);
879: MatDenseGetArray(MS,&S);
880: for (i=0;i<sr->sPres->nconv[0];i++) {
881: sr->eigr[aux[i]] = eigr[i];
882: sr->errest[aux[i]] = errest[i];
883: BVGetColumn(pep->V,i,&w);
884: BVInsertVec(sr->V,aux[i],w);
885: BVRestoreColumn(pep->V,i,&w);
886: idx = sr->ld*d*aux[i];
887: PetscArrayzero(sr->S+idx,sr->ld*d);
888: PetscArraycpy(sr->S+idx,S+i*(ld*d),sr->sPres->nconv[0]);
889: PetscArraycpy(sr->S+idx+sr->ld,S+i*(ld*d)+ld,sr->sPres->nconv[0]);
890: PetscFree(sr->qinfo[aux[i]].q);
891: PetscMalloc1(sr->sPres->nconv[0],&sr->qinfo[aux[i]].q);
892: PetscArraycpy(sr->qinfo[aux[i]].q,aux,sr->sPres->nconv[0]);
893: sr->qinfo[aux[i]].nq = sr->sPres->nconv[0];
894: }
895: MatDenseRestoreArray(MS,&S);
896: BVTensorRestoreFactors(ctx->V,NULL,&MS);
897: }
899: if (sr->sPres->nconv[1]) {
900: if (divide) {
901: BVTensorGetFactors(ctx->V,NULL,&MS);
902: MatDenseGetArray(MS,&S);
903: for (i=0;i<sr->sPres->nconv[1];i++) {
904: PetscArraycpy(S+i*(d*ld),tS+(sr->sPres->nconv[0]+i)*nconv*d,count);
905: PetscArraycpy(S+i*(d*ld)+ld,tS+(sr->sPres->nconv[0]+i)*nconv*d+nconv,count);
906: }
907: MatDenseRestoreArray(MS,&S);
908: BVTensorRestoreFactors(ctx->V,NULL,&MS);
909: BVSetActiveColumns(pep->V,0,count);
910: BVCopy(tV,pep->V);
911: BVSetActiveColumns(ctx->V,0,sr->sPres->nconv[1]);
912: BVTensorCompress(ctx->V,sr->sPres->nconv[1]);
913: }
914: for (i=0;i<sr->ndef1;i++) aux[i] = sr->idxDef1[i];
915: for (i=sr->ndef1;i<sr->sPres->nconv[1];i++) aux[i] = sr->indexEig+sr->sPres->nconv[0]-sr->ndef0+i-sr->ndef1;
916: BVTensorGetFactors(ctx->V,NULL,&MS);
917: MatDenseGetArray(MS,&S);
918: for (i=0;i<sr->sPres->nconv[1];i++) {
919: sr->eigr[aux[i]] = eigr[sr->sPres->nconv[0]+i];
920: sr->errest[aux[i]] = errest[sr->sPres->nconv[0]+i];
921: BVGetColumn(pep->V,i,&w);
922: BVInsertVec(sr->V,aux[i],w);
923: BVRestoreColumn(pep->V,i,&w);
924: idx = sr->ld*d*aux[i];
925: PetscArrayzero(sr->S+idx,sr->ld*d);
926: PetscArraycpy(sr->S+idx,S+i*(ld*d),sr->sPres->nconv[1]);
927: PetscArraycpy(sr->S+idx+sr->ld,S+i*(ld*d)+ld,sr->sPres->nconv[1]);
928: PetscFree(sr->qinfo[aux[i]].q);
929: PetscMalloc1(sr->sPres->nconv[1],&sr->qinfo[aux[i]].q);
930: PetscArraycpy(sr->qinfo[aux[i]].q,aux,sr->sPres->nconv[1]);
931: sr->qinfo[aux[i]].nq = sr->sPres->nconv[1];
932: }
933: MatDenseRestoreArray(MS,&S);
934: BVTensorRestoreFactors(ctx->V,NULL,&MS);
935: }
936: sPres->neigs = count-sr->ndef0-sr->ndef1;
937: sr->indexEig += sPres->neigs;
938: sPres->nconv[0]-= sr->ndef0;
939: sPres->nconv[1]-= sr->ndef1;
940: PetscFree4(eigr,errest,tS,aux);
941: } else {
942: sPres->neigs = 0;
943: sPres->nconv[0]= 0;
944: sPres->nconv[1]= 0;
945: }
946: /* Global ordering array updating */
947: sortRealEigenvalues(sr->eigr,sr->perm,sr->indexEig,PETSC_FALSE,sr->dir);
948: /* Check for completion */
949: sPres->comp[0] = PetscNot(sPres->nconv[0] < sPres->nsch[0]);
950: sPres->comp[1] = PetscNot(sPres->nconv[1] < sPres->nsch[1]);
951: 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");
952: if (divide) { BVDestroy(&tV); }
953: return(0);
954: }
956: static PetscErrorCode PEPLookForDeflation(PEP pep)
957: {
958: PetscReal val;
959: PetscInt i,count0=0,count1=0;
960: PEP_shift sPres;
961: PetscInt ini,fin;
962: PEP_SR sr;
963: PEP_STOAR *ctx=(PEP_STOAR*)pep->data;
966: sr = ctx->sr;
967: sPres = sr->sPres;
969: if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
970: else ini = 0;
971: fin = sr->indexEig;
972: /* Selection of ends for searching new values */
973: if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
974: else sPres->ext[0] = sPres->neighb[0]->value;
975: if (!sPres->neighb[1]) {
976: if (sr->hasEnd) sPres->ext[1] = sr->int1;
977: else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
978: } else sPres->ext[1] = sPres->neighb[1]->value;
979: /* Selection of values between right and left ends */
980: for (i=ini;i<fin;i++) {
981: val=PetscRealPart(sr->eigr[sr->perm[i]]);
982: /* Values to the right of left shift */
983: if ((sr->dir)*(val - sPres->ext[1]) < 0) {
984: if ((sr->dir)*(val - sPres->value) < 0) count0++;
985: else count1++;
986: } else break;
987: }
988: /* The number of values on each side are found */
989: if (sPres->neighb[0]) {
990: sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
991: if (sPres->nsch[0]<0) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Mismatch between number of values found and information from inertia");
992: } else sPres->nsch[0] = 0;
994: if (sPres->neighb[1]) {
995: sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
996: if (sPres->nsch[1]<0) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Mismatch between number of values found and information from inertia");
997: } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);
999: /* Completing vector of indexes for deflation */
1000: for (i=0;i<count0;i++) sr->idxDef0[i] = sr->perm[ini+i];
1001: sr->ndef0 = count0;
1002: for (i=0;i<count1;i++) sr->idxDef1[i] = sr->perm[ini+count0+i];
1003: sr->ndef1 = count1;
1004: return(0);
1005: }
1007: /*
1008: Compute a run of Lanczos iterations
1009: */
1010: static PetscErrorCode PEPSTOARrun_QSlice(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,Vec *t_)
1011: {
1013: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
1014: PetscInt i,j,m=*M,l,lock;
1015: PetscInt lds,d,ld,offq,nqt,ldds;
1016: Vec v=t_[0],t=t_[1],q=t_[2];
1017: PetscReal norm,sym=0.0,fro=0.0,*f;
1018: PetscScalar *y,*S,sigma;
1019: PetscBLASInt j_,one=1;
1020: PetscBool lindep;
1021: Mat MS;
1024: PetscMalloc1(*M,&y);
1025: BVGetSizes(pep->V,NULL,NULL,&ld);
1026: BVTensorGetDegree(ctx->V,&d);
1027: BVGetActiveColumns(pep->V,&lock,&nqt);
1028: lds = d*ld;
1029: offq = ld;
1030: DSGetLeadingDimension(pep->ds,&ldds);
1032: *breakdown = PETSC_FALSE; /* ----- */
1033: STGetShift(pep->st,&sigma);
1034: DSGetDimensions(pep->ds,NULL,NULL,&l,NULL,NULL);
1035: BVSetActiveColumns(ctx->V,0,m);
1036: BVSetActiveColumns(pep->V,0,nqt);
1037: for (j=k;j<m;j++) {
1038: /* apply operator */
1039: BVTensorGetFactors(ctx->V,NULL,&MS);
1040: MatDenseGetArray(MS,&S);
1041: BVGetColumn(pep->V,nqt,&t);
1042: BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
1043: MatMult(pep->A[1],v,q);
1044: MatMult(pep->A[2],v,t);
1045: VecAXPY(q,sigma*pep->sfactor,t);
1046: VecScale(q,pep->sfactor);
1047: BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
1048: MatMult(pep->A[2],v,t);
1049: VecAXPY(q,pep->sfactor*pep->sfactor,t);
1050: STMatSolve(pep->st,q,t);
1051: VecScale(t,-1.0);
1052: BVRestoreColumn(pep->V,nqt,&t);
1054: /* orthogonalize */
1055: BVOrthogonalizeColumn(pep->V,nqt,S+(j+1)*lds,&norm,&lindep);
1056: if (!lindep) {
1057: *(S+(j+1)*lds+nqt) = norm;
1058: BVScaleColumn(pep->V,nqt,1.0/norm);
1059: nqt++;
1060: }
1061: for (i=0;i<nqt;i++) *(S+(j+1)*lds+offq+i) = *(S+j*lds+i)+sigma*(*(S+(j+1)*lds+i));
1062: BVSetActiveColumns(pep->V,0,nqt);
1063: MatDenseRestoreArray(MS,&S);
1064: BVTensorRestoreFactors(ctx->V,NULL,&MS);
1066: /* level-2 orthogonalization */
1067: BVOrthogonalizeColumn(ctx->V,j+1,y,&norm,&lindep);
1068: a[j] = PetscRealPart(y[j]);
1069: omega[j+1] = (norm > 0)?1.0:-1.0;
1070: BVScaleColumn(ctx->V,j+1,1.0/norm);
1071: b[j] = PetscAbsReal(norm);
1073: /* check symmetry */
1074: DSGetArrayReal(pep->ds,DS_MAT_T,&f);
1075: if (j==k) {
1076: for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ldds+i]);
1077: for (i=0;i<l;i++) y[i] = 0.0;
1078: }
1079: DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
1080: if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
1081: PetscBLASIntCast(j,&j_);
1082: sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
1083: fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
1084: if (j>0) fro = SlepcAbs(fro,b[j-1]);
1085: if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
1086: *symmlost = PETSC_TRUE;
1087: *M=j;
1088: break;
1089: }
1090: }
1091: BVSetActiveColumns(pep->V,lock,nqt);
1092: BVSetActiveColumns(ctx->V,0,*M);
1093: PetscFree(y);
1094: return(0);
1095: }
1097: static PetscErrorCode PEPSTOAR_QSlice(PEP pep,Mat B)
1098: {
1100: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
1101: PetscInt j,k,l,nv=0,ld,ldds,t,nq=0,idx;
1102: PetscInt nconv=0,deg=pep->nmat-1,count0=0,count1=0;
1103: PetscScalar *Q,*om,sigma,*back,*S,*pQ;
1104: PetscReal beta,norm=1.0,*omega,*a,*b,*r,eta,lambda;
1105: PetscBool breakdown,symmlost=PETSC_FALSE,sinv,falselock=PETSC_TRUE;
1106: Mat MS,MQ;
1107: Vec v,vomega;
1108: PEP_SR sr;
1109: BVOrthogType otype;
1110: BVOrthogBlockType obtype;
1113: /* Resize if needed for deflating vectors */
1114: sr = ctx->sr;
1115: sigma = sr->sPres->value;
1116: k = sr->ndef0+sr->ndef1;
1117: pep->ncv = ctx->ncv+k;
1118: pep->nev = ctx->nev+k;
1119: PEPAllocateSolution(pep,3);
1120: BVDestroy(&ctx->V);
1121: BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
1122: BVGetOrthogonalization(pep->V,&otype,NULL,&eta,&obtype);
1123: BVSetOrthogonalization(ctx->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
1124: DSAllocate(pep->ds,pep->ncv+2);
1125: PetscMalloc1(pep->ncv,&back);
1126: DSGetLeadingDimension(pep->ds,&ldds);
1127: BVSetMatrix(ctx->V,B,PETSC_TRUE);
1128: if (ctx->lock) {
1129: /* undocumented option to use a cheaper locking instead of the true locking */
1130: PetscOptionsGetBool(NULL,NULL,"-pep_stoar_falselocking",&falselock,NULL);
1131: } else SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"A locking variant is needed for spectrum slicing");
1132: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
1133: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
1134: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
1136: /* Get the starting Arnoldi vector */
1137: BVSetActiveColumns(pep->V,0,1);
1138: BVTensorBuildFirstColumn(ctx->V,pep->nini);
1139: BVSetActiveColumns(ctx->V,0,1);
1140: if (k) {
1141: /* Insert deflated vectors */
1142: BVSetActiveColumns(pep->V,0,0);
1143: idx = sr->ndef0?sr->idxDef0[0]:sr->idxDef1[0];
1144: for (j=0;j<k;j++) {
1145: BVGetColumn(pep->V,j,&v);
1146: BVCopyVec(sr->V,sr->qinfo[idx].q[j],v);
1147: BVRestoreColumn(pep->V,j,&v);
1148: }
1149: /* Update innerproduct matrix */
1150: BVSetActiveColumns(ctx->V,0,0);
1151: BVTensorGetFactors(ctx->V,NULL,&MS);
1152: BVSetActiveColumns(pep->V,0,k);
1153: BVTensorRestoreFactors(ctx->V,NULL,&MS);
1155: BVGetSizes(pep->V,NULL,NULL,&ld);
1156: BVTensorGetFactors(ctx->V,NULL,&MS);
1157: MatDenseGetArray(MS,&S);
1158: for (j=0;j<sr->ndef0;j++) {
1159: PetscArrayzero(S+j*ld*deg,ld*deg);
1160: PetscArraycpy(S+j*ld*deg,sr->S+sr->idxDef0[j]*sr->ld*deg,k);
1161: PetscArraycpy(S+j*ld*deg+ld,sr->S+sr->idxDef0[j]*sr->ld*deg+sr->ld,k);
1162: pep->eigr[j] = sr->eigr[sr->idxDef0[j]];
1163: pep->errest[j] = sr->errest[sr->idxDef0[j]];
1164: }
1165: for (j=0;j<sr->ndef1;j++) {
1166: PetscArrayzero(S+(j+sr->ndef0)*ld*deg,ld*deg);
1167: PetscArraycpy(S+(j+sr->ndef0)*ld*deg,sr->S+sr->idxDef1[j]*sr->ld*deg,k);
1168: PetscArraycpy(S+(j+sr->ndef0)*ld*deg+ld,sr->S+sr->idxDef1[j]*sr->ld*deg+sr->ld,k);
1169: pep->eigr[j+sr->ndef0] = sr->eigr[sr->idxDef1[j]];
1170: pep->errest[j+sr->ndef0] = sr->errest[sr->idxDef1[j]];
1171: }
1172: MatDenseRestoreArray(MS,&S);
1173: BVTensorRestoreFactors(ctx->V,NULL,&MS);
1174: BVSetActiveColumns(ctx->V,0,k+1);
1175: VecCreateSeq(PETSC_COMM_SELF,k+1,&vomega);
1176: VecGetArray(vomega,&om);
1177: for (j=0;j<k;j++) {
1178: BVOrthogonalizeColumn(ctx->V,j,NULL,&norm,NULL);
1179: BVScaleColumn(ctx->V,j,1/norm);
1180: om[j] = (norm>=0.0)?1.0:-1.0;
1181: }
1182: BVTensorGetFactors(ctx->V,NULL,&MS);
1183: MatDenseGetArray(MS,&S);
1184: for (j=0;j<deg;j++) {
1185: BVSetRandomColumn(pep->V,k+j);
1186: BVOrthogonalizeColumn(pep->V,k+j,S+k*ld*deg+j*ld,&norm,NULL);
1187: BVScaleColumn(pep->V,k+j,1.0/norm);
1188: S[k*ld*deg+j*ld+k+j] = norm;
1189: }
1190: MatDenseRestoreArray(MS,&S);
1191: BVSetActiveColumns(pep->V,0,k+deg);
1192: BVTensorRestoreFactors(ctx->V,NULL,&MS);
1193: BVOrthogonalizeColumn(ctx->V,k,NULL,&norm,NULL);
1194: BVScaleColumn(ctx->V,k,1.0/norm);
1195: om[k] = (norm>=0.0)?1.0:-1.0;
1196: VecRestoreArray(vomega,&om);
1197: BVSetSignature(ctx->V,vomega);
1198: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
1199: VecGetArray(vomega,&om);
1200: for (j=0;j<k;j++) a[j] = PetscRealPart(om[j]/(pep->eigr[j]-sigma));
1201: VecRestoreArray(vomega,&om);
1202: VecDestroy(&vomega);
1203: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
1204: DSGetArray(pep->ds,DS_MAT_Q,&pQ);
1205: PetscArrayzero(pQ,ldds*k);
1206: for (j=0;j<k;j++) pQ[j+j*ldds] = 1.0;
1207: DSRestoreArray(pep->ds,DS_MAT_Q,&pQ);
1208: }
1209: BVSetActiveColumns(ctx->V,0,k+1);
1210: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
1211: VecCreateSeq(PETSC_COMM_SELF,k+1,&vomega);
1212: BVGetSignature(ctx->V,vomega);
1213: VecGetArray(vomega,&om);
1214: for (j=0;j<k+1;j++) omega[j] = PetscRealPart(om[j]);
1215: VecRestoreArray(vomega,&om);
1216: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
1217: VecDestroy(&vomega);
1219: 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]);
1221: /* Restart loop */
1222: l = 0;
1223: pep->nconv = k;
1224: while (pep->reason == PEP_CONVERGED_ITERATING) {
1225: pep->its++;
1226: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
1227: b = a+ldds;
1228: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
1230: /* Compute an nv-step Lanczos factorization */
1231: nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
1232: PEPSTOARrun_QSlice(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,pep->work);
1233: beta = b[nv-1];
1234: if (symmlost && nv==pep->nconv+l) {
1235: pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
1236: pep->nconv = nconv;
1237: PetscInfo2(pep,"Symmetry lost in STOAR sigma=%g nconv=%D\n",(double)sr->sPres->value,nconv);
1238: if (falselock || !ctx->lock) {
1239: BVSetActiveColumns(ctx->V,0,pep->nconv);
1240: BVTensorCompress(ctx->V,0);
1241: }
1242: break;
1243: }
1244: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
1245: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
1246: DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
1247: if (l==0) {
1248: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
1249: } else {
1250: DSSetState(pep->ds,DS_STATE_RAW);
1251: }
1253: /* Solve projected problem */
1254: DSSolve(pep->ds,pep->eigr,pep->eigi);
1255: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
1256: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
1258: /* Check convergence */
1259: /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
1260: norm = 1.0;
1261: DSGetDimensions(pep->ds,NULL,NULL,NULL,NULL,&t);
1262: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
1263: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
1264: for (j=0;j<k;j++) back[j] = pep->eigr[j];
1265: STBackTransform(pep->st,k,back,pep->eigi);
1266: count0=count1=0;
1267: for (j=0;j<k;j++) {
1268: lambda = PetscRealPart(back[j]);
1269: if (((sr->dir)*(sr->sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sr->sPres->ext[0]) > 0)) count0++;
1270: if (((sr->dir)*(lambda - sr->sPres->value) > 0) && ((sr->dir)*(sr->sPres->ext[1] - lambda) > 0)) count1++;
1271: }
1272: if ((count0-sr->ndef0 >= sr->sPres->nsch[0]) && (count1-sr->ndef1 >= sr->sPres->nsch[1])) pep->reason = PEP_CONVERGED_TOL;
1273: /* Update l */
1274: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
1275: else {
1276: l = PetscMax(1,(PetscInt)((nv-k)/2));
1277: l = PetscMin(l,t);
1278: if (!breakdown) {
1279: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
1280: if (*(a+ldds+k+l-1)!=0) {
1281: if (k+l<nv-1) l = l+1;
1282: else l = l-1;
1283: }
1284: /* Prepare the Rayleigh quotient for restart */
1285: DSGetArray(pep->ds,DS_MAT_Q,&Q);
1286: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
1287: r = a + 2*ldds;
1288: for (j=k;j<k+l;j++) {
1289: r[j] = PetscRealPart(Q[nv-1+j*ldds]*beta);
1290: }
1291: b = a+ldds;
1292: b[k+l-1] = r[k+l-1];
1293: omega[k+l] = omega[nv];
1294: DSRestoreArray(pep->ds,DS_MAT_Q,&Q);
1295: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
1296: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
1297: }
1298: }
1299: nconv = k;
1300: if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
1302: /* Update S */
1303: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
1304: BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
1305: MatDestroy(&MQ);
1307: /* Copy last column of S */
1308: BVCopyColumn(ctx->V,nv,k+l);
1309: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
1310: VecCreateSeq(PETSC_COMM_SELF,k+l,&vomega);
1311: VecGetArray(vomega,&om);
1312: for (j=0;j<k+l;j++) om[j] = omega[j];
1313: VecRestoreArray(vomega,&om);
1314: BVSetActiveColumns(ctx->V,0,k+l);
1315: BVSetSignature(ctx->V,vomega);
1316: VecDestroy(&vomega);
1317: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
1319: if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
1320: /* stop if breakdown */
1321: PetscInfo2(pep,"Breakdown TOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
1322: pep->reason = PEP_DIVERGED_BREAKDOWN;
1323: }
1324: if (pep->reason != PEP_CONVERGED_ITERATING) l--;
1325: BVGetActiveColumns(pep->V,NULL,&nq);
1326: if (k+l+deg<=nq) {
1327: BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
1328: if (!falselock && ctx->lock) {
1329: BVTensorCompress(ctx->V,k-pep->nconv);
1330: } else {
1331: BVTensorCompress(ctx->V,0);
1332: }
1333: }
1334: pep->nconv = k;
1335: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
1336: }
1337: sr->itsKs += pep->its;
1338: if (pep->nconv>0) {
1339: BVSetActiveColumns(ctx->V,0,pep->nconv);
1340: BVGetActiveColumns(pep->V,NULL,&nq);
1341: BVSetActiveColumns(pep->V,0,nq);
1342: if (nq>pep->nconv) {
1343: BVTensorCompress(ctx->V,pep->nconv);
1344: BVSetActiveColumns(pep->V,0,pep->nconv);
1345: }
1346: for (j=0;j<pep->nconv;j++) {
1347: pep->eigr[j] *= pep->sfactor;
1348: pep->eigi[j] *= pep->sfactor;
1349: }
1350: }
1351: 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);
1352: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
1353: RGPopScale(pep->rg);
1355: 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);
1356: if (pep->reason == PEP_DIVERGED_SYMMETRY_LOST && nconv==sr->ndef0+sr->ndef1) {
1357: if (++sr->symmlost>10) SETERRQ1(PetscObjectComm((PetscObject)pep),1,"Symmetry lost at sigma=%g",(double)sr->sPres->value);
1358: } else sr->symmlost = 0;
1360: /* truncate Schur decomposition and change the state to raw so that
1361: DSVectors() computes eigenvectors from scratch */
1362: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
1363: DSSetState(pep->ds,DS_STATE_RAW);
1364: PetscFree(back);
1365: return(0);
1366: }
1368: #define SWAP(a,b,t) {t=a;a=b;b=t;}
1370: static PetscErrorCode PEPQSliceGetInertias(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
1371: {
1372: PetscErrorCode ierr;
1373: PEP_STOAR *ctx=(PEP_STOAR*)pep->data;
1374: PEP_SR sr=ctx->sr;
1375: PetscInt i=0,j,tmpi;
1376: PetscReal v,tmpr;
1377: PEP_shift s;
1380: if (!pep->state) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSetUp() first");
1381: if (!sr) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see PEPSetInterval()");
1382: if (!sr->s0) { /* PEPSolve not called yet */
1383: *n = 2;
1384: } else {
1385: *n = 1;
1386: s = sr->s0;
1387: while (s) {
1388: (*n)++;
1389: s = s->neighb[1];
1390: }
1391: }
1392: PetscMalloc1(*n,shifts);
1393: PetscMalloc1(*n,inertias);
1394: if (!sr->s0) { /* PEPSolve not called yet */
1395: (*shifts)[0] = sr->int0;
1396: (*shifts)[1] = sr->int1;
1397: (*inertias)[0] = sr->inertia0;
1398: (*inertias)[1] = sr->inertia1;
1399: } else {
1400: s = sr->s0;
1401: while (s) {
1402: (*shifts)[i] = s->value;
1403: (*inertias)[i++] = s->inertia;
1404: s = s->neighb[1];
1405: }
1406: (*shifts)[i] = sr->int1;
1407: (*inertias)[i] = sr->inertia1;
1408: }
1409: /* remove possible duplicate in last position */
1410: if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
1411: /* sort result */
1412: for (i=0;i<*n;i++) {
1413: v = (*shifts)[i];
1414: for (j=i+1;j<*n;j++) {
1415: if (v > (*shifts)[j]) {
1416: SWAP((*shifts)[i],(*shifts)[j],tmpr);
1417: SWAP((*inertias)[i],(*inertias)[j],tmpi);
1418: v = (*shifts)[i];
1419: }
1420: }
1421: }
1422: return(0);
1423: }
1425: PetscErrorCode PEPSolve_STOAR_QSlice(PEP pep)
1426: {
1428: PetscInt i,j,ti,deg=pep->nmat-1;
1429: PetscReal newS;
1430: PEP_STOAR *ctx=(PEP_STOAR*)pep->data;
1431: PEP_SR sr=ctx->sr;
1432: Mat S,B;
1433: PetscScalar *pS;
1436: PetscCitationsRegister(citation,&cited);
1438: /* Only with eigenvalues present in the interval ...*/
1439: if (sr->numEigs==0) {
1440: pep->reason = PEP_CONVERGED_TOL;
1441: return(0);
1442: }
1444: /* Inner product matrix */
1445: PEPSTOARSetUpInnerMatrix(pep,&B);
1447: /* Array of pending shifts */
1448: sr->maxPend = 100; /* Initial size */
1449: sr->nPend = 0;
1450: PetscMalloc1(sr->maxPend,&sr->pending);
1451: PetscLogObjectMemory((PetscObject)pep,(sr->maxPend)*sizeof(PEP_shift));
1452: PEPCreateShift(pep,sr->int0,NULL,NULL);
1453: /* extract first shift */
1454: sr->sPrev = NULL;
1455: sr->sPres = sr->pending[--sr->nPend];
1456: sr->sPres->inertia = sr->inertia0;
1457: pep->target = sr->sPres->value;
1458: sr->s0 = sr->sPres;
1459: sr->indexEig = 0;
1461: /* Memory reservation for auxiliary variables */
1462: PetscLogObjectMemory((PetscObject)pep,(sr->numEigs+2*pep->ncv)*sizeof(PetscScalar));
1463: for (i=0;i<sr->numEigs;i++) {
1464: sr->eigr[i] = 0.0;
1465: sr->eigi[i] = 0.0;
1466: sr->errest[i] = 0.0;
1467: sr->perm[i] = i;
1468: }
1469: /* Vectors for deflation */
1470: PetscMalloc2(sr->numEigs,&sr->idxDef0,sr->numEigs,&sr->idxDef1);
1471: PetscLogObjectMemory((PetscObject)pep,sr->numEigs*sizeof(PetscInt));
1472: sr->indexEig = 0;
1473: while (sr->sPres) {
1474: /* Search for deflation */
1475: PEPLookForDeflation(pep);
1476: /* KrylovSchur */
1477: PEPSTOAR_QSlice(pep,B);
1479: PEPStoreEigenpairs(pep);
1480: /* Select new shift */
1481: if (!sr->sPres->comp[1]) {
1482: PEPGetNewShiftValue(pep,1,&newS);
1483: PEPCreateShift(pep,newS,sr->sPres,sr->sPres->neighb[1]);
1484: }
1485: if (!sr->sPres->comp[0]) {
1486: /* Completing earlier interval */
1487: PEPGetNewShiftValue(pep,0,&newS);
1488: PEPCreateShift(pep,newS,sr->sPres->neighb[0],sr->sPres);
1489: }
1490: /* Preparing for a new search of values */
1491: PEPExtractShift(pep);
1492: }
1494: /* Updating pep values prior to exit */
1495: PetscFree2(sr->idxDef0,sr->idxDef1);
1496: PetscFree(sr->pending);
1497: pep->nconv = sr->indexEig;
1498: pep->reason = PEP_CONVERGED_TOL;
1499: pep->its = sr->itsKs;
1500: pep->nev = sr->indexEig;
1501: MatCreateSeqDense(PETSC_COMM_SELF,pep->nconv,pep->nconv,NULL,&S);
1502: MatDenseGetArray(S,&pS);
1503: for (i=0;i<pep->nconv;i++) {
1504: 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);
1505: }
1506: MatDenseRestoreArray(S,&pS);
1507: BVSetActiveColumns(sr->V,0,pep->nconv);
1508: BVMultInPlace(sr->V,S,0,pep->nconv);
1509: MatDestroy(&S);
1510: BVDestroy(&pep->V);
1511: pep->V = sr->V;
1512: PetscFree4(pep->eigr,pep->eigi,pep->errest,pep->perm);
1513: pep->eigr = sr->eigr;
1514: pep->eigi = sr->eigi;
1515: pep->perm = sr->perm;
1516: pep->errest = sr->errest;
1517: if (sr->dir<0) {
1518: for (i=0;i<pep->nconv/2;i++) {
1519: ti = sr->perm[i]; sr->perm[i] = sr->perm[pep->nconv-1-i]; sr->perm[pep->nconv-1-i] = ti;
1520: }
1521: }
1522: PetscFree(ctx->inertias);
1523: PetscFree(ctx->shifts);
1524: MatDestroy(&B);
1525: PEPQSliceGetInertias(pep,&ctx->nshifts,&ctx->shifts,&ctx->inertias);
1526: return(0);
1527: }