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;
54: if (qep->problem_type != QEP_GENERAL) SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of qep->problem_type");
56: QEPAllocateSolution(qep);
57: QEPDefaultGetWork(qep,4);
59: DSSetType(qep->ds,DSNHEP);
60: DSSetExtraRow(qep->ds,PETSC_TRUE);
61: DSAllocate(qep->ds,qep->ncv+1);
63: KSPSetOperators(ctx->ksp,qep->M,qep->M,DIFFERENT_NONZERO_PATTERN);
64: KSPSetUp(ctx->ksp);
65: return(0);
66: }
70: /*
71: Compute a step of Classical Gram-Schmidt orthogonalization
72: */
73: 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)
74: {
76: PetscBLASInt ione = 1,j_1 = j+1;
77: PetscReal x,y;
78: PetscScalar dot,one = 1.0,zero = 0.0;
81: /* compute norm of v and w */
82: if (onorm) {
83: VecNorm(v,NORM_2,&x);
84: VecNorm(w,NORM_2,&y);
85: *onorm = PetscSqrtReal(x*x+y*y);
86: }
88: /* orthogonalize: compute h */
89: VecMDot(v,j_1,V,h);
90: VecMDot(w,j_1,V,work);
91: if (j>0)
92: BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione);
93: VecDot(w,t,&dot);
94: h[j] += dot;
96: /* orthogonalize: update v and w */
97: SlepcVecMAXPBY(v,1.0,-1.0,j_1,h,V);
98: if (j>0) {
99: BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione);
100: SlepcVecMAXPBY(w,1.0,-1.0,j_1,work,V);
101: }
102: VecAXPY(w,-h[j],t);
103:
104: /* compute norm of v and w */
105: if (norm) {
106: VecNorm(v,NORM_2,&x);
107: VecNorm(w,NORM_2,&y);
108: *norm = PetscSqrtReal(x*x+y*y);
109: }
110: return(0);
111: }
115: /*
116: Compute a run of Q-Arnoldi iterations
117: */
118: PetscErrorCode QEPQArnoldi(QEP qep,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
119: {
120: PetscErrorCode ierr;
121: PetscInt i,j,l,m = *M;
122: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
123: Vec t = qep->work[2],u = qep->work[3];
124: IPOrthogRefineType refinement;
125: PetscReal norm,onorm,eta;
126: PetscScalar *c = work + m;
129: IPGetOrthogonalization(qep->ip,PETSC_NULL,&refinement,&eta);
130: VecCopy(v,qep->V[k]);
131:
132: for (j=k;j<m;j++) {
133: /* apply operator */
134: VecCopy(w,t);
135: MatMult(qep->K,v,u);
136: MatMult(qep->C,t,w);
137: VecAXPY(u,qep->sfactor,w);
138: KSPSolve(ctx->ksp,u,w);
139: VecScale(w,-1.0/(qep->sfactor*qep->sfactor));
140: VecCopy(t,v);
142: /* orthogonalize */
143: switch (refinement) {
144: case IP_ORTHOG_REFINE_NEVER:
145: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,PETSC_NULL,&norm,work);
146: *breakdown = PETSC_FALSE;
147: break;
148: case IP_ORTHOG_REFINE_ALWAYS:
149: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,PETSC_NULL,PETSC_NULL,work);
150: QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,&onorm,&norm,work);
151: for (i=0;i<j;i++) H[ldh*j+i] += c[i];
152: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
153: else *breakdown = PETSC_FALSE;
154: break;
155: case IP_ORTHOG_REFINE_IFNEEDED:
156: QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,&onorm,&norm,work);
157: /* ||q|| < eta ||h|| */
158: l = 1;
159: while (l<3 && norm < eta * onorm) {
160: l++;
161: onorm = norm;
162: QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,PETSC_NULL,&norm,work);
163: for (i=0;i<j;i++) H[ldh*j+i] += c[i];
164: }
165: if (norm < eta * onorm) *breakdown = PETSC_TRUE;
166: else *breakdown = PETSC_FALSE;
167: break;
168: default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of ip->orth_ref");
169: }
170: VecScale(v,1.0/norm);
171: VecScale(w,1.0/norm);
172:
173: H[j+1+ldh*j] = norm;
174: if (j<m-1) {
175: VecCopy(v,V[j+1]);
176: }
177: }
178: *beta = norm;
179: return(0);
180: }
184: PetscErrorCode QEPSolve_QArnoldi(QEP qep)
185: {
187: PetscInt j,k,l,lwork,nv,ld,newn;
188: Vec v=qep->work[0],w=qep->work[1];
189: PetscScalar *S,*Q,*work;
190: PetscReal beta,norm,x,y;
191: PetscBool breakdown;
194: DSGetLeadingDimension(qep->ds,&ld);
195: lwork = 7*qep->ncv;
196: PetscMalloc(lwork*sizeof(PetscScalar),&work);
198: /* Get the starting Arnoldi vector */
199: if (qep->nini>0) {
200: VecCopy(qep->V[0],v);
201: } else {
202: SlepcVecSetRandom(v,qep->rand);
203: }
204: /* w is always a random vector */
205: SlepcVecSetRandom(w,qep->rand);
206: VecNorm(v,NORM_2,&x);
207: VecNorm(w,NORM_2,&y);
208: norm = PetscSqrtReal(x*x+y*y);
209: VecScale(v,1.0/norm);
210: VecScale(w,1.0/norm);
211:
212: /* Restart loop */
213: l = 0;
214: while (qep->reason == QEP_CONVERGED_ITERATING) {
215: qep->its++;
217: /* Compute an nv-step Arnoldi factorization */
218: nv = PetscMin(qep->nconv+qep->mpd,qep->ncv);
219: DSGetArray(qep->ds,DS_MAT_A,&S);
220: QEPQArnoldi(qep,S,ld,qep->V,qep->nconv+l,&nv,v,w,&beta,&breakdown,work);
221: DSRestoreArray(qep->ds,DS_MAT_A,&S);
222: DSSetDimensions(qep->ds,nv,PETSC_IGNORE,qep->nconv,qep->nconv+l);
223: if (l==0) {
224: DSSetState(qep->ds,DS_STATE_INTERMEDIATE);
225: } else {
226: DSSetState(qep->ds,DS_STATE_RAW);
227: }
229: /* Solve projected problem */
230: DSSolve(qep->ds,qep->eigr,qep->eigi);
231: DSSort(qep->ds,qep->eigr,qep->eigi,PETSC_NULL,PETSC_NULL,PETSC_NULL);
232: DSUpdateExtraRow(qep->ds);
234: /* Check convergence */
235: QEPKrylovConvergence(qep,PETSC_FALSE,qep->nconv,nv-qep->nconv,nv,beta,&k);
236: if (qep->its >= qep->max_it) qep->reason = QEP_DIVERGED_ITS;
237: if (k >= qep->nev) qep->reason = QEP_CONVERGED_TOL;
238:
239: /* Update l */
240: if (qep->reason != QEP_CONVERGED_ITERATING || breakdown) l = 0;
241: else l = (nv-k)/2;
243: if (qep->reason == QEP_CONVERGED_ITERATING) {
244: if (breakdown) {
245: /* Stop if breakdown */
246: PetscInfo2(qep,"Breakdown Quadratic Arnoldi method (it=%D norm=%G)\n",qep->its,beta);
247: qep->reason = QEP_DIVERGED_BREAKDOWN;
248: } else {
249: /* Prepare the Rayleigh quotient for restart */
250: DSTruncate(qep->ds,k+l);
251: DSGetDimensions(qep->ds,&newn,PETSC_NULL,PETSC_NULL,PETSC_NULL);
252: l = newn-k;
253: }
254: }
255: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
256: DSGetArray(qep->ds,DS_MAT_Q,&Q);
257: SlepcUpdateVectors(nv,qep->V,qep->nconv,k+l,Q,ld,PETSC_FALSE);
258: DSRestoreArray(qep->ds,DS_MAT_Q,&Q);
260: qep->nconv = k;
261: QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,nv);
262: }
264: for (j=0;j<qep->nconv;j++) {
265: qep->eigr[j] *= qep->sfactor;
266: qep->eigi[j] *= qep->sfactor;
267: }
269: /* truncate Schur decomposition and change the state to raw so that
270: PSVectors() computes eigenvectors from scratch */
271: DSSetDimensions(qep->ds,qep->nconv,PETSC_IGNORE,0,0);
272: DSSetState(qep->ds,DS_STATE_RAW);
274: /* Compute eigenvectors */
275: if (qep->nconv > 0) {
276: QEPComputeVectors_Schur(qep);
277: }
278: PetscFree(work);
279: return(0);
280: }
284: PetscErrorCode QEPSetFromOptions_QArnoldi(QEP qep)
285: {
287: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
288:
290: KSPSetFromOptions(ctx->ksp);
291: return(0);
292: }
296: PetscErrorCode QEPView_QArnoldi(QEP qep,PetscViewer viewer)
297: {
299: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
302: PetscViewerASCIIPushTab(viewer);
303: KSPView(ctx->ksp,viewer);
304: PetscViewerASCIIPopTab(viewer);
305: return(0);
306: }
310: PetscErrorCode QEPReset_QArnoldi(QEP qep)
311: {
313: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
316: KSPReset(ctx->ksp);
317: QEPDefaultFreeWork(qep);
318: QEPFreeSolution(qep);
319: return(0);
320: }
324: PetscErrorCode QEPDestroy_QArnoldi(QEP qep)
325: {
327: QEP_QARNOLDI *ctx = (QEP_QARNOLDI*)qep->data;
330: KSPDestroy(&ctx->ksp);
331: PetscFree(qep->data);
332: return(0);
333: }
335: EXTERN_C_BEGIN
338: PetscErrorCode QEPCreate_QArnoldi(QEP qep)
339: {
341: QEP_QARNOLDI *ctx;
344: PetscNewLog(qep,QEP_QARNOLDI,&ctx);
345: qep->data = ctx;
346: qep->ops->solve = QEPSolve_QArnoldi;
347: qep->ops->setup = QEPSetUp_QArnoldi;
348: qep->ops->setfromoptions = QEPSetFromOptions_QArnoldi;
349: qep->ops->destroy = QEPDestroy_QArnoldi;
350: qep->ops->reset = QEPReset_QArnoldi;
351: qep->ops->view = QEPView_QArnoldi;
352: KSPCreate(((PetscObject)qep)->comm,&ctx->ksp);
353: KSPSetOptionsPrefix(ctx->ksp,((PetscObject)qep)->prefix);
354: KSPAppendOptionsPrefix(ctx->ksp,"qep_");
355: PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)qep,1);
356: PetscLogObjectParent(qep,ctx->ksp);
357: return(0);
358: }
359: EXTERN_C_END