xref: /petsc/src/snes/impls/qn/qn.c (revision 5ba6227bf8f603f1c4b79384d5830ef7c2ec2598)
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 */
124b11644fSPeter Brune } QNContext;
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 
204b11644fSPeter Brune   QNContext * qn = (QNContext *)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;
744b11644fSPeter Brune   QNContext * qn = (QNContext*) 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 */
1343af51624SPeter Brune   if (snes->usegs && snes->ops->computegs) {
1353af51624SPeter Brune     ierr = VecCopy(X, D);CHKERRQ(ierr);
1363af51624SPeter Brune     ierr = SNESComputeGS(snes, B, D);CHKERRQ(ierr);
1373af51624SPeter Brune     ierr = VecAYPX(D, -1.0, X);CHKERRQ(ierr);
1383af51624SPeter Brune   } else if (snes->pc) {
1393af51624SPeter Brune     ierr = VecCopy(X, D);CHKERRQ(ierr);
1403af51624SPeter Brune     ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr);
1413af51624SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
1423af51624SPeter Brune     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
1433af51624SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1443af51624SPeter Brune       PetscFunctionReturn(0);
1453af51624SPeter Brune     }
1463af51624SPeter Brune     ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr);
1473af51624SPeter Brune   } else {
1483af51624SPeter Brune     ierr = VecCopy(F, D);CHKERRQ(ierr);
1493af51624SPeter Brune   }
1503af51624SPeter Brune   ierr = VecCopy(D, Y);CHKERRQ(ierr);
1513af51624SPeter Brune 
1525ba6227bSPeter Brune   for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
15370d3b23bSPeter Brune     /* line search for lambda */
15470d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
1553af51624SPeter Brune     ierr = VecCopy(D, Dold);CHKERRQ(ierr);
15615f5eeeaSPeter Brune     ierr = VecCopy(X, Xold);CHKERRQ(ierr);
1579f3a0142SPeter Brune     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1583af51624SPeter Brune 
1599f3a0142SPeter 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);
1609f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1619f3a0142SPeter Brune     if (snes->domainerror) {
1629f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1639f3a0142SPeter Brune       PetscFunctionReturn(0);
1649f3a0142SPeter Brune     }
1659f3a0142SPeter Brune     if (!lssucceed) {
1669f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
1679f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
1689f3a0142SPeter Brune         break;
1699f3a0142SPeter Brune       }
1709f3a0142SPeter Brune     }
1719f3a0142SPeter Brune     /* Update function and solution vectors */
1729f3a0142SPeter Brune     fnorm = gnorm;
1739f3a0142SPeter Brune     ierr = VecCopy(G,F);CHKERRQ(ierr);
17415f5eeeaSPeter Brune     ierr = VecCopy(W,X);CHKERRQ(ierr);
1753af51624SPeter Brune 
1764b11644fSPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1778409ca45SMatthew G Knepley     snes->iter = i + 1;
1784b11644fSPeter Brune     snes->norm = fnorm;
1794b11644fSPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1808409ca45SMatthew G Knepley     SNESLogConvHistory(snes,snes->norm,snes->iter);
1818409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1824b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
1835ba6227bSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1844b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
1854b11644fSPeter Brune 
1865ba6227bSPeter Brune     /* create the new direction */
1873af51624SPeter Brune     if (snes->usegs && snes->ops->computegs) {
1883af51624SPeter Brune       ierr = VecCopy(X, D);CHKERRQ(ierr);
1893af51624SPeter Brune       ierr = SNESComputeGS(snes, B, D);CHKERRQ(ierr);
1903af51624SPeter Brune       ierr = VecAYPX(D, -1.0, X);CHKERRQ(ierr);
1913af51624SPeter Brune     } else if (snes->pc) {
1923af51624SPeter Brune       ierr = VecCopy(X, D);CHKERRQ(ierr);
1933af51624SPeter Brune       ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr);
1943af51624SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
1953af51624SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
1963af51624SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
1973af51624SPeter Brune         PetscFunctionReturn(0);
1983af51624SPeter Brune       }
1993af51624SPeter Brune       ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr);
2003af51624SPeter Brune     } else {
2013af51624SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
2023af51624SPeter Brune     }
2033af51624SPeter Brune 
2045ba6227bSPeter Brune     /* check restart by Powell's Criterion: |D^T H_0 Dold| > 0.2 * |Dold^T H_0 Dold| */
2055ba6227bSPeter Brune     ierr = VecDot(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
2065ba6227bSPeter Brune     ierr = VecDot(Dold, D, &DolddotD);CHKERRQ(ierr);
2075ba6227bSPeter Brune     if (PetscAbs(PetscRealPart(DolddotD)) > qn->gamma*PetscAbs(PetscRealPart(DolddotDold))) {
2085ba6227bSPeter Brune       if (qn->monitor) {
2095ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2105ba6227bSPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e|\n", k, PetscRealPart(DolddotD), qn->gamma, PetscRealPart(DolddotDold));CHKERRQ(ierr);
2115ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2125ba6227bSPeter Brune       }
2135ba6227bSPeter Brune       i_r = -1;
2145ba6227bSPeter Brune       /* general purpose update */
2155ba6227bSPeter Brune       if (snes->ops->update) {
2165ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2175ba6227bSPeter Brune       }
2185ba6227bSPeter Brune       ierr = VecCopy(D, Y);CHKERRQ(ierr);
2195ba6227bSPeter Brune     } else {
220bd052dfeSPeter Brune       /* set the differences */
2215ba6227bSPeter Brune       k = i_r % m;
2223af51624SPeter Brune       ierr = VecCopy(D, dD[k]);CHKERRQ(ierr);
2233af51624SPeter Brune       ierr = VecAXPY(dD[k], -1.0, Dold);CHKERRQ(ierr);
22415f5eeeaSPeter Brune       ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
22515f5eeeaSPeter Brune       ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
2263af51624SPeter Brune       ierr = VecDot(dX[k], dD[k], &rhosc);CHKERRQ(ierr);
227bd052dfeSPeter Brune       rho[k] = 1. / rhosc;
22870d3b23bSPeter Brune 
22970d3b23bSPeter Brune       /* general purpose update */
23070d3b23bSPeter Brune       if (snes->ops->update) {
23170d3b23bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
23270d3b23bSPeter Brune       }
23370d3b23bSPeter Brune       /* apply the current iteration of the approximate jacobian in order to get the next search direction*/
2345ba6227bSPeter Brune       ierr = LBGFSApplyJinv_Private(snes, i_r+1, D, Y);CHKERRQ(ierr);
2355ba6227bSPeter Brune     }
2364b11644fSPeter Brune }
2374b11644fSPeter Brune   if (i == snes->max_its) {
2384b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
2394b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
2404b11644fSPeter Brune   }
2414b11644fSPeter Brune   PetscFunctionReturn(0);
2424b11644fSPeter Brune }
2434b11644fSPeter Brune 
2444b11644fSPeter Brune 
2454b11644fSPeter Brune #undef __FUNCT__
2464b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
2474b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
2484b11644fSPeter Brune {
2494b11644fSPeter Brune   QNContext * qn = (QNContext *)snes->data;
2504b11644fSPeter Brune   PetscErrorCode ierr;
2514b11644fSPeter Brune   PetscFunctionBegin;
2524b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
2533af51624SPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dD);CHKERRQ(ierr);
254bd052dfeSPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr);
2553af51624SPeter Brune   ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr);
2564b11644fSPeter Brune   PetscFunctionReturn(0);
2574b11644fSPeter Brune }
2584b11644fSPeter Brune 
2594b11644fSPeter Brune #undef __FUNCT__
2604b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
2614b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
2624b11644fSPeter Brune {
2634b11644fSPeter Brune   PetscErrorCode ierr;
2644b11644fSPeter Brune   QNContext * qn;
2654b11644fSPeter Brune   PetscFunctionBegin;
2664b11644fSPeter Brune   if (snes->data) {
2674b11644fSPeter Brune     qn = (QNContext *)snes->data;
2684b11644fSPeter Brune     if (qn->dX) {
2694b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
2704b11644fSPeter Brune     }
2713af51624SPeter Brune     if (qn->dD) {
2723af51624SPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dD);CHKERRQ(ierr);
2734b11644fSPeter Brune     }
274bd052dfeSPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr);
2754b11644fSPeter Brune   }
2764b11644fSPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
2774b11644fSPeter Brune   PetscFunctionReturn(0);
2784b11644fSPeter Brune }
2794b11644fSPeter Brune 
2804b11644fSPeter Brune #undef __FUNCT__
2814b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
2824b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
2834b11644fSPeter Brune {
2844b11644fSPeter Brune   PetscErrorCode ierr;
2854b11644fSPeter Brune   PetscFunctionBegin;
2864b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
2874b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2884b11644fSPeter Brune   PetscFunctionReturn(0);
2894b11644fSPeter Brune }
2904b11644fSPeter Brune 
2914b11644fSPeter Brune #undef __FUNCT__
2924b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
2934b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
2944b11644fSPeter Brune {
2954b11644fSPeter Brune 
2964b11644fSPeter Brune   PetscErrorCode ierr;
2974b11644fSPeter Brune   QNContext * qn;
29844f7e39eSPeter Brune   PetscBool monflg = PETSC_FALSE;
2994b11644fSPeter Brune   PetscFunctionBegin;
3004b11644fSPeter Brune 
3014b11644fSPeter Brune   qn = (QNContext *)snes->data;
3024b11644fSPeter Brune 
3034b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
3044b11644fSPeter Brune   ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
3055ba6227bSPeter Brune   ierr = PetscOptionsReal("-snes_qn_gamma", "Restart condition for L-Broyden methods", "SNES", qn->gamma, &qn->gamma, PETSC_NULL);CHKERRQ(ierr);
30644f7e39eSPeter Brune   ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNES", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr);
3074b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
30844f7e39eSPeter Brune   if (monflg) {
30944f7e39eSPeter Brune     qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
31044f7e39eSPeter Brune   }
3114b11644fSPeter Brune   PetscFunctionReturn(0);
3124b11644fSPeter Brune }
3134b11644fSPeter Brune 
31492c02d66SPeter Brune EXTERN_C_BEGIN
31592c02d66SPeter Brune #undef __FUNCT__
31692c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_QN"
31792c02d66SPeter Brune PetscErrorCode  SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type)
31892c02d66SPeter Brune {
31992c02d66SPeter Brune   PetscErrorCode ierr;
32092c02d66SPeter Brune   PetscFunctionBegin;
32192c02d66SPeter Brune 
32292c02d66SPeter Brune   switch (type) {
32392c02d66SPeter Brune   case SNES_LS_BASIC:
32492c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
32592c02d66SPeter Brune     break;
32692c02d66SPeter Brune   case SNES_LS_BASIC_NONORMS:
32792c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
32892c02d66SPeter Brune     break;
32992c02d66SPeter Brune   case SNES_LS_QUADRATIC:
330*af60355fSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
33192c02d66SPeter Brune     break;
332cf9edfecSPeter Brune   case SNES_LS_SECANT:
333cf9edfecSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr);
334cf9edfecSPeter Brune     break;
33592c02d66SPeter Brune   default:
33692c02d66SPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
33792c02d66SPeter Brune     break;
33892c02d66SPeter Brune   }
33992c02d66SPeter Brune   snes->ls_type = type;
34092c02d66SPeter Brune   PetscFunctionReturn(0);
34192c02d66SPeter Brune }
34292c02d66SPeter Brune EXTERN_C_END
34392c02d66SPeter Brune 
34492c02d66SPeter Brune 
3454b11644fSPeter Brune /* -------------------------------------------------------------------------- */
3464b11644fSPeter Brune /*MC
3474b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
3484b11644fSPeter Brune 
3496cc8130cSPeter Brune       Options Database:
3506cc8130cSPeter Brune 
3516cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
35244f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
3536cc8130cSPeter Brune 
35444f7e39eSPeter Brune       Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to
3556cc8130cSPeter Brune       form the approximate inverse Jacobian using a series of multiplicative rank-one updates.  This will eventually be
3566cc8130cSPeter Brune       generalized to implement several limited-memory Broyden methods.
3576cc8130cSPeter Brune 
3586cc8130cSPeter Brune       References:
3596cc8130cSPeter Brune 
3606cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
3616cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
3626cc8130cSPeter Brune 
3634b11644fSPeter Brune 
3644b11644fSPeter Brune       Level: beginner
3654b11644fSPeter Brune 
3664b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
3676cc8130cSPeter Brune 
3684b11644fSPeter Brune M*/
3694b11644fSPeter Brune EXTERN_C_BEGIN
3704b11644fSPeter Brune #undef __FUNCT__
3714b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
3724b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
3734b11644fSPeter Brune {
3744b11644fSPeter Brune 
3754b11644fSPeter Brune   PetscErrorCode ierr;
3764b11644fSPeter Brune   QNContext * qn;
3774b11644fSPeter Brune 
3784b11644fSPeter Brune   PetscFunctionBegin;
3794b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
3804b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
3814b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
3824b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
3834b11644fSPeter Brune   snes->ops->view            = 0;
3844b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
3854b11644fSPeter Brune 
38642f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
38742f4f86dSBarry Smith   snes->usesksp         = PETSC_FALSE;
38842f4f86dSBarry Smith 
3894b11644fSPeter Brune   ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr);
3904b11644fSPeter Brune   snes->data = (void *) qn;
39170d3b23bSPeter Brune   qn->m       = 10;
3924b11644fSPeter Brune   qn->dX      = PETSC_NULL;
3933af51624SPeter Brune   qn->dD      = PETSC_NULL;
39444f7e39eSPeter Brune   qn->monitor = PETSC_NULL;
3955ba6227bSPeter Brune   qn->gamma   = 0.999;
396ea630c6eSPeter Brune 
39792c02d66SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr);
39870d3b23bSPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
3999f3a0142SPeter Brune 
4004b11644fSPeter Brune   PetscFunctionReturn(0);
4014b11644fSPeter Brune }
4024b11644fSPeter Brune EXTERN_C_END
403