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__ 664b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 674b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 684b11644fSPeter Brune { 694b11644fSPeter Brune 704b11644fSPeter Brune PetscErrorCode ierr; 714b11644fSPeter Brune QNContext * qn = (QNContext*) snes->data; 724b11644fSPeter Brune 7315f5eeeaSPeter Brune Vec X, Xold; 749f3a0142SPeter Brune Vec F, Fold, G; 7515f5eeeaSPeter Brune Vec W, Y; 764b11644fSPeter Brune 77bd052dfeSPeter Brune PetscInt i, k; 784b11644fSPeter Brune 799f3a0142SPeter Brune PetscReal fnorm, xnorm = 0, ynorm, gnorm; 804b11644fSPeter Brune PetscInt m = qn->m; 819f3a0142SPeter Brune PetscBool lssucceed; 824b11644fSPeter Brune 83bd052dfeSPeter Brune PetscScalar rhosc; 84bd052dfeSPeter Brune 85bd052dfeSPeter Brune Vec * dX = qn->dX; 86bd052dfeSPeter Brune Vec * dF = qn->dF; 87bd052dfeSPeter Brune PetscScalar * rho = qn->rho; 884b11644fSPeter Brune 894b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 904b11644fSPeter Brune PetscFunctionBegin; 914b11644fSPeter Brune 929f3a0142SPeter Brune X = snes->vec_sol; /* solution vector */ 939f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 9470d3b23bSPeter Brune Y = snes->vec_sol_update; /* search direction */ 9570d3b23bSPeter Brune G = snes->work[0]; 9670d3b23bSPeter Brune W = snes->work[1]; 9770d3b23bSPeter Brune Xold = snes->work[2]; 9870d3b23bSPeter Brune Fold = snes->work[3]; 994b11644fSPeter Brune 1004b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1014b11644fSPeter Brune 1024b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1034b11644fSPeter Brune snes->iter = 0; 1044b11644fSPeter Brune snes->norm = 0.; 1054b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 10615f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1074b11644fSPeter Brune if (snes->domainerror) { 1084b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1094b11644fSPeter Brune PetscFunctionReturn(0); 1104b11644fSPeter Brune } 11115f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1124b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1134b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1144b11644fSPeter Brune snes->norm = fnorm; 1154b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1164b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 1174b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1184b11644fSPeter Brune 1194b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1204b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 1214b11644fSPeter Brune /* test convergence */ 1224b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1234b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 12470d3b23bSPeter Brune 12570d3b23bSPeter Brune /* initialize the search direction as steepest descent */ 12670d3b23bSPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 12770d3b23bSPeter Brune for(i = 0; i < snes->max_its; i++) { 12870d3b23bSPeter Brune /* line search for lambda */ 12970d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 13015f5eeeaSPeter Brune ierr = VecCopy(F, Fold);CHKERRQ(ierr); 13115f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 13270d3b23bSPeter Brune ierr = VecScale(Y, -1.0);CHKERRQ(ierr); 1339f3a0142SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1349f3a0142SPeter 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); 1359f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1369f3a0142SPeter Brune if (snes->domainerror) { 1379f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1389f3a0142SPeter Brune PetscFunctionReturn(0); 1399f3a0142SPeter Brune } 1409f3a0142SPeter Brune if (!lssucceed) { 1419f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1429f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 1439f3a0142SPeter Brune break; 1449f3a0142SPeter Brune } 1459f3a0142SPeter Brune } 1469f3a0142SPeter Brune /* Update function and solution vectors */ 1479f3a0142SPeter Brune fnorm = gnorm; 1489f3a0142SPeter Brune ierr = VecCopy(G,F);CHKERRQ(ierr); 14915f5eeeaSPeter Brune ierr = VecCopy(W,X);CHKERRQ(ierr); 1509f3a0142SPeter Brune 1514b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1528409ca45SMatthew G Knepley snes->iter = i+1; 1534b11644fSPeter Brune snes->norm = fnorm; 1544b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1558409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 1568409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1574b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1584b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1594b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 1604b11644fSPeter Brune 161bd052dfeSPeter Brune /* set the differences */ 162bd052dfeSPeter Brune k = i % m; 16315f5eeeaSPeter Brune ierr = VecCopy(F, dF[k]);CHKERRQ(ierr); 16415f5eeeaSPeter Brune ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr); 16515f5eeeaSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 16615f5eeeaSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 167bd052dfeSPeter Brune ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr); 168bd052dfeSPeter Brune rho[k] = 1. / rhosc; 16970d3b23bSPeter Brune 17070d3b23bSPeter Brune /* general purpose update */ 17170d3b23bSPeter Brune if (snes->ops->update) { 17270d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 17370d3b23bSPeter Brune } 17470d3b23bSPeter Brune 17570d3b23bSPeter Brune /* apply the current iteration of the approximate jacobian in order to get the next search direction*/ 17670d3b23bSPeter Brune ierr = LBGFSApplyJinv_Private(snes, i, F, Y);CHKERRQ(ierr); 1774b11644fSPeter Brune } 1784b11644fSPeter Brune if (i == snes->max_its) { 1794b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 1804b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1814b11644fSPeter Brune } 1824b11644fSPeter Brune PetscFunctionReturn(0); 1834b11644fSPeter Brune } 1844b11644fSPeter Brune 1854b11644fSPeter Brune 1864b11644fSPeter Brune #undef __FUNCT__ 1874b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 1884b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 1894b11644fSPeter Brune { 1904b11644fSPeter Brune QNContext * qn = (QNContext *)snes->data; 1914b11644fSPeter Brune PetscErrorCode ierr; 1924b11644fSPeter Brune PetscFunctionBegin; 1934b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 1944b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 195bd052dfeSPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr); 19670d3b23bSPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 1974b11644fSPeter Brune PetscFunctionReturn(0); 1984b11644fSPeter Brune } 1994b11644fSPeter Brune 2004b11644fSPeter Brune #undef __FUNCT__ 2014b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 2024b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 2034b11644fSPeter Brune { 2044b11644fSPeter Brune PetscErrorCode ierr; 2054b11644fSPeter Brune QNContext * qn; 2064b11644fSPeter Brune PetscFunctionBegin; 2074b11644fSPeter Brune if (snes->data) { 2084b11644fSPeter Brune qn = (QNContext *)snes->data; 2094b11644fSPeter Brune if (qn->dX) { 2104b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 2114b11644fSPeter Brune } 2124b11644fSPeter Brune if (qn->dF) { 2134b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 2144b11644fSPeter Brune } 215bd052dfeSPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr); 2164b11644fSPeter Brune } 2174b11644fSPeter Brune if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 2184b11644fSPeter Brune PetscFunctionReturn(0); 2194b11644fSPeter Brune } 2204b11644fSPeter Brune 2214b11644fSPeter Brune #undef __FUNCT__ 2224b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 2234b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 2244b11644fSPeter Brune { 2254b11644fSPeter Brune PetscErrorCode ierr; 2264b11644fSPeter Brune PetscFunctionBegin; 2274b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 2284b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 2294b11644fSPeter Brune PetscFunctionReturn(0); 2304b11644fSPeter Brune } 2314b11644fSPeter Brune 2324b11644fSPeter Brune #undef __FUNCT__ 2334b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 2344b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 2354b11644fSPeter Brune { 2364b11644fSPeter Brune 2374b11644fSPeter Brune PetscErrorCode ierr; 2384b11644fSPeter Brune QNContext * qn; 2394b11644fSPeter Brune 2404b11644fSPeter Brune PetscFunctionBegin; 2414b11644fSPeter Brune 2424b11644fSPeter Brune qn = (QNContext *)snes->data; 2434b11644fSPeter Brune 2444b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 2454b11644fSPeter Brune ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 2464b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 2474b11644fSPeter Brune PetscFunctionReturn(0); 2484b11644fSPeter Brune } 2494b11644fSPeter Brune 25092c02d66SPeter Brune EXTERN_C_BEGIN 25192c02d66SPeter Brune #undef __FUNCT__ 25292c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_QN" 25392c02d66SPeter Brune PetscErrorCode SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type) 25492c02d66SPeter Brune { 25592c02d66SPeter Brune PetscErrorCode ierr; 25692c02d66SPeter Brune PetscFunctionBegin; 25792c02d66SPeter Brune 25892c02d66SPeter Brune switch (type) { 25992c02d66SPeter Brune case SNES_LS_BASIC: 26092c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 26192c02d66SPeter Brune break; 26292c02d66SPeter Brune case SNES_LS_BASIC_NONORMS: 26392c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 26492c02d66SPeter Brune break; 26592c02d66SPeter Brune case SNES_LS_QUADRATIC: 266*af60355fSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 26792c02d66SPeter Brune break; 268cf9edfecSPeter Brune case SNES_LS_SECANT: 269cf9edfecSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr); 270cf9edfecSPeter Brune break; 27192c02d66SPeter Brune default: 27292c02d66SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 27392c02d66SPeter Brune break; 27492c02d66SPeter Brune } 27592c02d66SPeter Brune snes->ls_type = type; 27692c02d66SPeter Brune PetscFunctionReturn(0); 27792c02d66SPeter Brune } 27892c02d66SPeter Brune EXTERN_C_END 27992c02d66SPeter Brune 28092c02d66SPeter Brune 2814b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 2824b11644fSPeter Brune /*MC 2834b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 2844b11644fSPeter Brune 2856cc8130cSPeter Brune Options Database: 2866cc8130cSPeter Brune 2876cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 2886cc8130cSPeter Brune 2896cc8130cSPeter Brune Notes: This implements the L-BFGS algorithm for the solution of F(x) = 0 using previous change in F(x) and x to 2906cc8130cSPeter Brune form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 2916cc8130cSPeter Brune generalized to implement several limited-memory Broyden methods. 2926cc8130cSPeter Brune 2936cc8130cSPeter Brune References: 2946cc8130cSPeter Brune 2956cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 2966cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 2976cc8130cSPeter Brune 2984b11644fSPeter Brune 2994b11644fSPeter Brune Level: beginner 3004b11644fSPeter Brune 3014b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 3026cc8130cSPeter Brune 3034b11644fSPeter Brune M*/ 3044b11644fSPeter Brune EXTERN_C_BEGIN 3054b11644fSPeter Brune #undef __FUNCT__ 3064b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 3074b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 3084b11644fSPeter Brune { 3094b11644fSPeter Brune 3104b11644fSPeter Brune PetscErrorCode ierr; 3114b11644fSPeter Brune QNContext * qn; 3124b11644fSPeter Brune 3134b11644fSPeter Brune PetscFunctionBegin; 3144b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 3154b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 3164b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 3174b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 3184b11644fSPeter Brune snes->ops->view = 0; 3194b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 3204b11644fSPeter Brune 32142f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 32242f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 32342f4f86dSBarry Smith 3244b11644fSPeter Brune ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr); 3254b11644fSPeter Brune snes->data = (void *) qn; 32670d3b23bSPeter Brune qn->m = 10; 3274b11644fSPeter Brune qn->dX = PETSC_NULL; 3284b11644fSPeter Brune qn->dF = PETSC_NULL; 329ea630c6eSPeter Brune 33092c02d66SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr); 33170d3b23bSPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 3329f3a0142SPeter Brune 3334b11644fSPeter Brune PetscFunctionReturn(0); 3344b11644fSPeter Brune } 3354b11644fSPeter Brune EXTERN_C_END 336