Actual source code: arnoldi.c
1: /*
3: SLEPc eigensolver: "arnoldi"
5: Method: Explicitly Restarted Arnoldi
7: Description:
9: This solver implements the Arnoldi method with explicit restart
10: and deflation.
12: Algorithm:
14: The implemented algorithm builds an Arnoldi factorization of order
15: ncv. Converged eigenpairs are locked and the iteration is restarted
16: with the rest of the columns being the active columns for the next
17: Arnoldi factorization. Currently, no filtering is applied to the
18: vector used for restarting.
20: References:
22: [1] Z. Bai et al. (eds.), "Templates for the Solution of Algebraic
23: Eigenvalue Problems: A Practical Guide", SIAM (2000), pp 161-165.
25: [2] Y. Saad, "Numerical Methods for Large Eigenvalue Problems",
26: John Wiley (1992), pp 172-183.
28: Last update: June 2004
30: */
31: #include src/eps/epsimpl.h
32: #include slepcblaslapack.h
36: PetscErrorCode EPSSetUp_ARNOLDI(EPS eps)
37: {
38: PetscErrorCode ierr;
39: int N;
42: VecGetSize(eps->vec_initial,&N);
43: if (eps->ncv) {
44: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
45: }
46: else eps->ncv = PetscMax(2*eps->nev,eps->nev+15);
47: if (!eps->max_it) eps->max_it = PetscMax(100,N);
48: if (!eps->tol) eps->tol = 1.e-7;
49: if (eps->which!=EPS_LARGEST_MAGNITUDE)
50: SETERRQ(1,"Wrong value of eps->which");
51: EPSAllocateSolution(eps);
52: if (eps->T) { PetscFree(eps->T); }
53: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
54: EPSDefaultGetWork(eps,eps->ncv+1);
55: return(0);
56: }
60: /*
61: EPSBasicArnoldi - Computes an m-step Arnoldi factorization. The first k
62: columns are assumed to be locked and therefore they are not modified. On
63: exit, the following relation is satisfied:
65: OP * V - V * H = f * e_m^T
67: where the columns of V are the Arnoldi vectors (which are B-orthonormal),
68: H is an upper Hessenberg matrix, f is the residual vector and e_m is
69: the m-th vector of the canonical basis. The vector f is B-orthogonal to
70: the columns of V. On exit, beta contains the B-norm of f and the next
71: Arnoldi vector can be computed as v_{m+1} = f / beta.
72: */
73: static PetscErrorCode EPSBasicArnoldi(EPS eps,PetscScalar *H,Vec *V,int k,int m,Vec f,PetscReal *beta)
74: {
75: PetscErrorCode ierr;
76: int j;
77: PetscReal norm;
78: PetscScalar t;
79: PetscTruth breakdown;
82: for (j=k;j<m-1;j++) {
83: STApply(eps->OP,V[j],V[j+1]);
84: EPSOrthogonalize(eps,eps->nds,eps->DS,V[j+1],PETSC_NULL,PETSC_NULL,PETSC_NULL);
85: EPSOrthogonalize(eps,j+1,V,V[j+1],H+m*j,&norm,&breakdown);
86: H[(m+1)*j+1] = norm;
87: if (breakdown) {
88: PetscLogInfo(eps,"Breakdown in Arnoldi method (norm=%g)\n",norm);
89: SlepcVecSetRandom(V[j+1]);
90: STNorm(eps->OP,V[j+1],&norm);
91: }
92: t = 1 / norm;
93: VecScale(&t,V[j+1]);
94: }
95: STApply(eps->OP,V[m-1],f);
96: EPSOrthogonalize(eps,m,V,f,H+m*(m-1),beta,PETSC_NULL);
97: return(0);
98: }
100: #define SWAP(a,b,t) {t=a;a=b;b=t;}
104: PetscErrorCode EPSSolve_ARNOLDI(EPS eps)
105: {
106: PetscErrorCode ierr;
107: int i,k,mout,info,ifst,ilst,ncv=eps->ncv;
108: Vec f=eps->work[ncv];
109: PetscScalar *H=eps->T,*U,*Y,*work,ts;
110: PetscReal norm,beta,tr;
111: #if defined(PETSC_USE_COMPLEX)
112: PetscReal *rwork;
113: #endif
116: #if defined(PETSC_BLASLAPACK_ESSL_ONLY)
117: SETERRQ(PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
118: #endif
119: PetscMemzero(H,ncv*ncv*sizeof(PetscScalar));
120: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&U);
121: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Y);
122: PetscMalloc(3*ncv*sizeof(PetscScalar),&work);
123: #if defined(PETSC_USE_COMPLEX)
124: PetscMalloc(ncv*sizeof(PetscReal),&rwork);
125: #endif
127: /* The first Arnoldi vector is the normalized initial vector */
128: VecCopy(eps->vec_initial,eps->V[0]);
129: STNorm(eps->OP,eps->V[0],&norm);
130: ts = 1 / norm;
131: VecScale(&ts,eps->V[0]);
132:
133: eps->nconv = 0;
134: eps->its = 0;
135: /* Restart loop */
136: while (eps->its<eps->max_it) {
137: /* Compute an ncv-step Arnoldi factorization */
138: EPSBasicArnoldi(eps,H,eps->V,eps->nconv,ncv,f,&beta);
139: eps->its = eps->its + ncv - eps->nconv;
141: /* At this point, H has the following structure
143: | * * | * * * * |
144: | * | * * * * |
145: | ------|-------------- |
146: H = | | * * * * |
147: | | * * * * |
148: | | * * * |
149: | | * * |
151: that is, an upper Hessenberg matrix of order ncv whose principal
152: submatrix of order nconv is (quasi-)triangular. */
154: /* Reduce H to (quasi-)triangular form, H <- U H U' */
155: PetscMemzero(U,ncv*ncv*sizeof(PetscScalar));
156: for (i=0;i<ncv;i++) { U[i*(ncv+1)] = 1.0; }
157: EPSDenseSchur(ncv,eps->nconv,H,U,eps->eigr,eps->eigi);
159: /* Sort the remaining columns of the Schur form */
160: EPSSortDenseSchur(ncv,eps->nconv,H,U,eps->eigr,eps->eigi);
162: /* Compute eigenvectors Y of H */
163: PetscMemcpy(Y,U,ncv*ncv*sizeof(PetscScalar));
164: #if !defined(PETSC_USE_COMPLEX)
165: LAtrevc_("R","B",PETSC_NULL,&ncv,H,&ncv,PETSC_NULL,&ncv,Y,&ncv,&ncv,&mout,work,&info,1,1);
166: #else
167: LAtrevc_("R","B",PETSC_NULL,&ncv,H,&ncv,PETSC_NULL,&ncv,Y,&ncv,&ncv,&mout,work,rwork,&info,1,1);
168: #endif
169: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
171: /* Compute residual norm estimates as beta*abs(Y(m,:)) */
172: for (i=eps->nconv;i<ncv;i++) {
173: eps->errest[i] = beta*PetscAbsScalar(Y[i*ncv+ncv-1]);
174: }
176: /* Look for converged eigenpairs. If necessary, reorder the Arnoldi
177: factorization so that all converged eigenvalues are first */
178: k = eps->nconv;
179: while (k<ncv && eps->errest[k]<eps->tol) { k++; }
180: for (i=k;i<ncv;i++) {
181: if (eps->errest[i]<eps->tol) {
182: ifst = i + 1;
183: ilst = k + 1;
184: #if !defined(PETSC_USE_COMPLEX)
185: LAtrexc_("V",&ncv,H,&ncv,U,&ncv,&ifst,&ilst,work,&info,1);
186: SWAP(eps->eigr[k],eps->eigr[i],ts);
187: SWAP(eps->eigi[k],eps->eigi[i],ts);
188: SWAP(eps->errest[k],eps->errest[i],tr);
189: if (eps->eigi[i] != 0) {
190: SWAP(eps->eigr[k+1],eps->eigr[i+1],ts);
191: SWAP(eps->eigi[k+1],eps->eigi[i+1],ts);
192: SWAP(eps->errest[k+1],eps->errest[i+1],tr);
193: k++;
194: }
195: #else
196: LAtrexc_("V",&ncv,H,&ncv,U,&ncv,&ifst,&ilst,&info,1);
197: SWAP(eps->eigr[k],eps->eigr[i],ts);
198: SWAP(eps->errest[k],eps->errest[i],tr);
199: #endif
200: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
201: k++;
202: }
203: #if !defined(PETSC_USE_COMPLEX)
204: if (eps->eigi[i] != 0) i++;
205: #endif
206: }
208: /* Update V(:,idx) = V*U(:,idx) */
209: EPSReverseProjection(eps,eps->V,U,eps->nconv,ncv,eps->work);
210: eps->nconv = k;
211: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,ncv);
212: if (eps->nconv >= eps->nev) break;
213: }
214:
215: if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
216: else eps->reason = EPS_DIVERGED_ITS;
217: #if defined(PETSC_USE_COMPLEX)
218: for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;
219: #endif
221: PetscFree(U);
222: PetscFree(Y);
223: PetscFree(work);
224: #if defined(PETSC_USE_COMPLEX)
225: PetscFree(rwork);
226: #endif
227: return(0);
228: }
230: EXTERN_C_BEGIN
233: PetscErrorCode EPSCreate_ARNOLDI(EPS eps)
234: {
236: eps->data = (void *) 0;
237: eps->ops->solve = EPSSolve_ARNOLDI;
238: eps->ops->setup = EPSSetUp_ARNOLDI;
239: eps->ops->destroy = EPSDestroy_Default;
240: eps->ops->backtransform = EPSBackTransform_Default;
241: eps->ops->computevectors = EPSComputeVectors_Schur;
242: return(0);
243: }
244: EXTERN_C_END