Actual source code: qarnoldi.c
1: /*
3: Q-Arnoldi method for quadratic eigenproblems.
5: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6: SLEPc - Scalable Library for Eigenvalue Problem Computations
7: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
9: This file is part of SLEPc.
10:
11: SLEPc is free software: you can redistribute it and/or modify it under the
12: terms of version 3 of the GNU Lesser General Public License as published by
13: the Free Software Foundation.
15: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
16: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
18: more details.
20: You should have received a copy of the GNU Lesser General Public License
21: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: */
25: #include private/qepimpl.h
26: #include private/epsimpl.h
27: #include "petscblaslapack.h"
29: typedef struct {
30: KSP ksp;
31: } QEP_QARNOLDI;
35: PetscErrorCode QEPSetUp_QARNOLDI(QEP qep)
36: {
38: QEP_QARNOLDI *ctx = (QEP_QARNOLDI *)qep->data;
39:
42: if (qep->ncv) { /* ncv set */
43: if (qep->ncv<qep->nev) SETERRQ(1,"The value of ncv must be at least nev");
44: }
45: else if (qep->mpd) { /* mpd set */
46: qep->ncv = PetscMin(qep->n,qep->nev+qep->mpd);
47: }
48: else { /* neither set: defaults depend on nev being small or large */
49: if (qep->nev<500) qep->ncv = PetscMin(qep->n,PetscMax(2*qep->nev,qep->nev+15));
50: else { qep->mpd = 500; qep->ncv = PetscMin(qep->n,qep->nev+qep->mpd); }
51: }
52: if (!qep->mpd) qep->mpd = qep->ncv;
53: if (qep->ncv>qep->nev+qep->mpd) SETERRQ(1,"The value of ncv must not be larger than nev+mpd");
54: if (!qep->max_it) qep->max_it = PetscMax(100,2*qep->n/qep->ncv);
55: if (!qep->which) qep->which = QEP_LARGEST_MAGNITUDE;
56: if (qep->problem_type != QEP_GENERAL)
57: SETERRQ(1,"Wrong value of qep->problem_type");
59: PetscFree(qep->T);
60: PetscMalloc(qep->ncv*qep->ncv*sizeof(PetscScalar),&qep->T);
61: QEPDefaultGetWork(qep,4);
63: KSPSetOperators(ctx->ksp,qep->M,qep->M,DIFFERENT_NONZERO_PATTERN);
64: KSPSetUp(ctx->ksp);
66: return(0);
67: }
71: /*
72: Compute a step of Classical Gram-Schmidt orthogonalization
73: */
74: 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)
75: {
77: PetscBLASInt ione = 1, j_1 = j+1;
78: PetscReal x, y;
79: PetscScalar dot, one = 1.0, zero = 0.0;
82: /* compute norm of v and w */
83: if (onorm) {
84: VecNorm(v,NORM_2,&x);
85: VecNorm(w,NORM_2,&y);
86: *onorm = sqrt(x*x+y*y);
87: }
89: /* orthogonalize: compute h */
90: VecMDot(v,j_1,V,h);
91: VecMDot(w,j_1,V,work);
92: if (j>0)
93: BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione);
94: VecDot(t,w,&dot);
95: h[j] += dot;
97: /* orthogonalize: update v and w */
98: SlepcVecMAXPBY(v,1.0,-1.0,j_1,h,V);
99: if (j>0) {
100: BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione);
101: SlepcVecMAXPBY(w,1.0,-1.0,j_1,work,V);
102: }
103: VecAXPY(w,-h[j],t);
104:
105: /* compute norm of v and w */
106: if (norm) {
107: VecNorm(v,NORM_2,&x);
108: VecNorm(w,NORM_2,&y);
109: *norm = sqrt(x*x+y*y);
110: }
111: return(0);
112: }
116: /*
117: Compute a run of Q-Arnoldi iterations
118: */
119: PetscErrorCode QEPQArnoldi(QEP qep,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscTruth *breakdown,PetscScalar *work)
120: {
122: PetscInt i,j,l,m = *M;
123: QEP_QARNOLDI *ctx = (QEP_QARNOLDI *)qep->data;
124: Vec t = qep->work[2], u = qep->work[3];
125: IPOrthogonalizationRefinementType refinement;
126: PetscReal norm,onorm,eta;
127: PetscScalar *c = work + m;
130: IPGetOrthogonalization(qep->ip,PETSC_NULL,&refinement,&eta);
131: VecCopy(v,qep->V[k]);
132:
133: for (j=k;j<m;j++) {
134: /* apply operator */
135: VecCopy(w,t);
136: MatMult(qep->K,v,u);
137: MatMult(qep->C,t,w);
138: VecAXPY(u,qep->sfactor,w);
139: KSPSolve(ctx->ksp,u,w);
140: VecScale(w,-1.0/(qep->sfactor*qep->sfactor));
141: VecCopy(t,v);
143: /* orthogonalize */
144: switch (refinement) {
145: case IP_ORTH_REFINE_NEVER:
146: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,PETSC_NULL,&norm,work);
147: *breakdown = PETSC_FALSE;
148: break;
149: case IP_ORTH_REFINE_ALWAYS:
150: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,PETSC_NULL,PETSC_NULL,work);
151: QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,&onorm,&norm,work);
152: for (i=0;i<j;i++) H[ldh*j+i] += c[i];
153: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
154: else *breakdown = PETSC_FALSE;
155: break;
156: case IP_ORTH_REFINE_IFNEEDED:
157: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,&onorm,&norm,work);
158: /* ||q|| < eta ||h|| */
159: l = 1;
160: while (l<3 && norm < eta * onorm) {
161: l++;
162: onorm = norm;
163: QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,PETSC_NULL,&norm,work);
164: for (i=0;i<j;i++) H[ldh*j+i] += c[i];
165: }
166: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
167: else *breakdown = PETSC_FALSE;
168: break;
169: default: SETERRQ(1,"Wrong value of ip->orth_ref");
170: }
171: VecScale(v,1.0/norm);
172: VecScale(w,1.0/norm);
173:
174: if (j<m-1) {
175: H[j+1+ldh*j] = norm;
176: VecCopy(v,V[j+1]);
177: }
178: }
179: *beta = norm;
180: return(0);
181: }
186: /*
187: QEPProjectedKSNonsym - Solves the projected eigenproblem in the Krylov-Schur
188: method (non-symmetric case).
190: On input:
191: l is the number of vectors kept in previous restart (0 means first restart)
192: S is the projected matrix (leading dimension is lds)
194: On output:
195: S has (real) Schur form with diagonal blocks sorted appropriately
196: Q contains the corresponding Schur vectors (order n, leading dimension n)
197: */
198: PetscErrorCode QEPProjectedKSNonsym(QEP qep,PetscInt l,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt n)
199: {
201: PetscInt i;
204: if (l==0) {
205: PetscMemzero(Q,n*n*sizeof(PetscScalar));
206: for (i=0;i<n;i++)
207: Q[i*(n+1)] = 1.0;
208: } else {
209: /* Reduce S to Hessenberg form, S <- Q S Q' */
210: EPSDenseHessenberg(n,qep->nconv,S,lds,Q);
211: }
212: /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
213: EPSDenseSchur(n,qep->nconv,S,lds,Q,qep->eigr,qep->eigi);
214: /* Sort the remaining columns of the Schur form */
215: QEPSortDenseSchur(qep,n,qep->nconv,S,lds,Q,qep->eigr,qep->eigi);
216: return(0);
217: }
221: PetscErrorCode QEPSolve_QARNOLDI(QEP qep)
222: {
224: PetscInt i,j,k,l,lwork,nv;
225: Vec v=qep->work[0],w=qep->work[1];
226: PetscScalar *S=qep->T,*Q,*work;
227: PetscReal beta,norm,x,y;
228: PetscTruth breakdown;
232: PetscMemzero(S,qep->ncv*qep->ncv*sizeof(PetscScalar));
233: PetscMalloc(qep->ncv*qep->ncv*sizeof(PetscScalar),&Q);
234: lwork = 7*qep->ncv;
235: PetscMalloc(lwork*sizeof(PetscScalar),&work);
237: /* Get the starting Arnoldi vector */
238: if (qep->nini>0) {
239: VecCopy(qep->V[0],v);
240: } else {
241: SlepcVecSetRandom(v,qep->rand);
242: }
243: /* w is always a random vector */
244: SlepcVecSetRandom(w,qep->rand);
245: VecNorm(v,NORM_2,&x);
246: VecNorm(w,NORM_2,&y);
247: norm = sqrt(x*x+y*y);
248: VecScale(v,1.0/norm);
249: VecScale(w,1.0/norm);
250:
251: /* Restart loop */
252: l = 0;
253: while (qep->reason == QEP_CONVERGED_ITERATING) {
254: qep->its++;
256: /* Compute an nv-step Arnoldi factorization */
257: nv = PetscMin(qep->nconv+qep->mpd,qep->ncv);
258: QEPQArnoldi(qep,S,qep->ncv,qep->V,qep->nconv+l,&nv,v,w,&beta,&breakdown,work);
260: /* Solve projected problem */
261: QEPProjectedKSNonsym(qep,l,S,qep->ncv,Q,nv);
263: /* Check convergence */
264: QEPKrylovConvergence(qep,qep->nconv,nv-qep->nconv,S,qep->ncv,Q,nv,beta,&k,work);
265: if (qep->its >= qep->max_it) qep->reason = QEP_DIVERGED_ITS;
266: if (k >= qep->nev) qep->reason = QEP_CONVERGED_TOL;
267:
268: /* Update l */
269: if (qep->reason != QEP_CONVERGED_ITERATING || breakdown) l = 0;
270: else {
271: l = (nv-k)/2;
272: #if !defined(PETSC_USE_COMPLEX)
273: if (S[(k+l-1)*(qep->ncv+1)+1] != 0.0) {
274: if (k+l<nv-1) l = l+1;
275: else l = l-1;
276: }
277: #endif
278: }
280: if (qep->reason == QEP_CONVERGED_ITERATING) {
281: if (breakdown) {
282: /* Stop if breakdown */
283: PetscInfo2(qep,"Breakdown Quadratic Arnoldi method (it=%i norm=%g)\n",qep->its,beta);
284: qep->reason = QEP_DIVERGED_BREAKDOWN;
285: } else {
286: /* Prepare the Rayleigh quotient for restart */
287: for (i=k;i<k+l;i++) {
288: S[i*qep->ncv+k+l] = Q[(i+1)*nv-1]*beta;
289: }
290: }
291: }
292: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
293: SlepcUpdateVectors(nv,qep->V,qep->nconv,k+l,Q,nv,PETSC_FALSE);
295: qep->nconv = k;
297: QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,nv);
298:
299: }
301: for (j=0;j<qep->nconv;j++) {
302: qep->eigr[j] *= qep->sfactor;
303: qep->eigi[j] *= qep->sfactor;
304: }
306: /* Compute eigenvectors */
307: if (qep->nconv > 0) {
308: QEPComputeVectors_Schur(qep);
309: }
311: PetscFree(Q);
312: PetscFree(work);
313: return(0);
314: }
318: PetscErrorCode QEPSetFromOptions_QARNOLDI(QEP qep)
319: {
321: QEP_QARNOLDI *ctx = (QEP_QARNOLDI *)qep->data;
322:
324: KSPSetFromOptions(ctx->ksp);
325: return(0);
326: }
330: PetscErrorCode QEPView_QARNOLDI(QEP qep,PetscViewer viewer)
331: {
333: QEP_QARNOLDI *ctx = (QEP_QARNOLDI *)qep->data;
336: KSPView(ctx->ksp,viewer);
337: return(0);
338: }
342: PetscErrorCode QEPDestroy_QARNOLDI(QEP qep)
343: {
345: QEP_QARNOLDI *ctx = (QEP_QARNOLDI *)qep->data;
348: KSPDestroy(ctx->ksp);
349: QEPDestroy_Default(qep);
350: return(0);
351: }
356: PetscErrorCode QEPCreate_QARNOLDI(QEP qep)
357: {
359: QEP_QARNOLDI *ctx;
362: PetscNew(QEP_QARNOLDI,&ctx);
363: PetscLogObjectMemory(qep,sizeof(QEP_QARNOLDI));
364: qep->data = ctx;
365: qep->ops->solve = QEPSolve_QARNOLDI;
366: qep->ops->setup = QEPSetUp_QARNOLDI;
367: qep->ops->setfromoptions = QEPSetFromOptions_QARNOLDI;
368: qep->ops->destroy = QEPDestroy_QARNOLDI;
369: qep->ops->view = QEPView_QARNOLDI;
371: KSPCreate(((PetscObject)qep)->comm,&ctx->ksp);
372: KSPSetOptionsPrefix(ctx->ksp,((PetscObject)qep)->prefix);
373: KSPAppendOptionsPrefix(ctx->ksp,"qep_");
374: PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)qep,1);
375: PetscLogObjectParent(qep,ctx->ksp);
376: return(0);
377: }