Actual source code: qarnoldi.c
slepc-3.7.0 2016-05-16
1: /*
3: SLEPc quadratic eigensolver: "qarnoldi"
5: Method: Q-Arnoldi
7: Algorithm:
9: Quadratic Arnoldi with Krylov-Schur type restart.
11: References:
13: [1] K. Meerbergen, "The Quadratic Arnoldi method for the solution
14: of the quadratic eigenvalue problem", SIAM J. Matrix Anal.
15: Appl. 30(4):1462-1482, 2008.
17: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18: SLEPc - Scalable Library for Eigenvalue Problem Computations
19: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
21: This file is part of SLEPc.
23: SLEPc is free software: you can redistribute it and/or modify it under the
24: terms of version 3 of the GNU Lesser General Public License as published by
25: the Free Software Foundation.
27: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
28: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
29: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
30: more details.
32: You should have received a copy of the GNU Lesser General Public License
33: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
34: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: */
37: #include <slepc/private/pepimpl.h> /*I "slepcpep.h" I*/
38: #include <petscblaslapack.h>
40: typedef struct {
41: PetscReal keep; /* restart parameter */
42: PetscBool lock; /* locking/non-locking variant */
43: } PEP_QARNOLDI;
47: PetscErrorCode PEPSetUp_QArnoldi(PEP pep)
48: {
50: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
51: PetscBool sinv,flg;
54: pep->lineariz = PETSC_TRUE;
55: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
56: if (!ctx->lock && pep->mpd<pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
57: if (!pep->max_it) pep->max_it = PetscMax(100,4*pep->n/pep->ncv);
58: if (!pep->which) {
59: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
60: if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
61: else pep->which = PEP_LARGEST_MAGNITUDE;
62: }
64: if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver only available for quadratic problems");
65: if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
66: STGetTransform(pep->st,&flg);
67: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag set, see STSetTransform()");
69: /* set default extraction */
70: if (!pep->extract) {
71: pep->extract = PEP_EXTRACT_NONE;
72: }
73: if (pep->extract!=PEP_EXTRACT_NONE) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver does not support requested extraction");
74:
75: if (!ctx->keep) ctx->keep = 0.5;
77: PEPAllocateSolution(pep,0);
78: PEPSetWorkVecs(pep,4);
80: DSSetType(pep->ds,DSNHEP);
81: DSSetExtraRow(pep->ds,PETSC_TRUE);
82: DSAllocate(pep->ds,pep->ncv+1);
84: /* process starting vector */
85: if (pep->nini>-2) {
86: BVSetRandomColumn(pep->V,0);
87: BVSetRandomColumn(pep->V,1);
88: } else {
89: BVInsertVec(pep->V,0,pep->IS[0]);
90: BVInsertVec(pep->V,1,pep->IS[1]);
91: }
92: if (pep->nini<0) {
93: SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
94: }
95: return(0);
96: }
100: PetscErrorCode PEPExtractVectors_QArnoldi(PEP pep)
101: {
103: PetscInt i,k=pep->nconv,ldds;
104: PetscScalar *X,*pX0;
105: Mat X0;
108: if (pep->nconv==0) return(0);
109: DSGetLeadingDimension(pep->ds,&ldds);
110: DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
111: DSGetArray(pep->ds,DS_MAT_X,&X);
113: /* update vectors V = V*X */
114: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&X0);
115: MatDenseGetArray(X0,&pX0);
116: for (i=0;i<k;i++) {
117: PetscMemcpy(pX0+i*k,X+i*ldds,k*sizeof(PetscScalar));
118: }
119: MatDenseRestoreArray(X0,&pX0);
120: BVSetActiveColumns(pep->V,0,k);
121: BVMultInPlace(pep->V,X0,0,k);
122: MatDestroy(&X0);
123: BVSetActiveColumns(pep->V,0,k);
124: DSRestoreArray(pep->ds,DS_MAT_X,&X);
125: return(0);
126: }
130: /*
131: Compute a step of Classical Gram-Schmidt orthogonalization
132: */
133: 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)
134: {
136: PetscBLASInt ione = 1,j_1 = j+1;
137: PetscReal x,y;
138: PetscScalar dot,one = 1.0,zero = 0.0;
141: /* compute norm of v and w */
142: if (onorm) {
143: VecNorm(v,NORM_2,&x);
144: VecNorm(w,NORM_2,&y);
145: *onorm = PetscSqrtReal(x*x+y*y);
146: }
148: /* orthogonalize: compute h */
149: BVDotVec(V,v,h);
150: BVDotVec(V,w,work);
151: if (j>0)
152: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione));
153: VecDot(w,t,&dot);
154: h[j] += dot;
156: /* orthogonalize: update v and w */
157: BVMultVec(V,-1.0,1.0,v,h);
158: if (j>0) {
159: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione));
160: BVMultVec(V,-1.0,1.0,w,work);
161: }
162: VecAXPY(w,-h[j],t);
164: /* compute norm of v and w */
165: if (norm) {
166: VecNorm(v,NORM_2,&x);
167: VecNorm(w,NORM_2,&y);
168: *norm = PetscSqrtReal(x*x+y*y);
169: }
170: return(0);
171: }
175: /*
176: Compute a run of Q-Arnoldi iterations
177: */
178: static PetscErrorCode PEPQArnoldi(PEP pep,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
179: {
180: PetscErrorCode ierr;
181: PetscInt i,j,l,m = *M;
182: Vec t = pep->work[2],u = pep->work[3];
183: BVOrthogRefineType refinement;
184: PetscReal norm=0.0,onorm,eta;
185: PetscScalar *c = work + m;
188: BVGetOrthogonalization(pep->V,NULL,&refinement,&eta,NULL);
189: BVInsertVec(pep->V,k,v);
190: for (j=k;j<m;j++) {
191: /* apply operator */
192: VecCopy(w,t);
193: if (pep->Dr) {
194: VecPointwiseMult(v,v,pep->Dr);
195: }
196: STMatMult(pep->st,0,v,u);
197: VecCopy(t,v);
198: if (pep->Dr) {
199: VecPointwiseMult(t,t,pep->Dr);
200: }
201: STMatMult(pep->st,1,t,w);
202: VecAXPY(u,pep->sfactor,w);
203: STMatSolve(pep->st,u,w);
204: VecScale(w,-1.0/(pep->sfactor*pep->sfactor));
205: if (pep->Dr) {
206: VecPointwiseDivide(w,w,pep->Dr);
207: }
208: VecCopy(v,t);
209: BVSetActiveColumns(pep->V,0,j+1);
211: /* orthogonalize */
212: switch (refinement) {
213: case BV_ORTHOG_REFINE_NEVER:
214: PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,&norm,work);
215: *breakdown = PETSC_FALSE;
216: break;
217: case BV_ORTHOG_REFINE_ALWAYS:
218: PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,NULL,work);
219: PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,&onorm,&norm,work);
220: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
221: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
222: else *breakdown = PETSC_FALSE;
223: break;
224: case BV_ORTHOG_REFINE_IFNEEDED:
225: PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,&onorm,&norm,work);
226: /* ||q|| < eta ||h|| */
227: l = 1;
228: while (l<3 && norm < eta * onorm) {
229: l++;
230: onorm = norm;
231: PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,NULL,&norm,work);
232: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
233: }
234: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
235: else *breakdown = PETSC_FALSE;
236: break;
237: default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of ip->orth_ref");
238: }
239: VecScale(v,1.0/norm);
240: VecScale(w,1.0/norm);
242: H[j+1+ldh*j] = norm;
243: if (j<m-1) {
244: BVInsertVec(pep->V,j+1,v);
245: }
246: }
247: *beta = norm;
248: return(0);
249: }
253: PetscErrorCode PEPSolve_QArnoldi(PEP pep)
254: {
256: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
257: PetscInt j,k,l,lwork,nv,ld,newn,nconv;
258: Vec v=pep->work[0],w=pep->work[1];
259: Mat Q;
260: PetscScalar *S,*work;
261: PetscReal beta=0.0,norm,x,y;
262: PetscBool breakdown=PETSC_FALSE,sinv;
265: DSGetLeadingDimension(pep->ds,&ld);
266: lwork = 7*pep->ncv;
267: PetscMalloc1(lwork,&work);
268: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
269: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
270: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
272: /* Get the starting Arnoldi vector */
273: BVCopyVec(pep->V,0,v);
274: BVCopyVec(pep->V,1,w);
275: VecNorm(v,NORM_2,&x);
276: VecNorm(w,NORM_2,&y);
277: norm = PetscSqrtReal(x*x+y*y);
278: VecScale(v,1.0/norm);
279: VecScale(w,1.0/norm);
281: /* Restart loop */
282: l = 0;
283: while (pep->reason == PEP_CONVERGED_ITERATING) {
284: pep->its++;
286: /* Compute an nv-step Arnoldi factorization */
287: nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
288: DSGetArray(pep->ds,DS_MAT_A,&S);
289: PEPQArnoldi(pep,S,ld,pep->nconv+l,&nv,v,w,&beta,&breakdown,work);
290: DSRestoreArray(pep->ds,DS_MAT_A,&S);
291: DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
292: if (l==0) {
293: DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
294: } else {
295: DSSetState(pep->ds,DS_STATE_RAW);
296: }
297: BVSetActiveColumns(pep->V,pep->nconv,nv);
299: /* Solve projected problem */
300: DSSolve(pep->ds,pep->eigr,pep->eigi);
301: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
302: DSUpdateExtraRow(pep->ds);
304: /* Check convergence */
305: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
306: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
307: nconv = k;
309: /* Update l */
310: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
311: else l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
312: if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
314: if (pep->reason == PEP_CONVERGED_ITERATING) {
315: if (breakdown) {
316: /* Stop if breakdown */
317: PetscInfo2(pep,"Breakdown Quadratic Arnoldi method (it=%D norm=%g)\n",pep->its,(double)beta);
318: pep->reason = PEP_DIVERGED_BREAKDOWN;
319: } else {
320: /* Prepare the Rayleigh quotient for restart */
321: DSTruncate(pep->ds,k+l);
322: DSGetDimensions(pep->ds,&newn,NULL,NULL,NULL,NULL);
323: l = newn-k;
324: }
325: }
326: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
327: DSGetMat(pep->ds,DS_MAT_Q,&Q);
328: BVMultInPlace(pep->V,Q,pep->nconv,k+l);
329: MatDestroy(&Q);
331: pep->nconv = k;
332: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
333: }
335: for (j=0;j<pep->nconv;j++) {
336: pep->eigr[j] *= pep->sfactor;
337: pep->eigi[j] *= pep->sfactor;
338: }
340: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
341: RGPopScale(pep->rg);
343: /* truncate Schur decomposition and change the state to raw so that
344: DSVectors() computes eigenvectors from scratch */
345: DSSetDimensions(pep->ds,pep->nconv,0,0,0);
346: DSSetState(pep->ds,DS_STATE_RAW);
347: PetscFree(work);
348: return(0);
349: }
353: static PetscErrorCode PEPQArnoldiSetRestart_QArnoldi(PEP pep,PetscReal keep)
354: {
355: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
358: if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
359: else {
360: 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]");
361: ctx->keep = keep;
362: }
363: return(0);
364: }
368: /*@
369: PEPQArnoldiSetRestart - Sets the restart parameter for the Q-Arnoldi
370: method, in particular the proportion of basis vectors that must be kept
371: after restart.
373: Logically Collective on PEP
375: Input Parameters:
376: + pep - the eigenproblem solver context
377: - keep - the number of vectors to be kept at restart
379: Options Database Key:
380: . -pep_qarnoldi_restart - Sets the restart parameter
382: Notes:
383: Allowed values are in the range [0.1,0.9]. The default is 0.5.
385: Level: advanced
387: .seealso: PEPQArnoldiGetRestart()
388: @*/
389: PetscErrorCode PEPQArnoldiSetRestart(PEP pep,PetscReal keep)
390: {
396: PetscTryMethod(pep,"PEPQArnoldiSetRestart_C",(PEP,PetscReal),(pep,keep));
397: return(0);
398: }
402: static PetscErrorCode PEPQArnoldiGetRestart_QArnoldi(PEP pep,PetscReal *keep)
403: {
404: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
407: *keep = ctx->keep;
408: return(0);
409: }
413: /*@
414: PEPQArnoldiGetRestart - Gets the restart parameter used in the Q-Arnoldi method.
416: Not Collective
418: Input Parameter:
419: . pep - the eigenproblem solver context
421: Output Parameter:
422: . keep - the restart parameter
424: Level: advanced
426: .seealso: PEPQArnoldiSetRestart()
427: @*/
428: PetscErrorCode PEPQArnoldiGetRestart(PEP pep,PetscReal *keep)
429: {
435: PetscUseMethod(pep,"PEPQArnoldiGetRestart_C",(PEP,PetscReal*),(pep,keep));
436: return(0);
437: }
441: static PetscErrorCode PEPQArnoldiSetLocking_QArnoldi(PEP pep,PetscBool lock)
442: {
443: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
446: ctx->lock = lock;
447: return(0);
448: }
452: /*@
453: PEPQArnoldiSetLocking - Choose between locking and non-locking variants of
454: the Q-Arnoldi method.
456: Logically Collective on PEP
458: Input Parameters:
459: + pep - the eigenproblem solver context
460: - lock - true if the locking variant must be selected
462: Options Database Key:
463: . -pep_qarnoldi_locking - Sets the locking flag
465: Notes:
466: The default is to keep all directions in the working subspace even if
467: already converged to working accuracy (the non-locking variant).
468: This behaviour can be changed so that converged eigenpairs are locked
469: when the method restarts.
471: Note that the default behaviour is the opposite to Krylov solvers in EPS.
473: Level: advanced
475: .seealso: PEPQArnoldiGetLocking()
476: @*/
477: PetscErrorCode PEPQArnoldiSetLocking(PEP pep,PetscBool lock)
478: {
484: PetscTryMethod(pep,"PEPQArnoldiSetLocking_C",(PEP,PetscBool),(pep,lock));
485: return(0);
486: }
490: static PetscErrorCode PEPQArnoldiGetLocking_QArnoldi(PEP pep,PetscBool *lock)
491: {
492: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
495: *lock = ctx->lock;
496: return(0);
497: }
501: /*@
502: PEPQArnoldiGetLocking - Gets the locking flag used in the Q-Arnoldi method.
504: Not Collective
506: Input Parameter:
507: . pep - the eigenproblem solver context
509: Output Parameter:
510: . lock - the locking flag
512: Level: advanced
514: .seealso: PEPQArnoldiSetLocking()
515: @*/
516: PetscErrorCode PEPQArnoldiGetLocking(PEP pep,PetscBool *lock)
517: {
523: PetscUseMethod(pep,"PEPQArnoldiGetLocking_C",(PEP,PetscBool*),(pep,lock));
524: return(0);
525: }
529: PetscErrorCode PEPSetFromOptions_QArnoldi(PetscOptionItems *PetscOptionsObject,PEP pep)
530: {
532: PetscBool flg,lock;
533: PetscReal keep;
536: PetscOptionsHead(PetscOptionsObject,"PEP Q-Arnoldi Options");
537: PetscOptionsReal("-pep_qarnoldi_restart","Proportion of vectors kept after restart","PEPQArnoldiSetRestart",0.5,&keep,&flg);
538: if (flg) {
539: PEPQArnoldiSetRestart(pep,keep);
540: }
541: PetscOptionsBool("-pep_qarnoldi_locking","Choose between locking and non-locking variants","PEPQArnoldiSetLocking",PETSC_FALSE,&lock,&flg);
542: if (flg) {
543: PEPQArnoldiSetLocking(pep,lock);
544: }
545: PetscOptionsTail();
546: return(0);
547: }
551: PetscErrorCode PEPView_QArnoldi(PEP pep,PetscViewer viewer)
552: {
554: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
555: PetscBool isascii;
558: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
559: if (isascii) {
560: PetscViewerASCIIPrintf(viewer," Q-Arnoldi: %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
561: PetscViewerASCIIPrintf(viewer," Q-Arnoldi: using the %slocking variant\n",ctx->lock?"":"non-");
562: }
563: return(0);
564: }
568: PetscErrorCode PEPDestroy_QArnoldi(PEP pep)
569: {
573: PetscFree(pep->data);
574: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",NULL);
575: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",NULL);
576: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",NULL);
577: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",NULL);
578: return(0);
579: }
583: PETSC_EXTERN PetscErrorCode PEPCreate_QArnoldi(PEP pep)
584: {
585: PEP_QARNOLDI *ctx;
589: PetscNewLog(pep,&ctx);
590: pep->data = (void*)ctx;
591: ctx->lock = PETSC_TRUE;
593: pep->ops->solve = PEPSolve_QArnoldi;
594: pep->ops->setup = PEPSetUp_QArnoldi;
595: pep->ops->setfromoptions = PEPSetFromOptions_QArnoldi;
596: pep->ops->destroy = PEPDestroy_QArnoldi;
597: pep->ops->view = PEPView_QArnoldi;
598: pep->ops->backtransform = PEPBackTransform_Default;
599: pep->ops->computevectors = PEPComputeVectors_Default;
600: pep->ops->computevectors = PEPExtractVectors_QArnoldi;
601: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",PEPQArnoldiSetRestart_QArnoldi);
602: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",PEPQArnoldiGetRestart_QArnoldi);
603: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",PEPQArnoldiSetLocking_QArnoldi);
604: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",PEPQArnoldiGetLocking_QArnoldi);
605: return(0);
606: }