Actual source code: anderson.c
petsc-3.9.4 2018-09-11
1: #include <../src/snes/impls/ngmres/snesngmres.h>
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: static PetscErrorCode SNESSetFromOptions_Anderson(PetscOptionItems *PetscOptionsObject,SNES snes)
9: {
10: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
12: PetscBool monitor = PETSC_FALSE;
13: SNESLineSearch linesearch;
16: PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");
17: PetscOptionsInt("-snes_anderson_m", "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);
18: PetscOptionsReal("-snes_anderson_beta", "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);
19: PetscOptionsInt("-snes_anderson_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);
20: PetscOptionsInt("-snes_anderson_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);
21: PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,(PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);
22: PetscOptionsBool("-snes_anderson_monitor", "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&monitor,NULL);
23: if (monitor) {
24: ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
25: }
26: PetscOptionsTail();
27: /* set the default type of the line search if the user hasn't already. */
28: if (!snes->linesearch) {
29: SNESGetLineSearch(snes,&linesearch);
30: SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
31: }
32: return(0);
33: }
35: static PetscErrorCode SNESSolve_Anderson(SNES snes)
36: {
37: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
38: /* present solution, residual, and preconditioned residual */
39: Vec X,F,B,D;
40: /* candidate linear combination answers */
41: Vec XA,FA,XM,FM;
43: /* coefficients and RHS to the minimization problem */
44: PetscReal fnorm,fMnorm,fAnorm;
45: PetscReal xnorm,ynorm;
46: PetscReal dnorm,dminnorm=0.0,fminnorm;
47: PetscInt restart_count=0;
48: PetscInt k,k_restart,l,ivec;
49: PetscBool selectRestart;
50: SNESConvergedReason reason;
51: PetscErrorCode ierr;
54: 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);
56: PetscCitationsRegister(SNESCitation,&SNEScite);
57: /* variable initialization */
58: snes->reason = SNES_CONVERGED_ITERATING;
59: X = snes->vec_sol;
60: F = snes->vec_func;
61: B = snes->vec_rhs;
62: XA = snes->vec_sol_update;
63: FA = snes->work[0];
64: D = snes->work[1];
66: /* work for the line search */
67: XM = snes->work[3];
68: FM = snes->work[4];
70: PetscObjectSAWsTakeAccess((PetscObject)snes);
71: snes->iter = 0;
72: snes->norm = 0.;
73: PetscObjectSAWsGrantAccess((PetscObject)snes);
75: /* initialization */
77: /* r = F(x) */
79: if (snes->npc && snes->npcside== PC_LEFT) {
80: SNESApplyNPC(snes,X,NULL,F);
81: SNESGetConvergedReason(snes->npc,&reason);
82: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
83: snes->reason = SNES_DIVERGED_INNER;
84: return(0);
85: }
86: VecNorm(F,NORM_2,&fnorm);
87: } else {
88: if (!snes->vec_func_init_set) {
89: SNESComputeFunction(snes,X,F);
90: } else snes->vec_func_init_set = PETSC_FALSE;
92: VecNorm(F,NORM_2,&fnorm);
93: SNESCheckFunctionNorm(snes,fnorm);
94: }
95: fminnorm = fnorm;
97: PetscObjectSAWsTakeAccess((PetscObject)snes);
98: snes->norm = fnorm;
99: PetscObjectSAWsGrantAccess((PetscObject)snes);
100: SNESLogConvergenceHistory(snes,fnorm,0);
101: SNESMonitor(snes,0,fnorm);
102: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
103: if (snes->reason) return(0);
105: k_restart = 0;
106: l = 0;
107: ivec = 0;
108: for (k=1; k < snes->max_its+1; k++) {
109: /* select which vector of the stored subspace will be updated */
110: if (snes->npc && snes->npcside== PC_RIGHT) {
111: VecCopy(X,XM);
112: SNESSetInitialFunction(snes->npc,F);
114: PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0);
115: SNESSolve(snes->npc,B,XM);
116: PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0);
118: SNESGetConvergedReason(snes->npc,&reason);
119: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
120: snes->reason = SNES_DIVERGED_INNER;
121: return(0);
122: }
123: SNESGetNPCFunction(snes,FM,&fMnorm);
124: if (ngmres->andersonBeta != 1.0) {
125: VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);
126: }
127: } else {
128: VecCopy(F,FM);
129: VecCopy(X,XM);
130: VecAXPY(XM,-ngmres->andersonBeta,FM);
131: fMnorm = fnorm;
132: }
134: SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);
135: ivec = k_restart % ngmres->msize;
136: if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
137: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
138: SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);
139: /* if the restart conditions persist for more than restart_it iterations, restart. */
140: if (selectRestart) restart_count++;
141: else restart_count = 0;
142: } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
143: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
144: if (k_restart > ngmres->restart_periodic) {
145: if (ngmres->monitor) PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);
146: restart_count = ngmres->restart_it;
147: }
148: } else {
149: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);
150: }
151: /* restart after restart conditions have persisted for a fixed number of iterations */
152: if (restart_count >= ngmres->restart_it) {
153: if (ngmres->monitor) {
154: PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);
155: }
156: restart_count = 0;
157: k_restart = 0;
158: l = 0;
159: ivec = 0;
160: } else {
161: if (l < ngmres->msize) l++;
162: k_restart++;
163: SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);
164: }
166: fnorm = fAnorm;
167: if (fminnorm > fnorm) fminnorm = fnorm;
169: VecCopy(XA,X);
170: VecCopy(FA,F);
172: PetscObjectSAWsTakeAccess((PetscObject)snes);
173: snes->iter = k;
174: snes->norm = fnorm;
175: PetscObjectSAWsGrantAccess((PetscObject)snes);
176: SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
177: SNESMonitor(snes,snes->iter,snes->norm);
178: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
179: if (snes->reason) return(0);
180: }
181: snes->reason = SNES_DIVERGED_MAX_IT;
182: return(0);
183: }
185: /*MC
186: SNESANDERSON - Anderson Mixing method.
188: Level: beginner
190: Options Database:
191: + -snes_anderson_m - Number of stored previous solutions and residuals
192: . -snes_anderson_beta - Anderson mixing parameter
193: . -snes_anderson_restart_type - Type of restart (see SNESNGMRES)
194: . -snes_anderson_restart_it - Number of iterations of restart conditions before restart
195: . -snes_anderson_restart - Number of iterations before periodic restart
196: - -snes_anderson_monitor - Prints relevant information about the ngmres iteration
198: Notes:
200: The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
201: optimization problem at each iteration.
203: Very similar to the SNESNGMRES algorithm.
205: References:
206: + 1. - D. G. Anderson. Iterative procedures for nonlinear integral equations.
207: J. Assoc. Comput. Mach., 12, 1965."
208: - 2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
209: SIAM Review, 57(4), 2015
212: .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
213: M*/
215: PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
216: {
217: SNES_NGMRES *ngmres;
221: snes->ops->destroy = SNESDestroy_NGMRES;
222: snes->ops->setup = SNESSetUp_NGMRES;
223: snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
224: snes->ops->view = SNESView_NGMRES;
225: snes->ops->solve = SNESSolve_Anderson;
226: snes->ops->reset = SNESReset_NGMRES;
228: snes->usesnpc = PETSC_TRUE;
229: snes->usesksp = PETSC_FALSE;
230: snes->npcside = PC_RIGHT;
232: snes->alwayscomputesfinalresidual = PETSC_TRUE;
234: PetscNewLog(snes,&ngmres);
235: snes->data = (void*) ngmres;
236: ngmres->msize = 30;
238: if (!snes->tolerancesset) {
239: snes->max_funcs = 30000;
240: snes->max_its = 10000;
241: }
243: ngmres->additive_linesearch = NULL;
244: ngmres->approxfunc = PETSC_FALSE;
245: ngmres->restart_type = SNES_NGMRES_RESTART_NONE;
246: ngmres->restart_it = 2;
247: ngmres->restart_periodic = 30;
248: ngmres->gammaA = 2.0;
249: ngmres->gammaC = 2.0;
250: ngmres->deltaB = 0.9;
251: ngmres->epsilonB = 0.1;
253: ngmres->andersonBeta = 1.0;
254: return(0);
255: }