Actual source code: narnoldi.c

  1: /*

  3:    SLEPc nonlinear eigensolver: "narnoldi"

  5:    Method: Nonlinear Arnoldi

  7:    Algorithm:

  9:        Arnoldi for nonlinear eigenproblems.

 11:    References:

 13:        [1] H. Voss, "An Arnoldi method for nonlinear eigenvalue problems",
 14:            BIT 44:387-401, 2004.

 16:    Last update: Mar 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_NARNOLDI(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->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
 61:   if (!nep->max_funcs) nep->max_funcs = nep->max_it;
 62:   if (!nep->split) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NARNOLDI only available for split operator");

 64:   NEPAllocateSolution(nep);
 65:   NEPSetWorkVecs(nep,3);

 67:   /* set-up DS and transfer split operator functions */
 68:   DSSetType(nep->ds,DSNEP);
 69:   DSSetFN(nep->ds,nep->nt,nep->f);
 70:   DSAllocate(nep->ds,nep->ncv);
 71:   return(0);
 72: }

 76: PetscErrorCode NEPSolve_NARNOLDI(NEP nep)
 77: {
 78:   PetscErrorCode     ierr;
 79:   Mat                T=nep->function,Tsigma;
 80:   Vec                f,u=nep->V[0],r=nep->work[0],x=nep->work[1],w=nep->work[2];
 81:   PetscScalar        *X,lambda;
 82:   PetscReal          beta,resnorm=0.0;
 83:   PetscInt           n;
 84:   PetscBool          breakdown;
 85:   MatStructure       mats;
 86:   KSPConvergedReason kspreason;

 89:   /* get initial space and shift */
 90:   NEPGetDefaultShift(nep,&lambda);
 91:   if (!nep->nini) {
 92:     SlepcVecSetRandom(u,nep->rand);
 93:     VecNormalize(u,NULL);
 94:     n = 1;
 95:   } else n = nep->nini;

 97:   /* build projected matrices for initial space */
 98:   NEPProjectOperator(nep,0,n,r);

100:   /* prepare linear solver */
101:   NEPComputeFunction(nep,lambda,&T,&T,&mats);
102:   MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
103:   KSPSetOperators(nep->ksp,Tsigma,Tsigma,DIFFERENT_NONZERO_PATTERN);

105:   /* Restart loop */
106:   while (nep->reason == NEP_CONVERGED_ITERATING) {
107:     nep->its++;

109:     /* solve projected problem */
110:     DSSetDimensions(nep->ds,n,0,0,0);
111:     DSSetState(nep->ds,DS_STATE_RAW);
112:     DSSolve(nep->ds,nep->eig,NULL);
113:     lambda = nep->eig[0];

115:     /* compute Ritz vector, x = V*s */
116:     DSGetArray(nep->ds,DS_MAT_X,&X);
117:     SlepcVecMAXPBY(x,0.0,1.0,n,X,nep->V);
118:     DSRestoreArray(nep->ds,DS_MAT_X,&X);

120:     /* compute the residual, r = T(lambda)*x */
121:     NEPApplyFunction(nep,lambda,x,w,r,NULL,NULL,NULL);

123:     /* convergence test */
124:     VecNorm(r,NORM_2,&resnorm);
125:     nep->errest[nep->nconv] = resnorm;
126:     if (resnorm<=nep->rtol) {
127:       VecCopy(x,nep->V[nep->nconv]);
128:       nep->nconv = nep->nconv + 1;
129:       nep->reason = NEP_CONVERGED_FNORM_RELATIVE;
130:     }
131:     NEPMonitor(nep,nep->its,nep->nconv,nep->eig,nep->errest,1);

133:     if (nep->reason == NEP_CONVERGED_ITERATING) {

135:       /* continuation vector: f = T(sigma)\r */
136:       f = nep->V[n];
137:       NEP_KSPSolve(nep,r,f);
138:       KSPGetConvergedReason(nep->ksp,&kspreason);
139:       if (kspreason<0) {
140:         PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
141:         nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
142:         break;
143:       }

145:       /* orthonormalize */
146:       IPOrthogonalize(nep->ip,0,NULL,n,NULL,nep->V,f,NULL,&beta,&breakdown);
147:       if (breakdown || beta==0.0) {
148:         PetscInfo1(nep,"iter=%D, orthogonalization failed, stopping solve\n",nep->its);
149:         nep->reason = NEP_DIVERGED_BREAKDOWN;
150:         break;
151:       }
152:       VecScale(f,1.0/beta);

154:       /* update projected matrices */
155:       NEPProjectOperator(nep,n,n+1,r);
156:       n++;
157:     }
158:     if (nep->its >= nep->max_it) nep->reason = NEP_DIVERGED_MAX_IT;
159:   }
160:   MatDestroy(&Tsigma);
161:   return(0);
162: }

166: PETSC_EXTERN PetscErrorCode NEPCreate_NARNOLDI(NEP nep)
167: {
169:   nep->ops->solve          = NEPSolve_NARNOLDI;
170:   nep->ops->setup          = NEPSetUp_NARNOLDI;
171:   nep->ops->reset          = NEPReset_Default;
172:   return(0);
173: }