14b11644fSPeter Brune #include <private/snesimpl.h> 24b11644fSPeter Brune 34b11644fSPeter Brune typedef struct { 44b11644fSPeter Brune Vec *dX; /* The change in X */ 53af51624SPeter Brune Vec *dD; /* 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; 115ba6227bSPeter Brune PetscReal gamma; /* Powell restart constant */ 12b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 139f83bee8SJed Brown } SNES_QN; 144b11644fSPeter Brune 154b11644fSPeter Brune #undef __FUNCT__ 164b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private" 173af51624SPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) { 184b11644fSPeter Brune 194b11644fSPeter Brune PetscErrorCode ierr; 204b11644fSPeter Brune 219f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 224b11644fSPeter Brune 234b11644fSPeter Brune Vec *dX = qn->dX; 243af51624SPeter Brune Vec *dD = qn->dD; 254b11644fSPeter Brune 26bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 27bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 28bd052dfeSPeter Brune PetscScalar *rho = qn->rho; 29bd052dfeSPeter Brune 304b11644fSPeter Brune PetscInt k, i; 314b11644fSPeter Brune PetscInt m = qn->m; 324b11644fSPeter Brune PetscScalar t; 334b11644fSPeter Brune PetscInt l = m; 344b11644fSPeter Brune 354b11644fSPeter Brune PetscFunctionBegin; 364b11644fSPeter Brune 373af51624SPeter Brune ierr = VecCopy(D, Y);CHKERRQ(ierr); 384b11644fSPeter Brune 395ba6227bSPeter Brune if (it < m) l = it; 404b11644fSPeter Brune 414b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 424b11644fSPeter Brune for (i = 0; i < l; i++) { 43b21d5a53SPeter Brune k = (it - i - 1) % l; 443af51624SPeter Brune ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr); 45bd052dfeSPeter Brune alpha[k] = t*rho[k]; 4644f7e39eSPeter Brune if (qn->monitor) { 473af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 485ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 493af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 5044f7e39eSPeter Brune } 513af51624SPeter Brune ierr = VecAXPY(Y, -alpha[k], dD[k]);CHKERRQ(ierr); 524b11644fSPeter Brune } 534b11644fSPeter Brune 54b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 55b21d5a53SPeter Brune 564b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 57bd052dfeSPeter Brune for (i = 0; i < l; i++) { 58b21d5a53SPeter Brune k = (it + i - l) % l; 593af51624SPeter Brune ierr = VecDot(dD[k], Y, &t);CHKERRQ(ierr); 604b11644fSPeter Brune beta[k] = rho[k]*t; 613af51624SPeter Brune ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]); 6244f7e39eSPeter Brune if (qn->monitor) { 633af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 645ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 653af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 6644f7e39eSPeter Brune } 674b11644fSPeter Brune } 684b11644fSPeter Brune PetscFunctionReturn(0); 694b11644fSPeter Brune } 704b11644fSPeter Brune 714b11644fSPeter Brune #undef __FUNCT__ 724b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 734b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 744b11644fSPeter Brune { 754b11644fSPeter Brune 764b11644fSPeter Brune PetscErrorCode ierr; 779f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 784b11644fSPeter Brune 7915f5eeeaSPeter Brune Vec X, Xold; 803af51624SPeter Brune Vec F, G, B; 813af51624SPeter Brune Vec W, Y, D, Dold; 823af51624SPeter Brune SNESConvergedReason reason; 835ba6227bSPeter Brune PetscInt i, i_r, k; 844b11644fSPeter Brune 859f3a0142SPeter Brune PetscReal fnorm, xnorm = 0, ynorm, gnorm; 864b11644fSPeter Brune PetscInt m = qn->m; 879f3a0142SPeter Brune PetscBool lssucceed; 884b11644fSPeter Brune 89bd052dfeSPeter Brune PetscScalar rhosc; 90bd052dfeSPeter Brune 91bd052dfeSPeter Brune Vec *dX = qn->dX; 923af51624SPeter Brune Vec *dD = qn->dD; 93bd052dfeSPeter Brune PetscScalar *rho = qn->rho; 94b21d5a53SPeter Brune PetscScalar DolddotD, DolddotDold, a, fdotf; 954b11644fSPeter Brune 964b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 974b11644fSPeter Brune PetscFunctionBegin; 984b11644fSPeter Brune 999f3a0142SPeter Brune X = snes->vec_sol; /* solution vector */ 1009f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 1013af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 1023af51624SPeter Brune B = snes->vec_rhs; 10370d3b23bSPeter Brune G = snes->work[0]; 10470d3b23bSPeter Brune W = snes->work[1]; 10570d3b23bSPeter Brune Xold = snes->work[2]; 1063af51624SPeter Brune 1073af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 1083af51624SPeter Brune D = snes->work[3]; 1093af51624SPeter Brune Dold = snes->work[4]; 1104b11644fSPeter Brune 1114b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1124b11644fSPeter Brune 1134b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1144b11644fSPeter Brune snes->iter = 0; 1154b11644fSPeter Brune snes->norm = 0.; 1164b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 11715f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1184b11644fSPeter Brune if (snes->domainerror) { 1194b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1204b11644fSPeter Brune PetscFunctionReturn(0); 1214b11644fSPeter Brune } 12215f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1234b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1244b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1254b11644fSPeter Brune snes->norm = fnorm; 1264b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1274b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 1284b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1294b11644fSPeter Brune 1304b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1314b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 1324b11644fSPeter Brune /* test convergence */ 1334b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1344b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 13570d3b23bSPeter Brune 136b21d5a53SPeter Brune /* initialize the scaling */ 137b21d5a53SPeter Brune ierr = VecDot(X, F, &fdotf);CHKERRQ(ierr); 138*1b2f85ebSPeter Brune qn->scaling = PetscRealPart(fdotf); 139b21d5a53SPeter Brune ierr = VecDot(F, F, &fdotf);CHKERRQ(ierr); 140*1b2f85ebSPeter Brune qn->scaling = qn->scaling / PetscRealPart(fdotf); 141b21d5a53SPeter Brune 14270d3b23bSPeter Brune /* initialize the search direction as steepest descent */ 1433af51624SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 1443af51624SPeter Brune ierr = VecCopy(D, Y);CHKERRQ(ierr); 145b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 1463af51624SPeter Brune 1475ba6227bSPeter Brune for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 14870d3b23bSPeter Brune /* line search for lambda */ 14970d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 1503af51624SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 15115f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 1529f3a0142SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1539f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 1549f3a0142SPeter Brune if (snes->domainerror) { 1559f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1569f3a0142SPeter Brune PetscFunctionReturn(0); 1579f3a0142SPeter Brune } 1589f3a0142SPeter Brune if (!lssucceed) { 1599f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1609f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 1619f3a0142SPeter Brune break; 1629f3a0142SPeter Brune } 1639f3a0142SPeter Brune } 164b21d5a53SPeter Brune 1659f3a0142SPeter Brune /* Update function and solution vectors */ 1669f3a0142SPeter Brune fnorm = gnorm; 1679f3a0142SPeter Brune ierr = VecCopy(G,F);CHKERRQ(ierr); 16815f5eeeaSPeter Brune ierr = VecCopy(W,X);CHKERRQ(ierr); 1693af51624SPeter Brune 170b21d5a53SPeter Brune if (snes->pc) { 171b21d5a53SPeter Brune ierr = SNESSolve(snes->pc, snes->vec_rhs, X);CHKERRQ(ierr); 172b21d5a53SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 173b21d5a53SPeter Brune ierr = VecNorm(F, NORM_2, &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 1885ba6227bSPeter Brune /* create the new direction */ 189b21d5a53SPeter Brune if (0 && snes->pc) { 1903af51624SPeter Brune ierr = VecCopy(X, D);CHKERRQ(ierr); 1913af51624SPeter Brune ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr); 1923af51624SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 1933af51624SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 1943af51624SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1953af51624SPeter Brune PetscFunctionReturn(0); 1963af51624SPeter Brune } 1973af51624SPeter Brune ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr); 1983af51624SPeter Brune } else { 1993af51624SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 2003af51624SPeter Brune } 2013af51624SPeter Brune 2025ba6227bSPeter Brune /* check restart by Powell's Criterion: |D^T H_0 Dold| > 0.2 * |Dold^T H_0 Dold| */ 2035ba6227bSPeter Brune ierr = VecDot(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 2045ba6227bSPeter Brune ierr = VecDot(Dold, D, &DolddotD);CHKERRQ(ierr); 205b21d5a53SPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->gamma*PetscAbs(PetscRealPart(DolddotDold)) && i_r > 2) { 2065ba6227bSPeter Brune if (qn->monitor) { 2075ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2085ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e|\n", k, PetscRealPart(DolddotD), qn->gamma, PetscRealPart(DolddotDold));CHKERRQ(ierr); 2095ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2105ba6227bSPeter Brune } 2115ba6227bSPeter Brune i_r = -1; 2125ba6227bSPeter Brune /* general purpose update */ 2135ba6227bSPeter Brune if (snes->ops->update) { 2145ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 2155ba6227bSPeter Brune } 2165ba6227bSPeter Brune ierr = VecCopy(D, Y);CHKERRQ(ierr); 2175ba6227bSPeter Brune } else { 218bd052dfeSPeter Brune /* set the differences */ 2195ba6227bSPeter Brune k = i_r % m; 2203af51624SPeter Brune ierr = VecCopy(D, dD[k]);CHKERRQ(ierr); 2213af51624SPeter Brune ierr = VecAXPY(dD[k], -1.0, Dold);CHKERRQ(ierr); 22215f5eeeaSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 22315f5eeeaSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 2243af51624SPeter Brune ierr = VecDot(dX[k], dD[k], &rhosc);CHKERRQ(ierr); 225bd052dfeSPeter Brune rho[k] = 1. / rhosc; 226b21d5a53SPeter Brune ierr = VecDot(dD[k], dD[k], &a);CHKERRQ(ierr); 227*1b2f85ebSPeter Brune qn->scaling = PetscRealPart(rhosc) / PetscRealPart(a); 22870d3b23bSPeter Brune 22970d3b23bSPeter Brune /* general purpose update */ 23070d3b23bSPeter Brune if (snes->ops->update) { 23170d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 23270d3b23bSPeter Brune } 23370d3b23bSPeter Brune /* apply the current iteration of the approximate jacobian in order to get the next search direction*/ 2345ba6227bSPeter Brune ierr = LBGFSApplyJinv_Private(snes, i_r+1, D, Y);CHKERRQ(ierr); 2355ba6227bSPeter Brune } 2364b11644fSPeter Brune } 2374b11644fSPeter Brune if (i == snes->max_its) { 2384b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 2394b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 2404b11644fSPeter Brune } 2414b11644fSPeter Brune PetscFunctionReturn(0); 2424b11644fSPeter Brune } 2434b11644fSPeter Brune 2444b11644fSPeter Brune 2454b11644fSPeter Brune #undef __FUNCT__ 2464b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 2474b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 2484b11644fSPeter Brune { 2499f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 2504b11644fSPeter Brune PetscErrorCode ierr; 2514b11644fSPeter Brune PetscFunctionBegin; 2524b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 2533af51624SPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dD);CHKERRQ(ierr); 254bd052dfeSPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr); 2553af51624SPeter Brune ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr); 2564b11644fSPeter Brune PetscFunctionReturn(0); 2574b11644fSPeter Brune } 2584b11644fSPeter Brune 2594b11644fSPeter Brune #undef __FUNCT__ 2604b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 2614b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 2624b11644fSPeter Brune { 2634b11644fSPeter Brune PetscErrorCode ierr; 2649f83bee8SJed Brown SNES_QN *qn; 2654b11644fSPeter Brune PetscFunctionBegin; 2664b11644fSPeter Brune if (snes->data) { 2679f83bee8SJed Brown qn = (SNES_QN*)snes->data; 2684b11644fSPeter Brune if (qn->dX) { 2694b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 2704b11644fSPeter Brune } 2713af51624SPeter Brune if (qn->dD) { 2723af51624SPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dD);CHKERRQ(ierr); 2734b11644fSPeter Brune } 274bd052dfeSPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr); 2754b11644fSPeter Brune } 2764b11644fSPeter Brune PetscFunctionReturn(0); 2774b11644fSPeter Brune } 2784b11644fSPeter Brune 2794b11644fSPeter Brune #undef __FUNCT__ 2804b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 2814b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 2824b11644fSPeter Brune { 2834b11644fSPeter Brune PetscErrorCode ierr; 2844b11644fSPeter Brune PetscFunctionBegin; 2854b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 2864b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 2879c05af5eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 2884b11644fSPeter Brune PetscFunctionReturn(0); 2894b11644fSPeter Brune } 2904b11644fSPeter Brune 2914b11644fSPeter Brune #undef __FUNCT__ 2924b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 2934b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 2944b11644fSPeter Brune { 2954b11644fSPeter Brune 2964b11644fSPeter Brune PetscErrorCode ierr; 2979f83bee8SJed Brown SNES_QN *qn; 29844f7e39eSPeter Brune PetscBool monflg = PETSC_FALSE; 2994b11644fSPeter Brune PetscFunctionBegin; 3004b11644fSPeter Brune 3019f83bee8SJed Brown qn = (SNES_QN*)snes->data; 3024b11644fSPeter Brune 3034b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 3044b11644fSPeter Brune ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 3055ba6227bSPeter Brune ierr = PetscOptionsReal("-snes_qn_gamma", "Restart condition for L-Broyden methods", "SNES", qn->gamma, &qn->gamma, PETSC_NULL);CHKERRQ(ierr); 30644f7e39eSPeter Brune ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNES", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 3074b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 30844f7e39eSPeter Brune if (monflg) { 30944f7e39eSPeter Brune qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 31044f7e39eSPeter Brune } 3114b11644fSPeter Brune PetscFunctionReturn(0); 3124b11644fSPeter Brune } 3134b11644fSPeter Brune 31492c02d66SPeter Brune EXTERN_C_BEGIN 31592c02d66SPeter Brune #undef __FUNCT__ 31692c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_QN" 31792c02d66SPeter Brune PetscErrorCode SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type) 31892c02d66SPeter Brune { 31992c02d66SPeter Brune PetscErrorCode ierr; 32092c02d66SPeter Brune PetscFunctionBegin; 32192c02d66SPeter Brune 32292c02d66SPeter Brune switch (type) { 32392c02d66SPeter Brune case SNES_LS_BASIC: 32492c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 32592c02d66SPeter Brune break; 32692c02d66SPeter Brune case SNES_LS_BASIC_NONORMS: 32792c02d66SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 32892c02d66SPeter Brune break; 32992c02d66SPeter Brune case SNES_LS_QUADRATIC: 330af60355fSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 33192c02d66SPeter Brune break; 332cf9edfecSPeter Brune case SNES_LS_SECANT: 333cf9edfecSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr); 334cf9edfecSPeter Brune break; 33592c02d66SPeter Brune default: 33692c02d66SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 33792c02d66SPeter Brune break; 33892c02d66SPeter Brune } 33992c02d66SPeter Brune snes->ls_type = type; 34092c02d66SPeter Brune PetscFunctionReturn(0); 34192c02d66SPeter Brune } 34292c02d66SPeter Brune EXTERN_C_END 34392c02d66SPeter Brune 34492c02d66SPeter Brune 3454b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 3464b11644fSPeter Brune /*MC 3474b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 3484b11644fSPeter Brune 3496cc8130cSPeter Brune Options Database: 3506cc8130cSPeter Brune 3516cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 35244f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 3536cc8130cSPeter Brune 35444f7e39eSPeter Brune Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to 3556cc8130cSPeter Brune form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 3566cc8130cSPeter Brune generalized to implement several limited-memory Broyden methods. 3576cc8130cSPeter Brune 3586cc8130cSPeter Brune References: 3596cc8130cSPeter Brune 3606cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 3616cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 3626cc8130cSPeter Brune 3634b11644fSPeter Brune 3644b11644fSPeter Brune Level: beginner 3654b11644fSPeter Brune 3664b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 3676cc8130cSPeter Brune 3684b11644fSPeter Brune M*/ 3694b11644fSPeter Brune EXTERN_C_BEGIN 3704b11644fSPeter Brune #undef __FUNCT__ 3714b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 3724b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 3734b11644fSPeter Brune { 3744b11644fSPeter Brune 3754b11644fSPeter Brune PetscErrorCode ierr; 3769f83bee8SJed Brown SNES_QN *qn; 3774b11644fSPeter Brune 3784b11644fSPeter Brune PetscFunctionBegin; 3794b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 3804b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 3814b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 3824b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 3834b11644fSPeter Brune snes->ops->view = 0; 3844b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 3854b11644fSPeter Brune 38642f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 38742f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 38842f4f86dSBarry Smith 3890e444f03SPeter Brune snes->max_funcs = 30000; 3900e444f03SPeter Brune snes->max_its = 10000; 3910e444f03SPeter Brune 3929f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 3934b11644fSPeter Brune snes->data = (void *) qn; 39470d3b23bSPeter Brune qn->m = 10; 395b21d5a53SPeter Brune qn->scaling = 1.0; 3964b11644fSPeter Brune qn->dX = PETSC_NULL; 3973af51624SPeter Brune qn->dD = PETSC_NULL; 39844f7e39eSPeter Brune qn->monitor = PETSC_NULL; 399b21d5a53SPeter Brune qn->gamma = 0.5; 400ea630c6eSPeter Brune 40192c02d66SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr); 402b21d5a53SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_SECANT);CHKERRQ(ierr); 4039f3a0142SPeter Brune 4044b11644fSPeter Brune PetscFunctionReturn(0); 4054b11644fSPeter Brune } 4064b11644fSPeter Brune EXTERN_C_END 407