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,PETSC_NULL);CHKERRQ(ierr); 22 ierr = PetscOptionsReal("-snes_anderson_beta","Number of directions","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,PETSC_NULL);CHKERRQ(ierr); 23 ierr = PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes, 24 (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,PETSC_NULL);CHKERRQ(ierr); 25 ierr = PetscOptionsBool("-snes_anderson_monitor","Monitor actions of NGMRES","SNES",ngmres->monitor ? PETSC_TRUE: PETSC_FALSE,&debug,PETSC_NULL);CHKERRQ(ierr); 26 ierr = PetscOptionsInt("-snes_anderson_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,PETSC_NULL);CHKERRQ(ierr); 27 ierr = PetscOptionsInt("-snes_anderson_restart_it","Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,PETSC_NULL);CHKERRQ(ierr); 28 if (debug) { 29 ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);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 = SNESGetSNESLineSearch(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,FPC; 50 51 /* coefficients and RHS to the minimization problem */ 52 PetscReal fnorm,fMnorm; 53 PetscReal dnorm,dminnorm=0.0,fminnorm; 54 PetscInt restart_count=0; 55 PetscInt k,k_restart,l,ivec; 56 57 PetscBool selectRestart; 58 59 SNESConvergedReason reason; 60 PetscErrorCode ierr; 61 62 PetscFunctionBegin; 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 = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 77 snes->iter = 0; 78 snes->norm = 0.; 79 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 80 81 /* initialization */ 82 83 /* r = F(x) */ 84 if (!snes->vec_func_init_set) { 85 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 86 if (snes->domainerror) { 87 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 88 PetscFunctionReturn(0); 89 } 90 } else { 91 snes->vec_func_init_set = PETSC_FALSE; 92 } 93 94 if (!snes->norm_init_set) { 95 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 96 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in function evaluation"); 97 } else { 98 fnorm = snes->norm_init; 99 snes->norm_init_set = PETSC_FALSE; 100 } 101 fminnorm = fnorm; 102 103 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 104 snes->norm = fnorm; 105 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 106 SNESLogConvHistory(snes,fnorm,0); 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 for (k=1; k < snes->max_its+1; k++) { 114 /* select which vector of the stored subspace will be updated */ 115 ivec = k_restart % ngmres->msize; 116 if (snes->pc && snes->pcside == PC_RIGHT) { 117 ierr = VecCopy(X,XM);CHKERRQ(ierr); 118 ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr); 119 ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr); 120 121 ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 122 ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr); 123 ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 124 125 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 126 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 127 snes->reason = SNES_DIVERGED_INNER; 128 PetscFunctionReturn(0); 129 } 130 ierr = SNESGetFunction(snes->pc,&FPC,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 131 ierr = VecCopy(FPC,FM);CHKERRQ(ierr); 132 ierr = SNESGetFunctionNorm(snes->pc,&fMnorm);CHKERRQ(ierr); 133 if (ngmres->andersonBeta != 1.0) { 134 VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr); 135 } 136 } else { 137 ierr = VecCopy(F,FM);CHKERRQ(ierr); 138 ierr = VecCopy(X,XM);CHKERRQ(ierr); 139 ierr = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr); 140 fMnorm = fnorm; 141 } 142 143 ierr = SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr); 144 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 145 ierr = SNESNGMRESCalculateDifferences_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,PETSC_NULL);CHKERRQ(ierr); 146 ierr = SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 147 /* if the restart conditions persist for more than restart_it iterations, restart. */ 148 if (selectRestart) { 149 restart_count++; 150 } else { 151 restart_count = 0; 152 } 153 } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 154 if (k_restart > ngmres->restart_periodic) { 155 if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr); 156 restart_count = ngmres->restart_it; 157 } 158 } 159 /* restart after restart conditions have persisted for a fixed number of iterations */ 160 if (restart_count >= ngmres->restart_it) { 161 if (ngmres->monitor) { 162 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr); 163 } 164 restart_count = 0; 165 k_restart = 0; 166 l = 0; 167 } else { 168 if (l < ngmres->msize)l++; 169 k_restart++; 170 ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr); 171 } 172 173 ierr = VecNorm(FA,NORM_2,&fnorm);CHKERRQ(ierr); 174 if (fminnorm > fnorm) fminnorm = fnorm; 175 176 ierr = VecCopy(XA,X);CHKERRQ(ierr); 177 ierr = VecCopy(FA,F);CHKERRQ(ierr); 178 179 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 180 snes->iter = k; 181 snes->norm = fnorm; 182 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 183 SNESLogConvHistory(snes,snes->norm,snes->iter); 184 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 185 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 186 if (snes->reason) PetscFunctionReturn(0); 187 } 188 snes->reason = SNES_DIVERGED_MAX_IT; 189 PetscFunctionReturn(0); 190 } 191 192 /*MC 193 SNESAnderson - Anderson Mixing method. 194 195 Level: beginner 196 197 Options Database: 198 + -snes_anderson_m - Number of stored previous solutions and residuals 199 . -snes_anderson_beta - Relaxation parameter; X_{update} = X + \beta F 200 - -snes_anderson_monitor - Prints relevant information about the ngmres iteration 201 202 Notes: 203 204 The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized 205 optimization problem at each iteration. 206 207 References: 208 209 "D. G. Anderson. Iterative procedures for nonlinear integral equations. 210 J. Assoc. Comput. Mach., 12:547–560, 1965." 211 212 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 213 M*/ 214 215 EXTERN_C_BEGIN 216 #undef __FUNCT__ 217 #define __FUNCT__ "SNESCreate_Anderson" 218 PetscErrorCode SNESCreate_Anderson(SNES snes) 219 { 220 SNES_NGMRES *ngmres; 221 PetscErrorCode ierr; 222 223 PetscFunctionBegin; 224 snes->ops->destroy = SNESDestroy_NGMRES; 225 snes->ops->setup = SNESSetUp_NGMRES; 226 snes->ops->setfromoptions = SNESSetFromOptions_Anderson; 227 snes->ops->view = SNESView_NGMRES; 228 snes->ops->solve = SNESSolve_Anderson; 229 snes->ops->reset = SNESReset_NGMRES; 230 231 snes->usespc = PETSC_TRUE; 232 snes->usesksp = PETSC_FALSE; 233 234 ierr = PetscNewLog(snes,SNES_NGMRES,&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 = PETSC_NULL; 244 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 EXTERN_C_END 257