Actual source code: cheby.c
petsc-3.14.2 2020-12-03
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;
10: if (cheb->kspest) {
11: KSPReset(cheb->kspest);
12: }
13: return(0);
14: }
16: PETSC_STATIC_INLINE PetscScalar chebyhash(PetscInt xx)
17: {
18: unsigned int x = xx;
19: x = ((x >> 16) ^ x) * 0x45d9f3b;
20: x = ((x >> 16) ^ x) * 0x45d9f3b;
21: x = ((x >> 16) ^ x);
22: return (PetscScalar)((PetscInt64)x-2147483648)*5.e-10; /* center around zero, scaled about -1. to 1.*/
23: }
25: /*
26: * Must be passed a KSP solver that has "converged", with KSPSetComputeEigenvalues() called before the solve
27: */
28: static PetscErrorCode KSPChebyshevComputeExtremeEigenvalues_Private(KSP kspest,PetscReal *emin,PetscReal *emax)
29: {
31: PetscInt n,neig;
32: PetscReal *re,*im,min,max;
35: KSPGetIterationNumber(kspest,&n);
36: PetscMalloc2(n,&re,n,&im);
37: KSPComputeEigenvalues(kspest,n,re,im,&neig);
38: min = PETSC_MAX_REAL;
39: max = PETSC_MIN_REAL;
40: for (n=0; n<neig; n++) {
41: min = PetscMin(min,re[n]);
42: max = PetscMax(max,re[n]);
43: }
44: PetscFree2(re,im);
45: *emax = max;
46: *emin = min;
47: return(0);
48: }
50: static PetscErrorCode KSPSetUp_Chebyshev(KSP ksp)
51: {
52: KSP_Chebyshev *cheb = (KSP_Chebyshev*)ksp->data;
53: PetscErrorCode ierr;
54: PetscBool flg;
55: Mat Pmat,Amat;
56: PetscObjectId amatid, pmatid;
57: PetscObjectState amatstate, pmatstate;
59: KSPSetWorkVecs(ksp,3);
60: if ((cheb->emin == 0. || cheb->emax == 0.) && !cheb->kspest) { /* We need to estimate eigenvalues */
61: KSPChebyshevEstEigSet(ksp,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);
62: }
63: if (cheb->kspest) {
64: KSPGetOperators(ksp,&Amat,&Pmat);
65: MatGetOption(Pmat, MAT_SPD, &flg);
66: if (flg) {
67: const char *prefix;
68: KSPGetOptionsPrefix(cheb->kspest,&prefix);
69: PetscOptionsHasName(NULL,prefix,"-ksp_type",&flg);
70: if (!flg) {
71: KSPSetType(cheb->kspest, KSPCG);
72: }
73: }
74: PetscObjectGetId((PetscObject)Amat,&amatid);
75: PetscObjectGetId((PetscObject)Pmat,&pmatid);
76: PetscObjectStateGet((PetscObject)Amat,&amatstate);
77: PetscObjectStateGet((PetscObject)Pmat,&pmatstate);
78: if (amatid != cheb->amatid || pmatid != cheb->pmatid || amatstate != cheb->amatstate || pmatstate != cheb->pmatstate) {
79: PetscReal max=0.0,min=0.0;
80: Vec B;
81: KSPConvergedReason reason;
82: KSPSetPC(cheb->kspest,ksp->pc);
83: if (cheb->usenoisy) {
84: PetscInt n,i,istart;
85: PetscScalar *xx;
87: B = ksp->work[1];
88: VecGetOwnershipRange(B,&istart,NULL);
89: VecGetLocalSize(B,&n);
90: VecGetArrayWrite(B,&xx);
91: for (i=0; i<n; i++) xx[i] = chebyhash(i+istart);
92: VecRestoreArrayWrite(B,&xx);
93: } else {
94: PetscBool change;
96: PCPreSolveChangeRHS(ksp->pc,&change);
97: if (change) {
98: B = ksp->work[1];
99: VecCopy(ksp->vec_rhs,B);
100: } else B = ksp->vec_rhs;
101: }
102: KSPSolve(cheb->kspest,B,ksp->work[0]);
103: KSPGetConvergedReason(cheb->kspest,&reason);
104: if (reason == KSP_DIVERGED_ITS) {
105: PetscInfo(ksp,"Eigen estimator ran for prescribed number of iterations\n");
106: } else if (reason == KSP_DIVERGED_PC_FAILED) {
107: PetscInt its;
108: PCFailedReason pcreason;
110: KSPGetIterationNumber(cheb->kspest,&its);
111: if (ksp->normtype == KSP_NORM_NONE) {
112: PetscInt sendbuf,recvbuf;
113: PCGetFailedReasonRank(ksp->pc,&pcreason);
114: sendbuf = (PetscInt)pcreason;
115: MPI_Allreduce(&sendbuf,&recvbuf,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ksp));
116: PCSetFailedReason(ksp->pc,(PCFailedReason) recvbuf);
117: }
118: PCGetFailedReason(ksp->pc,&pcreason);
119: ksp->reason = KSP_DIVERGED_PC_FAILED;
120: PetscInfo3(ksp,"Eigen estimator failed: %s %s at iteration %D",KSPConvergedReasons[reason],PCFailedReasons[pcreason],its);
121: return(0);
122: } else if (reason == KSP_CONVERGED_RTOL || reason == KSP_CONVERGED_ATOL) {
123: PetscInfo(ksp,"Eigen estimator converged prematurely. Should not happen except for small or low rank problem\n");
124: } else if (reason < 0) {
125: PetscInfo1(ksp,"Eigen estimator failed %s, using estimates anyway\n",KSPConvergedReasons[reason]);
126: }
128: KSPChebyshevComputeExtremeEigenvalues_Private(cheb->kspest,&min,&max);
129: KSPSetPC(cheb->kspest,NULL);
131: cheb->emin_computed = min;
132: cheb->emax_computed = max;
133: cheb->emin = cheb->tform[0]*min + cheb->tform[1]*max;
134: cheb->emax = cheb->tform[2]*min + cheb->tform[3]*max;
136: cheb->amatid = amatid;
137: cheb->pmatid = pmatid;
138: cheb->amatstate = amatstate;
139: cheb->pmatstate = pmatstate;
140: }
141: }
142: return(0);
143: }
145: static PetscErrorCode KSPChebyshevSetEigenvalues_Chebyshev(KSP ksp,PetscReal emax,PetscReal emin)
146: {
147: KSP_Chebyshev *chebyshevP = (KSP_Chebyshev*)ksp->data;
151: if (emax <= emin) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %g",(double)emax,(double)emin);
152: if (emax*emin <= 0.0) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %g min %g",(double)emax,(double)emin);
153: chebyshevP->emax = emax;
154: chebyshevP->emin = emin;
156: KSPChebyshevEstEigSet(ksp,0.,0.,0.,0.); /* Destroy any estimation setup */
157: return(0);
158: }
160: static PetscErrorCode KSPChebyshevEstEigSet_Chebyshev(KSP ksp,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
161: {
162: KSP_Chebyshev *cheb = (KSP_Chebyshev*)ksp->data;
166: if (a != 0.0 || b != 0.0 || c != 0.0 || d != 0.0) {
167: if (!cheb->kspest) { /* should this block of code be moved to KSPSetUp_Chebyshev()? */
168: KSPCreate(PetscObjectComm((PetscObject)ksp),&cheb->kspest);
169: KSPSetErrorIfNotConverged(cheb->kspest,ksp->errorifnotconverged);
170: PetscObjectIncrementTabLevel((PetscObject)cheb->kspest,(PetscObject)ksp,1);
171: /* use PetscObjectSet/AppendOptionsPrefix() instead of KSPSet/AppendOptionsPrefix() so that the PC prefix is not changed */
172: PetscObjectSetOptionsPrefix((PetscObject)cheb->kspest,((PetscObject)ksp)->prefix);
173: PetscObjectAppendOptionsPrefix((PetscObject)cheb->kspest,"esteig_");
174: KSPSetSkipPCSetFromOptions(cheb->kspest,PETSC_TRUE);
176: KSPSetComputeEigenvalues(cheb->kspest,PETSC_TRUE);
178: /* We cannot turn off convergence testing because GMRES will break down if you attempt to keep iterating after a zero norm is obtained */
179: KSPSetTolerances(cheb->kspest,1.e-12,PETSC_DEFAULT,PETSC_DEFAULT,cheb->eststeps);
180: }
181: if (a >= 0) cheb->tform[0] = a;
182: if (b >= 0) cheb->tform[1] = b;
183: if (c >= 0) cheb->tform[2] = c;
184: if (d >= 0) cheb->tform[3] = d;
185: cheb->amatid = 0;
186: cheb->pmatid = 0;
187: cheb->amatstate = -1;
188: cheb->pmatstate = -1;
189: } else {
190: KSPDestroy(&cheb->kspest);
191: }
192: return(0);
193: }
195: static PetscErrorCode KSPChebyshevEstEigSetUseNoisy_Chebyshev(KSP ksp,PetscBool use)
196: {
197: KSP_Chebyshev *cheb = (KSP_Chebyshev*)ksp->data;
200: cheb->usenoisy = use;
201: return(0);
202: }
204: /*@
205: KSPChebyshevSetEigenvalues - Sets estimates for the extreme eigenvalues
206: of the preconditioned problem.
208: Logically Collective on ksp
210: Input Parameters:
211: + ksp - the Krylov space context
212: - emax, emin - the eigenvalue estimates
214: Options Database:
215: . -ksp_chebyshev_eigenvalues emin,emax
217: Note: Call KSPChebyshevEstEigSet() or use the option -ksp_chebyshev_esteig a,b,c,d to have the KSP
218: estimate the eigenvalues and use these estimated values automatically
220: Level: intermediate
222: @*/
223: PetscErrorCode KSPChebyshevSetEigenvalues(KSP ksp,PetscReal emax,PetscReal emin)
224: {
231: PetscTryMethod(ksp,"KSPChebyshevSetEigenvalues_C",(KSP,PetscReal,PetscReal),(ksp,emax,emin));
232: return(0);
233: }
235: /*@
236: KSPChebyshevEstEigSet - Automatically estimate the eigenvalues to use for Chebyshev
238: Logically Collective on ksp
240: Input Parameters:
241: + ksp - the Krylov space context
242: . a - multiple of min eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
243: . b - multiple of max eigenvalue estimate to use for min Chebyshev bound (or PETSC_DECIDE)
244: . c - multiple of min eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)
245: - d - multiple of max eigenvalue estimate to use for max Chebyshev bound (or PETSC_DECIDE)
247: Options Database:
248: . -ksp_chebyshev_esteig a,b,c,d
250: Notes:
251: The Chebyshev bounds are set using
252: .vb
253: minbound = a*minest + b*maxest
254: maxbound = c*minest + d*maxest
255: .ve
256: The default configuration targets the upper part of the spectrum for use as a multigrid smoother, so only the maximum eigenvalue estimate is used.
257: The minimum eigenvalue estimate obtained by Krylov iteration is typically not accurate until the method has converged.
259: If 0.0 is passed for all transform arguments (a,b,c,d), eigenvalue estimation is disabled.
261: The default transform is (0,0.1; 0,1.1) which targets the "upper" part of the spectrum, as desirable for use with multigrid.
263: The eigenvalues are estimated using the Lanczo (KSPCG) or Arnoldi (KSPGMRES) process using a noisy right hand side vector.
265: Level: intermediate
267: @*/
268: PetscErrorCode KSPChebyshevEstEigSet(KSP ksp,PetscReal a,PetscReal b,PetscReal c,PetscReal d)
269: {
278: PetscTryMethod(ksp,"KSPChebyshevEstEigSet_C",(KSP,PetscReal,PetscReal,PetscReal,PetscReal),(ksp,a,b,c,d));
279: return(0);
280: }
282: /*@
283: KSPChebyshevEstEigSetUseNoisy - use a noisy right hand side in order to do the estimate instead of the given right hand side
285: Logically Collective
287: Input Arguments:
288: + ksp - linear solver context
289: - use - PETSC_TRUE to use noisy
291: Options Database:
292: . -ksp_chebyshev_esteig_noisy <true,false>
294: Notes:
295: This alledgely works better for multigrid smoothers
297: Level: intermediate
299: .seealso: KSPChebyshevEstEigSet()
300: @*/
301: PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP ksp,PetscBool use)
302: {
307: PetscTryMethod(ksp,"KSPChebyshevEstEigSetUseNoisy_C",(KSP,PetscBool),(ksp,use));
308: return(0);
309: }
311: /*@
312: KSPChebyshevEstEigGetKSP - Get the Krylov method context used to estimate eigenvalues for the Chebyshev method. If
313: a Krylov method is not being used for this purpose, NULL is returned. The reference count of the returned KSP is
314: not incremented: it should not be destroyed by the user.
316: Input Parameters:
317: . ksp - the Krylov space context
319: Output Parameters:
320: . kspest the eigenvalue estimation Krylov space context
322: Level: intermediate
324: .seealso: KSPChebyshevEstEigSet()
325: @*/
326: PetscErrorCode KSPChebyshevEstEigGetKSP(KSP ksp, KSP *kspest)
327: {
333: *kspest = NULL;
334: PetscTryMethod(ksp,"KSPChebyshevEstEigGetKSP_C",(KSP,KSP*),(ksp,kspest));
335: return(0);
336: }
338: static PetscErrorCode KSPChebyshevEstEigGetKSP_Chebyshev(KSP ksp, KSP *kspest)
339: {
340: KSP_Chebyshev *cheb = (KSP_Chebyshev*)ksp->data;
343: *kspest = cheb->kspest;
344: return(0);
345: }
347: static PetscErrorCode KSPSetFromOptions_Chebyshev(PetscOptionItems *PetscOptionsObject,KSP ksp)
348: {
349: KSP_Chebyshev *cheb = (KSP_Chebyshev*)ksp->data;
351: PetscInt neigarg = 2, nestarg = 4;
352: PetscReal eminmax[2] = {0., 0.};
353: PetscReal tform[4] = {PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE};
354: PetscBool flgeig, flgest;
357: PetscOptionsHead(PetscOptionsObject,"KSP Chebyshev Options");
358: PetscOptionsInt("-ksp_chebyshev_esteig_steps","Number of est steps in Chebyshev","",cheb->eststeps,&cheb->eststeps,NULL);
359: PetscOptionsRealArray("-ksp_chebyshev_eigenvalues","extreme eigenvalues","KSPChebyshevSetEigenvalues",eminmax,&neigarg,&flgeig);
360: if (flgeig) {
361: if (neigarg != 2) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"-ksp_chebyshev_eigenvalues: must specify 2 parameters, min and max eigenvalues");
362: KSPChebyshevSetEigenvalues(ksp, eminmax[1], eminmax[0]);
363: }
364: PetscOptionsRealArray("-ksp_chebyshev_esteig","estimate eigenvalues using a Krylov method, then use this transform for Chebyshev eigenvalue bounds","KSPChebyshevEstEigSet",tform,&nestarg,&flgest);
365: if (flgest) {
366: switch (nestarg) {
367: case 0:
368: KSPChebyshevEstEigSet(ksp,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);
369: break;
370: case 2: /* Base everything on the max eigenvalues */
371: KSPChebyshevEstEigSet(ksp,PETSC_DECIDE,tform[0],PETSC_DECIDE,tform[1]);
372: break;
373: case 4: /* Use the full 2x2 linear transformation */
374: KSPChebyshevEstEigSet(ksp,tform[0],tform[1],tform[2],tform[3]);
375: break;
376: default: SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"Must specify either 0, 2, or 4 parameters for eigenvalue estimation");
377: }
378: }
380: /* We need to estimate eigenvalues; need to set this here so that KSPSetFromOptions() is called on the estimator */
381: if ((cheb->emin == 0. || cheb->emax == 0.) && !cheb->kspest) {
382: KSPChebyshevEstEigSet(ksp,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);
383: }
385: if (cheb->kspest) {
386: PetscOptionsBool("-ksp_chebyshev_esteig_noisy","Use noisy right hand side for estimate","KSPChebyshevEstEigSetUseNoisy",cheb->usenoisy,&cheb->usenoisy,NULL);
387: KSPSetFromOptions(cheb->kspest);
388: }
389: PetscOptionsTail();
390: return(0);
391: }
393: static PetscErrorCode KSPSolve_Chebyshev(KSP ksp)
394: {
395: KSP_Chebyshev *cheb = (KSP_Chebyshev*)ksp->data;
397: PetscInt k,kp1,km1,ktmp,i;
398: PetscScalar alpha,omegaprod,mu,omega,Gamma,c[3],scale;
399: PetscReal rnorm = 0.0;
400: Vec sol_orig,b,p[3],r;
401: Mat Amat,Pmat;
402: PetscBool diagonalscale;
405: PCGetDiagonalScale(ksp->pc,&diagonalscale);
406: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
408: PCGetOperators(ksp->pc,&Amat,&Pmat);
409: PetscObjectSAWsTakeAccess((PetscObject)ksp);
410: ksp->its = 0;
411: PetscObjectSAWsGrantAccess((PetscObject)ksp);
412: /* These three point to the three active solutions, we
413: rotate these three at each solution update */
414: km1 = 0; k = 1; kp1 = 2;
415: sol_orig = ksp->vec_sol; /* ksp->vec_sol will be asigned to rotating vector p[k], thus save its address */
416: b = ksp->vec_rhs;
417: p[km1] = sol_orig;
418: p[k] = ksp->work[0];
419: p[kp1] = ksp->work[1];
420: r = ksp->work[2];
422: /* use scale*B as our preconditioner */
423: scale = 2.0/(cheb->emax + cheb->emin);
425: /* -alpha <= scale*lambda(B^{-1}A) <= alpha */
426: alpha = 1.0 - scale*(cheb->emin);
427: Gamma = 1.0;
428: mu = 1.0/alpha;
429: omegaprod = 2.0/alpha;
431: c[km1] = 1.0;
432: c[k] = mu;
434: if (!ksp->guess_zero) {
435: KSP_MatMult(ksp,Amat,sol_orig,r); /* r = b - A*p[km1] */
436: VecAYPX(r,-1.0,b);
437: } else {
438: VecCopy(b,r);
439: }
441: /* calculate residual norm if requested, we have done one iteration */
442: if (ksp->normtype) {
443: switch (ksp->normtype) {
444: case KSP_NORM_PRECONDITIONED:
445: KSP_PCApply(ksp,r,p[k]); /* p[k] = B^{-1}r */
446: VecNorm(p[k],NORM_2,&rnorm);
447: break;
448: case KSP_NORM_UNPRECONDITIONED:
449: case KSP_NORM_NATURAL:
450: VecNorm(r,NORM_2,&rnorm);
451: break;
452: default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
453: }
454: PetscObjectSAWsTakeAccess((PetscObject)ksp);
455: ksp->rnorm = rnorm;
456: PetscObjectSAWsGrantAccess((PetscObject)ksp);
457: KSPLogResidualHistory(ksp,rnorm);
458: KSPMonitor(ksp,0,rnorm);
459: (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
460: } else ksp->reason = KSP_CONVERGED_ITERATING;
461: if (ksp->reason || ksp->max_it==0) {
462: if (ksp->max_it==0) ksp->reason = KSP_DIVERGED_ITS; /* This for a V(0,x) cycle */
463: return(0);
464: }
465: if (ksp->normtype != KSP_NORM_PRECONDITIONED) {
466: KSP_PCApply(ksp,r,p[k]); /* p[k] = B^{-1}r */
467: }
468: VecAYPX(p[k],scale,p[km1]); /* p[k] = scale B^{-1}r + p[km1] */
469: PetscObjectSAWsTakeAccess((PetscObject)ksp);
470: ksp->its = 1;
471: PetscObjectSAWsGrantAccess((PetscObject)ksp);
473: for (i=1; i<ksp->max_it; i++) {
474: PetscObjectSAWsTakeAccess((PetscObject)ksp);
475: ksp->its++;
476: PetscObjectSAWsGrantAccess((PetscObject)ksp);
478: KSP_MatMult(ksp,Amat,p[k],r); /* r = b - Ap[k] */
479: VecAYPX(r,-1.0,b);
480: /* calculate residual norm if requested */
481: if (ksp->normtype) {
482: switch (ksp->normtype) {
483: case KSP_NORM_PRECONDITIONED:
484: KSP_PCApply(ksp,r,p[kp1]); /* p[kp1] = B^{-1}r */
485: VecNorm(p[kp1],NORM_2,&rnorm);
486: break;
487: case KSP_NORM_UNPRECONDITIONED:
488: case KSP_NORM_NATURAL:
489: VecNorm(r,NORM_2,&rnorm);
490: break;
491: default:
492: rnorm = 0.0;
493: break;
494: }
495: KSPCheckNorm(ksp,rnorm);
496: PetscObjectSAWsTakeAccess((PetscObject)ksp);
497: ksp->rnorm = rnorm;
498: PetscObjectSAWsGrantAccess((PetscObject)ksp);
499: KSPLogResidualHistory(ksp,rnorm);
500: KSPMonitor(ksp,i,rnorm);
501: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
502: if (ksp->reason) break;
503: if (ksp->normtype != KSP_NORM_PRECONDITIONED) {
504: KSP_PCApply(ksp,r,p[kp1]); /* p[kp1] = B^{-1}r */
505: }
506: } else {
507: KSP_PCApply(ksp,r,p[kp1]); /* p[kp1] = B^{-1}r */
508: }
509: ksp->vec_sol = p[k];
511: c[kp1] = 2.0*mu*c[k] - c[km1];
512: omega = omegaprod*c[k]/c[kp1];
514: /* y^{k+1} = omega(y^{k} - y^{k-1} + Gamma*r^{k}) + y^{k-1} */
515: VecAXPBYPCZ(p[kp1],1.0-omega,omega,omega*Gamma*scale,p[km1],p[k]);
517: ktmp = km1;
518: km1 = k;
519: k = kp1;
520: kp1 = ktmp;
521: }
522: if (!ksp->reason) {
523: if (ksp->normtype) {
524: KSP_MatMult(ksp,Amat,p[k],r); /* r = b - Ap[k] */
525: VecAYPX(r,-1.0,b);
526: switch (ksp->normtype) {
527: case KSP_NORM_PRECONDITIONED:
528: KSP_PCApply(ksp,r,p[kp1]); /* p[kp1] = B^{-1}r */
529: VecNorm(p[kp1],NORM_2,&rnorm);
530: break;
531: case KSP_NORM_UNPRECONDITIONED:
532: case KSP_NORM_NATURAL:
533: VecNorm(r,NORM_2,&rnorm);
534: break;
535: default:
536: rnorm = 0.0;
537: break;
538: }
539: KSPCheckNorm(ksp,rnorm);
540: PetscObjectSAWsTakeAccess((PetscObject)ksp);
541: ksp->rnorm = rnorm;
542: PetscObjectSAWsGrantAccess((PetscObject)ksp);
543: KSPLogResidualHistory(ksp,rnorm);
544: KSPMonitor(ksp,i,rnorm);
545: }
546: if (ksp->its >= ksp->max_it) {
547: if (ksp->normtype != KSP_NORM_NONE) {
548: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
549: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
550: } else ksp->reason = KSP_CONVERGED_ITS;
551: }
552: }
554: /* make sure solution is in vector x */
555: ksp->vec_sol = sol_orig;
556: if (k) {
557: VecCopy(p[k],sol_orig);
558: }
559: return(0);
560: }
562: static PetscErrorCode KSPView_Chebyshev(KSP ksp,PetscViewer viewer)
563: {
564: KSP_Chebyshev *cheb = (KSP_Chebyshev*)ksp->data;
566: PetscBool iascii;
569: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
570: if (iascii) {
571: PetscViewerASCIIPrintf(viewer," eigenvalue estimates used: min = %g, max = %g\n",(double)cheb->emin,(double)cheb->emax);
572: if (cheb->kspest) {
573: PetscViewerASCIIPrintf(viewer," eigenvalues estimate via %s min %g, max %g\n",((PetscObject)(cheb->kspest))->type_name,(double)cheb->emin_computed,(double)cheb->emax_computed);
574: PetscViewerASCIIPrintf(viewer," eigenvalues estimated using %s with translations [%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]);
575: PetscViewerASCIIPushTab(viewer);
576: KSPView(cheb->kspest,viewer);
577: PetscViewerASCIIPopTab(viewer);
578: if (cheb->usenoisy) {
579: PetscViewerASCIIPrintf(viewer," estimating eigenvalues using noisy right hand side\n");
580: }
581: }
582: }
583: return(0);
584: }
586: static PetscErrorCode KSPDestroy_Chebyshev(KSP ksp)
587: {
588: KSP_Chebyshev *cheb = (KSP_Chebyshev*)ksp->data;
592: KSPDestroy(&cheb->kspest);
593: PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevSetEigenvalues_C",NULL);
594: PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigSet_C",NULL);
595: PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigSetUseNoisy_C",NULL);
596: PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigGetKSP_C",NULL);
597: KSPDestroyDefault(ksp);
598: return(0);
599: }
601: /*MC
602: KSPCHEBYSHEV - The preconditioned Chebyshev iterative method
604: Options Database Keys:
605: + -ksp_chebyshev_eigenvalues <emin,emax> - set approximations to the smallest and largest eigenvalues
606: of the preconditioned operator. If these are accurate you will get much faster convergence.
607: . -ksp_chebyshev_esteig <a,b,c,d> - estimate eigenvalues using a Krylov method, then use this
608: transform for Chebyshev eigenvalue bounds (KSPChebyshevEstEigSet())
609: . -ksp_chebyshev_esteig_steps - number of estimation steps
610: - -ksp_chebyshev_esteig_noisy - use noisy number generator to create right hand side for eigenvalue estimator
612: Level: beginner
614: Notes:
615: The Chebyshev method requires both the matrix and preconditioner to
616: be symmetric positive (semi) definite.
617: Only support for left preconditioning.
619: Chebyshev is configured as a smoother by default, targetting the "upper" part of the spectrum.
620: The user should call KSPChebyshevSetEigenvalues() if they have eigenvalue estimates.
622: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
623: KSPChebyshevSetEigenvalues(), KSPChebyshevEstEigSet(), KSPChebyshevEstEigSetUseNoisy()
624: KSPRICHARDSON, KSPCG, PCMG
626: M*/
628: PETSC_EXTERN PetscErrorCode KSPCreate_Chebyshev(KSP ksp)
629: {
631: KSP_Chebyshev *chebyshevP;
634: PetscNewLog(ksp,&chebyshevP);
636: ksp->data = (void*)chebyshevP;
637: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
638: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
639: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
640: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
642: chebyshevP->emin = 0.;
643: chebyshevP->emax = 0.;
645: chebyshevP->tform[0] = 0.0;
646: chebyshevP->tform[1] = 0.1;
647: chebyshevP->tform[2] = 0;
648: chebyshevP->tform[3] = 1.1;
649: chebyshevP->eststeps = 10;
650: chebyshevP->usenoisy = PETSC_TRUE;
651: ksp->setupnewmatrix = PETSC_TRUE;
653: ksp->ops->setup = KSPSetUp_Chebyshev;
654: ksp->ops->solve = KSPSolve_Chebyshev;
655: ksp->ops->destroy = KSPDestroy_Chebyshev;
656: ksp->ops->buildsolution = KSPBuildSolutionDefault;
657: ksp->ops->buildresidual = KSPBuildResidualDefault;
658: ksp->ops->setfromoptions = KSPSetFromOptions_Chebyshev;
659: ksp->ops->view = KSPView_Chebyshev;
660: ksp->ops->reset = KSPReset_Chebyshev;
662: PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevSetEigenvalues_C",KSPChebyshevSetEigenvalues_Chebyshev);
663: PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigSet_C",KSPChebyshevEstEigSet_Chebyshev);
664: PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigSetUseNoisy_C",KSPChebyshevEstEigSetUseNoisy_Chebyshev);
665: PetscObjectComposeFunction((PetscObject)ksp,"KSPChebyshevEstEigGetKSP_C",KSPChebyshevEstEigGetKSP_Chebyshev);
666: return(0);
667: }