Actual source code: rii.c
1: /*
3: SLEPc nonlinear eigensolver: "rii"
5: Method: Residual inverse iteration
7: Algorithm:
9: Simple residual inverse iteration with varying shift.
11: References:
13: [1] A. Neumaier, "Residual inverse iteration for the nonlinear
14: eigenvalue problem", SIAM J. Numer. Anal. 22(5):914-923, 1985.
16: Last update: Feb 2013
18: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
19: SLEPc - Scalable Library for Eigenvalue Problem Computations
20: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
22: This file is part of SLEPc.
24: SLEPc is free software: you can redistribute it and/or modify it under the
25: terms of version 3 of the GNU Lesser General Public License as published by
26: the Free Software Foundation.
28: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
29: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
30: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
31: more details.
33: You should have received a copy of the GNU Lesser General Public License
34: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: */
38: #include <slepc-private/nepimpl.h> /*I "slepcnep.h" I*/
42: PetscErrorCode NEPSetUp_RII(NEP nep)
43: {
47: if (nep->ncv) { /* ncv set */
48: if (nep->ncv<nep->nev) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must be at least nev");
49: } else if (nep->mpd) { /* mpd set */
50: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
51: } else { /* neither set: defaults depend on nev being small or large */
52: if (nep->nev<500) nep->ncv = PetscMin(nep->n,PetscMax(2*nep->nev,nep->nev+15));
53: else {
54: nep->mpd = 500;
55: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
56: }
57: }
58: if (!nep->mpd) nep->mpd = nep->ncv;
59: if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
60: if (nep->nev>1) { PetscInfo(nep,"Warning: requested more than one eigenpair but RII can only compute one\n"); }
61: if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
62: if (!nep->max_funcs) nep->max_funcs = nep->max_it;
64: NEPAllocateSolution(nep);
65: NEPSetWorkVecs(nep,2);
66: return(0);
67: }
71: PetscErrorCode NEPSolve_RII(NEP nep)
72: {
73: PetscErrorCode ierr;
74: Mat T=nep->function,Tp=nep->jacobian,Tsigma;
75: Vec u=nep->V[0],r=nep->work[0],delta=nep->work[1];
76: PetscScalar lambda,a1,a2;
77: PetscReal relerr;
78: PetscBool hascopy;
79: MatStructure mats;
80: KSPConvergedReason kspreason;
83: /* get initial approximation of eigenvalue and eigenvector */
84: NEPGetDefaultShift(nep,&lambda);
85: if (!nep->nini) {
86: SlepcVecSetRandom(u,nep->rand);
87: }
89: /* correct eigenvalue approximation: lambda = lambda - (u'*T*u)/(u'*Tp*u) */
90: NEPComputeFunction(nep,lambda,&T,&T,&mats);
91: MatMult(T,u,r);
92: VecDot(u,r,&a1);
93: NEPApplyJacobian(nep,lambda,u,delta,r,&Tp,&mats);
94: VecDot(u,r,&a2);
95: lambda = lambda - a1/a2;
97: /* prepare linear solver */
98: MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
99: KSPSetOperators(nep->ksp,Tsigma,Tsigma,DIFFERENT_NONZERO_PATTERN);
101: /* Restart loop */
102: while (nep->reason == NEP_CONVERGED_ITERATING) {
103: nep->its++;
105: /* update preconditioner and set adaptive tolerance */
106: if (nep->lag && !(nep->its%nep->lag) && nep->its>2*nep->lag && relerr<1e-2) {
107: MatHasOperation(T,MATOP_COPY,&hascopy);
108: if (hascopy) {
109: MatCopy(T,Tsigma,mats);
110: } else {
111: MatDestroy(&Tsigma);
112: MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
113: }
114: KSPSetOperators(nep->ksp,Tsigma,Tsigma,mats);
115: }
116: if (!nep->cctol) {
117: nep->ktol = PetscMax(nep->ktol/2.0,PETSC_MACHINE_EPSILON*10.0);
118: KSPSetTolerances(nep->ksp,nep->ktol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
119: }
121: /* form residual, r = T(lambda)*u */
122: NEPApplyFunction(nep,lambda,u,delta,r,&T,&T,&mats);
124: /* convergence test */
125: VecNorm(r,NORM_2,&relerr);
126: nep->errest[nep->nconv] = relerr;
127: nep->eig[nep->nconv] = lambda;
128: if (relerr<=nep->rtol) {
129: nep->nconv = nep->nconv + 1;
130: nep->reason = NEP_CONVERGED_FNORM_RELATIVE;
131: }
132: NEPMonitor(nep,nep->its,nep->nconv,nep->eig,nep->errest,1);
134: if (!nep->nconv) {
135: /* eigenvector correction: delta = T(sigma)\r */
136: NEP_KSPSolve(nep,r,delta);
137: KSPGetConvergedReason(nep->ksp,&kspreason);
138: if (kspreason<0) {
139: PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
140: nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
141: break;
142: }
144: /* update eigenvector: u = u - delta */
145: VecAXPY(u,-1.0,delta);
147: /* normalize eigenvector */
148: VecNormalize(u,NULL);
150: /* correct eigenvalue: lambda = lambda - (u'*T*u)/(u'*Tp*u) */
151: NEPApplyFunction(nep,lambda,u,delta,r,&T,&T,&mats);
152: VecDot(u,r,&a1);
153: NEPApplyJacobian(nep,lambda,u,delta,r,&Tp,&mats);
154: VecDot(u,r,&a2);
155: lambda = lambda - a1/a2;
156: }
157: if (nep->its >= nep->max_it) nep->reason = NEP_DIVERGED_MAX_IT;
158: }
159: MatDestroy(&Tsigma);
160: return(0);
161: }
165: PETSC_EXTERN PetscErrorCode NEPCreate_RII(NEP nep)
166: {
168: nep->ops->solve = NEPSolve_RII;
169: nep->ops->setup = NEPSetUp_RII;
170: nep->ops->reset = NEPReset_Default;
171: return(0);
172: }