xref: /petsc/src/snes/impls/qn/qn.c (revision 335fb97548d572aa4df905bb6590c5360b52bc0f)
14b11644fSPeter Brune #include <private/snesimpl.h>
2*335fb975SPeter Brune #include <petsclinesearch.h>
34b11644fSPeter Brune 
488d374b2SPeter Brune typedef enum {SNES_QN_SEQUENTIAL, SNES_QN_COMPOSED} SNESQNCompositionType;
50ecc9e1dSPeter Brune typedef enum {SNES_QN_SHANNOSCALE, SNES_QN_LSSCALE, SNES_QN_JACOBIANSCALE} SNESQNScalingType;
66bf1b2e5SPeter Brune 
74b11644fSPeter Brune typedef struct {
84b11644fSPeter Brune   Vec          *dX;              /* The change in X */
96bf1b2e5SPeter Brune   Vec          *dF;              /* The change in F */
104b11644fSPeter Brune   PetscInt     m;                /* the number of kept previous steps */
11bd052dfeSPeter Brune   PetscScalar  *alpha;
12bd052dfeSPeter Brune   PetscScalar  *beta;
13bd052dfeSPeter Brune   PetscScalar  *rho;
1444f7e39eSPeter Brune   PetscViewer  monitor;
156bf1b2e5SPeter Brune   PetscReal    powell_gamma;     /* Powell angle restart condition */
166bf1b2e5SPeter Brune   PetscReal    powell_downhill;  /* Powell descent restart condition */
17b21d5a53SPeter Brune   PetscReal    scaling;          /* scaling of H0 */
180ecc9e1dSPeter Brune   PetscInt     n_restart;        /* the maximum iterations between restart */
1988d374b2SPeter Brune 
20*335fb975SPeter Brune   LineSearch   linesearch;       /* line search */
21*335fb975SPeter 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 
35*335fb975SPeter 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;
42bd052dfeSPeter Brune   PetscScalar *rho = qn->rho;
43bd052dfeSPeter Brune 
440ecc9e1dSPeter Brune   /* ksp thing for jacobian scaling */
450ecc9e1dSPeter Brune   KSPConvergedReason kspreason;
460ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
470ecc9e1dSPeter Brune 
480ecc9e1dSPeter Brune   PetscInt k, i, lits;
494b11644fSPeter Brune   PetscInt m = qn->m;
504b11644fSPeter Brune   PetscScalar t;
514b11644fSPeter Brune   PetscInt l = m;
524b11644fSPeter Brune 
534b11644fSPeter Brune   PetscFunctionBegin;
544b11644fSPeter Brune 
553af51624SPeter Brune   ierr = VecCopy(D, Y);CHKERRQ(ierr);
564b11644fSPeter Brune 
575ba6227bSPeter Brune   if (it < m) l = it;
584b11644fSPeter Brune 
594b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
604b11644fSPeter Brune   for (i = 0; i < l; i++) {
61b21d5a53SPeter Brune     k = (it - i - 1) % l;
623af51624SPeter Brune     ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr);
63bd052dfeSPeter Brune     alpha[k] = t*rho[k];
6444f7e39eSPeter Brune     if (qn->monitor) {
653af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
665ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
673af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
6844f7e39eSPeter Brune     }
696bf1b2e5SPeter Brune     ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr);
704b11644fSPeter Brune   }
714b11644fSPeter Brune 
720ecc9e1dSPeter Brune   if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
730ecc9e1dSPeter Brune     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
740ecc9e1dSPeter Brune     ierr = SNES_KSPSolve(snes,snes->ksp,Y,Yin);CHKERRQ(ierr);
750ecc9e1dSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
760ecc9e1dSPeter Brune     if (kspreason < 0) {
770ecc9e1dSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
780ecc9e1dSPeter 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);
790ecc9e1dSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
800ecc9e1dSPeter Brune         PetscFunctionReturn(0);
810ecc9e1dSPeter Brune       }
820ecc9e1dSPeter Brune     }
830ecc9e1dSPeter Brune     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
840ecc9e1dSPeter Brune     snes->linear_its += lits;
850ecc9e1dSPeter Brune     ierr = VecCopy(Yin, Y);CHKERRQ(ierr);
860ecc9e1dSPeter Brune   } else {
87b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
880ecc9e1dSPeter Brune   }
894b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
90bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
91b21d5a53SPeter Brune     k = (it + i - l) % l;
926bf1b2e5SPeter Brune     ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
934b11644fSPeter Brune     beta[k] = rho[k]*t;
943af51624SPeter Brune     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
9544f7e39eSPeter Brune     if (qn->monitor) {
963af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
975ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
983af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
9944f7e39eSPeter Brune     }
1004b11644fSPeter Brune   }
1014b11644fSPeter Brune   PetscFunctionReturn(0);
1024b11644fSPeter Brune }
1034b11644fSPeter Brune 
1044b11644fSPeter Brune #undef __FUNCT__
1054b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
1064b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
1074b11644fSPeter Brune {
1084b11644fSPeter Brune 
1094b11644fSPeter Brune   PetscErrorCode ierr;
1109f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*) snes->data;
1114b11644fSPeter Brune 
11215f5eeeaSPeter Brune   Vec X, Xold;
113*335fb975SPeter Brune   Vec F, B;
114*335fb975SPeter Brune   Vec Y, FPC, D, Dold;
1153af51624SPeter Brune   SNESConvergedReason reason;
1165ba6227bSPeter Brune   PetscInt i, i_r, k;
1174b11644fSPeter Brune 
118*335fb975SPeter Brune   PetscReal fnorm, xnorm, ynorm, gnorm;
1194b11644fSPeter Brune   PetscInt m = qn->m;
120*335fb975SPeter Brune   PetscBool lssucceed;
1214b11644fSPeter Brune 
122bd052dfeSPeter Brune   PetscScalar rhosc;
123bd052dfeSPeter Brune 
124bd052dfeSPeter Brune   Vec *dX = qn->dX;
1256bf1b2e5SPeter Brune   Vec *dF = qn->dF;
126bd052dfeSPeter Brune   PetscScalar *rho = qn->rho;
12788d374b2SPeter Brune   PetscScalar DolddotD, DolddotDold, DdotD, YdotD, a;
1284b11644fSPeter Brune 
1290ecc9e1dSPeter Brune   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
1300ecc9e1dSPeter Brune 
1314b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
1324b11644fSPeter Brune   PetscFunctionBegin;
1334b11644fSPeter Brune 
1349f3a0142SPeter Brune   X             = snes->vec_sol;        /* solution vector */
1359f3a0142SPeter Brune   F             = snes->vec_func;       /* residual vector */
1363af51624SPeter Brune   Y             = snes->vec_sol_update; /* search direction generated by J^-1D*/
1373af51624SPeter Brune   B             = snes->vec_rhs;
138*335fb975SPeter Brune   Xold          = snes->work[0];
1393af51624SPeter Brune 
1403af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
141*335fb975SPeter Brune   D             = snes->work[1];
142*335fb975SPeter Brune   Dold          = snes->work[2];
1434b11644fSPeter Brune 
1444b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1454b11644fSPeter Brune 
1464b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1474b11644fSPeter Brune   snes->iter = 0;
1484b11644fSPeter Brune   snes->norm = 0.;
1494b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15015f5eeeaSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1514b11644fSPeter Brune   if (snes->domainerror) {
1524b11644fSPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1534b11644fSPeter Brune     PetscFunctionReturn(0);
1544b11644fSPeter Brune   }
15515f5eeeaSPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1564b11644fSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1574b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1584b11644fSPeter Brune   snes->norm = fnorm;
1594b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1604b11644fSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1614b11644fSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1624b11644fSPeter Brune 
1634b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
1644b11644fSPeter Brune    snes->ttol = fnorm*snes->rtol;
1654b11644fSPeter Brune   /* test convergence */
1664b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1674b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
16870d3b23bSPeter Brune 
16988d374b2SPeter Brune   /* composed solve -- either sequential or composed */
17088d374b2SPeter Brune   if (snes->pc) {
17188d374b2SPeter Brune     if (qn->compositiontype == SNES_QN_SEQUENTIAL) {
17288d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
17388d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
17488d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
17588d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
17688d374b2SPeter Brune         PetscFunctionReturn(0);
17788d374b2SPeter Brune       }
17888d374b2SPeter Brune       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
17988d374b2SPeter Brune       ierr = VecCopy(FPC, F);CHKERRQ(ierr);
18088d374b2SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
1816bf1b2e5SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
18288d374b2SPeter Brune     } else {
18388d374b2SPeter Brune       ierr = VecCopy(X, Y);CHKERRQ(ierr);
18488d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, Y);CHKERRQ(ierr);
18588d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
18688d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
18788d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
18888d374b2SPeter Brune         PetscFunctionReturn(0);
18988d374b2SPeter Brune       }
19088d374b2SPeter Brune       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
19188d374b2SPeter Brune     }
19288d374b2SPeter Brune   } else {
19388d374b2SPeter Brune     ierr = VecCopy(F, Y);CHKERRQ(ierr);
19488d374b2SPeter Brune   }
19588d374b2SPeter Brune   ierr = VecCopy(Y, D);CHKERRQ(ierr);
1963af51624SPeter Brune 
197f8e15203SPeter Brune   /* scale the initial update */
1980ecc9e1dSPeter Brune   if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
1990ecc9e1dSPeter Brune     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
2000ecc9e1dSPeter Brune   }
2010ecc9e1dSPeter Brune 
2025ba6227bSPeter Brune   for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
203f8e15203SPeter Brune     ierr = LBGFSApplyJinv_Private(snes, i_r, D, Y);CHKERRQ(ierr);
20470d3b23bSPeter Brune     /* line search for lambda */
20570d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
20688d374b2SPeter Brune     ierr = VecCopy(D, Dold);CHKERRQ(ierr);
20715f5eeeaSPeter Brune     ierr = VecCopy(X, Xold);CHKERRQ(ierr);
208*335fb975SPeter Brune     ierr = LineSearchApply(qn->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
2099f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2109f3a0142SPeter Brune     if (snes->domainerror) {
2119f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2129f3a0142SPeter Brune       PetscFunctionReturn(0);
2139f3a0142SPeter Brune       }
214*335fb975SPeter Brune     ierr = LineSearchGetSuccess(qn->linesearch, &lssucceed);CHKERRQ(ierr);
2159f3a0142SPeter Brune     if (!lssucceed) {
2169f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
2179f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2189f3a0142SPeter Brune         break;
2199f3a0142SPeter Brune       }
2209f3a0142SPeter Brune     }
221*335fb975SPeter Brune     ierr = LineSearchGetNorms(qn->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
2220ecc9e1dSPeter Brune     if (qn->scalingtype == SNES_QN_LSSCALE) {
223*335fb975SPeter Brune       ierr = LineSearchGetLambda(qn->linesearch, &qn->scaling);CHKERRQ(ierr);
2240ecc9e1dSPeter Brune     }
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);
2700ecc9e1dSPeter 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);
2730ecc9e1dSPeter 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       }
2810ecc9e1dSPeter Brune       if (qn->scalingtype == SNES_QN_JACOBIANSCALE) {
2820ecc9e1dSPeter Brune         ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
2830ecc9e1dSPeter Brune       }
2845ba6227bSPeter Brune     } else {
285bd052dfeSPeter Brune       /* set the differences */
2865ba6227bSPeter Brune       k = i_r % m;
28788d374b2SPeter Brune       ierr = VecCopy(D, dF[k]);CHKERRQ(ierr);
28888d374b2SPeter Brune       ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
28915f5eeeaSPeter Brune       ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
29015f5eeeaSPeter Brune       ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
2916bf1b2e5SPeter Brune       ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr);
2926bf1b2e5SPeter Brune 
2936bf1b2e5SPeter Brune       /* set scaling to be shanno scaling */
294bd052dfeSPeter Brune       rho[k] = 1. / rhosc;
2950ecc9e1dSPeter Brune       if (qn->scalingtype == SNES_QN_SHANNOSCALE) {
2966bf1b2e5SPeter Brune         ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr);
2971b2f85ebSPeter Brune         qn->scaling = PetscRealPart(rhosc) / PetscRealPart(a);
2980ecc9e1dSPeter Brune       }
29970d3b23bSPeter Brune       /* general purpose update */
30070d3b23bSPeter Brune       if (snes->ops->update) {
30170d3b23bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
30270d3b23bSPeter Brune       }
3035ba6227bSPeter Brune     }
3044b11644fSPeter Brune   }
3054b11644fSPeter Brune   if (i == snes->max_its) {
3064b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
3074b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
3084b11644fSPeter Brune   }
3094b11644fSPeter Brune   PetscFunctionReturn(0);
3104b11644fSPeter Brune }
3114b11644fSPeter Brune 
3124b11644fSPeter Brune 
3134b11644fSPeter Brune #undef __FUNCT__
3144b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
3154b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
3164b11644fSPeter Brune {
3179f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
3184b11644fSPeter Brune   PetscErrorCode ierr;
319*335fb975SPeter Brune   const char     *optionsprefix;
320*335fb975SPeter Brune 
3214b11644fSPeter Brune   PetscFunctionBegin;
3224b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
3236bf1b2e5SPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr);
324bd052dfeSPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr);
325*335fb975SPeter Brune   ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr);
326*335fb975SPeter Brune 
327*335fb975SPeter Brune   /* set up the line search */
328*335fb975SPeter Brune   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
329*335fb975SPeter Brune   ierr = LineSearchCreate(((PetscObject)snes)->comm, &qn->linesearch);CHKERRQ(ierr);
330*335fb975SPeter Brune   ierr = LineSearchSetSNES(qn->linesearch, snes);CHKERRQ(ierr);
331*335fb975SPeter Brune   if (!snes->pc || qn->compositiontype == SNES_QN_SEQUENTIAL) {
332*335fb975SPeter Brune     ierr = LineSearchSetType(qn->linesearch, LINESEARCHCP);CHKERRQ(ierr);
333*335fb975SPeter Brune   } else {
334*335fb975SPeter Brune     ierr = LineSearchSetType(qn->linesearch, LINESEARCHL2);CHKERRQ(ierr);
335*335fb975SPeter Brune   }
336*335fb975SPeter Brune   ierr = LineSearchAppendOptionsPrefix(qn->linesearch, optionsprefix);CHKERRQ(ierr);
337*335fb975SPeter Brune   ierr = LineSearchSetFromOptions(qn->linesearch);CHKERRQ(ierr);
3380ecc9e1dSPeter Brune 
3394b11644fSPeter Brune   PetscFunctionReturn(0);
3404b11644fSPeter Brune }
3414b11644fSPeter Brune 
3424b11644fSPeter Brune #undef __FUNCT__
3434b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
3444b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
3454b11644fSPeter Brune {
3464b11644fSPeter Brune   PetscErrorCode ierr;
3479f83bee8SJed Brown   SNES_QN *qn;
3484b11644fSPeter Brune   PetscFunctionBegin;
3494b11644fSPeter Brune   if (snes->data) {
3509f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
3514b11644fSPeter Brune     if (qn->dX) {
3524b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
3534b11644fSPeter Brune     }
3546bf1b2e5SPeter Brune     if (qn->dF) {
3556bf1b2e5SPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr);
3564b11644fSPeter Brune     }
357bd052dfeSPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr);
3584b11644fSPeter Brune   }
3594b11644fSPeter Brune   PetscFunctionReturn(0);
3604b11644fSPeter Brune }
3614b11644fSPeter Brune 
3624b11644fSPeter Brune #undef __FUNCT__
3634b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
3644b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
3654b11644fSPeter Brune {
3664b11644fSPeter Brune   PetscErrorCode ierr;
3674b11644fSPeter Brune   PetscFunctionBegin;
3684b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
3694b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3709c05af5eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr);
3714b11644fSPeter Brune   PetscFunctionReturn(0);
3724b11644fSPeter Brune }
3734b11644fSPeter Brune 
3744b11644fSPeter Brune #undef __FUNCT__
3754b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
3764b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
3774b11644fSPeter Brune {
3784b11644fSPeter Brune 
3794b11644fSPeter Brune   PetscErrorCode ierr;
3809f83bee8SJed Brown   SNES_QN    *qn;
3815b4627e8SPeter Brune   const char *compositions[] = {"sequential", "composed"};
3820ecc9e1dSPeter Brune   const char *scalings[]     = {"shanno", "ls", "jacobian"};
38388d374b2SPeter Brune   PetscInt   indx = 0;
38488d374b2SPeter Brune   PetscBool  flg;
38544f7e39eSPeter Brune   PetscBool  monflg = PETSC_FALSE;
3864b11644fSPeter Brune   PetscFunctionBegin;
3874b11644fSPeter Brune 
3889f83bee8SJed Brown   qn = (SNES_QN*)snes->data;
3894b11644fSPeter Brune 
3904b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
3915b4627e8SPeter Brune   ierr = PetscOptionsInt("-snes_qn_m",                "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
3920ecc9e1dSPeter Brune   ierr = PetscOptionsInt("-snes_qn_restart",                "Maximum number of iterations between restarts", "SNESQN", qn->n_restart, &qn->n_restart, PETSC_NULL);CHKERRQ(ierr);
3935b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_gamma",    "Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr);
3945b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_downhill", "Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr);
3955b4627e8SPeter Brune   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr);
3965b4627e8SPeter Brune   ierr = PetscOptionsEList("-snes_qn_composition",    "Composition type",                "SNESQN",compositions,2,"sequential",&indx,&flg);CHKERRQ(ierr);
39788d374b2SPeter Brune   if (flg) {
39888d374b2SPeter Brune     switch (indx) {
3995b4627e8SPeter Brune     case 0: qn->compositiontype = SNES_QN_SEQUENTIAL;
4005b4627e8SPeter Brune       break;
4015b4627e8SPeter Brune     case 1: qn->compositiontype = SNES_QN_COMPOSED;
4025b4627e8SPeter Brune       break;
40388d374b2SPeter Brune     }
40488d374b2SPeter Brune   }
4050ecc9e1dSPeter Brune   ierr = PetscOptionsEList("-snes_qn_scaling",        "Scaling type",                    "SNESQN",scalings,3,"shanno",&indx,&flg);CHKERRQ(ierr);
4060ecc9e1dSPeter Brune   if (flg) {
4070ecc9e1dSPeter Brune     switch (indx) {
4080ecc9e1dSPeter Brune     case 0: qn->scalingtype = SNES_QN_SHANNOSCALE;
4090ecc9e1dSPeter Brune       break;
4100ecc9e1dSPeter Brune     case 1: qn->scalingtype = SNES_QN_LSSCALE;
4110ecc9e1dSPeter Brune       break;
4120ecc9e1dSPeter Brune     case 2: qn->scalingtype = SNES_QN_JACOBIANSCALE;
4130ecc9e1dSPeter Brune       snes->usesksp = PETSC_TRUE;
4140ecc9e1dSPeter Brune       break;
4150ecc9e1dSPeter Brune     }
4160ecc9e1dSPeter Brune   }
4170ecc9e1dSPeter Brune 
4184b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
41944f7e39eSPeter Brune   if (monflg) {
42044f7e39eSPeter Brune     qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
42144f7e39eSPeter Brune   }
4224b11644fSPeter Brune   PetscFunctionReturn(0);
4234b11644fSPeter Brune }
4244b11644fSPeter Brune 
42592c02d66SPeter Brune 
4264b11644fSPeter Brune /* -------------------------------------------------------------------------- */
4274b11644fSPeter Brune /*MC
4284b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
4294b11644fSPeter Brune 
4306cc8130cSPeter Brune       Options Database:
4316cc8130cSPeter Brune 
4326cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
4331867fe5bSPeter Brune .     -snes_qn_powell_angle - Angle condition for restart.
4341867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
4351867fe5bSPeter Brune .     -snes_qn_composition <sequential, composed>- Type of composition.
436f3aaf736SPeter Brune .     -snes_ls <basic, basicnonorms, quadratic, critical> - Type of line search.
43744f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
4386cc8130cSPeter Brune 
43944f7e39eSPeter Brune       Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to
4406cc8130cSPeter Brune       form the approximate inverse Jacobian using a series of multiplicative rank-one updates.  This will eventually be
4416cc8130cSPeter Brune       generalized to implement several limited-memory Broyden methods.
4426cc8130cSPeter Brune 
4431867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
4441867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
4451867fe5bSPeter Brune       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
4461867fe5bSPeter Brune       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
4471867fe5bSPeter Brune 
4486cc8130cSPeter Brune       References:
4496cc8130cSPeter Brune 
4506cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
4516cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
4526cc8130cSPeter Brune 
4534b11644fSPeter Brune 
4544b11644fSPeter Brune       Level: beginner
4554b11644fSPeter Brune 
4564b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
4576cc8130cSPeter Brune 
4584b11644fSPeter Brune M*/
4594b11644fSPeter Brune EXTERN_C_BEGIN
4604b11644fSPeter Brune #undef __FUNCT__
4614b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
4624b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
4634b11644fSPeter Brune {
4644b11644fSPeter Brune 
4654b11644fSPeter Brune   PetscErrorCode ierr;
4669f83bee8SJed Brown   SNES_QN *qn;
4674b11644fSPeter Brune 
4684b11644fSPeter Brune   PetscFunctionBegin;
4694b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
4704b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
4714b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
4724b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
4734b11644fSPeter Brune   snes->ops->view            = 0;
4744b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
4754b11644fSPeter Brune 
47642f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
47742f4f86dSBarry Smith   snes->usesksp         = PETSC_FALSE;
47842f4f86dSBarry Smith 
4790e444f03SPeter Brune   snes->max_funcs = 30000;
4800e444f03SPeter Brune   snes->max_its   = 10000;
4810e444f03SPeter Brune 
4829f83bee8SJed Brown   ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
4834b11644fSPeter Brune   snes->data = (void *) qn;
4840ecc9e1dSPeter Brune   qn->m               = 10;
485b21d5a53SPeter Brune   qn->scaling         = 1.0;
4864b11644fSPeter Brune   qn->dX              = PETSC_NULL;
4876bf1b2e5SPeter Brune   qn->dF              = PETSC_NULL;
48844f7e39eSPeter Brune   qn->monitor         = PETSC_NULL;
48988d374b2SPeter Brune   qn->powell_gamma    = 0.9;
49088d374b2SPeter Brune   qn->powell_downhill = 0.2;
49188d374b2SPeter Brune   qn->compositiontype = SNES_QN_SEQUENTIAL;
4920ecc9e1dSPeter Brune   qn->scalingtype     = SNES_QN_SHANNOSCALE;
4930ecc9e1dSPeter Brune   qn->n_restart       = -1;
494ea630c6eSPeter Brune 
4954b11644fSPeter Brune   PetscFunctionReturn(0);
4964b11644fSPeter Brune }
4974b11644fSPeter Brune EXTERN_C_END
498