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: int KSPSetUp_Richardson(KSP ksp)
  9: {

 13:   if (ksp->pc_side == PC_RIGHT) {SETERRQ(2,"no right preconditioning for KSPRICHARDSON");}
 14:   else if (ksp->pc_side == PC_SYMMETRIC) {SETERRQ(2,"no symmetric preconditioning for KSPRICHARDSON");}
 15:   KSPDefaultGetWork(ksp,2);
 16:   return(0);
 17: }

 19: int  KSPSolve_Richardson(KSP ksp,int *its)
 20: {
 21:   int             i,maxit,ierr;
 22:   MatStructure    pflag;
 23:   PetscReal       rnorm = 0.0;
 24:   PetscScalar     scale,mone = -1.0;
 25:   Vec             x,b,r,z;
 26:   Mat             Amat,Pmat;
 27:   KSP_Richardson  *richardsonP = (KSP_Richardson*)ksp->data;
 28:   PetscTruth      exists,diagonalscale;

 31:   ierr    = PCDiagonalScale(ksp->B,&diagonalscale);
 32:   if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);

 34:   ksp->its = 0;

 36:   ierr    = PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
 37:   x       = ksp->vec_sol;
 38:   b       = ksp->vec_rhs;
 39:   r       = ksp->work[0];
 40:   maxit   = ksp->max_it;

 42:   /* if user has provided fast Richardson code use that */
 43:   PCApplyRichardsonExists(ksp->B,&exists);
 44:   if (exists && !ksp->numbermonitors && !ksp->transpose_solve) {
 45:     if (its) *its = maxit;
 46:     PCApplyRichardson(ksp->B,b,x,r,ksp->rtol,ksp->atol,ksp->divtol,maxit);
 47:     ksp->reason = KSP_DIVERGED_ITS; /* what should we really put here? */
 48:     return(0);
 49:   }

 51:   z       = ksp->work[1];
 52:   scale   = richardsonP->scale;

 54:   if (!ksp->guess_zero) {                          /*   r <- b - A x     */
 55:     KSP_MatMult(ksp,Amat,x,r);
 56:     VecAYPX(&mone,b,r);
 57:   } else {
 58:     VecCopy(b,r);
 59:   }

 61:   for (i=0; i<maxit; i++) {
 62:     PetscObjectTakeAccess(ksp);
 63:     ksp->its++;
 64:     PetscObjectGrantAccess(ksp);

 66:     if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
 67:       VecNorm(r,NORM_2,&rnorm); /*   rnorm <- r'*r     */
 68:       KSPMonitor(ksp,i,rnorm);
 69:     }

 71:     KSP_PCApply(ksp,ksp->B,r,z);    /*   z <- B r          */

 73:     if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
 74:       VecNorm(z,NORM_2,&rnorm); /*   rnorm <- z'*z     */
 75:       KSPMonitor(ksp,i,rnorm);
 76:     }

 78:     VecAXPY(&scale,z,x);    /*   x  <- x + scale z */
 79:     if (ksp->normtype != KSP_NO_NORM) {
 80:       ierr       = PetscObjectTakeAccess(ksp);
 81:       ksp->rnorm = rnorm;
 82:       ierr       = PetscObjectGrantAccess(ksp);
 83:       KSPLogResidualHistory(ksp,rnorm);

 85:       (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
 86:       if (ksp->reason) break;
 87:     }
 88: 
 89:     KSP_MatMult(ksp,Amat,x,r);      /*   r  <- b - Ax      */
 90:     VecAYPX(&mone,b,r);
 91:   }
 92:   if (!ksp->reason) {
 93:     ksp->reason = KSP_DIVERGED_ITS;
 94:     if (ksp->normtype != KSP_NO_NORM) {
 95:       if (ksp->normtype == KSP_UNPRECONDITIONED_NORM){
 96:         VecNorm(r,NORM_2,&rnorm);     /*   rnorm <- r'*r     */
 97:       } else {
 98:         KSP_PCApply(ksp,ksp->B,r,z);   /*   z <- B r          */
 99:         VecNorm(z,NORM_2,&rnorm);     /*   rnorm <- z'*z     */
100:       }
101:     }
102:     PetscObjectTakeAccess(ksp);
103:     ksp->rnorm = rnorm;
104:     PetscObjectGrantAccess(ksp);
105:     KSPLogResidualHistory(ksp,rnorm);
106:     KSPMonitor(ksp,i,rnorm);
107:     i--;
108:   } else if (!ksp->reason) {
109:     ksp->reason = KSP_DIVERGED_ITS;
110:     i--;
111:   }

113:   if (its) *its = i+1;
114:   return(0);
115: }

117: int KSPView_Richardson(KSP ksp,PetscViewer viewer)
118: {
119:   KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
120:   int            ierr;
121:   PetscTruth     isascii;

124:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
125:   if (isascii) {
126:     PetscViewerASCIIPrintf(viewer,"  Richardson: damping factor=%gn",richardsonP->scale);
127:   } else {
128:     SETERRQ1(1,"Viewer type %s not supported for KSP Richardson",((PetscObject)viewer)->type_name);
129:   }
130:   return(0);
131: }

133: int KSPSetFromOptions_Richardson(KSP ksp)
134: {
135:   KSP_Richardson *rich = (KSP_Richardson*)ksp->data;
136:   int            ierr;
137:   PetscReal      tmp;
138:   PetscTruth     flg;

141:   PetscOptionsHead("KSP Richardson Options");
142:     PetscOptionsReal("-ksp_richardson_scale","damping factor","KSPRichardsonSetScale",rich->scale,&tmp,&flg);
143:     if (flg) { KSPRichardsonSetScale(ksp,tmp); }
144:   PetscOptionsTail();
145:   return(0);
146: }

148: EXTERN_C_BEGIN
149: int KSPRichardsonSetScale_Richardson(KSP ksp,PetscReal scale)
150: {
151:   KSP_Richardson *richardsonP;

154:   richardsonP = (KSP_Richardson*)ksp->data;
155:   richardsonP->scale = scale;
156:   return(0);
157: }
158: EXTERN_C_END

160: EXTERN_C_BEGIN
161: int KSPCreate_Richardson(KSP ksp)
162: {
163:   int            ierr;
164:   KSP_Richardson *richardsonP;

167:   PetscNew(KSP_Richardson,&richardsonP);
168:   PetscLogObjectMemory(ksp,sizeof(KSP_Richardson));
169:   ksp->data                        = (void*)richardsonP;
170:   richardsonP->scale               = 1.0;
171:   ksp->ops->setup                  = KSPSetUp_Richardson;
172:   ksp->ops->solve                  = KSPSolve_Richardson;
173:   ksp->ops->destroy                = KSPDefaultDestroy;
174:   ksp->ops->buildsolution          = KSPDefaultBuildSolution;
175:   ksp->ops->buildresidual          = KSPDefaultBuildResidual;
176:   ksp->ops->view                   = KSPView_Richardson;
177:   ksp->ops->setfromoptions         = KSPSetFromOptions_Richardson;

179:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetScale_C",
180:                                     "KSPRichardsonSetScale_Richardson",
181:                                     KSPRichardsonSetScale_Richardson);
182:   return(0);
183: }
184: EXTERN_C_END