xref: /petsc/src/snes/impls/qn/qn.c (revision bd89dbf26d8a5efecb980364933175da61864cd7)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2fced5a79SAsbjørn Nilsen Riseth #include <petscdm.h>
34b11644fSPeter Brune 
49e5d0892SLisandro Dalcin const char *const SNESQNScaleTypes[]   = {"DEFAULT", "NONE", "SCALAR", "DIAGONAL", "JACOBIAN", "SNESQNScaleType", "SNES_QN_SCALING_", NULL};
59e5d0892SLisandro Dalcin const char *const SNESQNRestartTypes[] = {"DEFAULT", "NONE", "POWELL", "PERIODIC", "SNESQNRestartType", "SNES_QN_RESTART_", NULL};
69e5d0892SLisandro Dalcin const char *const SNESQNTypes[]        = {"LBFGS", "BROYDEN", "BADBROYDEN", "SNESQNType", "SNES_QN_", NULL};
7b8840d6bSPeter Brune 
84b11644fSPeter Brune typedef struct {
992f76d53SAlp Dener   Mat               B;      /* Quasi-Newton approximation Matrix (MATLMVM) */
108e231d97SPeter Brune   PetscInt          m;      /* The number of kept previous steps */
115c7a0a03SPeter Brune   PetscReal        *lambda; /* The line search history of the method */
1292f76d53SAlp Dener   PetscBool         monflg;
1344f7e39eSPeter Brune   PetscViewer       monitor;
146bf1b2e5SPeter Brune   PetscReal         powell_gamma; /* Powell angle restart condition */
15b21d5a53SPeter Brune   PetscReal         scaling;      /* scaling of H0 */
16b8840d6bSPeter Brune   SNESQNType        type;         /* the type of quasi-newton method used */
1788f769c5SPeter Brune   SNESQNScaleType   scale_type;   /* the type of scaling used */
180c777b0cSPeter Brune   SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */
199f83bee8SJed Brown } SNES_QN;
204b11644fSPeter Brune 
SNESQNGetMatrix_Private(SNES snes,Mat * B)2106594da7SStefano Zampini static PetscErrorCode SNESQNGetMatrix_Private(SNES snes, Mat *B)
2206594da7SStefano Zampini {
2306594da7SStefano Zampini   SNES_QN    *qn = (SNES_QN *)snes->data;
2406594da7SStefano Zampini   const char *optionsprefix;
2506594da7SStefano Zampini 
2606594da7SStefano Zampini   PetscFunctionBegin;
2706594da7SStefano Zampini   if (!qn->B) {
2806594da7SStefano Zampini     PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B));
2906594da7SStefano Zampini     PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
3006594da7SStefano Zampini     PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix));
3106594da7SStefano Zampini     PetscCall(MatAppendOptionsPrefix(qn->B, "qn_"));
3206594da7SStefano Zampini     switch (qn->type) {
3306594da7SStefano Zampini     case SNES_QN_BROYDEN:
3406594da7SStefano Zampini       PetscCall(MatSetType(qn->B, MATLMVMBROYDEN));
3506594da7SStefano Zampini       qn->scale_type = SNES_QN_SCALE_NONE;
3606594da7SStefano Zampini       break;
3706594da7SStefano Zampini     case SNES_QN_BADBROYDEN:
3806594da7SStefano Zampini       PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN));
3906594da7SStefano Zampini       qn->scale_type = SNES_QN_SCALE_NONE;
4006594da7SStefano Zampini       break;
4106594da7SStefano Zampini     default:
4206594da7SStefano Zampini       PetscCall(MatSetType(qn->B, MATLMVMBFGS));
4306594da7SStefano Zampini       switch (qn->scale_type) {
4406594da7SStefano Zampini       case SNES_QN_SCALE_NONE:
455839f191SToby Isaac       case SNES_QN_SCALE_JACOBIAN:
4606594da7SStefano Zampini         PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE));
4706594da7SStefano Zampini         break;
4806594da7SStefano Zampini       case SNES_QN_SCALE_SCALAR:
4906594da7SStefano Zampini       case SNES_QN_SCALE_DEFAULT:
5006594da7SStefano Zampini         PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR));
5106594da7SStefano Zampini         break;
5206594da7SStefano Zampini       default:
5306594da7SStefano Zampini         break;
5406594da7SStefano Zampini       }
5506594da7SStefano Zampini       break;
5606594da7SStefano Zampini     }
5706594da7SStefano Zampini   }
5806594da7SStefano Zampini   *B = qn->B;
5906594da7SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
6006594da7SStefano Zampini }
6106594da7SStefano Zampini 
SNESSolve_QN(SNES snes)62d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_QN(SNES snes)
63d71ae5a4SJacob Faibussowitsch {
649f83bee8SJed Brown   SNES_QN             *qn = (SNES_QN *)snes->data;
6506594da7SStefano Zampini   Vec                  X, F, W;
6647f26062SPeter Brune   Vec                  Y, D, Dold;
67b8840d6bSPeter Brune   PetscInt             i, i_r;
684279555eSSatish Balay   PetscReal            fnorm, xnorm, ynorm;
69*76c63389SBarry Smith   SNESLineSearchReason lsreason;
70*76c63389SBarry Smith   PetscBool            powell, periodic, restart;
711c6523dcSPeter Brune   PetscScalar          DolddotD, DolddotDold;
72b7281c8aSPeter Brune   SNESConvergedReason  reason;
734279555eSSatish Balay #if defined(PETSC_USE_INFO)
744279555eSSatish Balay   PetscReal gnorm;
754279555eSSatish Balay #endif
760ecc9e1dSPeter Brune 
776e111a19SKarl Rupp   PetscFunctionBegin;
780b121fc5SBarry 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);
79c579b300SPatrick Farrell 
809566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
819f3a0142SPeter Brune   F = snes->vec_func;       /* residual vector */
823af51624SPeter Brune   Y = snes->vec_sol_update; /* search direction generated by J^-1D*/
83b8840d6bSPeter Brune   X = snes->vec_sol;        /* solution vector */
843af51624SPeter Brune 
8506594da7SStefano Zampini   /* directions */
8606594da7SStefano Zampini   W    = snes->work[0];
87335fb975SPeter Brune   D    = snes->work[1];
88335fb975SPeter Brune   Dold = snes->work[2];
894b11644fSPeter Brune 
904b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
914b11644fSPeter Brune 
929566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
934b11644fSPeter Brune   snes->iter = 0;
944b11644fSPeter Brune   snes->norm = 0.;
959566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
9647f26062SPeter Brune 
97efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
9806594da7SStefano Zampini     /* we need to sample the preconditioned function only once here,
9906594da7SStefano Zampini        since linesearch will always use the preconditioned function */
1009566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
1019566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
10206594da7SStefano Zampini     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
103b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1043ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
105b7281c8aSPeter Brune     }
10647f26062SPeter Brune   } else {
10706594da7SStefano Zampini     if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F));
10806594da7SStefano Zampini     else snes->vec_func_init_set = PETSC_FALSE;
10906594da7SStefano Zampini   }
1109566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm));
111*76c63389SBarry Smith   SNESCheckFunctionDomainError(snes, fnorm);
112b28a06ddSPeter Brune 
1139566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1144b11644fSPeter Brune   snes->norm = fnorm;
1159566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1169566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
1179566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
1184b11644fSPeter Brune 
1194b11644fSPeter Brune   /* test convergence */
1202d157150SStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
1213ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
12270d3b23bSPeter Brune 
12306594da7SStefano Zampini   restart = PETSC_TRUE;
12406594da7SStefano Zampini   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
12506594da7SStefano Zampini     PetscScalar gdx;
12606594da7SStefano Zampini 
12706594da7SStefano Zampini     /* general purpose update */
12806594da7SStefano Zampini     PetscTryTypeMethod(snes, update, snes->iter);
12906594da7SStefano Zampini 
13006594da7SStefano Zampini     if (snes->npc) {
13106594da7SStefano Zampini       if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
13206594da7SStefano Zampini         PetscCall(SNESApplyNPC(snes, X, F, D));
13306594da7SStefano Zampini         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
13406594da7SStefano Zampini         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
13506594da7SStefano Zampini           snes->reason = SNES_DIVERGED_INNER;
13606594da7SStefano Zampini           PetscFunctionReturn(PETSC_SUCCESS);
13706594da7SStefano Zampini         }
13806594da7SStefano Zampini       } else if (snes->npcside == PC_RIGHT) {
1399566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0));
1409566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1419566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0));
1429566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
14306594da7SStefano Zampini         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
144ddd40ce5SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
1453ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
146ddd40ce5SPeter Brune         }
1479566063dSJacob Faibussowitsch         PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
1489566063dSJacob Faibussowitsch         PetscCall(VecCopy(F, D));
14906594da7SStefano Zampini       } else {
15006594da7SStefano Zampini         PetscCall(VecCopy(F, D));
15106594da7SStefano Zampini       }
15206594da7SStefano Zampini     } else {
15306594da7SStefano Zampini       PetscCall(VecCopy(F, D));
1543cf07b75SPeter Brune     }
1553cf07b75SPeter Brune 
156f8e15203SPeter Brune     /* scale the initial update */
15706594da7SStefano Zampini     if (qn->scale_type == SNES_QN_SCALE_JACOBIAN && restart) {
1589566063dSJacob Faibussowitsch       PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
159*76c63389SBarry Smith       SNESCheckJacobianDomainError(snes);
1609566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
1619566063dSJacob Faibussowitsch       PetscCall(MatLMVMSetJ0KSP(qn->B, snes->ksp));
1620ecc9e1dSPeter Brune     }
1630ecc9e1dSPeter Brune 
16492f76d53SAlp Dener     /* update QN approx and calculate step */
1659566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(qn->B, X, D));
1669566063dSJacob Faibussowitsch     PetscCall(MatSolve(qn->B, D, Y));
16706594da7SStefano Zampini     PetscCall(VecDot(F, Y, &gdx));
16806594da7SStefano Zampini     if (PetscAbsScalar(gdx) <= 0) {
16906594da7SStefano Zampini       /* Step is not descent or solve was not successful
17006594da7SStefano Zampini          Use steepest descent direction */
17106594da7SStefano Zampini       PetscCall(VecCopy(F, Y));
17206594da7SStefano Zampini       PetscCall(MatLMVMReset(qn->B, PETSC_FALSE));
17306594da7SStefano Zampini     }
17492f76d53SAlp Dener 
17570d3b23bSPeter Brune     /* line search for lambda */
1769371c9d4SSatish Balay     ynorm = 1;
1774279555eSSatish Balay #if defined(PETSC_USE_INFO)
1789371c9d4SSatish Balay     gnorm = fnorm;
1794279555eSSatish Balay #endif
1809566063dSJacob Faibussowitsch     PetscCall(VecCopy(D, Dold));
1819566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
182*76c63389SBarry Smith     if (snes->reason) break;
183*76c63389SBarry Smith     SNESCheckLineSearchFailure(snes);
184*76c63389SBarry Smith     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsreason));
1853af51624SPeter Brune 
186*76c63389SBarry Smith     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
18788d374b2SPeter Brune     /* convergence monitoring */
188*76c63389SBarry Smith     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lsreason));
189b21d5a53SPeter Brune 
19006594da7SStefano Zampini     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
19106594da7SStefano Zampini     snes->iter  = i + 1;
19271dbe336SPeter Brune     snes->norm  = fnorm;
193c1e67a49SFande Kong     snes->xnorm = xnorm;
194c1e67a49SFande Kong     snes->ynorm = ynorm;
19506594da7SStefano Zampini     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
196360c497dSPeter Brune 
1979566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
19892f76d53SAlp Dener 
1994b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
2002d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
2012d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
2023ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
20306594da7SStefano Zampini 
204efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2059566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, X, F, D));
2069566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
20706594da7SStefano Zampini       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
208b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
2093ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
210b7281c8aSPeter Brune       }
21188d374b2SPeter Brune     } else {
2129566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, D));
21388d374b2SPeter Brune     }
21401fe78eaSStefano Zampini 
21592f76d53SAlp Dener     /* restart conditions */
2160c777b0cSPeter Brune     powell = PETSC_FALSE;
2176bdcc836SBarry Smith     if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) {
2186bdcc836SBarry Smith       /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */
21923c918e7SPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
220574a8f54SStefano Zampini         PetscCall(MatMult(snes->jacobian, Dold, W));
22123c918e7SPeter Brune       } else {
2229566063dSJacob Faibussowitsch         PetscCall(VecCopy(Dold, W));
22323c918e7SPeter Brune       }
2249566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(W, Dold, &DolddotDold));
2259566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(W, D, &DolddotD));
2269566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(W, Dold, &DolddotDold));
2279566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(W, D, &DolddotD));
2280c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma * PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
2290c777b0cSPeter Brune     }
2300c777b0cSPeter Brune     periodic = PETSC_FALSE;
231b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
232b8840d6bSPeter Brune       if (i_r > qn->m - 1) periodic = PETSC_TRUE;
2330c777b0cSPeter Brune     }
23406594da7SStefano Zampini 
2350c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
23606594da7SStefano Zampini     restart = PETSC_FALSE;
237*76c63389SBarry Smith     if (lsreason || powell || periodic) {
23806594da7SStefano Zampini       restart = PETSC_TRUE;
23992f76d53SAlp Dener       if (qn->monflg) {
2409566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
2416bdcc836SBarry Smith         if (powell) {
24263a3b9bcSJacob 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));
2436bdcc836SBarry Smith         } else {
24463a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %" PetscInt_FMT "\n", i_r));
2456bdcc836SBarry Smith         }
2469566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(qn->monitor, ((PetscObject)snes)->tablevel + 2));
2475ba6227bSPeter Brune       }
2485ba6227bSPeter Brune       i_r = -1;
2499566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(qn->B, PETSC_FALSE));
2500ecc9e1dSPeter Brune     }
2515ba6227bSPeter Brune   }
2523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2534b11644fSPeter Brune }
2544b11644fSPeter Brune 
SNESSetUp_QN(SNES snes)255d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_QN(SNES snes)
256d71ae5a4SJacob Faibussowitsch {
2579f83bee8SJed Brown   SNES_QN *qn = (SNES_QN *)snes->data;
258fced5a79SAsbjørn Nilsen Riseth   DM       dm;
25992f76d53SAlp Dener   PetscInt n, N;
260335fb975SPeter Brune 
2614b11644fSPeter Brune   PetscFunctionBegin;
262fced5a79SAsbjørn Nilsen Riseth   if (!snes->vec_sol) {
2639566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dm));
2649566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol));
265fced5a79SAsbjørn Nilsen Riseth   }
26606594da7SStefano Zampini   PetscCall(SNESSetWorkVecs(snes, 3));
26706594da7SStefano Zampini   PetscCall(SNESQNGetMatrix_Private(snes, &qn->B));
26892f76d53SAlp Dener 
26948a46eb9SPierre Jolivet   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) PetscCall(SNESSetUpMatrices(snes));
27006594da7SStefano Zampini   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
27192f76d53SAlp Dener 
27260c08b40SPeter Brune   /* set method defaults */
2731efc8c45SPeter Brune   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
27460c08b40SPeter Brune     if (qn->type == SNES_QN_BADBROYDEN) {
27560c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_NONE;
27660c08b40SPeter Brune     } else {
27792f76d53SAlp Dener       qn->scale_type = SNES_QN_SCALE_SCALAR;
27860c08b40SPeter Brune     }
27960c08b40SPeter Brune   }
2801efc8c45SPeter Brune   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
28160c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
28260c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_POWELL;
28360c08b40SPeter Brune     } else {
28460c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_PERIODIC;
28560c08b40SPeter Brune     }
28660c08b40SPeter Brune   }
28792f76d53SAlp Dener   /* Set up the LMVM matrix */
28806594da7SStefano Zampini   PetscCall(SNESQNGetMatrix_Private(snes, &qn->B));
2899566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
2909566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol, &N));
2919566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(qn->B, n, n, N, N));
2929566063dSJacob Faibussowitsch   PetscCall(MatSetUp(qn->B));
2939566063dSJacob Faibussowitsch   PetscCall(MatLMVMReset(qn->B, PETSC_TRUE));
2949566063dSJacob Faibussowitsch   PetscCall(MatLMVMSetHistorySize(qn->B, qn->m));
2959566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func));
2963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2974b11644fSPeter Brune }
2984b11644fSPeter Brune 
SNESReset_QN(SNES snes)299d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_QN(SNES snes)
300d71ae5a4SJacob Faibussowitsch {
30106594da7SStefano Zampini   SNES_QN *qn = (SNES_QN *)snes->data;
3020adebc6cSBarry Smith 
3034b11644fSPeter Brune   PetscFunctionBegin;
3049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&qn->B));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3064b11644fSPeter Brune }
3074b11644fSPeter Brune 
SNESDestroy_QN(SNES snes)308d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_QN(SNES snes)
309d71ae5a4SJacob Faibussowitsch {
3104b11644fSPeter Brune   PetscFunctionBegin;
3119566063dSJacob Faibussowitsch   PetscCall(SNESReset_QN(snes));
3129566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3132e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", NULL));
3142e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", NULL));
3152e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", NULL));
3163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3174b11644fSPeter Brune }
3184b11644fSPeter Brune 
SNESSetFromOptions_QN(SNES snes,PetscOptionItems PetscOptionsObject)319ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_QN(SNES snes, PetscOptionItems PetscOptionsObject)
320d71ae5a4SJacob Faibussowitsch {
3212150357eSBarry Smith   SNES_QN          *qn = (SNES_QN *)snes->data;
32292f76d53SAlp Dener   PetscBool         flg;
323f1c6b773SPeter Brune   SNESLineSearch    linesearch;
3242150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
3252150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
3261efc8c45SPeter Brune   SNESQNType        qtype = qn->type;
3272150357eSBarry Smith 
3284b11644fSPeter Brune   PetscFunctionBegin;
329d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES QN options");
3309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, NULL));
3319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL));
3329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", qn->monflg, &qn->monflg, NULL));
3339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_scale_type", "Scaling type", "SNESQNSetScaleType", SNESQNScaleTypes, (PetscEnum)stype, (PetscEnum *)&stype, &flg));
3349566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetScaleType(snes, stype));
33588f769c5SPeter Brune 
3369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_restart_type", "Restart type", "SNESQNSetRestartType", SNESQNRestartTypes, (PetscEnum)rtype, (PetscEnum *)&rtype, &flg));
3379566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetRestartType(snes, rtype));
33888f769c5SPeter Brune 
3399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_type", "Quasi-Newton update type", "", SNESQNTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, &flg));
3409566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetType(snes, qtype));
341d0609cedSBarry Smith   PetscOptionsHeadEnd();
34206594da7SStefano Zampini   PetscCall(SNESQNGetMatrix_Private(snes, &qn->B));
34306594da7SStefano Zampini   PetscCall(MatSetFromOptions(qn->B));
3449e764e56SPeter Brune   if (!snes->linesearch) {
3459566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
346ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
34760c08b40SPeter Brune       if (qn->type == SNES_QN_LBFGS) {
3489566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
34960c08b40SPeter Brune       } else if (qn->type == SNES_QN_BROYDEN) {
3509566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
35160c08b40SPeter Brune       } else {
352a99ef635SJonas Heinzmann         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSECANT));
35360c08b40SPeter Brune       }
3549e764e56SPeter Brune     }
355ec786807SJed Brown   }
35648a46eb9SPierre Jolivet   if (qn->monflg) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor));
3573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3584b11644fSPeter Brune }
3594b11644fSPeter Brune 
SNESView_QN(SNES snes,PetscViewer viewer)360d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
361d71ae5a4SJacob Faibussowitsch {
3625cd83419SPeter Brune   SNES_QN  *qn = (SNES_QN *)snes->data;
3639f196a02SMartin Diehl   PetscBool isascii;
3645cd83419SPeter Brune 
3655cd83419SPeter Brune   PetscFunctionBegin;
3669f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3679f196a02SMartin Diehl   if (isascii) {
3689566063dSJacob 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]));
36963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Stored subspace size: %" PetscInt_FMT "\n", qn->m));
37006594da7SStefano Zampini     PetscCall(MatView(qn->B, viewer));
3715cd83419SPeter Brune   }
3723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3735cd83419SPeter Brune }
37492c02d66SPeter Brune 
3750c777b0cSPeter Brune /*@
376f6dfbefdSBarry Smith   SNESQNSetRestartType - Sets the restart type for `SNESQN`.
3770c777b0cSPeter Brune 
378c3339decSBarry Smith   Logically Collective
3790c777b0cSPeter Brune 
3800c777b0cSPeter Brune   Input Parameters:
3810c777b0cSPeter Brune + snes  - the iterative context
382ceaaa498SBarry Smith - rtype - restart type, see `SNESQNRestartType`
3830c777b0cSPeter Brune 
384f6dfbefdSBarry Smith   Options Database Keys:
3850c777b0cSPeter Brune + -snes_qn_restart_type <powell,periodic,none> - set the restart type
38684c577daSBarry Smith - -snes_qn_m <m>                               - sets the number of stored updates and the restart period for periodic
3870c777b0cSPeter Brune 
3880c777b0cSPeter Brune   Level: intermediate
3890c777b0cSPeter Brune 
390420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESQN`, `SNESQNRestartType`, `SNES_QN_RESTART_NONE`, `SNES_QN_RESTART_POWELL`, `SNES_QN_RESTART_PERIODIC`,
391ceaaa498SBarry Smith           `SNESQNType`, `SNESQNScaleType`
3920c777b0cSPeter Brune @*/
SNESQNSetRestartType(SNES snes,SNESQNRestartType rtype)393d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
394d71ae5a4SJacob Faibussowitsch {
3950c777b0cSPeter Brune   PetscFunctionBegin;
3960c777b0cSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
397cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetRestartType_C", (SNES, SNESQNRestartType), (snes, rtype));
3983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3990c777b0cSPeter Brune }
4000c777b0cSPeter Brune 
4010c777b0cSPeter Brune /*@
402f6dfbefdSBarry Smith   SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in `SNESQN`.
4030c777b0cSPeter Brune 
404c3339decSBarry Smith   Logically Collective
4050c777b0cSPeter Brune 
4060c777b0cSPeter Brune   Input Parameters:
407f6dfbefdSBarry Smith + snes  - the nonlinear solver context
408ceaaa498SBarry Smith - stype - scale type, see `SNESQNScaleType`
4090c777b0cSPeter Brune 
4103c7db156SBarry Smith   Options Database Key:
41167b8a455SSatish Balay . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
4120c777b0cSPeter Brune 
4130c777b0cSPeter Brune   Level: intermediate
4140c777b0cSPeter Brune 
415420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()`, `SNESQNType`, `SNESQNRestartType`
4160c777b0cSPeter Brune @*/
SNESQNSetScaleType(SNES snes,SNESQNScaleType stype)417d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
418d71ae5a4SJacob Faibussowitsch {
4190c777b0cSPeter Brune   PetscFunctionBegin;
4200c777b0cSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
421cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetScaleType_C", (SNES, SNESQNScaleType), (snes, stype));
4223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4230c777b0cSPeter Brune }
4240c777b0cSPeter Brune 
SNESQNSetScaleType_QN(SNES snes,SNESQNScaleType stype)42566976f2fSJacob Faibussowitsch static PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
426d71ae5a4SJacob Faibussowitsch {
4270c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4286e111a19SKarl Rupp 
4290c777b0cSPeter Brune   PetscFunctionBegin;
4300c777b0cSPeter Brune   qn->scale_type = stype;
4317044b745SBarry Smith   if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE;
4323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4330c777b0cSPeter Brune }
4340c777b0cSPeter Brune 
SNESQNSetRestartType_QN(SNES snes,SNESQNRestartType rtype)43566976f2fSJacob Faibussowitsch static PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
436d71ae5a4SJacob Faibussowitsch {
4370c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4386e111a19SKarl Rupp 
4390c777b0cSPeter Brune   PetscFunctionBegin;
4400c777b0cSPeter Brune   qn->restart_type = rtype;
4413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4420c777b0cSPeter Brune }
4430c777b0cSPeter Brune 
4441efc8c45SPeter Brune /*@
445f6dfbefdSBarry Smith   SNESQNSetType - Sets the quasi-Newton variant to be used in `SNESQN`.
4461efc8c45SPeter Brune 
447c3339decSBarry Smith   Logically Collective
4481efc8c45SPeter Brune 
4491efc8c45SPeter Brune   Input Parameters:
4501efc8c45SPeter Brune + snes  - the iterative context
451ceaaa498SBarry Smith - qtype - variant type, see `SNESQNType`
4521efc8c45SPeter Brune 
4533c7db156SBarry Smith   Options Database Key:
45467b8a455SSatish Balay . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
4551efc8c45SPeter Brune 
456ceaaa498SBarry Smith   Level: intermediate
4571efc8c45SPeter Brune 
458420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESQN`, `SNES_QN_LBFGS`, `SNES_QN_BROYDEN`, `SNES_QN_BADBROYDEN`, `SNESQNType`,  `SNESQNScaleType`, `TAOLMVM`, `TAOBLMVM`
4591efc8c45SPeter Brune @*/
SNESQNSetType(SNES snes,SNESQNType qtype)460d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
461d71ae5a4SJacob Faibussowitsch {
4621efc8c45SPeter Brune   PetscFunctionBegin;
4631efc8c45SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
464cac4c232SBarry Smith   PetscTryMethod(snes, "SNESQNSetType_C", (SNES, SNESQNType), (snes, qtype));
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4661efc8c45SPeter Brune }
4671efc8c45SPeter Brune 
SNESQNSetType_QN(SNES snes,SNESQNType qtype)46866976f2fSJacob Faibussowitsch static PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
469d71ae5a4SJacob Faibussowitsch {
4701efc8c45SPeter Brune   SNES_QN *qn = (SNES_QN *)snes->data;
4711efc8c45SPeter Brune 
4721efc8c45SPeter Brune   PetscFunctionBegin;
4731efc8c45SPeter Brune   qn->type = qtype;
4743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4751efc8c45SPeter Brune }
4761efc8c45SPeter Brune 
4774b11644fSPeter Brune /*MC
4784b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
4794b11644fSPeter Brune 
480f6dfbefdSBarry Smith       Options Database Keys:
48184c577daSBarry Smith +     -snes_qn_m <m>                                      - Number of past states saved for the L-Broyden methods.
482f6dfbefdSBarry Smith .     -snes_qn_restart_type <powell,periodic,none>        - set the restart type
4836bdcc836SBarry Smith .     -snes_qn_powell_gamma                               - Angle condition for restart.
4841867fe5bSPeter Brune .     -snes_qn_powell_descent                             - Descent condition for restart.
48584c577daSBarry Smith .     -snes_qn_type <lbfgs,broyden,badbroyden>            - QN type
48692f76d53SAlp Dener .     -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian
48772835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic>               - Type of line search.
48884c577daSBarry Smith -     -snes_qn_monitor                                    - Monitors the quasi-newton Jacobian.
4896cc8130cSPeter Brune 
4904b11644fSPeter Brune       Level: beginner
4914b11644fSPeter Brune 
492f6dfbefdSBarry Smith       Notes:
4931d27aa22SBarry Smith       This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of $F(x) = b$ using
4941d27aa22SBarry Smith       previous change in $F(x)$ and $x$ to form the approximate inverse Jacobian using a series of multiplicative rank-one
495f6dfbefdSBarry Smith       updates.
4966cc8130cSPeter Brune 
497f6dfbefdSBarry Smith       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
498f6dfbefdSBarry Smith       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
499f6dfbefdSBarry Smith       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
5001d27aa22SBarry Smith       perturbs the problem the Jacobian represents to be $P(x, b) - x = 0 $, where $P(x, b)$ is the preconditioner.
501f6dfbefdSBarry Smith 
502f6dfbefdSBarry Smith       Uses left nonlinear preconditioning by default.
503f6dfbefdSBarry Smith 
5041d27aa22SBarry Smith       See {cite}`byrd1994representations`, {cite}`kelley95`, {cite}`brown1985experiments`, {cite}`griewank2012broyden`, {cite}`gilbert1989some`,
5051d27aa22SBarry Smith       {cite}`dener2019accelerating`, and {cite}`bruneknepleysmithtu15`
506420bcc1bSBarry Smith 
507420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESQNRestartType`, `SNESQNSetRestartType()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`,
508420bcc1bSBarry Smith           `SNESQNScaleType`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`
5094b11644fSPeter Brune M*/
SNESCreate_QN(SNES snes)510d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
511d71ae5a4SJacob Faibussowitsch {
5129f83bee8SJed Brown   SNES_QN *qn;
5134b11644fSPeter Brune 
5144b11644fSPeter Brune   PetscFunctionBegin;
5154b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
5164b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
5174b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
5184b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
5195cd83419SPeter Brune   snes->ops->view           = SNESView_QN;
5204b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
5214b11644fSPeter Brune 
522efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
52347f26062SPeter Brune 
524efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
52542f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
52642f4f86dSBarry Smith 
5274fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
5284fc747eaSLawrence Mitchell 
52977e5a1f9SBarry Smith   PetscCall(SNESParametersInitialize(snes));
53077e5a1f9SBarry Smith   PetscObjectParameterSetDefault(snes, max_funcs, 30000);
53177e5a1f9SBarry Smith   PetscObjectParameterSetDefault(snes, max_its, 10000);
5320e444f03SPeter Brune 
5334dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&qn));
5344b11644fSPeter Brune   snes->data       = (void *)qn;
5350ecc9e1dSPeter Brune   qn->m            = 10;
536b21d5a53SPeter Brune   qn->scaling      = 1.0;
5370298fd71SBarry Smith   qn->monitor      = NULL;
53892f76d53SAlp Dener   qn->monflg       = PETSC_FALSE;
539b8840d6bSPeter Brune   qn->powell_gamma = 0.9999;
54060c08b40SPeter Brune   qn->scale_type   = SNES_QN_SCALE_DEFAULT;
54160c08b40SPeter Brune   qn->restart_type = SNES_QN_RESTART_DEFAULT;
542b8840d6bSPeter Brune   qn->type         = SNES_QN_LBFGS;
543ea630c6eSPeter Brune 
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", SNESQNSetScaleType_QN));
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", SNESQNSetRestartType_QN));
5469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", SNESQNSetType_QN));
5473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5484b11644fSPeter Brune }
549