Actual source code: cheby.c
2: #include <../src/ksp/ksp/impls/cheby/chebyshevimpl.h>
4: static PetscErrorCode KSPReset_Chebyshev(KSP ksp)
5: {
6: KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
8: if (cheb->kspest) KSPReset(cheb->kspest);
9: return 0;
10: }
12: /*
13: * Must be passed a KSP solver that has "converged", with KSPSetComputeEigenvalues() called before the solve
14: */
15: static PetscErrorCode KSPChebyshevComputeExtremeEigenvalues_Private(KSP kspest, PetscReal *emin, PetscReal *emax)
16: {
17: PetscInt n, neig;
18: PetscReal *re, *im, min, max;
20: KSPGetIterationNumber(kspest, &n);
21: PetscMalloc2(n, &re, n, &im);
22: KSPComputeEigenvalues(kspest, n, re, im, &neig);
23: min = PETSC_MAX_REAL;
24: max = PETSC_MIN_REAL;
25: for (n = 0; n < neig; n++) {
26: min = PetscMin(min, re[n]);
27: max = PetscMax(max, re[n]);
28: }
29: PetscFree2(re, im);
30: *emax = max;
31: *emin = min;
32: return 0;
33: }
35: static PetscErrorCode KSPSetUp_Chebyshev(KSP ksp)
36: {
37: KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
38: PetscBool isset, flg;
39: Mat Pmat, Amat;
40: PetscObjectId amatid, pmatid;
41: PetscObjectState amatstate, pmatstate;
43: KSPSetWorkVecs(ksp, 3);
44: if (cheb->emin == 0. || cheb->emax == 0.) { // User did not specify eigenvalues
45: PC pc;
46: KSPGetPC(ksp, &pc);
47: PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &flg);
48: if (!flg) { // Provided estimates are only relevant for Jacobi
49: cheb->emax_provided = 0;
50: cheb->emin_provided = 0;
51: }
52: if (!cheb->kspest) { /* We need to estimate eigenvalues */
53: KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE);
54: }
55: }
56: if (cheb->kspest) {
57: KSPGetOperators(ksp, &Amat, &Pmat);
58: MatIsSPDKnown(Pmat, &isset, &flg);
59: if (isset && flg) {
60: const char *prefix;
61: KSPGetOptionsPrefix(cheb->kspest, &prefix);
62: PetscOptionsHasName(NULL, prefix, "-ksp_type", &flg);
63: if (!flg) KSPSetType(cheb->kspest, KSPCG);
64: }
65: PetscObjectGetId((PetscObject)Amat, &amatid);
66: PetscObjectGetId((PetscObject)Pmat, &pmatid);
67: PetscObjectStateGet((PetscObject)Amat, &amatstate);
68: PetscObjectStateGet((PetscObject)Pmat, &pmatstate);
69: if (amatid != cheb->amatid || pmatid != cheb->pmatid || amatstate != cheb->amatstate || pmatstate != cheb->pmatstate) {
70: PetscReal max = 0.0, min = 0.0;
71: Vec B;
72: KSPConvergedReason reason;
73: KSPSetPC(cheb->kspest, ksp->pc);
74: if (cheb->usenoisy) {
75: B = ksp->work[1];
76: KSPSetNoisy_Private(B);
77: } else {
78: PetscBool change;
81: PCPreSolveChangeRHS(ksp->pc, &change);
82: if (change) {
83: B = ksp->work[1];
84: VecCopy(ksp->vec_rhs, B);
85: } else B = ksp->vec_rhs;
86: }
87: KSPSolve(cheb->kspest, B, ksp->work[0]);
88: KSPGetConvergedReason(cheb->kspest, &reason);
89: if (reason == KSP_DIVERGED_ITS) {
90: PetscInfo(ksp, "Eigen estimator ran for prescribed number of iterations\n");
91: } else if (reason == KSP_DIVERGED_PC_FAILED) {
92: PetscInt its;
93: PCFailedReason pcreason;
95: KSPGetIterationNumber(cheb->kspest, &its);
96: if (ksp->normtype == KSP_NORM_NONE) {
97: PetscInt sendbuf, recvbuf;
98: PCGetFailedReasonRank(ksp->pc, &pcreason);
99: sendbuf = (PetscInt)pcreason;
100: MPI_Allreduce(&sendbuf, &recvbuf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ksp));
101: PCSetFailedReason(ksp->pc, (PCFailedReason)recvbuf);
102: }
103: PCGetFailedReason(ksp->pc, &pcreason);
104: ksp->reason = KSP_DIVERGED_PC_FAILED;
105: PetscInfo(ksp, "Eigen estimator failed: %s %s at iteration %" PetscInt_FMT, KSPConvergedReasons[reason], PCFailedReasons[pcreason], its);
106: return 0;
107: } else if (reason == KSP_CONVERGED_RTOL || reason == KSP_CONVERGED_ATOL) {
108: PetscInfo(ksp, "Eigen estimator converged prematurely. Should not happen except for small or low rank problem\n");
109: } else if (reason < 0) {
110: PetscInfo(ksp, "Eigen estimator failed %s, using estimates anyway\n", KSPConvergedReasons[reason]);
111: }
113: KSPChebyshevComputeExtremeEigenvalues_Private(cheb->kspest, &min, &max);
114: KSPSetPC(cheb->kspest, NULL);
116: cheb->emin_computed = min;
117: cheb->emax_computed = max;
119: cheb->amatid = amatid;
120: cheb->pmatid = pmatid;
121: cheb->amatstate = amatstate;
122: cheb->pmatstate = pmatstate;
123: }
124: }
125: return 0;
126: }
128: static PetscErrorCode KSPChebyshevGetEigenvalues_Chebyshev(KSP ksp, PetscReal *emax, PetscReal *emin)
129: {
130: KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
132: *emax = 0;
133: *emin = 0;
134: if (cheb->emax != 0.) {
135: *emax = cheb->emax;
136: } else if (cheb->emax_computed != 0.) {
137: *emax = cheb->tform[2] * cheb->emin_computed + cheb->tform[3] * cheb->emax_computed;
138: } else if (cheb->emax_provided != 0.) {
139: *emax = cheb->tform[2] * cheb->emin_provided + cheb->tform[3] * cheb->emax_provided;
140: }
141: if (cheb->emin != 0.) {
142: *emin = cheb->emin;
143: } else if (cheb->emin_computed != 0.) {
144: *emin = cheb->tform[0] * cheb->emin_computed + cheb->tform[1] * cheb->emax_computed;
145: } else if (cheb->emin_provided != 0.) {
146: *emin = cheb->tform[0] * cheb->emin_provided + cheb->tform[1] * cheb->emax_provided;
147: }
148: return 0;
149: }
151: static PetscErrorCode KSPChebyshevSetEigenvalues_Chebyshev(KSP ksp, PetscReal emax, PetscReal emin)
152: {
153: KSP_Chebyshev *chebyshevP = (KSP_Chebyshev *)ksp->data;
157: chebyshevP->emax = emax;
158: chebyshevP->emin = emin;
160: KSPChebyshevEstEigSet(ksp, 0., 0., 0., 0.); /* Destroy any estimation setup */
161: return 0;
162: }
164: static PetscErrorCode KSPChebyshevEstEigSet_Chebyshev(KSP ksp, PetscReal a, PetscReal b, PetscReal c, PetscReal d)
165: {
166: KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
168: if (a != 0.0 || b != 0.0 || c != 0.0 || d != 0.0) {
169: if ((cheb->emin_provided == 0. || cheb->emax_provided == 0.) && !cheb->kspest) { /* should this block of code be moved to KSPSetUp_Chebyshev()? */
170: KSPCreate(PetscObjectComm((PetscObject)ksp), &cheb->kspest);
171: KSPSetErrorIfNotConverged(cheb->kspest, ksp->errorifnotconverged);
172: PetscObjectIncrementTabLevel((PetscObject)cheb->kspest, (PetscObject)ksp, 1);
173: /* use PetscObjectSet/AppendOptionsPrefix() instead of KSPSet/AppendOptionsPrefix() so that the PC prefix is not changed */
174: PetscObjectSetOptionsPrefix((PetscObject)cheb->kspest, ((PetscObject)ksp)->prefix);
175: PetscObjectAppendOptionsPrefix((PetscObject)cheb->kspest, "esteig_");
176: KSPSetSkipPCSetFromOptions(cheb->kspest, PETSC_TRUE);
178: KSPSetComputeEigenvalues(cheb->kspest, PETSC_TRUE);
180: /* We cannot turn off convergence testing because GMRES will break down if you attempt to keep iterating after a zero norm is obtained */
181: KSPSetTolerances(cheb->kspest, 1.e-12, PETSC_DEFAULT, PETSC_DEFAULT, cheb->eststeps);
182: }
183: if (a >= 0) cheb->tform[0] = a;
184: if (b >= 0) cheb->tform[1] = b;
185: if (c >= 0) cheb->tform[2] = c;
186: if (d >= 0) cheb->tform[3] = d;
187: cheb->amatid = 0;
188: cheb->pmatid = 0;
189: cheb->amatstate = -1;
190: cheb->pmatstate = -1;
191: } else {
192: KSPDestroy(&cheb->kspest);
193: }
194: return 0;
195: }
197: static PetscErrorCode KSPChebyshevEstEigSetUseNoisy_Chebyshev(KSP ksp, PetscBool use)
198: {
199: KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
201: cheb->usenoisy = use;
202: return 0;
203: }
205: /*@
206: KSPChebyshevSetEigenvalues - Sets estimates for the extreme eigenvalues
207: of the preconditioned problem.
209: Logically Collective on ksp
211: Input Parameters:
212: + ksp - the Krylov space context
213: - emax, emin - the eigenvalue estimates
215: Options Database:
216: . -ksp_chebyshev_eigenvalues emin,emax - extreme eigenvalues
218: Note:
219: Call KSPChebyshevEstEigSet() or use the option -ksp_chebyshev_esteig a,b,c,d to have the KSP
220: estimate the eigenvalues and use these estimated values automatically.
222: When KSPCHEBYSHEV is used as a smoother, one often wants to target a portion of the spectrum rather than the entire
223: spectrum. This function takes the range of target eigenvalues for Chebyshev, which will often slightly over-estimate
224: the largest eigenvalue of the actual operator (for safety) and greatly overestimate the smallest eigenvalue to
225: improve the smoothing properties of Chebyshev iteration on the higher frequencies in the spectrum.
227: Level: intermediate
229: .seealso: `KSPChebyshevEstEigSet()`
230: @*/
231: PetscErrorCode KSPChebyshevSetEigenvalues(KSP ksp, PetscReal emax, PetscReal emin)
232: {
236: PetscTryMethod(ksp, "KSPChebyshevSetEigenvalues_C", (KSP, PetscReal, PetscReal), (ksp, emax, emin));
237: return 0;
238: }
240: /*@
241: KSPChebyshevEstEigSet - Automatically estimate the eigenvalues to use for Chebyshev
243: Logically Collective on ksp
245: Input Parameters:
246: + ksp - the Krylov space context
247: . a - multiple of min eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
248: . b - multiple of max eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
249: . c - multiple of min eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)
250: - d - multiple of max eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)
252: Options Database:
253: . -ksp_chebyshev_esteig a,b,c,d - estimate eigenvalues using a Krylov method, then use this transform for Chebyshev eigenvalue bounds
255: Notes:
256: The Chebyshev bounds are set using
257: .vb
258: minbound = a*minest + b*maxest
259: maxbound = c*minest + d*maxest
260: .ve
261: The default configuration targets the upper part of the spectrum for use as a multigrid smoother, so only the maximum eigenvalue estimate is used.
262: The minimum eigenvalue estimate obtained by Krylov iteration is typically not accurate until the method has converged.
264: If 0.0 is passed for all transform arguments (a,b,c,d), eigenvalue estimation is disabled.
266: The default transform is (0,0.1; 0,1.1) which targets the "upper" part of the spectrum, as desirable for use with multigrid.
268: The eigenvalues are estimated using the Lanczo (KSPCG) or Arnoldi (KSPGMRES) process using a noisy right hand side vector.
270: Level: intermediate
272: @*/
273: PetscErrorCode KSPChebyshevEstEigSet(KSP ksp, PetscReal a, PetscReal b, PetscReal c, PetscReal d)
274: {
280: PetscTryMethod(ksp, "KSPChebyshevEstEigSet_C", (KSP, PetscReal, PetscReal, PetscReal, PetscReal), (ksp, a, b, c, d));
281: return 0;
282: }
284: /*@
285: KSPChebyshevEstEigSetUseNoisy - use a noisy right hand side in order to do the estimate instead of the given right hand side
287: Logically Collective
289: Input Parameters:
290: + ksp - linear solver context
291: - use - PETSC_TRUE to use noisy
293: Options Database:
294: . -ksp_chebyshev_esteig_noisy <true,false> - Use noisy right hand side for estimate
296: Notes:
297: This alledgely works better for multigrid smoothers
299: Level: intermediate
301: .seealso: `KSPChebyshevEstEigSet()`
302: @*/
303: PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP ksp, PetscBool use)
304: {
306: PetscTryMethod(ksp, "KSPChebyshevEstEigSetUseNoisy_C", (KSP, PetscBool), (ksp, use));
307: return 0;
308: }
310: /*@
311: KSPChebyshevEstEigGetKSP - Get the Krylov method context used to estimate eigenvalues for the Chebyshev method. If
312: a Krylov method is not being used for this purpose, NULL is returned. The reference count of the returned KSP is
313: not incremented: it should not be destroyed by the user.
315: Input Parameters:
316: . ksp - the Krylov space context
318: Output Parameters:
319: . kspest - the eigenvalue estimation Krylov space context
321: Level: intermediate
323: .seealso: `KSPChebyshevEstEigSet()`
324: @*/
325: PetscErrorCode KSPChebyshevEstEigGetKSP(KSP ksp, KSP *kspest)
326: {
329: *kspest = NULL;
330: PetscTryMethod(ksp, "KSPChebyshevEstEigGetKSP_C", (KSP, KSP *), (ksp, kspest));
331: return 0;
332: }
334: static PetscErrorCode KSPChebyshevEstEigGetKSP_Chebyshev(KSP ksp, KSP *kspest)
335: {
336: KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
338: *kspest = cheb->kspest;
339: return 0;
340: }
342: static PetscErrorCode KSPSetFromOptions_Chebyshev(KSP ksp, PetscOptionItems *PetscOptionsObject)
343: {
344: KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
345: PetscInt neigarg = 2, nestarg = 4;
346: PetscReal eminmax[2] = {0., 0.};
347: PetscReal tform[4] = {PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE};
348: PetscBool flgeig, flgest;
350: PetscOptionsHeadBegin(PetscOptionsObject, "KSP Chebyshev Options");
351: PetscOptionsInt("-ksp_chebyshev_esteig_steps", "Number of est steps in Chebyshev", "", cheb->eststeps, &cheb->eststeps, NULL);
352: PetscOptionsRealArray("-ksp_chebyshev_eigenvalues", "extreme eigenvalues", "KSPChebyshevSetEigenvalues", eminmax, &neigarg, &flgeig);
353: if (flgeig) {
355: KSPChebyshevSetEigenvalues(ksp, eminmax[1], eminmax[0]);
356: }
357: PetscOptionsRealArray("-ksp_chebyshev_esteig", "estimate eigenvalues using a Krylov method, then use this transform for Chebyshev eigenvalue bounds", "KSPChebyshevEstEigSet", tform, &nestarg, &flgest);
358: if (flgest) {
359: switch (nestarg) {
360: case 0:
361: KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE);
362: break;
363: case 2: /* Base everything on the max eigenvalues */
364: KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, tform[0], PETSC_DECIDE, tform[1]);
365: break;
366: case 4: /* Use the full 2x2 linear transformation */
367: KSPChebyshevEstEigSet(ksp, tform[0], tform[1], tform[2], tform[3]);
368: break;
369: default:
370: SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "Must specify either 0, 2, or 4 parameters for eigenvalue estimation");
371: }
372: }
374: /* We need to estimate eigenvalues; need to set this here so that KSPSetFromOptions() is called on the estimator */
375: if ((cheb->emin == 0. || cheb->emax == 0.) && !cheb->kspest) KSPChebyshevEstEigSet(ksp, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE);
377: if (cheb->kspest) {
378: PetscOptionsBool("-ksp_chebyshev_esteig_noisy", "Use noisy right hand side for estimate", "KSPChebyshevEstEigSetUseNoisy", cheb->usenoisy, &cheb->usenoisy, NULL);
379: KSPSetFromOptions(cheb->kspest);
380: }
381: PetscOptionsHeadEnd();
382: return 0;
383: }
385: static PetscErrorCode KSPSolve_Chebyshev(KSP ksp)
386: {
387: PetscInt k, kp1, km1, ktmp, i;
388: PetscScalar alpha, omegaprod, mu, omega, Gamma, c[3], scale;
389: PetscReal rnorm = 0.0, emax, emin;
390: Vec sol_orig, b, p[3], r;
391: Mat Amat, Pmat;
392: PetscBool diagonalscale;
394: PCGetDiagonalScale(ksp->pc, &diagonalscale);
397: PCGetOperators(ksp->pc, &Amat, &Pmat);
398: PetscObjectSAWsTakeAccess((PetscObject)ksp);
399: ksp->its = 0;
400: PetscObjectSAWsGrantAccess((PetscObject)ksp);
401: /* These three point to the three active solutions, we
402: rotate these three at each solution update */
403: km1 = 0;
404: k = 1;
405: kp1 = 2;
406: sol_orig = ksp->vec_sol; /* ksp->vec_sol will be assigned to rotating vector p[k], thus save its address */
407: b = ksp->vec_rhs;
408: p[km1] = sol_orig;
409: p[k] = ksp->work[0];
410: p[kp1] = ksp->work[1];
411: r = ksp->work[2];
413: KSPChebyshevGetEigenvalues_Chebyshev(ksp, &emax, &emin);
414: /* use scale*B as our preconditioner */
415: scale = 2.0 / (emax + emin);
417: /* -alpha <= scale*lambda(B^{-1}A) <= alpha */
418: alpha = 1.0 - scale * emin;
419: Gamma = 1.0;
420: mu = 1.0 / alpha;
421: omegaprod = 2.0 / alpha;
423: c[km1] = 1.0;
424: c[k] = mu;
426: if (!ksp->guess_zero) {
427: KSP_MatMult(ksp, Amat, sol_orig, r); /* r = b - A*p[km1] */
428: VecAYPX(r, -1.0, b);
429: } else {
430: VecCopy(b, r);
431: }
433: /* calculate residual norm if requested, we have done one iteration */
434: if (ksp->normtype) {
435: switch (ksp->normtype) {
436: case KSP_NORM_PRECONDITIONED:
437: KSP_PCApply(ksp, r, p[k]); /* p[k] = B^{-1}r */
438: VecNorm(p[k], NORM_2, &rnorm);
439: break;
440: case KSP_NORM_UNPRECONDITIONED:
441: case KSP_NORM_NATURAL:
442: VecNorm(r, NORM_2, &rnorm);
443: break;
444: default:
445: SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
446: }
447: PetscObjectSAWsTakeAccess((PetscObject)ksp);
448: ksp->rnorm = rnorm;
449: PetscObjectSAWsGrantAccess((PetscObject)ksp);
450: KSPLogResidualHistory(ksp, rnorm);
451: KSPLogErrorHistory(ksp);
452: KSPMonitor(ksp, 0, rnorm);
453: (*ksp->converged)(ksp, 0, rnorm, &ksp->reason, ksp->cnvP);
454: } else ksp->reason = KSP_CONVERGED_ITERATING;
455: if (ksp->reason || ksp->max_it == 0) {
456: if (ksp->max_it == 0) ksp->reason = KSP_DIVERGED_ITS; /* This for a V(0,x) cycle */
457: return 0;
458: }
459: if (ksp->normtype != KSP_NORM_PRECONDITIONED) { KSP_PCApply(ksp, r, p[k]); /* p[k] = B^{-1}r */ }
460: VecAYPX(p[k], scale, p[km1]); /* p[k] = scale B^{-1}r + p[km1] */
461: PetscObjectSAWsTakeAccess((PetscObject)ksp);
462: ksp->its = 1;
463: PetscObjectSAWsGrantAccess((PetscObject)ksp);
465: for (i = 1; i < ksp->max_it; i++) {
466: PetscObjectSAWsTakeAccess((PetscObject)ksp);
467: ksp->its++;
468: PetscObjectSAWsGrantAccess((PetscObject)ksp);
470: KSP_MatMult(ksp, Amat, p[k], r); /* r = b - Ap[k] */
471: VecAYPX(r, -1.0, b);
472: /* calculate residual norm if requested */
473: if (ksp->normtype) {
474: switch (ksp->normtype) {
475: case KSP_NORM_PRECONDITIONED:
476: KSP_PCApply(ksp, r, p[kp1]); /* p[kp1] = B^{-1}r */
477: VecNorm(p[kp1], NORM_2, &rnorm);
478: break;
479: case KSP_NORM_UNPRECONDITIONED:
480: case KSP_NORM_NATURAL:
481: VecNorm(r, NORM_2, &rnorm);
482: break;
483: default:
484: rnorm = 0.0;
485: break;
486: }
487: KSPCheckNorm(ksp, rnorm);
488: PetscObjectSAWsTakeAccess((PetscObject)ksp);
489: ksp->rnorm = rnorm;
490: PetscObjectSAWsGrantAccess((PetscObject)ksp);
491: KSPLogResidualHistory(ksp, rnorm);
492: KSPMonitor(ksp, i, rnorm);
493: (*ksp->converged)(ksp, i, rnorm, &ksp->reason, ksp->cnvP);
494: if (ksp->reason) break;
495: if (ksp->normtype != KSP_NORM_PRECONDITIONED) { KSP_PCApply(ksp, r, p[kp1]); /* p[kp1] = B^{-1}r */ }
496: } else {
497: KSP_PCApply(ksp, r, p[kp1]); /* p[kp1] = B^{-1}r */
498: }
499: ksp->vec_sol = p[k];
500: KSPLogErrorHistory(ksp);
502: c[kp1] = 2.0 * mu * c[k] - c[km1];
503: omega = omegaprod * c[k] / c[kp1];
505: /* y^{k+1} = omega(y^{k} - y^{k-1} + Gamma*r^{k}) + y^{k-1} */
506: VecAXPBYPCZ(p[kp1], 1.0 - omega, omega, omega * Gamma * scale, p[km1], p[k]);
508: ktmp = km1;
509: km1 = k;
510: k = kp1;
511: kp1 = ktmp;
512: }
513: if (!ksp->reason) {
514: if (ksp->normtype) {
515: KSP_MatMult(ksp, Amat, p[k], r); /* r = b - Ap[k] */
516: VecAYPX(r, -1.0, b);
517: switch (ksp->normtype) {
518: case KSP_NORM_PRECONDITIONED:
519: KSP_PCApply(ksp, r, p[kp1]); /* p[kp1] = B^{-1}r */
520: VecNorm(p[kp1], NORM_2, &rnorm);
521: break;
522: case KSP_NORM_UNPRECONDITIONED:
523: case KSP_NORM_NATURAL:
524: VecNorm(r, NORM_2, &rnorm);
525: break;
526: default:
527: rnorm = 0.0;
528: break;
529: }
530: KSPCheckNorm(ksp, rnorm);
531: PetscObjectSAWsTakeAccess((PetscObject)ksp);
532: ksp->rnorm = rnorm;
533: PetscObjectSAWsGrantAccess((PetscObject)ksp);
534: KSPLogResidualHistory(ksp, rnorm);
535: KSPMonitor(ksp, i, rnorm);
536: }
537: if (ksp->its >= ksp->max_it) {
538: if (ksp->normtype != KSP_NORM_NONE) {
539: (*ksp->converged)(ksp, i, rnorm, &ksp->reason, ksp->cnvP);
540: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
541: } else ksp->reason = KSP_CONVERGED_ITS;
542: }
543: }
545: /* make sure solution is in vector x */
546: ksp->vec_sol = sol_orig;
547: if (k) VecCopy(p[k], sol_orig);
548: if (ksp->reason == KSP_CONVERGED_ITS) KSPLogErrorHistory(ksp);
549: return 0;
550: }
552: static PetscErrorCode KSPView_Chebyshev(KSP ksp, PetscViewer viewer)
553: {
554: KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
555: PetscBool iascii;
557: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
558: if (iascii) {
559: PetscReal emax, emin;
560: KSPChebyshevGetEigenvalues_Chebyshev(ksp, &emax, &emin);
561: PetscViewerASCIIPrintf(viewer, " eigenvalue targets used: min %g, max %g\n", (double)emin, (double)emax);
562: if (cheb->kspest) {
563: PetscViewerASCIIPrintf(viewer, " eigenvalues estimated via %s: min %g, max %g\n", ((PetscObject)(cheb->kspest))->type_name, (double)cheb->emin_computed, (double)cheb->emax_computed);
564: PetscViewerASCIIPrintf(viewer, " eigenvalues estimated using %s with transform: [%g %g; %g %g]\n", ((PetscObject)cheb->kspest)->type_name, (double)cheb->tform[0], (double)cheb->tform[1], (double)cheb->tform[2], (double)cheb->tform[3]);
565: PetscViewerASCIIPushTab(viewer);
566: KSPView(cheb->kspest, viewer);
567: PetscViewerASCIIPopTab(viewer);
568: if (cheb->usenoisy) PetscViewerASCIIPrintf(viewer, " estimating eigenvalues using noisy right hand side\n");
569: } else if (cheb->emax_provided != 0.) {
570: PetscCall(PetscViewerASCIIPrintf(viewer, " eigenvalues provided (min %g, max %g) with transform: [%g %g; %g %g]\n", (double)cheb->emin_provided, (double)cheb->emax_provided, (double)cheb->tform[0], (double)cheb->tform[1], (double)cheb->tform[2],
571: (double)cheb->tform[3]));
572: }
573: }
574: return 0;
575: }
577: static PetscErrorCode KSPDestroy_Chebyshev(KSP ksp)
578: {
579: KSP_Chebyshev *cheb = (KSP_Chebyshev *)ksp->data;
581: KSPDestroy(&cheb->kspest);
582: PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevSetEigenvalues_C", NULL);
583: PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSet_C", NULL);
584: PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSetUseNoisy_C", NULL);
585: PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigGetKSP_C", NULL);
586: KSPDestroyDefault(ksp);
587: return 0;
588: }
590: /*MC
591: KSPCHEBYSHEV - The preconditioned Chebyshev iterative method
593: Options Database Keys:
594: + -ksp_chebyshev_eigenvalues <emin,emax> - set approximations to the smallest and largest eigenvalues
595: of the preconditioned operator. If these are accurate you will get much faster convergence.
596: . -ksp_chebyshev_esteig <a,b,c,d> - estimate eigenvalues using a Krylov method, then use this
597: transform for Chebyshev eigenvalue bounds (KSPChebyshevEstEigSet())
598: . -ksp_chebyshev_esteig_steps - number of estimation steps
599: - -ksp_chebyshev_esteig_noisy - use noisy number generator to create right hand side for eigenvalue estimator
601: Level: beginner
603: Notes:
604: The Chebyshev method requires both the matrix and preconditioner to
605: be symmetric positive (semi) definite.
606: Only support for left preconditioning.
608: Chebyshev is configured as a smoother by default, targetting the "upper" part of the spectrum.
609: The user should call KSPChebyshevSetEigenvalues() if they have eigenvalue estimates.
611: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`,
612: `KSPChebyshevSetEigenvalues()`, `KSPChebyshevEstEigSet()`, `KSPChebyshevEstEigSetUseNoisy()`
613: `KSPRICHARDSON`, `KSPCG`, `PCMG`
615: M*/
617: PETSC_EXTERN PetscErrorCode KSPCreate_Chebyshev(KSP ksp)
618: {
619: KSP_Chebyshev *chebyshevP;
621: PetscNew(&chebyshevP);
623: ksp->data = (void *)chebyshevP;
624: KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
625: KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2);
626: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1);
627: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);
629: chebyshevP->emin = 0.;
630: chebyshevP->emax = 0.;
632: chebyshevP->tform[0] = 0.0;
633: chebyshevP->tform[1] = 0.1;
634: chebyshevP->tform[2] = 0;
635: chebyshevP->tform[3] = 1.1;
636: chebyshevP->eststeps = 10;
637: chebyshevP->usenoisy = PETSC_TRUE;
638: ksp->setupnewmatrix = PETSC_TRUE;
640: ksp->ops->setup = KSPSetUp_Chebyshev;
641: ksp->ops->solve = KSPSolve_Chebyshev;
642: ksp->ops->destroy = KSPDestroy_Chebyshev;
643: ksp->ops->buildsolution = KSPBuildSolutionDefault;
644: ksp->ops->buildresidual = KSPBuildResidualDefault;
645: ksp->ops->setfromoptions = KSPSetFromOptions_Chebyshev;
646: ksp->ops->view = KSPView_Chebyshev;
647: ksp->ops->reset = KSPReset_Chebyshev;
649: PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevSetEigenvalues_C", KSPChebyshevSetEigenvalues_Chebyshev);
650: PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSet_C", KSPChebyshevEstEigSet_Chebyshev);
651: PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigSetUseNoisy_C", KSPChebyshevEstEigSetUseNoisy_Chebyshev);
652: PetscObjectComposeFunction((PetscObject)ksp, "KSPChebyshevEstEigGetKSP_C", KSPChebyshevEstEigGetKSP_Chebyshev);
653: return 0;
654: }