113a62661SPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 219653cdaSPeter Brune #include <petscblaslapack.h> 3fced5a79SAsbjørn Nilsen Riseth #include <petscdm.h> 4a312c225SMatthew G Knepley 56a6fc655SJed Brown const char *const SNESNGMRESRestartTypes[] = {"NONE","PERIODIC","DIFFERENCE","SNESNGMRESRestartType","SNES_NGMRES_RESTART_",0}; 66a6fc655SJed Brown const char *const SNESNGMRESSelectTypes[] = {"NONE","DIFFERENCE","LINESEARCH","SNESNGMRESSelectType","SNES_NGMRES_SELECT_",0}; 713a62661SPeter Brune 8a312c225SMatthew G Knepley #undef __FUNCT__ 9a312c225SMatthew G Knepley #define __FUNCT__ "SNESReset_NGMRES" 10a312c225SMatthew G Knepley PetscErrorCode SNESReset_NGMRES(SNES snes) 11a312c225SMatthew G Knepley { 12a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 13a312c225SMatthew G Knepley PetscErrorCode ierr; 14a312c225SMatthew G Knepley 15a312c225SMatthew G Knepley PetscFunctionBegin; 16f109b39eSPeter Brune ierr = VecDestroyVecs(ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr); 17f109b39eSPeter Brune ierr = VecDestroyVecs(ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr); 18f1c6b773SPeter Brune ierr = SNESLineSearchDestroy(&ngmres->additive_linesearch);CHKERRQ(ierr); 19a312c225SMatthew G Knepley PetscFunctionReturn(0); 20a312c225SMatthew G Knepley } 21a312c225SMatthew G Knepley 22a312c225SMatthew G Knepley #undef __FUNCT__ 23a312c225SMatthew G Knepley #define __FUNCT__ "SNESDestroy_NGMRES" 24a312c225SMatthew G Knepley PetscErrorCode SNESDestroy_NGMRES(SNES snes) 25a312c225SMatthew G Knepley { 26a312c225SMatthew G Knepley PetscErrorCode ierr; 2778440776SJed Brown SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 28a312c225SMatthew G Knepley 29a312c225SMatthew G Knepley PetscFunctionBegin; 30a312c225SMatthew G Knepley ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr); 31f109b39eSPeter Brune ierr = PetscFree5(ngmres->h,ngmres->beta,ngmres->xi,ngmres->fnorms,ngmres->q);CHKERRQ(ierr); 3219653cdaSPeter Brune ierr = PetscFree(ngmres->s);CHKERRQ(ierr); 3318aa0c0cSPeter Brune ierr = PetscFree(ngmres->xnorms);CHKERRQ(ierr); 3419653cdaSPeter Brune #if PETSC_USE_COMPLEX 3522d28d08SBarry Smith ierr = PetscFree(ngmres->rwork);CHKERRQ(ierr); 3619653cdaSPeter Brune #endif 3722d28d08SBarry Smith ierr = PetscFree(ngmres->work);CHKERRQ(ierr); 3822d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 39a312c225SMatthew G Knepley PetscFunctionReturn(0); 40a312c225SMatthew G Knepley } 41a312c225SMatthew G Knepley 42a312c225SMatthew G Knepley #undef __FUNCT__ 43a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetUp_NGMRES" 44a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes) 45a312c225SMatthew G Knepley { 46a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 47e7058c64SPeter Brune const char *optionsprefix; 4819653cdaSPeter Brune PetscInt msize,hsize; 49a312c225SMatthew G Knepley PetscErrorCode ierr; 50fced5a79SAsbjørn Nilsen Riseth DM dm; 51a312c225SMatthew G Knepley 52a312c225SMatthew G Knepley PetscFunctionBegin; 5346159c86SPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 5446159c86SPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNGMRES does not support left preconditioning with unpreconditioned function"); 5546159c86SPeter Brune } 566c67d002SPeter Brune if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 57fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,5);CHKERRQ(ierr); 58fced5a79SAsbjørn Nilsen Riseth 59fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) { 60fced5a79SAsbjørn Nilsen Riseth ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 61fced5a79SAsbjørn Nilsen Riseth ierr = DMCreateGlobalVector(dm,&snes->vec_sol);CHKERRQ(ierr); 62fced5a79SAsbjørn Nilsen Riseth } 63fced5a79SAsbjørn Nilsen Riseth 6478440776SJed Brown if (!ngmres->Xdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);} 6578440776SJed Brown if (!ngmres->Fdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);} 6678440776SJed Brown if (!ngmres->setup_called) { 67087dfb9eSxuemin msize = ngmres->msize; /* restart size */ 6819653cdaSPeter Brune hsize = msize * msize; 69087dfb9eSxuemin 7098b3e84cSPeter Brune /* explicit least squares minimization solve */ 71dcca6d9dSJed Brown ierr = PetscMalloc5(hsize,&ngmres->h, msize,&ngmres->beta, msize,&ngmres->xi, msize,&ngmres->fnorms, hsize,&ngmres->q);CHKERRQ(ierr); 72785e854fSJed Brown ierr = PetscMalloc1(msize,&ngmres->xnorms);CHKERRQ(ierr); 7319653cdaSPeter Brune ngmres->nrhs = 1; 7419653cdaSPeter Brune ngmres->lda = msize; 7519653cdaSPeter Brune ngmres->ldb = msize; 76785e854fSJed Brown ierr = PetscMalloc1(msize,&ngmres->s);CHKERRQ(ierr); 7719653cdaSPeter Brune ierr = PetscMemzero(ngmres->h, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 7819653cdaSPeter Brune ierr = PetscMemzero(ngmres->q, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 7919653cdaSPeter Brune ierr = PetscMemzero(ngmres->xi, msize*sizeof(PetscScalar));CHKERRQ(ierr); 8019653cdaSPeter Brune ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr); 8119653cdaSPeter Brune ngmres->lwork = 12*msize; 8219653cdaSPeter Brune #if PETSC_USE_COMPLEX 83854ce69bSBarry Smith ierr = PetscMalloc1(ngmres->lwork,&ngmres->rwork);CHKERRQ(ierr); 8419653cdaSPeter Brune #endif 85854ce69bSBarry Smith ierr = PetscMalloc1(ngmres->lwork,&ngmres->work);CHKERRQ(ierr); 8678440776SJed Brown } 87e7058c64SPeter Brune 88e7058c64SPeter Brune /* linesearch setup */ 89e7058c64SPeter Brune ierr = SNESGetOptionsPrefix(snes,&optionsprefix);CHKERRQ(ierr); 90e7058c64SPeter Brune 9113a62661SPeter Brune if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 92ce94432eSBarry Smith ierr = SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch);CHKERRQ(ierr); 93f1c6b773SPeter Brune ierr = SNESLineSearchSetSNES(ngmres->additive_linesearch,snes);CHKERRQ(ierr); 941a4f838cSPeter Brune ierr = SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2);CHKERRQ(ierr); 95f1c6b773SPeter Brune ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_");CHKERRQ(ierr); 96f1c6b773SPeter Brune ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix);CHKERRQ(ierr); 97f1c6b773SPeter Brune ierr = SNESLineSearchSetFromOptions(ngmres->additive_linesearch);CHKERRQ(ierr); 98e7058c64SPeter Brune } 99e7058c64SPeter Brune 10078440776SJed Brown ngmres->setup_called = PETSC_TRUE; 101a312c225SMatthew G Knepley PetscFunctionReturn(0); 102a312c225SMatthew G Knepley } 103a312c225SMatthew G Knepley 104a312c225SMatthew G Knepley #undef __FUNCT__ 105a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetFromOptions_NGMRES" 1064416b707SBarry Smith PetscErrorCode SNESSetFromOptions_NGMRES(PetscOptionItems *PetscOptionsObject,SNES snes) 107a312c225SMatthew G Knepley { 108a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 109a312c225SMatthew G Knepley PetscErrorCode ierr; 11094ae4db5SBarry Smith PetscBool debug = PETSC_FALSE; 111f1c6b773SPeter Brune SNESLineSearch linesearch; 1120adebc6cSBarry Smith 113a312c225SMatthew G Knepley PetscFunctionBegin; 114e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");CHKERRQ(ierr); 11513a62661SPeter Brune ierr = PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes, 1160298fd71SBarry Smith (PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL);CHKERRQ(ierr); 11713a62661SPeter Brune ierr = PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes, 1180298fd71SBarry Smith (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr); 1190298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,NULL);CHKERRQ(ierr); 120077c4231SPeter Brune ierr = PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL);CHKERRQ(ierr); 1210298fd71SBarry Smith ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr); 1220298fd71SBarry Smith ierr = PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr); 1230298fd71SBarry Smith ierr = PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr); 1240298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);CHKERRQ(ierr); 125dfbf837cSBarry Smith if (debug) { 126ce94432eSBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 127dfbf837cSBarry Smith } 1280298fd71SBarry Smith ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,NULL);CHKERRQ(ierr); 1290298fd71SBarry Smith ierr = PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,NULL);CHKERRQ(ierr); 1300298fd71SBarry Smith ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL);CHKERRQ(ierr); 1310298fd71SBarry Smith ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL);CHKERRQ(ierr); 1320298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL);CHKERRQ(ierr); 13323b3e82cSAsbjørn Nilsen Riseth 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); 134a312c225SMatthew G Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 1356a7cf640SPeter Brune if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 1369e764e56SPeter Brune 1379e764e56SPeter Brune /* set the default type of the line search if the user hasn't already. */ 1389e764e56SPeter Brune if (!snes->linesearch) { 1397601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 1401a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 1419e764e56SPeter Brune } 142a312c225SMatthew G Knepley PetscFunctionReturn(0); 143a312c225SMatthew G Knepley } 144a312c225SMatthew G Knepley 145a312c225SMatthew G Knepley #undef __FUNCT__ 146a312c225SMatthew G Knepley #define __FUNCT__ "SNESView_NGMRES" 147a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer) 148a312c225SMatthew G Knepley { 149a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 150a312c225SMatthew G Knepley PetscBool iascii; 151a312c225SMatthew G Knepley PetscErrorCode ierr; 152a312c225SMatthew G Knepley 153a312c225SMatthew G Knepley PetscFunctionBegin; 154251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 155a312c225SMatthew G Knepley if (iascii) { 156f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr); 157f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);CHKERRQ(ierr); 158f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);CHKERRQ(ierr); 15923b3e82cSAsbjørn Nilsen Riseth ierr = PetscViewerASCIIPrintf(viewer," Restart on F_M residual increase: %s\n",ngmres->restart_fm_rise?"TRUE":"FALSE");CHKERRQ(ierr); 160a312c225SMatthew G Knepley } 161a312c225SMatthew G Knepley PetscFunctionReturn(0); 162a312c225SMatthew G Knepley } 163a312c225SMatthew G Knepley 164a312c225SMatthew G Knepley #undef __FUNCT__ 165a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES" 166a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes) 167a312c225SMatthew G Knepley { 168087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 16998b3e84cSPeter Brune /* present solution, residual, and preconditioned residual */ 1709f425c49SPeter Brune Vec X,F,B,D,Y; 171f109b39eSPeter Brune 172f109b39eSPeter Brune /* candidate linear combination answers */ 173ddd40ce5SPeter Brune Vec XA,FA,XM,FM; 17419653cdaSPeter Brune 17598b3e84cSPeter Brune /* coefficients and RHS to the minimization problem */ 17618aa0c0cSPeter Brune PetscReal fnorm,fMnorm,fAnorm; 177b3c6a99cSPeter Brune PetscReal xnorm,xMnorm,xAnorm; 178b3c6a99cSPeter Brune PetscReal ynorm,yMnorm,yAnorm; 17938774f0aSPeter Brune PetscInt k,k_restart,l,ivec,restart_count = 0; 18019653cdaSPeter Brune 18198b3e84cSPeter Brune /* solution selection data */ 18238774f0aSPeter Brune PetscBool selectRestart; 18361ba4676SBarry Smith /* 18461ba4676SBarry Smith These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed 18561ba4676SBarry Smith to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case 18661ba4676SBarry Smith so the code is correct as written. 18761ba4676SBarry Smith */ 18861ba4676SBarry Smith PetscReal dnorm = 0.0,dminnorm = 0.0; 189b3c6a99cSPeter Brune PetscReal fminnorm; 19019653cdaSPeter Brune 1911e633543SBarry Smith SNESConvergedReason reason; 192422a814eSBarry Smith SNESLineSearchReason lssucceed; 193a312c225SMatthew G Knepley PetscErrorCode ierr; 194a312c225SMatthew G Knepley 195a312c225SMatthew G Knepley PetscFunctionBegin; 1966c4ed002SBarry Smith 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); 197c579b300SPatrick Farrell 198fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 19998b3e84cSPeter Brune /* variable initialization */ 200a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING; 201f109b39eSPeter Brune X = snes->vec_sol; 202f109b39eSPeter Brune F = snes->vec_func; 203f109b39eSPeter Brune B = snes->vec_rhs; 204f109b39eSPeter Brune XA = snes->vec_sol_update; 205f109b39eSPeter Brune FA = snes->work[0]; 206f109b39eSPeter Brune D = snes->work[1]; 207f109b39eSPeter Brune 208f109b39eSPeter Brune /* work for the line search */ 209f109b39eSPeter Brune Y = snes->work[2]; 2109f425c49SPeter Brune XM = snes->work[3]; 2119f425c49SPeter Brune FM = snes->work[4]; 212a312c225SMatthew G Knepley 213e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 214a312c225SMatthew G Knepley snes->iter = 0; 215a312c225SMatthew G Knepley snes->norm = 0.; 216e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 21719653cdaSPeter Brune 21898b3e84cSPeter Brune /* initialization */ 21919653cdaSPeter Brune 2203a2ae377SPeter Brune if (snes->pc && snes->pcside == PC_LEFT) { 221be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 2223a2ae377SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 2233a2ae377SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2243a2ae377SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2253a2ae377SPeter Brune PetscFunctionReturn(0); 2263a2ae377SPeter Brune } 2273a2ae377SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 2283a2ae377SPeter Brune } else { 229e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 230f109b39eSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 2311aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 232c1c75074SPeter Brune 233f109b39eSPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 234422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2353a2ae377SPeter Brune } 236e4ed7901SPeter Brune fminnorm = fnorm; 23719653cdaSPeter Brune 238e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 239f109b39eSPeter Brune snes->norm = fnorm; 240e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 241a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 242f109b39eSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 243f109b39eSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 244a312c225SMatthew G Knepley if (snes->reason) PetscFunctionReturn(0); 245b3c6a99cSPeter Brune SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X); 246a312c225SMatthew G Knepley 24719653cdaSPeter Brune k_restart = 1; 24819653cdaSPeter Brune l = 1; 249b3c6a99cSPeter Brune ivec = 0; 25009c08436SPeter Brune for (k=1; k < snes->max_its+1; k++) { 25198b3e84cSPeter Brune /* Computation of x^M */ 252c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 2539f425c49SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 254e4ed7901SPeter Brune ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr); 25563e7833aSPeter Brune 25663e7833aSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 2579f425c49SPeter Brune ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr); 25863e7833aSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 25963e7833aSPeter Brune 2608cc86e31SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 2618cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2628cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2638cc86e31SPeter Brune PetscFunctionReturn(0); 2648cc86e31SPeter Brune } 265be95d8f1SBarry Smith ierr = SNESGetNPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr); 2668cc86e31SPeter Brune } else { 267f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */ 268f109b39eSPeter Brune ierr = VecCopy(F,Y);CHKERRQ(ierr); 269e7058c64SPeter Brune ierr = VecCopy(F,FM);CHKERRQ(ierr); 270e7058c64SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 2711aa26658SKarl Rupp 272e7058c64SPeter Brune fMnorm = fnorm; 2731aa26658SKarl Rupp 274f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr); 275422a814eSBarry Smith ierr = SNESLineSearchGetReason(snes->linesearch,&lssucceed);CHKERRQ(ierr); 276422a814eSBarry Smith if (lssucceed) { 277f109b39eSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 278f109b39eSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 279f109b39eSPeter Brune PetscFunctionReturn(0); 280f109b39eSPeter Brune } 281f109b39eSPeter Brune } 2826634f59bSPeter Brune } 28323b3e82cSAsbjørn Nilsen Riseth 284b3c6a99cSPeter Brune ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr); 28598b3e84cSPeter Brune /* r = F(x) */ 2869f425c49SPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 28719653cdaSPeter Brune 2889f425c49SPeter Brune /* differences for selection and restart */ 28913a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 290b3c6a99cSPeter Brune ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);CHKERRQ(ierr); 29113a62661SPeter Brune } else { 292b3c6a99cSPeter Brune ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);CHKERRQ(ierr); 29313a62661SPeter Brune } 294422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2951aa26658SKarl Rupp 2969f425c49SPeter Brune /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 297b3c6a99cSPeter Brune 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); 29819653cdaSPeter Brune selectRestart = PETSC_FALSE; 29923b3e82cSAsbjørn Nilsen Riseth 30013a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 30123b3e82cSAsbjørn Nilsen Riseth ierr = SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 30223b3e82cSAsbjørn Nilsen Riseth 30328ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 3041aa26658SKarl Rupp if (selectRestart) restart_count++; 3051aa26658SKarl Rupp else restart_count = 0; 30613a62661SPeter Brune } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 30713a62661SPeter Brune if (k_restart > ngmres->restart_periodic) { 30813a62661SPeter Brune if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr); 30913a62661SPeter Brune restart_count = ngmres->restart_it; 31013a62661SPeter Brune } 31113a62661SPeter Brune } 31223b3e82cSAsbjørn Nilsen Riseth 313b3c6a99cSPeter Brune ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 31423b3e82cSAsbjørn Nilsen Riseth 31528ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 31628ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) { 317dfbf837cSBarry Smith if (ngmres->monitor) { 318dfbf837cSBarry Smith ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr); 319dfbf837cSBarry Smith } 32028ed4a04SPeter Brune restart_count = 0; 32119653cdaSPeter Brune k_restart = 1; 32219653cdaSPeter Brune l = 1; 323b3c6a99cSPeter Brune ivec = 0; 32498b3e84cSPeter Brune /* q_{00} = nu */ 325fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM);CHKERRQ(ierr); 326d2e16ddcSPeter Brune } else { 32798b3e84cSPeter Brune /* select the current size of the subspace */ 3281e633543SBarry Smith if (l < ngmres->msize) l++; 32919653cdaSPeter Brune k_restart++; 33098b3e84cSPeter Brune /* place the current entry in the list of previous entries */ 33138774f0aSPeter Brune if (ngmres->candidate) { 33238774f0aSPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; 333fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM);CHKERRQ(ierr); 334d2e16ddcSPeter Brune } else { 33538774f0aSPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; 336fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X);CHKERRQ(ierr); 33719653cdaSPeter Brune } 338d2e16ddcSPeter Brune } 33919653cdaSPeter Brune 340e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 341087dfb9eSxuemin snes->iter = k; 342f109b39eSPeter Brune snes->norm = fnorm; 343e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 344a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 3458409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 346b3c6a99cSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0,0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 347087dfb9eSxuemin if (snes->reason) PetscFunctionReturn(0); 348a312c225SMatthew G Knepley } 349a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 350a312c225SMatthew G Knepley PetscFunctionReturn(0); 351a312c225SMatthew G Knepley } 352a312c225SMatthew G Knepley 35313a62661SPeter Brune #undef __FUNCT__ 35423b3e82cSAsbjørn Nilsen Riseth #define __FUNCT__ "SNESNGMRESSetRestartFmRise" 35523b3e82cSAsbjørn Nilsen Riseth /*@ 35623b3e82cSAsbjørn Nilsen Riseth SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M 35723b3e82cSAsbjørn Nilsen Riseth 35823b3e82cSAsbjørn Nilsen Riseth Input Parameters: 35923b3e82cSAsbjørn Nilsen Riseth + snes - the SNES context. 36023b3e82cSAsbjørn Nilsen Riseth - flg - boolean value deciding whether to use the option or not 36123b3e82cSAsbjørn Nilsen Riseth 36223b3e82cSAsbjørn Nilsen Riseth Options Database: 36323b3e82cSAsbjørn Nilsen Riseth + -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M 36423b3e82cSAsbjørn Nilsen Riseth 36523b3e82cSAsbjørn Nilsen Riseth Level: intermediate 36623b3e82cSAsbjørn Nilsen Riseth 36723b3e82cSAsbjørn Nilsen Riseth Notes: 36823b3e82cSAsbjørn Nilsen Riseth If the proposed step x_M increases the residual F_M, it might be trying to get out of a stagnation area. 36923b3e82cSAsbjørn Nilsen Riseth To help the solver do that, reset the Krylov subspace whenever F_M increases. 37023b3e82cSAsbjørn Nilsen Riseth 37123b3e82cSAsbjørn Nilsen Riseth This option must be used with SNES_NGMRES_RESTART_DIFFERENCE 37223b3e82cSAsbjørn Nilsen Riseth 37323b3e82cSAsbjørn Nilsen Riseth The default is FALSE. 37423b3e82cSAsbjørn Nilsen Riseth .seealso: SNES_NGMRES_RESTART_DIFFERENCE 37523b3e82cSAsbjørn Nilsen Riseth @*/ 37623b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes,PetscBool flg) 37723b3e82cSAsbjørn Nilsen Riseth { 37823b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES,PetscBool); 37923b3e82cSAsbjørn Nilsen Riseth PetscErrorCode ierr; 38023b3e82cSAsbjørn Nilsen Riseth 38123b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 38223b3e82cSAsbjørn Nilsen Riseth ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",&f);CHKERRQ(ierr); 38323b3e82cSAsbjørn Nilsen Riseth if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 38423b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 38523b3e82cSAsbjørn Nilsen Riseth } 38623b3e82cSAsbjørn Nilsen Riseth 38723b3e82cSAsbjørn Nilsen Riseth #undef __FUNCT__ 38823b3e82cSAsbjørn Nilsen Riseth #define __FUNCT__ "SNESNGMRESSetRestartFmRise_NGMRES" 38923b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes,PetscBool flg) 39023b3e82cSAsbjørn Nilsen Riseth { 39123b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 39223b3e82cSAsbjørn Nilsen Riseth 39323b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 39423b3e82cSAsbjørn Nilsen Riseth ngmres->restart_fm_rise = flg; 39523b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 39623b3e82cSAsbjørn Nilsen Riseth } 39723b3e82cSAsbjørn Nilsen Riseth 39823b3e82cSAsbjørn Nilsen Riseth #undef __FUNCT__ 39923b3e82cSAsbjørn Nilsen Riseth #define __FUNCT__ "SNESNGMRESGetRestartFmRise" 40023b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes,PetscBool *flg) 40123b3e82cSAsbjørn Nilsen Riseth { 40223b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES,PetscBool*); 40323b3e82cSAsbjørn Nilsen Riseth PetscErrorCode ierr; 40423b3e82cSAsbjørn Nilsen Riseth 40523b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 40623b3e82cSAsbjørn Nilsen Riseth ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",&f);CHKERRQ(ierr); 40723b3e82cSAsbjørn Nilsen Riseth if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 40823b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 40923b3e82cSAsbjørn Nilsen Riseth } 41023b3e82cSAsbjørn Nilsen Riseth 41123b3e82cSAsbjørn Nilsen Riseth #undef __FUNCT__ 41223b3e82cSAsbjørn Nilsen Riseth #define __FUNCT__ "SNESNGMRESGetRestartFmRise_NGMRES" 41323b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes,PetscBool *flg) 41423b3e82cSAsbjørn Nilsen Riseth { 41523b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 41623b3e82cSAsbjørn Nilsen Riseth 41723b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 41823b3e82cSAsbjørn Nilsen Riseth *flg = ngmres->restart_fm_rise; 41923b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 42023b3e82cSAsbjørn Nilsen Riseth } 42123b3e82cSAsbjørn Nilsen Riseth 42223b3e82cSAsbjørn Nilsen Riseth 42323b3e82cSAsbjørn Nilsen Riseth #undef __FUNCT__ 42413a62661SPeter Brune #define __FUNCT__ "SNESNGMRESSetRestartType" 42513a62661SPeter Brune /*@ 42613a62661SPeter Brune SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES. 42713a62661SPeter Brune 42813a62661SPeter Brune Logically Collective on SNES 42913a62661SPeter Brune 43013a62661SPeter Brune Input Parameters: 43113a62661SPeter Brune + snes - the iterative context 43213a62661SPeter Brune - rtype - restart type 43313a62661SPeter Brune 43413a62661SPeter Brune Options Database: 43513a62661SPeter Brune + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type 4360c777b0cSPeter Brune - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic 43713a62661SPeter Brune 43813a62661SPeter Brune Level: intermediate 43913a62661SPeter Brune 44013a62661SPeter Brune SNESNGMRESRestartTypes: 44113a62661SPeter Brune + SNES_NGMRES_RESTART_NONE - never restart 44213a62661SPeter Brune . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria 44313a62661SPeter Brune - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations 44413a62661SPeter Brune 44513a62661SPeter Brune Notes: 44613a62661SPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 44713a62661SPeter Brune 44813a62661SPeter Brune .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch 44913a62661SPeter Brune @*/ 4500adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype) 4510adebc6cSBarry Smith { 45213a62661SPeter Brune PetscErrorCode ierr; 4535fd66863SKarl Rupp 45413a62661SPeter Brune PetscFunctionBegin; 45513a62661SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 45613a62661SPeter Brune ierr = PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));CHKERRQ(ierr); 45713a62661SPeter Brune PetscFunctionReturn(0); 45813a62661SPeter Brune } 45913a62661SPeter Brune 46013a62661SPeter Brune #undef __FUNCT__ 46113a62661SPeter Brune #define __FUNCT__ "SNESNGMRESSetSelectType" 46213a62661SPeter Brune /*@ 46313a62661SPeter Brune SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and 46413a62661SPeter Brune combined solution are used to create the next iterate. 46513a62661SPeter Brune 46613a62661SPeter Brune Logically Collective on SNES 46713a62661SPeter Brune 46813a62661SPeter Brune Input Parameters: 46913a62661SPeter Brune + snes - the iterative context 47013a62661SPeter Brune - stype - selection type 47113a62661SPeter Brune 47213a62661SPeter Brune Options Database: 47313a62661SPeter Brune . -snes_ngmres_select_type<difference,none,linesearch> 47413a62661SPeter Brune 47513a62661SPeter Brune Level: intermediate 47613a62661SPeter Brune 47713a62661SPeter Brune SNESNGMRESSelectTypes: 47813a62661SPeter Brune + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time 47913a62661SPeter Brune . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria 48013a62661SPeter Brune - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination 48113a62661SPeter Brune 48213a62661SPeter Brune Notes: 48313a62661SPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 48413a62661SPeter Brune 48513a62661SPeter Brune .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch 48613a62661SPeter Brune @*/ 4870adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype) 4880adebc6cSBarry Smith { 48913a62661SPeter Brune PetscErrorCode ierr; 4905fd66863SKarl Rupp 49113a62661SPeter Brune PetscFunctionBegin; 49213a62661SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 49313a62661SPeter Brune ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr); 49413a62661SPeter Brune PetscFunctionReturn(0); 49513a62661SPeter Brune } 49613a62661SPeter Brune 49713a62661SPeter Brune #undef __FUNCT__ 49813a62661SPeter Brune #define __FUNCT__ "SNESNGMRESSetSelectType_NGMRES" 4990adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype) 5000adebc6cSBarry Smith { 50113a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 5025fd66863SKarl Rupp 50313a62661SPeter Brune PetscFunctionBegin; 50413a62661SPeter Brune ngmres->select_type = stype; 50513a62661SPeter Brune PetscFunctionReturn(0); 50613a62661SPeter Brune } 50713a62661SPeter Brune 50813a62661SPeter Brune #undef __FUNCT__ 50913a62661SPeter Brune #define __FUNCT__ "SNESNGMRESSetRestartType_NGMRES" 5100adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype) 5110adebc6cSBarry Smith { 51213a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 5135fd66863SKarl Rupp 51413a62661SPeter Brune PetscFunctionBegin; 51513a62661SPeter Brune ngmres->restart_type = rtype; 51613a62661SPeter Brune PetscFunctionReturn(0); 51713a62661SPeter Brune } 51813a62661SPeter Brune 519dfbf837cSBarry Smith /*MC 5201867fe5bSPeter Brune SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 521a312c225SMatthew G Knepley 522dfbf837cSBarry Smith Level: beginner 523dfbf837cSBarry Smith 5241867fe5bSPeter Brune Options Database: 52513a62661SPeter Brune + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution 52638774f0aSPeter Brune . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions 52738774f0aSPeter Brune . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions 52813a62661SPeter Brune . -snes_ngmres_m - Number of stored previous solutions and residuals 52913a62661SPeter Brune . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart 53013a62661SPeter Brune . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination 53113a62661SPeter Brune . -snes_ngmres_gammaC - Residual tolerance for restart 53213a62661SPeter Brune . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart 53313a62661SPeter Brune . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart 53423b3e82cSAsbjørn Nilsen Riseth . -snes_ngmres_restart_fm_rise - Restart on residual rise from x_M step 53513a62661SPeter Brune . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration 5365c3e6ab7SPeter Brune . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother 53713a62661SPeter Brune - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type 5381867fe5bSPeter Brune 5391867fe5bSPeter Brune Notes: 5401867fe5bSPeter Brune 5411867fe5bSPeter Brune The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 5421867fe5bSPeter Brune optimization problem at each iteration. 5431867fe5bSPeter Brune 5444f02bc6aSBarry Smith Very similar to the SNESANDERSON algorithm. 5454f02bc6aSBarry Smith 5461867fe5bSPeter Brune References: 54796a0c994SBarry Smith + 1. - C. W. Oosterlee and T. Washio, "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", 548dfbf837cSBarry Smith SIAM Journal on Scientific Computing, 21(5), 2000. 54996a0c994SBarry Smith - 2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 5504f02bc6aSBarry Smith SIAM Review, 57(4), 2015 5514f02bc6aSBarry Smith 5524f02bc6aSBarry Smith 553dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 554dfbf837cSBarry Smith M*/ 555a312c225SMatthew G Knepley 556a312c225SMatthew G Knepley #undef __FUNCT__ 557a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES" 5588cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes) 559a312c225SMatthew G Knepley { 560a312c225SMatthew G Knepley SNES_NGMRES *ngmres; 561a312c225SMatthew G Knepley PetscErrorCode ierr; 562a312c225SMatthew G Knepley 563a312c225SMatthew G Knepley PetscFunctionBegin; 564a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES; 565a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES; 566a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 567a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES; 568a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES; 569a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES; 570a312c225SMatthew G Knepley 57142f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 5722c155ee1SBarry Smith snes->usesksp = PETSC_FALSE; 57346159c86SPeter Brune snes->pcside = PC_RIGHT; 5742c155ee1SBarry Smith 575*4fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 576*4fc747eaSLawrence Mitchell 577b00a9115SJed Brown ierr = PetscNewLog(snes,&ngmres);CHKERRQ(ierr); 578a312c225SMatthew G Knepley snes->data = (void*) ngmres; 579d2e16ddcSPeter Brune ngmres->msize = 30; 58019653cdaSPeter Brune 58188976e71SPeter Brune if (!snes->tolerancesset) { 5820e444f03SPeter Brune snes->max_funcs = 30000; 5830e444f03SPeter Brune snes->max_its = 10000; 58488976e71SPeter Brune } 5850e444f03SPeter Brune 58638774f0aSPeter Brune ngmres->candidate = PETSC_FALSE; 587d2e16ddcSPeter Brune 5880298fd71SBarry Smith ngmres->additive_linesearch = NULL; 589077c4231SPeter Brune ngmres->approxfunc = PETSC_FALSE; 59028ed4a04SPeter Brune ngmres->restart_it = 2; 59113a62661SPeter Brune ngmres->restart_periodic = 30; 592f109b39eSPeter Brune ngmres->gammaA = 2.0; 593f109b39eSPeter Brune ngmres->gammaC = 2.0; 594cac108bcSPeter Brune ngmres->deltaB = 0.9; 595cac108bcSPeter Brune ngmres->epsilonB = 0.1; 59623b3e82cSAsbjørn Nilsen Riseth ngmres->restart_fm_rise = PETSC_FALSE; 597e7058c64SPeter Brune 59813a62661SPeter Brune ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE; 59913a62661SPeter Brune ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE; 60013a62661SPeter Brune 601bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr); 602bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr); 60323b3e82cSAsbjørn Nilsen Riseth ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",SNESNGMRESSetRestartFmRise_NGMRES);CHKERRQ(ierr); 60423b3e82cSAsbjørn Nilsen Riseth ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",SNESNGMRESGetRestartFmRise_NGMRES);CHKERRQ(ierr); 605a312c225SMatthew G Knepley PetscFunctionReturn(0); 606a312c225SMatthew G Knepley } 607