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