Actual source code: anderson.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  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(PetscOptions *PetscOptionsObject,SNES snes)
 13: {
 14:   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
 16:   PetscBool      debug;
 17:   SNESLineSearch linesearch;

 20:   PetscOptionsHead(PetscOptionsObject,"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;
 47:   /* candidate linear combination answers */
 48:   Vec                 XA,FA,XM,FM;

 50:   /* coefficients and RHS to the minimization problem */
 51:   PetscReal           fnorm,fMnorm,fAnorm;
 52:   PetscReal           xnorm,ynorm;
 53:   PetscReal           dnorm,dminnorm=0.0,fminnorm;
 54:   PetscInt            restart_count=0;
 55:   PetscInt            k,k_restart,l,ivec;
 56:   PetscBool           selectRestart;
 57:   SNESConvergedReason reason;
 58:   PetscErrorCode      ierr;


 62:   if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
 63:     SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
 64:   }

 66:   PetscCitationsRegister(SNESCitation,&SNEScite);
 67:   /* variable initialization */
 68:   snes->reason = SNES_CONVERGED_ITERATING;
 69:   X            = snes->vec_sol;
 70:   F            = snes->vec_func;
 71:   B            = snes->vec_rhs;
 72:   XA           = snes->vec_sol_update;
 73:   FA           = snes->work[0];
 74:   D            = snes->work[1];

 76:   /* work for the line search */
 77:   XM = snes->work[3];
 78:   FM = snes->work[4];

 80:   PetscObjectSAWsTakeAccess((PetscObject)snes);
 81:   snes->iter = 0;
 82:   snes->norm = 0.;
 83:   PetscObjectSAWsGrantAccess((PetscObject)snes);

 85:   /* initialization */

 87:   /* r = F(x) */

 89:   if (snes->pc && snes->pcside == PC_LEFT) {
 90:     SNESApplyNPC(snes,X,NULL,F);
 91:     SNESGetConvergedReason(snes->pc,&reason);
 92:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
 93:       snes->reason = SNES_DIVERGED_INNER;
 94:       return(0);
 95:     }
 96:     VecNorm(F,NORM_2,&fnorm);
 97:   } else {
 98:     if (!snes->vec_func_init_set) {
 99:       SNESComputeFunction(snes,X,F);
100:     } else snes->vec_func_init_set = PETSC_FALSE;

102:     VecNorm(F,NORM_2,&fnorm);
103:     SNESCheckFunctionNorm(snes,fnorm);
104:   }
105:   fminnorm = fnorm;

107:   PetscObjectSAWsTakeAccess((PetscObject)snes);
108:   snes->norm = fnorm;
109:   PetscObjectSAWsGrantAccess((PetscObject)snes);
110:   SNESLogConvergenceHistory(snes,fnorm,0);
111:   SNESMonitor(snes,0,fnorm);
112:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
113:   if (snes->reason) return(0);

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

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

128:       SNESGetConvergedReason(snes->pc,&reason);
129:       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
130:         snes->reason = SNES_DIVERGED_INNER;
131:         return(0);
132:       }
133:       SNESGetNPCFunction(snes,FM,&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,ivec,l,XM,FM,fMnorm,X,XA,FA);
145:     ivec = k_restart % ngmres->msize;
146:     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
147:       SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
148:       SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);
149:       /* if the restart conditions persist for more than restart_it iterations, restart. */
150:       if (selectRestart) restart_count++;
151:       else restart_count = 0;
152:     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
153:       SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
154:       if (k_restart > ngmres->restart_periodic) {
155:         if (ngmres->monitor) PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);
156:         restart_count = ngmres->restart_it;
157:       }
158:     } else {
159:       SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
160:     }
161:     /* restart after restart conditions have persisted for a fixed number of iterations */
162:     if (restart_count >= ngmres->restart_it) {
163:       if (ngmres->monitor) {
164:         PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);
165:       }
166:       restart_count = 0;
167:       k_restart     = 0;
168:       l             = 0;
169:       ivec          = 0;
170:     } else {
171:       if (l < ngmres->msize) l++;
172:       k_restart++;
173:       SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);
174:     }

176:     fnorm = fAnorm;
177:     if (fminnorm > fnorm) fminnorm = fnorm;

179:     VecCopy(XA,X);
180:     VecCopy(FA,F);

182:     PetscObjectSAWsTakeAccess((PetscObject)snes);
183:     snes->iter = k;
184:     snes->norm = fnorm;
185:     PetscObjectSAWsGrantAccess((PetscObject)snes);
186:     SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
187:     SNESMonitor(snes,snes->iter,snes->norm);
188:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
189:     if (snes->reason) return(0);
190:   }
191:   snes->reason = SNES_DIVERGED_MAX_IT;
192:   return(0);
193: }

195: /*MC
196:   SNESAnderson - Anderson Mixing method.

198:    Level: beginner

200:    Options Database:
201: +  -snes_anderson_m                - Number of stored previous solutions and residuals
202: .  -snes_anderson_beta             - Relaxation parameter; X_{update} = X + \beta F
203: .  -snes_anderson_restart_type     - Type of restart (see SNESNGMRES)
204: .  -snes_anderson_restart_it       - Number of iterations of restart conditions before restart
205: .  -snes_anderson_restart          - Number of iterations before periodic restart
206: -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration

208:    Notes:

210:    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
211:    optimization problem at each iteration.

213:    References:

215:     "D. G. Anderson. Iterative procedures for nonlinear integral equations.
216:     J. Assoc. Comput. Mach., 12:547-560, 1965."

218: .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
219: M*/

223: PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
224: {
225:   SNES_NGMRES    *ngmres;

229:   snes->ops->destroy        = SNESDestroy_NGMRES;
230:   snes->ops->setup          = SNESSetUp_NGMRES;
231:   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
232:   snes->ops->view           = SNESView_NGMRES;
233:   snes->ops->solve          = SNESSolve_Anderson;
234:   snes->ops->reset          = SNESReset_NGMRES;

236:   snes->usespc  = PETSC_TRUE;
237:   snes->usesksp = PETSC_FALSE;
238:   snes->pcside  = PC_RIGHT;

240:   PetscNewLog(snes,&ngmres);
241:   snes->data    = (void*) ngmres;
242:   ngmres->msize = 30;

244:   if (!snes->tolerancesset) {
245:     snes->max_funcs = 30000;
246:     snes->max_its   = 10000;
247:   }

249:   ngmres->additive_linesearch = NULL;
250:   ngmres->approxfunc       = PETSC_FALSE;
251:   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
252:   ngmres->restart_it       = 2;
253:   ngmres->restart_periodic = 30;
254:   ngmres->gammaA           = 2.0;
255:   ngmres->gammaC           = 2.0;
256:   ngmres->deltaB           = 0.9;
257:   ngmres->epsilonB         = 0.1;

259:   ngmres->andersonBeta = 1.0;
260:   return(0);
261: }