10c777b0cSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 24b11644fSPeter Brune 38e231d97SPeter Brune #define H(i,j) qn->dXdFmat[i*qn->m + j] 48e231d97SPeter Brune 56a6fc655SJed Brown const char *const SNESQNScaleTypes[] = {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0}; 66a6fc655SJed Brown const char *const SNESQNRestartTypes[] = {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0}; 7b8840d6bSPeter Brune const char *const SNESQNTypes[] = {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_RESTART_",0}; 8b8840d6bSPeter Brune 9b8840d6bSPeter Brune typedef enum {SNES_QN_LBFGS = 0, 10b8840d6bSPeter Brune SNES_QN_BROYDEN = 1, 110ad7597dSKarl Rupp SNES_QN_BADBROYDEN = 2 120ad7597dSKarl Rupp } SNESQNType; 136bf1b2e5SPeter Brune 144b11644fSPeter Brune typedef struct { 15b8840d6bSPeter Brune Vec *U; /* Stored past states (vary from method to method) */ 16b8840d6bSPeter Brune Vec *V; /* Stored past states (vary from method to method) */ 178e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 185c7a0a03SPeter Brune PetscReal *lambda; /* The line search history of the method */ 195c7a0a03SPeter Brune PetscReal *norm; /* norms of the steps */ 208e231d97SPeter Brune PetscScalar *alpha, *beta; 218e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 2218aa0c0cSPeter Brune PetscBool singlereduction; /* Aggregated reduction implementation */ 238e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 2444f7e39eSPeter Brune PetscViewer monitor; 256bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 266bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 27b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 28b8840d6bSPeter Brune SNESQNType type; /* the type of quasi-newton method used */ 2988f769c5SPeter Brune SNESQNScaleType scale_type; /* the type of scaling used */ 300c777b0cSPeter Brune SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 319f83bee8SJed Brown } SNES_QN; 324b11644fSPeter Brune 334b11644fSPeter Brune #undef __FUNCT__ 34b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden" 355051edb1SPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D) 362150357eSBarry Smith { 374b11644fSPeter Brune PetscErrorCode ierr; 389f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 39b8840d6bSPeter Brune Vec W = snes->work[3]; 40b8840d6bSPeter Brune Vec *U = qn->U; 41b8840d6bSPeter Brune KSPConvergedReason kspreason; 42b8840d6bSPeter Brune PetscInt m = qn->m; 435c7a0a03SPeter Brune PetscInt k,i,j,lits,l = m; 445c7a0a03SPeter Brune PetscReal unorm,a,b; 455c7a0a03SPeter Brune PetscReal *lambda=qn->lambda; 46b8840d6bSPeter Brune PetscScalar gdot; 4759e7931aSPeter Brune PetscReal udot; 482150357eSBarry Smith 49b8840d6bSPeter Brune PetscFunctionBegin; 50b8840d6bSPeter Brune if (it < m) l = it; 515c7a0a03SPeter Brune if (it > 0) { 525c7a0a03SPeter Brune k = (it-1)%l; 535c7a0a03SPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr); 54*a49fa0d8SPeter Brune ierr = VecCopy(Xold,U[k]);CHKERRQ(ierr); 555051edb1SPeter Brune ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr); 565c7a0a03SPeter Brune if (qn->monitor) { 575c7a0a03SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 585c7a0a03SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %d of %d by lambda: %14.12e \n",k,l,lambda[k]);CHKERRQ(ierr); 595c7a0a03SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 605c7a0a03SPeter Brune } 615c7a0a03SPeter Brune } 62b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 63d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr); 64b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 65b8840d6bSPeter Brune if (kspreason < 0) { 66b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 67b8840d6bSPeter 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); 68b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 69b8840d6bSPeter Brune PetscFunctionReturn(0); 70b8840d6bSPeter Brune } 71b8840d6bSPeter Brune } 72b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 73b8840d6bSPeter Brune snes->linear_its += lits; 74b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 75b8840d6bSPeter Brune } else { 76b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 77b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 78b8840d6bSPeter Brune } 79b8840d6bSPeter Brune 80b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 81b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 825c7a0a03SPeter Brune j = (it+i-l)%l; 835c7a0a03SPeter Brune k = (it+i-l+1)%l; 845c7a0a03SPeter Brune ierr = VecNorm(U[j],NORM_2,&unorm);CHKERRQ(ierr); 855c7a0a03SPeter Brune ierr = VecDot(U[j],Y,&gdot);CHKERRQ(ierr); 8659e7931aSPeter Brune unorm *= unorm; 8759e7931aSPeter Brune udot = PetscRealPart(gdot); 885c7a0a03SPeter Brune a = (lambda[j]/lambda[k]); 895c7a0a03SPeter Brune b = -(1.-lambda[j]); 9059e7931aSPeter Brune a *= udot/unorm; 9159e7931aSPeter Brune b *= udot/unorm; 925c7a0a03SPeter Brune ierr = VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);CHKERRQ(ierr); 935c7a0a03SPeter Brune 94b8840d6bSPeter Brune if (qn->monitor) { 95b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 965c7a0a03SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d and %d, gdot: %14.12e\n",k,j,gdot);CHKERRQ(ierr); 97b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 98b8840d6bSPeter Brune } 99b8840d6bSPeter Brune } 100b8840d6bSPeter Brune if (l > 0) { 101b8840d6bSPeter Brune k = (it-1)%l; 102b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 1035c7a0a03SPeter Brune ierr = VecNorm(U[k],NORM_2,&unorm);CHKERRQ(ierr); 10459e7931aSPeter Brune unorm *= unorm; 10559e7931aSPeter Brune udot = PetscRealPart(gdot); 10659e7931aSPeter Brune a = unorm/(unorm-lambda[k]*udot); 10759e7931aSPeter Brune b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot); 108b8840d6bSPeter Brune if (qn->monitor) { 109b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 11059e7931aSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d: a: %14.12e b: %14.12e \n",k,a,b);CHKERRQ(ierr); 111b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 112b8840d6bSPeter Brune } 1135c7a0a03SPeter Brune ierr = VecAXPBY(Y,b,a,U[k]);CHKERRQ(ierr); 114b8840d6bSPeter Brune } 1155c7a0a03SPeter Brune l = m; 1165c7a0a03SPeter Brune if (it+1<m)l=it+1; 1175c7a0a03SPeter Brune k = it%l; 1185c7a0a03SPeter Brune if (qn->monitor) { 1195c7a0a03SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1205c7a0a03SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "setting vector %d of %d\n",k,l);CHKERRQ(ierr); 1215c7a0a03SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1225c7a0a03SPeter Brune } 123b8840d6bSPeter Brune PetscFunctionReturn(0); 124b8840d6bSPeter Brune } 125b8840d6bSPeter Brune 126b8840d6bSPeter Brune #undef __FUNCT__ 127b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden" 1280adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1290adebc6cSBarry Smith { 130b8840d6bSPeter Brune PetscErrorCode ierr; 131b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 132b8840d6bSPeter Brune Vec W = snes->work[3]; 133b8840d6bSPeter Brune Vec *U = qn->U; 134b8840d6bSPeter Brune Vec *T = qn->V; 135b8840d6bSPeter Brune 136b8840d6bSPeter Brune /* ksp thing for jacobian scaling */ 137b8840d6bSPeter Brune KSPConvergedReason kspreason; 138c9e59058SPeter Brune PetscInt h,k,j,i,lits; 139b8840d6bSPeter Brune PetscInt m = qn->m; 1405051edb1SPeter Brune PetscScalar gdot,udot; 141b8840d6bSPeter Brune PetscInt l = m; 1420adebc6cSBarry Smith 143b8840d6bSPeter Brune PetscFunctionBegin; 144b8840d6bSPeter Brune if (it < m) l = it; 145b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 146b8840d6bSPeter Brune if (l > 0) { 147b8840d6bSPeter Brune k = (it-1)%l; 148c9e59058SPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);CHKERRQ(ierr); 149b8840d6bSPeter Brune ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 150b8840d6bSPeter Brune ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 1515051edb1SPeter Brune ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr); 1525051edb1SPeter Brune ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr); 153b8840d6bSPeter Brune } 154b8840d6bSPeter Brune 155b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 156d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 157b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 158b8840d6bSPeter Brune if (kspreason < 0) { 159b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 160b8840d6bSPeter 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); 161b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 162b8840d6bSPeter Brune PetscFunctionReturn(0); 163b8840d6bSPeter Brune } 164b8840d6bSPeter Brune } 165b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 166b8840d6bSPeter Brune snes->linear_its += lits; 167b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 168b8840d6bSPeter Brune } else { 169b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 170b8840d6bSPeter Brune } 1715051edb1SPeter Brune 1725051edb1SPeter Brune /* inward recursion starting at the first update and working forward */ 1735051edb1SPeter Brune if (l > 0) { 1745051edb1SPeter Brune for (i = 0; i < l-1; i++) { 175c9e59058SPeter Brune j = (it+i-l)%l; 176c9e59058SPeter Brune k = (it+i-l+1)%l; 177c9e59058SPeter Brune h = (it-1)%l; 178c9e59058SPeter Brune ierr = VecDotBegin(U[j],U[h],&gdot);CHKERRQ(ierr); 179c9e59058SPeter Brune ierr = VecDotBegin(U[j],U[j],&udot);CHKERRQ(ierr); 180c9e59058SPeter Brune ierr = VecDotEnd(U[j],U[h],&gdot);CHKERRQ(ierr); 181c9e59058SPeter Brune ierr = VecDotEnd(U[j],U[j],&udot);CHKERRQ(ierr); 1825051edb1SPeter Brune ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr); 183c9e59058SPeter Brune ierr = VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);CHKERRQ(ierr); 1845051edb1SPeter Brune if (qn->monitor) { 1855051edb1SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1865051edb1SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 1875051edb1SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1885051edb1SPeter Brune } 1895051edb1SPeter Brune } 1905051edb1SPeter Brune } 191b8840d6bSPeter Brune PetscFunctionReturn(0); 192b8840d6bSPeter Brune } 193b8840d6bSPeter Brune 194b8840d6bSPeter Brune #undef __FUNCT__ 195b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 1960adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1970adebc6cSBarry Smith { 198b8840d6bSPeter Brune PetscErrorCode ierr; 199b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 200b8840d6bSPeter Brune Vec W = snes->work[3]; 201b8840d6bSPeter Brune Vec *dX = qn->U; 202b8840d6bSPeter Brune Vec *dF = qn->V; 203bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 204bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 2058e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 206b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 2078e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 208bd052dfeSPeter Brune 2090ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 2100ecc9e1dSPeter Brune KSPConvergedReason kspreason; 2118e231d97SPeter Brune PetscInt k,i,j,g,lits; 2124b11644fSPeter Brune PetscInt m = qn->m; 2134b11644fSPeter Brune PetscScalar t; 2144b11644fSPeter Brune PetscInt l = m; 2158e231d97SPeter Brune 2164b11644fSPeter Brune PetscFunctionBegin; 2175ba6227bSPeter Brune if (it < m) l = it; 2181c6523dcSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 219b8840d6bSPeter Brune if (it > 0) { 220b8840d6bSPeter Brune k = (it - 1) % l; 221b8840d6bSPeter Brune ierr = VecCopy(D,dF[k]);CHKERRQ(ierr); 222b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 223b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 224b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 22518aa0c0cSPeter Brune if (qn->singlereduction) { 2261c6523dcSPeter Brune ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2271c6523dcSPeter Brune ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2281c6523dcSPeter Brune ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr); 2291c6523dcSPeter Brune ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2301c6523dcSPeter Brune ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2311c6523dcSPeter Brune ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr); 232b8840d6bSPeter Brune for (j = 0; j < l; j++) { 233b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 234b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 235b8840d6bSPeter Brune } 236b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 2371aa26658SKarl Rupp for (j = 0; j < l; j++) dXtdF[j] = H(j, j); 238b8840d6bSPeter Brune } else { 239b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 24023d44fbcSPeter Brune } 24123d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 242b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 243b8840d6bSPeter Brune } 244b8840d6bSPeter Brune } 245b8840d6bSPeter Brune 2464b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2474b11644fSPeter Brune for (i=0; i<l; i++) { 248b21d5a53SPeter Brune k = (it-i-1)%l; 24918aa0c0cSPeter Brune if (qn->singlereduction) { 2508e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2518e231d97SPeter Brune t = YtdX[k]; 2528e231d97SPeter Brune for (j=0; j<i; j++) { 2538e231d97SPeter Brune g = (it-j-1)%l; 254b3cfccb2SPeter Brune t -= alpha[g]*H(k, g); 2558e231d97SPeter Brune } 2568e231d97SPeter Brune alpha[k] = t/H(k,k); 2578e231d97SPeter Brune } else { 2583af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2598e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2608e231d97SPeter Brune } 26144f7e39eSPeter Brune if (qn->monitor) { 2623af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2635ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2643af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 26544f7e39eSPeter Brune } 2666bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2674b11644fSPeter Brune } 2684b11644fSPeter Brune 2690c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 270d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 2710ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2720ecc9e1dSPeter Brune if (kspreason < 0) { 2730ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 2740ecc9e1dSPeter 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); 2750ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 2760ecc9e1dSPeter Brune PetscFunctionReturn(0); 2770ecc9e1dSPeter Brune } 2780ecc9e1dSPeter Brune } 2790ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2800ecc9e1dSPeter Brune snes->linear_its += lits; 281b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2820ecc9e1dSPeter Brune } else { 283b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2840ecc9e1dSPeter Brune } 28518aa0c0cSPeter Brune if (qn->singlereduction) { 286b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 2878e231d97SPeter Brune } 2884b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 289bd052dfeSPeter Brune for (i = 0; i < l; i++) { 290b21d5a53SPeter Brune k = (it + i - l) % l; 29118aa0c0cSPeter Brune if (qn->singlereduction) { 2928e231d97SPeter Brune t = YtdX[k]; 2938e231d97SPeter Brune for (j = 0; j < i; j++) { 2948e231d97SPeter Brune g = (it + j - l) % l; 295b3cfccb2SPeter Brune t += (alpha[g] - beta[g])*H(g, k); 2968e231d97SPeter Brune } 2978e231d97SPeter Brune beta[k] = t / H(k, k); 2988e231d97SPeter Brune } else { 2996bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 3008e231d97SPeter Brune beta[k] = t / dXtdF[k]; 3018e231d97SPeter Brune } 30222d28d08SBarry Smith ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 30344f7e39eSPeter Brune if (qn->monitor) { 3043af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3055ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 3063af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 30744f7e39eSPeter Brune } 3084b11644fSPeter Brune } 3094b11644fSPeter Brune PetscFunctionReturn(0); 3104b11644fSPeter Brune } 3114b11644fSPeter Brune 3124b11644fSPeter Brune #undef __FUNCT__ 3134b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 3144b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 3154b11644fSPeter Brune { 3164b11644fSPeter Brune PetscErrorCode ierr; 3179f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 31815f5eeeaSPeter Brune Vec X,Xold; 3190a3905e1SPeter Brune Vec F,W; 32047f26062SPeter Brune Vec Y,D,Dold; 321b8840d6bSPeter Brune PetscInt i, i_r; 322335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 3230c777b0cSPeter Brune PetscBool lssucceed,powell,periodic; 3241c6523dcSPeter Brune PetscScalar DolddotD,DolddotDold; 3250ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 326b7281c8aSPeter Brune SNESConvergedReason reason; 3270ecc9e1dSPeter Brune 3284b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 3294b11644fSPeter Brune 3306e111a19SKarl Rupp PetscFunctionBegin; 3319f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3323af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 3330a3905e1SPeter Brune W = snes->work[3]; 334b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 335335fb975SPeter Brune Xold = snes->work[0]; 3363af51624SPeter Brune 3373af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 338335fb975SPeter Brune D = snes->work[1]; 339335fb975SPeter Brune Dold = snes->work[2]; 3404b11644fSPeter Brune 3414b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3424b11644fSPeter Brune 343ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3444b11644fSPeter Brune snes->iter = 0; 3454b11644fSPeter Brune snes->norm = 0.; 346ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 34747f26062SPeter Brune 348b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 34947f26062SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 350b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 351b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 352b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 353b7281c8aSPeter Brune PetscFunctionReturn(0); 354b7281c8aSPeter Brune } 35547f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 35647f26062SPeter Brune } else { 357e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 35815f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3594b11644fSPeter Brune if (snes->domainerror) { 3604b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3614b11644fSPeter Brune PetscFunctionReturn(0); 3624b11644fSPeter Brune } 3631aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 364c1c75074SPeter Brune 36547f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 366189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 367189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 368189a9710SBarry Smith PetscFunctionReturn(0); 369189a9710SBarry Smith } 37047f26062SPeter Brune } 371b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 37247f26062SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 373b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 374b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 375b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 376b7281c8aSPeter Brune PetscFunctionReturn(0); 377b7281c8aSPeter Brune } 37847f26062SPeter Brune } else { 37947f26062SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 38047f26062SPeter Brune } 381b28a06ddSPeter Brune 382ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3834b11644fSPeter Brune snes->norm = fnorm; 384ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 385a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 3864b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3874b11644fSPeter Brune 3884b11644fSPeter Brune /* test convergence */ 3894b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3904b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 39170d3b23bSPeter Brune 3923cf07b75SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 3933cf07b75SPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 3943cf07b75SPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 3953cf07b75SPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 396ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 397ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 398ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 399ddd40ce5SPeter Brune PetscFunctionReturn(0); 400ddd40ce5SPeter Brune } 401ddd40ce5SPeter Brune ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 4023cf07b75SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 4033cf07b75SPeter Brune } 4043cf07b75SPeter Brune 405f8e15203SPeter Brune /* scale the initial update */ 4060c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4070ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4088cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 4090ecc9e1dSPeter Brune } 4100ecc9e1dSPeter Brune 4115ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 4120a3905e1SPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) { 4130a3905e1SPeter Brune PetscScalar ff,xf; 4140a3905e1SPeter Brune ierr = VecCopy(Dold,Y);CHKERRQ(ierr); 4150a3905e1SPeter Brune ierr = VecCopy(Xold,W);CHKERRQ(ierr); 4160a3905e1SPeter Brune ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr); 4170a3905e1SPeter Brune ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 4180a3905e1SPeter Brune ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr); 4190a3905e1SPeter Brune ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr); 4200a3905e1SPeter Brune ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr); 4210a3905e1SPeter Brune ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr); 4220a3905e1SPeter Brune qn->scaling = PetscRealPart(xf)/PetscRealPart(ff); 4230a3905e1SPeter Brune } 424b8840d6bSPeter Brune switch (qn->type) { 425b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 426b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 427b8840d6bSPeter Brune break; 428b8840d6bSPeter Brune case SNES_QN_BROYDEN: 4295051edb1SPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr); 430b8840d6bSPeter Brune break; 431b8840d6bSPeter Brune case SNES_QN_LBFGS: 432b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 433b8840d6bSPeter Brune break; 434b8840d6bSPeter Brune } 43570d3b23bSPeter Brune /* line search for lambda */ 43670d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 4370a3905e1SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 43815f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 439f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4409f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 4419f3a0142SPeter Brune if (snes->domainerror) { 4429f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 4439f3a0142SPeter Brune PetscFunctionReturn(0); 4449f3a0142SPeter Brune } 445f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 4469f3a0142SPeter Brune if (!lssucceed) { 4479f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4489f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4499f3a0142SPeter Brune break; 4509f3a0142SPeter Brune } 4519f3a0142SPeter Brune } 452f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 4530c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 45405ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4550ecc9e1dSPeter Brune } 4563af51624SPeter Brune 45788d374b2SPeter Brune /* convergence monitoring */ 458b21d5a53SPeter 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); 459b21d5a53SPeter Brune 460b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 461b28a06ddSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 462b28a06ddSPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 463b28a06ddSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 464ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 465ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 466ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 467ddd40ce5SPeter Brune PetscFunctionReturn(0); 468ddd40ce5SPeter Brune } 469ddd40ce5SPeter Brune ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 470b28a06ddSPeter Brune } 471b28a06ddSPeter Brune 472360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 473360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 474360c497dSPeter Brune 475a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 4768409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4774b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 478d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4794b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 480b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 48147f26062SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 482b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 483b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 484b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 485b7281c8aSPeter Brune PetscFunctionReturn(0); 486b7281c8aSPeter Brune } 48788d374b2SPeter Brune } else { 48888d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 48988d374b2SPeter Brune } 4900c777b0cSPeter Brune powell = PETSC_FALSE; 4910c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 4926bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 49323c918e7SPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 49423c918e7SPeter Brune ierr = MatMult(snes->jacobian_pre,Dold,W);CHKERRQ(ierr); 49523c918e7SPeter Brune } else { 49623c918e7SPeter Brune ierr = VecCopy(Dold,W);CHKERRQ(ierr); 49723c918e7SPeter Brune } 49823c918e7SPeter Brune ierr = VecDotBegin(W, Dold, &DolddotDold);CHKERRQ(ierr); 49923c918e7SPeter Brune ierr = VecDotBegin(W, D, &DolddotD);CHKERRQ(ierr); 50023c918e7SPeter Brune ierr = VecDotEnd(W, Dold, &DolddotDold);CHKERRQ(ierr); 50123c918e7SPeter Brune ierr = VecDotEnd(W, D, &DolddotD);CHKERRQ(ierr); 5020c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 5030c777b0cSPeter Brune } 5040c777b0cSPeter Brune periodic = PETSC_FALSE; 505b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 506b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 5070c777b0cSPeter Brune } 5080c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 5090c777b0cSPeter Brune if (powell || periodic) { 5105ba6227bSPeter Brune if (qn->monitor) { 5115ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 512516fe3c3SPeter 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); 5135ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 5145ba6227bSPeter Brune } 5155ba6227bSPeter Brune i_r = -1; 5165ba6227bSPeter Brune /* general purpose update */ 5175ba6227bSPeter Brune if (snes->ops->update) { 5185ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 5195ba6227bSPeter Brune } 5200c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5210ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 5228cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 5230ecc9e1dSPeter Brune } 5240ecc9e1dSPeter Brune } 52570d3b23bSPeter Brune /* general purpose update */ 52670d3b23bSPeter Brune if (snes->ops->update) { 52770d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 52870d3b23bSPeter Brune } 5295ba6227bSPeter Brune } 5304b11644fSPeter Brune if (i == snes->max_its) { 5314b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 5324b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 5334b11644fSPeter Brune } 5344b11644fSPeter Brune PetscFunctionReturn(0); 5354b11644fSPeter Brune } 5364b11644fSPeter Brune 5374b11644fSPeter Brune #undef __FUNCT__ 5384b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5394b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5404b11644fSPeter Brune { 5419f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5424b11644fSPeter Brune PetscErrorCode ierr; 543335fb975SPeter Brune 5444b11644fSPeter Brune PetscFunctionBegin; 545b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 54659e7931aSPeter Brune if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 5475c7a0a03SPeter Brune ierr = PetscMalloc4(qn->m, PetscScalar, &qn->alpha, 5488e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 5495c7a0a03SPeter Brune qn->m, PetscScalar, &qn->dXtdF, 5505c7a0a03SPeter Brune qn->m, PetscReal, &qn->lambda);CHKERRQ(ierr); 5518e231d97SPeter Brune 55218aa0c0cSPeter Brune if (qn->singlereduction) { 5538e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 5548e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 5558e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 5568e231d97SPeter Brune } 557fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 558335fb975SPeter Brune /* set up the line search */ 5590c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5608e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5618e231d97SPeter Brune } 5626c67d002SPeter Brune if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;} 5634b11644fSPeter Brune PetscFunctionReturn(0); 5644b11644fSPeter Brune } 5654b11644fSPeter Brune 5664b11644fSPeter Brune #undef __FUNCT__ 5674b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5684b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5694b11644fSPeter Brune { 5704b11644fSPeter Brune PetscErrorCode ierr; 5719f83bee8SJed Brown SNES_QN *qn; 5720adebc6cSBarry Smith 5734b11644fSPeter Brune PetscFunctionBegin; 5744b11644fSPeter Brune if (snes->data) { 5759f83bee8SJed Brown qn = (SNES_QN*)snes->data; 576b8840d6bSPeter Brune if (qn->U) { 577b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5784b11644fSPeter Brune } 579b8840d6bSPeter Brune if (qn->V) { 580b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5814b11644fSPeter Brune } 58218aa0c0cSPeter Brune if (qn->singlereduction) { 5838e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5848e231d97SPeter Brune } 5855c7a0a03SPeter Brune ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr); 5864b11644fSPeter Brune } 5874b11644fSPeter Brune PetscFunctionReturn(0); 5884b11644fSPeter Brune } 5894b11644fSPeter Brune 5904b11644fSPeter Brune #undef __FUNCT__ 5914b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5924b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5934b11644fSPeter Brune { 5944b11644fSPeter Brune PetscErrorCode ierr; 5956e111a19SKarl Rupp 5964b11644fSPeter Brune PetscFunctionBegin; 5974b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 5984b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 599bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr); 6004b11644fSPeter Brune PetscFunctionReturn(0); 6014b11644fSPeter Brune } 6024b11644fSPeter Brune 6034b11644fSPeter Brune #undef __FUNCT__ 6044b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 6054b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 6064b11644fSPeter Brune { 6074b11644fSPeter Brune 6084b11644fSPeter Brune PetscErrorCode ierr; 6092150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 61088f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 611f1c6b773SPeter Brune SNESLineSearch linesearch; 6122150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 6132150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 6142150357eSBarry Smith 6154b11644fSPeter Brune PetscFunctionBegin; 6164b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 6170298fd71SBarry Smith ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr); 6180298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr); 6190298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr); 6200298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr); 6210298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr); 62288f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 62388f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 62488f769c5SPeter Brune 62588f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 62688f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 62788f769c5SPeter Brune 628b8840d6bSPeter Brune ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 6290298fd71SBarry Smith (PetscEnum)qn->type,(PetscEnum*)&qn->type,NULL);CHKERRQ(ierr); 6304b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6319e764e56SPeter Brune if (!snes->linesearch) { 6327601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 6331a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 6349e764e56SPeter Brune } 63544f7e39eSPeter Brune if (monflg) { 636ce94432eSBarry Smith qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 63744f7e39eSPeter Brune } 6384b11644fSPeter Brune PetscFunctionReturn(0); 6394b11644fSPeter Brune } 6404b11644fSPeter Brune 64192c02d66SPeter Brune 6420c777b0cSPeter Brune #undef __FUNCT__ 6430c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6440c777b0cSPeter Brune /*@ 6450c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6460c777b0cSPeter Brune 6470c777b0cSPeter Brune Logically Collective on SNES 6480c777b0cSPeter Brune 6490c777b0cSPeter Brune Input Parameters: 6500c777b0cSPeter Brune + snes - the iterative context 6510c777b0cSPeter Brune - rtype - restart type 6520c777b0cSPeter Brune 6530c777b0cSPeter Brune Options Database: 6540c777b0cSPeter Brune + -snes_qn_restart_type<powell,periodic,none> - set the restart type 655b8840d6bSPeter Brune - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 6560c777b0cSPeter Brune 6570c777b0cSPeter Brune Level: intermediate 6580c777b0cSPeter Brune 6590c777b0cSPeter Brune SNESQNRestartTypes: 6600c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6610c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6620c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6630c777b0cSPeter Brune 6640c777b0cSPeter Brune Notes: 6650c777b0cSPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 6660c777b0cSPeter Brune 6670c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 6680c777b0cSPeter Brune @*/ 6692150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6702150357eSBarry Smith { 6710c777b0cSPeter Brune PetscErrorCode ierr; 6726e111a19SKarl Rupp 6730c777b0cSPeter Brune PetscFunctionBegin; 6740c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6750c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6760c777b0cSPeter Brune PetscFunctionReturn(0); 6770c777b0cSPeter Brune } 6780c777b0cSPeter Brune 6790c777b0cSPeter Brune #undef __FUNCT__ 6800c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6810c777b0cSPeter Brune /*@ 6820c777b0cSPeter Brune SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 6830c777b0cSPeter Brune 6840c777b0cSPeter Brune Logically Collective on SNES 6850c777b0cSPeter Brune 6860c777b0cSPeter Brune Input Parameters: 6870c777b0cSPeter Brune + snes - the iterative context 6880c777b0cSPeter Brune - stype - scale type 6890c777b0cSPeter Brune 6900c777b0cSPeter Brune Options Database: 6910c777b0cSPeter Brune . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 6920c777b0cSPeter Brune 6930c777b0cSPeter Brune Level: intermediate 6940c777b0cSPeter Brune 6950c777b0cSPeter Brune SNESQNSelectTypes: 6960c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 6970c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 6980c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 6990c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 7000c777b0cSPeter Brune 7010c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 7020c777b0cSPeter Brune @*/ 7030c777b0cSPeter Brune 7042150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 7052150357eSBarry Smith { 7060c777b0cSPeter Brune PetscErrorCode ierr; 7076e111a19SKarl Rupp 7080c777b0cSPeter Brune PetscFunctionBegin; 7090c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7100c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 7110c777b0cSPeter Brune PetscFunctionReturn(0); 7120c777b0cSPeter Brune } 7130c777b0cSPeter Brune 7140c777b0cSPeter Brune #undef __FUNCT__ 7150c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 7160adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 7170adebc6cSBarry Smith { 7180c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7196e111a19SKarl Rupp 7200c777b0cSPeter Brune PetscFunctionBegin; 7210c777b0cSPeter Brune qn->scale_type = stype; 7220c777b0cSPeter Brune PetscFunctionReturn(0); 7230c777b0cSPeter Brune } 7240c777b0cSPeter Brune 7250c777b0cSPeter Brune #undef __FUNCT__ 7260c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7270adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 7280adebc6cSBarry Smith { 7290c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7306e111a19SKarl Rupp 7310c777b0cSPeter Brune PetscFunctionBegin; 7320c777b0cSPeter Brune qn->restart_type = rtype; 7330c777b0cSPeter Brune PetscFunctionReturn(0); 7340c777b0cSPeter Brune } 7350c777b0cSPeter Brune 7364b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 7374b11644fSPeter Brune /*MC 7384b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 7394b11644fSPeter Brune 7406cc8130cSPeter Brune Options Database: 7416cc8130cSPeter Brune 7426cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 7431867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 7441867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 74572835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 74644f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 7476cc8130cSPeter Brune 748b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 749b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 750b8840d6bSPeter Brune updates. 7516cc8130cSPeter Brune 7521867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 7531867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 7541867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 7551867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 7561867fe5bSPeter Brune 7576cc8130cSPeter Brune References: 7586cc8130cSPeter Brune 7590a8e8854SPeter Brune Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995. 7606cc8130cSPeter Brune 7610a8e8854SPeter Brune R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi-Newton Matrices and their use in Limited Memory Methods, 7620a8e8854SPeter Brune Technical Report, Northwestern University, June 1992. 7630a8e8854SPeter Brune 7640a8e8854SPeter Brune Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE 7650a8e8854SPeter Brune Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 766b8840d6bSPeter Brune 7674b11644fSPeter Brune 7684b11644fSPeter Brune Level: beginner 7694b11644fSPeter Brune 77004d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 7716cc8130cSPeter Brune 7724b11644fSPeter Brune M*/ 7734b11644fSPeter Brune #undef __FUNCT__ 7744b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 7758cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 7764b11644fSPeter Brune { 7774b11644fSPeter Brune PetscErrorCode ierr; 7789f83bee8SJed Brown SNES_QN *qn; 7794b11644fSPeter Brune 7804b11644fSPeter Brune PetscFunctionBegin; 7814b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 7824b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 7834b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 7844b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 7854b11644fSPeter Brune snes->ops->view = 0; 7864b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 7874b11644fSPeter Brune 78847f26062SPeter Brune snes->pcside = PC_LEFT; 78947f26062SPeter Brune 79042f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 79142f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 79242f4f86dSBarry Smith 79388976e71SPeter Brune if (!snes->tolerancesset) { 7940e444f03SPeter Brune snes->max_funcs = 30000; 7950e444f03SPeter Brune snes->max_its = 10000; 79688976e71SPeter Brune } 7970e444f03SPeter Brune 7989f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 7994b11644fSPeter Brune snes->data = (void*) qn; 8000ecc9e1dSPeter Brune qn->m = 10; 801b21d5a53SPeter Brune qn->scaling = 1.0; 8020298fd71SBarry Smith qn->U = NULL; 8030298fd71SBarry Smith qn->V = NULL; 8040298fd71SBarry Smith qn->dXtdF = NULL; 8050298fd71SBarry Smith qn->dFtdX = NULL; 8060298fd71SBarry Smith qn->dXdFmat = NULL; 8070298fd71SBarry Smith qn->monitor = NULL; 8081c6523dcSPeter Brune qn->singlereduction = PETSC_TRUE; 809b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 8100c777b0cSPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 8110c777b0cSPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 812b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 813ea630c6eSPeter Brune 814bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr); 815bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr); 8164b11644fSPeter Brune PetscFunctionReturn(0); 8174b11644fSPeter Brune } 818