1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 2 #include <petscblaslapack.h> 3 4 const char *const SNESNGMRESRestartTypes[] = {"NONE","PERIODIC","DIFFERENCE","SNESNGMRESRestartType","SNES_NGMRES_RESTART_",0}; 5 const char *const SNESNGMRESSelectTypes[] = {"NONE","DIFFERENCE","LINESEARCH","SNESNGMRESSelectType","SNES_NGMRES_SELECT_",0}; 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "SNESReset_NGMRES" 9 PetscErrorCode SNESReset_NGMRES(SNES snes) 10 { 11 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 12 PetscErrorCode ierr; 13 14 PetscFunctionBegin; 15 ierr = VecDestroyVecs(ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr); 16 ierr = VecDestroyVecs(ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr); 17 ierr = SNESLineSearchDestroy(&ngmres->additive_linesearch);CHKERRQ(ierr); 18 PetscFunctionReturn(0); 19 } 20 21 #undef __FUNCT__ 22 #define __FUNCT__ "SNESDestroy_NGMRES" 23 PetscErrorCode SNESDestroy_NGMRES(SNES snes) 24 { 25 PetscErrorCode ierr; 26 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 27 28 PetscFunctionBegin; 29 ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr); 30 ierr = PetscFree5(ngmres->h,ngmres->beta,ngmres->xi,ngmres->fnorms,ngmres->q);CHKERRQ(ierr); 31 ierr = PetscFree(ngmres->s);CHKERRQ(ierr); 32 ierr = PetscFree(ngmres->xnorms);CHKERRQ(ierr); 33 #if PETSC_USE_COMPLEX 34 ierr = PetscFree(ngmres->rwork);CHKERRQ(ierr); 35 #endif 36 ierr = PetscFree(ngmres->work);CHKERRQ(ierr); 37 ierr = PetscFree(snes->data);CHKERRQ(ierr); 38 PetscFunctionReturn(0); 39 } 40 41 #undef __FUNCT__ 42 #define __FUNCT__ "SNESSetUp_NGMRES" 43 PetscErrorCode SNESSetUp_NGMRES(SNES snes) 44 { 45 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 46 const char *optionsprefix; 47 PetscInt msize,hsize; 48 PetscErrorCode ierr; 49 50 PetscFunctionBegin; 51 ierr = SNESSetWorkVecs(snes,5);CHKERRQ(ierr); 52 if (!ngmres->Xdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);} 53 if (!ngmres->Fdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);} 54 if (!ngmres->setup_called) { 55 msize = ngmres->msize; /* restart size */ 56 hsize = msize * msize; 57 58 /* explicit least squares minimization solve */ 59 ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h, 60 msize,PetscScalar,&ngmres->beta, 61 msize,PetscScalar,&ngmres->xi, 62 msize,PetscReal, &ngmres->fnorms, 63 hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr); 64 if (ngmres->singlereduction) { 65 ierr = PetscMalloc(msize*sizeof(PetscReal),&ngmres->xnorms);CHKERRQ(ierr); 66 } 67 ngmres->nrhs = 1; 68 ngmres->lda = msize; 69 ngmres->ldb = msize; 70 ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr); 71 ierr = PetscMemzero(ngmres->h, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 72 ierr = PetscMemzero(ngmres->q, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 73 ierr = PetscMemzero(ngmres->xi, msize*sizeof(PetscScalar));CHKERRQ(ierr); 74 ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr); 75 ngmres->lwork = 12*msize; 76 #if PETSC_USE_COMPLEX 77 ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork);CHKERRQ(ierr); 78 #endif 79 ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work);CHKERRQ(ierr); 80 } 81 82 /* linesearch setup */ 83 ierr = SNESGetOptionsPrefix(snes,&optionsprefix);CHKERRQ(ierr); 84 85 if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 86 ierr = SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch);CHKERRQ(ierr); 87 ierr = SNESLineSearchSetSNES(ngmres->additive_linesearch,snes);CHKERRQ(ierr); 88 ierr = SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2);CHKERRQ(ierr); 89 ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_");CHKERRQ(ierr); 90 ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix);CHKERRQ(ierr); 91 ierr = SNESLineSearchSetFromOptions(ngmres->additive_linesearch);CHKERRQ(ierr); 92 } 93 94 ngmres->setup_called = PETSC_TRUE; 95 PetscFunctionReturn(0); 96 } 97 98 #undef __FUNCT__ 99 #define __FUNCT__ "SNESSetFromOptions_NGMRES" 100 PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes) 101 { 102 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 103 PetscErrorCode ierr; 104 PetscBool debug; 105 SNESLineSearch linesearch; 106 107 PetscFunctionBegin; 108 ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr); 109 ierr = PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes, 110 (PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL);CHKERRQ(ierr); 111 ierr = PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes, 112 (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr); 113 ierr = PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,NULL);CHKERRQ(ierr); 114 ierr = PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL);CHKERRQ(ierr); 115 ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr); 116 ierr = PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr); 117 ierr = PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr); 118 ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);CHKERRQ(ierr); 119 if (debug) { 120 ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 121 } 122 ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,NULL);CHKERRQ(ierr); 123 ierr = PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,NULL);CHKERRQ(ierr); 124 ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL);CHKERRQ(ierr); 125 ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL);CHKERRQ(ierr); 126 ierr = PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL);CHKERRQ(ierr); 127 ierr = PetscOptionsTail();CHKERRQ(ierr); 128 if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 129 130 /* set the default type of the line search if the user hasn't already. */ 131 if (!snes->linesearch) { 132 ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 133 ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 134 } 135 PetscFunctionReturn(0); 136 } 137 138 #undef __FUNCT__ 139 #define __FUNCT__ "SNESView_NGMRES" 140 PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer) 141 { 142 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 143 PetscBool iascii; 144 PetscErrorCode ierr; 145 146 PetscFunctionBegin; 147 ierr = PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 148 if (iascii) { 149 ierr = PetscViewerASCIIPrintf(viewer," Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr); 150 ierr = PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);CHKERRQ(ierr); 151 ierr = PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);CHKERRQ(ierr); 152 } 153 PetscFunctionReturn(0); 154 } 155 156 #undef __FUNCT__ 157 #define __FUNCT__ "SNESSolve_NGMRES" 158 PetscErrorCode SNESSolve_NGMRES(SNES snes) 159 { 160 161 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 162 /* present solution, residual, and preconditioned residual */ 163 Vec X,F,B,D,Y; 164 165 /* candidate linear combination answers */ 166 Vec XA,FA,XM,FM,FPC; 167 168 /* coefficients and RHS to the minimization problem */ 169 PetscReal fnorm,fMnorm,fAnorm; 170 PetscInt k,k_restart,l,ivec,restart_count = 0; 171 172 /* solution selection data */ 173 PetscBool selectRestart; 174 PetscReal dnorm,dminnorm = 0.0; 175 PetscReal fminnorm,xnorm,ynorm; 176 177 SNESConvergedReason reason; 178 PetscBool lssucceed; 179 PCSide npcside; 180 PetscErrorCode ierr; 181 182 PetscFunctionBegin; 183 /* variable initialization */ 184 snes->reason = SNES_CONVERGED_ITERATING; 185 X = snes->vec_sol; 186 F = snes->vec_func; 187 B = snes->vec_rhs; 188 XA = snes->vec_sol_update; 189 FA = snes->work[0]; 190 D = snes->work[1]; 191 192 /* work for the line search */ 193 Y = snes->work[2]; 194 XM = snes->work[3]; 195 FM = snes->work[4]; 196 197 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 198 snes->iter = 0; 199 snes->norm = 0.; 200 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 201 202 /* initialization */ 203 204 /* r = F(x) */ 205 if (!snes->vec_func_init_set) { 206 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 207 if (snes->domainerror) { 208 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 209 PetscFunctionReturn(0); 210 } 211 } else snes->vec_func_init_set = PETSC_FALSE; 212 213 if (!snes->norm_init_set) { 214 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 215 if (PetscIsInfOrNanReal(fnorm)) { 216 snes->reason = SNES_DIVERGED_FNORM_NAN; 217 PetscFunctionReturn(0); 218 } 219 } else { 220 fnorm = snes->norm_init; 221 snes->norm_init_set = PETSC_FALSE; 222 } 223 fminnorm = fnorm; 224 225 /* q_{00} = nu */ 226 SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X); 227 228 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 229 snes->norm = fnorm; 230 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 231 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 232 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 233 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 234 if (snes->reason) PetscFunctionReturn(0); 235 236 k_restart = 1; 237 l = 1; 238 for (k=1; k < snes->max_its+1; k++) { 239 /* select which vector of the stored subspace will be updated */ 240 ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 241 242 /* Computation of x^M */ 243 if (snes->pc && snes->pcside == PC_RIGHT) { 244 ierr = VecCopy(X,XM);CHKERRQ(ierr); 245 ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr); 246 ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr); 247 248 ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 249 ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr); 250 ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 251 252 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 253 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 254 snes->reason = SNES_DIVERGED_INNER; 255 PetscFunctionReturn(0); 256 } 257 ierr = SNESGetPCSide(snes->pc,&npcside);CHKERRQ(ierr); 258 if (npcside == PC_RIGHT) { 259 ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr); 260 ierr = VecCopy(FPC,FM);CHKERRQ(ierr); 261 ierr = SNESGetFunctionNorm(snes->pc,&fMnorm);CHKERRQ(ierr); 262 } else { 263 ierr = SNESComputeFunction(snes,XM,FM);CHKERRQ(ierr); 264 ierr = VecNorm(FM,NORM_2,&fMnorm);CHKERRQ(ierr); 265 } 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 = SNESLineSearchGetSuccess(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 ierr = SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr); 284 /* r = F(x) */ 285 if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 286 287 /* differences for selection and restart */ 288 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 289 ierr = SNESNGMRESCalculateDifferences_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&fAnorm);CHKERRQ(ierr); 290 } else { 291 ierr = VecNorm(FA,NORM_2,&fAnorm);CHKERRQ(ierr); 292 } 293 if (PetscIsInfOrNanReal(fAnorm)) { 294 snes->reason = SNES_DIVERGED_FNORM_NAN; 295 PetscFunctionReturn(0); 296 } 297 298 /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 299 ierr = SNESNGMRESSelect_Private(snes,k_restart,XM,FM,fMnorm,XA,FA,fAnorm,dnorm,fminnorm,dminnorm,X,F,Y,&fnorm);CHKERRQ(ierr); 300 selectRestart = PETSC_FALSE; 301 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 302 ierr = SNESNGMRESSelectRestart_Private(snes,l,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 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 /* restart after restart conditions have persisted for a fixed number of iterations */ 313 if (restart_count >= ngmres->restart_it) { 314 if (ngmres->monitor) { 315 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr); 316 } 317 restart_count = 0; 318 k_restart = 1; 319 l = 1; 320 /* q_{00} = nu */ 321 if (ngmres->candidate) { 322 ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM);CHKERRQ(ierr); 323 } else { 324 ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fMnorm,X);CHKERRQ(ierr); 325 } 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 = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 341 snes->iter = k; 342 snes->norm = fnorm; 343 ierr = PetscObjectAMSGrantAccess((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 = VecNormBegin(Y,NORM_2,&ynorm);CHKERRQ(ierr); 347 ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); 348 ierr = VecNormEnd(Y,NORM_2,&ynorm);CHKERRQ(ierr); 349 ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 350 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 351 if (snes->reason) PetscFunctionReturn(0); 352 } 353 snes->reason = SNES_DIVERGED_MAX_IT; 354 PetscFunctionReturn(0); 355 } 356 357 #undef __FUNCT__ 358 #define __FUNCT__ "SNESNGMRESSetRestartType" 359 /*@ 360 SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES. 361 362 Logically Collective on SNES 363 364 Input Parameters: 365 + snes - the iterative context 366 - rtype - restart type 367 368 Options Database: 369 + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type 370 - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic 371 372 Level: intermediate 373 374 SNESNGMRESRestartTypes: 375 + SNES_NGMRES_RESTART_NONE - never restart 376 . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria 377 - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations 378 379 Notes: 380 The default line search used is the L2 line search and it requires two additional function evaluations. 381 382 .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch 383 @*/ 384 PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype) 385 { 386 PetscErrorCode ierr; 387 388 PetscFunctionBegin; 389 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 390 ierr = PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));CHKERRQ(ierr); 391 PetscFunctionReturn(0); 392 } 393 394 #undef __FUNCT__ 395 #define __FUNCT__ "SNESNGMRESSetSelectType" 396 /*@ 397 SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and 398 combined solution are used to create the next iterate. 399 400 Logically Collective on SNES 401 402 Input Parameters: 403 + snes - the iterative context 404 - stype - selection type 405 406 Options Database: 407 . -snes_ngmres_select_type<difference,none,linesearch> 408 409 Level: intermediate 410 411 SNESNGMRESSelectTypes: 412 + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time 413 . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria 414 - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination 415 416 Notes: 417 The default line search used is the L2 line search and it requires two additional function evaluations. 418 419 .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch 420 @*/ 421 PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype) 422 { 423 PetscErrorCode ierr; 424 425 PetscFunctionBegin; 426 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 427 ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr); 428 PetscFunctionReturn(0); 429 } 430 431 #undef __FUNCT__ 432 #define __FUNCT__ "SNESNGMRESSetSelectType_NGMRES" 433 PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype) 434 { 435 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 436 437 PetscFunctionBegin; 438 ngmres->select_type = stype; 439 PetscFunctionReturn(0); 440 } 441 442 #undef __FUNCT__ 443 #define __FUNCT__ "SNESNGMRESSetRestartType_NGMRES" 444 PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype) 445 { 446 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 447 448 PetscFunctionBegin; 449 ngmres->restart_type = rtype; 450 PetscFunctionReturn(0); 451 } 452 453 /*MC 454 SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 455 456 Level: beginner 457 458 Options Database: 459 + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution 460 . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions 461 . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions 462 . -snes_ngmres_m - Number of stored previous solutions and residuals 463 . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart 464 . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination 465 . -snes_ngmres_gammaC - Residual tolerance for restart 466 . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart 467 . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart 468 . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration 469 . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother 470 - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type 471 472 Notes: 473 474 The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 475 optimization problem at each iteration. 476 477 References: 478 479 "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 480 SIAM Journal on Scientific Computing, 21(5), 2000. 481 482 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 483 M*/ 484 485 #undef __FUNCT__ 486 #define __FUNCT__ "SNESCreate_NGMRES" 487 PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes) 488 { 489 SNES_NGMRES *ngmres; 490 PetscErrorCode ierr; 491 492 PetscFunctionBegin; 493 snes->ops->destroy = SNESDestroy_NGMRES; 494 snes->ops->setup = SNESSetUp_NGMRES; 495 snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 496 snes->ops->view = SNESView_NGMRES; 497 snes->ops->solve = SNESSolve_NGMRES; 498 snes->ops->reset = SNESReset_NGMRES; 499 500 snes->usespc = PETSC_TRUE; 501 snes->usesksp = PETSC_FALSE; 502 503 ierr = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr); 504 snes->data = (void*) ngmres; 505 ngmres->msize = 30; 506 507 if (!snes->tolerancesset) { 508 snes->max_funcs = 30000; 509 snes->max_its = 10000; 510 } 511 512 ngmres->candidate = PETSC_FALSE; 513 514 ngmres->additive_linesearch = NULL; 515 ngmres->approxfunc = PETSC_FALSE; 516 ngmres->restart_it = 2; 517 ngmres->restart_periodic = 30; 518 ngmres->gammaA = 2.0; 519 ngmres->gammaC = 2.0; 520 ngmres->deltaB = 0.9; 521 ngmres->epsilonB = 0.1; 522 523 ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE; 524 ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE; 525 526 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C","SNESNGMRESSetSelectType_NGMRES",SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr); 527 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C","SNESNGMRESSetRestartType_NGMRES",SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr); 528 PetscFunctionReturn(0); 529 } 530 531