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 */
112:     if (ksp->normtype != KSP_NORM_NO) {
113:       PetscObjectTakeAccess(ksp);
114:       ksp->rnorm = rnorm;
115:       PetscObjectGrantAccess(ksp);
116:       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:     ksp->reason = KSP_DIVERGED_ITS;
127:     if (ksp->normtype != KSP_NORM_NO) {
128:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED){
129:         VecNorm(r,NORM_2,&rnorm);     /*   rnorm <- r'*r     */
130:       } else {
131:         KSP_PCApply(ksp,r,z);   /*   z <- B r          */
132:         VecNorm(z,NORM_2,&rnorm);     /*   rnorm <- z'*z     */
133:       }
134:     }
135:     PetscObjectTakeAccess(ksp);
136:     ksp->rnorm = rnorm;
137:     PetscObjectGrantAccess(ksp);
138:     KSPLogResidualHistory(ksp,rnorm);
139:     KSPMonitor(ksp,i,rnorm);
140:     i--;
141:   } else if (!ksp->reason) {
142:     ksp->reason = KSP_DIVERGED_ITS;
143:   }

145:   return(0);
146: }

150: PetscErrorCode KSPView_Richardson(KSP ksp,PetscViewer viewer)
151: {
152:   KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
154:   PetscTruth     iascii;

157:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
158:   if (iascii) {
159:     PetscViewerASCIIPrintf(viewer,"  Richardson: damping factor=%G\n",richardsonP->scale);
160:   } else {
161:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP Richardson",((PetscObject)viewer)->type_name);
162:   }
163:   return(0);
164: }

168: PetscErrorCode KSPSetFromOptions_Richardson(KSP ksp)
169: {
170:   KSP_Richardson *rich = (KSP_Richardson*)ksp->data;
172:   PetscReal      tmp;
173:   PetscTruth     flg;

176:   PetscOptionsHead("KSP Richardson Options");
177:     PetscOptionsReal("-ksp_richardson_scale","damping factor","KSPRichardsonSetScale",rich->scale,&tmp,&flg);
178:     if (flg) { KSPRichardsonSetScale(ksp,tmp); }
179:   PetscOptionsTail();
180:   return(0);
181: }

185: PetscErrorCode KSPDestroy_Richardson(KSP ksp)
186: {

190:   KSPDefaultDestroy(ksp);
191:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetScale_C","",PETSC_NULL);
192:   return(0);
193: }

198: PetscErrorCode  KSPRichardsonSetScale_Richardson(KSP ksp,PetscReal scale)
199: {
200:   KSP_Richardson *richardsonP;

203:   richardsonP = (KSP_Richardson*)ksp->data;
204:   richardsonP->scale = scale;
205:   return(0);
206: }

209: /*MC
210:      KSPRICHARDSON - The preconditioned Richardson iterative method

212:    Options Database Keys:
213: .   -ksp_richardson_scale - damping factor on the correction (defaults to 1.0)

215:    Level: beginner

217:    Notes: x^{n+1} = x^{n} + scale*B(b - A x^{n})

219:           This method often (usually) will not converge unless scale is very small

221: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
222:            KSPRichardsonSetScale()

224: M*/

229: PetscErrorCode  KSPCreate_Richardson(KSP ksp)
230: {
232:   KSP_Richardson *richardsonP;

235:   PetscNew(KSP_Richardson,&richardsonP);
236:   PetscLogObjectMemory(ksp,sizeof(KSP_Richardson));
237:   ksp->data                        = (void*)richardsonP;

239:   ksp->normtype                    = KSP_NORM_PRECONDITIONED;
240:   ksp->pc_side                     = PC_LEFT;

242:   ksp->ops->setup                  = KSPSetUp_Richardson;
243:   ksp->ops->solve                  = KSPSolve_Richardson;
244:   ksp->ops->destroy                = KSPDestroy_Richardson;
245:   ksp->ops->buildsolution          = KSPDefaultBuildSolution;
246:   ksp->ops->buildresidual          = KSPDefaultBuildResidual;
247:   ksp->ops->view                   = KSPView_Richardson;
248:   ksp->ops->setfromoptions         = KSPSetFromOptions_Richardson;

250:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetScale_C",
251:                                     "KSPRichardsonSetScale_Richardson",
252:                                     KSPRichardsonSetScale_Richardson);
253:   richardsonP->scale               = 1.0;
254:   return(0);
255: }