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