14b11644fSPeter Brune #include <private/snesimpl.h> 24b11644fSPeter Brune 3*6bf1b2e5SPeter Brune 4*6bf1b2e5SPeter Brune 54b11644fSPeter Brune typedef struct { 64b11644fSPeter Brune Vec *dX; /* The change in X */ 7*6bf1b2e5SPeter 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; 13*6bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 14*6bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 15b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 169f83bee8SJed Brown } SNES_QN; 174b11644fSPeter Brune 184b11644fSPeter Brune #undef __FUNCT__ 194b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private" 203af51624SPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) { 214b11644fSPeter Brune 224b11644fSPeter Brune PetscErrorCode ierr; 234b11644fSPeter Brune 249f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 254b11644fSPeter Brune 264b11644fSPeter Brune Vec *dX = qn->dX; 27*6bf1b2e5SPeter Brune Vec *dF = qn->dF; 284b11644fSPeter Brune 29bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 30bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 31bd052dfeSPeter Brune PetscScalar *rho = qn->rho; 32bd052dfeSPeter Brune 334b11644fSPeter Brune PetscInt k, i; 344b11644fSPeter Brune PetscInt m = qn->m; 354b11644fSPeter Brune PetscScalar t; 364b11644fSPeter Brune PetscInt l = m; 374b11644fSPeter Brune 384b11644fSPeter Brune PetscFunctionBegin; 394b11644fSPeter Brune 403af51624SPeter Brune ierr = VecCopy(D, Y);CHKERRQ(ierr); 414b11644fSPeter Brune 425ba6227bSPeter Brune if (it < m) l = it; 434b11644fSPeter Brune 444b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 454b11644fSPeter Brune for (i = 0; i < l; i++) { 46b21d5a53SPeter Brune k = (it - i - 1) % l; 473af51624SPeter Brune ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr); 48bd052dfeSPeter Brune alpha[k] = t*rho[k]; 4944f7e39eSPeter Brune if (qn->monitor) { 503af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 515ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 523af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 5344f7e39eSPeter Brune } 54*6bf1b2e5SPeter Brune ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr); 554b11644fSPeter Brune } 564b11644fSPeter Brune 57b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 58b21d5a53SPeter Brune 594b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 60bd052dfeSPeter Brune for (i = 0; i < l; i++) { 61b21d5a53SPeter Brune k = (it + i - l) % l; 62*6bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 634b11644fSPeter Brune beta[k] = rho[k]*t; 643af51624SPeter Brune ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]); 6544f7e39eSPeter Brune if (qn->monitor) { 663af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 675ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 683af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 6944f7e39eSPeter Brune } 704b11644fSPeter Brune } 714b11644fSPeter Brune PetscFunctionReturn(0); 724b11644fSPeter Brune } 734b11644fSPeter Brune 744b11644fSPeter Brune #undef __FUNCT__ 754b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 764b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 774b11644fSPeter Brune { 784b11644fSPeter Brune 794b11644fSPeter Brune PetscErrorCode ierr; 809f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 814b11644fSPeter Brune 8215f5eeeaSPeter Brune Vec X, Xold; 833af51624SPeter Brune Vec F, G, B; 84*6bf1b2e5SPeter Brune Vec W, Y, FPC, Fold; 853af51624SPeter Brune SNESConvergedReason reason; 865ba6227bSPeter Brune PetscInt i, i_r, k; 874b11644fSPeter Brune 889f3a0142SPeter Brune PetscReal fnorm, xnorm = 0, ynorm, gnorm; 894b11644fSPeter Brune PetscInt m = qn->m; 909f3a0142SPeter Brune PetscBool lssucceed; 914b11644fSPeter Brune 92bd052dfeSPeter Brune PetscScalar rhosc; 93bd052dfeSPeter Brune 94bd052dfeSPeter Brune Vec *dX = qn->dX; 95*6bf1b2e5SPeter Brune Vec *dF = qn->dF; 96bd052dfeSPeter Brune PetscScalar *rho = qn->rho; 97*6bf1b2e5SPeter Brune PetscScalar FolddotF, FolddotFold, FdotF, YdotF, a; 984b11644fSPeter Brune 994b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 1004b11644fSPeter Brune PetscFunctionBegin; 1014b11644fSPeter Brune 1029f3a0142SPeter Brune X = snes->vec_sol; /* solution vector */ 1039f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 1043af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 1053af51624SPeter Brune B = snes->vec_rhs; 10670d3b23bSPeter Brune G = snes->work[0]; 10770d3b23bSPeter Brune W = snes->work[1]; 10870d3b23bSPeter Brune Xold = snes->work[2]; 1093af51624SPeter Brune 1103af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 111*6bf1b2e5SPeter Brune Fold = snes->work[3]; 1124b11644fSPeter Brune 1134b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1144b11644fSPeter Brune 1154b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1164b11644fSPeter Brune snes->iter = 0; 1174b11644fSPeter Brune snes->norm = 0.; 1184b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 11915f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1204b11644fSPeter Brune if (snes->domainerror) { 1214b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1224b11644fSPeter Brune PetscFunctionReturn(0); 1234b11644fSPeter Brune } 12415f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1254b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1264b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1274b11644fSPeter Brune snes->norm = fnorm; 1284b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1294b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 1304b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1314b11644fSPeter Brune 1324b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1334b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 1344b11644fSPeter Brune /* test convergence */ 1354b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1364b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 13770d3b23bSPeter Brune 13870d3b23bSPeter Brune /* initialize the search direction as steepest descent */ 139*6bf1b2e5SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 1403af51624SPeter Brune 1415ba6227bSPeter Brune for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 14270d3b23bSPeter Brune /* line search for lambda */ 14370d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 144*6bf1b2e5SPeter Brune ierr = VecCopy(F, Fold);CHKERRQ(ierr); 14515f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 1469f3a0142SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1479f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1489f3a0142SPeter Brune if (snes->domainerror) { 1499f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1509f3a0142SPeter Brune PetscFunctionReturn(0); 1519f3a0142SPeter Brune } 1529f3a0142SPeter Brune if (!lssucceed) { 1539f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1549f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 1559f3a0142SPeter Brune break; 1569f3a0142SPeter Brune } 1579f3a0142SPeter Brune } 158b21d5a53SPeter Brune 1599f3a0142SPeter Brune /* Update function and solution vectors */ 1609f3a0142SPeter Brune fnorm = gnorm; 1619f3a0142SPeter Brune ierr = VecCopy(G,F);CHKERRQ(ierr); 16215f5eeeaSPeter Brune ierr = VecCopy(W,X);CHKERRQ(ierr); 1633af51624SPeter Brune 164b21d5a53SPeter Brune if (snes->pc) { 165*6bf1b2e5SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 166*6bf1b2e5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 167*6bf1b2e5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 168*6bf1b2e5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 169*6bf1b2e5SPeter Brune PetscFunctionReturn(0); 170*6bf1b2e5SPeter Brune } 171*6bf1b2e5SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 172*6bf1b2e5SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 173*6bf1b2e5SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 174b21d5a53SPeter Brune } 175b21d5a53SPeter Brune 176b21d5a53SPeter 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); 177b21d5a53SPeter Brune 1784b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1798409ca45SMatthew G Knepley snes->iter = i + 1; 1804b11644fSPeter Brune snes->norm = fnorm; 1814b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1828409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 1838409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1844b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1855ba6227bSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1864b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 1874b11644fSPeter Brune 188*6bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 189*6bf1b2e5SPeter Brune ierr = VecDot(Fold, Fold, &FolddotFold);CHKERRQ(ierr); 190*6bf1b2e5SPeter Brune ierr = VecDot(Fold, F, &FolddotF);CHKERRQ(ierr); 191*6bf1b2e5SPeter Brune ierr = VecDot(F, F, &FdotF);CHKERRQ(ierr); 192*6bf1b2e5SPeter Brune ierr = VecDot(Y, F, &YdotF);CHKERRQ(ierr); 193*6bf1b2e5SPeter Brune if (PetscAbs(PetscRealPart(FolddotF)) > qn->powell_gamma*PetscAbs(PetscRealPart(FolddotFold))) { 1945ba6227bSPeter Brune if (qn->monitor) { 1955ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 196*6bf1b2e5SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e|\n", k, PetscRealPart(FolddotF), qn->powell_gamma, PetscRealPart(FolddotFold));CHKERRQ(ierr); 1975ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1985ba6227bSPeter Brune } 1995ba6227bSPeter Brune i_r = -1; 2005ba6227bSPeter Brune /* general purpose update */ 2015ba6227bSPeter Brune if (snes->ops->update) { 2025ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 2035ba6227bSPeter Brune } 204*6bf1b2e5SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 2055ba6227bSPeter Brune } else { 206bd052dfeSPeter Brune /* set the differences */ 2075ba6227bSPeter Brune k = i_r % m; 208*6bf1b2e5SPeter Brune ierr = VecCopy(F, dF[k]);CHKERRQ(ierr); 209*6bf1b2e5SPeter Brune ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr); 21015f5eeeaSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 21115f5eeeaSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 212*6bf1b2e5SPeter Brune ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr); 213*6bf1b2e5SPeter Brune 214*6bf1b2e5SPeter Brune /* set scaling to be shanno scaling */ 215bd052dfeSPeter Brune rho[k] = 1. / rhosc; 216*6bf1b2e5SPeter Brune ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr); 2171b2f85ebSPeter Brune qn->scaling = PetscRealPart(rhosc) / PetscRealPart(a); 21870d3b23bSPeter Brune 21970d3b23bSPeter Brune /* general purpose update */ 22070d3b23bSPeter Brune if (snes->ops->update) { 22170d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 22270d3b23bSPeter Brune } 22370d3b23bSPeter Brune /* apply the current iteration of the approximate jacobian in order to get the next search direction*/ 224*6bf1b2e5SPeter Brune ierr = LBGFSApplyJinv_Private(snes, i_r+1, F, Y);CHKERRQ(ierr); 2255ba6227bSPeter Brune } 2264b11644fSPeter Brune } 2274b11644fSPeter Brune if (i == snes->max_its) { 2284b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 2294b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 2304b11644fSPeter Brune } 2314b11644fSPeter Brune PetscFunctionReturn(0); 2324b11644fSPeter Brune } 2334b11644fSPeter Brune 2344b11644fSPeter Brune 2354b11644fSPeter Brune #undef __FUNCT__ 2364b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 2374b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 2384b11644fSPeter Brune { 2399f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 2404b11644fSPeter Brune PetscErrorCode ierr; 2414b11644fSPeter Brune PetscFunctionBegin; 2424b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 243*6bf1b2e5SPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 244bd052dfeSPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr); 245*6bf1b2e5SPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 2464b11644fSPeter Brune PetscFunctionReturn(0); 2474b11644fSPeter Brune } 2484b11644fSPeter Brune 2494b11644fSPeter Brune #undef __FUNCT__ 2504b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 2514b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 2524b11644fSPeter Brune { 2534b11644fSPeter Brune PetscErrorCode ierr; 2549f83bee8SJed Brown SNES_QN *qn; 2554b11644fSPeter Brune PetscFunctionBegin; 2564b11644fSPeter Brune if (snes->data) { 2579f83bee8SJed Brown qn = (SNES_QN*)snes->data; 2584b11644fSPeter Brune if (qn->dX) { 2594b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 2604b11644fSPeter Brune } 261*6bf1b2e5SPeter Brune if (qn->dF) { 262*6bf1b2e5SPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 2634b11644fSPeter Brune } 264bd052dfeSPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr); 2654b11644fSPeter Brune } 2664b11644fSPeter Brune PetscFunctionReturn(0); 2674b11644fSPeter Brune } 2684b11644fSPeter Brune 2694b11644fSPeter Brune #undef __FUNCT__ 2704b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 2714b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 2724b11644fSPeter Brune { 2734b11644fSPeter Brune PetscErrorCode ierr; 2744b11644fSPeter Brune PetscFunctionBegin; 2754b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 2764b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 2779c05af5eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 2784b11644fSPeter Brune PetscFunctionReturn(0); 2794b11644fSPeter Brune } 2804b11644fSPeter Brune 2814b11644fSPeter Brune #undef __FUNCT__ 2824b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 2834b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 2844b11644fSPeter Brune { 2854b11644fSPeter Brune 2864b11644fSPeter Brune PetscErrorCode ierr; 2879f83bee8SJed Brown SNES_QN *qn; 28844f7e39eSPeter Brune PetscBool monflg = PETSC_FALSE; 2894b11644fSPeter Brune PetscFunctionBegin; 2904b11644fSPeter Brune 2919f83bee8SJed Brown qn = (SNES_QN*)snes->data; 2924b11644fSPeter Brune 2934b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 294*6bf1b2e5SPeter Brune ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 295*6bf1b2e5SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNES", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr); 296*6bf1b2e5SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_downhill", "Powell descent tolerance", "SNES", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr); 29744f7e39eSPeter Brune ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNES", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 2984b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 29944f7e39eSPeter Brune if (monflg) { 30044f7e39eSPeter Brune qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 30144f7e39eSPeter Brune } 3024b11644fSPeter Brune PetscFunctionReturn(0); 3034b11644fSPeter Brune } 3044b11644fSPeter Brune 30592c02d66SPeter Brune EXTERN_C_BEGIN 30692c02d66SPeter Brune #undef __FUNCT__ 30792c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_QN" 30892c02d66SPeter Brune PetscErrorCode SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type) 30992c02d66SPeter Brune { 31092c02d66SPeter Brune PetscErrorCode ierr; 31192c02d66SPeter Brune PetscFunctionBegin; 31292c02d66SPeter Brune 31392c02d66SPeter Brune switch (type) { 31492c02d66SPeter Brune case SNES_LS_BASIC: 31592c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 31692c02d66SPeter Brune break; 31792c02d66SPeter Brune case SNES_LS_BASIC_NONORMS: 31892c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 31992c02d66SPeter Brune break; 32092c02d66SPeter Brune case SNES_LS_QUADRATIC: 321af60355fSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 32292c02d66SPeter Brune break; 323cf9edfecSPeter Brune case SNES_LS_SECANT: 324cf9edfecSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr); 325cf9edfecSPeter Brune break; 32692c02d66SPeter Brune default: 32792c02d66SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 32892c02d66SPeter Brune break; 32992c02d66SPeter Brune } 33092c02d66SPeter Brune snes->ls_type = type; 33192c02d66SPeter Brune PetscFunctionReturn(0); 33292c02d66SPeter Brune } 33392c02d66SPeter Brune EXTERN_C_END 33492c02d66SPeter Brune 33592c02d66SPeter Brune 3364b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 3374b11644fSPeter Brune /*MC 3384b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 3394b11644fSPeter Brune 3406cc8130cSPeter Brune Options Database: 3416cc8130cSPeter Brune 3426cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 34344f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 3446cc8130cSPeter Brune 34544f7e39eSPeter Brune Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to 3466cc8130cSPeter Brune form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 3476cc8130cSPeter Brune generalized to implement several limited-memory Broyden methods. 3486cc8130cSPeter Brune 3496cc8130cSPeter Brune References: 3506cc8130cSPeter Brune 3516cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 3526cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 3536cc8130cSPeter Brune 3544b11644fSPeter Brune 3554b11644fSPeter Brune Level: beginner 3564b11644fSPeter Brune 3574b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 3586cc8130cSPeter Brune 3594b11644fSPeter Brune M*/ 3604b11644fSPeter Brune EXTERN_C_BEGIN 3614b11644fSPeter Brune #undef __FUNCT__ 3624b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 3634b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 3644b11644fSPeter Brune { 3654b11644fSPeter Brune 3664b11644fSPeter Brune PetscErrorCode ierr; 3679f83bee8SJed Brown SNES_QN *qn; 3684b11644fSPeter Brune 3694b11644fSPeter Brune PetscFunctionBegin; 3704b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 3714b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 3724b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 3734b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 3744b11644fSPeter Brune snes->ops->view = 0; 3754b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 3764b11644fSPeter Brune 37742f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 37842f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 37942f4f86dSBarry Smith 3800e444f03SPeter Brune snes->max_funcs = 30000; 3810e444f03SPeter Brune snes->max_its = 10000; 3820e444f03SPeter Brune 3839f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 3844b11644fSPeter Brune snes->data = (void *) qn; 38570d3b23bSPeter Brune qn->m = 10; 386b21d5a53SPeter Brune qn->scaling = 1.0; 3874b11644fSPeter Brune qn->dX = PETSC_NULL; 388*6bf1b2e5SPeter Brune qn->dF = PETSC_NULL; 38944f7e39eSPeter Brune qn->monitor = PETSC_NULL; 390*6bf1b2e5SPeter Brune qn->powell_gamma = 0.8; 391*6bf1b2e5SPeter Brune qn->powell_downhill = 0.4; 392ea630c6eSPeter Brune 39392c02d66SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr); 394b21d5a53SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_SECANT);CHKERRQ(ierr); 3959f3a0142SPeter Brune 3964b11644fSPeter Brune PetscFunctionReturn(0); 3974b11644fSPeter Brune } 3984b11644fSPeter Brune EXTERN_C_END 399