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