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: #include src/eps/epsimpl.h
24: #include slepcblaslapack.h
28: PetscErrorCode EPSSetUp_KRYLOVSCHUR(EPS eps)
29: {
31: PetscInt N;
34: VecGetSize(eps->vec_initial,&N);
35: if (eps->ncv) {
36: if (eps->ncv<eps->nev+1) SETERRQ(1,"The value of ncv must be at least nev+1");
37: }
38: else eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
39: if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
40: if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
41: SETERRQ(1,"Wrong value of eps->which");
43: EPSAllocateSolution(eps);
44: PetscFree(eps->T);
45: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
46: EPSDefaultGetWork(eps,1);
47: return(0);
48: }
52: PetscErrorCode EPSSolve_KRYLOVSCHUR(EPS eps)
53: {
55: int i,j,k,l,n,lwork,*perm;
56: Vec u=eps->work[0];
57: PetscScalar *S=eps->T,*Q,*work,*b;
58: PetscReal beta,*ritz;
59: PetscTruth breakdown;
62: PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));
63: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);
64: PetscMalloc(eps->ncv*sizeof(PetscScalar),&b);
65: lwork = (eps->ncv+4)*eps->ncv;
66: if (!eps->ishermitian) {
67: PetscMalloc(lwork*sizeof(PetscScalar),&work);
68: } else {
69: PetscMalloc(eps->ncv*sizeof(PetscReal),&ritz);
70: PetscMalloc(eps->ncv*sizeof(int),&perm);
71: }
73: /* Get the starting Arnoldi vector */
74: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
75: l = 0;
76:
77: /* Restart loop */
78: while (eps->reason == EPS_CONVERGED_ITERATING) {
79: eps->its++;
81: /* Compute an nv-step Arnoldi factorization */
82: eps->nv = eps->ncv;
83: EPSBasicArnoldi(eps,PETSC_FALSE,S,eps->V,eps->nconv+l,&eps->nv,u,&beta,&breakdown);
84: VecScale(u,1.0/beta);
86: if (!eps->ishermitian) {
87: n = eps->nv; /* size of Q */
88: if (l==0) {
89: PetscMemzero(Q,n*n*sizeof(PetscScalar));
90: for (i=0;i<n;i++)
91: Q[i*(n+1)] = 1.0;
92: } else {
93: /* Reduce S to Hessenberg form, S <- Q S Q' */
94: EPSDenseHessenberg(n,eps->nconv,S,eps->ncv,Q);
95: }
96: /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
97: EPSDenseSchur(n,eps->nconv,S,eps->ncv,Q,eps->eigr,eps->eigi);
98: /* Sort the remaining columns of the Schur form */
99: EPSSortDenseSchur(n,eps->nconv,S,eps->ncv,Q,eps->eigr,eps->eigi,eps->which);
100: /* Compute residual norm estimates */
101: ArnoldiResiduals(S,eps->ncv,Q,beta,eps->nconv,n,eps->eigr,eps->eigi,eps->errest,work);
102: } else {
103: n = eps->nv-eps->nconv; /* size of Q */
104: /* Reduce S to diagonal form, S <- Q S Q' */
105: if (l==0) {
106: EPSDenseTridiagonal(n,S+eps->nconv*(eps->ncv+1),eps->ncv,ritz,Q+eps->nconv*n);
107: } else {
108: EPSDenseHEP(n,S+eps->nconv*(eps->ncv+1),eps->ncv,ritz,Q+eps->nconv*n);
109: }
110: /* Sort the remaining columns of the Schur form */
111: if (eps->which == EPS_SMALLEST_REAL) {
112: for (i=0;i<n;i++)
113: eps->eigr[i+eps->nconv] = ritz[i];
114: } else {
115: #ifdef PETSC_USE_COMPLEX
116: for (i=0;i<n;i++)
117: eps->eigr[i+eps->nconv] = ritz[i];
118: EPSSortEigenvalues(n,eps->eigr+eps->nconv,eps->eigi,eps->which,n,perm);
119: #else
120: EPSSortEigenvalues(n,ritz,eps->eigi+eps->nconv,eps->which,n,perm);
121: #endif
122: for (i=0;i<n;i++)
123: eps->eigr[i+eps->nconv] = ritz[perm[i]];
124: PetscMemcpy(S,Q+eps->nconv*n,n*n*sizeof(PetscScalar));
125: for (j=0;j<n;j++)
126: for (i=0;i<n;i++)
127: Q[(j+eps->nconv)*n+i] = S[perm[j]*n+i];
128: }
129: /* rebuild S from eigr */
130: for (i=eps->nconv;i<eps->nv;i++) {
131: S[i*(eps->ncv+1)] = eps->eigr[i];
132: for (j=i+1;j<eps->ncv;j++)
133: S[i*eps->ncv+j] = 0.0;
134: }
135: /* Compute residual norm estimates */
136: for (i=eps->nconv;i<eps->nv;i++)
137: eps->errest[i] = beta*PetscAbsScalar(Q[(i+1)*n-1]) / PetscAbsScalar(eps->eigr[i]);
138: }
140: /* Check convergence */
141: k = eps->nconv;
142: while (k<eps->nv && eps->errest[k]<eps->tol) k++;
143: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
144: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
145:
146: /* Update l */
147: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
148: else {
149: l = (eps->nv-k)/2;
150: #if !defined(PETSC_USE_COMPLEX)
151: if (S[(k+l-1)*(eps->ncv+1)+1] != 0.0) {
152: if (k+l<eps->nv-1) l = l+1;
153: else l = l-1;
154: }
155: #endif
156: }
157:
158: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
159: for (i=eps->nconv;i<k+l;i++) {
160: VecSet(eps->AV[i],0.0);
161: if (!eps->ishermitian) {
162: VecMAXPY(eps->AV[i],n,Q+i*n,eps->V);
163: } else {
164: VecMAXPY(eps->AV[i],n,Q+i*n,eps->V+eps->nconv);
165: }
166: }
167: for (i=eps->nconv;i<k+l;i++) {
168: VecCopy(eps->AV[i],eps->V[i]);
169: }
170: eps->nconv = k;
172: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nv);
173:
174: if (eps->reason == EPS_CONVERGED_ITERATING) {
175: if (breakdown) {
176: /* start a new Arnoldi factorization */
177: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%i norm=%g)\n",eps->its,beta);
178: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
179: if (breakdown) {
180: eps->reason = EPS_DIVERGED_BREAKDOWN;
181: PetscInfo(eps,"Unable to generate more start vectors\n");
182: }
183: } else {
184: /* update the Arnoldi-Schur decomposition */
185: for (i=k;i<k+l;i++) {
186: S[i*eps->ncv+k+l] = Q[(i+1)*n-1]*beta;
187: }
188: VecCopy(u,eps->V[k+l]);
189: }
190: }
191: }
193: PetscFree(Q);
194: PetscFree(b);
195: if (!eps->ishermitian) {
196: PetscFree(work);
197: } else {
198: PetscFree(ritz);
199: PetscFree(perm);
200: }
201: return(0);
202: }
207: PetscErrorCode EPSCreate_KRYLOVSCHUR(EPS eps)
208: {
210: eps->data = PETSC_NULL;
211: eps->ops->solve = EPSSolve_KRYLOVSCHUR;
212: eps->ops->solvets = PETSC_NULL;
213: eps->ops->setup = EPSSetUp_KRYLOVSCHUR;
214: eps->ops->setfromoptions = PETSC_NULL;
215: eps->ops->destroy = EPSDestroy_Default;
216: eps->ops->view = PETSC_NULL;
217: eps->ops->backtransform = EPSBackTransform_Default;
218: eps->ops->computevectors = EPSComputeVectors_Schur;
219: return(0);
220: }