xref: /petsc/src/snes/impls/qn/qn.c (revision 0ecc9e1ddfa6d6389c5c5f547ce092d58e8ad8b9)
14b11644fSPeter Brune #include <private/snesimpl.h>
24b11644fSPeter Brune 
388d374b2SPeter Brune typedef enum {SNES_QN_SEQUENTIAL, SNES_QN_COMPOSED} SNESQNCompositionType;
4*0ecc9e1dSPeter Brune typedef enum {SNES_QN_SHANNOSCALE, SNES_QN_LSSCALE, SNES_QN_JACOBIANSCALE} SNESQNScalingType;
56bf1b2e5SPeter Brune 
64b11644fSPeter Brune typedef struct {
74b11644fSPeter Brune   Vec         *dX;              /* The change in X */
86bf1b2e5SPeter Brune   Vec         *dF;              /* The change in F */
94b11644fSPeter Brune   PetscInt    m;                /* the number of kept previous steps */
10bd052dfeSPeter Brune   PetscScalar *alpha;
11bd052dfeSPeter Brune   PetscScalar *beta;
12bd052dfeSPeter Brune   PetscScalar *rho;
1344f7e39eSPeter Brune   PetscViewer monitor;
146bf1b2e5SPeter Brune   PetscReal   powell_gamma;     /* Powell angle restart condition */
156bf1b2e5SPeter Brune   PetscReal   powell_downhill;  /* Powell descent restart condition */
16b21d5a53SPeter Brune   PetscReal   scaling;          /* scaling of H0 */
17*0ecc9e1dSPeter Brune   PetscInt    n_restart;        /* the maximum iterations between restart */
1888d374b2SPeter Brune 
1988d374b2SPeter Brune   SNESQNCompositionType compositiontype; /* determine if the composition is done sequentially or as a composition */
20*0ecc9e1dSPeter Brune   SNESQNScalingType scalingtype; /* determine if the composition is done sequentially or as a composition */
2188d374b2SPeter Brune 
229f83bee8SJed Brown } SNES_QN;
234b11644fSPeter Brune 
244b11644fSPeter Brune #undef __FUNCT__
254b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private"
263af51624SPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) {
274b11644fSPeter Brune 
284b11644fSPeter Brune   PetscErrorCode ierr;
294b11644fSPeter Brune 
309f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*)snes->data;
314b11644fSPeter Brune 
32*0ecc9e1dSPeter Brune   Vec Yin = snes->work[0];
33*0ecc9e1dSPeter Brune 
344b11644fSPeter Brune   Vec *dX = qn->dX;
356bf1b2e5SPeter Brune   Vec *dF = qn->dF;
364b11644fSPeter Brune 
37bd052dfeSPeter Brune   PetscScalar *alpha = qn->alpha;
38bd052dfeSPeter Brune   PetscScalar *beta = qn->beta;
39bd052dfeSPeter Brune   PetscScalar *rho = qn->rho;
40bd052dfeSPeter Brune 
41*0ecc9e1dSPeter Brune   /* ksp thing for jacobian scaling */
42*0ecc9e1dSPeter Brune   KSPConvergedReason kspreason;
43*0ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
44*0ecc9e1dSPeter Brune 
45*0ecc9e1dSPeter Brune   PetscInt k, i, lits;
464b11644fSPeter Brune   PetscInt m = qn->m;
474b11644fSPeter Brune   PetscScalar t;
484b11644fSPeter Brune   PetscInt l = m;
494b11644fSPeter Brune 
504b11644fSPeter Brune   PetscFunctionBegin;
514b11644fSPeter Brune 
523af51624SPeter Brune   ierr = VecCopy(D, Y);CHKERRQ(ierr);
534b11644fSPeter Brune 
545ba6227bSPeter Brune   if (it < m) l = it;
554b11644fSPeter Brune 
564b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
574b11644fSPeter Brune   for (i = 0; i < l; i++) {
58b21d5a53SPeter Brune     k = (it - i - 1) % l;
593af51624SPeter Brune     ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr);
60bd052dfeSPeter Brune     alpha[k] = t*rho[k];
6144f7e39eSPeter Brune     if (qn->monitor) {
623af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
635ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
643af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
6544f7e39eSPeter Brune     }
666bf1b2e5SPeter Brune     ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr);
674b11644fSPeter Brune   }
684b11644fSPeter Brune 
69*0ecc9e1dSPeter Brune   if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
70*0ecc9e1dSPeter Brune     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
71*0ecc9e1dSPeter Brune     ierr = SNES_KSPSolve(snes,snes->ksp,Y,Yin);CHKERRQ(ierr);
72*0ecc9e1dSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
73*0ecc9e1dSPeter Brune     if (kspreason < 0) {
74*0ecc9e1dSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
75*0ecc9e1dSPeter Brune         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
76*0ecc9e1dSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
77*0ecc9e1dSPeter Brune         PetscFunctionReturn(0);
78*0ecc9e1dSPeter Brune       }
79*0ecc9e1dSPeter Brune     }
80*0ecc9e1dSPeter Brune     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
81*0ecc9e1dSPeter Brune     snes->linear_its += lits;
82*0ecc9e1dSPeter Brune     ierr = VecCopy(Yin, Y);CHKERRQ(ierr);
83*0ecc9e1dSPeter Brune   } else {
84b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
85*0ecc9e1dSPeter Brune   }
864b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
87bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
88b21d5a53SPeter Brune     k = (it + i - l) % l;
896bf1b2e5SPeter Brune     ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
904b11644fSPeter Brune     beta[k] = rho[k]*t;
913af51624SPeter Brune     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
9244f7e39eSPeter Brune     if (qn->monitor) {
933af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
945ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
953af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
9644f7e39eSPeter Brune     }
974b11644fSPeter Brune   }
984b11644fSPeter Brune   PetscFunctionReturn(0);
994b11644fSPeter Brune }
1004b11644fSPeter Brune 
1014b11644fSPeter Brune #undef __FUNCT__
1024b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
1034b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
1044b11644fSPeter Brune {
1054b11644fSPeter Brune 
1064b11644fSPeter Brune   PetscErrorCode ierr;
1079f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*) snes->data;
1084b11644fSPeter Brune 
10915f5eeeaSPeter Brune   Vec X, Xold;
1103af51624SPeter Brune   Vec F, G, B;
11188d374b2SPeter Brune   Vec W, Y, FPC, D, Dold;
1123af51624SPeter Brune   SNESConvergedReason reason;
1135ba6227bSPeter Brune   PetscInt i, i_r, k;
1144b11644fSPeter Brune 
1159f3a0142SPeter Brune   PetscReal fnorm, xnorm = 0, ynorm, gnorm;
1164b11644fSPeter Brune   PetscInt m = qn->m;
11705b53524SPeter Brune   PetscBool lssucceed, changed;
1184b11644fSPeter Brune 
119bd052dfeSPeter Brune   PetscScalar rhosc;
120bd052dfeSPeter Brune 
121bd052dfeSPeter Brune   Vec *dX = qn->dX;
1226bf1b2e5SPeter Brune   Vec *dF = qn->dF;
123bd052dfeSPeter Brune   PetscScalar *rho = qn->rho;
12488d374b2SPeter Brune   PetscScalar DolddotD, DolddotDold, DdotD, YdotD, a;
1254b11644fSPeter Brune 
126*0ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
127*0ecc9e1dSPeter Brune 
1284b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
1294b11644fSPeter Brune   PetscFunctionBegin;
1304b11644fSPeter Brune 
1319f3a0142SPeter Brune   X             = snes->vec_sol;        /* solution vector */
1329f3a0142SPeter Brune   F             = snes->vec_func;       /* residual vector */
1333af51624SPeter Brune   Y             = snes->vec_sol_update; /* search direction generated by J^-1D*/
1343af51624SPeter Brune   B             = snes->vec_rhs;
13570d3b23bSPeter Brune   G             = snes->work[0];
13670d3b23bSPeter Brune   W             = snes->work[1];
13770d3b23bSPeter Brune   Xold          = snes->work[2];
1383af51624SPeter Brune 
1393af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
14088d374b2SPeter Brune   D             = snes->work[3];
14188d374b2SPeter Brune   Dold          = snes->work[4];
1424b11644fSPeter Brune 
1434b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1444b11644fSPeter Brune 
1454b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1464b11644fSPeter Brune   snes->iter = 0;
1474b11644fSPeter Brune   snes->norm = 0.;
1484b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
14915f5eeeaSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1504b11644fSPeter Brune   if (snes->domainerror) {
1514b11644fSPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1524b11644fSPeter Brune     PetscFunctionReturn(0);
1534b11644fSPeter Brune   }
15415f5eeeaSPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1554b11644fSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1564b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1574b11644fSPeter Brune   snes->norm = fnorm;
1584b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1594b11644fSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1604b11644fSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1614b11644fSPeter Brune 
1624b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
1634b11644fSPeter Brune    snes->ttol = fnorm*snes->rtol;
1644b11644fSPeter Brune   /* test convergence */
1654b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1664b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
16770d3b23bSPeter Brune 
16888d374b2SPeter Brune   /* composed solve -- either sequential or composed */
16988d374b2SPeter Brune   if (snes->pc) {
17088d374b2SPeter Brune     if (qn->compositiontype == SNES_QN_SEQUENTIAL) {
17188d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
17288d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
17388d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
17488d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
17588d374b2SPeter Brune         PetscFunctionReturn(0);
17688d374b2SPeter Brune       }
17788d374b2SPeter Brune       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
17888d374b2SPeter Brune       ierr = VecCopy(FPC, F);CHKERRQ(ierr);
17988d374b2SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
1806bf1b2e5SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
18188d374b2SPeter Brune     } else {
18288d374b2SPeter Brune       ierr = VecCopy(X, Y);CHKERRQ(ierr);
18388d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, Y);CHKERRQ(ierr);
18488d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
18588d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
18688d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
18788d374b2SPeter Brune         PetscFunctionReturn(0);
18888d374b2SPeter Brune       }
18988d374b2SPeter Brune       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
19088d374b2SPeter Brune     }
19188d374b2SPeter Brune   } else {
19288d374b2SPeter Brune     ierr = VecCopy(F, Y);CHKERRQ(ierr);
19388d374b2SPeter Brune   }
19488d374b2SPeter Brune   ierr = VecCopy(Y, D);CHKERRQ(ierr);
1953af51624SPeter Brune 
196*0ecc9e1dSPeter Brune   if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
197*0ecc9e1dSPeter Brune     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
198*0ecc9e1dSPeter Brune   }
199*0ecc9e1dSPeter Brune 
2005ba6227bSPeter Brune   for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
20170d3b23bSPeter Brune     /* line search for lambda */
20270d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
20388d374b2SPeter Brune     ierr = VecCopy(D, Dold);CHKERRQ(ierr);
20415f5eeeaSPeter Brune     ierr = VecCopy(X, Xold);CHKERRQ(ierr);
20505b53524SPeter Brune     ierr = SNESLineSearchPreCheckApply(snes,X,Y,&changed);CHKERRQ(ierr);
20605b53524SPeter Brune     ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
2079f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2089f3a0142SPeter Brune     if (snes->domainerror) {
2099f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2109f3a0142SPeter Brune       PetscFunctionReturn(0);
2119f3a0142SPeter Brune       }
2129f3a0142SPeter Brune     if (!lssucceed) {
2139f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
2149f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2159f3a0142SPeter Brune         break;
2169f3a0142SPeter Brune       }
2179f3a0142SPeter Brune     }
218*0ecc9e1dSPeter Brune     if (qn->scalingtype == SNES_QN_LSSCALE) {
219*0ecc9e1dSPeter Brune       qn->scaling = ynorm;
220*0ecc9e1dSPeter Brune     }
2219f3a0142SPeter Brune     /* Update function and solution vectors */
2229f3a0142SPeter Brune     fnorm = gnorm;
2239f3a0142SPeter Brune     ierr = VecCopy(G,F);CHKERRQ(ierr);
22415f5eeeaSPeter Brune     ierr = VecCopy(W,X);CHKERRQ(ierr);
2253af51624SPeter Brune 
22688d374b2SPeter Brune     /* convergence monitoring */
227b21d5a53SPeter Brune     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
228b21d5a53SPeter Brune 
229360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
230360c497dSPeter Brune     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
231360c497dSPeter Brune 
2328409ca45SMatthew G Knepley     SNESLogConvHistory(snes,snes->norm,snes->iter);
2338409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
2344b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
2355ba6227bSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2364b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
2374b11644fSPeter Brune 
23888d374b2SPeter Brune 
23988d374b2SPeter Brune     if (snes->pc) {
24088d374b2SPeter Brune       if (qn->compositiontype == SNES_QN_SEQUENTIAL) {
24188d374b2SPeter Brune         ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
24288d374b2SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
24388d374b2SPeter Brune         if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
24488d374b2SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
24588d374b2SPeter Brune           PetscFunctionReturn(0);
24688d374b2SPeter Brune         }
24788d374b2SPeter Brune         ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
24888d374b2SPeter Brune         ierr = VecCopy(FPC, F);CHKERRQ(ierr);
24988d374b2SPeter Brune         ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
25088d374b2SPeter Brune         ierr = VecCopy(F, D);CHKERRQ(ierr);
25188d374b2SPeter Brune       } else {
25288d374b2SPeter Brune         ierr = VecCopy(X, D);CHKERRQ(ierr);
25388d374b2SPeter Brune         ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr);
25488d374b2SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
25588d374b2SPeter Brune         if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
25688d374b2SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
25788d374b2SPeter Brune           PetscFunctionReturn(0);
25888d374b2SPeter Brune         }
25988d374b2SPeter Brune         ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr);
26088d374b2SPeter Brune       }
26188d374b2SPeter Brune     } else {
26288d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
26388d374b2SPeter Brune     }
26488d374b2SPeter Brune 
2656bf1b2e5SPeter Brune     /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
26688d374b2SPeter Brune     ierr = VecDot(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
26788d374b2SPeter Brune     ierr = VecDot(Dold, D, &DolddotD);CHKERRQ(ierr);
26888d374b2SPeter Brune     ierr = VecDot(D, D, &DdotD);CHKERRQ(ierr);
26988d374b2SPeter Brune     ierr = VecDot(Y, D, &YdotD);CHKERRQ(ierr);
270*0ecc9e1dSPeter Brune     if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold)) || (i_r > qn->n_restart - 1 && qn->n_restart > 0)) {
2715ba6227bSPeter Brune       if (qn->monitor) {
2725ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
273*0ecc9e1dSPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", k, PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);CHKERRQ(ierr);
2745ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2755ba6227bSPeter Brune       }
2765ba6227bSPeter Brune       i_r = -1;
2775ba6227bSPeter Brune       /* general purpose update */
2785ba6227bSPeter Brune       if (snes->ops->update) {
2795ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2805ba6227bSPeter Brune       }
281*0ecc9e1dSPeter Brune       if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
282*0ecc9e1dSPeter Brune         ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
283*0ecc9e1dSPeter Brune       }
28488d374b2SPeter Brune       ierr = VecCopy(D, Y);CHKERRQ(ierr);
2855ba6227bSPeter Brune     } else {
286bd052dfeSPeter Brune       /* set the differences */
2875ba6227bSPeter Brune       k = i_r % m;
28888d374b2SPeter Brune       ierr = VecCopy(D, dF[k]);CHKERRQ(ierr);
28988d374b2SPeter Brune       ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
29015f5eeeaSPeter Brune       ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
29115f5eeeaSPeter Brune       ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
2926bf1b2e5SPeter Brune       ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr);
2936bf1b2e5SPeter Brune 
2946bf1b2e5SPeter Brune       /* set scaling to be shanno scaling */
295bd052dfeSPeter Brune       rho[k] = 1. / rhosc;
296*0ecc9e1dSPeter Brune       if (qn->scalingtype == SNES_QN_SHANNOSCALE) {
2976bf1b2e5SPeter Brune         ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr);
2981b2f85ebSPeter Brune         qn->scaling = PetscRealPart(rhosc) / PetscRealPart(a);
299*0ecc9e1dSPeter Brune       }
30070d3b23bSPeter Brune       /* general purpose update */
30170d3b23bSPeter Brune       if (snes->ops->update) {
30270d3b23bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
30370d3b23bSPeter Brune       }
30470d3b23bSPeter Brune       /* apply the current iteration of the approximate jacobian in order to get the next search direction*/
30588d374b2SPeter Brune       ierr = LBGFSApplyJinv_Private(snes, i_r+1, D, Y);CHKERRQ(ierr);
3065ba6227bSPeter Brune     }
3074b11644fSPeter Brune }
3084b11644fSPeter Brune   if (i == snes->max_its) {
3094b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
3104b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
3114b11644fSPeter Brune   }
3124b11644fSPeter Brune   PetscFunctionReturn(0);
3134b11644fSPeter Brune }
3144b11644fSPeter Brune 
3154b11644fSPeter Brune 
3164b11644fSPeter Brune #undef __FUNCT__
3174b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
3184b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
3194b11644fSPeter Brune {
3209f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*)snes->data;
3214b11644fSPeter Brune   PetscErrorCode ierr;
3224b11644fSPeter Brune   PetscFunctionBegin;
3234b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
3246bf1b2e5SPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr);
325bd052dfeSPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr);
32688d374b2SPeter Brune   ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr);
327*0ecc9e1dSPeter Brune 
3284b11644fSPeter Brune   PetscFunctionReturn(0);
3294b11644fSPeter Brune }
3304b11644fSPeter Brune 
3314b11644fSPeter Brune #undef __FUNCT__
3324b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
3334b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
3344b11644fSPeter Brune {
3354b11644fSPeter Brune   PetscErrorCode ierr;
3369f83bee8SJed Brown   SNES_QN *qn;
3374b11644fSPeter Brune   PetscFunctionBegin;
3384b11644fSPeter Brune   if (snes->data) {
3399f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
3404b11644fSPeter Brune     if (qn->dX) {
3414b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
3424b11644fSPeter Brune     }
3436bf1b2e5SPeter Brune     if (qn->dF) {
3446bf1b2e5SPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr);
3454b11644fSPeter Brune     }
346bd052dfeSPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr);
3474b11644fSPeter Brune   }
3484b11644fSPeter Brune   PetscFunctionReturn(0);
3494b11644fSPeter Brune }
3504b11644fSPeter Brune 
3514b11644fSPeter Brune #undef __FUNCT__
3524b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
3534b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
3544b11644fSPeter Brune {
3554b11644fSPeter Brune   PetscErrorCode ierr;
3564b11644fSPeter Brune   PetscFunctionBegin;
3574b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
3584b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3599c05af5eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr);
3604b11644fSPeter Brune   PetscFunctionReturn(0);
3614b11644fSPeter Brune }
3624b11644fSPeter Brune 
3634b11644fSPeter Brune #undef __FUNCT__
3644b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
3654b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
3664b11644fSPeter Brune {
3674b11644fSPeter Brune 
3684b11644fSPeter Brune   PetscErrorCode ierr;
3699f83bee8SJed Brown   SNES_QN    *qn;
3705b4627e8SPeter Brune   const char *compositions[] = {"sequential", "composed"};
371*0ecc9e1dSPeter Brune   const char *scalings[]     = {"shanno", "ls", "jacobian"};
37288d374b2SPeter Brune   PetscInt   indx = 0;
37388d374b2SPeter Brune   PetscBool  flg;
37444f7e39eSPeter Brune   PetscBool  monflg = PETSC_FALSE;
3754b11644fSPeter Brune   PetscFunctionBegin;
3764b11644fSPeter Brune 
3779f83bee8SJed Brown   qn = (SNES_QN*)snes->data;
3784b11644fSPeter Brune 
3794b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
3805b4627e8SPeter Brune   ierr = PetscOptionsInt("-snes_qn_m",                "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
381*0ecc9e1dSPeter Brune   ierr = PetscOptionsInt("-snes_qn_restart",                "Maximum number of iterations between restarts", "SNESQN", qn->n_restart, &qn->n_restart, PETSC_NULL);CHKERRQ(ierr);
3825b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_gamma",    "Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr);
3835b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_downhill", "Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr);
3845b4627e8SPeter Brune   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr);
3855b4627e8SPeter Brune   ierr = PetscOptionsEList("-snes_qn_composition",    "Composition type",                "SNESQN",compositions,2,"sequential",&indx,&flg);CHKERRQ(ierr);
38688d374b2SPeter Brune   if (flg) {
38788d374b2SPeter Brune     switch (indx) {
3885b4627e8SPeter Brune     case 0: qn->compositiontype = SNES_QN_SEQUENTIAL;
3895b4627e8SPeter Brune       break;
3905b4627e8SPeter Brune     case 1: qn->compositiontype = SNES_QN_COMPOSED;
3915b4627e8SPeter Brune       break;
39288d374b2SPeter Brune     }
39388d374b2SPeter Brune   }
394*0ecc9e1dSPeter Brune   ierr = PetscOptionsEList("-snes_qn_scaling",        "Scaling type",                    "SNESQN",scalings,3,"shanno",&indx,&flg);CHKERRQ(ierr);
395*0ecc9e1dSPeter Brune   if (flg) {
396*0ecc9e1dSPeter Brune     switch (indx) {
397*0ecc9e1dSPeter Brune     case 0: qn->scalingtype = SNES_QN_SHANNOSCALE;
398*0ecc9e1dSPeter Brune       break;
399*0ecc9e1dSPeter Brune     case 1: qn->scalingtype = SNES_QN_LSSCALE;
400*0ecc9e1dSPeter Brune       break;
401*0ecc9e1dSPeter Brune     case 2: qn->scalingtype = SNES_QN_JACOBIANSCALE;
402*0ecc9e1dSPeter Brune       snes->usesksp = PETSC_TRUE;
403*0ecc9e1dSPeter Brune       break;
404*0ecc9e1dSPeter Brune     }
405*0ecc9e1dSPeter Brune   }
406*0ecc9e1dSPeter Brune 
4074b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
40844f7e39eSPeter Brune   if (monflg) {
40944f7e39eSPeter Brune     qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
41044f7e39eSPeter Brune   }
4114b11644fSPeter Brune   PetscFunctionReturn(0);
4124b11644fSPeter Brune }
4134b11644fSPeter Brune 
41492c02d66SPeter Brune EXTERN_C_BEGIN
41592c02d66SPeter Brune #undef __FUNCT__
41692c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_QN"
41792c02d66SPeter Brune PetscErrorCode  SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type)
41892c02d66SPeter Brune {
41992c02d66SPeter Brune   PetscErrorCode ierr;
42092c02d66SPeter Brune   PetscFunctionBegin;
42192c02d66SPeter Brune 
42292c02d66SPeter Brune   switch (type) {
42392c02d66SPeter Brune   case SNES_LS_BASIC:
42492c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
42592c02d66SPeter Brune     break;
42692c02d66SPeter Brune   case SNES_LS_BASIC_NONORMS:
42792c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
42892c02d66SPeter Brune     break;
42992c02d66SPeter Brune   case SNES_LS_QUADRATIC:
430af60355fSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
43192c02d66SPeter Brune     break;
432f3aaf736SPeter Brune   case SNES_LS_CRITICAL:
433f3aaf736SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchCriticalSecant,PETSC_NULL);CHKERRQ(ierr);
434cf9edfecSPeter Brune     break;
43592c02d66SPeter Brune   default:
43692c02d66SPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
43792c02d66SPeter Brune     break;
43892c02d66SPeter Brune   }
43992c02d66SPeter Brune   snes->ls_type = type;
44092c02d66SPeter Brune   PetscFunctionReturn(0);
44192c02d66SPeter Brune }
44292c02d66SPeter Brune EXTERN_C_END
44392c02d66SPeter Brune 
44492c02d66SPeter Brune 
4454b11644fSPeter Brune /* -------------------------------------------------------------------------- */
4464b11644fSPeter Brune /*MC
4474b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
4484b11644fSPeter Brune 
4496cc8130cSPeter Brune       Options Database:
4506cc8130cSPeter Brune 
4516cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
4521867fe5bSPeter Brune .     -snes_qn_powell_angle - Angle condition for restart.
4531867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
4541867fe5bSPeter Brune .     -snes_qn_composition <sequential, composed>- Type of composition.
455f3aaf736SPeter Brune .     -snes_ls <basic, basicnonorms, quadratic, critical> - Type of line search.
45644f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
4576cc8130cSPeter Brune 
45844f7e39eSPeter Brune       Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to
4596cc8130cSPeter Brune       form the approximate inverse Jacobian using a series of multiplicative rank-one updates.  This will eventually be
4606cc8130cSPeter Brune       generalized to implement several limited-memory Broyden methods.
4616cc8130cSPeter Brune 
4621867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
4631867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
4641867fe5bSPeter Brune       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
4651867fe5bSPeter Brune       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
4661867fe5bSPeter Brune 
4676cc8130cSPeter Brune       References:
4686cc8130cSPeter Brune 
4696cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
4706cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
4716cc8130cSPeter Brune 
4724b11644fSPeter Brune 
4734b11644fSPeter Brune       Level: beginner
4744b11644fSPeter Brune 
4754b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
4766cc8130cSPeter Brune 
4774b11644fSPeter Brune M*/
4784b11644fSPeter Brune EXTERN_C_BEGIN
4794b11644fSPeter Brune #undef __FUNCT__
4804b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
4814b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
4824b11644fSPeter Brune {
4834b11644fSPeter Brune 
4844b11644fSPeter Brune   PetscErrorCode ierr;
4859f83bee8SJed Brown   SNES_QN *qn;
4864b11644fSPeter Brune 
4874b11644fSPeter Brune   PetscFunctionBegin;
4884b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
4894b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
4904b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
4914b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
4924b11644fSPeter Brune   snes->ops->view            = 0;
4934b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
4944b11644fSPeter Brune 
49542f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
49642f4f86dSBarry Smith   snes->usesksp         = PETSC_FALSE;
49742f4f86dSBarry Smith 
4980e444f03SPeter Brune   snes->max_funcs = 30000;
4990e444f03SPeter Brune   snes->max_its   = 10000;
5000e444f03SPeter Brune 
5019f83bee8SJed Brown   ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
5024b11644fSPeter Brune   snes->data = (void *) qn;
503*0ecc9e1dSPeter Brune   qn->m               = 10;
504b21d5a53SPeter Brune   qn->scaling         = 1.0;
5054b11644fSPeter Brune   qn->dX              = PETSC_NULL;
5066bf1b2e5SPeter Brune   qn->dF              = PETSC_NULL;
50744f7e39eSPeter Brune   qn->monitor         = PETSC_NULL;
50888d374b2SPeter Brune   qn->powell_gamma    = 0.9;
50988d374b2SPeter Brune   qn->powell_downhill = 0.2;
51088d374b2SPeter Brune   qn->compositiontype = SNES_QN_SEQUENTIAL;
511*0ecc9e1dSPeter Brune   qn->scalingtype     = SNES_QN_SHANNOSCALE;
512*0ecc9e1dSPeter Brune   qn->n_restart       = -1;
513ea630c6eSPeter Brune 
51492c02d66SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr);
515f3aaf736SPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_CRITICAL);CHKERRQ(ierr);
5169f3a0142SPeter Brune 
5174b11644fSPeter Brune   PetscFunctionReturn(0);
5184b11644fSPeter Brune }
5194b11644fSPeter Brune EXTERN_C_END
520