xref: /petsc/src/snes/impls/qn/qn.c (revision 88976e7148bb891b38778b9a5d5ae4c2b2cf44bf)
14b11644fSPeter Brune #include <private/snesimpl.h>
24b11644fSPeter Brune 
38e231d97SPeter Brune #define H(i,j)  qn->dXdFmat[i*qn->m + j]
48e231d97SPeter Brune 
588d374b2SPeter Brune typedef enum {SNES_QN_SEQUENTIAL, SNES_QN_COMPOSED} SNESQNCompositionType;
60ecc9e1dSPeter Brune typedef enum {SNES_QN_SHANNOSCALE, SNES_QN_LSSCALE, SNES_QN_JACOBIANSCALE} SNESQNScalingType;
76bf1b2e5SPeter Brune 
84b11644fSPeter Brune typedef struct {
94b11644fSPeter Brune   Vec          *dX;              /* The change in X */
106bf1b2e5SPeter Brune   Vec          *dF;              /* The change in F */
118e231d97SPeter Brune   PetscInt     m;                /* The number of kept previous steps */
128e231d97SPeter Brune   PetscScalar  *alpha, *beta;
138e231d97SPeter Brune   PetscScalar  *dXtdF, *dFtdX, *YtdX;
148e231d97SPeter Brune   PetscBool    aggreduct;        /* Aggregated reduction implementation */
158e231d97SPeter Brune   PetscScalar  *dXdFmat;         /* A matrix of values for dX_i dot dF_j */
1644f7e39eSPeter Brune   PetscViewer  monitor;
176bf1b2e5SPeter Brune   PetscReal    powell_gamma;     /* Powell angle restart condition */
186bf1b2e5SPeter Brune   PetscReal    powell_downhill;  /* Powell descent restart condition */
19b21d5a53SPeter Brune   PetscReal    scaling;          /* scaling of H0 */
200ecc9e1dSPeter Brune   PetscInt     n_restart;        /* the maximum iterations between restart */
2188d374b2SPeter Brune 
2288d374b2SPeter Brune   SNESQNCompositionType compositiontype; /* determine if the composition is done sequentially or as a composition */
230ecc9e1dSPeter Brune   SNESQNScalingType scalingtype; /* determine if the composition is done sequentially or as a composition */
2488d374b2SPeter Brune 
259f83bee8SJed Brown } SNES_QN;
264b11644fSPeter Brune 
274b11644fSPeter Brune #undef __FUNCT__
284b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private"
293af51624SPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) {
304b11644fSPeter Brune 
314b11644fSPeter Brune   PetscErrorCode ierr;
324b11644fSPeter Brune 
339f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*)snes->data;
344b11644fSPeter Brune 
35335fb975SPeter Brune   Vec Yin = snes->work[3];
360ecc9e1dSPeter Brune 
374b11644fSPeter Brune   Vec *dX = qn->dX;
386bf1b2e5SPeter Brune   Vec *dF = qn->dF;
394b11644fSPeter Brune 
40bd052dfeSPeter Brune   PetscScalar *alpha    = qn->alpha;
41bd052dfeSPeter Brune   PetscScalar *beta     = qn->beta;
428e231d97SPeter Brune   PetscScalar *dXtdF    = qn->dXtdF;
438e231d97SPeter Brune   PetscScalar *YtdX     = qn->YtdX;
44bd052dfeSPeter Brune 
450ecc9e1dSPeter Brune   /* ksp thing for jacobian scaling */
460ecc9e1dSPeter Brune   KSPConvergedReason kspreason;
470ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
480ecc9e1dSPeter Brune 
498e231d97SPeter Brune   PetscInt k, i, j, g, lits;
504b11644fSPeter Brune   PetscInt m = qn->m;
514b11644fSPeter Brune   PetscScalar t;
524b11644fSPeter Brune   PetscInt l = m;
534b11644fSPeter Brune 
548e231d97SPeter Brune   Mat jac, jac_pre;
558e231d97SPeter Brune 
564b11644fSPeter Brune   PetscFunctionBegin;
574b11644fSPeter Brune 
583af51624SPeter Brune   ierr = VecCopy(D, Y);CHKERRQ(ierr);
594b11644fSPeter Brune 
605ba6227bSPeter Brune   if (it < m) l = it;
614b11644fSPeter Brune 
628e231d97SPeter Brune   if (qn->aggreduct) {
638e231d97SPeter Brune     ierr = VecMDot(Y, l, qn->dX, YtdX);CHKERRQ(ierr);
648e231d97SPeter Brune   }
654b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
664b11644fSPeter Brune   for (i = 0; i < l; i++) {
67b21d5a53SPeter Brune     k = (it - i - 1) % l;
688e231d97SPeter Brune     if (qn->aggreduct) {
698e231d97SPeter Brune       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
708e231d97SPeter Brune       t = YtdX[k];
718e231d97SPeter Brune       for (j = 0; j < i; j++) {
728e231d97SPeter Brune         g = (it - j - 1) % l;
738e231d97SPeter Brune         t += -alpha[g]*H(g, k);
748e231d97SPeter Brune       }
758e231d97SPeter Brune       alpha[k] = t / H(k, k);
768e231d97SPeter Brune     } else {
773af51624SPeter Brune       ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr);
788e231d97SPeter Brune       alpha[k] = t / dXtdF[k];
798e231d97SPeter Brune     }
8044f7e39eSPeter Brune     if (qn->monitor) {
813af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
825ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
833af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
8444f7e39eSPeter Brune     }
856bf1b2e5SPeter Brune     ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr);
864b11644fSPeter Brune   }
874b11644fSPeter Brune 
880ecc9e1dSPeter Brune   if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
898e231d97SPeter Brune     ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
908e231d97SPeter Brune     ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr);
910ecc9e1dSPeter Brune     ierr = SNES_KSPSolve(snes,snes->ksp,Y,Yin);CHKERRQ(ierr);
920ecc9e1dSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
930ecc9e1dSPeter Brune     if (kspreason < 0) {
940ecc9e1dSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
950ecc9e1dSPeter 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);
960ecc9e1dSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
970ecc9e1dSPeter Brune         PetscFunctionReturn(0);
980ecc9e1dSPeter Brune       }
990ecc9e1dSPeter Brune     }
1000ecc9e1dSPeter Brune     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
1010ecc9e1dSPeter Brune     snes->linear_its += lits;
1020ecc9e1dSPeter Brune     ierr = VecCopy(Yin, Y);CHKERRQ(ierr);
1030ecc9e1dSPeter Brune   } else {
104b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
1050ecc9e1dSPeter Brune   }
1068e231d97SPeter Brune   if (qn->aggreduct) {
1078e231d97SPeter Brune     ierr = VecMDot(Y, l, qn->dF, YtdX);CHKERRQ(ierr);
1088e231d97SPeter Brune   }
1094b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
110bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
111b21d5a53SPeter Brune     k = (it + i - l) % l;
1128e231d97SPeter Brune     if (qn->aggreduct) {
1138e231d97SPeter Brune       t = YtdX[k];
1148e231d97SPeter Brune       for (j = 0; j < i; j++) {
1158e231d97SPeter Brune         g = (it + j - l) % l;
1168e231d97SPeter Brune         t += (alpha[g] - beta[g])*H(k, g);
1178e231d97SPeter Brune       }
1188e231d97SPeter Brune       beta[k] = t / H(k, k);
1198e231d97SPeter Brune     } else {
1206bf1b2e5SPeter Brune       ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
1218e231d97SPeter Brune       beta[k] = t / dXtdF[k];
1228e231d97SPeter Brune     }
1233af51624SPeter Brune     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
12444f7e39eSPeter Brune     if (qn->monitor) {
1253af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1265ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
1273af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
12844f7e39eSPeter Brune     }
1294b11644fSPeter Brune   }
1304b11644fSPeter Brune   PetscFunctionReturn(0);
1314b11644fSPeter Brune }
1324b11644fSPeter Brune 
1334b11644fSPeter Brune #undef __FUNCT__
1344b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
1354b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
1364b11644fSPeter Brune {
1374b11644fSPeter Brune 
1384b11644fSPeter Brune   PetscErrorCode ierr;
1399f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*) snes->data;
1404b11644fSPeter Brune 
14115f5eeeaSPeter Brune   Vec X, Xold;
142335fb975SPeter Brune   Vec F, B;
143335fb975SPeter Brune   Vec Y, FPC, D, Dold;
1443af51624SPeter Brune   SNESConvergedReason reason;
1458e231d97SPeter Brune   PetscInt i, i_r, k, l, j;
1464b11644fSPeter Brune 
147335fb975SPeter Brune   PetscReal fnorm, xnorm, ynorm, gnorm;
1484b11644fSPeter Brune   PetscInt m = qn->m;
149335fb975SPeter Brune   PetscBool lssucceed;
1504b11644fSPeter Brune 
151bd052dfeSPeter Brune   Vec *dX = qn->dX;
1526bf1b2e5SPeter Brune   Vec *dF = qn->dF;
1538e231d97SPeter Brune   PetscScalar *dXtdF = qn->dXtdF;
1548e231d97SPeter Brune   PetscScalar *dFtdX = qn->dFtdX;
15588d374b2SPeter Brune   PetscScalar DolddotD, DolddotDold, DdotD, YdotD, a;
1564b11644fSPeter Brune 
1570ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
1580ecc9e1dSPeter Brune 
1594b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
1604b11644fSPeter Brune   PetscFunctionBegin;
1614b11644fSPeter Brune 
1629f3a0142SPeter Brune   X             = snes->vec_sol;        /* solution vector */
1639f3a0142SPeter Brune   F             = snes->vec_func;       /* residual vector */
1643af51624SPeter Brune   Y             = snes->vec_sol_update; /* search direction generated by J^-1D*/
1653af51624SPeter Brune   B             = snes->vec_rhs;
166335fb975SPeter Brune   Xold          = snes->work[0];
1673af51624SPeter Brune 
1683af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
169335fb975SPeter Brune   D             = snes->work[1];
170335fb975SPeter Brune   Dold          = snes->work[2];
1714b11644fSPeter Brune 
1724b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1734b11644fSPeter Brune 
1744b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1754b11644fSPeter Brune   snes->iter = 0;
1764b11644fSPeter Brune   snes->norm = 0.;
1774b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
17815f5eeeaSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1794b11644fSPeter Brune   if (snes->domainerror) {
1804b11644fSPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1814b11644fSPeter Brune     PetscFunctionReturn(0);
1824b11644fSPeter Brune   }
18315f5eeeaSPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1844b11644fSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1854b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1864b11644fSPeter Brune   snes->norm = fnorm;
1874b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1884b11644fSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1894b11644fSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1904b11644fSPeter Brune 
1914b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
1924b11644fSPeter Brune    snes->ttol = fnorm*snes->rtol;
1934b11644fSPeter Brune   /* test convergence */
1944b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1954b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
19670d3b23bSPeter Brune 
19788d374b2SPeter Brune   /* composed solve -- either sequential or composed */
19888d374b2SPeter Brune   if (snes->pc) {
19988d374b2SPeter Brune     if (qn->compositiontype == SNES_QN_SEQUENTIAL) {
20088d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
20188d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
20288d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
20388d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
20488d374b2SPeter Brune         PetscFunctionReturn(0);
20588d374b2SPeter Brune       }
20688d374b2SPeter Brune       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
20788d374b2SPeter Brune       ierr = VecCopy(FPC, F);CHKERRQ(ierr);
20888d374b2SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
2096bf1b2e5SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
21088d374b2SPeter Brune     } else {
21188d374b2SPeter Brune       ierr = VecCopy(X, Y);CHKERRQ(ierr);
21288d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, Y);CHKERRQ(ierr);
21388d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
21488d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
21588d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
21688d374b2SPeter Brune         PetscFunctionReturn(0);
21788d374b2SPeter Brune       }
21888d374b2SPeter Brune       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
21988d374b2SPeter Brune     }
22088d374b2SPeter Brune   } else {
22188d374b2SPeter Brune     ierr = VecCopy(F, Y);CHKERRQ(ierr);
22288d374b2SPeter Brune   }
22388d374b2SPeter Brune   ierr = VecCopy(Y, D);CHKERRQ(ierr);
2243af51624SPeter Brune 
225f8e15203SPeter Brune   /* scale the initial update */
2260ecc9e1dSPeter Brune   if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
2270ecc9e1dSPeter Brune     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
2280ecc9e1dSPeter Brune   }
2290ecc9e1dSPeter Brune 
2305ba6227bSPeter Brune   for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
231f8e15203SPeter Brune     ierr = LBGFSApplyJinv_Private(snes, i_r, D, Y);CHKERRQ(ierr);
23270d3b23bSPeter Brune     /* line search for lambda */
23370d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
23488d374b2SPeter Brune     ierr = VecCopy(D, Dold);CHKERRQ(ierr);
23515f5eeeaSPeter Brune     ierr = VecCopy(X, Xold);CHKERRQ(ierr);
236f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
2379f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2389f3a0142SPeter Brune     if (snes->domainerror) {
2399f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2409f3a0142SPeter Brune       PetscFunctionReturn(0);
2419f3a0142SPeter Brune       }
242f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr);
2439f3a0142SPeter Brune     if (!lssucceed) {
2449f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
2459f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2469f3a0142SPeter Brune         break;
2479f3a0142SPeter Brune       }
2489f3a0142SPeter Brune     }
249f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
2500ecc9e1dSPeter Brune     if (qn->scalingtype == SNES_QN_LSSCALE) {
251f1c6b773SPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
2520ecc9e1dSPeter Brune     }
2533af51624SPeter Brune 
25488d374b2SPeter Brune     /* convergence monitoring */
255b21d5a53SPeter 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);
256b21d5a53SPeter Brune 
257360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
258360c497dSPeter Brune     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
259360c497dSPeter Brune 
2608409ca45SMatthew G Knepley     SNESLogConvHistory(snes,snes->norm,snes->iter);
2618409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
2624b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
2635ba6227bSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2644b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
2654b11644fSPeter Brune 
26688d374b2SPeter Brune 
26788d374b2SPeter Brune     if (snes->pc) {
26888d374b2SPeter Brune       if (qn->compositiontype == SNES_QN_SEQUENTIAL) {
26988d374b2SPeter Brune         ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
27088d374b2SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
27188d374b2SPeter Brune         if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
27288d374b2SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
27388d374b2SPeter Brune           PetscFunctionReturn(0);
27488d374b2SPeter Brune         }
27588d374b2SPeter Brune         ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
27688d374b2SPeter Brune         ierr = VecCopy(FPC, F);CHKERRQ(ierr);
27788d374b2SPeter Brune         ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
27888d374b2SPeter Brune         ierr = VecCopy(F, D);CHKERRQ(ierr);
27988d374b2SPeter Brune       } else {
28088d374b2SPeter Brune         ierr = VecCopy(X, D);CHKERRQ(ierr);
28188d374b2SPeter Brune         ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr);
28288d374b2SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
28388d374b2SPeter Brune         if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
28488d374b2SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
28588d374b2SPeter Brune           PetscFunctionReturn(0);
28688d374b2SPeter Brune         }
28788d374b2SPeter Brune         ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr);
28888d374b2SPeter Brune       }
28988d374b2SPeter Brune     } else {
29088d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
29188d374b2SPeter Brune     }
29288d374b2SPeter Brune 
2936bf1b2e5SPeter Brune     /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
2948e231d97SPeter Brune     ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
2958e231d97SPeter Brune     ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr);
2968e231d97SPeter Brune     ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr);
2978e231d97SPeter Brune     ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr);
2988e231d97SPeter Brune     ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
2998e231d97SPeter Brune     ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr);
3008e231d97SPeter Brune     ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr);
3018e231d97SPeter Brune     ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr);
3020ecc9e1dSPeter Brune     if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold)) || (i_r > qn->n_restart - 1 && qn->n_restart > 0)) {
3035ba6227bSPeter Brune       if (qn->monitor) {
3045ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
305516fe3c3SPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);CHKERRQ(ierr);
3065ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
3075ba6227bSPeter Brune       }
3085ba6227bSPeter Brune       i_r = -1;
3095ba6227bSPeter Brune       /* general purpose update */
3105ba6227bSPeter Brune       if (snes->ops->update) {
3115ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
3125ba6227bSPeter Brune       }
3130ecc9e1dSPeter Brune       if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
3140ecc9e1dSPeter Brune         ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
3150ecc9e1dSPeter Brune       }
3165ba6227bSPeter Brune     } else {
317bd052dfeSPeter Brune       /* set the differences */
3185ba6227bSPeter Brune       k = i_r % m;
3198e231d97SPeter Brune       l = m;
3208e231d97SPeter Brune       if (i_r + 1 < m) l = i_r + 1;
32188d374b2SPeter Brune       ierr = VecCopy(D, dF[k]);CHKERRQ(ierr);
32288d374b2SPeter Brune       ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
32315f5eeeaSPeter Brune       ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
32415f5eeeaSPeter Brune       ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
3258e231d97SPeter Brune       if (qn->aggreduct) {
3268e231d97SPeter Brune         ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr);
3278e231d97SPeter Brune         ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr);
3288e231d97SPeter Brune         for (j = 0; j < l; j++) {
3298e231d97SPeter Brune           H(k, j) = dFtdX[j];
3308e231d97SPeter Brune           H(j, k) = dXtdF[j];
3318e231d97SPeter Brune         }
3328e231d97SPeter Brune         /* copy back over to make the computation of alpha and beta easier */
3338e231d97SPeter Brune         for (j = 0; j < l; j++) {
3348e231d97SPeter Brune           dXtdF[j] = H(j, j);
3358e231d97SPeter Brune         }
3368e231d97SPeter Brune       } else {
3378e231d97SPeter Brune         ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr);
3388e231d97SPeter Brune       }
3396bf1b2e5SPeter Brune       /* set scaling to be shanno scaling */
3400ecc9e1dSPeter Brune       if (qn->scalingtype == SNES_QN_SHANNOSCALE) {
3416bf1b2e5SPeter Brune         ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr);
3428e231d97SPeter Brune         qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(a);
3430ecc9e1dSPeter Brune       }
34470d3b23bSPeter Brune       /* general purpose update */
34570d3b23bSPeter Brune       if (snes->ops->update) {
34670d3b23bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
34770d3b23bSPeter Brune       }
3485ba6227bSPeter Brune     }
3494b11644fSPeter Brune   }
3504b11644fSPeter Brune   if (i == snes->max_its) {
3514b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
3524b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
3534b11644fSPeter Brune   }
3544b11644fSPeter Brune   PetscFunctionReturn(0);
3554b11644fSPeter Brune }
3564b11644fSPeter Brune 
3574b11644fSPeter Brune 
3584b11644fSPeter Brune #undef __FUNCT__
3594b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
3604b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
3614b11644fSPeter Brune {
3629f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
3634b11644fSPeter Brune   PetscErrorCode ierr;
364335fb975SPeter Brune 
3654b11644fSPeter Brune   PetscFunctionBegin;
3664b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
3676bf1b2e5SPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr);
3688e231d97SPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha,
3698e231d97SPeter Brune                       qn->m, PetscScalar, &qn->beta,
3708e231d97SPeter Brune                       qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr);
3718e231d97SPeter Brune 
3728e231d97SPeter Brune   if (qn->aggreduct) {
3738e231d97SPeter Brune     ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat,
3748e231d97SPeter Brune                         qn->m, PetscScalar, &qn->dFtdX,
3758e231d97SPeter Brune                         qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr);
3768e231d97SPeter Brune   }
377335fb975SPeter Brune   ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr);
378335fb975SPeter Brune 
379335fb975SPeter Brune   /* set up the line search */
3808e231d97SPeter Brune   if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
3818e231d97SPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
3828e231d97SPeter Brune   }
3834b11644fSPeter Brune   PetscFunctionReturn(0);
3844b11644fSPeter Brune }
3854b11644fSPeter Brune 
3864b11644fSPeter Brune #undef __FUNCT__
3874b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
3884b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
3894b11644fSPeter Brune {
3904b11644fSPeter Brune   PetscErrorCode ierr;
3919f83bee8SJed Brown   SNES_QN *qn;
3924b11644fSPeter Brune   PetscFunctionBegin;
3934b11644fSPeter Brune   if (snes->data) {
3949f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
3954b11644fSPeter Brune     if (qn->dX) {
3964b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
3974b11644fSPeter Brune     }
3986bf1b2e5SPeter Brune     if (qn->dF) {
3996bf1b2e5SPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr);
4004b11644fSPeter Brune     }
4018e231d97SPeter Brune     if (qn->aggreduct) {
4028e231d97SPeter Brune       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
4038e231d97SPeter Brune     }
4048e231d97SPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr);
4054b11644fSPeter Brune   }
4064b11644fSPeter Brune   PetscFunctionReturn(0);
4074b11644fSPeter Brune }
4084b11644fSPeter Brune 
4094b11644fSPeter Brune #undef __FUNCT__
4104b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
4114b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
4124b11644fSPeter Brune {
4134b11644fSPeter Brune   PetscErrorCode ierr;
4144b11644fSPeter Brune   PetscFunctionBegin;
4154b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
4164b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
4179c05af5eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr);
4184b11644fSPeter Brune   PetscFunctionReturn(0);
4194b11644fSPeter Brune }
4204b11644fSPeter Brune 
4214b11644fSPeter Brune #undef __FUNCT__
4224b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
4234b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
4244b11644fSPeter Brune {
4254b11644fSPeter Brune 
4264b11644fSPeter Brune   PetscErrorCode ierr;
4279f83bee8SJed Brown   SNES_QN    *qn;
4285b4627e8SPeter Brune   const char *compositions[] = {"sequential", "composed"};
4290ecc9e1dSPeter Brune   const char *scalings[]     = {"shanno", "ls", "jacobian"};
43088d374b2SPeter Brune   PetscInt   indx = 0;
43188d374b2SPeter Brune   PetscBool  flg;
43244f7e39eSPeter Brune   PetscBool  monflg = PETSC_FALSE;
433f1c6b773SPeter Brune   SNESLineSearch linesearch;
4344b11644fSPeter Brune   PetscFunctionBegin;
4354b11644fSPeter Brune 
4369f83bee8SJed Brown   qn = (SNES_QN*)snes->data;
4374b11644fSPeter Brune 
4384b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
4395b4627e8SPeter Brune   ierr = PetscOptionsInt("-snes_qn_m",                "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
4400ecc9e1dSPeter Brune   ierr = PetscOptionsInt("-snes_qn_restart",                "Maximum number of iterations between restarts", "SNESQN", qn->n_restart, &qn->n_restart, PETSC_NULL);CHKERRQ(ierr);
4415b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_gamma",    "Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr);
4425b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_downhill", "Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr);
4435b4627e8SPeter Brune   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr);
4448e231d97SPeter Brune   ierr = PetscOptionsBool("-snes_qn_aggreduct",       "Aggregate reductions",            "SNESQN", qn->aggreduct, &qn->aggreduct, PETSC_NULL);CHKERRQ(ierr);
4455b4627e8SPeter Brune   ierr = PetscOptionsEList("-snes_qn_composition",    "Composition type",                "SNESQN",compositions,2,"sequential",&indx,&flg);CHKERRQ(ierr);
44688d374b2SPeter Brune   if (flg) {
44788d374b2SPeter Brune     switch (indx) {
4485b4627e8SPeter Brune     case 0: qn->compositiontype = SNES_QN_SEQUENTIAL;
4495b4627e8SPeter Brune       break;
4505b4627e8SPeter Brune     case 1: qn->compositiontype = SNES_QN_COMPOSED;
4515b4627e8SPeter Brune       break;
45288d374b2SPeter Brune     }
45388d374b2SPeter Brune   }
4540ecc9e1dSPeter Brune   ierr = PetscOptionsEList("-snes_qn_scaling",        "Scaling type",                    "SNESQN",scalings,3,"shanno",&indx,&flg);CHKERRQ(ierr);
4550ecc9e1dSPeter Brune   if (flg) {
4560ecc9e1dSPeter Brune     switch (indx) {
4570ecc9e1dSPeter Brune     case 0: qn->scalingtype = SNES_QN_SHANNOSCALE;
4580ecc9e1dSPeter Brune       break;
4590ecc9e1dSPeter Brune     case 1: qn->scalingtype = SNES_QN_LSSCALE;
4600ecc9e1dSPeter Brune       break;
4610ecc9e1dSPeter Brune     case 2: qn->scalingtype = SNES_QN_JACOBIANSCALE;
4620ecc9e1dSPeter Brune       snes->usesksp = PETSC_TRUE;
4630ecc9e1dSPeter Brune       break;
4640ecc9e1dSPeter Brune     }
4650ecc9e1dSPeter Brune   }
4660ecc9e1dSPeter Brune 
4674b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
4689e764e56SPeter Brune   if (!snes->linesearch) {
469f1c6b773SPeter Brune     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
4709e764e56SPeter Brune     if (!snes->pc || qn->compositiontype == SNES_QN_SEQUENTIAL) {
471f1c6b773SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_CP);CHKERRQ(ierr);
4729e764e56SPeter Brune     } else {
473f1c6b773SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr);
4749e764e56SPeter Brune     }
4759e764e56SPeter Brune   }
47644f7e39eSPeter Brune   if (monflg) {
47744f7e39eSPeter Brune     qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
47844f7e39eSPeter Brune   }
4794b11644fSPeter Brune   PetscFunctionReturn(0);
4804b11644fSPeter Brune }
4814b11644fSPeter Brune 
48292c02d66SPeter Brune 
4834b11644fSPeter Brune /* -------------------------------------------------------------------------- */
4844b11644fSPeter Brune /*MC
4854b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
4864b11644fSPeter Brune 
4876cc8130cSPeter Brune       Options Database:
4886cc8130cSPeter Brune 
4896cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
4901867fe5bSPeter Brune .     -snes_qn_powell_angle - Angle condition for restart.
4911867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
4921867fe5bSPeter Brune .     -snes_qn_composition <sequential, composed>- Type of composition.
493f3aaf736SPeter Brune .     -snes_ls <basic, basicnonorms, quadratic, critical> - Type of line search.
49444f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
4956cc8130cSPeter Brune 
49644f7e39eSPeter Brune       Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to
4976cc8130cSPeter Brune       form the approximate inverse Jacobian using a series of multiplicative rank-one updates.  This will eventually be
4986cc8130cSPeter Brune       generalized to implement several limited-memory Broyden methods.
4996cc8130cSPeter Brune 
5001867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
5011867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
5021867fe5bSPeter Brune       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
5031867fe5bSPeter Brune       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
5041867fe5bSPeter Brune 
5056cc8130cSPeter Brune       References:
5066cc8130cSPeter Brune 
5076cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
5086cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
5096cc8130cSPeter Brune 
5104b11644fSPeter Brune 
5114b11644fSPeter Brune       Level: beginner
5124b11644fSPeter Brune 
5134b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
5146cc8130cSPeter Brune 
5154b11644fSPeter Brune M*/
5164b11644fSPeter Brune EXTERN_C_BEGIN
5174b11644fSPeter Brune #undef __FUNCT__
5184b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
5194b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
5204b11644fSPeter Brune {
5214b11644fSPeter Brune 
5224b11644fSPeter Brune   PetscErrorCode ierr;
5239f83bee8SJed Brown   SNES_QN *qn;
5244b11644fSPeter Brune 
5254b11644fSPeter Brune   PetscFunctionBegin;
5264b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
5274b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
5284b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
5294b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
5304b11644fSPeter Brune   snes->ops->view            = 0;
5314b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
5324b11644fSPeter Brune 
53342f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
53442f4f86dSBarry Smith   snes->usesksp         = PETSC_FALSE;
53542f4f86dSBarry Smith 
536*88976e71SPeter Brune   if (!snes->tolerancesset) {
5370e444f03SPeter Brune     snes->max_funcs = 30000;
5380e444f03SPeter Brune     snes->max_its   = 10000;
539*88976e71SPeter Brune   }
5400e444f03SPeter Brune 
5419f83bee8SJed Brown   ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
5424b11644fSPeter Brune   snes->data = (void *) qn;
5430ecc9e1dSPeter Brune   qn->m               = 10;
544b21d5a53SPeter Brune   qn->scaling         = 1.0;
5454b11644fSPeter Brune   qn->dX              = PETSC_NULL;
5466bf1b2e5SPeter Brune   qn->dF              = PETSC_NULL;
5478e231d97SPeter Brune   qn->dXtdF           = PETSC_NULL;
5488e231d97SPeter Brune   qn->dFtdX           = PETSC_NULL;
5498e231d97SPeter Brune   qn->dXdFmat         = PETSC_NULL;
55044f7e39eSPeter Brune   qn->monitor         = PETSC_NULL;
5518e231d97SPeter Brune   qn->aggreduct       = PETSC_FALSE;
55288d374b2SPeter Brune   qn->powell_gamma    = 0.9;
55388d374b2SPeter Brune   qn->powell_downhill = 0.2;
55488d374b2SPeter Brune   qn->compositiontype = SNES_QN_SEQUENTIAL;
5550ecc9e1dSPeter Brune   qn->scalingtype     = SNES_QN_SHANNOSCALE;
5560ecc9e1dSPeter Brune   qn->n_restart       = -1;
557ea630c6eSPeter Brune 
5584b11644fSPeter Brune   PetscFunctionReturn(0);
5594b11644fSPeter Brune }
5608e231d97SPeter Brune 
5614b11644fSPeter Brune EXTERN_C_END
562