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 59e5d0892SLisandro Dalcin const char *const SNESNGMRESRestartTypes[] = {"NONE","PERIODIC","DIFFERENCE","SNESNGMRESRestartType","SNES_NGMRES_RESTART_",NULL}; 69e5d0892SLisandro Dalcin const char *const SNESNGMRESSelectTypes[] = {"NONE","DIFFERENCE","LINESEARCH","SNESNGMRESSelectType","SNES_NGMRES_SELECT_",NULL}; 713a62661SPeter Brune 8a312c225SMatthew G Knepley PetscErrorCode SNESReset_NGMRES(SNES snes) 9a312c225SMatthew G Knepley { 10a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 11a312c225SMatthew G Knepley 12a312c225SMatthew G Knepley PetscFunctionBegin; 139566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ngmres->msize,&ngmres->Fdot)); 149566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ngmres->msize,&ngmres->Xdot)); 159566063dSJacob Faibussowitsch PetscCall(SNESLineSearchDestroy(&ngmres->additive_linesearch)); 16a312c225SMatthew G Knepley PetscFunctionReturn(0); 17a312c225SMatthew G Knepley } 18a312c225SMatthew G Knepley 19a312c225SMatthew G Knepley PetscErrorCode SNESDestroy_NGMRES(SNES snes) 20a312c225SMatthew G Knepley { 2178440776SJed Brown SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 22a312c225SMatthew G Knepley 23a312c225SMatthew G Knepley PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(SNESReset_NGMRES(snes)); 259566063dSJacob Faibussowitsch PetscCall(PetscFree4(ngmres->h,ngmres->beta,ngmres->xi,ngmres->q)); 269566063dSJacob Faibussowitsch PetscCall(PetscFree3(ngmres->xnorms,ngmres->fnorms,ngmres->s)); 27dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX) 289566063dSJacob Faibussowitsch PetscCall(PetscFree(ngmres->rwork)); 2919653cdaSPeter Brune #endif 309566063dSJacob Faibussowitsch PetscCall(PetscFree(ngmres->work)); 319566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 32a312c225SMatthew G Knepley PetscFunctionReturn(0); 33a312c225SMatthew G Knepley } 34a312c225SMatthew G Knepley 35a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes) 36a312c225SMatthew G Knepley { 37a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 38e7058c64SPeter Brune const char *optionsprefix; 3919653cdaSPeter Brune PetscInt msize,hsize; 40fced5a79SAsbjørn Nilsen Riseth DM dm; 41a312c225SMatthew G Knepley 42a312c225SMatthew G Knepley PetscFunctionBegin; 43efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 4446159c86SPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNGMRES does not support left preconditioning with unpreconditioned function"); 4546159c86SPeter Brune } 46efd4aadfSBarry Smith if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 479566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes,5)); 48fced5a79SAsbjørn Nilsen Riseth 49fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) { 509566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 519566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm,&snes->vec_sol)); 52fced5a79SAsbjørn Nilsen Riseth } 53fced5a79SAsbjørn Nilsen Riseth 549566063dSJacob Faibussowitsch if (!ngmres->Xdot) PetscCall(VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot)); 559566063dSJacob Faibussowitsch if (!ngmres->Fdot) PetscCall(VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot)); 5678440776SJed Brown if (!ngmres->setup_called) { 57087dfb9eSxuemin msize = ngmres->msize; /* restart size */ 5819653cdaSPeter Brune hsize = msize * msize; 59087dfb9eSxuemin 6098b3e84cSPeter Brune /* explicit least squares minimization solve */ 619566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(hsize,&ngmres->h, msize,&ngmres->beta, msize,&ngmres->xi, hsize,&ngmres->q)); 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(msize,&ngmres->xnorms,msize,&ngmres->fnorms,msize,&ngmres->s)); 6319653cdaSPeter Brune ngmres->nrhs = 1; 6419653cdaSPeter Brune ngmres->lda = msize; 6519653cdaSPeter Brune ngmres->ldb = msize; 6619653cdaSPeter Brune ngmres->lwork = 12*msize; 67dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX) 689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ngmres->lwork,&ngmres->rwork)); 6919653cdaSPeter Brune #endif 709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ngmres->lwork,&ngmres->work)); 7178440776SJed Brown } 72e7058c64SPeter Brune 73e7058c64SPeter Brune /* linesearch setup */ 749566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes,&optionsprefix)); 75e7058c64SPeter Brune 7613a62661SPeter Brune if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 779566063dSJacob Faibussowitsch PetscCall(SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch)); 789566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetSNES(ngmres->additive_linesearch,snes)); 79ec786807SJed Brown if (!((PetscObject)ngmres->additive_linesearch)->type_name) { 809566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2)); 81ec786807SJed Brown } 829566063dSJacob Faibussowitsch PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_")); 839566063dSJacob Faibussowitsch PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix)); 849566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetFromOptions(ngmres->additive_linesearch)); 85e7058c64SPeter Brune } 86e7058c64SPeter Brune 8778440776SJed Brown ngmres->setup_called = PETSC_TRUE; 88a312c225SMatthew G Knepley PetscFunctionReturn(0); 89a312c225SMatthew G Knepley } 90a312c225SMatthew G Knepley 914416b707SBarry Smith PetscErrorCode SNESSetFromOptions_NGMRES(PetscOptionItems *PetscOptionsObject,SNES snes) 92a312c225SMatthew G Knepley { 93a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 94a312c225SMatthew G Knepley PetscErrorCode ierr; 9594ae4db5SBarry Smith PetscBool debug = PETSC_FALSE; 960adebc6cSBarry Smith 97a312c225SMatthew G Knepley PetscFunctionBegin; 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options")); 9913a62661SPeter Brune ierr = PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes, 1009566063dSJacob Faibussowitsch (PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL);PetscCall(ierr); 10113a62661SPeter Brune ierr = PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes, 1029566063dSJacob Faibussowitsch (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);PetscCall(ierr); 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL)); 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,NULL)); 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL)); 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL)); 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL)); 109dfbf837cSBarry Smith if (debug) { 1102f613bf5SBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); 111dfbf837cSBarry Smith } 1129566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,NULL)); 1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,NULL)); 1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL)); 1159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL)); 1169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL)); 1179566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise",ngmres->restart_fm_rise,&ngmres->restart_fm_rise,NULL)); 1189566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 1196a7cf640SPeter Brune if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 120a312c225SMatthew G Knepley PetscFunctionReturn(0); 121a312c225SMatthew G Knepley } 122a312c225SMatthew G Knepley 123a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer) 124a312c225SMatthew G Knepley { 125a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 126a312c225SMatthew G Knepley PetscBool iascii; 127a312c225SMatthew G Knepley 128a312c225SMatthew G Knepley PetscFunctionBegin; 1299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii)); 130a312c225SMatthew G Knepley if (iascii) { 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Number of stored past updates: %d\n", ngmres->msize)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB)); 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Restart on F_M residual increase: %s\n",ngmres->restart_fm_rise?"TRUE":"FALSE")); 135a312c225SMatthew G Knepley } 136a312c225SMatthew G Knepley PetscFunctionReturn(0); 137a312c225SMatthew G Knepley } 138a312c225SMatthew G Knepley 139a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes) 140a312c225SMatthew G Knepley { 141087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 14298b3e84cSPeter Brune /* present solution, residual, and preconditioned residual */ 1439f425c49SPeter Brune Vec X,F,B,D,Y; 144f109b39eSPeter Brune 145f109b39eSPeter Brune /* candidate linear combination answers */ 146ddd40ce5SPeter Brune Vec XA,FA,XM,FM; 14719653cdaSPeter Brune 14898b3e84cSPeter Brune /* coefficients and RHS to the minimization problem */ 14918aa0c0cSPeter Brune PetscReal fnorm,fMnorm,fAnorm; 150b3c6a99cSPeter Brune PetscReal xnorm,xMnorm,xAnorm; 151b3c6a99cSPeter Brune PetscReal ynorm,yMnorm,yAnorm; 15238774f0aSPeter Brune PetscInt k,k_restart,l,ivec,restart_count = 0; 15319653cdaSPeter Brune 15498b3e84cSPeter Brune /* solution selection data */ 15538774f0aSPeter Brune PetscBool selectRestart; 15661ba4676SBarry Smith /* 15761ba4676SBarry Smith These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed 15861ba4676SBarry Smith to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case 15961ba4676SBarry Smith so the code is correct as written. 16061ba4676SBarry Smith */ 16161ba4676SBarry Smith PetscReal dnorm = 0.0,dminnorm = 0.0; 162b3c6a99cSPeter Brune PetscReal fminnorm; 16319653cdaSPeter Brune 1641e633543SBarry Smith SNESConvergedReason reason; 165422a814eSBarry Smith SNESLineSearchReason lssucceed; 166a312c225SMatthew G Knepley 167a312c225SMatthew G Knepley PetscFunctionBegin; 1682c71b3e2SJacob Faibussowitsch PetscCheckFalse(snes->xl || snes->xu || snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 169c579b300SPatrick Farrell 1709566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 17198b3e84cSPeter Brune /* variable initialization */ 172a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING; 173f109b39eSPeter Brune X = snes->vec_sol; 174f109b39eSPeter Brune F = snes->vec_func; 175f109b39eSPeter Brune B = snes->vec_rhs; 176f109b39eSPeter Brune XA = snes->vec_sol_update; 177f109b39eSPeter Brune FA = snes->work[0]; 178f109b39eSPeter Brune D = snes->work[1]; 179f109b39eSPeter Brune 180f109b39eSPeter Brune /* work for the line search */ 181f109b39eSPeter Brune Y = snes->work[2]; 1829f425c49SPeter Brune XM = snes->work[3]; 1839f425c49SPeter Brune FM = snes->work[4]; 184a312c225SMatthew G Knepley 1859566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 186a312c225SMatthew G Knepley snes->iter = 0; 187a312c225SMatthew G Knepley snes->norm = 0.; 1889566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 18919653cdaSPeter Brune 19098b3e84cSPeter Brune /* initialization */ 19119653cdaSPeter Brune 192efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT) { 1939566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes,X,NULL,F)); 1949566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 1953a2ae377SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1963a2ae377SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1973a2ae377SPeter Brune PetscFunctionReturn(0); 1983a2ae377SPeter Brune } 1999566063dSJacob Faibussowitsch PetscCall(VecNorm(F,NORM_2,&fnorm)); 2003a2ae377SPeter Brune } else { 201e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2029566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 2031aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 204c1c75074SPeter Brune 2059566063dSJacob Faibussowitsch PetscCall(VecNorm(F,NORM_2,&fnorm)); 206422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2073a2ae377SPeter Brune } 208e4ed7901SPeter Brune fminnorm = fnorm; 20919653cdaSPeter Brune 2109566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 211f109b39eSPeter Brune snes->norm = fnorm; 2129566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2139566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,fnorm,0)); 2149566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,0,fnorm)); 2159566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 216a312c225SMatthew G Knepley if (snes->reason) PetscFunctionReturn(0); 217b3c6a99cSPeter Brune SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X); 218a312c225SMatthew G Knepley 21919653cdaSPeter Brune k_restart = 1; 22019653cdaSPeter Brune l = 1; 221b3c6a99cSPeter Brune ivec = 0; 22209c08436SPeter Brune for (k=1; k < snes->max_its+1; k++) { 22398b3e84cSPeter Brune /* Computation of x^M */ 224efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_RIGHT) { 2259566063dSJacob Faibussowitsch PetscCall(VecCopy(X,XM)); 2269566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(snes->npc,F)); 22763e7833aSPeter Brune 2289566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0)); 2299566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc,B,XM)); 2309566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0)); 23163e7833aSPeter Brune 2329566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 2338cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2348cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2358cc86e31SPeter Brune PetscFunctionReturn(0); 2368cc86e31SPeter Brune } 2379566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes,FM,&fMnorm)); 2388cc86e31SPeter Brune } else { 239f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */ 2409566063dSJacob Faibussowitsch PetscCall(VecCopy(F,Y)); 2419566063dSJacob Faibussowitsch PetscCall(VecCopy(F,FM)); 2429566063dSJacob Faibussowitsch PetscCall(VecCopy(X,XM)); 2431aa26658SKarl Rupp 244e7058c64SPeter Brune fMnorm = fnorm; 2451aa26658SKarl Rupp 2469566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y)); 2479566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch,&lssucceed)); 248422a814eSBarry Smith if (lssucceed) { 249f109b39eSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 250f109b39eSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 251f109b39eSPeter Brune PetscFunctionReturn(0); 252f109b39eSPeter Brune } 253f109b39eSPeter Brune } 2546634f59bSPeter Brune } 25523b3e82cSAsbjørn Nilsen Riseth 2569566063dSJacob Faibussowitsch PetscCall(SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA)); 25798b3e84cSPeter Brune /* r = F(x) */ 2589f425c49SPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 25919653cdaSPeter Brune 2609f425c49SPeter Brune /* differences for selection and restart */ 26113a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 2629566063dSJacob Faibussowitsch PetscCall(SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm)); 26313a62661SPeter Brune } else { 2649566063dSJacob Faibussowitsch PetscCall(SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm)); 26513a62661SPeter Brune } 266422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2671aa26658SKarl Rupp 2689f425c49SPeter Brune /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 2699566063dSJacob Faibussowitsch PetscCall(SNESNGMRESSelect_Private(snes,k_restart,XM,FM,xMnorm,fMnorm,yMnorm,XA,FA,xAnorm,fAnorm,yAnorm,dnorm,fminnorm,dminnorm,X,F,Y,&xnorm,&fnorm,&ynorm)); 27019653cdaSPeter Brune selectRestart = PETSC_FALSE; 27123b3e82cSAsbjørn Nilsen Riseth 27213a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 2739566063dSJacob Faibussowitsch PetscCall(SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart)); 27423b3e82cSAsbjørn Nilsen Riseth 27528ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 2761aa26658SKarl Rupp if (selectRestart) restart_count++; 2771aa26658SKarl Rupp else restart_count = 0; 27813a62661SPeter Brune } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 27913a62661SPeter Brune if (k_restart > ngmres->restart_periodic) { 2809566063dSJacob Faibussowitsch if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart)); 28113a62661SPeter Brune restart_count = ngmres->restart_it; 28213a62661SPeter Brune } 28313a62661SPeter Brune } 28423b3e82cSAsbjørn Nilsen Riseth 285b3c6a99cSPeter Brune ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 28623b3e82cSAsbjørn Nilsen Riseth 28728ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 28828ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) { 289dfbf837cSBarry Smith if (ngmres->monitor) { 2909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart)); 291dfbf837cSBarry Smith } 29228ed4a04SPeter Brune restart_count = 0; 29319653cdaSPeter Brune k_restart = 1; 29419653cdaSPeter Brune l = 1; 295b3c6a99cSPeter Brune ivec = 0; 29698b3e84cSPeter Brune /* q_{00} = nu */ 2979566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM)); 298d2e16ddcSPeter Brune } else { 29998b3e84cSPeter Brune /* select the current size of the subspace */ 3001e633543SBarry Smith if (l < ngmres->msize) l++; 30119653cdaSPeter Brune k_restart++; 30298b3e84cSPeter Brune /* place the current entry in the list of previous entries */ 30338774f0aSPeter Brune if (ngmres->candidate) { 30438774f0aSPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; 3059566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM)); 306d2e16ddcSPeter Brune } else { 30738774f0aSPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; 3089566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X)); 30919653cdaSPeter Brune } 310d2e16ddcSPeter Brune } 31119653cdaSPeter Brune 3129566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 313087dfb9eSxuemin snes->iter = k; 314f109b39eSPeter Brune snes->norm = fnorm; 3159566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3169566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,snes->iter)); 3179566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 3189566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,snes->iter,0,0,fnorm,&snes->reason,snes->cnvP)); 319087dfb9eSxuemin if (snes->reason) PetscFunctionReturn(0); 320a312c225SMatthew G Knepley } 321a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 322a312c225SMatthew G Knepley PetscFunctionReturn(0); 323a312c225SMatthew G Knepley } 324a312c225SMatthew G Knepley 32523b3e82cSAsbjørn Nilsen Riseth /*@ 32623b3e82cSAsbjørn Nilsen Riseth SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M 32723b3e82cSAsbjørn Nilsen Riseth 32823b3e82cSAsbjørn Nilsen Riseth Input Parameters: 32923b3e82cSAsbjørn Nilsen Riseth + snes - the SNES context. 33023b3e82cSAsbjørn Nilsen Riseth - flg - boolean value deciding whether to use the option or not 33123b3e82cSAsbjørn Nilsen Riseth 33223b3e82cSAsbjørn Nilsen Riseth Options Database: 33323b3e82cSAsbjørn Nilsen Riseth + -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M 33423b3e82cSAsbjørn Nilsen Riseth 33523b3e82cSAsbjørn Nilsen Riseth Level: intermediate 33623b3e82cSAsbjørn Nilsen Riseth 33723b3e82cSAsbjørn Nilsen Riseth Notes: 33823b3e82cSAsbjø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. 33923b3e82cSAsbjørn Nilsen Riseth To help the solver do that, reset the Krylov subspace whenever F_M increases. 34023b3e82cSAsbjørn Nilsen Riseth 34123b3e82cSAsbjørn Nilsen Riseth This option must be used with SNES_NGMRES_RESTART_DIFFERENCE 34223b3e82cSAsbjørn Nilsen Riseth 34323b3e82cSAsbjørn Nilsen Riseth The default is FALSE. 34423b3e82cSAsbjørn Nilsen Riseth .seealso: SNES_NGMRES_RESTART_DIFFERENCE 34523b3e82cSAsbjørn Nilsen Riseth @*/ 34623b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes,PetscBool flg) 34723b3e82cSAsbjørn Nilsen Riseth { 34823b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES,PetscBool); 34923b3e82cSAsbjørn Nilsen Riseth 35023b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 3519566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",&f)); 3529566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,flg)); 35323b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 35423b3e82cSAsbjørn Nilsen Riseth } 35523b3e82cSAsbjørn Nilsen Riseth 35623b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes,PetscBool flg) 35723b3e82cSAsbjørn Nilsen Riseth { 35823b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 35923b3e82cSAsbjørn Nilsen Riseth 36023b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 36123b3e82cSAsbjørn Nilsen Riseth ngmres->restart_fm_rise = flg; 36223b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 36323b3e82cSAsbjørn Nilsen Riseth } 36423b3e82cSAsbjørn Nilsen Riseth 36523b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes,PetscBool *flg) 36623b3e82cSAsbjørn Nilsen Riseth { 36723b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES,PetscBool*); 36823b3e82cSAsbjørn Nilsen Riseth 36923b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 3709566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",&f)); 3719566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,flg)); 37223b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 37323b3e82cSAsbjørn Nilsen Riseth } 37423b3e82cSAsbjørn Nilsen Riseth 37523b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes,PetscBool *flg) 37623b3e82cSAsbjørn Nilsen Riseth { 37723b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 37823b3e82cSAsbjørn Nilsen Riseth 37923b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 38023b3e82cSAsbjørn Nilsen Riseth *flg = ngmres->restart_fm_rise; 38123b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 38223b3e82cSAsbjørn Nilsen Riseth } 38323b3e82cSAsbjørn Nilsen Riseth 38413a62661SPeter Brune /*@ 38513a62661SPeter Brune SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES. 38613a62661SPeter Brune 38713a62661SPeter Brune Logically Collective on SNES 38813a62661SPeter Brune 38913a62661SPeter Brune Input Parameters: 39013a62661SPeter Brune + snes - the iterative context 39113a62661SPeter Brune - rtype - restart type 39213a62661SPeter Brune 39313a62661SPeter Brune Options Database: 39413a62661SPeter Brune + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type 3950c777b0cSPeter Brune - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic 39613a62661SPeter Brune 39713a62661SPeter Brune Level: intermediate 39813a62661SPeter Brune 39913a62661SPeter Brune SNESNGMRESRestartTypes: 40013a62661SPeter Brune + SNES_NGMRES_RESTART_NONE - never restart 40113a62661SPeter Brune . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria 40213a62661SPeter Brune - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations 40313a62661SPeter Brune 40413a62661SPeter Brune Notes: 40513a62661SPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 40613a62661SPeter Brune 40713a62661SPeter Brune @*/ 4080adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype) 4090adebc6cSBarry Smith { 41013a62661SPeter Brune PetscFunctionBegin; 41113a62661SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 412*cac4c232SBarry Smith PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype)); 41313a62661SPeter Brune PetscFunctionReturn(0); 41413a62661SPeter Brune } 41513a62661SPeter Brune 41613a62661SPeter Brune /*@ 41713a62661SPeter Brune SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and 41813a62661SPeter Brune combined solution are used to create the next iterate. 41913a62661SPeter Brune 42013a62661SPeter Brune Logically Collective on SNES 42113a62661SPeter Brune 42213a62661SPeter Brune Input Parameters: 42313a62661SPeter Brune + snes - the iterative context 42413a62661SPeter Brune - stype - selection type 42513a62661SPeter Brune 42613a62661SPeter Brune Options Database: 42767b8a455SSatish Balay . -snes_ngmres_select_type<difference,none,linesearch> - select type 42813a62661SPeter Brune 42913a62661SPeter Brune Level: intermediate 43013a62661SPeter Brune 43113a62661SPeter Brune SNESNGMRESSelectTypes: 43213a62661SPeter Brune + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time 43313a62661SPeter Brune . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria 43413a62661SPeter Brune - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination 43513a62661SPeter Brune 43613a62661SPeter Brune Notes: 43713a62661SPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 43813a62661SPeter Brune 43913a62661SPeter Brune @*/ 4400adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype) 4410adebc6cSBarry Smith { 44213a62661SPeter Brune PetscFunctionBegin; 44313a62661SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 444*cac4c232SBarry Smith PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype)); 44513a62661SPeter Brune PetscFunctionReturn(0); 44613a62661SPeter Brune } 44713a62661SPeter Brune 4480adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype) 4490adebc6cSBarry Smith { 45013a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 4515fd66863SKarl Rupp 45213a62661SPeter Brune PetscFunctionBegin; 45313a62661SPeter Brune ngmres->select_type = stype; 45413a62661SPeter Brune PetscFunctionReturn(0); 45513a62661SPeter Brune } 45613a62661SPeter Brune 4570adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype) 4580adebc6cSBarry Smith { 45913a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 4605fd66863SKarl Rupp 46113a62661SPeter Brune PetscFunctionBegin; 46213a62661SPeter Brune ngmres->restart_type = rtype; 46313a62661SPeter Brune PetscFunctionReturn(0); 46413a62661SPeter Brune } 46513a62661SPeter Brune 466dfbf837cSBarry Smith /*MC 4671867fe5bSPeter Brune SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 468a312c225SMatthew G Knepley 469dfbf837cSBarry Smith Level: beginner 470dfbf837cSBarry Smith 4711867fe5bSPeter Brune Options Database: 47213a62661SPeter Brune + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution 47338774f0aSPeter Brune . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions 47438774f0aSPeter Brune . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions 47513a62661SPeter Brune . -snes_ngmres_m - Number of stored previous solutions and residuals 47613a62661SPeter Brune . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart 47713a62661SPeter Brune . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination 47813a62661SPeter Brune . -snes_ngmres_gammaC - Residual tolerance for restart 47913a62661SPeter Brune . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart 48013a62661SPeter Brune . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart 48123b3e82cSAsbjørn Nilsen Riseth . -snes_ngmres_restart_fm_rise - Restart on residual rise from x_M step 48213a62661SPeter Brune . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration 4835c3e6ab7SPeter Brune . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother 48413a62661SPeter Brune - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type 4851867fe5bSPeter Brune 4861867fe5bSPeter Brune Notes: 4871867fe5bSPeter Brune 4881867fe5bSPeter Brune The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 4891867fe5bSPeter Brune optimization problem at each iteration. 4901867fe5bSPeter Brune 4914f02bc6aSBarry Smith Very similar to the SNESANDERSON algorithm. 4924f02bc6aSBarry Smith 4931867fe5bSPeter Brune References: 494606c0280SSatish Balay + * - C. W. Oosterlee and T. Washio, "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", 495dfbf837cSBarry Smith SIAM Journal on Scientific Computing, 21(5), 2000. 496606c0280SSatish Balay - * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 4974f02bc6aSBarry Smith SIAM Review, 57(4), 2015 4984f02bc6aSBarry Smith 499dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 500dfbf837cSBarry Smith M*/ 501a312c225SMatthew G Knepley 5028cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes) 503a312c225SMatthew G Knepley { 504a312c225SMatthew G Knepley SNES_NGMRES *ngmres; 505d8d34be6SBarry Smith SNESLineSearch linesearch; 506a312c225SMatthew G Knepley 507a312c225SMatthew G Knepley PetscFunctionBegin; 508a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES; 509a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES; 510a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 511a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES; 512a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES; 513a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES; 514a312c225SMatthew G Knepley 515efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 5162c155ee1SBarry Smith snes->usesksp = PETSC_FALSE; 517efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 5182c155ee1SBarry Smith 5194fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 5204fc747eaSLawrence Mitchell 5219566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes,&ngmres)); 522a312c225SMatthew G Knepley snes->data = (void*) ngmres; 523d2e16ddcSPeter Brune ngmres->msize = 30; 52419653cdaSPeter Brune 52588976e71SPeter Brune if (!snes->tolerancesset) { 5260e444f03SPeter Brune snes->max_funcs = 30000; 5270e444f03SPeter Brune snes->max_its = 10000; 52888976e71SPeter Brune } 5290e444f03SPeter Brune 53038774f0aSPeter Brune ngmres->candidate = PETSC_FALSE; 531d2e16ddcSPeter Brune 5329566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes,&linesearch)); 533ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 5349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC)); 535ec786807SJed Brown } 536d8d34be6SBarry Smith 5370298fd71SBarry Smith ngmres->additive_linesearch = NULL; 538077c4231SPeter Brune ngmres->approxfunc = PETSC_FALSE; 53928ed4a04SPeter Brune ngmres->restart_it = 2; 54013a62661SPeter Brune ngmres->restart_periodic = 30; 541f109b39eSPeter Brune ngmres->gammaA = 2.0; 542f109b39eSPeter Brune ngmres->gammaC = 2.0; 543cac108bcSPeter Brune ngmres->deltaB = 0.9; 544cac108bcSPeter Brune ngmres->epsilonB = 0.1; 54523b3e82cSAsbjørn Nilsen Riseth ngmres->restart_fm_rise = PETSC_FALSE; 546e7058c64SPeter Brune 54713a62661SPeter Brune ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE; 54813a62661SPeter Brune ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE; 54913a62661SPeter Brune 5509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES)); 5519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES)); 5529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",SNESNGMRESSetRestartFmRise_NGMRES)); 5539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",SNESNGMRESGetRestartFmRise_NGMRES)); 554a312c225SMatthew G Knepley PetscFunctionReturn(0); 555a312c225SMatthew G Knepley } 556