Actual source code: anderson.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  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;
  8:   SNESLineSearch linesearch;

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

 30: static PetscErrorCode SNESSolve_Anderson(SNES snes)
 31: {
 32:   SNES_NGMRES         *ngmres = (SNES_NGMRES*) snes->data;
 33:   /* present solution, residual, and preconditioned residual */
 34:   Vec                 X,F,B,D;
 35:   /* candidate linear combination answers */
 36:   Vec                 XA,FA,XM,FM;

 38:   /* coefficients and RHS to the minimization problem */
 39:   PetscReal           fnorm,fMnorm,fAnorm;
 40:   PetscReal           xnorm,ynorm;
 41:   PetscReal           dnorm,dminnorm=0.0,fminnorm;
 42:   PetscInt            restart_count=0;
 43:   PetscInt            k,k_restart,l,ivec;
 44:   PetscBool           selectRestart;
 45:   SNESConvergedReason reason;
 46:   PetscErrorCode      ierr;

 49:   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);

 51:   PetscCitationsRegister(SNESCitation,&SNEScite);
 52:   /* variable initialization */
 53:   snes->reason = SNES_CONVERGED_ITERATING;
 54:   X            = snes->vec_sol;
 55:   F            = snes->vec_func;
 56:   B            = snes->vec_rhs;
 57:   XA           = snes->vec_sol_update;
 58:   FA           = snes->work[0];
 59:   D            = snes->work[1];

 61:   /* work for the line search */
 62:   XM = snes->work[3];
 63:   FM = snes->work[4];

 65:   PetscObjectSAWsTakeAccess((PetscObject)snes);
 66:   snes->iter = 0;
 67:   snes->norm = 0.;
 68:   PetscObjectSAWsGrantAccess((PetscObject)snes);

 70:   /* initialization */

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

 74:   if (snes->npc && snes->npcside== PC_LEFT) {
 75:     SNESApplyNPC(snes,X,NULL,F);
 76:     SNESGetConvergedReason(snes->npc,&reason);
 77:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
 78:       snes->reason = SNES_DIVERGED_INNER;
 79:       return(0);
 80:     }
 81:     VecNorm(F,NORM_2,&fnorm);
 82:   } else {
 83:     if (!snes->vec_func_init_set) {
 84:       SNESComputeFunction(snes,X,F);
 85:     } else snes->vec_func_init_set = PETSC_FALSE;

 87:     VecNorm(F,NORM_2,&fnorm);
 88:     SNESCheckFunctionNorm(snes,fnorm);
 89:   }
 90:   fminnorm = fnorm;

 92:   PetscObjectSAWsTakeAccess((PetscObject)snes);
 93:   snes->norm = fnorm;
 94:   PetscObjectSAWsGrantAccess((PetscObject)snes);
 95:   SNESLogConvergenceHistory(snes,fnorm,0);
 96:   SNESMonitor(snes,0,fnorm);
 97:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
 98:   if (snes->reason) return(0);

100:   k_restart = 0;
101:   l         = 0;
102:   ivec      = 0;
103:   for (k=1; k < snes->max_its+1; k++) {
104:     /* select which vector of the stored subspace will be updated */
105:     if (snes->npc && snes->npcside== PC_RIGHT) {
106:       VecCopy(X,XM);
107:       SNESSetInitialFunction(snes->npc,F);

109:       PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0);
110:       SNESSolve(snes->npc,B,XM);
111:       PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0);

113:       SNESGetConvergedReason(snes->npc,&reason);
114:       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
115:         snes->reason = SNES_DIVERGED_INNER;
116:         return(0);
117:       }
118:       SNESGetNPCFunction(snes,FM,&fMnorm);
119:       if (ngmres->andersonBeta != 1.0) {
120:         VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);
121:       }
122:     } else {
123:       VecCopy(F,FM);
124:       VecCopy(X,XM);
125:       VecAXPY(XM,-ngmres->andersonBeta,FM);
126:       fMnorm = fnorm;
127:     }

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

161:     fnorm = fAnorm;
162:     if (fminnorm > fnorm) fminnorm = fnorm;

164:     VecCopy(XA,X);
165:     VecCopy(FA,F);

167:     PetscObjectSAWsTakeAccess((PetscObject)snes);
168:     snes->iter = k;
169:     snes->norm = fnorm;
170:     PetscObjectSAWsGrantAccess((PetscObject)snes);
171:     SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
172:     SNESMonitor(snes,snes->iter,snes->norm);
173:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
174:     if (snes->reason) return(0);
175:   }
176:   snes->reason = SNES_DIVERGED_MAX_IT;
177:   return(0);
178: }

180: /*MC
181:   SNESANDERSON - Anderson Mixing method.

183:    Level: beginner

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

193:    Notes:

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

198:    Very similar to the SNESNGMRES algorithm.

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


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

210: PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
211: {
212:   SNES_NGMRES    *ngmres;

216:   snes->ops->destroy        = SNESDestroy_NGMRES;
217:   snes->ops->setup          = SNESSetUp_NGMRES;
218:   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
219:   snes->ops->view           = SNESView_NGMRES;
220:   snes->ops->solve          = SNESSolve_Anderson;
221:   snes->ops->reset          = SNESReset_NGMRES;

223:   snes->usesnpc = PETSC_TRUE;
224:   snes->usesksp = PETSC_FALSE;
225:   snes->npcside = PC_RIGHT;

227:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

229:   PetscNewLog(snes,&ngmres);
230:   snes->data    = (void*) ngmres;
231:   ngmres->msize = 30;

233:   if (!snes->tolerancesset) {
234:     snes->max_funcs = 30000;
235:     snes->max_its   = 10000;
236:   }

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

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