1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 2 3 static PetscErrorCode SNESSetFromOptions_Anderson(SNES snes, PetscOptionItems *PetscOptionsObject) 4 { 5 SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 6 PetscBool monitor = PETSC_FALSE; 7 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 } 20 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; 28 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; 37 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); 40 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->vec_sol_update; 48 FA = snes->work[0]; 49 D = snes->work[1]; 50 51 /* work for the line search */ 52 XM = snes->work[3]; 53 FM = snes->work[4]; 54 55 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 56 snes->iter = 0; 57 snes->norm = 0.; 58 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 59 60 /* initialization */ 61 62 /* r = F(x) */ 63 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; 76 77 PetscCall(VecNorm(F, NORM_2, &fnorm)); 78 SNESCheckFunctionNorm(snes, fnorm); 79 } 80 fminnorm = fnorm; 81 82 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 83 snes->norm = fnorm; 84 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 85 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 86 PetscCall(SNESMonitor(snes, 0, fnorm)); 87 PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 88 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 89 90 k_restart = 0; 91 l = 0; 92 ivec = 0; 93 for (k = 1; k < snes->max_its + 1; k++) { 94 /* select which vector of the stored subspace will be updated */ 95 if (snes->npc && snes->npcside == PC_RIGHT) { 96 PetscCall(VecCopy(X, XM)); 97 PetscCall(SNESSetInitialFunction(snes->npc, F)); 98 99 PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, XM, B, 0)); 100 PetscCall(SNESSolve(snes->npc, B, XM)); 101 PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, XM, B, 0)); 102 103 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 104 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 105 snes->reason = SNES_DIVERGED_INNER; 106 PetscFunctionReturn(PETSC_SUCCESS); 107 } 108 PetscCall(SNESGetNPCFunction(snes, FM, &fMnorm)); 109 if (ngmres->andersonBeta != 1.0) PetscCall(VecAXPBY(XM, (1.0 - ngmres->andersonBeta), ngmres->andersonBeta, X)); 110 } else { 111 PetscCall(VecCopy(F, FM)); 112 PetscCall(VecCopy(X, XM)); 113 PetscCall(VecAXPY(XM, -ngmres->andersonBeta, FM)); 114 fMnorm = fnorm; 115 } 116 117 PetscCall(SNESNGMRESFormCombinedSolution_Private(snes, ivec, l, XM, FM, fMnorm, X, XA, FA)); 118 ivec = k_restart % ngmres->msize; 119 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 120 PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, &dnorm, &dminnorm, NULL, NULL, NULL, &xnorm, &fAnorm, &ynorm)); 121 PetscCall(SNESNGMRESSelectRestart_Private(snes, l, fMnorm, fnorm, dnorm, fminnorm, dminnorm, &selectRestart)); 122 /* if the restart conditions persist for more than restart_it iterations, restart. */ 123 if (selectRestart) restart_count++; 124 else restart_count = 0; 125 } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 126 PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, NULL, NULL, NULL, &xnorm, &fAnorm, &ynorm)); 127 if (k_restart > ngmres->restart_periodic) { 128 if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %" PetscInt_FMT " iterations\n", k_restart)); 129 restart_count = ngmres->restart_it; 130 } 131 } else { 132 PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, NULL, NULL, NULL, &xnorm, &fAnorm, &ynorm)); 133 } 134 /* restart after restart conditions have persisted for a fixed number of iterations */ 135 if (restart_count >= ngmres->restart_it) { 136 if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %" PetscInt_FMT "\n", k_restart)); 137 restart_count = 0; 138 k_restart = 0; 139 l = 0; 140 ivec = 0; 141 } else { 142 if (l < ngmres->msize) l++; 143 k_restart++; 144 PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, FM, fnorm, XM)); 145 } 146 147 fnorm = fAnorm; 148 if (fminnorm > fnorm) fminnorm = fnorm; 149 150 PetscCall(VecCopy(XA, X)); 151 PetscCall(VecCopy(FA, F)); 152 153 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 154 snes->iter = k; 155 snes->norm = fnorm; 156 snes->xnorm = xnorm; 157 snes->ynorm = ynorm; 158 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 159 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter)); 160 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 161 PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 162 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 163 } 164 snes->reason = SNES_DIVERGED_MAX_IT; 165 PetscFunctionReturn(PETSC_SUCCESS); 166 } 167 168 /*MC 169 SNESANDERSON - Anderson Mixing nonlinear solver 170 171 Level: beginner 172 173 Options Database Keys: 174 + -snes_anderson_m - Number of stored previous solutions and residuals 175 . -snes_anderson_beta - Anderson mixing parameter 176 . -snes_anderson_restart_type - Type of restart (see SNESNGMRES) 177 . -snes_anderson_restart_it - Number of iterations of restart conditions before restart 178 . -snes_anderson_restart - Number of iterations before periodic restart 179 - -snes_anderson_monitor - Prints relevant information about the ngmres iteration 180 181 Notes: 182 The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized 183 optimization problem at each iteration. 184 185 Very similar to the `SNESNGMRES` algorithm. 186 187 References: 188 + * - D. G. Anderson. Iterative procedures for nonlinear integral equations. 189 J. Assoc. Comput. Mach., 12, 1965." 190 - * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers", 191 SIAM Review, 57(4), 2015 192 193 .seealso: `SNESNGMRES`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType` 194 M*/ 195 196 PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes) 197 { 198 SNES_NGMRES *ngmres; 199 SNESLineSearch linesearch; 200 201 PetscFunctionBegin; 202 snes->ops->destroy = SNESDestroy_NGMRES; 203 snes->ops->setup = SNESSetUp_NGMRES; 204 snes->ops->setfromoptions = SNESSetFromOptions_Anderson; 205 snes->ops->view = SNESView_NGMRES; 206 snes->ops->solve = SNESSolve_Anderson; 207 snes->ops->reset = SNESReset_NGMRES; 208 209 snes->usesnpc = PETSC_TRUE; 210 snes->usesksp = PETSC_FALSE; 211 snes->npcside = PC_RIGHT; 212 213 snes->alwayscomputesfinalresidual = PETSC_TRUE; 214 215 PetscCall(PetscNew(&ngmres)); 216 snes->data = (void *)ngmres; 217 ngmres->msize = 30; 218 219 if (!snes->tolerancesset) { 220 snes->max_funcs = 30000; 221 snes->max_its = 10000; 222 } 223 224 PetscCall(SNESGetLineSearch(snes, &linesearch)); 225 if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 226 227 ngmres->additive_linesearch = NULL; 228 ngmres->approxfunc = PETSC_FALSE; 229 ngmres->restart_type = SNES_NGMRES_RESTART_NONE; 230 ngmres->restart_it = 2; 231 ngmres->restart_periodic = 30; 232 ngmres->gammaA = 2.0; 233 ngmres->gammaC = 2.0; 234 ngmres->deltaB = 0.9; 235 ngmres->epsilonB = 0.1; 236 237 ngmres->andersonBeta = 1.0; 238 PetscFunctionReturn(PETSC_SUCCESS); 239 } 240