Actual source code: qarnoldi.c
slepc-3.8.0 2017-10-20
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, 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 quadratic eigensolver: "qarnoldi"
13: Method: Q-Arnoldi
15: Algorithm:
17: Quadratic Arnoldi with Krylov-Schur type restart.
19: References:
21: [1] K. Meerbergen, "The Quadratic Arnoldi method for the solution
22: of the quadratic eigenvalue problem", SIAM J. Matrix Anal.
23: Appl. 30(4):1462-1482, 2008.
24: */
26: #include <slepc/private/pepimpl.h> /*I "slepcpep.h" I*/
27: #include <petscblaslapack.h>
29: typedef struct {
30: PetscReal keep; /* restart parameter */
31: PetscBool lock; /* locking/non-locking variant */
32: } PEP_QARNOLDI;
34: PetscErrorCode PEPSetUp_QArnoldi(PEP pep)
35: {
37: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
38: PetscBool shift,sinv,flg;
41: pep->lineariz = PETSC_TRUE;
42: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
43: if (!ctx->lock && pep->mpd<pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
44: if (!pep->max_it) pep->max_it = PetscMax(100,4*pep->n/pep->ncv);
46: PetscObjectTypeCompare((PetscObject)pep->st,STSHIFT,&shift);
47: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
48: if (!shift && !sinv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Only STSHIFT and STSINVERT spectral transformations can be used");
49: if (!pep->which) {
50: if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
51: else pep->which = PEP_LARGEST_MAGNITUDE;
52: }
54: if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver only available for quadratic problems");
55: if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
56: STGetTransform(pep->st,&flg);
57: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag set, see STSetTransform()");
59: /* set default extraction */
60: if (!pep->extract) {
61: pep->extract = PEP_EXTRACT_NONE;
62: }
63: if (pep->extract!=PEP_EXTRACT_NONE) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver does not support requested extraction");
65: if (!ctx->keep) ctx->keep = 0.5;
67: PEPAllocateSolution(pep,0);
68: PEPSetWorkVecs(pep,4);
70: DSSetType(pep->ds,DSNHEP);
71: DSSetExtraRow(pep->ds,PETSC_TRUE);
72: DSAllocate(pep->ds,pep->ncv+1);
74: return(0);
75: }
77: PetscErrorCode PEPExtractVectors_QArnoldi(PEP pep)
78: {
80: PetscInt i,k=pep->nconv,ldds;
81: PetscScalar *X,*pX0;
82: Mat X0;
85: if (pep->nconv==0) return(0);
86: DSGetLeadingDimension(pep->ds,&ldds);
87: DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
88: DSGetArray(pep->ds,DS_MAT_X,&X);
90: /* update vectors V = V*X */
91: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&X0);
92: MatDenseGetArray(X0,&pX0);
93: for (i=0;i<k;i++) {
94: PetscMemcpy(pX0+i*k,X+i*ldds,k*sizeof(PetscScalar));
95: }
96: MatDenseRestoreArray(X0,&pX0);
97: BVMultInPlace(pep->V,X0,0,k);
98: MatDestroy(&X0);
99: DSRestoreArray(pep->ds,DS_MAT_X,&X);
100: return(0);
101: }
103: /*
104: Compute a step of Classical Gram-Schmidt orthogonalization
105: */
106: static PetscErrorCode PEPQArnoldiCGS(PEP pep,PetscScalar *H,PetscBLASInt ldh,PetscScalar *h,PetscBLASInt j,BV V,Vec t,Vec v,Vec w,PetscReal *onorm,PetscReal *norm,PetscScalar *work)
107: {
109: PetscBLASInt ione = 1,j_1 = j+1;
110: PetscReal x,y;
111: PetscScalar dot,one = 1.0,zero = 0.0;
114: /* compute norm of v and w */
115: if (onorm) {
116: VecNorm(v,NORM_2,&x);
117: VecNorm(w,NORM_2,&y);
118: *onorm = PetscSqrtReal(x*x+y*y);
119: }
121: /* orthogonalize: compute h */
122: BVDotVec(V,v,h);
123: BVDotVec(V,w,work);
124: if (j>0)
125: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione));
126: VecDot(w,t,&dot);
127: h[j] += dot;
129: /* orthogonalize: update v and w */
130: BVMultVec(V,-1.0,1.0,v,h);
131: if (j>0) {
132: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione));
133: BVMultVec(V,-1.0,1.0,w,work);
134: }
135: VecAXPY(w,-h[j],t);
137: /* compute norm of v and w */
138: if (norm) {
139: VecNorm(v,NORM_2,&x);
140: VecNorm(w,NORM_2,&y);
141: *norm = PetscSqrtReal(x*x+y*y);
142: }
143: return(0);
144: }
146: /*
147: Compute a run of Q-Arnoldi iterations
148: */
149: static PetscErrorCode PEPQArnoldi(PEP pep,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
150: {
151: PetscErrorCode ierr;
152: PetscInt i,j,l,m = *M;
153: Vec t = pep->work[2],u = pep->work[3];
154: BVOrthogRefineType refinement;
155: PetscReal norm=0.0,onorm,eta;
156: PetscScalar *c = work + m;
159: BVGetOrthogonalization(pep->V,NULL,&refinement,&eta,NULL);
160: BVInsertVec(pep->V,k,v);
161: for (j=k;j<m;j++) {
162: /* apply operator */
163: VecCopy(w,t);
164: if (pep->Dr) {
165: VecPointwiseMult(v,v,pep->Dr);
166: }
167: STMatMult(pep->st,0,v,u);
168: VecCopy(t,v);
169: if (pep->Dr) {
170: VecPointwiseMult(t,t,pep->Dr);
171: }
172: STMatMult(pep->st,1,t,w);
173: VecAXPY(u,pep->sfactor,w);
174: STMatSolve(pep->st,u,w);
175: VecScale(w,-1.0/(pep->sfactor*pep->sfactor));
176: if (pep->Dr) {
177: VecPointwiseDivide(w,w,pep->Dr);
178: }
179: VecCopy(v,t);
180: BVSetActiveColumns(pep->V,0,j+1);
182: /* orthogonalize */
183: switch (refinement) {
184: case BV_ORTHOG_REFINE_NEVER:
185: PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,&norm,work);
186: *breakdown = PETSC_FALSE;
187: break;
188: case BV_ORTHOG_REFINE_ALWAYS:
189: PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,NULL,work);
190: PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,&onorm,&norm,work);
191: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
192: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
193: else *breakdown = PETSC_FALSE;
194: break;
195: case BV_ORTHOG_REFINE_IFNEEDED:
196: PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,&onorm,&norm,work);
197: /* ||q|| < eta ||h|| */
198: l = 1;
199: while (l<3 && norm < eta * onorm) {
200: l++;
201: onorm = norm;
202: PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,NULL,&norm,work);
203: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
204: }
205: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
206: else *breakdown = PETSC_FALSE;
207: break;
208: }
209: VecScale(v,1.0/norm);
210: VecScale(w,1.0/norm);
212: H[j+1+ldh*j] = norm;
213: if (j<m-1) {
214: BVInsertVec(pep->V,j+1,v);
215: }
216: }
217: *beta = norm;
218: return(0);
219: }
221: PetscErrorCode PEPSolve_QArnoldi(PEP pep)
222: {
224: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
225: PetscInt j,k,l,lwork,nv,ld,newn,nconv;
226: Vec v=pep->work[0],w=pep->work[1];
227: Mat Q;
228: PetscScalar *S,*work;
229: PetscReal beta=0.0,norm,x,y;
230: PetscBool breakdown=PETSC_FALSE,sinv;
233: DSGetLeadingDimension(pep->ds,&ld);
234: lwork = 7*pep->ncv;
235: PetscMalloc1(lwork,&work);
236: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
237: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
238: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
240: /* Get the starting Arnoldi vector */
241: for (j=0;j<2;j++) {
242: if (j>=pep->nini) {
243: BVSetRandomColumn(pep->V,j);
244: }
245: }
246: BVCopyVec(pep->V,0,v);
247: BVCopyVec(pep->V,1,w);
248: VecNorm(v,NORM_2,&x);
249: VecNorm(w,NORM_2,&y);
250: norm = PetscSqrtReal(x*x+y*y);
251: VecScale(v,1.0/norm);
252: VecScale(w,1.0/norm);
254: /* clean projected matrix (including the extra-arrow) */
255: DSGetArray(pep->ds,DS_MAT_A,&S);
256: PetscMemzero(S,ld*ld*sizeof(PetscScalar));
257: DSRestoreArray(pep->ds,DS_MAT_A,&S);
258:
259: /* Restart loop */
260: l = 0;
261: while (pep->reason == PEP_CONVERGED_ITERATING) {
262: pep->its++;
264: /* Compute an nv-step Arnoldi factorization */
265: nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
266: DSGetArray(pep->ds,DS_MAT_A,&S);
267: PEPQArnoldi(pep,S,ld,pep->nconv+l,&nv,v,w,&beta,&breakdown,work);
268: DSRestoreArray(pep->ds,DS_MAT_A,&S);
269: DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
270: if (l==0) {
271: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
272: } else {
273: DSSetState(pep->ds,DS_STATE_RAW);
274: }
275: BVSetActiveColumns(pep->V,pep->nconv,nv);
277: /* Solve projected problem */
278: DSSolve(pep->ds,pep->eigr,pep->eigi);
279: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
280: DSUpdateExtraRow(pep->ds);
281: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
283: /* Check convergence */
284: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
285: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
286: nconv = k;
288: /* Update l */
289: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
290: else l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
291: if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
293: if (pep->reason == PEP_CONVERGED_ITERATING) {
294: if (breakdown) {
295: /* Stop if breakdown */
296: PetscInfo2(pep,"Breakdown Quadratic Arnoldi method (it=%D norm=%g)\n",pep->its,(double)beta);
297: pep->reason = PEP_DIVERGED_BREAKDOWN;
298: } else {
299: /* Prepare the Rayleigh quotient for restart */
300: DSTruncate(pep->ds,k+l);
301: DSGetDimensions(pep->ds,&newn,NULL,NULL,NULL,NULL);
302: l = newn-k;
303: }
304: }
305: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
306: DSGetMat(pep->ds,DS_MAT_Q,&Q);
307: BVMultInPlace(pep->V,Q,pep->nconv,k+l);
308: MatDestroy(&Q);
310: pep->nconv = k;
311: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
312: }
314: for (j=0;j<pep->nconv;j++) {
315: pep->eigr[j] *= pep->sfactor;
316: pep->eigi[j] *= pep->sfactor;
317: }
319: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
320: RGPopScale(pep->rg);
322: /* truncate Schur decomposition and change the state to raw so that
323: DSVectors() computes eigenvectors from scratch */
324: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
325: DSSetState(pep->ds,DS_STATE_RAW);
326: PetscFree(work);
327: return(0);
328: }
330: static PetscErrorCode PEPQArnoldiSetRestart_QArnoldi(PEP pep,PetscReal keep)
331: {
332: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
335: if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
336: else {
337: if (keep<0.1 || keep>0.9) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
338: ctx->keep = keep;
339: }
340: return(0);
341: }
343: /*@
344: PEPQArnoldiSetRestart - Sets the restart parameter for the Q-Arnoldi
345: method, in particular the proportion of basis vectors that must be kept
346: after restart.
348: Logically Collective on PEP
350: Input Parameters:
351: + pep - the eigenproblem solver context
352: - keep - the number of vectors to be kept at restart
354: Options Database Key:
355: . -pep_qarnoldi_restart - Sets the restart parameter
357: Notes:
358: Allowed values are in the range [0.1,0.9]. The default is 0.5.
360: Level: advanced
362: .seealso: PEPQArnoldiGetRestart()
363: @*/
364: PetscErrorCode PEPQArnoldiSetRestart(PEP pep,PetscReal keep)
365: {
371: PetscTryMethod(pep,"PEPQArnoldiSetRestart_C",(PEP,PetscReal),(pep,keep));
372: return(0);
373: }
375: static PetscErrorCode PEPQArnoldiGetRestart_QArnoldi(PEP pep,PetscReal *keep)
376: {
377: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
380: *keep = ctx->keep;
381: return(0);
382: }
384: /*@
385: PEPQArnoldiGetRestart - Gets the restart parameter used in the Q-Arnoldi method.
387: Not Collective
389: Input Parameter:
390: . pep - the eigenproblem solver context
392: Output Parameter:
393: . keep - the restart parameter
395: Level: advanced
397: .seealso: PEPQArnoldiSetRestart()
398: @*/
399: PetscErrorCode PEPQArnoldiGetRestart(PEP pep,PetscReal *keep)
400: {
406: PetscUseMethod(pep,"PEPQArnoldiGetRestart_C",(PEP,PetscReal*),(pep,keep));
407: return(0);
408: }
410: static PetscErrorCode PEPQArnoldiSetLocking_QArnoldi(PEP pep,PetscBool lock)
411: {
412: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
415: ctx->lock = lock;
416: return(0);
417: }
419: /*@
420: PEPQArnoldiSetLocking - Choose between locking and non-locking variants of
421: the Q-Arnoldi method.
423: Logically Collective on PEP
425: Input Parameters:
426: + pep - the eigenproblem solver context
427: - lock - true if the locking variant must be selected
429: Options Database Key:
430: . -pep_qarnoldi_locking - Sets the locking flag
432: Notes:
433: The default is to keep all directions in the working subspace even if
434: already converged to working accuracy (the non-locking variant).
435: This behaviour can be changed so that converged eigenpairs are locked
436: when the method restarts.
438: Note that the default behaviour is the opposite to Krylov solvers in EPS.
440: Level: advanced
442: .seealso: PEPQArnoldiGetLocking()
443: @*/
444: PetscErrorCode PEPQArnoldiSetLocking(PEP pep,PetscBool lock)
445: {
451: PetscTryMethod(pep,"PEPQArnoldiSetLocking_C",(PEP,PetscBool),(pep,lock));
452: return(0);
453: }
455: static PetscErrorCode PEPQArnoldiGetLocking_QArnoldi(PEP pep,PetscBool *lock)
456: {
457: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
460: *lock = ctx->lock;
461: return(0);
462: }
464: /*@
465: PEPQArnoldiGetLocking - Gets the locking flag used in the Q-Arnoldi method.
467: Not Collective
469: Input Parameter:
470: . pep - the eigenproblem solver context
472: Output Parameter:
473: . lock - the locking flag
475: Level: advanced
477: .seealso: PEPQArnoldiSetLocking()
478: @*/
479: PetscErrorCode PEPQArnoldiGetLocking(PEP pep,PetscBool *lock)
480: {
486: PetscUseMethod(pep,"PEPQArnoldiGetLocking_C",(PEP,PetscBool*),(pep,lock));
487: return(0);
488: }
490: PetscErrorCode PEPSetFromOptions_QArnoldi(PetscOptionItems *PetscOptionsObject,PEP pep)
491: {
493: PetscBool flg,lock;
494: PetscReal keep;
497: PetscOptionsHead(PetscOptionsObject,"PEP Q-Arnoldi Options");
499: PetscOptionsReal("-pep_qarnoldi_restart","Proportion of vectors kept after restart","PEPQArnoldiSetRestart",0.5,&keep,&flg);
500: if (flg) { PEPQArnoldiSetRestart(pep,keep); }
502: PetscOptionsBool("-pep_qarnoldi_locking","Choose between locking and non-locking variants","PEPQArnoldiSetLocking",PETSC_FALSE,&lock,&flg);
503: if (flg) { PEPQArnoldiSetLocking(pep,lock); }
505: PetscOptionsTail();
506: return(0);
507: }
509: PetscErrorCode PEPView_QArnoldi(PEP pep,PetscViewer viewer)
510: {
512: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
513: PetscBool isascii;
516: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
517: if (isascii) {
518: PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
519: PetscViewerASCIIPrintf(viewer," using the %slocking variant\n",ctx->lock?"":"non-");
520: }
521: return(0);
522: }
524: PetscErrorCode PEPDestroy_QArnoldi(PEP pep)
525: {
529: PetscFree(pep->data);
530: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",NULL);
531: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",NULL);
532: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",NULL);
533: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",NULL);
534: return(0);
535: }
537: PETSC_EXTERN PetscErrorCode PEPCreate_QArnoldi(PEP pep)
538: {
539: PEP_QARNOLDI *ctx;
543: PetscNewLog(pep,&ctx);
544: pep->data = (void*)ctx;
545: ctx->lock = PETSC_TRUE;
547: pep->ops->solve = PEPSolve_QArnoldi;
548: pep->ops->setup = PEPSetUp_QArnoldi;
549: pep->ops->setfromoptions = PEPSetFromOptions_QArnoldi;
550: pep->ops->destroy = PEPDestroy_QArnoldi;
551: pep->ops->view = PEPView_QArnoldi;
552: pep->ops->backtransform = PEPBackTransform_Default;
553: pep->ops->computevectors = PEPComputeVectors_Default;
554: pep->ops->extractvectors = PEPExtractVectors_QArnoldi;
555: pep->ops->setdefaultst = PEPSetDefaultST_Transform;
557: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",PEPQArnoldiSetRestart_QArnoldi);
558: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",PEPQArnoldiGetRestart_QArnoldi);
559: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",PEPQArnoldiSetLocking_QArnoldi);
560: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",PEPQArnoldiGetLocking_QArnoldi);
561: return(0);
562: }