Actual source code: rii.c
slepc-3.6.1 2015-09-03
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: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
17: SLEPc - Scalable Library for Eigenvalue Problem Computations
18: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
20: This file is part of SLEPc.
22: SLEPc is free software: you can redistribute it and/or modify it under the
23: terms of version 3 of the GNU Lesser General Public License as published by
24: the Free Software Foundation.
26: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
27: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
28: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
29: more details.
31: You should have received a copy of the GNU Lesser General Public License
32: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
33: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: */
36: #include <slepc/private/nepimpl.h>
40: PetscErrorCode NEPSetUp_RII(NEP nep)
41: {
43: PetscBool istrivial;
46: if (nep->ncv) { /* ncv set */
47: if (nep->ncv<nep->nev) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must be at least nev");
48: } else if (nep->mpd) { /* mpd set */
49: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
50: } else { /* neither set: defaults depend on nev being small or large */
51: if (nep->nev<500) nep->ncv = PetscMin(nep->n,PetscMax(2*nep->nev,nep->nev+15));
52: else {
53: nep->mpd = 500;
54: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
55: }
56: }
57: if (!nep->mpd) nep->mpd = nep->ncv;
58: if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
59: if (nep->nev>1) { PetscInfo(nep,"Warning: requested more than one eigenpair but RII can only compute one\n"); }
60: if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
61: if (!nep->max_funcs) nep->max_funcs = nep->max_it;
63: RGIsTrivial(nep->rg,&istrivial);
64: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support region filtering");
66: NEPAllocateSolution(nep,0);
67: NEPSetWorkVecs(nep,2);
68: return(0);
69: }
73: PetscErrorCode NEPSolve_RII(NEP nep)
74: {
75: PetscErrorCode ierr;
76: Mat T=nep->function,Tp=nep->jacobian,Tsigma;
77: Vec u,r=nep->work[0],delta=nep->work[1];
78: PetscScalar lambda,a1,a2;
79: PetscReal relerr;
80: PetscBool hascopy;
81: KSPConvergedReason kspreason;
84: /* get initial approximation of eigenvalue and eigenvector */
85: NEPGetDefaultShift(nep,&lambda);
86: if (!nep->nini) {
87: BVSetRandomColumn(nep->V,0,nep->rand);
88: }
89: BVGetColumn(nep->V,0,&u);
91: /* correct eigenvalue approximation: lambda = lambda - (u'*T*u)/(u'*Tp*u) */
92: NEPComputeFunction(nep,lambda,T,T);
93: MatMult(T,u,r);
94: VecDot(u,r,&a1);
95: NEPApplyJacobian(nep,lambda,u,delta,r,Tp);
96: VecDot(u,r,&a2);
97: lambda = lambda - a1/a2;
99: /* prepare linear solver */
100: MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
101: KSPSetOperators(nep->ksp,Tsigma,Tsigma);
103: /* Restart loop */
104: while (nep->reason == NEP_CONVERGED_ITERATING) {
105: nep->its++;
107: /* update preconditioner and set adaptive tolerance */
108: if (nep->lag && !(nep->its%nep->lag) && nep->its>2*nep->lag && relerr<1e-2) {
109: MatHasOperation(T,MATOP_COPY,&hascopy);
110: if (hascopy) {
111: MatCopy(T,Tsigma,DIFFERENT_NONZERO_PATTERN);
112: } else {
113: MatDestroy(&Tsigma);
114: MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
115: }
116: KSPSetOperators(nep->ksp,Tsigma,Tsigma);
117: }
118: if (!nep->cctol) {
119: nep->ktol = PetscMax(nep->ktol/2.0,PETSC_MACHINE_EPSILON*10.0);
120: KSPSetTolerances(nep->ksp,nep->ktol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
121: }
123: /* form residual, r = T(lambda)*u */
124: NEPApplyFunction(nep,lambda,u,delta,r,T,T);
126: /* convergence test */
127: VecNorm(r,NORM_2,&relerr);
128: nep->errest[nep->nconv] = relerr;
129: nep->eigr[nep->nconv] = lambda;
130: if (relerr<=nep->rtol) {
131: nep->nconv = nep->nconv + 1;
132: nep->reason = NEP_CONVERGED_FNORM_RELATIVE;
133: }
134: NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->errest,1);
136: if (!nep->nconv) {
137: /* eigenvector correction: delta = T(sigma)\r */
138: NEP_KSPSolve(nep,r,delta);
139: KSPGetConvergedReason(nep->ksp,&kspreason);
140: if (kspreason<0) {
141: PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
142: nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
143: break;
144: }
146: /* update eigenvector: u = u - delta */
147: VecAXPY(u,-1.0,delta);
149: /* normalize eigenvector */
150: VecNormalize(u,NULL);
152: /* correct eigenvalue: lambda = lambda - (u'*T*u)/(u'*Tp*u) */
153: NEPApplyFunction(nep,lambda,u,delta,r,T,T);
154: VecDot(u,r,&a1);
155: NEPApplyJacobian(nep,lambda,u,delta,r,Tp);
156: VecDot(u,r,&a2);
157: lambda = lambda - a1/a2;
158: }
159: if (nep->its >= nep->max_it) nep->reason = NEP_DIVERGED_MAX_IT;
160: }
161: MatDestroy(&Tsigma);
162: BVRestoreColumn(nep->V,0,&u);
163: return(0);
164: }
168: PETSC_EXTERN PetscErrorCode NEPCreate_RII(NEP nep)
169: {
171: nep->ops->solve = NEPSolve_RII;
172: nep->ops->setup = NEPSetUp_RII;
173: return(0);
174: }