1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 2 3 static PetscErrorCode SNESSetFromOptions_Anderson(PetscOptionItems *PetscOptionsObject,SNES snes) 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) { 17 ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); 18 } 19 PetscOptionsHeadEnd(); 20 PetscFunctionReturn(0); 21 } 22 23 static PetscErrorCode SNESSolve_Anderson(SNES snes) 24 { 25 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 26 /* present solution, residual, and preconditioned residual */ 27 Vec X,F,B,D; 28 /* candidate linear combination answers */ 29 Vec XA,FA,XM,FM; 30 31 /* coefficients and RHS to the minimization problem */ 32 PetscReal fnorm,fMnorm,fAnorm; 33 PetscReal xnorm,ynorm; 34 PetscReal dnorm,dminnorm=0.0,fminnorm; 35 PetscInt restart_count=0; 36 PetscInt k,k_restart,l,ivec; 37 PetscBool selectRestart; 38 SNESConvergedReason reason; 39 40 PetscFunctionBegin; 41 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); 42 43 PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 44 /* variable initialization */ 45 snes->reason = SNES_CONVERGED_ITERATING; 46 X = snes->vec_sol; 47 F = snes->vec_func; 48 B = snes->vec_rhs; 49 XA = snes->vec_sol_update; 50 FA = snes->work[0]; 51 D = snes->work[1]; 52 53 /* work for the line search */ 54 XM = snes->work[3]; 55 FM = snes->work[4]; 56 57 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 58 snes->iter = 0; 59 snes->norm = 0.; 60 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 61 62 /* initialization */ 63 64 /* r = F(x) */ 65 66 if (snes->npc && snes->npcside== PC_LEFT) { 67 PetscCall(SNESApplyNPC(snes,X,NULL,F)); 68 PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 69 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 70 snes->reason = SNES_DIVERGED_INNER; 71 PetscFunctionReturn(0); 72 } 73 PetscCall(VecNorm(F,NORM_2,&fnorm)); 74 } else { 75 if (!snes->vec_func_init_set) { 76 PetscCall(SNESComputeFunction(snes,X,F)); 77 } else snes->vec_func_init_set = PETSC_FALSE; 78 79 PetscCall(VecNorm(F,NORM_2,&fnorm)); 80 SNESCheckFunctionNorm(snes,fnorm); 81 } 82 fminnorm = fnorm; 83 84 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 85 snes->norm = fnorm; 86 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 87 PetscCall(SNESLogConvergenceHistory(snes,fnorm,0)); 88 PetscCall(SNESMonitor(snes,0,fnorm)); 89 PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 90 if (snes->reason) PetscFunctionReturn(0); 91 92 k_restart = 0; 93 l = 0; 94 ivec = 0; 95 for (k=1; k < snes->max_its+1; k++) { 96 /* select which vector of the stored subspace will be updated */ 97 if (snes->npc && snes->npcside== PC_RIGHT) { 98 PetscCall(VecCopy(X,XM)); 99 PetscCall(SNESSetInitialFunction(snes->npc,F)); 100 101 PetscCall(PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0)); 102 PetscCall(SNESSolve(snes->npc,B,XM)); 103 PetscCall(PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0)); 104 105 PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 106 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 107 snes->reason = SNES_DIVERGED_INNER; 108 PetscFunctionReturn(0); 109 } 110 PetscCall(SNESGetNPCFunction(snes,FM,&fMnorm)); 111 if (ngmres->andersonBeta != 1.0) { 112 PetscCall(VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X)); 113 } 114 } else { 115 PetscCall(VecCopy(F,FM)); 116 PetscCall(VecCopy(X,XM)); 117 PetscCall(VecAXPY(XM,-ngmres->andersonBeta,FM)); 118 fMnorm = fnorm; 119 } 120 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,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) { 141 PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %" PetscInt_FMT "\n",k_restart)); 142 } 143 restart_count = 0; 144 k_restart = 0; 145 l = 0; 146 ivec = 0; 147 } else { 148 if (l < ngmres->msize) l++; 149 k_restart++; 150 PetscCall(SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM)); 151 } 152 153 fnorm = fAnorm; 154 if (fminnorm > fnorm) fminnorm = fnorm; 155 156 PetscCall(VecCopy(XA,X)); 157 PetscCall(VecCopy(FA,F)); 158 159 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 160 snes->iter = k; 161 snes->norm = fnorm; 162 snes->xnorm = xnorm; 163 snes->ynorm = ynorm; 164 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 165 PetscCall(SNESLogConvergenceHistory(snes,snes->norm,snes->iter)); 166 PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 167 PetscCall((*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP)); 168 if (snes->reason) PetscFunctionReturn(0); 169 } 170 snes->reason = SNES_DIVERGED_MAX_IT; 171 PetscFunctionReturn(0); 172 } 173 174 /*MC 175 SNESANDERSON - Anderson Mixing method. 176 177 Level: beginner 178 179 Options Database: 180 + -snes_anderson_m - Number of stored previous solutions and residuals 181 . -snes_anderson_beta - Anderson mixing parameter 182 . -snes_anderson_restart_type - Type of restart (see SNESNGMRES) 183 . -snes_anderson_restart_it - Number of iterations of restart conditions before restart 184 . -snes_anderson_restart - Number of iterations before periodic restart 185 - -snes_anderson_monitor - Prints relevant information about the ngmres iteration 186 187 Notes: 188 189 The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized 190 optimization problem at each iteration. 191 192 Very similar to the SNESNGMRES algorithm. 193 194 References: 195 + * - D. G. Anderson. Iterative procedures for nonlinear integral equations. 196 J. Assoc. Comput. Mach., 12, 1965." 197 - * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers", 198 SIAM Review, 57(4), 2015 199 200 .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 201 M*/ 202 203 PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes) 204 { 205 SNES_NGMRES *ngmres; 206 SNESLineSearch linesearch; 207 208 PetscFunctionBegin; 209 snes->ops->destroy = SNESDestroy_NGMRES; 210 snes->ops->setup = SNESSetUp_NGMRES; 211 snes->ops->setfromoptions = SNESSetFromOptions_Anderson; 212 snes->ops->view = SNESView_NGMRES; 213 snes->ops->solve = SNESSolve_Anderson; 214 snes->ops->reset = SNESReset_NGMRES; 215 216 snes->usesnpc = PETSC_TRUE; 217 snes->usesksp = PETSC_FALSE; 218 snes->npcside = PC_RIGHT; 219 220 snes->alwayscomputesfinalresidual = PETSC_TRUE; 221 222 PetscCall(PetscNewLog(snes,&ngmres)); 223 snes->data = (void*) ngmres; 224 ngmres->msize = 30; 225 226 if (!snes->tolerancesset) { 227 snes->max_funcs = 30000; 228 snes->max_its = 10000; 229 } 230 231 PetscCall(SNESGetLineSearch(snes,&linesearch)); 232 if (!((PetscObject)linesearch)->type_name) { 233 PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC)); 234 } 235 236 ngmres->additive_linesearch = NULL; 237 ngmres->approxfunc = PETSC_FALSE; 238 ngmres->restart_type = SNES_NGMRES_RESTART_NONE; 239 ngmres->restart_it = 2; 240 ngmres->restart_periodic = 30; 241 ngmres->gammaA = 2.0; 242 ngmres->gammaC = 2.0; 243 ngmres->deltaB = 0.9; 244 ngmres->epsilonB = 0.1; 245 246 ngmres->andersonBeta = 1.0; 247 PetscFunctionReturn(0); 248 } 249