Actual source code: qarnoldi.c
slepc-3.17.0 2022-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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>
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: {
36: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
37: PetscBool flg;
39: PEPCheckQuadratic(pep);
40: PEPCheckShiftSinvert(pep);
41: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
43: if (pep->max_it==PETSC_DEFAULT) pep->max_it = PetscMax(100,4*pep->n/pep->ncv);
44: if (!pep->which) PEPSetWhichEigenpairs_Default(pep);
47: STGetTransform(pep->st,&flg);
50: /* set default extraction */
51: if (!pep->extract) pep->extract = PEP_EXTRACT_NONE;
52: PEPCheckUnsupported(pep,PEP_FEATURE_NONMONOMIAL | PEP_FEATURE_EXTRACT);
54: if (!ctx->keep) ctx->keep = 0.5;
56: PEPAllocateSolution(pep,0);
57: PEPSetWorkVecs(pep,4);
59: DSSetType(pep->ds,DSNHEP);
60: DSSetExtraRow(pep->ds,PETSC_TRUE);
61: DSAllocate(pep->ds,pep->ncv+1);
63: PetscFunctionReturn(0);
64: }
66: PetscErrorCode PEPExtractVectors_QArnoldi(PEP pep)
67: {
68: PetscInt i,k=pep->nconv,ldds;
69: PetscScalar *X,*pX0;
70: Mat X0;
72: if (pep->nconv==0) PetscFunctionReturn(0);
73: DSGetLeadingDimension(pep->ds,&ldds);
74: DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
75: DSGetArray(pep->ds,DS_MAT_X,&X);
77: /* update vectors V = V*X */
78: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&X0);
79: MatDenseGetArrayWrite(X0,&pX0);
80: for (i=0;i<k;i++) PetscArraycpy(pX0+i*k,X+i*ldds,k);
81: MatDenseRestoreArrayWrite(X0,&pX0);
82: BVMultInPlace(pep->V,X0,0,k);
83: MatDestroy(&X0);
84: DSRestoreArray(pep->ds,DS_MAT_X,&X);
85: PetscFunctionReturn(0);
86: }
88: /*
89: Compute a step of Classical Gram-Schmidt orthogonalization
90: */
91: 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)
92: {
93: PetscBLASInt ione = 1,j_1 = j+1;
94: PetscReal x,y;
95: PetscScalar dot,one = 1.0,zero = 0.0;
97: /* compute norm of v and w */
98: if (onorm) {
99: VecNorm(v,NORM_2,&x);
100: VecNorm(w,NORM_2,&y);
101: *onorm = PetscSqrtReal(x*x+y*y);
102: }
104: /* orthogonalize: compute h */
105: BVDotVec(V,v,h);
106: BVDotVec(V,w,work);
107: if (j>0)
108: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione));
109: VecDot(w,t,&dot);
110: h[j] += dot;
112: /* orthogonalize: update v and w */
113: BVMultVec(V,-1.0,1.0,v,h);
114: if (j>0) {
115: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione));
116: BVMultVec(V,-1.0,1.0,w,work);
117: }
118: VecAXPY(w,-h[j],t);
120: /* compute norm of v and w */
121: if (norm) {
122: VecNorm(v,NORM_2,&x);
123: VecNorm(w,NORM_2,&y);
124: *norm = PetscSqrtReal(x*x+y*y);
125: }
126: PetscFunctionReturn(0);
127: }
129: /*
130: Compute a run of Q-Arnoldi iterations
131: */
132: static PetscErrorCode PEPQArnoldi(PEP pep,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
133: {
134: PetscInt i,j,l,m = *M;
135: Vec t = pep->work[2],u = pep->work[3];
136: BVOrthogRefineType refinement;
137: PetscReal norm=0.0,onorm,eta;
138: PetscScalar *c = work + m;
140: BVGetOrthogonalization(pep->V,NULL,&refinement,&eta,NULL);
141: BVInsertVec(pep->V,k,v);
142: for (j=k;j<m;j++) {
143: /* apply operator */
144: VecCopy(w,t);
145: if (pep->Dr) VecPointwiseMult(v,v,pep->Dr);
146: STMatMult(pep->st,0,v,u);
147: VecCopy(t,v);
148: if (pep->Dr) VecPointwiseMult(t,t,pep->Dr);
149: STMatMult(pep->st,1,t,w);
150: VecAXPY(u,pep->sfactor,w);
151: STMatSolve(pep->st,u,w);
152: VecScale(w,-1.0/(pep->sfactor*pep->sfactor));
153: if (pep->Dr) VecPointwiseDivide(w,w,pep->Dr);
154: VecCopy(v,t);
155: BVSetActiveColumns(pep->V,0,j+1);
157: /* orthogonalize */
158: switch (refinement) {
159: case BV_ORTHOG_REFINE_NEVER:
160: PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,&norm,work);
161: *breakdown = PETSC_FALSE;
162: break;
163: case BV_ORTHOG_REFINE_ALWAYS:
164: PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,NULL,work);
165: PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,&onorm,&norm,work);
166: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
167: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
168: else *breakdown = PETSC_FALSE;
169: break;
170: case BV_ORTHOG_REFINE_IFNEEDED:
171: PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,&onorm,&norm,work);
172: /* ||q|| < eta ||h|| */
173: l = 1;
174: while (l<3 && norm < eta * onorm) {
175: l++;
176: onorm = norm;
177: PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,NULL,&norm,work);
178: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
179: }
180: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
181: else *breakdown = PETSC_FALSE;
182: break;
183: }
184: VecScale(v,1.0/norm);
185: VecScale(w,1.0/norm);
187: H[j+1+ldh*j] = norm;
188: if (j<m-1) BVInsertVec(pep->V,j+1,v);
189: }
190: *beta = norm;
191: PetscFunctionReturn(0);
192: }
194: PetscErrorCode PEPSolve_QArnoldi(PEP pep)
195: {
196: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
197: PetscInt j,k,l,lwork,nv,ld,nconv;
198: Vec v=pep->work[0],w=pep->work[1];
199: Mat Q;
200: PetscScalar *S,*work;
201: PetscReal beta=0.0,norm,x,y;
202: PetscBool breakdown=PETSC_FALSE,sinv;
204: DSGetLeadingDimension(pep->ds,&ld);
205: lwork = 7*pep->ncv;
206: PetscMalloc1(lwork,&work);
207: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
208: RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor);
209: STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor);
211: /* Get the starting Arnoldi vector */
212: for (j=0;j<2;j++) {
213: if (j>=pep->nini) BVSetRandomColumn(pep->V,j);
214: }
215: BVCopyVec(pep->V,0,v);
216: BVCopyVec(pep->V,1,w);
217: VecNorm(v,NORM_2,&x);
218: VecNorm(w,NORM_2,&y);
219: norm = PetscSqrtReal(x*x+y*y);
220: VecScale(v,1.0/norm);
221: VecScale(w,1.0/norm);
223: /* clean projected matrix (including the extra-arrow) */
224: DSGetArray(pep->ds,DS_MAT_A,&S);
225: PetscArrayzero(S,ld*ld);
226: DSRestoreArray(pep->ds,DS_MAT_A,&S);
228: /* Restart loop */
229: l = 0;
230: while (pep->reason == PEP_CONVERGED_ITERATING) {
231: pep->its++;
233: /* Compute an nv-step Arnoldi factorization */
234: nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
235: DSGetArray(pep->ds,DS_MAT_A,&S);
236: PEPQArnoldi(pep,S,ld,pep->nconv+l,&nv,v,w,&beta,&breakdown,work);
237: DSRestoreArray(pep->ds,DS_MAT_A,&S);
238: DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l);
239: DSSetState(pep->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE);
240: BVSetActiveColumns(pep->V,pep->nconv,nv);
242: /* Solve projected problem */
243: DSSolve(pep->ds,pep->eigr,pep->eigi);
244: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
245: DSUpdateExtraRow(pep->ds);
246: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
248: /* Check convergence */
249: PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
250: (*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx);
251: nconv = k;
253: /* Update l */
254: if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
255: else {
256: l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
257: DSGetTruncateSize(pep->ds,k,nv,&l);
258: }
259: if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
260: if (l) PetscInfo(pep,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l);
262: if (pep->reason == PEP_CONVERGED_ITERATING) {
263: if (PetscUnlikely(breakdown)) {
264: /* Stop if breakdown */
265: PetscInfo(pep,"Breakdown Quadratic Arnoldi method (it=%" PetscInt_FMT " norm=%g)\n",pep->its,(double)beta);
266: pep->reason = PEP_DIVERGED_BREAKDOWN;
267: } else {
268: /* Prepare the Rayleigh quotient for restart */
269: DSTruncate(pep->ds,k+l,PETSC_FALSE);
270: }
271: }
272: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
273: DSGetMat(pep->ds,DS_MAT_Q,&Q);
274: BVMultInPlace(pep->V,Q,pep->nconv,k+l);
275: MatDestroy(&Q);
277: pep->nconv = k;
278: PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv);
279: }
280: BVSetActiveColumns(pep->V,0,pep->nconv);
281: for (j=0;j<pep->nconv;j++) {
282: pep->eigr[j] *= pep->sfactor;
283: pep->eigi[j] *= pep->sfactor;
284: }
286: STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor);
287: RGPopScale(pep->rg);
289: DSTruncate(pep->ds,pep->nconv,PETSC_TRUE);
290: PetscFree(work);
291: PetscFunctionReturn(0);
292: }
294: static PetscErrorCode PEPQArnoldiSetRestart_QArnoldi(PEP pep,PetscReal keep)
295: {
296: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
298: if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
299: else {
301: ctx->keep = keep;
302: }
303: PetscFunctionReturn(0);
304: }
306: /*@
307: PEPQArnoldiSetRestart - Sets the restart parameter for the Q-Arnoldi
308: method, in particular the proportion of basis vectors that must be kept
309: after restart.
311: Logically Collective on pep
313: Input Parameters:
314: + pep - the eigenproblem solver context
315: - keep - the number of vectors to be kept at restart
317: Options Database Key:
318: . -pep_qarnoldi_restart - Sets the restart parameter
320: Notes:
321: Allowed values are in the range [0.1,0.9]. The default is 0.5.
323: Level: advanced
325: .seealso: PEPQArnoldiGetRestart()
326: @*/
327: PetscErrorCode PEPQArnoldiSetRestart(PEP pep,PetscReal keep)
328: {
331: PetscTryMethod(pep,"PEPQArnoldiSetRestart_C",(PEP,PetscReal),(pep,keep));
332: PetscFunctionReturn(0);
333: }
335: static PetscErrorCode PEPQArnoldiGetRestart_QArnoldi(PEP pep,PetscReal *keep)
336: {
337: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
339: *keep = ctx->keep;
340: PetscFunctionReturn(0);
341: }
343: /*@
344: PEPQArnoldiGetRestart - Gets the restart parameter used in the Q-Arnoldi method.
346: Not Collective
348: Input Parameter:
349: . pep - the eigenproblem solver context
351: Output Parameter:
352: . keep - the restart parameter
354: Level: advanced
356: .seealso: PEPQArnoldiSetRestart()
357: @*/
358: PetscErrorCode PEPQArnoldiGetRestart(PEP pep,PetscReal *keep)
359: {
362: PetscUseMethod(pep,"PEPQArnoldiGetRestart_C",(PEP,PetscReal*),(pep,keep));
363: PetscFunctionReturn(0);
364: }
366: static PetscErrorCode PEPQArnoldiSetLocking_QArnoldi(PEP pep,PetscBool lock)
367: {
368: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
370: ctx->lock = lock;
371: PetscFunctionReturn(0);
372: }
374: /*@
375: PEPQArnoldiSetLocking - Choose between locking and non-locking variants of
376: the Q-Arnoldi method.
378: Logically Collective on pep
380: Input Parameters:
381: + pep - the eigenproblem solver context
382: - lock - true if the locking variant must be selected
384: Options Database Key:
385: . -pep_qarnoldi_locking - Sets the locking flag
387: Notes:
388: The default is to lock converged eigenpairs when the method restarts.
389: This behaviour can be changed so that all directions are kept in the
390: working subspace even if already converged to working accuracy (the
391: non-locking variant).
393: Level: advanced
395: .seealso: PEPQArnoldiGetLocking()
396: @*/
397: PetscErrorCode PEPQArnoldiSetLocking(PEP pep,PetscBool lock)
398: {
401: PetscTryMethod(pep,"PEPQArnoldiSetLocking_C",(PEP,PetscBool),(pep,lock));
402: PetscFunctionReturn(0);
403: }
405: static PetscErrorCode PEPQArnoldiGetLocking_QArnoldi(PEP pep,PetscBool *lock)
406: {
407: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
409: *lock = ctx->lock;
410: PetscFunctionReturn(0);
411: }
413: /*@
414: PEPQArnoldiGetLocking - Gets the locking flag used in the Q-Arnoldi method.
416: Not Collective
418: Input Parameter:
419: . pep - the eigenproblem solver context
421: Output Parameter:
422: . lock - the locking flag
424: Level: advanced
426: .seealso: PEPQArnoldiSetLocking()
427: @*/
428: PetscErrorCode PEPQArnoldiGetLocking(PEP pep,PetscBool *lock)
429: {
432: PetscUseMethod(pep,"PEPQArnoldiGetLocking_C",(PEP,PetscBool*),(pep,lock));
433: PetscFunctionReturn(0);
434: }
436: PetscErrorCode PEPSetFromOptions_QArnoldi(PetscOptionItems *PetscOptionsObject,PEP pep)
437: {
438: PetscBool flg,lock;
439: PetscReal keep;
441: PetscOptionsHead(PetscOptionsObject,"PEP Q-Arnoldi Options");
443: PetscOptionsReal("-pep_qarnoldi_restart","Proportion of vectors kept after restart","PEPQArnoldiSetRestart",0.5,&keep,&flg);
444: if (flg) PEPQArnoldiSetRestart(pep,keep);
446: PetscOptionsBool("-pep_qarnoldi_locking","Choose between locking and non-locking variants","PEPQArnoldiSetLocking",PETSC_FALSE,&lock,&flg);
447: if (flg) PEPQArnoldiSetLocking(pep,lock);
449: PetscOptionsTail();
450: PetscFunctionReturn(0);
451: }
453: PetscErrorCode PEPView_QArnoldi(PEP pep,PetscViewer viewer)
454: {
455: PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;
456: PetscBool isascii;
458: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
459: if (isascii) {
460: PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
461: PetscViewerASCIIPrintf(viewer," using the %slocking variant\n",ctx->lock?"":"non-");
462: }
463: PetscFunctionReturn(0);
464: }
466: PetscErrorCode PEPDestroy_QArnoldi(PEP pep)
467: {
468: PetscFree(pep->data);
469: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",NULL);
470: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",NULL);
471: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",NULL);
472: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",NULL);
473: PetscFunctionReturn(0);
474: }
476: SLEPC_EXTERN PetscErrorCode PEPCreate_QArnoldi(PEP pep)
477: {
478: PEP_QARNOLDI *ctx;
480: PetscNewLog(pep,&ctx);
481: pep->data = (void*)ctx;
483: pep->lineariz = PETSC_TRUE;
484: ctx->lock = PETSC_TRUE;
486: pep->ops->solve = PEPSolve_QArnoldi;
487: pep->ops->setup = PEPSetUp_QArnoldi;
488: pep->ops->setfromoptions = PEPSetFromOptions_QArnoldi;
489: pep->ops->destroy = PEPDestroy_QArnoldi;
490: pep->ops->view = PEPView_QArnoldi;
491: pep->ops->backtransform = PEPBackTransform_Default;
492: pep->ops->computevectors = PEPComputeVectors_Default;
493: pep->ops->extractvectors = PEPExtractVectors_QArnoldi;
494: pep->ops->setdefaultst = PEPSetDefaultST_Transform;
496: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",PEPQArnoldiSetRestart_QArnoldi);
497: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",PEPQArnoldiGetRestart_QArnoldi);
498: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetLocking_C",PEPQArnoldiSetLocking_QArnoldi);
499: PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetLocking_C",PEPQArnoldiGetLocking_QArnoldi);
500: PetscFunctionReturn(0);
501: }