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