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