10c777b0cSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 24b11644fSPeter Brune 38e231d97SPeter Brune #define H(i,j) qn->dXdFmat[i*qn->m + j] 48e231d97SPeter Brune 56a6fc655SJed Brown const char *const SNESQNScaleTypes[] = {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0}; 66a6fc655SJed Brown const char *const SNESQNRestartTypes[] = {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0}; 7b8840d6bSPeter Brune const char *const SNESQNTypes[] = {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_RESTART_",0}; 8b8840d6bSPeter Brune 9b8840d6bSPeter Brune typedef enum {SNES_QN_LBFGS = 0, 10b8840d6bSPeter Brune SNES_QN_BROYDEN = 1, 11b8840d6bSPeter Brune SNES_QN_BADBROYDEN = 2} SNESQNType; 126bf1b2e5SPeter Brune 134b11644fSPeter Brune typedef struct { 14b8840d6bSPeter Brune Vec *U; /* Stored past states (vary from method to method) */ 15b8840d6bSPeter Brune Vec *V; /* Stored past states (vary from method to method) */ 168e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 178e231d97SPeter Brune PetscScalar *alpha, *beta; 188e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 1918aa0c0cSPeter Brune PetscBool singlereduction; /* Aggregated reduction implementation */ 208e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 2144f7e39eSPeter Brune PetscViewer monitor; 226bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 236bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 24b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 2588d374b2SPeter Brune 26b8840d6bSPeter Brune SNESQNType type; /* the type of quasi-newton method used */ 2788f769c5SPeter Brune SNESQNScaleType scale_type; /* the type of scaling used */ 280c777b0cSPeter Brune SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 299f83bee8SJed Brown } SNES_QN; 304b11644fSPeter Brune 314b11644fSPeter Brune #undef __FUNCT__ 32b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden" 332150357eSBarry Smith PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold, Vec D, Vec Dold) 342150357eSBarry Smith { 354b11644fSPeter Brune PetscErrorCode ierr; 369f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 37b8840d6bSPeter Brune Vec W = snes->work[3]; 38b8840d6bSPeter Brune Vec *U = qn->U; 39b8840d6bSPeter Brune Vec *V = qn->V; 40b8840d6bSPeter Brune KSPConvergedReason kspreason; 41b8840d6bSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 42b8840d6bSPeter Brune PetscInt k,i,lits; 43b8840d6bSPeter Brune PetscInt m = qn->m; 44b8840d6bSPeter Brune PetscScalar gdot; 45b8840d6bSPeter Brune PetscInt l = m; 46b8840d6bSPeter Brune Mat jac,jac_pre; 472150357eSBarry Smith 48b8840d6bSPeter Brune PetscFunctionBegin; 49b8840d6bSPeter Brune if (it < m) l = it; 50b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 51b8840d6bSPeter Brune ierr = SNESGetJacobian(snes,&jac,&jac_pre,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 52b8840d6bSPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 53b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,D,W);CHKERRQ(ierr); 54b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 55b8840d6bSPeter Brune if (kspreason < 0) { 56b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 57b8840d6bSPeter 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); 58b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 59b8840d6bSPeter Brune PetscFunctionReturn(0); 60b8840d6bSPeter Brune } 61b8840d6bSPeter Brune } 62b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 63b8840d6bSPeter Brune snes->linear_its += lits; 64b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 65b8840d6bSPeter Brune } else { 66b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 67b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 68b8840d6bSPeter Brune } 69b8840d6bSPeter Brune 70b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 71b8840d6bSPeter Brune if (it > 1) { 72b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 73b8840d6bSPeter Brune k = (it+i-l)%l; 74b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 75b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 76b8840d6bSPeter Brune if (qn->monitor) { 77b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 78b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 79b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 80b8840d6bSPeter Brune } 81b8840d6bSPeter Brune } 82b8840d6bSPeter Brune } 83b8840d6bSPeter Brune if (it < m) l = it; 84b8840d6bSPeter Brune 85b8840d6bSPeter Brune /* set W to be the last step, post-linesearch */ 86b8840d6bSPeter Brune ierr = VecCopy(Xold,W);CHKERRQ(ierr); 87b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 88b8840d6bSPeter Brune 89b8840d6bSPeter Brune if (l > 0) { 90b8840d6bSPeter Brune k = (it-1)%l; 91b8840d6bSPeter Brune ierr = VecCopy(W,U[k]);CHKERRQ(ierr); 92b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,Y);CHKERRQ(ierr); 93b8840d6bSPeter Brune ierr = VecDot(U[k],W,&gdot);CHKERRQ(ierr); 94b8840d6bSPeter Brune ierr = VecCopy(Y,V[k]);CHKERRQ(ierr); 95b8840d6bSPeter Brune ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr); 96b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 97b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 98b8840d6bSPeter Brune if (qn->monitor) { 99b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 100b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 101b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 102b8840d6bSPeter Brune } 103b8840d6bSPeter Brune } 104b8840d6bSPeter Brune PetscFunctionReturn(0); 105b8840d6bSPeter Brune } 106b8840d6bSPeter Brune 107b8840d6bSPeter Brune #undef __FUNCT__ 108b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden" 1090adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1100adebc6cSBarry Smith { 111b8840d6bSPeter Brune PetscErrorCode ierr; 112b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 113b8840d6bSPeter Brune Vec W = snes->work[3]; 114b8840d6bSPeter Brune Vec *U = qn->U; 115b8840d6bSPeter Brune Vec *T = qn->V; 116b8840d6bSPeter Brune 117b8840d6bSPeter Brune /* ksp thing for jacobian scaling */ 118b8840d6bSPeter Brune KSPConvergedReason kspreason; 119b8840d6bSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 120b8840d6bSPeter Brune PetscInt k, i, lits; 121b8840d6bSPeter Brune PetscInt m = qn->m; 122b8840d6bSPeter Brune PetscScalar gdot; 123b8840d6bSPeter Brune PetscInt l = m; 124b8840d6bSPeter Brune Mat jac, jac_pre; 1250adebc6cSBarry Smith 126b8840d6bSPeter Brune PetscFunctionBegin; 127b8840d6bSPeter Brune if (it < m) l = it; 128b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 129b8840d6bSPeter Brune if (l > 0) { 130b8840d6bSPeter Brune k = (it-1)%l; 131b8840d6bSPeter Brune ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 132b8840d6bSPeter Brune ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 133b8840d6bSPeter Brune ierr = VecDot(U[k],U[k],&gdot);CHKERRQ(ierr); 134b8840d6bSPeter Brune ierr = VecCopy(D,T[k]);CHKERRQ(ierr); 135b8840d6bSPeter Brune ierr = VecScale(T[k],1.0/gdot);CHKERRQ(ierr); 136b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 137b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 138b8840d6bSPeter Brune if (qn->monitor) { 139b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 140b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 141b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 142b8840d6bSPeter Brune } 143b8840d6bSPeter Brune } 144b8840d6bSPeter Brune 145b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 146b8840d6bSPeter Brune if (it > 1) { 147b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 148b8840d6bSPeter Brune k = (it+i-l)%l; 149b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 150b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 151b8840d6bSPeter Brune if (qn->monitor) { 152b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 153b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 154b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 155b8840d6bSPeter Brune } 156b8840d6bSPeter Brune } 157b8840d6bSPeter Brune } 158b8840d6bSPeter Brune 159b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 160b8840d6bSPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 161b8840d6bSPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 162b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr); 163b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 164b8840d6bSPeter Brune if (kspreason < 0) { 165b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 166b8840d6bSPeter 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); 167b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 168b8840d6bSPeter Brune PetscFunctionReturn(0); 169b8840d6bSPeter Brune } 170b8840d6bSPeter Brune } 171b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 172b8840d6bSPeter Brune snes->linear_its += lits; 173b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 174b8840d6bSPeter Brune } else { 175b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 176b8840d6bSPeter Brune } 177b8840d6bSPeter Brune PetscFunctionReturn(0); 178b8840d6bSPeter Brune } 179b8840d6bSPeter Brune 180b8840d6bSPeter Brune #undef __FUNCT__ 181b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 1820adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1830adebc6cSBarry Smith { 184b8840d6bSPeter Brune PetscErrorCode ierr; 185b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 186b8840d6bSPeter Brune Vec W = snes->work[3]; 187b8840d6bSPeter Brune Vec *dX = qn->U; 188b8840d6bSPeter Brune Vec *dF = qn->V; 189bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 190bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 1918e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 192b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 1938e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 194bd052dfeSPeter Brune 1950ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 1960ecc9e1dSPeter Brune KSPConvergedReason kspreason; 1970ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 1988e231d97SPeter Brune PetscInt k,i,j,g,lits; 1994b11644fSPeter Brune PetscInt m = qn->m; 2004b11644fSPeter Brune PetscScalar t; 2014b11644fSPeter Brune PetscInt l = m; 2028e231d97SPeter Brune Mat jac,jac_pre; 2038e231d97SPeter Brune 2044b11644fSPeter Brune PetscFunctionBegin; 2055ba6227bSPeter Brune if (it < m) l = it; 206b8840d6bSPeter Brune if (it > 0) { 207b8840d6bSPeter Brune k = (it - 1) % l; 208b8840d6bSPeter Brune ierr = VecCopy(D, dF[k]);CHKERRQ(ierr); 209b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 210b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 211b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 21218aa0c0cSPeter Brune if (qn->singlereduction) { 213b8840d6bSPeter Brune ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr); 214b8840d6bSPeter Brune ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr); 215b8840d6bSPeter Brune for (j = 0; j < l; j++) { 216b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 217b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 218b8840d6bSPeter Brune } 219b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 220b8840d6bSPeter Brune for (j = 0; j < l; j++) { 221b8840d6bSPeter Brune dXtdF[j] = H(j, j); 222b8840d6bSPeter Brune } 223b8840d6bSPeter Brune } else { 224b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 225b8840d6bSPeter Brune } 226b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 22767392de3SBarry Smith PetscReal dFtdF; 22867392de3SBarry Smith ierr = VecDotRealPart(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 22967392de3SBarry Smith qn->scaling = PetscRealPart(dXtdF[k])/dFtdF; 230b8840d6bSPeter Brune } else if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 231b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 232b8840d6bSPeter Brune } 233b8840d6bSPeter Brune } 234b8840d6bSPeter Brune 235b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 236b8840d6bSPeter Brune 237b8840d6bSPeter Brune if (qn->singlereduction) { 238b8840d6bSPeter Brune ierr = VecMDot(Y,l,dX,YtdX);CHKERRQ(ierr); 2398e231d97SPeter Brune } 2404b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2414b11644fSPeter Brune for (i=0;i<l;i++) { 242b21d5a53SPeter Brune k = (it-i-1)%l; 24318aa0c0cSPeter Brune if (qn->singlereduction) { 2448e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2458e231d97SPeter Brune t = YtdX[k]; 2468e231d97SPeter Brune for (j=0;j<i;j++) { 2478e231d97SPeter Brune g = (it-j-1)%l; 2488e231d97SPeter Brune t += -alpha[g]*H(g, k); 2498e231d97SPeter Brune } 2508e231d97SPeter Brune alpha[k] = t/H(k,k); 2518e231d97SPeter Brune } else { 2523af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2538e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2548e231d97SPeter Brune } 25544f7e39eSPeter Brune if (qn->monitor) { 2563af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2575ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2583af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 25944f7e39eSPeter Brune } 2606bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2614b11644fSPeter Brune } 2624b11644fSPeter Brune 2630c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 2648e231d97SPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 2658e231d97SPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 266b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr); 2670ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2680ecc9e1dSPeter Brune if (kspreason < 0) { 2690ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 2700ecc9e1dSPeter 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); 2710ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 2720ecc9e1dSPeter Brune PetscFunctionReturn(0); 2730ecc9e1dSPeter Brune } 2740ecc9e1dSPeter Brune } 2750ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2760ecc9e1dSPeter Brune snes->linear_its += lits; 277b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2780ecc9e1dSPeter Brune } else { 279b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2800ecc9e1dSPeter Brune } 28118aa0c0cSPeter Brune if (qn->singlereduction) { 282b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 2838e231d97SPeter Brune } 2844b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 285bd052dfeSPeter Brune for (i = 0; i < l; i++) { 286b21d5a53SPeter Brune k = (it + i - l) % l; 28718aa0c0cSPeter Brune if (qn->singlereduction) { 2888e231d97SPeter Brune t = YtdX[k]; 2898e231d97SPeter Brune for (j = 0; j < i; j++) { 2908e231d97SPeter Brune g = (it + j - l) % l; 2918e231d97SPeter Brune t += (alpha[g] - beta[g])*H(k, g); 2928e231d97SPeter Brune } 2938e231d97SPeter Brune beta[k] = t / H(k, k); 2948e231d97SPeter Brune } else { 2956bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 2968e231d97SPeter Brune beta[k] = t / dXtdF[k]; 2978e231d97SPeter Brune } 29822d28d08SBarry Smith ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 29944f7e39eSPeter Brune if (qn->monitor) { 3003af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3015ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 3023af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 30344f7e39eSPeter Brune } 3044b11644fSPeter Brune } 3054b11644fSPeter Brune PetscFunctionReturn(0); 3064b11644fSPeter Brune } 3074b11644fSPeter Brune 3084b11644fSPeter Brune #undef __FUNCT__ 3094b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 3104b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 3114b11644fSPeter Brune { 3124b11644fSPeter Brune PetscErrorCode ierr; 3139f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 31415f5eeeaSPeter Brune Vec X,Xold; 315335fb975SPeter Brune Vec F,B; 316335fb975SPeter Brune Vec Y,FPC,D,Dold; 3173af51624SPeter Brune SNESConvergedReason reason; 318b8840d6bSPeter Brune PetscInt i, i_r; 319335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 3200c777b0cSPeter Brune PetscBool lssucceed,powell,periodic; 321b8840d6bSPeter Brune PetscScalar DolddotD,DolddotDold,DdotD,YdotD; 3220ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 3230ecc9e1dSPeter Brune 3244b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 3254b11644fSPeter Brune 3266e111a19SKarl Rupp PetscFunctionBegin; 3279f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3283af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 3293af51624SPeter Brune B = snes->vec_rhs; 330b8840d6bSPeter Brune 331b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 332335fb975SPeter Brune Xold = snes->work[0]; 3333af51624SPeter Brune 3343af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 335335fb975SPeter Brune D = snes->work[1]; 336335fb975SPeter Brune Dold = snes->work[2]; 3374b11644fSPeter Brune 3384b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3394b11644fSPeter Brune 3404b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3414b11644fSPeter Brune snes->iter = 0; 3424b11644fSPeter Brune snes->norm = 0.; 3434b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 344e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 34515f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3464b11644fSPeter Brune if (snes->domainerror) { 3474b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3484b11644fSPeter Brune PetscFunctionReturn(0); 3494b11644fSPeter Brune } 350e4ed7901SPeter Brune } else { 351e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 352e4ed7901SPeter Brune } 353e4ed7901SPeter Brune 354e4ed7901SPeter Brune if (!snes->norm_init_set) { 35515f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 3564b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 357e4ed7901SPeter Brune } else { 358e4ed7901SPeter Brune fnorm = snes->norm_init; 359e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 360e4ed7901SPeter Brune } 361e4ed7901SPeter Brune 3624b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3634b11644fSPeter Brune snes->norm = fnorm; 3644b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 3654b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 3664b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3674b11644fSPeter Brune 3684b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 3694b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 3704b11644fSPeter Brune /* test convergence */ 3714b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3724b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 37370d3b23bSPeter Brune 37488f769c5SPeter Brune /* composed solve */ 375c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 376217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 377217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 378*63e7833aSPeter Brune 379*63e7833aSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,B,0);CHKERRQ(ierr); 38088d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 381*63e7833aSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,B,0);CHKERRQ(ierr); 382*63e7833aSPeter Brune 38388d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 38488d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 38588d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 38688d374b2SPeter Brune PetscFunctionReturn(0); 38788d374b2SPeter Brune } 38888d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 38988d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 39088d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 3916bf1b2e5SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 39288d374b2SPeter Brune } else { 39388d374b2SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 39488d374b2SPeter Brune } 39588d374b2SPeter Brune ierr = VecCopy(Y, D);CHKERRQ(ierr); 3963af51624SPeter Brune 397f8e15203SPeter Brune /* scale the initial update */ 3980c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 3990ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4000ecc9e1dSPeter Brune } 4010ecc9e1dSPeter Brune 4025ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 403b8840d6bSPeter Brune switch (qn->type) { 404b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 405b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 406b8840d6bSPeter Brune break; 407b8840d6bSPeter Brune case SNES_QN_BROYDEN: 408b8840d6bSPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 409b8840d6bSPeter Brune break; 410b8840d6bSPeter Brune case SNES_QN_LBFGS: 411b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 412b8840d6bSPeter Brune break; 413b8840d6bSPeter Brune } 41470d3b23bSPeter Brune /* line search for lambda */ 41570d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 41688d374b2SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 41715f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 418f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4199f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 4209f3a0142SPeter Brune if (snes->domainerror) { 4219f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 4229f3a0142SPeter Brune PetscFunctionReturn(0); 4239f3a0142SPeter Brune } 424f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 4259f3a0142SPeter Brune if (!lssucceed) { 4269f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4279f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4289f3a0142SPeter Brune break; 4299f3a0142SPeter Brune } 4309f3a0142SPeter Brune } 431f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 4320c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 43305ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4340ecc9e1dSPeter Brune } 4353af51624SPeter Brune 43688d374b2SPeter Brune /* convergence monitoring */ 437b21d5a53SPeter 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); 438b21d5a53SPeter Brune 439360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 440360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 441360c497dSPeter Brune 4428409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 4438409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4444b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 445d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4464b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 4474b11644fSPeter Brune 448c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 449217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 450217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 45188d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 45288d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 45388d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 45488d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 45588d374b2SPeter Brune PetscFunctionReturn(0); 45688d374b2SPeter Brune } 45788d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 45888d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 45988d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 46088d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 46188d374b2SPeter Brune } else { 46288d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 46388d374b2SPeter Brune } 46488d374b2SPeter Brune 4650c777b0cSPeter Brune powell = PETSC_FALSE; 4660c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 4676bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 4688e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4698e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 4708e231d97SPeter Brune ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr); 4718e231d97SPeter Brune ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr); 4728e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4738e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 4748e231d97SPeter Brune ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr); 4758e231d97SPeter Brune ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr); 4760c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 4770c777b0cSPeter Brune } 4780c777b0cSPeter Brune periodic = PETSC_FALSE; 479b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 480b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 4810c777b0cSPeter Brune } 4820c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 4830c777b0cSPeter Brune if (powell || periodic) { 4845ba6227bSPeter Brune if (qn->monitor) { 4855ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 486516fe3c3SPeter 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); 4875ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 4885ba6227bSPeter Brune } 4895ba6227bSPeter Brune i_r = -1; 4905ba6227bSPeter Brune /* general purpose update */ 4915ba6227bSPeter Brune if (snes->ops->update) { 4925ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 4935ba6227bSPeter Brune } 4940c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4950ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4960ecc9e1dSPeter Brune } 4970ecc9e1dSPeter Brune } 49870d3b23bSPeter Brune /* general purpose update */ 49970d3b23bSPeter Brune if (snes->ops->update) { 50070d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 50170d3b23bSPeter Brune } 5025ba6227bSPeter Brune } 5034b11644fSPeter Brune if (i == snes->max_its) { 5044b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 5054b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 5064b11644fSPeter Brune } 5074b11644fSPeter Brune PetscFunctionReturn(0); 5084b11644fSPeter Brune } 5094b11644fSPeter Brune 5104b11644fSPeter Brune #undef __FUNCT__ 5114b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5124b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5134b11644fSPeter Brune { 5149f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5154b11644fSPeter Brune PetscErrorCode ierr; 516335fb975SPeter Brune 5174b11644fSPeter Brune PetscFunctionBegin; 518b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 519b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 5208e231d97SPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 5218e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 5228e231d97SPeter Brune qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 5238e231d97SPeter Brune 52418aa0c0cSPeter Brune if (qn->singlereduction) { 5258e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 5268e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 5278e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 5288e231d97SPeter Brune } 529335fb975SPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 530335fb975SPeter Brune 531335fb975SPeter Brune /* set up the line search */ 5320c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5338e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5348e231d97SPeter Brune } 5354b11644fSPeter Brune PetscFunctionReturn(0); 5364b11644fSPeter Brune } 5374b11644fSPeter Brune 5384b11644fSPeter Brune #undef __FUNCT__ 5394b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5404b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5414b11644fSPeter Brune { 5424b11644fSPeter Brune PetscErrorCode ierr; 5439f83bee8SJed Brown SNES_QN *qn; 5440adebc6cSBarry Smith 5454b11644fSPeter Brune PetscFunctionBegin; 5464b11644fSPeter Brune if (snes->data) { 5479f83bee8SJed Brown qn = (SNES_QN*)snes->data; 548b8840d6bSPeter Brune if (qn->U) { 549b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5504b11644fSPeter Brune } 551b8840d6bSPeter Brune if (qn->V) { 552b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5534b11644fSPeter Brune } 55418aa0c0cSPeter Brune if (qn->singlereduction) { 5558e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5568e231d97SPeter Brune } 5578e231d97SPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 5584b11644fSPeter Brune } 5594b11644fSPeter Brune PetscFunctionReturn(0); 5604b11644fSPeter Brune } 5614b11644fSPeter Brune 5624b11644fSPeter Brune #undef __FUNCT__ 5634b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5644b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5654b11644fSPeter Brune { 5664b11644fSPeter Brune PetscErrorCode ierr; 5676e111a19SKarl Rupp 5684b11644fSPeter Brune PetscFunctionBegin; 5694b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 5704b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 5719c05af5eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 5724b11644fSPeter Brune PetscFunctionReturn(0); 5734b11644fSPeter Brune } 5744b11644fSPeter Brune 5754b11644fSPeter Brune #undef __FUNCT__ 5764b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 5774b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 5784b11644fSPeter Brune { 5794b11644fSPeter Brune 5804b11644fSPeter Brune PetscErrorCode ierr; 5812150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 58288f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 583f1c6b773SPeter Brune SNESLineSearch linesearch; 5842150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 5852150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 5862150357eSBarry Smith 5874b11644fSPeter Brune PetscFunctionBegin; 5884b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 5895b4627e8SPeter Brune ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,PETSC_NULL);CHKERRQ(ierr); 5905b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr); 5915b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr); 5925b4627e8SPeter Brune ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 5934d116c7dSPeter Brune ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, PETSC_NULL);CHKERRQ(ierr); 59488f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 59588f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 59688f769c5SPeter Brune 59788f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 59888f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 59988f769c5SPeter Brune 600b8840d6bSPeter Brune ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 601b8840d6bSPeter Brune (PetscEnum)qn->type,(PetscEnum*)&qn->type,PETSC_NULL);CHKERRQ(ierr); 6024b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6039e764e56SPeter Brune if (!snes->linesearch) { 604f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 6051a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 6069e764e56SPeter Brune } 60744f7e39eSPeter Brune if (monflg) { 60844f7e39eSPeter Brune qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 60944f7e39eSPeter Brune } 6104b11644fSPeter Brune PetscFunctionReturn(0); 6114b11644fSPeter Brune } 6124b11644fSPeter Brune 61392c02d66SPeter Brune 6140c777b0cSPeter Brune #undef __FUNCT__ 6150c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6160c777b0cSPeter Brune /*@ 6170c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6180c777b0cSPeter Brune 6190c777b0cSPeter Brune Logically Collective on SNES 6200c777b0cSPeter Brune 6210c777b0cSPeter Brune Input Parameters: 6220c777b0cSPeter Brune + snes - the iterative context 6230c777b0cSPeter Brune - rtype - restart type 6240c777b0cSPeter Brune 6250c777b0cSPeter Brune Options Database: 6260c777b0cSPeter Brune + -snes_qn_restart_type<powell,periodic,none> - set the restart type 627b8840d6bSPeter Brune - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 6280c777b0cSPeter Brune 6290c777b0cSPeter Brune Level: intermediate 6300c777b0cSPeter Brune 6310c777b0cSPeter Brune SNESQNRestartTypes: 6320c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6330c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6340c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6350c777b0cSPeter Brune 6360c777b0cSPeter Brune Notes: 6370c777b0cSPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 6380c777b0cSPeter Brune 6390c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 6400c777b0cSPeter Brune @*/ 6412150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6422150357eSBarry Smith { 6430c777b0cSPeter Brune PetscErrorCode ierr; 6446e111a19SKarl Rupp 6450c777b0cSPeter Brune PetscFunctionBegin; 6460c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6470c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6480c777b0cSPeter Brune PetscFunctionReturn(0); 6490c777b0cSPeter Brune } 6500c777b0cSPeter Brune 6510c777b0cSPeter Brune #undef __FUNCT__ 6520c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6530c777b0cSPeter Brune /*@ 6540c777b0cSPeter Brune SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 6550c777b0cSPeter Brune 6560c777b0cSPeter Brune Logically Collective on SNES 6570c777b0cSPeter Brune 6580c777b0cSPeter Brune Input Parameters: 6590c777b0cSPeter Brune + snes - the iterative context 6600c777b0cSPeter Brune - stype - scale type 6610c777b0cSPeter Brune 6620c777b0cSPeter Brune Options Database: 6630c777b0cSPeter Brune . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 6640c777b0cSPeter Brune 6650c777b0cSPeter Brune Level: intermediate 6660c777b0cSPeter Brune 6670c777b0cSPeter Brune SNESQNSelectTypes: 6680c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 6690c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 6700c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 6710c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 6720c777b0cSPeter Brune 6730c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 6740c777b0cSPeter Brune @*/ 6750c777b0cSPeter Brune 6762150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 6772150357eSBarry Smith { 6780c777b0cSPeter Brune PetscErrorCode ierr; 6796e111a19SKarl Rupp 6800c777b0cSPeter Brune PetscFunctionBegin; 6810c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6820c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 6830c777b0cSPeter Brune PetscFunctionReturn(0); 6840c777b0cSPeter Brune } 6850c777b0cSPeter Brune 6860c777b0cSPeter Brune EXTERN_C_BEGIN 6870c777b0cSPeter Brune #undef __FUNCT__ 6880c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 6890adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 6900adebc6cSBarry Smith { 6910c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 6926e111a19SKarl Rupp 6930c777b0cSPeter Brune PetscFunctionBegin; 6940c777b0cSPeter Brune qn->scale_type = stype; 6950c777b0cSPeter Brune PetscFunctionReturn(0); 6960c777b0cSPeter Brune } 6970c777b0cSPeter Brune 6980c777b0cSPeter Brune #undef __FUNCT__ 6990c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7000adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 7010adebc6cSBarry Smith { 7020c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 7036e111a19SKarl Rupp 7040c777b0cSPeter Brune PetscFunctionBegin; 7050c777b0cSPeter Brune qn->restart_type = rtype; 7060c777b0cSPeter Brune PetscFunctionReturn(0); 7070c777b0cSPeter Brune } 7080c777b0cSPeter Brune EXTERN_C_END 7090c777b0cSPeter Brune 7104b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 7114b11644fSPeter Brune /*MC 7124b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 7134b11644fSPeter Brune 7146cc8130cSPeter Brune Options Database: 7156cc8130cSPeter Brune 7166cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 7171867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 7181867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 71972835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 72044f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 7216cc8130cSPeter Brune 722b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 723b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 724b8840d6bSPeter Brune updates. 7256cc8130cSPeter Brune 7261867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 7271867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 7281867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 7291867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 7301867fe5bSPeter Brune 7316cc8130cSPeter Brune References: 7326cc8130cSPeter Brune 7336cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 7346cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 7356cc8130cSPeter Brune 736b8840d6bSPeter Brune Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker 737b8840d6bSPeter Brune SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 738b8840d6bSPeter Brune 7394b11644fSPeter Brune 7404b11644fSPeter Brune Level: beginner 7414b11644fSPeter Brune 74204d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 7436cc8130cSPeter Brune 7444b11644fSPeter Brune M*/ 7454b11644fSPeter Brune EXTERN_C_BEGIN 7464b11644fSPeter Brune #undef __FUNCT__ 7474b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 7484b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 7494b11644fSPeter Brune { 7504b11644fSPeter Brune PetscErrorCode ierr; 7519f83bee8SJed Brown SNES_QN *qn; 7524b11644fSPeter Brune 7534b11644fSPeter Brune PetscFunctionBegin; 7544b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 7554b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 7564b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 7574b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 7584b11644fSPeter Brune snes->ops->view = 0; 7594b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 7604b11644fSPeter Brune 76142f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 76242f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 76342f4f86dSBarry Smith 76488976e71SPeter Brune if (!snes->tolerancesset) { 7650e444f03SPeter Brune snes->max_funcs = 30000; 7660e444f03SPeter Brune snes->max_its = 10000; 76788976e71SPeter Brune } 7680e444f03SPeter Brune 7699f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 7704b11644fSPeter Brune snes->data = (void *) qn; 7710ecc9e1dSPeter Brune qn->m = 10; 772b21d5a53SPeter Brune qn->scaling = 1.0; 773b8840d6bSPeter Brune qn->U = PETSC_NULL; 774b8840d6bSPeter Brune qn->V = PETSC_NULL; 7758e231d97SPeter Brune qn->dXtdF = PETSC_NULL; 7768e231d97SPeter Brune qn->dFtdX = PETSC_NULL; 7778e231d97SPeter Brune qn->dXdFmat = PETSC_NULL; 77844f7e39eSPeter Brune qn->monitor = PETSC_NULL; 77918aa0c0cSPeter Brune qn->singlereduction = PETSC_FALSE; 780b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 7810c777b0cSPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 7820c777b0cSPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 783b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 784ea630c6eSPeter Brune 78588f769c5SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetScaleType_C","SNESQNSetScaleType_QN",SNESQNSetScaleType_QN);CHKERRQ(ierr); 78688f769c5SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetRestartType_C","SNESQNSetRestartType_QN",SNESQNSetRestartType_QN);CHKERRQ(ierr); 7874b11644fSPeter Brune PetscFunctionReturn(0); 7884b11644fSPeter Brune } 7898e231d97SPeter Brune 7904b11644fSPeter Brune EXTERN_C_END 791