xref: /petsc/src/snes/impls/qn/qn.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
234b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
244b11644fSPeter Brune {
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 */
93*dbbe0bcdSBarry 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 */
110*dbbe0bcdSBarry 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 */
12670d3b23bSPeter Brune     ynorm = 1; 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 */
166*dbbe0bcdSBarry 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 */
180*dbbe0bcdSBarry 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 
2274b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
2284b11644fSPeter Brune {
2299f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
230fced5a79SAsbjørn Nilsen Riseth   DM             dm;
23192f76d53SAlp Dener   PetscInt       n, N;
232335fb975SPeter Brune 
2334b11644fSPeter Brune   PetscFunctionBegin;
234fced5a79SAsbjørn Nilsen Riseth 
235fced5a79SAsbjørn Nilsen Riseth   if (!snes->vec_sol) {
2369566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes,&dm));
2379566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm,&snes->vec_sol));
238fced5a79SAsbjørn Nilsen Riseth   }
2399566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes,4));
24092f76d53SAlp Dener 
24192f76d53SAlp Dener   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
2429566063dSJacob Faibussowitsch     PetscCall(SNESSetUpMatrices(snes));
24392f76d53SAlp Dener   }
24492f76d53SAlp Dener   if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;}
24592f76d53SAlp Dener 
24660c08b40SPeter Brune   /* set method defaults */
2471efc8c45SPeter Brune   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
24860c08b40SPeter Brune     if (qn->type == SNES_QN_BADBROYDEN) {
24960c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_NONE;
25060c08b40SPeter Brune     } else {
25192f76d53SAlp Dener       qn->scale_type = SNES_QN_SCALE_SCALAR;
25260c08b40SPeter Brune     }
25360c08b40SPeter Brune   }
2541efc8c45SPeter Brune   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
25560c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
25660c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_POWELL;
25760c08b40SPeter Brune     } else {
25860c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_PERIODIC;
25960c08b40SPeter Brune     }
26060c08b40SPeter Brune   }
26160c08b40SPeter Brune 
26292f76d53SAlp Dener   /* Set up the LMVM matrix */
26392f76d53SAlp Dener   switch (qn->type) {
26492f76d53SAlp Dener     case SNES_QN_BROYDEN:
2659566063dSJacob Faibussowitsch       PetscCall(MatSetType(qn->B, MATLMVMBROYDEN));
26692f76d53SAlp Dener       qn->scale_type = SNES_QN_SCALE_NONE;
26792f76d53SAlp Dener       break;
26892f76d53SAlp Dener     case SNES_QN_BADBROYDEN:
2699566063dSJacob Faibussowitsch       PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN));
27092f76d53SAlp Dener       qn->scale_type = SNES_QN_SCALE_NONE;
27192f76d53SAlp Dener       break;
27292f76d53SAlp Dener     default:
2739566063dSJacob Faibussowitsch       PetscCall(MatSetType(qn->B, MATLMVMBFGS));
27492f76d53SAlp Dener       switch (qn->scale_type) {
27592f76d53SAlp Dener         case SNES_QN_SCALE_NONE:
2769566063dSJacob Faibussowitsch           PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE));
27792f76d53SAlp Dener           break;
27892f76d53SAlp Dener         case SNES_QN_SCALE_SCALAR:
2799566063dSJacob Faibussowitsch           PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR));
28092f76d53SAlp Dener           break;
28192f76d53SAlp Dener         case SNES_QN_SCALE_JACOBIAN:
2829566063dSJacob Faibussowitsch           PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_USER));
28392f76d53SAlp Dener           break;
28492f76d53SAlp Dener         case SNES_QN_SCALE_DIAGONAL:
28592f76d53SAlp Dener         case SNES_QN_SCALE_DEFAULT:
28692f76d53SAlp Dener         default:
28792f76d53SAlp Dener           break;
2888e231d97SPeter Brune       }
28992f76d53SAlp Dener       break;
29092f76d53SAlp Dener   }
2919566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(snes->vec_sol, &n));
2929566063dSJacob Faibussowitsch   PetscCall(VecGetSize(snes->vec_sol, &N));
2939566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(qn->B, n, n, N, N));
2949566063dSJacob Faibussowitsch   PetscCall(MatSetUp(qn->B));
2959566063dSJacob Faibussowitsch   PetscCall(MatLMVMReset(qn->B, PETSC_TRUE));
2969566063dSJacob Faibussowitsch   PetscCall(MatLMVMSetHistorySize(qn->B, qn->m));
2979566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func));
2984b11644fSPeter Brune   PetscFunctionReturn(0);
2994b11644fSPeter Brune }
3004b11644fSPeter Brune 
3014b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
3024b11644fSPeter Brune {
3039f83bee8SJed Brown   SNES_QN        *qn;
3040adebc6cSBarry Smith 
3054b11644fSPeter Brune   PetscFunctionBegin;
3064b11644fSPeter Brune   if (snes->data) {
3079f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
3089566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&qn->B));
3094b11644fSPeter Brune   }
3104b11644fSPeter Brune   PetscFunctionReturn(0);
3114b11644fSPeter Brune }
3124b11644fSPeter Brune 
3134b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
3144b11644fSPeter Brune {
3154b11644fSPeter Brune   PetscFunctionBegin;
3169566063dSJacob Faibussowitsch   PetscCall(SNESReset_QN(snes));
3179566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3182e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",NULL));
3192e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",NULL));
3202e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",NULL));
3214b11644fSPeter Brune   PetscFunctionReturn(0);
3224b11644fSPeter Brune }
3234b11644fSPeter Brune 
324*dbbe0bcdSBarry Smith static PetscErrorCode SNESSetFromOptions_QN(SNES snes,PetscOptionItems *PetscOptionsObject)
3254b11644fSPeter Brune {
3264b11644fSPeter Brune 
3272150357eSBarry Smith   SNES_QN           *qn    = (SNES_QN*)snes->data;
32892f76d53SAlp Dener   PetscBool         flg;
329f1c6b773SPeter Brune   SNESLineSearch    linesearch;
3302150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
3312150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
3321efc8c45SPeter Brune   SNESQNType        qtype = qn->type;
3332150357eSBarry Smith 
3344b11644fSPeter Brune   PetscFunctionBegin;
335d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"SNES QN options");
3369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL));
3379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL));
3389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", qn->monflg, &qn->monflg, NULL));
3399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg));
3409566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetScaleType(snes,stype));
34188f769c5SPeter Brune 
3429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg));
3439566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetRestartType(snes,rtype));
34488f769c5SPeter Brune 
3459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg));
3469566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESQNSetType(snes,qtype));
3479566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(qn->B));
348d0609cedSBarry Smith   PetscOptionsHeadEnd();
3499e764e56SPeter Brune   if (!snes->linesearch) {
3509566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
351ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
35260c08b40SPeter Brune       if (qn->type == SNES_QN_LBFGS) {
3539566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
35460c08b40SPeter Brune       } else if (qn->type == SNES_QN_BROYDEN) {
3559566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
35660c08b40SPeter Brune       } else {
3579566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
35860c08b40SPeter Brune       }
3599e764e56SPeter Brune     }
360ec786807SJed Brown   }
36192f76d53SAlp Dener   if (qn->monflg) {
3629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor));
36344f7e39eSPeter Brune   }
3644b11644fSPeter Brune   PetscFunctionReturn(0);
3654b11644fSPeter Brune }
3664b11644fSPeter Brune 
3675cd83419SPeter Brune static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
3685cd83419SPeter Brune {
3695cd83419SPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
3705cd83419SPeter Brune   PetscBool      iascii;
3715cd83419SPeter Brune 
3725cd83419SPeter Brune   PetscFunctionBegin;
3739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
3745cd83419SPeter Brune   if (iascii) {
3759566063dSJacob 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]));
37663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Stored subspace size: %" PetscInt_FMT "\n", qn->m));
3775cd83419SPeter Brune   }
3785cd83419SPeter Brune   PetscFunctionReturn(0);
3795cd83419SPeter Brune }
38092c02d66SPeter Brune 
3810c777b0cSPeter Brune /*@
3820c777b0cSPeter Brune     SNESQNSetRestartType - Sets the restart type for SNESQN.
3830c777b0cSPeter Brune 
3840c777b0cSPeter Brune     Logically Collective on SNES
3850c777b0cSPeter Brune 
3860c777b0cSPeter Brune     Input Parameters:
3870c777b0cSPeter Brune +   snes - the iterative context
3880c777b0cSPeter Brune -   rtype - restart type
3890c777b0cSPeter Brune 
3900c777b0cSPeter Brune     Options Database:
3910c777b0cSPeter Brune +   -snes_qn_restart_type <powell,periodic,none> - set the restart type
39284c577daSBarry Smith -   -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
3930c777b0cSPeter Brune 
3940c777b0cSPeter Brune     Level: intermediate
3950c777b0cSPeter Brune 
3960c777b0cSPeter Brune     SNESQNRestartTypes:
3970c777b0cSPeter Brune +   SNES_QN_RESTART_NONE - never restart
3980c777b0cSPeter Brune .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
3990c777b0cSPeter Brune -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
4000c777b0cSPeter Brune 
4010c777b0cSPeter Brune @*/
4022150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
4032150357eSBarry Smith {
4040c777b0cSPeter Brune   PetscFunctionBegin;
4050c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
406cac4c232SBarry Smith   PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));
4070c777b0cSPeter Brune   PetscFunctionReturn(0);
4080c777b0cSPeter Brune }
4090c777b0cSPeter Brune 
4100c777b0cSPeter Brune /*@
41184c577daSBarry Smith     SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN.
4120c777b0cSPeter Brune 
4130c777b0cSPeter Brune     Logically Collective on SNES
4140c777b0cSPeter Brune 
4150c777b0cSPeter Brune     Input Parameters:
4160c777b0cSPeter Brune +   snes - the iterative context
4170c777b0cSPeter Brune -   stype - scale type
4180c777b0cSPeter Brune 
4190c777b0cSPeter Brune     Options Database:
42067b8a455SSatish Balay .   -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
4210c777b0cSPeter Brune 
4220c777b0cSPeter Brune     Level: intermediate
4230c777b0cSPeter Brune 
42484c577daSBarry Smith     SNESQNScaleTypes:
4250c777b0cSPeter Brune +   SNES_QN_SCALE_NONE - don't scale the problem
4267044b745SBarry Smith .   SNES_QN_SCALE_SCALAR - use Shanno scaling
42792f76d53SAlp Dener .   SNES_QN_SCALE_DIAGONAL - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available
428a01a0525SBarry Smith -   SNES_QN_SCALE_JACOBIAN - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration
429a01a0525SBarry Smith                              of QN and at ever restart.
4300c777b0cSPeter Brune 
431db781477SPatrick Sanan .seealso: `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()`
4320c777b0cSPeter Brune @*/
4330c777b0cSPeter Brune 
4342150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
4352150357eSBarry Smith {
4360c777b0cSPeter Brune   PetscFunctionBegin;
4370c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
438cac4c232SBarry Smith   PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));
4390c777b0cSPeter Brune   PetscFunctionReturn(0);
4400c777b0cSPeter Brune }
4410c777b0cSPeter Brune 
4420adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
4430adebc6cSBarry Smith {
4440c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
4456e111a19SKarl Rupp 
4460c777b0cSPeter Brune   PetscFunctionBegin;
4470c777b0cSPeter Brune   qn->scale_type = stype;
4487044b745SBarry Smith   if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE;
4490c777b0cSPeter Brune   PetscFunctionReturn(0);
4500c777b0cSPeter Brune }
4510c777b0cSPeter Brune 
4520adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
4530adebc6cSBarry Smith {
4540c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
4556e111a19SKarl Rupp 
4560c777b0cSPeter Brune   PetscFunctionBegin;
4570c777b0cSPeter Brune   qn->restart_type = rtype;
4580c777b0cSPeter Brune   PetscFunctionReturn(0);
4590c777b0cSPeter Brune }
4600c777b0cSPeter Brune 
4611efc8c45SPeter Brune /*@
4621efc8c45SPeter Brune     SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN.
4631efc8c45SPeter Brune 
4641efc8c45SPeter Brune     Logically Collective on SNES
4651efc8c45SPeter Brune 
4661efc8c45SPeter Brune     Input Parameters:
4671efc8c45SPeter Brune +   snes - the iterative context
4681efc8c45SPeter Brune -   qtype - variant type
4691efc8c45SPeter Brune 
4701efc8c45SPeter Brune     Options Database:
47167b8a455SSatish Balay .   -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
4721efc8c45SPeter Brune 
4731efc8c45SPeter Brune     Level: beginner
4741efc8c45SPeter Brune 
4751efc8c45SPeter Brune     SNESQNTypes:
4761efc8c45SPeter Brune +   SNES_QN_LBFGS - LBFGS variant
4771efc8c45SPeter Brune .   SNES_QN_BROYDEN - Broyden variant
4781efc8c45SPeter Brune -   SNES_QN_BADBROYDEN - Bad Broyden variant
4791efc8c45SPeter Brune 
4801efc8c45SPeter Brune @*/
4811efc8c45SPeter Brune 
4821efc8c45SPeter Brune PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
4831efc8c45SPeter Brune {
4841efc8c45SPeter Brune   PetscFunctionBegin;
4851efc8c45SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
486cac4c232SBarry Smith   PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));
4871efc8c45SPeter Brune   PetscFunctionReturn(0);
4881efc8c45SPeter Brune }
4891efc8c45SPeter Brune 
4901efc8c45SPeter Brune PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
4911efc8c45SPeter Brune {
4921efc8c45SPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
4931efc8c45SPeter Brune 
4941efc8c45SPeter Brune   PetscFunctionBegin;
4951efc8c45SPeter Brune   qn->type = qtype;
4961efc8c45SPeter Brune   PetscFunctionReturn(0);
4971efc8c45SPeter Brune }
4981efc8c45SPeter Brune 
4994b11644fSPeter Brune /* -------------------------------------------------------------------------- */
5004b11644fSPeter Brune /*MC
5014b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
5024b11644fSPeter Brune 
5036cc8130cSPeter Brune       Options Database:
5046cc8130cSPeter Brune 
50584c577daSBarry Smith +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
50684c577daSBarry Smith +     -snes_qn_restart_type <powell,periodic,none> - set the restart type
5076bdcc836SBarry Smith .     -snes_qn_powell_gamma - Angle condition for restart.
5081867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
50984c577daSBarry Smith .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
51092f76d53SAlp Dener .     -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian
51172835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
51284c577daSBarry Smith -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.
5136cc8130cSPeter Brune 
51495452b02SPatrick Sanan       Notes:
51595452b02SPatrick Sanan     This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
516b8840d6bSPeter Brune       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
517b8840d6bSPeter Brune       updates.
5186cc8130cSPeter Brune 
5191867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
5201867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
52184c577daSBarry Smith       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
52284c577daSBarry Smith       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
5231867fe5bSPeter Brune 
5242d547940SBarry Smith       Uses left nonlinear preconditioning by default.
5252d547940SBarry Smith 
5266cc8130cSPeter Brune       References:
527606c0280SSatish Balay +   * -   Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
528606c0280SSatish Balay .   * -   R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods,
5290a8e8854SPeter Brune       Technical Report, Northwestern University, June 1992.
530606c0280SSatish Balay .   * -   Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
5310a8e8854SPeter Brune       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
532606c0280SSatish Balay .   * -   Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
5334f02bc6aSBarry Smith        SIAM Review, 57(4), 2015
534606c0280SSatish Balay .   * -   Griewank, Andreas. "Broyden updating, the good and the bad!." Doc. Math (2012): 301-315.
535606c0280SSatish Balay .   * -   Gilbert, Jean Charles, and Claude Lemar{\'e}chal. "Some numerical experiments with variable-storage quasi-Newton algorithms."
53692f76d53SAlp Dener       Mathematical programming 45.1-3 (1989): 407-435.
537606c0280SSatish Balay -   * -   Dener A., Munson T. "Accelerating Limited-Memory Quasi-Newton Convergence for Large-Scale Optimization"
538110fc3b0SBarry Smith       Computational Science - ICCS 2019. ICCS 2019. Lecture Notes in Computer Science, vol 11538. Springer, Cham
5394b11644fSPeter Brune 
5404b11644fSPeter Brune       Level: beginner
5414b11644fSPeter Brune 
542db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`
5436cc8130cSPeter Brune 
5444b11644fSPeter Brune M*/
5458cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
5464b11644fSPeter Brune {
5479f83bee8SJed Brown   SNES_QN        *qn;
54892f76d53SAlp Dener   const char     *optionsprefix;
5494b11644fSPeter Brune 
5504b11644fSPeter Brune   PetscFunctionBegin;
5514b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
5524b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
5534b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
5544b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
5555cd83419SPeter Brune   snes->ops->view           = SNESView_QN;
5564b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
5574b11644fSPeter Brune 
558efd4aadfSBarry Smith   snes->npcside= PC_LEFT;
55947f26062SPeter Brune 
560efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
56142f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
56242f4f86dSBarry Smith 
5634fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
5644fc747eaSLawrence Mitchell 
56588976e71SPeter Brune   if (!snes->tolerancesset) {
5660e444f03SPeter Brune     snes->max_funcs = 30000;
5670e444f03SPeter Brune     snes->max_its   = 10000;
56888976e71SPeter Brune   }
5690e444f03SPeter Brune 
5709566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes,&qn));
5714b11644fSPeter Brune   snes->data          = (void*) qn;
5720ecc9e1dSPeter Brune   qn->m               = 10;
573b21d5a53SPeter Brune   qn->scaling         = 1.0;
5740298fd71SBarry Smith   qn->monitor         = NULL;
57592f76d53SAlp Dener   qn->monflg          = PETSC_FALSE;
576b8840d6bSPeter Brune   qn->powell_gamma    = 0.9999;
57760c08b40SPeter Brune   qn->scale_type      = SNES_QN_SCALE_DEFAULT;
57860c08b40SPeter Brune   qn->restart_type    = SNES_QN_RESTART_DEFAULT;
579b8840d6bSPeter Brune   qn->type            = SNES_QN_LBFGS;
580ea630c6eSPeter Brune 
5819566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B));
5829566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
5839566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix));
58492f76d53SAlp Dener 
5859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN));
5869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN));
5879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN));
5884b11644fSPeter Brune   PetscFunctionReturn(0);
5894b11644fSPeter Brune }
590