Actual source code: rii.c
slepc-3.14.0 2020-09-30
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>
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: {
42: if (nep->ncv!=PETSC_DEFAULT) { PetscInfo(nep,"Setting ncv = nev, ignoring user-provided value\n"); }
43: nep->ncv = nep->nev;
44: if (nep->mpd!=PETSC_DEFAULT) { PetscInfo(nep,"Setting mpd = nev, ignoring user-provided value\n"); }
45: nep->mpd = nep->nev;
46: if (nep->max_it==PETSC_DEFAULT) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
47: if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
48: if (nep->which!=NEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver supports only target magnitude eigenvalues");
49: NEPCheckUnsupported(nep,NEP_FEATURE_REGION | NEP_FEATURE_TWOSIDED);
50: NEPAllocateSolution(nep,0);
51: NEPSetWorkVecs(nep,2);
52: return(0);
53: }
55: PetscErrorCode NEPSolve_RII(NEP nep)
56: {
57: PetscErrorCode ierr;
58: NEP_RII *ctx = (NEP_RII*)nep->data;
59: Mat T,Tp,H;
60: Vec uu,u,r,delta,t;
61: PetscScalar lambda,lambda2,sigma,a1,a2,corr,*Hp,*Ap;
62: PetscReal nrm,resnorm=1.0,ktol=0.1,perr,rtol;
63: PetscBool skip=PETSC_FALSE,lock=PETSC_FALSE;
64: PetscInt inner_its,its=0,ldh,ldds,i,j;
65: NEP_EXT_OP extop=NULL;
66: KSPConvergedReason kspreason;
69: /* get initial approximation of eigenvalue and eigenvector */
70: NEPGetDefaultShift(nep,&sigma);
71: lambda = sigma;
72: if (!nep->nini) {
73: BVSetRandomColumn(nep->V,0);
74: BVNormColumn(nep->V,0,NORM_2,&nrm);
75: BVScaleColumn(nep->V,0,1.0/nrm);
76: }
77: if (!ctx->ksp) { NEPRIIGetKSP(nep,&ctx->ksp); }
78: NEPDeflationInitialize(nep,nep->V,ctx->ksp,PETSC_FALSE,nep->nev,&extop);
79: NEPDeflationCreateVec(extop,&u);
80: VecDuplicate(u,&r);
81: VecDuplicate(u,&delta);
82: BVGetColumn(nep->V,0,&uu);
83: NEPDeflationCopyToExtendedVec(extop,uu,NULL,u,PETSC_FALSE);
84: BVRestoreColumn(nep->V,0,&uu);
86: /* prepare linear solver */
87: NEPDeflationSolveSetUp(extop,sigma);
88: KSPGetTolerances(ctx->ksp,&rtol,NULL,NULL,NULL);
90: VecCopy(u,r);
91: NEPDeflationFunctionSolve(extop,r,u);
92: VecNorm(u,NORM_2,&nrm);
93: VecScale(u,1.0/nrm);
95: /* Restart loop */
96: while (nep->reason == NEP_CONVERGED_ITERATING) {
97: its++;
99: /* Use Newton's method to compute nonlinear Rayleigh functional. Current eigenvalue
100: estimate as starting value. */
101: inner_its=0;
102: lambda2 = lambda;
103: do {
104: NEPDeflationComputeFunction(extop,lambda,&T);
105: MatMult(T,u,r);
106: if (!ctx->herm) {
107: NEPDeflationFunctionSolve(extop,r,delta);
108: KSPGetConvergedReason(ctx->ksp,&kspreason);
109: if (kspreason<0) {
110: PetscInfo1(nep,"iter=%D, linear solve failed\n",nep->its);
111: }
112: t = delta;
113: } else t = r;
114: VecDot(t,u,&a1);
115: NEPDeflationComputeJacobian(extop,lambda,&Tp);
116: MatMult(Tp,u,r);
117: if (!ctx->herm) {
118: NEPDeflationFunctionSolve(extop,r,delta);
119: KSPGetConvergedReason(ctx->ksp,&kspreason);
120: if (kspreason<0) {
121: PetscInfo1(nep,"iter=%D, linear solve failed\n",nep->its);
122: }
123: t = delta;
124: } else t = r;
125: VecDot(t,u,&a2);
126: corr = a1/a2;
127: lambda = lambda - corr;
128: inner_its++;
129: } while (PetscAbsScalar(corr)/PetscAbsScalar(lambda)>PETSC_SQRT_MACHINE_EPSILON && inner_its<ctx->max_inner_it);
131: /* form residual, r = T(lambda)*u */
132: NEPDeflationComputeFunction(extop,lambda,&T);
133: MatMult(T,u,r);
135: /* convergence test */
136: perr = nep->errest[nep->nconv];
137: VecNorm(r,NORM_2,&resnorm);
138: (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
139: nep->eigr[nep->nconv] = lambda;
140: if (its>1 && (nep->errest[nep->nconv]<=nep->tol || nep->errest[nep->nconv]<=ctx->deftol)) {
141: if (nep->errest[nep->nconv]<=ctx->deftol && !extop->ref && nep->nconv) {
142: NEPDeflationExtractEigenpair(extop,nep->nconv,u,lambda,nep->ds);
143: NEPDeflationSetRefine(extop,PETSC_TRUE);
144: MatMult(T,u,r);
145: VecNorm(r,NORM_2,&resnorm);
146: (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
147: if (nep->errest[nep->nconv]<=nep->tol) lock = PETSC_TRUE;
148: } else if (nep->errest[nep->nconv]<=nep->tol) lock = PETSC_TRUE;
149: }
150: if (lock) {
151: NEPDeflationSetRefine(extop,PETSC_FALSE);
152: nep->nconv = nep->nconv + 1;
153: NEPDeflationLocking(extop,u,lambda);
154: nep->its += its;
155: skip = PETSC_TRUE;
156: lock = PETSC_FALSE;
157: }
158: (*nep->stopping)(nep,nep->its+its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx);
159: if (!skip || nep->reason>0) {
160: NEPMonitor(nep,nep->its+its,nep->nconv,nep->eigr,nep->eigi,nep->errest,(nep->reason>0)?nep->nconv:nep->nconv+1);
161: }
163: if (nep->reason == NEP_CONVERGED_ITERATING) {
164: if (!skip) {
165: /* update preconditioner and set adaptive tolerance */
166: if (ctx->lag && !(its%ctx->lag) && its>=2*ctx->lag && perr && nep->errest[nep->nconv]>.5*perr) {
167: NEPDeflationSolveSetUp(extop,lambda2);
168: }
169: if (!ctx->cctol) {
170: ktol = PetscMax(ktol/2.0,rtol);
171: KSPSetTolerances(ctx->ksp,ktol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
172: }
174: /* eigenvector correction: delta = T(sigma)\r */
175: NEPDeflationFunctionSolve(extop,r,delta);
176: KSPGetConvergedReason(ctx->ksp,&kspreason);
177: if (kspreason<0) {
178: PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
179: nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
180: break;
181: }
183: /* update eigenvector: u = u - delta */
184: VecAXPY(u,-1.0,delta);
186: /* normalize eigenvector */
187: VecNormalize(u,NULL);
188: } else {
189: its = -1;
190: NEPDeflationSetRandomVec(extop,u);
191: NEPDeflationSolveSetUp(extop,sigma);
192: VecCopy(u,r);
193: NEPDeflationFunctionSolve(extop,r,u);
194: VecNorm(u,NORM_2,&nrm);
195: VecScale(u,1.0/nrm);
196: lambda = sigma;
197: skip = PETSC_FALSE;
198: }
199: }
200: }
201: NEPDeflationGetInvariantPair(extop,NULL,&H);
202: MatGetSize(H,NULL,&ldh);
203: DSSetType(nep->ds,DSNHEP);
204: DSReset(nep->ds);
205: DSAllocate(nep->ds,PetscMax(nep->nconv,1));
206: DSGetLeadingDimension(nep->ds,&ldds);
207: MatDenseGetArray(H,&Hp);
208: DSGetArray(nep->ds,DS_MAT_A,&Ap);
209: for (j=0;j<nep->nconv;j++)
210: for (i=0;i<nep->nconv;i++) Ap[j*ldds+i] = Hp[j*ldh+i];
211: DSRestoreArray(nep->ds,DS_MAT_A,&Ap);
212: MatDenseRestoreArray(H,&Hp);
213: MatDestroy(&H);
214: DSSetDimensions(nep->ds,nep->nconv,0,0,nep->nconv);
215: DSSolve(nep->ds,nep->eigr,nep->eigi);
216: NEPDeflationReset(extop);
217: VecDestroy(&u);
218: VecDestroy(&r);
219: VecDestroy(&delta);
220: return(0);
221: }
223: PetscErrorCode NEPSetFromOptions_RII(PetscOptionItems *PetscOptionsObject,NEP nep)
224: {
226: NEP_RII *ctx = (NEP_RII*)nep->data;
227: PetscBool flg;
228: PetscInt i;
229: PetscReal r;
232: PetscOptionsHead(PetscOptionsObject,"NEP RII Options");
234: i = 0;
235: PetscOptionsInt("-nep_rii_max_it","Maximum number of Newton iterations for updating Rayleigh functional","NEPRIISetMaximumIterations",ctx->max_inner_it,&i,&flg);
236: if (flg) { NEPRIISetMaximumIterations(nep,i); }
238: PetscOptionsBool("-nep_rii_const_correction_tol","Constant correction tolerance for the linear solver","NEPRIISetConstCorrectionTol",ctx->cctol,&ctx->cctol,NULL);
240: PetscOptionsBool("-nep_rii_hermitian","Use Hermitian version of the scalar nonlinear equation","NEPRIISetHermitian",ctx->herm,&ctx->herm,NULL);
242: i = 0;
243: PetscOptionsInt("-nep_rii_lag_preconditioner","Interval to rebuild preconditioner","NEPRIISetLagPreconditioner",ctx->lag,&i,&flg);
244: if (flg) { NEPRIISetLagPreconditioner(nep,i); }
246: r = 0.0;
247: PetscOptionsReal("-nep_rii_deflation_threshold","Tolerance used as a threshold for including deflated eigenpairs","NEPRIISetDeflationThreshold",ctx->deftol,&r,&flg);
248: if (flg) { NEPRIISetDeflationThreshold(nep,r); }
250: PetscOptionsTail();
252: if (!ctx->ksp) { NEPRIIGetKSP(nep,&ctx->ksp); }
253: KSPSetFromOptions(ctx->ksp);
254: return(0);
255: }
257: static PetscErrorCode NEPRIISetMaximumIterations_RII(NEP nep,PetscInt its)
258: {
259: NEP_RII *ctx = (NEP_RII*)nep->data;
262: if (its==PETSC_DEFAULT) ctx->max_inner_it = 10;
263: else {
264: if (its<=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations must be >0");
265: ctx->max_inner_it = its;
266: }
267: return(0);
268: }
270: /*@
271: NEPRIISetMaximumIterations - Sets the maximum number of inner iterations to be
272: used in the RII solver. These are the Newton iterations related to the computation
273: of the nonlinear Rayleigh functional.
275: Logically Collective on nep
277: Input Parameters:
278: + nep - nonlinear eigenvalue solver
279: - its - maximum inner iterations
281: Level: advanced
283: .seealso: NEPRIIGetMaximumIterations()
284: @*/
285: PetscErrorCode NEPRIISetMaximumIterations(NEP nep,PetscInt its)
286: {
292: PetscTryMethod(nep,"NEPRIISetMaximumIterations_C",(NEP,PetscInt),(nep,its));
293: return(0);
294: }
296: static PetscErrorCode NEPRIIGetMaximumIterations_RII(NEP nep,PetscInt *its)
297: {
298: NEP_RII *ctx = (NEP_RII*)nep->data;
301: *its = ctx->max_inner_it;
302: return(0);
303: }
305: /*@
306: NEPRIIGetMaximumIterations - Gets the maximum number of inner iterations of RII.
308: Not Collective
310: Input Parameter:
311: . nep - nonlinear eigenvalue solver
313: Output Parameter:
314: . its - maximum inner iterations
316: Level: advanced
318: .seealso: NEPRIISetMaximumIterations()
319: @*/
320: PetscErrorCode NEPRIIGetMaximumIterations(NEP nep,PetscInt *its)
321: {
327: PetscUseMethod(nep,"NEPRIIGetMaximumIterations_C",(NEP,PetscInt*),(nep,its));
328: return(0);
329: }
331: static PetscErrorCode NEPRIISetLagPreconditioner_RII(NEP nep,PetscInt lag)
332: {
333: NEP_RII *ctx = (NEP_RII*)nep->data;
336: if (lag<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be non-negative");
337: ctx->lag = lag;
338: return(0);
339: }
341: /*@
342: NEPRIISetLagPreconditioner - Determines when the preconditioner is rebuilt in the
343: nonlinear solve.
345: Logically Collective on nep
347: Input Parameters:
348: + nep - nonlinear eigenvalue solver
349: - lag - 0 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is
350: computed within the nonlinear iteration, 2 means every second time
351: the Jacobian is built, etc.
353: Options Database Keys:
354: . -nep_rii_lag_preconditioner <lag>
356: Notes:
357: The default is 1.
358: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve.
360: Level: intermediate
362: .seealso: NEPRIIGetLagPreconditioner()
363: @*/
364: PetscErrorCode NEPRIISetLagPreconditioner(NEP nep,PetscInt lag)
365: {
371: PetscTryMethod(nep,"NEPRIISetLagPreconditioner_C",(NEP,PetscInt),(nep,lag));
372: return(0);
373: }
375: static PetscErrorCode NEPRIIGetLagPreconditioner_RII(NEP nep,PetscInt *lag)
376: {
377: NEP_RII *ctx = (NEP_RII*)nep->data;
380: *lag = ctx->lag;
381: return(0);
382: }
384: /*@
385: NEPRIIGetLagPreconditioner - Indicates how often the preconditioner is rebuilt.
387: Not Collective
389: Input Parameter:
390: . nep - nonlinear eigenvalue solver
392: Output Parameter:
393: . lag - the lag parameter
395: Level: intermediate
397: .seealso: NEPRIISetLagPreconditioner()
398: @*/
399: PetscErrorCode NEPRIIGetLagPreconditioner(NEP nep,PetscInt *lag)
400: {
406: PetscUseMethod(nep,"NEPRIIGetLagPreconditioner_C",(NEP,PetscInt*),(nep,lag));
407: return(0);
408: }
410: static PetscErrorCode NEPRIISetConstCorrectionTol_RII(NEP nep,PetscBool cct)
411: {
412: NEP_RII *ctx = (NEP_RII*)nep->data;
415: ctx->cctol = cct;
416: return(0);
417: }
419: /*@
420: NEPRIISetConstCorrectionTol - Sets a flag to keep the tolerance used
421: in the linear solver constant.
423: Logically Collective on nep
425: Input Parameters:
426: + nep - nonlinear eigenvalue solver
427: - cct - a boolean value
429: Options Database Keys:
430: . -nep_rii_const_correction_tol <bool> - set the boolean flag
432: Notes:
433: By default, an exponentially decreasing tolerance is set in the KSP used
434: within the nonlinear iteration, so that each Newton iteration requests
435: better accuracy than the previous one. The constant correction tolerance
436: flag stops this behaviour.
438: Level: intermediate
440: .seealso: NEPRIIGetConstCorrectionTol()
441: @*/
442: PetscErrorCode NEPRIISetConstCorrectionTol(NEP nep,PetscBool cct)
443: {
449: PetscTryMethod(nep,"NEPRIISetConstCorrectionTol_C",(NEP,PetscBool),(nep,cct));
450: return(0);
451: }
453: static PetscErrorCode NEPRIIGetConstCorrectionTol_RII(NEP nep,PetscBool *cct)
454: {
455: NEP_RII *ctx = (NEP_RII*)nep->data;
458: *cct = ctx->cctol;
459: return(0);
460: }
462: /*@
463: NEPRIIGetConstCorrectionTol - Returns the constant tolerance flag.
465: Not Collective
467: Input Parameter:
468: . nep - nonlinear eigenvalue solver
470: Output Parameter:
471: . cct - the value of the constant tolerance flag
473: Level: intermediate
475: .seealso: NEPRIISetConstCorrectionTol()
476: @*/
477: PetscErrorCode NEPRIIGetConstCorrectionTol(NEP nep,PetscBool *cct)
478: {
484: PetscUseMethod(nep,"NEPRIIGetConstCorrectionTol_C",(NEP,PetscBool*),(nep,cct));
485: return(0);
486: }
488: static PetscErrorCode NEPRIISetHermitian_RII(NEP nep,PetscBool herm)
489: {
490: NEP_RII *ctx = (NEP_RII*)nep->data;
493: ctx->herm = herm;
494: return(0);
495: }
497: /*@
498: NEPRIISetHermitian - Sets a flag to indicate if the Hermitian version of the
499: scalar nonlinear equation must be used by the solver.
501: Logically Collective on nep
503: Input Parameters:
504: + nep - nonlinear eigenvalue solver
505: - herm - a boolean value
507: Options Database Keys:
508: . -nep_rii_hermitian <bool> - set the boolean flag
510: Notes:
511: By default, the scalar nonlinear equation x'*inv(T(sigma))*T(z)*x=0 is solved
512: at each step of the nonlinear iteration. When this flag is set the simpler
513: form x'*T(z)*x=0 is used, which is supposed to be valid only for Hermitian
514: problems.
516: Level: intermediate
518: .seealso: NEPRIIGetHermitian()
519: @*/
520: PetscErrorCode NEPRIISetHermitian(NEP nep,PetscBool herm)
521: {
527: PetscTryMethod(nep,"NEPRIISetHermitian_C",(NEP,PetscBool),(nep,herm));
528: return(0);
529: }
531: static PetscErrorCode NEPRIIGetHermitian_RII(NEP nep,PetscBool *herm)
532: {
533: NEP_RII *ctx = (NEP_RII*)nep->data;
536: *herm = ctx->herm;
537: return(0);
538: }
540: /*@
541: NEPRIIGetHermitian - Returns the flag about using the Hermitian version of
542: the scalar nonlinear equation.
544: Not Collective
546: Input Parameter:
547: . nep - nonlinear eigenvalue solver
549: Output Parameter:
550: . herm - the value of the hermitian flag
552: Level: intermediate
554: .seealso: NEPRIISetHermitian()
555: @*/
556: PetscErrorCode NEPRIIGetHermitian(NEP nep,PetscBool *herm)
557: {
563: PetscUseMethod(nep,"NEPRIIGetHermitian_C",(NEP,PetscBool*),(nep,herm));
564: return(0);
565: }
567: static PetscErrorCode NEPRIISetDeflationThreshold_RII(NEP nep,PetscReal deftol)
568: {
569: NEP_RII *ctx = (NEP_RII*)nep->data;
572: ctx->deftol = deftol;
573: return(0);
574: }
576: /*@
577: NEPRIISetDeflationThreshold - Sets the threshold value used to switch between
578: deflated and non-deflated iteration.
580: Logically Collective on nep
582: Input Parameters:
583: + nep - nonlinear eigenvalue solver
584: - deftol - the threshold value
586: Options Database Keys:
587: . -nep_rii_deflation_threshold <deftol> - set the threshold
589: Notes:
590: Normally, the solver iterates on the extended problem in order to deflate
591: previously converged eigenpairs. If this threshold is set to a nonzero value,
592: then once the residual error is below this threshold the solver will
593: continue the iteration without deflation. The intention is to be able to
594: improve the current eigenpair further, despite having previous eigenpairs
595: with somewhat bad precision.
597: Level: advanced
599: .seealso: NEPRIIGetDeflationThreshold()
600: @*/
601: PetscErrorCode NEPRIISetDeflationThreshold(NEP nep,PetscReal deftol)
602: {
608: PetscTryMethod(nep,"NEPRIISetDeflationThreshold_C",(NEP,PetscReal),(nep,deftol));
609: return(0);
610: }
612: static PetscErrorCode NEPRIIGetDeflationThreshold_RII(NEP nep,PetscReal *deftol)
613: {
614: NEP_RII *ctx = (NEP_RII*)nep->data;
617: *deftol = ctx->deftol;
618: return(0);
619: }
621: /*@
622: NEPRIIGetDeflationThreshold - Returns the threshold value that controls deflation.
624: Not Collective
626: Input Parameter:
627: . nep - nonlinear eigenvalue solver
629: Output Parameter:
630: . deftol - the threshold
632: Level: advanced
634: .seealso: NEPRIISetDeflationThreshold()
635: @*/
636: PetscErrorCode NEPRIIGetDeflationThreshold(NEP nep,PetscReal *deftol)
637: {
643: PetscUseMethod(nep,"NEPRIIGetDeflationThreshold_C",(NEP,PetscReal*),(nep,deftol));
644: return(0);
645: }
647: static PetscErrorCode NEPRIISetKSP_RII(NEP nep,KSP ksp)
648: {
650: NEP_RII *ctx = (NEP_RII*)nep->data;
653: PetscObjectReference((PetscObject)ksp);
654: KSPDestroy(&ctx->ksp);
655: ctx->ksp = ksp;
656: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
657: nep->state = NEP_STATE_INITIAL;
658: return(0);
659: }
661: /*@
662: NEPRIISetKSP - Associate a linear solver object (KSP) to the nonlinear
663: eigenvalue solver.
665: Collective on nep
667: Input Parameters:
668: + nep - eigenvalue solver
669: - ksp - the linear solver object
671: Level: advanced
673: .seealso: NEPRIIGetKSP()
674: @*/
675: PetscErrorCode NEPRIISetKSP(NEP nep,KSP ksp)
676: {
683: PetscTryMethod(nep,"NEPRIISetKSP_C",(NEP,KSP),(nep,ksp));
684: return(0);
685: }
687: static PetscErrorCode NEPRIIGetKSP_RII(NEP nep,KSP *ksp)
688: {
690: NEP_RII *ctx = (NEP_RII*)nep->data;
693: if (!ctx->ksp) {
694: KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp);
695: PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)nep,1);
696: KSPSetOptionsPrefix(ctx->ksp,((PetscObject)nep)->prefix);
697: KSPAppendOptionsPrefix(ctx->ksp,"nep_rii_");
698: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
699: PetscObjectSetOptions((PetscObject)ctx->ksp,((PetscObject)nep)->options);
700: KSPSetErrorIfNotConverged(ctx->ksp,PETSC_TRUE);
701: KSPSetTolerances(ctx->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
702: }
703: *ksp = ctx->ksp;
704: return(0);
705: }
707: /*@
708: NEPRIIGetKSP - Retrieve the linear solver object (KSP) associated with
709: the nonlinear eigenvalue solver.
711: Not Collective
713: Input Parameter:
714: . nep - nonlinear eigenvalue solver
716: Output Parameter:
717: . ksp - the linear solver object
719: Level: advanced
721: .seealso: NEPRIISetKSP()
722: @*/
723: PetscErrorCode NEPRIIGetKSP(NEP nep,KSP *ksp)
724: {
730: PetscUseMethod(nep,"NEPRIIGetKSP_C",(NEP,KSP*),(nep,ksp));
731: return(0);
732: }
734: PetscErrorCode NEPView_RII(NEP nep,PetscViewer viewer)
735: {
737: NEP_RII *ctx = (NEP_RII*)nep->data;
738: PetscBool isascii;
741: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
742: if (isascii) {
743: PetscViewerASCIIPrintf(viewer," maximum number of inner iterations: %D\n",ctx->max_inner_it);
744: if (ctx->cctol) {
745: PetscViewerASCIIPrintf(viewer," using a constant tolerance for the linear solver\n");
746: }
747: if (ctx->herm) {
748: PetscViewerASCIIPrintf(viewer," using the Hermitian version of the scalar nonlinear equation\n");
749: }
750: if (ctx->lag) {
751: PetscViewerASCIIPrintf(viewer," updating the preconditioner every %D iterations\n",ctx->lag);
752: }
753: if (ctx->deftol) {
754: PetscViewerASCIIPrintf(viewer," deflation threshold: %g\n",(double)ctx->deftol);
755: }
756: if (!ctx->ksp) { NEPRIIGetKSP(nep,&ctx->ksp); }
757: PetscViewerASCIIPushTab(viewer);
758: KSPView(ctx->ksp,viewer);
759: PetscViewerASCIIPopTab(viewer);
760: }
761: return(0);
762: }
764: PetscErrorCode NEPReset_RII(NEP nep)
765: {
767: NEP_RII *ctx = (NEP_RII*)nep->data;
770: KSPReset(ctx->ksp);
771: return(0);
772: }
774: PetscErrorCode NEPDestroy_RII(NEP nep)
775: {
777: NEP_RII *ctx = (NEP_RII*)nep->data;
780: KSPDestroy(&ctx->ksp);
781: PetscFree(nep->data);
782: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetMaximumIterations_C",NULL);
783: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetMaximumIterations_C",NULL);
784: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetLagPreconditioner_C",NULL);
785: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetLagPreconditioner_C",NULL);
786: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetConstCorrectionTol_C",NULL);
787: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetConstCorrectionTol_C",NULL);
788: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetHermitian_C",NULL);
789: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetHermitian_C",NULL);
790: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetDeflationThreshold_C",NULL);
791: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetDeflationThreshold_C",NULL);
792: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetKSP_C",NULL);
793: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetKSP_C",NULL);
794: return(0);
795: }
797: SLEPC_EXTERN PetscErrorCode NEPCreate_RII(NEP nep)
798: {
800: NEP_RII *ctx;
803: PetscNewLog(nep,&ctx);
804: nep->data = (void*)ctx;
805: ctx->max_inner_it = 10;
806: ctx->lag = 1;
807: ctx->cctol = PETSC_FALSE;
808: ctx->herm = PETSC_FALSE;
809: ctx->deftol = 0.0;
811: nep->useds = PETSC_TRUE;
813: nep->ops->solve = NEPSolve_RII;
814: nep->ops->setup = NEPSetUp_RII;
815: nep->ops->setfromoptions = NEPSetFromOptions_RII;
816: nep->ops->reset = NEPReset_RII;
817: nep->ops->destroy = NEPDestroy_RII;
818: nep->ops->view = NEPView_RII;
819: nep->ops->computevectors = NEPComputeVectors_Schur;
821: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetMaximumIterations_C",NEPRIISetMaximumIterations_RII);
822: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetMaximumIterations_C",NEPRIIGetMaximumIterations_RII);
823: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetLagPreconditioner_C",NEPRIISetLagPreconditioner_RII);
824: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetLagPreconditioner_C",NEPRIIGetLagPreconditioner_RII);
825: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetConstCorrectionTol_C",NEPRIISetConstCorrectionTol_RII);
826: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetConstCorrectionTol_C",NEPRIIGetConstCorrectionTol_RII);
827: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetHermitian_C",NEPRIISetHermitian_RII);
828: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetHermitian_C",NEPRIIGetHermitian_RII);
829: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetDeflationThreshold_C",NEPRIISetDeflationThreshold_RII);
830: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetDeflationThreshold_C",NEPRIIGetDeflationThreshold_RII);
831: PetscObjectComposeFunction((PetscObject)nep,"NEPRIISetKSP_C",NEPRIISetKSP_RII);
832: PetscObjectComposeFunction((PetscObject)nep,"NEPRIIGetKSP_C",NEPRIIGetKSP_RII);
833: return(0);
834: }