Actual source code: stoar.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: SLEPc polynomial eigensolver: "stoar"
13: Method: S-TOAR
15: Algorithm:
17: Symmetric Two-Level Orthogonal Arnoldi.
19: References:
21: [1] C. Campos and J.E. Roman, "Restarted Q-Arnoldi-type methods
22: exploiting symmetry in quadratic eigenvalue problems", BIT
23: Numer. Math. 56(4):1213-1236, 2016.
24: */
26: #include <slepc/private/pepimpl.h>
27: #include "../src/pep/impls/krylov/pepkrylov.h"
28: #include <slepcblaslapack.h>
30: static PetscBool cited = PETSC_FALSE;
31: static const char citation[] =
32: "@Article{slepc-stoar,\n"
33: " author = \"C. Campos and J. E. Roman\",\n"
34: " title = \"Restarted {Q-Arnoldi-type} methods exploiting symmetry in quadratic eigenvalue problems\",\n"
35: " journal = \"{BIT} Numer. Math.\",\n"
36: " volume = \"56\",\n"
37: " number = \"4\",\n"
38: " pages = \"1213--1236\",\n"
39: " year = \"2016,\"\n"
40: " doi = \"https://doi.org/10.1007/s10543-016-0601-5\"\n"
41: "}\n";
43: typedef struct {
44: PetscReal scal[2];
45: Mat A[2];
46: Vec t;
47: } PEP_STOAR_MATSHELL;
49: static PetscErrorCode MatMult_STOAR(Mat A,Vec x,Vec y)
50: {
51: PetscErrorCode ierr;
52: PEP_STOAR_MATSHELL *ctx;
55: MatShellGetContext(A,&ctx);
56: MatMult(ctx->A[0],x,y);
57: VecScale(y,ctx->scal[0]);
58: if (ctx->scal[1]) {
59: MatMult(ctx->A[1],x,ctx->t);
60: VecAXPY(y,ctx->scal[1],ctx->t);
61: }
62: return(0);
63: }
65: static PetscErrorCode MatDestroy_STOAR(Mat A)
66: {
67: PEP_STOAR_MATSHELL *ctx;
68: PetscErrorCode ierr;
71: MatShellGetContext(A,&ctx);
72: VecDestroy(&ctx->t);
73: PetscFree(ctx);
74: return(0);
75: }
77: PetscErrorCode PEPSTOARSetUpInnerMatrix(PEP pep,Mat *B)
78: {
79: Mat pB[4],Bs[3],D[3];
80: PetscInt i,j,n,m;
81: PEP_STOAR_MATSHELL *ctxMat[3];
82: PEP_STOAR *ctx=(PEP_STOAR*)pep->data;
83: PetscErrorCode ierr;
86: for (i=0;i<3;i++) {
87: STGetMatrixTransformed(pep->st,i,&D[i]); /* D[2] = M */
88: }
89: MatGetLocalSize(D[2],&m,&n);
91: for (j=0;j<3;j++) {
92: PetscNew(ctxMat+j);
93: MatCreateShell(PetscObjectComm((PetscObject)pep),m,n,PETSC_DETERMINE,PETSC_DETERMINE,ctxMat[j],&Bs[j]);
94: MatShellSetOperation(Bs[j],MATOP_MULT,(void(*)(void))MatMult_STOAR);
95: MatShellSetOperation(Bs[j],MATOP_DESTROY,(void(*)(void))MatDestroy_STOAR);
96: }
97: for (i=0;i<4;i++) pB[i] = NULL;
98: if (ctx->alpha) {
99: ctxMat[0]->A[0] = D[0]; ctxMat[0]->scal[0] = ctx->alpha; ctxMat[0]->scal[1] = 0.0;
100: ctxMat[2]->A[0] = D[2]; ctxMat[2]->scal[0] = -ctx->alpha*pep->sfactor*pep->sfactor; ctxMat[2]->scal[1] = 0.0;
101: pB[0] = Bs[0]; pB[3] = Bs[2];
102: }
103: if (ctx->beta) {
104: i = (ctx->alpha)?1:0;
105: ctxMat[0]->scal[1] = 0.0;
106: ctxMat[0]->A[i] = D[1]; ctxMat[0]->scal[i] = -ctx->beta*pep->sfactor;
107: ctxMat[1]->A[0] = D[2]; ctxMat[1]->scal[0] = -ctx->beta*pep->sfactor*pep->sfactor; ctxMat[1]->scal[1] = 0.0;
108: pB[0] = Bs[0]; pB[1] = pB[2] = Bs[1];
109: }
110: BVCreateVec(pep->V,&ctxMat[0]->t);
111: MatCreateNest(PetscObjectComm((PetscObject)pep),2,NULL,2,NULL,pB,B);
112: for (j=0;j<3;j++) { MatDestroy(&Bs[j]); }
113: return(0);
114: }
116: PetscErrorCode PEPSetUp_STOAR(PEP pep)
117: {
118: PetscErrorCode ierr;
119: PetscBool sinv,flg;
120: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
121: PetscInt ld,i;
122: PetscReal eta;
123: BVOrthogType otype;
124: BVOrthogBlockType obtype;
127: PEPCheckHermitian(pep);
128: PEPCheckQuadratic(pep);
129: PEPCheckShiftSinvert(pep);
130: /* spectrum slicing requires special treatment of default values */
131: if (pep->which==PEP_ALL) {
132: pep->ops->solve = PEPSolve_STOAR_QSlice;
133: pep->ops->extractvectors = NULL;
134: pep->ops->setdefaultst = NULL;
135: PEPSetUp_STOAR_QSlice(pep);
136: } else {
137: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
138: if (!ctx->lock && pep->mpd<pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
139: if (pep->max_it==PETSC_DEFAULT) pep->max_it = PetscMax(100,2*(pep->nmat-1)*pep->n/pep->ncv);
140: pep->ops->solve = PEPSolve_STOAR;
141: ld = pep->ncv+2;
142: DSSetType(pep->ds,DSGHIEP);
143: DSSetCompact(pep->ds,PETSC_TRUE);
144: DSSetExtraRow(pep->ds,PETSC_TRUE);
145: DSAllocate(pep->ds,ld);
146: PEPBasisCoefficients(pep,pep->pbc);
147: STGetTransform(pep->st,&flg);
148: if (!flg) {
149: PetscFree(pep->solvematcoeffs);
150: PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
151: PetscLogObjectMemory((PetscObject)pep,pep->nmat*sizeof(PetscScalar));
152: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
153: if (sinv) {
154: PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL);
155: } else {
156: for (i=0;i<pep->nmat-1;i++) pep->solvematcoeffs[i] = 0.0;
157: pep->solvematcoeffs[pep->nmat-1] = 1.0;
158: }
159: }
160: }
161: if (!pep->which) { PEPSetWhichEigenpairs_Default(pep); }
162: PEPCheckUnsupported(pep,PEP_FEATURE_NONMONOMIAL | PEP_FEATURE_REGION);
164: PEPAllocateSolution(pep,2);
165: PEPSetWorkVecs(pep,4);
166: BVDestroy(&ctx->V);
167: BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
168: BVGetOrthogonalization(pep->V,&otype,NULL,&eta,&obtype);
169: BVSetOrthogonalization(ctx->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
170: return(0);
171: }
173: /*
174: Compute a run of Lanczos iterations. dim(work)=(ctx->ld)*4
175: */
176: static PetscErrorCode PEPSTOARrun(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,Vec *t_)
177: {
179: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
180: PetscInt i,j,m=*M,l,lock;
181: PetscInt lds,d,ld,offq,nqt,ldds;
182: Vec v=t_[0],t=t_[1],q=t_[2];
183: PetscReal norm,sym=0.0,fro=0.0,*f;
184: PetscScalar *y,*S,*x,sigma;
185: PetscBLASInt j_,one=1;
186: PetscBool lindep,flg,sinvert=PETSC_FALSE;
187: Mat MS;
190: PetscMalloc1(*M,&y);
191: BVGetSizes(pep->V,NULL,NULL,&ld);
192: BVTensorGetDegree(ctx->V,&d);
193: BVGetActiveColumns(pep->V,&lock,&nqt);
194: lds = d*ld;
195: offq = ld;
196: DSGetLeadingDimension(pep->ds,&ldds);
197: *breakdown = PETSC_FALSE; /* ----- */
198: DSGetDimensions(pep->ds,NULL,&l,NULL,NULL);
199: BVSetActiveColumns(ctx->V,0,m);
200: BVSetActiveColumns(pep->V,0,nqt);
201: STGetTransform(pep->st,&flg);
202: if (!flg) {
203: /* spectral transformation handled by the solver */
204: PetscObjectTypeCompareAny((PetscObject)pep->st,&flg,STSINVERT,STSHIFT,"");
205: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"ST type not supported for TOAR without transforming matrices");
206: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinvert);
207: STGetShift(pep->st,&sigma);
208: }
209: for (j=k;j<m;j++) {
210: /* apply operator */
211: BVTensorGetFactors(ctx->V,NULL,&MS);
212: MatDenseGetArray(MS,&S);
213: BVGetColumn(pep->V,nqt,&t);
214: BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
215: if (!sinvert) {
216: STMatMult(pep->st,0,v,q);
217: BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
218: STMatMult(pep->st,1,v,t);
219: VecAXPY(q,pep->sfactor,t);
220: if (ctx->beta && ctx->alpha) {
221: STMatMult(pep->st,2,v,t);
222: VecAXPY(q,-pep->sfactor*pep->sfactor*ctx->beta/ctx->alpha,t);
223: }
224: STMatSolve(pep->st,q,t);
225: VecScale(t,-1.0/(pep->sfactor*pep->sfactor));
226: } else {
227: STMatMult(pep->st,1,v,q);
228: STMatMult(pep->st,2,v,t);
229: VecAXPY(q,sigma*pep->sfactor,t);
230: VecScale(q,pep->sfactor);
231: BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
232: STMatMult(pep->st,2,v,t);
233: VecAXPY(q,pep->sfactor*pep->sfactor,t);
234: STMatSolve(pep->st,q,t);
235: VecScale(t,-1.0);
236: }
237: BVRestoreColumn(pep->V,nqt,&t);
239: /* orthogonalize */
240: if (!sinvert) x = S+offq+(j+1)*lds;
241: else x = S+(j+1)*lds;
242: BVOrthogonalizeColumn(pep->V,nqt,x,&norm,&lindep);
244: if (!lindep) {
245: if (!sinvert) *(S+offq+(j+1)*lds+nqt) = norm;
246: else *(S+(j+1)*lds+nqt) = norm;
247: BVScaleColumn(pep->V,nqt,1.0/norm);
248: nqt++;
249: }
250: if (!sinvert) {
251: for (i=0;i<=nqt-1;i++) *(S+(j+1)*lds+i) = *(S+offq+j*lds+i);
252: if (ctx->beta && ctx->alpha) {
253: for (i=0;i<=nqt-1;i++) *(S+(j+1)*lds+offq+i) -= *(S+(j+1)*lds+i)*ctx->beta/ctx->alpha;
254: }
255: } else for (i=0;i<nqt;i++) *(S+(j+1)*lds+offq+i) = *(S+j*lds+i)+sigma*(*(S+(j+1)*lds+i));
256: BVSetActiveColumns(pep->V,0,nqt);
257: MatDenseRestoreArray(MS,&S);
258: BVTensorRestoreFactors(ctx->V,NULL,&MS);
260: /* level-2 orthogonalization */
261: BVOrthogonalizeColumn(ctx->V,j+1,y,&norm,&lindep);
262: a[j] = PetscRealPart(y[j]);
263: omega[j+1] = (norm > 0)?1.0:-1.0;
264: BVScaleColumn(ctx->V,j+1,1.0/norm);
265: b[j] = PetscAbsReal(norm);
267: /* check symmetry */
268: DSGetArrayReal(pep->ds,DS_MAT_T,&f);
269: if (j==k) {
270: for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ldds+i]);
271: for (i=0;i<l;i++) y[i] = 0.0;
272: }
273: DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
274: if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
275: PetscBLASIntCast(j,&j_);
276: sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
277: fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
278: if (j>0) fro = SlepcAbs(fro,b[j-1]);
279: if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
280: *symmlost = PETSC_TRUE;
281: *M=j;
282: break;
283: }
284: }
285: BVSetActiveColumns(pep->V,lock,nqt);
286: BVSetActiveColumns(ctx->V,0,*M);
287: PetscFree(y);
288: return(0);
289: }
291: #if 0
292: static PetscErrorCode PEPSTOARpreKConvergence(PEP pep,PetscInt nv,PetscReal *norm,Vec *w)
293: {
295: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
296: PetscBLASInt n_,one=1;
297: PetscInt lds=2*ctx->ld;
298: PetscReal t1,t2;
299: PetscScalar *S=ctx->S;
302: PetscBLASIntCast(nv+2,&n_);
303: t1 = BLASnrm2_(&n_,S+nv*2*ctx->ld,&one);
304: t2 = BLASnrm2_(&n_,S+(nv*2+1)*ctx->ld,&one);
305: *norm = SlepcAbs(t1,t2);
306: BVSetActiveColumns(pep->V,0,nv+2);
307: BVMultVec(pep->V,1.0,0.0,w[1],S+nv*lds);
308: STMatMult(pep->st,0,w[1],w[2]);
309: VecNorm(w[2],NORM_2,&t1);
310: BVMultVec(pep->V,1.0,0.0,w[1],S+ctx->ld+nv*lds);
311: STMatMult(pep->st,2,w[1],w[2]);
312: VecNorm(w[2],NORM_2,&t2);
313: t2 *= pep->sfactor*pep->sfactor;
314: *norm = PetscMax(*norm,SlepcAbs(t1,t2));
315: return(0);
316: }
317: #endif
319: PetscErrorCode PEPSolve_STOAR(PEP pep)
320: {
322: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
323: PetscInt j,k,l,nv=0,ld,ldds,t,nq=0;
324: PetscInt nconv=0,deg=pep->nmat-1;
325: PetscScalar *om,sigma;
326: PetscReal beta,norm=1.0,*omega,*a,*b;
327: PetscBool breakdown,symmlost=PETSC_FALSE,sinv=PETSC_FALSE,falselock=PETSC_TRUE,flg;
328: Mat MQ,A;
329: Vec vomega;
332: PetscCitationsRegister(citation,&cited);
333: PEPSTOARSetUpInnerMatrix(pep,&A);
334: BVSetMatrix(ctx->V,A,PETSC_TRUE);
335: MatDestroy(&A);
336: if (ctx->lock) {
337: /* undocumented option to use a cheaper locking instead of the true locking */
338: PetscOptionsGetBool(NULL,NULL,"-pep_stoar_falselocking",&falselock,NULL);
339: }
340: BVGetSizes(pep->V,NULL,NULL,&ld);
341: STGetShift(pep->st,&sigma);
342: STGetTransform(pep->st,&flg);
343: if (pep->sfactor!=1.0) {
344: if (!flg) {
345: pep->target /= pep->sfactor;
346: RGPushScale(pep->rg,1.0/pep->sfactor);
347: STScaleShift(pep->st,1.0/pep->sfactor);
348: sigma /= pep->sfactor;
349: } else {
350: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
351: pep->target = sinv?pep->target*pep->sfactor:pep->target/pep->sfactor;
352: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
353: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
354: }
355: }
356: if (flg) sigma = 0.0;
358: /* Get the starting Arnoldi vector */
359: BVTensorBuildFirstColumn(ctx->V,pep->nini);
360: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
361: VecCreateSeq(PETSC_COMM_SELF,1,&vomega);
362: BVSetActiveColumns(ctx->V,0,1);
363: BVGetSignature(ctx->V,vomega);
364: VecGetArray(vomega,&om);
365: omega[0] = PetscRealPart(om[0]);
366: VecRestoreArray(vomega,&om);
367: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
368: VecDestroy(&vomega);
370: /* Restart loop */
371: l = 0;
372: DSGetLeadingDimension(pep->ds,&ldds);
373: while (pep->reason == PEP_CONVERGED_ITERATING) {
374: pep->its++;
375: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
376: b = a+ldds;
377: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
379: /* Compute an nv-step Lanczos factorization */
380: nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
381: PEPSTOARrun(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,pep->work);
382: beta = b[nv-1];
383: if (symmlost && nv==pep->nconv+l) {
384: pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
385: pep->nconv = nconv;
386: if (falselock || !ctx->lock) {
387: BVSetActiveColumns(ctx->V,0,pep->nconv);
388: BVTensorCompress(ctx->V,0);
389: }
390: break;
391: }
392: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
393: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
394: DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l);
395: if (l==0) {
396: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
397: } else {
398: DSSetState(pep->ds,DS_STATE_RAW);
399: }
401: /* Solve projected problem */
402: DSSolve(pep->ds,pep->eigr,pep->eigi);
403: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
404: DSUpdateExtraRow(pep->ds);
405: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
407: /* Check convergence */
408: /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
409: norm = 1.0;
410: DSGetDimensions(pep->ds,NULL,NULL,NULL,&t);
411: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
412: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
414: /* Update l */
415: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
416: else {
417: l = PetscMax(1,(PetscInt)((nv-k)/2));
418: l = PetscMin(l,t);
419: DSGetTruncateSize(pep->ds,k,t,&l);
420: if (!breakdown) {
421: /* Prepare the Rayleigh quotient for restart */
422: DSTruncate(pep->ds,k+l,PETSC_FALSE);
423: }
424: }
425: nconv = k;
426: if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
427: if (l) { PetscInfo1(pep,"Preparing to restart keeping l=%D vectors\n",l); }
429: /* Update S */
430: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
431: BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
432: MatDestroy(&MQ);
434: /* Copy last column of S */
435: BVCopyColumn(ctx->V,nv,k+l);
436: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
437: VecCreateSeq(PETSC_COMM_SELF,k+l,&vomega);
438: VecGetArray(vomega,&om);
439: for (j=0;j<k+l;j++) om[j] = omega[j];
440: VecRestoreArray(vomega,&om);
441: BVSetActiveColumns(ctx->V,0,k+l);
442: BVSetSignature(ctx->V,vomega);
443: VecDestroy(&vomega);
444: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
446: if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
447: /* stop if breakdown */
448: PetscInfo2(pep,"Breakdown TOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
449: pep->reason = PEP_DIVERGED_BREAKDOWN;
450: }
451: if (pep->reason != PEP_CONVERGED_ITERATING) l--;
452: BVGetActiveColumns(pep->V,NULL,&nq);
453: if (k+l+deg<=nq) {
454: BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
455: if (!falselock && ctx->lock) {
456: BVTensorCompress(ctx->V,k-pep->nconv);
457: } else {
458: BVTensorCompress(ctx->V,0);
459: }
460: }
461: pep->nconv = k;
462: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
463: }
465: if (pep->nconv>0) {
466: BVSetActiveColumns(ctx->V,0,pep->nconv);
467: BVGetActiveColumns(pep->V,NULL,&nq);
468: BVSetActiveColumns(pep->V,0,nq);
469: if (nq>pep->nconv) {
470: BVTensorCompress(ctx->V,pep->nconv);
471: BVSetActiveColumns(pep->V,0,pep->nconv);
472: }
473: }
474: STGetTransform(pep->st,&flg);
475: if (!flg && pep->ops->backtransform) {
476: (*pep->ops->backtransform)(pep);
477: }
478: if (pep->sfactor!=1.0) {
479: for (j=0;j<pep->nconv;j++) {
480: pep->eigr[j] *= pep->sfactor;
481: pep->eigi[j] *= pep->sfactor;
482: }
483: }
484: /* restore original values */
485: if (!flg) {
486: pep->target *= pep->sfactor;
487: STScaleShift(pep->st,pep->sfactor);
488: } else {
489: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
490: pep->target = (sinv)?pep->target/pep->sfactor:pep->target*pep->sfactor;
491: }
492: if (pep->sfactor!=1.0) { RGPopScale(pep->rg); }
494: DSTruncate(pep->ds,pep->nconv,PETSC_TRUE);
495: return(0);
496: }
498: PetscErrorCode PEPSetFromOptions_STOAR(PetscOptionItems *PetscOptionsObject,PEP pep)
499: {
501: PetscBool flg,lock,b,f1,f2,f3;
502: PetscInt i,j,k;
503: PetscReal array[2]={0,0};
504: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
507: PetscOptionsHead(PetscOptionsObject,"PEP STOAR Options");
509: PetscOptionsBool("-pep_stoar_locking","Choose between locking and non-locking variants","PEPSTOARSetLocking",PETSC_FALSE,&lock,&flg);
510: if (flg) { PEPSTOARSetLocking(pep,lock); }
512: b = ctx->detect;
513: PetscOptionsBool("-pep_stoar_detect_zeros","Check zeros during factorizations at interval boundaries","PEPSTOARSetDetectZeros",ctx->detect,&b,&flg);
514: if (flg) { PEPSTOARSetDetectZeros(pep,b); }
516: i = 1;
517: j = k = PETSC_DECIDE;
518: PetscOptionsInt("-pep_stoar_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","PEPSTOARSetDimensions",20,&i,&f1);
519: PetscOptionsInt("-pep_stoar_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","PEPSTOARSetDimensions",40,&j,&f2);
520: PetscOptionsInt("-pep_stoar_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","PEPSTOARSetDimensions",40,&k,&f3);
521: if (f1 || f2 || f3) { PEPSTOARSetDimensions(pep,i,j,k); }
523: k = 2;
524: PetscOptionsRealArray("-pep_stoar_linearization","Parameters of the linearization","PEPSTOARSetLinearization",array,&k,&flg);
525: if (flg) {
526: PEPSTOARSetLinearization(pep,array[0],array[1]);
527: }
529: b = ctx->checket;
530: PetscOptionsBool("-pep_stoar_check_eigenvalue_type","Check eigenvalue type during spectrum slicing","PEPSTOARSetCheckEigenvalueType",ctx->checket,&b,&flg);
531: if (flg) { PEPSTOARSetCheckEigenvalueType(pep,b); }
533: PetscOptionsTail();
534: return(0);
535: }
537: static PetscErrorCode PEPSTOARSetLocking_STOAR(PEP pep,PetscBool lock)
538: {
539: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
542: ctx->lock = lock;
543: return(0);
544: }
546: /*@
547: PEPSTOARSetLocking - Choose between locking and non-locking variants of
548: the STOAR method.
550: Logically Collective on pep
552: Input Parameters:
553: + pep - the eigenproblem solver context
554: - lock - true if the locking variant must be selected
556: Options Database Key:
557: . -pep_stoar_locking - Sets the locking flag
559: Notes:
560: The default is to lock converged eigenpairs when the method restarts.
561: This behaviour can be changed so that all directions are kept in the
562: working subspace even if already converged to working accuracy (the
563: non-locking variant).
565: Level: advanced
567: .seealso: PEPSTOARGetLocking()
568: @*/
569: PetscErrorCode PEPSTOARSetLocking(PEP pep,PetscBool lock)
570: {
576: PetscTryMethod(pep,"PEPSTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
577: return(0);
578: }
580: static PetscErrorCode PEPSTOARGetLocking_STOAR(PEP pep,PetscBool *lock)
581: {
582: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
585: *lock = ctx->lock;
586: return(0);
587: }
589: /*@
590: PEPSTOARGetLocking - Gets the locking flag used in the STOAR method.
592: Not Collective
594: Input Parameter:
595: . pep - the eigenproblem solver context
597: Output Parameter:
598: . lock - the locking flag
600: Level: advanced
602: .seealso: PEPSTOARSetLocking()
603: @*/
604: PetscErrorCode PEPSTOARGetLocking(PEP pep,PetscBool *lock)
605: {
611: PetscUseMethod(pep,"PEPSTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
612: return(0);
613: }
615: static PetscErrorCode PEPSTOARGetInertias_STOAR(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
616: {
618: PetscInt i,numsh;
619: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
620: PEP_SR sr = ctx->sr;
623: if (!pep->state) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSetUp() first");
624: if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see PEPSetInterval()");
625: switch (pep->state) {
626: case PEP_STATE_INITIAL:
627: break;
628: case PEP_STATE_SETUP:
629: if (n) *n = 2;
630: if (shifts) {
631: PetscMalloc1(2,shifts);
632: (*shifts)[0] = pep->inta;
633: (*shifts)[1] = pep->intb;
634: }
635: if (inertias) {
636: PetscMalloc1(2,inertias);
637: (*inertias)[0] = (sr->dir==1)?sr->inertia0:sr->inertia1;
638: (*inertias)[1] = (sr->dir==1)?sr->inertia1:sr->inertia0;
639: }
640: break;
641: case PEP_STATE_SOLVED:
642: case PEP_STATE_EIGENVECTORS:
643: numsh = ctx->nshifts;
644: if (n) *n = numsh;
645: if (shifts) {
646: PetscMalloc1(numsh,shifts);
647: for (i=0;i<numsh;i++) (*shifts)[i] = ctx->shifts[i];
648: }
649: if (inertias) {
650: PetscMalloc1(numsh,inertias);
651: for (i=0;i<numsh;i++) (*inertias)[i] = ctx->inertias[i];
652: }
653: break;
654: }
655: return(0);
656: }
658: /*@C
659: PEPSTOARGetInertias - Gets the values of the shifts and their
660: corresponding inertias in case of doing spectrum slicing for a
661: computational interval.
663: Not Collective
665: Input Parameter:
666: . pep - the eigenproblem solver context
668: Output Parameters:
669: + n - number of shifts, including the endpoints of the interval
670: . shifts - the values of the shifts used internally in the solver
671: - inertias - the values of the inertia in each shift
673: Notes:
674: If called after PEPSolve(), all shifts used internally by the solver are
675: returned (including both endpoints and any intermediate ones). If called
676: before PEPSolve() and after PEPSetUp() then only the information of the
677: endpoints of subintervals is available.
679: This function is only available for spectrum slicing runs.
681: The returned arrays should be freed by the user. Can pass NULL in any of
682: the two arrays if not required.
684: Fortran Notes:
685: The calling sequence from Fortran is
686: .vb
687: PEPSTOARGetInertias(pep,n,shifts,inertias,ierr)
688: integer n
689: double precision shifts(*)
690: integer inertias(*)
691: .ve
692: The arrays should be at least of length n. The value of n can be determined
693: by an initial call
694: .vb
695: PEPSTOARGetInertias(pep,n,PETSC_NULL_REAL,PETSC_NULL_INTEGER,ierr)
696: .ve
698: Level: advanced
700: .seealso: PEPSetInterval()
701: @*/
702: PetscErrorCode PEPSTOARGetInertias(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
703: {
709: PetscUseMethod(pep,"PEPSTOARGetInertias_C",(PEP,PetscInt*,PetscReal**,PetscInt**),(pep,n,shifts,inertias));
710: return(0);
711: }
713: static PetscErrorCode PEPSTOARSetDetectZeros_STOAR(PEP pep,PetscBool detect)
714: {
715: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
718: ctx->detect = detect;
719: pep->state = PEP_STATE_INITIAL;
720: return(0);
721: }
723: /*@
724: PEPSTOARSetDetectZeros - Sets a flag to enforce detection of
725: zeros during the factorizations throughout the spectrum slicing computation.
727: Logically Collective on pep
729: Input Parameters:
730: + pep - the eigenproblem solver context
731: - detect - check for zeros
733: Options Database Key:
734: . -pep_stoar_detect_zeros - Check for zeros; this takes an optional
735: bool value (0/1/no/yes/true/false)
737: Notes:
738: A zero in the factorization indicates that a shift coincides with an eigenvalue.
740: This flag is turned off by default, and may be necessary in some cases.
741: This feature currently requires an external package for factorizations
742: with support for zero detection, e.g. MUMPS.
744: Level: advanced
746: .seealso: PEPSetInterval()
747: @*/
748: PetscErrorCode PEPSTOARSetDetectZeros(PEP pep,PetscBool detect)
749: {
755: PetscTryMethod(pep,"PEPSTOARSetDetectZeros_C",(PEP,PetscBool),(pep,detect));
756: return(0);
757: }
759: static PetscErrorCode PEPSTOARGetDetectZeros_STOAR(PEP pep,PetscBool *detect)
760: {
761: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
764: *detect = ctx->detect;
765: return(0);
766: }
768: /*@
769: PEPSTOARGetDetectZeros - Gets the flag that enforces zero detection
770: in spectrum slicing.
772: Not Collective
774: Input Parameter:
775: . pep - the eigenproblem solver context
777: Output Parameter:
778: . detect - whether zeros detection is enforced during factorizations
780: Level: advanced
782: .seealso: PEPSTOARSetDetectZeros()
783: @*/
784: PetscErrorCode PEPSTOARGetDetectZeros(PEP pep,PetscBool *detect)
785: {
791: PetscUseMethod(pep,"PEPSTOARGetDetectZeros_C",(PEP,PetscBool*),(pep,detect));
792: return(0);
793: }
795: static PetscErrorCode PEPSTOARSetLinearization_STOAR(PEP pep,PetscReal alpha,PetscReal beta)
796: {
797: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
800: if (beta==0.0 && alpha==0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Parameters alpha and beta cannot be zero simultaneously");
801: ctx->alpha = alpha;
802: ctx->beta = beta;
803: return(0);
804: }
806: /*@
807: PEPSTOARSetLinearization - Set the coefficients that define
808: the linearization of a quadratic eigenproblem.
810: Logically Collective on pep
812: Input Parameters:
813: + pep - polynomial eigenvalue solver
814: . alpha - first parameter of the linearization
815: - beta - second parameter of the linearization
817: Options Database Key:
818: . -pep_stoar_linearization <alpha,beta> - Sets the coefficients
820: Notes:
821: Cannot pass zero for both alpha and beta. The default values are
822: alpha=1 and beta=0.
824: Level: advanced
826: .seealso: PEPSTOARGetLinearization()
827: @*/
828: PetscErrorCode PEPSTOARSetLinearization(PEP pep,PetscReal alpha,PetscReal beta)
829: {
836: PetscTryMethod(pep,"PEPSTOARSetLinearization_C",(PEP,PetscReal,PetscReal),(pep,alpha,beta));
837: return(0);
838: }
840: static PetscErrorCode PEPSTOARGetLinearization_STOAR(PEP pep,PetscReal *alpha,PetscReal *beta)
841: {
842: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
845: if (alpha) *alpha = ctx->alpha;
846: if (beta) *beta = ctx->beta;
847: return(0);
848: }
850: /*@
851: PEPSTOARGetLinearization - Returns the coefficients that define
852: the linearization of a quadratic eigenproblem.
854: Not Collective
856: Input Parameter:
857: . pep - polynomial eigenvalue solver
859: Output Parameters:
860: + alpha - the first parameter of the linearization
861: - beta - the second parameter of the linearization
863: Level: advanced
865: .seealso: PEPSTOARSetLinearization()
866: @*/
867: PetscErrorCode PEPSTOARGetLinearization(PEP pep,PetscReal *alpha,PetscReal *beta)
868: {
873: PetscUseMethod(pep,"PEPSTOARGetLinearization_C",(PEP,PetscReal*,PetscReal*),(pep,alpha,beta));
874: return(0);
875: }
877: static PetscErrorCode PEPSTOARSetDimensions_STOAR(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
878: {
879: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
882: if (nev<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
883: ctx->nev = nev;
884: if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
885: ctx->ncv = PETSC_DEFAULT;
886: } else {
887: if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
888: ctx->ncv = ncv;
889: }
890: if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
891: ctx->mpd = PETSC_DEFAULT;
892: } else {
893: if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
894: ctx->mpd = mpd;
895: }
896: pep->state = PEP_STATE_INITIAL;
897: return(0);
898: }
900: /*@
901: PEPSTOARSetDimensions - Sets the dimensions used for each subsolve
902: step in case of doing spectrum slicing for a computational interval.
903: The meaning of the parameters is the same as in PEPSetDimensions().
905: Logically Collective on pep
907: Input Parameters:
908: + pep - the eigenproblem solver context
909: . nev - number of eigenvalues to compute
910: . ncv - the maximum dimension of the subspace to be used by the subsolve
911: - mpd - the maximum dimension allowed for the projected problem
913: Options Database Key:
914: + -eps_stoar_nev <nev> - Sets the number of eigenvalues
915: . -eps_stoar_ncv <ncv> - Sets the dimension of the subspace
916: - -eps_stoar_mpd <mpd> - Sets the maximum projected dimension
918: Level: advanced
920: .seealso: PEPSTOARGetDimensions(), PEPSetDimensions(), PEPSetInterval()
921: @*/
922: PetscErrorCode PEPSTOARSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
923: {
931: PetscTryMethod(pep,"PEPSTOARSetDimensions_C",(PEP,PetscInt,PetscInt,PetscInt),(pep,nev,ncv,mpd));
932: return(0);
933: }
935: static PetscErrorCode PEPSTOARGetDimensions_STOAR(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
936: {
937: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
940: if (nev) *nev = ctx->nev;
941: if (ncv) *ncv = ctx->ncv;
942: if (mpd) *mpd = ctx->mpd;
943: return(0);
944: }
946: /*@
947: PEPSTOARGetDimensions - Gets the dimensions used for each subsolve
948: step in case of doing spectrum slicing for a computational interval.
950: Not Collective
952: Input Parameter:
953: . pep - the eigenproblem solver context
955: Output Parameters:
956: + nev - number of eigenvalues to compute
957: . ncv - the maximum dimension of the subspace to be used by the subsolve
958: - mpd - the maximum dimension allowed for the projected problem
960: Level: advanced
962: .seealso: PEPSTOARSetDimensions()
963: @*/
964: PetscErrorCode PEPSTOARGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
965: {
970: PetscUseMethod(pep,"PEPSTOARGetDimensions_C",(PEP,PetscInt*,PetscInt*,PetscInt*),(pep,nev,ncv,mpd));
971: return(0);
972: }
974: static PetscErrorCode PEPSTOARSetCheckEigenvalueType_STOAR(PEP pep,PetscBool checket)
975: {
976: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
979: ctx->checket = checket;
980: pep->state = PEP_STATE_INITIAL;
981: return(0);
982: }
984: /*@
985: PEPSTOARSetCheckEigenvalueType - Sets a flag to check that all the eigenvalues
986: obtained throughout the spectrum slicing computation have the same definite type.
988: Logically Collective on pep
990: Input Parameters:
991: + pep - the eigenproblem solver context
992: - checket - check eigenvalue type
994: Options Database Key:
995: . -pep_stoar_check_eigenvalue_type - Check eigenvalue type; this takes an optional
996: bool value (0/1/no/yes/true/false)
998: Notes:
999: This option is relevant only for spectrum slicing computations, but it is
1000: ignored if the problem type is PEP_HYPERBOLIC.
1002: This flag is turned on by default, to guarantee that the computed eigenvalues
1003: have the same type (otherwise the computed solution might be wrong). But since
1004: the check is computationally quite expensive, the check may be turned off if
1005: the user knows for sure that all eigenvalues in the requested interval have
1006: the same type.
1008: Level: advanced
1010: .seealso: PEPSetProblemType(), PEPSetInterval()
1011: @*/
1012: PetscErrorCode PEPSTOARSetCheckEigenvalueType(PEP pep,PetscBool checket)
1013: {
1019: PetscTryMethod(pep,"PEPSTOARSetCheckEigenvalueType_C",(PEP,PetscBool),(pep,checket));
1020: return(0);
1021: }
1023: static PetscErrorCode PEPSTOARGetCheckEigenvalueType_STOAR(PEP pep,PetscBool *checket)
1024: {
1025: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
1028: *checket = ctx->checket;
1029: return(0);
1030: }
1032: /*@
1033: PEPSTOARGetCheckEigenvalueType - Gets the flag for the eigenvalue type
1034: check in spectrum slicing.
1036: Not Collective
1038: Input Parameter:
1039: . pep - the eigenproblem solver context
1041: Output Parameter:
1042: . checket - whether eigenvalue type must be checked during spectrum slcing
1044: Level: advanced
1046: .seealso: PEPSTOARSetCheckEigenvalueType()
1047: @*/
1048: PetscErrorCode PEPSTOARGetCheckEigenvalueType(PEP pep,PetscBool *checket)
1049: {
1055: PetscUseMethod(pep,"PEPSTOARGetCheckEigenvalueType_C",(PEP,PetscBool*),(pep,checket));
1056: return(0);
1057: }
1059: PetscErrorCode PEPView_STOAR(PEP pep,PetscViewer viewer)
1060: {
1062: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
1063: PetscBool isascii;
1066: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1067: if (isascii) {
1068: PetscViewerASCIIPrintf(viewer," using the %slocking variant\n",ctx->lock?"":"non-");
1069: PetscViewerASCIIPrintf(viewer," linearization parameters: alpha=%g beta=%g\n",(double)ctx->alpha,(double)ctx->beta);
1070: if (pep->which==PEP_ALL && !ctx->hyperbolic) {
1071: PetscViewerASCIIPrintf(viewer," checking eigenvalue type: %s\n",ctx->checket?"enabled":"disabled");
1072: }
1073: }
1074: return(0);
1075: }
1077: PetscErrorCode PEPReset_STOAR(PEP pep)
1078: {
1082: if (pep->which==PEP_ALL) {
1083: PEPReset_STOAR_QSlice(pep);
1084: }
1085: return(0);
1086: }
1088: PetscErrorCode PEPDestroy_STOAR(PEP pep)
1089: {
1091: PEP_STOAR *ctx = (PEP_STOAR*)pep->data;
1094: BVDestroy(&ctx->V);
1095: PetscFree(pep->data);
1096: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",NULL);
1097: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",NULL);
1098: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDetectZeros_C",NULL);
1099: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDetectZeros_C",NULL);
1100: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetInertias_C",NULL);
1101: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDimensions_C",NULL);
1102: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDimensions_C",NULL);
1103: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLinearization_C",NULL);
1104: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLinearization_C",NULL);
1105: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetCheckEigenvalueType_C",NULL);
1106: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetCheckEigenvalueType_C",NULL);
1107: return(0);
1108: }
1110: SLEPC_EXTERN PetscErrorCode PEPCreate_STOAR(PEP pep)
1111: {
1113: PEP_STOAR *ctx;
1116: PetscNewLog(pep,&ctx);
1117: pep->data = (void*)ctx;
1119: pep->lineariz = PETSC_TRUE;
1120: ctx->lock = PETSC_TRUE;
1121: ctx->nev = 1;
1122: ctx->ncv = PETSC_DEFAULT;
1123: ctx->mpd = PETSC_DEFAULT;
1124: ctx->alpha = 1.0;
1125: ctx->beta = 0.0;
1126: ctx->checket = PETSC_TRUE;
1128: pep->ops->setup = PEPSetUp_STOAR;
1129: pep->ops->setfromoptions = PEPSetFromOptions_STOAR;
1130: pep->ops->destroy = PEPDestroy_STOAR;
1131: pep->ops->view = PEPView_STOAR;
1132: pep->ops->backtransform = PEPBackTransform_Default;
1133: pep->ops->computevectors = PEPComputeVectors_Default;
1134: pep->ops->extractvectors = PEPExtractVectors_TOAR;
1135: pep->ops->reset = PEPReset_STOAR;
1137: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",PEPSTOARSetLocking_STOAR);
1138: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",PEPSTOARGetLocking_STOAR);
1139: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDetectZeros_C",PEPSTOARSetDetectZeros_STOAR);
1140: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDetectZeros_C",PEPSTOARGetDetectZeros_STOAR);
1141: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetInertias_C",PEPSTOARGetInertias_STOAR);
1142: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDimensions_C",PEPSTOARGetDimensions_STOAR);
1143: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDimensions_C",PEPSTOARSetDimensions_STOAR);
1144: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLinearization_C",PEPSTOARSetLinearization_STOAR);
1145: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLinearization_C",PEPSTOARGetLinearization_STOAR);
1146: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetCheckEigenvalueType_C",PEPSTOARSetCheckEigenvalueType_STOAR);
1147: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetCheckEigenvalueType_C",PEPSTOARGetCheckEigenvalueType_STOAR);
1148: return(0);
1149: }