Actual source code: anderson.c

  1: #include <../src/snes/impls/ngmres/snesngmres.h>

  3: static PetscErrorCode SNESSetFromOptions_Anderson(PetscOptionItems *PetscOptionsObject,SNES snes)
  4: {
  5:   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
  7:   PetscBool      monitor = PETSC_FALSE;

 10:   PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");
 11:   PetscOptionsInt("-snes_anderson_m",            "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);
 12:   PetscOptionsReal("-snes_anderson_beta",        "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);
 13:   PetscOptionsInt("-snes_anderson_restart",      "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);
 14:   PetscOptionsInt("-snes_anderson_restart_it",   "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);
 15:   PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,(PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);
 16:   PetscOptionsBool("-snes_anderson_monitor",     "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&monitor,NULL);
 17:   if (monitor) {
 18:     ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
 19:   }
 20:   PetscOptionsTail();
 21:   return(0);
 22: }

 24: static PetscErrorCode SNESSolve_Anderson(SNES snes)
 25: {
 26:   SNES_NGMRES         *ngmres = (SNES_NGMRES*) snes->data;
 27:   /* present solution, residual, and preconditioned residual */
 28:   Vec                 X,F,B,D;
 29:   /* candidate linear combination answers */
 30:   Vec                 XA,FA,XM,FM;

 32:   /* coefficients and RHS to the minimization problem */
 33:   PetscReal           fnorm,fMnorm,fAnorm;
 34:   PetscReal           xnorm,ynorm;
 35:   PetscReal           dnorm,dminnorm=0.0,fminnorm;
 36:   PetscInt            restart_count=0;
 37:   PetscInt            k,k_restart,l,ivec;
 38:   PetscBool           selectRestart;
 39:   SNESConvergedReason reason;
 40:   PetscErrorCode      ierr;

 43:   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);

 45:   PetscCitationsRegister(SNESCitation,&SNEScite);
 46:   /* variable initialization */
 47:   snes->reason = SNES_CONVERGED_ITERATING;
 48:   X            = snes->vec_sol;
 49:   F            = snes->vec_func;
 50:   B            = snes->vec_rhs;
 51:   XA           = snes->vec_sol_update;
 52:   FA           = snes->work[0];
 53:   D            = snes->work[1];

 55:   /* work for the line search */
 56:   XM = snes->work[3];
 57:   FM = snes->work[4];

 59:   PetscObjectSAWsTakeAccess((PetscObject)snes);
 60:   snes->iter = 0;
 61:   snes->norm = 0.;
 62:   PetscObjectSAWsGrantAccess((PetscObject)snes);

 64:   /* initialization */

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

 68:   if (snes->npc && snes->npcside== PC_LEFT) {
 69:     SNESApplyNPC(snes,X,NULL,F);
 70:     SNESGetConvergedReason(snes->npc,&reason);
 71:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
 72:       snes->reason = SNES_DIVERGED_INNER;
 73:       return(0);
 74:     }
 75:     VecNorm(F,NORM_2,&fnorm);
 76:   } else {
 77:     if (!snes->vec_func_init_set) {
 78:       SNESComputeFunction(snes,X,F);
 79:     } else snes->vec_func_init_set = PETSC_FALSE;

 81:     VecNorm(F,NORM_2,&fnorm);
 82:     SNESCheckFunctionNorm(snes,fnorm);
 83:   }
 84:   fminnorm = fnorm;

 86:   PetscObjectSAWsTakeAccess((PetscObject)snes);
 87:   snes->norm = fnorm;
 88:   PetscObjectSAWsGrantAccess((PetscObject)snes);
 89:   SNESLogConvergenceHistory(snes,fnorm,0);
 90:   SNESMonitor(snes,0,fnorm);
 91:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
 92:   if (snes->reason) return(0);

 94:   k_restart = 0;
 95:   l         = 0;
 96:   ivec      = 0;
 97:   for (k=1; k < snes->max_its+1; k++) {
 98:     /* select which vector of the stored subspace will be updated */
 99:     if (snes->npc && snes->npcside== PC_RIGHT) {
100:       VecCopy(X,XM);
101:       SNESSetInitialFunction(snes->npc,F);

103:       PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0);
104:       SNESSolve(snes->npc,B,XM);
105:       PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0);

107:       SNESGetConvergedReason(snes->npc,&reason);
108:       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
109:         snes->reason = SNES_DIVERGED_INNER;
110:         return(0);
111:       }
112:       SNESGetNPCFunction(snes,FM,&fMnorm);
113:       if (ngmres->andersonBeta != 1.0) {
114:         VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);
115:       }
116:     } else {
117:       VecCopy(F,FM);
118:       VecCopy(X,XM);
119:       VecAXPY(XM,-ngmres->andersonBeta,FM);
120:       fMnorm = fnorm;
121:     }

123:     SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);
124:     ivec = k_restart % ngmres->msize;
125:     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
126:       SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
127:       SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);
128:       /* if the restart conditions persist for more than restart_it iterations, restart. */
129:       if (selectRestart) restart_count++;
130:       else restart_count = 0;
131:     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
132:       SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
133:       if (k_restart > ngmres->restart_periodic) {
134:         if (ngmres->monitor) {PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);}
135:         restart_count = ngmres->restart_it;
136:       }
137:     } else {
138:       SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
139:     }
140:     /* restart after restart conditions have persisted for a fixed number of iterations */
141:     if (restart_count >= ngmres->restart_it) {
142:       if (ngmres->monitor) {
143:         PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);
144:       }
145:       restart_count = 0;
146:       k_restart     = 0;
147:       l             = 0;
148:       ivec          = 0;
149:     } else {
150:       if (l < ngmres->msize) l++;
151:       k_restart++;
152:       SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);
153:     }

