xref: /petsc/src/snes/impls/qn/qn.c (revision f3aaf7366bbc7a33966ee24f401bc15bc56f5905)
14b11644fSPeter Brune #include <private/snesimpl.h>
24b11644fSPeter Brune 
388d374b2SPeter Brune typedef enum {SNES_QN_SEQUENTIAL, SNES_QN_COMPOSED} SNESQNCompositionType;
46bf1b2e5SPeter Brune 
54b11644fSPeter Brune typedef struct {
64b11644fSPeter Brune   Vec         *dX;              /* The change in X */
76bf1b2e5SPeter Brune   Vec         *dF;              /* The change in F */
84b11644fSPeter Brune   PetscInt    m;                /* the number of kept previous steps */
9bd052dfeSPeter Brune   PetscScalar *alpha;
10bd052dfeSPeter Brune   PetscScalar *beta;
11bd052dfeSPeter Brune   PetscScalar *rho;
1244f7e39eSPeter Brune   PetscViewer monitor;
136bf1b2e5SPeter Brune   PetscReal   powell_gamma;     /* Powell angle restart condition */
146bf1b2e5SPeter Brune   PetscReal   powell_downhill;  /* Powell descent restart condition */
15b21d5a53SPeter Brune   PetscReal   scaling;          /* scaling of H0 */
1688d374b2SPeter Brune 
1788d374b2SPeter Brune   SNESQNCompositionType compositiontype; /* determine if the composition is done sequentially or as a composition */
1888d374b2SPeter Brune 
199f83bee8SJed Brown } SNES_QN;
204b11644fSPeter Brune 
214b11644fSPeter Brune #undef __FUNCT__
224b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private"
233af51624SPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) {
244b11644fSPeter Brune 
254b11644fSPeter Brune   PetscErrorCode ierr;
264b11644fSPeter Brune 
279f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*)snes->data;
284b11644fSPeter Brune 
294b11644fSPeter Brune   Vec *dX = qn->dX;
306bf1b2e5SPeter Brune   Vec *dF = qn->dF;
314b11644fSPeter Brune 
32bd052dfeSPeter Brune   PetscScalar *alpha = qn->alpha;
33bd052dfeSPeter Brune   PetscScalar *beta = qn->beta;
34bd052dfeSPeter Brune   PetscScalar *rho = qn->rho;
35bd052dfeSPeter Brune 
364b11644fSPeter Brune   PetscInt k, i;
374b11644fSPeter Brune   PetscInt m = qn->m;
384b11644fSPeter Brune   PetscScalar t;
394b11644fSPeter Brune   PetscInt l = m;
404b11644fSPeter Brune 
414b11644fSPeter Brune   PetscFunctionBegin;
424b11644fSPeter Brune 
433af51624SPeter Brune   ierr = VecCopy(D, Y);CHKERRQ(ierr);
444b11644fSPeter Brune 
455ba6227bSPeter Brune   if (it < m) l = it;
464b11644fSPeter Brune 
474b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
484b11644fSPeter Brune   for (i = 0; i < l; i++) {
49b21d5a53SPeter Brune     k = (it - i - 1) % l;
503af51624SPeter Brune     ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr);
51bd052dfeSPeter Brune     alpha[k] = t*rho[k];
5244f7e39eSPeter Brune     if (qn->monitor) {
533af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
545ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
553af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
5644f7e39eSPeter Brune     }
576bf1b2e5SPeter Brune     ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr);
584b11644fSPeter Brune   }
594b11644fSPeter Brune 
60b21d5a53SPeter Brune   ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
61b21d5a53SPeter Brune 
624b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
63bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
64b21d5a53SPeter Brune     k = (it + i - l) % l;
656bf1b2e5SPeter Brune     ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
664b11644fSPeter Brune     beta[k] = rho[k]*t;
673af51624SPeter Brune     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
6844f7e39eSPeter Brune     if (qn->monitor) {
693af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
705ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
713af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
7244f7e39eSPeter Brune     }
734b11644fSPeter Brune   }
744b11644fSPeter Brune   PetscFunctionReturn(0);
754b11644fSPeter Brune }
764b11644fSPeter Brune 
774b11644fSPeter Brune #undef __FUNCT__
784b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
794b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
804b11644fSPeter Brune {
814b11644fSPeter Brune 
824b11644fSPeter Brune   PetscErrorCode ierr;
839f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*) snes->data;
844b11644fSPeter Brune 
8515f5eeeaSPeter Brune   Vec X, Xold;
863af51624SPeter Brune   Vec F, G, B;
8788d374b2SPeter Brune   Vec W, Y, FPC, D, Dold;
883af51624SPeter Brune   SNESConvergedReason reason;
895ba6227bSPeter Brune   PetscInt i, i_r, k;
904b11644fSPeter Brune 
919f3a0142SPeter Brune   PetscReal fnorm, xnorm = 0, ynorm, gnorm;
924b11644fSPeter Brune   PetscInt m = qn->m;
9305b53524SPeter Brune   PetscBool lssucceed, changed;
944b11644fSPeter Brune 
95bd052dfeSPeter Brune   PetscScalar rhosc;
96bd052dfeSPeter Brune 
97bd052dfeSPeter Brune   Vec *dX = qn->dX;
986bf1b2e5SPeter Brune   Vec *dF = qn->dF;
99bd052dfeSPeter Brune   PetscScalar *rho = qn->rho;
10088d374b2SPeter Brune   PetscScalar DolddotD, DolddotDold, DdotD, YdotD, a;
1014b11644fSPeter Brune 
1024b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
1034b11644fSPeter Brune   PetscFunctionBegin;
1044b11644fSPeter Brune 
1059f3a0142SPeter Brune   X             = snes->vec_sol;        /* solution vector */
1069f3a0142SPeter Brune   F             = snes->vec_func;       /* residual vector */
1073af51624SPeter Brune   Y             = snes->vec_sol_update; /* search direction generated by J^-1D*/
1083af51624SPeter Brune   B             = snes->vec_rhs;
10970d3b23bSPeter Brune   G             = snes->work[0];
11070d3b23bSPeter Brune   W             = snes->work[1];
11170d3b23bSPeter Brune   Xold          = snes->work[2];
1123af51624SPeter Brune 
1133af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
11488d374b2SPeter Brune   D             = snes->work[3];
11588d374b2SPeter Brune   Dold          = snes->work[4];
1164b11644fSPeter Brune 
1174b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1184b11644fSPeter Brune 
1194b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1204b11644fSPeter Brune   snes->iter = 0;
1214b11644fSPeter Brune   snes->norm = 0.;
1224b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
12315f5eeeaSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1244b11644fSPeter Brune   if (snes->domainerror) {
1254b11644fSPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1264b11644fSPeter Brune     PetscFunctionReturn(0);
1274b11644fSPeter Brune   }
12815f5eeeaSPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1294b11644fSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1304b11644fSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1314b11644fSPeter Brune   snes->norm = fnorm;
1324b11644fSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1334b11644fSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1344b11644fSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1354b11644fSPeter Brune 
1364b11644fSPeter Brune   /* set parameter for default relative tolerance convergence test */
1374b11644fSPeter Brune    snes->ttol = fnorm*snes->rtol;
1384b11644fSPeter Brune   /* test convergence */
1394b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1404b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
14170d3b23bSPeter Brune 
14288d374b2SPeter Brune   /* composed solve -- either sequential or composed */
14388d374b2SPeter Brune   if (snes->pc) {
14488d374b2SPeter Brune     if (qn->compositiontype == SNES_QN_SEQUENTIAL) {
14588d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
14688d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
14788d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
14888d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
14988d374b2SPeter Brune         PetscFunctionReturn(0);
15088d374b2SPeter Brune       }
15188d374b2SPeter Brune       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
15288d374b2SPeter Brune       ierr = VecCopy(FPC, F);CHKERRQ(ierr);
15388d374b2SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
1546bf1b2e5SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
15588d374b2SPeter Brune     } else {
15688d374b2SPeter Brune       ierr = VecCopy(X, Y);CHKERRQ(ierr);
15788d374b2SPeter Brune       ierr = SNESSolve(snes->pc, B, Y);CHKERRQ(ierr);
15888d374b2SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
15988d374b2SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
16088d374b2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
16188d374b2SPeter Brune         PetscFunctionReturn(0);
16288d374b2SPeter Brune       }
16388d374b2SPeter Brune       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
16488d374b2SPeter Brune     }
16588d374b2SPeter Brune   } else {
16688d374b2SPeter Brune     ierr = VecCopy(F, Y);CHKERRQ(ierr);
16788d374b2SPeter Brune   }
16888d374b2SPeter Brune   ierr = VecCopy(Y, D);CHKERRQ(ierr);
1693af51624SPeter Brune 
1705ba6227bSPeter Brune   for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
17170d3b23bSPeter Brune     /* line search for lambda */
17270d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
17388d374b2SPeter Brune     ierr = VecCopy(D, Dold);CHKERRQ(ierr);
17415f5eeeaSPeter Brune     ierr = VecCopy(X, Xold);CHKERRQ(ierr);
17505b53524SPeter Brune     ierr = SNESLineSearchPreCheckApply(snes,X,Y,&changed);CHKERRQ(ierr);
17605b53524SPeter Brune     ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1779f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
1789f3a0142SPeter Brune     if (snes->domainerror) {
1799f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1809f3a0142SPeter Brune       PetscFunctionReturn(0);
1819f3a0142SPeter Brune       }
1829f3a0142SPeter Brune     if (!lssucceed) {
1839f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
1849f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
1859f3a0142SPeter Brune         break;
1869f3a0142SPeter Brune       }
1879f3a0142SPeter Brune     }
188b21d5a53SPeter Brune 
1899f3a0142SPeter Brune     /* Update function and solution vectors */
1909f3a0142SPeter Brune     fnorm = gnorm;
1919f3a0142SPeter Brune     ierr = VecCopy(G,F);CHKERRQ(ierr);
19215f5eeeaSPeter Brune     ierr = VecCopy(W,X);CHKERRQ(ierr);
1933af51624SPeter Brune 
19488d374b2SPeter Brune     /* convergence monitoring */
195b21d5a53SPeter 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);
196b21d5a53SPeter Brune 
197360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
198360c497dSPeter Brune     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
199360c497dSPeter Brune 
2008409ca45SMatthew G Knepley     SNESLogConvHistory(snes,snes->norm,snes->iter);
2018409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
2024b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
2035ba6227bSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2044b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
2054b11644fSPeter Brune 
20688d374b2SPeter Brune 
20788d374b2SPeter Brune     if (snes->pc) {
20888d374b2SPeter Brune       if (qn->compositiontype == SNES_QN_SEQUENTIAL) {
20988d374b2SPeter Brune         ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
21088d374b2SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
21188d374b2SPeter Brune         if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
21288d374b2SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
21388d374b2SPeter Brune           PetscFunctionReturn(0);
21488d374b2SPeter Brune         }
21588d374b2SPeter Brune         ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
21688d374b2SPeter Brune         ierr = VecCopy(FPC, F);CHKERRQ(ierr);
21788d374b2SPeter Brune         ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
21888d374b2SPeter Brune         ierr = VecCopy(F, D);CHKERRQ(ierr);
21988d374b2SPeter Brune       } else {
22088d374b2SPeter Brune         ierr = VecCopy(X, D);CHKERRQ(ierr);
22188d374b2SPeter Brune         ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr);
22288d374b2SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
22388d374b2SPeter Brune         if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
22488d374b2SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
22588d374b2SPeter Brune           PetscFunctionReturn(0);
22688d374b2SPeter Brune         }
22788d374b2SPeter Brune         ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr);
22888d374b2SPeter Brune       }
22988d374b2SPeter Brune     } else {
23088d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
23188d374b2SPeter Brune     }
23288d374b2SPeter Brune 
2336bf1b2e5SPeter Brune     /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
23488d374b2SPeter Brune     ierr = VecDot(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
23588d374b2SPeter Brune     ierr = VecDot(Dold, D, &DolddotD);CHKERRQ(ierr);
23688d374b2SPeter Brune     ierr = VecDot(D, D, &DdotD);CHKERRQ(ierr);
23788d374b2SPeter Brune     ierr = VecDot(Y, D, &YdotD);CHKERRQ(ierr);
23888d374b2SPeter Brune     if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) {
2395ba6227bSPeter Brune       if (qn->monitor) {
2405ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
24188d374b2SPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e|\n", k, PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold));CHKERRQ(ierr);
2425ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2435ba6227bSPeter Brune       }
2445ba6227bSPeter Brune       i_r = -1;
2455ba6227bSPeter Brune       /* general purpose update */
2465ba6227bSPeter Brune       if (snes->ops->update) {
2475ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2485ba6227bSPeter Brune       }
24988d374b2SPeter Brune       ierr = VecCopy(D, Y);CHKERRQ(ierr);
2505ba6227bSPeter Brune     } else {
251bd052dfeSPeter Brune       /* set the differences */
2525ba6227bSPeter Brune       k = i_r % m;
25388d374b2SPeter Brune       ierr = VecCopy(D, dF[k]);CHKERRQ(ierr);
25488d374b2SPeter Brune       ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
25515f5eeeaSPeter Brune       ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
25615f5eeeaSPeter Brune       ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
2576bf1b2e5SPeter Brune       ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr);
2586bf1b2e5SPeter Brune 
2596bf1b2e5SPeter Brune       /* set scaling to be shanno scaling */
260bd052dfeSPeter Brune       rho[k] = 1. / rhosc;
2616bf1b2e5SPeter Brune       ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr);
2621b2f85ebSPeter Brune       qn->scaling = PetscRealPart(rhosc) / PetscRealPart(a);
26370d3b23bSPeter Brune 
26470d3b23bSPeter Brune       /* general purpose update */
26570d3b23bSPeter Brune       if (snes->ops->update) {
26670d3b23bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
26770d3b23bSPeter Brune       }
26870d3b23bSPeter Brune       /* apply the current iteration of the approximate jacobian in order to get the next search direction*/
26988d374b2SPeter Brune       ierr = LBGFSApplyJinv_Private(snes, i_r+1, D, Y);CHKERRQ(ierr);
2705ba6227bSPeter Brune     }
2714b11644fSPeter Brune }
2724b11644fSPeter Brune   if (i == snes->max_its) {
2734b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
2744b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
2754b11644fSPeter Brune   }
2764b11644fSPeter Brune   PetscFunctionReturn(0);
2774b11644fSPeter Brune }
2784b11644fSPeter Brune 
2794b11644fSPeter Brune 
2804b11644fSPeter Brune #undef __FUNCT__
2814b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
2824b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
2834b11644fSPeter Brune {
2849f83bee8SJed Brown   SNES_QN *qn = (SNES_QN*)snes->data;
2854b11644fSPeter Brune   PetscErrorCode ierr;
2864b11644fSPeter Brune   PetscFunctionBegin;
2874b11644fSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
2886bf1b2e5SPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr);
289bd052dfeSPeter Brune   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr);
29088d374b2SPeter Brune   ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr);
2914b11644fSPeter Brune   PetscFunctionReturn(0);
2924b11644fSPeter Brune }
2934b11644fSPeter Brune 
2944b11644fSPeter Brune #undef __FUNCT__
2954b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
2964b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
2974b11644fSPeter Brune {
2984b11644fSPeter Brune   PetscErrorCode ierr;
2999f83bee8SJed Brown   SNES_QN *qn;
3004b11644fSPeter Brune   PetscFunctionBegin;
3014b11644fSPeter Brune   if (snes->data) {
3029f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
3034b11644fSPeter Brune     if (qn->dX) {
3044b11644fSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
3054b11644fSPeter Brune     }
3066bf1b2e5SPeter Brune     if (qn->dF) {
3076bf1b2e5SPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr);
3084b11644fSPeter Brune     }
309bd052dfeSPeter Brune     ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr);
3104b11644fSPeter Brune   }
3114b11644fSPeter Brune   PetscFunctionReturn(0);
3124b11644fSPeter Brune }
3134b11644fSPeter Brune 
3144b11644fSPeter Brune #undef __FUNCT__
3154b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
3164b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
3174b11644fSPeter Brune {
3184b11644fSPeter Brune   PetscErrorCode ierr;
3194b11644fSPeter Brune   PetscFunctionBegin;
3204b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
3214b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3229c05af5eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr);
3234b11644fSPeter Brune   PetscFunctionReturn(0);
3244b11644fSPeter Brune }
3254b11644fSPeter Brune 
3264b11644fSPeter Brune #undef __FUNCT__
3274b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
3284b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
3294b11644fSPeter Brune {
3304b11644fSPeter Brune 
3314b11644fSPeter Brune   PetscErrorCode ierr;
3329f83bee8SJed Brown   SNES_QN    *qn;
3335b4627e8SPeter Brune   const char *compositions[] = {"sequential", "composed"};
33488d374b2SPeter Brune   PetscInt   indx = 0;
33588d374b2SPeter Brune   PetscBool  flg;
33644f7e39eSPeter Brune   PetscBool  monflg = PETSC_FALSE;
3374b11644fSPeter Brune   PetscFunctionBegin;
3384b11644fSPeter Brune 
3399f83bee8SJed Brown   qn = (SNES_QN*)snes->data;
3404b11644fSPeter Brune 
3414b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
3425b4627e8SPeter Brune   ierr = PetscOptionsInt("-snes_qn_m",                "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
3435b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_gamma",    "Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr);
3445b4627e8SPeter Brune   ierr = PetscOptionsReal("-snes_qn_powell_downhill", "Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr);
3455b4627e8SPeter Brune   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr);
3465b4627e8SPeter Brune   ierr = PetscOptionsEList("-snes_qn_composition",    "Composition type",                "SNESQN",compositions,2,"sequential",&indx,&flg);CHKERRQ(ierr);
34788d374b2SPeter Brune   if (flg) {
34888d374b2SPeter Brune     switch (indx) {
3495b4627e8SPeter Brune     case 0: qn->compositiontype = SNES_QN_SEQUENTIAL;
3505b4627e8SPeter Brune       break;
3515b4627e8SPeter Brune     case 1: qn->compositiontype = SNES_QN_COMPOSED;
3525b4627e8SPeter Brune       break;
35388d374b2SPeter Brune     }
35488d374b2SPeter Brune   }
3554b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
35644f7e39eSPeter Brune   if (monflg) {
35744f7e39eSPeter Brune     qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
35844f7e39eSPeter Brune   }
3594b11644fSPeter Brune   PetscFunctionReturn(0);
3604b11644fSPeter Brune }
3614b11644fSPeter Brune 
36292c02d66SPeter Brune EXTERN_C_BEGIN
36392c02d66SPeter Brune #undef __FUNCT__
36492c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_QN"
36592c02d66SPeter Brune PetscErrorCode  SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type)
36692c02d66SPeter Brune {
36792c02d66SPeter Brune   PetscErrorCode ierr;
36892c02d66SPeter Brune   PetscFunctionBegin;
36992c02d66SPeter Brune 
37092c02d66SPeter Brune   switch (type) {
37192c02d66SPeter Brune   case SNES_LS_BASIC:
37292c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
37392c02d66SPeter Brune     break;
37492c02d66SPeter Brune   case SNES_LS_BASIC_NONORMS:
37592c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
37692c02d66SPeter Brune     break;
37792c02d66SPeter Brune   case SNES_LS_QUADRATIC:
378af60355fSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
37992c02d66SPeter Brune     break;
380*f3aaf736SPeter Brune   case SNES_LS_CRITICAL:
381*f3aaf736SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchCriticalSecant,PETSC_NULL);CHKERRQ(ierr);
382cf9edfecSPeter Brune     break;
38392c02d66SPeter Brune   default:
38492c02d66SPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
38592c02d66SPeter Brune     break;
38692c02d66SPeter Brune   }
38792c02d66SPeter Brune   snes->ls_type = type;
38892c02d66SPeter Brune   PetscFunctionReturn(0);
38992c02d66SPeter Brune }
39092c02d66SPeter Brune EXTERN_C_END
39192c02d66SPeter Brune 
39292c02d66SPeter Brune 
3934b11644fSPeter Brune /* -------------------------------------------------------------------------- */
3944b11644fSPeter Brune /*MC
3954b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
3964b11644fSPeter Brune 
3976cc8130cSPeter Brune       Options Database:
3986cc8130cSPeter Brune 
3996cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
4001867fe5bSPeter Brune .     -snes_qn_powell_angle - Angle condition for restart.
4011867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
4021867fe5bSPeter Brune .     -snes_qn_composition <sequential, composed>- Type of composition.
403*f3aaf736SPeter Brune .     -snes_ls <basic, basicnonorms, quadratic, critical> - Type of line search.
40444f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
4056cc8130cSPeter Brune 
40644f7e39eSPeter Brune       Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to
4076cc8130cSPeter Brune       form the approximate inverse Jacobian using a series of multiplicative rank-one updates.  This will eventually be
4086cc8130cSPeter Brune       generalized to implement several limited-memory Broyden methods.
4096cc8130cSPeter Brune 
4101867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
4111867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
4121867fe5bSPeter Brune       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
4131867fe5bSPeter Brune       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
4141867fe5bSPeter Brune 
4156cc8130cSPeter Brune       References:
4166cc8130cSPeter Brune 
4176cc8130cSPeter Brune       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
4186cc8130cSPeter Brune       International Journal of Computer Mathematics, vol. 86, 2009.
4196cc8130cSPeter Brune 
4204b11644fSPeter Brune 
4214b11644fSPeter Brune       Level: beginner
4224b11644fSPeter Brune 
4234b11644fSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
4246cc8130cSPeter Brune 
4254b11644fSPeter Brune M*/
4264b11644fSPeter Brune EXTERN_C_BEGIN
4274b11644fSPeter Brune #undef __FUNCT__
4284b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
4294b11644fSPeter Brune PetscErrorCode  SNESCreate_QN(SNES snes)
4304b11644fSPeter Brune {
4314b11644fSPeter Brune 
4324b11644fSPeter Brune   PetscErrorCode ierr;
4339f83bee8SJed Brown   SNES_QN *qn;
4344b11644fSPeter Brune 
4354b11644fSPeter Brune   PetscFunctionBegin;
4364b11644fSPeter Brune   snes->ops->setup           = SNESSetUp_QN;
4374b11644fSPeter Brune   snes->ops->solve           = SNESSolve_QN;
4384b11644fSPeter Brune   snes->ops->destroy         = SNESDestroy_QN;
4394b11644fSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
4404b11644fSPeter Brune   snes->ops->view            = 0;
4414b11644fSPeter Brune   snes->ops->reset           = SNESReset_QN;
4424b11644fSPeter Brune 
44342f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
44442f4f86dSBarry Smith   snes->usesksp         = PETSC_FALSE;
44542f4f86dSBarry Smith 
4460e444f03SPeter Brune   snes->max_funcs = 30000;
4470e444f03SPeter Brune   snes->max_its   = 10000;
4480e444f03SPeter Brune 
4499f83bee8SJed Brown   ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
4504b11644fSPeter Brune   snes->data = (void *) qn;
45188d374b2SPeter Brune   qn->m               = 30;
452b21d5a53SPeter Brune   qn->scaling         = 1.0;
4534b11644fSPeter Brune   qn->dX              = PETSC_NULL;
4546bf1b2e5SPeter Brune   qn->dF              = PETSC_NULL;
45544f7e39eSPeter Brune   qn->monitor         = PETSC_NULL;
45688d374b2SPeter Brune   qn->powell_gamma    = 0.9;
45788d374b2SPeter Brune   qn->powell_downhill = 0.2;
45888d374b2SPeter Brune   qn->compositiontype = SNES_QN_SEQUENTIAL;
459ea630c6eSPeter Brune 
46092c02d66SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr);
461*f3aaf736SPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_CRITICAL);CHKERRQ(ierr);
4629f3a0142SPeter Brune 
4634b11644fSPeter Brune   PetscFunctionReturn(0);
4644b11644fSPeter Brune }
4654b11644fSPeter Brune EXTERN_C_END
466