Actual source code: anderson.c

petsc-3.4.5 2014-06-29
  1: #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/

  3: extern PetscErrorCode SNESDestroy_NGMRES(SNES);
  4: extern PetscErrorCode SNESReset_NGMRES(SNES);
  5: extern PetscErrorCode SNESSetUp_NGMRES(SNES);
  6: extern PetscErrorCode SNESView_NGMRES(SNES,PetscViewer);

  8: PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];

 12: PetscErrorCode SNESSetFromOptions_Anderson(SNES snes)
 13: {
 14:   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
 16:   PetscBool      debug;
 17:   SNESLineSearch linesearch;

 20:   PetscOptionsHead("SNES NGMRES options");
 21:   PetscOptionsInt("-snes_anderson_m",            "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);
 22:   PetscOptionsReal("-snes_anderson_beta",        "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);
 23:   PetscOptionsBool("-snes_anderson_monitor",     "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);
 24:   PetscOptionsInt("-snes_anderson_restart",      "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);
 25:   PetscOptionsInt("-snes_anderson_restart_it",   "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);
 26:   PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,
 27:                           (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);
 28:   if (debug) {
 29:     ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
 30:   }
 31:   PetscOptionsTail();
 32:   /* set the default type of the line search if the user hasn't already. */
 33:   if (!snes->linesearch) {
 34:     SNESGetLineSearch(snes,&linesearch);
 35:     SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
 36:   }
 37:   return(0);
 38: }

 42: PetscErrorCode SNESSolve_Anderson(SNES snes)
 43: {
 44:   SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
 45:   /* present solution, residual, and preconditioned residual */
 46:   Vec X,F,B,D;

 48:   /* candidate linear combination answers */
 49:   Vec XA,FA,XM,FM,FPC;

 51:   /* coefficients and RHS to the minimization problem */
 52:   PetscReal fnorm,fMnorm;
 53:   PetscReal dnorm,dminnorm=0.0,fminnorm;
 54:   PetscInt  restart_count=0;
 55:   PetscInt  k,k_restart,l,ivec;

 57:   PetscBool selectRestart;

 59:   SNESConvergedReason reason;
 60:   PetscErrorCode      ierr;

 63:   /* variable initialization */
 64:   snes->reason = SNES_CONVERGED_ITERATING;
 65:   X            = snes->vec_sol;
 66:   F            = snes->vec_func;
 67:   B            = snes->vec_rhs;
 68:   XA           = snes->vec_sol_update;
 69:   FA           = snes->work[0];
 70:   D            = snes->work[1];

 72:   /* work for the line search */
 73:   XM = snes->work[3];
 74:   FM = snes->work[4];

 76:   PetscObjectAMSTakeAccess((PetscObject)snes);
 77:   snes->iter = 0;
 78:   snes->norm = 0.;
 79:   PetscObjectAMSGrantAccess((PetscObject)snes);

 81:   /* initialization */

 83:   /* r = F(x) */
 84:   if (!snes->vec_func_init_set) {
 85:     SNESComputeFunction(snes,X,F);
 86:     if (snes->domainerror) {
 87:       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
 88:       return(0);
 89:     }
 90:   } else snes->vec_func_init_set = PETSC_FALSE;

 92:   if (!snes->norm_init_set) {
 93:     VecNorm(F,NORM_2,&fnorm);
 94:     if (PetscIsInfOrNanReal(fnorm)) {
 95:       snes->reason = SNES_DIVERGED_FNORM_NAN;
 96:       return(0);
 97:     }
 98:   } else {
 99:     fnorm               = snes->norm_init;
100:     snes->norm_init_set = PETSC_FALSE;
101:   }
102:   fminnorm = fnorm;

104:   PetscObjectAMSTakeAccess((PetscObject)snes);
105:   snes->norm = fnorm;
106:   PetscObjectAMSGrantAccess((PetscObject)snes);
107:   SNESLogConvergenceHistory(snes,fnorm,0);
108:   SNESMonitor(snes,0,fnorm);
109:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
110:   if (snes->reason) return(0);

112:   k_restart = 0;
113:   l         = 0;
114:   for (k=1; k < snes->max_its+1; k++) {
115:     /* select which vector of the stored subspace will be updated */
116:     ivec = k_restart % ngmres->msize;
117:     if (snes->pc && snes->pcside == PC_RIGHT) {
118:       VecCopy(X,XM);
119:       SNESSetInitialFunction(snes->pc,F);
120:       SNESSetInitialFunctionNorm(snes->pc,fnorm);

122:       PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);
123:       SNESSolve(snes->pc,B,XM);
124:       PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);

126:       SNESGetConvergedReason(snes->pc,&reason);
127:       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
128:         snes->reason = SNES_DIVERGED_INNER;
129:         return(0);
130:       }
131:       SNESGetFunction(snes->pc,&FPC,NULL,NULL);
132:       VecCopy(FPC,FM);
133:       SNESGetFunctionNorm(snes->pc,&fMnorm);
134:       if (ngmres->andersonBeta != 1.0) {
135:         VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);
136:       }
137:     } else {
138:       VecCopy(F,FM);
139:       VecCopy(X,XM);
140:       VecAXPY(XM,-ngmres->andersonBeta,FM);
141:       fMnorm = fnorm;
142:     }

144:     SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA);
145:     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
146:       SNESNGMRESCalculateDifferences_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL);
147:       SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);
148:       /* if the restart conditions persist for more than restart_it iterations, restart. */
149:       if (selectRestart) restart_count++;
150:       else restart_count = 0;
151:     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
152:       if (k_restart > ngmres->restart_periodic) {
153:         if (ngmres->monitor) PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);
154:         restart_count = ngmres->restart_it;
155:       }
156:     }
157:     /* restart after restart conditions have persisted for a fixed number of iterations */
158:     if (restart_count >= ngmres->restart_it) {
159:       if (ngmres->monitor) {
160:         PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);
161:       }
162:       restart_count = 0;
163:       k_restart     = 0;
164:       l             = 0;
165:     } else {
166:       if (l < ngmres->msize) l++;
167:       k_restart++;
168:       SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);
169:     }

171:     VecNorm(FA,NORM_2,&fnorm);
172:     if (fminnorm > fnorm) fminnorm = fnorm;

174:     VecCopy(XA,X);
175:     VecCopy(FA,F);

177:     PetscObjectAMSTakeAccess((PetscObject)snes);
178:     snes->iter = k;
179:     snes->norm = fnorm;
180:     PetscObjectAMSGrantAccess((PetscObject)snes);
181:     SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
182:     SNESMonitor(snes,snes->iter,snes->norm);
183:     (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
184:     if (snes->reason) return(0);
185:   }
186:   snes->reason = SNES_DIVERGED_MAX_IT;
187:   return(0);
188: }

190: /*MC
191:   SNESAnderson - Anderson Mixing method.

193:    Level: beginner