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 821687c63SPeter Brune PETSC_EXTERN const char *const SNESNGMRESRestartTypes[]; 9f31c9d25SPeter Brune 10f31c9d25SPeter Brune #undef __FUNCT__ 11f31c9d25SPeter Brune #define __FUNCT__ "SNESSetFromOptions_Anderson" 1240244768SBarry Smith static PetscErrorCode SNESSetFromOptions_Anderson(PetscOptionItems *PetscOptionsObject,SNES snes) 13f31c9d25SPeter Brune { 14f31c9d25SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 15f31c9d25SPeter Brune PetscErrorCode ierr; 1653efadd4SBarry Smith PetscBool monitor = PETSC_FALSE; 17f31c9d25SPeter Brune SNESLineSearch linesearch; 18f31c9d25SPeter Brune 19f31c9d25SPeter Brune PetscFunctionBegin; 20e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");CHKERRQ(ierr); 210298fd71SBarry Smith ierr = PetscOptionsInt("-snes_anderson_m", "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr); 225c3e6ab7SPeter Brune ierr = PetscOptionsReal("-snes_anderson_beta", "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);CHKERRQ(ierr); 230298fd71SBarry Smith ierr = PetscOptionsInt("-snes_anderson_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr); 240298fd71SBarry Smith ierr = PetscOptionsInt("-snes_anderson_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr); 2589eaf18bSBarry Smith ierr = PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,(PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr); 2653efadd4SBarry Smith ierr = PetscOptionsBool("-snes_anderson_monitor", "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&monitor,NULL);CHKERRQ(ierr); 2753efadd4SBarry Smith if (monitor) { 28ce94432eSBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 29f31c9d25SPeter Brune } 30f31c9d25SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 31f31c9d25SPeter Brune /* set the default type of the line search if the user hasn't already. */ 32f31c9d25SPeter Brune if (!snes->linesearch) { 337601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 34f31c9d25SPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 35f31c9d25SPeter Brune } 36f31c9d25SPeter Brune PetscFunctionReturn(0); 37f31c9d25SPeter Brune } 38f31c9d25SPeter Brune 39f31c9d25SPeter Brune #undef __FUNCT__ 40f31c9d25SPeter Brune #define __FUNCT__ "SNESSolve_Anderson" 4140244768SBarry Smith static PetscErrorCode SNESSolve_Anderson(SNES snes) 42f31c9d25SPeter Brune { 43f31c9d25SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 44f31c9d25SPeter Brune /* present solution, residual, and preconditioned residual */ 4521687c63SPeter Brune Vec X,F,B,D; 46f31c9d25SPeter Brune /* candidate linear combination answers */ 47ddd40ce5SPeter Brune Vec XA,FA,XM,FM; 48f31c9d25SPeter Brune 49f31c9d25SPeter Brune /* coefficients and RHS to the minimization problem */ 50b3c6a99cSPeter Brune PetscReal fnorm,fMnorm,fAnorm; 51b3c6a99cSPeter Brune PetscReal xnorm,ynorm; 5221687c63SPeter Brune PetscReal dnorm,dminnorm=0.0,fminnorm; 5321687c63SPeter Brune PetscInt restart_count=0; 54f31c9d25SPeter Brune PetscInt k,k_restart,l,ivec; 5521687c63SPeter Brune PetscBool selectRestart; 56f31c9d25SPeter Brune SNESConvergedReason reason; 57f31c9d25SPeter Brune PetscErrorCode ierr; 58f31c9d25SPeter Brune 59f31c9d25SPeter Brune PetscFunctionBegin; 606c4ed002SBarry 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); 61c579b300SPatrick Farrell 62fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 63f31c9d25SPeter Brune /* variable initialization */ 64f31c9d25SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 65f31c9d25SPeter Brune X = snes->vec_sol; 66f31c9d25SPeter Brune F = snes->vec_func; 67f31c9d25SPeter Brune B = snes->vec_rhs; 68f31c9d25SPeter Brune XA = snes->vec_sol_update; 69f31c9d25SPeter Brune FA = snes->work[0]; 7021687c63SPeter Brune D = snes->work[1]; 71f31c9d25SPeter Brune 72f31c9d25SPeter Brune /* work for the line search */ 73f31c9d25SPeter Brune XM = snes->work[3]; 74f31c9d25SPeter Brune FM = snes->work[4]; 75f31c9d25SPeter Brune 76e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 77f31c9d25SPeter Brune snes->iter = 0; 78f31c9d25SPeter Brune snes->norm = 0.; 79e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 80f31c9d25SPeter Brune 81f31c9d25SPeter Brune /* initialization */ 82f31c9d25SPeter Brune 83f31c9d25SPeter Brune /* r = F(x) */ 843a2ae377SPeter Brune 853a2ae377SPeter Brune if (snes->pc && snes->pcside == PC_LEFT) { 86be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 873a2ae377SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 883a2ae377SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 893a2ae377SPeter Brune snes->reason = SNES_DIVERGED_INNER; 903a2ae377SPeter Brune PetscFunctionReturn(0); 913a2ae377SPeter Brune } 923a2ae377SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 933a2ae377SPeter Brune } else { 94f31c9d25SPeter Brune if (!snes->vec_func_init_set) { 95f31c9d25SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 961aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 97c1c75074SPeter Brune 98f31c9d25SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 99422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 1003a2ae377SPeter Brune } 10121687c63SPeter Brune fminnorm = fnorm; 10221687c63SPeter Brune 103e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 104f31c9d25SPeter Brune snes->norm = fnorm; 105e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 106a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 107f31c9d25SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 108f31c9d25SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 109f31c9d25SPeter Brune if (snes->reason) PetscFunctionReturn(0); 110f31c9d25SPeter Brune 111f31c9d25SPeter Brune k_restart = 0; 112f31c9d25SPeter Brune l = 0; 113b3c6a99cSPeter Brune ivec = 0; 114f31c9d25SPeter Brune for (k=1; k < snes->max_its+1; k++) { 115f31c9d25SPeter Brune /* select which vector of the stored subspace will be updated */ 116f31c9d25SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 117f31c9d25SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 118f31c9d25SPeter Brune ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr); 119f31c9d25SPeter Brune 120f31c9d25SPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 121f31c9d25SPeter Brune ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr); 122f31c9d25SPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 123f31c9d25SPeter Brune 124f31c9d25SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 125f31c9d25SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 126f31c9d25SPeter Brune snes->reason = SNES_DIVERGED_INNER; 127f31c9d25SPeter Brune PetscFunctionReturn(0); 128f31c9d25SPeter Brune } 129be95d8f1SBarry Smith ierr = SNESGetNPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr); 130f31c9d25SPeter Brune if (ngmres->andersonBeta != 1.0) { 131f31c9d25SPeter Brune VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr); 132f31c9d25SPeter Brune } 133f31c9d25SPeter Brune } else { 134f31c9d25SPeter Brune ierr = VecCopy(F,FM);CHKERRQ(ierr); 135f31c9d25SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 136f31c9d25SPeter Brune ierr = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr); 137f31c9d25SPeter Brune fMnorm = fnorm; 138f31c9d25SPeter Brune } 139f31c9d25SPeter Brune 140b3c6a99cSPeter Brune ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr); 141b3c6a99cSPeter Brune ivec = k_restart % ngmres->msize; 14221687c63SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 143b3c6a99cSPeter Brune ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr); 14423b3e82cSAsbjørn Nilsen Riseth ierr = SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 14521687c63SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 1461aa26658SKarl Rupp if (selectRestart) restart_count++; 1471aa26658SKarl Rupp else restart_count = 0; 14821687c63SPeter Brune } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 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 if (k_restart > ngmres->restart_periodic) { 15121687c63SPeter Brune if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr); 15221687c63SPeter Brune restart_count = ngmres->restart_it; 15321687c63SPeter Brune } 154b3c6a99cSPeter Brune } else { 155b3c6a99cSPeter Brune ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr); 15621687c63SPeter Brune } 15721687c63SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 15821687c63SPeter Brune if (restart_count >= ngmres->restart_it) { 15921687c63SPeter Brune if (ngmres->monitor) { 16021687c63SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr); 16121687c63SPeter Brune } 16221687c63SPeter Brune restart_count = 0; 16321687c63SPeter Brune k_restart = 0; 16421687c63SPeter Brune l = 0; 165b3c6a99cSPeter Brune ivec = 0; 16621687c63SPeter Brune } else { 167f31c9d25SPeter Brune if (l < ngmres->msize) l++; 168f31c9d25SPeter Brune k_restart++; 169fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr); 17021687c63SPeter Brune } 171f31c9d25SPeter Brune 172b3c6a99cSPeter Brune fnorm = fAnorm; 17321687c63SPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; 17421687c63SPeter Brune 175f31c9d25SPeter Brune ierr = VecCopy(XA,X);CHKERRQ(ierr); 176f31c9d25SPeter Brune ierr = VecCopy(FA,F);CHKERRQ(ierr); 177f31c9d25SPeter Brune 178e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 179f31c9d25SPeter Brune snes->iter = k; 180f31c9d25SPeter Brune snes->norm = fnorm; 181e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 182a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 183f31c9d25SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 184b3c6a99cSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 185f31c9d25SPeter Brune if (snes->reason) PetscFunctionReturn(0); 186f31c9d25SPeter Brune } 187f31c9d25SPeter Brune snes->reason = SNES_DIVERGED_MAX_IT; 188f31c9d25SPeter Brune PetscFunctionReturn(0); 189f31c9d25SPeter Brune } 190f31c9d25SPeter Brune 191f31c9d25SPeter Brune /*MC 1924f02bc6aSBarry Smith SNESANDERSON - Anderson Mixing method. 193f31c9d25SPeter Brune 194f31c9d25SPeter Brune Level: beginner 195f31c9d25SPeter Brune 196f31c9d25SPeter Brune Options Database: 197f31c9d25SPeter Brune + -snes_anderson_m - Number of stored previous solutions and residuals 198f31c9d25SPeter Brune . -snes_anderson_beta - Relaxation parameter; X_{update} = X + \beta F 1995c3e6ab7SPeter Brune . -snes_anderson_restart_type - Type of restart (see SNESNGMRES) 2005c3e6ab7SPeter Brune . -snes_anderson_restart_it - Number of iterations of restart conditions before restart 2015c3e6ab7SPeter Brune . -snes_anderson_restart - Number of iterations before periodic restart 202f31c9d25SPeter Brune - -snes_anderson_monitor - Prints relevant information about the ngmres iteration 203f31c9d25SPeter Brune 204f31c9d25SPeter Brune Notes: 205f31c9d25SPeter Brune 206f31c9d25SPeter Brune The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized 207f31c9d25SPeter Brune optimization problem at each iteration. 208f31c9d25SPeter Brune 2094f02bc6aSBarry Smith Very similar to the SNESNGMRES algorithm. 2104f02bc6aSBarry Smith 211f31c9d25SPeter Brune References: 21296a0c994SBarry Smith + 1. - D. G. Anderson. Iterative procedures for nonlinear integral equations. 21396a0c994SBarry Smith J. Assoc. Comput. Mach., 12, 1965." 21496a0c994SBarry Smith - 2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers", 2154f02bc6aSBarry Smith SIAM Review, 57(4), 2015 2164f02bc6aSBarry Smith 2174f02bc6aSBarry Smith 2185c3e6ab7SPeter Brune .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 219f31c9d25SPeter Brune M*/ 220f31c9d25SPeter Brune 221f31c9d25SPeter Brune #undef __FUNCT__ 222f31c9d25SPeter Brune #define __FUNCT__ "SNESCreate_Anderson" 2238cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes) 224f31c9d25SPeter Brune { 225f31c9d25SPeter Brune SNES_NGMRES *ngmres; 226f31c9d25SPeter Brune PetscErrorCode ierr; 227f31c9d25SPeter Brune 228f31c9d25SPeter Brune PetscFunctionBegin; 229f31c9d25SPeter Brune snes->ops->destroy = SNESDestroy_NGMRES; 230f31c9d25SPeter Brune snes->ops->setup = SNESSetUp_NGMRES; 231f31c9d25SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Anderson; 232f31c9d25SPeter Brune snes->ops->view = SNESView_NGMRES; 233f31c9d25SPeter Brune snes->ops->solve = SNESSolve_Anderson; 234f31c9d25SPeter Brune snes->ops->reset = SNESReset_NGMRES; 235f31c9d25SPeter Brune 236f31c9d25SPeter Brune snes->usespc = PETSC_TRUE; 237f31c9d25SPeter Brune snes->usesksp = PETSC_FALSE; 2386c67d002SPeter Brune snes->pcside = PC_RIGHT; 239f31c9d25SPeter Brune 240*4fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 241*4fc747eaSLawrence Mitchell 242b00a9115SJed Brown ierr = PetscNewLog(snes,&ngmres);CHKERRQ(ierr); 243f31c9d25SPeter Brune snes->data = (void*) ngmres; 244f31c9d25SPeter Brune ngmres->msize = 30; 245f31c9d25SPeter Brune 246f31c9d25SPeter Brune if (!snes->tolerancesset) { 247f31c9d25SPeter Brune snes->max_funcs = 30000; 248f31c9d25SPeter Brune snes->max_its = 10000; 249f31c9d25SPeter Brune } 250f31c9d25SPeter Brune 2510298fd71SBarry Smith ngmres->additive_linesearch = NULL; 252077c4231SPeter Brune ngmres->approxfunc = PETSC_FALSE; 25321687c63SPeter Brune ngmres->restart_type = SNES_NGMRES_RESTART_NONE; 254f31c9d25SPeter Brune ngmres->restart_it = 2; 255f31c9d25SPeter Brune ngmres->restart_periodic = 30; 256f31c9d25SPeter Brune ngmres->gammaA = 2.0; 257f31c9d25SPeter Brune ngmres->gammaC = 2.0; 258f31c9d25SPeter Brune ngmres->deltaB = 0.9; 259f31c9d25SPeter Brune ngmres->epsilonB = 0.1; 260f31c9d25SPeter Brune 261f31c9d25SPeter Brune ngmres->andersonBeta = 1.0; 262f31c9d25SPeter Brune PetscFunctionReturn(0); 263f31c9d25SPeter Brune } 264