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