Actual source code: arnoldi-ts.c
1: /*
3: SLEPc eigensolver: "arnoldi"
5: Method: Explicitly Restarted Arnoldi (two-sided)
7: */
8: #include src/eps/epsimpl.h
9: #include slepcblaslapack.h
13: PetscErrorCode EPSSolve_TS_ARNOLDI(EPS eps)
14: {
16: int i,k,ncv=eps->ncv;
17: Vec fr=eps->work[0];
18: Vec fl=eps->work[1];
19: Vec *Qr=eps->V, *Ql=eps->W;
20: PetscScalar *Hr=eps->T,*Ur,*work;
21: PetscScalar *Hl=eps->Tl,*Ul;
22: PetscReal beta,g;
23: PetscScalar *eigr,*eigi,*aux;
26: PetscMemzero(Hr,ncv*ncv*sizeof(PetscScalar));
27: PetscMemzero(Hl,ncv*ncv*sizeof(PetscScalar));
28: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Ur);
29: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Ul);
30: PetscMalloc((ncv+4)*ncv*sizeof(PetscScalar),&work);
31: PetscMalloc(ncv*sizeof(PetscScalar),&eigr);
32: PetscMalloc(ncv*sizeof(PetscScalar),&eigi);
33: PetscMalloc(ncv*sizeof(PetscScalar),&aux);
35: eps->nconv = 0;
36: eps->its = 0;
37: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,ncv);
39: /* Get the starting Arnoldi vector */
40: EPSGetStartVector(eps,eps->its,Qr[0],PETSC_NULL);
41: EPSGetLeftStartVector(eps,eps->its,Ql[0]);
42:
43: /* Restart loop */
44: while (eps->its<eps->max_it) {
46: /* Compute an ncv-step Arnoldi factorization for both A and A' */
47: EPSBasicArnoldi(eps,PETSC_FALSE,Hr,Qr,eps->nconv,&ncv,fr,&beta);
48: EPSBasicArnoldi(eps,PETSC_TRUE,Hl,Ql,eps->nconv,&ncv,fl,&g);
50: EPSBiOrthogonalize(eps,ncv,Qr,Ql,fr,aux,PETSC_NULL);
51: for (i=0;i<ncv;i++) {
52: Hr[ncv*(ncv-1)+i] += beta * aux[i];
53: }
54: EPSBiOrthogonalize(eps,ncv,Ql,Qr,fl,aux,PETSC_NULL);
55: for (i=0;i<ncv;i++) {
56: Hl[ncv*(ncv-1)+i] += g * aux[i];
57: }
59: /* Reduce H to (quasi-)triangular form, H <- U H U' */
60: PetscMemzero(Ur,ncv*ncv*sizeof(PetscScalar));
61: for (i=0;i<ncv;i++) { Ur[i*(ncv+1)] = 1.0; }
62: EPSDenseSchur(ncv,eps->nconv,Hr,ncv,Ur,eps->eigr,eps->eigi);
64: PetscMemzero(Ul,ncv*ncv*sizeof(PetscScalar));
65: for (i=0;i<ncv;i++) { Ul[i*(ncv+1)] = 1.0; }
66: EPSDenseSchur(ncv,eps->nconv,Hl,ncv,Ul,eigr,eigi);
68: /* Sort the remaining columns of the Schur form */
69: EPSSortDenseSchur(ncv,eps->nconv,Hr,ncv,Ur,eps->eigr,eps->eigi);
70: EPSSortDenseSchur(ncv,eps->nconv,Hl,ncv,Ul,eigr,eigi);
72: /* Compute residual norm estimates */
73: ArnoldiResiduals(Hr,ncv,Ur,beta,eps->nconv,ncv,eps->eigr,eps->eigi,eps->errest,work);
74: ArnoldiResiduals(Hl,ncv,Ul,g,eps->nconv,ncv,eigr,eigi,eps->errest_left,work);
76: /* Lock converged eigenpairs and update the corresponding vectors,
77: including the restart vector: V(:,idx) = V*U(:,idx) */
78: k = eps->nconv;
79: while (k<ncv && eps->errest[k]<eps->tol && eps->errest_left[k]<eps->tol) k++;
80: for (i=eps->nconv;i<=k && i<ncv;i++) {
81: VecSet(eps->AV[i],0.0);
82: VecMAXPY(eps->AV[i],ncv,Ur+ncv*i,Qr);
83: }
84: for (i=eps->nconv;i<=k && i<ncv;i++) {
85: VecCopy(eps->AV[i],Qr[i]);
86: }
87: for (i=eps->nconv;i<=k && i<ncv;i++) {
88: VecSet(eps->AW[i],0.0);
89: VecMAXPY(eps->AW[i],ncv,Ul+ncv*i,Ql);
90: }
91: for (i=eps->nconv;i<=k && i<ncv;i++) {
92: VecCopy(eps->AW[i],Ql[i]);
93: }
94: eps->nconv = k;
96: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,ncv);
97: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest_left,ncv);
98: if (eps->nconv >= eps->nev) break;
99: }
100:
101: if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
102: else eps->reason = EPS_DIVERGED_ITS;
103: #if defined(PETSC_USE_COMPLEX)
104: for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;
105: #endif
107: PetscFree(Ur);
108: PetscFree(Ul);
109: PetscFree(work);
110: PetscFree(eigr);
111: PetscFree(eigi);
112: PetscFree(aux);
113: return(0);
114: }