155:     fnorm = fAnorm;
156:     if (fminnorm > fnorm) fminnorm = fnorm;

158:     VecCopy(XA,X);
159:     VecCopy(FA,F);

161:     PetscObjectSAWsTakeAccess((PetscObject)snes);
162:     snes->iter = k;
163:     snes->norm = fnorm;
164:     snes->xnorm = xnorm;
165:     snes->ynorm = ynorm;
166:     PetscObjectSAWsGrantAccess((PetscObject)snes);
167:     SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
168:     SNESMonitor(snes,snes->iter,snes->norm);
169:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
170:     if (snes->reason) return(0);
171:   }
172:   snes->reason = SNES_DIVERGED_MAX_IT;
173:   return(0);
174: }

176: /*MC
177:   SNESANDERSON - Anderson Mixing method.

179:    Level: beginner

181:    Options Database:
182: +  -snes_anderson_m                - Number of stored previous solutions and residuals
183: .  -snes_anderson_beta             - Anderson mixing parameter
184: .  -snes_anderson_restart_type     - Type of restart (see SNESNGMRES)
185: .  -snes_anderson_restart_it       - Number of iterations of restart conditions before restart
186: .  -snes_anderson_restart          - Number of iterations before periodic restart
187: -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration

189:    Notes:

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

194:    Very similar to the SNESNGMRES algorithm.

196:    References:
197: +  1. -  D. G. Anderson. Iterative procedures for nonlinear integral equations.
198:     J. Assoc. Comput. Mach., 12, 1965."
199: -  2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
200:    SIAM Review, 57(4), 2015

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

205: PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
206: {
207:   SNES_NGMRES    *ngmres;
209:   SNESLineSearch linesearch;

212:   snes->ops->destroy        = SNESDestroy_NGMRES;
213:   snes->ops->setup          = SNESSetUp_NGMRES;
214:   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
215:   snes->ops->view           = SNESView_NGMRES;
216:   snes->ops->solve          = SNESSolve_Anderson;
217:   snes->ops->reset          = SNESReset_NGMRES;

219:   snes->usesnpc = PETSC_TRUE;
220:   snes->usesksp = PETSC_FALSE;
221:   snes->npcside = PC_RIGHT;

223:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

225:   PetscNewLog(snes,&ngmres);
226:   snes->data    = (void*) ngmres;
227:   ngmres->msize = 30;

229:   if (!snes->tolerancesset) {
230:     snes->max_funcs = 30000;
231:     snes->max_its   = 10000;
232:   }

234:   SNESGetLineSearch(snes,&linesearch);
235:   if (!((PetscObject)linesearch)->type_name) {
236:     SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
237:   }

239:   ngmres->additive_linesearch = NULL;
240:   ngmres->approxfunc       = PETSC_FALSE;
241:   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
242:   ngmres->restart_it       = 2;
243:   ngmres->restart_periodic = 30;
244:   ngmres->gammaA           = 2.0;
245:   ngmres->gammaC           = 2.0;
246:   ngmres->deltaB           = 0.9;
247:   ngmres->epsilonB         = 0.1;

249:   ngmres->andersonBeta = 1.0;
250:   return(0);
251: }