1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2*fced5a79SAsbjørn Nilsen Riseth #include <petscdm.h> 34b11644fSPeter Brune 48e231d97SPeter Brune #define H(i,j) qn->dXdFmat[i*qn->m + j] 58e231d97SPeter Brune 61efc8c45SPeter Brune const char *const SNESQNScaleTypes[] = {"DEFAULT","NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0}; 71efc8c45SPeter Brune const char *const SNESQNRestartTypes[] = {"DEFAULT","NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0}; 860c08b40SPeter Brune const char *const SNESQNTypes[] = {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_",0}; 9b8840d6bSPeter Brune 104b11644fSPeter Brune typedef struct { 11b8840d6bSPeter Brune Vec *U; /* Stored past states (vary from method to method) */ 12b8840d6bSPeter Brune Vec *V; /* Stored past states (vary from method to method) */ 138e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 145c7a0a03SPeter Brune PetscReal *lambda; /* The line search history of the method */ 155c7a0a03SPeter Brune PetscReal *norm; /* norms of the steps */ 168e231d97SPeter Brune PetscScalar *alpha, *beta; 178e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 1818aa0c0cSPeter Brune PetscBool singlereduction; /* Aggregated reduction implementation */ 198e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 2044f7e39eSPeter Brune PetscViewer monitor; 216bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 226bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 23b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 24b8840d6bSPeter Brune SNESQNType type; /* the type of quasi-newton method used */ 2588f769c5SPeter Brune SNESQNScaleType scale_type; /* the type of scaling used */ 260c777b0cSPeter Brune SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 279f83bee8SJed Brown } SNES_QN; 284b11644fSPeter Brune 294b11644fSPeter Brune #undef __FUNCT__ 30b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden" 315051edb1SPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D) 322150357eSBarry Smith { 334b11644fSPeter Brune PetscErrorCode ierr; 349f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 35b8840d6bSPeter Brune Vec W = snes->work[3]; 36b8840d6bSPeter Brune Vec *U = qn->U; 37b8840d6bSPeter Brune PetscInt m = qn->m; 385c7a0a03SPeter Brune PetscInt k,i,j,lits,l = m; 395c7a0a03SPeter Brune PetscReal unorm,a,b; 405c7a0a03SPeter Brune PetscReal *lambda=qn->lambda; 41b8840d6bSPeter Brune PetscScalar gdot; 4259e7931aSPeter Brune PetscReal udot; 432150357eSBarry Smith 44b8840d6bSPeter Brune PetscFunctionBegin; 45b8840d6bSPeter Brune if (it < m) l = it; 465c7a0a03SPeter Brune if (it > 0) { 475c7a0a03SPeter Brune k = (it-1)%l; 485c7a0a03SPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr); 49a49fa0d8SPeter Brune ierr = VecCopy(Xold,U[k]);CHKERRQ(ierr); 505051edb1SPeter Brune ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr); 515c7a0a03SPeter Brune if (qn->monitor) { 525c7a0a03SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 535c7a0a03SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %d of %d by lambda: %14.12e \n",k,l,lambda[k]);CHKERRQ(ierr); 545c7a0a03SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 555c7a0a03SPeter Brune } 565c7a0a03SPeter Brune } 57b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 58d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr); 59422a814eSBarry Smith SNESCheckKSPSolve(snes); 60b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 61b8840d6bSPeter Brune snes->linear_its += lits; 62b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 63b8840d6bSPeter Brune } else { 64b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 65b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 66b8840d6bSPeter Brune } 67b8840d6bSPeter Brune 68b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 69b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 705c7a0a03SPeter Brune j = (it+i-l)%l; 715c7a0a03SPeter Brune k = (it+i-l+1)%l; 725c7a0a03SPeter Brune ierr = VecNorm(U[j],NORM_2,&unorm);CHKERRQ(ierr); 735c7a0a03SPeter Brune ierr = VecDot(U[j],Y,&gdot);CHKERRQ(ierr); 7459e7931aSPeter Brune unorm *= unorm; 7559e7931aSPeter Brune udot = PetscRealPart(gdot); 765c7a0a03SPeter Brune a = (lambda[j]/lambda[k]); 775c7a0a03SPeter Brune b = -(1.-lambda[j]); 7859e7931aSPeter Brune a *= udot/unorm; 7959e7931aSPeter Brune b *= udot/unorm; 805c7a0a03SPeter Brune ierr = VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);CHKERRQ(ierr); 815c7a0a03SPeter Brune 82b8840d6bSPeter Brune if (qn->monitor) { 83b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 84e639b251SJed Brown ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d and %d, gdot: %14.12e\n",k,j,PetscRealPart(gdot));CHKERRQ(ierr); 85b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 86b8840d6bSPeter Brune } 87b8840d6bSPeter Brune } 88b8840d6bSPeter Brune if (l > 0) { 89b8840d6bSPeter Brune k = (it-1)%l; 90b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 915c7a0a03SPeter Brune ierr = VecNorm(U[k],NORM_2,&unorm);CHKERRQ(ierr); 9259e7931aSPeter Brune unorm *= unorm; 9359e7931aSPeter Brune udot = PetscRealPart(gdot); 9459e7931aSPeter Brune a = unorm/(unorm-lambda[k]*udot); 9559e7931aSPeter Brune b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot); 96b8840d6bSPeter Brune if (qn->monitor) { 97b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 9859e7931aSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d: a: %14.12e b: %14.12e \n",k,a,b);CHKERRQ(ierr); 99b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 100b8840d6bSPeter Brune } 1015c7a0a03SPeter Brune ierr = VecAXPBY(Y,b,a,U[k]);CHKERRQ(ierr); 102b8840d6bSPeter Brune } 1035c7a0a03SPeter Brune l = m; 1045c7a0a03SPeter Brune if (it+1<m)l=it+1; 1055c7a0a03SPeter Brune k = it%l; 1065c7a0a03SPeter Brune if (qn->monitor) { 1075c7a0a03SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1085c7a0a03SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "setting vector %d of %d\n",k,l);CHKERRQ(ierr); 1095c7a0a03SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1105c7a0a03SPeter Brune } 111b8840d6bSPeter Brune PetscFunctionReturn(0); 112b8840d6bSPeter Brune } 113b8840d6bSPeter Brune 114b8840d6bSPeter Brune #undef __FUNCT__ 115b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden" 1160adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1170adebc6cSBarry Smith { 118b8840d6bSPeter Brune PetscErrorCode ierr; 119b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 120b8840d6bSPeter Brune Vec W = snes->work[3]; 121b8840d6bSPeter Brune Vec *U = qn->U; 122b8840d6bSPeter Brune Vec *T = qn->V; 123b8840d6bSPeter Brune 12484c577daSBarry Smith /* ksp thing for Jacobian scaling */ 125c9e59058SPeter Brune PetscInt h,k,j,i,lits; 126b8840d6bSPeter Brune PetscInt m = qn->m; 1275051edb1SPeter Brune PetscScalar gdot,udot; 128b8840d6bSPeter Brune PetscInt l = m; 1290adebc6cSBarry Smith 130b8840d6bSPeter Brune PetscFunctionBegin; 131b8840d6bSPeter Brune if (it < m) l = it; 132b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 133b8840d6bSPeter Brune if (l > 0) { 134b8840d6bSPeter Brune k = (it-1)%l; 135c9e59058SPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);CHKERRQ(ierr); 136b8840d6bSPeter Brune ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 137b8840d6bSPeter Brune ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 1385051edb1SPeter Brune ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr); 1395051edb1SPeter Brune ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr); 140b8840d6bSPeter Brune } 141b8840d6bSPeter Brune 142b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 143d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 144422a814eSBarry Smith SNESCheckKSPSolve(snes); 145b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 146b8840d6bSPeter Brune snes->linear_its += lits; 147b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 148b8840d6bSPeter Brune } else { 149b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 150b8840d6bSPeter Brune } 1515051edb1SPeter Brune 1525051edb1SPeter Brune /* inward recursion starting at the first update and working forward */ 1535051edb1SPeter Brune if (l > 0) { 1545051edb1SPeter Brune for (i = 0; i < l-1; i++) { 155c9e59058SPeter Brune j = (it+i-l)%l; 156c9e59058SPeter Brune k = (it+i-l+1)%l; 157c9e59058SPeter Brune h = (it-1)%l; 158c9e59058SPeter Brune ierr = VecDotBegin(U[j],U[h],&gdot);CHKERRQ(ierr); 159c9e59058SPeter Brune ierr = VecDotBegin(U[j],U[j],&udot);CHKERRQ(ierr); 160c9e59058SPeter Brune ierr = VecDotEnd(U[j],U[h],&gdot);CHKERRQ(ierr); 161c9e59058SPeter Brune ierr = VecDotEnd(U[j],U[j],&udot);CHKERRQ(ierr); 1625051edb1SPeter Brune ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr); 163c9e59058SPeter Brune ierr = VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);CHKERRQ(ierr); 1645051edb1SPeter Brune if (qn->monitor) { 1655051edb1SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1665051edb1SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 1675051edb1SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1685051edb1SPeter Brune } 1695051edb1SPeter Brune } 1705051edb1SPeter Brune } 171b8840d6bSPeter Brune PetscFunctionReturn(0); 172b8840d6bSPeter Brune } 173b8840d6bSPeter Brune 174b8840d6bSPeter Brune #undef __FUNCT__ 175b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 1760adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1770adebc6cSBarry Smith { 178b8840d6bSPeter Brune PetscErrorCode ierr; 179b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 180b8840d6bSPeter Brune Vec W = snes->work[3]; 181b8840d6bSPeter Brune Vec *dX = qn->U; 182b8840d6bSPeter Brune Vec *dF = qn->V; 183bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 184bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 1858e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 186b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 1878e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 188bd052dfeSPeter Brune 18984c577daSBarry Smith /* ksp thing for Jacobian scaling */ 1908e231d97SPeter Brune PetscInt k,i,j,g,lits; 1914b11644fSPeter Brune PetscInt m = qn->m; 1924b11644fSPeter Brune PetscScalar t; 1934b11644fSPeter Brune PetscInt l = m; 1948e231d97SPeter Brune 1954b11644fSPeter Brune PetscFunctionBegin; 1965ba6227bSPeter Brune if (it < m) l = it; 1971c6523dcSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 198b8840d6bSPeter Brune if (it > 0) { 199b8840d6bSPeter Brune k = (it - 1) % l; 200b8840d6bSPeter Brune ierr = VecCopy(D,dF[k]);CHKERRQ(ierr); 201b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 202b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 203b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 20418aa0c0cSPeter Brune if (qn->singlereduction) { 2051c6523dcSPeter Brune ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2061c6523dcSPeter Brune ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2071c6523dcSPeter Brune ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr); 2081c6523dcSPeter Brune ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2091c6523dcSPeter Brune ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2101c6523dcSPeter Brune ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr); 211b8840d6bSPeter Brune for (j = 0; j < l; j++) { 212b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 213b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 214b8840d6bSPeter Brune } 215b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 2161aa26658SKarl Rupp for (j = 0; j < l; j++) dXtdF[j] = H(j, j); 217b8840d6bSPeter Brune } else { 218b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 21923d44fbcSPeter Brune } 22023d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 221b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 222b8840d6bSPeter Brune } 223b8840d6bSPeter Brune } 224b8840d6bSPeter Brune 2254b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2264b11644fSPeter Brune for (i=0; i<l; i++) { 227b21d5a53SPeter Brune k = (it-i-1)%l; 22818aa0c0cSPeter Brune if (qn->singlereduction) { 2298e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2308e231d97SPeter Brune t = YtdX[k]; 2318e231d97SPeter Brune for (j=0; j<i; j++) { 2328e231d97SPeter Brune g = (it-j-1)%l; 233b3cfccb2SPeter Brune t -= alpha[g]*H(k, g); 2348e231d97SPeter Brune } 2358e231d97SPeter Brune alpha[k] = t/H(k,k); 2368e231d97SPeter Brune } else { 2373af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2388e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2398e231d97SPeter Brune } 24044f7e39eSPeter Brune if (qn->monitor) { 2413af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2425ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2433af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 24444f7e39eSPeter Brune } 2456bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2464b11644fSPeter Brune } 2474b11644fSPeter Brune 2480c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 249d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 250422a814eSBarry Smith SNESCheckKSPSolve(snes); 2510ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2520ecc9e1dSPeter Brune snes->linear_its += lits; 253b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2540ecc9e1dSPeter Brune } else { 255b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2560ecc9e1dSPeter Brune } 25718aa0c0cSPeter Brune if (qn->singlereduction) { 258b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 2598e231d97SPeter Brune } 2604b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 261bd052dfeSPeter Brune for (i = 0; i < l; i++) { 262b21d5a53SPeter Brune k = (it + i - l) % l; 26318aa0c0cSPeter Brune if (qn->singlereduction) { 2648e231d97SPeter Brune t = YtdX[k]; 2658e231d97SPeter Brune for (j = 0; j < i; j++) { 2668e231d97SPeter Brune g = (it + j - l) % l; 267b3cfccb2SPeter Brune t += (alpha[g] - beta[g])*H(g, k); 2688e231d97SPeter Brune } 2698e231d97SPeter Brune beta[k] = t / H(k, k); 2708e231d97SPeter Brune } else { 2716bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 2728e231d97SPeter Brune beta[k] = t / dXtdF[k]; 2738e231d97SPeter Brune } 27422d28d08SBarry Smith ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 27544f7e39eSPeter Brune if (qn->monitor) { 2763af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2775ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 2783af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 27944f7e39eSPeter Brune } 2804b11644fSPeter Brune } 2814b11644fSPeter Brune PetscFunctionReturn(0); 2824b11644fSPeter Brune } 2834b11644fSPeter Brune 2844b11644fSPeter Brune #undef __FUNCT__ 2854b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 2864b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 2874b11644fSPeter Brune { 2884b11644fSPeter Brune PetscErrorCode ierr; 2899f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 29015f5eeeaSPeter Brune Vec X,Xold; 2910a3905e1SPeter Brune Vec F,W; 29247f26062SPeter Brune Vec Y,D,Dold; 293b8840d6bSPeter Brune PetscInt i, i_r; 294335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 295422a814eSBarry Smith SNESLineSearchReason lssucceed; 296422a814eSBarry Smith PetscBool powell,periodic; 2971c6523dcSPeter Brune PetscScalar DolddotD,DolddotDold; 298b7281c8aSPeter Brune SNESConvergedReason reason; 2990ecc9e1dSPeter Brune 30084c577daSBarry Smith /* basically just a regular newton's method except for the application of the Jacobian */ 3014b11644fSPeter Brune 3026e111a19SKarl Rupp PetscFunctionBegin; 303c579b300SPatrick Farrell 304c579b300SPatrick Farrell if (snes->xl || snes->xu || snes->ops->computevariablebounds) { 305c579b300SPatrick Farrell SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 306c579b300SPatrick Farrell } 307c579b300SPatrick Farrell 308fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 3099f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3103af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 3110a3905e1SPeter Brune W = snes->work[3]; 312b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 313335fb975SPeter Brune Xold = snes->work[0]; 3143af51624SPeter Brune 3153af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 316335fb975SPeter Brune D = snes->work[1]; 317335fb975SPeter Brune Dold = snes->work[2]; 3184b11644fSPeter Brune 3194b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3204b11644fSPeter Brune 321e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3224b11644fSPeter Brune snes->iter = 0; 3234b11644fSPeter Brune snes->norm = 0.; 324e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 32547f26062SPeter Brune 326b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 327be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 328b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 329b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 330b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 331b7281c8aSPeter Brune PetscFunctionReturn(0); 332b7281c8aSPeter Brune } 33347f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 33447f26062SPeter Brune } else { 335e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 33615f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3371aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 338c1c75074SPeter Brune 33947f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 340422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 34147f26062SPeter Brune } 342b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 343be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr); 344b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 345b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 346b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 347b7281c8aSPeter Brune PetscFunctionReturn(0); 348b7281c8aSPeter Brune } 34947f26062SPeter Brune } else { 35047f26062SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 35147f26062SPeter Brune } 352b28a06ddSPeter Brune 353e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3544b11644fSPeter Brune snes->norm = fnorm; 355e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 356a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 3574b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3584b11644fSPeter Brune 3594b11644fSPeter Brune /* test convergence */ 3604b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3614b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 36270d3b23bSPeter Brune 3633cf07b75SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 3643cf07b75SPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 3653cf07b75SPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 3663cf07b75SPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 367ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 368ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 369ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 370ddd40ce5SPeter Brune PetscFunctionReturn(0); 371ddd40ce5SPeter Brune } 372be95d8f1SBarry Smith ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 3733cf07b75SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 3743cf07b75SPeter Brune } 3753cf07b75SPeter Brune 376f8e15203SPeter Brune /* scale the initial update */ 3770c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 378d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 3790ecc9e1dSPeter Brune } 3800ecc9e1dSPeter Brune 3815ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 3820a3905e1SPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) { 3830a3905e1SPeter Brune PetscScalar ff,xf; 3840a3905e1SPeter Brune ierr = VecCopy(Dold,Y);CHKERRQ(ierr); 3850a3905e1SPeter Brune ierr = VecCopy(Xold,W);CHKERRQ(ierr); 3860a3905e1SPeter Brune ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr); 3870a3905e1SPeter Brune ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 3880a3905e1SPeter Brune ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr); 3890a3905e1SPeter Brune ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr); 3900a3905e1SPeter Brune ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr); 3910a3905e1SPeter Brune ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr); 3920a3905e1SPeter Brune qn->scaling = PetscRealPart(xf)/PetscRealPart(ff); 3930a3905e1SPeter Brune } 394b8840d6bSPeter Brune switch (qn->type) { 395b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 396b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 397b8840d6bSPeter Brune break; 398b8840d6bSPeter Brune case SNES_QN_BROYDEN: 3995051edb1SPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr); 400b8840d6bSPeter Brune break; 401b8840d6bSPeter Brune case SNES_QN_LBFGS: 402b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 403b8840d6bSPeter Brune break; 404b8840d6bSPeter Brune } 40570d3b23bSPeter Brune /* line search for lambda */ 40670d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 4070a3905e1SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 40815f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 409f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4109f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 411422a814eSBarry Smith ierr = SNESLineSearchGetReason(snes->linesearch, &lssucceed);CHKERRQ(ierr); 412422a814eSBarry Smith ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 413422a814eSBarry Smith if (lssucceed) { 4149f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4159f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4169f3a0142SPeter Brune break; 4179f3a0142SPeter Brune } 4189f3a0142SPeter Brune } 4190c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 42005ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4210ecc9e1dSPeter Brune } 4223af51624SPeter Brune 42388d374b2SPeter Brune /* convergence monitoring */ 424b21d5a53SPeter 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); 425b21d5a53SPeter Brune 426b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 427b28a06ddSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 428b28a06ddSPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 429b28a06ddSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 430ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 431ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 432ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 433ddd40ce5SPeter Brune PetscFunctionReturn(0); 434ddd40ce5SPeter Brune } 435be95d8f1SBarry Smith ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 436b28a06ddSPeter Brune } 437b28a06ddSPeter Brune 438360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 43971dbe336SPeter Brune snes->norm = fnorm; 440360c497dSPeter Brune 441a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 4428409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4434b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 444d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4454b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 446b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 447be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr); 448b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 449b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 450b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 451b7281c8aSPeter Brune PetscFunctionReturn(0); 452b7281c8aSPeter Brune } 45388d374b2SPeter Brune } else { 45488d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 45588d374b2SPeter Brune } 4560c777b0cSPeter Brune powell = PETSC_FALSE; 4570c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 4586bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 45923c918e7SPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 46023c918e7SPeter Brune ierr = MatMult(snes->jacobian_pre,Dold,W);CHKERRQ(ierr); 46123c918e7SPeter Brune } else { 46223c918e7SPeter Brune ierr = VecCopy(Dold,W);CHKERRQ(ierr); 46323c918e7SPeter Brune } 46423c918e7SPeter Brune ierr = VecDotBegin(W, Dold, &DolddotDold);CHKERRQ(ierr); 46523c918e7SPeter Brune ierr = VecDotBegin(W, D, &DolddotD);CHKERRQ(ierr); 46623c918e7SPeter Brune ierr = VecDotEnd(W, Dold, &DolddotDold);CHKERRQ(ierr); 46723c918e7SPeter Brune ierr = VecDotEnd(W, D, &DolddotD);CHKERRQ(ierr); 4680c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 4690c777b0cSPeter Brune } 4700c777b0cSPeter Brune periodic = PETSC_FALSE; 471b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 472b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 4730c777b0cSPeter Brune } 4740c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 4750c777b0cSPeter Brune if (powell || periodic) { 4765ba6227bSPeter Brune if (qn->monitor) { 4775ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 478516fe3c3SPeter 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); 4795ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 4805ba6227bSPeter Brune } 4815ba6227bSPeter Brune i_r = -1; 4825ba6227bSPeter Brune /* general purpose update */ 4835ba6227bSPeter Brune if (snes->ops->update) { 4845ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 4855ba6227bSPeter Brune } 4860c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 487d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 4880ecc9e1dSPeter Brune } 4890ecc9e1dSPeter Brune } 49070d3b23bSPeter Brune /* general purpose update */ 49170d3b23bSPeter Brune if (snes->ops->update) { 49270d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 49370d3b23bSPeter Brune } 4945ba6227bSPeter Brune } 4954b11644fSPeter Brune if (i == snes->max_its) { 4964b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 4974b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 4984b11644fSPeter Brune } 4994b11644fSPeter Brune PetscFunctionReturn(0); 5004b11644fSPeter Brune } 5014b11644fSPeter Brune 5024b11644fSPeter Brune #undef __FUNCT__ 5034b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5044b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5054b11644fSPeter Brune { 5069f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5074b11644fSPeter Brune PetscErrorCode ierr; 508*fced5a79SAsbjørn Nilsen Riseth DM dm; 509335fb975SPeter Brune 5104b11644fSPeter Brune PetscFunctionBegin; 511*fced5a79SAsbjørn Nilsen Riseth 512*fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) { 513*fced5a79SAsbjørn Nilsen Riseth ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 514*fced5a79SAsbjørn Nilsen Riseth ierr = DMCreateGlobalVector(dm,&snes->vec_sol);CHKERRQ(ierr); 515*fced5a79SAsbjørn Nilsen Riseth } 516*fced5a79SAsbjørn Nilsen Riseth 517b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 51859e7931aSPeter Brune if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 519dcca6d9dSJed Brown ierr = PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);CHKERRQ(ierr); 5208e231d97SPeter Brune 52118aa0c0cSPeter Brune if (qn->singlereduction) { 522dcca6d9dSJed Brown ierr = PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);CHKERRQ(ierr); 5238e231d97SPeter Brune } 524fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 52560c08b40SPeter Brune /* set method defaults */ 5261efc8c45SPeter Brune if (qn->scale_type == SNES_QN_SCALE_DEFAULT) { 52760c08b40SPeter Brune if (qn->type == SNES_QN_BADBROYDEN) { 52860c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_NONE; 52960c08b40SPeter Brune } else { 53060c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 53160c08b40SPeter Brune } 53260c08b40SPeter Brune } 5331efc8c45SPeter Brune if (qn->restart_type == SNES_QN_RESTART_DEFAULT) { 53460c08b40SPeter Brune if (qn->type == SNES_QN_LBFGS) { 53560c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 53660c08b40SPeter Brune } else { 53760c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_PERIODIC; 53860c08b40SPeter Brune } 53960c08b40SPeter Brune } 54060c08b40SPeter Brune 5410c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5428e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5438e231d97SPeter Brune } 5446c67d002SPeter Brune if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;} 5454b11644fSPeter Brune PetscFunctionReturn(0); 5464b11644fSPeter Brune } 5474b11644fSPeter Brune 5484b11644fSPeter Brune #undef __FUNCT__ 5494b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5504b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5514b11644fSPeter Brune { 5524b11644fSPeter Brune PetscErrorCode ierr; 5539f83bee8SJed Brown SNES_QN *qn; 5540adebc6cSBarry Smith 5554b11644fSPeter Brune PetscFunctionBegin; 5564b11644fSPeter Brune if (snes->data) { 5579f83bee8SJed Brown qn = (SNES_QN*)snes->data; 558b8840d6bSPeter Brune if (qn->U) { 559b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5604b11644fSPeter Brune } 561b8840d6bSPeter Brune if (qn->V) { 562b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5634b11644fSPeter Brune } 56418aa0c0cSPeter Brune if (qn->singlereduction) { 5658e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5668e231d97SPeter Brune } 5675c7a0a03SPeter Brune ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr); 5684b11644fSPeter Brune } 5694b11644fSPeter Brune PetscFunctionReturn(0); 5704b11644fSPeter Brune } 5714b11644fSPeter Brune 5724b11644fSPeter Brune #undef __FUNCT__ 5734b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5744b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5754b11644fSPeter Brune { 5764b11644fSPeter Brune PetscErrorCode ierr; 5776e111a19SKarl Rupp 5784b11644fSPeter Brune PetscFunctionBegin; 5794b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 5804b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 581bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr); 5824b11644fSPeter Brune PetscFunctionReturn(0); 5834b11644fSPeter Brune } 5844b11644fSPeter Brune 5854b11644fSPeter Brune #undef __FUNCT__ 5864b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 5878c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_QN(PetscOptions *PetscOptionsObject,SNES snes) 5884b11644fSPeter Brune { 5894b11644fSPeter Brune 5904b11644fSPeter Brune PetscErrorCode ierr; 5912150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 59288f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 593f1c6b773SPeter Brune SNESLineSearch linesearch; 5942150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 5952150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 5961efc8c45SPeter Brune SNESQNType qtype = qn->type; 5972150357eSBarry Smith 5984b11644fSPeter Brune PetscFunctionBegin; 599e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES QN options");CHKERRQ(ierr); 6000298fd71SBarry Smith ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr); 6010298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr); 6020298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr); 6030298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr); 6040298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr); 60588f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 60688f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 60788f769c5SPeter Brune 60888f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 60988f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 61088f769c5SPeter Brune 6110fdccdaeSBarry Smith ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg);CHKERRQ(ierr); 6121efc8c45SPeter Brune if (flg) {ierr = SNESQNSetType(snes,qtype);CHKERRQ(ierr);} 6134b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6149e764e56SPeter Brune if (!snes->linesearch) { 6157601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 61660c08b40SPeter Brune if (qn->type == SNES_QN_LBFGS) { 6171a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 61860c08b40SPeter Brune } else if (qn->type == SNES_QN_BROYDEN) { 61960c08b40SPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);CHKERRQ(ierr); 62060c08b40SPeter Brune } else { 62160c08b40SPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 62260c08b40SPeter Brune } 6239e764e56SPeter Brune } 62444f7e39eSPeter Brune if (monflg) { 625ce94432eSBarry Smith qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 62644f7e39eSPeter Brune } 6274b11644fSPeter Brune PetscFunctionReturn(0); 6284b11644fSPeter Brune } 6294b11644fSPeter Brune 6305cd83419SPeter Brune #undef __FUNCT__ 6315cd83419SPeter Brune #define __FUNCT__ "SNESView_QN" 6325cd83419SPeter Brune static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer) 6335cd83419SPeter Brune { 6345cd83419SPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 6355cd83419SPeter Brune PetscBool iascii; 6365cd83419SPeter Brune PetscErrorCode ierr; 6375cd83419SPeter Brune 6385cd83419SPeter Brune PetscFunctionBegin; 6395cd83419SPeter Brune ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 6405cd83419SPeter Brune if (iascii) { 6415cd83419SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," QN type is %s, restart type is %s, scale type is %s\n",SNESQNTypes[qn->type],SNESQNRestartTypes[qn->restart_type],SNESQNScaleTypes[qn->scale_type]);CHKERRQ(ierr); 6425cd83419SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Stored subspace size: %d\n", qn->m);CHKERRQ(ierr); 6435cd83419SPeter Brune if (qn->singlereduction) { 6445cd83419SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using the single reduction variant.\n");CHKERRQ(ierr); 6455cd83419SPeter Brune } 6465cd83419SPeter Brune } 6475cd83419SPeter Brune PetscFunctionReturn(0); 6485cd83419SPeter Brune } 64992c02d66SPeter Brune 6500c777b0cSPeter Brune #undef __FUNCT__ 6510c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6520c777b0cSPeter Brune /*@ 6530c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6540c777b0cSPeter Brune 6550c777b0cSPeter Brune Logically Collective on SNES 6560c777b0cSPeter Brune 6570c777b0cSPeter Brune Input Parameters: 6580c777b0cSPeter Brune + snes - the iterative context 6590c777b0cSPeter Brune - rtype - restart type 6600c777b0cSPeter Brune 6610c777b0cSPeter Brune Options Database: 6620c777b0cSPeter Brune + -snes_qn_restart_type <powell,periodic,none> - set the restart type 66384c577daSBarry Smith - -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic 6640c777b0cSPeter Brune 6650c777b0cSPeter Brune Level: intermediate 6660c777b0cSPeter Brune 6670c777b0cSPeter Brune SNESQNRestartTypes: 6680c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6690c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6700c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6710c777b0cSPeter Brune 67284c577daSBarry Smith .keywords: SNES, SNESQN, restart, type, set SNESLineSearch, SNESQNRestartType 6730c777b0cSPeter Brune @*/ 6742150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6752150357eSBarry Smith { 6760c777b0cSPeter Brune PetscErrorCode ierr; 6776e111a19SKarl Rupp 6780c777b0cSPeter Brune PetscFunctionBegin; 6790c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6800c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6810c777b0cSPeter Brune PetscFunctionReturn(0); 6820c777b0cSPeter Brune } 6830c777b0cSPeter Brune 6840c777b0cSPeter Brune #undef __FUNCT__ 6850c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6860c777b0cSPeter Brune /*@ 68784c577daSBarry Smith SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN. 6880c777b0cSPeter Brune 6890c777b0cSPeter Brune Logically Collective on SNES 6900c777b0cSPeter Brune 6910c777b0cSPeter Brune Input Parameters: 6920c777b0cSPeter Brune + snes - the iterative context 6930c777b0cSPeter Brune - stype - scale type 6940c777b0cSPeter Brune 6950c777b0cSPeter Brune Options Database: 6960c777b0cSPeter Brune . -snes_qn_scale_type <shanno,none,linesearch,jacobian> 6970c777b0cSPeter Brune 6980c777b0cSPeter Brune Level: intermediate 6990c777b0cSPeter Brune 70084c577daSBarry Smith SNESQNScaleTypes: 7010c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 7020c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 7030c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 7040c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 7050c777b0cSPeter Brune 70684c577daSBarry Smith .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch, SNESQNScaleType 7070c777b0cSPeter Brune @*/ 7080c777b0cSPeter Brune 7092150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 7102150357eSBarry Smith { 7110c777b0cSPeter Brune PetscErrorCode ierr; 7126e111a19SKarl Rupp 7130c777b0cSPeter Brune PetscFunctionBegin; 7140c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7150c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 7160c777b0cSPeter Brune PetscFunctionReturn(0); 7170c777b0cSPeter Brune } 7180c777b0cSPeter Brune 7190c777b0cSPeter Brune #undef __FUNCT__ 7200c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 7210adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 7220adebc6cSBarry Smith { 7230c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7246e111a19SKarl Rupp 7250c777b0cSPeter Brune PetscFunctionBegin; 7260c777b0cSPeter Brune qn->scale_type = stype; 7270c777b0cSPeter Brune PetscFunctionReturn(0); 7280c777b0cSPeter Brune } 7290c777b0cSPeter Brune 7300c777b0cSPeter Brune #undef __FUNCT__ 7310c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7320adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 7330adebc6cSBarry Smith { 7340c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7356e111a19SKarl Rupp 7360c777b0cSPeter Brune PetscFunctionBegin; 7370c777b0cSPeter Brune qn->restart_type = rtype; 7380c777b0cSPeter Brune PetscFunctionReturn(0); 7390c777b0cSPeter Brune } 7400c777b0cSPeter Brune 7411efc8c45SPeter Brune #undef __FUNCT__ 7421efc8c45SPeter Brune #define __FUNCT__ "SNESQNSetType" 7431efc8c45SPeter Brune /*@ 7441efc8c45SPeter Brune SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN. 7451efc8c45SPeter Brune 7461efc8c45SPeter Brune Logically Collective on SNES 7471efc8c45SPeter Brune 7481efc8c45SPeter Brune Input Parameters: 7491efc8c45SPeter Brune + snes - the iterative context 7501efc8c45SPeter Brune - qtype - variant type 7511efc8c45SPeter Brune 7521efc8c45SPeter Brune Options Database: 75384c577daSBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> 7541efc8c45SPeter Brune 7551efc8c45SPeter Brune Level: beginner 7561efc8c45SPeter Brune 7571efc8c45SPeter Brune SNESQNTypes: 7581efc8c45SPeter Brune + SNES_QN_LBFGS - LBFGS variant 7591efc8c45SPeter Brune . SNES_QN_BROYDEN - Broyden variant 7601efc8c45SPeter Brune - SNES_QN_BADBROYDEN - Bad Broyden variant 7611efc8c45SPeter Brune 76284c577daSBarry Smith .keywords: SNES, SNESQN, type, set, SNESQNType 7631efc8c45SPeter Brune @*/ 7641efc8c45SPeter Brune 7651efc8c45SPeter Brune PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype) 7661efc8c45SPeter Brune { 7671efc8c45SPeter Brune PetscErrorCode ierr; 7681efc8c45SPeter Brune 7691efc8c45SPeter Brune PetscFunctionBegin; 7701efc8c45SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7711efc8c45SPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));CHKERRQ(ierr); 7721efc8c45SPeter Brune PetscFunctionReturn(0); 7731efc8c45SPeter Brune } 7741efc8c45SPeter Brune 7751efc8c45SPeter Brune #undef __FUNCT__ 7761efc8c45SPeter Brune #define __FUNCT__ "SNESQNSetType_QN" 7771efc8c45SPeter Brune PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype) 7781efc8c45SPeter Brune { 7791efc8c45SPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7801efc8c45SPeter Brune 7811efc8c45SPeter Brune PetscFunctionBegin; 7821efc8c45SPeter Brune qn->type = qtype; 7831efc8c45SPeter Brune PetscFunctionReturn(0); 7841efc8c45SPeter Brune } 7851efc8c45SPeter Brune 7864b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 7874b11644fSPeter Brune /*MC 7884b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 7894b11644fSPeter Brune 7906cc8130cSPeter Brune Options Database: 7916cc8130cSPeter Brune 79284c577daSBarry Smith + -snes_qn_m <m> - Number of past states saved for the L-Broyden methods. 79384c577daSBarry Smith + -snes_qn_restart_type <powell,periodic,none> - set the restart type 7941867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 7951867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 79684c577daSBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> - QN type 79784c577daSBarry Smith . -snes_qn_scale_type <shanno,none,linesearch,jacobian> - scaling performed on inner Jacobian 79872835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 79984c577daSBarry Smith - -snes_qn_monitor - Monitors the quasi-newton Jacobian. 8006cc8130cSPeter Brune 801b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 802b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 803b8840d6bSPeter Brune updates. 8046cc8130cSPeter Brune 8051867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 8061867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 80784c577daSBarry Smith iteration as the current iteration's values when constructing the approximate Jacobian. The second, composed, 80884c577daSBarry Smith perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 8091867fe5bSPeter Brune 8106cc8130cSPeter Brune References: 8116cc8130cSPeter Brune 8120a8e8854SPeter Brune Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995. 8136cc8130cSPeter Brune 8140a8e8854SPeter Brune R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi-Newton Matrices and their use in Limited Memory Methods, 8150a8e8854SPeter Brune Technical Report, Northwestern University, June 1992. 8160a8e8854SPeter Brune 8170a8e8854SPeter Brune Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE 8180a8e8854SPeter Brune Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 819b8840d6bSPeter Brune 8204b11644fSPeter Brune 8214b11644fSPeter Brune Level: beginner 8224b11644fSPeter Brune 82304d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 8246cc8130cSPeter Brune 8254b11644fSPeter Brune M*/ 8264b11644fSPeter Brune #undef __FUNCT__ 8274b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 8288cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 8294b11644fSPeter Brune { 8304b11644fSPeter Brune PetscErrorCode ierr; 8319f83bee8SJed Brown SNES_QN *qn; 8324b11644fSPeter Brune 8334b11644fSPeter Brune PetscFunctionBegin; 8344b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 8354b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 8364b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 8374b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 8385cd83419SPeter Brune snes->ops->view = SNESView_QN; 8394b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 8404b11644fSPeter Brune 84147f26062SPeter Brune snes->pcside = PC_LEFT; 84247f26062SPeter Brune 84342f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 84442f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 84542f4f86dSBarry Smith 84688976e71SPeter Brune if (!snes->tolerancesset) { 8470e444f03SPeter Brune snes->max_funcs = 30000; 8480e444f03SPeter Brune snes->max_its = 10000; 84988976e71SPeter Brune } 8500e444f03SPeter Brune 851b00a9115SJed Brown ierr = PetscNewLog(snes,&qn);CHKERRQ(ierr); 8524b11644fSPeter Brune snes->data = (void*) qn; 8530ecc9e1dSPeter Brune qn->m = 10; 854b21d5a53SPeter Brune qn->scaling = 1.0; 8550298fd71SBarry Smith qn->U = NULL; 8560298fd71SBarry Smith qn->V = NULL; 8570298fd71SBarry Smith qn->dXtdF = NULL; 8580298fd71SBarry Smith qn->dFtdX = NULL; 8590298fd71SBarry Smith qn->dXdFmat = NULL; 8600298fd71SBarry Smith qn->monitor = NULL; 8611c6523dcSPeter Brune qn->singlereduction = PETSC_TRUE; 862b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 86360c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_DEFAULT; 86460c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_DEFAULT; 865b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 866ea630c6eSPeter Brune 867bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr); 868bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr); 8691efc8c45SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);CHKERRQ(ierr); 8704b11644fSPeter Brune PetscFunctionReturn(0); 8714b11644fSPeter Brune } 872