xref: /petsc/src/snes/impls/qn/qn.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 
239371c9d4SSatish Balay static PetscErrorCode SNESSolve_QN(SNES snes) {
249f83bee8SJed Brown   SNES_QN             *qn = (SNES_QN *)snes->data;
2515f5eeeaSPeter Brune   Vec                  X, Xold;
260a3905e1SPeter Brune   Vec                  F, W;
2747f26062SPeter Brune   Vec                  Y, D, Dold;
28b8840d6bSPeter Brune   PetscInt             i, i_r;
29335fb975SPeter Brune   PetscReal            fnorm, xnorm, ynorm, gnorm;
30422a814eSBarry Smith   SNESLineSearchReason lssucceed;
3192f76d53SAlp Dener   PetscBool            badstep, powell, periodic;
321c6523dcSPeter Brune   PetscScalar          DolddotD, DolddotDold;
33b7281c8aSPeter Brune   SNESConvergedReason  reason;
340ecc9e1dSPeter Brune 
3584c577daSBarry Smith   /* basically just a regular newton's method except for the application of the Jacobian */
364b11644fSPeter Brune 
376e111a19SKarl Rupp   PetscFunctionBegin;
380b121fc5SBarry 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);
39c579b300SPatrick Farrell 
409566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
419f3a0142SPeter Brune   F    = snes->vec_func;       /* residual vector */
423af51624SPeter Brune   Y    = snes->vec_sol_update; /* search direction generated by J^-1D*/
430a3905e1SPeter Brune   W    = snes->work[3];
44b8840d6bSPeter Brune   X    = snes->vec_sol; /* solution vector */
45335fb975SPeter Brune   Xold = snes->work[0];
463af51624SPeter Brune 
473af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
48335fb975SPeter Brune   D    = snes->work[1];
49335fb975SPeter Brune   Dold = snes->work[2];
504b11644fSPeter Brune 
514b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
524b11644fSPeter Brune 
539566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
544b11644fSPeter Brune   snes->iter = 0;
554b11644fSPeter Brune   snes->norm = 0.;
569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
5747f26062SPeter Brune 
58efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
599566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
609566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
61b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
62b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
63b7281c8aSPeter Brune       PetscFunctionReturn(0);
64b7281c8aSPeter Brune     }
659566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
6647f26062SPeter Brune   } else {
67e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
689566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
691aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
70c1c75074SPeter Brune 
719566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
72422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
7347f26062SPeter Brune   }
74efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
759566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, F, D));
769566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
77b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
78b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
79b7281c8aSPeter Brune       PetscFunctionReturn(0);
80b7281c8aSPeter Brune     }
8147f26062SPeter Brune   } else {
829566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, D));
8347f26062SPeter Brune   }
84b28a06ddSPeter Brune 
859566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
864b11644fSPeter Brune   snes->norm = fnorm;
879566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
889566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
899566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
904b11644fSPeter Brune 
914b11644fSPeter Brune   /* test convergence */
92dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
934b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
9470d3b23bSPeter Brune 
95efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_RIGHT) {
969566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0));
979566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
989566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0));
999566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
100ddd40ce5SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
101ddd40ce5SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
102ddd40ce5SPeter Brune       PetscFunctionReturn(0);
103ddd40ce5SPeter Brune     }
1049566063dSJacob Faibussowitsch     PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
1059566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, D));
1063cf07b75SPeter Brune   }
1073cf07b75SPeter Brune 
10801fe78eaSStefano Zampini   /* general purpose update */
109dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
11001fe78eaSStefano Zampini 
111f8e15203SPeter Brune   /* scale the initial update */
1120c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
1139566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
11407b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
1159566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
1169566063dSJacob Faibussowitsch     PetscCall(MatLMVMSetJ0KSP(qn->B, snes->ksp));
1170ecc9e1dSPeter Brune   }
1180ecc9e1dSPeter Brune 
1195ba6227bSPeter Brune   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
12092f76d53SAlp Dener     /* update QN approx and calculate step */
1219566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(qn->B, X, D));
1229566063dSJacob Faibussowitsch     PetscCall(MatSolve(qn->B, D, Y));
12392f76d53SAlp Dener 
12470d3b23bSPeter Brune     /* line search for lambda */
1259371c9d4SSatish Balay     ynorm = 1;
1269371c9d4SSatish Balay     gnorm = fnorm;
1279566063dSJacob Faibussowitsch     PetscCall(VecCopy(D, Dold));
1289566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xold));
1299566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
1309f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1319566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
1329566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
13392f76d53SAlp Dener     badstep = PETSC_FALSE;
134422a814eSBarry Smith     if (lssucceed) {
1359f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
1369f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
1379f3a0142SPeter Brune         break;
1389f3a0142SPeter Brune       }
13992f76d53SAlp Dener       badstep = PETSC_TRUE;
1400ecc9e1dSPeter Brune     }
1413af51624SPeter Brune 
14288d374b2SPeter Brune     /* convergence monitoring */
1439566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed));
144b21d5a53SPeter Brune 
145efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_RIGHT) {
1469566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0));
1479566063dSJacob Faibussowitsch       PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1489566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0));
1499566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
150ddd40ce5SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
151ddd40ce5SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
152ddd40ce5SPeter Brune         PetscFunctionReturn(0);
153ddd40ce5SPeter Brune       }
1549566063dSJacob Faibussowitsch       PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
155b28a06ddSPeter Brune     }
156b28a06ddSPeter Brune 
1579566063dSJacob Faibussowitsch     PetscCall(SNESSetIterationNumber(snes, i + 1));
15871dbe336SPeter Brune     snes->norm  = fnorm;
159c1e67a49SFande Kong     snes->xnorm = xnorm;
160c1e67a49SFande Kong     snes->ynorm = ynorm;
161360c497dSPeter Brune 
1629566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
1639566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
16492f76d53SAlp Dener 
1654b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
166dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
1674b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
168efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
1699566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, X, F, D));
1709566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
171b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
172b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
173b7281c8aSPeter Brune         PetscFunctionReturn(0);
174b7281c8aSPeter Brune       }
17588d374b2SPeter Brune     } else {
1769566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, D));
17788d374b2SPeter Brune     }
17801fe78eaSStefano Zampini 
17901fe78eaSStefano Zampini     /* general purpose update */
180dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
18101fe78eaSStefano Zampini 
18292f76d53SAlp Dener     /* restart conditions */
1830c777b0cSPeter Brune     powell = PETSC_FALSE;
1846bdcc836SBarry Smith     if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) {
1856bdcc836SBarry Smith       /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */
18623c918e7SPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
1879566063dSJacob Faibussowitsch         PetscCall(MatMult(snes->jacobian_pre, Dold, W));
18823c918e7SPeter Brune       } else {
1899566063dSJacob Faibussowitsch         PetscCall(VecCopy(Dold, W));
19023c918e7SPeter Brune       }
1919566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(W, Dold, &DolddotDold));
1929566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(W, D, &DolddotD));
1939566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(W, Dold, &DolddotDold));
1949566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(W, D, &DolddotD));
1950c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma * PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
1960c777b0cSPeter Brune     }
1970c777b0cSPeter Brune     periodic = PETSC_FALSE;
198b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
199b8840d6bSPeter Brune       if (i_r > qn->m - 1) periodic = PETSC_TRUE;
2000c777b0cSPeter Brune     }
2010c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
20292f76d53SAlp Dener     if (badstep || powell || periodic) {
20392f76d53SAlp Dener       if (qn->monflg) {
2049566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
2056bdcc836SBarry Smith         if (powell) {
20663a3b9bcSJacob 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));
2076bdcc836SBarry Smith         } else {
20863a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %" PetscInt_FMT "\n", i_r));
2096bdcc836SBarry Smith         }
2109566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
2115ba6227bSPeter Brune       }
2125ba6227bSPeter Brune       i_r = -1;
2130c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
2149566063dSJacob Faibussowitsch         PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
21507b62357SFande Kong         SNESCheckJacobianDomainerror(snes);
2160ecc9e1dSPeter Brune       }
2179566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(qn->B, PETSC_FALSE));
2180ecc9e1dSPeter Brune     }
2195ba6227bSPeter Brune   }
2204b11644fSPeter Brune   if (i == snes->max_its) {
22163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its));
2224b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
2234b11644fSPeter Brune   }
2244b11644fSPeter Brune   PetscFunctionReturn(0);
2254b11644fSPeter Brune }
2264b11644fSPeter Brune 
2279371c9d4SSatish Balay static PetscErrorCode SNESSetUp_QN(SNES snes) {
2289f83bee8SJed Brown   SNES_QN *qn = (SNES_QN *)snes->data;
229fced5a79SAsbjørn Nilsen Riseth   DM       dm;
23092f76d53SAlp Dener   PetscInt n, N;
231335fb975SPeter Brune 
2324b11644fSPeter Brune   PetscFunctionBegin;
233fced5a79SAsbjørn Nilsen Riseth 
234fced5a79SAsbjørn Nilsen Riseth   if (!snes->vec_sol) {
2359566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dm));
2369566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol));
237fced5a79SAsbjørn Nilsen Riseth   }
2389566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 4));
23992f76d53SAlp Dener 
24048a46eb9SPierre Jolivet   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) PetscCall(SNESSetUpMatrices(snes));
241ad540459SPierre Jolivet   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
24292f76d53SAlp Dener 
24360c08b40SPeter Brune   /* set method defaults */
2441efc8c45SPeter Brune   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
24560c08b40SPeter Brune     if (qn->type == SNES_QN_BADBROYDEN) {
24660c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_NONE;
24760c08b40SPeter Brune     } else {
24892f76d53SAlp Dener       qn->scale_type = SNES_QN_SCALE_SCALAR;
24960c08b40SPeter Brune     }
25060c08b40SPeter Brune   }
2511efc8c45SPeter Brune   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
25260c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
25360c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_POWELL;
25460c08b40SPeter Brune     } else {
25560c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_PERIODIC;
25660c08b40SPeter Brune     }
25760c08b40SPeter Brune   }
25860c08b40SPeter Brune 
25992f76d53SAlp Dener   /* Set up the LMVM matrix */
26092f76d53SAlp Dener   switch (qn->type) {
26192f76d53SAlp Dener   case SNES_QN_BROYDEN:
2629566063dSJacob Faibussowitsch     PetscCall(MatSetType(qn->B, MATLMVMBROYDEN));
26392f76d53SAlp Dener     qn->scale_type = SNES_QN_SCALE_NONE;
26492f76d53SAlp Dener     break;
26592f76d53SAlp Dener   case SNES_QN_BADBROYDEN:
2669566063dSJacob Faibussowitsch     PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN));
26792f76d53SAlp Dener     qn->scale_type = SNES_QN_SCALE_NONE;
26892f76d53SAlp Dener     break;
26992f76d53SAlp Dener   default:
2709566063dSJacob Faibussowitsch     PetscCall(MatSetType(qn->B, MATLMVMBFGS));
27192f76d53SAlp Dener     switch (qn->scale_type) {
2729371c9d4SSatish Balay     case SNES_QN_SCALE_NONE: PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE)); break;
2739371c9d4SSatish Balay     case SNES_QN_SCALE_SCALAR: PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR)); break;
2749371c9d4SSatish Balay     case SNES_QN_SCALE_JACOBIAN: PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_USER)); break;
27592f76d53SAlp Dener     case SNES_QN_SCALE_DIAGONAL:
27692f76d53SAlp Dener     case SNES_QN_SCALE_DEFAULT:
2779371c9d4SSatish Balay     default: break;
2788e231d97SPeter Brune     }
27992f76d53SAlp Dener     break;
28092f76d53SAlp Dener   }
2819566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
2829566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol, &N));
2839566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(qn->B, n, n, N, N));
2849566063dSJacob Faibussowitsch   PetscCall(MatSetUp(qn->B));
2859566063dSJacob Faibussowitsch   PetscCall(MatLMVMReset(qn->B, PETSC_TRUE));
2869566063dSJacob Faibussowitsch   PetscCall(MatLMVMSetHistorySize(qn->B, qn->m));
2879566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func));
2884b11644fSPeter Brune   PetscFunctionReturn(0);
2894b11644fSPeter Brune }
2904b11644fSPeter Brune 
2919371c9d4SSatish Balay static PetscErrorCode SNESReset_QN(SNES snes) {
2929f83bee8SJed Brown   SNES_QN *qn;
2930adebc6cSBarry Smith 
2944b11644fSPeter Brune   PetscFunctionBegin;
2954b11644fSPeter Brune   if (snes->data) {
2969f83bee8SJed Brown     qn = (SNES_QN *)snes->data;
2979566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&qn->B));
2984b11644fSPeter Brune   }
2994b11644fSPeter Brune   PetscFunctionReturn(0);
3004b11644fSPeter Brune }
3014b11644fSPeter Brune 
3029371c9d4SSatish Balay static PetscErrorCode SNESDestroy_QN(SNES snes) {
3034b11644fSPeter Brune   PetscFunctionBegin;
3049566063dSJacob Faibussowitsch   PetscCall(SNESReset_QN(snes));
3059566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3062e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", NULL));
3072e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", NULL));
3082e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", NULL));
3094b11644fSPeter Brune   PetscFunctionReturn(0);
3104b11644fSPeter Brune }
3114b11644fSPeter Brune 
3129371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_QN(SNES snes, PetscOptionItems *PetscOptionsObject) {
3132150357eSBarry Smith   SNES_QN          *qn = (SNES_QN *)snes->data;
31492f76d53SAlp Dener   PetscBool         flg;
315f1c6b773SPeter Brune   SNESLineSearch    linesearch;
3162150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
3172150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
3181efc8c45SPeter Brune   SNESQNType        qtype = qn->type;
3192150357eSBarry Smith 
3204b11644fSPeter Brune   PetscFunctionBegin;
321d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES QN options");
3229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, NULL));
3239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL));
3249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", qn->monflg, &qn->monflg, NULL));
3259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_scale_type", "Scaling type", "SNESQNSetScaleType", SNESQNScaleTypes, (PetscEnum)stype, (PetscEnum *)&stype, &flg));
3269566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetScaleType(snes, stype));
32788f769c5SPeter Brune 
3289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_restart_type", "Restart type", "SNESQNSetRestartType", SNESQNRestartTypes, (PetscEnum)rtype, (PetscEnum *)&rtype, &flg));
3299566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetRestartType(snes, rtype));
33088f769c5SPeter Brune 
3319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_type", "Quasi-Newton update type", "", SNESQNTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, &flg));
3329566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetType(snes, qtype));
3339566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(qn->B));
334d0609cedSBarry Smith   PetscOptionsHeadEnd();
3359e764e56SPeter Brune   if (!snes->linesearch) {
3369566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
337ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
33860c08b40SPeter Brune       if (qn->type == SNES_QN_LBFGS) {
3399566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
34060c08b40SPeter Brune       } else if (qn->type == SNES_QN_BROYDEN) {
3419566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
34260c08b40SPeter Brune       } else {
3439566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
34460c08b40SPeter Brune       }
3459e764e56SPeter Brune     }
346ec786807SJed Brown   }
34748a46eb9SPierre Jolivet   if (qn->monflg) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor));
3484b11644fSPeter Brune   PetscFunctionReturn(0);
3494b11644fSPeter Brune }
3504b11644fSPeter Brune 
3519371c9d4SSatish Balay static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer) {
3525cd83419SPeter Brune   SNES_QN  *qn = (SNES_QN *)snes->data;
3535cd83419SPeter Brune   PetscBool iascii;
3545cd83419SPeter Brune 
3555cd83419SPeter Brune   PetscFunctionBegin;
3569566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3575cd83419SPeter Brune   if (iascii) {
3589566063dSJacob 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]));
35963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Stored subspace size: %" PetscInt_FMT "\n", qn->m));
3605cd83419SPeter Brune   }
3615cd83419SPeter Brune   PetscFunctionReturn(0);
3625cd83419SPeter Brune }
36392c02d66SPeter Brune 
3640c777b0cSPeter Brune /*@
365*f6dfbefdSBarry Smith     SNESQNSetRestartType - Sets the restart type for `SNESQN`.
3660c777b0cSPeter Brune 
367*f6dfbefdSBarry Smith     Logically Collective on snes
3680c777b0cSPeter Brune 
3690c777b0cSPeter Brune     Input Parameters:
3700c777b0cSPeter Brune +   snes - the iterative context
3710c777b0cSPeter Brune -   rtype - restart type
3720c777b0cSPeter Brune 
373*f6dfbefdSBarry Smith     Options Database Keys:
3740c777b0cSPeter Brune +   -snes_qn_restart_type <powell,periodic,none> - set the restart type
37584c577daSBarry Smith -   -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
3760c777b0cSPeter Brune 
3770c777b0cSPeter Brune     Level: intermediate
3780c777b0cSPeter Brune 
379*f6dfbefdSBarry Smith     `SNESQNRestartType`s:
380*f6dfbefdSBarry Smith +   `SNES_QN_RESTART_NONE` - never restart
381*f6dfbefdSBarry Smith .   `SNES_QN_RESTART_POWELL` - restart based upon descent criteria
382*f6dfbefdSBarry Smith -   `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations
3830c777b0cSPeter Brune 
384*f6dfbefdSBarry Smith .seealso: `SNESQN`, `SNESQNRestartType`, `SNES_QN_RESTART_NONE`, `SNES_QN_RESTART_POWELL`, `SNES_QN_RESTART_PERIODIC`
3850c777b0cSPeter Brune @*/
3869371c9d4SSatish Balay PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) {
3870c777b0cSPeter Brune   PetscFunctionBegin;
3880c777b0cSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
389cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetRestartType_C", (SNES, SNESQNRestartType), (snes, rtype));
3900c777b0cSPeter Brune   PetscFunctionReturn(0);
3910c777b0cSPeter Brune }
3920c777b0cSPeter Brune 
3930c777b0cSPeter Brune /*@
394*f6dfbefdSBarry Smith     SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in `SNESQN`.
3950c777b0cSPeter Brune 
396*f6dfbefdSBarry Smith     Logically Collective on snes
3970c777b0cSPeter Brune 
3980c777b0cSPeter Brune     Input Parameters:
399*f6dfbefdSBarry Smith +   snes - the nonlinear solver context
4000c777b0cSPeter Brune -   stype - scale type
4010c777b0cSPeter Brune 
4020c777b0cSPeter Brune     Options Database:
40367b8a455SSatish Balay .   -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
4040c777b0cSPeter Brune 
4050c777b0cSPeter Brune     Level: intermediate
4060c777b0cSPeter Brune 
407*f6dfbefdSBarry Smith     `SNESQNScaleType`s:
408*f6dfbefdSBarry Smith +   `SNES_QN_SCALE_NONE` - don't scale the problem
409*f6dfbefdSBarry Smith .   `SNES_QN_SCALE_SCALAR` - use Shanno scaling
410*f6dfbefdSBarry Smith .   `SNES_QN_SCALE_DIAGONAL` - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available
411*f6dfbefdSBarry Smith -   `SNES_QN_SCALE_JACOBIAN` - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration
412a01a0525SBarry Smith                              of QN and at ever restart.
4130c777b0cSPeter Brune 
414db781477SPatrick Sanan .seealso: `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()`
4150c777b0cSPeter Brune @*/
4160c777b0cSPeter Brune 
4179371c9d4SSatish Balay PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) {
4180c777b0cSPeter Brune   PetscFunctionBegin;
4190c777b0cSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
420cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetScaleType_C", (SNES, SNESQNScaleType), (snes, stype));
4210c777b0cSPeter Brune   PetscFunctionReturn(0);
4220c777b0cSPeter Brune }
4230c777b0cSPeter Brune 
4249371c9d4SSatish Balay PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) {
4250c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4266e111a19SKarl Rupp 
4270c777b0cSPeter Brune   PetscFunctionBegin;
4280c777b0cSPeter Brune   qn->scale_type = stype;
4297044b745SBarry Smith   if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE;
4300c777b0cSPeter Brune   PetscFunctionReturn(0);
4310c777b0cSPeter Brune }
4320c777b0cSPeter Brune 
4339371c9d4SSatish Balay PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) {
4340c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4356e111a19SKarl Rupp 
4360c777b0cSPeter Brune   PetscFunctionBegin;
4370c777b0cSPeter Brune   qn->restart_type = rtype;
4380c777b0cSPeter Brune   PetscFunctionReturn(0);
4390c777b0cSPeter Brune }
4400c777b0cSPeter Brune 
4411efc8c45SPeter Brune /*@
442*f6dfbefdSBarry Smith     SNESQNSetType - Sets the quasi-Newton variant to be used in `SNESQN`.
4431efc8c45SPeter Brune 
444*f6dfbefdSBarry Smith     Logically Collective on snes
4451efc8c45SPeter Brune 
4461efc8c45SPeter Brune     Input Parameters:
4471efc8c45SPeter Brune +   snes - the iterative context
4481efc8c45SPeter Brune -   qtype - variant type
4491efc8c45SPeter Brune 
4501efc8c45SPeter Brune     Options Database:
45167b8a455SSatish Balay .   -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
4521efc8c45SPeter Brune 
4531efc8c45SPeter Brune     Level: beginner
4541efc8c45SPeter Brune 
455*f6dfbefdSBarry Smith     `SNESQNType`s:
456*f6dfbefdSBarry Smith +   `SNES_QN_LBFGS` - LBFGS variant
457*f6dfbefdSBarry Smith .   `SNES_QN_BROYDEN` - Broyden variant
458*f6dfbefdSBarry Smith -   `SNES_QN_BADBROYDEN` - Bad Broyden variant
4591efc8c45SPeter Brune 
460*f6dfbefdSBarry Smith .seealso: `SNESQN`, `SNES_QN_LBFGS`, `SNES_QN_BROYDEN`, `SNES_QN_BADBROYDEN`, `SNESQNType`, `TAOLMVM`, `TAOBLMVM`
4611efc8c45SPeter Brune @*/
4621efc8c45SPeter Brune 
4639371c9d4SSatish Balay PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype) {
4641efc8c45SPeter Brune   PetscFunctionBegin;
4651efc8c45SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
466cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetType_C", (SNES, SNESQNType), (snes, qtype));
4671efc8c45SPeter Brune   PetscFunctionReturn(0);
4681efc8c45SPeter Brune }
4691efc8c45SPeter Brune 
4709371c9d4SSatish Balay PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype) {
4711efc8c45SPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4721efc8c45SPeter Brune 
4731efc8c45SPeter Brune   PetscFunctionBegin;
4741efc8c45SPeter Brune   qn->type = qtype;
4751efc8c45SPeter Brune   PetscFunctionReturn(0);
4761efc8c45SPeter Brune }
4771efc8c45SPeter Brune 
4784b11644fSPeter Brune /*MC
4794b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
4804b11644fSPeter Brune 
481*f6dfbefdSBarry Smith       Options Database Keys:
48284c577daSBarry Smith +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
483*f6dfbefdSBarry Smith .     -snes_qn_restart_type <powell,periodic,none> - set the restart type
4846bdcc836SBarry Smith .     -snes_qn_powell_gamma - Angle condition for restart.
4851867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
48684c577daSBarry Smith .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
48792f76d53SAlp Dener .     -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian
48872835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
48984c577daSBarry Smith -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.
4906cc8130cSPeter Brune 
4916cc8130cSPeter Brune       References:
492606c0280SSatish Balay +   * -   Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
493606c0280SSatish Balay .   * -   R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods,
4940a8e8854SPeter Brune       Technical Report, Northwestern University, June 1992.
495606c0280SSatish Balay .   * -   Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
4960a8e8854SPeter Brune       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
497606c0280SSatish Balay .   * -   Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
4984f02bc6aSBarry Smith        SIAM Review, 57(4), 2015
499606c0280SSatish Balay .   * -   Griewank, Andreas. "Broyden updating, the good and the bad!." Doc. Math (2012): 301-315.
500606c0280SSatish Balay .   * -   Gilbert, Jean Charles, and Claude Lemar{\'e}chal. "Some numerical experiments with variable-storage quasi-Newton algorithms."
50192f76d53SAlp Dener       Mathematical programming 45.1-3 (1989): 407-435.
502606c0280SSatish Balay -   * -   Dener A., Munson T. "Accelerating Limited-Memory Quasi-Newton Convergence for Large-Scale Optimization"
503110fc3b0SBarry Smith       Computational Science - ICCS 2019. ICCS 2019. Lecture Notes in Computer Science, vol 11538. Springer, Cham
5044b11644fSPeter Brune 
5054b11644fSPeter Brune       Level: beginner
5064b11644fSPeter Brune 
507*f6dfbefdSBarry Smith       Notes:
508*f6dfbefdSBarry Smith       This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
509*f6dfbefdSBarry Smith       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
510*f6dfbefdSBarry Smith       updates.
5116cc8130cSPeter Brune 
512*f6dfbefdSBarry Smith       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
513*f6dfbefdSBarry Smith       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
514*f6dfbefdSBarry Smith       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
515*f6dfbefdSBarry Smith       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
516*f6dfbefdSBarry Smith 
517*f6dfbefdSBarry Smith       Uses left nonlinear preconditioning by default.
518*f6dfbefdSBarry Smith 
519*f6dfbefdSBarry Smith .seealso: `SNESQNRestartType`, `SNESQNSetRestartType()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`,
520*f6dfbefdSBarry Smith           `SNESQNScaleType`, `SNESQNSetScaleType()`, `SNESQNSetType`, `SNESQNSetType ()`
5214b11644fSPeter Brune M*/
5229371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) {
5239f83bee8SJed Brown   SNES_QN    *qn;
52492f76d53SAlp Dener   const char *optionsprefix;
5254b11644fSPeter Brune 
5264b11644fSPeter Brune   PetscFunctionBegin;
5274b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
5284b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
5294b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
5304b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
5315cd83419SPeter Brune   snes->ops->view           = SNESView_QN;
5324b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
5334b11644fSPeter Brune 
534efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
53547f26062SPeter Brune 
536efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
53742f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
53842f4f86dSBarry Smith 
5394fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
5404fc747eaSLawrence Mitchell 
54188976e71SPeter Brune   if (!snes->tolerancesset) {
5420e444f03SPeter Brune     snes->max_funcs = 30000;
5430e444f03SPeter Brune     snes->max_its   = 10000;
54488976e71SPeter Brune   }
5450e444f03SPeter Brune 
5469566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &qn));
5474b11644fSPeter Brune   snes->data       = (void *)qn;
5480ecc9e1dSPeter Brune   qn->m            = 10;
549b21d5a53SPeter Brune   qn->scaling      = 1.0;
5500298fd71SBarry Smith   qn->monitor      = NULL;
55192f76d53SAlp Dener   qn->monflg       = PETSC_FALSE;
552b8840d6bSPeter Brune   qn->powell_gamma = 0.9999;
55360c08b40SPeter Brune   qn->scale_type   = SNES_QN_SCALE_DEFAULT;
55460c08b40SPeter Brune   qn->restart_type = SNES_QN_RESTART_DEFAULT;
555b8840d6bSPeter Brune   qn->type         = SNES_QN_LBFGS;
556ea630c6eSPeter Brune 
5579566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B));
5589566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
5599566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix));
56092f76d53SAlp Dener 
5619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", SNESQNSetScaleType_QN));
5629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", SNESQNSetRestartType_QN));
5639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", SNESQNSetType_QN));
5644b11644fSPeter Brune   PetscFunctionReturn(0);
5654b11644fSPeter Brune }
566