xref: /petsc/src/snes/impls/qn/qn.c (revision 70d3b23b377b60e976d343365f7c92b8b2f0f13f)
14b11644fSPeter Brune #include <private/snesimpl.h>
24b11644fSPeter Brune 
34b11644fSPeter Brune typedef struct {
44b11644fSPeter Brune   Vec * dX;           /* The change in X */
54b11644fSPeter Brune   Vec * dF;           /* The change in F */
64b11644fSPeter Brune   PetscInt m;         /* the number of kept previous steps */
7bd052dfeSPeter Brune   PetscScalar * alpha;
8bd052dfeSPeter Brune   PetscScalar * beta;
9bd052dfeSPeter Brune   PetscScalar * rho;
104b11644fSPeter Brune } QNContext;
114b11644fSPeter Brune 
124b11644fSPeter Brune #undef __FUNCT__
134b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private"
144b11644fSPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec g, Vec z) {
154b11644fSPeter Brune 
164b11644fSPeter Brune   PetscErrorCode ierr;
174b11644fSPeter Brune 
184b11644fSPeter Brune   QNContext * qn = (QNContext *)snes->data;
194b11644fSPeter Brune 
204b11644fSPeter Brune   Vec * dX = qn->dX;
214b11644fSPeter Brune   Vec * dF = qn->dF;
224b11644fSPeter Brune 
23bd052dfeSPeter Brune   PetscScalar * alpha = qn->alpha;
24bd052dfeSPeter Brune   PetscScalar * beta = qn->beta;
25bd052dfeSPeter Brune   PetscScalar * rho = qn->rho;
26bd052dfeSPeter Brune 
274b11644fSPeter Brune   PetscInt k, i;
284b11644fSPeter Brune   PetscInt m = qn->m;
294b11644fSPeter Brune   PetscScalar t;
304b11644fSPeter Brune   PetscInt l = m;
314b11644fSPeter Brune 
324b11644fSPeter Brune   PetscFunctionBegin;
334b11644fSPeter Brune 
34bd052dfeSPeter Brune   if (it < m) l = it;
354b11644fSPeter Brune 
364b11644fSPeter Brune   ierr = VecCopy(g, z);CHKERRQ(ierr);
374b11644fSPeter Brune 
384b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
394b11644fSPeter Brune   for (i = 0; i < l; i++) {
40*70d3b23bSPeter Brune     k = (it - i - 1) % l;
414b11644fSPeter Brune     ierr = VecDot(dX[k], z, &t);CHKERRQ(ierr);
42bd052dfeSPeter Brune     alpha[k] = t*rho[k];
43*70d3b23bSPeter Brune     ierr = PetscPrintf(PETSC_COMM_WORLD, "%d, alpha %g\n", k, alpha[k]);CHKERRQ(ierr);
444b11644fSPeter Brune     ierr = VecAXPY(z, -alpha[k], dF[k]);CHKERRQ(ierr);
454b11644fSPeter Brune   }
464b11644fSPeter Brune 
474b11644fSPeter Brune   /* inner application of the initial inverse jacobian approximation */
484b11644fSPeter Brune   /* right now it's just the identity. Nothing needs to go here. */
494b11644fSPeter Brune 
504b11644fSPeter Brune   /* inward recursion starting at the first update and working forward*/
51bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
52*70d3b23bSPeter Brune     k = (it + i - l) % l;
534b11644fSPeter Brune     ierr = VecDot(dF[k], z, &t);CHKERRQ(ierr);
544b11644fSPeter Brune     beta[k] = rho[k]*t;
554b11644fSPeter Brune     ierr = VecAXPY(z, (alpha[k] - beta[k]), dX[k]);
56*70d3b23bSPeter Brune     ierr = PetscPrintf(PETSC_COMM_WORLD, "%d, alpha - beta %g\n", k, alpha[k] - beta[k]);CHKERRQ(ierr);
574b11644fSPeter Brune   }
584b11644fSPeter Brune   PetscFunctionReturn(0);
594b11644fSPeter Brune }
604b11644fSPeter Brune 
614b11644fSPeter Brune #undef __FUNCT__
62*70d3b23bSPeter Brune #define __FUNCT__ "SNESLineSearchQuadratic_QN"
63*70d3b23bSPeter Brune PetscErrorCode  SNESLineSearchQuadratic_QN(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal xnorm,Vec G,Vec W,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
64*70d3b23bSPeter Brune {
65*70d3b23bSPeter Brune   PetscInt         i;
66*70d3b23bSPeter Brune   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
67*70d3b23bSPeter Brune   PetscReal        norms[3];
68*70d3b23bSPeter Brune   PetscReal        alpha,a,b;
69*70d3b23bSPeter Brune   PetscErrorCode   ierr;
70*70d3b23bSPeter Brune 
71*70d3b23bSPeter Brune   PetscFunctionBegin;
72*70d3b23bSPeter Brune   norms[0]  = fnorm;
73*70d3b23bSPeter Brune   for(i=1; i < 3; ++i) {
74*70d3b23bSPeter Brune     ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
75*70d3b23bSPeter Brune     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
76*70d3b23bSPeter Brune     ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr);
77*70d3b23bSPeter Brune   }
78*70d3b23bSPeter Brune   for(i = 0; i < 3; ++i) {
79*70d3b23bSPeter Brune     norms[i] = PetscSqr(norms[i]);
80*70d3b23bSPeter Brune   }
81*70d3b23bSPeter Brune   /* Fit a quadratic:
82*70d3b23bSPeter Brune        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
83*70d3b23bSPeter Brune        a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1)
84*70d3b23bSPeter Brune        b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1)
85*70d3b23bSPeter Brune        c = y_0
86*70d3b23bSPeter Brune        x_min = -b/2a
87*70d3b23bSPeter Brune 
88*70d3b23bSPeter Brune        If we let x_{0,1,2} = 0, 0.5, 1.0
89*70d3b23bSPeter Brune        a = 2 y_2 - 4 y_1 + 2 y_0
90*70d3b23bSPeter Brune        b =  -y_2 + 4 y_1 - 3 y_0
91*70d3b23bSPeter Brune        c =   y_0
92*70d3b23bSPeter Brune   */
93*70d3b23bSPeter Brune   a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1]));
94*70d3b23bSPeter Brune   b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]);
95*70d3b23bSPeter Brune   /* Check for positive a (concave up) */
96*70d3b23bSPeter Brune   if (a >= 0.0) {
97*70d3b23bSPeter Brune     alpha = -b/(2.0*a);
98*70d3b23bSPeter Brune     alpha = PetscMin(alpha, alphas[2]);
99*70d3b23bSPeter Brune     alpha = PetscMax(alpha, alphas[0]);
100*70d3b23bSPeter Brune   } else {
101*70d3b23bSPeter Brune     alpha = 1.0;
102*70d3b23bSPeter Brune   }
103*70d3b23bSPeter Brune   if (snes->ls_monitor) {
104*70d3b23bSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
105*70d3b23bSPeter Brune     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr);
106*70d3b23bSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
107*70d3b23bSPeter Brune   }
108*70d3b23bSPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
109*70d3b23bSPeter Brune   ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr);
110*70d3b23bSPeter Brune   if (alpha != 1.0) {
111*70d3b23bSPeter Brune     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
112*70d3b23bSPeter Brune     ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr);
113*70d3b23bSPeter Brune   } else {
114*70d3b23bSPeter Brune     *gnorm = PetscSqrtReal(norms[2]);
115*70d3b23bSPeter Brune   }
116*70d3b23bSPeter Brune   if (alpha == 0.0) *flag = PETSC_FALSE;
117*70d3b23bSPeter Brune   else              *flag = PETSC_TRUE;
118*70d3b23bSPeter Brune   PetscFunctionReturn(0);
119*70d3b23bSPeter Brune }
120*70d3b23bSPeter Brune 
121*70d3b23bSPeter Brune #undef __FUNCT__
1224b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
1234b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
1244b11644fSPeter Brune {
1254b11644fSPeter Brune 
1264b11644fSPeter Brune   PetscErrorCode ierr;
1274b11644fSPeter Brune   QNContext * qn = (QNContext*) snes->data;
1284b11644fSPeter Brune 
12915f5eeeaSPeter Brune   Vec X, Xold;
1309f3a0142SPeter Brune   Vec F, Fold, G;
13115f5eeeaSPeter Brune   Vec W, Y;
1324b11644fSPeter Brune 
133bd052dfeSPeter Brune   PetscInt i, k;
1344b11644fSPeter Brune 
1359f3a0142SPeter Brune   PetscReal fnorm, xnorm = 0, ynorm, gnorm;
1364b11644fSPeter Brune   PetscInt m = qn->m;
1379f3a0142SPeter Brune   PetscBool lssucceed;
1384b11644fSPeter Brune 
139bd052dfeSPeter Brune   PetscScalar rhosc;
140bd052dfeSPeter Brune 
141bd052dfeSPeter Brune   Vec * dX = qn->dX;
142bd052dfeSPeter Brune   Vec * dF = qn->dF;
143bd052dfeSPeter Brune   PetscScalar * rho = qn->rho;
1444b11644fSPeter Brune 
1454b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
1464b11644fSPeter Brune   PetscFunctionBegin;
1474b11644fSPeter Brune 
1489f3a0142SPeter Brune   X             = snes->vec_sol;        /* solution vector */
1499f3a0142SPeter Brune   F             = snes->vec_func;       /* residual vector */
150*70d3b23bSPeter Brune   Y             = snes->vec_sol_update; /* search direction */
151*70d3b23bSPeter Brune   G             = snes->work[0];
152*70d3b23bSPeter Brune   W             = snes->work[1];
153*70d3b23bSPeter Brune   Xold          = snes->work[2];
154*70d3b23bSPeter Brune   Fold          = snes->work[3];
1554b11644fSPeter Brune 
1564b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1574b11644fSPeter Brune 
1584b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1594b11644fSPeter Brune   snes->iter = 0;
1604b11644fSPeter Brune   snes->norm = 0.;
1614b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16215f5eeeaSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1634b11644fSPeter Brune   if (snes->domainerror) {
1644b11644fSPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1654b11644fSPeter Brune     PetscFunctionReturn(0);
1664b11644fSPeter Brune   }
16715f5eeeaSPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1684b11644fSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1694b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1704b11644fSPeter Brune   snes->norm = fnorm;
1714b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1724b11644fSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1734b11644fSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1744b11644fSPeter Brune 
1754b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
1764b11644fSPeter Brune    snes->ttol = fnorm*snes->rtol;
1774b11644fSPeter Brune   /* test convergence */
1784b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1794b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
180*70d3b23bSPeter Brune 
181*70d3b23bSPeter Brune   /* initialize the search direction as steepest descent */
182*70d3b23bSPeter Brune   ierr = VecCopy(F, Y);CHKERRQ(ierr);
183*70d3b23bSPeter Brune  for(i = 0; i < snes->max_its; i++) {
184*70d3b23bSPeter Brune    /* line search for lambda */
185*70d3b23bSPeter Brune    ynorm = 1; gnorm = fnorm;
18615f5eeeaSPeter Brune    ierr = VecCopy(F, Fold);CHKERRQ(ierr);
18715f5eeeaSPeter Brune    ierr = VecCopy(X, Xold);CHKERRQ(ierr);
188*70d3b23bSPeter Brune    ierr = VecScale(Y, -1.0);CHKERRQ(ierr);
1899f3a0142SPeter Brune    ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1909f3a0142SPeter 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);
1919f3a0142SPeter Brune    if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1929f3a0142SPeter Brune    if (snes->domainerror) {
1939f3a0142SPeter Brune      snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1949f3a0142SPeter Brune      PetscFunctionReturn(0);
1959f3a0142SPeter Brune    }
1969f3a0142SPeter Brune    if (!lssucceed) {
1979f3a0142SPeter Brune      if (++snes->numFailures >= snes->maxFailures) {
1989f3a0142SPeter Brune        snes->reason = SNES_DIVERGED_LINE_SEARCH;
1999f3a0142SPeter Brune        break;
2009f3a0142SPeter Brune      }
2019f3a0142SPeter Brune    }
2029f3a0142SPeter Brune    /* Update function and solution vectors */
2039f3a0142SPeter Brune    fnorm = gnorm;
2049f3a0142SPeter Brune    ierr = VecCopy(G,F);CHKERRQ(ierr);
20515f5eeeaSPeter Brune    ierr = VecCopy(W,X);CHKERRQ(ierr);
2069f3a0142SPeter Brune 
2074b11644fSPeter Brune    ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2088409ca45SMatthew G Knepley    snes->iter = i+1;
2094b11644fSPeter Brune    snes->norm = fnorm;
2104b11644fSPeter Brune    ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2118409ca45SMatthew G Knepley    SNESLogConvHistory(snes,snes->norm,snes->iter);
2128409ca45SMatthew G Knepley    ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
2134b11644fSPeter Brune    /* set parameter for default relative tolerance convergence test */
2144b11644fSPeter Brune    ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2154b11644fSPeter Brune    if (snes->reason) PetscFunctionReturn(0);
2164b11644fSPeter Brune 
217bd052dfeSPeter Brune    /* set the differences */
218bd052dfeSPeter Brune    k = i % m;
21915f5eeeaSPeter Brune    ierr = VecCopy(F, dF[k]);CHKERRQ(ierr);
22015f5eeeaSPeter Brune    ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr);
22115f5eeeaSPeter Brune    ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
22215f5eeeaSPeter Brune    ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
223bd052dfeSPeter Brune    ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr);
224bd052dfeSPeter Brune    rho[k] = 1. / rhosc;
225*70d3b23bSPeter Brune 
226*70d3b23bSPeter Brune    /* general purpose update */
227*70d3b23bSPeter Brune    if (snes->ops->update) {
228*70d3b23bSPeter Brune      ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
229*70d3b23bSPeter Brune    }
230*70d3b23bSPeter Brune 
231*70d3b23bSPeter Brune    /* apply the current iteration of the approximate jacobian in order to get the next search direction*/
232*70d3b23bSPeter Brune    ierr = LBGFSApplyJinv_Private(snes, i, F, Y);CHKERRQ(ierr);
2334b11644fSPeter Brune  }
2344b11644fSPeter Brune   if (i == snes->max_its) {
2354b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
2364b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
2374b11644fSPeter Brune   }
2384b11644fSPeter Brune   PetscFunctionReturn(0);
2394b11644fSPeter Brune }
2404b11644fSPeter Brune 
2414b11644fSPeter Brune 
2424b11644fSPeter Brune #undef __FUNCT__
2434b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
2444b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
2454b11644fSPeter Brune {
2464b11644fSPeter Brune   QNContext * qn = (QNContext *)snes->data;
2474b11644fSPeter Brune   PetscErrorCode ierr;
2484b11644fSPeter Brune   PetscFunctionBegin;
2494b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
2504b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr);
251bd052dfeSPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr);
252*70d3b23bSPeter Brune   ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr);
2534b11644fSPeter Brune   PetscFunctionReturn(0);
2544b11644fSPeter Brune }
2554b11644fSPeter Brune 
2564b11644fSPeter Brune #undef __FUNCT__
2574b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
2584b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
2594b11644fSPeter Brune {
2604b11644fSPeter Brune   PetscErrorCode ierr;
2614b11644fSPeter Brune   QNContext * qn;
2624b11644fSPeter Brune   PetscFunctionBegin;
2634b11644fSPeter Brune   if (snes->data) {
2644b11644fSPeter Brune     qn = (QNContext *)snes->data;
2654b11644fSPeter Brune     if (qn->dX) {
2664b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
2674b11644fSPeter Brune     }
2684b11644fSPeter Brune     if (qn->dF) {
2694b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr);
2704b11644fSPeter Brune     }
271bd052dfeSPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr);
2724b11644fSPeter Brune   }
2734b11644fSPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
2744b11644fSPeter Brune   PetscFunctionReturn(0);
2754b11644fSPeter Brune }
2764b11644fSPeter Brune 
2774b11644fSPeter Brune #undef __FUNCT__
2784b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
2794b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
2804b11644fSPeter Brune {
2814b11644fSPeter Brune   PetscErrorCode ierr;
2824b11644fSPeter Brune   PetscFunctionBegin;
2834b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
2844b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2854b11644fSPeter Brune   PetscFunctionReturn(0);
2864b11644fSPeter Brune }
2874b11644fSPeter Brune 
2884b11644fSPeter Brune #undef __FUNCT__
2894b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
2904b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
2914b11644fSPeter Brune {
2924b11644fSPeter Brune 
2934b11644fSPeter Brune   PetscErrorCode ierr;
2944b11644fSPeter Brune   QNContext * qn;
2954b11644fSPeter Brune 
2964b11644fSPeter Brune   PetscFunctionBegin;
2974b11644fSPeter Brune 
2984b11644fSPeter Brune   qn = (QNContext *)snes->data;
2994b11644fSPeter Brune 
3004b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
3014b11644fSPeter Brune   ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
3024b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
3034b11644fSPeter Brune   PetscFunctionReturn(0);
3044b11644fSPeter Brune }
3054b11644fSPeter Brune 
3064b11644fSPeter Brune /* -------------------------------------------------------------------------- */
3074b11644fSPeter Brune /*MC
3084b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
3094b11644fSPeter Brune 
3106cc8130cSPeter Brune       Options Database:
3116cc8130cSPeter Brune 
3126cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
3136cc8130cSPeter Brune 
3146cc8130cSPeter Brune       Notes: This implements the L-BFGS algorithm for the solution of F(x) = 0 using previous change in F(x) and x to
3156cc8130cSPeter Brune       form the approximate inverse Jacobian using a series of multiplicative rank-one updates.  This will eventually be
3166cc8130cSPeter Brune       generalized to implement several limited-memory Broyden methods.
3176cc8130cSPeter Brune 
3186cc8130cSPeter Brune       References:
3196cc8130cSPeter Brune 
3206cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
3216cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
3226cc8130cSPeter Brune 
3234b11644fSPeter Brune 
3244b11644fSPeter Brune       Level: beginner
3254b11644fSPeter Brune 
3264b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
3276cc8130cSPeter Brune 
3284b11644fSPeter Brune M*/
3294b11644fSPeter Brune EXTERN_C_BEGIN
3304b11644fSPeter Brune #undef __FUNCT__
3314b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
3324b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
3334b11644fSPeter Brune {
3344b11644fSPeter Brune 
3354b11644fSPeter Brune   PetscErrorCode ierr;
3364b11644fSPeter Brune   QNContext * qn;
3374b11644fSPeter Brune 
3384b11644fSPeter Brune   PetscFunctionBegin;
3394b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
3404b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
3414b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
3424b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
3434b11644fSPeter Brune   snes->ops->view            = 0;
3444b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
3454b11644fSPeter Brune 
34642f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
34742f4f86dSBarry Smith   snes->usesksp         = PETSC_FALSE;
34842f4f86dSBarry Smith 
3494b11644fSPeter Brune   ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr);
3504b11644fSPeter Brune   snes->data = (void *) qn;
351*70d3b23bSPeter Brune   qn->m = 10;
3524b11644fSPeter Brune   qn->dX = PETSC_NULL;
3534b11644fSPeter Brune   qn->dF = PETSC_NULL;
354ea630c6eSPeter Brune 
3559f3a0142SPeter Brune   snes->ops->linesearchcubic     = SNESLineSearchCubic;
356*70d3b23bSPeter Brune   snes->ops->linesearchquadratic = SNESLineSearchQuadratic_QN;
357ea630c6eSPeter Brune   snes->ops->linesearchno        = SNESLineSearchNo;
358ea630c6eSPeter Brune   snes->ops->linesearchnonorms   = SNESLineSearchNoNorms;
3599f3a0142SPeter Brune 
360*70d3b23bSPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
3619f3a0142SPeter Brune 
3624b11644fSPeter Brune   PetscFunctionReturn(0);
3634b11644fSPeter Brune }
3644b11644fSPeter Brune EXTERN_C_END
365