14b11644fSPeter Brune 24b11644fSPeter Brune #include <private/snesimpl.h> 34b11644fSPeter Brune 44b11644fSPeter Brune typedef struct { 5*40012acdSPeter Brune PetscReal lambda; /* The default step length for the update */ 64b11644fSPeter Brune Vec * dX; /* The change in X */ 74b11644fSPeter Brune Vec * dF; /* The change in F */ 84b11644fSPeter Brune PetscInt m; /* the number of kept previous steps */ 94b11644fSPeter Brune } QNContext; 104b11644fSPeter Brune 114b11644fSPeter Brune #undef __FUNCT__ 124b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private" 134b11644fSPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec g, Vec z) { 144b11644fSPeter Brune 154b11644fSPeter Brune PetscErrorCode ierr; 164b11644fSPeter Brune 174b11644fSPeter Brune QNContext * qn = (QNContext *)snes->data; 184b11644fSPeter Brune 194b11644fSPeter Brune Vec * dX = qn->dX; 204b11644fSPeter Brune Vec * dF = qn->dF; 214b11644fSPeter Brune 224b11644fSPeter Brune PetscInt k, i; 234b11644fSPeter Brune PetscInt m = qn->m; 244b11644fSPeter Brune PetscScalar * alpha = PETSC_NULL; 254b11644fSPeter Brune PetscScalar * beta = PETSC_NULL; 264b11644fSPeter Brune PetscScalar * rho = PETSC_NULL; 274b11644fSPeter Brune PetscScalar t; 284b11644fSPeter Brune PetscInt l = m; 294b11644fSPeter Brune 304b11644fSPeter Brune PetscFunctionBegin; 314b11644fSPeter Brune 324b11644fSPeter Brune if (it < m) l = it+1; 334b11644fSPeter Brune ierr = PetscMalloc3(m, PetscScalar, &alpha, m, PetscScalar, &beta, m, PetscScalar, &rho);CHKERRQ(ierr); 344b11644fSPeter Brune 354b11644fSPeter Brune /* precalculate alpha, beta, rho corresponding to the normal indices*/ 364b11644fSPeter Brune for (i = 0; i < l; i++) { 374b11644fSPeter Brune ierr = VecDot(dX[i], dF[i], &t);CHKERRQ(ierr); 384b11644fSPeter Brune rho[i] = 1. / t; 394b11644fSPeter Brune beta[i] = 0; 404b11644fSPeter Brune alpha[i] = 0; 414b11644fSPeter Brune } 424b11644fSPeter Brune 434b11644fSPeter Brune ierr = VecCopy(g, z);CHKERRQ(ierr); 444b11644fSPeter Brune 454b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 464b11644fSPeter Brune for (i = 0; i < l; i++) { 474b11644fSPeter Brune k = (it - i) % m; 484b11644fSPeter Brune ierr = VecDot(dX[k], z, &t);CHKERRQ(ierr); 494b11644fSPeter Brune alpha[k] = t / rho[k]; 504b11644fSPeter Brune ierr = VecAXPY(z, -alpha[k], dF[k]);CHKERRQ(ierr); 514b11644fSPeter Brune } 524b11644fSPeter Brune 534b11644fSPeter Brune /* inner application of the initial inverse jacobian approximation */ 544b11644fSPeter Brune /* right now it's just the identity. Nothing needs to go here. */ 554b11644fSPeter Brune 564b11644fSPeter Brune /* inward recursion starting at the first update and working forward*/ 574b11644fSPeter Brune for (i = l - 1; i >= 0; i--) { 584b11644fSPeter Brune k = (it - i) % m; 594b11644fSPeter Brune ierr = VecDot(dF[k], z, &t);CHKERRQ(ierr); 604b11644fSPeter Brune beta[k] = rho[k]*t; 614b11644fSPeter Brune ierr = VecAXPY(z, (alpha[k] - beta[k]), dX[k]); 624b11644fSPeter Brune } 634b11644fSPeter Brune ierr = PetscFree3(alpha, beta, rho);CHKERRQ(ierr); 644b11644fSPeter Brune 654b11644fSPeter Brune PetscFunctionReturn(0); 664b11644fSPeter Brune } 674b11644fSPeter Brune 684b11644fSPeter Brune #undef __FUNCT__ 694b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 704b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 714b11644fSPeter Brune { 724b11644fSPeter Brune 734b11644fSPeter Brune PetscErrorCode ierr; 744b11644fSPeter Brune QNContext * qn = (QNContext*) snes->data; 754b11644fSPeter Brune 764b11644fSPeter Brune Vec x; 774b11644fSPeter Brune Vec f; 784b11644fSPeter Brune Vec p, pold; 794b11644fSPeter Brune 804b11644fSPeter Brune PetscInt i, j, k, l; 814b11644fSPeter Brune 824b11644fSPeter Brune PetscReal fnorm; 834b11644fSPeter Brune PetscScalar gdot; 844b11644fSPeter Brune PetscInt m = qn->m; 854b11644fSPeter Brune 864b11644fSPeter Brune Vec * W = qn->dX; 874b11644fSPeter Brune Vec * V = qn->dF; 884b11644fSPeter Brune 894b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 904b11644fSPeter Brune PetscFunctionBegin; 914b11644fSPeter Brune 924b11644fSPeter Brune x = snes->vec_sol; 934b11644fSPeter Brune f = snes->vec_func; 944b11644fSPeter Brune p = snes->vec_sol_update; 954b11644fSPeter Brune pold = snes->work[0]; 964b11644fSPeter Brune 974b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 984b11644fSPeter Brune 994b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1004b11644fSPeter Brune snes->iter = 0; 1014b11644fSPeter Brune snes->norm = 0.; 1024b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1034b11644fSPeter Brune ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); 1044b11644fSPeter Brune if (snes->domainerror) { 1054b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1064b11644fSPeter Brune PetscFunctionReturn(0); 1074b11644fSPeter Brune } 1084b11644fSPeter Brune ierr = VecNorm(f, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1094b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1104b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1114b11644fSPeter Brune snes->norm = fnorm; 1124b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1134b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 1144b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1154b11644fSPeter Brune 1164b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1174b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 1184b11644fSPeter Brune /* test convergence */ 1194b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1204b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 1214b11644fSPeter Brune ierr = VecCopy(f, pold);CHKERRQ(ierr); 1224b11644fSPeter Brune ierr = VecAXPY(x, -1.0, pold);CHKERRQ(ierr); 1234b11644fSPeter Brune for(i = 0; i < snes->max_its; i++) { 1244b11644fSPeter Brune 1254b11644fSPeter Brune ierr = SNESComputeFunction(snes, x, f);CHKERRQ(ierr); 1264b11644fSPeter Brune ierr = VecNorm(f, NORM_2, &fnorm);CHKERRQ(ierr); 1274b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1284b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1294b11644fSPeter Brune snes->norm = fnorm; 1304b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1314b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,i+1); 1324b11644fSPeter Brune ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 1334b11644fSPeter Brune 1344b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1354b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 1364b11644fSPeter Brune /* test convergence */ 1374b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1384b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 1394b11644fSPeter Brune k = (i) % m; 1404b11644fSPeter Brune l = m; 1414b11644fSPeter Brune if (i < l) l = i; 1424b11644fSPeter Brune ierr = VecCopy(f, p);CHKERRQ(ierr); 1434b11644fSPeter Brune for (j=0; j<k; j++) { /* p = product_{j<i} [I+v(j)w(j)^T]*p */ 1444b11644fSPeter Brune ierr = VecDot(W[j],p,&gdot);CHKERRQ(ierr); 1454b11644fSPeter Brune ierr = VecAXPY(p,gdot,V[j]);CHKERRQ(ierr); 1464b11644fSPeter Brune } 1474b11644fSPeter Brune ierr = VecCopy(pold,W[k]);CHKERRQ(ierr); /* w[i] = pold */ 1484b11644fSPeter Brune ierr = VecAXPY(pold,-1.0,p);CHKERRQ(ierr); /* v[i] = p */ 1494b11644fSPeter Brune ierr = VecDot(W[k],pold,&gdot);CHKERRQ(ierr); /* ----------------- */ 1504b11644fSPeter Brune ierr = VecCopy(p,V[k]);CHKERRQ(ierr); /* w[i]'*(Pold - p) */ 1514b11644fSPeter Brune ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr); 1524b11644fSPeter Brune 1534b11644fSPeter Brune ierr = VecDot(W[k],p,&gdot);CHKERRQ(ierr); /* p = (I + v[i]*w[i]')*p */ 1544b11644fSPeter Brune ierr = VecAXPY(p,gdot,V[k]);CHKERRQ(ierr); 1554b11644fSPeter Brune ierr = VecCopy(p,pold);CHKERRQ(ierr); 1564b11644fSPeter Brune 1574b11644fSPeter Brune ierr = VecAXPY(x, -1.0, p);CHKERRQ(ierr); 1584b11644fSPeter Brune } 1594b11644fSPeter Brune if (i == snes->max_its) { 1604b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 1614b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1624b11644fSPeter Brune } 1634b11644fSPeter Brune PetscFunctionReturn(0); 1644b11644fSPeter Brune } 1654b11644fSPeter Brune 1664b11644fSPeter Brune 1674b11644fSPeter Brune #undef __FUNCT__ 1684b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 1694b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 1704b11644fSPeter Brune { 1714b11644fSPeter Brune QNContext * qn = (QNContext *)snes->data; 1724b11644fSPeter Brune PetscErrorCode ierr; 1734b11644fSPeter Brune PetscFunctionBegin; 1744b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 1754b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 17635f34c5cSPeter Brune ierr = SNESDefaultGetWork(snes,1);CHKERRQ(ierr); 1774b11644fSPeter Brune PetscFunctionReturn(0); 1784b11644fSPeter Brune } 1794b11644fSPeter Brune 1804b11644fSPeter Brune #undef __FUNCT__ 1814b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 1824b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 1834b11644fSPeter Brune { 1844b11644fSPeter Brune PetscErrorCode ierr; 1854b11644fSPeter Brune QNContext * qn; 1864b11644fSPeter Brune PetscFunctionBegin; 1874b11644fSPeter Brune if (snes->data) { 1884b11644fSPeter Brune qn = (QNContext *)snes->data; 1894b11644fSPeter Brune if (qn->dX) { 1904b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 1914b11644fSPeter Brune } 1924b11644fSPeter Brune if (qn->dF) { 1934b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 1944b11644fSPeter Brune } 1954b11644fSPeter Brune } 1964b11644fSPeter Brune if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 1974b11644fSPeter Brune PetscFunctionReturn(0); 1984b11644fSPeter Brune } 1994b11644fSPeter Brune 2004b11644fSPeter Brune #undef __FUNCT__ 2014b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 2024b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 2034b11644fSPeter Brune { 2044b11644fSPeter Brune PetscErrorCode ierr; 2054b11644fSPeter Brune PetscFunctionBegin; 2064b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 2074b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 2084b11644fSPeter Brune PetscFunctionReturn(0); 2094b11644fSPeter Brune } 2104b11644fSPeter Brune 2114b11644fSPeter Brune #undef __FUNCT__ 2124b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 2134b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 2144b11644fSPeter Brune { 2154b11644fSPeter Brune 2164b11644fSPeter Brune PetscErrorCode ierr; 2174b11644fSPeter Brune QNContext * qn; 2184b11644fSPeter Brune 2194b11644fSPeter Brune PetscFunctionBegin; 2204b11644fSPeter Brune 2214b11644fSPeter Brune qn = (QNContext *)snes->data; 2224b11644fSPeter Brune 2234b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 2244b11644fSPeter Brune ierr = PetscOptionsReal("-snes_qn_lambda", "SOR mixing parameter", "SNES", qn->lambda, &qn->lambda, PETSC_NULL);CHKERRQ(ierr); 2254b11644fSPeter Brune ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 2264b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 2274b11644fSPeter Brune PetscFunctionReturn(0); 2284b11644fSPeter Brune } 2294b11644fSPeter Brune 2304b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 2314b11644fSPeter Brune /*MC 2324b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 2334b11644fSPeter Brune 23435f34c5cSPeter Brune Implements a limited-memory "good" Broyden update method. 2354b11644fSPeter Brune 2364b11644fSPeter Brune Level: beginner 2374b11644fSPeter Brune 2384b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 2394b11644fSPeter Brune M*/ 2404b11644fSPeter Brune EXTERN_C_BEGIN 2414b11644fSPeter Brune #undef __FUNCT__ 2424b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 2434b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 2444b11644fSPeter Brune { 2454b11644fSPeter Brune 2464b11644fSPeter Brune PetscErrorCode ierr; 2474b11644fSPeter Brune QNContext * qn; 2484b11644fSPeter Brune 2494b11644fSPeter Brune PetscFunctionBegin; 2504b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 2514b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 2524b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 2534b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 2544b11644fSPeter Brune snes->ops->view = 0; 2554b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 2564b11644fSPeter Brune 2574b11644fSPeter Brune ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr); 2584b11644fSPeter Brune snes->data = (void *) qn; 2594b11644fSPeter Brune qn->m = 10; 2604b11644fSPeter Brune qn->lambda = 1.; 2614b11644fSPeter Brune qn->dX = PETSC_NULL; 2624b11644fSPeter Brune qn->dF = PETSC_NULL; 2634b11644fSPeter Brune PetscFunctionReturn(0); 2644b11644fSPeter Brune } 2654b11644fSPeter Brune EXTERN_C_END 266