Actual source code: qarnoldi.c
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: Last update: Nov 2012
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: SLEPc - Scalable Library for Eigenvalue Problem Computations
21: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
23: This file is part of SLEPc.
25: SLEPc is free software: you can redistribute it and/or modify it under the
26: terms of version 3 of the GNU Lesser General Public License as published by
27: the Free Software Foundation.
29: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
30: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
31: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
32: more details.
34: You should have received a copy of the GNU Lesser General Public License
35: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
36: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: */
39: #include <slepc-private/qepimpl.h> /*I "slepcqep.h" I*/
40: #include <petscblaslapack.h>
44: PetscErrorCode QEPSetUp_QArnoldi(QEP qep)
45: {
47: PetscBool sinv;
50: if (qep->ncv) { /* ncv set */
51: if (qep->ncv<qep->nev) SETERRQ(PetscObjectComm((PetscObject)qep),1,"The value of ncv must be at least nev");
52: } else if (qep->mpd) { /* mpd set */
53: qep->ncv = PetscMin(qep->n,qep->nev+qep->mpd);
54: } else { /* neither set: defaults depend on nev being small or large */
55: if (qep->nev<500) qep->ncv = PetscMin(qep->n,PetscMax(2*qep->nev,qep->nev+15));
56: else {
57: qep->mpd = 500;
58: qep->ncv = PetscMin(qep->n,qep->nev+qep->mpd);
59: }
60: }
61: if (!qep->mpd) qep->mpd = qep->ncv;
62: if (qep->ncv>qep->nev+qep->mpd) SETERRQ(PetscObjectComm((PetscObject)qep),1,"The value of ncv must not be larger than nev+mpd");
63: if (!qep->max_it) qep->max_it = PetscMax(100,2*qep->n/qep->ncv);
64: if (!qep->which) {
65: PetscObjectTypeCompare((PetscObject)qep->st,STSINVERT,&sinv);
66: if (sinv) qep->which = QEP_TARGET_MAGNITUDE;
67: else qep->which = QEP_LARGEST_MAGNITUDE;
68: }
70: QEPAllocateSolution(qep);
71: QEPSetWorkVecs(qep,4);
73: DSSetType(qep->ds,DSNHEP);
74: DSSetExtraRow(qep->ds,PETSC_TRUE);
75: DSAllocate(qep->ds,qep->ncv+1);
76: return(0);
77: }
81: /*
82: Compute a step of Classical Gram-Schmidt orthogonalization
83: */
84: static PetscErrorCode QEPQArnoldiCGS(QEP qep,PetscScalar *H,PetscBLASInt ldh,PetscScalar *h,PetscBLASInt j,Vec *V,Vec t,Vec v,Vec w,PetscReal *onorm,PetscReal *norm,PetscScalar *work)
85: {
87: PetscBLASInt ione = 1,j_1 = j+1;
88: PetscReal x,y;
89: PetscScalar dot,one = 1.0,zero = 0.0;
92: /* compute norm of v and w */
93: if (onorm) {
94: VecNorm(v,NORM_2,&x);
95: VecNorm(w,NORM_2,&y);
96: *onorm = PetscSqrtReal(x*x+y*y);
97: }
99: /* orthogonalize: compute h */
100: VecMDot(v,j_1,V,h);
101: VecMDot(w,j_1,V,work);
102: if (j>0)
103: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione));
104: VecDot(w,t,&dot);
105: h[j] += dot;
107: /* orthogonalize: update v and w */
108: SlepcVecMAXPBY(v,1.0,-1.0,j_1,h,V);
109: if (j>0) {
110: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione));
111: SlepcVecMAXPBY(w,1.0,-1.0,j_1,work,V);
112: }
113: VecAXPY(w,-h[j],t);
115: /* compute norm of v and w */
116: if (norm) {
117: VecNorm(v,NORM_2,&x);
118: VecNorm(w,NORM_2,&y);
119: *norm = PetscSqrtReal(x*x+y*y);
120: }
121: return(0);
122: }
126: /*
127: Compute a run of Q-Arnoldi iterations
128: */
129: static PetscErrorCode QEPQArnoldi(QEP qep,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
130: {
131: PetscErrorCode ierr;
132: PetscInt i,j,l,m = *M;
133: Vec t = qep->work[2],u = qep->work[3];
134: IPOrthogRefineType refinement;
135: PetscReal norm,onorm,eta;
136: PetscScalar *c = work + m;
139: IPGetOrthogonalization(qep->ip,NULL,&refinement,&eta);
140: VecCopy(v,qep->V[k]);
142: for (j=k;j<m;j++) {
143: /* apply operator */
144: VecCopy(w,t);
145: STMatMult(qep->st,0,v,u);
146: STMatMult(qep->st,1,t,w);
147: VecAXPY(u,qep->sfactor,w);
148: STMatSolve(qep->st,2,u,w);
149: VecScale(w,-1.0/(qep->sfactor*qep->sfactor));
150: VecCopy(t,v);
152: /* orthogonalize */
153: switch (refinement) {
154: case IP_ORTHOG_REFINE_NEVER:
155: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,NULL,&norm,work);
156: *breakdown = PETSC_FALSE;
157: break;
158: case IP_ORTHOG_REFINE_ALWAYS:
159: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,NULL,NULL,work);
160: QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,&onorm,&norm,work);
161: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
162: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
163: else *breakdown = PETSC_FALSE;
164: break;
165: case IP_ORTHOG_REFINE_IFNEEDED:
166: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,&onorm,&norm,work);
167: /* ||q|| < eta ||h|| */
168: l = 1;
169: while (l<3 && norm < eta * onorm) {
170: l++;
171: onorm = norm;
172: QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,NULL,&norm,work);
173: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
174: }
175: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
176: else *breakdown = PETSC_FALSE;
177: break;
178: default: SETERRQ(PetscObjectComm((PetscObject)qep),1,"Wrong value of ip->orth_ref");
179: }
180: VecScale(v,1.0/norm);
181: VecScale(w,1.0/norm);
183: H[j+1+ldh*j] = norm;
184: if (j<m-1) {
185: VecCopy(v,V[j+1]);
186: }
187: }
188: *beta = norm;
189: return(0);
190: }
194: PetscErrorCode QEPSolve_QArnoldi(QEP qep)
195: {
197: PetscInt j,k,l,lwork,nv,ld,newn;
198: Vec v=qep->work[0],w=qep->work[1],v_=qep->work[2],w_=qep->work[3];
199: PetscScalar *S,*Q,*work,r,s;
200: PetscReal beta=0.0,norm,x,y,t;
201: PetscBool breakdown=PETSC_FALSE,issinv;
204: DSGetLeadingDimension(qep->ds,&ld);
205: lwork = 7*qep->ncv;
206: PetscMalloc(lwork*sizeof(PetscScalar),&work);
208: /* Get the starting Arnoldi vector */
209: if (qep->nini>0) {
210: VecCopy(qep->V[0],v);
211: } else {
212: SlepcVecSetRandom(v,qep->rand);
213: }
214: /* w is always a random vector */
215: SlepcVecSetRandom(w,qep->rand);
216: VecNorm(v,NORM_2,&x);
217: VecNorm(w,NORM_2,&y);
218: norm = PetscSqrtReal(x*x+y*y);
219: VecScale(v,1.0/norm);
220: VecScale(w,1.0/norm);
222: /* Compute scaling factor if not set by user */
223: PetscObjectTypeCompare((PetscObject)qep->st,STSINVERT,&issinv);
224: if (issinv && !qep->sfactor_set) {
225: STMatMult(qep->st,1,w,w_);
226: STMatMult(qep->st,0,v,v_);
227: VecAXPY(v_,1.0,w_);
228: STMatSolve(qep->st,2,v_,w_);
229: VecScale(w_,-1.0);
230: VecCopy(w,v_);
231: VecDot(v_,v,&r);
232: VecDot(w_,w,&s);
233: t = PetscAbsScalar(r+s);
234: qep->sfactor = 1.0;
235: while (t > 1.0) {
236: qep->sfactor *=10.0;
237: t /= 10.0;
238: }
239: }
240: /* Restart loop */
241: l = 0;
242: while (qep->reason == QEP_CONVERGED_ITERATING) {
243: qep->its++;
245: /* Compute an nv-step Arnoldi factorization */
246: nv = PetscMin(qep->nconv+qep->mpd,qep->ncv);
247: DSGetArray(qep->ds,DS_MAT_A,&S);
248: QEPQArnoldi(qep,S,ld,qep->V,qep->nconv+l,&nv,v,w,&beta,&breakdown,work);
249: DSRestoreArray(qep->ds,DS_MAT_A,&S);
250: DSSetDimensions(qep->ds,nv,0,qep->nconv,qep->nconv+l);
251: if (l==0) {
252: DSSetState(qep->ds,DS_STATE_INTERMEDIATE);
253: } else {
254: DSSetState(qep->ds,DS_STATE_RAW);
255: }
257: /* Solve projected problem */
258: DSSolve(qep->ds,qep->eigr,qep->eigi);
259: DSSort(qep->ds,qep->eigr,qep->eigi,NULL,NULL,NULL);
260: DSUpdateExtraRow(qep->ds);
262: /* Check convergence */
263: QEPKrylovConvergence(qep,PETSC_FALSE,qep->nconv,nv-qep->nconv,nv,beta,&k);
264: if (qep->its >= qep->max_it) qep->reason = QEP_DIVERGED_ITS;
265: if (k >= qep->nev) qep->reason = QEP_CONVERGED_TOL;
267: /* Update l */
268: if (qep->reason != QEP_CONVERGED_ITERATING || breakdown) l = 0;
269: else l = (nv-k)/2;
271: if (qep->reason == QEP_CONVERGED_ITERATING) {
272: if (breakdown) {
273: /* Stop if breakdown */
274: PetscInfo2(qep,"Breakdown Quadratic Arnoldi method (it=%D norm=%G)\n",qep->its,beta);
275: qep->reason = QEP_DIVERGED_BREAKDOWN;
276: } else {
277: /* Prepare the Rayleigh quotient for restart */
278: DSTruncate(qep->ds,k+l);
279: DSGetDimensions(qep->ds,&newn,NULL,NULL,NULL,NULL);
280: l = newn-k;
281: }
282: }
283: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
284: DSGetArray(qep->ds,DS_MAT_Q,&Q);
285: SlepcUpdateVectors(nv,qep->V,qep->nconv,k+l,Q,ld,PETSC_FALSE);
286: DSRestoreArray(qep->ds,DS_MAT_Q,&Q);
288: qep->nconv = k;
289: QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,nv);
290: }
292: for (j=0;j<qep->nconv;j++) {
293: qep->eigr[j] *= qep->sfactor;
294: qep->eigi[j] *= qep->sfactor;
295: }
297: /* truncate Schur decomposition and change the state to raw so that
298: DSVectors() computes eigenvectors from scratch */
299: DSSetDimensions(qep->ds,qep->nconv,0,0,0);
300: DSSetState(qep->ds,DS_STATE_RAW);
302: /* Compute eigenvectors */
303: if (qep->nconv > 0) {
304: QEPComputeVectors_Schur(qep);
305: }
306: PetscFree(work);
307: return(0);
308: }
312: PETSC_EXTERN PetscErrorCode QEPCreate_QArnoldi(QEP qep)
313: {
315: qep->ops->solve = QEPSolve_QArnoldi;
316: qep->ops->setup = QEPSetUp_QArnoldi;
317: qep->ops->reset = QEPReset_Default;
318: return(0);
319: }