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" 12f31c9d25SPeter Brune PetscErrorCode SNESSetFromOptions_Anderson(SNES snes) 13f31c9d25SPeter Brune { 14f31c9d25SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 15f31c9d25SPeter Brune PetscErrorCode ierr; 16f31c9d25SPeter Brune PetscBool debug; 17f31c9d25SPeter Brune SNESLineSearch linesearch; 18f31c9d25SPeter Brune 19f31c9d25SPeter Brune PetscFunctionBegin; 20f31c9d25SPeter Brune ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr); 210298fd71SBarry Smith ierr = PetscOptionsInt("-snes_anderson_m", "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr); 220298fd71SBarry Smith ierr = PetscOptionsReal("-snes_anderson_beta", "Number of directions","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);CHKERRQ(ierr); 230298fd71SBarry Smith ierr = PetscOptionsBool("-snes_anderson_monitor", "Monitor actions of NGMRES","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);CHKERRQ(ierr); 240298fd71SBarry Smith ierr = PetscOptionsInt("-snes_anderson_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr); 250298fd71SBarry Smith ierr = PetscOptionsInt("-snes_anderson_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr); 26077c4231SPeter Brune ierr = PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes, 27077c4231SPeter Brune (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr); 28f31c9d25SPeter Brune if (debug) { 29ce94432eSBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 30f31c9d25SPeter Brune } 31f31c9d25SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 32f31c9d25SPeter Brune /* set the default type of the line search if the user hasn't already. */ 33f31c9d25SPeter Brune if (!snes->linesearch) { 34f31c9d25SPeter Brune ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 35f31c9d25SPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 36f31c9d25SPeter Brune } 37f31c9d25SPeter Brune PetscFunctionReturn(0); 38f31c9d25SPeter Brune } 39f31c9d25SPeter Brune 40f31c9d25SPeter Brune #undef __FUNCT__ 41f31c9d25SPeter Brune #define __FUNCT__ "SNESSolve_Anderson" 42f31c9d25SPeter Brune PetscErrorCode SNESSolve_Anderson(SNES snes) 43f31c9d25SPeter Brune { 44f31c9d25SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 45f31c9d25SPeter Brune /* present solution, residual, and preconditioned residual */ 4621687c63SPeter Brune Vec X,F,B,D; 47f31c9d25SPeter Brune 48f31c9d25SPeter Brune /* candidate linear combination answers */ 49f31c9d25SPeter Brune Vec XA,FA,XM,FM,FPC; 50f31c9d25SPeter Brune 51f31c9d25SPeter Brune /* coefficients and RHS to the minimization problem */ 52f31c9d25SPeter Brune PetscReal fnorm,fMnorm; 5321687c63SPeter Brune PetscReal dnorm,dminnorm=0.0,fminnorm; 5421687c63SPeter Brune PetscInt restart_count=0; 55f31c9d25SPeter Brune PetscInt k,k_restart,l,ivec; 56f31c9d25SPeter Brune 5721687c63SPeter Brune PetscBool selectRestart; 5821687c63SPeter Brune 59f31c9d25SPeter Brune SNESConvergedReason reason; 60f31c9d25SPeter Brune PetscErrorCode ierr; 61f31c9d25SPeter Brune 62f31c9d25SPeter Brune PetscFunctionBegin; 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 76f31c9d25SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 77f31c9d25SPeter Brune snes->iter = 0; 78f31c9d25SPeter Brune snes->norm = 0.; 79f31c9d25SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 80f31c9d25SPeter Brune 81f31c9d25SPeter Brune /* initialization */ 82f31c9d25SPeter Brune 83f31c9d25SPeter Brune /* r = F(x) */ 84f31c9d25SPeter Brune if (!snes->vec_func_init_set) { 85f31c9d25SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 86f31c9d25SPeter Brune if (snes->domainerror) { 87f31c9d25SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 88f31c9d25SPeter Brune PetscFunctionReturn(0); 89f31c9d25SPeter Brune } 901aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 91f31c9d25SPeter Brune 92f31c9d25SPeter Brune if (!snes->norm_init_set) { 93f31c9d25SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 94189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 95189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 96189a9710SBarry Smith PetscFunctionReturn(0); 97189a9710SBarry Smith } 98f31c9d25SPeter Brune } else { 99f31c9d25SPeter Brune fnorm = snes->norm_init; 100f31c9d25SPeter Brune snes->norm_init_set = PETSC_FALSE; 101f31c9d25SPeter Brune } 10221687c63SPeter Brune fminnorm = fnorm; 10321687c63SPeter Brune 104f31c9d25SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 105f31c9d25SPeter Brune snes->norm = fnorm; 106f31c9d25SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 107a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 108f31c9d25SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 109f31c9d25SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 110f31c9d25SPeter Brune if (snes->reason) PetscFunctionReturn(0); 111f31c9d25SPeter Brune 112f31c9d25SPeter Brune k_restart = 0; 113f31c9d25SPeter Brune l = 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 ivec = k_restart % ngmres->msize; 117f31c9d25SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 118f31c9d25SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 119f31c9d25SPeter Brune ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr); 120f31c9d25SPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr); 121f31c9d25SPeter Brune 122f31c9d25SPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 123f31c9d25SPeter Brune ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr); 124f31c9d25SPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 125f31c9d25SPeter Brune 126f31c9d25SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 127f31c9d25SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 128f31c9d25SPeter Brune snes->reason = SNES_DIVERGED_INNER; 129f31c9d25SPeter Brune PetscFunctionReturn(0); 130f31c9d25SPeter Brune } 1310298fd71SBarry Smith ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr); 132f31c9d25SPeter Brune ierr = VecCopy(FPC,FM);CHKERRQ(ierr); 133f31c9d25SPeter Brune ierr = SNESGetFunctionNorm(snes->pc,&fMnorm);CHKERRQ(ierr); 134f31c9d25SPeter Brune if (ngmres->andersonBeta != 1.0) { 135f31c9d25SPeter Brune VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr); 136f31c9d25SPeter Brune } 137f31c9d25SPeter Brune } else { 138f31c9d25SPeter Brune ierr = VecCopy(F,FM);CHKERRQ(ierr); 139f31c9d25SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 140f31c9d25SPeter Brune ierr = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr); 141f31c9d25SPeter Brune fMnorm = fnorm; 142f31c9d25SPeter Brune } 143f31c9d25SPeter Brune 144fa8c639aSPeter Brune ierr = SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr); 14521687c63SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 1460298fd71SBarry Smith ierr = SNESNGMRESCalculateDifferences_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL);CHKERRQ(ierr); 14721687c63SPeter Brune ierr = SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 14821687c63SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 1491aa26658SKarl Rupp if (selectRestart) restart_count++; 1501aa26658SKarl Rupp else restart_count = 0; 15121687c63SPeter Brune } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 15221687c63SPeter Brune if (k_restart > ngmres->restart_periodic) { 15321687c63SPeter Brune if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr); 15421687c63SPeter Brune restart_count = ngmres->restart_it; 15521687c63SPeter Brune } 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; 16521687c63SPeter Brune } else { 166f31c9d25SPeter Brune if (l < ngmres->msize) l++; 167f31c9d25SPeter Brune k_restart++; 168fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr); 16921687c63SPeter Brune } 170f31c9d25SPeter Brune 171f31c9d25SPeter Brune ierr = VecNorm(FA,NORM_2,&fnorm);CHKERRQ(ierr); 17221687c63SPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; 17321687c63SPeter Brune 174f31c9d25SPeter Brune ierr = VecCopy(XA,X);CHKERRQ(ierr); 175f31c9d25SPeter Brune ierr = VecCopy(FA,F);CHKERRQ(ierr); 176f31c9d25SPeter Brune 177f31c9d25SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 178f31c9d25SPeter Brune snes->iter = k; 179f31c9d25SPeter Brune snes->norm = fnorm; 180f31c9d25SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 181a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 182f31c9d25SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 183f31c9d25SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 184f31c9d25SPeter Brune if (snes->reason) PetscFunctionReturn(0); 185f31c9d25SPeter Brune } 186f31c9d25SPeter Brune snes->reason = SNES_DIVERGED_MAX_IT; 187f31c9d25SPeter Brune PetscFunctionReturn(0); 188f31c9d25SPeter Brune } 189f31c9d25SPeter Brune 190f31c9d25SPeter Brune /*MC 191f31c9d25SPeter Brune SNESAnderson - Anderson Mixing method. 192f31c9d25SPeter Brune 193f31c9d25SPeter Brune Level: beginner 194f31c9d25SPeter Brune 195f31c9d25SPeter Brune Options Database: 196f31c9d25SPeter Brune + -snes_anderson_m - Number of stored previous solutions and residuals 197f31c9d25SPeter Brune . -snes_anderson_beta - Relaxation parameter; X_{update} = X + \beta F 198f31c9d25SPeter Brune - -snes_anderson_monitor - Prints relevant information about the ngmres iteration 199f31c9d25SPeter Brune 200f31c9d25SPeter Brune Notes: 201f31c9d25SPeter Brune 202f31c9d25SPeter Brune The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized 203f31c9d25SPeter Brune optimization problem at each iteration. 204f31c9d25SPeter Brune 205f31c9d25SPeter Brune References: 206f31c9d25SPeter Brune 207f31c9d25SPeter Brune "D. G. Anderson. Iterative procedures for nonlinear integral equations. 208f31c9d25SPeter Brune J. Assoc. Comput. Mach., 12:547–560, 1965." 209f31c9d25SPeter Brune 210f31c9d25SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 211f31c9d25SPeter Brune M*/ 212f31c9d25SPeter Brune 213f31c9d25SPeter Brune #undef __FUNCT__ 214f31c9d25SPeter Brune #define __FUNCT__ "SNESCreate_Anderson" 215*8cc058d9SJed 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 228f31c9d25SPeter Brune snes->usespc = PETSC_TRUE; 229f31c9d25SPeter Brune snes->usesksp = PETSC_FALSE; 230f31c9d25SPeter Brune 231f31c9d25SPeter Brune ierr = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr); 232f31c9d25SPeter Brune snes->data = (void*) ngmres; 233f31c9d25SPeter Brune ngmres->msize = 30; 234f31c9d25SPeter Brune 235f31c9d25SPeter Brune if (!snes->tolerancesset) { 236f31c9d25SPeter Brune snes->max_funcs = 30000; 237f31c9d25SPeter Brune snes->max_its = 10000; 238f31c9d25SPeter Brune } 239f31c9d25SPeter Brune 2400298fd71SBarry Smith ngmres->additive_linesearch = NULL; 241077c4231SPeter Brune ngmres->approxfunc = PETSC_FALSE; 24221687c63SPeter Brune ngmres->restart_type = SNES_NGMRES_RESTART_NONE; 243f31c9d25SPeter Brune ngmres->restart_it = 2; 244f31c9d25SPeter Brune ngmres->restart_periodic = 30; 245f31c9d25SPeter Brune ngmres->gammaA = 2.0; 246f31c9d25SPeter Brune ngmres->gammaC = 2.0; 247f31c9d25SPeter Brune ngmres->deltaB = 0.9; 248f31c9d25SPeter Brune ngmres->epsilonB = 0.1; 249f31c9d25SPeter Brune 250f31c9d25SPeter Brune ngmres->andersonBeta = 1.0; 251f31c9d25SPeter Brune PetscFunctionReturn(0); 252f31c9d25SPeter Brune } 25399e0435eSBarry Smith 254