1b45d2f2cSJed Brown #include <petsc-private/snesimpl.h> 24b11644fSPeter Brune 38e231d97SPeter Brune #define H(i,j) qn->dXdFmat[i*qn->m + j] 48e231d97SPeter Brune 588d374b2SPeter Brune typedef enum {SNES_QN_SEQUENTIAL, SNES_QN_COMPOSED} SNESQNCompositionType; 60ecc9e1dSPeter Brune typedef enum {SNES_QN_SHANNOSCALE, SNES_QN_LSSCALE, SNES_QN_JACOBIANSCALE} SNESQNScalingType; 76bf1b2e5SPeter Brune 84b11644fSPeter Brune typedef struct { 94b11644fSPeter Brune Vec *dX; /* The change in X */ 106bf1b2e5SPeter Brune Vec *dF; /* The change in F */ 118e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 128e231d97SPeter Brune PetscScalar *alpha, *beta; 138e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 14*18aa0c0cSPeter Brune PetscBool singlereduction; /* Aggregated reduction implementation */ 158e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 1644f7e39eSPeter Brune PetscViewer monitor; 176bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 186bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 19b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 200ecc9e1dSPeter Brune PetscInt n_restart; /* the maximum iterations between restart */ 2188d374b2SPeter Brune 2288d374b2SPeter Brune SNESQNCompositionType compositiontype; /* determine if the composition is done sequentially or as a composition */ 230ecc9e1dSPeter Brune SNESQNScalingType scalingtype; /* determine if the composition is done sequentially or as a composition */ 2488d374b2SPeter Brune 259f83bee8SJed Brown } SNES_QN; 264b11644fSPeter Brune 274b11644fSPeter Brune #undef __FUNCT__ 288c40d5fbSBarry Smith #define __FUNCT__ "SNESQNApplyJinv_Private" 298c40d5fbSBarry Smith PetscErrorCode SNESQNApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) { 304b11644fSPeter Brune 314b11644fSPeter Brune PetscErrorCode ierr; 324b11644fSPeter Brune 339f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 344b11644fSPeter Brune 35335fb975SPeter Brune Vec Yin = snes->work[3]; 360ecc9e1dSPeter Brune 374b11644fSPeter Brune Vec *dX = qn->dX; 386bf1b2e5SPeter Brune Vec *dF = qn->dF; 394b11644fSPeter Brune 40bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 41bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 428e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 438e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 44bd052dfeSPeter Brune 450ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 460ecc9e1dSPeter Brune KSPConvergedReason kspreason; 470ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 480ecc9e1dSPeter Brune 498e231d97SPeter Brune PetscInt k, i, j, g, lits; 504b11644fSPeter Brune PetscInt m = qn->m; 514b11644fSPeter Brune PetscScalar t; 524b11644fSPeter Brune PetscInt l = m; 534b11644fSPeter Brune 548e231d97SPeter Brune Mat jac, jac_pre; 558e231d97SPeter Brune 564b11644fSPeter Brune PetscFunctionBegin; 574b11644fSPeter Brune 583af51624SPeter Brune ierr = VecCopy(D, Y);CHKERRQ(ierr); 594b11644fSPeter Brune 605ba6227bSPeter Brune if (it < m) l = it; 614b11644fSPeter Brune 62*18aa0c0cSPeter Brune if (qn->singlereduction) { 638e231d97SPeter Brune ierr = VecMDot(Y, l, qn->dX, YtdX);CHKERRQ(ierr); 648e231d97SPeter Brune } 654b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 664b11644fSPeter Brune for (i = 0; i < l; i++) { 67b21d5a53SPeter Brune k = (it - i - 1) % l; 68*18aa0c0cSPeter Brune if (qn->singlereduction) { 698e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 708e231d97SPeter Brune t = YtdX[k]; 718e231d97SPeter Brune for (j = 0; j < i; j++) { 728e231d97SPeter Brune g = (it - j - 1) % l; 738e231d97SPeter Brune t += -alpha[g]*H(g, k); 748e231d97SPeter Brune } 758e231d97SPeter Brune alpha[k] = t / H(k, k); 768e231d97SPeter Brune } else { 773af51624SPeter Brune ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr); 788e231d97SPeter Brune alpha[k] = t / dXtdF[k]; 798e231d97SPeter Brune } 8044f7e39eSPeter Brune if (qn->monitor) { 813af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 825ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 833af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 8444f7e39eSPeter Brune } 856bf1b2e5SPeter Brune ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr); 864b11644fSPeter Brune } 874b11644fSPeter Brune 880ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_JACOBIANSCALE) { 898e231d97SPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 908e231d97SPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 910ecc9e1dSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,Yin);CHKERRQ(ierr); 920ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 930ecc9e1dSPeter Brune if (kspreason < 0) { 940ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 950ecc9e1dSPeter 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); 960ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 970ecc9e1dSPeter Brune PetscFunctionReturn(0); 980ecc9e1dSPeter Brune } 990ecc9e1dSPeter Brune } 1000ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 1010ecc9e1dSPeter Brune snes->linear_its += lits; 1020ecc9e1dSPeter Brune ierr = VecCopy(Yin, Y);CHKERRQ(ierr); 1030ecc9e1dSPeter Brune } else { 104b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 1050ecc9e1dSPeter Brune } 106*18aa0c0cSPeter Brune if (qn->singlereduction) { 1078e231d97SPeter Brune ierr = VecMDot(Y, l, qn->dF, YtdX);CHKERRQ(ierr); 1088e231d97SPeter Brune } 1094b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 110bd052dfeSPeter Brune for (i = 0; i < l; i++) { 111b21d5a53SPeter Brune k = (it + i - l) % l; 112*18aa0c0cSPeter Brune if (qn->singlereduction) { 1138e231d97SPeter Brune t = YtdX[k]; 1148e231d97SPeter Brune for (j = 0; j < i; j++) { 1158e231d97SPeter Brune g = (it + j - l) % l; 1168e231d97SPeter Brune t += (alpha[g] - beta[g])*H(k, g); 1178e231d97SPeter Brune } 1188e231d97SPeter Brune beta[k] = t / H(k, k); 1198e231d97SPeter Brune } else { 1206bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 1218e231d97SPeter Brune beta[k] = t / dXtdF[k]; 1228e231d97SPeter Brune } 1233af51624SPeter Brune ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]); 12444f7e39eSPeter Brune if (qn->monitor) { 1253af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1265ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 1273af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 12844f7e39eSPeter Brune } 1294b11644fSPeter Brune } 1304b11644fSPeter Brune PetscFunctionReturn(0); 1314b11644fSPeter Brune } 1324b11644fSPeter Brune 1334b11644fSPeter Brune #undef __FUNCT__ 1344b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 1354b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 1364b11644fSPeter Brune { 1374b11644fSPeter Brune 1384b11644fSPeter Brune PetscErrorCode ierr; 1399f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 1404b11644fSPeter Brune 14115f5eeeaSPeter Brune Vec X, Xold; 142335fb975SPeter Brune Vec F, B; 143335fb975SPeter Brune Vec Y, FPC, D, Dold; 1443af51624SPeter Brune SNESConvergedReason reason; 1458e231d97SPeter Brune PetscInt i, i_r, k, l, j; 1464b11644fSPeter Brune 147335fb975SPeter Brune PetscReal fnorm, xnorm, ynorm, gnorm; 1484b11644fSPeter Brune PetscInt m = qn->m; 149335fb975SPeter Brune PetscBool lssucceed; 1504b11644fSPeter Brune 151bd052dfeSPeter Brune Vec *dX = qn->dX; 1526bf1b2e5SPeter Brune Vec *dF = qn->dF; 1538e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 1548e231d97SPeter Brune PetscScalar *dFtdX = qn->dFtdX; 15588d374b2SPeter Brune PetscScalar DolddotD, DolddotDold, DdotD, YdotD, a; 1564b11644fSPeter Brune 1570ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 1580ecc9e1dSPeter Brune 1594b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 1604b11644fSPeter Brune PetscFunctionBegin; 1614b11644fSPeter Brune 1629f3a0142SPeter Brune X = snes->vec_sol; /* solution vector */ 1639f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 1643af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 1653af51624SPeter Brune B = snes->vec_rhs; 166335fb975SPeter Brune Xold = snes->work[0]; 1673af51624SPeter Brune 1683af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 169335fb975SPeter Brune D = snes->work[1]; 170335fb975SPeter Brune Dold = snes->work[2]; 1714b11644fSPeter Brune 1724b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1734b11644fSPeter Brune 1744b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1754b11644fSPeter Brune snes->iter = 0; 1764b11644fSPeter Brune snes->norm = 0.; 1774b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 178e4ed7901SPeter Brune if (!snes->vec_func_init_set){ 17915f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1804b11644fSPeter Brune if (snes->domainerror) { 1814b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1824b11644fSPeter Brune PetscFunctionReturn(0); 1834b11644fSPeter Brune } 184e4ed7901SPeter Brune } else { 185e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 186e4ed7901SPeter Brune } 187e4ed7901SPeter Brune 188e4ed7901SPeter Brune if (!snes->norm_init_set) { 18915f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1904b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 191e4ed7901SPeter Brune } else { 192e4ed7901SPeter Brune fnorm = snes->norm_init; 193e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 194e4ed7901SPeter Brune } 195e4ed7901SPeter Brune 1964b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1974b11644fSPeter Brune snes->norm = fnorm; 1984b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1994b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 2004b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 2014b11644fSPeter Brune 2024b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 2034b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 2044b11644fSPeter Brune /* test convergence */ 2054b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2064b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 20770d3b23bSPeter Brune 20888d374b2SPeter Brune /* composed solve -- either sequential or composed */ 20988d374b2SPeter Brune if (snes->pc) { 21088d374b2SPeter Brune if (qn->compositiontype == SNES_QN_SEQUENTIAL) { 211217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 212217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 21388d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 21488d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 21588d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 21688d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 21788d374b2SPeter Brune PetscFunctionReturn(0); 21888d374b2SPeter Brune } 21988d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 22088d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 22188d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 2226bf1b2e5SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 22388d374b2SPeter Brune } else { 22488d374b2SPeter Brune ierr = VecCopy(X, Y);CHKERRQ(ierr); 225217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 226217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 22788d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, Y);CHKERRQ(ierr); 22888d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 22988d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 23088d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 23188d374b2SPeter Brune PetscFunctionReturn(0); 23288d374b2SPeter Brune } 23388d374b2SPeter Brune ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 23488d374b2SPeter Brune } 23588d374b2SPeter Brune } else { 23688d374b2SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 23788d374b2SPeter Brune } 23888d374b2SPeter Brune ierr = VecCopy(Y, D);CHKERRQ(ierr); 2393af51624SPeter Brune 240f8e15203SPeter Brune /* scale the initial update */ 2410ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_JACOBIANSCALE) { 2420ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 2430ecc9e1dSPeter Brune } 2440ecc9e1dSPeter Brune 2455ba6227bSPeter Brune for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 2468c40d5fbSBarry Smith ierr = SNESQNApplyJinv_Private(snes, i_r, D, Y);CHKERRQ(ierr); 24770d3b23bSPeter Brune /* line search for lambda */ 24870d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 24988d374b2SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 25015f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 251f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 2529f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 2539f3a0142SPeter Brune if (snes->domainerror) { 2549f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2559f3a0142SPeter Brune PetscFunctionReturn(0); 2569f3a0142SPeter Brune } 257f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 2589f3a0142SPeter Brune if (!lssucceed) { 2599f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 2609f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 2619f3a0142SPeter Brune break; 2629f3a0142SPeter Brune } 2639f3a0142SPeter Brune } 264f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 2650ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_LSSCALE) { 266f1c6b773SPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 2670ecc9e1dSPeter Brune } 2683af51624SPeter Brune 26988d374b2SPeter Brune /* convergence monitoring */ 270b21d5a53SPeter 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); 271b21d5a53SPeter Brune 272360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 273360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 274360c497dSPeter Brune 2758409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 2768409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 2774b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 2785ba6227bSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2794b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 2804b11644fSPeter Brune 28188d374b2SPeter Brune 28288d374b2SPeter Brune if (snes->pc) { 28388d374b2SPeter Brune if (qn->compositiontype == SNES_QN_SEQUENTIAL) { 284217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 285217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 28688d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 28788d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 28888d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 28988d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 29088d374b2SPeter Brune PetscFunctionReturn(0); 29188d374b2SPeter Brune } 29288d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 29388d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 29488d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 29588d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 29688d374b2SPeter Brune } else { 29788d374b2SPeter Brune ierr = VecCopy(X, D);CHKERRQ(ierr); 298217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 299217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 30088d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr); 30188d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 30288d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 30388d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 30488d374b2SPeter Brune PetscFunctionReturn(0); 30588d374b2SPeter Brune } 30688d374b2SPeter Brune ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr); 30788d374b2SPeter Brune } 30888d374b2SPeter Brune } else { 30988d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 31088d374b2SPeter Brune } 31188d374b2SPeter Brune 3126bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 3138e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 3148e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 3158e231d97SPeter Brune ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr); 3168e231d97SPeter Brune ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr); 3178e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 3188e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 3198e231d97SPeter Brune ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr); 3208e231d97SPeter Brune ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr); 3210ecc9e1dSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold)) || (i_r > qn->n_restart - 1 && qn->n_restart > 0)) { 3225ba6227bSPeter Brune if (qn->monitor) { 3235ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 324516fe3c3SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);CHKERRQ(ierr); 3255ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3265ba6227bSPeter Brune } 3275ba6227bSPeter Brune i_r = -1; 3285ba6227bSPeter Brune /* general purpose update */ 3295ba6227bSPeter Brune if (snes->ops->update) { 3305ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 3315ba6227bSPeter Brune } 3320ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_JACOBIANSCALE) { 3330ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 3340ecc9e1dSPeter Brune } 3355ba6227bSPeter Brune } else { 336bd052dfeSPeter Brune /* set the differences */ 3375ba6227bSPeter Brune k = i_r % m; 3388e231d97SPeter Brune l = m; 3398e231d97SPeter Brune if (i_r + 1 < m) l = i_r + 1; 34088d374b2SPeter Brune ierr = VecCopy(D, dF[k]);CHKERRQ(ierr); 34188d374b2SPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 34215f5eeeaSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 34315f5eeeaSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 344*18aa0c0cSPeter Brune if (qn->singlereduction) { 3458e231d97SPeter Brune ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr); 3468e231d97SPeter Brune ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr); 3478e231d97SPeter Brune for (j = 0; j < l; j++) { 3488e231d97SPeter Brune H(k, j) = dFtdX[j]; 3498e231d97SPeter Brune H(j, k) = dXtdF[j]; 3508e231d97SPeter Brune } 3518e231d97SPeter Brune /* copy back over to make the computation of alpha and beta easier */ 3528e231d97SPeter Brune for (j = 0; j < l; j++) { 3538e231d97SPeter Brune dXtdF[j] = H(j, j); 3548e231d97SPeter Brune } 3558e231d97SPeter Brune } else { 3568e231d97SPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 3578e231d97SPeter Brune } 3586bf1b2e5SPeter Brune /* set scaling to be shanno scaling */ 3590ecc9e1dSPeter Brune if (qn->scalingtype == SNES_QN_SHANNOSCALE) { 3606bf1b2e5SPeter Brune ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr); 3618e231d97SPeter Brune qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(a); 3620ecc9e1dSPeter Brune } 36370d3b23bSPeter Brune /* general purpose update */ 36470d3b23bSPeter Brune if (snes->ops->update) { 36570d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 36670d3b23bSPeter Brune } 3675ba6227bSPeter Brune } 3684b11644fSPeter Brune } 3694b11644fSPeter Brune if (i == snes->max_its) { 3704b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 3714b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 3724b11644fSPeter Brune } 3734b11644fSPeter Brune PetscFunctionReturn(0); 3744b11644fSPeter Brune } 3754b11644fSPeter Brune 3764b11644fSPeter Brune 3774b11644fSPeter Brune #undef __FUNCT__ 3784b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 3794b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 3804b11644fSPeter Brune { 3819f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 3824b11644fSPeter Brune PetscErrorCode ierr; 383335fb975SPeter Brune 3844b11644fSPeter Brune PetscFunctionBegin; 3854b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 3866bf1b2e5SPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 3878e231d97SPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 3888e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 3898e231d97SPeter Brune qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 3908e231d97SPeter Brune 391*18aa0c0cSPeter Brune if (qn->singlereduction) { 3928e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 3938e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 3948e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 3958e231d97SPeter Brune } 396335fb975SPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 397335fb975SPeter Brune 398335fb975SPeter Brune /* set up the line search */ 3998e231d97SPeter Brune if (qn->scalingtype == SNES_QN_JACOBIANSCALE) { 4008e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 4018e231d97SPeter Brune } 4024b11644fSPeter Brune PetscFunctionReturn(0); 4034b11644fSPeter Brune } 4044b11644fSPeter Brune 4054b11644fSPeter Brune #undef __FUNCT__ 4064b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 4074b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 4084b11644fSPeter Brune { 4094b11644fSPeter Brune PetscErrorCode ierr; 4109f83bee8SJed Brown SNES_QN *qn; 4114b11644fSPeter Brune PetscFunctionBegin; 4124b11644fSPeter Brune if (snes->data) { 4139f83bee8SJed Brown qn = (SNES_QN*)snes->data; 4144b11644fSPeter Brune if (qn->dX) { 4154b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 4164b11644fSPeter Brune } 4176bf1b2e5SPeter Brune if (qn->dF) { 4186bf1b2e5SPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 4194b11644fSPeter Brune } 420*18aa0c0cSPeter Brune if (qn->singlereduction) { 4218e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 4228e231d97SPeter Brune } 4238e231d97SPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 4244b11644fSPeter Brune } 4254b11644fSPeter Brune PetscFunctionReturn(0); 4264b11644fSPeter Brune } 4274b11644fSPeter Brune 4284b11644fSPeter Brune #undef __FUNCT__ 4294b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 4304b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 4314b11644fSPeter Brune { 4324b11644fSPeter Brune PetscErrorCode ierr; 4334b11644fSPeter Brune PetscFunctionBegin; 4344b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 4354b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 4369c05af5eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 4374b11644fSPeter Brune PetscFunctionReturn(0); 4384b11644fSPeter Brune } 4394b11644fSPeter Brune 4404b11644fSPeter Brune #undef __FUNCT__ 4414b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 4424b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 4434b11644fSPeter Brune { 4444b11644fSPeter Brune 4454b11644fSPeter Brune PetscErrorCode ierr; 4469f83bee8SJed Brown SNES_QN *qn; 4475b4627e8SPeter Brune const char *compositions[] = {"sequential", "composed"}; 4480ecc9e1dSPeter Brune const char *scalings[] = {"shanno", "ls", "jacobian"}; 44988d374b2SPeter Brune PetscInt indx = 0; 45088d374b2SPeter Brune PetscBool flg; 45144f7e39eSPeter Brune PetscBool monflg = PETSC_FALSE; 452f1c6b773SPeter Brune SNESLineSearch linesearch; 4534b11644fSPeter Brune PetscFunctionBegin; 4544b11644fSPeter Brune 4559f83bee8SJed Brown qn = (SNES_QN*)snes->data; 4564b11644fSPeter Brune 4574b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 4585b4627e8SPeter Brune ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 4590ecc9e1dSPeter Brune ierr = PetscOptionsInt("-snes_qn_restart", "Maximum number of iterations between restarts", "SNESQN", qn->n_restart, &qn->n_restart, PETSC_NULL);CHKERRQ(ierr); 4605b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr); 4615b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_downhill", "Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr); 4625b4627e8SPeter Brune ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 463*18aa0c0cSPeter Brune ierr = PetscOptionsBool("-snes_qn_singlereduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, PETSC_NULL);CHKERRQ(ierr); 4645b4627e8SPeter Brune ierr = PetscOptionsEList("-snes_qn_composition", "Composition type", "SNESQN",compositions,2,"sequential",&indx,&flg);CHKERRQ(ierr); 46588d374b2SPeter Brune if (flg) { 46688d374b2SPeter Brune switch (indx) { 4675b4627e8SPeter Brune case 0: qn->compositiontype = SNES_QN_SEQUENTIAL; 4685b4627e8SPeter Brune break; 4695b4627e8SPeter Brune case 1: qn->compositiontype = SNES_QN_COMPOSED; 4705b4627e8SPeter Brune break; 47188d374b2SPeter Brune } 47288d374b2SPeter Brune } 4730ecc9e1dSPeter Brune ierr = PetscOptionsEList("-snes_qn_scaling", "Scaling type", "SNESQN",scalings,3,"shanno",&indx,&flg);CHKERRQ(ierr); 4740ecc9e1dSPeter Brune if (flg) { 4750ecc9e1dSPeter Brune switch (indx) { 4760ecc9e1dSPeter Brune case 0: qn->scalingtype = SNES_QN_SHANNOSCALE; 4770ecc9e1dSPeter Brune break; 4780ecc9e1dSPeter Brune case 1: qn->scalingtype = SNES_QN_LSSCALE; 4790ecc9e1dSPeter Brune break; 4800ecc9e1dSPeter Brune case 2: qn->scalingtype = SNES_QN_JACOBIANSCALE; 4810ecc9e1dSPeter Brune snes->usesksp = PETSC_TRUE; 4820ecc9e1dSPeter Brune break; 4830ecc9e1dSPeter Brune } 4840ecc9e1dSPeter Brune } 4850ecc9e1dSPeter Brune 4864b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 4879e764e56SPeter Brune if (!snes->linesearch) { 488f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 4899e764e56SPeter Brune if (!snes->pc || qn->compositiontype == SNES_QN_SEQUENTIAL) { 4901a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 4919e764e56SPeter Brune } else { 4921a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 4939e764e56SPeter Brune } 4949e764e56SPeter Brune } 49544f7e39eSPeter Brune if (monflg) { 49644f7e39eSPeter Brune qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 49744f7e39eSPeter Brune } 4984b11644fSPeter Brune PetscFunctionReturn(0); 4994b11644fSPeter Brune } 5004b11644fSPeter Brune 50192c02d66SPeter Brune 5024b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 5034b11644fSPeter Brune /*MC 5044b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 5054b11644fSPeter Brune 5066cc8130cSPeter Brune Options Database: 5076cc8130cSPeter Brune 5086cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 5091867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 5101867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 5111867fe5bSPeter Brune . -snes_qn_composition <sequential, composed>- Type of composition. 512f3aaf736SPeter Brune . -snes_ls <basic, basicnonorms, quadratic, critical> - Type of line search. 51344f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 5146cc8130cSPeter Brune 51544f7e39eSPeter Brune Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to 5166cc8130cSPeter Brune form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 5176cc8130cSPeter Brune generalized to implement several limited-memory Broyden methods. 5186cc8130cSPeter Brune 5191867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 5201867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 5211867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 5221867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 5231867fe5bSPeter Brune 5246cc8130cSPeter Brune References: 5256cc8130cSPeter Brune 5266cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 5276cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 5286cc8130cSPeter Brune 5294b11644fSPeter Brune 5304b11644fSPeter Brune Level: beginner 5314b11644fSPeter Brune 5324b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 5336cc8130cSPeter Brune 5344b11644fSPeter Brune M*/ 5354b11644fSPeter Brune EXTERN_C_BEGIN 5364b11644fSPeter Brune #undef __FUNCT__ 5374b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 5384b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 5394b11644fSPeter Brune { 5404b11644fSPeter Brune 5414b11644fSPeter Brune PetscErrorCode ierr; 5429f83bee8SJed Brown SNES_QN *qn; 5434b11644fSPeter Brune 5444b11644fSPeter Brune PetscFunctionBegin; 5454b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 5464b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 5474b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 5484b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 5494b11644fSPeter Brune snes->ops->view = 0; 5504b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 5514b11644fSPeter Brune 55242f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 55342f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 55442f4f86dSBarry Smith 55588976e71SPeter Brune if (!snes->tolerancesset) { 5560e444f03SPeter Brune snes->max_funcs = 30000; 5570e444f03SPeter Brune snes->max_its = 10000; 55888976e71SPeter Brune } 5590e444f03SPeter Brune 5609f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 5614b11644fSPeter Brune snes->data = (void *) qn; 5620ecc9e1dSPeter Brune qn->m = 10; 563b21d5a53SPeter Brune qn->scaling = 1.0; 5644b11644fSPeter Brune qn->dX = PETSC_NULL; 5656bf1b2e5SPeter Brune qn->dF = PETSC_NULL; 5668e231d97SPeter Brune qn->dXtdF = PETSC_NULL; 5678e231d97SPeter Brune qn->dFtdX = PETSC_NULL; 5688e231d97SPeter Brune qn->dXdFmat = PETSC_NULL; 56944f7e39eSPeter Brune qn->monitor = PETSC_NULL; 570*18aa0c0cSPeter Brune qn->singlereduction = PETSC_FALSE; 57188d374b2SPeter Brune qn->powell_gamma = 0.9; 57288d374b2SPeter Brune qn->powell_downhill = 0.2; 57388d374b2SPeter Brune qn->compositiontype = SNES_QN_SEQUENTIAL; 5740ecc9e1dSPeter Brune qn->scalingtype = SNES_QN_SHANNOSCALE; 5750ecc9e1dSPeter Brune qn->n_restart = -1; 576ea630c6eSPeter Brune 5774b11644fSPeter Brune PetscFunctionReturn(0); 5784b11644fSPeter Brune } 5798e231d97SPeter Brune 5804b11644fSPeter Brune EXTERN_C_END 581