1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 2 #include <petscblaslapack.h> 3 #include <petscdm.h> 4 5 const char *const SNESNGMRESRestartTypes[] = {"NONE","PERIODIC","DIFFERENCE","SNESNGMRESRestartType","SNES_NGMRES_RESTART_",0}; 6 const char *const SNESNGMRESSelectTypes[] = {"NONE","DIFFERENCE","LINESEARCH","SNESNGMRESSelectType","SNES_NGMRES_SELECT_",0}; 7 8 PetscErrorCode SNESReset_NGMRES(SNES snes) 9 { 10 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 11 PetscErrorCode ierr; 12 13 PetscFunctionBegin; 14 ierr = VecDestroyVecs(ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr); 15 ierr = VecDestroyVecs(ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr); 16 ierr = SNESLineSearchDestroy(&ngmres->additive_linesearch);CHKERRQ(ierr); 17 PetscFunctionReturn(0); 18 } 19 20 PetscErrorCode SNESDestroy_NGMRES(SNES snes) 21 { 22 PetscErrorCode ierr; 23 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 24 25 PetscFunctionBegin; 26 ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr); 27 ierr = PetscFree4(ngmres->h,ngmres->beta,ngmres->xi,ngmres->q);CHKERRQ(ierr); 28 ierr = PetscFree3(ngmres->xnorms,ngmres->fnorms,ngmres->s);CHKERRQ(ierr); 29 #if defined(PETSC_USE_COMPLEX) 30 ierr = PetscFree(ngmres->rwork);CHKERRQ(ierr); 31 #endif 32 ierr = PetscFree(ngmres->work);CHKERRQ(ierr); 33 ierr = PetscFree(snes->data);CHKERRQ(ierr); 34 PetscFunctionReturn(0); 35 } 36 37 PetscErrorCode SNESSetUp_NGMRES(SNES snes) 38 { 39 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 40 const char *optionsprefix; 41 PetscInt msize,hsize; 42 PetscErrorCode ierr; 43 DM dm; 44 45 PetscFunctionBegin; 46 if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 47 SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNGMRES does not support left preconditioning with unpreconditioned function"); 48 } 49 if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 50 ierr = SNESSetWorkVecs(snes,5);CHKERRQ(ierr); 51 52 if (!snes->vec_sol) { 53 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 54 ierr = DMCreateGlobalVector(dm,&snes->vec_sol);CHKERRQ(ierr); 55 } 56 57 if (!ngmres->Xdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);} 58 if (!ngmres->Fdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);} 59 if (!ngmres->setup_called) { 60 msize = ngmres->msize; /* restart size */ 61 hsize = msize * msize; 62 63 /* explicit least squares minimization solve */ 64 ierr = PetscCalloc4(hsize,&ngmres->h, msize,&ngmres->beta, msize,&ngmres->xi, hsize,&ngmres->q);CHKERRQ(ierr); 65 ierr = PetscMalloc3(msize,&ngmres->xnorms,msize,&ngmres->fnorms,msize,&ngmres->s);CHKERRQ(ierr); 66 ngmres->nrhs = 1; 67 ngmres->lda = msize; 68 ngmres->ldb = msize; 69 ngmres->lwork = 12*msize; 70 #if defined(PETSC_USE_COMPLEX) 71 ierr = PetscMalloc1(ngmres->lwork,&ngmres->rwork);CHKERRQ(ierr); 72 #endif 73 ierr = PetscMalloc1(ngmres->lwork,&ngmres->work);CHKERRQ(ierr); 74 } 75 76 /* linesearch setup */ 77 ierr = SNESGetOptionsPrefix(snes,&optionsprefix);CHKERRQ(ierr); 78 79 if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 80 ierr = SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch);CHKERRQ(ierr); 81 ierr = SNESLineSearchSetSNES(ngmres->additive_linesearch,snes);CHKERRQ(ierr); 82 ierr = SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2);CHKERRQ(ierr); 83 ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_");CHKERRQ(ierr); 84 ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix);CHKERRQ(ierr); 85 ierr = SNESLineSearchSetFromOptions(ngmres->additive_linesearch);CHKERRQ(ierr); 86 } 87 88 ngmres->setup_called = PETSC_TRUE; 89 PetscFunctionReturn(0); 90 } 91 92 PetscErrorCode SNESSetFromOptions_NGMRES(PetscOptionItems *PetscOptionsObject,SNES snes) 93 { 94 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 95 PetscErrorCode ierr; 96 PetscBool debug = PETSC_FALSE; 97 98 PetscFunctionBegin; 99 ierr = PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");CHKERRQ(ierr); 100 ierr = PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes, 101 (PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL);CHKERRQ(ierr); 102 ierr = PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes, 103 (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr); 104 ierr = PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,NULL);CHKERRQ(ierr); 105 ierr = PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL);CHKERRQ(ierr); 106 ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr); 107 ierr = PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr); 108 ierr = PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr); 109 ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);CHKERRQ(ierr); 110 if (debug) { 111 ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 112 } 113 ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,NULL);CHKERRQ(ierr); 114 ierr = PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,NULL);CHKERRQ(ierr); 115 ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL);CHKERRQ(ierr); 116 ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL);CHKERRQ(ierr); 117 ierr = PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL);CHKERRQ(ierr); 118 ierr = PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise",ngmres->restart_fm_rise,&ngmres->restart_fm_rise,NULL);CHKERRQ(ierr); 119 ierr = PetscOptionsTail();CHKERRQ(ierr); 120 if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 121 PetscFunctionReturn(0); 122 } 123 124 PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer) 125 { 126 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 127 PetscBool iascii; 128 PetscErrorCode ierr; 129 130 PetscFunctionBegin; 131 ierr = PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 132 if (iascii) { 133 ierr = PetscViewerASCIIPrintf(viewer," Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr); 134 ierr = PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);CHKERRQ(ierr); 135 ierr = PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);CHKERRQ(ierr); 136 ierr = PetscViewerASCIIPrintf(viewer," Restart on F_M residual increase: %s\n",ngmres->restart_fm_rise?"TRUE":"FALSE");CHKERRQ(ierr); 137 } 138 PetscFunctionReturn(0); 139 } 140 141 PetscErrorCode SNESSolve_NGMRES(SNES snes) 142 { 143 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 144 /* present solution, residual, and preconditioned residual */ 145 Vec X,F,B,D,Y; 146 147 /* candidate linear combination answers */ 148 Vec XA,FA,XM,FM; 149 150 /* coefficients and RHS to the minimization problem */ 151 PetscReal fnorm,fMnorm,fAnorm; 152 PetscReal xnorm,xMnorm,xAnorm; 153 PetscReal ynorm,yMnorm,yAnorm; 154 PetscInt k,k_restart,l,ivec,restart_count = 0; 155 156 /* solution selection data */ 157 PetscBool selectRestart; 158 /* 159 These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed 160 to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case 161 so the code is correct as written. 162 */ 163 PetscReal dnorm = 0.0,dminnorm = 0.0; 164 PetscReal fminnorm; 165 166 SNESConvergedReason reason; 167 SNESLineSearchReason lssucceed; 168 PetscErrorCode ierr; 169 170 PetscFunctionBegin; 171 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); 172 173 ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 174 /* variable initialization */ 175 snes->reason = SNES_CONVERGED_ITERATING; 176 X = snes->vec_sol; 177 F = snes->vec_func; 178 B = snes->vec_rhs; 179 XA = snes->vec_sol_update; 180 FA = snes->work[0]; 181 D = snes->work[1]; 182 183 /* work for the line search */ 184 Y = snes->work[2]; 185 XM = snes->work[3]; 186 FM = snes->work[4]; 187 188 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 189 snes->iter = 0; 190 snes->norm = 0.; 191 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 192 193 /* initialization */ 194 195 if (snes->npc && snes->npcside== PC_LEFT) { 196 ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 197 ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 198 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 199 snes->reason = SNES_DIVERGED_INNER; 200 PetscFunctionReturn(0); 201 } 202 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 203 } else { 204 if (!snes->vec_func_init_set) { 205 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 206 } else snes->vec_func_init_set = PETSC_FALSE; 207 208 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 209 SNESCheckFunctionNorm(snes,fnorm); 210 } 211 fminnorm = fnorm; 212 213 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 214 snes->norm = fnorm; 215 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 216 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 217 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 218 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 219 if (snes->reason) PetscFunctionReturn(0); 220 SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X); 221 222 k_restart = 1; 223 l = 1; 224 ivec = 0; 225 for (k=1; k < snes->max_its+1; k++) { 226 /* Computation of x^M */ 227 if (snes->npc && snes->npcside== PC_RIGHT) { 228 ierr = VecCopy(X,XM);CHKERRQ(ierr); 229 ierr = SNESSetInitialFunction(snes->npc,F);CHKERRQ(ierr); 230 231 ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0);CHKERRQ(ierr); 232 ierr = SNESSolve(snes->npc,B,XM);CHKERRQ(ierr); 233 ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0);CHKERRQ(ierr); 234 235 ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 236 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 237 snes->reason = SNES_DIVERGED_INNER; 238 PetscFunctionReturn(0); 239 } 240 ierr = SNESGetNPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr); 241 } else { 242 /* no preconditioner -- just take gradient descent with line search */ 243 ierr = VecCopy(F,Y);CHKERRQ(ierr); 244 ierr = VecCopy(F,FM);CHKERRQ(ierr); 245 ierr = VecCopy(X,XM);CHKERRQ(ierr); 246 247 fMnorm = fnorm; 248 249 ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr); 250 ierr = SNESLineSearchGetReason(snes->linesearch,&lssucceed);CHKERRQ(ierr); 251 if (lssucceed) { 252 if (++snes->numFailures >= snes->maxFailures) { 253 snes->reason = SNES_DIVERGED_LINE_SEARCH; 254 PetscFunctionReturn(0); 255 } 256 } 257 } 258 259 ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr); 260 /* r = F(x) */ 261 if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 262 263 /* differences for selection and restart */ 264 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 265 ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);CHKERRQ(ierr); 266 } else { 267 ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);CHKERRQ(ierr); 268 } 269 SNESCheckFunctionNorm(snes,fnorm); 270 271 /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 272 ierr = SNESNGMRESSelect_Private(snes,k_restart,XM,FM,xMnorm,fMnorm,yMnorm,XA,FA,xAnorm,fAnorm,yAnorm,dnorm,fminnorm,dminnorm,X,F,Y,&xnorm,&fnorm,&ynorm);CHKERRQ(ierr); 273 selectRestart = PETSC_FALSE; 274 275 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 276 ierr = SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 277 278 /* if the restart conditions persist for more than restart_it iterations, restart. */ 279 if (selectRestart) restart_count++; 280 else restart_count = 0; 281 } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 282 if (k_restart > ngmres->restart_periodic) { 283 if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr); 284 restart_count = ngmres->restart_it; 285 } 286 } 287 288 ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 289 290 /* restart after restart conditions have persisted for a fixed number of iterations */ 291 if (restart_count >= ngmres->restart_it) { 292 if (ngmres->monitor) { 293 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr); 294 } 295 restart_count = 0; 296 k_restart = 1; 297 l = 1; 298 ivec = 0; 299 /* q_{00} = nu */ 300 ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM);CHKERRQ(ierr); 301 } else { 302 /* select the current size of the subspace */ 303 if (l < ngmres->msize) l++; 304 k_restart++; 305 /* place the current entry in the list of previous entries */ 306 if (ngmres->candidate) { 307 if (fminnorm > fMnorm) fminnorm = fMnorm; 308 ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM);CHKERRQ(ierr); 309 } else { 310 if (fminnorm > fnorm) fminnorm = fnorm; 311 ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X);CHKERRQ(ierr); 312 } 313 } 314 315 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 316 snes->iter = k; 317 snes->norm = fnorm; 318 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 319 ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 320 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 321 ierr = (*snes->ops->converged)(snes,snes->iter,0,0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 322 if (snes->reason) PetscFunctionReturn(0); 323 } 324 snes->reason = SNES_DIVERGED_MAX_IT; 325 PetscFunctionReturn(0); 326 } 327 328 /*@ 329 SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M 330 331 Input Parameters: 332 + snes - the SNES context. 333 - flg - boolean value deciding whether to use the option or not 334 335 Options Database: 336 + -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M 337 338 Level: intermediate 339 340 Notes: 341 If the proposed step x_M increases the residual F_M, it might be trying to get out of a stagnation area. 342 To help the solver do that, reset the Krylov subspace whenever F_M increases. 343 344 This option must be used with SNES_NGMRES_RESTART_DIFFERENCE 345 346 The default is FALSE. 347 .seealso: SNES_NGMRES_RESTART_DIFFERENCE 348 @*/ 349 PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes,PetscBool flg) 350 { 351 PetscErrorCode (*f)(SNES,PetscBool); 352 PetscErrorCode ierr; 353 354 PetscFunctionBegin; 355 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",&f);CHKERRQ(ierr); 356 if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 357 PetscFunctionReturn(0); 358 } 359 360 PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes,PetscBool flg) 361 { 362 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 363 364 PetscFunctionBegin; 365 ngmres->restart_fm_rise = flg; 366 PetscFunctionReturn(0); 367 } 368 369 PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes,PetscBool *flg) 370 { 371 PetscErrorCode (*f)(SNES,PetscBool*); 372 PetscErrorCode ierr; 373 374 PetscFunctionBegin; 375 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",&f);CHKERRQ(ierr); 376 if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 377 PetscFunctionReturn(0); 378 } 379 380 PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes,PetscBool *flg) 381 { 382 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 383 384 PetscFunctionBegin; 385 *flg = ngmres->restart_fm_rise; 386 PetscFunctionReturn(0); 387 } 388 389 390 /*@ 391 SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES. 392 393 Logically Collective on SNES 394 395 Input Parameters: 396 + snes - the iterative context 397 - rtype - restart type 398 399 Options Database: 400 + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type 401 - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic 402 403 Level: intermediate 404 405 SNESNGMRESRestartTypes: 406 + SNES_NGMRES_RESTART_NONE - never restart 407 . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria 408 - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations 409 410 Notes: 411 The default line search used is the L2 line search and it requires two additional function evaluations. 412 413 @*/ 414 PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype) 415 { 416 PetscErrorCode ierr; 417 418 PetscFunctionBegin; 419 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 420 ierr = PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));CHKERRQ(ierr); 421 PetscFunctionReturn(0); 422 } 423 424 /*@ 425 SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and 426 combined solution are used to create the next iterate. 427 428 Logically Collective on SNES 429 430 Input Parameters: 431 + snes - the iterative context 432 - stype - selection type 433 434 Options Database: 435 . -snes_ngmres_select_type<difference,none,linesearch> 436 437 Level: intermediate 438 439 SNESNGMRESSelectTypes: 440 + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time 441 . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria 442 - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination 443 444 Notes: 445 The default line search used is the L2 line search and it requires two additional function evaluations. 446 447 @*/ 448 PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype) 449 { 450 PetscErrorCode ierr; 451 452 PetscFunctionBegin; 453 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 454 ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr); 455 PetscFunctionReturn(0); 456 } 457 458 PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype) 459 { 460 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 461 462 PetscFunctionBegin; 463 ngmres->select_type = stype; 464 PetscFunctionReturn(0); 465 } 466 467 PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype) 468 { 469 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 470 471 PetscFunctionBegin; 472 ngmres->restart_type = rtype; 473 PetscFunctionReturn(0); 474 } 475 476 /*MC 477 SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 478 479 Level: beginner 480 481 Options Database: 482 + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution 483 . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions 484 . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions 485 . -snes_ngmres_m - Number of stored previous solutions and residuals 486 . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart 487 . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination 488 . -snes_ngmres_gammaC - Residual tolerance for restart 489 . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart 490 . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart 491 . -snes_ngmres_restart_fm_rise - Restart on residual rise from x_M step 492 . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration 493 . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother 494 - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type 495 496 Notes: 497 498 The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 499 optimization problem at each iteration. 500 501 Very similar to the SNESANDERSON algorithm. 502 503 References: 504 + 1. - C. W. Oosterlee and T. Washio, "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", 505 SIAM Journal on Scientific Computing, 21(5), 2000. 506 - 2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 507 SIAM Review, 57(4), 2015 508 509 510 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 511 M*/ 512 513 PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes) 514 { 515 SNES_NGMRES *ngmres; 516 PetscErrorCode ierr; 517 SNESLineSearch linesearch; 518 519 PetscFunctionBegin; 520 snes->ops->destroy = SNESDestroy_NGMRES; 521 snes->ops->setup = SNESSetUp_NGMRES; 522 snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 523 snes->ops->view = SNESView_NGMRES; 524 snes->ops->solve = SNESSolve_NGMRES; 525 snes->ops->reset = SNESReset_NGMRES; 526 527 snes->usesnpc = PETSC_TRUE; 528 snes->usesksp = PETSC_FALSE; 529 snes->npcside = PC_RIGHT; 530 531 snes->alwayscomputesfinalresidual = PETSC_TRUE; 532 533 ierr = PetscNewLog(snes,&ngmres);CHKERRQ(ierr); 534 snes->data = (void*) ngmres; 535 ngmres->msize = 30; 536 537 if (!snes->tolerancesset) { 538 snes->max_funcs = 30000; 539 snes->max_its = 10000; 540 } 541 542 ngmres->candidate = PETSC_FALSE; 543 544 ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 545 ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 546 547 ngmres->additive_linesearch = NULL; 548 ngmres->approxfunc = PETSC_FALSE; 549 ngmres->restart_it = 2; 550 ngmres->restart_periodic = 30; 551 ngmres->gammaA = 2.0; 552 ngmres->gammaC = 2.0; 553 ngmres->deltaB = 0.9; 554 ngmres->epsilonB = 0.1; 555 ngmres->restart_fm_rise = PETSC_FALSE; 556 557 ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE; 558 ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE; 559 560 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr); 561 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr); 562 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",SNESNGMRESSetRestartFmRise_NGMRES);CHKERRQ(ierr); 563 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",SNESNGMRESGetRestartFmRise_NGMRES);CHKERRQ(ierr); 564 PetscFunctionReturn(0); 565 } 566