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