1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 2 3 PETSC_EXTERN PetscErrorCode SNESDestroy_NGMRES(SNES); 4 PETSC_EXTERN PetscErrorCode SNESReset_NGMRES(SNES); 5 PETSC_EXTERN PetscErrorCode SNESSetUp_NGMRES(SNES); 6 PETSC_EXTERN PetscErrorCode SNESView_NGMRES(SNES,PetscViewer); 7 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "SNESSetFromOptions_Anderson" 11 PetscErrorCode SNESSetFromOptions_Anderson(SNES snes) 12 { 13 SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 14 PetscErrorCode ierr; 15 PetscBool debug; 16 SNESLineSearch linesearch; 17 18 PetscFunctionBegin; 19 ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr); 20 ierr = PetscOptionsInt("-snes_anderson_m","Number of directions","SNES",ngmres->msize,&ngmres->msize,PETSC_NULL);CHKERRQ(ierr); 21 ierr = PetscOptionsReal("-snes_anderson_beta","Number of directions","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,PETSC_NULL);CHKERRQ(ierr); 22 ierr = PetscOptionsBool("-snes_anderson_monitor","Monitor actions of NGMRES","SNES",ngmres->monitor ? PETSC_TRUE: PETSC_FALSE,&debug,PETSC_NULL);CHKERRQ(ierr); 23 if (debug) { 24 ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 25 } 26 ierr = PetscOptionsTail();CHKERRQ(ierr); 27 /* set the default type of the line search if the user hasn't already. */ 28 if (!snes->linesearch) { 29 ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 30 ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 31 } 32 PetscFunctionReturn(0); 33 } 34 35 #undef __FUNCT__ 36 #define __FUNCT__ "SNESSolve_Anderson" 37 PetscErrorCode SNESSolve_Anderson(SNES snes) 38 { 39 SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 40 /* present solution, residual, and preconditioned residual */ 41 Vec X,F,B; 42 43 /* candidate linear combination answers */ 44 Vec XA,FA,XM,FM,FPC; 45 46 /* coefficients and RHS to the minimization problem */ 47 PetscReal fnorm,fMnorm; 48 PetscInt k,k_restart,l,ivec; 49 50 SNESConvergedReason reason; 51 PetscErrorCode ierr; 52 53 PetscFunctionBegin; 54 /* variable initialization */ 55 snes->reason = SNES_CONVERGED_ITERATING; 56 X = snes->vec_sol; 57 F = snes->vec_func; 58 B = snes->vec_rhs; 59 XA = snes->vec_sol_update; 60 FA = snes->work[0]; 61 62 /* work for the line search */ 63 XM = snes->work[3]; 64 FM = snes->work[4]; 65 66 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 67 snes->iter = 0; 68 snes->norm = 0.; 69 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 70 71 /* initialization */ 72 73 /* r = F(x) */ 74 if (!snes->vec_func_init_set) { 75 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 76 if (snes->domainerror) { 77 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 78 PetscFunctionReturn(0); 79 } 80 } else { 81 snes->vec_func_init_set = PETSC_FALSE; 82 } 83 84 if (!snes->norm_init_set) { 85 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 86 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in function evaluation"); 87 } else { 88 fnorm = snes->norm_init; 89 snes->norm_init_set = PETSC_FALSE; 90 } 91 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 92 snes->norm = fnorm; 93 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 94 SNESLogConvHistory(snes,fnorm,0); 95 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 96 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 97 if (snes->reason) PetscFunctionReturn(0); 98 99 k_restart = 0; 100 l = 0; 101 for (k=1; k < snes->max_its+1; k++) { 102 /* select which vector of the stored subspace will be updated */ 103 ivec = k_restart % ngmres->msize; 104 if (snes->pc && snes->pcside == PC_RIGHT) { 105 ierr = VecCopy(X,XM);CHKERRQ(ierr); 106 ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr); 107 ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr); 108 109 ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 110 ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr); 111 ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 112 113 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 114 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 115 snes->reason = SNES_DIVERGED_INNER; 116 PetscFunctionReturn(0); 117 } 118 ierr = SNESGetFunction(snes->pc,&FPC,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 119 ierr = VecCopy(FPC,FM);CHKERRQ(ierr); 120 ierr = SNESGetFunctionNorm(snes->pc,&fMnorm);CHKERRQ(ierr); 121 if (ngmres->andersonBeta != 1.0) { 122 VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr); 123 } 124 } else { 125 ierr = VecCopy(F,FM);CHKERRQ(ierr); 126 ierr = VecCopy(X,XM);CHKERRQ(ierr); 127 ierr = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr); 128 fMnorm = fnorm; 129 } 130 131 ierr = SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA); 132 133 if (l < ngmres->msize)l++; 134 k_restart++; 135 ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM); 136 137 ierr = VecNorm(FA,NORM_2,&fnorm);CHKERRQ(ierr); 138 ierr = VecCopy(XA,X);CHKERRQ(ierr); 139 ierr = VecCopy(FA,F);CHKERRQ(ierr); 140 141 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 142 snes->iter = k; 143 snes->norm = fnorm; 144 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 145 SNESLogConvHistory(snes,snes->norm,snes->iter); 146 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 147 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 148 if (snes->reason) PetscFunctionReturn(0); 149 } 150 snes->reason = SNES_DIVERGED_MAX_IT; 151 PetscFunctionReturn(0); 152 } 153 154 /*MC 155 SNESAnderson - Anderson Mixing method. 156 157 Level: beginner 158 159 Options Database: 160 + -snes_anderson_m - Number of stored previous solutions and residuals 161 . -snes_anderson_beta - Relaxation parameter; X_{update} = X + \beta F 162 - -snes_anderson_monitor - Prints relevant information about the ngmres iteration 163 164 Notes: 165 166 The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized 167 optimization problem at each iteration. 168 169 References: 170 171 "D. G. Anderson. Iterative procedures for nonlinear integral equations. 172 J. Assoc. Comput. Mach., 12:547–560, 1965." 173 174 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 175 M*/ 176 177 EXTERN_C_BEGIN 178 #undef __FUNCT__ 179 #define __FUNCT__ "SNESCreate_Anderson" 180 PetscErrorCode SNESCreate_Anderson(SNES snes) 181 { 182 SNES_NGMRES *ngmres; 183 PetscErrorCode ierr; 184 185 PetscFunctionBegin; 186 snes->ops->destroy = SNESDestroy_NGMRES; 187 snes->ops->setup = SNESSetUp_NGMRES; 188 snes->ops->setfromoptions = SNESSetFromOptions_Anderson; 189 snes->ops->view = SNESView_NGMRES; 190 snes->ops->solve = SNESSolve_Anderson; 191 snes->ops->reset = SNESReset_NGMRES; 192 193 snes->usespc = PETSC_TRUE; 194 snes->usesksp = PETSC_FALSE; 195 196 ierr = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr); 197 snes->data = (void*) ngmres; 198 ngmres->msize = 30; 199 200 if (!snes->tolerancesset) { 201 snes->max_funcs = 30000; 202 snes->max_its = 10000; 203 } 204 205 ngmres->additive_linesearch = PETSC_NULL; 206 207 ngmres->restart_it = 2; 208 ngmres->restart_periodic = 30; 209 ngmres->gammaA = 2.0; 210 ngmres->gammaC = 2.0; 211 ngmres->deltaB = 0.9; 212 ngmres->epsilonB = 0.1; 213 214 ngmres->andersonBeta = 1.0; 215 PetscFunctionReturn(0); 216 } 217 EXTERN_C_END 218