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