Actual source code: qslice.c
slepc-3.9.0 2018-04-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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", in preparation,
23: 2018.
24: */
26: #include <slepc/private/pepimpl.h> /*I "slepcpep.h" I*/
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 = \"In preparation\",\n"
36: " volume = \"xx\",\n"
37: " number = \"x\",\n"
38: " pages = \"xx--xx\",\n"
39: " year = \"2018,\"\n"
40: " doi = \"https://doi.org/10.1007/xxx\"\n"
41: "}\n";
43: #define SLICE_PTOL PETSC_SQRT_MACHINE_EPSILON
45: static PetscErrorCode PEPQSliceResetSR(PEP pep)
46: {
48: PEP_TOAR *ctx=(PEP_TOAR*)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_TOAR *ctx=(PEP_TOAR*)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_TOAR *ctx=(PEP_TOAR*)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 PEPQSliceGetInertia(PEP pep,PetscReal shift,PetscInt *inertia,PetscInt *zeros,PetscInt correction)
140: {
142: KSP ksp;
143: PC pc;
144: Mat F,P;
145: PetscReal nzshift=0.0;
146: PetscScalar dot;
147: PetscRandom rand;
148: PetscInt nconv;
149: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
150: PEP_SR sr=ctx->sr;
153: if (shift >= PETSC_MAX_REAL) { /* Right-open interval */
154: *inertia = 0;
155: } else if (shift <= PETSC_MIN_REAL) {
156: *inertia = 0;
157: if (zeros) *zeros = 0;
158: } else {
159: /* If the shift is zero, perturb it to a very small positive value.
160: The goal is that the nonzero pattern is the same in all cases and reuse
161: the symbolic factorizations */
162: nzshift = (shift==0.0)? 10.0/PETSC_MAX_REAL: shift;
163: STSetShift(pep->st,nzshift);
164: PEPEvaluateBasis(pep,nzshift,0,pep->solvematcoeffs,NULL);
165: STSetUp(pep->st);
166: STMatSetUp(pep->st,pep->sfactor,pep->solvematcoeffs);
167: STGetKSP(pep->st,&ksp);
168: KSPGetPC(ksp,&pc);
169: PCFactorGetMatrix(pc,&F);
170: MatGetInertia(F,inertia,zeros,NULL);
171: }
172: if (!correction) {
173: if (shift >= PETSC_MAX_REAL) *inertia = 2*pep->n;
174: else if (shift>PETSC_MIN_REAL) {
175: KSPGetOperators(ksp,&P,NULL);
176: if (*inertia!=pep->n && !sr->v[0]) {
177: MatCreateVecs(P,&sr->v[0],NULL);
178: VecDuplicate(sr->v[0],&sr->v[1]);
179: VecDuplicate(sr->v[0],&sr->v[2]);
180: BVGetRandomContext(pep->V,&rand);
181: VecSetRandom(sr->v[0],rand);
182: }
183: if (*inertia<pep->n && *inertia>0) {
184: if (!sr->eps) {
185: EPSCreate(PetscObjectComm((PetscObject)pep),&sr->eps);
186: EPSSetProblemType(sr->eps,EPS_HEP);
187: EPSSetWhichEigenpairs(sr->eps,EPS_LARGEST_REAL);
188: EPSSetConvergenceTestFunction(sr->eps,ConvergedPositive,NULL,NULL);
189: }
190: EPSSetOperators(sr->eps,P,NULL);
191: EPSSolve(sr->eps);
192: EPSGetConverged(sr->eps,&nconv);
193: if (!nconv) SETERRQ1(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",nzshift);
194: EPSGetEigenpair(sr->eps,0,NULL,NULL,sr->v[0],sr->v[1]);
195: }
196: if (*inertia!=pep->n) {
197: MatMult(pep->A[1],sr->v[0],sr->v[1]);
198: MatMult(pep->A[2],sr->v[0],sr->v[2]);
199: VecAXPY(sr->v[1],2*nzshift,sr->v[2]);
200: VecDot(sr->v[1],sr->v[0],&dot);
201: if (PetscRealPart(dot)>0.0) *inertia = 2*pep->n-*inertia;
202: }
203: }
204: } else if (correction<0) *inertia = 2*pep->n-*inertia;
205: return(0);
206: }
208: /*
209: Dummy backtransform operation
210: */
211: static PetscErrorCode PEPBackTransform_Skip(PEP pep)
212: {
214: return(0);
215: }
217: PetscErrorCode PEPSetUp_STOAR_QSlice(PEP pep)
218: {
220: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
221: PEP_SR sr;
222: PetscInt ld,i,zeros=0;
223: SlepcSC sc;
224: PetscBool issinv;
225: PetscReal r;
228: 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");
229: PetscObjectTypeCompareAny((PetscObject)pep->st,&issinv,STSINVERT,STCAYLEY,"");
230: if (!issinv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
231: if (pep->tol==PETSC_DEFAULT) pep->tol = SLEPC_DEFAULT_TOL*1e-2; /* use tighter tolerance */
232: if (ctx->nev==0) ctx->nev = PetscMin(20,pep->n); /* nev not set, use default value */
233: 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");
234: pep->ops->backtransform = PEPBackTransform_Skip;
236: /* create spectrum slicing context and initialize it */
237: PEPQSliceResetSR(pep);
238: PetscNewLog(pep,&sr);
239: ctx->sr = sr;
240: sr->itsKs = 0;
241: sr->nleap = 0;
242: sr->sPres = NULL;
243: /* check presence of ends and finding direction */
244: if (pep->inta > PETSC_MIN_REAL || pep->intb >= PETSC_MAX_REAL) {
245: sr->int0 = pep->inta;
246: sr->int1 = pep->intb;
247: sr->dir = 1;
248: if (pep->intb >= PETSC_MAX_REAL) { /* Right-open interval */
249: sr->hasEnd = PETSC_FALSE;
250: } else sr->hasEnd = PETSC_TRUE;
251: } else {
252: sr->int0 = pep->intb;
253: sr->int1 = pep->inta;
254: sr->dir = -1;
255: sr->hasEnd = PetscNot(pep->inta <= PETSC_MIN_REAL);
256: }
258: PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
259: if (!pep->st) {PEPGetST(pep,&pep->st);}
260: STSetTransform(pep->st,PETSC_FALSE);
261: STSetUp(pep->st);
263: /* compute inertia0 */
264: ctx->hyperbolic = (pep->problem_type==PEP_HYPERBOLIC)? PETSC_TRUE: PETSC_FALSE;
265: PEPQSliceGetInertia(pep,sr->int0,&sr->inertia0,ctx->detect?&zeros:NULL,ctx->hyperbolic?0:1);
266: 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");
268: /* compute inertia1 */
269: PEPQSliceGetInertia(pep,sr->int1,&sr->inertia1,ctx->detect?&zeros:NULL,ctx->hyperbolic?0:1);
270: if (zeros) SETERRQ(((PetscObject)pep)->comm,PETSC_ERR_USER,"Found singular matrix for the transformed problem in an interval endpoint defined by user");
271: if (!ctx->hyperbolic) {
272: if (sr->dir*(sr->inertia1-sr->inertia0)<0) {
273: sr->intcorr = -1;
274: sr->inertia0 = 2*pep->n-sr->inertia0;
275: sr->inertia1 = 2*pep->n-sr->inertia1;
276: } else sr->intcorr = 1;
277: } else {
278: if (sr->inertia0<=pep->n && sr->inertia1<=pep->n) sr->intcorr = 1;
279: else if (sr->inertia0>=pep->n && sr->inertia1>=pep->n) sr->intcorr = -1;
280: }
281:
282: if (sr->hasEnd) {
283: sr->dir = -sr->dir; r = sr->int0; sr->int0 = sr->int1; sr->int1 = r;
284: i = sr->inertia0; sr->inertia0 = sr->inertia1; sr->inertia1 = i;
285: }
287: /* number of eigenvalues in interval */
288: sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
289: if (sr->numEigs) {
290: PEPQSliceAllocateSolution(pep);
291: PEPSetDimensions_Default(pep,ctx->nev,&ctx->ncv,&ctx->mpd);
292: pep->nev = ctx->nev; pep->ncv = ctx->ncv; pep->mpd = ctx->mpd;
293: if (!pep->max_it) pep->max_it = PetscMax(100,2*(pep->nmat-1)*pep->n/pep->ncv);
294: ld = ctx->ncv+2;
295: DSSetType(pep->ds,DSGHIEP);
296: DSSetCompact(pep->ds,PETSC_TRUE);
297: DSAllocate(pep->ds,ld);
298: DSGetSlepcSC(pep->ds,&sc);
299: sc->rg = NULL;
300: sc->comparison = SlepcCompareLargestMagnitude;
301: sc->comparisonctx = NULL;
302: sc->map = NULL;
303: sc->mapobj = NULL;
304: }
305: return(0);
306: }
308: /*
309: Fills the fields of a shift structure
310: */
311: static PetscErrorCode PEPCreateShift(PEP pep,PetscReal val,PEP_shift neighb0,PEP_shift neighb1)
312: {
314: PEP_shift s,*pending2;
315: PetscInt i;
316: PEP_SR sr;
317: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
320: sr = ctx->sr;
321: PetscNewLog(pep,&s);
322: s->value = val;
323: s->neighb[0] = neighb0;
324: if (neighb0) neighb0->neighb[1] = s;
325: s->neighb[1] = neighb1;
326: if (neighb1) neighb1->neighb[0] = s;
327: s->comp[0] = PETSC_FALSE;
328: s->comp[1] = PETSC_FALSE;
329: s->index = -1;
330: s->neigs = 0;
331: s->nconv[0] = s->nconv[1] = 0;
332: s->nsch[0] = s->nsch[1]=0;
333: /* Inserts in the stack of pending shifts */
334: /* If needed, the array is resized */
335: if (sr->nPend >= sr->maxPend) {
336: sr->maxPend *= 2;
337: PetscMalloc1(sr->maxPend,&pending2);
338: PetscLogObjectMemory((PetscObject)pep,sizeof(PEP_shift));
339: for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
340: PetscFree(sr->pending);
341: sr->pending = pending2;
342: }
343: sr->pending[sr->nPend++]=s;
344: return(0);
345: }
347: /* Provides next shift to be computed */
348: static PetscErrorCode PEPExtractShift(PEP pep)
349: {
351: PetscInt iner,zeros=0;
352: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
353: PEP_SR sr;
354: PetscReal newShift;
355: PEP_shift sPres;
358: sr = ctx->sr;
359: if (sr->nPend > 0) {
360: sr->sPrev = sr->sPres;
361: sr->sPres = sr->pending[--sr->nPend];
362: sPres = sr->sPres;
363: PEPQSliceGetInertia(pep,sPres->value,&iner,ctx->detect?&zeros:NULL,sr->intcorr);
364: if (zeros) {
365: newShift = sPres->value*(1.0+SLICE_PTOL);
366: if (sr->dir*(sPres->neighb[0] && newShift-sPres->neighb[0]->value) < 0) newShift = (sPres->value+sPres->neighb[0]->value)/2;
367: else if (sPres->neighb[1] && sr->dir*(sPres->neighb[1]->value-newShift) < 0) newShift = (sPres->value+sPres->neighb[1]->value)/2;
368: PEPQSliceGetInertia(pep,newShift,&iner,&zeros,sr->intcorr);
369: if (zeros) SETERRQ1(((PetscObject)pep)->comm,PETSC_ERR_CONV_FAILED,"Inertia computation fails in %g",newShift);
370: sPres->value = newShift;
371: }
372: sr->sPres->inertia = iner;
373: pep->target = sr->sPres->value;
374: pep->reason = PEP_CONVERGED_ITERATING;
375: pep->its = 0;
376: } else sr->sPres = NULL;
377: return(0);
378: }
380: /*
381: Obtains value of subsequent shift
382: */
383: static PetscErrorCode PEPGetNewShiftValue(PEP pep,PetscInt side,PetscReal *newS)
384: {
385: PetscReal lambda,d_prev;
386: PetscInt i,idxP;
387: PEP_SR sr;
388: PEP_shift sPres,s;
389: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
392: sr = ctx->sr;
393: sPres = sr->sPres;
394: if (sPres->neighb[side]) {
395: /* Completing a previous interval */
396: if (!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0) { /* One of the ends might be too far from eigenvalues */
397: if (side) *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[sr->indexEig-1]]))/2;
398: else *newS = (sPres->value + PetscRealPart(sr->eigr[sr->perm[0]]))/2;
399: } else *newS=(sPres->value + sPres->neighb[side]->value)/2;
400: } else { /* (Only for side=1). Creating a new interval. */
401: if (sPres->neigs==0) {/* No value has been accepted*/
402: if (sPres->neighb[0]) {
403: /* Multiplying by 10 the previous distance */
404: *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
405: sr->nleap++;
406: /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
407: if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Unable to compute the wanted eigenvalues with open interval");
408: } else { /* First shift */
409: if (pep->nconv != 0) {
410: /* Unaccepted values give information for next shift */
411: idxP=0;/* Number of values left from shift */
412: for (i=0;i<pep->nconv;i++) {
413: lambda = PetscRealPart(pep->eigr[i]);
414: if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
415: else break;
416: }
417: /* Avoiding subtraction of eigenvalues (might be the same).*/
418: if (idxP>0) {
419: d_prev = PetscAbsReal(sPres->value - PetscRealPart(pep->eigr[0]))/(idxP+0.3);
420: } else {
421: d_prev = PetscAbsReal(sPres->value - PetscRealPart(pep->eigr[pep->nconv-1]))/(pep->nconv+0.3);
422: }
423: *newS = sPres->value + ((sr->dir)*d_prev*pep->nev)/2;
424: } else { /* No values found, no information for next shift */
425: SETERRQ(PetscObjectComm((PetscObject)pep),1,"First shift renders no information");
426: }
427: }
428: } else { /* Accepted values found */
429: sr->nleap = 0;
430: /* Average distance of values in previous subinterval */
431: s = sPres->neighb[0];
432: while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
433: s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
434: }
435: if (s) {
436: d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
437: } else { /* First shift. Average distance obtained with values in this shift */
438: /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
439: 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)) {
440: d_prev = PetscAbsReal((PetscRealPart(sr->eigr[sr->indexEig-1]) - PetscRealPart(sr->eigr[0])))/(sPres->neigs+0.3);
441: } else {
442: d_prev = PetscAbsReal(PetscRealPart(sr->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
443: }
444: }
445: /* Average distance is used for next shift by adding it to value on the right or to shift */
446: if ((sr->dir)*(PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
447: *newS = PetscRealPart(sr->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(pep->nev))/2;
448: } else { /* Last accepted value is on the left of shift. Adding to shift */
449: *newS = sPres->value + ((sr->dir)*d_prev*(pep->nev))/2;
450: }
451: }
452: /* End of interval can not be surpassed */
453: if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
454: }/* of neighb[side]==null */
455: return(0);
456: }
458: /*
459: Function for sorting an array of real values
460: */
461: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
462: {
463: PetscReal re;
464: PetscInt i,j,tmp;
467: if (!prev) for (i=0;i<nr;i++) perm[i] = i;
468: /* Insertion sort */
469: for (i=1;i<nr;i++) {
470: re = PetscRealPart(r[perm[i]]);
471: j = i-1;
472: while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
473: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
474: }
475: }
476: return(0);
477: }
479: /* Stores the pairs obtained since the last shift in the global arrays */
480: static PetscErrorCode PEPStoreEigenpairs(PEP pep)
481: {
483: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
484: PetscReal lambda,err,*errest;
485: PetscInt i,*aux,count=0,ndef,ld,nconv=pep->nconv,d=pep->nmat-1,idx;
486: PetscBool iscayley,divide=PETSC_FALSE;
487: PEP_SR sr = ctx->sr;
488: PEP_shift sPres;
489: Vec w;
490: Mat MS;
491: BV tV;
492: PetscScalar *S,*eigr,*tS;
495: sPres = sr->sPres;
496: sPres->index = sr->indexEig;
498: /* Back-transform */
499: STBackTransform(pep->st,nconv,pep->eigr,pep->eigi);
501: PetscObjectTypeCompare((PetscObject)pep->st,STCAYLEY,&iscayley);
502: /* Sort eigenvalues */
503: sortRealEigenvalues(pep->eigr,pep->perm,nconv,PETSC_FALSE,sr->dir);
504: BVGetSizes(pep->V,NULL,NULL,&ld);
505: BVTensorGetFactors(ctx->V,NULL,&MS);
506: MatDenseGetArray(MS,&S);
507: /* Values stored in global array */
508: PetscCalloc4(nconv,&eigr,nconv,&errest,nconv*nconv*d,&tS,nconv,&aux);
509: ndef = sr->ndef0+sr->ndef1;
510: for (i=0;i<nconv;i++) {
511: lambda = PetscRealPart(pep->eigr[pep->perm[i]]);
512: err = pep->errest[pep->perm[i]];
513: if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
514: if (sr->indexEig+count-ndef>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Unexpected error in Spectrum Slicing");
515: eigr[count] = lambda;
516: errest[count] = err;
517: if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) sPres->nconv[0]++;
518: if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) sPres->nconv[1]++;
519: PetscMemcpy(tS+count*(d*nconv),S+pep->perm[i]*(d*ld),nconv*sizeof(PetscScalar));
520: PetscMemcpy(tS+count*(d*nconv)+nconv,S+pep->perm[i]*(d*ld)+ld,nconv*sizeof(PetscScalar));
521: count++;
522: }
523: }
524: for (i=0;i<count;i++) {
525: PetscMemcpy(S+i*(d*ld),tS+i*nconv*d,nconv*sizeof(PetscScalar));
526: PetscMemcpy(S+i*(d*ld)+ld,tS+i*nconv*d+nconv,nconv*sizeof(PetscScalar));
527: }
528: MatDenseRestoreArray(MS,&S);
529: BVTensorRestoreFactors(ctx->V,NULL,&MS);
530: BVSetActiveColumns(ctx->V,0,count);
531: BVTensorCompress(ctx->V,count);
532: if (sr->sPres->nconv[0] && sr->sPres->nconv[1]) {
533: divide = PETSC_TRUE;
534: BVTensorGetFactors(ctx->V,NULL,&MS);
535: MatDenseGetArray(MS,&S);
536: PetscMemzero(tS,nconv*nconv*d*sizeof(PetscScalar));
537: for (i=0;i<count;i++) {
538: PetscMemcpy(tS+i*nconv*d,S+i*(d*ld),count*sizeof(PetscScalar));
539: PetscMemcpy(tS+i*nconv*d+nconv,S+i*(d*ld)+ld,count*sizeof(PetscScalar));
540: }
541: MatDenseRestoreArray(MS,&S);
542: BVTensorRestoreFactors(ctx->V,NULL,&MS);
543: BVSetActiveColumns(pep->V,0,count);
544: BVDuplicateResize(pep->V,count,&tV);
545: BVCopy(pep->V,tV);
546: }
547: if (sr->sPres->nconv[0]) {
548: if (divide) {
549: BVSetActiveColumns(ctx->V,0,sr->sPres->nconv[0]);
550: BVTensorCompress(ctx->V,sr->sPres->nconv[0]);
551: }
552: for (i=0;i<sr->ndef0;i++) aux[i] = sr->idxDef0[i];
553: for (i=sr->ndef0;i<sr->sPres->nconv[0];i++) aux[i] = sr->indexEig+i-sr->ndef0;
554: BVTensorGetFactors(ctx->V,NULL,&MS);
555: MatDenseGetArray(MS,&S);
556: for (i=0;i<sr->sPres->nconv[0];i++) {
557: sr->eigr[aux[i]] = eigr[i];
558: sr->errest[aux[i]] = errest[i];
559: BVGetColumn(pep->V,i,&w);
560: BVInsertVec(sr->V,aux[i],w);
561: BVRestoreColumn(pep->V,i,&w);
562: idx = sr->ld*d*aux[i];
563: PetscMemzero(sr->S+idx,sr->ld*d*sizeof(PetscScalar));
564: PetscMemcpy(sr->S+idx,S+i*(ld*d),sr->sPres->nconv[0]*sizeof(PetscScalar));
565: PetscMemcpy(sr->S+idx+sr->ld,S+i*(ld*d)+ld,sr->sPres->nconv[0]*sizeof(PetscScalar));
566: PetscFree(sr->qinfo[aux[i]].q);
567: PetscMalloc1(sr->sPres->nconv[0],&sr->qinfo[aux[i]].q);
568: PetscMemcpy(sr->qinfo[aux[i]].q,aux,sr->sPres->nconv[0]*sizeof(PetscInt));
569: sr->qinfo[aux[i]].nq = sr->sPres->nconv[0];
570: }
571: MatDenseRestoreArray(MS,&S);
572: BVTensorRestoreFactors(ctx->V,NULL,&MS);
573: }
575: if (sr->sPres->nconv[1]) {
576: if (divide) {
577: BVTensorGetFactors(ctx->V,NULL,&MS);
578: MatDenseGetArray(MS,&S);
579: for (i=0;i<sr->sPres->nconv[1];i++) {
580: PetscMemcpy(S+i*(d*ld),tS+(sr->sPres->nconv[0]+i)*nconv*d,count*sizeof(PetscScalar));
581: PetscMemcpy(S+i*(d*ld)+ld,tS+(sr->sPres->nconv[0]+i)*nconv*d+nconv,count*sizeof(PetscScalar));
582: }
583: MatDenseRestoreArray(MS,&S);
584: BVTensorRestoreFactors(ctx->V,NULL,&MS);
585: BVSetActiveColumns(pep->V,0,count);
586: BVCopy(tV,pep->V);
587: BVSetActiveColumns(ctx->V,0,sr->sPres->nconv[1]);
588: BVTensorCompress(ctx->V,sr->sPres->nconv[1]);
589: }
590: for (i=0;i<sr->ndef1;i++) aux[i] = sr->idxDef1[i];
591: for (i=sr->ndef1;i<sr->sPres->nconv[1];i++) aux[i] = sr->indexEig+sr->sPres->nconv[0]-sr->ndef0+i-sr->ndef1;
592: BVTensorGetFactors(ctx->V,NULL,&MS);
593: MatDenseGetArray(MS,&S);
594: for (i=0;i<sr->sPres->nconv[1];i++) {
595: sr->eigr[aux[i]] = eigr[sr->sPres->nconv[0]+i];
596: sr->errest[aux[i]] = errest[sr->sPres->nconv[0]+i];
597: BVGetColumn(pep->V,i,&w);
598: BVInsertVec(sr->V,aux[i],w);
599: BVRestoreColumn(pep->V,i,&w);
600: idx = sr->ld*d*aux[i];
601: PetscMemzero(sr->S+idx,sr->ld*d*sizeof(PetscScalar));
602: PetscMemcpy(sr->S+idx,S+i*(ld*d),sr->sPres->nconv[1]*sizeof(PetscScalar));
603: PetscMemcpy(sr->S+idx+sr->ld,S+i*(ld*d)+ld,sr->sPres->nconv[1]*sizeof(PetscScalar));
604: PetscFree(sr->qinfo[aux[i]].q);
605: PetscMalloc1(sr->sPres->nconv[1],&sr->qinfo[aux[i]].q);
606: PetscMemcpy(sr->qinfo[aux[i]].q,aux,sr->sPres->nconv[1]*sizeof(PetscInt));
607: sr->qinfo[aux[i]].nq = sr->sPres->nconv[1];
608: }
609: MatDenseRestoreArray(MS,&S);
610: BVTensorRestoreFactors(ctx->V,NULL,&MS);
611: }
612: sPres->neigs = count-sr->ndef0-sr->ndef1;
613: sr->indexEig += sPres->neigs;
614: sPres->nconv[0]-= sr->ndef0;
615: sPres->nconv[1]-= sr->ndef1;
616: /* Global ordering array updating */
617: sortRealEigenvalues(sr->eigr,sr->perm,sr->indexEig,PETSC_FALSE,sr->dir);
618: /* Check for completion */
619: sPres->comp[0] = PetscNot(sPres->nconv[0] < sPres->nsch[0]);
620: sPres->comp[1] = PetscNot(sPres->nconv[1] < sPres->nsch[1]);
621: 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, consider using PEPKrylovSchurSetDetectZeros()");
622: PetscFree4(eigr,errest,tS,aux);
623: if (divide) { BVDestroy(&tV); }
624: return(0);
625: }
627: static PetscErrorCode PEPLookForDeflation(PEP pep)
628: {
629: PetscReal val;
630: PetscInt i,count0=0,count1=0;
631: PEP_shift sPres;
632: PetscInt ini,fin;
633: PEP_SR sr;
634: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
637: sr = ctx->sr;
638: sPres = sr->sPres;
640: if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
641: else ini = 0;
642: fin = sr->indexEig;
643: /* Selection of ends for searching new values */
644: if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
645: else sPres->ext[0] = sPres->neighb[0]->value;
646: if (!sPres->neighb[1]) {
647: if (sr->hasEnd) sPres->ext[1] = sr->int1;
648: else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
649: } else sPres->ext[1] = sPres->neighb[1]->value;
650: /* Selection of values between right and left ends */
651: for (i=ini;i<fin;i++) {
652: val=PetscRealPart(sr->eigr[sr->perm[i]]);
653: /* Values to the right of left shift */
654: if ((sr->dir)*(val - sPres->ext[1]) < 0) {
655: if ((sr->dir)*(val - sPres->value) < 0) count0++;
656: else count1++;
657: } else break;
658: }
659: /* The number of values on each side are found */
660: if (sPres->neighb[0]) {
661: sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
662: if (sPres->nsch[0]<0) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Mismatch between number of values found and information from inertia, consider using PEPSTOARSetDetectZeros()");
663: } else sPres->nsch[0] = 0;
665: if (sPres->neighb[1]) {
666: sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
667: if (sPres->nsch[1]<0) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Mismatch between number of values found and information from inertia, consider using PEPSTOARSetDetectZeros()");
668: } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);
670: /* Completing vector of indexes for deflation */
671: for (i=0;i<count0;i++) sr->idxDef0[i] = sr->perm[ini+i];
672: sr->ndef0 = count0;
673: for (i=0;i<count1;i++) sr->idxDef1[i] = sr->perm[ini+count0+i];
674: sr->ndef1 = count1;
675: return(0);
676: }
678: /*
679: Compute a run of Lanczos iterations
680: */
681: static PetscErrorCode PEPSTOARrun_QSlice(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,Vec *t_)
682: {
684: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
685: PetscInt i,j,m=*M,l,lock;
686: PetscInt lds,d,ld,offq,nqt;
687: Vec v=t_[0],t=t_[1],q=t_[2];
688: PetscReal norm,sym=0.0,fro=0.0,*f;
689: PetscScalar *y,*S,sigma;
690: PetscBLASInt j_,one=1;
691: PetscBool lindep;
692: Mat MS;
695: PetscMalloc1(*M,&y);
696: BVGetSizes(pep->V,NULL,NULL,&ld);
697: BVTensorGetDegree(ctx->V,&d);
698: BVGetActiveColumns(pep->V,&lock,&nqt);
699: lds = d*ld;
700: offq = ld;
702: *breakdown = PETSC_FALSE; /* ----- */
703: STGetShift(pep->st,&sigma);
704: DSGetDimensions(pep->ds,NULL,NULL,&l,NULL,NULL);
705: BVSetActiveColumns(ctx->V,0,m);
706: BVSetActiveColumns(pep->V,0,nqt);
707: for (j=k;j<m;j++) {
708: /* apply operator */
709: BVTensorGetFactors(ctx->V,NULL,&MS);
710: MatDenseGetArray(MS,&S);
711: BVGetColumn(pep->V,nqt,&t);
712: BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
713: MatMult(pep->A[1],v,q);
714: MatMult(pep->A[2],v,t);
715: VecAXPY(q,sigma*pep->sfactor,t);
716: VecScale(q,pep->sfactor);
717: BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
718: MatMult(pep->A[2],v,t);
719: VecAXPY(q,pep->sfactor*pep->sfactor,t);
720: STMatSolve(pep->st,q,t);
721: VecScale(t,-1.0);
722: BVRestoreColumn(pep->V,nqt,&t);
724: /* orthogonalize */
725: BVOrthogonalizeColumn(pep->V,nqt,S+(j+1)*lds,&norm,&lindep);
726: if (!lindep) {
727: *(S+(j+1)*lds+nqt) = norm;
728: BVScaleColumn(pep->V,nqt,1.0/norm);
729: nqt++;
730: }
731: for (i=0;i<nqt;i++) *(S+(j+1)*lds+offq+i) = *(S+j*lds+i)+sigma*(*(S+(j+1)*lds+i));
732: BVSetActiveColumns(pep->V,0,nqt);
733: MatDenseRestoreArray(MS,&S);
734: BVTensorRestoreFactors(ctx->V,NULL,&MS);
736: /* level-2 orthogonalization */
737: BVOrthogonalizeColumn(ctx->V,j+1,y,&norm,&lindep);
738: a[j] = PetscRealPart(y[j]);
739: omega[j+1] = (norm > 0)?1.0:-1.0;
740: BVScaleColumn(ctx->V,j+1,1.0/norm);
741: b[j] = PetscAbsReal(norm);
743: /* check symmetry */
744: DSGetArrayReal(pep->ds,DS_MAT_T,&f);
745: if (j==k) {
746: for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ld+i]);
747: for (i=0;i<l;i++) y[i] = 0.0;
748: }
749: DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
750: if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
751: PetscBLASIntCast(j,&j_);
752: sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
753: fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
754: if (j>0) fro = SlepcAbs(fro,b[j-1]);
755: if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
756: *symmlost = PETSC_TRUE;
757: *M=j;
758: break;
759: }
760: }
761: BVSetActiveColumns(pep->V,lock,nqt);
762: BVSetActiveColumns(ctx->V,0,*M);
763: PetscFree(y);
764: return(0);
765: }
767: static PetscErrorCode PEPSTOAR_QSlice(PEP pep)
768: {
770: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
771: PetscInt j,k,l,nv=0,ld,ldds,t,nq=0,m,n,idx;
772: PetscInt nconv=0,deg=pep->nmat-1,count0=0,count1=0;
773: PetscScalar *Q,*om,scal[2],sigma,*back,*S,*pQ;
774: PetscReal beta,norm=1.0,*omega,*a,*b,*r,eta,lambda;
775: PetscBool breakdown,symmlost=PETSC_FALSE,sinv,falselock=PETSC_TRUE;
776: Mat MS,MQ,A,pA[4],As[2],D[2];
777: Vec v,vomega;
778: ShellMatCtx *ctxMat[2];
779: PEP_SR sr;
780: BVOrthogType otype;
781: BVOrthogBlockType obtype;
784: PetscCitationsRegister(citation,&cited);
786: /* Resize if needed for deflating vectors */
787: sr = ctx->sr;
788: sigma = sr->sPres->value;
789: k = sr->ndef0+sr->ndef1;
790: pep->ncv = ctx->ncv+k;
791: pep->nev = ctx->nev+k;
792: PEPAllocateSolution(pep,2);
793: BVDestroy(&ctx->V);
794: BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
795: BVGetOrthogonalization(pep->V,&otype,NULL,&eta,&obtype);
796: BVSetOrthogonalization(ctx->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
797: DSAllocate(pep->ds,pep->ncv+2);
798: PetscMalloc1(pep->ncv,&back);
799: DSGetLeadingDimension(pep->ds,&ldds);
801: STGetMatrixTransformed(pep->st,2,&D[1]); /* M */
802: MatGetLocalSize(D[1],&m,&n);
803: STGetMatrixTransformed(pep->st,0,&D[0]); /* K */
804: scal[0] = -1.0; scal[1] = pep->sfactor*pep->sfactor;
805: for (j=0;j<2;j++) {
806: PetscNew(ctxMat+j);
807: (ctxMat[j])->A = D[j]; (ctxMat[j])->scal = scal[j];
808: MatCreateShell(PetscObjectComm((PetscObject)pep),m,n,PETSC_DETERMINE,PETSC_DETERMINE,ctxMat[j],&As[j]);
809: MatShellSetOperation(As[j],MATOP_MULT,(void(*)())MatMult_Func);
810: MatShellSetOperation(As[j],MATOP_DESTROY,(void(*)())MatDestroy_Func);
811: }
812: pA[0] = As[0]; pA[1] = pA[2] = NULL; pA[3] = As[1];
813: MatCreateNest(PetscObjectComm((PetscObject)pep),2,NULL,2,NULL,pA,&A);
814: for (j=0;j<2;j++) { MatDestroy(&As[j]); }
815: BVSetMatrix(ctx->V,A,PETSC_TRUE);
816: MatDestroy(&A);
817: if (ctx->lock) {
818: PetscOptionsGetBool(NULL,NULL,"-pep_stoar_falselocking",&falselock,NULL);
819: } else SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"A locking variant is needed for spectrum slicing");
820: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
821: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
822: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
824: /* Get the starting Arnoldi vector */
825: BVTensorBuildFirstColumn(ctx->V,pep->nini);
826: BVSetActiveColumns(ctx->V,0,1);
827: if (k) {
828: /* Insert deflated vectors */
829: BVSetActiveColumns(pep->V,0,0);
830: idx = sr->ndef0?sr->idxDef0[0]:sr->idxDef1[0];
831: for (j=0;j<k;j++) {
832: BVGetColumn(pep->V,j,&v);
833: BVCopyVec(sr->V,sr->qinfo[idx].q[j],v);
834: BVRestoreColumn(pep->V,j,&v);
835: }
836: /* Update innerproduct matrix */
837: BVSetActiveColumns(ctx->V,0,0);
838: BVTensorGetFactors(ctx->V,NULL,&MS);
839: BVSetActiveColumns(pep->V,0,k);
840: BVTensorRestoreFactors(ctx->V,NULL,&MS);
842: BVGetSizes(pep->V,NULL,NULL,&ld);
843: BVTensorGetFactors(ctx->V,NULL,&MS);
844: MatDenseGetArray(MS,&S);
845: for (j=0;j<sr->ndef0;j++) {
846: PetscMemzero(S+j*ld*deg,ld*deg*sizeof(PetscScalar));
847: PetscMemcpy(S+j*ld*deg,sr->S+sr->idxDef0[j]*sr->ld*deg,k*sizeof(PetscScalar));
848: PetscMemcpy(S+j*ld*deg+ld,sr->S+sr->idxDef0[j]*sr->ld*deg+sr->ld,k*sizeof(PetscScalar));
849: pep->eigr[j] = sr->eigr[sr->idxDef0[j]];
850: pep->errest[j] = sr->errest[sr->idxDef0[j]];
851: }
852: for (j=0;j<sr->ndef1;j++) {
853: PetscMemzero(S+(j+sr->ndef0)*ld*deg,ld*deg*sizeof(PetscScalar));
854: PetscMemcpy(S+(j+sr->ndef0)*ld*deg,sr->S+sr->idxDef1[j]*sr->ld*deg,k*sizeof(PetscScalar));
855: PetscMemcpy(S+(j+sr->ndef0)*ld*deg+ld,sr->S+sr->idxDef1[j]*sr->ld*deg+sr->ld,k*sizeof(PetscScalar));
856: pep->eigr[j+sr->ndef0] = sr->eigr[sr->idxDef1[j]];
857: pep->errest[j+sr->ndef0] = sr->errest[sr->idxDef1[j]];
858: }
859: MatDenseRestoreArray(MS,&S);
860: BVTensorRestoreFactors(ctx->V,NULL,&MS);
861: BVSetActiveColumns(ctx->V,0,k+1);
862: VecCreateSeq(PETSC_COMM_SELF,k+1,&vomega);
863: VecGetArray(vomega,&om);
864: for (j=0;j<k;j++) {
865: BVOrthogonalizeColumn(ctx->V,j,NULL,&norm,NULL);
866: BVScaleColumn(ctx->V,j,1/norm);
867: om[j] = (norm>=0.0)?1.0:-1.0;
868: }
869: BVTensorGetFactors(ctx->V,NULL,&MS);
870: MatDenseGetArray(MS,&S);
871: for (j=0;j<deg;j++) {
872: BVSetRandomColumn(pep->V,k+j);
873: BVOrthogonalizeColumn(pep->V,k+j,S+k*ld*deg+j*ld,&norm,NULL);
874: BVScaleColumn(pep->V,k+j,1.0/norm);
875: S[k*ld*deg+j*ld+k+j] = norm;
876: }
877: MatDenseRestoreArray(MS,&S);
878: BVSetActiveColumns(pep->V,0,k+deg);
879: BVTensorRestoreFactors(ctx->V,NULL,&MS);
880: BVOrthogonalizeColumn(ctx->V,k,NULL,&norm,NULL);
881: BVScaleColumn(ctx->V,k,1.0/norm);
882: om[k] = (norm>=0.0)?1.0:-1.0;
883: VecRestoreArray(vomega,&om);
884: BVSetSignature(ctx->V,vomega);
885: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
886: VecGetArray(vomega,&om);
887: for (j=0;j<k;j++) a[j] = PetscRealPart(om[j]/(pep->eigr[j]-sigma));
888: VecRestoreArray(vomega,&om);
889: VecDestroy(&vomega);
890: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
891: DSGetArray(pep->ds,DS_MAT_Q,&pQ);
892: PetscMemzero(pQ,ldds*k*sizeof(PetscScalar));
893: for (j=0;j<k;j++) pQ[j+j*ldds] = 1.0;
894: DSRestoreArray(pep->ds,DS_MAT_Q,&pQ);
895: }
896: BVSetActiveColumns(ctx->V,0,k+1);
897: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
898: VecCreateSeq(PETSC_COMM_SELF,k+1,&vomega);
899: BVGetSignature(ctx->V,vomega);
900: VecGetArray(vomega,&om);
901: for (j=0;j<k+1;j++) omega[j] = PetscRealPart(om[j]);
902: VecRestoreArray(vomega,&om);
903: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
904: VecDestroy(&vomega);
906: 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]);
908: /* Restart loop */
909: l = 0;
910: pep->nconv = k;
911: while (pep->reason == PEP_CONVERGED_ITERATING) {
912: pep->its++;
913: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
914: b = a+ldds;
915: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
917: /* Compute an nv-step Lanczos factorization */
918: nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
919: PEPSTOARrun_QSlice(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,pep->work);
920: beta = b[nv-1];
921: if (symmlost && nv==pep->nconv+l) {
922: pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
923: pep->nconv = nconv;
924: if (falselock || !ctx->lock) {
925: BVSetActiveColumns(ctx->V,0,pep->nconv);
926: BVTensorCompress(ctx->V,0);
927: }
928: break;
929: }
930: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
931: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
932: DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
933: if (l==0) {
934: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
935: } else {
936: DSSetState(pep->ds,DS_STATE_RAW);
937: }
939: /* Solve projected problem */
940: DSSolve(pep->ds,pep->eigr,pep->eigi);
941: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
942: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
944: /* Check convergence */
945: /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
946: norm = 1.0;
947: DSGetDimensions(pep->ds,NULL,NULL,NULL,NULL,&t);
948: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
949: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
950: for (j=0;j<k;j++) back[j] = pep->eigr[j];
951: STBackTransform(pep->st,k,back,pep->eigi);
952: count0=count1=0;
953: for (j=0;j<k;j++) {
954: lambda = PetscRealPart(back[j]);
955: if (((sr->dir)*(sr->sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sr->sPres->ext[0]) > 0)) count0++;
956: if (((sr->dir)*(lambda - sr->sPres->value) > 0) && ((sr->dir)*(sr->sPres->ext[1] - lambda) > 0)) count1++;
957: }
958: if ((count0-sr->ndef0 >= sr->sPres->nsch[0]) && (count1-sr->ndef1 >= sr->sPres->nsch[1])) pep->reason = PEP_CONVERGED_TOL;
959: /* Update l */
960: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
961: else {
962: l = PetscMax(1,(PetscInt)((nv-k)/2));
963: l = PetscMin(l,t);
964: if (!breakdown) {
965: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
966: if (*(a+ldds+k+l-1)!=0) {
967: if (k+l<nv-1) l = l+1;
968: else l = l-1;
969: }
970: /* Prepare the Rayleigh quotient for restart */
971: DSGetArray(pep->ds,DS_MAT_Q,&Q);
972: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
973: r = a + 2*ldds;
974: for (j=k;j<k+l;j++) {
975: r[j] = PetscRealPart(Q[nv-1+j*ldds]*beta);
976: }
977: b = a+ldds;
978: b[k+l-1] = r[k+l-1];
979: omega[k+l] = omega[nv];
980: DSRestoreArray(pep->ds,DS_MAT_Q,&Q);
981: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
982: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
983: }
984: }
985: nconv = k;
986: if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
988: /* Update S */
989: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
990: BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
991: MatDestroy(&MQ);
993: /* Copy last column of S */
994: BVCopyColumn(ctx->V,nv,k+l);
995: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
996: VecCreateSeq(PETSC_COMM_SELF,k+l,&vomega);
997: VecGetArray(vomega,&om);
998: for (j=0;j<k+l;j++) om[j] = omega[j];
999: VecRestoreArray(vomega,&om);
1000: BVSetActiveColumns(ctx->V,0,k+l);
1001: BVSetSignature(ctx->V,vomega);
1002: VecDestroy(&vomega);
1003: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
1005: if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
1006: /* stop if breakdown */
1007: PetscInfo2(pep,"Breakdown TOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
1008: pep->reason = PEP_DIVERGED_BREAKDOWN;
1009: }
1010: if (pep->reason != PEP_CONVERGED_ITERATING) l--;
1011: BVGetActiveColumns(pep->V,NULL,&nq);
1012: if (k+l+deg<=nq) {
1013: BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
1014: if (!falselock && ctx->lock) {
1015: BVTensorCompress(ctx->V,k-pep->nconv);
1016: } else {
1017: BVTensorCompress(ctx->V,0);
1018: }
1019: }
1020: pep->nconv = k;
1021: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
1022: }
1023: sr->itsKs += pep->its;
1024: if (pep->nconv>0) {
1025: BVSetActiveColumns(ctx->V,0,pep->nconv);
1026: BVGetActiveColumns(pep->V,NULL,&nq);
1027: BVSetActiveColumns(pep->V,0,nq);
1028: if (nq>pep->nconv) {
1029: BVTensorCompress(ctx->V,pep->nconv);
1030: BVSetActiveColumns(pep->V,0,pep->nconv);
1031: }
1032: for (j=0;j<pep->nconv;j++) {
1033: pep->eigr[j] *= pep->sfactor;
1034: pep->eigi[j] *= pep->sfactor;
1035: }
1036: }
1037: 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);
1038: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
1039: RGPopScale(pep->rg);
1041: /* truncate Schur decomposition and change the state to raw so that
1042: DSVectors() computes eigenvectors from scratch */
1043: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
1044: DSSetState(pep->ds,DS_STATE_RAW);
1045: PetscFree(back);
1046: return(0);
1047: }
1049: #define SWAP(a,b,t) {t=a;a=b;b=t;}
1051: static PetscErrorCode PEPQSliceGetInertias(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
1052: {
1053: PetscErrorCode ierr;
1054: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
1055: PEP_SR sr=ctx->sr;
1056: PetscInt i=0,j,tmpi;
1057: PetscReal v,tmpr;
1058: PEP_shift s;
1061: if (!pep->state) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSetUp() first");
1062: if (!sr) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see PEPSetInterval()");
1063: if (!sr->s0) { /* PEPSolve not called yet */
1064: *n = 2;
1065: } else {
1066: *n = 1;
1067: s = sr->s0;
1068: while (s) {
1069: (*n)++;
1070: s = s->neighb[1];
1071: }
1072: }
1073: PetscMalloc1(*n,shifts);
1074: PetscMalloc1(*n,inertias);
1075: if (!sr->s0) { /* PEPSolve not called yet */
1076: (*shifts)[0] = sr->int0;
1077: (*shifts)[1] = sr->int1;
1078: (*inertias)[0] = sr->inertia0;
1079: (*inertias)[1] = sr->inertia1;
1080: } else {
1081: s = sr->s0;
1082: while (s) {
1083: (*shifts)[i] = s->value;
1084: (*inertias)[i++] = s->inertia;
1085: s = s->neighb[1];
1086: }
1087: (*shifts)[i] = sr->int1;
1088: (*inertias)[i] = sr->inertia1;
1089: }
1090: /* remove possible duplicate in last position */
1091: if ((*shifts)[(*n)-1]==(*shifts)[(*n)-2]) (*n)--;
1092: /* sort result */
1093: for (i=0;i<*n;i++) {
1094: v = (*shifts)[i];
1095: for (j=i+1;j<*n;j++) {
1096: if (v > (*shifts)[j]) {
1097: SWAP((*shifts)[i],(*shifts)[j],tmpr);
1098: SWAP((*inertias)[i],(*inertias)[j],tmpi);
1099: v = (*shifts)[i];
1100: }
1101: }
1102: }
1103: return(0);
1104: }
1106: PetscErrorCode PEPSolve_STOAR_QSlice(PEP pep)
1107: {
1109: PetscInt i,j,ti,deg=pep->nmat-1;
1110: PetscReal newS;
1111: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
1112: PEP_SR sr=ctx->sr;
1113: Mat S;
1114: PetscScalar *pS;
1117: /* Only with eigenvalues present in the interval ...*/
1118: if (sr->numEigs==0) {
1119: pep->reason = PEP_CONVERGED_TOL;
1120: return(0);
1121: }
1122: /* Array of pending shifts */
1123: sr->maxPend = 100; /* Initial size */
1124: sr->nPend = 0;
1125: PetscMalloc1(sr->maxPend,&sr->pending);
1126: PetscLogObjectMemory((PetscObject)pep,(sr->maxPend)*sizeof(PEP_shift));
1127: PEPCreateShift(pep,sr->int0,NULL,NULL);
1128: /* extract first shift */
1129: sr->sPrev = NULL;
1130: sr->sPres = sr->pending[--sr->nPend];
1131: sr->sPres->inertia = sr->inertia0;
1132: pep->target = sr->sPres->value;
1133: sr->s0 = sr->sPres;
1134: sr->indexEig = 0;
1135: /* Memory reservation for auxiliary variables */
1136: PetscLogObjectMemory((PetscObject)pep,(sr->numEigs+2*pep->ncv)*sizeof(PetscScalar));
1137: for (i=0;i<sr->numEigs;i++) {
1138: sr->eigr[i] = 0.0;
1139: sr->eigi[i] = 0.0;
1140: sr->errest[i] = 0.0;
1141: sr->perm[i] = i;
1142: }
1143: /* Vectors for deflation */
1144: PetscMalloc2(sr->numEigs,&sr->idxDef0,sr->numEigs,&sr->idxDef1);
1145: PetscLogObjectMemory((PetscObject)pep,sr->numEigs*sizeof(PetscInt));
1146: sr->indexEig = 0;
1147: while (sr->sPres) {
1148: /* Search for deflation */
1149: PEPLookForDeflation(pep);
1150: /* KrylovSchur */
1151: PEPSTOAR_QSlice(pep);
1153: PEPStoreEigenpairs(pep);
1154: /* Select new shift */
1155: if (!sr->sPres->comp[1]) {
1156: PEPGetNewShiftValue(pep,1,&newS);
1157: PEPCreateShift(pep,newS,sr->sPres,sr->sPres->neighb[1]);
1158: }
1159: if (!sr->sPres->comp[0]) {
1160: /* Completing earlier interval */
1161: PEPGetNewShiftValue(pep,0,&newS);
1162: PEPCreateShift(pep,newS,sr->sPres->neighb[0],sr->sPres);
1163: }
1164: /* Preparing for a new search of values */
1165: PEPExtractShift(pep);
1166: }
1168: /* Updating pep values prior to exit */
1169: PetscFree2(sr->idxDef0,sr->idxDef1);
1170: PetscFree(sr->pending);
1171: pep->nconv = sr->indexEig;
1172: pep->reason = PEP_CONVERGED_TOL;
1173: pep->its = sr->itsKs;
1174: pep->nev = sr->indexEig;
1175: MatCreateSeqDense(PETSC_COMM_SELF,pep->nconv,pep->nconv,NULL,&S);
1176: MatDenseGetArray(S,&pS);
1177: for (i=0;i<pep->nconv;i++) {
1178: 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);
1179: }
1180: MatDenseRestoreArray(S,&pS);
1181: BVSetActiveColumns(sr->V,0,pep->nconv);
1182: BVMultInPlace(sr->V,S,0,pep->nconv);
1183: MatDestroy(&S);
1184: BVDestroy(&pep->V);
1185: pep->V = sr->V;
1186: PetscFree4(pep->eigr,pep->eigi,pep->errest,pep->perm);
1187: pep->eigr = sr->eigr;
1188: pep->eigi = sr->eigi;
1189: pep->perm = sr->perm;
1190: pep->errest = sr->errest;
1191: if (sr->dir<0) {
1192: for (i=0;i<pep->nconv/2;i++) {
1193: ti = sr->perm[i]; sr->perm[i] = sr->perm[pep->nconv-1-i]; sr->perm[pep->nconv-1-i] = ti;
1194: }
1195: }
1196: PetscFree(ctx->inertias);
1197: PetscFree(ctx->shifts);
1198: PEPQSliceGetInertias(pep,&ctx->nshifts,&ctx->shifts,&ctx->inertias);
1199: return(0);
1200: }