Actual source code: qlanczos.c
1: /*
3: SLEPc quadratic eigensolver: "qlanczos"
5: Method: Q-Lanczos
7: Algorithm:
9: Quadratic Lanczos with indefinite B-orthogonalization and
10: thick restart.
12: References:
14: [1] C. Campos and J.E. Roman, "A thick-restart Q-Lanczos method
15: for quadratic eigenvalue problems", in preparation, 2013.
17: Last update: Mar 2013
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_QLanczos(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: }
69: if (qep->problem_type!=QEP_HERMITIAN) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
71: QEPAllocateSolution(qep);
72: QEPSetWorkVecs(qep,4);
74: DSSetType(qep->ds,DSGHIEP);
75: DSSetCompact(qep->ds,PETSC_TRUE);
76: DSAllocate(qep->ds,qep->ncv+1);
77: return(0);
78: }
83: /*
84: Compute B-norm of v=[v1;v2] whith B=diag(-qep->T[0],qep->T[2])
85: */
86: PetscErrorCode QEPQLanczosNorm_private(QEP qep,Vec v1,Vec v2,PetscReal *norm,Vec vw)
87: {
89: PetscScalar p1,p2;
90:
92: STMatMult(qep->st,0,v1,vw);
93: VecDot(vw,v1,&p1);
94: STMatMult(qep->st,2,v2,vw);
95: VecDot(vw,v2,&p2);
96: *norm = PetscRealPart(-p1+qep->sfactor*qep->sfactor*p2);
97: *norm = (*norm>0.0)?PetscSqrtReal(*norm):-PetscSqrtReal(-*norm);
98: return(0);
99: }
104: /*
105: Compute a step of Classical Gram-Schmidt orthogonalization
106: */
107: static PetscErrorCode QEPQLanczosCGS(QEP qep,PetscScalar *H,PetscBLASInt ldh,PetscReal *omega,PetscScalar *h,PetscBLASInt j,Vec *V,Vec t,Vec v,Vec w,PetscReal *onorm,PetscReal *norm,PetscScalar *work,Vec vw)
108: {
110: PetscBLASInt ione = 1,j_1 = j+1,i;
111: PetscScalar dot,one = 1.0,zero = 0.0;
114: /* compute norm of [v;w] */
115: if (onorm) {
116: QEPQLanczosNorm_private(qep,v,w,onorm,vw);
117: }
119: /* orthogonalize: compute h */
120: STMatMult(qep->st,0,v,vw);
121: VecMDot(vw,j_1,V,h);
122: for (i=0;i<=j;i++) h[i] = -h[i];
123: STMatMult(qep->st,2,w,vw);
124: if (j>0) {
125: VecMDot(vw,j_1,V,work);
126: for (i=0;i<j_1;i++) work[i] *= qep->sfactor*qep->sfactor;
127: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione));
128: }
129: VecDot(vw,t,&dot);
130: h[j] += dot*qep->sfactor*qep->sfactor;
131: /* apply inverse of signature */
132: for (i=0;i<=j;i++) h[i] /= omega[i];
133: /* orthogonalize: update v and w */
134: SlepcVecMAXPBY(v,1.0,-1.0,j_1,h,V);
135: if (j>0) {
136: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione));
137: SlepcVecMAXPBY(w,1.0,-1.0,j_1,work,V);
138: }
139: VecAXPY(w,-h[j],t);
140: /* H pseudosymmetric, signature is not revert*/
142: /* compute norm of v and w */
143: if (norm) {
144: QEPQLanczosNorm_private(qep,v,w,norm,vw);
145: }
146: return(0);
147: }
151: /*
152: Compute a run of Q-Lanczos iterations
153: */
154: static PetscErrorCode QEPQLanczos(QEP qep,PetscScalar *H,PetscInt ldh,PetscReal *omega,Vec *V,PetscInt k,PetscInt *M,Vec v,Vec w,PetscBool *breakdown,PetscScalar *work,Vec *t_)
155: {
156: PetscErrorCode ierr;
157: PetscInt i,j,m = *M,nwu = 0,l;
158: Vec t = t_[0],u = t_[1];
159: IPOrthogRefineType refinement;
160: PetscReal onorm,norm,eta;
161: PetscScalar *c;
164: IPGetOrthogonalization(qep->ip,NULL,&refinement,&eta);
165: nwu += m+1;
166: c = work+nwu;
168: for (j=k;j<m;j++) {
169: /* apply operator */
170: VecCopy(w,t);
171: STMatMult(qep->st,0,v,u);
172: STMatMult(qep->st,1,t,w);
173: VecAXPY(u,qep->sfactor,w);
174: STMatSolve(qep->st,2,u,w);
175: VecScale(w,-1.0/(qep->sfactor*qep->sfactor));
176: VecCopy(t,v);
178: /* orthogonalize */
179: switch (refinement) {
180: case IP_ORTHOG_REFINE_NEVER:
181: QEPQLanczosCGS(qep,H,ldh,omega,H+ldh*j,j,V,t,v,w,NULL,&norm,work,u);
182: break;
183: case IP_ORTHOG_REFINE_ALWAYS:
184: QEPQLanczosCGS(qep,H,ldh,omega,H+ldh*j,j,V,t,v,w,NULL,NULL,work,u);
185: QEPQLanczosCGS(qep,H,ldh,omega,c,j,V,t,v,w,&onorm,&norm,work,u);
186: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
187: break;
188: case IP_ORTHOG_REFINE_IFNEEDED: /* -- */
189: QEPQLanczosCGS(qep,H,ldh,omega,H+ldh*j,j,V,t,v,w,&onorm,&norm,work,u);
190: /* ||q|| < eta ||h|| */
191: l = 1;
192: while (l<3 && PetscAbsReal(norm) < eta * PetscAbsReal(onorm)) {
193: l++;
194: onorm = norm;
195: QEPQLanczosCGS(qep,H,ldh,omega,c,j,V,t,v,w,NULL,&norm,work,u);
196: for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
197: }
198: break;
199: default: SETERRQ(PetscObjectComm((PetscObject)qep),1,"Wrong value of ip->orth_ref");
200: }
201: if (breakdown) *breakdown = PETSC_FALSE; /* -- */
202: VecScale(v,1.0/norm);
203: VecScale(w,1.0/norm);
205: H[j+1+ldh*j] = norm;
206: omega[j+1] = (norm<0.0)?-1.0:1.0;
207: if (j<m-1) {
208: VecCopy(v,V[j+1]);
209: }
210: }
211: return(0);
212: }
216: PetscErrorCode QEPSolve_QLanczos(QEP qep)
217: {
219: PetscInt j,k,l,lwork,nv,ld,ti,t;
220: Vec v=qep->work[0],w=qep->work[1],w_=qep->work[2]; /* v_=qep->work[3]; */
221: PetscScalar *S,*Q,*work;
222: PetscReal beta,norm,*omega,*a,*b,*r;
223: PetscBool breakdown;
226: DSGetLeadingDimension(qep->ds,&ld);
227: lwork = 7*qep->ncv; /* -- */
228: PetscMalloc(lwork*sizeof(PetscScalar),&work);
230: /* Get the starting Lanczos vector */
231: if (qep->nini==0) {
232: SlepcVecSetRandom(qep->V[0],qep->rand);
233: }
234: /* w is always a random vector */
235: SlepcVecSetRandom(w,qep->rand);
236: QEPQLanczosNorm_private(qep,qep->V[0],w,&norm,w_);
237: if (norm==0.0) SETERRQ(((PetscObject)qep)->comm,1,"Initial vector is zero or belongs to the deflation space");
238: VecScale(qep->V[0],1.0/norm);
239: VecScale(w,1.0/norm);
240: VecCopy(qep->V[0],v);
241: DSGetArrayReal(qep->ds,DS_MAT_D,&omega);
242: omega[0] = (norm > 0)?1.0:-1.0;
243: DSRestoreArrayReal(qep->ds,DS_MAT_D,&omega);
245: /* Restart loop */
246: l = 0;
247: while (qep->reason == QEP_CONVERGED_ITERATING) {
248: qep->its++;
250: /* Form projected problem matrix on S */
251: DSGetArray(qep->ds,DS_MAT_A,&S);
252: DSGetArrayReal(qep->ds,DS_MAT_T,&a);
253: b = a+ld;
254: r = b+ld;
255: DSGetArrayReal(qep->ds,DS_MAT_D,&omega);
256: PetscMemzero(S,ld*ld*sizeof(PetscScalar));
257: ti = qep->nconv+l;
258: for (j=0;j<ti;j++) {
259: S[j+j*ld] = a[j]*omega[j];
260: S[j+1+j*ld] = b[j]*omega[j+1];
261: S[j+(j+1)*ld] = b[j]*omega[j];
262: }
263: for (j=qep->nconv;j<ti;j++) {
264: S[ti+j*ld] = r[j]*omega[ti];
265: }
267: /* Compute an nv-step Lanczos factorization */
268: nv = PetscMin(qep->nconv+qep->mpd,qep->ncv);
269: QEPQLanczos(qep,S,ld,omega,qep->V,qep->nconv+l,&nv,v,w,&breakdown,work,&qep->work[2]);
270: for (j=ti;j<nv-1;j++) {
271: a[j] = PetscRealPart(S[j+j*ld]*omega[j]);
272: b[j] = PetscAbsReal(PetscRealPart(S[j+1+j*ld]));
273: }
274: a[nv-1] = PetscRealPart(S[nv-1+(nv-1)*ld]*omega[nv-1]);
275: beta = PetscAbsReal(PetscRealPart(S[nv+(nv-1)*ld]));
276: DSRestoreArray(qep->ds,DS_MAT_A,&S);
277: DSRestoreArrayReal(qep->ds,DS_MAT_T,&a);
278: DSRestoreArrayReal(qep->ds,DS_MAT_D,&omega);
279: DSSetDimensions(qep->ds,nv,0,qep->nconv,qep->nconv+l);
280: if (l==0) {
281: DSSetState(qep->ds,DS_STATE_INTERMEDIATE);
282: } else {
283: DSSetState(qep->ds,DS_STATE_RAW);
284: }
286: /* Solve projected problem */
287: DSSolve(qep->ds,qep->eigr,qep->eigi);
288: DSSort(qep->ds,qep->eigr,qep->eigi,NULL,NULL,NULL);
290: /* Check convergence */
291: DSGetDimensions(qep->ds,NULL,NULL,NULL,NULL,&t);
292: QEPKrylovConvergence(qep,PETSC_FALSE,qep->nconv,t-qep->nconv,nv,beta,&k);
293: if (qep->its >= qep->max_it) qep->reason = QEP_DIVERGED_ITS;
294: if (k >= qep->nev) qep->reason = QEP_CONVERGED_TOL;
296: /* Update l */
297: if (qep->reason != QEP_CONVERGED_ITERATING || breakdown) l = 0;
298: else {
299: l = PetscMax(1,(PetscInt)((nv-k)/2));
300: l = PetscMin(l,t);
301: DSGetArrayReal(qep->ds,DS_MAT_T,&a);
302: if (*(a+ld+k+l-1)!=0) {
303: if (k+l<nv-1) l = l+1;
304: else l = l-1;
305: }
306: DSRestoreArrayReal(qep->ds,DS_MAT_T,&a);
307: }
309: if (qep->reason == QEP_CONVERGED_ITERATING) {
310: if (breakdown) {
311: /* Stop if breakdown */
312: PetscInfo2(qep,"Breakdown Quadratic Lanczos method (it=%D norm=%G)\n",qep->its,beta);
313: qep->reason = QEP_DIVERGED_BREAKDOWN;
314: } else {
315: /* Prepare the Rayleigh quotient for restart */
316: DSGetArray(qep->ds,DS_MAT_Q,&Q);
317: DSGetArrayReal(qep->ds,DS_MAT_T,&a);
318: DSGetArrayReal(qep->ds,DS_MAT_D,&omega);
319: r = a + 2*ld;
320: for (j=k;j<k+l;j++) {
321: r[j] = PetscRealPart(Q[nv-1+j*ld]*beta);
322: }
323: b = a+ld;
324: b[k+l-1] = r[k+l-1];
325: omega[k+l] = omega[nv];
326: DSRestoreArray(qep->ds,DS_MAT_Q,&Q);
327: DSRestoreArrayReal(qep->ds,DS_MAT_T,&a);
328: DSRestoreArrayReal(qep->ds,DS_MAT_D,&omega);
329: }
330: }
331: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
332: DSGetArray(qep->ds,DS_MAT_Q,&Q);
333: SlepcUpdateVectors(nv,qep->V,qep->nconv,k+l,Q,ld,PETSC_FALSE);
334: DSRestoreArray(qep->ds,DS_MAT_Q,&Q);
335:
336: /* Append v to V*/
337: if (qep->reason == QEP_CONVERGED_ITERATING) {
338: VecCopy(v,qep->V[k+l]);
339: }
340: qep->nconv = k;
341: QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,nv);
342: }
344: for (j=0;j<qep->nconv;j++) {
345: qep->eigr[j] *= qep->sfactor;
346: qep->eigi[j] *= qep->sfactor;
347: }
349: /* truncate Schur decomposition and change the state to raw so that
350: DSVectors() computes eigenvectors from scratch */
351: DSSetDimensions(qep->ds,qep->nconv,0,0,0);
352: DSSetState(qep->ds,DS_STATE_RAW);
354: /* Compute eigenvectors */
355: if (qep->nconv > 0) {
356: QEPComputeVectors_Schur(qep);
357: }
358: PetscFree(work);
359: return(0);
360: }
364: PETSC_EXTERN PetscErrorCode QEPCreate_QLanczos(QEP qep)
365: {
367: qep->ops->solve = QEPSolve_QLanczos;
368: qep->ops->setup = QEPSetUp_QLanczos;
369: qep->ops->reset = QEPReset_Default;
370: return(0);
371: }