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/sles/ksp/kspimpl.h
6: #include src/sles/ksp/impls/rich/richctx.h
8: #undef __FUNCT__
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: }
21: #undef __FUNCT__
23: int KSPSolve_Richardson(KSP ksp,int *its)
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: ierr = 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: ierr = 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: if (its) *its = maxit;
50: PCApplyRichardson(ksp->B,b,x,r,ksp->rtol,ksp->atol,ksp->divtol,maxit);
51: ksp->reason = KSP_DIVERGED_ITS; /* what should we really put here? */
52: return(0);
53: }
55: z = ksp->work[1];
56: scale = richardsonP->scale;
58: if (!ksp->guess_zero) { /* r <- b - A x */
59: KSP_MatMult(ksp,Amat,x,r);
60: VecAYPX(&mone,b,r);
61: } else {
62: VecCopy(b,r);
63: }
65: for (i=0; i<maxit; i++) {
66: PetscObjectTakeAccess(ksp);
67: ksp->its++;
68: PetscObjectGrantAccess(ksp);
70: if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
71: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
72: KSPMonitor(ksp,i,rnorm);
73: }
75: KSP_PCApply(ksp,ksp->B,r,z); /* z <- B r */
77: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
78: VecNorm(z,NORM_2,&rnorm); /* rnorm <- z'*z */
79: KSPMonitor(ksp,i,rnorm);
80: }
82: VecAXPY(&scale,z,x); /* x <- x + scale z */
83: if (ksp->normtype != KSP_NO_NORM) {
84: ierr = PetscObjectTakeAccess(ksp);
85: ksp->rnorm = rnorm;
86: ierr = PetscObjectGrantAccess(ksp);
87: KSPLogResidualHistory(ksp,rnorm);
89: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
90: if (ksp->reason) break;
91: }
92:
93: KSP_MatMult(ksp,Amat,x,r); /* r <- b - Ax */
94: VecAYPX(&mone,b,r);
95: }
96: if (!ksp->reason) {
97: ksp->reason = KSP_DIVERGED_ITS;
98: if (ksp->normtype != KSP_NO_NORM) {
99: if (ksp->normtype == KSP_UNPRECONDITIONED_NORM){
100: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
101: } else {
102: KSP_PCApply(ksp,ksp->B,r,z); /* z <- B r */
103: VecNorm(z,NORM_2,&rnorm); /* rnorm <- z'*z */
104: }
105: }
106: PetscObjectTakeAccess(ksp);
107: ksp->rnorm = rnorm;
108: PetscObjectGrantAccess(ksp);
109: KSPLogResidualHistory(ksp,rnorm);
110: KSPMonitor(ksp,i,rnorm);
111: i--;
112: } else if (!ksp->reason) {
113: ksp->reason = KSP_DIVERGED_ITS;
114: }
116: if (its) *its = ksp->its;
117: return(0);
118: }
120: #undef __FUNCT__
122: int KSPView_Richardson(KSP ksp,PetscViewer viewer)
123: {
124: KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
125: int ierr;
126: PetscTruth isascii;
129: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
130: if (isascii) {
131: PetscViewerASCIIPrintf(viewer," Richardson: damping factor=%gn",richardsonP->scale);
132: } else {
133: SETERRQ1(1,"Viewer type %s not supported for KSP Richardson",((PetscObject)viewer)->type_name);
134: }
135: return(0);
136: }
138: #undef __FUNCT__
140: int KSPSetFromOptions_Richardson(KSP ksp)
141: {
142: KSP_Richardson *rich = (KSP_Richardson*)ksp->data;
143: int ierr;
144: PetscReal tmp;
145: PetscTruth flg;
148: PetscOptionsHead("KSP Richardson Options");
149: PetscOptionsReal("-ksp_richardson_scale","damping factor","KSPRichardsonSetScale",rich->scale,&tmp,&flg);
150: if (flg) { KSPRichardsonSetScale(ksp,tmp); }
151: PetscOptionsTail();
152: return(0);
153: }
155: EXTERN_C_BEGIN
156: #undef __FUNCT__
158: int KSPRichardsonSetScale_Richardson(KSP ksp,PetscReal scale)
159: {
160: KSP_Richardson *richardsonP;
163: richardsonP = (KSP_Richardson*)ksp->data;
164: richardsonP->scale = scale;
165: return(0);
166: }
167: EXTERN_C_END
169: EXTERN_C_BEGIN
170: #undef __FUNCT__
172: int KSPCreate_Richardson(KSP ksp)
173: {
174: int ierr;
175: KSP_Richardson *richardsonP;
178: PetscNew(KSP_Richardson,&richardsonP);
179: PetscLogObjectMemory(ksp,sizeof(KSP_Richardson));
180: ksp->data = (void*)richardsonP;
181: richardsonP->scale = 1.0;
182: ksp->ops->setup = KSPSetUp_Richardson;
183: ksp->ops->solve = KSPSolve_Richardson;
184: ksp->ops->destroy = KSPDefaultDestroy;
185: ksp->ops->buildsolution = KSPDefaultBuildSolution;
186: ksp->ops->buildresidual = KSPDefaultBuildResidual;
187: ksp->ops->view = KSPView_Richardson;
188: ksp->ops->setfromoptions = KSPSetFromOptions_Richardson;
190: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetScale_C",
191: "KSPRichardsonSetScale_Richardson",
192: KSPRichardsonSetScale_Richardson);
193: return(0);
194: }
195: EXTERN_C_END