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, 110ad7597dSKarl Rupp SNES_QN_BADBROYDEN = 2 120ad7597dSKarl 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 */ 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 PetscInt k,i,lits; 42b8840d6bSPeter Brune PetscInt m = qn->m; 43b8840d6bSPeter Brune PetscScalar gdot; 44b8840d6bSPeter Brune PetscInt l = m; 452150357eSBarry Smith 46b8840d6bSPeter Brune PetscFunctionBegin; 47b8840d6bSPeter Brune if (it < m) l = it; 48b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 49d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr); 50b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 51b8840d6bSPeter Brune if (kspreason < 0) { 52b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 53b8840d6bSPeter 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); 54b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 55b8840d6bSPeter Brune PetscFunctionReturn(0); 56b8840d6bSPeter Brune } 57b8840d6bSPeter Brune } 58b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 59b8840d6bSPeter Brune snes->linear_its += lits; 60b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 61b8840d6bSPeter Brune } else { 62b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 63b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 64b8840d6bSPeter Brune } 65b8840d6bSPeter Brune 66b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 67b8840d6bSPeter Brune if (it > 1) { 68b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 69b8840d6bSPeter Brune k = (it+i-l)%l; 70b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 71b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 72b8840d6bSPeter Brune if (qn->monitor) { 73b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 74b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 75b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 76b8840d6bSPeter Brune } 77b8840d6bSPeter Brune } 78b8840d6bSPeter Brune } 79b8840d6bSPeter Brune if (it < m) l = it; 80b8840d6bSPeter Brune 81b8840d6bSPeter Brune /* set W to be the last step, post-linesearch */ 82b8840d6bSPeter Brune ierr = VecCopy(Xold,W);CHKERRQ(ierr); 83b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 84b8840d6bSPeter Brune 85b8840d6bSPeter Brune if (l > 0) { 86b8840d6bSPeter Brune k = (it-1)%l; 87b8840d6bSPeter Brune ierr = VecCopy(W,U[k]);CHKERRQ(ierr); 88b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,Y);CHKERRQ(ierr); 89b8840d6bSPeter Brune ierr = VecDot(U[k],W,&gdot);CHKERRQ(ierr); 90b8840d6bSPeter Brune ierr = VecCopy(Y,V[k]);CHKERRQ(ierr); 91b8840d6bSPeter Brune ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr); 92b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 93b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 94b8840d6bSPeter Brune if (qn->monitor) { 95b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 96b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 97b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 98b8840d6bSPeter Brune } 99b8840d6bSPeter Brune } 100b8840d6bSPeter Brune PetscFunctionReturn(0); 101b8840d6bSPeter Brune } 102b8840d6bSPeter Brune 103b8840d6bSPeter Brune #undef __FUNCT__ 104b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden" 1050adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1060adebc6cSBarry Smith { 107b8840d6bSPeter Brune PetscErrorCode ierr; 108b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 109b8840d6bSPeter Brune Vec W = snes->work[3]; 110b8840d6bSPeter Brune Vec *U = qn->U; 111b8840d6bSPeter Brune Vec *T = qn->V; 112b8840d6bSPeter Brune 113b8840d6bSPeter Brune /* ksp thing for jacobian scaling */ 114b8840d6bSPeter Brune KSPConvergedReason kspreason; 115b8840d6bSPeter Brune PetscInt k, i, lits; 116b8840d6bSPeter Brune PetscInt m = qn->m; 117b8840d6bSPeter Brune PetscScalar gdot; 118b8840d6bSPeter Brune PetscInt l = m; 1190adebc6cSBarry Smith 120b8840d6bSPeter Brune PetscFunctionBegin; 121b8840d6bSPeter Brune if (it < m) l = it; 122b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 123b8840d6bSPeter Brune if (l > 0) { 124b8840d6bSPeter Brune k = (it-1)%l; 125b8840d6bSPeter Brune ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 126b8840d6bSPeter Brune ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 127b8840d6bSPeter Brune ierr = VecDot(U[k],U[k],&gdot);CHKERRQ(ierr); 128b8840d6bSPeter Brune ierr = VecCopy(D,T[k]);CHKERRQ(ierr); 129b8840d6bSPeter Brune ierr = VecScale(T[k],1.0/gdot);CHKERRQ(ierr); 130b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 131b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 132b8840d6bSPeter Brune if (qn->monitor) { 133b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 134b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 135b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 136b8840d6bSPeter Brune } 137b8840d6bSPeter Brune } 138b8840d6bSPeter Brune 139b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 140b8840d6bSPeter Brune if (it > 1) { 141b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 142b8840d6bSPeter Brune k = (it+i-l)%l; 143b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 144b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 145b8840d6bSPeter Brune if (qn->monitor) { 146b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 147b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 148b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 149b8840d6bSPeter Brune } 150b8840d6bSPeter Brune } 151b8840d6bSPeter Brune } 152b8840d6bSPeter Brune 153b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 154d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 155b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 156b8840d6bSPeter Brune if (kspreason < 0) { 157b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 158b8840d6bSPeter 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); 159b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 160b8840d6bSPeter Brune PetscFunctionReturn(0); 161b8840d6bSPeter Brune } 162b8840d6bSPeter Brune } 163b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 164b8840d6bSPeter Brune snes->linear_its += lits; 165b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 166b8840d6bSPeter Brune } else { 167b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 168b8840d6bSPeter Brune } 169b8840d6bSPeter Brune PetscFunctionReturn(0); 170b8840d6bSPeter Brune } 171b8840d6bSPeter Brune 172b8840d6bSPeter Brune #undef __FUNCT__ 173b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 1740adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1750adebc6cSBarry Smith { 176b8840d6bSPeter Brune PetscErrorCode ierr; 177b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 178b8840d6bSPeter Brune Vec W = snes->work[3]; 179b8840d6bSPeter Brune Vec *dX = qn->U; 180b8840d6bSPeter Brune Vec *dF = qn->V; 181bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 182bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 1838e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 184b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 1858e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 186bd052dfeSPeter Brune 1870ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 1880ecc9e1dSPeter Brune KSPConvergedReason kspreason; 1898e231d97SPeter Brune PetscInt k,i,j,g,lits; 1904b11644fSPeter Brune PetscInt m = qn->m; 1914b11644fSPeter Brune PetscScalar t; 1924b11644fSPeter Brune PetscInt l = m; 1938e231d97SPeter Brune 1944b11644fSPeter Brune PetscFunctionBegin; 1955ba6227bSPeter Brune if (it < m) l = it; 1961c6523dcSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 197b8840d6bSPeter Brune if (it > 0) { 198b8840d6bSPeter Brune k = (it - 1) % l; 199b8840d6bSPeter Brune ierr = VecCopy(D,dF[k]);CHKERRQ(ierr); 200b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 201b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 202b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 20318aa0c0cSPeter Brune if (qn->singlereduction) { 204*23d44fbcSPeter Brune PetscScalar dFtdF; 2051c6523dcSPeter Brune ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2061c6523dcSPeter Brune ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2071c6523dcSPeter Brune ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr); 208*23d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) {ierr = VecDotBegin(dF[k],dF[k],&dFtdF);CHKERRQ(ierr);} 2091c6523dcSPeter Brune ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2101c6523dcSPeter Brune ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2111c6523dcSPeter Brune ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr); 212*23d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 213*23d44fbcSPeter Brune ierr = VecDotEnd(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 214*23d44fbcSPeter Brune qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(dFtdF); 215*23d44fbcSPeter Brune } 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 */ 2211aa26658SKarl Rupp for (j = 0; j < l; j++) dXtdF[j] = H(j, j); 222b8840d6bSPeter Brune } else { 223b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 224b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 22567392de3SBarry Smith PetscReal dFtdF; 22667392de3SBarry Smith ierr = VecDotRealPart(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 22767392de3SBarry Smith qn->scaling = PetscRealPart(dXtdF[k])/dFtdF; 228*23d44fbcSPeter Brune } 229*23d44fbcSPeter Brune } 230*23d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 231b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 232b8840d6bSPeter Brune } 233b8840d6bSPeter Brune } 234b8840d6bSPeter Brune 2354b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2364b11644fSPeter Brune for (i=0; i<l; i++) { 237b21d5a53SPeter Brune k = (it-i-1)%l; 23818aa0c0cSPeter Brune if (qn->singlereduction) { 2398e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2408e231d97SPeter Brune t = YtdX[k]; 2418e231d97SPeter Brune for (j=0; j<i; j++) { 2428e231d97SPeter Brune g = (it-j-1)%l; 2438e231d97SPeter Brune t += -alpha[g]*H(g, k); 2448e231d97SPeter Brune } 2458e231d97SPeter Brune alpha[k] = t/H(k,k); 2468e231d97SPeter Brune } else { 2473af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2488e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2498e231d97SPeter Brune } 25044f7e39eSPeter Brune if (qn->monitor) { 2513af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2525ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2533af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 25444f7e39eSPeter Brune } 2556bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2564b11644fSPeter Brune } 2574b11644fSPeter Brune 2580c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 259d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 2600ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2610ecc9e1dSPeter Brune if (kspreason < 0) { 2620ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 2630ecc9e1dSPeter 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); 2640ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 2650ecc9e1dSPeter Brune PetscFunctionReturn(0); 2660ecc9e1dSPeter Brune } 2670ecc9e1dSPeter Brune } 2680ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2690ecc9e1dSPeter Brune snes->linear_its += lits; 270b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2710ecc9e1dSPeter Brune } else { 272b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2730ecc9e1dSPeter Brune } 27418aa0c0cSPeter Brune if (qn->singlereduction) { 275b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 2768e231d97SPeter Brune } 2774b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 278bd052dfeSPeter Brune for (i = 0; i < l; i++) { 279b21d5a53SPeter Brune k = (it + i - l) % l; 28018aa0c0cSPeter Brune if (qn->singlereduction) { 2818e231d97SPeter Brune t = YtdX[k]; 2828e231d97SPeter Brune for (j = 0; j < i; j++) { 2838e231d97SPeter Brune g = (it + j - l) % l; 2848e231d97SPeter Brune t += (alpha[g] - beta[g])*H(k, g); 2858e231d97SPeter Brune } 2868e231d97SPeter Brune beta[k] = t / H(k, k); 2878e231d97SPeter Brune } else { 2886bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 2898e231d97SPeter Brune beta[k] = t / dXtdF[k]; 2908e231d97SPeter Brune } 29122d28d08SBarry Smith ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 29244f7e39eSPeter Brune if (qn->monitor) { 2933af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2945ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 2953af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 29644f7e39eSPeter Brune } 2974b11644fSPeter Brune } 2984b11644fSPeter Brune PetscFunctionReturn(0); 2994b11644fSPeter Brune } 3004b11644fSPeter Brune 3014b11644fSPeter Brune #undef __FUNCT__ 3024b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 3034b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 3044b11644fSPeter Brune { 3054b11644fSPeter Brune PetscErrorCode ierr; 3069f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 30715f5eeeaSPeter Brune Vec X,Xold; 30847f26062SPeter Brune Vec F; 30947f26062SPeter Brune Vec Y,D,Dold; 310b8840d6bSPeter Brune PetscInt i, i_r; 311335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 3120c777b0cSPeter Brune PetscBool lssucceed,powell,periodic; 3131c6523dcSPeter Brune PetscScalar DolddotD,DolddotDold; 3140ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 315b7281c8aSPeter Brune SNESConvergedReason reason; 3160ecc9e1dSPeter Brune 3174b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 3184b11644fSPeter Brune 3196e111a19SKarl Rupp PetscFunctionBegin; 3209f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3213af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 322b8840d6bSPeter Brune 323b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 324335fb975SPeter Brune Xold = snes->work[0]; 3253af51624SPeter Brune 3263af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 327335fb975SPeter Brune D = snes->work[1]; 328335fb975SPeter Brune Dold = snes->work[2]; 3294b11644fSPeter Brune 3304b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3314b11644fSPeter Brune 332ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3334b11644fSPeter Brune snes->iter = 0; 3344b11644fSPeter Brune snes->norm = 0.; 335ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 33647f26062SPeter Brune 337b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 33847f26062SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 339b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 340b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 341b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 342b7281c8aSPeter Brune PetscFunctionReturn(0); 343b7281c8aSPeter Brune } 34447f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 34547f26062SPeter Brune } else { 346e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 34715f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3484b11644fSPeter Brune if (snes->domainerror) { 3494b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3504b11644fSPeter Brune PetscFunctionReturn(0); 3514b11644fSPeter Brune } 3521aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 353e4ed7901SPeter Brune if (!snes->norm_init_set) { 35447f26062SPeter Brune /* convergence test */ 35547f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 356189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 357189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 358189a9710SBarry Smith PetscFunctionReturn(0); 359189a9710SBarry Smith } 360e4ed7901SPeter Brune } else { 361e4ed7901SPeter Brune fnorm = snes->norm_init; 362e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 363e4ed7901SPeter Brune } 36447f26062SPeter Brune } 365b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 36647f26062SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 367b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 368b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 369b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 370b7281c8aSPeter Brune PetscFunctionReturn(0); 371b7281c8aSPeter Brune } 37247f26062SPeter Brune } else { 37347f26062SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 37447f26062SPeter Brune } 375b28a06ddSPeter Brune 376ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3774b11644fSPeter Brune snes->norm = fnorm; 378ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 379a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 3804b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3814b11644fSPeter Brune 3824b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 3834b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 3844b11644fSPeter Brune /* test convergence */ 3854b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3864b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 38770d3b23bSPeter Brune 3883cf07b75SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 3893cf07b75SPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 3903cf07b75SPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 3913cf07b75SPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 392ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 393ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 394ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 395ddd40ce5SPeter Brune PetscFunctionReturn(0); 396ddd40ce5SPeter Brune } 397ddd40ce5SPeter Brune ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 3983cf07b75SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 3993cf07b75SPeter Brune } 4003cf07b75SPeter Brune 401f8e15203SPeter Brune /* scale the initial update */ 4020c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4030ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4048cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 4050ecc9e1dSPeter Brune } 4060ecc9e1dSPeter Brune 4075ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 408b8840d6bSPeter Brune switch (qn->type) { 409b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 410b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 411b8840d6bSPeter Brune break; 412b8840d6bSPeter Brune case SNES_QN_BROYDEN: 413b8840d6bSPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 414b8840d6bSPeter Brune break; 415b8840d6bSPeter Brune case SNES_QN_LBFGS: 416b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 417b8840d6bSPeter Brune break; 418b8840d6bSPeter Brune } 41970d3b23bSPeter Brune /* line search for lambda */ 42070d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 42188d374b2SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 42215f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 423f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4249f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 4259f3a0142SPeter Brune if (snes->domainerror) { 4269f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 4279f3a0142SPeter Brune PetscFunctionReturn(0); 4289f3a0142SPeter Brune } 429f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 4309f3a0142SPeter Brune if (!lssucceed) { 4319f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4329f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4339f3a0142SPeter Brune break; 4349f3a0142SPeter Brune } 4359f3a0142SPeter Brune } 436f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 4370c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 43805ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4390ecc9e1dSPeter Brune } 4403af51624SPeter Brune 44188d374b2SPeter Brune /* convergence monitoring */ 442b21d5a53SPeter 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); 443b21d5a53SPeter Brune 444b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 445b28a06ddSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 446b28a06ddSPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 447b28a06ddSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 448ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 449ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 450ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 451ddd40ce5SPeter Brune PetscFunctionReturn(0); 452ddd40ce5SPeter Brune } 453ddd40ce5SPeter Brune ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 454b28a06ddSPeter Brune } 455b28a06ddSPeter Brune 456360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 457360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 458360c497dSPeter Brune 459a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 4608409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4614b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 462d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4634b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 464b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 46547f26062SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 466b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 467b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 468b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 469b7281c8aSPeter Brune PetscFunctionReturn(0); 470b7281c8aSPeter Brune } 47188d374b2SPeter Brune } else { 47288d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 47388d374b2SPeter Brune } 47488d374b2SPeter Brune 4750c777b0cSPeter Brune powell = PETSC_FALSE; 4760c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 4776bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 4788e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4798e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 4808e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4818e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 4820c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 4830c777b0cSPeter Brune } 4840c777b0cSPeter Brune periodic = PETSC_FALSE; 485b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 486b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 4870c777b0cSPeter Brune } 4880c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 4890c777b0cSPeter Brune if (powell || periodic) { 4905ba6227bSPeter Brune if (qn->monitor) { 4915ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 492516fe3c3SPeter 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); 4935ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 4945ba6227bSPeter Brune } 4955ba6227bSPeter Brune i_r = -1; 4965ba6227bSPeter Brune /* general purpose update */ 4975ba6227bSPeter Brune if (snes->ops->update) { 4985ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 4995ba6227bSPeter Brune } 5000c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5010ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 5028cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 5030ecc9e1dSPeter Brune } 5040ecc9e1dSPeter Brune } 50570d3b23bSPeter Brune /* general purpose update */ 50670d3b23bSPeter Brune if (snes->ops->update) { 50770d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 50870d3b23bSPeter Brune } 5095ba6227bSPeter Brune } 5104b11644fSPeter Brune if (i == snes->max_its) { 5114b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 5124b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 5134b11644fSPeter Brune } 5144b11644fSPeter Brune PetscFunctionReturn(0); 5154b11644fSPeter Brune } 5164b11644fSPeter Brune 5174b11644fSPeter Brune #undef __FUNCT__ 5184b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5194b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5204b11644fSPeter Brune { 5219f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5224b11644fSPeter Brune PetscErrorCode ierr; 523335fb975SPeter Brune 5244b11644fSPeter Brune PetscFunctionBegin; 525b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 526b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 5278e231d97SPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 5288e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 5298e231d97SPeter Brune qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 5308e231d97SPeter Brune 53118aa0c0cSPeter Brune if (qn->singlereduction) { 5328e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 5338e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 5348e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 5358e231d97SPeter Brune } 536fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 537335fb975SPeter Brune 538335fb975SPeter Brune /* set up the line search */ 5390c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5408e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5418e231d97SPeter Brune } 5426c67d002SPeter Brune 5436c67d002SPeter Brune if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;} 5446c67d002SPeter Brune 5454b11644fSPeter Brune PetscFunctionReturn(0); 5464b11644fSPeter Brune } 5474b11644fSPeter Brune 5484b11644fSPeter Brune #undef __FUNCT__ 5494b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5504b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5514b11644fSPeter Brune { 5524b11644fSPeter Brune PetscErrorCode ierr; 5539f83bee8SJed Brown SNES_QN *qn; 5540adebc6cSBarry Smith 5554b11644fSPeter Brune PetscFunctionBegin; 5564b11644fSPeter Brune if (snes->data) { 5579f83bee8SJed Brown qn = (SNES_QN*)snes->data; 558b8840d6bSPeter Brune if (qn->U) { 559b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5604b11644fSPeter Brune } 561b8840d6bSPeter Brune if (qn->V) { 562b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5634b11644fSPeter Brune } 56418aa0c0cSPeter Brune if (qn->singlereduction) { 5658e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5668e231d97SPeter Brune } 5678e231d97SPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 5684b11644fSPeter Brune } 5694b11644fSPeter Brune PetscFunctionReturn(0); 5704b11644fSPeter Brune } 5714b11644fSPeter Brune 5724b11644fSPeter Brune #undef __FUNCT__ 5734b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5744b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5754b11644fSPeter Brune { 5764b11644fSPeter Brune PetscErrorCode ierr; 5776e111a19SKarl Rupp 5784b11644fSPeter Brune PetscFunctionBegin; 5794b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 5804b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 581bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr); 5824b11644fSPeter Brune PetscFunctionReturn(0); 5834b11644fSPeter Brune } 5844b11644fSPeter Brune 5854b11644fSPeter Brune #undef __FUNCT__ 5864b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 5874b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 5884b11644fSPeter Brune { 5894b11644fSPeter Brune 5904b11644fSPeter Brune PetscErrorCode ierr; 5912150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 59288f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 593f1c6b773SPeter Brune SNESLineSearch linesearch; 5942150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 5952150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 5962150357eSBarry Smith 5974b11644fSPeter Brune PetscFunctionBegin; 5984b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 5990298fd71SBarry Smith ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr); 6000298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr); 6010298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr); 6020298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr); 6030298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr); 60488f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 60588f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 60688f769c5SPeter Brune 60788f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 60888f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 60988f769c5SPeter Brune 610b8840d6bSPeter Brune ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 6110298fd71SBarry Smith (PetscEnum)qn->type,(PetscEnum*)&qn->type,NULL);CHKERRQ(ierr); 6124b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6139e764e56SPeter Brune if (!snes->linesearch) { 6147601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 6151a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 6169e764e56SPeter Brune } 61744f7e39eSPeter Brune if (monflg) { 618ce94432eSBarry Smith qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 61944f7e39eSPeter Brune } 6204b11644fSPeter Brune PetscFunctionReturn(0); 6214b11644fSPeter Brune } 6224b11644fSPeter Brune 62392c02d66SPeter Brune 6240c777b0cSPeter Brune #undef __FUNCT__ 6250c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6260c777b0cSPeter Brune /*@ 6270c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6280c777b0cSPeter Brune 6290c777b0cSPeter Brune Logically Collective on SNES 6300c777b0cSPeter Brune 6310c777b0cSPeter Brune Input Parameters: 6320c777b0cSPeter Brune + snes - the iterative context 6330c777b0cSPeter Brune - rtype - restart type 6340c777b0cSPeter Brune 6350c777b0cSPeter Brune Options Database: 6360c777b0cSPeter Brune + -snes_qn_restart_type<powell,periodic,none> - set the restart type 637b8840d6bSPeter Brune - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 6380c777b0cSPeter Brune 6390c777b0cSPeter Brune Level: intermediate 6400c777b0cSPeter Brune 6410c777b0cSPeter Brune SNESQNRestartTypes: 6420c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6430c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6440c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6450c777b0cSPeter Brune 6460c777b0cSPeter Brune Notes: 6470c777b0cSPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 6480c777b0cSPeter Brune 6490c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 6500c777b0cSPeter Brune @*/ 6512150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6522150357eSBarry Smith { 6530c777b0cSPeter Brune PetscErrorCode ierr; 6546e111a19SKarl Rupp 6550c777b0cSPeter Brune PetscFunctionBegin; 6560c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6570c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6580c777b0cSPeter Brune PetscFunctionReturn(0); 6590c777b0cSPeter Brune } 6600c777b0cSPeter Brune 6610c777b0cSPeter Brune #undef __FUNCT__ 6620c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6630c777b0cSPeter Brune /*@ 6640c777b0cSPeter Brune SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 6650c777b0cSPeter Brune 6660c777b0cSPeter Brune Logically Collective on SNES 6670c777b0cSPeter Brune 6680c777b0cSPeter Brune Input Parameters: 6690c777b0cSPeter Brune + snes - the iterative context 6700c777b0cSPeter Brune - stype - scale type 6710c777b0cSPeter Brune 6720c777b0cSPeter Brune Options Database: 6730c777b0cSPeter Brune . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 6740c777b0cSPeter Brune 6750c777b0cSPeter Brune Level: intermediate 6760c777b0cSPeter Brune 6770c777b0cSPeter Brune SNESQNSelectTypes: 6780c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 6790c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 6800c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 6810c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 6820c777b0cSPeter Brune 6830c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 6840c777b0cSPeter Brune @*/ 6850c777b0cSPeter Brune 6862150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 6872150357eSBarry Smith { 6880c777b0cSPeter Brune PetscErrorCode ierr; 6896e111a19SKarl Rupp 6900c777b0cSPeter Brune PetscFunctionBegin; 6910c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6920c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 6930c777b0cSPeter Brune PetscFunctionReturn(0); 6940c777b0cSPeter Brune } 6950c777b0cSPeter Brune 6960c777b0cSPeter Brune #undef __FUNCT__ 6970c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 6980adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 6990adebc6cSBarry Smith { 7000c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7016e111a19SKarl Rupp 7020c777b0cSPeter Brune PetscFunctionBegin; 7030c777b0cSPeter Brune qn->scale_type = stype; 7040c777b0cSPeter Brune PetscFunctionReturn(0); 7050c777b0cSPeter Brune } 7060c777b0cSPeter Brune 7070c777b0cSPeter Brune #undef __FUNCT__ 7080c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7090adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 7100adebc6cSBarry Smith { 7110c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7126e111a19SKarl Rupp 7130c777b0cSPeter Brune PetscFunctionBegin; 7140c777b0cSPeter Brune qn->restart_type = rtype; 7150c777b0cSPeter Brune PetscFunctionReturn(0); 7160c777b0cSPeter Brune } 7170c777b0cSPeter Brune 7184b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 7194b11644fSPeter Brune /*MC 7204b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 7214b11644fSPeter Brune 7226cc8130cSPeter Brune Options Database: 7236cc8130cSPeter Brune 7246cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 7251867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 7261867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 72772835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 72844f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 7296cc8130cSPeter Brune 730b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 731b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 732b8840d6bSPeter Brune updates. 7336cc8130cSPeter Brune 7341867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 7351867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 7361867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 7371867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 7381867fe5bSPeter Brune 7396cc8130cSPeter Brune References: 7406cc8130cSPeter Brune 7416cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 7426cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 7436cc8130cSPeter Brune 744b8840d6bSPeter Brune Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker 745b8840d6bSPeter Brune SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 746b8840d6bSPeter Brune 7474b11644fSPeter Brune 7484b11644fSPeter Brune Level: beginner 7494b11644fSPeter Brune 75004d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 7516cc8130cSPeter Brune 7524b11644fSPeter Brune M*/ 7534b11644fSPeter Brune #undef __FUNCT__ 7544b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 7558cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 7564b11644fSPeter Brune { 7574b11644fSPeter Brune PetscErrorCode ierr; 7589f83bee8SJed Brown SNES_QN *qn; 7594b11644fSPeter Brune 7604b11644fSPeter Brune PetscFunctionBegin; 7614b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 7624b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 7634b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 7644b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 7654b11644fSPeter Brune snes->ops->view = 0; 7664b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 7674b11644fSPeter Brune 76847f26062SPeter Brune snes->pcside = PC_LEFT; 76947f26062SPeter Brune 77042f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 77142f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 77242f4f86dSBarry Smith 77388976e71SPeter Brune if (!snes->tolerancesset) { 7740e444f03SPeter Brune snes->max_funcs = 30000; 7750e444f03SPeter Brune snes->max_its = 10000; 77688976e71SPeter Brune } 7770e444f03SPeter Brune 7789f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 7794b11644fSPeter Brune snes->data = (void*) qn; 7800ecc9e1dSPeter Brune qn->m = 10; 781b21d5a53SPeter Brune qn->scaling = 1.0; 7820298fd71SBarry Smith qn->U = NULL; 7830298fd71SBarry Smith qn->V = NULL; 7840298fd71SBarry Smith qn->dXtdF = NULL; 7850298fd71SBarry Smith qn->dFtdX = NULL; 7860298fd71SBarry Smith qn->dXdFmat = NULL; 7870298fd71SBarry Smith qn->monitor = NULL; 7881c6523dcSPeter Brune qn->singlereduction = PETSC_TRUE; 789b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 7900c777b0cSPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 7910c777b0cSPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 792b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 793ea630c6eSPeter Brune 794bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr); 795bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr); 7964b11644fSPeter Brune PetscFunctionReturn(0); 7974b11644fSPeter Brune } 798