14b11644fSPeter Brune #include <private/snesimpl.h> 24b11644fSPeter Brune 34b11644fSPeter Brune typedef struct { 44b11644fSPeter Brune Vec * dX; /* The change in X */ 54b11644fSPeter Brune Vec * dF; /* The change in F */ 64b11644fSPeter Brune PetscInt m; /* the number of kept previous steps */ 7bd052dfeSPeter Brune PetscScalar * alpha; 8bd052dfeSPeter Brune PetscScalar * beta; 9bd052dfeSPeter Brune PetscScalar * rho; 104b11644fSPeter Brune } QNContext; 114b11644fSPeter Brune 124b11644fSPeter Brune #undef __FUNCT__ 134b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private" 144b11644fSPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec g, Vec z) { 154b11644fSPeter Brune 164b11644fSPeter Brune PetscErrorCode ierr; 174b11644fSPeter Brune 184b11644fSPeter Brune QNContext * qn = (QNContext *)snes->data; 194b11644fSPeter Brune 204b11644fSPeter Brune Vec * dX = qn->dX; 214b11644fSPeter Brune Vec * dF = qn->dF; 224b11644fSPeter Brune 23bd052dfeSPeter Brune PetscScalar * alpha = qn->alpha; 24bd052dfeSPeter Brune PetscScalar * beta = qn->beta; 25bd052dfeSPeter Brune PetscScalar * rho = qn->rho; 26bd052dfeSPeter Brune 274b11644fSPeter Brune PetscInt k, i; 284b11644fSPeter Brune PetscInt m = qn->m; 294b11644fSPeter Brune PetscScalar t; 304b11644fSPeter Brune PetscInt l = m; 314b11644fSPeter Brune 324b11644fSPeter Brune PetscFunctionBegin; 334b11644fSPeter Brune 34bd052dfeSPeter Brune if (it < m) l = it; 354b11644fSPeter Brune 364b11644fSPeter Brune ierr = VecCopy(g, z);CHKERRQ(ierr); 374b11644fSPeter Brune 384b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 394b11644fSPeter Brune for (i = 0; i < l; i++) { 4070d3b23bSPeter Brune k = (it - i - 1) % l; 414b11644fSPeter Brune ierr = VecDot(dX[k], z, &t);CHKERRQ(ierr); 42bd052dfeSPeter Brune alpha[k] = t*rho[k]; 4392c02d66SPeter Brune /* 4470d3b23bSPeter Brune ierr = PetscPrintf(PETSC_COMM_WORLD, "%d, alpha %g\n", k, alpha[k]);CHKERRQ(ierr); 4592c02d66SPeter Brune */ 464b11644fSPeter Brune ierr = VecAXPY(z, -alpha[k], dF[k]);CHKERRQ(ierr); 474b11644fSPeter Brune } 484b11644fSPeter Brune 494b11644fSPeter Brune /* inner application of the initial inverse jacobian approximation */ 504b11644fSPeter Brune /* right now it's just the identity. Nothing needs to go here. */ 514b11644fSPeter Brune 524b11644fSPeter Brune /* inward recursion starting at the first update and working forward*/ 53bd052dfeSPeter Brune for (i = 0; i < l; i++) { 5470d3b23bSPeter Brune k = (it + i - l) % l; 554b11644fSPeter Brune ierr = VecDot(dF[k], z, &t);CHKERRQ(ierr); 564b11644fSPeter Brune beta[k] = rho[k]*t; 574b11644fSPeter Brune ierr = VecAXPY(z, (alpha[k] - beta[k]), dX[k]); 5892c02d66SPeter Brune /* 5970d3b23bSPeter Brune ierr = PetscPrintf(PETSC_COMM_WORLD, "%d, alpha - beta %g\n", k, alpha[k] - beta[k]);CHKERRQ(ierr); 6092c02d66SPeter Brune */ 614b11644fSPeter Brune } 624b11644fSPeter Brune PetscFunctionReturn(0); 634b11644fSPeter Brune } 644b11644fSPeter Brune 654b11644fSPeter Brune #undef __FUNCT__ 6670d3b23bSPeter Brune #define __FUNCT__ "SNESLineSearchQuadratic_QN" 6770d3b23bSPeter Brune PetscErrorCode SNESLineSearchQuadratic_QN(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal xnorm,Vec G,Vec W,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 6870d3b23bSPeter Brune { 6970d3b23bSPeter Brune PetscInt i; 7070d3b23bSPeter Brune PetscReal alphas[3] = {0.0, 0.5, 1.0}; 7170d3b23bSPeter Brune PetscReal norms[3]; 7270d3b23bSPeter Brune PetscReal alpha,a,b; 7370d3b23bSPeter Brune PetscErrorCode ierr; 7470d3b23bSPeter Brune 7570d3b23bSPeter Brune PetscFunctionBegin; 7670d3b23bSPeter Brune norms[0] = fnorm; 7770d3b23bSPeter Brune for(i=1; i < 3; ++i) { 7870d3b23bSPeter Brune ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr); /* W = X^n - \alpha Y */ 7970d3b23bSPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 8070d3b23bSPeter Brune ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr); 8170d3b23bSPeter Brune } 8270d3b23bSPeter Brune for(i = 0; i < 3; ++i) { 8370d3b23bSPeter Brune norms[i] = PetscSqr(norms[i]); 8470d3b23bSPeter Brune } 8570d3b23bSPeter Brune /* Fit a quadratic: 8670d3b23bSPeter Brune If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2} 8770d3b23bSPeter Brune a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1) 8870d3b23bSPeter Brune b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1) 8970d3b23bSPeter Brune c = y_0 9070d3b23bSPeter Brune x_min = -b/2a 9170d3b23bSPeter Brune 9270d3b23bSPeter Brune If we let x_{0,1,2} = 0, 0.5, 1.0 9370d3b23bSPeter Brune a = 2 y_2 - 4 y_1 + 2 y_0 9470d3b23bSPeter Brune b = -y_2 + 4 y_1 - 3 y_0 9570d3b23bSPeter Brune c = y_0 9670d3b23bSPeter Brune */ 9770d3b23bSPeter Brune a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1])); 9870d3b23bSPeter Brune b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]); 9970d3b23bSPeter Brune /* Check for positive a (concave up) */ 10070d3b23bSPeter Brune if (a >= 0.0) { 10170d3b23bSPeter Brune alpha = -b/(2.0*a); 10270d3b23bSPeter Brune alpha = PetscMin(alpha, alphas[2]); 10370d3b23bSPeter Brune alpha = PetscMax(alpha, alphas[0]); 10470d3b23bSPeter Brune } else { 10570d3b23bSPeter Brune alpha = 1.0; 10670d3b23bSPeter Brune } 10770d3b23bSPeter Brune if (snes->ls_monitor) { 10870d3b23bSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 10970d3b23bSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr); 11070d3b23bSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 11170d3b23bSPeter Brune } 11270d3b23bSPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 11370d3b23bSPeter Brune ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr); 11470d3b23bSPeter Brune if (alpha != 1.0) { 11570d3b23bSPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 11670d3b23bSPeter Brune ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr); 11770d3b23bSPeter Brune } else { 11870d3b23bSPeter Brune *gnorm = PetscSqrtReal(norms[2]); 11970d3b23bSPeter Brune } 12070d3b23bSPeter Brune if (alpha == 0.0) *flag = PETSC_FALSE; 12170d3b23bSPeter Brune else *flag = PETSC_TRUE; 12270d3b23bSPeter Brune PetscFunctionReturn(0); 12370d3b23bSPeter Brune } 12470d3b23bSPeter Brune 12570d3b23bSPeter Brune #undef __FUNCT__ 1264b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 1274b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 1284b11644fSPeter Brune { 1294b11644fSPeter Brune 1304b11644fSPeter Brune PetscErrorCode ierr; 1314b11644fSPeter Brune QNContext * qn = (QNContext*) snes->data; 1324b11644fSPeter Brune 13315f5eeeaSPeter Brune Vec X, Xold; 1349f3a0142SPeter Brune Vec F, Fold, G; 13515f5eeeaSPeter Brune Vec W, Y; 1364b11644fSPeter Brune 137bd052dfeSPeter Brune PetscInt i, k; 1384b11644fSPeter Brune 1399f3a0142SPeter Brune PetscReal fnorm, xnorm = 0, ynorm, gnorm; 1404b11644fSPeter Brune PetscInt m = qn->m; 1419f3a0142SPeter Brune PetscBool lssucceed; 1424b11644fSPeter Brune 143bd052dfeSPeter Brune PetscScalar rhosc; 144bd052dfeSPeter Brune 145bd052dfeSPeter Brune Vec * dX = qn->dX; 146bd052dfeSPeter Brune Vec * dF = qn->dF; 147bd052dfeSPeter Brune PetscScalar * rho = qn->rho; 1484b11644fSPeter Brune 1494b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 1504b11644fSPeter Brune PetscFunctionBegin; 1514b11644fSPeter Brune 1529f3a0142SPeter Brune X = snes->vec_sol; /* solution vector */ 1539f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 15470d3b23bSPeter Brune Y = snes->vec_sol_update; /* search direction */ 15570d3b23bSPeter Brune G = snes->work[0]; 15670d3b23bSPeter Brune W = snes->work[1]; 15770d3b23bSPeter Brune Xold = snes->work[2]; 15870d3b23bSPeter Brune Fold = snes->work[3]; 1594b11644fSPeter Brune 1604b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1614b11644fSPeter Brune 1624b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1634b11644fSPeter Brune snes->iter = 0; 1644b11644fSPeter Brune snes->norm = 0.; 1654b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16615f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1674b11644fSPeter Brune if (snes->domainerror) { 1684b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1694b11644fSPeter Brune PetscFunctionReturn(0); 1704b11644fSPeter Brune } 17115f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1724b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1734b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1744b11644fSPeter Brune snes->norm = fnorm; 1754b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1764b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 1774b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1784b11644fSPeter Brune 1794b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1804b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 1814b11644fSPeter Brune /* test convergence */ 1824b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1834b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 18470d3b23bSPeter Brune 18570d3b23bSPeter Brune /* initialize the search direction as steepest descent */ 18670d3b23bSPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 18770d3b23bSPeter Brune for(i = 0; i < snes->max_its; i++) { 18870d3b23bSPeter Brune /* line search for lambda */ 18970d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 19015f5eeeaSPeter Brune ierr = VecCopy(F, Fold);CHKERRQ(ierr); 19115f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 19270d3b23bSPeter Brune ierr = VecScale(Y, -1.0);CHKERRQ(ierr); 1939f3a0142SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1949f3a0142SPeter 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); 1959f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1969f3a0142SPeter Brune if (snes->domainerror) { 1979f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1989f3a0142SPeter Brune PetscFunctionReturn(0); 1999f3a0142SPeter Brune } 2009f3a0142SPeter Brune if (!lssucceed) { 2019f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 2029f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 2039f3a0142SPeter Brune break; 2049f3a0142SPeter Brune } 2059f3a0142SPeter Brune } 2069f3a0142SPeter Brune /* Update function and solution vectors */ 2079f3a0142SPeter Brune fnorm = gnorm; 2089f3a0142SPeter Brune ierr = VecCopy(G,F);CHKERRQ(ierr); 20915f5eeeaSPeter Brune ierr = VecCopy(W,X);CHKERRQ(ierr); 2109f3a0142SPeter Brune 2114b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2128409ca45SMatthew G Knepley snes->iter = i+1; 2134b11644fSPeter Brune snes->norm = fnorm; 2144b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2158409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 2168409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 2174b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 2184b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2194b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 2204b11644fSPeter Brune 221bd052dfeSPeter Brune /* set the differences */ 222bd052dfeSPeter Brune k = i % m; 22315f5eeeaSPeter Brune ierr = VecCopy(F, dF[k]);CHKERRQ(ierr); 22415f5eeeaSPeter Brune ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr); 22515f5eeeaSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 22615f5eeeaSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 227bd052dfeSPeter Brune ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr); 228bd052dfeSPeter Brune rho[k] = 1. / rhosc; 22970d3b23bSPeter Brune 23070d3b23bSPeter Brune /* general purpose update */ 23170d3b23bSPeter Brune if (snes->ops->update) { 23270d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 23370d3b23bSPeter Brune } 23470d3b23bSPeter Brune 23570d3b23bSPeter Brune /* apply the current iteration of the approximate jacobian in order to get the next search direction*/ 23670d3b23bSPeter Brune ierr = LBGFSApplyJinv_Private(snes, i, F, Y);CHKERRQ(ierr); 2374b11644fSPeter Brune } 2384b11644fSPeter Brune if (i == snes->max_its) { 2394b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 2404b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 2414b11644fSPeter Brune } 2424b11644fSPeter Brune PetscFunctionReturn(0); 2434b11644fSPeter Brune } 2444b11644fSPeter Brune 2454b11644fSPeter Brune 2464b11644fSPeter Brune #undef __FUNCT__ 2474b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 2484b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 2494b11644fSPeter Brune { 2504b11644fSPeter Brune QNContext * qn = (QNContext *)snes->data; 2514b11644fSPeter Brune PetscErrorCode ierr; 2524b11644fSPeter Brune PetscFunctionBegin; 2534b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 2544b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 255bd052dfeSPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr); 25670d3b23bSPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 2574b11644fSPeter Brune PetscFunctionReturn(0); 2584b11644fSPeter Brune } 2594b11644fSPeter Brune 2604b11644fSPeter Brune #undef __FUNCT__ 2614b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 2624b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 2634b11644fSPeter Brune { 2644b11644fSPeter Brune PetscErrorCode ierr; 2654b11644fSPeter Brune QNContext * qn; 2664b11644fSPeter Brune PetscFunctionBegin; 2674b11644fSPeter Brune if (snes->data) { 2684b11644fSPeter Brune qn = (QNContext *)snes->data; 2694b11644fSPeter Brune if (qn->dX) { 2704b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 2714b11644fSPeter Brune } 2724b11644fSPeter Brune if (qn->dF) { 2734b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 2744b11644fSPeter Brune } 275bd052dfeSPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr); 2764b11644fSPeter Brune } 2774b11644fSPeter Brune if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 2784b11644fSPeter Brune PetscFunctionReturn(0); 2794b11644fSPeter Brune } 2804b11644fSPeter Brune 2814b11644fSPeter Brune #undef __FUNCT__ 2824b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 2834b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 2844b11644fSPeter Brune { 2854b11644fSPeter Brune PetscErrorCode ierr; 2864b11644fSPeter Brune PetscFunctionBegin; 2874b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 2884b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 2894b11644fSPeter Brune PetscFunctionReturn(0); 2904b11644fSPeter Brune } 2914b11644fSPeter Brune 2924b11644fSPeter Brune #undef __FUNCT__ 2934b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 2944b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 2954b11644fSPeter Brune { 2964b11644fSPeter Brune 2974b11644fSPeter Brune PetscErrorCode ierr; 2984b11644fSPeter Brune QNContext * qn; 2994b11644fSPeter Brune 3004b11644fSPeter Brune PetscFunctionBegin; 3014b11644fSPeter Brune 3024b11644fSPeter Brune qn = (QNContext *)snes->data; 3034b11644fSPeter Brune 3044b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 3054b11644fSPeter Brune ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 3064b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3074b11644fSPeter Brune PetscFunctionReturn(0); 3084b11644fSPeter Brune } 3094b11644fSPeter Brune 31092c02d66SPeter Brune EXTERN_C_BEGIN 31192c02d66SPeter Brune #undef __FUNCT__ 31292c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_QN" 31392c02d66SPeter Brune PetscErrorCode SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type) 31492c02d66SPeter Brune { 31592c02d66SPeter Brune PetscErrorCode ierr; 31692c02d66SPeter Brune PetscFunctionBegin; 31792c02d66SPeter Brune 31892c02d66SPeter Brune switch (type) { 31992c02d66SPeter Brune case SNES_LS_BASIC: 32092c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 32192c02d66SPeter Brune break; 32292c02d66SPeter Brune case SNES_LS_BASIC_NONORMS: 32392c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 32492c02d66SPeter Brune break; 32592c02d66SPeter Brune case SNES_LS_QUADRATIC: 32692c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_QN,PETSC_NULL);CHKERRQ(ierr); 32792c02d66SPeter Brune break; 328*cf9edfecSPeter Brune case SNES_LS_SECANT: 329*cf9edfecSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr); 330*cf9edfecSPeter Brune break; 33192c02d66SPeter Brune default: 33292c02d66SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 33392c02d66SPeter Brune break; 33492c02d66SPeter Brune } 33592c02d66SPeter Brune snes->ls_type = type; 33692c02d66SPeter Brune PetscFunctionReturn(0); 33792c02d66SPeter Brune } 33892c02d66SPeter Brune EXTERN_C_END 33992c02d66SPeter Brune 34092c02d66SPeter Brune 3414b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 3424b11644fSPeter Brune /*MC 3434b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 3444b11644fSPeter Brune 3456cc8130cSPeter Brune Options Database: 3466cc8130cSPeter Brune 3476cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 3486cc8130cSPeter Brune 3496cc8130cSPeter Brune Notes: This implements the L-BFGS algorithm for the solution of F(x) = 0 using previous change in F(x) and x to 3506cc8130cSPeter Brune form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 3516cc8130cSPeter Brune generalized to implement several limited-memory Broyden methods. 3526cc8130cSPeter Brune 3536cc8130cSPeter Brune References: 3546cc8130cSPeter Brune 3556cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 3566cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 3576cc8130cSPeter Brune 3584b11644fSPeter Brune 3594b11644fSPeter Brune Level: beginner 3604b11644fSPeter Brune 3614b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 3626cc8130cSPeter Brune 3634b11644fSPeter Brune M*/ 3644b11644fSPeter Brune EXTERN_C_BEGIN 3654b11644fSPeter Brune #undef __FUNCT__ 3664b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 3674b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 3684b11644fSPeter Brune { 3694b11644fSPeter Brune 3704b11644fSPeter Brune PetscErrorCode ierr; 3714b11644fSPeter Brune QNContext * qn; 3724b11644fSPeter Brune 3734b11644fSPeter Brune PetscFunctionBegin; 3744b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 3754b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 3764b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 3774b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 3784b11644fSPeter Brune snes->ops->view = 0; 3794b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 3804b11644fSPeter Brune 38142f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 38242f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 38342f4f86dSBarry Smith 3844b11644fSPeter Brune ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr); 3854b11644fSPeter Brune snes->data = (void *) qn; 38670d3b23bSPeter Brune qn->m = 10; 3874b11644fSPeter Brune qn->dX = PETSC_NULL; 3884b11644fSPeter Brune qn->dF = PETSC_NULL; 389ea630c6eSPeter Brune 39092c02d66SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr); 39170d3b23bSPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 3929f3a0142SPeter Brune 3934b11644fSPeter Brune PetscFunctionReturn(0); 3944b11644fSPeter Brune } 3954b11644fSPeter Brune EXTERN_C_END 396