Actual source code: rich.c
petsc-3.10.5 2019-03-28
2: /*
3: This implements Richardson Iteration.
4: */
5: #include <../src/ksp/ksp/impls/rich/richardsonimpl.h>
7: PetscErrorCode KSPSetUp_Richardson(KSP ksp)
8: {
10: KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
13: if (richardsonP->selfscale) {
14: KSPSetWorkVecs(ksp,4);
15: } else {
16: KSPSetWorkVecs(ksp,2);
17: }
18: return(0);
19: }
21: PetscErrorCode KSPSolve_Richardson(KSP ksp)
22: {
24: PetscInt i,maxit;
25: PetscReal rnorm = 0.0,abr;
26: PetscScalar scale,rdot;
27: Vec x,b,r,z,w = NULL,y = NULL;
28: PetscInt xs, ws;
29: Mat Amat,Pmat;
30: KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
31: PetscBool exists,diagonalscale;
32: MatNullSpace nullsp;
35: PCGetDiagonalScale(ksp->pc,&diagonalscale);
36: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
38: ksp->its = 0;
40: PCGetOperators(ksp->pc,&Amat,&Pmat);
41: x = ksp->vec_sol;
42: b = ksp->vec_rhs;
43: VecGetSize(x,&xs);
44: VecGetSize(ksp->work[0],&ws);
45: if (xs != ws) {
46: if (richardsonP->selfscale) {
47: KSPSetWorkVecs(ksp,4);
48: } else {
49: KSPSetWorkVecs(ksp,2);
50: }
51: }
52: r = ksp->work[0];
53: z = ksp->work[1];
54: if (richardsonP->selfscale) {
55: w = ksp->work[2];
56: y = ksp->work[3];
57: }
58: maxit = ksp->max_it;
60: /* if user has provided fast Richardson code use that */
61: PCApplyRichardsonExists(ksp->pc,&exists);
62: MatGetNullSpace(Pmat,&nullsp);
63: if (exists && maxit > 0 && richardsonP->scale == 1.0 && !ksp->numbermonitors && !ksp->transpose_solve && !nullsp) {
64: PCRichardsonConvergedReason reason;
65: PCApplyRichardson(ksp->pc,b,x,r,ksp->rtol,ksp->abstol,ksp->divtol,maxit,ksp->guess_zero,&ksp->its,&reason);
66: ksp->reason = (KSPConvergedReason)reason;
67: return(0);
68: } else if (exists && maxit > 0 && !ksp->numbermonitors && !ksp->transpose_solve && !nullsp) {
69: PetscInfo(ksp,"KSPSolve_Richardson: Warning, skipping optimized PCApplyRichardson() because scale factor is not 1.0\n");
70: }
72: if (!ksp->guess_zero) { /* r <- b - A x */
73: KSP_MatMult(ksp,Amat,x,r);
74: VecAYPX(r,-1.0,b);
75: } else {
76: VecCopy(b,r);
77: }
79: ksp->its = 0;
80: if (richardsonP->selfscale) {
81: KSP_PCApply(ksp,r,z); /* z <- B r */
82: for (i=0; i<maxit; i++) {
84: if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
85: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
86: KSPMonitor(ksp,i,rnorm);
87: ksp->rnorm = rnorm;
88: KSPLogResidualHistory(ksp,rnorm);
89: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
90: if (ksp->reason) break;
91: } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
92: VecNorm(z,NORM_2,&rnorm); /* rnorm <- z'*z */
93: KSPMonitor(ksp,i,rnorm);
94: ksp->rnorm = rnorm;
95: KSPLogResidualHistory(ksp,rnorm);
96: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
97: if (ksp->reason) break;
98: }
99: KSP_PCApplyBAorAB(ksp,z,y,w); /* y = BAz = BABr */
100: VecDotNorm2(z,y,&rdot,&abr); /* rdot = (Br)^T(BABR); abr = (BABr)^T (BABr) */
101: scale = rdot/abr;
102: PetscInfo1(ksp,"Self-scale factor %g\n",(double)PetscRealPart(scale));
103: VecAXPY(x,scale,z); /* x <- x + scale z */
104: VecAXPY(r,-scale,w); /* r <- r - scale*Az */
105: VecAXPY(z,-scale,y); /* z <- z - scale*y */
106: ksp->its++;
107: }
108: } else {
109: for (i=0; i<maxit; i++) {
111: if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
112: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
113: KSPMonitor(ksp,i,rnorm);
114: ksp->rnorm = rnorm;
115: KSPLogResidualHistory(ksp,rnorm);
116: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
117: if (ksp->reason) break;
118: }
120: KSP_PCApply(ksp,r,z); /* z <- B r */
122: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
123: VecNorm(z,NORM_2,&rnorm); /* rnorm <- z'*z */
124: KSPMonitor(ksp,i,rnorm);
125: ksp->rnorm = rnorm;
126: KSPLogResidualHistory(ksp,rnorm);
127: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
128: if (ksp->reason) break;
129: }
131: VecAXPY(x,richardsonP->scale,z); /* x <- x + scale z */
132: ksp->its++;
134: if (i+1 < maxit || ksp->normtype != KSP_NORM_NONE) {
135: KSP_MatMult(ksp,Amat,x,r); /* r <- b - Ax */
136: VecAYPX(r,-1.0,b);
137: }
138: }
139: }
140: if (!ksp->reason) {
141: if (ksp->normtype != KSP_NORM_NONE) {
142: if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
143: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
144: } else {
145: KSP_PCApply(ksp,r,z); /* z <- B r */
146: VecNorm(z,NORM_2,&rnorm); /* rnorm <- z'*z */
147: }
148: ksp->rnorm = rnorm;
149: KSPLogResidualHistory(ksp,rnorm);
150: KSPMonitor(ksp,i,rnorm);
151: }
152: if (ksp->its >= ksp->max_it) {
153: if (ksp->normtype != KSP_NORM_NONE) {
154: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
155: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
156: } else {
157: ksp->reason = KSP_CONVERGED_ITS;
158: }
159: }
160: }
161: return(0);
162: }
164: PetscErrorCode KSPView_Richardson(KSP ksp,PetscViewer viewer)
165: {
166: KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
168: PetscBool iascii;
171: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
172: if (iascii) {
173: if (richardsonP->selfscale) {
174: PetscViewerASCIIPrintf(viewer," using self-scale best computed damping factor\n");
175: } else {
176: PetscViewerASCIIPrintf(viewer," damping factor=%g\n",(double)richardsonP->scale);
177: }
178: }
179: return(0);
180: }
182: PetscErrorCode KSPSetFromOptions_Richardson(PetscOptionItems *PetscOptionsObject,KSP ksp)
183: {
184: KSP_Richardson *rich = (KSP_Richardson*)ksp->data;
186: PetscReal tmp;
187: PetscBool flg,flg2;
190: PetscOptionsHead(PetscOptionsObject,"KSP Richardson Options");
191: PetscOptionsReal("-ksp_richardson_scale","damping factor","KSPRichardsonSetScale",rich->scale,&tmp,&flg);
192: if (flg) { KSPRichardsonSetScale(ksp,tmp); }
193: PetscOptionsBool("-ksp_richardson_self_scale","dynamically determine optimal damping factor","KSPRichardsonSetSelfScale",rich->selfscale,&flg2,&flg);
194: if (flg) { KSPRichardsonSetSelfScale(ksp,flg2); }
195: PetscOptionsTail();
196: return(0);
197: }
199: PetscErrorCode KSPDestroy_Richardson(KSP ksp)
200: {
204: PetscObjectComposeFunction((PetscObject)ksp,"KSPRichardsonSetScale_C",NULL);
205: KSPDestroyDefault(ksp);
206: return(0);
207: }
209: static PetscErrorCode KSPRichardsonSetScale_Richardson(KSP ksp,PetscReal scale)
210: {
211: KSP_Richardson *richardsonP;
214: richardsonP = (KSP_Richardson*)ksp->data;
215: richardsonP->scale = scale;
216: return(0);
217: }
219: static PetscErrorCode KSPRichardsonSetSelfScale_Richardson(KSP ksp,PetscBool selfscale)
220: {
221: KSP_Richardson *richardsonP;
224: richardsonP = (KSP_Richardson*)ksp->data;
225: richardsonP->selfscale = selfscale;
226: return(0);
227: }
229: /*MC
230: KSPRICHARDSON - The preconditioned Richardson iterative method
232: Options Database Keys:
233: . -ksp_richardson_scale - damping factor on the correction (defaults to 1.0)
235: Level: beginner
237: Notes:
238: x^{n+1} = x^{n} + scale*B(b - A x^{n})
240: Here B is the application of the preconditioner
242: This method often (usually) will not converge unless scale is very small.
244: Notes:
245: For some preconditioners, currently SOR, the convergence test is skipped to improve speed,
246: thus it always iterates the maximum number of iterations you've selected. When -ksp_monitor
247: (or any other monitor) is turned on, the norm is computed at each iteration and so the convergence test is run unless
248: you specifically call KSPSetNormType(ksp,KSP_NORM_NONE);
250: For some preconditioners, currently PCMG and PCHYPRE with BoomerAMG if -ksp_monitor (and also
251: any other monitor) is not turned on then the convergence test is done by the preconditioner itself and
252: so the solver may run more or fewer iterations then if -ksp_monitor is selected.
254: Supports only left preconditioning
256: References:
257: . 1. - L. F. Richardson, "The Approximate Arithmetical Solution by Finite Differences of Physical Problems Involving
258: Differential Equations, with an Application to the Stresses in a Masonry Dam",
259: Philosophical Transactions of the Royal Society of London. Series A,
260: Containing Papers of a Mathematical or Physical Character, Vol. 210, 1911 (1911).
262: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
263: KSPRichardsonSetScale()
265: M*/
267: PETSC_EXTERN PetscErrorCode KSPCreate_Richardson(KSP ksp)
268: {
270: KSP_Richardson *richardsonP;
273: PetscNewLog(ksp,&richardsonP);
274: ksp->data = (void*)richardsonP;
276: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
277: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
278: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
280: ksp->ops->setup = KSPSetUp_Richardson;
281: ksp->ops->solve = KSPSolve_Richardson;
282: ksp->ops->destroy = KSPDestroy_Richardson;
283: ksp->ops->buildsolution = KSPBuildSolutionDefault;
284: ksp->ops->buildresidual = KSPBuildResidualDefault;
285: ksp->ops->view = KSPView_Richardson;
286: ksp->ops->setfromoptions = KSPSetFromOptions_Richardson;
288: PetscObjectComposeFunction((PetscObject)ksp,"KSPRichardsonSetScale_C",KSPRichardsonSetScale_Richardson);
289: PetscObjectComposeFunction((PetscObject)ksp,"KSPRichardsonSetSelfScale_C",KSPRichardsonSetSelfScale_Richardson);
291: richardsonP->scale = 1.0;
292: return(0);
293: }