Actual source code: krylovschur.c
1: /*
3: SLEPc eigensolver: "krylovschur"
5: Method: Krylov-Schur
7: Algorithm:
9: Single-vector Krylov-Schur method for both symmetric and non-symmetric
10: problems.
12: References:
14: [1] "Krylov-Schur Methods in SLEPc", SLEPc Technical Report STR-7,
15: available at http://www.grycap.upv.es/slepc.
17: [2] G.W. Stewart, "A Krylov-Schur Algorithm for Large Eigenproblems",
18: SIAM J. Matrix Analysis and App., 23(3), pp. 601-614, 2001.
20: Last update: Oct 2006
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: SLEPc - Scalable Library for Eigenvalue Problem Computations
24: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
26: This file is part of SLEPc. See the README file for conditions of use
27: and additional information.
28: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
29: */
31: #include src/eps/epsimpl.h
32: #include slepcblaslapack.h
36: PetscErrorCode EPSSetUp_KRYLOVSCHUR(EPS eps)
37: {
39: PetscInt N;
42: VecGetSize(eps->vec_initial,&N);
43: if (eps->ncv) {
44: if (eps->ncv<eps->nev+1) SETERRQ(1,"The value of ncv must be at least nev+1");
45: }
46: else eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
47: if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
48: if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
49: SETERRQ(1,"Wrong value of eps->which");
51: EPSAllocateSolution(eps);
52: PetscFree(eps->T);
53: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
54: EPSDefaultGetWork(eps,2);
55: return(0);
56: }
60: PetscErrorCode EPSSolve_KRYLOVSCHUR(EPS eps)
61: {
63: int i,j,k,l,n,lwork,*perm;
64: Vec u=eps->work[1];
65: PetscScalar *S=eps->T,*Q,*work,*b;
66: PetscReal beta,*ritz;
67: PetscTruth breakdown;
70: PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));
71: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);
72: PetscMalloc(eps->ncv*sizeof(PetscScalar),&b);
73: lwork = (eps->ncv+4)*eps->ncv;
74: if (!eps->ishermitian) {
75: PetscMalloc(lwork*sizeof(PetscScalar),&work);
76: } else {
77: PetscMalloc(eps->ncv*sizeof(PetscReal),&ritz);
78: PetscMalloc(eps->ncv*sizeof(int),&perm);
79: }
81: /* Get the starting Arnoldi vector */
82: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
83: l = 0;
84:
85: /* Restart loop */
86: while (eps->reason == EPS_CONVERGED_ITERATING) {
87: eps->its++;
89: /* Compute an nv-step Arnoldi factorization */
90: eps->nv = eps->ncv;
91: EPSBasicArnoldi(eps,PETSC_FALSE,S,eps->V,eps->nconv+l,&eps->nv,u,&beta,&breakdown);
92: VecScale(u,1.0/beta);
94: if (!eps->ishermitian) {
95: n = eps->nv; /* size of Q */
96: if (l==0) {
97: PetscMemzero(Q,n*n*sizeof(PetscScalar));
98: for (i=0;i<n;i++)
99: Q[i*(n+1)] = 1.0;
100: } else {
101: /* Reduce S to Hessenberg form, S <- Q S Q' */
102: EPSDenseHessenberg(n,eps->nconv,S,eps->ncv,Q);
103: }
104: /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
105: EPSDenseSchur(n,eps->nconv,S,eps->ncv,Q,eps->eigr,eps->eigi);
106: /* Sort the remaining columns of the Schur form */
107: EPSSortDenseSchur(n,eps->nconv,S,eps->ncv,Q,eps->eigr,eps->eigi,eps->which);
108: /* Compute residual norm estimates */
109: ArnoldiResiduals(S,eps->ncv,Q,beta,eps->nconv,n,eps->eigr,eps->eigi,eps->errest,work);
110: } else {
111: n = eps->nv-eps->nconv; /* size of Q */
112: /* Reduce S to diagonal form, S <- Q S Q' */
113: if (l==0) {
114: EPSDenseTridiagonal(n,S+eps->nconv*(eps->ncv+1),eps->ncv,ritz,Q+eps->nconv*n);
115: } else {
116: EPSDenseHEP(n,S+eps->nconv*(eps->ncv+1),eps->ncv,ritz,Q+eps->nconv*n);
117: }
118: /* Sort the remaining columns of the Schur form */
119: if (eps->which == EPS_SMALLEST_REAL) {
120: for (i=0;i<n;i++)
121: eps->eigr[i+eps->nconv] = ritz[i];
122: } else {
123: #ifdef PETSC_USE_COMPLEX
124: for (i=0;i<n;i++)
125: eps->eigr[i+eps->nconv] = ritz[i];
126: EPSSortEigenvalues(n,eps->eigr+eps->nconv,eps->eigi,eps->which,n,perm);
127: #else
128: EPSSortEigenvalues(n,ritz,eps->eigi+eps->nconv,eps->which,n,perm);
129: #endif
130: for (i=0;i<n;i++)
131: eps->eigr[i+eps->nconv] = ritz[perm[i]];
132: PetscMemcpy(S,Q+eps->nconv*n,n*n*sizeof(PetscScalar));
133: for (j=0;j<n;j++)
134: for (i=0;i<n;i++)
135: Q[(j+eps->nconv)*n+i] = S[perm[j]*n+i];
136: }
137: /* rebuild S from eigr */
138: for (i=eps->nconv;i<eps->nv;i++) {
139: S[i*(eps->ncv+1)] = eps->eigr[i];
140: for (j=i+1;j<eps->ncv;j++)
141: S[i*eps->ncv+j] = 0.0;
142: }
143: /* Compute residual norm estimates */
144: for (i=eps->nconv;i<eps->nv;i++)
145: eps->errest[i] = beta*PetscAbsScalar(Q[(i+1)*n-1]) / PetscAbsScalar(eps->eigr[i]);
146: }
148: /* Check convergence */
149: k = eps->nconv;
150: while (k<eps->nv && eps->errest[k]<eps->tol) k++;
151: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
152: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
153:
154: /* Update l */
155: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
156: else {
157: l = (eps->nv-k)/2;
158: #if !defined(PETSC_USE_COMPLEX)
159: if (S[(k+l-1)*(eps->ncv+1)+1] != 0.0) {
160: if (k+l<eps->nv-1) l = l+1;
161: else l = l-1;
162: }
163: #endif
164: }
165:
166: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
167: for (i=eps->nconv;i<k+l;i++) {
168: VecSet(eps->AV[i],0.0);
169: if (!eps->ishermitian) {
170: VecMAXPY(eps->AV[i],n,Q+i*n,eps->V);
171: } else {
172: VecMAXPY(eps->AV[i],n,Q+i*n,eps->V+eps->nconv);
173: }
174: }
175: for (i=eps->nconv;i<k+l;i++) {
176: VecCopy(eps->AV[i],eps->V[i]);
177: }
178: eps->nconv = k;
180: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nv);
181:
182: if (eps->reason == EPS_CONVERGED_ITERATING) {
183: if (breakdown) {
184: /* start a new Arnoldi factorization */
185: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%i norm=%g)\n",eps->its,beta);
186: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
187: if (breakdown) {
188: eps->reason = EPS_DIVERGED_BREAKDOWN;
189: PetscInfo(eps,"Unable to generate more start vectors\n");
190: }
191: } else {
192: /* update the Arnoldi-Schur decomposition */
193: for (i=k;i<k+l;i++) {
194: S[i*eps->ncv+k+l] = Q[(i+1)*n-1]*beta;
195: }
196: VecCopy(u,eps->V[k+l]);
197: }
198: }
199: }
201: PetscFree(Q);
202: PetscFree(b);
203: if (!eps->ishermitian) {
204: PetscFree(work);
205: } else {
206: PetscFree(ritz);
207: PetscFree(perm);
208: }
209: return(0);
210: }
215: PetscErrorCode EPSCreate_KRYLOVSCHUR(EPS eps)
216: {
218: eps->data = PETSC_NULL;
219: eps->ops->solve = EPSSolve_KRYLOVSCHUR;
220: eps->ops->solvets = PETSC_NULL;
221: eps->ops->setup = EPSSetUp_KRYLOVSCHUR;
222: eps->ops->setfromoptions = PETSC_NULL;
223: eps->ops->destroy = EPSDestroy_Default;
224: eps->ops->view = PETSC_NULL;
225: eps->ops->backtransform = EPSBackTransform_Default;
226: eps->ops->computevectors = EPSComputeVectors_Schur;
227: return(0);
228: }