xref: /petsc/src/snes/impls/qn/qn.c (revision 9f83bee822802db247e6a60c934d2076658eefbe)
14b11644fSPeter Brune #include <private/snesimpl.h>
24b11644fSPeter Brune 
34b11644fSPeter Brune typedef struct {
44b11644fSPeter Brune   Vec *dX;           /* The change in X */
53af51624SPeter Brune   Vec *dD;           /* 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;
1044f7e39eSPeter Brune   PetscViewer monitor;
115ba6227bSPeter Brune   PetscReal gamma;    /* Powell restart constant */
12*9f83bee8SJed Brown } SNES_QN;
134b11644fSPeter Brune 
144b11644fSPeter Brune #undef __FUNCT__
154b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private"
163af51624SPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) {
174b11644fSPeter Brune 
184b11644fSPeter Brune   PetscErrorCode ierr;
194b11644fSPeter Brune 
20*9f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*)snes->data;
214b11644fSPeter Brune 
224b11644fSPeter Brune   Vec *dX = qn->dX;
233af51624SPeter Brune   Vec *dD = qn->dD;
244b11644fSPeter Brune 
25bd052dfeSPeter Brune   PetscScalar *alpha = qn->alpha;
26bd052dfeSPeter Brune   PetscScalar *beta = qn->beta;
27bd052dfeSPeter Brune   PetscScalar *rho = qn->rho;
28bd052dfeSPeter Brune 
294b11644fSPeter Brune   PetscInt k, i;
304b11644fSPeter Brune   PetscInt m = qn->m;
314b11644fSPeter Brune   PetscScalar t;
324b11644fSPeter Brune   PetscInt l = m;
334b11644fSPeter Brune 
344b11644fSPeter Brune   PetscFunctionBegin;
354b11644fSPeter Brune 
363af51624SPeter Brune   ierr = VecCopy(D, Y);CHKERRQ(ierr);
374b11644fSPeter Brune 
385ba6227bSPeter Brune   if (it < m) l = it;
394b11644fSPeter Brune 
404b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
414b11644fSPeter Brune   for (i = 0; i < l; i++) {
423af51624SPeter Brune     k = (it - i) % l;
433af51624SPeter Brune     ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr);
44bd052dfeSPeter Brune     alpha[k] = t*rho[k];
4544f7e39eSPeter Brune     if (qn->monitor) {
463af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
475ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
483af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
4944f7e39eSPeter Brune     }
503af51624SPeter Brune     ierr = VecAXPY(Y, -alpha[k], dD[k]);CHKERRQ(ierr);
514b11644fSPeter Brune   }
524b11644fSPeter Brune 
534b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
54bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
553af51624SPeter Brune     k = (it + i - l + 1) % l;
563af51624SPeter Brune     ierr = VecDot(dD[k], Y, &t);CHKERRQ(ierr);
574b11644fSPeter Brune     beta[k] = rho[k]*t;
583af51624SPeter Brune     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
5944f7e39eSPeter Brune     if (qn->monitor) {
603af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
615ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
623af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
6344f7e39eSPeter Brune     }
644b11644fSPeter Brune   }
654b11644fSPeter Brune   PetscFunctionReturn(0);
664b11644fSPeter Brune }
674b11644fSPeter Brune 
684b11644fSPeter Brune #undef __FUNCT__
694b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
704b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
714b11644fSPeter Brune {
724b11644fSPeter Brune 
734b11644fSPeter Brune   PetscErrorCode ierr;
74*9f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*) snes->data;
754b11644fSPeter Brune 
7615f5eeeaSPeter Brune   Vec X, Xold;
773af51624SPeter Brune   Vec F, G, B;
783af51624SPeter Brune   Vec W, Y, D, Dold;
793af51624SPeter Brune   SNESConvergedReason reason;
805ba6227bSPeter Brune   PetscInt i, i_r, k;
814b11644fSPeter Brune 
829f3a0142SPeter Brune   PetscReal fnorm, xnorm = 0, ynorm, gnorm;
834b11644fSPeter Brune   PetscInt m = qn->m;
849f3a0142SPeter Brune   PetscBool lssucceed;
854b11644fSPeter Brune 
86bd052dfeSPeter Brune   PetscScalar rhosc;
87bd052dfeSPeter Brune 
88bd052dfeSPeter Brune   Vec *dX = qn->dX;
893af51624SPeter Brune   Vec *dD = qn->dD;
90bd052dfeSPeter Brune   PetscScalar *rho = qn->rho;
915ba6227bSPeter Brune   PetscScalar DolddotD, DolddotDold;
924b11644fSPeter Brune 
934b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
944b11644fSPeter Brune   PetscFunctionBegin;
954b11644fSPeter Brune 
969f3a0142SPeter Brune   X             = snes->vec_sol;        /* solution vector */
979f3a0142SPeter Brune   F             = snes->vec_func;       /* residual vector */
983af51624SPeter Brune   Y             = snes->vec_sol_update; /* search direction generated by J^-1D*/
993af51624SPeter Brune   B             = snes->vec_rhs;
10070d3b23bSPeter Brune   G             = snes->work[0];
10170d3b23bSPeter Brune   W             = snes->work[1];
10270d3b23bSPeter Brune   Xold          = snes->work[2];
1033af51624SPeter Brune 
1043af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
1053af51624SPeter Brune   D             = snes->work[3];
1063af51624SPeter Brune   Dold          = snes->work[4];
1074b11644fSPeter Brune 
1084b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1094b11644fSPeter Brune 
1104b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1114b11644fSPeter Brune   snes->iter = 0;
1124b11644fSPeter Brune   snes->norm = 0.;
1134b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
11415f5eeeaSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1154b11644fSPeter Brune   if (snes->domainerror) {
1164b11644fSPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1174b11644fSPeter Brune     PetscFunctionReturn(0);
1184b11644fSPeter Brune   }
11915f5eeeaSPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1204b11644fSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1214b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1224b11644fSPeter Brune   snes->norm = fnorm;
1234b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1244b11644fSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1254b11644fSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1264b11644fSPeter Brune 
1274b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
1284b11644fSPeter Brune    snes->ttol = fnorm*snes->rtol;
1294b11644fSPeter Brune   /* test convergence */
1304b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1314b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
13270d3b23bSPeter Brune 
13370d3b23bSPeter Brune   /* initialize the search direction as steepest descent */
13483947a81SPeter Brune   if (snes->pc) {
1353af51624SPeter Brune     ierr = VecCopy(X, D);CHKERRQ(ierr);
1363af51624SPeter Brune     ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr);
1373af51624SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
1383af51624SPeter Brune     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
1393af51624SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1403af51624SPeter Brune       PetscFunctionReturn(0);
1413af51624SPeter Brune     }
1423af51624SPeter Brune     ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr);
1433af51624SPeter Brune   } else {
1443af51624SPeter Brune     ierr = VecCopy(F, D);CHKERRQ(ierr);
1453af51624SPeter Brune   }
1463af51624SPeter Brune   ierr = VecCopy(D, Y);CHKERRQ(ierr);
1473af51624SPeter Brune 
1485ba6227bSPeter Brune   for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
14970d3b23bSPeter Brune     /* line search for lambda */
15070d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
1513af51624SPeter Brune     ierr = VecCopy(D, Dold);CHKERRQ(ierr);
15215f5eeeaSPeter Brune     ierr = VecCopy(X, Xold);CHKERRQ(ierr);
1539f3a0142SPeter Brune     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1543af51624SPeter Brune 
1559f3a0142SPeter 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);
1569f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1579f3a0142SPeter Brune     if (snes->domainerror) {
1589f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1599f3a0142SPeter Brune       PetscFunctionReturn(0);
1609f3a0142SPeter Brune     }
1619f3a0142SPeter Brune     if (!lssucceed) {
1629f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
1639f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
1649f3a0142SPeter Brune         break;
1659f3a0142SPeter Brune       }
1669f3a0142SPeter Brune     }
1679f3a0142SPeter Brune     /* Update function and solution vectors */
1689f3a0142SPeter Brune     fnorm = gnorm;
1699f3a0142SPeter Brune     ierr = VecCopy(G,F);CHKERRQ(ierr);
17015f5eeeaSPeter Brune     ierr = VecCopy(W,X);CHKERRQ(ierr);
1713af51624SPeter Brune 
1724b11644fSPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1738409ca45SMatthew G Knepley     snes->iter = i + 1;
1744b11644fSPeter Brune     snes->norm = fnorm;
1754b11644fSPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1768409ca45SMatthew G Knepley     SNESLogConvHistory(snes,snes->norm,snes->iter);
1778409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1784b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
1795ba6227bSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1804b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
1814b11644fSPeter Brune 
1825ba6227bSPeter Brune     /* create the new direction */
18383947a81SPeter Brune     if (snes->pc) {
1843af51624SPeter Brune       ierr = VecCopy(X, D);CHKERRQ(ierr);
1853af51624SPeter Brune       ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr);
1863af51624SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
1873af51624SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
1883af51624SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
1893af51624SPeter Brune         PetscFunctionReturn(0);
1903af51624SPeter Brune       }
1913af51624SPeter Brune       ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr);
1923af51624SPeter Brune     } else {
1933af51624SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
1943af51624SPeter Brune     }
1953af51624SPeter Brune 
1965ba6227bSPeter Brune     /* check restart by Powell's Criterion: |D^T H_0 Dold| > 0.2 * |Dold^T H_0 Dold| */
1975ba6227bSPeter Brune     ierr = VecDot(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
1985ba6227bSPeter Brune     ierr = VecDot(Dold, D, &DolddotD);CHKERRQ(ierr);
1995ba6227bSPeter Brune     if (PetscAbs(PetscRealPart(DolddotD)) > qn->gamma*PetscAbs(PetscRealPart(DolddotDold))) {
2005ba6227bSPeter Brune       if (qn->monitor) {
2015ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2025ba6227bSPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e|\n", k, PetscRealPart(DolddotD), qn->gamma, PetscRealPart(DolddotDold));CHKERRQ(ierr);
2035ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2045ba6227bSPeter Brune       }
2055ba6227bSPeter Brune       i_r = -1;
2065ba6227bSPeter Brune       /* general purpose update */
2075ba6227bSPeter Brune       if (snes->ops->update) {
2085ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2095ba6227bSPeter Brune       }
2105ba6227bSPeter Brune       ierr = VecCopy(D, Y);CHKERRQ(ierr);
2115ba6227bSPeter Brune     } else {
212bd052dfeSPeter Brune       /* set the differences */
2135ba6227bSPeter Brune       k = i_r % m;
2143af51624SPeter Brune       ierr = VecCopy(D, dD[k]);CHKERRQ(ierr);
2153af51624SPeter Brune       ierr = VecAXPY(dD[k], -1.0, Dold);CHKERRQ(ierr);
21615f5eeeaSPeter Brune       ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
21715f5eeeaSPeter Brune       ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
2183af51624SPeter Brune       ierr = VecDot(dX[k], dD[k], &rhosc);CHKERRQ(ierr);
219bd052dfeSPeter Brune       rho[k] = 1. / rhosc;
22070d3b23bSPeter Brune 
22170d3b23bSPeter Brune       /* general purpose update */
22270d3b23bSPeter Brune       if (snes->ops->update) {
22370d3b23bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
22470d3b23bSPeter Brune       }
22570d3b23bSPeter Brune       /* apply the current iteration of the approximate jacobian in order to get the next search direction*/
2265ba6227bSPeter Brune       ierr = LBGFSApplyJinv_Private(snes, i_r+1, D, Y);CHKERRQ(ierr);
2275ba6227bSPeter Brune     }
2284b11644fSPeter Brune }
2294b11644fSPeter Brune   if (i == snes->max_its) {
2304b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
2314b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
2324b11644fSPeter Brune   }
2334b11644fSPeter Brune   PetscFunctionReturn(0);
2344b11644fSPeter Brune }
2354b11644fSPeter Brune 
2364b11644fSPeter Brune 
2374b11644fSPeter Brune #undef __FUNCT__
2384b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
2394b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
2404b11644fSPeter Brune {
241*9f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*)snes->data;
2424b11644fSPeter Brune   PetscErrorCode ierr;
2434b11644fSPeter Brune   PetscFunctionBegin;
2444b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
2453af51624SPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dD);CHKERRQ(ierr);
246bd052dfeSPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr);
2473af51624SPeter Brune   ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr);
2484b11644fSPeter Brune   PetscFunctionReturn(0);
2494b11644fSPeter Brune }
2504b11644fSPeter Brune 
2514b11644fSPeter Brune #undef __FUNCT__
2524b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
2534b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
2544b11644fSPeter Brune {
2554b11644fSPeter Brune   PetscErrorCode ierr;
256*9f83bee8SJed Brown   SNES_QN *qn;
2574b11644fSPeter Brune   PetscFunctionBegin;
2584b11644fSPeter Brune   if (snes->data) {
259*9f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
2604b11644fSPeter Brune     if (qn->dX) {
2614b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
2624b11644fSPeter Brune     }
2633af51624SPeter Brune     if (qn->dD) {
2643af51624SPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dD);CHKERRQ(ierr);
2654b11644fSPeter Brune     }
266bd052dfeSPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr);
2674b11644fSPeter Brune   }
2684b11644fSPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
2694b11644fSPeter Brune   PetscFunctionReturn(0);
2704b11644fSPeter Brune }
2714b11644fSPeter Brune 
2724b11644fSPeter Brune #undef __FUNCT__
2734b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
2744b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
2754b11644fSPeter Brune {
2764b11644fSPeter Brune   PetscErrorCode ierr;
2774b11644fSPeter Brune   PetscFunctionBegin;
2784b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
2794b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2804b11644fSPeter Brune   PetscFunctionReturn(0);
2814b11644fSPeter Brune }
2824b11644fSPeter Brune 
2834b11644fSPeter Brune #undef __FUNCT__
2844b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
2854b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
2864b11644fSPeter Brune {
2874b11644fSPeter Brune 
2884b11644fSPeter Brune   PetscErrorCode ierr;
289*9f83bee8SJed Brown   SNES_QN *qn;
29044f7e39eSPeter Brune   PetscBool monflg = PETSC_FALSE;
2914b11644fSPeter Brune   PetscFunctionBegin;
2924b11644fSPeter Brune 
293*9f83bee8SJed Brown   qn = (SNES_QN*)snes->data;
2944b11644fSPeter Brune 
2954b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
2964b11644fSPeter Brune   ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
2975ba6227bSPeter Brune   ierr = PetscOptionsReal("-snes_qn_gamma", "Restart condition for L-Broyden methods", "SNES", qn->gamma, &qn->gamma, PETSC_NULL);CHKERRQ(ierr);
29844f7e39eSPeter Brune   ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNES", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr);
2994b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
30044f7e39eSPeter Brune   if (monflg) {
30144f7e39eSPeter Brune     qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
30244f7e39eSPeter Brune   }
3034b11644fSPeter Brune   PetscFunctionReturn(0);
3044b11644fSPeter Brune }
3054b11644fSPeter Brune 
30692c02d66SPeter Brune EXTERN_C_BEGIN
30792c02d66SPeter Brune #undef __FUNCT__
30892c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_QN"
30992c02d66SPeter Brune PetscErrorCode  SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type)
31092c02d66SPeter Brune {
31192c02d66SPeter Brune   PetscErrorCode ierr;
31292c02d66SPeter Brune   PetscFunctionBegin;
31392c02d66SPeter Brune 
31492c02d66SPeter Brune   switch (type) {
31592c02d66SPeter Brune   case SNES_LS_BASIC:
31692c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
31792c02d66SPeter Brune     break;
31892c02d66SPeter Brune   case SNES_LS_BASIC_NONORMS:
31992c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
32092c02d66SPeter Brune     break;
32192c02d66SPeter Brune   case SNES_LS_QUADRATIC:
322af60355fSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
32392c02d66SPeter Brune     break;
324cf9edfecSPeter Brune   case SNES_LS_SECANT:
325cf9edfecSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr);
326cf9edfecSPeter Brune     break;
32792c02d66SPeter Brune   default:
32892c02d66SPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
32992c02d66SPeter Brune     break;
33092c02d66SPeter Brune   }
33192c02d66SPeter Brune   snes->ls_type = type;
33292c02d66SPeter Brune   PetscFunctionReturn(0);
33392c02d66SPeter Brune }
33492c02d66SPeter Brune EXTERN_C_END
33592c02d66SPeter Brune 
33692c02d66SPeter Brune 
3374b11644fSPeter Brune /* -------------------------------------------------------------------------- */
3384b11644fSPeter Brune /*MC
3394b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
3404b11644fSPeter Brune 
3416cc8130cSPeter Brune       Options Database:
3426cc8130cSPeter Brune 
3436cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
34444f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
3456cc8130cSPeter Brune 
34644f7e39eSPeter Brune       Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to
3476cc8130cSPeter Brune       form the approximate inverse Jacobian using a series of multiplicative rank-one updates.  This will eventually be
3486cc8130cSPeter Brune       generalized to implement several limited-memory Broyden methods.
3496cc8130cSPeter Brune 
3506cc8130cSPeter Brune       References:
3516cc8130cSPeter Brune 
3526cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
3536cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
3546cc8130cSPeter Brune 
3554b11644fSPeter Brune 
3564b11644fSPeter Brune       Level: beginner
3574b11644fSPeter Brune 
3584b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
3596cc8130cSPeter Brune 
3604b11644fSPeter Brune M*/
3614b11644fSPeter Brune EXTERN_C_BEGIN
3624b11644fSPeter Brune #undef __FUNCT__
3634b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
3644b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
3654b11644fSPeter Brune {
3664b11644fSPeter Brune 
3674b11644fSPeter Brune   PetscErrorCode ierr;
368*9f83bee8SJed Brown   SNES_QN *qn;
3694b11644fSPeter Brune 
3704b11644fSPeter Brune   PetscFunctionBegin;
3714b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
3724b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
3734b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
3744b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
3754b11644fSPeter Brune   snes->ops->view            = 0;
3764b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
3774b11644fSPeter Brune 
37842f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
37942f4f86dSBarry Smith   snes->usesksp         = PETSC_FALSE;
38042f4f86dSBarry Smith 
381*9f83bee8SJed Brown   ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
3824b11644fSPeter Brune   snes->data = (void *) qn;
38370d3b23bSPeter Brune   qn->m       = 10;
3844b11644fSPeter Brune   qn->dX      = PETSC_NULL;
3853af51624SPeter Brune   qn->dD      = PETSC_NULL;
38644f7e39eSPeter Brune   qn->monitor = PETSC_NULL;
3875ba6227bSPeter Brune   qn->gamma   = 0.999;
388ea630c6eSPeter Brune 
38992c02d66SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr);
39070d3b23bSPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
3919f3a0142SPeter Brune 
3924b11644fSPeter Brune   PetscFunctionReturn(0);
3934b11644fSPeter Brune }
3944b11644fSPeter Brune EXTERN_C_END
395