xref: /petsc/src/snes/impls/qn/qn.c (revision 574a8f54e195df5ade05daee45527e64c2280ff7)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2fced5a79SAsbjørn Nilsen Riseth #include <petscdm.h>
34b11644fSPeter Brune 
48e231d97SPeter Brune #define H(i, j) qn->dXdFmat[i * qn->m + j]
58e231d97SPeter Brune 
69e5d0892SLisandro Dalcin const char *const SNESQNScaleTypes[]   = {"DEFAULT", "NONE", "SCALAR", "DIAGONAL", "JACOBIAN", "SNESQNScaleType", "SNES_QN_SCALING_", NULL};
79e5d0892SLisandro Dalcin const char *const SNESQNRestartTypes[] = {"DEFAULT", "NONE", "POWELL", "PERIODIC", "SNESQNRestartType", "SNES_QN_RESTART_", NULL};
89e5d0892SLisandro Dalcin const char *const SNESQNTypes[]        = {"LBFGS", "BROYDEN", "BADBROYDEN", "SNESQNType", "SNES_QN_", NULL};
9b8840d6bSPeter Brune 
104b11644fSPeter Brune typedef struct {
1192f76d53SAlp Dener   Mat               B;      /* Quasi-Newton approximation Matrix (MATLMVM) */
128e231d97SPeter Brune   PetscInt          m;      /* The number of kept previous steps */
135c7a0a03SPeter Brune   PetscReal        *lambda; /* The line search history of the method */
1492f76d53SAlp Dener   PetscBool         monflg;
1544f7e39eSPeter Brune   PetscViewer       monitor;
166bf1b2e5SPeter Brune   PetscReal         powell_gamma; /* Powell angle restart condition */
17b21d5a53SPeter Brune   PetscReal         scaling;      /* scaling of H0 */
18b8840d6bSPeter Brune   SNESQNType        type;         /* the type of quasi-newton method used */
1988f769c5SPeter Brune   SNESQNScaleType   scale_type;   /* the type of scaling used */
200c777b0cSPeter Brune   SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */
219f83bee8SJed Brown } SNES_QN;
224b11644fSPeter Brune 
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_QN(SNES snes)
24d71ae5a4SJacob Faibussowitsch {
259f83bee8SJed Brown   SNES_QN             *qn = (SNES_QN *)snes->data;
2615f5eeeaSPeter Brune   Vec                  X, Xold;
270a3905e1SPeter Brune   Vec                  F, W;
2847f26062SPeter Brune   Vec                  Y, D, Dold;
29b8840d6bSPeter Brune   PetscInt             i, i_r;
30335fb975SPeter Brune   PetscReal            fnorm, xnorm, ynorm, gnorm;
31422a814eSBarry Smith   SNESLineSearchReason lssucceed;
3292f76d53SAlp Dener   PetscBool            badstep, powell, periodic;
331c6523dcSPeter Brune   PetscScalar          DolddotD, DolddotDold;
34b7281c8aSPeter Brune   SNESConvergedReason  reason;
350ecc9e1dSPeter Brune 
3684c577daSBarry Smith   /* basically just a regular newton's method except for the application of the Jacobian */
374b11644fSPeter Brune 
386e111a19SKarl Rupp   PetscFunctionBegin;
390b121fc5SBarry Smith   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);
40c579b300SPatrick Farrell 
419566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
429f3a0142SPeter Brune   F    = snes->vec_func;       /* residual vector */
433af51624SPeter Brune   Y    = snes->vec_sol_update; /* search direction generated by J^-1D*/
440a3905e1SPeter Brune   W    = snes->work[3];
45b8840d6bSPeter Brune   X    = snes->vec_sol; /* solution vector */
46335fb975SPeter Brune   Xold = snes->work[0];
473af51624SPeter Brune 
483af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
49335fb975SPeter Brune   D    = snes->work[1];
50335fb975SPeter Brune   Dold = snes->work[2];
514b11644fSPeter Brune 
524b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
534b11644fSPeter Brune 
549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
554b11644fSPeter Brune   snes->iter = 0;
564b11644fSPeter Brune   snes->norm = 0.;
579566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
5847f26062SPeter Brune 
59efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
609566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
619566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
62b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
63b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
64b7281c8aSPeter Brune       PetscFunctionReturn(0);
65b7281c8aSPeter Brune     }
669566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
6747f26062SPeter Brune   } else {
68e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
699566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
701aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
71c1c75074SPeter Brune 
729566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
73422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
7447f26062SPeter Brune   }
75efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
769566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, F, D));
779566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
78b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
79b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
80b7281c8aSPeter Brune       PetscFunctionReturn(0);
81b7281c8aSPeter Brune     }
8247f26062SPeter Brune   } else {
839566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, D));
8447f26062SPeter Brune   }
85b28a06ddSPeter Brune 
869566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
874b11644fSPeter Brune   snes->norm = fnorm;
889566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
899566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
909566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
914b11644fSPeter Brune 
924b11644fSPeter Brune   /* test convergence */
93dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
944b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
9570d3b23bSPeter Brune 
96efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_RIGHT) {
979566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0));
989566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
999566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0));
1009566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
101ddd40ce5SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
102ddd40ce5SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
103ddd40ce5SPeter Brune       PetscFunctionReturn(0);
104ddd40ce5SPeter Brune     }
1059566063dSJacob Faibussowitsch     PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
1069566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, D));
1073cf07b75SPeter Brune   }
1083cf07b75SPeter Brune 
10901fe78eaSStefano Zampini   /* general purpose update */
110dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
11101fe78eaSStefano Zampini 
112f8e15203SPeter Brune   /* scale the initial update */
1130c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
1149566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
11507b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
1169566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
1179566063dSJacob Faibussowitsch     PetscCall(MatLMVMSetJ0KSP(qn->B, snes->ksp));
1180ecc9e1dSPeter Brune   }
1190ecc9e1dSPeter Brune 
1205ba6227bSPeter Brune   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
12192f76d53SAlp Dener     /* update QN approx and calculate step */
1229566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(qn->B, X, D));
1239566063dSJacob Faibussowitsch     PetscCall(MatSolve(qn->B, D, Y));
12492f76d53SAlp Dener 
12570d3b23bSPeter Brune     /* line search for lambda */
1269371c9d4SSatish Balay     ynorm = 1;
1279371c9d4SSatish Balay     gnorm = fnorm;
1289566063dSJacob Faibussowitsch     PetscCall(VecCopy(D, Dold));
1299566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xold));
1309566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
1319f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1329566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
1339566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
13492f76d53SAlp Dener     badstep = PETSC_FALSE;
135422a814eSBarry Smith     if (lssucceed) {
1369f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
1379f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
1389f3a0142SPeter Brune         break;
1399f3a0142SPeter Brune       }
14092f76d53SAlp Dener       badstep = PETSC_TRUE;
1410ecc9e1dSPeter Brune     }
1423af51624SPeter Brune 
14388d374b2SPeter Brune     /* convergence monitoring */
1449566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed));
145b21d5a53SPeter Brune 
146efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_RIGHT) {
1479566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0));
1489566063dSJacob Faibussowitsch       PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1499566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0));
1509566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
151ddd40ce5SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
152ddd40ce5SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
153ddd40ce5SPeter Brune         PetscFunctionReturn(0);
154ddd40ce5SPeter Brune       }
1559566063dSJacob Faibussowitsch       PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
156b28a06ddSPeter Brune     }
157b28a06ddSPeter Brune 
1589566063dSJacob Faibussowitsch     PetscCall(SNESSetIterationNumber(snes, i + 1));
15971dbe336SPeter Brune     snes->norm  = fnorm;
160c1e67a49SFande Kong     snes->xnorm = xnorm;
161c1e67a49SFande Kong     snes->ynorm = ynorm;
162360c497dSPeter Brune 
1639566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
1649566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
16592f76d53SAlp Dener 
1664b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
167dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
1684b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
169efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
1709566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, X, F, D));
1719566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
172b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
173b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
174b7281c8aSPeter Brune         PetscFunctionReturn(0);
175b7281c8aSPeter Brune       }
17688d374b2SPeter Brune     } else {
1779566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, D));
17888d374b2SPeter Brune     }
17901fe78eaSStefano Zampini 
18001fe78eaSStefano Zampini     /* general purpose update */
181dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
18201fe78eaSStefano Zampini 
18392f76d53SAlp Dener     /* restart conditions */
1840c777b0cSPeter Brune     powell = PETSC_FALSE;
1856bdcc836SBarry Smith     if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) {
1866bdcc836SBarry Smith       /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */
18723c918e7SPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
188*574a8f54SStefano Zampini         PetscCall(MatMult(snes->jacobian, Dold, W));
18923c918e7SPeter Brune       } else {
1909566063dSJacob Faibussowitsch         PetscCall(VecCopy(Dold, W));
19123c918e7SPeter Brune       }
1929566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(W, Dold, &DolddotDold));
1939566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(W, D, &DolddotD));
1949566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(W, Dold, &DolddotDold));
1959566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(W, D, &DolddotD));
1960c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma * PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
1970c777b0cSPeter Brune     }
1980c777b0cSPeter Brune     periodic = PETSC_FALSE;
199b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
200b8840d6bSPeter Brune       if (i_r > qn->m - 1) periodic = PETSC_TRUE;
2010c777b0cSPeter Brune     }
2020c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
20392f76d53SAlp Dener     if (badstep || powell || periodic) {
20492f76d53SAlp Dener       if (qn->monflg) {
2059566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
2066bdcc836SBarry Smith         if (powell) {
20763a3b9bcSJacob Faibussowitsch           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));
2086bdcc836SBarry Smith         } else {
20963a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %" PetscInt_FMT "\n", i_r));
2106bdcc836SBarry Smith         }
2119566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
2125ba6227bSPeter Brune       }
2135ba6227bSPeter Brune       i_r = -1;
2140c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
2159566063dSJacob Faibussowitsch         PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
21607b62357SFande Kong         SNESCheckJacobianDomainerror(snes);
2170ecc9e1dSPeter Brune       }
2189566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(qn->B, PETSC_FALSE));
2190ecc9e1dSPeter Brune     }
2205ba6227bSPeter Brune   }
2214b11644fSPeter Brune   if (i == snes->max_its) {
22263a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its));
2234b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
2244b11644fSPeter Brune   }
2254b11644fSPeter Brune   PetscFunctionReturn(0);
2264b11644fSPeter Brune }
2274b11644fSPeter Brune 
228d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_QN(SNES snes)
229d71ae5a4SJacob Faibussowitsch {
2309f83bee8SJed Brown   SNES_QN *qn = (SNES_QN *)snes->data;
231fced5a79SAsbjørn Nilsen Riseth   DM       dm;
23292f76d53SAlp Dener   PetscInt n, N;
233335fb975SPeter Brune 
2344b11644fSPeter Brune   PetscFunctionBegin;
235fced5a79SAsbjørn Nilsen Riseth 
236fced5a79SAsbjørn Nilsen Riseth   if (!snes->vec_sol) {
2379566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dm));
2389566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol));
239fced5a79SAsbjørn Nilsen Riseth   }
2409566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 4));
24192f76d53SAlp Dener 
24248a46eb9SPierre Jolivet   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) PetscCall(SNESSetUpMatrices(snes));
243ad540459SPierre Jolivet   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
24492f76d53SAlp Dener 
24560c08b40SPeter Brune   /* set method defaults */
2461efc8c45SPeter Brune   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
24760c08b40SPeter Brune     if (qn->type == SNES_QN_BADBROYDEN) {
24860c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_NONE;
24960c08b40SPeter Brune     } else {
25092f76d53SAlp Dener       qn->scale_type = SNES_QN_SCALE_SCALAR;
25160c08b40SPeter Brune     }
25260c08b40SPeter Brune   }
2531efc8c45SPeter Brune   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
25460c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
25560c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_POWELL;
25660c08b40SPeter Brune     } else {
25760c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_PERIODIC;
25860c08b40SPeter Brune     }
25960c08b40SPeter Brune   }
26060c08b40SPeter Brune 
26192f76d53SAlp Dener   /* Set up the LMVM matrix */
26292f76d53SAlp Dener   switch (qn->type) {
26392f76d53SAlp Dener   case SNES_QN_BROYDEN:
2649566063dSJacob Faibussowitsch     PetscCall(MatSetType(qn->B, MATLMVMBROYDEN));
26592f76d53SAlp Dener     qn->scale_type = SNES_QN_SCALE_NONE;
26692f76d53SAlp Dener     break;
26792f76d53SAlp Dener   case SNES_QN_BADBROYDEN:
2689566063dSJacob Faibussowitsch     PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN));
26992f76d53SAlp Dener     qn->scale_type = SNES_QN_SCALE_NONE;
27092f76d53SAlp Dener     break;
27192f76d53SAlp Dener   default:
2729566063dSJacob Faibussowitsch     PetscCall(MatSetType(qn->B, MATLMVMBFGS));
27392f76d53SAlp Dener     switch (qn->scale_type) {
274d71ae5a4SJacob Faibussowitsch     case SNES_QN_SCALE_NONE:
275d71ae5a4SJacob Faibussowitsch       PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE));
276d71ae5a4SJacob Faibussowitsch       break;
277d71ae5a4SJacob Faibussowitsch     case SNES_QN_SCALE_SCALAR:
278d71ae5a4SJacob Faibussowitsch       PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR));
279d71ae5a4SJacob Faibussowitsch       break;
280d71ae5a4SJacob Faibussowitsch     case SNES_QN_SCALE_JACOBIAN:
281d71ae5a4SJacob Faibussowitsch       PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_USER));
282d71ae5a4SJacob Faibussowitsch       break;
28392f76d53SAlp Dener     case SNES_QN_SCALE_DIAGONAL:
28492f76d53SAlp Dener     case SNES_QN_SCALE_DEFAULT:
285d71ae5a4SJacob Faibussowitsch     default:
286d71ae5a4SJacob Faibussowitsch       break;
2878e231d97SPeter Brune     }
28892f76d53SAlp Dener     break;
28992f76d53SAlp Dener   }
2909566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
2919566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol, &N));
2929566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(qn->B, n, n, N, N));
2939566063dSJacob Faibussowitsch   PetscCall(MatSetUp(qn->B));
2949566063dSJacob Faibussowitsch   PetscCall(MatLMVMReset(qn->B, PETSC_TRUE));
2959566063dSJacob Faibussowitsch   PetscCall(MatLMVMSetHistorySize(qn->B, qn->m));
2969566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func));
2974b11644fSPeter Brune   PetscFunctionReturn(0);
2984b11644fSPeter Brune }
2994b11644fSPeter Brune 
300d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_QN(SNES snes)
301d71ae5a4SJacob Faibussowitsch {
3029f83bee8SJed Brown   SNES_QN *qn;
3030adebc6cSBarry Smith 
3044b11644fSPeter Brune   PetscFunctionBegin;
3054b11644fSPeter Brune   if (snes->data) {
3069f83bee8SJed Brown     qn = (SNES_QN *)snes->data;
3079566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&qn->B));
3084b11644fSPeter Brune   }
3094b11644fSPeter Brune   PetscFunctionReturn(0);
3104b11644fSPeter Brune }
3114b11644fSPeter Brune 
312d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_QN(SNES snes)
313d71ae5a4SJacob Faibussowitsch {
3144b11644fSPeter Brune   PetscFunctionBegin;
3159566063dSJacob Faibussowitsch   PetscCall(SNESReset_QN(snes));
3169566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3172e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", NULL));
3182e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", NULL));
3192e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", NULL));
3204b11644fSPeter Brune   PetscFunctionReturn(0);
3214b11644fSPeter Brune }
3224b11644fSPeter Brune 
323d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_QN(SNES snes, PetscOptionItems *PetscOptionsObject)
324d71ae5a4SJacob Faibussowitsch {
3252150357eSBarry Smith   SNES_QN          *qn = (SNES_QN *)snes->data;
32692f76d53SAlp Dener   PetscBool         flg;
327f1c6b773SPeter Brune   SNESLineSearch    linesearch;
3282150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
3292150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
3301efc8c45SPeter Brune   SNESQNType        qtype = qn->type;
3312150357eSBarry Smith 
3324b11644fSPeter Brune   PetscFunctionBegin;
333d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES QN options");
3349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, NULL));
3359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL));
3369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", qn->monflg, &qn->monflg, NULL));
3379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_scale_type", "Scaling type", "SNESQNSetScaleType", SNESQNScaleTypes, (PetscEnum)stype, (PetscEnum *)&stype, &flg));
3389566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetScaleType(snes, stype));
33988f769c5SPeter Brune 
3409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_restart_type", "Restart type", "SNESQNSetRestartType", SNESQNRestartTypes, (PetscEnum)rtype, (PetscEnum *)&rtype, &flg));
3419566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetRestartType(snes, rtype));
34288f769c5SPeter Brune 
3439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_type", "Quasi-Newton update type", "", SNESQNTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, &flg));
3449566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetType(snes, qtype));
3459566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(qn->B));
346d0609cedSBarry Smith   PetscOptionsHeadEnd();
3479e764e56SPeter Brune   if (!snes->linesearch) {
3489566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
349ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
35060c08b40SPeter Brune       if (qn->type == SNES_QN_LBFGS) {
3519566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
35260c08b40SPeter Brune       } else if (qn->type == SNES_QN_BROYDEN) {
3539566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
35460c08b40SPeter Brune       } else {
3559566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
35660c08b40SPeter Brune       }
3579e764e56SPeter Brune     }
358ec786807SJed Brown   }
35948a46eb9SPierre Jolivet   if (qn->monflg) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor));
3604b11644fSPeter Brune   PetscFunctionReturn(0);
3614b11644fSPeter Brune }
3624b11644fSPeter Brune 
363d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
364d71ae5a4SJacob Faibussowitsch {
3655cd83419SPeter Brune   SNES_QN  *qn = (SNES_QN *)snes->data;
3665cd83419SPeter Brune   PetscBool iascii;
3675cd83419SPeter Brune 
3685cd83419SPeter Brune   PetscFunctionBegin;
3699566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3705cd83419SPeter Brune   if (iascii) {
3719566063dSJacob Faibussowitsch     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]));
37263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Stored subspace size: %" PetscInt_FMT "\n", qn->m));
3735cd83419SPeter Brune   }
3745cd83419SPeter Brune   PetscFunctionReturn(0);
3755cd83419SPeter Brune }
37692c02d66SPeter Brune 
3770c777b0cSPeter Brune /*@
378f6dfbefdSBarry Smith     SNESQNSetRestartType - Sets the restart type for `SNESQN`.
3790c777b0cSPeter Brune 
380f6dfbefdSBarry Smith     Logically Collective on snes
3810c777b0cSPeter Brune 
3820c777b0cSPeter Brune     Input Parameters:
3830c777b0cSPeter Brune +   snes - the iterative context
3840c777b0cSPeter Brune -   rtype - restart type
3850c777b0cSPeter Brune 
386f6dfbefdSBarry Smith     Options Database Keys:
3870c777b0cSPeter Brune +   -snes_qn_restart_type <powell,periodic,none> - set the restart type
38884c577daSBarry Smith -   -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
3890c777b0cSPeter Brune 
3900c777b0cSPeter Brune     Level: intermediate
3910c777b0cSPeter Brune 
392f6dfbefdSBarry Smith     `SNESQNRestartType`s:
393f6dfbefdSBarry Smith +   `SNES_QN_RESTART_NONE` - never restart
394f6dfbefdSBarry Smith .   `SNES_QN_RESTART_POWELL` - restart based upon descent criteria
395f6dfbefdSBarry Smith -   `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations
3960c777b0cSPeter Brune 
397f6dfbefdSBarry Smith .seealso: `SNESQN`, `SNESQNRestartType`, `SNES_QN_RESTART_NONE`, `SNES_QN_RESTART_POWELL`, `SNES_QN_RESTART_PERIODIC`
3980c777b0cSPeter Brune @*/
399d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
400d71ae5a4SJacob Faibussowitsch {
4010c777b0cSPeter Brune   PetscFunctionBegin;
4020c777b0cSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
403cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetRestartType_C", (SNES, SNESQNRestartType), (snes, rtype));
4040c777b0cSPeter Brune   PetscFunctionReturn(0);
4050c777b0cSPeter Brune }
4060c777b0cSPeter Brune 
4070c777b0cSPeter Brune /*@
408f6dfbefdSBarry Smith     SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in `SNESQN`.
4090c777b0cSPeter Brune 
410f6dfbefdSBarry Smith     Logically Collective on snes
4110c777b0cSPeter Brune 
4120c777b0cSPeter Brune     Input Parameters:
413f6dfbefdSBarry Smith +   snes - the nonlinear solver context
4140c777b0cSPeter Brune -   stype - scale type
4150c777b0cSPeter Brune 
4163c7db156SBarry Smith     Options Database Key:
41767b8a455SSatish Balay .   -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
4180c777b0cSPeter Brune 
4190c777b0cSPeter Brune     Level: intermediate
4200c777b0cSPeter Brune 
421f6dfbefdSBarry Smith     `SNESQNScaleType`s:
422f6dfbefdSBarry Smith +   `SNES_QN_SCALE_NONE` - don't scale the problem
423f6dfbefdSBarry Smith .   `SNES_QN_SCALE_SCALAR` - use Shanno scaling
424f6dfbefdSBarry Smith .   `SNES_QN_SCALE_DIAGONAL` - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available
425f6dfbefdSBarry Smith -   `SNES_QN_SCALE_JACOBIAN` - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration
426a01a0525SBarry Smith                              of QN and at ever restart.
4270c777b0cSPeter Brune 
428db781477SPatrick Sanan .seealso: `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()`
4290c777b0cSPeter Brune @*/
4300c777b0cSPeter Brune 
431d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
432d71ae5a4SJacob Faibussowitsch {
4330c777b0cSPeter Brune   PetscFunctionBegin;
4340c777b0cSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
435cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetScaleType_C", (SNES, SNESQNScaleType), (snes, stype));
4360c777b0cSPeter Brune   PetscFunctionReturn(0);
4370c777b0cSPeter Brune }
4380c777b0cSPeter Brune 
439d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
440d71ae5a4SJacob Faibussowitsch {
4410c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4426e111a19SKarl Rupp 
4430c777b0cSPeter Brune   PetscFunctionBegin;
4440c777b0cSPeter Brune   qn->scale_type = stype;
4457044b745SBarry Smith   if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE;
4460c777b0cSPeter Brune   PetscFunctionReturn(0);
4470c777b0cSPeter Brune }
4480c777b0cSPeter Brune 
449d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
450d71ae5a4SJacob Faibussowitsch {
4510c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4526e111a19SKarl Rupp 
4530c777b0cSPeter Brune   PetscFunctionBegin;
4540c777b0cSPeter Brune   qn->restart_type = rtype;
4550c777b0cSPeter Brune   PetscFunctionReturn(0);
4560c777b0cSPeter Brune }
4570c777b0cSPeter Brune 
4581efc8c45SPeter Brune /*@
459f6dfbefdSBarry Smith     SNESQNSetType - Sets the quasi-Newton variant to be used in `SNESQN`.
4601efc8c45SPeter Brune 
461f6dfbefdSBarry Smith     Logically Collective on snes
4621efc8c45SPeter Brune 
4631efc8c45SPeter Brune     Input Parameters:
4641efc8c45SPeter Brune +   snes - the iterative context
4651efc8c45SPeter Brune -   qtype - variant type
4661efc8c45SPeter Brune 
4673c7db156SBarry Smith     Options Database Key:
46867b8a455SSatish Balay .   -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
4691efc8c45SPeter Brune 
4701efc8c45SPeter Brune     Level: beginner
4711efc8c45SPeter Brune 
472f6dfbefdSBarry Smith     `SNESQNType`s:
473f6dfbefdSBarry Smith +   `SNES_QN_LBFGS` - LBFGS variant
474f6dfbefdSBarry Smith .   `SNES_QN_BROYDEN` - Broyden variant
475f6dfbefdSBarry Smith -   `SNES_QN_BADBROYDEN` - Bad Broyden variant
4761efc8c45SPeter Brune 
477f6dfbefdSBarry Smith .seealso: `SNESQN`, `SNES_QN_LBFGS`, `SNES_QN_BROYDEN`, `SNES_QN_BADBROYDEN`, `SNESQNType`, `TAOLMVM`, `TAOBLMVM`
4781efc8c45SPeter Brune @*/
4791efc8c45SPeter Brune 
480d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
481d71ae5a4SJacob Faibussowitsch {
4821efc8c45SPeter Brune   PetscFunctionBegin;
4831efc8c45SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
484cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetType_C", (SNES, SNESQNType), (snes, qtype));
4851efc8c45SPeter Brune   PetscFunctionReturn(0);
4861efc8c45SPeter Brune }
4871efc8c45SPeter Brune 
488d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
489d71ae5a4SJacob Faibussowitsch {
4901efc8c45SPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4911efc8c45SPeter Brune 
4921efc8c45SPeter Brune   PetscFunctionBegin;
4931efc8c45SPeter Brune   qn->type = qtype;
4941efc8c45SPeter Brune   PetscFunctionReturn(0);
4951efc8c45SPeter Brune }
4961efc8c45SPeter Brune 
4974b11644fSPeter Brune /*MC
4984b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
4994b11644fSPeter Brune 
500f6dfbefdSBarry Smith       Options Database Keys:
50184c577daSBarry Smith +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
502f6dfbefdSBarry Smith .     -snes_qn_restart_type <powell,periodic,none> - set the restart type
5036bdcc836SBarry Smith .     -snes_qn_powell_gamma - Angle condition for restart.
5041867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
50584c577daSBarry Smith .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
50692f76d53SAlp Dener .     -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian
50772835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
50884c577daSBarry Smith -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.
5096cc8130cSPeter Brune 
5106cc8130cSPeter Brune       References:
511606c0280SSatish Balay +   * -   Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
512606c0280SSatish Balay .   * -   R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods,
5130a8e8854SPeter Brune       Technical Report, Northwestern University, June 1992.
514606c0280SSatish Balay .   * -   Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
5150a8e8854SPeter Brune       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
516606c0280SSatish Balay .   * -   Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
5174f02bc6aSBarry Smith        SIAM Review, 57(4), 2015
518606c0280SSatish Balay .   * -   Griewank, Andreas. "Broyden updating, the good and the bad!." Doc. Math (2012): 301-315.
519606c0280SSatish Balay .   * -   Gilbert, Jean Charles, and Claude Lemar{\'e}chal. "Some numerical experiments with variable-storage quasi-Newton algorithms."
52092f76d53SAlp Dener       Mathematical programming 45.1-3 (1989): 407-435.
521606c0280SSatish Balay -   * -   Dener A., Munson T. "Accelerating Limited-Memory Quasi-Newton Convergence for Large-Scale Optimization"
522110fc3b0SBarry Smith       Computational Science - ICCS 2019. ICCS 2019. Lecture Notes in Computer Science, vol 11538. Springer, Cham
5234b11644fSPeter Brune 
5244b11644fSPeter Brune       Level: beginner
5254b11644fSPeter Brune 
526f6dfbefdSBarry Smith       Notes:
527f6dfbefdSBarry Smith       This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
528f6dfbefdSBarry Smith       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
529f6dfbefdSBarry Smith       updates.
5306cc8130cSPeter Brune 
531f6dfbefdSBarry Smith       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
532f6dfbefdSBarry Smith       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
533f6dfbefdSBarry Smith       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
534f6dfbefdSBarry Smith       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
535f6dfbefdSBarry Smith 
536f6dfbefdSBarry Smith       Uses left nonlinear preconditioning by default.
537f6dfbefdSBarry Smith 
538f6dfbefdSBarry Smith .seealso: `SNESQNRestartType`, `SNESQNSetRestartType()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`,
539f6dfbefdSBarry Smith           `SNESQNScaleType`, `SNESQNSetScaleType()`, `SNESQNSetType`, `SNESQNSetType ()`
5404b11644fSPeter Brune M*/
541d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
542d71ae5a4SJacob Faibussowitsch {
5439f83bee8SJed Brown   SNES_QN    *qn;
54492f76d53SAlp Dener   const char *optionsprefix;
5454b11644fSPeter Brune 
5464b11644fSPeter Brune   PetscFunctionBegin;
5474b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
5484b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
5494b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
5504b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
5515cd83419SPeter Brune   snes->ops->view           = SNESView_QN;
5524b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
5534b11644fSPeter Brune 
554efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
55547f26062SPeter Brune 
556efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
55742f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
55842f4f86dSBarry Smith 
5594fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
5604fc747eaSLawrence Mitchell 
56188976e71SPeter Brune   if (!snes->tolerancesset) {
5620e444f03SPeter Brune     snes->max_funcs = 30000;
5630e444f03SPeter Brune     snes->max_its   = 10000;
56488976e71SPeter Brune   }
5650e444f03SPeter Brune 
5664dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&qn));
5674b11644fSPeter Brune   snes->data       = (void *)qn;
5680ecc9e1dSPeter Brune   qn->m            = 10;
569b21d5a53SPeter Brune   qn->scaling      = 1.0;
5700298fd71SBarry Smith   qn->monitor      = NULL;
57192f76d53SAlp Dener   qn->monflg       = PETSC_FALSE;
572b8840d6bSPeter Brune   qn->powell_gamma = 0.9999;
57360c08b40SPeter Brune   qn->scale_type   = SNES_QN_SCALE_DEFAULT;
57460c08b40SPeter Brune   qn->restart_type = SNES_QN_RESTART_DEFAULT;
575b8840d6bSPeter Brune   qn->type         = SNES_QN_LBFGS;
576ea630c6eSPeter Brune 
5779566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B));
5789566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
5799566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix));
58092f76d53SAlp Dener 
5819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", SNESQNSetScaleType_QN));
5829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", SNESQNSetRestartType_QN));
5839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", SNESQNSetType_QN));
5844b11644fSPeter Brune   PetscFunctionReturn(0);
5854b11644fSPeter Brune }
586