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-2012, Universitat 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 <slepc-private/qepimpl.h> /*I "slepcqep.h" I*/
26: #include <petscblaslapack.h>
28: typedef struct {
29: KSP ksp;
30: } QEP_QARNOLDI;
34: PetscErrorCode QEPSetUp_QArnoldi(QEP qep)
35: {
37: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
38:
40: if (qep->ncv) { /* ncv set */
41: if (qep->ncv<qep->nev) SETERRQ(((PetscObject)qep)->comm,1,"The value of ncv must be at least nev");
42: }
43: else if (qep->mpd) { /* mpd set */
44: qep->ncv = PetscMin(qep->n,qep->nev+qep->mpd);
45: }
46: else { /* neither set: defaults depend on nev being small or large */
47: if (qep->nev<500) qep->ncv = PetscMin(qep->n,PetscMax(2*qep->nev,qep->nev+15));
48: else { qep->mpd = 500; qep->ncv = PetscMin(qep->n,qep->nev+qep->mpd); }
49: }
50: if (!qep->mpd) qep->mpd = qep->ncv;
51: if (qep->ncv>qep->nev+qep->mpd) SETERRQ(((PetscObject)qep)->comm,1,"The value of ncv must not be larger than nev+mpd");
52: if (!qep->max_it) qep->max_it = PetscMax(100,2*qep->n/qep->ncv);
53: if (!qep->which) qep->which = QEP_LARGEST_MAGNITUDE;
55: QEPAllocateSolution(qep);
56: QEPDefaultGetWork(qep,4);
58: DSSetType(qep->ds,DSNHEP);
59: DSSetExtraRow(qep->ds,PETSC_TRUE);
60: DSAllocate(qep->ds,qep->ncv+1);
62: KSPSetOperators(ctx->ksp,qep->M,qep->M,DIFFERENT_NONZERO_PATTERN);
63: KSPSetUp(ctx->ksp);
64: return(0);
65: }
69: /*
70: Compute a step of Classical Gram-Schmidt orthogonalization
71: */
72: 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)
73: {
75: PetscBLASInt ione = 1,j_1 = j+1;
76: PetscReal x,y;
77: PetscScalar dot,one = 1.0,zero = 0.0;
80: /* compute norm of v and w */
81: if (onorm) {
82: VecNorm(v,NORM_2,&x);
83: VecNorm(w,NORM_2,&y);
84: *onorm = PetscSqrtReal(x*x+y*y);
85: }
87: /* orthogonalize: compute h */
88: VecMDot(v,j_1,V,h);
89: VecMDot(w,j_1,V,work);
90: if (j>0)
91: BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione);
92: VecDot(w,t,&dot);
93: h[j] += dot;
95: /* orthogonalize: update v and w */
96: SlepcVecMAXPBY(v,1.0,-1.0,j_1,h,V);
97: if (j>0) {
98: BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione);
99: SlepcVecMAXPBY(w,1.0,-1.0,j_1,work,V);
100: }
101: VecAXPY(w,-h[j],t);
102:
103: /* compute norm of v and w */
104: if (norm) {
105: VecNorm(v,NORM_2,&x);
106: VecNorm(w,NORM_2,&y);
107: *norm = PetscSqrtReal(x*x+y*y);
108: }
109: return(0);
110: }
114: /*
115: Compute a run of Q-Arnoldi iterations
116: */
117: PetscErrorCode QEPQArnoldi(QEP qep,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
118: {
119: PetscErrorCode ierr;
120: PetscInt i,j,l,m = *M;
121: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
122: Vec t = qep->work[2],u = qep->work[3];
123: IPOrthogRefineType refinement;
124: PetscReal norm,onorm,eta;
125: PetscScalar *c = work + m;
128: IPGetOrthogonalization(qep->ip,PETSC_NULL,&refinement,&eta);
129: VecCopy(v,qep->V[k]);
130:
131: for (j=k;j<m;j++) {
132: /* apply operator */
133: VecCopy(w,t);
134: MatMult(qep->K,v,u);
135: MatMult(qep->C,t,w);
136: VecAXPY(u,qep->sfactor,w);
137: KSPSolve(ctx->ksp,u,w);
138: VecScale(w,-1.0/(qep->sfactor*qep->sfactor));
139: VecCopy(t,v);
141: /* orthogonalize */
142: switch (refinement) {
143: case IP_ORTHOG_REFINE_NEVER:
144: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,PETSC_NULL,&norm,work);
145: *breakdown = PETSC_FALSE;
146: break;
147: case IP_ORTHOG_REFINE_ALWAYS:
148: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,PETSC_NULL,PETSC_NULL,work);
149: QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,&onorm,&norm,work);
150: for (i=0;i<j;i++) H[ldh*j+i] += c[i];
151: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
152: else *breakdown = PETSC_FALSE;
153: break;
154: case IP_ORTHOG_REFINE_IFNEEDED:
155: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,&onorm,&norm,work);
156: /* ||q|| < eta ||h|| */
157: l = 1;
158: while (l<3 && norm < eta * onorm) {
159: l++;
160: onorm = norm;
161: QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,PETSC_NULL,&norm,work);
162: for (i=0;i<j;i++) H[ldh*j+i] += c[i];
163: }
164: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
165: else *breakdown = PETSC_FALSE;
166: break;
167: default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of ip->orth_ref");
168: }
169: VecScale(v,1.0/norm);
170: VecScale(w,1.0/norm);
171:
172: H[j+1+ldh*j] = norm;
173: if (j<m-1) {
174: VecCopy(v,V[j+1]);
175: }
176: }
177: *beta = norm;
178: return(0);
179: }
183: PetscErrorCode QEPSolve_QArnoldi(QEP qep)
184: {
186: PetscInt j,k,l,lwork,nv,ld,newn;
187: Vec v=qep->work[0],w=qep->work[1];
188: PetscScalar *S,*Q,*work;
189: PetscReal beta,norm,x,y;
190: PetscBool breakdown;
193: DSGetLeadingDimension(qep->ds,&ld);
194: lwork = 7*qep->ncv;
195: PetscMalloc(lwork*sizeof(PetscScalar),&work);
197: /* Get the starting Arnoldi vector */
198: if (qep->nini>0) {
199: VecCopy(qep->V[0],v);
200: } else {
201: SlepcVecSetRandom(v,qep->rand);
202: }
203: /* w is always a random vector */
204: SlepcVecSetRandom(w,qep->rand);
205: VecNorm(v,NORM_2,&x);
206: VecNorm(w,NORM_2,&y);
207: norm = PetscSqrtReal(x*x+y*y);
208: VecScale(v,1.0/norm);
209: VecScale(w,1.0/norm);
210:
211: /* Restart loop */
212: l = 0;
213: while (qep->reason == QEP_CONVERGED_ITERATING) {
214: qep->its++;
216: /* Compute an nv-step Arnoldi factorization */
217: nv = PetscMin(qep->nconv+qep->mpd,qep->ncv);
218: DSGetArray(qep->ds,DS_MAT_A,&S);
219: QEPQArnoldi(qep,S,ld,qep->V,qep->nconv+l,&nv,v,w,&beta,&breakdown,work);
220: DSRestoreArray(qep->ds,DS_MAT_A,&S);
221: DSSetDimensions(qep->ds,nv,PETSC_IGNORE,qep->nconv,qep->nconv+l);
222: if (l==0) {
223: DSSetState(qep->ds,DS_STATE_INTERMEDIATE);
224: } else {
225: DSSetState(qep->ds,DS_STATE_RAW);
226: }
228: /* Solve projected problem */
229: DSSolve(qep->ds,qep->eigr,qep->eigi);
230: DSSort(qep->ds,qep->eigr,qep->eigi,PETSC_NULL,PETSC_NULL,PETSC_NULL);
231: DSUpdateExtraRow(qep->ds);
233: /* Check convergence */
234: QEPKrylovConvergence(qep,PETSC_FALSE,qep->nconv,nv-qep->nconv,nv,beta,&k);
235: if (qep->its >= qep->max_it) qep->reason = QEP_DIVERGED_ITS;
236: if (k >= qep->nev) qep->reason = QEP_CONVERGED_TOL;
237:
238: /* Update l */
239: if (qep->reason != QEP_CONVERGED_ITERATING || breakdown) l = 0;
240: else l = (nv-k)/2;
242: if (qep->reason == QEP_CONVERGED_ITERATING) {
243: if (breakdown) {
244: /* Stop if breakdown */
245: PetscInfo2(qep,"Breakdown Quadratic Arnoldi method (it=%D norm=%G)\n",qep->its,beta);
246: qep->reason = QEP_DIVERGED_BREAKDOWN;
247: } else {
248: /* Prepare the Rayleigh quotient for restart */
249: DSTruncate(qep->ds,k+l);
250: DSGetDimensions(qep->ds,&newn,PETSC_NULL,PETSC_NULL,PETSC_NULL);
251: l = newn-k;
252: }
253: }
254: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
255: DSGetArray(qep->ds,DS_MAT_Q,&Q);
256: SlepcUpdateVectors(nv,qep->V,qep->nconv,k+l,Q,ld,PETSC_FALSE);
257: DSRestoreArray(qep->ds,DS_MAT_Q,&Q);
259: qep->nconv = k;
260: QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,nv);
261: }
263: for (j=0;j<qep->nconv;j++) {
264: qep->eigr[j] *= qep->sfactor;
265: qep->eigi[j] *= qep->sfactor;
266: }
268: /* truncate Schur decomposition and change the state to raw so that
269: DSVectors() computes eigenvectors from scratch */
270: DSSetDimensions(qep->ds,qep->nconv,PETSC_IGNORE,0,0);
271: DSSetState(qep->ds,DS_STATE_RAW);
273: /* Compute eigenvectors */
274: if (qep->nconv > 0) {
275: QEPComputeVectors_Schur(qep);
276: }
277: PetscFree(work);
278: return(0);
279: }
283: PetscErrorCode QEPSetFromOptions_QArnoldi(QEP qep)
284: {
286: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
287:
289: KSPSetFromOptions(ctx->ksp);
290: return(0);
291: }
295: PetscErrorCode QEPView_QArnoldi(QEP qep,PetscViewer viewer)
296: {
298: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
301: PetscViewerASCIIPushTab(viewer);
302: KSPView(ctx->ksp,viewer);
303: PetscViewerASCIIPopTab(viewer);
304: return(0);
305: }
309: PetscErrorCode QEPReset_QArnoldi(QEP qep)
310: {
312: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
315: KSPReset(ctx->ksp);
316: QEPDefaultFreeWork(qep);
317: QEPFreeSolution(qep);
318: return(0);
319: }
323: PetscErrorCode QEPDestroy_QArnoldi(QEP qep)
324: {
326: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
329: KSPDestroy(&ctx->ksp);
330: PetscFree(qep->data);
331: return(0);
332: }
334: EXTERN_C_BEGIN
337: PetscErrorCode QEPCreate_QArnoldi(QEP qep)
338: {
340: QEP_QARNOLDI *ctx;
343: PetscNewLog(qep,QEP_QARNOLDI,&ctx);
344: qep->data = ctx;
345: qep->ops->solve = QEPSolve_QArnoldi;
346: qep->ops->setup = QEPSetUp_QArnoldi;
347: qep->ops->setfromoptions = QEPSetFromOptions_QArnoldi;
348: qep->ops->destroy = QEPDestroy_QArnoldi;
349: qep->ops->reset = QEPReset_QArnoldi;
350: qep->ops->view = QEPView_QArnoldi;
351: KSPCreate(((PetscObject)qep)->comm,&ctx->ksp);
352: KSPSetOptionsPrefix(ctx->ksp,((PetscObject)qep)->prefix);
353: KSPAppendOptionsPrefix(ctx->ksp,"qep_");
354: PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)qep,1);
355: PetscLogObjectParent(qep,ctx->ksp);
356: return(0);
357: }
358: EXTERN_C_END