xref: /petsc/src/snes/impls/qn/qn.c (revision 6a4a1270853b5b1995c94de98f44d28bec57470a)
1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2 #include <petscdm.h>
3 
4 #define H(i, j) qn->dXdFmat[i * qn->m + j]
5 
6 const char *const SNESQNScaleTypes[]   = {"DEFAULT", "NONE", "SCALAR", "DIAGONAL", "JACOBIAN", "SNESQNScaleType", "SNES_QN_SCALING_", NULL};
7 const char *const SNESQNRestartTypes[] = {"DEFAULT", "NONE", "POWELL", "PERIODIC", "SNESQNRestartType", "SNES_QN_RESTART_", NULL};
8 const char *const SNESQNTypes[]        = {"LBFGS", "BROYDEN", "BADBROYDEN", "SNESQNType", "SNES_QN_", NULL};
9 
10 typedef struct {
11   Mat               B;      /* Quasi-Newton approximation Matrix (MATLMVM) */
12   PetscInt          m;      /* The number of kept previous steps */
13   PetscReal        *lambda; /* The line search history of the method */
14   PetscBool         monflg;
15   PetscViewer       monitor;
16   PetscReal         powell_gamma; /* Powell angle restart condition */
17   PetscReal         scaling;      /* scaling of H0 */
18   SNESQNType        type;         /* the type of quasi-newton method used */
19   SNESQNScaleType   scale_type;   /* the type of scaling used */
20   SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */
21 } SNES_QN;
22 
23 static PetscErrorCode SNESSolve_QN(SNES snes) {
24   SNES_QN             *qn = (SNES_QN *)snes->data;
25   Vec                  X, Xold;
26   Vec                  F, W;
27   Vec                  Y, D, Dold;
28   PetscInt             i, i_r;
29   PetscReal            fnorm, xnorm, ynorm, gnorm;
30   SNESLineSearchReason lssucceed;
31   PetscBool            badstep, powell, periodic;
32   PetscScalar          DolddotD, DolddotDold;
33   SNESConvergedReason  reason;
34 
35   /* basically just a regular newton's method except for the application of the Jacobian */
36 
37   PetscFunctionBegin;
38   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);
39 
40   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
41   F    = snes->vec_func;       /* residual vector */
42   Y    = snes->vec_sol_update; /* search direction generated by J^-1D*/
43   W    = snes->work[3];
44   X    = snes->vec_sol; /* solution vector */
45   Xold = snes->work[0];
46 
47   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
48   D    = snes->work[1];
49   Dold = snes->work[2];
50 
51   snes->reason = SNES_CONVERGED_ITERATING;
52 
53   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
54   snes->iter = 0;
55   snes->norm = 0.;
56   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
57 
58   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
59     PetscCall(SNESApplyNPC(snes, X, NULL, F));
60     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
61     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
62       snes->reason = SNES_DIVERGED_INNER;
63       PetscFunctionReturn(0);
64     }
65     PetscCall(VecNorm(F, NORM_2, &fnorm));
66   } else {
67     if (!snes->vec_func_init_set) {
68       PetscCall(SNESComputeFunction(snes, X, F));
69     } else snes->vec_func_init_set = PETSC_FALSE;
70 
71     PetscCall(VecNorm(F, NORM_2, &fnorm));
72     SNESCheckFunctionNorm(snes, fnorm);
73   }
74   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
75     PetscCall(SNESApplyNPC(snes, X, F, D));
76     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
77     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
78       snes->reason = SNES_DIVERGED_INNER;
79       PetscFunctionReturn(0);
80     }
81   } else {
82     PetscCall(VecCopy(F, D));
83   }
84 
85   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
86   snes->norm = fnorm;
87   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
88   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
89   PetscCall(SNESMonitor(snes, 0, fnorm));
90 
91   /* test convergence */
92   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
93   if (snes->reason) PetscFunctionReturn(0);
94 
95   if (snes->npc && snes->npcside == PC_RIGHT) {
96     PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0));
97     PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
98     PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0));
99     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
100     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
101       snes->reason = SNES_DIVERGED_INNER;
102       PetscFunctionReturn(0);
103     }
104     PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
105     PetscCall(VecCopy(F, D));
106   }
107 
108   /* general purpose update */
109   PetscTryTypeMethod(snes, update, snes->iter);
110 
111   /* scale the initial update */
112   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
113     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
114     SNESCheckJacobianDomainerror(snes);
115     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
116     PetscCall(MatLMVMSetJ0KSP(qn->B, snes->ksp));
117   }
118 
119   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
120     /* update QN approx and calculate step */
121     PetscCall(MatLMVMUpdate(qn->B, X, D));
122     PetscCall(MatSolve(qn->B, D, Y));
123 
124     /* line search for lambda */
125     ynorm = 1;
126     gnorm = fnorm;
127     PetscCall(VecCopy(D, Dold));
128     PetscCall(VecCopy(X, Xold));
129     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
130     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
131     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
132     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
133     badstep = PETSC_FALSE;
134     if (lssucceed) {
135       if (++snes->numFailures >= snes->maxFailures) {
136         snes->reason = SNES_DIVERGED_LINE_SEARCH;
137         break;
138       }
139       badstep = PETSC_TRUE;
140     }
141 
142     /* convergence monitoring */
143     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed));
144 
145     if (snes->npc && snes->npcside == PC_RIGHT) {
146       PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0));
147       PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
148       PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0));
149       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
150       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
151         snes->reason = SNES_DIVERGED_INNER;
152         PetscFunctionReturn(0);
153       }
154       PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
155     }
156 
157     PetscCall(SNESSetIterationNumber(snes, i + 1));
158     snes->norm  = fnorm;
159     snes->xnorm = xnorm;
160     snes->ynorm = ynorm;
161 
162     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
163     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
164 
165     /* set parameter for default relative tolerance convergence test */
166     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
167     if (snes->reason) PetscFunctionReturn(0);
168     if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
169       PetscCall(SNESApplyNPC(snes, X, F, D));
170       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
171       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
172         snes->reason = SNES_DIVERGED_INNER;
173         PetscFunctionReturn(0);
174       }
175     } else {
176       PetscCall(VecCopy(F, D));
177     }
178 
179     /* general purpose update */
180     PetscTryTypeMethod(snes, update, snes->iter);
181 
182     /* restart conditions */
183     powell = PETSC_FALSE;
184     if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) {
185       /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */
186       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
187         PetscCall(MatMult(snes->jacobian_pre, Dold, W));
188       } else {
189         PetscCall(VecCopy(Dold, W));
190       }
191       PetscCall(VecDotBegin(W, Dold, &DolddotDold));
192       PetscCall(VecDotBegin(W, D, &DolddotD));
193       PetscCall(VecDotEnd(W, Dold, &DolddotDold));
194       PetscCall(VecDotEnd(W, D, &DolddotD));
195       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma * PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
196     }
197     periodic = PETSC_FALSE;
198     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
199       if (i_r > qn->m - 1) periodic = PETSC_TRUE;
200     }
201     /* restart if either powell or periodic restart is satisfied. */
202     if (badstep || powell || periodic) {
203       if (qn->monflg) {
204         PetscCall(PetscViewerASCIIAddTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
205         if (powell) {
206           PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Powell restart! |%14.12e| > %6.4f*|%14.12e| i_r = %" PetscInt_FMT "\n", (double)PetscRealPart(DolddotD), (double)qn->powell_gamma, (double)PetscRealPart(DolddotDold), i_r));
207         } else {
208           PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %" PetscInt_FMT "\n", i_r));
209         }
210         PetscCall(PetscViewerASCIISubtractTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
211       }
212       i_r = -1;
213       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
214         PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
215         SNESCheckJacobianDomainerror(snes);
216       }
217       PetscCall(MatLMVMReset(qn->B, PETSC_FALSE));
218     }
219   }
220   if (i == snes->max_its) {
221     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its));
222     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
223   }
224   PetscFunctionReturn(0);
225 }
226 
227 static PetscErrorCode SNESSetUp_QN(SNES snes) {
228   SNES_QN *qn = (SNES_QN *)snes->data;
229   DM       dm;
230   PetscInt n, N;
231 
232   PetscFunctionBegin;
233 
234   if (!snes->vec_sol) {
235     PetscCall(SNESGetDM(snes, &dm));
236     PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol));
237   }
238   PetscCall(SNESSetWorkVecs(snes, 4));
239 
240   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) PetscCall(SNESSetUpMatrices(snes));
241   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
242 
243   /* set method defaults */
244   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
245     if (qn->type == SNES_QN_BADBROYDEN) {
246       qn->scale_type = SNES_QN_SCALE_NONE;
247     } else {
248       qn->scale_type = SNES_QN_SCALE_SCALAR;
249     }
250   }
251   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
252     if (qn->type == SNES_QN_LBFGS) {
253       qn->restart_type = SNES_QN_RESTART_POWELL;
254     } else {
255       qn->restart_type = SNES_QN_RESTART_PERIODIC;
256     }
257   }
258 
259   /* Set up the LMVM matrix */
260   switch (qn->type) {
261   case SNES_QN_BROYDEN:
262     PetscCall(MatSetType(qn->B, MATLMVMBROYDEN));
263     qn->scale_type = SNES_QN_SCALE_NONE;
264     break;
265   case SNES_QN_BADBROYDEN:
266     PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN));
267     qn->scale_type = SNES_QN_SCALE_NONE;
268     break;
269   default:
270     PetscCall(MatSetType(qn->B, MATLMVMBFGS));
271     switch (qn->scale_type) {
272     case SNES_QN_SCALE_NONE: PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE)); break;
273     case SNES_QN_SCALE_SCALAR: PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR)); break;
274     case SNES_QN_SCALE_JACOBIAN: PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_USER)); break;
275     case SNES_QN_SCALE_DIAGONAL:
276     case SNES_QN_SCALE_DEFAULT:
277     default: break;
278     }
279     break;
280   }
281   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
282   PetscCall(VecGetSize(snes->vec_sol, &N));
283   PetscCall(MatSetSizes(qn->B, n, n, N, N));
284   PetscCall(MatSetUp(qn->B));
285   PetscCall(MatLMVMReset(qn->B, PETSC_TRUE));
286   PetscCall(MatLMVMSetHistorySize(qn->B, qn->m));
287   PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func));
288   PetscFunctionReturn(0);
289 }
290 
291 static PetscErrorCode SNESReset_QN(SNES snes) {
292   SNES_QN *qn;
293 
294   PetscFunctionBegin;
295   if (snes->data) {
296     qn = (SNES_QN *)snes->data;
297     PetscCall(MatDestroy(&qn->B));
298   }
299   PetscFunctionReturn(0);
300 }
301 
302 static PetscErrorCode SNESDestroy_QN(SNES snes) {
303   PetscFunctionBegin;
304   PetscCall(SNESReset_QN(snes));
305   PetscCall(PetscFree(snes->data));
306   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", NULL));
307   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", NULL));
308   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", NULL));
309   PetscFunctionReturn(0);
310 }
311 
312 static PetscErrorCode SNESSetFromOptions_QN(SNES snes, PetscOptionItems *PetscOptionsObject) {
313   SNES_QN          *qn = (SNES_QN *)snes->data;
314   PetscBool         flg;
315   SNESLineSearch    linesearch;
316   SNESQNRestartType rtype = qn->restart_type;
317   SNESQNScaleType   stype = qn->scale_type;
318   SNESQNType        qtype = qn->type;
319 
320   PetscFunctionBegin;
321   PetscOptionsHeadBegin(PetscOptionsObject, "SNES QN options");
322   PetscCall(PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, NULL));
323   PetscCall(PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL));
324   PetscCall(PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", qn->monflg, &qn->monflg, NULL));
325   PetscCall(PetscOptionsEnum("-snes_qn_scale_type", "Scaling type", "SNESQNSetScaleType", SNESQNScaleTypes, (PetscEnum)stype, (PetscEnum *)&stype, &flg));
326   if (flg) PetscCall(SNESQNSetScaleType(snes, stype));
327 
328   PetscCall(PetscOptionsEnum("-snes_qn_restart_type", "Restart type", "SNESQNSetRestartType", SNESQNRestartTypes, (PetscEnum)rtype, (PetscEnum *)&rtype, &flg));
329   if (flg) PetscCall(SNESQNSetRestartType(snes, rtype));
330 
331   PetscCall(PetscOptionsEnum("-snes_qn_type", "Quasi-Newton update type", "", SNESQNTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, &flg));
332   if (flg) PetscCall(SNESQNSetType(snes, qtype));
333   PetscCall(MatSetFromOptions(qn->B));
334   PetscOptionsHeadEnd();
335   if (!snes->linesearch) {
336     PetscCall(SNESGetLineSearch(snes, &linesearch));
337     if (!((PetscObject)linesearch)->type_name) {
338       if (qn->type == SNES_QN_LBFGS) {
339         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
340       } else if (qn->type == SNES_QN_BROYDEN) {
341         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
342       } else {
343         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
344       }
345     }
346   }
347   if (qn->monflg) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor));
348   PetscFunctionReturn(0);
349 }
350 
351 static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer) {
352   SNES_QN  *qn = (SNES_QN *)snes->data;
353   PetscBool iascii;
354 
355   PetscFunctionBegin;
356   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
357   if (iascii) {
358     PetscCall(PetscViewerASCIIPrintf(viewer, "  type is %s, restart type is %s, scale type is %s\n", SNESQNTypes[qn->type], SNESQNRestartTypes[qn->restart_type], SNESQNScaleTypes[qn->scale_type]));
359     PetscCall(PetscViewerASCIIPrintf(viewer, "  Stored subspace size: %" PetscInt_FMT "\n", qn->m));
360   }
361   PetscFunctionReturn(0);
362 }
363 
364 /*@
365     SNESQNSetRestartType - Sets the restart type for `SNESQN`.
366 
367     Logically Collective on snes
368 
369     Input Parameters:
370 +   snes - the iterative context
371 -   rtype - restart type
372 
373     Options Database Keys:
374 +   -snes_qn_restart_type <powell,periodic,none> - set the restart type
375 -   -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
376 
377     Level: intermediate
378 
379     `SNESQNRestartType`s:
380 +   `SNES_QN_RESTART_NONE` - never restart
381 .   `SNES_QN_RESTART_POWELL` - restart based upon descent criteria
382 -   `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations
383 
384 .seealso: `SNESQN`, `SNESQNRestartType`, `SNES_QN_RESTART_NONE`, `SNES_QN_RESTART_POWELL`, `SNES_QN_RESTART_PERIODIC`
385 @*/
386 PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) {
387   PetscFunctionBegin;
388   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
389   PetscTryMethod(snes, "SNESQNSetRestartType_C", (SNES, SNESQNRestartType), (snes, rtype));
390   PetscFunctionReturn(0);
391 }
392 
393 /*@
394     SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in `SNESQN`.
395 
396     Logically Collective on snes
397 
398     Input Parameters:
399 +   snes - the nonlinear solver context
400 -   stype - scale type
401 
402     Options Database:
403 .   -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
404 
405     Level: intermediate
406 
407     `SNESQNScaleType`s:
408 +   `SNES_QN_SCALE_NONE` - don't scale the problem
409 .   `SNES_QN_SCALE_SCALAR` - use Shanno scaling
410 .   `SNES_QN_SCALE_DIAGONAL` - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available
411 -   `SNES_QN_SCALE_JACOBIAN` - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration
412                              of QN and at ever restart.
413 
414 .seealso: `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()`
415 @*/
416 
417 PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) {
418   PetscFunctionBegin;
419   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
420   PetscTryMethod(snes, "SNESQNSetScaleType_C", (SNES, SNESQNScaleType), (snes, stype));
421   PetscFunctionReturn(0);
422 }
423 
424 PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) {
425   SNES_QN *qn = (SNES_QN *)snes->data;
426 
427   PetscFunctionBegin;
428   qn->scale_type = stype;
429   if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE;
430   PetscFunctionReturn(0);
431 }
432 
433 PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) {
434   SNES_QN *qn = (SNES_QN *)snes->data;
435 
436   PetscFunctionBegin;
437   qn->restart_type = rtype;
438   PetscFunctionReturn(0);
439 }
440 
441 /*@
442     SNESQNSetType - Sets the quasi-Newton variant to be used in `SNESQN`.
443 
444     Logically Collective on snes
445 
446     Input Parameters:
447 +   snes - the iterative context
448 -   qtype - variant type
449 
450     Options Database:
451 .   -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
452 
453     Level: beginner
454 
455     `SNESQNType`s:
456 +   `SNES_QN_LBFGS` - LBFGS variant
457 .   `SNES_QN_BROYDEN` - Broyden variant
458 -   `SNES_QN_BADBROYDEN` - Bad Broyden variant
459 
460 .seealso: `SNESQN`, `SNES_QN_LBFGS`, `SNES_QN_BROYDEN`, `SNES_QN_BADBROYDEN`, `SNESQNType`, `TAOLMVM`, `TAOBLMVM`
461 @*/
462 
463 PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype) {
464   PetscFunctionBegin;
465   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
466   PetscTryMethod(snes, "SNESQNSetType_C", (SNES, SNESQNType), (snes, qtype));
467   PetscFunctionReturn(0);
468 }
469 
470 PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype) {
471   SNES_QN *qn = (SNES_QN *)snes->data;
472 
473   PetscFunctionBegin;
474   qn->type = qtype;
475   PetscFunctionReturn(0);
476 }
477 
478 /*MC
479       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
480 
481       Options Database Keys:
482 +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
483 .     -snes_qn_restart_type <powell,periodic,none> - set the restart type
484 .     -snes_qn_powell_gamma - Angle condition for restart.
485 .     -snes_qn_powell_descent - Descent condition for restart.
486 .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
487 .     -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian
488 .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
489 -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.
490 
491       References:
492 +   * -   Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
493 .   * -   R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods,
494       Technical Report, Northwestern University, June 1992.
495 .   * -   Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
496       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
497 .   * -   Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
498        SIAM Review, 57(4), 2015
499 .   * -   Griewank, Andreas. "Broyden updating, the good and the bad!." Doc. Math (2012): 301-315.
500 .   * -   Gilbert, Jean Charles, and Claude Lemar{\'e}chal. "Some numerical experiments with variable-storage quasi-Newton algorithms."
501       Mathematical programming 45.1-3 (1989): 407-435.
502 -   * -   Dener A., Munson T. "Accelerating Limited-Memory Quasi-Newton Convergence for Large-Scale Optimization"
503       Computational Science - ICCS 2019. ICCS 2019. Lecture Notes in Computer Science, vol 11538. Springer, Cham
504 
505       Level: beginner
506 
507       Notes:
508       This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
509       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
510       updates.
511 
512       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
513       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
514       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
515       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
516 
517       Uses left nonlinear preconditioning by default.
518 
519 .seealso: `SNESQNRestartType`, `SNESQNSetRestartType()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`,
520           `SNESQNScaleType`, `SNESQNSetScaleType()`, `SNESQNSetType`, `SNESQNSetType ()`
521 M*/
522 PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) {
523   SNES_QN    *qn;
524   const char *optionsprefix;
525 
526   PetscFunctionBegin;
527   snes->ops->setup          = SNESSetUp_QN;
528   snes->ops->solve          = SNESSolve_QN;
529   snes->ops->destroy        = SNESDestroy_QN;
530   snes->ops->setfromoptions = SNESSetFromOptions_QN;
531   snes->ops->view           = SNESView_QN;
532   snes->ops->reset          = SNESReset_QN;
533 
534   snes->npcside = PC_LEFT;
535 
536   snes->usesnpc = PETSC_TRUE;
537   snes->usesksp = PETSC_FALSE;
538 
539   snes->alwayscomputesfinalresidual = PETSC_TRUE;
540 
541   if (!snes->tolerancesset) {
542     snes->max_funcs = 30000;
543     snes->max_its   = 10000;
544   }
545 
546   PetscCall(PetscNewLog(snes, &qn));
547   snes->data       = (void *)qn;
548   qn->m            = 10;
549   qn->scaling      = 1.0;
550   qn->monitor      = NULL;
551   qn->monflg       = PETSC_FALSE;
552   qn->powell_gamma = 0.9999;
553   qn->scale_type   = SNES_QN_SCALE_DEFAULT;
554   qn->restart_type = SNES_QN_RESTART_DEFAULT;
555   qn->type         = SNES_QN_LBFGS;
556 
557   PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B));
558   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
559   PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix));
560 
561   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", SNESQNSetScaleType_QN));
562   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", SNESQNSetRestartType_QN));
563   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", SNESQNSetType_QN));
564   PetscFunctionReturn(0);
565 }
566