Actual source code: anderson.c

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/

  3: extern PetscErrorCode SNESDestroy_NGMRES(SNES);
  4: extern PetscErrorCode SNESReset_NGMRES(SNES);
  5: extern PetscErrorCode SNESSetUp_NGMRES(SNES);
  6: extern PetscErrorCode SNESView_NGMRES(SNES,PetscViewer);

  8: PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];

 12: PetscErrorCode SNESSetFromOptions_Anderson(SNES snes)
 13: {
 14:   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
 16:   PetscBool      debug;
 17:   SNESLineSearch linesearch;

 20:   PetscOptionsHead("SNES NGMRES options");
 21:   PetscOptionsInt("-snes_anderson_m",            "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);
 22:   PetscOptionsReal("-snes_anderson_beta",        "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);
 23:   PetscOptionsBool("-snes_anderson_monitor",     "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);
 24:   PetscOptionsInt("-snes_anderson_restart",      "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);
 25:   PetscOptionsInt("-snes_anderson_restart_it",   "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);
 26:   PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,
 27:                           (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);
 28:   if (debug) {
 29:     ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
 30:   }
 31:   PetscOptionsTail();
 32:   /* set the default type of the line search if the user hasn't already. */
 33:   if (!snes->linesearch) {
 34:     SNESGetLineSearch(snes,&linesearch);
 35:     SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
 36:   }
 37:   return(0);
 38: }

 42: PetscErrorCode SNESSolve_Anderson(SNES snes)
 43: {
 44:   SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
 45:   /* present solution, residual, and preconditioned residual */
 46:   Vec X,F,B,D;

 48:   /* candidate linear combination answers */
 49:   Vec XA,FA,XM,FM;

 51:   /* coefficients and RHS to the minimization problem */
 52:   PetscReal fnorm,fMnorm,fAnorm;
 53:   PetscReal  xnorm,ynorm;
 54:   PetscReal dnorm,dminnorm=0.0,fminnorm;
 55:   PetscInt  restart_count=0;
 56:   PetscInt  k,k_restart,l,ivec;

 58:   PetscBool selectRestart;

 60:   SNESConvergedReason reason;
 61:   PetscErrorCode      ierr;

 64:   /* variable initialization */
 65:   snes->reason = SNES_CONVERGED_ITERATING;
 66:   X            = snes->vec_sol;
 67:   F            = snes->vec_func;
 68:   B            = snes->vec_rhs;
 69:   XA           = snes->vec_sol_update;
 70:   FA           = snes->work[0];
 71:   D            = snes->work[1];

 73:   /* work for the line search */
 74:   XM = snes->work[3];
 75:   FM = snes->work[4];

 77:   PetscObjectSAWsTakeAccess((PetscObject)snes);
 78:   snes->iter = 0;
 79:   snes->norm = 0.;
 80:   PetscObjectSAWsGrantAccess((PetscObject)snes);

 82:   /* initialization */

 84:   /* r = F(x) */

 86:   if (snes->pc && snes->pcside == PC_LEFT) {
 87:     SNESApplyPC(snes,X,NULL,NULL,F);
 88:     SNESGetConvergedReason(snes->pc,&reason);
 89:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
 90:       snes->reason = SNES_DIVERGED_INNER;
 91:       return(0);
 92:     }
 93:     VecNorm(F,NORM_2,&fnorm);
 94:   } else {
 95:     if (!snes->vec_func_init_set) {
 96:       SNESComputeFunction(snes,X,F);
 97:       if (snes->domainerror) {
 98:         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
 99:         return(0);
100:       }
101:     } else snes->vec_func_init_set = PETSC_FALSE;

103:     VecNorm(F,NORM_2,&fnorm);
104:     if (PetscIsInfOrNanReal(fnorm)) {
105:       snes->reason = SNES_DIVERGED_FNORM_NAN;
106:       return(0);
107:     }
108:   }
109:   fminnorm = fnorm;

111:   PetscObjectSAWsTakeAccess((PetscObject)snes);
112:   snes->norm = fnorm;
113:   PetscObjectSAWsGrantAccess((PetscObject)snes);
114:   SNESLogConvergenceHistory(snes,fnorm,0);
115:   SNESMonitor(snes,0,fnorm);
116:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
117:   if (snes->reason) return(0);

119:   k_restart = 0;
120:   l         = 0;
121:   ivec      = 0;
122:   for (k=1; k < snes->max_its+1; k++) {
123:     /* select which vector of the stored subspace will be updated */
124:     if (snes->pc && snes->pcside == PC_RIGHT) {
125:       VecCopy(X,XM);
126:       SNESSetInitialFunction(snes->pc,F);

128:       PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);
129:       SNESSolve(snes->pc,B,XM);
130:       PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);

132:       SNESGetConvergedReason(snes->pc,&reason);
133:       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
134:         snes->reason = SNES_DIVERGED_INNER;
135:         return(0);
136:       }
137:       SNESGetPCFunction(snes,FM,&fMnorm);
138:       if (ngmres->andersonBeta != 1.0) {
139:         VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);
140:       }
141:     } else {
142:       VecCopy(F,FM);
143:       VecCopy(X,XM);
144:       VecAXPY(XM,-ngmres->andersonBeta,FM);
145:       fMnorm = fnorm;
146:     }

148:     SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);
149:     ivec = k_restart % ngmres->msize;
150:     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
151:       SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
152:       SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);
153:       /* if the restart conditions persist for more than restart_it iterations, restart. */
154:       if (selectRestart) restart_count++;
155:       else restart_count = 0;
156:     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
157:       SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
158:       if (k_restart > ngmres->restart_periodic) {
159:         if (ngmres->monitor) PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);
160:         restart_count = ngmres->restart_it;
161:       }
162:     } else {
163:       SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
164:     }
165:     /* restart after restart conditions have persisted for a fixed number of iterations */
166:     if (restart_count >= ngmres->restart_it) {
167:       if (ngmres->monitor) {
168:         PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);
169:       }
170:       restart_count = 0;
171:       k_restart     = 0;
172:       l             = 0;
173:       ivec          = 0;
174:     } else {
175:       if (l < ngmres->msize) l++;
176:       k_restart++;
177:       SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);
178:     }

