Actual source code: arnoldi-ts.c
1: /*
3: SLEPc eigensolver: "arnoldi"
5: Method: Explicitly Restarted Arnoldi (two-sided)
7: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
8: SLEPc - Scalable Library for Eigenvalue Problem Computations
9: Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
11: This file is part of SLEPc.
12:
13: SLEPc is free software: you can redistribute it and/or modify it under the
14: terms of version 3 of the GNU Lesser General Public License as published by
15: the Free Software Foundation.
17: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
18: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
20: more details.
22: You should have received a copy of the GNU Lesser General Public License
23: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
24: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25: */
27: #include private/epsimpl.h
28: #include slepcblaslapack.h
32: PetscErrorCode EPSSolve_TS_ARNOLDI(EPS eps)
33: {
35: PetscInt i,k,ncv=eps->ncv;
36: Vec fr=eps->work[0];
37: Vec fl=eps->work[1];
38: Vec *Qr=eps->V, *Ql=eps->W;
39: PetscScalar *Hr=eps->T,*Ur,*work;
40: PetscScalar *Hl=eps->Tl,*Ul;
41: PetscReal beta,g;
42: PetscScalar *eigr,*eigi,*aux;
43: PetscTruth breakdown;
46: PetscMemzero(Hr,ncv*ncv*sizeof(PetscScalar));
47: PetscMemzero(Hl,ncv*ncv*sizeof(PetscScalar));
48: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Ur);
49: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Ul);
50: PetscMalloc((ncv+4)*ncv*sizeof(PetscScalar),&work);
51: PetscMalloc(ncv*sizeof(PetscScalar),&eigr);
52: PetscMalloc(ncv*sizeof(PetscScalar),&eigi);
53: PetscMalloc(ncv*sizeof(PetscScalar),&aux);
55: /* Get the starting Arnoldi vector */
56: EPSGetStartVector(eps,eps->its,Qr[0],PETSC_NULL);
57: EPSGetLeftStartVector(eps,eps->its,Ql[0]);
58:
59: /* Restart loop */
60: while (eps->its<eps->max_it) {
61: eps->its++;
63: /* Compute an ncv-step Arnoldi factorization for both A and A' */
64: EPSBasicArnoldi(eps,PETSC_FALSE,Hr,ncv,Qr,eps->nconv,&ncv,fr,&beta,&breakdown);
65: EPSBasicArnoldi(eps,PETSC_TRUE,Hl,ncv,Ql,eps->nconv,&ncv,fl,&g,&breakdown);
67: IPBiOrthogonalize(eps->ip,ncv,Qr,Ql,fr,aux,PETSC_NULL);
68: for (i=0;i<ncv;i++) {
69: Hr[ncv*(ncv-1)+i] += beta * aux[i];
70: }
71: IPBiOrthogonalize(eps->ip,ncv,Ql,Qr,fl,aux,PETSC_NULL);
72: for (i=0;i<ncv;i++) {
73: Hl[ncv*(ncv-1)+i] += g * aux[i];
74: }
76: /* Reduce H to (quasi-)triangular form, H <- U H U' */
77: PetscMemzero(Ur,ncv*ncv*sizeof(PetscScalar));
78: for (i=0;i<ncv;i++) { Ur[i*(ncv+1)] = 1.0; }
79: EPSDenseSchur(ncv,eps->nconv,Hr,ncv,Ur,eps->eigr,eps->eigi);
81: PetscMemzero(Ul,ncv*ncv*sizeof(PetscScalar));
82: for (i=0;i<ncv;i++) { Ul[i*(ncv+1)] = 1.0; }
83: EPSDenseSchur(ncv,eps->nconv,Hl,ncv,Ul,eigr,eigi);
85: /* Sort the remaining columns of the Schur form */
86: EPSSortDenseSchur(ncv,eps->nconv,Hr,ncv,Ur,eps->eigr,eps->eigi,eps->which);
87: EPSSortDenseSchur(ncv,eps->nconv,Hl,ncv,Ul,eigr,eigi,eps->which);
89: /* Compute residual norm estimates */
90: ArnoldiResiduals(Hr,ncv,Ur,beta,eps->nconv,ncv,eps->eigr,eps->eigi,eps->errest,work);
91: ArnoldiResiduals(Hl,ncv,Ul,g,eps->nconv,ncv,eigr,eigi,eps->errest_left,work);
93: /* Lock converged eigenpairs and update the corresponding vectors,
94: including the restart vector: V(:,idx) = V*U(:,idx) */
95: k = eps->nconv;
96: while (k<ncv && eps->errest[k]<eps->tol && eps->errest_left[k]<eps->tol) k++;
97: SlepcUpdateVectors(ncv,Qr,eps->nconv,PetscMin(k+1,ncv),Ur,ncv,PETSC_FALSE);
98: SlepcUpdateVectors(ncv,Ql,eps->nconv,PetscMin(k+1,ncv),Ul,ncv,PETSC_FALSE);
99: eps->nconv = k;
101: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,ncv);
102: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest_left,ncv);
103: if (eps->nconv >= eps->nev) break;
104: }
105:
106: if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
107: else eps->reason = EPS_DIVERGED_ITS;
109: PetscFree(Ur);
110: PetscFree(Ul);
111: PetscFree(work);
112: PetscFree(eigr);
113: PetscFree(eigi);
114: PetscFree(aux);
115: return(0);
116: }