1f31c9d25SPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 2f31c9d25SPeter Brune 39df56c09SSatish Balay extern PetscErrorCode SNESDestroy_NGMRES(SNES); 49df56c09SSatish Balay extern PetscErrorCode SNESReset_NGMRES(SNES); 59df56c09SSatish Balay extern PetscErrorCode SNESSetUp_NGMRES(SNES); 69df56c09SSatish Balay extern PetscErrorCode SNESView_NGMRES(SNES,PetscViewer); 7f31c9d25SPeter Brune 840244768SBarry Smith static PetscErrorCode SNESSetFromOptions_Anderson(PetscOptionItems *PetscOptionsObject,SNES snes) 9f31c9d25SPeter Brune { 10f31c9d25SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 11f31c9d25SPeter Brune PetscErrorCode ierr; 1253efadd4SBarry Smith PetscBool monitor = PETSC_FALSE; 13f31c9d25SPeter Brune SNESLineSearch linesearch; 14f31c9d25SPeter Brune 15f31c9d25SPeter Brune PetscFunctionBegin; 16e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");CHKERRQ(ierr); 170298fd71SBarry Smith ierr = PetscOptionsInt("-snes_anderson_m", "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr); 185c3e6ab7SPeter Brune ierr = PetscOptionsReal("-snes_anderson_beta", "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);CHKERRQ(ierr); 190298fd71SBarry Smith ierr = PetscOptionsInt("-snes_anderson_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr); 200298fd71SBarry Smith ierr = PetscOptionsInt("-snes_anderson_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr); 2189eaf18bSBarry Smith ierr = PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,(PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr); 2253efadd4SBarry Smith ierr = PetscOptionsBool("-snes_anderson_monitor", "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&monitor,NULL);CHKERRQ(ierr); 2353efadd4SBarry Smith if (monitor) { 24ce94432eSBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 25f31c9d25SPeter Brune } 26f31c9d25SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 27f31c9d25SPeter Brune /* set the default type of the line search if the user hasn't already. */ 28f31c9d25SPeter Brune if (!snes->linesearch) { 297601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 30f31c9d25SPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 31f31c9d25SPeter Brune } 32f31c9d25SPeter Brune PetscFunctionReturn(0); 33f31c9d25SPeter Brune } 34f31c9d25SPeter Brune 3540244768SBarry Smith static PetscErrorCode SNESSolve_Anderson(SNES snes) 36f31c9d25SPeter Brune { 37f31c9d25SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 38f31c9d25SPeter Brune /* present solution, residual, and preconditioned residual */ 3921687c63SPeter Brune Vec X,F,B,D; 40f31c9d25SPeter Brune /* candidate linear combination answers */ 41ddd40ce5SPeter Brune Vec XA,FA,XM,FM; 42f31c9d25SPeter Brune 43f31c9d25SPeter Brune /* coefficients and RHS to the minimization problem */ 44b3c6a99cSPeter Brune PetscReal fnorm,fMnorm,fAnorm; 45b3c6a99cSPeter Brune PetscReal xnorm,ynorm; 4621687c63SPeter Brune PetscReal dnorm,dminnorm=0.0,fminnorm; 4721687c63SPeter Brune PetscInt restart_count=0; 48f31c9d25SPeter Brune PetscInt k,k_restart,l,ivec; 4921687c63SPeter Brune PetscBool selectRestart; 50f31c9d25SPeter Brune SNESConvergedReason reason; 51f31c9d25SPeter Brune PetscErrorCode ierr; 52f31c9d25SPeter Brune 53f31c9d25SPeter Brune PetscFunctionBegin; 546c4ed002SBarry 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); 55c579b300SPatrick Farrell 56fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 57f31c9d25SPeter Brune /* variable initialization */ 58f31c9d25SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 59f31c9d25SPeter Brune X = snes->vec_sol; 60f31c9d25SPeter Brune F = snes->vec_func; 61f31c9d25SPeter Brune B = snes->vec_rhs; 62f31c9d25SPeter Brune XA = snes->vec_sol_update; 63f31c9d25SPeter Brune FA = snes->work[0]; 6421687c63SPeter Brune D = snes->work[1]; 65f31c9d25SPeter Brune 66f31c9d25SPeter Brune /* work for the line search */ 67f31c9d25SPeter Brune XM = snes->work[3]; 68f31c9d25SPeter Brune FM = snes->work[4]; 69f31c9d25SPeter Brune 70e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 71f31c9d25SPeter Brune snes->iter = 0; 72f31c9d25SPeter Brune snes->norm = 0.; 73e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 74f31c9d25SPeter Brune 75f31c9d25SPeter Brune /* initialization */ 76f31c9d25SPeter Brune 77f31c9d25SPeter Brune /* r = F(x) */ 783a2ae377SPeter Brune 79*efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT) { 80be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 81*efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 823a2ae377SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 833a2ae377SPeter Brune snes->reason = SNES_DIVERGED_INNER; 843a2ae377SPeter Brune PetscFunctionReturn(0); 853a2ae377SPeter Brune } 863a2ae377SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 873a2ae377SPeter Brune } else { 88f31c9d25SPeter Brune if (!snes->vec_func_init_set) { 89f31c9d25SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 901aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 91c1c75074SPeter Brune 92f31c9d25SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 93422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 943a2ae377SPeter Brune } 9521687c63SPeter Brune fminnorm = fnorm; 9621687c63SPeter Brune 97e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 98f31c9d25SPeter Brune snes->norm = fnorm; 99e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 100a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 101f31c9d25SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 102f31c9d25SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 103f31c9d25SPeter Brune if (snes->reason) PetscFunctionReturn(0); 104f31c9d25SPeter Brune 105f31c9d25SPeter Brune k_restart = 0; 106f31c9d25SPeter Brune l = 0; 107b3c6a99cSPeter Brune ivec = 0; 108f31c9d25SPeter Brune for (k=1; k < snes->max_its+1; k++) { 109f31c9d25SPeter Brune /* select which vector of the stored subspace will be updated */ 110*efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_RIGHT) { 111f31c9d25SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 112*efd4aadfSBarry Smith ierr = SNESSetInitialFunction(snes->npc,F);CHKERRQ(ierr); 113f31c9d25SPeter Brune 114*efd4aadfSBarry Smith ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0);CHKERRQ(ierr); 115*efd4aadfSBarry Smith ierr = SNESSolve(snes->npc,B,XM);CHKERRQ(ierr); 116*efd4aadfSBarry Smith ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0);CHKERRQ(ierr); 117f31c9d25SPeter Brune 118*efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 119f31c9d25SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 120f31c9d25SPeter Brune snes->reason = SNES_DIVERGED_INNER; 121f31c9d25SPeter Brune PetscFunctionReturn(0); 122f31c9d25SPeter Brune } 123be95d8f1SBarry Smith ierr = SNESGetNPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr); 124f31c9d25SPeter Brune if (ngmres->andersonBeta != 1.0) { 12535f895b4SBarry Smith ierr = VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr); 126f31c9d25SPeter Brune } 127f31c9d25SPeter Brune } else { 128f31c9d25SPeter Brune ierr = VecCopy(F,FM);CHKERRQ(ierr); 129f31c9d25SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 130f31c9d25SPeter Brune ierr = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr); 131f31c9d25SPeter Brune fMnorm = fnorm; 132f31c9d25SPeter Brune } 133f31c9d25SPeter Brune 134b3c6a99cSPeter Brune ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr); 135b3c6a99cSPeter Brune ivec = k_restart % ngmres->msize; 13621687c63SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 137b3c6a99cSPeter Brune ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr); 13823b3e82cSAsbjørn Nilsen Riseth ierr = SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 13921687c63SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 1401aa26658SKarl Rupp if (selectRestart) restart_count++; 1411aa26658SKarl Rupp else restart_count = 0; 14221687c63SPeter Brune } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 143b3c6a99cSPeter Brune ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr); 14421687c63SPeter Brune if (k_restart > ngmres->restart_periodic) { 14521687c63SPeter Brune if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr); 14621687c63SPeter Brune restart_count = ngmres->restart_it; 14721687c63SPeter Brune } 148b3c6a99cSPeter Brune } else { 149b3c6a99cSPeter Brune ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr); 15021687c63SPeter Brune } 15121687c63SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 15221687c63SPeter Brune if (restart_count >= ngmres->restart_it) { 15321687c63SPeter Brune if (ngmres->monitor) { 15421687c63SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr); 15521687c63SPeter Brune } 15621687c63SPeter Brune restart_count = 0; 15721687c63SPeter Brune k_restart = 0; 15821687c63SPeter Brune l = 0; 159b3c6a99cSPeter Brune ivec = 0; 16021687c63SPeter Brune } else { 161f31c9d25SPeter Brune if (l < ngmres->msize) l++; 162f31c9d25SPeter Brune k_restart++; 163fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr); 16421687c63SPeter Brune } 165f31c9d25SPeter Brune 166b3c6a99cSPeter Brune fnorm = fAnorm; 16721687c63SPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; 16821687c63SPeter Brune 169f31c9d25SPeter Brune ierr = VecCopy(XA,X);CHKERRQ(ierr); 170f31c9d25SPeter Brune ierr = VecCopy(FA,F);CHKERRQ(ierr); 171f31c9d25SPeter Brune 172e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 173f31c9d25SPeter Brune snes->iter = k; 174f31c9d25SPeter Brune snes->norm = fnorm; 175e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 176a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 177f31c9d25SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 178b3c6a99cSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 179f31c9d25SPeter Brune if (snes->reason) PetscFunctionReturn(0); 180f31c9d25SPeter Brune } 181f31c9d25SPeter Brune snes->reason = SNES_DIVERGED_MAX_IT; 182f31c9d25SPeter Brune PetscFunctionReturn(0); 183f31c9d25SPeter Brune } 184f31c9d25SPeter Brune 185f31c9d25SPeter Brune /*MC 1864f02bc6aSBarry Smith SNESANDERSON - Anderson Mixing method. 187f31c9d25SPeter Brune 188f31c9d25SPeter Brune Level: beginner 189f31c9d25SPeter Brune 190f31c9d25SPeter Brune Options Database: 191f31c9d25SPeter Brune + -snes_anderson_m - Number of stored previous solutions and residuals 19235f895b4SBarry Smith . -snes_anderson_beta - Anderson mixing parameter 1935c3e6ab7SPeter Brune . -snes_anderson_restart_type - Type of restart (see SNESNGMRES) 1945c3e6ab7SPeter Brune . -snes_anderson_restart_it - Number of iterations of restart conditions before restart 1955c3e6ab7SPeter Brune . -snes_anderson_restart - Number of iterations before periodic restart 196f31c9d25SPeter Brune - -snes_anderson_monitor - Prints relevant information about the ngmres iteration 197f31c9d25SPeter Brune 198f31c9d25SPeter Brune Notes: 199f31c9d25SPeter Brune 200f31c9d25SPeter Brune The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized 201f31c9d25SPeter Brune optimization problem at each iteration. 202f31c9d25SPeter Brune 2034f02bc6aSBarry Smith Very similar to the SNESNGMRES algorithm. 2044f02bc6aSBarry Smith 205f31c9d25SPeter Brune References: 20696a0c994SBarry Smith + 1. - D. G. Anderson. Iterative procedures for nonlinear integral equations. 20796a0c994SBarry Smith J. Assoc. Comput. Mach., 12, 1965." 20896a0c994SBarry Smith - 2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers", 2094f02bc6aSBarry Smith SIAM Review, 57(4), 2015 2104f02bc6aSBarry Smith 2114f02bc6aSBarry Smith 2125c3e6ab7SPeter Brune .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 213f31c9d25SPeter Brune M*/ 214f31c9d25SPeter Brune 2158cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes) 216f31c9d25SPeter Brune { 217f31c9d25SPeter Brune SNES_NGMRES *ngmres; 218f31c9d25SPeter Brune PetscErrorCode ierr; 219f31c9d25SPeter Brune 220f31c9d25SPeter Brune PetscFunctionBegin; 221f31c9d25SPeter Brune snes->ops->destroy = SNESDestroy_NGMRES; 222f31c9d25SPeter Brune snes->ops->setup = SNESSetUp_NGMRES; 223f31c9d25SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Anderson; 224f31c9d25SPeter Brune snes->ops->view = SNESView_NGMRES; 225f31c9d25SPeter Brune snes->ops->solve = SNESSolve_Anderson; 226f31c9d25SPeter Brune snes->ops->reset = SNESReset_NGMRES; 227f31c9d25SPeter Brune 228*efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 229f31c9d25SPeter Brune snes->usesksp = PETSC_FALSE; 230*efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 231f31c9d25SPeter Brune 2324fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 2334fc747eaSLawrence Mitchell 234b00a9115SJed Brown ierr = PetscNewLog(snes,&ngmres);CHKERRQ(ierr); 235f31c9d25SPeter Brune snes->data = (void*) ngmres; 236f31c9d25SPeter Brune ngmres->msize = 30; 237f31c9d25SPeter Brune 238f31c9d25SPeter Brune if (!snes->tolerancesset) { 239f31c9d25SPeter Brune snes->max_funcs = 30000; 240f31c9d25SPeter Brune snes->max_its = 10000; 241f31c9d25SPeter Brune } 242f31c9d25SPeter Brune 2430298fd71SBarry Smith ngmres->additive_linesearch = NULL; 244077c4231SPeter Brune ngmres->approxfunc = PETSC_FALSE; 24521687c63SPeter Brune ngmres->restart_type = SNES_NGMRES_RESTART_NONE; 246f31c9d25SPeter Brune ngmres->restart_it = 2; 247f31c9d25SPeter Brune ngmres->restart_periodic = 30; 248f31c9d25SPeter Brune ngmres->gammaA = 2.0; 249f31c9d25SPeter Brune ngmres->gammaC = 2.0; 250f31c9d25SPeter Brune ngmres->deltaB = 0.9; 251f31c9d25SPeter Brune ngmres->epsilonB = 0.1; 252f31c9d25SPeter Brune 253f31c9d25SPeter Brune ngmres->andersonBeta = 1.0; 254f31c9d25SPeter Brune PetscFunctionReturn(0); 255f31c9d25SPeter Brune } 256