Actual source code: rii.c
slepc-3.7.3 2016-09-29
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-2016, 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> /*I "slepcnep.h" I*/
38: typedef struct {
39: PetscInt max_inner_it; /* maximum number of Newton iterations */
40: PetscInt lag; /* interval to rebuild preconditioner */
41: PetscBool cctol; /* constant correction tolerance */
42: KSP ksp; /* linear solver object */
43: } NEP_RII;
47: PETSC_STATIC_INLINE PetscErrorCode NEPRII_KSPSolve(NEP nep,Vec b,Vec x)
48: {
50: PetscInt lits;
51: NEP_RII *ctx = (NEP_RII*)nep->data;
54: KSPSolve(ctx->ksp,b,x);
55: KSPGetIterationNumber(ctx->ksp,&lits);
56: PetscInfo2(nep,"iter=%D, linear solve iterations=%D\n",nep->its,lits);
57: return(0);
58: }
62: PetscErrorCode NEPSetUp_RII(NEP nep)
63: {
65: PetscBool istrivial;
68: if (nep->nev>1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Requested several eigenpairs but this solver can compute only one");
69: if (nep->ncv) { PetscInfo(nep,"Setting ncv = 1, ignoring user-provided value\n"); }
70: nep->ncv = 1;
71: if (nep->mpd) { PetscInfo(nep,"Setting mpd = 1, ignoring user-provided value\n"); }
72: nep->mpd = 1;
73: if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
74: if (nep->which && nep->which!=NEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Wrong value of which");
76: RGIsTrivial(nep->rg,&istrivial);
77: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support region filtering");
79: NEPAllocateSolution(nep,0);
80: NEPSetWorkVecs(nep,2);
81: return(0);
82: }
86: PetscErrorCode NEPSolve_RII(NEP nep)
87: {
88: PetscErrorCode ierr;
89: NEP_RII *ctx = (NEP_RII*)nep->data;
90: Mat T=nep->function,Tp=nep->jacobian,Tsigma;
91: Vec u,r=nep->work[0],delta=nep->work[1];
92: PetscScalar lambda,a1,a2,corr;
93: PetscReal resnorm=1.0,ktol=0.1;
94: PetscBool hascopy;
95: PetscInt inner_its;
96: KSPConvergedReason kspreason;
99: /* get initial approximation of eigenvalue and eigenvector */
100: NEPGetDefaultShift(nep,&lambda);
101: if (!nep->nini) {
102: BVSetRandomColumn(nep->V,0);
103: }
104: BVGetColumn(nep->V,0,&u);
105: NEPComputeFunction(nep,lambda,T,T);
107: /* prepare linear solver */
108: if (!ctx->ksp) { NEPRIIGetKSP(nep,&ctx->ksp); }
109: MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
110: KSPSetOperators(ctx->ksp,Tsigma,Tsigma);
112: /* Restart loop */
113: while (nep->reason == NEP_CONVERGED_ITERATING) {
114: nep->its++;
116: /* Use Newton's method to compute nonlinear Rayleigh functional. Current eigenvalue
117: estimate as starting value. */
118: inner_its=0;
119: do {
120: NEPApplyFunction(nep,lambda,u,delta,r,T,T);
121: VecDot(r,u,&a1);
122: NEPApplyJacobian(nep,lambda,u,delta,r,Tp);
123: VecDot(r,u,&a2);
124: corr = a1/a2;
125: lambda = lambda - corr;
126: inner_its++;
127: } while (PetscAbsScalar(corr)>PETSC_SQRT_MACHINE_EPSILON && inner_its<ctx->max_inner_it);
129: /* update preconditioner and set adaptive tolerance */
130: if (ctx->lag && !(nep->its%ctx->lag) && nep->its>2*ctx->lag && resnorm<1e-2) {
131: MatHasOperation(T,MATOP_COPY,&hascopy);
132: if (hascopy) {
133: MatCopy(T,Tsigma,DIFFERENT_NONZERO_PATTERN);
134: } else {
135: MatDestroy(&Tsigma);
136: MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
137: }
138: KSPSetOperators(ctx->ksp,Tsigma,Tsigma);
139: }
140: if (!ctx->cctol) {
141: ktol = PetscMax(ktol/2.0,PETSC_MACHINE_EPSILON*10.0);
142: KSPSetTolerances(ctx->ksp,ktol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
143: }
145: /* form residual, r = T(lambda)*u */
146: NEPApplyFunction(nep,lambda,u,delta,r,T,T);
148: /* convergence test */
149: VecNorm(r,NORM_2,&resnorm);
150: (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
151: nep->eigr[nep->nconv] = lambda;
152: if (nep->errest[nep->nconv]<=nep->tol) {
153: nep->nconv = nep->nconv + 1;
154: }
155: (*nep->stopping)(nep,nep->its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx);
156: NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->eigi,nep->errest,1);
158: if (nep->reason == NEP_CONVERGED_ITERATING) {
159: /* eigenvector correction: delta = T(sigma)\r */
160: NEPRII_KSPSolve(nep,r,delta);
161: KSPGetConvergedReason(ctx->ksp,&kspreason);
162: if (kspreason<0) {
163: PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
164: nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
165: break;
166: }
168: /* update eigenvector: u = u - delta */
169: VecAXPY(u,-1.0,delta);
171: /* normalize eigenvector */
172: VecNormalize(u,NULL);
173: }
174: }
175: MatDestroy(&Tsigma);
176: BVRestoreColumn(nep->V,0,&u);
177: return(0);
178: }
182: PetscErrorCode NEPSetFromOptions_RII(PetscOptionItems *PetscOptionsObject,NEP nep)
183: {
185: NEP_RII *ctx = (NEP_RII*)nep->data;
186: PetscBool flg;
187: PetscInt i;
190: if (!ctx->ksp) { NEPRIIGetKSP(nep,&ctx->ksp); }
191: KSPSetOperators(ctx->ksp,nep->function,nep->function_pre);
192: KSPSetFromOptions(ctx->ksp);
193: PetscOptionsHead(PetscOptionsObject,"NEP RII Options");
194: PetscOptionsInt("-nep_rii_max_it","Maximum number of Newton iterations for updating Rayleigh functional","NEPRIISetMaximumIterations",ctx->max_inner_it,&ctx->max_inner_it,NULL);
195: PetscOptionsBool("-nep_rii_const_correction_tol","Constant correction tolerance for the linear solver","NEPRIISetConstCorrectionTol",ctx->cctol,&ctx->cctol,NULL);
196: i = 0;
197: PetscOptionsInt("-nep_rii_lag_preconditioner","Interval to rebuild preconditioner","NEPRIISetLagPreconditioner",ctx->lag,&i,&flg);
198: if (flg) { NEPRIISetLagPreconditioner(nep,i); }
199: PetscOptionsTail();
200: return(0);
201: }
205: static PetscErrorCode NEPRIISetMaximumIterations_RII(NEP nep,PetscInt its)
206: {
207: NEP_RII *ctx = (NEP_RII*)nep->data;
210: ctx->max_inner_it = its;
211: return(0);
212: }
216: /*@
217: NEPRIISetMaximumIterations - Sets the maximum number of inner iterations to be
218: used in the RII solver. These are the Newton iterations related to the computation
219: of the nonlinear Rayleigh functional.
221: Logically Collective on NEP
223: Input Parameters:
224: + nep - nonlinear eigenvalue solver
225: - its - maximum inner iterations
227: Level: advanced
229: .seealso: NEPRIIGetMaximumIterations()
230: @*/
231: PetscErrorCode NEPRIISetMaximumIterations(NEP nep,PetscInt its)
232: {
238: PetscTryMethod(nep,"NEPRIISetMaximumIterations_C",(NEP,PetscInt),(nep,its));
239: return(0);
240: }
244: static PetscErrorCode NEPRIIGetMaximumIterations_RII(NEP nep,PetscInt *its)
245: {
246: NEP_RII *ctx = (NEP_RII*)nep->data;
249: *its = ctx->max_inner_it;
250: return(0);
251: }
255: /*@
256: NEPRIIGetMaximumIterations - Gets the maximum number of inner iterations of RII.
258: Not Collective
260: Input Parameter:
261: . nep - nonlinear eigenvalue solver
263: Output Parameter:
264: . its - maximum inner iterations
266: Level: advanced
268: .seealso: NEPRIISetMaximumIterations()
269: @*/
270: PetscErrorCode NEPRIIGetMaximumIterations(NEP nep,PetscInt *its)
271: {
277: PetscUseMethod(nep,"NEPRIIGetMaximumIterations_C",(NEP,PetscInt*),(nep,its));
278: return(0);
279: }
283: static PetscErrorCode NEPRIISetLagPreconditioner_RII(NEP nep,PetscInt lag)
284: {
285: NEP_RII *ctx = (NEP_RII*)nep->data;
288: if (lag<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be non-negative");
289: ctx->lag = lag;
290: return(0);
291: }
295: /*@
296: NEPRIISetLagPreconditioner - Determines when the preconditioner is rebuilt in the
297: nonlinear solve.
299: Logically Collective on NEP
301: Input Parameters:
302: + nep - nonlinear eigenvalue solver
303: - lag - 0 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is
304: computed within the nonlinear iteration, 2 means every second time
305: the Jacobian is built, etc.
307: Options Database Keys:
308: . -nep_rii_lag_preconditioner <lag>
310: Notes:
311: The default is 1.
312: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve.
314: Level: intermediate
316: .seealso: NEPRIIGetLagPreconditioner()
317: @*/
318: PetscErrorCode NEPRIISetLagPreconditioner(NEP nep,PetscInt lag)
319: {
325: PetscTryMethod(nep,"NEPRIISetLagPreconditioner_C",(NEP,PetscInt),(nep,lag));
326: return(0);
327: }
331: static PetscErrorCode NEPRIIGetLagPreconditioner_RII(NEP nep,PetscInt *lag)
332: {
333: NEP_RII *ctx = (NEP_RII*)nep->data;
336: *lag = ctx->lag;
337: return(0);
338: }
342: /*@
343: NEPRIIGetLagPreconditioner - Indicates how often the preconditioner is rebuilt.
345: Not Collective
347: Input Parameter:
348: . nep - nonlinear eigenvalue solver
350: Output Parameter:
351: . lag - the lag parameter
353: Level: intermediate
355: .seealso: NEPRIISetLagPreconditioner()
356: @*/
357: PetscErrorCode NEPRIIGetLagPreconditioner(NEP nep,PetscInt *lag)
358: {
364: PetscUseMethod(nep,"NEPRIIGetLagPreconditioner_C",(NEP,PetscInt*),(nep,lag));
365: return(0);
366: }
370: static PetscErrorCode NEPRIISetConstCorrectionTol_RII(NEP nep,PetscBool cct)
371: {
372: NEP_RII *ctx = (NEP_RII*)nep->data;
375: ctx->cctol = cct;
376: return(0);
377: }
381: /*@
382: NEPRIISetConstCorrectionTol - Sets a flag to keep the tolerance used
383: in the linear solver constant.
385: Logically Collective on NEP
387: Input Parameters:
388: + nep - nonlinear eigenvalue solver
389: - cct - a boolean value
391: Options Database Keys:
392: . -nep_rii_const_correction_tol <bool> - set the boolean flag
394: Notes:
395: By default, an exponentially decreasing tolerance is set in the KSP used
396: within the nonlinear iteration, so that each Newton iteration requests
397: better accuracy than the previous one. The constant correction tolerance
398: flag stops this behaviour.
400: Level: intermediate
402: .seealso: NEPRIIGetConstCorrectionTol()
403: @*/
404: PetscErrorCode NEPRIISetConstCorrectionTol(NEP nep,PetscBool cct)
405: {
411: PetscTryMethod(nep,"NEPRIISetConstCorrectionTol_C",(NEP,PetscBool),(nep,cct));
412: return(0);
413: }
417: static PetscErrorCode NEPRIIGetConstCorrectionTol_RII(NEP nep,PetscBool *cct)
418: {
419: NEP_RII *ctx = (NEP_RII*)nep->data;
422: *cct = ctx->cctol;
423: return(0);
424: }
428: /*@
429: NEPRIIGetConstCorrectionTol - Returns the constant tolerance flag.
431: Not Collective
433: Input Parameter:
434: . nep - nonlinear eigenvalue solver
436: Output Parameter:
437: . cct - the value of the constant tolerance flag
439: Level: intermediate
441: .seealso: NEPRIISetConstCorrectionTol()
442: @*/
443: PetscErrorCode NEPRIIGetConstCorrectionTol(NEP nep,PetscBool *cct)
444: {
450: PetscUseMethod(nep,"NEPRIIGetConstCorrectionTol_C",(NEP,PetscBool*),(nep,cct));
451: return(0);
452: }
456: static PetscErrorCode NEPRIISetKSP_RII(NEP nep,KSP ksp)
457: {
459: NEP_RII *ctx = (NEP_RII*)nep->data;
462: PetscObjectReference((PetscObject)ksp);
463: KSPDestroy(&ctx->ksp);
464: ctx->ksp = ksp;
465: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
466: nep->state = NEP_STATE_INITIAL;
467: return(0);
468: }
472: /*@
473: NEPRIISetKSP - Associate a linear solver object (KSP) to the nonlinear
474: eigenvalue solver.
476: Collective on NEP
478: Input Parameters:
479: + nep - eigenvalue solver
480: - ksp - the linear solver object
482: Level: advanced
484: .seealso: NEPRIIGetKSP()
485: @*/
486: PetscErrorCode NEPRIISetKSP(NEP nep,KSP ksp)
487: {
494: PetscTryMethod(nep,"NEPRIISetKSP_C",(NEP,KSP),(nep,ksp));
495: return(0);
496: }
500: static PetscErrorCode NEPRIIGetKSP_RII(NEP nep,KSP *ksp)
501: {
503: NEP_RII *ctx = (NEP_RII*)nep->data;
506: if (!ctx->ksp) {
507: KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp);
508: KSPSetOptionsPrefix(ctx->ksp,((PetscObject)nep)->prefix);
509: KSPAppendOptionsPrefix(ctx->ksp,"nep_rii_");
510: PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)nep,1);
511: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
512: KSPSetErrorIfNotConverged(ctx->ksp,PETSC_TRUE);
513: }
514: *ksp = ctx->ksp;
515: return(0);
516: }
520: /*@
521: NEPRIIGetKSP - Retrieve the linear solver object (KSP) associated with
522: the nonlinear eigenvalue solver.
524: Not Collective
526: Input Parameter:
527: . nep - nonlinear eigenvalue solver
529: Output Parameter:
530: . ksp - the linear solver object
532: Level: advanced
534: .seealso: NEPRIISetKSP()
535: @*/
536: PetscErrorCode NEPRIIGetKSP(NEP nep,KSP *ksp)
537: {
543: PetscUseMethod(nep,"NEPRIIGetKSP_C",(NEP,KSP*),(nep,ksp));
544: return(0);
545: }
549: PetscErrorCode NEPView_RII(NEP nep,PetscViewer viewer)
550: {
552: NEP_RII *ctx = (NEP_RII*)nep->data;
553: PetscBool isascii;
556: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
557: if (isascii) {
558: if (!ctx->ksp) { NEPRIIGetKSP(nep,&ctx->ksp); }
559: PetscViewerASCIIPrintf(viewer," RII: maximum number of inner iterations: %D\n",ctx->max_inner_it);
560: if (ctx->cctol) {
561: PetscViewerASCIIPrintf(viewer," RII: using a constant tolerance for the linear solver\n");
562: }
563: if (ctx->lag) {
564: PetscViewerASCIIPrintf(viewer," RII: updating the preconditioner every %D iterations\n",ctx->lag);
565: }
566: PetscViewerASCIIPushTab(viewer);
567: KSPView(ctx->ksp,viewer);
568: PetscViewerASCIIPopTab(viewer);
569: }
570: return(0);
571: }
575: PetscErrorCode NEPDestroy_RII(NEP nep)
576: {
578: NEP_RII *ctx = (NEP_RII*)nep->data;
581: KSPDestroy(&ctx->ksp);
582: PetscFree(nep->data);
583: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetMaximumIterations_C",NULL);
584: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetMaximumIterations_C",NULL);
585: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetLagPreconditioner_C",NULL);
586: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetLagPreconditioner_C",NULL);
587: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetConstCorrectionTol_C",NULL);
588: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetConstCorrectionTol_C",NULL);
589: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetKSP_C",NULL);
590: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetKSP_C",NULL);
591: return(0);
592: }
596: PETSC_EXTERN PetscErrorCode NEPCreate_RII(NEP nep)
597: {
599: NEP_RII *ctx;
602: PetscNewLog(nep,&ctx);
603: ctx->max_inner_it = 10;
604: ctx->lag = 1;
605: ctx->cctol = PETSC_FALSE;
606: nep->data = (void*)ctx;
608: nep->ops->solve = NEPSolve_RII;
609: nep->ops->setup = NEPSetUp_RII;
610: nep->ops->setfromoptions = NEPSetFromOptions_RII;
611: nep->ops->destroy = NEPDestroy_RII;
612: nep->ops->view = NEPView_RII;
613: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetMaximumIterations_C",NEPRIISetMaximumIterations_RII);
614: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetMaximumIterations_C",NEPRIIGetMaximumIterations_RII);
615: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetLagPreconditioner_C",NEPRIISetLagPreconditioner_RII);
616: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetLagPreconditioner_C",NEPRIIGetLagPreconditioner_RII);
617: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetConstCorrectionTol_C",NEPRIISetConstCorrectionTol_RII);
618: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetConstCorrectionTol_C",NEPRIIGetConstCorrectionTol_RII);
619: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetKSP_C",NEPRIISetKSP_RII);
620: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetKSP_C",NEPRIIGetKSP_RII);
621: return(0);
622: }