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 SNESQNCompositionTypes[] = {"SEQUENTIAL","COMPOSED","SNESQNCompositionType","SNES_QN_",0}; 66a6fc655SJed Brown const char *const SNESQNScaleTypes[] = {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0}; 76a6fc655SJed Brown const char *const SNESQNRestartTypes[] = {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0}; 8*b8840d6bSPeter Brune const char *const SNESQNTypes[] = {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_RESTART_",0}; 9*b8840d6bSPeter Brune 10*b8840d6bSPeter Brune typedef enum {SNES_QN_LBFGS = 0, 11*b8840d6bSPeter Brune SNES_QN_BROYDEN = 1, 12*b8840d6bSPeter Brune SNES_QN_BADBROYDEN = 2} SNESQNType; 136bf1b2e5SPeter Brune 144b11644fSPeter Brune typedef struct { 15*b8840d6bSPeter Brune Vec *U; /* Stored past states (vary from method to method) */ 16*b8840d6bSPeter Brune Vec *V; /* Stored past states (vary from method to method) */ 178e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 188e231d97SPeter Brune PetscScalar *alpha, *beta; 198e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 2018aa0c0cSPeter Brune PetscBool singlereduction; /* Aggregated reduction implementation */ 218e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 2244f7e39eSPeter Brune PetscViewer monitor; 236bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 246bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 25b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 2688d374b2SPeter Brune 27*b8840d6bSPeter Brune SNESQNType type; /* the type of quasi-newton method used */ 280c777b0cSPeter Brune SNESQNCompositionType composition_type; /* determine if the composition is done sequentially or as a composition */ 290c777b0cSPeter Brune SNESQNScaleType scale_type; /* determine if the composition is done sequentially or as a composition */ 300c777b0cSPeter Brune SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 319f83bee8SJed Brown } SNES_QN; 324b11644fSPeter Brune 334b11644fSPeter Brune #undef __FUNCT__ 34*b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden" 35*b8840d6bSPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold, Vec D, Vec Dold) { 364b11644fSPeter Brune 374b11644fSPeter Brune PetscErrorCode ierr; 384b11644fSPeter Brune 399f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 404b11644fSPeter Brune 41*b8840d6bSPeter Brune Vec W = snes->work[3]; 420ecc9e1dSPeter Brune 43*b8840d6bSPeter Brune Vec *U = qn->U; 44*b8840d6bSPeter Brune Vec *V = qn->V; 45*b8840d6bSPeter Brune 46*b8840d6bSPeter Brune /* ksp thing for jacobian scaling */ 47*b8840d6bSPeter Brune KSPConvergedReason kspreason; 48*b8840d6bSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 49*b8840d6bSPeter Brune 50*b8840d6bSPeter Brune PetscInt k,i,lits; 51*b8840d6bSPeter Brune PetscInt m = qn->m; 52*b8840d6bSPeter Brune PetscScalar gdot; 53*b8840d6bSPeter Brune PetscInt l = m; 54*b8840d6bSPeter Brune 55*b8840d6bSPeter Brune Mat jac,jac_pre; 56*b8840d6bSPeter Brune PetscFunctionBegin; 57*b8840d6bSPeter Brune 58*b8840d6bSPeter Brune if (it < m) l = it; 59*b8840d6bSPeter Brune 60*b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 61*b8840d6bSPeter Brune ierr = SNESGetJacobian(snes,&jac,&jac_pre,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 62*b8840d6bSPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 63*b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,D,W);CHKERRQ(ierr); 64*b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 65*b8840d6bSPeter Brune if (kspreason < 0) { 66*b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 67*b8840d6bSPeter 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); 68*b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 69*b8840d6bSPeter Brune PetscFunctionReturn(0); 70*b8840d6bSPeter Brune } 71*b8840d6bSPeter Brune } 72*b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 73*b8840d6bSPeter Brune snes->linear_its += lits; 74*b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 75*b8840d6bSPeter Brune } else { 76*b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 77*b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 78*b8840d6bSPeter Brune } 79*b8840d6bSPeter Brune 80*b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 81*b8840d6bSPeter Brune if (it > 1) { 82*b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 83*b8840d6bSPeter Brune k = (it+i-l)%l; 84*b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 85*b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 86*b8840d6bSPeter Brune if (qn->monitor) { 87*b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 88*b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 89*b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 90*b8840d6bSPeter Brune } 91*b8840d6bSPeter Brune } 92*b8840d6bSPeter Brune } 93*b8840d6bSPeter Brune if (it < m) l = it; 94*b8840d6bSPeter Brune 95*b8840d6bSPeter Brune /* set W to be the last step, post-linesearch */ 96*b8840d6bSPeter Brune ierr = VecCopy(Xold,W);CHKERRQ(ierr); 97*b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 98*b8840d6bSPeter Brune 99*b8840d6bSPeter Brune if (l > 0) { 100*b8840d6bSPeter Brune k = (it-1)%l; 101*b8840d6bSPeter Brune ierr = VecCopy(W,U[k]);CHKERRQ(ierr); 102*b8840d6bSPeter Brune ierr = VecAXPY(W,-1.0,Y);CHKERRQ(ierr); 103*b8840d6bSPeter Brune ierr = VecDot(U[k],W,&gdot);CHKERRQ(ierr); 104*b8840d6bSPeter Brune ierr = VecCopy(Y,V[k]);CHKERRQ(ierr); 105*b8840d6bSPeter Brune ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr); 106*b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 107*b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 108*b8840d6bSPeter Brune if (qn->monitor) { 109*b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 110*b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 111*b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 112*b8840d6bSPeter Brune } 113*b8840d6bSPeter Brune } 114*b8840d6bSPeter Brune PetscFunctionReturn(0); 115*b8840d6bSPeter Brune } 116*b8840d6bSPeter Brune 117*b8840d6bSPeter Brune #undef __FUNCT__ 118*b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden" 119*b8840d6bSPeter Brune PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) { 120*b8840d6bSPeter Brune 121*b8840d6bSPeter Brune PetscErrorCode ierr; 122*b8840d6bSPeter Brune 123*b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 124*b8840d6bSPeter Brune 125*b8840d6bSPeter Brune Vec W = snes->work[3]; 126*b8840d6bSPeter Brune 127*b8840d6bSPeter Brune Vec *U = qn->U; 128*b8840d6bSPeter Brune Vec *T = qn->V; 129*b8840d6bSPeter Brune 130*b8840d6bSPeter Brune /* ksp thing for jacobian scaling */ 131*b8840d6bSPeter Brune KSPConvergedReason kspreason; 132*b8840d6bSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 133*b8840d6bSPeter Brune 134*b8840d6bSPeter Brune PetscInt k, i, lits; 135*b8840d6bSPeter Brune PetscInt m = qn->m; 136*b8840d6bSPeter Brune PetscScalar gdot; 137*b8840d6bSPeter Brune PetscInt l = m; 138*b8840d6bSPeter Brune 139*b8840d6bSPeter Brune Mat jac, jac_pre; 140*b8840d6bSPeter Brune PetscFunctionBegin; 141*b8840d6bSPeter Brune 142*b8840d6bSPeter Brune if (it < m) l = it; 143*b8840d6bSPeter Brune 144*b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 145*b8840d6bSPeter Brune 146*b8840d6bSPeter Brune if (l > 0) { 147*b8840d6bSPeter Brune k = (it-1)%l; 148*b8840d6bSPeter Brune ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 149*b8840d6bSPeter Brune ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 150*b8840d6bSPeter Brune ierr = VecDot(U[k],U[k],&gdot);CHKERRQ(ierr); 151*b8840d6bSPeter Brune ierr = VecCopy(D,T[k]);CHKERRQ(ierr); 152*b8840d6bSPeter Brune ierr = VecScale(T[k],1.0/gdot);CHKERRQ(ierr); 153*b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 154*b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 155*b8840d6bSPeter Brune if (qn->monitor) { 156*b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 157*b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 158*b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 159*b8840d6bSPeter Brune } 160*b8840d6bSPeter Brune } 161*b8840d6bSPeter Brune 162*b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 163*b8840d6bSPeter Brune if (it > 1) { 164*b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 165*b8840d6bSPeter Brune k = (it+i-l)%l; 166*b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 167*b8840d6bSPeter Brune ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 168*b8840d6bSPeter Brune if (qn->monitor) { 169*b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 170*b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 171*b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 172*b8840d6bSPeter Brune } 173*b8840d6bSPeter Brune } 174*b8840d6bSPeter Brune } 175*b8840d6bSPeter Brune 176*b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 177*b8840d6bSPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 178*b8840d6bSPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 179*b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr); 180*b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 181*b8840d6bSPeter Brune if (kspreason < 0) { 182*b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 183*b8840d6bSPeter 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); 184*b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 185*b8840d6bSPeter Brune PetscFunctionReturn(0); 186*b8840d6bSPeter Brune } 187*b8840d6bSPeter Brune } 188*b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 189*b8840d6bSPeter Brune snes->linear_its += lits; 190*b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 191*b8840d6bSPeter Brune } else { 192*b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 193*b8840d6bSPeter Brune } 194*b8840d6bSPeter Brune PetscFunctionReturn(0); 195*b8840d6bSPeter Brune } 196*b8840d6bSPeter Brune 197*b8840d6bSPeter Brune #undef __FUNCT__ 198*b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 199*b8840d6bSPeter Brune PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) { 200*b8840d6bSPeter Brune 201*b8840d6bSPeter Brune PetscErrorCode ierr; 202*b8840d6bSPeter Brune 203*b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 204*b8840d6bSPeter Brune 205*b8840d6bSPeter Brune Vec W = snes->work[3]; 206*b8840d6bSPeter Brune 207*b8840d6bSPeter Brune Vec *dX = qn->U; 208*b8840d6bSPeter Brune Vec *dF = qn->V; 2094b11644fSPeter Brune 210bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 211bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 2128e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 213*b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 2148e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 215bd052dfeSPeter Brune 2160ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 2170ecc9e1dSPeter Brune KSPConvergedReason kspreason; 2180ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 2190ecc9e1dSPeter Brune 2208e231d97SPeter Brune PetscInt k,i,j,g,lits; 2214b11644fSPeter Brune PetscInt m = qn->m; 2224b11644fSPeter Brune PetscScalar t; 2234b11644fSPeter Brune PetscInt l = m; 2244b11644fSPeter Brune 2258e231d97SPeter Brune Mat jac,jac_pre; 2268e231d97SPeter Brune 2274b11644fSPeter Brune PetscFunctionBegin; 2284b11644fSPeter Brune 2295ba6227bSPeter Brune if (it < m) l = it; 2304b11644fSPeter Brune 231*b8840d6bSPeter Brune if (it > 0) { 232*b8840d6bSPeter Brune k = (it - 1) % l; 233*b8840d6bSPeter Brune ierr = VecCopy(D, dF[k]);CHKERRQ(ierr); 234*b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 235*b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 236*b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 23718aa0c0cSPeter Brune if (qn->singlereduction) { 238*b8840d6bSPeter Brune ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr); 239*b8840d6bSPeter Brune ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr); 240*b8840d6bSPeter Brune for (j = 0; j < l; j++) { 241*b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 242*b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 243*b8840d6bSPeter Brune } 244*b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 245*b8840d6bSPeter Brune for (j = 0; j < l; j++) { 246*b8840d6bSPeter Brune dXtdF[j] = H(j, j); 247*b8840d6bSPeter Brune } 248*b8840d6bSPeter Brune } else { 249*b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 250*b8840d6bSPeter Brune } 251*b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 252*b8840d6bSPeter Brune PetscScalar dFtdF; 253*b8840d6bSPeter Brune ierr = VecDot(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 254*b8840d6bSPeter Brune qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(dFtdF); 255*b8840d6bSPeter Brune } else if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 256*b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 257*b8840d6bSPeter Brune } 258*b8840d6bSPeter Brune } 259*b8840d6bSPeter Brune 260*b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 261*b8840d6bSPeter Brune 262*b8840d6bSPeter Brune if (qn->singlereduction) { 263*b8840d6bSPeter Brune ierr = VecMDot(Y,l,dX,YtdX);CHKERRQ(ierr); 2648e231d97SPeter Brune } 2654b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2664b11644fSPeter Brune for (i=0;i<l;i++) { 267b21d5a53SPeter Brune k = (it-i-1)%l; 26818aa0c0cSPeter Brune if (qn->singlereduction) { 2698e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2708e231d97SPeter Brune t = YtdX[k]; 2718e231d97SPeter Brune for (j=0;j<i;j++) { 2728e231d97SPeter Brune g = (it-j-1)%l; 2738e231d97SPeter Brune t += -alpha[g]*H(g, k); 2748e231d97SPeter Brune } 2758e231d97SPeter Brune alpha[k] = t/H(k,k); 2768e231d97SPeter Brune } else { 2773af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2788e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2798e231d97SPeter Brune } 28044f7e39eSPeter Brune if (qn->monitor) { 2813af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2825ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2833af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 28444f7e39eSPeter Brune } 2856bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2864b11644fSPeter Brune } 2874b11644fSPeter Brune 2880c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 2898e231d97SPeter Brune ierr = SNESGetJacobian(snes, &jac, &jac_pre, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 2908e231d97SPeter Brune ierr = KSPSetOperators(snes->ksp,jac,jac_pre,flg);CHKERRQ(ierr); 291*b8840d6bSPeter Brune ierr = SNES_KSPSolve(snes,snes->ksp,Y,W);CHKERRQ(ierr); 2920ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2930ecc9e1dSPeter Brune if (kspreason < 0) { 2940ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 2950ecc9e1dSPeter 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); 2960ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 2970ecc9e1dSPeter Brune PetscFunctionReturn(0); 2980ecc9e1dSPeter Brune } 2990ecc9e1dSPeter Brune } 3000ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 3010ecc9e1dSPeter Brune snes->linear_its += lits; 302*b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 3030ecc9e1dSPeter Brune } else { 304b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 3050ecc9e1dSPeter Brune } 30618aa0c0cSPeter Brune if (qn->singlereduction) { 307*b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 3088e231d97SPeter Brune } 3094b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 310bd052dfeSPeter Brune for (i = 0; i < l; i++) { 311b21d5a53SPeter Brune k = (it + i - l) % l; 31218aa0c0cSPeter Brune if (qn->singlereduction) { 3138e231d97SPeter Brune t = YtdX[k]; 3148e231d97SPeter Brune for (j = 0; j < i; j++) { 3158e231d97SPeter Brune g = (it + j - l) % l; 3168e231d97SPeter Brune t += (alpha[g] - beta[g])*H(k, g); 3178e231d97SPeter Brune } 3188e231d97SPeter Brune beta[k] = t / H(k, k); 3198e231d97SPeter Brune } else { 3206bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 3218e231d97SPeter Brune beta[k] = t / dXtdF[k]; 3228e231d97SPeter Brune } 3233af51624SPeter Brune ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]); 32444f7e39eSPeter Brune if (qn->monitor) { 3253af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3265ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 3273af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 32844f7e39eSPeter Brune } 3294b11644fSPeter Brune } 3304b11644fSPeter Brune PetscFunctionReturn(0); 3314b11644fSPeter Brune } 3324b11644fSPeter Brune 3334b11644fSPeter Brune #undef __FUNCT__ 3344b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 3354b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 3364b11644fSPeter Brune { 3374b11644fSPeter Brune 3384b11644fSPeter Brune PetscErrorCode ierr; 3399f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 3404b11644fSPeter Brune 34115f5eeeaSPeter Brune Vec X,Xold; 342335fb975SPeter Brune Vec F,B; 343335fb975SPeter Brune Vec Y,FPC,D,Dold; 3443af51624SPeter Brune SNESConvergedReason reason; 345*b8840d6bSPeter Brune PetscInt i, i_r; 3464b11644fSPeter Brune 347335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 3480c777b0cSPeter Brune PetscBool lssucceed,powell,periodic; 3494b11644fSPeter Brune 350*b8840d6bSPeter Brune PetscScalar DolddotD,DolddotDold,DdotD,YdotD; 3514b11644fSPeter Brune 3520ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 3530ecc9e1dSPeter Brune 3544b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 3554b11644fSPeter Brune PetscFunctionBegin; 3564b11644fSPeter Brune 3579f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3583af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 3593af51624SPeter Brune B = snes->vec_rhs; 360*b8840d6bSPeter Brune 361*b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 362335fb975SPeter Brune Xold = snes->work[0]; 3633af51624SPeter Brune 3643af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 365335fb975SPeter Brune D = snes->work[1]; 366335fb975SPeter Brune Dold = snes->work[2]; 3674b11644fSPeter Brune 3684b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3694b11644fSPeter Brune 3704b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3714b11644fSPeter Brune snes->iter = 0; 3724b11644fSPeter Brune snes->norm = 0.; 3734b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 374e4ed7901SPeter Brune if (!snes->vec_func_init_set){ 37515f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3764b11644fSPeter Brune if (snes->domainerror) { 3774b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3784b11644fSPeter Brune PetscFunctionReturn(0); 3794b11644fSPeter Brune } 380e4ed7901SPeter Brune } else { 381e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 382e4ed7901SPeter Brune } 383e4ed7901SPeter Brune 384e4ed7901SPeter Brune if (!snes->norm_init_set) { 38515f5eeeaSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 3864b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 387e4ed7901SPeter Brune } else { 388e4ed7901SPeter Brune fnorm = snes->norm_init; 389e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 390e4ed7901SPeter Brune } 391e4ed7901SPeter Brune 3924b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3934b11644fSPeter Brune snes->norm = fnorm; 3944b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 3954b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 3964b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3974b11644fSPeter Brune 3984b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 3994b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 4004b11644fSPeter Brune /* test convergence */ 4014b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4024b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 40370d3b23bSPeter Brune 40488d374b2SPeter Brune /* composed solve -- either sequential or composed */ 40588d374b2SPeter Brune if (snes->pc) { 4060c777b0cSPeter Brune if (qn->composition_type == SNES_QN_SEQUENTIAL) { 407217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 408217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 40988d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 41088d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 41188d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 41288d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 41388d374b2SPeter Brune PetscFunctionReturn(0); 41488d374b2SPeter Brune } 41588d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 41688d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 41788d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 4186bf1b2e5SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 41988d374b2SPeter Brune } else { 42088d374b2SPeter Brune ierr = VecCopy(X, Y);CHKERRQ(ierr); 421217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 422217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 42388d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, Y);CHKERRQ(ierr); 42488d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 42588d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 42688d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 42788d374b2SPeter Brune PetscFunctionReturn(0); 42888d374b2SPeter Brune } 42988d374b2SPeter Brune ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 43088d374b2SPeter Brune } 43188d374b2SPeter Brune } else { 43288d374b2SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 43388d374b2SPeter Brune } 43488d374b2SPeter Brune ierr = VecCopy(Y, D);CHKERRQ(ierr); 4353af51624SPeter Brune 436f8e15203SPeter Brune /* scale the initial update */ 4370c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4380ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4390ecc9e1dSPeter Brune } 4400ecc9e1dSPeter Brune 4415ba6227bSPeter Brune for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 442*b8840d6bSPeter Brune switch(qn->type) { 443*b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 444*b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 445*b8840d6bSPeter Brune break; 446*b8840d6bSPeter Brune case SNES_QN_BROYDEN: 447*b8840d6bSPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 448*b8840d6bSPeter Brune break; 449*b8840d6bSPeter Brune case SNES_QN_LBFGS: 450*b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 451*b8840d6bSPeter Brune break; 452*b8840d6bSPeter Brune } 45370d3b23bSPeter Brune /* line search for lambda */ 45470d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 45588d374b2SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 45615f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 457f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4589f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 4599f3a0142SPeter Brune if (snes->domainerror) { 4609f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 4619f3a0142SPeter Brune PetscFunctionReturn(0); 4629f3a0142SPeter Brune } 463f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 4649f3a0142SPeter Brune if (!lssucceed) { 4659f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4669f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4679f3a0142SPeter Brune break; 4689f3a0142SPeter Brune } 4699f3a0142SPeter Brune } 470f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 4710c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 472*b8840d6bSPeter Brune PetscScalar scl; 473*b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &scl);CHKERRQ(ierr); 474*b8840d6bSPeter Brune qn->scaling = PetscRealPart(scl);CHKERRQ(ierr); 4750ecc9e1dSPeter Brune } 4763af51624SPeter Brune 47788d374b2SPeter Brune /* convergence monitoring */ 478b21d5a53SPeter 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); 479b21d5a53SPeter Brune 480360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 481360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 482360c497dSPeter Brune 4838409ca45SMatthew G Knepley SNESLogConvHistory(snes,snes->norm,snes->iter); 4848409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4854b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 486d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4874b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 4884b11644fSPeter Brune 48988d374b2SPeter Brune 49088d374b2SPeter Brune if (snes->pc) { 4910c777b0cSPeter Brune if (qn->composition_type == SNES_QN_SEQUENTIAL) { 492217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 493217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 49488d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 49588d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 49688d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 49788d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 49888d374b2SPeter Brune PetscFunctionReturn(0); 49988d374b2SPeter Brune } 50088d374b2SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 50188d374b2SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 50288d374b2SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 50388d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 50488d374b2SPeter Brune } else { 50588d374b2SPeter Brune ierr = VecCopy(X, D);CHKERRQ(ierr); 506217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 507217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 50888d374b2SPeter Brune ierr = SNESSolve(snes->pc, B, D);CHKERRQ(ierr); 50988d374b2SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 51088d374b2SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 51188d374b2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 51288d374b2SPeter Brune PetscFunctionReturn(0); 51388d374b2SPeter Brune } 51488d374b2SPeter Brune ierr = VecAYPX(D,-1.0,X);CHKERRQ(ierr); 51588d374b2SPeter Brune } 51688d374b2SPeter Brune } else { 51788d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 51888d374b2SPeter Brune } 51988d374b2SPeter Brune 5200c777b0cSPeter Brune powell = PETSC_FALSE; 5210c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 5226bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 5238e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 5248e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 5258e231d97SPeter Brune ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr); 5268e231d97SPeter Brune ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr); 5278e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 5288e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 5298e231d97SPeter Brune ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr); 5308e231d97SPeter Brune ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr); 5310c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 5320c777b0cSPeter Brune } 5330c777b0cSPeter Brune periodic = PETSC_FALSE; 534*b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 535*b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 5360c777b0cSPeter Brune } 5370c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 5380c777b0cSPeter Brune if (powell || periodic) { 5395ba6227bSPeter Brune if (qn->monitor) { 5405ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 541516fe3c3SPeter 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); 5425ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 5435ba6227bSPeter Brune } 5445ba6227bSPeter Brune i_r = -1; 5455ba6227bSPeter Brune /* general purpose update */ 5465ba6227bSPeter Brune if (snes->ops->update) { 5475ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 5485ba6227bSPeter Brune } 5490c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5500ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 5510ecc9e1dSPeter Brune } 5520ecc9e1dSPeter Brune } 55370d3b23bSPeter Brune /* general purpose update */ 55470d3b23bSPeter Brune if (snes->ops->update) { 55570d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 55670d3b23bSPeter Brune } 5575ba6227bSPeter Brune } 5584b11644fSPeter Brune if (i == snes->max_its) { 5594b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 5604b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 5614b11644fSPeter Brune } 5624b11644fSPeter Brune PetscFunctionReturn(0); 5634b11644fSPeter Brune } 5644b11644fSPeter Brune 5654b11644fSPeter Brune 5664b11644fSPeter Brune #undef __FUNCT__ 5674b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5684b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5694b11644fSPeter Brune { 5709f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5714b11644fSPeter Brune PetscErrorCode ierr; 572335fb975SPeter Brune 5734b11644fSPeter Brune PetscFunctionBegin; 574*b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 575*b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 5768e231d97SPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 5778e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 5788e231d97SPeter Brune qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 5798e231d97SPeter Brune 58018aa0c0cSPeter Brune if (qn->singlereduction) { 5818e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 5828e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 5838e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 5848e231d97SPeter Brune } 585335fb975SPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 586335fb975SPeter Brune 587335fb975SPeter Brune /* set up the line search */ 5880c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5898e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5908e231d97SPeter Brune } 5914b11644fSPeter Brune PetscFunctionReturn(0); 5924b11644fSPeter Brune } 5934b11644fSPeter Brune 5944b11644fSPeter Brune #undef __FUNCT__ 5954b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5964b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5974b11644fSPeter Brune { 5984b11644fSPeter Brune PetscErrorCode ierr; 5999f83bee8SJed Brown SNES_QN *qn; 6004b11644fSPeter Brune PetscFunctionBegin; 6014b11644fSPeter Brune if (snes->data) { 6029f83bee8SJed Brown qn = (SNES_QN*)snes->data; 603*b8840d6bSPeter Brune if (qn->U) { 604*b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 6054b11644fSPeter Brune } 606*b8840d6bSPeter Brune if (qn->V) { 607*b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 6084b11644fSPeter Brune } 60918aa0c0cSPeter Brune if (qn->singlereduction) { 6108e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 6118e231d97SPeter Brune } 6128e231d97SPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 6134b11644fSPeter Brune } 6144b11644fSPeter Brune PetscFunctionReturn(0); 6154b11644fSPeter Brune } 6164b11644fSPeter Brune 6174b11644fSPeter Brune #undef __FUNCT__ 6184b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 6194b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 6204b11644fSPeter Brune { 6214b11644fSPeter Brune PetscErrorCode ierr; 6224b11644fSPeter Brune PetscFunctionBegin; 6234b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 6244b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 6259c05af5eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 6264b11644fSPeter Brune PetscFunctionReturn(0); 6274b11644fSPeter Brune } 6284b11644fSPeter Brune 6294b11644fSPeter Brune #undef __FUNCT__ 6304b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 6314b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 6324b11644fSPeter Brune { 6334b11644fSPeter Brune 6344b11644fSPeter Brune PetscErrorCode ierr; 6359f83bee8SJed Brown SNES_QN *qn; 63644f7e39eSPeter Brune PetscBool monflg = PETSC_FALSE; 637f1c6b773SPeter Brune SNESLineSearch linesearch; 6384b11644fSPeter Brune PetscFunctionBegin; 6394b11644fSPeter Brune 6409f83bee8SJed Brown qn = (SNES_QN*)snes->data; 6414b11644fSPeter Brune 6424b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 6435b4627e8SPeter Brune ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,PETSC_NULL);CHKERRQ(ierr); 6445b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr); 6455b4627e8SPeter Brune ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr); 6465b4627e8SPeter Brune ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 6474d116c7dSPeter Brune ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, PETSC_NULL);CHKERRQ(ierr); 6480c777b0cSPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)qn->scale_type,(PetscEnum*)&qn->scale_type,PETSC_NULL);CHKERRQ(ierr); 6490c777b0cSPeter Brune ierr = PetscOptionsEnum("-snes_qn_composition_type","Composition type","SNESQNSetCompositionType",SNESQNCompositionTypes, 6500c777b0cSPeter Brune (PetscEnum)qn->composition_type,(PetscEnum*)&qn->composition_type,PETSC_NULL);CHKERRQ(ierr); 6510c777b0cSPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes, 6520c777b0cSPeter Brune (PetscEnum)qn->restart_type,(PetscEnum*)&qn->restart_type,PETSC_NULL);CHKERRQ(ierr); 653*b8840d6bSPeter Brune ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 654*b8840d6bSPeter Brune (PetscEnum)qn->type,(PetscEnum*)&qn->type,PETSC_NULL);CHKERRQ(ierr); 6554b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6569e764e56SPeter Brune if (!snes->linesearch) { 657f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 6580c777b0cSPeter Brune if (!snes->pc || qn->composition_type == SNES_QN_SEQUENTIAL) { 6591a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 6609e764e56SPeter Brune } else { 6611a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 6629e764e56SPeter Brune } 6639e764e56SPeter Brune } 66444f7e39eSPeter Brune if (monflg) { 66544f7e39eSPeter Brune qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 66644f7e39eSPeter Brune } 6674b11644fSPeter Brune PetscFunctionReturn(0); 6684b11644fSPeter Brune } 6694b11644fSPeter Brune 67092c02d66SPeter Brune 6710c777b0cSPeter Brune #undef __FUNCT__ 6720c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6730c777b0cSPeter Brune /*@ 6740c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6750c777b0cSPeter Brune 6760c777b0cSPeter Brune Logically Collective on SNES 6770c777b0cSPeter Brune 6780c777b0cSPeter Brune Input Parameters: 6790c777b0cSPeter Brune + snes - the iterative context 6800c777b0cSPeter Brune - rtype - restart type 6810c777b0cSPeter Brune 6820c777b0cSPeter Brune Options Database: 6830c777b0cSPeter Brune + -snes_qn_restart_type<powell,periodic,none> - set the restart type 684*b8840d6bSPeter Brune - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 6850c777b0cSPeter Brune 6860c777b0cSPeter Brune Level: intermediate 6870c777b0cSPeter Brune 6880c777b0cSPeter Brune SNESQNRestartTypes: 6890c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6900c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6910c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6920c777b0cSPeter Brune 6930c777b0cSPeter Brune Notes: 6940c777b0cSPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 6950c777b0cSPeter Brune 6960c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 6970c777b0cSPeter Brune @*/ 6980c777b0cSPeter Brune PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) { 6990c777b0cSPeter Brune PetscErrorCode ierr; 7000c777b0cSPeter Brune PetscFunctionBegin; 7010c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7020c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 7030c777b0cSPeter Brune PetscFunctionReturn(0); 7040c777b0cSPeter Brune } 7050c777b0cSPeter Brune 7060c777b0cSPeter Brune #undef __FUNCT__ 7070c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 7080c777b0cSPeter Brune /*@ 7090c777b0cSPeter Brune SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 7100c777b0cSPeter Brune 7110c777b0cSPeter Brune Logically Collective on SNES 7120c777b0cSPeter Brune 7130c777b0cSPeter Brune Input Parameters: 7140c777b0cSPeter Brune + snes - the iterative context 7150c777b0cSPeter Brune - stype - scale type 7160c777b0cSPeter Brune 7170c777b0cSPeter Brune Options Database: 7180c777b0cSPeter Brune . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 7190c777b0cSPeter Brune 7200c777b0cSPeter Brune Level: intermediate 7210c777b0cSPeter Brune 7220c777b0cSPeter Brune SNESQNSelectTypes: 7230c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 7240c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 7250c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 7260c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 7270c777b0cSPeter Brune 7280c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 7290c777b0cSPeter Brune @*/ 7300c777b0cSPeter Brune 7310c777b0cSPeter Brune PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) { 7320c777b0cSPeter Brune PetscErrorCode ierr; 7330c777b0cSPeter Brune PetscFunctionBegin; 7340c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7350c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 7360c777b0cSPeter Brune PetscFunctionReturn(0); 7370c777b0cSPeter Brune } 7380c777b0cSPeter Brune 7390c777b0cSPeter Brune #undef __FUNCT__ 7400c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetCompositionType" 7410c777b0cSPeter Brune /*@ 7420c777b0cSPeter Brune SNESQNSetCompositionType - Sets the composition type 7430c777b0cSPeter Brune 7440c777b0cSPeter Brune Logically Collective on SNES 7450c777b0cSPeter Brune 7460c777b0cSPeter Brune Input Parameters: 7470c777b0cSPeter Brune + snes - the iterative context 7480c777b0cSPeter Brune - stype - composition type 7490c777b0cSPeter Brune 7500c777b0cSPeter Brune Options Database: 7510c777b0cSPeter Brune . -snes_qn_composition_type<sequential, composed> 7520c777b0cSPeter Brune 7530c777b0cSPeter Brune Level: intermediate 7540c777b0cSPeter Brune 7550c777b0cSPeter Brune SNESQNSelectTypes: 7560c777b0cSPeter Brune + SNES_QN_COMPOSITION_SEQUENTIAL - Solve the system with X = PC(X) and D = F(PC(X)) 7570c777b0cSPeter Brune - SNES_QN_COMPOSITION_COMPOSED - solve the system with X = X and D = PC(X) - X 7580c777b0cSPeter Brune 7590c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 7600c777b0cSPeter Brune @*/ 7610c777b0cSPeter Brune 7620c777b0cSPeter Brune PetscErrorCode SNESQNSetCompositionType(SNES snes, SNESQNCompositionType ctype) { 7630c777b0cSPeter Brune PetscErrorCode ierr; 7640c777b0cSPeter Brune PetscFunctionBegin; 7650c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7660c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetCompositionType_C",(SNES,SNESQNCompositionType),(snes,ctype));CHKERRQ(ierr); 7670c777b0cSPeter Brune PetscFunctionReturn(0); 7680c777b0cSPeter Brune } 7690c777b0cSPeter Brune 7700c777b0cSPeter Brune EXTERN_C_BEGIN 7710c777b0cSPeter Brune #undef __FUNCT__ 7720c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 7730c777b0cSPeter Brune PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) { 7740c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 7750c777b0cSPeter Brune PetscFunctionBegin; 7760c777b0cSPeter Brune qn->scale_type = stype; 7770c777b0cSPeter Brune PetscFunctionReturn(0); 7780c777b0cSPeter Brune } 7790c777b0cSPeter Brune 7800c777b0cSPeter Brune #undef __FUNCT__ 7810c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7820c777b0cSPeter Brune PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) { 7830c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 7840c777b0cSPeter Brune PetscFunctionBegin; 7850c777b0cSPeter Brune qn->restart_type = rtype; 7860c777b0cSPeter Brune PetscFunctionReturn(0); 7870c777b0cSPeter Brune } 7880c777b0cSPeter Brune 7890c777b0cSPeter Brune #undef __FUNCT__ 7900c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetCompositionType_QN" 7910c777b0cSPeter Brune 7920c777b0cSPeter Brune PetscErrorCode SNESQNSetCompositionType_QN(SNES snes, SNESQNCompositionType ctype) { 7930c777b0cSPeter Brune SNES_QN *qn = (SNES_QN *)snes->data; 7940c777b0cSPeter Brune PetscFunctionBegin; 7950c777b0cSPeter Brune qn->composition_type = ctype; 7960c777b0cSPeter Brune PetscFunctionReturn(0); 7970c777b0cSPeter Brune } 7980c777b0cSPeter Brune EXTERN_C_END 7990c777b0cSPeter Brune 8004b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 8014b11644fSPeter Brune /*MC 8024b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 8034b11644fSPeter Brune 8046cc8130cSPeter Brune Options Database: 8056cc8130cSPeter Brune 8066cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 8071867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 8081867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 8091867fe5bSPeter Brune . -snes_qn_composition <sequential, composed>- Type of composition. 81072835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 81144f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 8126cc8130cSPeter Brune 813*b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 814*b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 815*b8840d6bSPeter Brune updates. 8166cc8130cSPeter Brune 8171867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 8181867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 8191867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 8201867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 8211867fe5bSPeter Brune 8226cc8130cSPeter Brune References: 8236cc8130cSPeter Brune 8246cc8130cSPeter Brune L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 8256cc8130cSPeter Brune International Journal of Computer Mathematics, vol. 86, 2009. 8266cc8130cSPeter Brune 827*b8840d6bSPeter Brune Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker 828*b8840d6bSPeter Brune SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 829*b8840d6bSPeter Brune 8304b11644fSPeter Brune 8314b11644fSPeter Brune Level: beginner 8324b11644fSPeter Brune 8334b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 8346cc8130cSPeter Brune 8354b11644fSPeter Brune M*/ 8364b11644fSPeter Brune EXTERN_C_BEGIN 8374b11644fSPeter Brune #undef __FUNCT__ 8384b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 8394b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 8404b11644fSPeter Brune { 8414b11644fSPeter Brune 8424b11644fSPeter Brune PetscErrorCode ierr; 8439f83bee8SJed Brown SNES_QN *qn; 8444b11644fSPeter Brune 8454b11644fSPeter Brune PetscFunctionBegin; 8464b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 8474b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 8484b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 8494b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 8504b11644fSPeter Brune snes->ops->view = 0; 8514b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 8524b11644fSPeter Brune 85342f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 85442f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 85542f4f86dSBarry Smith 85688976e71SPeter Brune if (!snes->tolerancesset) { 8570e444f03SPeter Brune snes->max_funcs = 30000; 8580e444f03SPeter Brune snes->max_its = 10000; 85988976e71SPeter Brune } 8600e444f03SPeter Brune 8619f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 8624b11644fSPeter Brune snes->data = (void *) qn; 8630ecc9e1dSPeter Brune qn->m = 10; 864b21d5a53SPeter Brune qn->scaling = 1.0; 865*b8840d6bSPeter Brune qn->U = PETSC_NULL; 866*b8840d6bSPeter Brune qn->V = PETSC_NULL; 8678e231d97SPeter Brune qn->dXtdF = PETSC_NULL; 8688e231d97SPeter Brune qn->dFtdX = PETSC_NULL; 8698e231d97SPeter Brune qn->dXdFmat = PETSC_NULL; 87044f7e39eSPeter Brune qn->monitor = PETSC_NULL; 87118aa0c0cSPeter Brune qn->singlereduction = PETSC_FALSE; 872*b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 8730c777b0cSPeter Brune qn->composition_type= SNES_QN_SEQUENTIAL; 8740c777b0cSPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 8750c777b0cSPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 876*b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 877ea630c6eSPeter Brune 8784b11644fSPeter Brune PetscFunctionReturn(0); 8794b11644fSPeter Brune } 8808e231d97SPeter Brune 8814b11644fSPeter Brune EXTERN_C_END 882