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, 11*0ad7597dSKarl Rupp SNES_QN_BADBROYDEN = 2 12*0ad7597dSKarl Rupp } SNESQNType; 136bf1b2e5SPeter Brune 144b11644fSPeter Brune typedef struct { 15b8840d6bSPeter Brune Vec *U; /* Stored past states (vary from method to method) */ 16b8840d6bSPeter Brune Vec *V; /* Stored past states (vary from method to method) */ 178e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 188e231d97SPeter Brune PetscScalar *alpha, *beta; 198e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 2018aa0c0cSPeter Brune PetscBool singlereduction; /* Aggregated reduction implementation */ 218e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 2244f7e39eSPeter Brune PetscViewer monitor; 236bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 246bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 25b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 2688d374b2SPeter Brune 27b8840d6bSPeter Brune SNESQNType type; /* the type of quasi-newton method used */ 2888f769c5SPeter Brune SNESQNScaleType scale_type; /* the type of scaling used */ 290c777b0cSPeter Brune SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 309f83bee8SJed Brown } SNES_QN; 314b11644fSPeter Brune 324b11644fSPeter Brune #undef __FUNCT__ 33b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden" 342150357eSBarry Smith PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold, Vec D, Vec Dold) 352150357eSBarry Smith { 364b11644fSPeter Brune PetscErrorCode ierr; 379f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 38b8840d6bSPeter Brune Vec W = snes->work[3]; 39b8840d6bSPeter Brune Vec *U = qn->U; 40b8840d6bSPeter Brune Vec *V = qn->V; 41b8840d6bSPeter Brune KSPConvergedReason kspreason; 42b8840d6bSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 43b8840d6bSPeter Brune PetscInt k,i,lits; 44b8840d6bSPeter Brune PetscInt m = qn->m; 45b8840d6bSPeter Brune PetscScalar gdot; 46b8840d6bSPeter Brune PetscInt l = m; 47b8840d6bSPeter Brune Mat jac,jac_pre; 482150357eSBarry Smith 49b8840d6bSPeter Brune PetscFunctionBegin; 50b8840d6bSPeter Brune if (it < m) l = it; 51b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 52b8840d6bSPeter Brune ierr = SNESGetJacobian(snes,&jac,&jac_pre,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 53b8840d6bSPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 54b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,D,W);CHKERRQ(ierr); 55b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 56b8840d6bSPeter Brune if (kspreason < 0) { 57b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 58b8840d6bSPeter 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); 59b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 60b8840d6bSPeter Brune PetscFunctionReturn(0); 61b8840d6bSPeter Brune } 62b8840d6bSPeter Brune } 63b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 64b8840d6bSPeter Brune snes->linear_its += lits; 65b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 66b8840d6bSPeter Brune } else { 67b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 68b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 69b8840d6bSPeter Brune } 70b8840d6bSPeter Brune 71b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 72b8840d6bSPeter Brune if (it > 1) { 73b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 74b8840d6bSPeter Brune k = (it+i-l)%l; 75b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 76b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 77b8840d6bSPeter Brune if (qn->monitor) { 78b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 79b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 80b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 81b8840d6bSPeter Brune } 82b8840d6bSPeter Brune } 83b8840d6bSPeter Brune } 84b8840d6bSPeter Brune if (it < m) l = it; 85b8840d6bSPeter Brune 86b8840d6bSPeter Brune /* set W to be the last step, post-linesearch */ 87b8840d6bSPeter Brune ierr = VecCopy(Xold,W);CHKERRQ(ierr); 88b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 89b8840d6bSPeter Brune 90b8840d6bSPeter Brune if (l > 0) { 91b8840d6bSPeter Brune k = (it-1)%l; 92b8840d6bSPeter Brune ierr = VecCopy(W,U[k]);CHKERRQ(ierr); 93b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,Y);CHKERRQ(ierr); 94b8840d6bSPeter Brune ierr = VecDot(U[k],W,&gdot);CHKERRQ(ierr); 95b8840d6bSPeter Brune ierr = VecCopy(Y,V[k]);CHKERRQ(ierr); 96b8840d6bSPeter Brune ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr); 97b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 98b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 99b8840d6bSPeter Brune if (qn->monitor) { 100b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 101b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 102b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 103b8840d6bSPeter Brune } 104b8840d6bSPeter Brune } 105b8840d6bSPeter Brune PetscFunctionReturn(0); 106b8840d6bSPeter Brune } 107b8840d6bSPeter Brune 108b8840d6bSPeter Brune #undef __FUNCT__ 109b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden" 1100adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1110adebc6cSBarry Smith { 112b8840d6bSPeter Brune PetscErrorCode ierr; 113b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 114b8840d6bSPeter Brune Vec W = snes->work[3]; 115b8840d6bSPeter Brune Vec *U = qn->U; 116b8840d6bSPeter Brune Vec *T = qn->V; 117b8840d6bSPeter Brune 118b8840d6bSPeter Brune /* ksp thing for jacobian scaling */ 119b8840d6bSPeter Brune KSPConvergedReason kspreason; 120b8840d6bSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 121b8840d6bSPeter Brune PetscInt k, i, lits; 122b8840d6bSPeter Brune PetscInt m = qn->m; 123b8840d6bSPeter Brune PetscScalar gdot; 124b8840d6bSPeter Brune PetscInt l = m; 125b8840d6bSPeter Brune Mat jac, jac_pre; 1260adebc6cSBarry Smith 127b8840d6bSPeter Brune PetscFunctionBegin; 128b8840d6bSPeter Brune if (it < m) l = it; 129b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 130b8840d6bSPeter Brune if (l > 0) { 131b8840d6bSPeter Brune k = (it-1)%l; 132b8840d6bSPeter Brune ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 133b8840d6bSPeter Brune ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 134b8840d6bSPeter Brune ierr = VecDot(U[k],U[k],&gdot);CHKERRQ(ierr); 135b8840d6bSPeter Brune ierr = VecCopy(D,T[k]);CHKERRQ(ierr); 136b8840d6bSPeter Brune ierr = VecScale(T[k],1.0/gdot);CHKERRQ(ierr); 137b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 138b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 139b8840d6bSPeter Brune if (qn->monitor) { 140b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 141b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 142b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 143b8840d6bSPeter Brune } 144b8840d6bSPeter Brune } 145b8840d6bSPeter Brune 146b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 147b8840d6bSPeter Brune if (it > 1) { 148b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 149b8840d6bSPeter Brune k = (it+i-l)%l; 150b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 151b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 152b8840d6bSPeter Brune if (qn->monitor) { 153b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 154b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 155b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 156b8840d6bSPeter Brune } 157b8840d6bSPeter Brune } 158b8840d6bSPeter Brune } 159b8840d6bSPeter Brune 160b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 161b8840d6bSPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 162b8840d6bSPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 163b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr); 164b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 165b8840d6bSPeter Brune if (kspreason < 0) { 166b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 167b8840d6bSPeter 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); 168b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 169b8840d6bSPeter Brune PetscFunctionReturn(0); 170b8840d6bSPeter Brune } 171b8840d6bSPeter Brune } 172b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 173b8840d6bSPeter Brune snes->linear_its += lits; 174b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 175b8840d6bSPeter Brune } else { 176b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 177b8840d6bSPeter Brune } 178b8840d6bSPeter Brune PetscFunctionReturn(0); 179b8840d6bSPeter Brune } 180b8840d6bSPeter Brune 181b8840d6bSPeter Brune #undef __FUNCT__ 182b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 1830adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1840adebc6cSBarry Smith { 185b8840d6bSPeter Brune PetscErrorCode ierr; 186b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 187b8840d6bSPeter Brune Vec W = snes->work[3]; 188b8840d6bSPeter Brune Vec *dX = qn->U; 189b8840d6bSPeter Brune Vec *dF = qn->V; 190bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 191bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 1928e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 193b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 1948e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 195bd052dfeSPeter Brune 1960ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 1970ecc9e1dSPeter Brune KSPConvergedReason kspreason; 1980ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 1998e231d97SPeter Brune PetscInt k,i,j,g,lits; 2004b11644fSPeter Brune PetscInt m = qn->m; 2014b11644fSPeter Brune PetscScalar t; 2024b11644fSPeter Brune PetscInt l = m; 2038e231d97SPeter Brune Mat jac,jac_pre; 2048e231d97SPeter Brune 2054b11644fSPeter Brune PetscFunctionBegin; 2065ba6227bSPeter Brune if (it < m) l = it; 207b8840d6bSPeter Brune if (it > 0) { 208b8840d6bSPeter Brune k = (it - 1) % l; 209b8840d6bSPeter Brune ierr = VecCopy(D, dF[k]);CHKERRQ(ierr); 210b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 211b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 212b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 21318aa0c0cSPeter Brune if (qn->singlereduction) { 214b8840d6bSPeter Brune ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr); 215b8840d6bSPeter Brune ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr); 216b8840d6bSPeter Brune for (j = 0; j < l; j++) { 217b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 218b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 219b8840d6bSPeter Brune } 220b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 221b8840d6bSPeter Brune for (j = 0; j < l; j++) { 222b8840d6bSPeter Brune dXtdF[j] = H(j, j); 223b8840d6bSPeter Brune } 224b8840d6bSPeter Brune } else { 225b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 226b8840d6bSPeter Brune } 227b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 22867392de3SBarry Smith PetscReal dFtdF; 22967392de3SBarry Smith ierr = VecDotRealPart(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 23067392de3SBarry Smith qn->scaling = PetscRealPart(dXtdF[k])/dFtdF; 231b8840d6bSPeter Brune } else if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 232b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 233b8840d6bSPeter Brune } 234b8840d6bSPeter Brune } 235b8840d6bSPeter Brune 236b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 237b8840d6bSPeter Brune 238b8840d6bSPeter Brune if (qn->singlereduction) { 239b8840d6bSPeter Brune ierr = VecMDot(Y,l,dX,YtdX);CHKERRQ(ierr); 2408e231d97SPeter Brune } 2414b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2424b11644fSPeter Brune for (i=0;i<l;i++) { 243b21d5a53SPeter Brune k = (it-i-1)%l; 24418aa0c0cSPeter Brune if (qn->singlereduction) { 2458e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2468e231d97SPeter Brune t = YtdX[k]; 2478e231d97SPeter Brune for (j=0;j<i;j++) { 2488e231d97SPeter Brune g = (it-j-1)%l; 2498e231d97SPeter Brune t += -alpha[g]*H(g, k); 2508e231d97SPeter Brune } 2518e231d97SPeter Brune alpha[k] = t/H(k,k); 2528e231d97SPeter Brune } else { 2533af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2548e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2558e231d97SPeter Brune } 25644f7e39eSPeter Brune if (qn->monitor) { 2573af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2585ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2593af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 26044f7e39eSPeter Brune } 2616bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2624b11644fSPeter Brune } 2634b11644fSPeter Brune 2640c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 2658e231d97SPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 2668e231d97SPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 267b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr); 2680ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2690ecc9e1dSPeter Brune if (kspreason < 0) { 2700ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 2710ecc9e1dSPeter 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); 2720ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 2730ecc9e1dSPeter Brune PetscFunctionReturn(0); 2740ecc9e1dSPeter Brune } 2750ecc9e1dSPeter Brune } 2760ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2770ecc9e1dSPeter Brune snes->linear_its += lits; 278b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2790ecc9e1dSPeter Brune } else { 280b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2810ecc9e1dSPeter Brune } 28218aa0c0cSPeter Brune if (qn->singlereduction) { 283b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 2848e231d97SPeter Brune } 2854b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 286bd052dfeSPeter Brune for (i = 0; i < l; i++) { 287b21d5a53SPeter Brune k = (it + i - l) % l; 28818aa0c0cSPeter Brune if (qn->singlereduction) { 2898e231d97SPeter Brune t = YtdX[k]; 2908e231d97SPeter Brune for (j = 0; j < i; j++) { 2918e231d97SPeter Brune g = (it + j - l) % l; 2928e231d97SPeter Brune t += (alpha[g] - beta[g])*H(k, g); 2938e231d97SPeter Brune } 2948e231d97SPeter Brune beta[k] = t / H(k, k); 2958e231d97SPeter Brune } else { 2966bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 2978e231d97SPeter Brune beta[k] = t / dXtdF[k]; 2988e231d97SPeter Brune } 29922d28d08SBarry Smith ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 30044f7e39eSPeter Brune if (qn->monitor) { 3013af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3025ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 3033af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 30444f7e39eSPeter Brune } 3054b11644fSPeter Brune } 3064b11644fSPeter Brune PetscFunctionReturn(0); 3074b11644fSPeter Brune } 3084b11644fSPeter Brune 3094b11644fSPeter Brune #undef __FUNCT__ 3104b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 3114b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 3124b11644fSPeter Brune { 3134b11644fSPeter Brune PetscErrorCode ierr; 3149f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 31515f5eeeaSPeter Brune Vec X,Xold; 316335fb975SPeter Brune Vec F,B; 317335fb975SPeter Brune Vec Y,FPC,D,Dold; 3183af51624SPeter Brune SNESConvergedReason reason; 319b8840d6bSPeter Brune PetscInt i, i_r; 320335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 3210c777b0cSPeter Brune PetscBool lssucceed,powell,periodic; 322b8840d6bSPeter Brune PetscScalar DolddotD,DolddotDold,DdotD,YdotD; 3230ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 3240ecc9e1dSPeter Brune 3254b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 3264b11644fSPeter Brune 3276e111a19SKarl Rupp PetscFunctionBegin; 3289f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3293af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 3303af51624SPeter Brune B = snes->vec_rhs; 331b8840d6bSPeter Brune 332b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 333335fb975SPeter Brune Xold = snes->work[0]; 3343af51624SPeter Brune 3353af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 336335fb975SPeter Brune D = snes->work[1]; 337335fb975SPeter Brune Dold = snes->work[2]; 3384b11644fSPeter Brune 3394b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3404b11644fSPeter Brune 3414b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3424b11644fSPeter Brune snes->iter = 0; 3434b11644fSPeter Brune snes->norm = 0.; 3444b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 345e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 34615f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3474b11644fSPeter Brune if (snes->domainerror) { 3484b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3494b11644fSPeter Brune PetscFunctionReturn(0); 3504b11644fSPeter Brune } 351e4ed7901SPeter Brune } else { 352e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 353e4ed7901SPeter Brune } 354e4ed7901SPeter Brune 355e4ed7901SPeter Brune if (!snes->norm_init_set) { 35615f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 3574b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 358e4ed7901SPeter Brune } else { 359e4ed7901SPeter Brune fnorm = snes->norm_init; 360e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 361e4ed7901SPeter Brune } 362e4ed7901SPeter Brune 3634b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3644b11644fSPeter Brune snes->norm = fnorm; 3654b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 3664b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 3674b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3684b11644fSPeter Brune 3694b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 3704b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 3714b11644fSPeter Brune /* test convergence */ 3724b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3734b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 37470d3b23bSPeter Brune 37588f769c5SPeter Brune /* composed solve */ 376c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 377217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 378217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 37963e7833aSPeter Brune 38063e7833aSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,B,0);CHKERRQ(ierr); 38188d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 38263e7833aSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,B,0);CHKERRQ(ierr); 38363e7833aSPeter Brune 38488d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 38588d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 38688d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 38788d374b2SPeter Brune PetscFunctionReturn(0); 38888d374b2SPeter Brune } 38988d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 39088d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 39188d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 3926bf1b2e5SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 39388d374b2SPeter Brune } else { 39488d374b2SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 39588d374b2SPeter Brune } 39688d374b2SPeter Brune ierr = VecCopy(Y, D);CHKERRQ(ierr); 3973af51624SPeter Brune 398f8e15203SPeter Brune /* scale the initial update */ 3990c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4000ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4010ecc9e1dSPeter Brune } 4020ecc9e1dSPeter Brune 4035ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 404b8840d6bSPeter Brune switch (qn->type) { 405b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 406b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 407b8840d6bSPeter Brune break; 408b8840d6bSPeter Brune case SNES_QN_BROYDEN: 409b8840d6bSPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 410b8840d6bSPeter Brune break; 411b8840d6bSPeter Brune case SNES_QN_LBFGS: 412b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 413b8840d6bSPeter Brune break; 414b8840d6bSPeter Brune } 41570d3b23bSPeter Brune /* line search for lambda */ 41670d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 41788d374b2SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 41815f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 419f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4209f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 4219f3a0142SPeter Brune if (snes->domainerror) { 4229f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 4239f3a0142SPeter Brune PetscFunctionReturn(0); 4249f3a0142SPeter Brune } 425f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 4269f3a0142SPeter Brune if (!lssucceed) { 4279f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4289f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4299f3a0142SPeter Brune break; 4309f3a0142SPeter Brune } 4319f3a0142SPeter Brune } 432f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 4330c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 43405ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4350ecc9e1dSPeter Brune } 4363af51624SPeter Brune 43788d374b2SPeter Brune /* convergence monitoring */ 438b21d5a53SPeter 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); 439b21d5a53SPeter Brune 440360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 441360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 442360c497dSPeter Brune 4438409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 4448409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4454b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 446d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4474b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 4484b11644fSPeter Brune 449c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 450217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 451217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 45288d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 45388d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 45488d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 45588d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 45688d374b2SPeter Brune PetscFunctionReturn(0); 45788d374b2SPeter Brune } 45888d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 45988d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 46088d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 46188d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 46288d374b2SPeter Brune } else { 46388d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 46488d374b2SPeter Brune } 46588d374b2SPeter Brune 4660c777b0cSPeter Brune powell = PETSC_FALSE; 4670c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 4686bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 4698e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4708e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 4718e231d97SPeter Brune ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr); 4728e231d97SPeter Brune ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr); 4738e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4748e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 4758e231d97SPeter Brune ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr); 4768e231d97SPeter Brune ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr); 4770c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 4780c777b0cSPeter Brune } 4790c777b0cSPeter Brune periodic = PETSC_FALSE; 480b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 481b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 4820c777b0cSPeter Brune } 4830c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 4840c777b0cSPeter Brune if (powell || periodic) { 4855ba6227bSPeter Brune if (qn->monitor) { 4865ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 487516fe3c3SPeter 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); 4885ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 4895ba6227bSPeter Brune } 4905ba6227bSPeter Brune i_r = -1; 4915ba6227bSPeter Brune /* general purpose update */ 4925ba6227bSPeter Brune if (snes->ops->update) { 4935ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 4945ba6227bSPeter Brune } 4950c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4960ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4970ecc9e1dSPeter Brune } 4980ecc9e1dSPeter Brune } 49970d3b23bSPeter Brune /* general purpose update */ 50070d3b23bSPeter Brune if (snes->ops->update) { 50170d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 50270d3b23bSPeter Brune } 5035ba6227bSPeter Brune } 5044b11644fSPeter Brune if (i == snes->max_its) { 5054b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 5064b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 5074b11644fSPeter Brune } 5084b11644fSPeter Brune PetscFunctionReturn(0); 5094b11644fSPeter Brune } 5104b11644fSPeter Brune 5114b11644fSPeter Brune #undef __FUNCT__ 5124b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5134b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5144b11644fSPeter Brune { 5159f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5164b11644fSPeter Brune PetscErrorCode ierr; 517335fb975SPeter Brune 5184b11644fSPeter Brune PetscFunctionBegin; 519b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 520b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 5218e231d97SPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 5228e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 5238e231d97SPeter Brune qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 5248e231d97SPeter Brune 52518aa0c0cSPeter Brune if (qn->singlereduction) { 5268e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 5278e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 5288e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 5298e231d97SPeter Brune } 530335fb975SPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 531335fb975SPeter Brune 532335fb975SPeter Brune /* set up the line search */ 5330c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5348e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5358e231d97SPeter Brune } 5364b11644fSPeter Brune PetscFunctionReturn(0); 5374b11644fSPeter Brune } 5384b11644fSPeter Brune 5394b11644fSPeter Brune #undef __FUNCT__ 5404b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5414b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5424b11644fSPeter Brune { 5434b11644fSPeter Brune PetscErrorCode ierr; 5449f83bee8SJed Brown SNES_QN *qn; 5450adebc6cSBarry Smith 5464b11644fSPeter Brune PetscFunctionBegin; 5474b11644fSPeter Brune if (snes->data) { 5489f83bee8SJed Brown qn = (SNES_QN*)snes->data; 549b8840d6bSPeter Brune if (qn->U) { 550b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5514b11644fSPeter Brune } 552b8840d6bSPeter Brune if (qn->V) { 553b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5544b11644fSPeter Brune } 55518aa0c0cSPeter Brune if (qn->singlereduction) { 5568e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5578e231d97SPeter Brune } 5588e231d97SPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 5594b11644fSPeter Brune } 5604b11644fSPeter Brune PetscFunctionReturn(0); 5614b11644fSPeter Brune } 5624b11644fSPeter Brune 5634b11644fSPeter Brune #undef __FUNCT__ 5644b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5654b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5664b11644fSPeter Brune { 5674b11644fSPeter Brune PetscErrorCode ierr; 5686e111a19SKarl Rupp 5694b11644fSPeter Brune PetscFunctionBegin; 5704b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 5714b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 5729c05af5eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 5734b11644fSPeter Brune PetscFunctionReturn(0); 5744b11644fSPeter Brune } 5754b11644fSPeter Brune 5764b11644fSPeter Brune #undef __FUNCT__ 5774b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 5784b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 5794b11644fSPeter Brune { 5804b11644fSPeter Brune 5814b11644fSPeter Brune PetscErrorCode ierr; 5822150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 58388f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 584f1c6b773SPeter Brune SNESLineSearch linesearch; 5852150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 5862150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 5872150357eSBarry Smith 5884b11644fSPeter Brune PetscFunctionBegin; 5894b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 5905b4627e8SPeter Brune ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,PETSC_NULL);CHKERRQ(ierr); 5915b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr); 5925b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr); 5935b4627e8SPeter Brune ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 5944d116c7dSPeter Brune ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, PETSC_NULL);CHKERRQ(ierr); 59588f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 59688f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 59788f769c5SPeter Brune 59888f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 59988f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 60088f769c5SPeter Brune 601b8840d6bSPeter Brune ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 602b8840d6bSPeter Brune (PetscEnum)qn->type,(PetscEnum*)&qn->type,PETSC_NULL);CHKERRQ(ierr); 6034b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6049e764e56SPeter Brune if (!snes->linesearch) { 605f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 6061a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 6079e764e56SPeter Brune } 60844f7e39eSPeter Brune if (monflg) { 60944f7e39eSPeter Brune qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 61044f7e39eSPeter Brune } 6114b11644fSPeter Brune PetscFunctionReturn(0); 6124b11644fSPeter Brune } 6134b11644fSPeter Brune 61492c02d66SPeter Brune 6150c777b0cSPeter Brune #undef __FUNCT__ 6160c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6170c777b0cSPeter Brune /*@ 6180c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6190c777b0cSPeter Brune 6200c777b0cSPeter Brune Logically Collective on SNES 6210c777b0cSPeter Brune 6220c777b0cSPeter Brune Input Parameters: 6230c777b0cSPeter Brune + snes - the iterative context 6240c777b0cSPeter Brune - rtype - restart type 6250c777b0cSPeter Brune 6260c777b0cSPeter Brune Options Database: 6270c777b0cSPeter Brune + -snes_qn_restart_type<powell,periodic,none> - set the restart type 628b8840d6bSPeter Brune - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 6290c777b0cSPeter Brune 6300c777b0cSPeter Brune Level: intermediate 6310c777b0cSPeter Brune 6320c777b0cSPeter Brune SNESQNRestartTypes: 6330c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6340c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6350c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6360c777b0cSPeter Brune 6370c777b0cSPeter Brune Notes: 6380c777b0cSPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 6390c777b0cSPeter Brune 6400c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 6410c777b0cSPeter Brune @*/ 6422150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6432150357eSBarry Smith { 6440c777b0cSPeter Brune PetscErrorCode ierr; 6456e111a19SKarl Rupp 6460c777b0cSPeter Brune PetscFunctionBegin; 6470c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6480c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6490c777b0cSPeter Brune PetscFunctionReturn(0); 6500c777b0cSPeter Brune } 6510c777b0cSPeter Brune 6520c777b0cSPeter Brune #undef __FUNCT__ 6530c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6540c777b0cSPeter Brune /*@ 6550c777b0cSPeter Brune SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 6560c777b0cSPeter Brune 6570c777b0cSPeter Brune Logically Collective on SNES 6580c777b0cSPeter Brune 6590c777b0cSPeter Brune Input Parameters: 6600c777b0cSPeter Brune + snes - the iterative context 6610c777b0cSPeter Brune - stype - scale type 6620c777b0cSPeter Brune 6630c777b0cSPeter Brune Options Database: 6640c777b0cSPeter Brune . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 6650c777b0cSPeter Brune 6660c777b0cSPeter Brune Level: intermediate 6670c777b0cSPeter Brune 6680c777b0cSPeter Brune SNESQNSelectTypes: 6690c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 6700c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 6710c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 6720c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 6730c777b0cSPeter Brune 6740c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 6750c777b0cSPeter Brune @*/ 6760c777b0cSPeter Brune 6772150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 6782150357eSBarry Smith { 6790c777b0cSPeter Brune PetscErrorCode ierr; 6806e111a19SKarl Rupp 6810c777b0cSPeter Brune PetscFunctionBegin; 6820c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6830c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 6840c777b0cSPeter Brune PetscFunctionReturn(0); 6850c777b0cSPeter Brune } 6860c777b0cSPeter Brune 6870c777b0cSPeter Brune EXTERN_C_BEGIN 6880c777b0cSPeter Brune #undef __FUNCT__ 6890c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 6900adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 6910adebc6cSBarry Smith { 6920c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 6936e111a19SKarl Rupp 6940c777b0cSPeter Brune PetscFunctionBegin; 6950c777b0cSPeter Brune qn->scale_type = stype; 6960c777b0cSPeter Brune PetscFunctionReturn(0); 6970c777b0cSPeter Brune } 6980c777b0cSPeter Brune 6990c777b0cSPeter Brune #undef __FUNCT__ 7000c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7010adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 7020adebc6cSBarry Smith { 7030c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 7046e111a19SKarl Rupp 7050c777b0cSPeter Brune PetscFunctionBegin; 7060c777b0cSPeter Brune qn->restart_type = rtype; 7070c777b0cSPeter Brune PetscFunctionReturn(0); 7080c777b0cSPeter Brune } 7090c777b0cSPeter Brune EXTERN_C_END 7100c777b0cSPeter Brune 7114b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 7124b11644fSPeter Brune /*MC 7134b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 7144b11644fSPeter Brune 7156cc8130cSPeter Brune Options Database: 7166cc8130cSPeter Brune 7176cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 7181867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 7191867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 72072835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 72144f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 7226cc8130cSPeter Brune 723b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 724b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 725b8840d6bSPeter Brune updates. 7266cc8130cSPeter Brune 7271867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 7281867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 7291867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 7301867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 7311867fe5bSPeter Brune 7326cc8130cSPeter Brune References: 7336cc8130cSPeter Brune 7346cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 7356cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 7366cc8130cSPeter Brune 737b8840d6bSPeter Brune Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker 738b8840d6bSPeter Brune SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 739b8840d6bSPeter Brune 7404b11644fSPeter Brune 7414b11644fSPeter Brune Level: beginner 7424b11644fSPeter Brune 74304d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 7446cc8130cSPeter Brune 7454b11644fSPeter Brune M*/ 7464b11644fSPeter Brune EXTERN_C_BEGIN 7474b11644fSPeter Brune #undef __FUNCT__ 7484b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 7494b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 7504b11644fSPeter Brune { 7514b11644fSPeter Brune PetscErrorCode ierr; 7529f83bee8SJed Brown SNES_QN *qn; 7534b11644fSPeter Brune 7544b11644fSPeter Brune PetscFunctionBegin; 7554b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 7564b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 7574b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 7584b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 7594b11644fSPeter Brune snes->ops->view = 0; 7604b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 7614b11644fSPeter Brune 76242f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 76342f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 76442f4f86dSBarry Smith 76588976e71SPeter Brune if (!snes->tolerancesset) { 7660e444f03SPeter Brune snes->max_funcs = 30000; 7670e444f03SPeter Brune snes->max_its = 10000; 76888976e71SPeter Brune } 7690e444f03SPeter Brune 7709f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 7714b11644fSPeter Brune snes->data = (void *) qn; 7720ecc9e1dSPeter Brune qn->m = 10; 773b21d5a53SPeter Brune qn->scaling = 1.0; 774b8840d6bSPeter Brune qn->U = PETSC_NULL; 775b8840d6bSPeter Brune qn->V = PETSC_NULL; 7768e231d97SPeter Brune qn->dXtdF = PETSC_NULL; 7778e231d97SPeter Brune qn->dFtdX = PETSC_NULL; 7788e231d97SPeter Brune qn->dXdFmat = PETSC_NULL; 77944f7e39eSPeter Brune qn->monitor = PETSC_NULL; 78018aa0c0cSPeter Brune qn->singlereduction = PETSC_FALSE; 781b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 7820c777b0cSPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 7830c777b0cSPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 784b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 785ea630c6eSPeter Brune 78688f769c5SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetScaleType_C","SNESQNSetScaleType_QN",SNESQNSetScaleType_QN);CHKERRQ(ierr); 78788f769c5SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetRestartType_C","SNESQNSetRestartType_QN",SNESQNSetRestartType_QN);CHKERRQ(ierr); 7884b11644fSPeter Brune PetscFunctionReturn(0); 7894b11644fSPeter Brune } 7908e231d97SPeter Brune 7914b11644fSPeter Brune EXTERN_C_END 792