Actual source code: anderson.c

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

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

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

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

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

 38:   PetscFunctionBegin;
 39:   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);

 41:   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
 42:   /* variable initialization */
 43:   snes->reason = SNES_CONVERGED_ITERATING;
 44:   X            = snes->vec_sol;
 45:   F            = snes->vec_func;
 46:   B            = snes->vec_rhs;
 47:   XA           = snes->work[2];
 48:   FA           = snes->work[0];
 49:   D            = snes->work[1];

 51:   /* work for the line search */
 52:   XM = snes->work[3];
 53:   FM = snes->work[4];

 55:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
 56:   snes->iter = 0;
 57:   snes->norm = 0.;
 58:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));

 60:   /* initialization */

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

 64:   if (snes->npc && snes->npcside == PC_LEFT) {
 65:     PetscCall(SNESApplyNPC(snes, X, NULL, F));
 66:     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
 67:     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
 68:       snes->reason = SNES_DIVERGED_INNER;
 69:       PetscFunctionReturn(PETSC_SUCCESS);
 70:     }
 71:     PetscCall(VecNorm(F, NORM_2, &fnorm));
 72:   } else {
 73:     if (!snes->vec_func_init_set) {
 74:       PetscCall(SNESComputeFunction(snes, X, F));
 75:     } else snes->vec_func_init_set = PETSC_FALSE;

 77:     PetscCall(VecNorm(F, NORM_2, &fnorm));
 78:     SNESCheckFunctionNorm(snes, fnorm);
 79:   }
 80:   fminnorm = fnorm;

 82:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
 83:   snes->norm = fnorm;
 84:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
 85:   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
 86:   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
 87:   PetscCall(SNESMonitor(snes, 0, fnorm));
 88:   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
 89:   PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, F, fnorm, X));

 91:   k_restart = 0;
 92:   l         = 0;
 93:   ivec      = 0;
 94:   for (k = 1; k < snes->max_its + 1; k++) {
 95:     /* Call general purpose update function */
 96:     PetscTryTypeMethod(snes, update, snes->iter);

 98:     /* select which vector of the stored subspace will be updated */
 99:     if (snes->npc && snes->npcside == PC_RIGHT) {
100:       PetscCall(VecCopy(X, XM));
101:       PetscCall(SNESSetInitialFunction(snes->npc, F));

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

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

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

151:     fnorm = fAnorm;
152:     if (fminnorm > fnorm) fminnorm = fnorm;

154:     PetscCall(VecCopy(XA, X));
155:     PetscCall(VecCopy(FA, F));

157:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
158:     snes->iter  = k;
159:     snes->norm  = fnorm;
160:     snes->xnorm = xnorm;
161:     snes->ynorm = ynorm;
162:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
163:     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
164:     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
165:     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
166:     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
167:   }
168:   snes->reason = SNES_DIVERGED_MAX_IT;
169:   PetscFunctionReturn(PETSC_SUCCESS);
170: }

172: /*MC
173:   SNESANDERSON - Anderson Mixing nonlinear solver {cite}`anderson1965`, {cite}`bruneknepleysmithtu15`

175:    Level: beginner

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

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

189:    Very similar to the `SNESNGMRES` algorithm.

191: .seealso: [](ch_snes), `SNESNGMRES`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`
192: M*/

194: PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
195: {
196:   SNES_NGMRES   *ngmres;
197:   SNESLineSearch linesearch;

199:   PetscFunctionBegin;
200:   snes->ops->destroy        = SNESDestroy_NGMRES;
201:   snes->ops->setup          = SNESSetUp_NGMRES;
202:   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
203:   snes->ops->view           = SNESView_NGMRES;
204:   snes->ops->solve          = SNESSolve_Anderson;
205:   snes->ops->reset          = SNESReset_NGMRES;

207:   snes->usesnpc = PETSC_TRUE;
208:   snes->usesksp = PETSC_FALSE;
209:   snes->npcside = PC_RIGHT;

211:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

213:   PetscCall(PetscNew(&ngmres));
214:   snes->data    = (void *)ngmres;
215:   ngmres->msize = 30;

217:   if (!snes->tolerancesset) {
218:     snes->max_funcs = 30000;
219:     snes->max_its   = 10000;
220:   }

222:   PetscCall(SNESGetLineSearch(snes, &linesearch));
223:   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));

225:   ngmres->additive_linesearch = NULL;
226:   ngmres->approxfunc          = PETSC_FALSE;
227:   ngmres->restart_type        = SNES_NGMRES_RESTART_NONE;
228:   ngmres->restart_it          = 2;
229:   ngmres->restart_periodic    = 30;
230:   ngmres->gammaA              = 2.0;
231:   ngmres->gammaC              = 2.0;
232:   ngmres->deltaB              = 0.9;
233:   ngmres->epsilonB            = 0.1;

235:   ngmres->andersonBeta = 1.0;
236:   PetscFunctionReturn(PETSC_SUCCESS);
237: }