Actual source code: narnoldi.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  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:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 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_NArnoldi(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->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
 60:   if (!nep->max_funcs) nep->max_funcs = nep->max_it;
 61:   if (!nep->split) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NARNOLDI only available for split operator");

 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,3);

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

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

 90:   /* get initial space and shift */
 91:   NEPGetDefaultShift(nep,&lambda);
 92:   if (!nep->nini) {
 93:     BVSetRandomColumn(nep->V,0,nep->rand);
 94:     BVNormColumn(nep->V,0,NORM_2,&nrm);
 95:     BVScaleColumn(nep->V,0,1.0/nrm);
 96:     n = 1;
 97:   } else n = nep->nini;

 99:   /* build projected matrices for initial space */
100:   DSSetDimensions(nep->ds,n,0,0,0);
101:   NEPProjectOperator(nep,0,n);

103:   /* prepare linear solver */
104:   NEPComputeFunction(nep,lambda,T,T);
105:   MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
106:   KSPSetOperators(nep->ksp,Tsigma,Tsigma);

108:   /* Restart loop */
109:   while (nep->reason == NEP_CONVERGED_ITERATING) {
110:     nep->its++;

112:     /* solve projected problem */
113:     DSSetDimensions(nep->ds,n,0,0,0);
114:     DSSetState(nep->ds,DS_STATE_RAW);
115:     DSSolve(nep->ds,nep->eigr,NULL);
116:     lambda = nep->eigr[0];

118:     /* compute Ritz vector, x = V*s */
119:     DSGetArray(nep->ds,DS_MAT_X,&X);
120:     BVSetActiveColumns(nep->V,0,n);
121:     BVMultVec(nep->V,1.0,0.0,x,X);
122:     DSRestoreArray(nep->ds,DS_MAT_X,&X);

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

127:     /* convergence test */
128:     VecNorm(r,NORM_2,&resnorm);
129:     nep->errest[nep->nconv] = resnorm;
130:     if (resnorm<=nep->rtol) {
131:       BVInsertVec(nep->V,nep->nconv,x);
132:       nep->nconv = nep->nconv + 1;
133:       nep->reason = NEP_CONVERGED_FNORM_RELATIVE;
134:     }
135:     NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->errest,1);

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

139:       /* continuation vector: f = T(sigma)\r */
140:       BVGetColumn(nep->V,n,&f);
141:       NEP_KSPSolve(nep,r,f);
142:       BVRestoreColumn(nep->V,n,&f);
143:       KSPGetConvergedReason(nep->ksp,&kspreason);
144:       if (kspreason<0) {
145:         PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
146:         nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
147:         break;
148:       }

150:       /* orthonormalize */
151:       BVOrthogonalizeColumn(nep->V,n,NULL,&beta,&breakdown);
152:       if (breakdown || beta==0.0) {
153:         PetscInfo1(nep,"iter=%D, orthogonalization failed, stopping solve\n",nep->its);
154:         nep->reason = NEP_DIVERGED_BREAKDOWN;
155:         break;
156:       }
157:       BVScaleColumn(nep->V,n,1.0/beta);

159:       /* update projected matrices */
160:       DSSetDimensions(nep->ds,n+1,0,0,0);
161:       NEPProjectOperator(nep,n,n+1);
162:       n++;
163:     }
164:     if (nep->its >= nep->max_it) nep->reason = NEP_DIVERGED_MAX_IT;
165:   }
166:   MatDestroy(&Tsigma);
167:   return(0);
168: }

172: PETSC_EXTERN PetscErrorCode NEPCreate_NArnoldi(NEP nep)
173: {
175:   nep->ops->solve          = NEPSolve_NArnoldi;
176:   nep->ops->setup          = NEPSetUp_NArnoldi;
177:   return(0);
178: }