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) { 49b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,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) { 154b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,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; 196b8840d6bSPeter Brune if (it > 0) { 197b8840d6bSPeter Brune k = (it - 1) % l; 198b8840d6bSPeter Brune ierr = VecCopy(D, dF[k]);CHKERRQ(ierr); 199b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 200b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 201b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 20218aa0c0cSPeter Brune if (qn->singlereduction) { 203b8840d6bSPeter Brune ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr); 204b8840d6bSPeter Brune ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr); 205b8840d6bSPeter Brune for (j = 0; j < l; j++) { 206b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 207b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 208b8840d6bSPeter Brune } 209b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 2101aa26658SKarl Rupp for (j = 0; j < l; j++) dXtdF[j] = H(j, j); 211b8840d6bSPeter Brune } else { 212b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 213b8840d6bSPeter Brune } 214b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 21567392de3SBarry Smith PetscReal dFtdF; 21667392de3SBarry Smith ierr = VecDotRealPart(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 21767392de3SBarry Smith qn->scaling = PetscRealPart(dXtdF[k])/dFtdF; 218b8840d6bSPeter Brune } else if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 219b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 220b8840d6bSPeter Brune } 221b8840d6bSPeter Brune } 222b8840d6bSPeter Brune 223b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 224b8840d6bSPeter Brune 225b8840d6bSPeter Brune if (qn->singlereduction) { 226b8840d6bSPeter Brune ierr = VecMDot(Y,l,dX,YtdX);CHKERRQ(ierr); 2278e231d97SPeter Brune } 2284b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2294b11644fSPeter Brune for (i=0; i<l; i++) { 230b21d5a53SPeter Brune k = (it-i-1)%l; 23118aa0c0cSPeter Brune if (qn->singlereduction) { 2328e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2338e231d97SPeter Brune t = YtdX[k]; 2348e231d97SPeter Brune for (j=0; j<i; j++) { 2358e231d97SPeter Brune g = (it-j-1)%l; 2368e231d97SPeter Brune t += -alpha[g]*H(g, k); 2378e231d97SPeter Brune } 2388e231d97SPeter Brune alpha[k] = t/H(k,k); 2398e231d97SPeter Brune } else { 2403af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2418e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2428e231d97SPeter Brune } 24344f7e39eSPeter Brune if (qn->monitor) { 2443af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2455ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2463af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 24744f7e39eSPeter Brune } 2486bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2494b11644fSPeter Brune } 2504b11644fSPeter Brune 2510c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 252b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr); 2530ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2540ecc9e1dSPeter Brune if (kspreason < 0) { 2550ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 2560ecc9e1dSPeter 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); 2570ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 2580ecc9e1dSPeter Brune PetscFunctionReturn(0); 2590ecc9e1dSPeter Brune } 2600ecc9e1dSPeter Brune } 2610ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2620ecc9e1dSPeter Brune snes->linear_its += lits; 263b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2640ecc9e1dSPeter Brune } else { 265b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2660ecc9e1dSPeter Brune } 26718aa0c0cSPeter Brune if (qn->singlereduction) { 268b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 2698e231d97SPeter Brune } 2704b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 271bd052dfeSPeter Brune for (i = 0; i < l; i++) { 272b21d5a53SPeter Brune k = (it + i - l) % l; 27318aa0c0cSPeter Brune if (qn->singlereduction) { 2748e231d97SPeter Brune t = YtdX[k]; 2758e231d97SPeter Brune for (j = 0; j < i; j++) { 2768e231d97SPeter Brune g = (it + j - l) % l; 2778e231d97SPeter Brune t += (alpha[g] - beta[g])*H(k, g); 2788e231d97SPeter Brune } 2798e231d97SPeter Brune beta[k] = t / H(k, k); 2808e231d97SPeter Brune } else { 2816bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 2828e231d97SPeter Brune beta[k] = t / dXtdF[k]; 2838e231d97SPeter Brune } 28422d28d08SBarry Smith ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 28544f7e39eSPeter Brune if (qn->monitor) { 2863af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2875ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 2883af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 28944f7e39eSPeter Brune } 2904b11644fSPeter Brune } 2914b11644fSPeter Brune PetscFunctionReturn(0); 2924b11644fSPeter Brune } 2934b11644fSPeter Brune 2944b11644fSPeter Brune #undef __FUNCT__ 2954b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 2964b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 2974b11644fSPeter Brune { 2984b11644fSPeter Brune PetscErrorCode ierr; 2999f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 30015f5eeeaSPeter Brune Vec X,Xold; 301335fb975SPeter Brune Vec F,B; 302335fb975SPeter Brune Vec Y,FPC,D,Dold; 3033af51624SPeter Brune SNESConvergedReason reason; 304b8840d6bSPeter Brune PetscInt i, i_r; 305335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 3060c777b0cSPeter Brune PetscBool lssucceed,powell,periodic; 307b8840d6bSPeter Brune PetscScalar DolddotD,DolddotDold,DdotD,YdotD; 3080ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 3090ecc9e1dSPeter Brune 3104b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 3114b11644fSPeter Brune 3126e111a19SKarl Rupp PetscFunctionBegin; 3139f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3143af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 3153af51624SPeter Brune B = snes->vec_rhs; 316b8840d6bSPeter Brune 317b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 318335fb975SPeter Brune Xold = snes->work[0]; 3193af51624SPeter Brune 3203af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 321335fb975SPeter Brune D = snes->work[1]; 322335fb975SPeter Brune Dold = snes->work[2]; 3234b11644fSPeter Brune 3244b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3254b11644fSPeter Brune 3264b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3274b11644fSPeter Brune snes->iter = 0; 3284b11644fSPeter Brune snes->norm = 0.; 3294b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 330e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 33115f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3324b11644fSPeter Brune if (snes->domainerror) { 3334b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3344b11644fSPeter Brune PetscFunctionReturn(0); 3354b11644fSPeter Brune } 3361aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 337e4ed7901SPeter Brune 338e4ed7901SPeter Brune if (!snes->norm_init_set) { 33915f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 3404b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 341e4ed7901SPeter Brune } else { 342e4ed7901SPeter Brune fnorm = snes->norm_init; 343e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 344e4ed7901SPeter Brune } 345e4ed7901SPeter Brune 3464b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3474b11644fSPeter Brune snes->norm = fnorm; 3484b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 349a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 3504b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3514b11644fSPeter Brune 3524b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 3534b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 3544b11644fSPeter Brune /* test convergence */ 3554b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3564b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 35770d3b23bSPeter Brune 35888f769c5SPeter Brune /* composed solve */ 359c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 360217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 361217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 36263e7833aSPeter Brune 36363e7833aSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,B,0);CHKERRQ(ierr); 36488d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 36563e7833aSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,B,0);CHKERRQ(ierr); 36663e7833aSPeter Brune 36788d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 36888d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 36988d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 37088d374b2SPeter Brune PetscFunctionReturn(0); 37188d374b2SPeter Brune } 3720298fd71SBarry Smith ierr = SNESGetFunction(snes->pc, &FPC, NULL, NULL);CHKERRQ(ierr); 37388d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 37488d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 3756bf1b2e5SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 37688d374b2SPeter Brune } else { 37788d374b2SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 37888d374b2SPeter Brune } 37988d374b2SPeter Brune ierr = VecCopy(Y, D);CHKERRQ(ierr); 3803af51624SPeter Brune 381f8e15203SPeter Brune /* scale the initial update */ 3820c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 3830ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 3848cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 3850ecc9e1dSPeter Brune } 3860ecc9e1dSPeter Brune 3875ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 388b8840d6bSPeter Brune switch (qn->type) { 389b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 390b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 391b8840d6bSPeter Brune break; 392b8840d6bSPeter Brune case SNES_QN_BROYDEN: 393b8840d6bSPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 394b8840d6bSPeter Brune break; 395b8840d6bSPeter Brune case SNES_QN_LBFGS: 396b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 397b8840d6bSPeter Brune break; 398b8840d6bSPeter Brune } 39970d3b23bSPeter Brune /* line search for lambda */ 40070d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 40188d374b2SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 40215f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 403f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4049f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 4059f3a0142SPeter Brune if (snes->domainerror) { 4069f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 4079f3a0142SPeter Brune PetscFunctionReturn(0); 4089f3a0142SPeter Brune } 409f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 4109f3a0142SPeter Brune if (!lssucceed) { 4119f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4129f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4139f3a0142SPeter Brune break; 4149f3a0142SPeter Brune } 4159f3a0142SPeter Brune } 416f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 4170c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 41805ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4190ecc9e1dSPeter Brune } 4203af51624SPeter Brune 42188d374b2SPeter Brune /* convergence monitoring */ 422b21d5a53SPeter 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); 423b21d5a53SPeter Brune 424360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 425360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 426360c497dSPeter Brune 427a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 4288409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4294b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 430d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4314b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 4324b11644fSPeter Brune 433c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 434217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 435217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 43688d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 43788d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 43888d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 43988d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 44088d374b2SPeter Brune PetscFunctionReturn(0); 44188d374b2SPeter Brune } 4420298fd71SBarry Smith ierr = SNESGetFunction(snes->pc, &FPC, NULL, NULL);CHKERRQ(ierr); 44388d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 44488d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 44588d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 44688d374b2SPeter Brune } else { 44788d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 44888d374b2SPeter Brune } 44988d374b2SPeter Brune 4500c777b0cSPeter Brune powell = PETSC_FALSE; 4510c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 4526bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 4538e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4548e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 4558e231d97SPeter Brune ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr); 4568e231d97SPeter Brune ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr); 4578e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4588e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 4598e231d97SPeter Brune ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr); 4608e231d97SPeter Brune ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr); 4610c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 4620c777b0cSPeter Brune } 4630c777b0cSPeter Brune periodic = PETSC_FALSE; 464b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 465b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 4660c777b0cSPeter Brune } 4670c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 4680c777b0cSPeter Brune if (powell || periodic) { 4695ba6227bSPeter Brune if (qn->monitor) { 4705ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 471516fe3c3SPeter 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); 4725ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 4735ba6227bSPeter Brune } 4745ba6227bSPeter Brune i_r = -1; 4755ba6227bSPeter Brune /* general purpose update */ 4765ba6227bSPeter Brune if (snes->ops->update) { 4775ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 4785ba6227bSPeter Brune } 4790c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4800ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4818cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 4820ecc9e1dSPeter Brune } 4830ecc9e1dSPeter Brune } 48470d3b23bSPeter Brune /* general purpose update */ 48570d3b23bSPeter Brune if (snes->ops->update) { 48670d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 48770d3b23bSPeter Brune } 4885ba6227bSPeter Brune } 4894b11644fSPeter Brune if (i == snes->max_its) { 4904b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 4914b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 4924b11644fSPeter Brune } 4934b11644fSPeter Brune PetscFunctionReturn(0); 4944b11644fSPeter Brune } 4954b11644fSPeter Brune 4964b11644fSPeter Brune #undef __FUNCT__ 4974b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 4984b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 4994b11644fSPeter Brune { 5009f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5014b11644fSPeter Brune PetscErrorCode ierr; 502335fb975SPeter Brune 5034b11644fSPeter Brune PetscFunctionBegin; 504b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 505b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 5068e231d97SPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 5078e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 5088e231d97SPeter Brune qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 5098e231d97SPeter Brune 51018aa0c0cSPeter Brune if (qn->singlereduction) { 5118e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 5128e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 5138e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 5148e231d97SPeter Brune } 515335fb975SPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 516335fb975SPeter Brune 517335fb975SPeter Brune /* set up the line search */ 5180c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5198e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5208e231d97SPeter Brune } 5214b11644fSPeter Brune PetscFunctionReturn(0); 5224b11644fSPeter Brune } 5234b11644fSPeter Brune 5244b11644fSPeter Brune #undef __FUNCT__ 5254b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5264b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5274b11644fSPeter Brune { 5284b11644fSPeter Brune PetscErrorCode ierr; 5299f83bee8SJed Brown SNES_QN *qn; 5300adebc6cSBarry Smith 5314b11644fSPeter Brune PetscFunctionBegin; 5324b11644fSPeter Brune if (snes->data) { 5339f83bee8SJed Brown qn = (SNES_QN*)snes->data; 534b8840d6bSPeter Brune if (qn->U) { 535b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5364b11644fSPeter Brune } 537b8840d6bSPeter Brune if (qn->V) { 538b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5394b11644fSPeter Brune } 54018aa0c0cSPeter Brune if (qn->singlereduction) { 5418e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5428e231d97SPeter Brune } 5438e231d97SPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 5444b11644fSPeter Brune } 5454b11644fSPeter Brune PetscFunctionReturn(0); 5464b11644fSPeter Brune } 5474b11644fSPeter Brune 5484b11644fSPeter Brune #undef __FUNCT__ 5494b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5504b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5514b11644fSPeter Brune { 5524b11644fSPeter Brune PetscErrorCode ierr; 5536e111a19SKarl Rupp 5544b11644fSPeter Brune PetscFunctionBegin; 5554b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 5564b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 5570298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",NULL);CHKERRQ(ierr); 5584b11644fSPeter Brune PetscFunctionReturn(0); 5594b11644fSPeter Brune } 5604b11644fSPeter Brune 5614b11644fSPeter Brune #undef __FUNCT__ 5624b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 5634b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 5644b11644fSPeter Brune { 5654b11644fSPeter Brune 5664b11644fSPeter Brune PetscErrorCode ierr; 5672150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 56888f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 569f1c6b773SPeter Brune SNESLineSearch linesearch; 5702150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 5712150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 5722150357eSBarry Smith 5734b11644fSPeter Brune PetscFunctionBegin; 5744b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 5750298fd71SBarry Smith ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr); 5760298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr); 5770298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr); 5780298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr); 5790298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr); 58088f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 58188f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 58288f769c5SPeter Brune 58388f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 58488f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 58588f769c5SPeter Brune 586b8840d6bSPeter Brune ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 5870298fd71SBarry Smith (PetscEnum)qn->type,(PetscEnum*)&qn->type,NULL);CHKERRQ(ierr); 5884b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 5899e764e56SPeter Brune if (!snes->linesearch) { 590f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 5911a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 5929e764e56SPeter Brune } 59344f7e39eSPeter Brune if (monflg) { 594*ce94432eSBarry Smith qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 59544f7e39eSPeter Brune } 5964b11644fSPeter Brune PetscFunctionReturn(0); 5974b11644fSPeter Brune } 5984b11644fSPeter Brune 59992c02d66SPeter Brune 6000c777b0cSPeter Brune #undef __FUNCT__ 6010c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6020c777b0cSPeter Brune /*@ 6030c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6040c777b0cSPeter Brune 6050c777b0cSPeter Brune Logically Collective on SNES 6060c777b0cSPeter Brune 6070c777b0cSPeter Brune Input Parameters: 6080c777b0cSPeter Brune + snes - the iterative context 6090c777b0cSPeter Brune - rtype - restart type 6100c777b0cSPeter Brune 6110c777b0cSPeter Brune Options Database: 6120c777b0cSPeter Brune + -snes_qn_restart_type<powell,periodic,none> - set the restart type 613b8840d6bSPeter Brune - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 6140c777b0cSPeter Brune 6150c777b0cSPeter Brune Level: intermediate 6160c777b0cSPeter Brune 6170c777b0cSPeter Brune SNESQNRestartTypes: 6180c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6190c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6200c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6210c777b0cSPeter Brune 6220c777b0cSPeter Brune Notes: 6230c777b0cSPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 6240c777b0cSPeter Brune 6250c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 6260c777b0cSPeter Brune @*/ 6272150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6282150357eSBarry Smith { 6290c777b0cSPeter Brune PetscErrorCode ierr; 6306e111a19SKarl Rupp 6310c777b0cSPeter Brune PetscFunctionBegin; 6320c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6330c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6340c777b0cSPeter Brune PetscFunctionReturn(0); 6350c777b0cSPeter Brune } 6360c777b0cSPeter Brune 6370c777b0cSPeter Brune #undef __FUNCT__ 6380c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6390c777b0cSPeter Brune /*@ 6400c777b0cSPeter Brune SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 6410c777b0cSPeter Brune 6420c777b0cSPeter Brune Logically Collective on SNES 6430c777b0cSPeter Brune 6440c777b0cSPeter Brune Input Parameters: 6450c777b0cSPeter Brune + snes - the iterative context 6460c777b0cSPeter Brune - stype - scale type 6470c777b0cSPeter Brune 6480c777b0cSPeter Brune Options Database: 6490c777b0cSPeter Brune . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 6500c777b0cSPeter Brune 6510c777b0cSPeter Brune Level: intermediate 6520c777b0cSPeter Brune 6530c777b0cSPeter Brune SNESQNSelectTypes: 6540c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 6550c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 6560c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 6570c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 6580c777b0cSPeter Brune 6590c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 6600c777b0cSPeter Brune @*/ 6610c777b0cSPeter Brune 6622150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 6632150357eSBarry Smith { 6640c777b0cSPeter Brune PetscErrorCode ierr; 6656e111a19SKarl Rupp 6660c777b0cSPeter Brune PetscFunctionBegin; 6670c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6680c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 6690c777b0cSPeter Brune PetscFunctionReturn(0); 6700c777b0cSPeter Brune } 6710c777b0cSPeter Brune 6720c777b0cSPeter Brune EXTERN_C_BEGIN 6730c777b0cSPeter Brune #undef __FUNCT__ 6740c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 6750adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 6760adebc6cSBarry Smith { 6770c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 6786e111a19SKarl Rupp 6790c777b0cSPeter Brune PetscFunctionBegin; 6800c777b0cSPeter Brune qn->scale_type = stype; 6810c777b0cSPeter Brune PetscFunctionReturn(0); 6820c777b0cSPeter Brune } 6830c777b0cSPeter Brune 6840c777b0cSPeter Brune #undef __FUNCT__ 6850c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 6860adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 6870adebc6cSBarry Smith { 6880c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 6896e111a19SKarl Rupp 6900c777b0cSPeter Brune PetscFunctionBegin; 6910c777b0cSPeter Brune qn->restart_type = rtype; 6920c777b0cSPeter Brune PetscFunctionReturn(0); 6930c777b0cSPeter Brune } 6940c777b0cSPeter Brune EXTERN_C_END 6950c777b0cSPeter Brune 6964b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 6974b11644fSPeter Brune /*MC 6984b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 6994b11644fSPeter Brune 7006cc8130cSPeter Brune Options Database: 7016cc8130cSPeter Brune 7026cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 7031867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 7041867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 70572835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 70644f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 7076cc8130cSPeter Brune 708b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 709b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 710b8840d6bSPeter Brune updates. 7116cc8130cSPeter Brune 7121867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 7131867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 7141867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 7151867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 7161867fe5bSPeter Brune 7176cc8130cSPeter Brune References: 7186cc8130cSPeter Brune 7196cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 7206cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 7216cc8130cSPeter Brune 722b8840d6bSPeter Brune Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker 723b8840d6bSPeter Brune SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 724b8840d6bSPeter Brune 7254b11644fSPeter Brune 7264b11644fSPeter Brune Level: beginner 7274b11644fSPeter Brune 72804d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 7296cc8130cSPeter Brune 7304b11644fSPeter Brune M*/ 7314b11644fSPeter Brune EXTERN_C_BEGIN 7324b11644fSPeter Brune #undef __FUNCT__ 7334b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 7344b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 7354b11644fSPeter Brune { 7364b11644fSPeter Brune PetscErrorCode ierr; 7379f83bee8SJed Brown SNES_QN *qn; 7384b11644fSPeter Brune 7394b11644fSPeter Brune PetscFunctionBegin; 7404b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 7414b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 7424b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 7434b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 7444b11644fSPeter Brune snes->ops->view = 0; 7454b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 7464b11644fSPeter Brune 74742f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 74842f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 74942f4f86dSBarry Smith 75088976e71SPeter Brune if (!snes->tolerancesset) { 7510e444f03SPeter Brune snes->max_funcs = 30000; 7520e444f03SPeter Brune snes->max_its = 10000; 75388976e71SPeter Brune } 7540e444f03SPeter Brune 7559f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 7564b11644fSPeter Brune snes->data = (void*) qn; 7570ecc9e1dSPeter Brune qn->m = 10; 758b21d5a53SPeter Brune qn->scaling = 1.0; 7590298fd71SBarry Smith qn->U = NULL; 7600298fd71SBarry Smith qn->V = NULL; 7610298fd71SBarry Smith qn->dXtdF = NULL; 7620298fd71SBarry Smith qn->dFtdX = NULL; 7630298fd71SBarry Smith qn->dXdFmat = NULL; 7640298fd71SBarry Smith qn->monitor = NULL; 76518aa0c0cSPeter Brune qn->singlereduction = PETSC_FALSE; 766b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 7670c777b0cSPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 7680c777b0cSPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 769b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 770ea630c6eSPeter Brune 77188f769c5SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetScaleType_C","SNESQNSetScaleType_QN",SNESQNSetScaleType_QN);CHKERRQ(ierr); 77288f769c5SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetRestartType_C","SNESQNSetRestartType_QN",SNESQNSetRestartType_QN);CHKERRQ(ierr); 7734b11644fSPeter Brune PetscFunctionReturn(0); 7744b11644fSPeter Brune } 7758e231d97SPeter Brune 7764b11644fSPeter Brune EXTERN_C_END 777