Actual source code: anderson.c
petsc-3.7.3 2016-08-01
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(PetscOptionItems *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;
61: 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);
63: PetscCitationsRegister(SNESCitation,&SNEScite);
64: /* variable initialization */
65: snes->reason = SNES_CONVERGED_ITERATING;
66: X = snes->vec_sol;
67: F = snes->vec_func;
68: B = snes->vec_rhs;
69: XA = snes->vec_sol_update;
70: FA = snes->work[0];
71: D = snes->work[1];
73: /* work for the line search */
74: XM = snes->work[3];
75: FM = snes->work[4];
77: PetscObjectSAWsTakeAccess((PetscObject)snes);
78: snes->iter = 0;
79: snes->norm = 0.;
80: PetscObjectSAWsGrantAccess((PetscObject)snes);
82: /* initialization */
84: /* r = F(x) */
86: if (snes->pc && snes->pcside == PC_LEFT) {
87: SNESApplyNPC(snes,X,NULL,F);
88: SNESGetConvergedReason(snes->pc,&reason);
89: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
90: snes->reason = SNES_DIVERGED_INNER;
91: return(0);
92: }
93: VecNorm(F,NORM_2,&fnorm);
94: } else {
95: if (!snes->vec_func_init_set) {
96: SNESComputeFunction(snes,X,F);
97: } else snes->vec_func_init_set = PETSC_FALSE;
99: VecNorm(F,NORM_2,&fnorm);
100: SNESCheckFunctionNorm(snes,fnorm);
101: }
102: fminnorm = fnorm;
104: PetscObjectSAWsTakeAccess((PetscObject)snes);
105: snes->norm = fnorm;
106: PetscObjectSAWsGrantAccess((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: ivec = 0;
115: for (k=1; k < snes->max_its+1; k++) {
116: /* select which vector of the stored subspace will be updated */
117: if (snes->pc && snes->pcside == PC_RIGHT) {
118: VecCopy(X,XM);
119: SNESSetInitialFunction(snes->pc,F);
121: PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);
122: SNESSolve(snes->pc,B,XM);
123: PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);
125: SNESGetConvergedReason(snes->pc,&reason);
126: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
127: snes->reason = SNES_DIVERGED_INNER;
128: return(0);
129: }
130: SNESGetNPCFunction(snes,FM,&fMnorm);
131: if (ngmres->andersonBeta != 1.0) {
132: VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);
133: }
134: } else {
135: VecCopy(F,FM);
136: VecCopy(X,XM);
137: VecAXPY(XM,-ngmres->andersonBeta,FM);
138: fMnorm = fnorm;
139: }
141: SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);
142: ivec = k_restart % ngmres->msize;
143: if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
144: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
145: SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);
146: /* if the restart conditions persist for more than restart_it iterations, restart. */
147: if (selectRestart) restart_count++;
148: else restart_count = 0;
149: } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
150: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
151: if (k_restart > ngmres->restart_periodic) {
152: if (ngmres->monitor) PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);
153: restart_count = ngmres->restart_it;
154: }
155: } else {
156: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
157: }
158: /* restart after restart conditions have persisted for a fixed number of iterations */
159: if (restart_count >= ngmres->restart_it) {
160: if (ngmres->monitor) {
161: PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);
162: }
163: restart_count = 0;
164: k_restart = 0;
165: l = 0;
166: ivec = 0;
167: } else {
168: if (l < ngmres->msize) l++;
169: k_restart++;
170: SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);
171: }
173: fnorm = fAnorm;
174: if (fminnorm > fnorm) fminnorm = fnorm;
176: VecCopy(XA,X);
177: VecCopy(FA,F);
179: PetscObjectSAWsTakeAccess((PetscObject)snes);
180: snes->iter = k;
181: snes->norm = fnorm;
182: PetscObjectSAWsGrantAccess((PetscObject)snes);
183: SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
184: SNESMonitor(snes,snes->iter,snes->norm);
185: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
186: if (snes->reason) return(0);
187: }
188: snes->reason = SNES_DIVERGED_MAX_IT;
189: return(0);
190: }
192: /*MC
193: SNESANDERSON - Anderson Mixing method.
195: Level: beginner
197: Options Database:
198: + -snes_anderson_m - Number of stored previous solutions and residuals
199: . -snes_anderson_beta - Relaxation parameter; X_{update} = X + \beta F
200: . -snes_anderson_restart_type - Type of restart (see SNESNGMRES)
201: . -snes_anderson_restart_it - Number of iterations of restart conditions before restart
202: . -snes_anderson_restart - Number of iterations before periodic restart
203: - -snes_anderson_monitor - Prints relevant information about the ngmres iteration
205: Notes:
207: The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
208: optimization problem at each iteration.
210: Very similar to the SNESNGMRES algorithm.
212: References:
213: + 1. - D. G. Anderson. Iterative procedures for nonlinear integral equations.
214: J. Assoc. Comput. Mach., 12, 1965."
215: - 2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
216: SIAM Review, 57(4), 2015
219: .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
220: M*/
224: PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
225: {
226: SNES_NGMRES *ngmres;
230: snes->ops->destroy = SNESDestroy_NGMRES;
231: snes->ops->setup = SNESSetUp_NGMRES;
232: snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
233: snes->ops->view = SNESView_NGMRES;
234: snes->ops->solve = SNESSolve_Anderson;
235: snes->ops->reset = SNESReset_NGMRES;
237: snes->usespc = PETSC_TRUE;
238: snes->usesksp = PETSC_FALSE;
239: snes->pcside = PC_RIGHT;
241: PetscNewLog(snes,&ngmres);
242: snes->data = (void*) ngmres;
243: ngmres->msize = 30;
245: if (!snes->tolerancesset) {
246: snes->max_funcs = 30000;
247: snes->max_its = 10000;
248: }
250: ngmres->additive_linesearch = NULL;
251: ngmres->approxfunc = PETSC_FALSE;
252: ngmres->restart_type = SNES_NGMRES_RESTART_NONE;
253: ngmres->restart_it = 2;
254: ngmres->restart_periodic = 30;
255: ngmres->gammaA = 2.0;
256: ngmres->gammaC = 2.0;
257: ngmres->deltaB = 0.9;
258: ngmres->epsilonB = 0.1;
260: ngmres->andersonBeta = 1.0;
261: return(0);
262: }