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(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("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 ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 62 /* variable initialization */ 63 snes->reason = SNES_CONVERGED_ITERATING; 64 X = snes->vec_sol; 65 F = snes->vec_func; 66 B = snes->vec_rhs; 67 XA = snes->vec_sol_update; 68 FA = snes->work[0]; 69 D = snes->work[1]; 70 71 /* work for the line search */ 72 XM = snes->work[3]; 73 FM = snes->work[4]; 74 75 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 76 snes->iter = 0; 77 snes->norm = 0.; 78 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 79 80 /* initialization */ 81 82 /* r = F(x) */ 83 84 if (snes->pc && snes->pcside == PC_LEFT) { 85 ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 86 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 87 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 88 snes->reason = SNES_DIVERGED_INNER; 89 PetscFunctionReturn(0); 90 } 91 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 92 } else { 93 if (!snes->vec_func_init_set) { 94 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 95 if (snes->domainerror) { 96 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 97 PetscFunctionReturn(0); 98 } 99 } else snes->vec_func_init_set = PETSC_FALSE; 100 101 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 102 if (PetscIsInfOrNanReal(fnorm)) { 103 snes->reason = SNES_DIVERGED_FNORM_NAN; 104 PetscFunctionReturn(0); 105 } 106 } 107 fminnorm = fnorm; 108 109 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 110 snes->norm = fnorm; 111 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 112 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 113 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 114 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 115 if (snes->reason) PetscFunctionReturn(0); 116 117 k_restart = 0; 118 l = 0; 119 ivec = 0; 120 for (k=1; k < snes->max_its+1; k++) { 121 /* select which vector of the stored subspace will be updated */ 122 if (snes->pc && snes->pcside == PC_RIGHT) { 123 ierr = VecCopy(X,XM);CHKERRQ(ierr); 124 ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr); 125 126 ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 127 ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr); 128 ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 129 130 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 131 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 132 snes->reason = SNES_DIVERGED_INNER; 133 PetscFunctionReturn(0); 134 } 135 ierr = SNESGetNPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr); 136 if (ngmres->andersonBeta != 1.0) { 137 VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr); 138 } 139 } else { 140 ierr = VecCopy(F,FM);CHKERRQ(ierr); 141 ierr = VecCopy(X,XM);CHKERRQ(ierr); 142 ierr = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr); 143 fMnorm = fnorm; 144 } 145 146 ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr); 147 ivec = k_restart % ngmres->msize; 148 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 149 ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr); 150 ierr = SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 151 /* if the restart conditions persist for more than restart_it iterations, restart. */ 152 if (selectRestart) restart_count++; 153 else restart_count = 0; 154 } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 155 ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr); 156 if (k_restart > ngmres->restart_periodic) { 157 if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr); 158 restart_count = ngmres->restart_it; 159 } 160 } else { 161 ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr); 162 } 163 /* restart after restart conditions have persisted for a fixed number of iterations */ 164 if (restart_count >= ngmres->restart_it) { 165 if (ngmres->monitor) { 166 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr); 167 } 168 restart_count = 0; 169 k_restart = 0; 170 l = 0; 171 ivec = 0; 172 } else { 173 if (l < ngmres->msize) l++; 174 k_restart++; 175 ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr); 176 } 177 178 fnorm = fAnorm; 179 if (fminnorm > fnorm) fminnorm = fnorm; 180 181 ierr = VecCopy(XA,X);CHKERRQ(ierr); 182 ierr = VecCopy(FA,F);CHKERRQ(ierr); 183 184 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 185 snes->iter = k; 186 snes->norm = fnorm; 187 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 188 ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 189 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 190 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 191 if (snes->reason) PetscFunctionReturn(0); 192 } 193 snes->reason = SNES_DIVERGED_MAX_IT; 194 PetscFunctionReturn(0); 195 } 196 197 /*MC 198 SNESAnderson - Anderson Mixing method. 199 200 Level: beginner 201 202 Options Database: 203 + -snes_anderson_m - Number of stored previous solutions and residuals 204 . -snes_anderson_beta - Relaxation parameter; X_{update} = X + \beta F 205 . -snes_anderson_restart_type - Type of restart (see SNESNGMRES) 206 . -snes_anderson_restart_it - Number of iterations of restart conditions before restart 207 . -snes_anderson_restart - Number of iterations before periodic restart 208 - -snes_anderson_monitor - Prints relevant information about the ngmres iteration 209 210 Notes: 211 212 The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized 213 optimization problem at each iteration. 214 215 References: 216 217 "D. G. Anderson. Iterative procedures for nonlinear integral equations. 218 J. Assoc. Comput. Mach., 12:547-560, 1965." 219 220 .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 221 M*/ 222 223 #undef __FUNCT__ 224 #define __FUNCT__ "SNESCreate_Anderson" 225 PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes) 226 { 227 SNES_NGMRES *ngmres; 228 PetscErrorCode ierr; 229 230 PetscFunctionBegin; 231 snes->ops->destroy = SNESDestroy_NGMRES; 232 snes->ops->setup = SNESSetUp_NGMRES; 233 snes->ops->setfromoptions = SNESSetFromOptions_Anderson; 234 snes->ops->view = SNESView_NGMRES; 235 snes->ops->solve = SNESSolve_Anderson; 236 snes->ops->reset = SNESReset_NGMRES; 237 238 snes->usespc = PETSC_TRUE; 239 snes->usesksp = PETSC_FALSE; 240 snes->pcside = PC_RIGHT; 241 242 ierr = PetscNewLog(snes,&ngmres);CHKERRQ(ierr); 243 snes->data = (void*) ngmres; 244 ngmres->msize = 30; 245 246 if (!snes->tolerancesset) { 247 snes->max_funcs = 30000; 248 snes->max_its = 10000; 249 } 250 251 ngmres->additive_linesearch = NULL; 252 ngmres->approxfunc = PETSC_FALSE; 253 ngmres->restart_type = SNES_NGMRES_RESTART_NONE; 254 ngmres->restart_it = 2; 255 ngmres->restart_periodic = 30; 256 ngmres->gammaA = 2.0; 257 ngmres->gammaC = 2.0; 258 ngmres->deltaB = 0.9; 259 ngmres->epsilonB = 0.1; 260 261 ngmres->andersonBeta = 1.0; 262 PetscFunctionReturn(0); 263 } 264 265