1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2
SNESSetFromOptions_Anderson(SNES snes,PetscOptionItems PetscOptionsObject)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
SNESSolve_Anderson(SNES snes)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->work[2];
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) PetscCall(SNESComputeFunction(snes, X, F));
74 else snes->vec_func_init_set = PETSC_FALSE;
75
76 PetscCall(VecNorm(F, NORM_2, &fnorm));
77 SNESCheckFunctionDomainError(snes, fnorm);
78 }
79 fminnorm = fnorm;
80
81 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
82 snes->norm = fnorm;
83 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
84 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
85 PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
86 PetscCall(SNESMonitor(snes, 0, fnorm));
87 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
88 PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, F, fnorm, X));
89
90 k_restart = 0;
91 l = 0;
92 ivec = 0;
93 for (k = 1; k < snes->max_its + 1; k++) {
94 /* Call general purpose update function */
95 PetscTryTypeMethod(snes, update, snes->iter);
96
97 /* select which vector of the stored subspace will be updated */
98 if (snes->npc && snes->npcside == PC_RIGHT) {
99 PetscCall(VecCopy(X, XM));
100 PetscCall(SNESSetInitialFunction(snes->npc, F));
101
102 PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, XM, B, 0));
103 PetscCall(SNESSolve(snes->npc, B, XM));
104 PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, XM, B, 0));
105
106 PetscCall(SNESGetConvergedReason(snes->npc, &reason));
107 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
108 snes->reason = SNES_DIVERGED_INNER;
109 PetscFunctionReturn(PETSC_SUCCESS);
110 }
111 PetscCall(SNESGetNPCFunction(snes, FM, &fMnorm));
112 PetscCall(VecAXPBY(XM, 1.0 - ngmres->andersonBeta, ngmres->andersonBeta, X));
113 } else {
114 PetscCall(VecCopy(F, FM));
115 PetscCall(VecWAXPY(XM, -ngmres->andersonBeta, FM, X));
116 fMnorm = fnorm;
117 }
118
119 PetscCall(SNESNGMRESFormCombinedSolution_Private(snes, ivec, l, XM, FM, fMnorm, X, XA, FA));
120 ivec = k_restart % ngmres->msize;
121 if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
122 PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, &dnorm, &dminnorm, NULL, NULL, NULL, &xnorm, &fAnorm, &ynorm));
123 PetscCall(SNESNGMRESSelectRestart_Private(snes, l, snes->norm, fMnorm, fnorm, dnorm, fminnorm, dminnorm, &selectRestart));
124 /* if the restart conditions persist for more than restart_it iterations, restart. */
125 if (selectRestart) restart_count++;
126 else restart_count = 0;
127 } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
128 PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, NULL, NULL, NULL, &xnorm, &fAnorm, &ynorm));
129 if (k_restart > ngmres->restart_periodic) {
130 if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %" PetscInt_FMT " iterations\n", k_restart));
131 restart_count = ngmres->restart_it;
132 }
133 } else {
134 PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, NULL, NULL, NULL, &xnorm, &fAnorm, &ynorm));
135 }
136 /* restart after restart conditions have persisted for a fixed number of iterations */
137 if (restart_count >= ngmres->restart_it) {
138 if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %" PetscInt_FMT "\n", k_restart));
139 restart_count = 0;
140 k_restart = 0;
141 l = 0;
142 ivec = 0;
143 } else {
144 if (l < ngmres->msize) l++;
145 k_restart++;
146 PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, FM, fnorm, XM));
147 }
148
149 fnorm = fAnorm;
150 if (fminnorm > fnorm) fminnorm = fnorm;
151
152 PetscCall(VecCopy(XA, X));
153 PetscCall(VecCopy(FA, F));
154
155 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
156 snes->iter = k;
157 snes->norm = fnorm;
158 snes->xnorm = xnorm;
159 snes->ynorm = ynorm;
160 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
161 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
162 PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
163 PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
164 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
165 }
166 snes->reason = SNES_DIVERGED_MAX_IT;
167 PetscFunctionReturn(PETSC_SUCCESS);
168 }
169
170 /*MC
171 SNESANDERSON - Implements the Anderson Mixing nonlinear solver {cite}`anderson1965`, {cite}`bruneknepleysmithtu15`
172
173 Level: beginner
174
175 Options Database Keys:
176 + -snes_anderson_m - Number of stored previous solutions and residuals
177 . -snes_anderson_beta - Anderson mixing parameter
178 . -snes_anderson_restart_type - Type of restart (see `SNESNGMRES`)
179 . -snes_anderson_restart_it - Number of iterations of restart conditions before restart
180 . -snes_anderson_restart - Number of iterations before periodic restart
181 - -snes_anderson_monitor - Prints relevant information about the Anderson mixing iteration
182
183 Notes:
184 The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
185 optimization problem at each iteration.
186
187 Very similar to the `SNESNGMRES` algorithm.
188
189 This algorithm ignores any Jacobian provided with `SNESSetJacobian()`
190
191 Only supports left non-linear preconditioning.
192
193 .seealso: [](ch_snes), `SNESNGMRES`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`
194 M*/
195
SNESCreate_Anderson(SNES snes)196 PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
197 {
198 SNES_NGMRES *ngmres;
199 SNESLineSearch linesearch;
200
201 PetscFunctionBegin;
202 snes->ops->destroy = SNESDestroy_NGMRES;
203 snes->ops->setup = SNESSetUp_NGMRES;
204 snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
205 snes->ops->view = SNESView_NGMRES;
206 snes->ops->solve = SNESSolve_Anderson;
207 snes->ops->reset = SNESReset_NGMRES;
208
209 snes->usesnpc = PETSC_TRUE;
210 snes->usesksp = PETSC_FALSE;
211 snes->npcside = PC_RIGHT;
212
213 snes->alwayscomputesfinalresidual = PETSC_TRUE;
214
215 PetscCall(PetscNew(&ngmres));
216 snes->data = (void *)ngmres;
217 ngmres->msize = 30;
218
219 PetscCall(SNESParametersInitialize(snes));
220 PetscObjectParameterSetDefault(snes, max_funcs, 30000);
221 PetscObjectParameterSetDefault(snes, max_its, 10000);
222
223 PetscCall(SNESGetLineSearch(snes, &linesearch));
224 if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
225
226 ngmres->additive_linesearch = NULL;
227 ngmres->approxfunc = PETSC_FALSE;
228 ngmres->restart_type = SNES_NGMRES_RESTART_NONE;
229 ngmres->restart_it = 2;
230 ngmres->restart_periodic = 30;
231 ngmres->gammaA = 2.0;
232 ngmres->gammaC = 2.0;
233 ngmres->deltaB = 0.9;
234 ngmres->epsilonB = 0.1;
235
236 ngmres->andersonBeta = 1.0;
237 PetscFunctionReturn(PETSC_SUCCESS);
238 }
239