Actual source code: rich.c
1: /*$Id: rich.c,v 1.104 2001/08/21 21:03:36 bsmith Exp $*/
2: /*
3: This implements Richardson Iteration.
4: */
5: #include src/ksp/ksp/kspimpl.h
6: #include src/ksp/ksp/impls/rich/richctx.h
10: int KSPSetUp_Richardson(KSP ksp)
11: {
15: if (ksp->pc_side == PC_RIGHT) {SETERRQ(2,"no right preconditioning for KSPRICHARDSON");}
16: else if (ksp->pc_side == PC_SYMMETRIC) {SETERRQ(2,"no symmetric preconditioning for KSPRICHARDSON");}
17: KSPDefaultGetWork(ksp,2);
18: return(0);
19: }
23: int KSPSolve_Richardson(KSP ksp)
24: {
25: int i,maxit,ierr;
26: MatStructure pflag;
27: PetscReal rnorm = 0.0;
28: PetscScalar scale,mone = -1.0;
29: Vec x,b,r,z;
30: Mat Amat,Pmat;
31: KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
32: PetscTruth exists,diagonalscale;
35: PCDiagonalScale(ksp->B,&diagonalscale);
36: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
38: ksp->its = 0;
40: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
41: x = ksp->vec_sol;
42: b = ksp->vec_rhs;
43: r = ksp->work[0];
44: maxit = ksp->max_it;
46: /* if user has provided fast Richardson code use that */
47: PCApplyRichardsonExists(ksp->B,&exists);
48: if (exists && !ksp->numbermonitors && !ksp->transpose_solve) {
49: PCApplyRichardson(ksp->B,b,x,r,ksp->rtol,ksp->atol,ksp->divtol,maxit);
50: ksp->reason = KSP_DIVERGED_ITS; /* what should we really put here? */
51: return(0);
52: }
54: z = ksp->work[1];
55: scale = richardsonP->scale;
57: if (!ksp->guess_zero) { /* r <- b - A x */
58: KSP_MatMult(ksp,Amat,x,r);
59: VecAYPX(&mone,b,r);
60: } else {
61: VecCopy(b,r);
62: }
64: for (i=0; i<maxit; i++) {
65: PetscObjectTakeAccess(ksp);
66: ksp->its++;
67: PetscObjectGrantAccess(ksp);
69: if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
70: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
71: KSPMonitor(ksp,i,rnorm);
72: }
74: KSP_PCApply(ksp,ksp->B,r,z); /* z <- B r */
76: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
77: VecNorm(z,NORM_2,&rnorm); /* rnorm <- z'*z */
78: KSPMonitor(ksp,i,rnorm);
79: }
81: VecAXPY(&scale,z,x); /* x <- x + scale z */
82: if (ksp->normtype != KSP_NO_NORM) {
83: PetscObjectTakeAccess(ksp);
84: ksp->rnorm = rnorm;
85: PetscObjectGrantAccess(ksp);
86: KSPLogResidualHistory(ksp,rnorm);
88: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
89: if (ksp->reason) break;
90: }
91:
92: KSP_MatMult(ksp,Amat,x,r); /* r <- b - Ax */
93: VecAYPX(&mone,b,r);
94: }
95: if (!ksp->reason) {
96: ksp->reason = KSP_DIVERGED_ITS;
97: if (ksp->normtype != KSP_NO_NORM) {
98: if (ksp->normtype == KSP_UNPRECONDITIONED_NORM){
99: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
100: } else {
101: KSP_PCApply(ksp,ksp->B,r,z); /* z <- B r */
102: VecNorm(z,NORM_2,&rnorm); /* rnorm <- z'*z */
103: }
104: }
105: PetscObjectTakeAccess(ksp);
106: ksp->rnorm = rnorm;
107: PetscObjectGrantAccess(ksp);
108: KSPLogResidualHistory(ksp,rnorm);
109: KSPMonitor(ksp,i,rnorm);
110: i--;
111: } else if (!ksp->reason) {
112: ksp->reason = KSP_DIVERGED_ITS;
113: }
115: return(0);
116: }
120: int KSPView_Richardson(KSP ksp,PetscViewer viewer)
121: {
122: KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
123: int ierr;
124: PetscTruth isascii;
127: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
128: if (isascii) {
129: PetscViewerASCIIPrintf(viewer," Richardson: damping factor=%g\n",richardsonP->scale);
130: } else {
131: SETERRQ1(1,"Viewer type %s not supported for KSP Richardson",((PetscObject)viewer)->type_name);
132: }
133: return(0);
134: }
138: int KSPSetFromOptions_Richardson(KSP ksp)
139: {
140: KSP_Richardson *rich = (KSP_Richardson*)ksp->data;
141: int ierr;
142: PetscReal tmp;
143: PetscTruth flg;
146: PetscOptionsHead("KSP Richardson Options");
147: PetscOptionsReal("-ksp_richardson_scale","damping factor","KSPRichardsonSetScale",rich->scale,&tmp,&flg);
148: if (flg) { KSPRichardsonSetScale(ksp,tmp); }
149: PetscOptionsTail();
150: return(0);
151: }
153: EXTERN_C_BEGIN
156: int KSPRichardsonSetScale_Richardson(KSP ksp,PetscReal scale)
157: {
158: KSP_Richardson *richardsonP;
161: richardsonP = (KSP_Richardson*)ksp->data;
162: richardsonP->scale = scale;
163: return(0);
164: }
165: EXTERN_C_END
167: /*MC
168: KSPRICHARDSON - The preconditioned Richardson iterative method
170: Options Database Keys:
171: . -ksp_richardson_scale - damping factor on the correction (defaults to 1.0)
173: Level: beginner
175: Notes: x^{n+1} = x^{n} + scale*B(b - A x^{n})
177: This method often (usually) will not converge unless scale is very small
179: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
180: KSPRichardsonSetScale()
182: M*/
184: EXTERN_C_BEGIN
187: int KSPCreate_Richardson(KSP ksp)
188: {
189: int ierr;
190: KSP_Richardson *richardsonP;
193: PetscNew(KSP_Richardson,&richardsonP);
194: PetscLogObjectMemory(ksp,sizeof(KSP_Richardson));
195: ksp->data = (void*)richardsonP;
196: richardsonP->scale = 1.0;
197: ksp->ops->setup = KSPSetUp_Richardson;
198: ksp->ops->solve = KSPSolve_Richardson;
199: ksp->ops->destroy = KSPDefaultDestroy;
200: ksp->ops->buildsolution = KSPDefaultBuildSolution;
201: ksp->ops->buildresidual = KSPDefaultBuildResidual;
202: ksp->ops->view = KSPView_Richardson;
203: ksp->ops->setfromoptions = KSPSetFromOptions_Richardson;
205: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetScale_C",
206: "KSPRichardsonSetScale_Richardson",
207: KSPRichardsonSetScale_Richardson);
208: return(0);
209: }
210: EXTERN_C_END