Actual source code: stoar.c
slepc-3.10.1 2018-10-23
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
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> /*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-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";
44: PetscErrorCode MatMult_Func(Mat A,Vec x,Vec y)
45: {
47: ShellMatCtx *ctx;
50: MatShellGetContext(A,&ctx);
51: MatMult(ctx->A[0],x,y);
52: VecScale(y,ctx->scal[0]);
53: if (ctx->scal[1]) {
54: MatMult(ctx->A[1],x,ctx->t);
55: VecAXPY(y,ctx->scal[1],ctx->t);
56: }
57: return(0);
58: }
60: PetscErrorCode MatDestroy_Func(Mat A)
61: {
62: ShellMatCtx *ctx;
66: MatShellGetContext(A,(void**)&ctx);
67: VecDestroy(&ctx->t);
68: PetscFree(ctx);
69: return(0);
70: }
72: PetscErrorCode PEPSTOARSetUpInnerMatrix(PEP pep,Mat *B)
73: {
74: Mat pB[4],Bs[3],D[3];
75: PetscInt i,j,n,m;
76: ShellMatCtx *ctxMat[3];
77: PEP_TOAR *ctx=(PEP_TOAR*)pep->data;
81: for (i=0;i<3;i++) {
82: STGetMatrixTransformed(pep->st,i,&D[i]); /* D[2] = M */
83: }
84: MatGetLocalSize(D[2],&m,&n);
85:
86: for (j=0;j<3;j++) {
87: PetscNew(ctxMat+j);
88: MatCreateShell(PetscObjectComm((PetscObject)pep),m,n,PETSC_DETERMINE,PETSC_DETERMINE,ctxMat[j],&Bs[j]);
89: MatShellSetOperation(Bs[j],MATOP_MULT,(void(*)())MatMult_Func);
90: MatShellSetOperation(Bs[j],MATOP_DESTROY,(void(*)())MatDestroy_Func);
91: }
92: for (i=0;i<4;i++) pB[i] = NULL;
93: if (ctx->alpha) {
94: ctxMat[0]->A[0] = D[0]; ctxMat[0]->scal[0] = ctx->alpha; ctxMat[0]->scal[1] = 0.0;
95: ctxMat[2]->A[0] = D[2]; ctxMat[2]->scal[0] = -ctx->alpha*pep->sfactor*pep->sfactor; ctxMat[2]->scal[1] = 0.0;
96: pB[0] = Bs[0]; pB[3] = Bs[2];
97: }
98: if (ctx->beta) {
99: i = (ctx->alpha)?1:0;
100: ctxMat[0]->scal[1] = 0.0;
101: ctxMat[0]->A[i] = D[1]; ctxMat[0]->scal[i] = -ctx->beta*pep->sfactor;
102: ctxMat[1]->A[0] = D[2]; ctxMat[1]->scal[0] = -ctx->beta*pep->sfactor*pep->sfactor; ctxMat[1]->scal[1] = 0.0;
103: pB[0] = Bs[0]; pB[1] = pB[2] = Bs[1];
104: }
105: BVCreateVec(pep->V,&ctxMat[0]->t);
106: MatCreateNest(PetscObjectComm((PetscObject)pep),2,NULL,2,NULL,pB,B);
107: for (j=0;j<3;j++) { MatDestroy(&Bs[j]); }
108: return(0);
109: }
111: PetscErrorCode PEPSetUp_STOAR(PEP pep)
112: {
113: PetscErrorCode ierr;
114: PetscBool shift,sinv,flg;
115: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
116: PetscInt ld;
117: PetscReal eta;
118: BVOrthogType otype;
119: BVOrthogBlockType obtype;
122: if (pep->problem_type!=PEP_HERMITIAN && pep->problem_type!=PEP_HYPERBOLIC) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Requested method is only available for Hermitian (or hyperbolic) problems");
123: if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver only available for quadratic problems");
124: if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
125: /* spectrum slicing requires special treatment of default values */
126: if (pep->which==PEP_ALL) {
127: if (pep->stopping!=PEPStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Spectrum slicing does not support user-defined stopping test");
128: pep->ops->solve = PEPSolve_STOAR_QSlice;
129: pep->ops->extractvectors = NULL;
130: pep->ops->setdefaultst = NULL;
131: PEPSetUp_STOAR_QSlice(pep);
132: } else {
133: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
134: if (!ctx->lock && pep->mpd<pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
135: if (!pep->max_it) pep->max_it = PetscMax(100,2*(pep->nmat-1)*pep->n/pep->ncv);
136: pep->ops->solve = PEPSolve_STOAR;
137: ld = pep->ncv+2;
138: DSSetType(pep->ds,DSGHIEP);
139: DSSetCompact(pep->ds,PETSC_TRUE);
140: DSAllocate(pep->ds,ld);
141: STGetTransform(pep->st,&flg);
142: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag set, see STSetTransform()");
143: }
145: pep->lineariz = PETSC_TRUE;
146: PetscObjectTypeCompare((PetscObject)pep->st,STSHIFT,&shift);
147: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
148: if (!shift && !sinv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Only STSHIFT and STSINVERT spectral transformations can be used");
149: if (!pep->which) {
150: if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
151: else pep->which = PEP_LARGEST_MAGNITUDE;
152: }
154: PEPAllocateSolution(pep,2);
155: PEPSetWorkVecs(pep,4);
156: BVDestroy(&ctx->V);
157: BVCreateTensor(pep->V,pep->nmat-1,&ctx->V);
158: BVGetOrthogonalization(pep->V,&otype,NULL,&eta,&obtype);
159: BVSetOrthogonalization(ctx->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
160: return(0);
161: }
163: /*
164: Compute a run of Lanczos iterations. dim(work)=(ctx->ld)*4
165: */
166: static PetscErrorCode PEPSTOARrun(PEP pep,PetscReal *a,PetscReal *b,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *breakdown,PetscBool *symmlost,Vec *t_)
167: {
169: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
170: PetscInt i,j,m=*M,l,lock;
171: PetscInt lds,d,ld,offq,nqt;
172: Vec v=t_[0],t=t_[1],q=t_[2];
173: PetscReal norm,sym=0.0,fro=0.0,*f;
174: PetscScalar *y,*S;
175: PetscBLASInt j_,one=1;
176: PetscBool lindep;
177: Mat MS;
180: PetscMalloc1(*M,&y);
181: BVGetSizes(pep->V,NULL,NULL,&ld);
182: BVTensorGetDegree(ctx->V,&d);
183: BVGetActiveColumns(pep->V,&lock,&nqt);
184: lds = d*ld;
185: offq = ld;
186: *breakdown = PETSC_FALSE; /* ----- */
187: DSGetDimensions(pep->ds,NULL,NULL,&l,NULL,NULL);
188: BVSetActiveColumns(ctx->V,0,m);
189: BVSetActiveColumns(pep->V,0,nqt);
190: for (j=k;j<m;j++) {
191: /* apply operator */
192: BVTensorGetFactors(ctx->V,NULL,&MS);
193: MatDenseGetArray(MS,&S);
194: BVGetColumn(pep->V,nqt,&t);
195: BVMultVec(pep->V,1.0,0.0,v,S+j*lds);
196: STMatMult(pep->st,0,v,q);
197: BVMultVec(pep->V,1.0,0.0,v,S+offq+j*lds);
198: STMatMult(pep->st,1,v,t);
199: VecAXPY(q,pep->sfactor,t);
200: if (ctx->beta && ctx->alpha) {
201: STMatMult(pep->st,2,v,t);
202: VecAXPY(q,-pep->sfactor*pep->sfactor*ctx->beta/ctx->alpha,t);
203: }
204: STMatSolve(pep->st,q,t);
205: VecScale(t,-1.0/(pep->sfactor*pep->sfactor));
206: BVRestoreColumn(pep->V,nqt,&t);
208: /* orthogonalize */
209: BVOrthogonalizeColumn(pep->V,nqt,S+offq+(j+1)*lds,&norm,&lindep);
210: if (!lindep) {
211: *(S+offq+(j+1)*lds+nqt) = norm;
212: BVScaleColumn(pep->V,nqt,1.0/norm);
213: nqt++;
214: }
215: for (i=0;i<=nqt-1;i++) *(S+(j+1)*lds+i) = *(S+offq+j*lds+i);
216: if (ctx->beta && ctx->alpha) {
217: for (i=0;i<=nqt-1;i++) *(S+(j+1)*lds+offq+i) -= *(S+(j+1)*lds+i)*ctx->beta/ctx->alpha;
218: }
219: BVSetActiveColumns(pep->V,0,nqt);
220: MatDenseRestoreArray(MS,&S);
221: BVTensorRestoreFactors(ctx->V,NULL,&MS);
223: /* level-2 orthogonalization */
224: BVOrthogonalizeColumn(ctx->V,j+1,y,&norm,&lindep);
225: a[j] = PetscRealPart(y[j]);
226: omega[j+1] = (norm > 0)?1.0:-1.0;
227: BVScaleColumn(ctx->V,j+1,1.0/norm);
228: b[j] = PetscAbsReal(norm);
230: /* check symmetry */
231: DSGetArrayReal(pep->ds,DS_MAT_T,&f);
232: if (j==k) {
233: for (i=l;i<j-1;i++) y[i] = PetscAbsScalar(y[i])-PetscAbsReal(f[2*ld+i]);
234: for (i=0;i<l;i++) y[i] = 0.0;
235: }
236: DSRestoreArrayReal(pep->ds,DS_MAT_T,&f);
237: if (j>0) y[j-1] = PetscAbsScalar(y[j-1])-PetscAbsReal(b[j-1]);
238: PetscBLASIntCast(j,&j_);
239: sym = SlepcAbs(BLASnrm2_(&j_,y,&one),sym);
240: fro = SlepcAbs(fro,SlepcAbs(a[j],b[j]));
241: if (j>0) fro = SlepcAbs(fro,b[j-1]);
242: if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,10*pep->tol)) {
243: *symmlost = PETSC_TRUE;
244: *M=j;
245: break;
246: }
247: }
248: BVSetActiveColumns(pep->V,lock,nqt);
249: BVSetActiveColumns(ctx->V,0,*M);
250: PetscFree(y);
251: return(0);
252: }
254: #if 0
255: static PetscErrorCode PEPSTOARpreKConvergence(PEP pep,PetscInt nv,PetscReal *norm,Vec *w)
256: {
258: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
259: PetscBLASInt n_,one=1;
260: PetscInt lds=2*ctx->ld;
261: PetscReal t1,t2;
262: PetscScalar *S=ctx->S;
265: PetscBLASIntCast(nv+2,&n_);
266: t1 = BLASnrm2_(&n_,S+nv*2*ctx->ld,&one);
267: t2 = BLASnrm2_(&n_,S+(nv*2+1)*ctx->ld,&one);
268: *norm = SlepcAbs(t1,t2);
269: BVSetActiveColumns(pep->V,0,nv+2);
270: BVMultVec(pep->V,1.0,0.0,w[1],S+nv*lds);
271: STMatMult(pep->st,0,w[1],w[2]);
272: VecNorm(w[2],NORM_2,&t1);
273: BVMultVec(pep->V,1.0,0.0,w[1],S+ctx->ld+nv*lds);
274: STMatMult(pep->st,2,w[1],w[2]);
275: VecNorm(w[2],NORM_2,&t2);
276: t2 *= pep->sfactor*pep->sfactor;
277: *norm = PetscMax(*norm,SlepcAbs(t1,t2));
278: return(0);
279: }
280: #endif
282: PetscErrorCode PEPSolve_STOAR(PEP pep)
283: {
285: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
286: PetscInt j,k,l,nv=0,ld,ldds,t,nq=0;
287: PetscInt nconv=0,deg=pep->nmat-1;
288: PetscScalar *Q,*om;
289: PetscReal beta,norm=1.0,*omega,*a,*b,*r;
290: PetscBool breakdown,symmlost=PETSC_FALSE,sinv,falselock=PETSC_TRUE;
291: Mat MQ,A;
292: Vec vomega;
295: PetscCitationsRegister(citation,&cited);
296: PEPSTOARSetUpInnerMatrix(pep,&A);
297: BVSetMatrix(ctx->V,A,PETSC_TRUE);
298: MatDestroy(&A);
299: if (ctx->lock) {
300: PetscOptionsGetBool(NULL,NULL,"-pep_stoar_falselocking",&falselock,NULL);
301: }
302: BVGetSizes(pep->V,NULL,NULL,&ld);
303: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
304: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
305: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
307: /* Get the starting Arnoldi vector */
308: BVTensorBuildFirstColumn(ctx->V,pep->nini);
309: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
310: VecCreateSeq(PETSC_COMM_SELF,1,&vomega);
311: BVSetActiveColumns(ctx->V,0,1);
312: BVGetSignature(ctx->V,vomega);
313: VecGetArray(vomega,&om);
314: omega[0] = PetscRealPart(om[0]);
315: VecRestoreArray(vomega,&om);
316: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
317: VecDestroy(&vomega);
319: /* Restart loop */
320: l = 0;
321: DSGetLeadingDimension(pep->ds,&ldds);
322: while (pep->reason == PEP_CONVERGED_ITERATING) {
323: pep->its++;
324: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
325: b = a+ldds;
326: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
328: /* Compute an nv-step Lanczos factorization */
329: nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
330: PEPSTOARrun(pep,a,b,omega,pep->nconv+l,&nv,&breakdown,&symmlost,pep->work);
331: beta = b[nv-1];
332: if (symmlost && nv==pep->nconv+l) {
333: pep->reason = PEP_DIVERGED_SYMMETRY_LOST;
334: pep->nconv = nconv;
335: if (falselock || !ctx->lock) {
336: BVSetActiveColumns(ctx->V,0,pep->nconv);
337: BVTensorCompress(ctx->V,0);
338: }
339: break;
340: }
341: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
342: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
343: DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
344: if (l==0) {
345: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
346: } else {
347: DSSetState(pep->ds,DS_STATE_RAW);
348: }
350: /* Solve projected problem */
351: DSSolve(pep->ds,pep->eigr,pep->eigi);
352: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
353: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
355: /* Check convergence */
356: /* PEPSTOARpreKConvergence(pep,nv,&norm,pep->work);*/
357: norm = 1.0;
358: DSGetDimensions(pep->ds,NULL,NULL,NULL,NULL,&t);
359: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,t-pep->nconv,PetscAbsReal(beta)*norm,&k);
360: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
362: /* Update l */
363: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
364: else {
365: l = PetscMax(1,(PetscInt)((nv-k)/2));
366: l = PetscMin(l,t);
367: if (!breakdown) {
368: DSGetArrayReal(pep->ds,DS_MAT_T,&a);
369: if (*(a+ldds+k+l-1)!=0) {
370: if (k+l<nv-1) l = l+1;
371: else l = l-1;
372: }
373: /* Prepare the Rayleigh quotient for restart */
374: DSGetArray(pep->ds,DS_MAT_Q,&Q);
375: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
376: r = a + 2*ldds;
377: for (j=k;j<k+l;j++) {
378: r[j] = PetscRealPart(Q[nv-1+j*ldds]*beta);
379: }
380: b = a+ldds;
381: b[k+l-1] = r[k+l-1];
382: omega[k+l] = omega[nv];
383: DSRestoreArray(pep->ds,DS_MAT_Q,&Q);
384: DSRestoreArrayReal(pep->ds,DS_MAT_T,&a);
385: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
386: }
387: }
388: nconv = k;
389: if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
391: /* Update S */
392: DSGetMat(pep->ds,DS_MAT_Q,&MQ);
393: BVMultInPlace(ctx->V,MQ,pep->nconv,k+l);
394: MatDestroy(&MQ);
396: /* Copy last column of S */
397: BVCopyColumn(ctx->V,nv,k+l);
398: DSGetArrayReal(pep->ds,DS_MAT_D,&omega);
399: VecCreateSeq(PETSC_COMM_SELF,k+l,&vomega);
400: VecGetArray(vomega,&om);
401: for (j=0;j<k+l;j++) om[j] = omega[j];
402: VecRestoreArray(vomega,&om);
403: BVSetActiveColumns(ctx->V,0,k+l);
404: BVSetSignature(ctx->V,vomega);
405: VecDestroy(&vomega);
406: DSRestoreArrayReal(pep->ds,DS_MAT_D,&omega);
408: if (breakdown && pep->reason == PEP_CONVERGED_ITERATING) {
409: /* stop if breakdown */
410: PetscInfo2(pep,"Breakdown TOAR method (it=%D norm=%g)\n",pep->its,(double)beta);
411: pep->reason = PEP_DIVERGED_BREAKDOWN;
412: }
413: if (pep->reason != PEP_CONVERGED_ITERATING) l--;
414: BVGetActiveColumns(pep->V,NULL,&nq);
415: if (k+l+deg<=nq) {
416: BVSetActiveColumns(ctx->V,pep->nconv,k+l+1);
417: if (!falselock && ctx->lock) {
418: BVTensorCompress(ctx->V,k-pep->nconv);
419: } else {
420: BVTensorCompress(ctx->V,0);
421: }
422: }
423: pep->nconv = k;
424: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
425: }
427: if (pep->nconv>0) {
428: BVSetActiveColumns(ctx->V,0,pep->nconv);
429: BVGetActiveColumns(pep->V,NULL,&nq);
430: BVSetActiveColumns(pep->V,0,nq);
431: if (nq>pep->nconv) {
432: BVTensorCompress(ctx->V,pep->nconv);
433: BVSetActiveColumns(pep->V,0,pep->nconv);
434: }
435: for (j=0;j<pep->nconv;j++) {
436: pep->eigr[j] *= pep->sfactor;
437: pep->eigi[j] *= pep->sfactor;
438: }
439: }
440: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
441: RGPopScale(pep->rg);
443: /* truncate Schur decomposition and change the state to raw so that
444: DSVectors() computes eigenvectors from scratch */
445: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
446: DSSetState(pep->ds,DS_STATE_RAW);
447: return(0);
448: }
450: PetscErrorCode PEPSetFromOptions_STOAR(PetscOptionItems *PetscOptionsObject,PEP pep)
451: {
453: PetscBool flg,lock,b,f1,f2,f3;
454: PetscInt i,j,k;
455: PetscReal array[2]={0,0};
456: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
459: PetscOptionsHead(PetscOptionsObject,"PEP STOAR Options");
461: PetscOptionsBool("-pep_stoar_locking","Choose between locking and non-locking variants","PEPSTOARSetLocking",PETSC_FALSE,&lock,&flg);
462: if (flg) { PEPSTOARSetLocking(pep,lock); }
464: b = ctx->detect;
465: PetscOptionsBool("-pep_stoar_detect_zeros","Check zeros during factorizations at interval boundaries","PEPSTOARSetDetectZeros",ctx->detect,&b,&flg);
466: if (flg) { PEPSTOARSetDetectZeros(pep,b); }
468: i = 1;
469: j = k = PETSC_DECIDE;
470: PetscOptionsInt("-pep_stoar_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","PEPSTOARSetDimensions",20,&i,&f1);
471: PetscOptionsInt("-pep_stoar_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","PEPSTOARSetDimensions",40,&j,&f2);
472: PetscOptionsInt("-pep_stoar_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","PEPSTOARSetDimensions",40,&k,&f3);
473: if (f1 || f2 || f3) { PEPSTOARSetDimensions(pep,i,j,k); }
475: k = 2;
476: PetscOptionsRealArray("-pep_stoar_linearization","Parameters of the linearization","PEPSTOARSetLinearization",array,&k,&flg);
477: if (flg) {
478: PEPSTOARSetLinearization(pep,array[0],array[1]);
479: }
481: PetscOptionsTail();
482: return(0);
483: }
485: static PetscErrorCode PEPSTOARSetLocking_STOAR(PEP pep,PetscBool lock)
486: {
487: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
490: ctx->lock = lock;
491: return(0);
492: }
494: /*@
495: PEPSTOARSetLocking - Choose between locking and non-locking variants of
496: the STOAR method.
498: Logically Collective on PEP
500: Input Parameters:
501: + pep - the eigenproblem solver context
502: - lock - true if the locking variant must be selected
504: Options Database Key:
505: . -pep_stoar_locking - Sets the locking flag
507: Notes:
508: The default is to lock converged eigenpairs when the method restarts.
509: This behaviour can be changed so that all directions are kept in the
510: working subspace even if already converged to working accuracy (the
511: non-locking variant).
513: Level: advanced
515: .seealso: PEPSTOARGetLocking()
516: @*/
517: PetscErrorCode PEPSTOARSetLocking(PEP pep,PetscBool lock)
518: {
524: PetscTryMethod(pep,"PEPSTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
525: return(0);
526: }
528: static PetscErrorCode PEPSTOARGetLocking_STOAR(PEP pep,PetscBool *lock)
529: {
530: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
533: *lock = ctx->lock;
534: return(0);
535: }
537: /*@
538: PEPSTOARGetLocking - Gets the locking flag used in the STOAR method.
540: Not Collective
542: Input Parameter:
543: . pep - the eigenproblem solver context
545: Output Parameter:
546: . lock - the locking flag
548: Level: advanced
550: .seealso: PEPSTOARSetLocking()
551: @*/
552: PetscErrorCode PEPSTOARGetLocking(PEP pep,PetscBool *lock)
553: {
559: PetscUseMethod(pep,"PEPSTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
560: return(0);
561: }
563: static PetscErrorCode PEPSTOARGetInertias_STOAR(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
564: {
566: PetscInt i,numsh;
567: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
568: PEP_SR sr = ctx->sr;
571: if (!pep->state) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSetUp() first");
572: if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see PEPSetInterval()");
573: switch (pep->state) {
574: case PEP_STATE_INITIAL:
575: break;
576: case PEP_STATE_SETUP:
577: if (n) *n = 2;
578: if (shifts) {
579: PetscMalloc1(2,shifts);
580: (*shifts)[0] = pep->inta;
581: (*shifts)[1] = pep->intb;
582: }
583: if (inertias) {
584: PetscMalloc1(2,inertias);
585: (*inertias)[0] = (sr->dir==1)?sr->inertia0:sr->inertia1;
586: (*inertias)[1] = (sr->dir==1)?sr->inertia1:sr->inertia0;
587: }
588: break;
589: case PEP_STATE_SOLVED:
590: case PEP_STATE_EIGENVECTORS:
591: numsh = ctx->nshifts;
592: if (n) *n = numsh;
593: if (shifts) {
594: PetscMalloc1(numsh,shifts);
595: for (i=0;i<numsh;i++) (*shifts)[i] = ctx->shifts[i];
596: }
597: if (inertias) {
598: PetscMalloc1(numsh,inertias);
599: for (i=0;i<numsh;i++) (*inertias)[i] = ctx->inertias[i];
600: }
601: break;
602: }
603: return(0);
604: }
606: /*@C
607: PEPSTOARGetInertias - Gets the values of the shifts and their
608: corresponding inertias in case of doing spectrum slicing for a
609: computational interval.
611: Not Collective
613: Input Parameter:
614: . pep - the eigenproblem solver context
616: Output Parameters:
617: + n - number of shifts, including the endpoints of the interval
618: . shifts - the values of the shifts used internally in the solver
619: - inertias - the values of the inertia in each shift
621: Notes:
622: If called after PEPSolve(), all shifts used internally by the solver are
623: returned (including both endpoints and any intermediate ones). If called
624: before PEPSolve() and after PEPSetUp() then only the information of the
625: endpoints of subintervals is available.
627: This function is only available for spectrum slicing runs.
629: The returned arrays should be freed by the user. Can pass NULL in any of
630: the two arrays if not required.
632: Fortran Notes:
633: The calling sequence from Fortran is
634: .vb
635: PEPSTOARGetInertias(pep,n,shifts,inertias,ierr)
636: integer n
637: double precision shifts(*)
638: integer inertias(*)
639: .ve
640: The arrays should be at least of length n. The value of n can be determined
641: by an initial call
642: .vb
643: PEPSTOARGetInertias(pep,n,PETSC_NULL_REAL,PETSC_NULL_INTEGER,ierr)
644: .ve
646: Level: advanced
648: .seealso: PEPSetInterval()
649: @*/
650: PetscErrorCode PEPSTOARGetInertias(PEP pep,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
651: {
657: PetscUseMethod(pep,"PEPSTOARGetInertias_C",(PEP,PetscInt*,PetscReal**,PetscInt**),(pep,n,shifts,inertias));
658: return(0);
659: }
661: static PetscErrorCode PEPSTOARSetDetectZeros_STOAR(PEP pep,PetscBool detect)
662: {
663: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
666: ctx->detect = detect;
667: pep->state = PEP_STATE_INITIAL;
668: return(0);
669: }
671: /*@
672: PEPSTOARSetDetectZeros - Sets a flag to enforce detection of
673: zeros during the factorizations throughout the spectrum slicing computation.
675: Logically Collective on PEP
677: Input Parameters:
678: + pep - the eigenproblem solver context
679: - detect - check for zeros
681: Options Database Key:
682: . -pep_stoar_detect_zeros - Check for zeros; this takes an optional
683: bool value (0/1/no/yes/true/false)
685: Notes:
686: A zero in the factorization indicates that a shift coincides with an eigenvalue.
688: This flag is turned off by default, and may be necessary in some cases.
689: This feature currently requires an external package for factorizations
690: with support for zero detection, e.g. MUMPS.
692: Level: advanced
694: .seealso: PEPSTOARSetPartitions(), PEPSetInterval()
695: @*/
696: PetscErrorCode PEPSTOARSetDetectZeros(PEP pep,PetscBool detect)
697: {
703: PetscTryMethod(pep,"PEPSTOARSetDetectZeros_C",(PEP,PetscBool),(pep,detect));
704: return(0);
705: }
707: static PetscErrorCode PEPSTOARGetDetectZeros_STOAR(PEP pep,PetscBool *detect)
708: {
709: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
712: *detect = ctx->detect;
713: return(0);
714: }
716: /*@
717: PEPSTOARGetDetectZeros - Gets the flag that enforces zero detection
718: in spectrum slicing.
720: Not Collective
722: Input Parameter:
723: . pep - the eigenproblem solver context
725: Output Parameter:
726: . detect - whether zeros detection is enforced during factorizations
728: Level: advanced
730: .seealso: PEPSTOARSetDetectZeros()
731: @*/
732: PetscErrorCode PEPSTOARGetDetectZeros(PEP pep,PetscBool *detect)
733: {
739: PetscUseMethod(pep,"PEPSTOARGetDetectZeros_C",(PEP,PetscBool*),(pep,detect));
740: return(0);
741: }
743: static PetscErrorCode PEPSTOARSetLinearization_STOAR(PEP pep,PetscReal alpha,PetscReal beta)
744: {
745: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
748: if (beta==0.0 && alpha==0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Parameters alpha and beta cannot be zero simultaneously");
749: ctx->alpha = alpha;
750: ctx->beta = beta;
751: return(0);
752: }
754: /*@
755: PEPSTOARSetLinearization - Set the coefficients that define
756: the linearization of a quadratic eigenproblem.
758: Logically Collective on PEP
760: Input Parameters:
761: + pep - polynomial eigenvalue solver
762: . alpha - first parameter of the linearization
763: - beta - second parameter of the linearization
765: Options Database Key:
766: . -pep_stoar_linearization <alpha,beta> - Sets the coefficients
768: Notes:
769: Cannot pass zero for both alpha and beta. The default values are
770: alpha=1 and beta=0.
772: Level: advanced
774: .seealso: PEPSTOARGetLinearization()
775: @*/
776: PetscErrorCode PEPSTOARSetLinearization(PEP pep,PetscReal alpha,PetscReal beta)
777: {
784: PetscTryMethod(pep,"PEPSTOARSetLinearization_C",(PEP,PetscReal,PetscReal),(pep,alpha,beta));
785: return(0);
786: }
788: static PetscErrorCode PEPSTOARGetLinearization_STOAR(PEP pep,PetscReal *alpha,PetscReal *beta)
789: {
790: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
793: if (alpha) *alpha = ctx->alpha;
794: if (beta) *beta = ctx->beta;
795: return(0);
796: }
798: /*@
799: PEPSTOARGetLinearization - Returns the coefficients that define
800: the linearization of a quadratic eigenproblem.
802: Not Collective
804: Input Parameter:
805: . pep - polynomial eigenvalue solver
807: Output Parameters:
808: + alpha - the first parameter of the linearization
809: - beta - the second parameter of the linearization
811: Level: advanced
813: .seealso: PEPSTOARSetLinearization()
814: @*/
815: PetscErrorCode PEPSTOARGetLinearization(PEP pep,PetscReal *alpha,PetscReal *beta)
816: {
821: PetscUseMethod(pep,"PEPSTOARGetLinearization_C",(PEP,PetscReal*,PetscReal*),(pep,alpha,beta));
822: return(0);
823: }
825: static PetscErrorCode PEPSTOARSetDimensions_STOAR(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
826: {
827: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
830: if (nev<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
831: ctx->nev = nev;
832: if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
833: ctx->ncv = 0;
834: } else {
835: if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
836: ctx->ncv = ncv;
837: }
838: if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
839: ctx->mpd = 0;
840: } else {
841: if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
842: ctx->mpd = mpd;
843: }
844: pep->state = PEP_STATE_INITIAL;
845: return(0);
846: }
848: /*@
849: PEPSTOARSetDimensions - Sets the dimensions used for each subsolve
850: step in case of doing spectrum slicing for a computational interval.
851: The meaning of the parameters is the same as in PEPSetDimensions().
853: Logically Collective on PEP
855: Input Parameters:
856: + pep - the eigenproblem solver context
857: . nev - number of eigenvalues to compute
858: . ncv - the maximum dimension of the subspace to be used by the subsolve
859: - mpd - the maximum dimension allowed for the projected problem
861: Options Database Key:
862: + -eps_stoar_nev <nev> - Sets the number of eigenvalues
863: . -eps_stoar_ncv <ncv> - Sets the dimension of the subspace
864: - -eps_stoar_mpd <mpd> - Sets the maximum projected dimension
866: Level: advanced
868: .seealso: PEPSTOARGetDimensions(), PEPSetDimensions(), PEPSetInterval()
869: @*/
870: PetscErrorCode PEPSTOARSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
871: {
879: PetscTryMethod(pep,"PEPSTOARSetDimensions_C",(PEP,PetscInt,PetscInt,PetscInt),(pep,nev,ncv,mpd));
880: return(0);
881: }
883: static PetscErrorCode PEPSTOARGetDimensions_STOAR(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
884: {
885: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
888: if (nev) *nev = ctx->nev;
889: if (ncv) *ncv = ctx->ncv;
890: if (mpd) *mpd = ctx->mpd;
891: return(0);
892: }
894: /*@
895: PEPSTOARGetDimensions - Gets the dimensions used for each subsolve
896: step in case of doing spectrum slicing for a computational interval.
898: Not Collective
900: Input Parameter:
901: . pep - the eigenproblem solver context
903: Output Parameters:
904: + nev - number of eigenvalues to compute
905: . ncv - the maximum dimension of the subspace to be used by the subsolve
906: - mpd - the maximum dimension allowed for the projected problem
908: Level: advanced
910: .seealso: PEPSTOARSetDimensions()
911: @*/
912: PetscErrorCode PEPSTOARGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
913: {
918: PetscUseMethod(pep,"PEPSTOARGetDimensions_C",(PEP,PetscInt*,PetscInt*,PetscInt*),(pep,nev,ncv,mpd));
919: return(0);
920: }
922: PetscErrorCode PEPView_STOAR(PEP pep,PetscViewer viewer)
923: {
925: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
926: PetscBool isascii;
929: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
930: if (isascii) {
931: PetscViewerASCIIPrintf(viewer," using the %slocking variant\n",ctx->lock?"":"non-");
932: PetscViewerASCIIPrintf(viewer," linearization parameters: alpha=%g beta=%g\n",(double)ctx->alpha,(double)ctx->beta);
933: }
934: return(0);
935: }
937: PetscErrorCode PEPReset_STOAR(PEP pep)
938: {
942: if (pep->which==PEP_ALL) {
943: PEPReset_STOAR_QSlice(pep);
944: }
945: return(0);
946: }
948: PetscErrorCode PEPDestroy_STOAR(PEP pep)
949: {
951: PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
954: BVDestroy(&ctx->V);
955: PetscFree(pep->data);
956: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",NULL);
957: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",NULL);
958: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDetectZeros_C",NULL);
959: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDetectZeros_C",NULL);
960: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetInertias_C",NULL);
961: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDimensions_C",NULL);
962: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDimensions_C",NULL);
963: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLinearization_C",NULL);
964: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLinearization_C",NULL);
965: return(0);
966: }
968: PETSC_EXTERN PetscErrorCode PEPCreate_STOAR(PEP pep)
969: {
971: PEP_TOAR *ctx;
974: PetscNewLog(pep,&ctx);
975: pep->data = (void*)ctx;
976: ctx->lock = PETSC_TRUE;
977: ctx->alpha = 1.0;
978: ctx->beta = 0.0;
980: pep->ops->setup = PEPSetUp_STOAR;
981: pep->ops->setfromoptions = PEPSetFromOptions_STOAR;
982: pep->ops->destroy = PEPDestroy_STOAR;
983: pep->ops->view = PEPView_STOAR;
984: pep->ops->backtransform = PEPBackTransform_Default;
985: pep->ops->computevectors = PEPComputeVectors_Default;
986: pep->ops->extractvectors = PEPExtractVectors_TOAR;
987: pep->ops->setdefaultst = PEPSetDefaultST_Transform;
988: pep->ops->reset = PEPReset_STOAR;
990: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLocking_C",PEPSTOARSetLocking_STOAR);
991: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLocking_C",PEPSTOARGetLocking_STOAR);
992: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDetectZeros_C",PEPSTOARSetDetectZeros_STOAR);
993: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDetectZeros_C",PEPSTOARGetDetectZeros_STOAR);
994: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetInertias_C",PEPSTOARGetInertias_STOAR);
995: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetDimensions_C",PEPSTOARGetDimensions_STOAR);
996: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetDimensions_C",PEPSTOARSetDimensions_STOAR);
997: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARSetLinearization_C",PEPSTOARSetLinearization_STOAR);
998: PetscObjectComposeFunction((PetscObject)pep,"PEPSTOARGetLinearization_C",PEPSTOARGetLinearization_STOAR);
999: return(0);
1000: }