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 = SNESDefaultGetWork(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(((PetscObject)snes)->comm,&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,PETSC_NULL);CHKERRQ(ierr); 111 ierr = PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes, 112 (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,PETSC_NULL);CHKERRQ(ierr); 113 ierr = PetscOptionsBool("-snes_ngmres_candidate","Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,PETSC_NULL);CHKERRQ(ierr); 114 ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,PETSC_NULL);CHKERRQ(ierr); 115 ierr = PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,PETSC_NULL);CHKERRQ(ierr); 116 ierr = PetscOptionsInt("-snes_ngmres_restart_it","Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,PETSC_NULL);CHKERRQ(ierr); 117 ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE: PETSC_FALSE,&debug,PETSC_NULL);CHKERRQ(ierr); 118 if (debug) { 119 ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 120 } 121 ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,PETSC_NULL);CHKERRQ(ierr); 122 ierr = PetscOptionsReal("-snes_ngmres_gammaC", " Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,PETSC_NULL);CHKERRQ(ierr); 123 ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,PETSC_NULL);CHKERRQ(ierr); 124 ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,PETSC_NULL);CHKERRQ(ierr); 125 ierr = PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,PETSC_NULL);CHKERRQ(ierr); 126 ierr = PetscOptionsTail();CHKERRQ(ierr); 127 if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 128 129 /* set the default type of the line search if the user hasn't already. */ 130 if (!snes->linesearch) { 131 ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 132 ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 133 } 134 PetscFunctionReturn(0); 135 } 136 137 #undef __FUNCT__ 138 #define __FUNCT__ "SNESView_NGMRES" 139 PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer) 140 { 141 SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 142 PetscBool iascii; 143 PetscErrorCode ierr; 144 145 PetscFunctionBegin; 146 ierr = PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 147 if (iascii) { 148 ierr = PetscViewerASCIIPrintf(viewer," Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr); 149 ierr = PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);CHKERRQ(ierr); 150 ierr = PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);CHKERRQ(ierr); 151 } 152 PetscFunctionReturn(0); 153 } 154 155 #undef __FUNCT__ 156 #define __FUNCT__ "SNESSolve_NGMRES" 157 PetscErrorCode SNESSolve_NGMRES(SNES snes) 158 { 159 160 SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 161 /* present solution, residual, and preconditioned residual */ 162 Vec X,F,B,D,Y; 163 164 /* candidate linear combination answers */ 165 Vec XA,FA,XM,FM,FPC; 166 167 /* coefficients and RHS to the minimization problem */ 168 PetscReal fnorm,fMnorm,fAnorm; 169 PetscInt k,k_restart,l,ivec,restart_count = 0; 170 171 /* solution selection data */ 172 PetscBool selectRestart; 173 PetscReal dnorm,dminnorm = 0.0; 174 PetscReal fminnorm,xnorm,ynorm; 175 176 SNESConvergedReason reason; 177 PetscBool lssucceed; 178 PetscErrorCode ierr; 179 180 PetscFunctionBegin; 181 /* variable initialization */ 182 snes->reason = SNES_CONVERGED_ITERATING; 183 X = snes->vec_sol; 184 F = snes->vec_func; 185 B = snes->vec_rhs; 186 XA = snes->vec_sol_update; 187 FA = snes->work[0]; 188 D = snes->work[1]; 189 190 /* work for the line search */ 191 Y = snes->work[2]; 192 XM = snes->work[3]; 193 FM = snes->work[4]; 194 195 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 196 snes->iter = 0; 197 snes->norm = 0.; 198 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 199 200 /* initialization */ 201 202 /* r = F(x) */ 203 if (!snes->vec_func_init_set) { 204 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 205 if (snes->domainerror) { 206 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 207 PetscFunctionReturn(0); 208 } 209 } else { 210 snes->vec_func_init_set = PETSC_FALSE; 211 } 212 213 if (!snes->norm_init_set) { 214 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 215 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in function evaluation"); 216 } else { 217 fnorm = snes->norm_init; 218 snes->norm_init_set = PETSC_FALSE; 219 } 220 fminnorm = fnorm; 221 222 /* q_{00} = nu */ 223 SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X); 224 225 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 226 snes->norm = fnorm; 227 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 228 SNESLogConvHistory(snes,fnorm,0); 229 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 230 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 231 if (snes->reason) PetscFunctionReturn(0); 232 233 k_restart = 1; 234 l = 1; 235 for (k=1; k < snes->max_its+1; k++) { 236 /* select which vector of the stored subspace will be updated */ 237 ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 238 239 /* Computation of x^M */ 240 if (snes->pc && snes->pcside == PC_RIGHT) { 241 ierr = VecCopy(X,XM);CHKERRQ(ierr); 242 ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr); 243 ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr); 244 245 ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 246 ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr); 247 ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 248 249 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 250 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 251 snes->reason = SNES_DIVERGED_INNER; 252 PetscFunctionReturn(0); 253 } 254 ierr = SNESGetFunction(snes->pc,&FPC,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 255 ierr = VecCopy(FPC,FM);CHKERRQ(ierr); 256 ierr = SNESGetFunctionNorm(snes->pc,&fMnorm);CHKERRQ(ierr); 257 } else { 258 /* no preconditioner -- just take gradient descent with line search */ 259 ierr = VecCopy(F,Y);CHKERRQ(ierr); 260 ierr = VecCopy(F,FM);CHKERRQ(ierr); 261 ierr = VecCopy(X,XM);CHKERRQ(ierr); 262 fMnorm = fnorm; 263 ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr); 264 ierr = SNESLineSearchGetSuccess(snes->linesearch,&lssucceed);CHKERRQ(ierr); 265 if (!lssucceed) { 266 if (++snes->numFailures >= snes->maxFailures) { 267 snes->reason = SNES_DIVERGED_LINE_SEARCH; 268 PetscFunctionReturn(0); 269 } 270 } 271 } 272 ierr = SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA); 273 /* r = F(x) */ 274 if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 275 276 /* differences for selection and restart */ 277 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 278 ierr = SNESNGMRESCalculateDifferences_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&fAnorm); 279 } else { 280 ierr = VecNorm(FA,NORM_2,&fAnorm);CHKERRQ(ierr); 281 } 282 if (PetscIsInfOrNanReal(fAnorm)) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in function evaluation"); 283 /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 284 ierr = SNESNGMRESSelect_Private(snes,k_restart,XM,FM,fMnorm,XA,FA,fAnorm,dnorm,fminnorm,dminnorm,X,F,Y,&fnorm); 285 selectRestart = PETSC_FALSE; 286 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 287 ierr = SNESNGMRESSelectRestart_Private(snes,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 288 /* if the restart conditions persist for more than restart_it iterations, restart. */ 289 if (selectRestart) { 290 restart_count++; 291 } else { 292 restart_count = 0; 293 } 294 } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 295 if (k_restart > ngmres->restart_periodic) { 296 if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr); 297 restart_count = ngmres->restart_it; 298 } 299 } 300 /* restart after restart conditions have persisted for a fixed number of iterations */ 301 if (restart_count >= ngmres->restart_it) { 302 if (ngmres->monitor) { 303 ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr); 304 } 305 restart_count = 0; 306 k_restart = 1; 307 l = 1; 308 /* q_{00} = nu */ 309 if (ngmres->candidate) { 310 ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM); 311 } else { 312 ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fMnorm,X); 313 } 314 } else { 315 /* select the current size of the subspace */ 316 if (l < ngmres->msize) l++; 317 k_restart++; 318 /* place the current entry in the list of previous entries */ 319 if (ngmres->candidate) { 320 if (fminnorm > fMnorm) fminnorm = fMnorm; 321 ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM); 322 } else { 323 if (fminnorm > fnorm) fminnorm = fnorm; 324 ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X); 325 } 326 } 327 328 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 329 snes->iter = k; 330 snes->norm = fnorm; 331 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 332 SNESLogConvHistory(snes,snes->norm,snes->iter); 333 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 334 ierr = VecNormBegin(Y,NORM_2,&ynorm);CHKERRQ(ierr); 335 ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); 336 ierr = VecNormEnd(Y,NORM_2,&ynorm);CHKERRQ(ierr); 337 ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 338 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 339 if (snes->reason) PetscFunctionReturn(0); 340 } 341 snes->reason = SNES_DIVERGED_MAX_IT; 342 PetscFunctionReturn(0); 343 } 344 345 #undef __FUNCT__ 346 #define __FUNCT__ "SNESNGMRESSetRestartType" 347 /*@ 348 SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES. 349 350 Logically Collective on SNES 351 352 Input Parameters: 353 + snes - the iterative context 354 - rtype - restart type 355 356 Options Database: 357 + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type 358 - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic 359 360 Level: intermediate 361 362 SNESNGMRESRestartTypes: 363 + SNES_NGMRES_RESTART_NONE - never restart 364 . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria 365 - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations 366 367 Notes: 368 The default line search used is the L2 line search and it requires two additional function evaluations. 369 370 .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch 371 @*/ 372 PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype) 373 { 374 PetscErrorCode ierr; 375 376 PetscFunctionBegin; 377 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 378 ierr = PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));CHKERRQ(ierr); 379 PetscFunctionReturn(0); 380 } 381 382 #undef __FUNCT__ 383 #define __FUNCT__ "SNESNGMRESSetSelectType" 384 /*@ 385 SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and 386 combined solution are used to create the next iterate. 387 388 Logically Collective on SNES 389 390 Input Parameters: 391 + snes - the iterative context 392 - stype - selection type 393 394 Options Database: 395 . -snes_ngmres_select_type<difference,none,linesearch> 396 397 Level: intermediate 398 399 SNESNGMRESSelectTypes: 400 + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time 401 . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria 402 - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination 403 404 Notes: 405 The default line search used is the L2 line search and it requires two additional function evaluations. 406 407 .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch 408 @*/ 409 PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype) 410 { 411 PetscErrorCode ierr; 412 413 PetscFunctionBegin; 414 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 415 ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr); 416 PetscFunctionReturn(0); 417 } 418 419 EXTERN_C_BEGIN 420 #undef __FUNCT__ 421 #define __FUNCT__ "SNESNGMRESSetSelectType_NGMRES" 422 PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype) 423 { 424 SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 425 426 PetscFunctionBegin; 427 ngmres->select_type = stype; 428 PetscFunctionReturn(0); 429 } 430 431 #undef __FUNCT__ 432 #define __FUNCT__ "SNESNGMRESSetRestartType_NGMRES" 433 PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype) 434 { 435 SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 436 437 PetscFunctionBegin; 438 ngmres->restart_type = rtype; 439 PetscFunctionReturn(0); 440 } 441 EXTERN_C_END 442 443 /*MC 444 SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 445 446 Level: beginner 447 448 Options Database: 449 + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution 450 . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions 451 . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions 452 . -snes_ngmres_m - Number of stored previous solutions and residuals 453 . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart 454 . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination 455 . -snes_ngmres_gammaC - Residual tolerance for restart 456 . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart 457 . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart 458 . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration 459 . -snes_linesearch_type <basic,basicnonorms,quadratic,critical> - Line search type used for the default smoother 460 - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type 461 462 Notes: 463 464 The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 465 optimization problem at each iteration. 466 467 References: 468 469 "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 470 SIAM Journal on Scientific Computing, 21(5), 2000. 471 472 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 473 M*/ 474 475 EXTERN_C_BEGIN 476 #undef __FUNCT__ 477 #define __FUNCT__ "SNESCreate_NGMRES" 478 PetscErrorCode SNESCreate_NGMRES(SNES snes) 479 { 480 SNES_NGMRES *ngmres; 481 PetscErrorCode ierr; 482 483 PetscFunctionBegin; 484 snes->ops->destroy = SNESDestroy_NGMRES; 485 snes->ops->setup = SNESSetUp_NGMRES; 486 snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 487 snes->ops->view = SNESView_NGMRES; 488 snes->ops->solve = SNESSolve_NGMRES; 489 snes->ops->reset = SNESReset_NGMRES; 490 491 snes->usespc = PETSC_TRUE; 492 snes->usesksp = PETSC_FALSE; 493 494 ierr = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr); 495 snes->data = (void*) ngmres; 496 ngmres->msize = 30; 497 498 if (!snes->tolerancesset) { 499 snes->max_funcs = 30000; 500 snes->max_its = 10000; 501 } 502 503 ngmres->candidate = PETSC_FALSE; 504 505 ngmres->additive_linesearch = PETSC_NULL; 506 507 ngmres->restart_it = 2; 508 ngmres->restart_periodic = 30; 509 ngmres->gammaA = 2.0; 510 ngmres->gammaC = 2.0; 511 ngmres->deltaB = 0.9; 512 ngmres->epsilonB = 0.1; 513 514 ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE; 515 ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE; 516 517 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNGMRESSetSelectType_C","SNESNGMRESSetSelectType_NGMRES",SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr); 518 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNGMRESSetRestartType_C","SNESNGMRESSetRestartType_NGMRES",SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr); 519 PetscFunctionReturn(0); 520 } 521 EXTERN_C_END 522