Actual source code: rich.c
1: #define PETSCKSP_DLL
3: /*
4: This implements Richardson Iteration.
5: */
6: #include include/private/kspimpl.h
7: #include src/ksp/ksp/impls/rich/richctx.h
11: PetscErrorCode KSPSetUp_Richardson(KSP ksp)
12: {
16: if (ksp->pc_side == PC_RIGHT) {SETERRQ(PETSC_ERR_SUP,"no right preconditioning for KSPRICHARDSON");}
17: else if (ksp->pc_side == PC_SYMMETRIC) {SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPRICHARDSON");}
18: KSPDefaultGetWork(ksp,2);
19: return(0);
20: }
24: PetscErrorCode KSPSolve_Richardson(KSP ksp)
25: {
27: PetscInt i,maxit;
28: MatStructure pflag;
29: PetscReal rnorm = 0.0;
30: PetscScalar scale,dt;
31: Vec x,b,r,z;
32: PetscInt xs, ws;
33: Mat Amat,Pmat;
34: KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
35: PetscTruth exists,diagonalscale;
38: PCDiagonalScale(ksp->pc,&diagonalscale);
39: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);
41: ksp->its = 0;
43: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
44: x = ksp->vec_sol;
45: b = ksp->vec_rhs;
46: VecGetSize(x,&xs);
47: VecGetSize(ksp->work[0],&ws);
48: if (xs != ws) {
49: KSPDefaultFreeWork(ksp);
50: KSPDefaultGetWork(ksp,2);
51: }
52: r = ksp->work[0];
53: z = ksp->work[1];
54: maxit = ksp->max_it;
56: /* if user has provided fast Richardson code use that */
57: PCApplyRichardsonExists(ksp->pc,&exists);
58: if (exists && !ksp->numbermonitors && !ksp->transpose_solve) {
59: ksp->normtype = KSP_NORM_NO;
60: PCApplyRichardson(ksp->pc,b,x,r,ksp->rtol,ksp->abstol,ksp->divtol,maxit);
61: if (ksp->normtype != KSP_NORM_NO) {
62: KSP_MatMult(ksp,Amat,x,r);
63: VecAYPX(r,-1.0,b);
64: if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->pc_side == PC_RIGHT) {
65: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
66: } else {
67: PetscExceptionTry1((KSP_PCApply(ksp,r,z)),PETSC_ERR_SUP);
68: if (PetscExceptionCaught(ierr,PETSC_ERR_SUP)) {
69: ksp->reason = KSP_CONVERGED_ITS;
70: return(0);
71: }
72: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
73: VecNorm(z,NORM_2,&rnorm); /* rnorm <- r'B'*Br */
74: } else {
75: VecDot(r,z,&dt);
76: rnorm = PetscAbsScalar(dt); /* rnorm <- z'*r = r'Br = e'*A'*B*A*e */
77: }
78: }
79: (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
80: }
81: if (!ksp->reason) ksp->reason = KSP_CONVERGED_ITS;
82: return(0);
83: }
85: scale = richardsonP->scale;
87: if (!ksp->guess_zero) { /* r <- b - A x */
88: KSP_MatMult(ksp,Amat,x,r);
89: VecAYPX(r,-1.0,b);
90: } else {
91: VecCopy(b,r);
92: }
94: for (i=0; i<maxit; i++) {
95: PetscObjectTakeAccess(ksp);
96: ksp->its++;
97: PetscObjectGrantAccess(ksp);
99: if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
100: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
101: KSPMonitor(ksp,i,rnorm);
102: }
104: KSP_PCApply(ksp,r,z); /* z <- B r */
106: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
107: VecNorm(z,NORM_2,&rnorm); /* rnorm <- z'*z */
108: KSPMonitor(ksp,i,rnorm);
109: }
111: VecAXPY(x,scale,z); /* x <- x + scale z */
113: if (ksp->normtype != KSP_NORM_NO) {
114: PetscObjectTakeAccess(ksp);
115: ksp->rnorm = rnorm;
116: PetscObjectGrantAccess(ksp);
117: KSPLogResidualHistory(ksp,rnorm);
118: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
119: if (ksp->reason) break;
120: }
121:
122: KSP_MatMult(ksp,Amat,x,r); /* r <- b - Ax */
123: VecAYPX(r,-1.0,b);
124: }
125: if (!ksp->reason) {
126: if (ksp->normtype != KSP_NORM_NO) {
127: if (ksp->normtype == KSP_NORM_UNPRECONDITIONED){
128: VecNorm(r,NORM_2,&rnorm); /* rnorm <- r'*r */
129: } else {
130: KSP_PCApply(ksp,r,z); /* z <- B r */
131: VecNorm(z,NORM_2,&rnorm); /* rnorm <- z'*z */
132: }
133: PetscObjectTakeAccess(ksp);
134: ksp->rnorm = rnorm;
135: PetscObjectGrantAccess(ksp);
136: KSPLogResidualHistory(ksp,rnorm);
137: KSPMonitor(ksp,i,rnorm);
138: }
139: if (ksp->its >= ksp->max_it) {
140: if (ksp->normtype != KSP_NORM_NO) {
141: (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
142: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
143: } else {
144: ksp->reason = KSP_CONVERGED_ITS;
145: }
146: }
147: }
148: return(0);
149: }
153: PetscErrorCode KSPView_Richardson(KSP ksp,PetscViewer viewer)
154: {
155: KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
157: PetscTruth iascii;
160: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
161: if (iascii) {
162: PetscViewerASCIIPrintf(viewer," Richardson: damping factor=%G\n",richardsonP->scale);
163: } else {
164: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP Richardson",((PetscObject)viewer)->type_name);
165: }
166: return(0);
167: }
171: PetscErrorCode KSPSetFromOptions_Richardson(KSP ksp)
172: {
173: KSP_Richardson *rich = (KSP_Richardson*)ksp->data;
175: PetscReal tmp;
176: PetscTruth flg;
179: PetscOptionsHead("KSP Richardson Options");
180: PetscOptionsReal("-ksp_richardson_scale","damping factor","KSPRichardsonSetScale",rich->scale,&tmp,&flg);
181: if (flg) { KSPRichardsonSetScale(ksp,tmp); }
182: PetscOptionsTail();
183: return(0);
184: }
188: PetscErrorCode KSPDestroy_Richardson(KSP ksp)
189: {
193: KSPDefaultDestroy(ksp);
194: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetScale_C","",PETSC_NULL);
195: return(0);
196: }
201: PetscErrorCode KSPRichardsonSetScale_Richardson(KSP ksp,PetscReal scale)
202: {
203: KSP_Richardson *richardsonP;
206: richardsonP = (KSP_Richardson*)ksp->data;
207: richardsonP->scale = scale;
208: return(0);
209: }
212: /*MC
213: KSPRICHARDSON - The preconditioned Richardson iterative method
215: Options Database Keys:
216: . -ksp_richardson_scale - damping factor on the correction (defaults to 1.0)
218: Level: beginner
220: Notes: x^{n+1} = x^{n} + scale*B(b - A x^{n})
222: This method often (usually) will not converge unless scale is very small
224: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
225: KSPRichardsonSetScale()
227: M*/
232: PetscErrorCode KSPCreate_Richardson(KSP ksp)
233: {
235: KSP_Richardson *richardsonP;
238: PetscNew(KSP_Richardson,&richardsonP);
239: PetscLogObjectMemory(ksp,sizeof(KSP_Richardson));
240: ksp->data = (void*)richardsonP;
242: ksp->normtype = KSP_NORM_PRECONDITIONED;
243: ksp->pc_side = PC_LEFT;
245: ksp->ops->setup = KSPSetUp_Richardson;
246: ksp->ops->solve = KSPSolve_Richardson;
247: ksp->ops->destroy = KSPDestroy_Richardson;
248: ksp->ops->buildsolution = KSPDefaultBuildSolution;
249: ksp->ops->buildresidual = KSPDefaultBuildResidual;
250: ksp->ops->view = KSPView_Richardson;
251: ksp->ops->setfromoptions = KSPSetFromOptions_Richardson;
253: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetScale_C",
254: "KSPRichardsonSetScale_Richardson",
255: KSPRichardsonSetScale_Richardson);
256: richardsonP->scale = 1.0;
257: return(0);
258: }