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