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