10c777b0cSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 24b11644fSPeter Brune 38e231d97SPeter Brune #define H(i,j) qn->dXdFmat[i*qn->m + j] 48e231d97SPeter Brune 56a6fc655SJed Brown const char *const SNESQNScaleTypes[] = {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0}; 66a6fc655SJed Brown const char *const SNESQNRestartTypes[] = {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0}; 7b8840d6bSPeter Brune const char *const SNESQNTypes[] = {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_RESTART_",0}; 8b8840d6bSPeter Brune 9b8840d6bSPeter Brune typedef enum {SNES_QN_LBFGS = 0, 10b8840d6bSPeter Brune SNES_QN_BROYDEN = 1, 11b8840d6bSPeter Brune SNES_QN_BADBROYDEN = 2} SNESQNType; 126bf1b2e5SPeter Brune 134b11644fSPeter Brune typedef struct { 14b8840d6bSPeter Brune Vec *U; /* Stored past states (vary from method to method) */ 15b8840d6bSPeter Brune Vec *V; /* Stored past states (vary from method to method) */ 168e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 178e231d97SPeter Brune PetscScalar *alpha, *beta; 188e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 1918aa0c0cSPeter Brune PetscBool singlereduction; /* Aggregated reduction implementation */ 208e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 2144f7e39eSPeter Brune PetscViewer monitor; 226bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 236bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 24b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 2588d374b2SPeter Brune 26b8840d6bSPeter Brune SNESQNType type; /* the type of quasi-newton method used */ 2788f769c5SPeter Brune SNESQNScaleType scale_type; /* the type of scaling used */ 280c777b0cSPeter Brune SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 299f83bee8SJed Brown } SNES_QN; 304b11644fSPeter Brune 314b11644fSPeter Brune #undef __FUNCT__ 32b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden" 332150357eSBarry Smith PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold, Vec D, Vec Dold) 342150357eSBarry Smith { 354b11644fSPeter Brune PetscErrorCode ierr; 369f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 37b8840d6bSPeter Brune Vec W = snes->work[3]; 38b8840d6bSPeter Brune Vec *U = qn->U; 39b8840d6bSPeter Brune Vec *V = qn->V; 40b8840d6bSPeter Brune KSPConvergedReason kspreason; 41b8840d6bSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 42b8840d6bSPeter Brune PetscInt k,i,lits; 43b8840d6bSPeter Brune PetscInt m = qn->m; 44b8840d6bSPeter Brune PetscScalar gdot; 45b8840d6bSPeter Brune PetscInt l = m; 46b8840d6bSPeter Brune Mat jac,jac_pre; 472150357eSBarry Smith 48b8840d6bSPeter Brune PetscFunctionBegin; 49b8840d6bSPeter Brune if (it < m) l = it; 50b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 51b8840d6bSPeter Brune ierr = SNESGetJacobian(snes,&jac,&jac_pre,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 52b8840d6bSPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 53b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,D,W);CHKERRQ(ierr); 54b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 55b8840d6bSPeter Brune if (kspreason < 0) { 56b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 57b8840d6bSPeter Brune ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 58b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 59b8840d6bSPeter Brune PetscFunctionReturn(0); 60b8840d6bSPeter Brune } 61b8840d6bSPeter Brune } 62b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 63b8840d6bSPeter Brune snes->linear_its += lits; 64b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 65b8840d6bSPeter Brune } else { 66b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 67b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 68b8840d6bSPeter Brune } 69b8840d6bSPeter Brune 70b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 71b8840d6bSPeter Brune if (it > 1) { 72b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 73b8840d6bSPeter Brune k = (it+i-l)%l; 74b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 75b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 76b8840d6bSPeter Brune if (qn->monitor) { 77b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 78b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 79b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 80b8840d6bSPeter Brune } 81b8840d6bSPeter Brune } 82b8840d6bSPeter Brune } 83b8840d6bSPeter Brune if (it < m) l = it; 84b8840d6bSPeter Brune 85b8840d6bSPeter Brune /* set W to be the last step, post-linesearch */ 86b8840d6bSPeter Brune ierr = VecCopy(Xold,W);CHKERRQ(ierr); 87b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 88b8840d6bSPeter Brune 89b8840d6bSPeter Brune if (l > 0) { 90b8840d6bSPeter Brune k = (it-1)%l; 91b8840d6bSPeter Brune ierr = VecCopy(W,U[k]);CHKERRQ(ierr); 92b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,Y);CHKERRQ(ierr); 93b8840d6bSPeter Brune ierr = VecDot(U[k],W,&gdot);CHKERRQ(ierr); 94b8840d6bSPeter Brune ierr = VecCopy(Y,V[k]);CHKERRQ(ierr); 95b8840d6bSPeter Brune ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr); 96b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 97b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 98b8840d6bSPeter Brune if (qn->monitor) { 99b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 100b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 101b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 102b8840d6bSPeter Brune } 103b8840d6bSPeter Brune } 104b8840d6bSPeter Brune PetscFunctionReturn(0); 105b8840d6bSPeter Brune } 106b8840d6bSPeter Brune 107b8840d6bSPeter Brune #undef __FUNCT__ 108b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden" 109b8840d6bSPeter Brune PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) { 110b8840d6bSPeter Brune 111b8840d6bSPeter Brune PetscErrorCode ierr; 112b8840d6bSPeter Brune 113b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 114b8840d6bSPeter Brune 115b8840d6bSPeter Brune Vec W = snes->work[3]; 116b8840d6bSPeter Brune 117b8840d6bSPeter Brune Vec *U = qn->U; 118b8840d6bSPeter Brune Vec *T = qn->V; 119b8840d6bSPeter Brune 120b8840d6bSPeter Brune /* ksp thing for jacobian scaling */ 121b8840d6bSPeter Brune KSPConvergedReason kspreason; 122b8840d6bSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 123b8840d6bSPeter Brune 124b8840d6bSPeter Brune PetscInt k, i, lits; 125b8840d6bSPeter Brune PetscInt m = qn->m; 126b8840d6bSPeter Brune PetscScalar gdot; 127b8840d6bSPeter Brune PetscInt l = m; 128b8840d6bSPeter Brune 129b8840d6bSPeter Brune Mat jac, jac_pre; 130b8840d6bSPeter Brune PetscFunctionBegin; 131b8840d6bSPeter Brune 132b8840d6bSPeter Brune if (it < m) l = it; 133b8840d6bSPeter Brune 134b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 135b8840d6bSPeter Brune 136b8840d6bSPeter Brune if (l > 0) { 137b8840d6bSPeter Brune k = (it-1)%l; 138b8840d6bSPeter Brune ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 139b8840d6bSPeter Brune ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 140b8840d6bSPeter Brune ierr = VecDot(U[k],U[k],&gdot);CHKERRQ(ierr); 141b8840d6bSPeter Brune ierr = VecCopy(D,T[k]);CHKERRQ(ierr); 142b8840d6bSPeter Brune ierr = VecScale(T[k],1.0/gdot);CHKERRQ(ierr); 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, "update: %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 /* inward recursion starting at the first update and working forward */ 153b8840d6bSPeter Brune if (it > 1) { 154b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 155b8840d6bSPeter Brune k = (it+i-l)%l; 156b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 157b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 158b8840d6bSPeter Brune if (qn->monitor) { 159b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 160b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 161b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 162b8840d6bSPeter Brune } 163b8840d6bSPeter Brune } 164b8840d6bSPeter Brune } 165b8840d6bSPeter Brune 166b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 167b8840d6bSPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 168b8840d6bSPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 169b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr); 170b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 171b8840d6bSPeter Brune if (kspreason < 0) { 172b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 173b8840d6bSPeter 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); 174b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 175b8840d6bSPeter Brune PetscFunctionReturn(0); 176b8840d6bSPeter Brune } 177b8840d6bSPeter Brune } 178b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 179b8840d6bSPeter Brune snes->linear_its += lits; 180b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 181b8840d6bSPeter Brune } else { 182b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 183b8840d6bSPeter Brune } 184b8840d6bSPeter Brune PetscFunctionReturn(0); 185b8840d6bSPeter Brune } 186b8840d6bSPeter Brune 187b8840d6bSPeter Brune #undef __FUNCT__ 188b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 189b8840d6bSPeter Brune PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) { 190b8840d6bSPeter Brune 191b8840d6bSPeter Brune PetscErrorCode ierr; 192b8840d6bSPeter Brune 193b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 194b8840d6bSPeter Brune 195b8840d6bSPeter Brune Vec W = snes->work[3]; 196b8840d6bSPeter Brune 197b8840d6bSPeter Brune Vec *dX = qn->U; 198b8840d6bSPeter Brune Vec *dF = qn->V; 1994b11644fSPeter Brune 200bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 201bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 2028e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 203b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 2048e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 205bd052dfeSPeter Brune 2060ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 2070ecc9e1dSPeter Brune KSPConvergedReason kspreason; 2080ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 2090ecc9e1dSPeter Brune 2108e231d97SPeter Brune PetscInt k,i,j,g,lits; 2114b11644fSPeter Brune PetscInt m = qn->m; 2124b11644fSPeter Brune PetscScalar t; 2134b11644fSPeter Brune PetscInt l = m; 2144b11644fSPeter Brune 2158e231d97SPeter Brune Mat jac,jac_pre; 2168e231d97SPeter Brune 2174b11644fSPeter Brune PetscFunctionBegin; 2184b11644fSPeter Brune 2195ba6227bSPeter Brune if (it < m) l = it; 2204b11644fSPeter Brune 221b8840d6bSPeter Brune if (it > 0) { 222b8840d6bSPeter Brune k = (it - 1) % l; 223b8840d6bSPeter Brune ierr = VecCopy(D, dF[k]);CHKERRQ(ierr); 224b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 225b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 226b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 22718aa0c0cSPeter Brune if (qn->singlereduction) { 228b8840d6bSPeter Brune ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr); 229b8840d6bSPeter Brune ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr); 230b8840d6bSPeter Brune for (j = 0; j < l; j++) { 231b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 232b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 233b8840d6bSPeter Brune } 234b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 235b8840d6bSPeter Brune for (j = 0; j < l; j++) { 236b8840d6bSPeter Brune dXtdF[j] = H(j, j); 237b8840d6bSPeter Brune } 238b8840d6bSPeter Brune } else { 239b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 240b8840d6bSPeter Brune } 241b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 242b8840d6bSPeter Brune PetscScalar dFtdF; 243b8840d6bSPeter Brune ierr = VecDot(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 244b8840d6bSPeter Brune qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(dFtdF); 245b8840d6bSPeter Brune } else if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 246b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 247b8840d6bSPeter Brune } 248b8840d6bSPeter Brune } 249b8840d6bSPeter Brune 250b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 251b8840d6bSPeter Brune 252b8840d6bSPeter Brune if (qn->singlereduction) { 253b8840d6bSPeter Brune ierr = VecMDot(Y,l,dX,YtdX);CHKERRQ(ierr); 2548e231d97SPeter Brune } 2554b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2564b11644fSPeter Brune for (i=0;i<l;i++) { 257b21d5a53SPeter Brune k = (it-i-1)%l; 25818aa0c0cSPeter Brune if (qn->singlereduction) { 2598e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2608e231d97SPeter Brune t = YtdX[k]; 2618e231d97SPeter Brune for (j=0;j<i;j++) { 2628e231d97SPeter Brune g = (it-j-1)%l; 2638e231d97SPeter Brune t += -alpha[g]*H(g, k); 2648e231d97SPeter Brune } 2658e231d97SPeter Brune alpha[k] = t/H(k,k); 2668e231d97SPeter Brune } else { 2673af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2688e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2698e231d97SPeter Brune } 27044f7e39eSPeter Brune if (qn->monitor) { 2713af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2725ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2733af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 27444f7e39eSPeter Brune } 2756bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2764b11644fSPeter Brune } 2774b11644fSPeter Brune 2780c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 2798e231d97SPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 2808e231d97SPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 281b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr); 2820ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2830ecc9e1dSPeter Brune if (kspreason < 0) { 2840ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 2850ecc9e1dSPeter 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); 2860ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 2870ecc9e1dSPeter Brune PetscFunctionReturn(0); 2880ecc9e1dSPeter Brune } 2890ecc9e1dSPeter Brune } 2900ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2910ecc9e1dSPeter Brune snes->linear_its += lits; 292b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2930ecc9e1dSPeter Brune } else { 294b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2950ecc9e1dSPeter Brune } 29618aa0c0cSPeter Brune if (qn->singlereduction) { 297b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 2988e231d97SPeter Brune } 2994b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 300bd052dfeSPeter Brune for (i = 0; i < l; i++) { 301b21d5a53SPeter Brune k = (it + i - l) % l; 30218aa0c0cSPeter Brune if (qn->singlereduction) { 3038e231d97SPeter Brune t = YtdX[k]; 3048e231d97SPeter Brune for (j = 0; j < i; j++) { 3058e231d97SPeter Brune g = (it + j - l) % l; 3068e231d97SPeter Brune t += (alpha[g] - beta[g])*H(k, g); 3078e231d97SPeter Brune } 3088e231d97SPeter Brune beta[k] = t / H(k, k); 3098e231d97SPeter Brune } else { 3106bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 3118e231d97SPeter Brune beta[k] = t / dXtdF[k]; 3128e231d97SPeter Brune } 3133af51624SPeter Brune ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]); 31444f7e39eSPeter Brune if (qn->monitor) { 3153af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3165ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 3173af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 31844f7e39eSPeter Brune } 3194b11644fSPeter Brune } 3204b11644fSPeter Brune PetscFunctionReturn(0); 3214b11644fSPeter Brune } 3224b11644fSPeter Brune 3234b11644fSPeter Brune #undef __FUNCT__ 3244b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 3254b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 3264b11644fSPeter Brune { 3274b11644fSPeter Brune 3284b11644fSPeter Brune PetscErrorCode ierr; 3299f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 3304b11644fSPeter Brune 33115f5eeeaSPeter Brune Vec X,Xold; 332335fb975SPeter Brune Vec F,B; 333335fb975SPeter Brune Vec Y,FPC,D,Dold; 3343af51624SPeter Brune SNESConvergedReason reason; 335b8840d6bSPeter Brune PetscInt i, i_r; 3364b11644fSPeter Brune 337335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 3380c777b0cSPeter Brune PetscBool lssucceed,powell,periodic; 3394b11644fSPeter Brune 340b8840d6bSPeter Brune PetscScalar DolddotD,DolddotDold,DdotD,YdotD; 3414b11644fSPeter Brune 3420ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 3430ecc9e1dSPeter Brune 3444b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 3454b11644fSPeter Brune PetscFunctionBegin; 3464b11644fSPeter Brune 3479f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3483af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 3493af51624SPeter Brune B = snes->vec_rhs; 350b8840d6bSPeter Brune 351b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 352335fb975SPeter Brune Xold = snes->work[0]; 3533af51624SPeter Brune 3543af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 355335fb975SPeter Brune D = snes->work[1]; 356335fb975SPeter Brune Dold = snes->work[2]; 3574b11644fSPeter Brune 3584b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3594b11644fSPeter Brune 3604b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3614b11644fSPeter Brune snes->iter = 0; 3624b11644fSPeter Brune snes->norm = 0.; 3634b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 364e4ed7901SPeter Brune if (!snes->vec_func_init_set){ 36515f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3664b11644fSPeter Brune if (snes->domainerror) { 3674b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3684b11644fSPeter Brune PetscFunctionReturn(0); 3694b11644fSPeter Brune } 370e4ed7901SPeter Brune } else { 371e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 372e4ed7901SPeter Brune } 373e4ed7901SPeter Brune 374e4ed7901SPeter Brune if (!snes->norm_init_set) { 37515f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 3764b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 377e4ed7901SPeter Brune } else { 378e4ed7901SPeter Brune fnorm = snes->norm_init; 379e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 380e4ed7901SPeter Brune } 381e4ed7901SPeter Brune 3824b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3834b11644fSPeter Brune snes->norm = fnorm; 3844b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 3854b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 3864b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3874b11644fSPeter Brune 3884b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 3894b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 3904b11644fSPeter Brune /* test convergence */ 3914b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3924b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 39370d3b23bSPeter Brune 39488f769c5SPeter Brune /* composed solve */ 395c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 396217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 397217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 39888d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 39988d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 40088d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 40188d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 40288d374b2SPeter Brune PetscFunctionReturn(0); 40388d374b2SPeter Brune } 40488d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 40588d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 40688d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 4076bf1b2e5SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 40888d374b2SPeter Brune } else { 40988d374b2SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 41088d374b2SPeter Brune } 41188d374b2SPeter Brune ierr = VecCopy(Y, D);CHKERRQ(ierr); 4123af51624SPeter Brune 413f8e15203SPeter Brune /* scale the initial update */ 4140c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4150ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4160ecc9e1dSPeter Brune } 4170ecc9e1dSPeter Brune 4185ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 419b8840d6bSPeter Brune switch(qn->type) { 420b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 421b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 422b8840d6bSPeter Brune break; 423b8840d6bSPeter Brune case SNES_QN_BROYDEN: 424b8840d6bSPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 425b8840d6bSPeter Brune break; 426b8840d6bSPeter Brune case SNES_QN_LBFGS: 427b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 428b8840d6bSPeter Brune break; 429b8840d6bSPeter Brune } 43070d3b23bSPeter Brune /* line search for lambda */ 43170d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 43288d374b2SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 43315f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 434f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4359f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 4369f3a0142SPeter Brune if (snes->domainerror) { 4379f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 4389f3a0142SPeter Brune PetscFunctionReturn(0); 4399f3a0142SPeter Brune } 440f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 4419f3a0142SPeter Brune if (!lssucceed) { 4429f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4439f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4449f3a0142SPeter Brune break; 4459f3a0142SPeter Brune } 4469f3a0142SPeter Brune } 447f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 4480c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 44905ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4500ecc9e1dSPeter Brune } 4513af51624SPeter Brune 45288d374b2SPeter Brune /* convergence monitoring */ 453b21d5a53SPeter 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); 454b21d5a53SPeter Brune 455360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 456360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 457360c497dSPeter Brune 4588409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 4598409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4604b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 461d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4624b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 4634b11644fSPeter Brune 464c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 465217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 466217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 46788d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 46888d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 46988d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 47088d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 47188d374b2SPeter Brune PetscFunctionReturn(0); 47288d374b2SPeter Brune } 47388d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 47488d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 47588d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 47688d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 47788d374b2SPeter Brune } else { 47888d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 47988d374b2SPeter Brune } 48088d374b2SPeter Brune 4810c777b0cSPeter Brune powell = PETSC_FALSE; 4820c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 4836bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 4848e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4858e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 4868e231d97SPeter Brune ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr); 4878e231d97SPeter Brune ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr); 4888e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4898e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 4908e231d97SPeter Brune ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr); 4918e231d97SPeter Brune ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr); 4920c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 4930c777b0cSPeter Brune } 4940c777b0cSPeter Brune periodic = PETSC_FALSE; 495b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 496b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 4970c777b0cSPeter Brune } 4980c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 4990c777b0cSPeter Brune if (powell || periodic) { 5005ba6227bSPeter Brune if (qn->monitor) { 5015ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 502516fe3c3SPeter 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); 5035ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 5045ba6227bSPeter Brune } 5055ba6227bSPeter Brune i_r = -1; 5065ba6227bSPeter Brune /* general purpose update */ 5075ba6227bSPeter Brune if (snes->ops->update) { 5085ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 5095ba6227bSPeter Brune } 5100c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5110ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 5120ecc9e1dSPeter Brune } 5130ecc9e1dSPeter Brune } 51470d3b23bSPeter Brune /* general purpose update */ 51570d3b23bSPeter Brune if (snes->ops->update) { 51670d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 51770d3b23bSPeter Brune } 5185ba6227bSPeter Brune } 5194b11644fSPeter Brune if (i == snes->max_its) { 5204b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 5214b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 5224b11644fSPeter Brune } 5234b11644fSPeter Brune PetscFunctionReturn(0); 5244b11644fSPeter Brune } 5254b11644fSPeter Brune 5264b11644fSPeter Brune 5274b11644fSPeter Brune #undef __FUNCT__ 5284b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5294b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5304b11644fSPeter Brune { 5319f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5324b11644fSPeter Brune PetscErrorCode ierr; 533335fb975SPeter Brune 5344b11644fSPeter Brune PetscFunctionBegin; 535b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 536b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 5378e231d97SPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 5388e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 5398e231d97SPeter Brune qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 5408e231d97SPeter Brune 54118aa0c0cSPeter Brune if (qn->singlereduction) { 5428e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 5438e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 5448e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 5458e231d97SPeter Brune } 546335fb975SPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 547335fb975SPeter Brune 548335fb975SPeter Brune /* set up the line search */ 5490c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5508e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5518e231d97SPeter Brune } 5524b11644fSPeter Brune PetscFunctionReturn(0); 5534b11644fSPeter Brune } 5544b11644fSPeter Brune 5554b11644fSPeter Brune #undef __FUNCT__ 5564b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5574b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5584b11644fSPeter Brune { 5594b11644fSPeter Brune PetscErrorCode ierr; 5609f83bee8SJed Brown SNES_QN *qn; 5614b11644fSPeter Brune PetscFunctionBegin; 5624b11644fSPeter Brune if (snes->data) { 5639f83bee8SJed Brown qn = (SNES_QN*)snes->data; 564b8840d6bSPeter Brune if (qn->U) { 565b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5664b11644fSPeter Brune } 567b8840d6bSPeter Brune if (qn->V) { 568b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5694b11644fSPeter Brune } 57018aa0c0cSPeter Brune if (qn->singlereduction) { 5718e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5728e231d97SPeter Brune } 5738e231d97SPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 5744b11644fSPeter Brune } 5754b11644fSPeter Brune PetscFunctionReturn(0); 5764b11644fSPeter Brune } 5774b11644fSPeter Brune 5784b11644fSPeter Brune #undef __FUNCT__ 5794b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5804b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5814b11644fSPeter Brune { 5824b11644fSPeter Brune PetscErrorCode ierr; 5834b11644fSPeter Brune PetscFunctionBegin; 5844b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 5854b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 5869c05af5eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 5874b11644fSPeter Brune PetscFunctionReturn(0); 5884b11644fSPeter Brune } 5894b11644fSPeter Brune 5904b11644fSPeter Brune #undef __FUNCT__ 5914b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 5924b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 5934b11644fSPeter Brune { 5944b11644fSPeter Brune 5954b11644fSPeter Brune PetscErrorCode ierr; 5962150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 59788f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 598f1c6b773SPeter Brune SNESLineSearch linesearch; 5992150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 6002150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 6012150357eSBarry Smith 6024b11644fSPeter Brune PetscFunctionBegin; 6034b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 6045b4627e8SPeter Brune ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,PETSC_NULL);CHKERRQ(ierr); 6055b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr); 6065b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr); 6075b4627e8SPeter Brune ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 6084d116c7dSPeter Brune ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, PETSC_NULL);CHKERRQ(ierr); 60988f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 61088f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 61188f769c5SPeter Brune 61288f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 61388f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 61488f769c5SPeter Brune 615b8840d6bSPeter Brune ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 616b8840d6bSPeter Brune (PetscEnum)qn->type,(PetscEnum*)&qn->type,PETSC_NULL);CHKERRQ(ierr); 6174b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6189e764e56SPeter Brune if (!snes->linesearch) { 619f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 6201a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 6219e764e56SPeter Brune } 62244f7e39eSPeter Brune if (monflg) { 62344f7e39eSPeter Brune qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 62444f7e39eSPeter Brune } 6254b11644fSPeter Brune PetscFunctionReturn(0); 6264b11644fSPeter Brune } 6274b11644fSPeter Brune 62892c02d66SPeter Brune 6290c777b0cSPeter Brune #undef __FUNCT__ 6300c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6310c777b0cSPeter Brune /*@ 6320c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6330c777b0cSPeter Brune 6340c777b0cSPeter Brune Logically Collective on SNES 6350c777b0cSPeter Brune 6360c777b0cSPeter Brune Input Parameters: 6370c777b0cSPeter Brune + snes - the iterative context 6380c777b0cSPeter Brune - rtype - restart type 6390c777b0cSPeter Brune 6400c777b0cSPeter Brune Options Database: 6410c777b0cSPeter Brune + -snes_qn_restart_type<powell,periodic,none> - set the restart type 642b8840d6bSPeter Brune - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 6430c777b0cSPeter Brune 6440c777b0cSPeter Brune Level: intermediate 6450c777b0cSPeter Brune 6460c777b0cSPeter Brune SNESQNRestartTypes: 6470c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6480c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6490c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6500c777b0cSPeter Brune 6510c777b0cSPeter Brune Notes: 6520c777b0cSPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 6530c777b0cSPeter Brune 6540c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 6550c777b0cSPeter Brune @*/ 6562150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6572150357eSBarry Smith { 6580c777b0cSPeter Brune PetscErrorCode ierr; 6590c777b0cSPeter Brune PetscFunctionBegin; 6600c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6610c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6620c777b0cSPeter Brune PetscFunctionReturn(0); 6630c777b0cSPeter Brune } 6640c777b0cSPeter Brune 6650c777b0cSPeter Brune #undef __FUNCT__ 6660c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6670c777b0cSPeter Brune /*@ 6680c777b0cSPeter Brune SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 6690c777b0cSPeter Brune 6700c777b0cSPeter Brune Logically Collective on SNES 6710c777b0cSPeter Brune 6720c777b0cSPeter Brune Input Parameters: 6730c777b0cSPeter Brune + snes - the iterative context 6740c777b0cSPeter Brune - stype - scale type 6750c777b0cSPeter Brune 6760c777b0cSPeter Brune Options Database: 6770c777b0cSPeter Brune . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 6780c777b0cSPeter Brune 6790c777b0cSPeter Brune Level: intermediate 6800c777b0cSPeter Brune 6810c777b0cSPeter Brune SNESQNSelectTypes: 6820c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 6830c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 6840c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 6850c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 6860c777b0cSPeter Brune 6870c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 6880c777b0cSPeter Brune @*/ 6890c777b0cSPeter Brune 6902150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 6912150357eSBarry Smith { 6920c777b0cSPeter Brune PetscErrorCode ierr; 6930c777b0cSPeter Brune PetscFunctionBegin; 6940c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6950c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 6960c777b0cSPeter Brune PetscFunctionReturn(0); 6970c777b0cSPeter Brune } 6980c777b0cSPeter Brune 6990c777b0cSPeter Brune EXTERN_C_BEGIN 7000c777b0cSPeter Brune #undef __FUNCT__ 7010c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 7020c777b0cSPeter Brune PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) { 7030c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 7040c777b0cSPeter Brune PetscFunctionBegin; 7050c777b0cSPeter Brune qn->scale_type = stype; 7060c777b0cSPeter Brune PetscFunctionReturn(0); 7070c777b0cSPeter Brune } 7080c777b0cSPeter Brune 7090c777b0cSPeter Brune #undef __FUNCT__ 7100c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7110c777b0cSPeter Brune PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) { 7120c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 7130c777b0cSPeter Brune PetscFunctionBegin; 7140c777b0cSPeter Brune qn->restart_type = rtype; 7150c777b0cSPeter Brune PetscFunctionReturn(0); 7160c777b0cSPeter Brune } 7170c777b0cSPeter Brune EXTERN_C_END 7180c777b0cSPeter Brune 7194b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 7204b11644fSPeter Brune /*MC 7214b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 7224b11644fSPeter Brune 7236cc8130cSPeter Brune Options Database: 7246cc8130cSPeter Brune 7256cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 7261867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 7271867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 72872835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 72944f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 7306cc8130cSPeter Brune 731b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 732b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 733b8840d6bSPeter Brune updates. 7346cc8130cSPeter Brune 7351867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 7361867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 7371867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 7381867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 7391867fe5bSPeter Brune 7406cc8130cSPeter Brune References: 7416cc8130cSPeter Brune 7426cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 7436cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 7446cc8130cSPeter Brune 745b8840d6bSPeter Brune Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker 746b8840d6bSPeter Brune SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 747b8840d6bSPeter Brune 7484b11644fSPeter Brune 7494b11644fSPeter Brune Level: beginner 7504b11644fSPeter Brune 751*04d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 7526cc8130cSPeter Brune 7534b11644fSPeter Brune M*/ 7544b11644fSPeter Brune EXTERN_C_BEGIN 7554b11644fSPeter Brune #undef __FUNCT__ 7564b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 7574b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 7584b11644fSPeter Brune { 7594b11644fSPeter Brune 7604b11644fSPeter Brune PetscErrorCode ierr; 7619f83bee8SJed Brown SNES_QN *qn; 7624b11644fSPeter Brune 7634b11644fSPeter Brune PetscFunctionBegin; 7644b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 7654b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 7664b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 7674b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 7684b11644fSPeter Brune snes->ops->view = 0; 7694b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 7704b11644fSPeter Brune 77142f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 77242f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 77342f4f86dSBarry Smith 77488976e71SPeter Brune if (!snes->tolerancesset) { 7750e444f03SPeter Brune snes->max_funcs = 30000; 7760e444f03SPeter Brune snes->max_its = 10000; 77788976e71SPeter Brune } 7780e444f03SPeter Brune 7799f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 7804b11644fSPeter Brune snes->data = (void *) qn; 7810ecc9e1dSPeter Brune qn->m = 10; 782b21d5a53SPeter Brune qn->scaling = 1.0; 783b8840d6bSPeter Brune qn->U = PETSC_NULL; 784b8840d6bSPeter Brune qn->V = PETSC_NULL; 7858e231d97SPeter Brune qn->dXtdF = PETSC_NULL; 7868e231d97SPeter Brune qn->dFtdX = PETSC_NULL; 7878e231d97SPeter Brune qn->dXdFmat = PETSC_NULL; 78844f7e39eSPeter Brune qn->monitor = PETSC_NULL; 78918aa0c0cSPeter Brune qn->singlereduction = PETSC_FALSE; 790b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 7910c777b0cSPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 7920c777b0cSPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 793b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 794ea630c6eSPeter Brune 79588f769c5SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetScaleType_C","SNESQNSetScaleType_QN",SNESQNSetScaleType_QN);CHKERRQ(ierr); 79688f769c5SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESQNSetRestartType_C","SNESQNSetRestartType_QN",SNESQNSetRestartType_QN);CHKERRQ(ierr); 79788f769c5SPeter Brune 7984b11644fSPeter Brune PetscFunctionReturn(0); 7994b11644fSPeter Brune } 8008e231d97SPeter Brune 8014b11644fSPeter Brune EXTERN_C_END 802