180:     fnorm = fAnorm;
181:     if (fminnorm > fnorm) fminnorm = fnorm;

183:     VecCopy(XA,X);
184:     VecCopy(FA,F);

186:     PetscObjectSAWsTakeAccess((PetscObject)snes);
187:     snes->iter = k;
188:     snes->norm = fnorm;
189:     PetscObjectSAWsGrantAccess((PetscObject)snes);
190:     SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
191:     SNESMonitor(snes,snes->iter,snes->norm);
192:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
193:     if (snes->reason) return(0);
194:   }
195:   snes->reason = SNES_DIVERGED_MAX_IT;
196:   return(0);
197: }

199: /*MC
200:   SNESAnderson - Anderson Mixing method.

202:    Level: beginner

204:    Options Database:
205: +  -snes_anderson_m                - Number of stored previous solutions and residuals
206: .  -snes_anderson_beta             - Relaxation parameter; X_{update} = X + \beta F
207: .  -snes_anderson_restart_type     - Type of restart (see SNESNGMRES)
208: .  -snes_anderson_restart_it       - Number of iterations of restart conditions before restart
209: .  -snes_anderson_restart          - Number of iterations before periodic restart
210: -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration

212:    Notes:

214:    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
215:    optimization problem at each iteration.

217:    References:

219:     "D. G. Anderson. Iterative procedures for nonlinear integral equations.
220:     J. Assoc. Comput. Mach., 12:547-560, 1965."

222: .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
223: M*/

227: PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
228: {
229:   SNES_NGMRES    *ngmres;

233:   snes->ops->destroy        = SNESDestroy_NGMRES;
234:   snes->ops->setup          = SNESSetUp_NGMRES;
235:   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
236:   snes->ops->view           = SNESView_NGMRES;
237:   snes->ops->solve          = SNESSolve_Anderson;
238:   snes->ops->reset          = SNESReset_NGMRES;

240:   snes->usespc  = PETSC_TRUE;
241:   snes->usesksp = PETSC_FALSE;
242:   snes->pcside  = PC_RIGHT;

244:   PetscNewLog(snes,&ngmres);
245:   snes->data    = (void*) ngmres;
246:   ngmres->msize = 30;

248:   if (!snes->tolerancesset) {
249:     snes->max_funcs = 30000;
250:     snes->max_its   = 10000;
251:   }

253:   ngmres->additive_linesearch = NULL;
254:   ngmres->approxfunc       = PETSC_FALSE;
255:   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
256:   ngmres->restart_it       = 2;
257:   ngmres->restart_periodic = 30;
258:   ngmres->gammaA           = 2.0;
259:   ngmres->gammaC           = 2.0;
260:   ngmres->deltaB           = 0.9;
261:   ngmres->epsilonB         = 0.1;

263:   ngmres->andersonBeta = 1.0;
264:   return(0);
265: }