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