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)); 31*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",NULL)); 32*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",NULL)); 33*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",NULL)); 34*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",NULL)); 359566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 36a312c225SMatthew G Knepley PetscFunctionReturn(0); 37a312c225SMatthew G Knepley } 38a312c225SMatthew G Knepley 39a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes) 40a312c225SMatthew G Knepley { 41a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 42e7058c64SPeter Brune const char *optionsprefix; 4319653cdaSPeter Brune PetscInt msize,hsize; 44fced5a79SAsbjørn Nilsen Riseth DM dm; 45a312c225SMatthew G Knepley 46a312c225SMatthew G Knepley PetscFunctionBegin; 47efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 4846159c86SPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNGMRES does not support left preconditioning with unpreconditioned function"); 4946159c86SPeter Brune } 50efd4aadfSBarry Smith if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 519566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes,5)); 52fced5a79SAsbjørn Nilsen Riseth 53fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) { 549566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 559566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm,&snes->vec_sol)); 56fced5a79SAsbjørn Nilsen Riseth } 57fced5a79SAsbjørn Nilsen Riseth 589566063dSJacob Faibussowitsch if (!ngmres->Xdot) PetscCall(VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot)); 599566063dSJacob Faibussowitsch if (!ngmres->Fdot) PetscCall(VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot)); 6078440776SJed Brown if (!ngmres->setup_called) { 61087dfb9eSxuemin msize = ngmres->msize; /* restart size */ 6219653cdaSPeter Brune hsize = msize * msize; 63087dfb9eSxuemin 6498b3e84cSPeter Brune /* explicit least squares minimization solve */ 659566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(hsize,&ngmres->h, msize,&ngmres->beta, msize,&ngmres->xi, hsize,&ngmres->q)); 669566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(msize,&ngmres->xnorms,msize,&ngmres->fnorms,msize,&ngmres->s)); 6719653cdaSPeter Brune ngmres->nrhs = 1; 6819653cdaSPeter Brune ngmres->lda = msize; 6919653cdaSPeter Brune ngmres->ldb = msize; 7019653cdaSPeter Brune ngmres->lwork = 12*msize; 71dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX) 729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ngmres->lwork,&ngmres->rwork)); 7319653cdaSPeter Brune #endif 749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ngmres->lwork,&ngmres->work)); 7578440776SJed Brown } 76e7058c64SPeter Brune 77e7058c64SPeter Brune /* linesearch setup */ 789566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes,&optionsprefix)); 79e7058c64SPeter Brune 8013a62661SPeter Brune if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 819566063dSJacob Faibussowitsch PetscCall(SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch)); 829566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetSNES(ngmres->additive_linesearch,snes)); 83ec786807SJed Brown if (!((PetscObject)ngmres->additive_linesearch)->type_name) { 849566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2)); 85ec786807SJed Brown } 869566063dSJacob Faibussowitsch PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_")); 879566063dSJacob Faibussowitsch PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix)); 889566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetFromOptions(ngmres->additive_linesearch)); 89e7058c64SPeter Brune } 90e7058c64SPeter Brune 9178440776SJed Brown ngmres->setup_called = PETSC_TRUE; 92a312c225SMatthew G Knepley PetscFunctionReturn(0); 93a312c225SMatthew G Knepley } 94a312c225SMatthew G Knepley 954416b707SBarry Smith PetscErrorCode SNESSetFromOptions_NGMRES(PetscOptionItems *PetscOptionsObject,SNES snes) 96a312c225SMatthew G Knepley { 97a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 9894ae4db5SBarry Smith PetscBool debug = PETSC_FALSE; 990adebc6cSBarry Smith 100a312c225SMatthew G Knepley PetscFunctionBegin; 101d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"SNES NGMRES options"); 102d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes,(PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL)); 103d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,(PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,NULL)); 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL)); 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,NULL)); 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL)); 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL)); 1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL)); 110dfbf837cSBarry Smith if (debug) { 1112f613bf5SBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); 112dfbf837cSBarry Smith } 1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,NULL)); 1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,NULL)); 1159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL)); 1169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL)); 1179566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL)); 1189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise",ngmres->restart_fm_rise,&ngmres->restart_fm_rise,NULL)); 119d0609cedSBarry Smith PetscOptionsHeadEnd(); 1206a7cf640SPeter Brune if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 121a312c225SMatthew G Knepley PetscFunctionReturn(0); 122a312c225SMatthew G Knepley } 123a312c225SMatthew G Knepley 124a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer) 125a312c225SMatthew G Knepley { 126a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 127a312c225SMatthew G Knepley PetscBool iascii; 128a312c225SMatthew G Knepley 129a312c225SMatthew G Knepley PetscFunctionBegin; 1309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii)); 131a312c225SMatthew G Knepley if (iascii) { 13263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Number of stored past updates: %" PetscInt_FMT "\n", ngmres->msize)); 13363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",(double)ngmres->gammaA,(double)ngmres->gammaC)); 13463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",(double)ngmres->epsilonB,(double)ngmres->deltaB)); 13563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Restart on F_M residual increase: %s\n",PetscBools[ngmres->restart_fm_rise])); 136a312c225SMatthew G Knepley } 137a312c225SMatthew G Knepley PetscFunctionReturn(0); 138a312c225SMatthew G Knepley } 139a312c225SMatthew G Knepley 140a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes) 141a312c225SMatthew G Knepley { 142087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 14398b3e84cSPeter Brune /* present solution, residual, and preconditioned residual */ 1449f425c49SPeter Brune Vec X,F,B,D,Y; 145f109b39eSPeter Brune 146f109b39eSPeter Brune /* candidate linear combination answers */ 147ddd40ce5SPeter Brune Vec XA,FA,XM,FM; 14819653cdaSPeter Brune 14998b3e84cSPeter Brune /* coefficients and RHS to the minimization problem */ 15018aa0c0cSPeter Brune PetscReal fnorm,fMnorm,fAnorm; 151b3c6a99cSPeter Brune PetscReal xnorm,xMnorm,xAnorm; 152b3c6a99cSPeter Brune PetscReal ynorm,yMnorm,yAnorm; 15338774f0aSPeter Brune PetscInt k,k_restart,l,ivec,restart_count = 0; 15419653cdaSPeter Brune 15598b3e84cSPeter Brune /* solution selection data */ 15638774f0aSPeter Brune PetscBool selectRestart; 15761ba4676SBarry Smith /* 15861ba4676SBarry Smith These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed 15961ba4676SBarry Smith to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case 16061ba4676SBarry Smith so the code is correct as written. 16161ba4676SBarry Smith */ 16261ba4676SBarry Smith PetscReal dnorm = 0.0,dminnorm = 0.0; 163b3c6a99cSPeter Brune PetscReal fminnorm; 16419653cdaSPeter Brune 1651e633543SBarry Smith SNESConvergedReason reason; 166422a814eSBarry Smith SNESLineSearchReason lssucceed; 167a312c225SMatthew G Knepley 168a312c225SMatthew G Knepley PetscFunctionBegin; 1690b121fc5SBarry Smith PetscCheck(!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); 170c579b300SPatrick Farrell 1719566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 17298b3e84cSPeter Brune /* variable initialization */ 173a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING; 174f109b39eSPeter Brune X = snes->vec_sol; 175f109b39eSPeter Brune F = snes->vec_func; 176f109b39eSPeter Brune B = snes->vec_rhs; 177f109b39eSPeter Brune XA = snes->vec_sol_update; 178f109b39eSPeter Brune FA = snes->work[0]; 179f109b39eSPeter Brune D = snes->work[1]; 180f109b39eSPeter Brune 181f109b39eSPeter Brune /* work for the line search */ 182f109b39eSPeter Brune Y = snes->work[2]; 1839f425c49SPeter Brune XM = snes->work[3]; 1849f425c49SPeter Brune FM = snes->work[4]; 185a312c225SMatthew G Knepley 1869566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 187a312c225SMatthew G Knepley snes->iter = 0; 188a312c225SMatthew G Knepley snes->norm = 0.; 1899566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 19019653cdaSPeter Brune 19198b3e84cSPeter Brune /* initialization */ 19219653cdaSPeter Brune 193efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT) { 1949566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes,X,NULL,F)); 1959566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 1963a2ae377SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1973a2ae377SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1983a2ae377SPeter Brune PetscFunctionReturn(0); 1993a2ae377SPeter Brune } 2009566063dSJacob Faibussowitsch PetscCall(VecNorm(F,NORM_2,&fnorm)); 2013a2ae377SPeter Brune } else { 202e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2039566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 2041aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 205c1c75074SPeter Brune 2069566063dSJacob Faibussowitsch PetscCall(VecNorm(F,NORM_2,&fnorm)); 207422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2083a2ae377SPeter Brune } 209e4ed7901SPeter Brune fminnorm = fnorm; 21019653cdaSPeter Brune 2119566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 212f109b39eSPeter Brune snes->norm = fnorm; 2139566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2149566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,fnorm,0)); 2159566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,0,fnorm)); 2169566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 217a312c225SMatthew G Knepley if (snes->reason) PetscFunctionReturn(0); 218b3c6a99cSPeter Brune SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X); 219a312c225SMatthew G Knepley 22019653cdaSPeter Brune k_restart = 1; 22119653cdaSPeter Brune l = 1; 222b3c6a99cSPeter Brune ivec = 0; 22309c08436SPeter Brune for (k=1; k < snes->max_its+1; k++) { 22498b3e84cSPeter Brune /* Computation of x^M */ 225efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_RIGHT) { 2269566063dSJacob Faibussowitsch PetscCall(VecCopy(X,XM)); 2279566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(snes->npc,F)); 22863e7833aSPeter Brune 2299566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0)); 2309566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc,B,XM)); 2319566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0)); 23263e7833aSPeter Brune 2339566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc,&reason)); 2348cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2358cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2368cc86e31SPeter Brune PetscFunctionReturn(0); 2378cc86e31SPeter Brune } 2389566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes,FM,&fMnorm)); 2398cc86e31SPeter Brune } else { 240f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */ 2419566063dSJacob Faibussowitsch PetscCall(VecCopy(F,Y)); 2429566063dSJacob Faibussowitsch PetscCall(VecCopy(F,FM)); 2439566063dSJacob Faibussowitsch PetscCall(VecCopy(X,XM)); 2441aa26658SKarl Rupp 245e7058c64SPeter Brune fMnorm = fnorm; 2461aa26658SKarl Rupp 2479566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y)); 2489566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch,&lssucceed)); 249422a814eSBarry Smith if (lssucceed) { 250f109b39eSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 251f109b39eSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 252f109b39eSPeter Brune PetscFunctionReturn(0); 253f109b39eSPeter Brune } 254f109b39eSPeter Brune } 2556634f59bSPeter Brune } 25623b3e82cSAsbjørn Nilsen Riseth 2579566063dSJacob Faibussowitsch PetscCall(SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA)); 25898b3e84cSPeter Brune /* r = F(x) */ 2599f425c49SPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 26019653cdaSPeter Brune 2619f425c49SPeter Brune /* differences for selection and restart */ 26213a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 2639566063dSJacob Faibussowitsch PetscCall(SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm)); 26413a62661SPeter Brune } else { 2659566063dSJacob Faibussowitsch PetscCall(SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm)); 26613a62661SPeter Brune } 267422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2681aa26658SKarl Rupp 2699f425c49SPeter Brune /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 2709566063dSJacob 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)); 27119653cdaSPeter Brune selectRestart = PETSC_FALSE; 27223b3e82cSAsbjørn Nilsen Riseth 27313a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 2749566063dSJacob Faibussowitsch PetscCall(SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart)); 27523b3e82cSAsbjørn Nilsen Riseth 27628ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 2771aa26658SKarl Rupp if (selectRestart) restart_count++; 2781aa26658SKarl Rupp else restart_count = 0; 27913a62661SPeter Brune } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 28013a62661SPeter Brune if (k_restart > ngmres->restart_periodic) { 28163a3b9bcSJacob Faibussowitsch if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %" PetscInt_FMT " iterations\n",k_restart)); 28213a62661SPeter Brune restart_count = ngmres->restart_it; 28313a62661SPeter Brune } 28413a62661SPeter Brune } 28523b3e82cSAsbjørn Nilsen Riseth 286b3c6a99cSPeter Brune ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 28723b3e82cSAsbjørn Nilsen Riseth 28828ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 28928ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) { 290dfbf837cSBarry Smith if (ngmres->monitor) { 29163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %" PetscInt_FMT "\n",k_restart)); 292dfbf837cSBarry Smith } 29328ed4a04SPeter Brune restart_count = 0; 29419653cdaSPeter Brune k_restart = 1; 29519653cdaSPeter Brune l = 1; 296b3c6a99cSPeter Brune ivec = 0; 29798b3e84cSPeter Brune /* q_{00} = nu */ 2989566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM)); 299d2e16ddcSPeter Brune } else { 30098b3e84cSPeter Brune /* select the current size of the subspace */ 3011e633543SBarry Smith if (l < ngmres->msize) l++; 30219653cdaSPeter Brune k_restart++; 30398b3e84cSPeter Brune /* place the current entry in the list of previous entries */ 30438774f0aSPeter Brune if (ngmres->candidate) { 30538774f0aSPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; 3069566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM)); 307d2e16ddcSPeter Brune } else { 30838774f0aSPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; 3099566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X)); 31019653cdaSPeter Brune } 311d2e16ddcSPeter Brune } 31219653cdaSPeter Brune 3139566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 314087dfb9eSxuemin snes->iter = k; 315f109b39eSPeter Brune snes->norm = fnorm; 3169566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3179566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,snes->iter)); 3189566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 3199566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,snes->iter,0,0,fnorm,&snes->reason,snes->cnvP)); 320087dfb9eSxuemin if (snes->reason) PetscFunctionReturn(0); 321a312c225SMatthew G Knepley } 322a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 323a312c225SMatthew G Knepley PetscFunctionReturn(0); 324a312c225SMatthew G Knepley } 325a312c225SMatthew G Knepley 32623b3e82cSAsbjørn Nilsen Riseth /*@ 32723b3e82cSAsbjørn Nilsen Riseth SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M 32823b3e82cSAsbjørn Nilsen Riseth 32923b3e82cSAsbjørn Nilsen Riseth Input Parameters: 33023b3e82cSAsbjørn Nilsen Riseth + snes - the SNES context. 33123b3e82cSAsbjørn Nilsen Riseth - flg - boolean value deciding whether to use the option or not 33223b3e82cSAsbjørn Nilsen Riseth 33323b3e82cSAsbjørn Nilsen Riseth Options Database: 33423b3e82cSAsbjørn Nilsen Riseth + -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M 33523b3e82cSAsbjørn Nilsen Riseth 33623b3e82cSAsbjørn Nilsen Riseth Level: intermediate 33723b3e82cSAsbjørn Nilsen Riseth 33823b3e82cSAsbjørn Nilsen Riseth Notes: 33923b3e82cSAsbjø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. 34023b3e82cSAsbjørn Nilsen Riseth To help the solver do that, reset the Krylov subspace whenever F_M increases. 34123b3e82cSAsbjørn Nilsen Riseth 34223b3e82cSAsbjørn Nilsen Riseth This option must be used with SNES_NGMRES_RESTART_DIFFERENCE 34323b3e82cSAsbjørn Nilsen Riseth 34423b3e82cSAsbjørn Nilsen Riseth The default is FALSE. 345db781477SPatrick Sanan .seealso: `SNES_NGMRES_RESTART_DIFFERENCE` 34623b3e82cSAsbjørn Nilsen Riseth @*/ 34723b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes,PetscBool flg) 34823b3e82cSAsbjørn Nilsen Riseth { 34923b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES,PetscBool); 35023b3e82cSAsbjørn Nilsen Riseth 35123b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 3529566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",&f)); 3539566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,flg)); 35423b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 35523b3e82cSAsbjørn Nilsen Riseth } 35623b3e82cSAsbjørn Nilsen Riseth 35723b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes,PetscBool flg) 35823b3e82cSAsbjørn Nilsen Riseth { 35923b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 36023b3e82cSAsbjørn Nilsen Riseth 36123b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 36223b3e82cSAsbjørn Nilsen Riseth ngmres->restart_fm_rise = flg; 36323b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 36423b3e82cSAsbjørn Nilsen Riseth } 36523b3e82cSAsbjørn Nilsen Riseth 36623b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes,PetscBool *flg) 36723b3e82cSAsbjørn Nilsen Riseth { 36823b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES,PetscBool*); 36923b3e82cSAsbjørn Nilsen Riseth 37023b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 3719566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",&f)); 3729566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,flg)); 37323b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 37423b3e82cSAsbjørn Nilsen Riseth } 37523b3e82cSAsbjørn Nilsen Riseth 37623b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes,PetscBool *flg) 37723b3e82cSAsbjørn Nilsen Riseth { 37823b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 37923b3e82cSAsbjørn Nilsen Riseth 38023b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 38123b3e82cSAsbjørn Nilsen Riseth *flg = ngmres->restart_fm_rise; 38223b3e82cSAsbjørn Nilsen Riseth PetscFunctionReturn(0); 38323b3e82cSAsbjørn Nilsen Riseth } 38423b3e82cSAsbjørn Nilsen Riseth 38513a62661SPeter Brune /*@ 38613a62661SPeter Brune SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES. 38713a62661SPeter Brune 38813a62661SPeter Brune Logically Collective on SNES 38913a62661SPeter Brune 39013a62661SPeter Brune Input Parameters: 39113a62661SPeter Brune + snes - the iterative context 39213a62661SPeter Brune - rtype - restart type 39313a62661SPeter Brune 39413a62661SPeter Brune Options Database: 39513a62661SPeter Brune + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type 3960c777b0cSPeter Brune - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic 39713a62661SPeter Brune 39813a62661SPeter Brune Level: intermediate 39913a62661SPeter Brune 40013a62661SPeter Brune SNESNGMRESRestartTypes: 40113a62661SPeter Brune + SNES_NGMRES_RESTART_NONE - never restart 40213a62661SPeter Brune . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria 40313a62661SPeter Brune - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations 40413a62661SPeter Brune 40513a62661SPeter Brune Notes: 40613a62661SPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 40713a62661SPeter Brune 40813a62661SPeter Brune @*/ 4090adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype) 4100adebc6cSBarry Smith { 41113a62661SPeter Brune PetscFunctionBegin; 41213a62661SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 413cac4c232SBarry Smith PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype)); 41413a62661SPeter Brune PetscFunctionReturn(0); 41513a62661SPeter Brune } 41613a62661SPeter Brune 41713a62661SPeter Brune /*@ 41813a62661SPeter Brune SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and 41913a62661SPeter Brune combined solution are used to create the next iterate. 42013a62661SPeter Brune 42113a62661SPeter Brune Logically Collective on SNES 42213a62661SPeter Brune 42313a62661SPeter Brune Input Parameters: 42413a62661SPeter Brune + snes - the iterative context 42513a62661SPeter Brune - stype - selection type 42613a62661SPeter Brune 42713a62661SPeter Brune Options Database: 42867b8a455SSatish Balay . -snes_ngmres_select_type<difference,none,linesearch> - select type 42913a62661SPeter Brune 43013a62661SPeter Brune Level: intermediate 43113a62661SPeter Brune 43213a62661SPeter Brune SNESNGMRESSelectTypes: 43313a62661SPeter Brune + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time 43413a62661SPeter Brune . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria 43513a62661SPeter Brune - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination 43613a62661SPeter Brune 43713a62661SPeter Brune Notes: 43813a62661SPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 43913a62661SPeter Brune 44013a62661SPeter Brune @*/ 4410adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype) 4420adebc6cSBarry Smith { 44313a62661SPeter Brune PetscFunctionBegin; 44413a62661SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 445cac4c232SBarry Smith PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype)); 44613a62661SPeter Brune PetscFunctionReturn(0); 44713a62661SPeter Brune } 44813a62661SPeter Brune 4490adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype) 4500adebc6cSBarry Smith { 45113a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 4525fd66863SKarl Rupp 45313a62661SPeter Brune PetscFunctionBegin; 45413a62661SPeter Brune ngmres->select_type = stype; 45513a62661SPeter Brune PetscFunctionReturn(0); 45613a62661SPeter Brune } 45713a62661SPeter Brune 4580adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype) 4590adebc6cSBarry Smith { 46013a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 4615fd66863SKarl Rupp 46213a62661SPeter Brune PetscFunctionBegin; 46313a62661SPeter Brune ngmres->restart_type = rtype; 46413a62661SPeter Brune PetscFunctionReturn(0); 46513a62661SPeter Brune } 46613a62661SPeter Brune 467dfbf837cSBarry Smith /*MC 4681867fe5bSPeter Brune SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 469a312c225SMatthew G Knepley 470dfbf837cSBarry Smith Level: beginner 471dfbf837cSBarry Smith 4721867fe5bSPeter Brune Options Database: 47313a62661SPeter Brune + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution 47438774f0aSPeter Brune . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions 47538774f0aSPeter Brune . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions 47613a62661SPeter Brune . -snes_ngmres_m - Number of stored previous solutions and residuals 47713a62661SPeter Brune . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart 47813a62661SPeter Brune . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination 47913a62661SPeter Brune . -snes_ngmres_gammaC - Residual tolerance for restart 48013a62661SPeter Brune . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart 48113a62661SPeter Brune . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart 48223b3e82cSAsbjørn Nilsen Riseth . -snes_ngmres_restart_fm_rise - Restart on residual rise from x_M step 48313a62661SPeter Brune . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration 4845c3e6ab7SPeter Brune . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother 48513a62661SPeter Brune - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type 4861867fe5bSPeter Brune 4871867fe5bSPeter Brune Notes: 4881867fe5bSPeter Brune 4891867fe5bSPeter Brune The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 4901867fe5bSPeter Brune optimization problem at each iteration. 4911867fe5bSPeter Brune 4924f02bc6aSBarry Smith Very similar to the SNESANDERSON algorithm. 4934f02bc6aSBarry Smith 4941867fe5bSPeter Brune References: 495606c0280SSatish Balay + * - C. W. Oosterlee and T. Washio, "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", 496dfbf837cSBarry Smith SIAM Journal on Scientific Computing, 21(5), 2000. 497606c0280SSatish Balay - * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 4984f02bc6aSBarry Smith SIAM Review, 57(4), 2015 4994f02bc6aSBarry Smith 500db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType` 501dfbf837cSBarry Smith M*/ 502a312c225SMatthew G Knepley 5038cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes) 504a312c225SMatthew G Knepley { 505a312c225SMatthew G Knepley SNES_NGMRES *ngmres; 506d8d34be6SBarry Smith SNESLineSearch linesearch; 507a312c225SMatthew G Knepley 508a312c225SMatthew G Knepley PetscFunctionBegin; 509a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES; 510a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES; 511a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 512a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES; 513a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES; 514a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES; 515a312c225SMatthew G Knepley 516efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 5172c155ee1SBarry Smith snes->usesksp = PETSC_FALSE; 518efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 5192c155ee1SBarry Smith 5204fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 5214fc747eaSLawrence Mitchell 5229566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes,&ngmres)); 523a312c225SMatthew G Knepley snes->data = (void*) ngmres; 524d2e16ddcSPeter Brune ngmres->msize = 30; 52519653cdaSPeter Brune 52688976e71SPeter Brune if (!snes->tolerancesset) { 5270e444f03SPeter Brune snes->max_funcs = 30000; 5280e444f03SPeter Brune snes->max_its = 10000; 52988976e71SPeter Brune } 5300e444f03SPeter Brune 53138774f0aSPeter Brune ngmres->candidate = PETSC_FALSE; 532d2e16ddcSPeter Brune 5339566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes,&linesearch)); 534ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 5359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC)); 536ec786807SJed Brown } 537d8d34be6SBarry Smith 5380298fd71SBarry Smith ngmres->additive_linesearch = NULL; 539077c4231SPeter Brune ngmres->approxfunc = PETSC_FALSE; 54028ed4a04SPeter Brune ngmres->restart_it = 2; 54113a62661SPeter Brune ngmres->restart_periodic = 30; 542f109b39eSPeter Brune ngmres->gammaA = 2.0; 543f109b39eSPeter Brune ngmres->gammaC = 2.0; 544cac108bcSPeter Brune ngmres->deltaB = 0.9; 545cac108bcSPeter Brune ngmres->epsilonB = 0.1; 54623b3e82cSAsbjørn Nilsen Riseth ngmres->restart_fm_rise = PETSC_FALSE; 547e7058c64SPeter Brune 54813a62661SPeter Brune ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE; 54913a62661SPeter Brune ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE; 55013a62661SPeter Brune 5519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES)); 5529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES)); 5539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",SNESNGMRESSetRestartFmRise_NGMRES)); 5549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",SNESNGMRESGetRestartFmRise_NGMRES)); 555a312c225SMatthew G Knepley PetscFunctionReturn(0); 556a312c225SMatthew G Knepley } 557