14b11644fSPeter Brune #include <private/snesimpl.h> 2335fb975SPeter Brune #include <petsclinesearch.h> 34b11644fSPeter Brune 4*8e231d97SPeter Brune #define H(i,j) qn->dXdFmat[i*qn->m + j] 5*8e231d97SPeter Brune 688d374b2SPeter Brune typedef enum {SNES_QN_SEQUENTIAL, SNES_QN_COMPOSED} SNESQNCompositionType; 70ecc9e1dSPeter Brune typedef enum {SNES_QN_SHANNOSCALE, SNES_QN_LSSCALE, SNES_QN_JACOBIANSCALE} SNESQNScalingType; 86bf1b2e5SPeter Brune 94b11644fSPeter Brune typedef struct { 104b11644fSPeter Brune Vec *dX; /* The change in X */ 116bf1b2e5SPeter Brune Vec *dF; /* The change in F */ 12*8e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 13*8e231d97SPeter Brune PetscScalar *alpha, *beta; 14*8e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 15*8e231d97SPeter Brune PetscBool aggreduct; /* Aggregated reduction implementation */ 16*8e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 1744f7e39eSPeter Brune PetscViewer monitor; 186bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 196bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 20b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 210ecc9e1dSPeter Brune PetscInt n_restart; /* the maximum iterations between restart */ 2288d374b2SPeter Brune 23335fb975SPeter Brune LineSearch linesearch; /* line search */ 24335fb975SPeter Brune 2588d374b2SPeter Brune SNESQNCompositionType compositiontype; /* determine if the composition is done sequentially or as a composition */ 260ecc9e1dSPeter Brune SNESQNScalingType scalingtype; /* determine if the composition is done sequentially or as a composition */ 2788d374b2SPeter Brune 289f83bee8SJed Brown } SNES_QN; 294b11644fSPeter Brune 304b11644fSPeter Brune #undef __FUNCT__ 314b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private" 323af51624SPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) { 334b11644fSPeter Brune 344b11644fSPeter Brune PetscErrorCode ierr; 354b11644fSPeter Brune 369f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 374b11644fSPeter Brune 38335fb975SPeter Brune Vec Yin = snes->work[3]; 390ecc9e1dSPeter Brune 404b11644fSPeter Brune Vec *dX = qn->dX; 416bf1b2e5SPeter Brune Vec *dF = qn->dF; 424b11644fSPeter Brune 43bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 44bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 45*8e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 46*8e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 47bd052dfeSPeter Brune 480ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 490ecc9e1dSPeter Brune KSPConvergedReason kspreason; 500ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 510ecc9e1dSPeter Brune 52*8e231d97SPeter Brune PetscInt k, i, j, g, lits; 534b11644fSPeter Brune PetscInt m = qn->m; 544b11644fSPeter Brune PetscScalar t; 554b11644fSPeter Brune PetscInt l = m; 564b11644fSPeter Brune 57*8e231d97SPeter Brune Mat jac, jac_pre; 58*8e231d97SPeter Brune 594b11644fSPeter Brune PetscFunctionBegin; 604b11644fSPeter Brune 613af51624SPeter Brune ierr = VecCopy(D, Y);CHKERRQ(ierr); 624b11644fSPeter Brune 635ba6227bSPeter Brune if (it < m) l = it; 644b11644fSPeter Brune 65*8e231d97SPeter Brune if (qn->aggreduct) { 66*8e231d97SPeter Brune ierr = VecMDot(Y, l, qn->dX, YtdX);CHKERRQ(ierr); 67*8e231d97SPeter Brune } 684b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 694b11644fSPeter Brune for (i = 0; i < l; i++) { 70b21d5a53SPeter Brune k = (it - i - 1) % l; 71*8e231d97SPeter Brune if (qn->aggreduct) { 72*8e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 73*8e231d97SPeter Brune t = YtdX[k]; 74*8e231d97SPeter Brune for (j = 0; j < i; j++) { 75*8e231d97SPeter Brune g = (it - j - 1) % l; 76*8e231d97SPeter Brune t += -alpha[g]*H(g, k); 77*8e231d97SPeter Brune } 78*8e231d97SPeter Brune alpha[k] = t / H(k, k); 79*8e231d97SPeter Brune } else { 803af51624SPeter Brune ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr); 81*8e231d97SPeter Brune alpha[k] = t / dXtdF[k]; 82*8e231d97SPeter Brune } 8344f7e39eSPeter Brune if (qn->monitor) { 843af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 855ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 863af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 8744f7e39eSPeter Brune } 886bf1b2e5SPeter Brune ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr); 894b11644fSPeter Brune } 904b11644fSPeter Brune 910ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_JACOBIANSCALE) { 92*8e231d97SPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 93*8e231d97SPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 940ecc9e1dSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,Yin);CHKERRQ(ierr); 950ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 960ecc9e1dSPeter Brune if (kspreason < 0) { 970ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 980ecc9e1dSPeter Brune ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 990ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 1000ecc9e1dSPeter Brune PetscFunctionReturn(0); 1010ecc9e1dSPeter Brune } 1020ecc9e1dSPeter Brune } 1030ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1040ecc9e1dSPeter Brune snes->linear_its += lits; 1050ecc9e1dSPeter Brune ierr = VecCopy(Yin, Y);CHKERRQ(ierr); 1060ecc9e1dSPeter Brune } else { 107b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 1080ecc9e1dSPeter Brune } 109*8e231d97SPeter Brune if (qn->aggreduct) { 110*8e231d97SPeter Brune ierr = VecMDot(Y, l, qn->dF, YtdX);CHKERRQ(ierr); 111*8e231d97SPeter Brune } 1124b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 113bd052dfeSPeter Brune for (i = 0; i < l; i++) { 114b21d5a53SPeter Brune k = (it + i - l) % l; 115*8e231d97SPeter Brune if (qn->aggreduct) { 116*8e231d97SPeter Brune t = YtdX[k]; 117*8e231d97SPeter Brune for (j = 0; j < i; j++) { 118*8e231d97SPeter Brune g = (it + j - l) % l; 119*8e231d97SPeter Brune t += (alpha[g] - beta[g])*H(k, g); 120*8e231d97SPeter Brune } 121*8e231d97SPeter Brune beta[k] = t / H(k, k); 122*8e231d97SPeter Brune } else { 1236bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 124*8e231d97SPeter Brune beta[k] = t / dXtdF[k]; 125*8e231d97SPeter Brune } 1263af51624SPeter Brune ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]); 12744f7e39eSPeter Brune if (qn->monitor) { 1283af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1295ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 1303af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 13144f7e39eSPeter Brune } 1324b11644fSPeter Brune } 1334b11644fSPeter Brune PetscFunctionReturn(0); 1344b11644fSPeter Brune } 1354b11644fSPeter Brune 1364b11644fSPeter Brune #undef __FUNCT__ 1374b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 1384b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 1394b11644fSPeter Brune { 1404b11644fSPeter Brune 1414b11644fSPeter Brune PetscErrorCode ierr; 1429f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 1434b11644fSPeter Brune 14415f5eeeaSPeter Brune Vec X, Xold; 145335fb975SPeter Brune Vec F, B; 146335fb975SPeter Brune Vec Y, FPC, D, Dold; 1473af51624SPeter Brune SNESConvergedReason reason; 148*8e231d97SPeter Brune PetscInt i, i_r, k, l, j; 1494b11644fSPeter Brune 150335fb975SPeter Brune PetscReal fnorm, xnorm, ynorm, gnorm; 1514b11644fSPeter Brune PetscInt m = qn->m; 152335fb975SPeter Brune PetscBool lssucceed; 1534b11644fSPeter Brune 154bd052dfeSPeter Brune Vec *dX = qn->dX; 1556bf1b2e5SPeter Brune Vec *dF = qn->dF; 156*8e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 157*8e231d97SPeter Brune PetscScalar *dFtdX = qn->dFtdX; 15888d374b2SPeter Brune PetscScalar DolddotD, DolddotDold, DdotD, YdotD, a; 1594b11644fSPeter Brune 1600ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 1610ecc9e1dSPeter Brune 1624b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 1634b11644fSPeter Brune PetscFunctionBegin; 1644b11644fSPeter Brune 1659f3a0142SPeter Brune X = snes->vec_sol; /* solution vector */ 1669f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 1673af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 1683af51624SPeter Brune B = snes->vec_rhs; 169335fb975SPeter Brune Xold = snes->work[0]; 1703af51624SPeter Brune 1713af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 172335fb975SPeter Brune D = snes->work[1]; 173335fb975SPeter Brune Dold = snes->work[2]; 1744b11644fSPeter Brune 1754b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1764b11644fSPeter Brune 1774b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1784b11644fSPeter Brune snes->iter = 0; 1794b11644fSPeter Brune snes->norm = 0.; 1804b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 18115f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1824b11644fSPeter Brune if (snes->domainerror) { 1834b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1844b11644fSPeter Brune PetscFunctionReturn(0); 1854b11644fSPeter Brune } 18615f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1874b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1884b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1894b11644fSPeter Brune snes->norm = fnorm; 1904b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1914b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 1924b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1934b11644fSPeter Brune 1944b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1954b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 1964b11644fSPeter Brune /* test convergence */ 1974b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1984b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 19970d3b23bSPeter Brune 20088d374b2SPeter Brune /* composed solve -- either sequential or composed */ 20188d374b2SPeter Brune if (snes->pc) { 20288d374b2SPeter Brune if (qn->compositiontype == SNES_QN_SEQUENTIAL) { 20388d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 20488d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 20588d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 20688d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 20788d374b2SPeter Brune PetscFunctionReturn(0); 20888d374b2SPeter Brune } 20988d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 21088d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 21188d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 2126bf1b2e5SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 21388d374b2SPeter Brune } else { 21488d374b2SPeter Brune ierr = VecCopy(X, Y);CHKERRQ(ierr); 21588d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, Y);CHKERRQ(ierr); 21688d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 21788d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 21888d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 21988d374b2SPeter Brune PetscFunctionReturn(0); 22088d374b2SPeter Brune } 22188d374b2SPeter Brune ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 22288d374b2SPeter Brune } 22388d374b2SPeter Brune } else { 22488d374b2SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 22588d374b2SPeter Brune } 22688d374b2SPeter Brune ierr = VecCopy(Y, D);CHKERRQ(ierr); 2273af51624SPeter Brune 228f8e15203SPeter Brune /* scale the initial update */ 2290ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_JACOBIANSCALE) { 2300ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 2310ecc9e1dSPeter Brune } 2320ecc9e1dSPeter Brune 2335ba6227bSPeter Brune for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 234f8e15203SPeter Brune ierr = LBGFSApplyJinv_Private(snes, i_r, D, Y);CHKERRQ(ierr); 23570d3b23bSPeter Brune /* line search for lambda */ 23670d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 23788d374b2SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 23815f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 239335fb975SPeter Brune ierr = LineSearchApply(qn->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 2409f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 2419f3a0142SPeter Brune if (snes->domainerror) { 2429f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2439f3a0142SPeter Brune PetscFunctionReturn(0); 2449f3a0142SPeter Brune } 245335fb975SPeter Brune ierr = LineSearchGetSuccess(qn->linesearch, &lssucceed);CHKERRQ(ierr); 2469f3a0142SPeter Brune if (!lssucceed) { 2479f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 2489f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 2499f3a0142SPeter Brune break; 2509f3a0142SPeter Brune } 2519f3a0142SPeter Brune } 252335fb975SPeter Brune ierr = LineSearchGetNorms(qn->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 2530ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_LSSCALE) { 254335fb975SPeter Brune ierr = LineSearchGetLambda(qn->linesearch, &qn->scaling);CHKERRQ(ierr); 2550ecc9e1dSPeter Brune } 2563af51624SPeter Brune 25788d374b2SPeter Brune /* convergence monitoring */ 258b21d5a53SPeter 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); 259b21d5a53SPeter Brune 260360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 261360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 262360c497dSPeter Brune 2638409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 2648409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 2654b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 2665ba6227bSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2674b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 2684b11644fSPeter Brune 26988d374b2SPeter Brune 27088d374b2SPeter Brune if (snes->pc) { 27188d374b2SPeter Brune if (qn->compositiontype == SNES_QN_SEQUENTIAL) { 27288d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 27388d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 27488d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 27588d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 27688d374b2SPeter Brune PetscFunctionReturn(0); 27788d374b2SPeter Brune } 27888d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 27988d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 28088d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 28188d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 28288d374b2SPeter Brune } else { 28388d374b2SPeter Brune ierr = VecCopy(X, D);CHKERRQ(ierr); 28488d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr); 28588d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 28688d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 28788d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 28888d374b2SPeter Brune PetscFunctionReturn(0); 28988d374b2SPeter Brune } 29088d374b2SPeter Brune ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr); 29188d374b2SPeter Brune } 29288d374b2SPeter Brune } else { 29388d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 29488d374b2SPeter Brune } 29588d374b2SPeter Brune 2966bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 297*8e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 298*8e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 299*8e231d97SPeter Brune ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr); 300*8e231d97SPeter Brune ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr); 301*8e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 302*8e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 303*8e231d97SPeter Brune ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr); 304*8e231d97SPeter Brune ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr); 3050ecc9e1dSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold)) || (i_r > qn->n_restart - 1 && qn->n_restart > 0)) { 3065ba6227bSPeter Brune if (qn->monitor) { 3075ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3080ecc9e1dSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", k, PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);CHKERRQ(ierr); 3095ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3105ba6227bSPeter Brune } 3115ba6227bSPeter Brune i_r = -1; 3125ba6227bSPeter Brune /* general purpose update */ 3135ba6227bSPeter Brune if (snes->ops->update) { 3145ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 3155ba6227bSPeter Brune } 3160ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_JACOBIANSCALE) { 3170ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 3180ecc9e1dSPeter Brune } 3195ba6227bSPeter Brune } else { 320bd052dfeSPeter Brune /* set the differences */ 3215ba6227bSPeter Brune k = i_r % m; 322*8e231d97SPeter Brune l = m; 323*8e231d97SPeter Brune if (i_r + 1 < m) l = i_r + 1; 32488d374b2SPeter Brune ierr = VecCopy(D, dF[k]);CHKERRQ(ierr); 32588d374b2SPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 32615f5eeeaSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 32715f5eeeaSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 328*8e231d97SPeter Brune if (qn->aggreduct) { 329*8e231d97SPeter Brune ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr); 330*8e231d97SPeter Brune ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr); 331*8e231d97SPeter Brune for (j = 0; j < l; j++) { 332*8e231d97SPeter Brune H(k, j) = dFtdX[j]; 333*8e231d97SPeter Brune H(j, k) = dXtdF[j]; 334*8e231d97SPeter Brune } 335*8e231d97SPeter Brune /* copy back over to make the computation of alpha and beta easier */ 336*8e231d97SPeter Brune for (j = 0; j < l; j++) { 337*8e231d97SPeter Brune dXtdF[j] = H(j, j); 338*8e231d97SPeter Brune } 339*8e231d97SPeter Brune } else { 340*8e231d97SPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 341*8e231d97SPeter Brune } 3426bf1b2e5SPeter Brune /* set scaling to be shanno scaling */ 3430ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_SHANNOSCALE) { 3446bf1b2e5SPeter Brune ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr); 345*8e231d97SPeter Brune qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(a); 3460ecc9e1dSPeter Brune } 34770d3b23bSPeter Brune /* general purpose update */ 34870d3b23bSPeter Brune if (snes->ops->update) { 34970d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 35070d3b23bSPeter Brune } 3515ba6227bSPeter Brune } 3524b11644fSPeter Brune } 3534b11644fSPeter Brune if (i == snes->max_its) { 3544b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 3554b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 3564b11644fSPeter Brune } 3574b11644fSPeter Brune PetscFunctionReturn(0); 3584b11644fSPeter Brune } 3594b11644fSPeter Brune 3604b11644fSPeter Brune 3614b11644fSPeter Brune #undef __FUNCT__ 3624b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 3634b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 3644b11644fSPeter Brune { 3659f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 3664b11644fSPeter Brune PetscErrorCode ierr; 367335fb975SPeter Brune const char *optionsprefix; 368335fb975SPeter Brune 3694b11644fSPeter Brune PetscFunctionBegin; 3704b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 3716bf1b2e5SPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 372*8e231d97SPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 373*8e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 374*8e231d97SPeter Brune qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 375*8e231d97SPeter Brune 376*8e231d97SPeter Brune if (qn->aggreduct) { 377*8e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 378*8e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 379*8e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 380*8e231d97SPeter Brune } 381335fb975SPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 382335fb975SPeter Brune 383335fb975SPeter Brune /* set up the line search */ 384335fb975SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 385335fb975SPeter Brune ierr = LineSearchCreate(((PetscObject)snes)->comm, &qn->linesearch);CHKERRQ(ierr); 386335fb975SPeter Brune ierr = LineSearchSetSNES(qn->linesearch, snes);CHKERRQ(ierr); 387335fb975SPeter Brune if (!snes->pc || qn->compositiontype == SNES_QN_SEQUENTIAL) { 388335fb975SPeter Brune ierr = LineSearchSetType(qn->linesearch, LINESEARCHCP);CHKERRQ(ierr); 389335fb975SPeter Brune } else { 390335fb975SPeter Brune ierr = LineSearchSetType(qn->linesearch, LINESEARCHL2);CHKERRQ(ierr); 391335fb975SPeter Brune } 392335fb975SPeter Brune ierr = LineSearchAppendOptionsPrefix(qn->linesearch, optionsprefix);CHKERRQ(ierr); 393335fb975SPeter Brune ierr = LineSearchSetFromOptions(qn->linesearch);CHKERRQ(ierr); 394*8e231d97SPeter Brune if (qn->scalingtype == SNES_QN_JACOBIANSCALE) { 395*8e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 396*8e231d97SPeter Brune } 3974b11644fSPeter Brune PetscFunctionReturn(0); 3984b11644fSPeter Brune } 3994b11644fSPeter Brune 4004b11644fSPeter Brune #undef __FUNCT__ 4014b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 4024b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 4034b11644fSPeter Brune { 4044b11644fSPeter Brune PetscErrorCode ierr; 4059f83bee8SJed Brown SNES_QN *qn; 4064b11644fSPeter Brune PetscFunctionBegin; 4074b11644fSPeter Brune if (snes->data) { 4089f83bee8SJed Brown qn = (SNES_QN*)snes->data; 4094b11644fSPeter Brune if (qn->dX) { 4104b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 4114b11644fSPeter Brune } 4126bf1b2e5SPeter Brune if (qn->dF) { 4136bf1b2e5SPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 4144b11644fSPeter Brune } 415*8e231d97SPeter Brune if (qn->aggreduct) { 416*8e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 417*8e231d97SPeter Brune } 418*8e231d97SPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 419*8e231d97SPeter Brune ierr = LineSearchDestroy(&qn->linesearch);CHKERRQ(ierr); 4204b11644fSPeter Brune } 4214b11644fSPeter Brune PetscFunctionReturn(0); 4224b11644fSPeter Brune } 4234b11644fSPeter Brune 4244b11644fSPeter Brune #undef __FUNCT__ 4254b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 4264b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 4274b11644fSPeter Brune { 4284b11644fSPeter Brune PetscErrorCode ierr; 4294b11644fSPeter Brune PetscFunctionBegin; 4304b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 4314b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 4329c05af5eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 4334b11644fSPeter Brune PetscFunctionReturn(0); 4344b11644fSPeter Brune } 4354b11644fSPeter Brune 4364b11644fSPeter Brune #undef __FUNCT__ 4374b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 4384b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 4394b11644fSPeter Brune { 4404b11644fSPeter Brune 4414b11644fSPeter Brune PetscErrorCode ierr; 4429f83bee8SJed Brown SNES_QN *qn; 4435b4627e8SPeter Brune const char *compositions[] = {"sequential", "composed"}; 4440ecc9e1dSPeter Brune const char *scalings[] = {"shanno", "ls", "jacobian"}; 44588d374b2SPeter Brune PetscInt indx = 0; 44688d374b2SPeter Brune PetscBool flg; 44744f7e39eSPeter Brune PetscBool monflg = PETSC_FALSE; 4484b11644fSPeter Brune PetscFunctionBegin; 4494b11644fSPeter Brune 4509f83bee8SJed Brown qn = (SNES_QN*)snes->data; 4514b11644fSPeter Brune 4524b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 4535b4627e8SPeter Brune ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 4540ecc9e1dSPeter Brune ierr = PetscOptionsInt("-snes_qn_restart", "Maximum number of iterations between restarts", "SNESQN", qn->n_restart, &qn->n_restart, PETSC_NULL);CHKERRQ(ierr); 4555b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr); 4565b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_downhill", "Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr); 4575b4627e8SPeter Brune ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 458*8e231d97SPeter Brune ierr = PetscOptionsBool("-snes_qn_aggreduct", "Aggregate reductions", "SNESQN", qn->aggreduct, &qn->aggreduct, PETSC_NULL);CHKERRQ(ierr); 4595b4627e8SPeter Brune ierr = PetscOptionsEList("-snes_qn_composition", "Composition type", "SNESQN",compositions,2,"sequential",&indx,&flg);CHKERRQ(ierr); 46088d374b2SPeter Brune if (flg) { 46188d374b2SPeter Brune switch (indx) { 4625b4627e8SPeter Brune case 0: qn->compositiontype = SNES_QN_SEQUENTIAL; 4635b4627e8SPeter Brune break; 4645b4627e8SPeter Brune case 1: qn->compositiontype = SNES_QN_COMPOSED; 4655b4627e8SPeter Brune break; 46688d374b2SPeter Brune } 46788d374b2SPeter Brune } 4680ecc9e1dSPeter Brune ierr = PetscOptionsEList("-snes_qn_scaling", "Scaling type", "SNESQN",scalings,3,"shanno",&indx,&flg);CHKERRQ(ierr); 4690ecc9e1dSPeter Brune if (flg) { 4700ecc9e1dSPeter Brune switch (indx) { 4710ecc9e1dSPeter Brune case 0: qn->scalingtype = SNES_QN_SHANNOSCALE; 4720ecc9e1dSPeter Brune break; 4730ecc9e1dSPeter Brune case 1: qn->scalingtype = SNES_QN_LSSCALE; 4740ecc9e1dSPeter Brune break; 4750ecc9e1dSPeter Brune case 2: qn->scalingtype = SNES_QN_JACOBIANSCALE; 4760ecc9e1dSPeter Brune snes->usesksp = PETSC_TRUE; 4770ecc9e1dSPeter Brune break; 4780ecc9e1dSPeter Brune } 4790ecc9e1dSPeter Brune } 4800ecc9e1dSPeter Brune 4814b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 48244f7e39eSPeter Brune if (monflg) { 48344f7e39eSPeter Brune qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 48444f7e39eSPeter Brune } 4854b11644fSPeter Brune PetscFunctionReturn(0); 4864b11644fSPeter Brune } 4874b11644fSPeter Brune 48892c02d66SPeter Brune 4894b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 4904b11644fSPeter Brune /*MC 4914b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 4924b11644fSPeter Brune 4936cc8130cSPeter Brune Options Database: 4946cc8130cSPeter Brune 4956cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 4961867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 4971867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 4981867fe5bSPeter Brune . -snes_qn_composition <sequential, composed>- Type of composition. 499f3aaf736SPeter Brune . -snes_ls <basic, basicnonorms, quadratic, critical> - Type of line search. 50044f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 5016cc8130cSPeter Brune 50244f7e39eSPeter Brune Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to 5036cc8130cSPeter Brune form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 5046cc8130cSPeter Brune generalized to implement several limited-memory Broyden methods. 5056cc8130cSPeter Brune 5061867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 5071867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 5081867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 5091867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 5101867fe5bSPeter Brune 5116cc8130cSPeter Brune References: 5126cc8130cSPeter Brune 5136cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 5146cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 5156cc8130cSPeter Brune 5164b11644fSPeter Brune 5174b11644fSPeter Brune Level: beginner 5184b11644fSPeter Brune 5194b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 5206cc8130cSPeter Brune 5214b11644fSPeter Brune M*/ 5224b11644fSPeter Brune EXTERN_C_BEGIN 5234b11644fSPeter Brune #undef __FUNCT__ 5244b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 5254b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 5264b11644fSPeter Brune { 5274b11644fSPeter Brune 5284b11644fSPeter Brune PetscErrorCode ierr; 5299f83bee8SJed Brown SNES_QN *qn; 5304b11644fSPeter Brune 5314b11644fSPeter Brune PetscFunctionBegin; 5324b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 5334b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 5344b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 5354b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 5364b11644fSPeter Brune snes->ops->view = 0; 5374b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 5384b11644fSPeter Brune 53942f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 54042f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 54142f4f86dSBarry Smith 5420e444f03SPeter Brune snes->max_funcs = 30000; 5430e444f03SPeter Brune snes->max_its = 10000; 5440e444f03SPeter Brune 5459f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 5464b11644fSPeter Brune snes->data = (void *) qn; 5470ecc9e1dSPeter Brune qn->m = 10; 548b21d5a53SPeter Brune qn->scaling = 1.0; 5494b11644fSPeter Brune qn->dX = PETSC_NULL; 5506bf1b2e5SPeter Brune qn->dF = PETSC_NULL; 551*8e231d97SPeter Brune qn->dXtdF = PETSC_NULL; 552*8e231d97SPeter Brune qn->dFtdX = PETSC_NULL; 553*8e231d97SPeter Brune qn->dXdFmat = PETSC_NULL; 55444f7e39eSPeter Brune qn->monitor = PETSC_NULL; 555*8e231d97SPeter Brune qn->aggreduct = PETSC_FALSE; 55688d374b2SPeter Brune qn->powell_gamma = 0.9; 55788d374b2SPeter Brune qn->powell_downhill = 0.2; 55888d374b2SPeter Brune qn->compositiontype = SNES_QN_SEQUENTIAL; 5590ecc9e1dSPeter Brune qn->scalingtype = SNES_QN_SHANNOSCALE; 5600ecc9e1dSPeter Brune qn->n_restart = -1; 561ea630c6eSPeter Brune 5624b11644fSPeter Brune PetscFunctionReturn(0); 5634b11644fSPeter Brune } 564*8e231d97SPeter Brune 5654b11644fSPeter Brune EXTERN_C_END 566