Actual source code: rii.c

slepc-3.13.1 2020-04-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc nonlinear eigensolver: "rii"

 13:    Method: Residual inverse iteration

 15:    Algorithm:

 17:        Simple residual inverse iteration with varying shift.

 19:    References:

 21:        [1] A. Neumaier, "Residual inverse iteration for the nonlinear
 22:            eigenvalue problem", SIAM J. Numer. Anal. 22(5):914-923, 1985.
 23: */

 25: #include <slepc/private/nepimpl.h>         /*I "slepcnep.h" I*/
 26: #include <../src/nep/impls/nepdefl.h>

 28: typedef struct {
 29:   PetscInt  max_inner_it;     /* maximum number of Newton iterations */
 30:   PetscInt  lag;              /* interval to rebuild preconditioner */
 31:   PetscBool cctol;            /* constant correction tolerance */
 32:   PetscBool herm;             /* whether the Hermitian version of the scalar equation must be used */
 33:   PetscReal deftol;           /* tolerance for the deflation (threshold) */
 34:   KSP       ksp;              /* linear solver object */
 35: } NEP_RII;

 37: PetscErrorCode NEPSetUp_RII(NEP nep)
 38: {
 40:   PetscBool      istrivial;

 43:   if (nep->ncv) { PetscInfo(nep,"Setting ncv = nev, ignoring user-provided value\n"); }
 44:   nep->ncv = nep->nev;
 45:   if (nep->mpd) { PetscInfo(nep,"Setting mpd = nev, ignoring user-provided value\n"); }
 46:   nep->mpd = nep->nev;
 47:   if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
 48:   if (nep->which && nep->which!=NEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Wrong value of which");

 50:   RGIsTrivial(nep->rg,&istrivial);
 51:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support region filtering");

 53:   NEPAllocateSolution(nep,0);
 54:   NEPSetWorkVecs(nep,2);
 55:   return(0);
 56: }

 58: PetscErrorCode NEPSolve_RII(NEP nep)
 59: {
 60:   PetscErrorCode     ierr;
 61:   NEP_RII            *ctx = (NEP_RII*)nep->data;
 62:   Mat                T,Tp,H;
 63:   Vec                uu,u,r,delta,t;
 64:   PetscScalar        lambda,lambda2,sigma,a1,a2,corr,*Hp,*Ap;
 65:   PetscReal          nrm,resnorm=1.0,ktol=0.1,perr,rtol;
 66:   PetscBool          skip=PETSC_FALSE,lock=PETSC_FALSE;
 67:   PetscInt           inner_its,its=0,ldh,ldds,i,j;
 68:   NEP_EXT_OP         extop=NULL;
 69:   KSPConvergedReason kspreason;

 72:   /* get initial approximation of eigenvalue and eigenvector */
 73:   NEPGetDefaultShift(nep,&sigma);
 74:   lambda = sigma;
 75:   if (!nep->nini) {
 76:     BVSetRandomColumn(nep->V,0);
 77:     BVNormColumn(nep->V,0,NORM_2,&nrm);
 78:     BVScaleColumn(nep->V,0,1.0/nrm);
 79:   }
 80:   if (!ctx->ksp) { NEPRIIGetKSP(nep,&ctx->ksp); }
 81:   NEPDeflationInitialize(nep,nep->V,ctx->ksp,PETSC_FALSE,nep->nev,&extop);
 82:   NEPDeflationCreateVec(extop,&u);
 83:   VecDuplicate(u,&r);
 84:   VecDuplicate(u,&delta);
 85:   BVGetColumn(nep->V,0,&uu);
 86:   NEPDeflationCopyToExtendedVec(extop,uu,NULL,u,PETSC_FALSE);
 87:   BVRestoreColumn(nep->V,0,&uu);

 89:   /* prepare linear solver */
 90:   NEPDeflationSolveSetUp(extop,sigma);
 91:   KSPGetTolerances(ctx->ksp,&rtol,NULL,NULL,NULL);

 93:   VecCopy(u,r);
 94:   NEPDeflationFunctionSolve(extop,r,u);
 95:   VecNorm(u,NORM_2,&nrm);
 96:   VecScale(u,1.0/nrm);

 98:   /* Restart loop */
 99:   while (nep->reason == NEP_CONVERGED_ITERATING) {
100:     its++;

102:     /* Use Newton's method to compute nonlinear Rayleigh functional. Current eigenvalue
103:        estimate as starting value. */
104:     inner_its=0;
105:     lambda2 = lambda;
106:     do {
107:       NEPDeflationComputeFunction(extop,lambda,&T);
108:       MatMult(T,u,r);
109:       if (!ctx->herm) {
110:         NEPDeflationFunctionSolve(extop,r,delta);
111:         KSPGetConvergedReason(ctx->ksp,&kspreason);
112:         if (kspreason<0) {
113:           PetscInfo1(nep,"iter=%D, linear solve failed\n",nep->its);
114:         }
115:         t = delta;
116:       } else t = r;
117:       VecDot(t,u,&a1);
118:       NEPDeflationComputeJacobian(extop,lambda,&Tp);
119:       MatMult(Tp,u,r);
120:       if (!ctx->herm) {
121:         NEPDeflationFunctionSolve(extop,r,delta);
122:         KSPGetConvergedReason(ctx->ksp,&kspreason);
123:         if (kspreason<0) {
124:           PetscInfo1(nep,"iter=%D, linear solve failed\n",nep->its);
125:         }
126:         t = delta;
127:       } else t = r;
128:       VecDot(t,u,&a2);
129:       corr = a1/a2;
130:       lambda = lambda - corr;
131:       inner_its++;
132:     } while (PetscAbsScalar(corr)/PetscAbsScalar(lambda)>PETSC_SQRT_MACHINE_EPSILON && inner_its<ctx->max_inner_it);

134:     /* form residual,  r = T(lambda)*u */
135:     NEPDeflationComputeFunction(extop,lambda,&T);
136:     MatMult(T,u,r);

138:     /* convergence test */
139:     perr = nep->errest[nep->nconv];
140:     VecNorm(r,NORM_2,&resnorm);
141:     (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
142:     nep->eigr[nep->nconv] = lambda;
143:     if (its>1 && (nep->errest[nep->nconv]<=nep->tol || nep->errest[nep->nconv]<=ctx->deftol)) {
144:       if (nep->errest[nep->nconv]<=ctx->deftol && !extop->ref && nep->nconv) {
145:         NEPDeflationExtractEigenpair(extop,nep->nconv,u,lambda,nep->ds);
146:         NEPDeflationSetRefine(extop,PETSC_TRUE);
147:         MatMult(T,u,r);
148:         VecNorm(r,NORM_2,&resnorm);
149:         (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
150:         if (nep->errest[nep->nconv]<=nep->tol) lock = PETSC_TRUE;
151:       } else if (nep->errest[nep->nconv]<=nep->tol) lock = PETSC_TRUE;
152:     }
153:     if (lock) {
154:       NEPDeflationSetRefine(extop,PETSC_FALSE);
155:       nep->nconv = nep->nconv + 1;
156:       NEPDeflationLocking(extop,u,lambda);
157:       nep->its += its;
158:       skip = PETSC_TRUE;
159:       lock = PETSC_FALSE;
160:     }
161:     (*nep->stopping)(nep,nep->its+its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx);
162:     if (!skip || nep->reason>0) {
163:       NEPMonitor(nep,nep->its+its,nep->nconv,nep->eigr,nep->eigi,nep->errest,(nep->reason>0)?nep->nconv:nep->nconv+1);
164:     }

166:     if (nep->reason == NEP_CONVERGED_ITERATING) {
167:       if (!skip) {
168:         /* update preconditioner and set adaptive tolerance */
169:         if (ctx->lag && !(its%ctx->lag) && its>=2*ctx->lag && perr && nep->errest[nep->nconv]>.5*perr) {
170:           NEPDeflationSolveSetUp(extop,lambda2);
171:         }
172:         if (!ctx->cctol) {
173:           ktol = PetscMax(ktol/2.0,rtol);
174:           KSPSetTolerances(ctx->ksp,ktol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
175:         }

177:         /* eigenvector correction: delta = T(sigma)\r */
178:         NEPDeflationFunctionSolve(extop,r,delta);
179:         KSPGetConvergedReason(ctx->ksp,&kspreason);
180:         if (kspreason<0) {
181:           PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
182:           nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
183:           break;
184:         }

186:         /* update eigenvector: u = u - delta */
187:         VecAXPY(u,-1.0,delta);

189:         /* normalize eigenvector */
190:         VecNormalize(u,NULL);
191:       } else {
192:         its = -1;
193:         NEPDeflationSetRandomVec(extop,u);
194:         NEPDeflationSolveSetUp(extop,sigma);
195:         VecCopy(u,r);
196:         NEPDeflationFunctionSolve(extop,r,u);
197:         VecNorm(u,NORM_2,&nrm);
198:         VecScale(u,1.0/nrm);
199:         lambda = sigma;
200:         skip = PETSC_FALSE;
201:       }
202:     }
203:   }
204:   NEPDeflationGetInvariantPair(extop,NULL,&H);
205:   MatGetSize(H,NULL,&ldh);
206:   DSSetType(nep->ds,DSNHEP);
207:   DSReset(nep->ds);
208:   DSAllocate(nep->ds,PetscMax(nep->nconv,1));
209:   DSGetLeadingDimension(nep->ds,&ldds);
210:   MatDenseGetArray(H,&Hp);
211:   DSGetArray(nep->ds,DS_MAT_A,&Ap);
212:   for (j=0;j<nep->nconv;j++)
213:     for (i=0;i<nep->nconv;i++) Ap[j*ldds+i] = Hp[j*ldh+i];
214:   DSRestoreArray(nep->ds,DS_MAT_A,&Ap);
215:   MatDenseRestoreArray(H,&Hp);
216:   MatDestroy(&H);
217:   DSSetDimensions(nep->ds,nep->nconv,0,0,nep->nconv);
218:   DSSolve(nep->ds,nep->eigr,nep->eigi);
219:   NEPDeflationReset(extop);
220:   VecDestroy(&u);
221:   VecDestroy(&r);
222:   VecDestroy(&delta);
223:   return(0);
224: }

226: PetscErrorCode NEPSetFromOptions_RII(PetscOptionItems *PetscOptionsObject,NEP nep)
227: {
229:   NEP_RII        *ctx = (NEP_RII*)nep->data;
230:   PetscBool      flg;
231:   PetscInt       i;
232:   PetscReal      r;

235:   PetscOptionsHead(PetscOptionsObject,"NEP RII Options");

237:     i = 0;
238:     PetscOptionsInt("-nep_rii_max_it","Maximum number of Newton iterations for updating Rayleigh functional","NEPRIISetMaximumIterations",ctx->max_inner_it,&i,&flg);
239:     if (flg) { NEPRIISetMaximumIterations(nep,i); }

241:     PetscOptionsBool("-nep_rii_const_correction_tol","Constant correction tolerance for the linear solver","NEPRIISetConstCorrectionTol",ctx->cctol,&ctx->cctol,NULL);

243:     PetscOptionsBool("-nep_rii_hermitian","Use Hermitian version of the scalar nonlinear equation","NEPRIISetHermitian",ctx->herm,&ctx->herm,NULL);

245:     i = 0;
246:     PetscOptionsInt("-nep_rii_lag_preconditioner","Interval to rebuild preconditioner","NEPRIISetLagPreconditioner",ctx->lag,&i,&flg);
247:     if (flg) { NEPRIISetLagPreconditioner(nep,i); }

249:     r = 0.0;
250:     PetscOptionsReal("-nep_rii_deflation_threshold","Tolerance used as a threshold for including deflated eigenpairs","NEPRIISetDeflationThreshold",ctx->deftol,&r,&flg);
251:     if (flg) { NEPRIISetDeflationThreshold(nep,r); }

253:   PetscOptionsTail();

255:   if (!ctx->ksp) { NEPRIIGetKSP(nep,&ctx->ksp); }
256:   KSPSetFromOptions(ctx->ksp);
257:   return(0);
258: }

260: static PetscErrorCode NEPRIISetMaximumIterations_RII(NEP nep,PetscInt its)
261: {
262:   NEP_RII *ctx = (NEP_RII*)nep->data;

265:   if (its==PETSC_DEFAULT) ctx->max_inner_it = 10;
266:   else {
267:     if (its<=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations must be >0");
268:     ctx->max_inner_it = its;
269:   }
270:   return(0);
271: }

273: /*@
274:    NEPRIISetMaximumIterations - Sets the maximum number of inner iterations to be
275:    used in the RII solver. These are the Newton iterations related to the computation
276:    of the nonlinear Rayleigh functional.

278:    Logically Collective on nep

280:    Input Parameters:
281: +  nep - nonlinear eigenvalue solver
282: -  its - maximum inner iterations

284:    Level: advanced

286: .seealso: NEPRIIGetMaximumIterations()
287: @*/
288: PetscErrorCode NEPRIISetMaximumIterations(NEP nep,PetscInt its)
289: {

295:   PetscTryMethod(nep,"NEPRIISetMaximumIterations_C",(NEP,PetscInt),(nep,its));
296:   return(0);
297: }

299: static PetscErrorCode NEPRIIGetMaximumIterations_RII(NEP nep,PetscInt *its)
300: {
301:   NEP_RII *ctx = (NEP_RII*)nep->data;

304:   *its = ctx->max_inner_it;
305:   return(0);
306: }

308: /*@
309:    NEPRIIGetMaximumIterations - Gets the maximum number of inner iterations of RII.

311:    Not Collective

313:    Input Parameter:
314: .  nep - nonlinear eigenvalue solver

316:    Output Parameter:
317: .  its - maximum inner iterations

319:    Level: advanced

321: .seealso: NEPRIISetMaximumIterations()
322: @*/
323: PetscErrorCode NEPRIIGetMaximumIterations(NEP nep,PetscInt *its)
324: {

330:   PetscUseMethod(nep,"NEPRIIGetMaximumIterations_C",(NEP,PetscInt*),(nep,its));
331:   return(0);
332: }

334: static PetscErrorCode NEPRIISetLagPreconditioner_RII(NEP nep,PetscInt lag)
335: {
336:   NEP_RII *ctx = (NEP_RII*)nep->data;

339:   if (lag<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be non-negative");
340:   ctx->lag = lag;
341:   return(0);
342: }

344: /*@
345:    NEPRIISetLagPreconditioner - Determines when the preconditioner is rebuilt in the
346:    nonlinear solve.

348:    Logically Collective on nep

350:    Input Parameters:
351: +  nep - nonlinear eigenvalue solver
352: -  lag - 0 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is
353:           computed within the nonlinear iteration, 2 means every second time
354:           the Jacobian is built, etc.

356:    Options Database Keys:
357: .  -nep_rii_lag_preconditioner <lag>

359:    Notes:
360:    The default is 1.
361:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve.

363:    Level: intermediate

365: .seealso: NEPRIIGetLagPreconditioner()
366: @*/
367: PetscErrorCode NEPRIISetLagPreconditioner(NEP nep,PetscInt lag)
368: {

374:   PetscTryMethod(nep,"NEPRIISetLagPreconditioner_C",(NEP,PetscInt),(nep,lag));
375:   return(0);
376: }

378: static PetscErrorCode NEPRIIGetLagPreconditioner_RII(NEP nep,PetscInt *lag)
379: {
380:   NEP_RII *ctx = (NEP_RII*)nep->data;

383:   *lag = ctx->lag;
384:   return(0);
385: }

387: /*@
388:    NEPRIIGetLagPreconditioner - Indicates how often the preconditioner is rebuilt.

390:    Not Collective

392:    Input Parameter:
393: .  nep - nonlinear eigenvalue solver

395:    Output Parameter:
396: .  lag - the lag parameter

398:    Level: intermediate

400: .seealso: NEPRIISetLagPreconditioner()
401: @*/
402: PetscErrorCode NEPRIIGetLagPreconditioner(NEP nep,PetscInt *lag)
403: {

409:   PetscUseMethod(nep,"NEPRIIGetLagPreconditioner_C",(NEP,PetscInt*),(nep,lag));
410:   return(0);
411: }

413: static PetscErrorCode NEPRIISetConstCorrectionTol_RII(NEP nep,PetscBool cct)
414: {
415:   NEP_RII *ctx = (NEP_RII*)nep->data;

418:   ctx->cctol = cct;
419:   return(0);
420: }

422: /*@
423:    NEPRIISetConstCorrectionTol - Sets a flag to keep the tolerance used
424:    in the linear solver constant.

426:    Logically Collective on nep

428:    Input Parameters:
429: +  nep - nonlinear eigenvalue solver
430: -  cct - a boolean value

432:    Options Database Keys:
433: .  -nep_rii_const_correction_tol <bool> - set the boolean flag

435:    Notes:
436:    By default, an exponentially decreasing tolerance is set in the KSP used
437:    within the nonlinear iteration, so that each Newton iteration requests
438:    better accuracy than the previous one. The constant correction tolerance
439:    flag stops this behaviour.

441:    Level: intermediate

443: .seealso: NEPRIIGetConstCorrectionTol()
444: @*/
445: PetscErrorCode NEPRIISetConstCorrectionTol(NEP nep,PetscBool cct)
446: {

452:   PetscTryMethod(nep,"NEPRIISetConstCorrectionTol_C",(NEP,PetscBool),(nep,cct));
453:   return(0);
454: }

456: static PetscErrorCode NEPRIIGetConstCorrectionTol_RII(NEP nep,PetscBool *cct)
457: {
458:   NEP_RII *ctx = (NEP_RII*)nep->data;

461:   *cct = ctx->cctol;
462:   return(0);
463: }

465: /*@
466:    NEPRIIGetConstCorrectionTol - Returns the constant tolerance flag.

468:    Not Collective

470:    Input Parameter:
471: .  nep - nonlinear eigenvalue solver

473:    Output Parameter:
474: .  cct - the value of the constant tolerance flag

476:    Level: intermediate

478: .seealso: NEPRIISetConstCorrectionTol()
479: @*/
480: PetscErrorCode NEPRIIGetConstCorrectionTol(NEP nep,PetscBool *cct)
481: {

487:   PetscUseMethod(nep,"NEPRIIGetConstCorrectionTol_C",(NEP,PetscBool*),(nep,cct));
488:   return(0);
489: }

491: static PetscErrorCode NEPRIISetHermitian_RII(NEP nep,PetscBool herm)
492: {
493:   NEP_RII *ctx = (NEP_RII*)nep->data;

496:   ctx->herm = herm;
497:   return(0);
498: }

500: /*@
501:    NEPRIISetHermitian - Sets a flag to indicate if the Hermitian version of the
502:    scalar nonlinear equation must be used by the solver.

504:    Logically Collective on nep

506:    Input Parameters:
507: +  nep  - nonlinear eigenvalue solver
508: -  herm - a boolean value

510:    Options Database Keys:
511: .  -nep_rii_hermitian <bool> - set the boolean flag

513:    Notes:
514:    By default, the scalar nonlinear equation x'*inv(T(sigma))*T(z)*x=0 is solved
515:    at each step of the nonlinear iteration. When this flag is set the simpler
516:    form x'*T(z)*x=0 is used, which is supposed to be valid only for Hermitian
517:    problems.

519:    Level: intermediate

521: .seealso: NEPRIIGetHermitian()
522: @*/
523: PetscErrorCode NEPRIISetHermitian(NEP nep,PetscBool herm)
524: {

530:   PetscTryMethod(nep,"NEPRIISetHermitian_C",(NEP,PetscBool),(nep,herm));
531:   return(0);
532: }

534: static PetscErrorCode NEPRIIGetHermitian_RII(NEP nep,PetscBool *herm)
535: {
536:   NEP_RII *ctx = (NEP_RII*)nep->data;

539:   *herm = ctx->herm;
540:   return(0);
541: }

543: /*@
544:    NEPRIIGetHermitian - Returns the flag about using the Hermitian version of
545:    the scalar nonlinear equation.

547:    Not Collective

549:    Input Parameter:
550: .  nep - nonlinear eigenvalue solver

552:    Output Parameter:
553: .  herm - the value of the hermitian flag

555:    Level: intermediate

557: .seealso: NEPRIISetHermitian()
558: @*/
559: PetscErrorCode NEPRIIGetHermitian(NEP nep,PetscBool *herm)
560: {

566:   PetscUseMethod(nep,"NEPRIIGetHermitian_C",(NEP,PetscBool*),(nep,herm));
567:   return(0);
568: }

570: static PetscErrorCode NEPRIISetDeflationThreshold_RII(NEP nep,PetscReal deftol)
571: {
572:   NEP_RII *ctx = (NEP_RII*)nep->data;

575:   ctx->deftol = deftol;
576:   return(0);
577: }

579: /*@
580:    NEPRIISetDeflationThreshold - Sets the threshold value used to switch between
581:    deflated and non-deflated iteration.

583:    Logically Collective on nep

585:    Input Parameters:
586: +  nep    - nonlinear eigenvalue solver
587: -  deftol - the threshold value

589:    Options Database Keys:
590: .  -nep_rii_deflation_threshold <deftol> - set the threshold

592:    Notes:
593:    Normally, the solver iterates on the extended problem in order to deflate
594:    previously converged eigenpairs. If this threshold is set to a nonzero value,
595:    then once the residual error is below this threshold the solver will
596:    continue the iteration without deflation. The intention is to be able to
597:    improve the current eigenpair further, despite having previous eigenpairs
598:    with somewhat bad precision.

600:    Level: advanced

602: .seealso: NEPRIIGetDeflationThreshold()
603: @*/
604: PetscErrorCode NEPRIISetDeflationThreshold(NEP nep,PetscReal deftol)
605: {

611:   PetscTryMethod(nep,"NEPRIISetDeflationThreshold_C",(NEP,PetscReal),(nep,deftol));
612:   return(0);
613: }

615: static PetscErrorCode NEPRIIGetDeflationThreshold_RII(NEP nep,PetscReal *deftol)
616: {
617:   NEP_RII *ctx = (NEP_RII*)nep->data;

620:   *deftol = ctx->deftol;
621:   return(0);
622: }

624: /*@
625:    NEPRIIGetDeflationThreshold - Returns the threshold value that controls deflation.

627:    Not Collective

629:    Input Parameter:
630: .  nep - nonlinear eigenvalue solver

632:    Output Parameter:
633: .  deftol - the threshold

635:    Level: advanced

637: .seealso: NEPRIISetDeflationThreshold()
638: @*/
639: PetscErrorCode NEPRIIGetDeflationThreshold(NEP nep,PetscReal *deftol)
640: {

646:   PetscUseMethod(nep,"NEPRIIGetDeflationThreshold_C",(NEP,PetscReal*),(nep,deftol));
647:   return(0);
648: }

650: static PetscErrorCode NEPRIISetKSP_RII(NEP nep,KSP ksp)
651: {
653:   NEP_RII        *ctx = (NEP_RII*)nep->data;

656:   PetscObjectReference((PetscObject)ksp);
657:   KSPDestroy(&ctx->ksp);
658:   ctx->ksp = ksp;
659:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
660:   nep->state = NEP_STATE_INITIAL;
661:   return(0);
662: }

664: /*@
665:    NEPRIISetKSP - Associate a linear solver object (KSP) to the nonlinear
666:    eigenvalue solver.

668:    Collective on nep

670:    Input Parameters:
671: +  nep - eigenvalue solver
672: -  ksp - the linear solver object

674:    Level: advanced

676: .seealso: NEPRIIGetKSP()
677: @*/
678: PetscErrorCode NEPRIISetKSP(NEP nep,KSP ksp)
679: {

686:   PetscTryMethod(nep,"NEPRIISetKSP_C",(NEP,KSP),(nep,ksp));
687:   return(0);
688: }

690: static PetscErrorCode NEPRIIGetKSP_RII(NEP nep,KSP *ksp)
691: {
693:   NEP_RII        *ctx = (NEP_RII*)nep->data;

696:   if (!ctx->ksp) {
697:     KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp);
698:     PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)nep,1);
699:     KSPSetOptionsPrefix(ctx->ksp,((PetscObject)nep)->prefix);
700:     KSPAppendOptionsPrefix(ctx->ksp,"nep_rii_");
701:     PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
702:     PetscObjectSetOptions((PetscObject)ctx->ksp,((PetscObject)nep)->options);
703:     KSPSetErrorIfNotConverged(ctx->ksp,PETSC_TRUE);
704:     KSPSetTolerances(ctx->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
705:   }
706:   *ksp = ctx->ksp;
707:   return(0);
708: }

710: /*@
711:    NEPRIIGetKSP - Retrieve the linear solver object (KSP) associated with
712:    the nonlinear eigenvalue solver.

714:    Not Collective

716:    Input Parameter:
717: .  nep - nonlinear eigenvalue solver

719:    Output Parameter:
720: .  ksp - the linear solver object

722:    Level: advanced

724: .seealso: NEPRIISetKSP()
725: @*/
726: PetscErrorCode NEPRIIGetKSP(NEP nep,KSP *ksp)
727: {

733:   PetscUseMethod(nep,"NEPRIIGetKSP_C",(NEP,KSP*),(nep,ksp));
734:   return(0);
735: }

737: PetscErrorCode NEPView_RII(NEP nep,PetscViewer viewer)
738: {
740:   NEP_RII        *ctx = (NEP_RII*)nep->data;
741:   PetscBool      isascii;

744:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
745:   if (isascii) {
746:     PetscViewerASCIIPrintf(viewer,"  maximum number of inner iterations: %D\n",ctx->max_inner_it);
747:     if (ctx->cctol) {
748:       PetscViewerASCIIPrintf(viewer,"  using a constant tolerance for the linear solver\n");
749:     }
750:     if (ctx->herm) {
751:       PetscViewerASCIIPrintf(viewer,"  using the Hermitian version of the scalar nonlinear equation\n");
752:     }
753:     if (ctx->lag) {
754:       PetscViewerASCIIPrintf(viewer,"  updating the preconditioner every %D iterations\n",ctx->lag);
755:     }
756:     if (ctx->deftol) {
757:       PetscViewerASCIIPrintf(viewer,"  deflation threshold: %g\n",(double)ctx->deftol);
758:     }
759:     if (!ctx->ksp) { NEPRIIGetKSP(nep,&ctx->ksp); }
760:     PetscViewerASCIIPushTab(viewer);
761:     KSPView(ctx->ksp,viewer);
762:     PetscViewerASCIIPopTab(viewer);
763:   }
764:   return(0);
765: }

767: PetscErrorCode NEPReset_RII(NEP nep)
768: {
770:   NEP_RII        *ctx = (NEP_RII*)nep->data;

773:   KSPReset(ctx->ksp);
774:   return(0);
775: }

777: PetscErrorCode NEPDestroy_RII(NEP nep)
778: {
780:   NEP_RII        *ctx = (NEP_RII*)nep->data;

783:   KSPDestroy(&ctx->ksp);
784:   PetscFree(nep->data);
785:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetMaximumIterations_C",NULL);
786:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetMaximumIterations_C",NULL);
787:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetLagPreconditioner_C",NULL);
788:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetLagPreconditioner_C",NULL);
789:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetConstCorrectionTol_C",NULL);
790:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetConstCorrectionTol_C",NULL);
791:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetHermitian_C",NULL);
792:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetHermitian_C",NULL);
793:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetDeflationThreshold_C",NULL);
794:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetDeflationThreshold_C",NULL);
795:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetKSP_C",NULL);
796:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetKSP_C",NULL);
797:   return(0);
798: }

800: SLEPC_EXTERN PetscErrorCode NEPCreate_RII(NEP nep)
801: {
803:   NEP_RII        *ctx;

806:   PetscNewLog(nep,&ctx);
807:   nep->data = (void*)ctx;
808:   ctx->max_inner_it = 10;
809:   ctx->lag          = 1;
810:   ctx->cctol        = PETSC_FALSE;
811:   ctx->herm         = PETSC_FALSE;
812:   ctx->deftol       = 0.0;

814:   nep->useds = PETSC_TRUE;

816:   nep->ops->solve          = NEPSolve_RII;
817:   nep->ops->setup          = NEPSetUp_RII;
818:   nep->ops->setfromoptions = NEPSetFromOptions_RII;
819:   nep->ops->reset          = NEPReset_RII;
820:   nep->ops->destroy        = NEPDestroy_RII;
821:   nep->ops->view           = NEPView_RII;
822:   nep->ops->computevectors = NEPComputeVectors_Schur;

824:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetMaximumIterations_C",NEPRIISetMaximumIterations_RII);
825:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetMaximumIterations_C",NEPRIIGetMaximumIterations_RII);
826:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetLagPreconditioner_C",NEPRIISetLagPreconditioner_RII);
827:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetLagPreconditioner_C",NEPRIIGetLagPreconditioner_RII);
828:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetConstCorrectionTol_C",NEPRIISetConstCorrectionTol_RII);
829:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetConstCorrectionTol_C",NEPRIIGetConstCorrectionTol_RII);
830:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetHermitian_C",NEPRIISetHermitian_RII);
831:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetHermitian_C",NEPRIIGetHermitian_RII);
832:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetDeflationThreshold_C",NEPRIISetDeflationThreshold_RII);
833:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetDeflationThreshold_C",NEPRIIGetDeflationThreshold_RII);
834:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetKSP_C",NEPRIISetKSP_RII);
835:   PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetKSP_C",NEPRIIGetKSP_RII);
836:   return(0);
837: }