Actual source code: rich.c

  1: /*$Id: rich.c,v 1.101 2001/03/23 23:23:40 balay 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:   Scalar          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,pres,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:     *its = maxit;
 46:     PCApplyRichardson(ksp->B,b,x,r,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;
 53:   pres    = ksp->use_pres;

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

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

 67:     if (ksp->calc_res && !ksp->avoidnorms && !pres) {
 68:       VecNorm(r,NORM_2,&rnorm); /*   rnorm <- r'*r     */
 69:       KSPMonitor(ksp,i,rnorm);
 70:     }

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

 74:     if (ksp->calc_res && !ksp->avoidnorms && pres) {
 75:       VecNorm(z,NORM_2,&rnorm); /*   rnorm <- z'*z     */
 76:       KSPMonitor(ksp,i,rnorm);
 77:     }

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

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

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

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

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

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

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

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

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

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

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

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