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; 1044f7e39eSPeter Brune PetscViewer monitor; 114b11644fSPeter Brune } QNContext; 124b11644fSPeter Brune 134b11644fSPeter Brune #undef __FUNCT__ 144b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private" 154b11644fSPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec g, Vec z) { 164b11644fSPeter Brune 174b11644fSPeter Brune PetscErrorCode ierr; 184b11644fSPeter Brune 194b11644fSPeter Brune QNContext * qn = (QNContext *)snes->data; 204b11644fSPeter Brune 214b11644fSPeter Brune Vec * dX = qn->dX; 224b11644fSPeter Brune Vec * dF = qn->dF; 234b11644fSPeter Brune 24bd052dfeSPeter Brune PetscScalar * alpha = qn->alpha; 25bd052dfeSPeter Brune PetscScalar * beta = qn->beta; 26bd052dfeSPeter Brune PetscScalar * rho = qn->rho; 27bd052dfeSPeter Brune 284b11644fSPeter Brune PetscInt k, i; 294b11644fSPeter Brune PetscInt m = qn->m; 304b11644fSPeter Brune PetscScalar t; 314b11644fSPeter Brune PetscInt l = m; 324b11644fSPeter Brune 334b11644fSPeter Brune PetscFunctionBegin; 344b11644fSPeter Brune 35bd052dfeSPeter Brune if (it < m) l = it; 364b11644fSPeter Brune 374b11644fSPeter Brune ierr = VecCopy(g, z);CHKERRQ(ierr); 384b11644fSPeter Brune 394b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 404b11644fSPeter Brune for (i = 0; i < l; i++) { 4170d3b23bSPeter Brune k = (it - i - 1) % l; 424b11644fSPeter Brune ierr = VecDot(dX[k], z, &t);CHKERRQ(ierr); 43bd052dfeSPeter Brune alpha[k] = t*rho[k]; 4444f7e39eSPeter Brune if (qn->monitor) { 4544f7e39eSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 4644f7e39eSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "k: %d alpha: %14.12e\n", k, alpha[k]);CHKERRQ(ierr); 4744f7e39eSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 4844f7e39eSPeter Brune } 494b11644fSPeter Brune ierr = VecAXPY(z, -alpha[k], dF[k]);CHKERRQ(ierr); 504b11644fSPeter Brune } 514b11644fSPeter Brune 524b11644fSPeter Brune /* inner application of the initial inverse jacobian approximation */ 534b11644fSPeter Brune /* right now it's just the identity. Nothing needs to go here. */ 544b11644fSPeter Brune 554b11644fSPeter Brune /* inward recursion starting at the first update and working forward*/ 56bd052dfeSPeter Brune for (i = 0; i < l; i++) { 5770d3b23bSPeter Brune k = (it + i - l) % l; 584b11644fSPeter Brune ierr = VecDot(dF[k], z, &t);CHKERRQ(ierr); 594b11644fSPeter Brune beta[k] = rho[k]*t; 604b11644fSPeter Brune ierr = VecAXPY(z, (alpha[k] - beta[k]), dX[k]); 6144f7e39eSPeter Brune if (qn->monitor) { 6244f7e39eSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 6344f7e39eSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "k: %d alpha - beta: %14.12e\n", k, alpha[k] - beta[k]);CHKERRQ(ierr); 6444f7e39eSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 6544f7e39eSPeter Brune } 664b11644fSPeter Brune } 674b11644fSPeter Brune PetscFunctionReturn(0); 684b11644fSPeter Brune } 694b11644fSPeter Brune 704b11644fSPeter Brune #undef __FUNCT__ 714b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 724b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 734b11644fSPeter Brune { 744b11644fSPeter Brune 754b11644fSPeter Brune PetscErrorCode ierr; 764b11644fSPeter Brune QNContext * qn = (QNContext*) snes->data; 774b11644fSPeter Brune 7815f5eeeaSPeter Brune Vec X, Xold; 799f3a0142SPeter Brune Vec F, Fold, G; 8015f5eeeaSPeter Brune Vec W, Y; 814b11644fSPeter Brune 82bd052dfeSPeter Brune PetscInt i, k; 834b11644fSPeter Brune 849f3a0142SPeter Brune PetscReal fnorm, xnorm = 0, ynorm, gnorm; 854b11644fSPeter Brune PetscInt m = qn->m; 869f3a0142SPeter Brune PetscBool lssucceed; 874b11644fSPeter Brune 88bd052dfeSPeter Brune PetscScalar rhosc; 89bd052dfeSPeter Brune 90bd052dfeSPeter Brune Vec * dX = qn->dX; 91bd052dfeSPeter Brune Vec * dF = qn->dF; 92bd052dfeSPeter Brune PetscScalar * rho = qn->rho; 934b11644fSPeter Brune 944b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 954b11644fSPeter Brune PetscFunctionBegin; 964b11644fSPeter Brune 979f3a0142SPeter Brune X = snes->vec_sol; /* solution vector */ 989f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 9970d3b23bSPeter Brune Y = snes->vec_sol_update; /* search direction */ 10070d3b23bSPeter Brune G = snes->work[0]; 10170d3b23bSPeter Brune W = snes->work[1]; 10270d3b23bSPeter Brune Xold = snes->work[2]; 10370d3b23bSPeter Brune Fold = snes->work[3]; 1044b11644fSPeter Brune 1054b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1064b11644fSPeter Brune 1074b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1084b11644fSPeter Brune snes->iter = 0; 1094b11644fSPeter Brune snes->norm = 0.; 1104b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 11115f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1124b11644fSPeter Brune if (snes->domainerror) { 1134b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1144b11644fSPeter Brune PetscFunctionReturn(0); 1154b11644fSPeter Brune } 11615f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1174b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1184b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1194b11644fSPeter Brune snes->norm = fnorm; 1204b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1214b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 1224b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1234b11644fSPeter Brune 1244b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1254b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 1264b11644fSPeter Brune /* test convergence */ 1274b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1284b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 12970d3b23bSPeter Brune 13070d3b23bSPeter Brune /* initialize the search direction as steepest descent */ 13170d3b23bSPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 13270d3b23bSPeter Brune for(i = 0; i < snes->max_its; i++) { 13370d3b23bSPeter Brune /* line search for lambda */ 13470d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 13515f5eeeaSPeter Brune ierr = VecCopy(F, Fold);CHKERRQ(ierr); 13615f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 13770d3b23bSPeter Brune ierr = VecScale(Y, -1.0);CHKERRQ(ierr); 1389f3a0142SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1399f3a0142SPeter 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); 1409f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1419f3a0142SPeter Brune if (snes->domainerror) { 1429f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1439f3a0142SPeter Brune PetscFunctionReturn(0); 1449f3a0142SPeter Brune } 1459f3a0142SPeter Brune if (!lssucceed) { 1469f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1479f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 1489f3a0142SPeter Brune break; 1499f3a0142SPeter Brune } 1509f3a0142SPeter Brune } 1519f3a0142SPeter Brune /* Update function and solution vectors */ 1529f3a0142SPeter Brune fnorm = gnorm; 1539f3a0142SPeter Brune ierr = VecCopy(G,F);CHKERRQ(ierr); 15415f5eeeaSPeter Brune ierr = VecCopy(W,X);CHKERRQ(ierr); 1554b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1568409ca45SMatthew G Knepley snes->iter = i+1; 1574b11644fSPeter Brune snes->norm = fnorm; 1584b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1598409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 1608409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1614b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1624b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1634b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 1644b11644fSPeter Brune 165bd052dfeSPeter Brune /* set the differences */ 166bd052dfeSPeter Brune k = i % m; 16715f5eeeaSPeter Brune ierr = VecCopy(F, dF[k]);CHKERRQ(ierr); 16815f5eeeaSPeter Brune ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr); 16915f5eeeaSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 17015f5eeeaSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 171bd052dfeSPeter Brune ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr); 172bd052dfeSPeter Brune rho[k] = 1. / rhosc; 17370d3b23bSPeter Brune 17470d3b23bSPeter Brune /* general purpose update */ 17570d3b23bSPeter Brune if (snes->ops->update) { 17670d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 17770d3b23bSPeter Brune } 17870d3b23bSPeter Brune /* apply the current iteration of the approximate jacobian in order to get the next search direction*/ 17970d3b23bSPeter Brune ierr = LBGFSApplyJinv_Private(snes, i, F, Y);CHKERRQ(ierr); 1804b11644fSPeter Brune } 1814b11644fSPeter Brune if (i == snes->max_its) { 1824b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 1834b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1844b11644fSPeter Brune } 1854b11644fSPeter Brune PetscFunctionReturn(0); 1864b11644fSPeter Brune } 1874b11644fSPeter Brune 1884b11644fSPeter Brune 1894b11644fSPeter Brune #undef __FUNCT__ 1904b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 1914b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 1924b11644fSPeter Brune { 1934b11644fSPeter Brune QNContext * qn = (QNContext *)snes->data; 1944b11644fSPeter Brune PetscErrorCode ierr; 1954b11644fSPeter Brune PetscFunctionBegin; 1964b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 1974b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 198bd052dfeSPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr); 19970d3b23bSPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 2004b11644fSPeter Brune PetscFunctionReturn(0); 2014b11644fSPeter Brune } 2024b11644fSPeter Brune 2034b11644fSPeter Brune #undef __FUNCT__ 2044b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 2054b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 2064b11644fSPeter Brune { 2074b11644fSPeter Brune PetscErrorCode ierr; 2084b11644fSPeter Brune QNContext * qn; 2094b11644fSPeter Brune PetscFunctionBegin; 2104b11644fSPeter Brune if (snes->data) { 2114b11644fSPeter Brune qn = (QNContext *)snes->data; 2124b11644fSPeter Brune if (qn->dX) { 2134b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 2144b11644fSPeter Brune } 2154b11644fSPeter Brune if (qn->dF) { 2164b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 2174b11644fSPeter Brune } 218bd052dfeSPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr); 2194b11644fSPeter Brune } 2204b11644fSPeter Brune if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 2214b11644fSPeter Brune PetscFunctionReturn(0); 2224b11644fSPeter Brune } 2234b11644fSPeter Brune 2244b11644fSPeter Brune #undef __FUNCT__ 2254b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 2264b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 2274b11644fSPeter Brune { 2284b11644fSPeter Brune PetscErrorCode ierr; 2294b11644fSPeter Brune PetscFunctionBegin; 2304b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 2314b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 2324b11644fSPeter Brune PetscFunctionReturn(0); 2334b11644fSPeter Brune } 2344b11644fSPeter Brune 2354b11644fSPeter Brune #undef __FUNCT__ 2364b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 2374b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 2384b11644fSPeter Brune { 2394b11644fSPeter Brune 2404b11644fSPeter Brune PetscErrorCode ierr; 2414b11644fSPeter Brune QNContext * qn; 24244f7e39eSPeter Brune PetscBool monflg = PETSC_FALSE; 2434b11644fSPeter Brune PetscFunctionBegin; 2444b11644fSPeter Brune 2454b11644fSPeter Brune qn = (QNContext *)snes->data; 2464b11644fSPeter Brune 2474b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 2484b11644fSPeter Brune ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 24944f7e39eSPeter Brune ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNES", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 2504b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 25144f7e39eSPeter Brune if (monflg) { 25244f7e39eSPeter Brune qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 25344f7e39eSPeter Brune } 2544b11644fSPeter Brune PetscFunctionReturn(0); 2554b11644fSPeter Brune } 2564b11644fSPeter Brune 25792c02d66SPeter Brune EXTERN_C_BEGIN 25892c02d66SPeter Brune #undef __FUNCT__ 25992c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_QN" 26092c02d66SPeter Brune PetscErrorCode SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type) 26192c02d66SPeter Brune { 26292c02d66SPeter Brune PetscErrorCode ierr; 26392c02d66SPeter Brune PetscFunctionBegin; 26492c02d66SPeter Brune 26592c02d66SPeter Brune switch (type) { 26692c02d66SPeter Brune case SNES_LS_BASIC: 26792c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 26892c02d66SPeter Brune break; 26992c02d66SPeter Brune case SNES_LS_BASIC_NONORMS: 27092c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 27192c02d66SPeter Brune break; 27292c02d66SPeter Brune case SNES_LS_QUADRATIC: 273*af60355fSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 27492c02d66SPeter Brune break; 275cf9edfecSPeter Brune case SNES_LS_SECANT: 276cf9edfecSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr); 277cf9edfecSPeter Brune break; 27892c02d66SPeter Brune default: 27992c02d66SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 28092c02d66SPeter Brune break; 28192c02d66SPeter Brune } 28292c02d66SPeter Brune snes->ls_type = type; 28392c02d66SPeter Brune PetscFunctionReturn(0); 28492c02d66SPeter Brune } 28592c02d66SPeter Brune EXTERN_C_END 28692c02d66SPeter Brune 28792c02d66SPeter Brune 2884b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 2894b11644fSPeter Brune /*MC 2904b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 2914b11644fSPeter Brune 2926cc8130cSPeter Brune Options Database: 2936cc8130cSPeter Brune 2946cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 29544f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 2966cc8130cSPeter Brune 29744f7e39eSPeter Brune Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to 2986cc8130cSPeter Brune form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 2996cc8130cSPeter Brune generalized to implement several limited-memory Broyden methods. 3006cc8130cSPeter Brune 3016cc8130cSPeter Brune References: 3026cc8130cSPeter Brune 3036cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 3046cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 3056cc8130cSPeter Brune 3064b11644fSPeter Brune 3074b11644fSPeter Brune Level: beginner 3084b11644fSPeter Brune 3094b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 3106cc8130cSPeter Brune 3114b11644fSPeter Brune M*/ 3124b11644fSPeter Brune EXTERN_C_BEGIN 3134b11644fSPeter Brune #undef __FUNCT__ 3144b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 3154b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 3164b11644fSPeter Brune { 3174b11644fSPeter Brune 3184b11644fSPeter Brune PetscErrorCode ierr; 3194b11644fSPeter Brune QNContext * qn; 3204b11644fSPeter Brune 3214b11644fSPeter Brune PetscFunctionBegin; 3224b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 3234b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 3244b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 3254b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 3264b11644fSPeter Brune snes->ops->view = 0; 3274b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 3284b11644fSPeter Brune 32942f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 33042f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 33142f4f86dSBarry Smith 3324b11644fSPeter Brune ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr); 3334b11644fSPeter Brune snes->data = (void *) qn; 33470d3b23bSPeter Brune qn->m = 10; 3354b11644fSPeter Brune qn->dX = PETSC_NULL; 3364b11644fSPeter Brune qn->dF = PETSC_NULL; 33744f7e39eSPeter Brune qn->monitor = PETSC_NULL; 338ea630c6eSPeter Brune 33992c02d66SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr); 34070d3b23bSPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 3419f3a0142SPeter Brune 3424b11644fSPeter Brune PetscFunctionReturn(0); 3434b11644fSPeter Brune } 3444b11644fSPeter Brune EXTERN_C_END 345