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