1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2fced5a79SAsbjø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 */ 22b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 23b8840d6bSPeter Brune SNESQNType type; /* the type of quasi-newton method used */ 2488f769c5SPeter Brune SNESQNScaleType scale_type; /* the type of scaling used */ 250c777b0cSPeter Brune SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 269f83bee8SJed Brown } SNES_QN; 274b11644fSPeter Brune 284b11644fSPeter Brune #undef __FUNCT__ 29b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden" 305051edb1SPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D) 312150357eSBarry Smith { 324b11644fSPeter Brune PetscErrorCode ierr; 339f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 34b8840d6bSPeter Brune Vec W = snes->work[3]; 35b8840d6bSPeter Brune Vec *U = qn->U; 36b8840d6bSPeter Brune PetscInt m = qn->m; 375c7a0a03SPeter Brune PetscInt k,i,j,lits,l = m; 385c7a0a03SPeter Brune PetscReal unorm,a,b; 395c7a0a03SPeter Brune PetscReal *lambda=qn->lambda; 40b8840d6bSPeter Brune PetscScalar gdot; 4159e7931aSPeter Brune PetscReal udot; 422150357eSBarry Smith 43b8840d6bSPeter Brune PetscFunctionBegin; 44b8840d6bSPeter Brune if (it < m) l = it; 455c7a0a03SPeter Brune if (it > 0) { 465c7a0a03SPeter Brune k = (it-1)%l; 475c7a0a03SPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr); 48a49fa0d8SPeter Brune ierr = VecCopy(Xold,U[k]);CHKERRQ(ierr); 495051edb1SPeter Brune ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr); 505c7a0a03SPeter Brune if (qn->monitor) { 515c7a0a03SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 52*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %D of %D by lambda: %14.12e \n",k,l,(double)lambda[k]);CHKERRQ(ierr); 535c7a0a03SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 545c7a0a03SPeter Brune } 555c7a0a03SPeter Brune } 56b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 57d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr); 58422a814eSBarry Smith SNESCheckKSPSolve(snes); 59b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 60b8840d6bSPeter Brune snes->linear_its += lits; 61b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 62b8840d6bSPeter Brune } else { 63b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 64b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 65b8840d6bSPeter Brune } 66b8840d6bSPeter Brune 67b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 68b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 695c7a0a03SPeter Brune j = (it+i-l)%l; 705c7a0a03SPeter Brune k = (it+i-l+1)%l; 715c7a0a03SPeter Brune ierr = VecNorm(U[j],NORM_2,&unorm);CHKERRQ(ierr); 725c7a0a03SPeter Brune ierr = VecDot(U[j],Y,&gdot);CHKERRQ(ierr); 7359e7931aSPeter Brune unorm *= unorm; 7459e7931aSPeter Brune udot = PetscRealPart(gdot); 755c7a0a03SPeter Brune a = (lambda[j]/lambda[k]); 765c7a0a03SPeter Brune b = -(1.-lambda[j]); 7759e7931aSPeter Brune a *= udot/unorm; 7859e7931aSPeter Brune b *= udot/unorm; 795c7a0a03SPeter Brune ierr = VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);CHKERRQ(ierr); 805c7a0a03SPeter Brune 81b8840d6bSPeter Brune if (qn->monitor) { 82b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 83*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %D and %D, gdot: %14.12e\n",k,j,(double)PetscRealPart(gdot));CHKERRQ(ierr); 84b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 85b8840d6bSPeter Brune } 86b8840d6bSPeter Brune } 87b8840d6bSPeter Brune if (l > 0) { 88b8840d6bSPeter Brune k = (it-1)%l; 89b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 905c7a0a03SPeter Brune ierr = VecNorm(U[k],NORM_2,&unorm);CHKERRQ(ierr); 9159e7931aSPeter Brune unorm *= unorm; 9259e7931aSPeter Brune udot = PetscRealPart(gdot); 9359e7931aSPeter Brune a = unorm/(unorm-lambda[k]*udot); 9459e7931aSPeter Brune b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot); 95b8840d6bSPeter Brune if (qn->monitor) { 96b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 97*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %D: a: %14.12e b: %14.12e \n",k,(double)a,(double)b);CHKERRQ(ierr); 98b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 99b8840d6bSPeter Brune } 1005c7a0a03SPeter Brune ierr = VecAXPBY(Y,b,a,U[k]);CHKERRQ(ierr); 101b8840d6bSPeter Brune } 1025c7a0a03SPeter Brune l = m; 1035c7a0a03SPeter Brune if (it+1<m)l=it+1; 1045c7a0a03SPeter Brune k = it%l; 1055c7a0a03SPeter Brune if (qn->monitor) { 1065c7a0a03SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 107*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "setting vector %D of %D\n",k,l);CHKERRQ(ierr); 1085c7a0a03SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1095c7a0a03SPeter Brune } 110b8840d6bSPeter Brune PetscFunctionReturn(0); 111b8840d6bSPeter Brune } 112b8840d6bSPeter Brune 113b8840d6bSPeter Brune #undef __FUNCT__ 114b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden" 1150adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1160adebc6cSBarry Smith { 117b8840d6bSPeter Brune PetscErrorCode ierr; 118b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 119b8840d6bSPeter Brune Vec W = snes->work[3]; 120b8840d6bSPeter Brune Vec *U = qn->U; 121b8840d6bSPeter Brune Vec *T = qn->V; 122b8840d6bSPeter Brune 12384c577daSBarry Smith /* ksp thing for Jacobian scaling */ 124c9e59058SPeter Brune PetscInt h,k,j,i,lits; 125b8840d6bSPeter Brune PetscInt m = qn->m; 1265051edb1SPeter Brune PetscScalar gdot,udot; 127b8840d6bSPeter Brune PetscInt l = m; 1280adebc6cSBarry Smith 129b8840d6bSPeter Brune PetscFunctionBegin; 130b8840d6bSPeter Brune if (it < m) l = it; 131b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 132b8840d6bSPeter Brune if (l > 0) { 133b8840d6bSPeter Brune k = (it-1)%l; 134c9e59058SPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);CHKERRQ(ierr); 135b8840d6bSPeter Brune ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 136b8840d6bSPeter Brune ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 1375051edb1SPeter Brune ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr); 1385051edb1SPeter Brune ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr); 139b8840d6bSPeter Brune } 140b8840d6bSPeter Brune 141b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 142d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 143422a814eSBarry Smith SNESCheckKSPSolve(snes); 144b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 145b8840d6bSPeter Brune snes->linear_its += lits; 146b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 147b8840d6bSPeter Brune } else { 148b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 149b8840d6bSPeter Brune } 1505051edb1SPeter Brune 1515051edb1SPeter Brune /* inward recursion starting at the first update and working forward */ 1525051edb1SPeter Brune if (l > 0) { 1535051edb1SPeter Brune for (i = 0; i < l-1; i++) { 154c9e59058SPeter Brune j = (it+i-l)%l; 155c9e59058SPeter Brune k = (it+i-l+1)%l; 156c9e59058SPeter Brune h = (it-1)%l; 157c9e59058SPeter Brune ierr = VecDotBegin(U[j],U[h],&gdot);CHKERRQ(ierr); 158c9e59058SPeter Brune ierr = VecDotBegin(U[j],U[j],&udot);CHKERRQ(ierr); 159c9e59058SPeter Brune ierr = VecDotEnd(U[j],U[h],&gdot);CHKERRQ(ierr); 160c9e59058SPeter Brune ierr = VecDotEnd(U[j],U[j],&udot);CHKERRQ(ierr); 1615051edb1SPeter Brune ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr); 162c9e59058SPeter Brune ierr = VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);CHKERRQ(ierr); 1635051edb1SPeter Brune if (qn->monitor) { 1645051edb1SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 165*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D gdot: %14.12e\n", it, k, (double)PetscRealPart(gdot));CHKERRQ(ierr); 1665051edb1SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1675051edb1SPeter Brune } 1685051edb1SPeter Brune } 1695051edb1SPeter Brune } 170b8840d6bSPeter Brune PetscFunctionReturn(0); 171b8840d6bSPeter Brune } 172b8840d6bSPeter Brune 173b8840d6bSPeter Brune #undef __FUNCT__ 174b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 1750adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1760adebc6cSBarry Smith { 177b8840d6bSPeter Brune PetscErrorCode ierr; 178b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 179b8840d6bSPeter Brune Vec W = snes->work[3]; 180b8840d6bSPeter Brune Vec *dX = qn->U; 181b8840d6bSPeter Brune Vec *dF = qn->V; 182bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 183bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 1848e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 185b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 1868e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 187bd052dfeSPeter Brune 18884c577daSBarry Smith /* ksp thing for Jacobian scaling */ 1898e231d97SPeter Brune PetscInt k,i,j,g,lits; 1904b11644fSPeter Brune PetscInt m = qn->m; 1914b11644fSPeter Brune PetscScalar t; 1924b11644fSPeter Brune PetscInt l = m; 1938e231d97SPeter Brune 1944b11644fSPeter Brune PetscFunctionBegin; 1955ba6227bSPeter Brune if (it < m) l = it; 1961c6523dcSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 197b8840d6bSPeter Brune if (it > 0) { 198b8840d6bSPeter Brune k = (it - 1) % l; 199b8840d6bSPeter Brune ierr = VecCopy(D,dF[k]);CHKERRQ(ierr); 200b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 201b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 202b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 20318aa0c0cSPeter Brune if (qn->singlereduction) { 2041c6523dcSPeter Brune ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2051c6523dcSPeter Brune ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2061c6523dcSPeter Brune ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr); 2071c6523dcSPeter Brune ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2081c6523dcSPeter Brune ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2091c6523dcSPeter Brune ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr); 210b8840d6bSPeter Brune for (j = 0; j < l; j++) { 211b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 212b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 213b8840d6bSPeter Brune } 214b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 2151aa26658SKarl Rupp for (j = 0; j < l; j++) dXtdF[j] = H(j, j); 216b8840d6bSPeter Brune } else { 217b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 21823d44fbcSPeter Brune } 21923d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 220b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 221b8840d6bSPeter Brune } 222b8840d6bSPeter Brune } 223b8840d6bSPeter Brune 2244b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2254b11644fSPeter Brune for (i=0; i<l; i++) { 226b21d5a53SPeter Brune k = (it-i-1)%l; 22718aa0c0cSPeter Brune if (qn->singlereduction) { 2288e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2298e231d97SPeter Brune t = YtdX[k]; 2308e231d97SPeter Brune for (j=0; j<i; j++) { 2318e231d97SPeter Brune g = (it-j-1)%l; 232b3cfccb2SPeter Brune t -= alpha[g]*H(k, g); 2338e231d97SPeter Brune } 2348e231d97SPeter Brune alpha[k] = t/H(k,k); 2358e231d97SPeter Brune } else { 2363af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2378e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2388e231d97SPeter Brune } 23944f7e39eSPeter Brune if (qn->monitor) { 2403af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 241*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D alpha: %14.12e\n", it, k, (double)PetscRealPart(alpha[k]));CHKERRQ(ierr); 2423af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 24344f7e39eSPeter Brune } 2446bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2454b11644fSPeter Brune } 2464b11644fSPeter Brune 2470c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 248d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 249422a814eSBarry Smith SNESCheckKSPSolve(snes); 2500ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2510ecc9e1dSPeter Brune snes->linear_its += lits; 252b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2530ecc9e1dSPeter Brune } else { 254b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2550ecc9e1dSPeter Brune } 25618aa0c0cSPeter Brune if (qn->singlereduction) { 257b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 2588e231d97SPeter Brune } 2594b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 260bd052dfeSPeter Brune for (i = 0; i < l; i++) { 261b21d5a53SPeter Brune k = (it + i - l) % l; 26218aa0c0cSPeter Brune if (qn->singlereduction) { 2638e231d97SPeter Brune t = YtdX[k]; 2648e231d97SPeter Brune for (j = 0; j < i; j++) { 2658e231d97SPeter Brune g = (it + j - l) % l; 266b3cfccb2SPeter Brune t += (alpha[g] - beta[g])*H(g, k); 2678e231d97SPeter Brune } 2688e231d97SPeter Brune beta[k] = t / H(k, k); 2698e231d97SPeter Brune } else { 2706bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 2718e231d97SPeter Brune beta[k] = t / dXtdF[k]; 2728e231d97SPeter Brune } 27322d28d08SBarry Smith ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 27444f7e39eSPeter Brune if (qn->monitor) { 2753af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 276*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D alpha - beta: %14.12e\n", it, k, (double)PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 2773af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 27844f7e39eSPeter Brune } 2794b11644fSPeter Brune } 2804b11644fSPeter Brune PetscFunctionReturn(0); 2814b11644fSPeter Brune } 2824b11644fSPeter Brune 2834b11644fSPeter Brune #undef __FUNCT__ 2844b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 2854b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 2864b11644fSPeter Brune { 2874b11644fSPeter Brune PetscErrorCode ierr; 2889f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 28915f5eeeaSPeter Brune Vec X,Xold; 2900a3905e1SPeter Brune Vec F,W; 29147f26062SPeter Brune Vec Y,D,Dold; 292b8840d6bSPeter Brune PetscInt i, i_r; 293335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 294422a814eSBarry Smith SNESLineSearchReason lssucceed; 295422a814eSBarry Smith PetscBool powell,periodic; 2961c6523dcSPeter Brune PetscScalar DolddotD,DolddotDold; 297b7281c8aSPeter Brune SNESConvergedReason reason; 2980ecc9e1dSPeter Brune 29984c577daSBarry Smith /* basically just a regular newton's method except for the application of the Jacobian */ 3004b11644fSPeter Brune 3016e111a19SKarl Rupp PetscFunctionBegin; 3026c4ed002SBarry Smith if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 303c579b300SPatrick Farrell 304fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 3059f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3063af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 3070a3905e1SPeter Brune W = snes->work[3]; 308b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 309335fb975SPeter Brune Xold = snes->work[0]; 3103af51624SPeter Brune 3113af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 312335fb975SPeter Brune D = snes->work[1]; 313335fb975SPeter Brune Dold = snes->work[2]; 3144b11644fSPeter Brune 3154b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3164b11644fSPeter Brune 317e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3184b11644fSPeter Brune snes->iter = 0; 3194b11644fSPeter Brune snes->norm = 0.; 320e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 32147f26062SPeter Brune 322b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 323be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 324b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 325b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 326b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 327b7281c8aSPeter Brune PetscFunctionReturn(0); 328b7281c8aSPeter Brune } 32947f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 33047f26062SPeter Brune } else { 331e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 33215f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3331aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 334c1c75074SPeter Brune 33547f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 336422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 33747f26062SPeter Brune } 338b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 339be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr); 340b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 341b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 342b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 343b7281c8aSPeter Brune PetscFunctionReturn(0); 344b7281c8aSPeter Brune } 34547f26062SPeter Brune } else { 34647f26062SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 34747f26062SPeter Brune } 348b28a06ddSPeter Brune 349e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3504b11644fSPeter Brune snes->norm = fnorm; 351e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 352a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 3534b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3544b11644fSPeter Brune 3554b11644fSPeter Brune /* test convergence */ 3564b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3574b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 35870d3b23bSPeter Brune 3593cf07b75SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 3603cf07b75SPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 3613cf07b75SPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 3623cf07b75SPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 363ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 364ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 365ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 366ddd40ce5SPeter Brune PetscFunctionReturn(0); 367ddd40ce5SPeter Brune } 368be95d8f1SBarry Smith ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 3693cf07b75SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 3703cf07b75SPeter Brune } 3713cf07b75SPeter Brune 372f8e15203SPeter Brune /* scale the initial update */ 3730c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 374d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 375aeb49b86SAsbjørn Nilsen Riseth ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 3760ecc9e1dSPeter Brune } 3770ecc9e1dSPeter Brune 3785ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 3790a3905e1SPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) { 3800a3905e1SPeter Brune PetscScalar ff,xf; 3810a3905e1SPeter Brune ierr = VecCopy(Dold,Y);CHKERRQ(ierr); 3820a3905e1SPeter Brune ierr = VecCopy(Xold,W);CHKERRQ(ierr); 3830a3905e1SPeter Brune ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr); 3840a3905e1SPeter Brune ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 3850a3905e1SPeter Brune ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr); 3860a3905e1SPeter Brune ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr); 3870a3905e1SPeter Brune ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr); 3880a3905e1SPeter Brune ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr); 3890a3905e1SPeter Brune qn->scaling = PetscRealPart(xf)/PetscRealPart(ff); 390*6bdcc836SBarry Smith if (qn->monitor) { 391*6bdcc836SBarry Smith ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 392*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "Shanno scaling %D %g\n", i,(double)qn->scaling);CHKERRQ(ierr); 393*6bdcc836SBarry Smith ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 394*6bdcc836SBarry Smith } 3950a3905e1SPeter Brune } 396b8840d6bSPeter Brune switch (qn->type) { 397b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 398b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 399b8840d6bSPeter Brune break; 400b8840d6bSPeter Brune case SNES_QN_BROYDEN: 4015051edb1SPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr); 402b8840d6bSPeter Brune break; 403b8840d6bSPeter Brune case SNES_QN_LBFGS: 404b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 405b8840d6bSPeter Brune break; 406b8840d6bSPeter Brune } 40770d3b23bSPeter Brune /* line search for lambda */ 40870d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 4090a3905e1SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 41015f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 411f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4129f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 413422a814eSBarry Smith ierr = SNESLineSearchGetReason(snes->linesearch, &lssucceed);CHKERRQ(ierr); 414422a814eSBarry Smith ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 415422a814eSBarry Smith if (lssucceed) { 4169f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4179f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4189f3a0142SPeter Brune break; 4199f3a0142SPeter Brune } 4209f3a0142SPeter Brune } 4210c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 42205ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4230ecc9e1dSPeter Brune } 4243af51624SPeter Brune 42588d374b2SPeter Brune /* convergence monitoring */ 426b21d5a53SPeter 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); 427b21d5a53SPeter Brune 428b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 429b28a06ddSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 430b28a06ddSPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 431b28a06ddSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 432ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 433ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 434ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 435ddd40ce5SPeter Brune PetscFunctionReturn(0); 436ddd40ce5SPeter Brune } 437be95d8f1SBarry Smith ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 438b28a06ddSPeter Brune } 439b28a06ddSPeter Brune 440360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 44171dbe336SPeter Brune snes->norm = fnorm; 442360c497dSPeter Brune 443a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 4448409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4454b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 446d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4474b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 448b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 449be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr); 450b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 451b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 452b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 453b7281c8aSPeter Brune PetscFunctionReturn(0); 454b7281c8aSPeter Brune } 45588d374b2SPeter Brune } else { 45688d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 45788d374b2SPeter Brune } 4580c777b0cSPeter Brune powell = PETSC_FALSE; 459*6bdcc836SBarry Smith if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) { 460*6bdcc836SBarry Smith /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */ 46123c918e7SPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 46223c918e7SPeter Brune ierr = MatMult(snes->jacobian_pre,Dold,W);CHKERRQ(ierr); 46323c918e7SPeter Brune } else { 46423c918e7SPeter Brune ierr = VecCopy(Dold,W);CHKERRQ(ierr); 46523c918e7SPeter Brune } 46623c918e7SPeter Brune ierr = VecDotBegin(W, Dold, &DolddotDold);CHKERRQ(ierr); 46723c918e7SPeter Brune ierr = VecDotBegin(W, D, &DolddotD);CHKERRQ(ierr); 46823c918e7SPeter Brune ierr = VecDotEnd(W, Dold, &DolddotDold);CHKERRQ(ierr); 46923c918e7SPeter Brune ierr = VecDotEnd(W, D, &DolddotD);CHKERRQ(ierr); 4700c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 4710c777b0cSPeter Brune } 4720c777b0cSPeter Brune periodic = PETSC_FALSE; 473b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 474b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 4750c777b0cSPeter Brune } 4760c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 4770c777b0cSPeter Brune if (powell || periodic) { 4785ba6227bSPeter Brune if (qn->monitor) { 4795ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 480*6bdcc836SBarry Smith if (powell) { 481*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "Powell restart! |%14.12e| > %6.4f*|%14.12e| i_r = %D\n", (double)PetscRealPart(DolddotD), (double)qn->powell_gamma, (double)PetscRealPart(DolddotDold),i_r);CHKERRQ(ierr); 482*6bdcc836SBarry Smith } else { 483*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %D\n", i_r);CHKERRQ(ierr); 484*6bdcc836SBarry Smith } 4855ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 4865ba6227bSPeter Brune } 4875ba6227bSPeter Brune i_r = -1; 4885ba6227bSPeter Brune /* general purpose update */ 4895ba6227bSPeter Brune if (snes->ops->update) { 4905ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 4915ba6227bSPeter Brune } 4920c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 493d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 4940ecc9e1dSPeter Brune } 4950ecc9e1dSPeter Brune } 49670d3b23bSPeter Brune /* general purpose update */ 49770d3b23bSPeter Brune if (snes->ops->update) { 49870d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 49970d3b23bSPeter Brune } 5005ba6227bSPeter Brune } 5014b11644fSPeter Brune if (i == snes->max_its) { 5024b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 5034b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 5044b11644fSPeter Brune } 5054b11644fSPeter Brune PetscFunctionReturn(0); 5064b11644fSPeter Brune } 5074b11644fSPeter Brune 5084b11644fSPeter Brune #undef __FUNCT__ 5094b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5104b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5114b11644fSPeter Brune { 5129f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5134b11644fSPeter Brune PetscErrorCode ierr; 514fced5a79SAsbjørn Nilsen Riseth DM dm; 515335fb975SPeter Brune 5164b11644fSPeter Brune PetscFunctionBegin; 517fced5a79SAsbjørn Nilsen Riseth 518fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) { 519fced5a79SAsbjørn Nilsen Riseth ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 520fced5a79SAsbjørn Nilsen Riseth ierr = DMCreateGlobalVector(dm,&snes->vec_sol);CHKERRQ(ierr); 521fced5a79SAsbjørn Nilsen Riseth } 522fced5a79SAsbjørn Nilsen Riseth 523b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 52459e7931aSPeter Brune if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 525dcca6d9dSJed Brown ierr = PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);CHKERRQ(ierr); 5268e231d97SPeter Brune 52718aa0c0cSPeter Brune if (qn->singlereduction) { 528dcca6d9dSJed Brown ierr = PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);CHKERRQ(ierr); 5298e231d97SPeter Brune } 530fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 53160c08b40SPeter Brune /* set method defaults */ 5321efc8c45SPeter Brune if (qn->scale_type == SNES_QN_SCALE_DEFAULT) { 53360c08b40SPeter Brune if (qn->type == SNES_QN_BADBROYDEN) { 53460c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_NONE; 53560c08b40SPeter Brune } else { 53660c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 53760c08b40SPeter Brune } 53860c08b40SPeter Brune } 5391efc8c45SPeter Brune if (qn->restart_type == SNES_QN_RESTART_DEFAULT) { 54060c08b40SPeter Brune if (qn->type == SNES_QN_LBFGS) { 54160c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 54260c08b40SPeter Brune } else { 54360c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_PERIODIC; 54460c08b40SPeter Brune } 54560c08b40SPeter Brune } 54660c08b40SPeter Brune 5470c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5488e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5498e231d97SPeter Brune } 5506c67d002SPeter Brune if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;} 5514b11644fSPeter Brune PetscFunctionReturn(0); 5524b11644fSPeter Brune } 5534b11644fSPeter Brune 5544b11644fSPeter Brune #undef __FUNCT__ 5554b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5564b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5574b11644fSPeter Brune { 5584b11644fSPeter Brune PetscErrorCode ierr; 5599f83bee8SJed Brown SNES_QN *qn; 5600adebc6cSBarry Smith 5614b11644fSPeter Brune PetscFunctionBegin; 5624b11644fSPeter Brune if (snes->data) { 5639f83bee8SJed Brown qn = (SNES_QN*)snes->data; 564b8840d6bSPeter Brune if (qn->U) { 565b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5664b11644fSPeter Brune } 567b8840d6bSPeter Brune if (qn->V) { 568b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5694b11644fSPeter Brune } 57018aa0c0cSPeter Brune if (qn->singlereduction) { 5718e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5728e231d97SPeter Brune } 5735c7a0a03SPeter Brune ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr); 5744b11644fSPeter Brune } 5754b11644fSPeter Brune PetscFunctionReturn(0); 5764b11644fSPeter Brune } 5774b11644fSPeter Brune 5784b11644fSPeter Brune #undef __FUNCT__ 5794b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5804b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5814b11644fSPeter Brune { 5824b11644fSPeter Brune PetscErrorCode ierr; 5836e111a19SKarl Rupp 5844b11644fSPeter Brune PetscFunctionBegin; 5854b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 5864b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 587bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr); 5884b11644fSPeter Brune PetscFunctionReturn(0); 5894b11644fSPeter Brune } 5904b11644fSPeter Brune 5914b11644fSPeter Brune #undef __FUNCT__ 5924b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 5934416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_QN(PetscOptionItems *PetscOptionsObject,SNES snes) 5944b11644fSPeter Brune { 5954b11644fSPeter Brune 5964b11644fSPeter Brune PetscErrorCode ierr; 5972150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 59888f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 599f1c6b773SPeter Brune SNESLineSearch linesearch; 6002150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 6012150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 6021efc8c45SPeter Brune SNESQNType qtype = qn->type; 6032150357eSBarry Smith 6044b11644fSPeter Brune PetscFunctionBegin; 605e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES QN options");CHKERRQ(ierr); 6060298fd71SBarry Smith ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr); 6070298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr); 6080298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr); 6090298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr); 61088f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 61188f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 61288f769c5SPeter Brune 61388f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 61488f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 61588f769c5SPeter Brune 6160fdccdaeSBarry Smith ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg);CHKERRQ(ierr); 6171efc8c45SPeter Brune if (flg) {ierr = SNESQNSetType(snes,qtype);CHKERRQ(ierr);} 6184b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6199e764e56SPeter Brune if (!snes->linesearch) { 6207601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 62160c08b40SPeter Brune if (qn->type == SNES_QN_LBFGS) { 6221a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 62360c08b40SPeter Brune } else if (qn->type == SNES_QN_BROYDEN) { 62460c08b40SPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);CHKERRQ(ierr); 62560c08b40SPeter Brune } else { 62660c08b40SPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 62760c08b40SPeter Brune } 6289e764e56SPeter Brune } 62944f7e39eSPeter Brune if (monflg) { 630ce94432eSBarry Smith qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 63144f7e39eSPeter Brune } 6324b11644fSPeter Brune PetscFunctionReturn(0); 6334b11644fSPeter Brune } 6344b11644fSPeter Brune 6355cd83419SPeter Brune #undef __FUNCT__ 6365cd83419SPeter Brune #define __FUNCT__ "SNESView_QN" 6375cd83419SPeter Brune static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer) 6385cd83419SPeter Brune { 6395cd83419SPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 6405cd83419SPeter Brune PetscBool iascii; 6415cd83419SPeter Brune PetscErrorCode ierr; 6425cd83419SPeter Brune 6435cd83419SPeter Brune PetscFunctionBegin; 6445cd83419SPeter Brune ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 6455cd83419SPeter Brune if (iascii) { 6465cd83419SPeter 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); 647*6bdcc836SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Stored subspace size: %D\n", qn->m);CHKERRQ(ierr); 6485cd83419SPeter Brune if (qn->singlereduction) { 6495cd83419SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using the single reduction variant.\n");CHKERRQ(ierr); 6505cd83419SPeter Brune } 6515cd83419SPeter Brune } 6525cd83419SPeter Brune PetscFunctionReturn(0); 6535cd83419SPeter Brune } 65492c02d66SPeter Brune 6550c777b0cSPeter Brune #undef __FUNCT__ 6560c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6570c777b0cSPeter Brune /*@ 6580c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6590c777b0cSPeter Brune 6600c777b0cSPeter Brune Logically Collective on SNES 6610c777b0cSPeter Brune 6620c777b0cSPeter Brune Input Parameters: 6630c777b0cSPeter Brune + snes - the iterative context 6640c777b0cSPeter Brune - rtype - restart type 6650c777b0cSPeter Brune 6660c777b0cSPeter Brune Options Database: 6670c777b0cSPeter Brune + -snes_qn_restart_type <powell,periodic,none> - set the restart type 66884c577daSBarry Smith - -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic 6690c777b0cSPeter Brune 6700c777b0cSPeter Brune Level: intermediate 6710c777b0cSPeter Brune 6720c777b0cSPeter Brune SNESQNRestartTypes: 6730c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6740c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6750c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6760c777b0cSPeter Brune 67784c577daSBarry Smith .keywords: SNES, SNESQN, restart, type, set SNESLineSearch, SNESQNRestartType 6780c777b0cSPeter Brune @*/ 6792150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6802150357eSBarry Smith { 6810c777b0cSPeter Brune PetscErrorCode ierr; 6826e111a19SKarl Rupp 6830c777b0cSPeter Brune PetscFunctionBegin; 6840c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6850c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6860c777b0cSPeter Brune PetscFunctionReturn(0); 6870c777b0cSPeter Brune } 6880c777b0cSPeter Brune 6890c777b0cSPeter Brune #undef __FUNCT__ 6900c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6910c777b0cSPeter Brune /*@ 69284c577daSBarry Smith SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN. 6930c777b0cSPeter Brune 6940c777b0cSPeter Brune Logically Collective on SNES 6950c777b0cSPeter Brune 6960c777b0cSPeter Brune Input Parameters: 6970c777b0cSPeter Brune + snes - the iterative context 6980c777b0cSPeter Brune - stype - scale type 6990c777b0cSPeter Brune 7000c777b0cSPeter Brune Options Database: 7010c777b0cSPeter Brune . -snes_qn_scale_type <shanno,none,linesearch,jacobian> 7020c777b0cSPeter Brune 7030c777b0cSPeter Brune Level: intermediate 7040c777b0cSPeter Brune 70584c577daSBarry Smith SNESQNScaleTypes: 7060c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 7070c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 7080c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 709a01a0525SBarry Smith - SNES_QN_SCALE_JACOBIAN - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration 710a01a0525SBarry Smith of QN and at ever restart. 7110c777b0cSPeter Brune 712a01a0525SBarry Smith .keywords: scaling type 713a01a0525SBarry Smith 714a01a0525SBarry Smith .seealso: SNES, SNESQN, SNESLineSearch, SNESQNScaleType, SNESetJacobian() 7150c777b0cSPeter Brune @*/ 7160c777b0cSPeter Brune 7172150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 7182150357eSBarry Smith { 7190c777b0cSPeter Brune PetscErrorCode ierr; 7206e111a19SKarl Rupp 7210c777b0cSPeter Brune PetscFunctionBegin; 7220c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7230c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 7240c777b0cSPeter Brune PetscFunctionReturn(0); 7250c777b0cSPeter Brune } 7260c777b0cSPeter Brune 7270c777b0cSPeter Brune #undef __FUNCT__ 7280c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 7290adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 7300adebc6cSBarry Smith { 7310c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7326e111a19SKarl Rupp 7330c777b0cSPeter Brune PetscFunctionBegin; 7340c777b0cSPeter Brune qn->scale_type = stype; 7350c777b0cSPeter Brune PetscFunctionReturn(0); 7360c777b0cSPeter Brune } 7370c777b0cSPeter Brune 7380c777b0cSPeter Brune #undef __FUNCT__ 7390c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7400adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 7410adebc6cSBarry Smith { 7420c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7436e111a19SKarl Rupp 7440c777b0cSPeter Brune PetscFunctionBegin; 7450c777b0cSPeter Brune qn->restart_type = rtype; 7460c777b0cSPeter Brune PetscFunctionReturn(0); 7470c777b0cSPeter Brune } 7480c777b0cSPeter Brune 7491efc8c45SPeter Brune #undef __FUNCT__ 7501efc8c45SPeter Brune #define __FUNCT__ "SNESQNSetType" 7511efc8c45SPeter Brune /*@ 7521efc8c45SPeter Brune SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN. 7531efc8c45SPeter Brune 7541efc8c45SPeter Brune Logically Collective on SNES 7551efc8c45SPeter Brune 7561efc8c45SPeter Brune Input Parameters: 7571efc8c45SPeter Brune + snes - the iterative context 7581efc8c45SPeter Brune - qtype - variant type 7591efc8c45SPeter Brune 7601efc8c45SPeter Brune Options Database: 76184c577daSBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> 7621efc8c45SPeter Brune 7631efc8c45SPeter Brune Level: beginner 7641efc8c45SPeter Brune 7651efc8c45SPeter Brune SNESQNTypes: 7661efc8c45SPeter Brune + SNES_QN_LBFGS - LBFGS variant 7671efc8c45SPeter Brune . SNES_QN_BROYDEN - Broyden variant 7681efc8c45SPeter Brune - SNES_QN_BADBROYDEN - Bad Broyden variant 7691efc8c45SPeter Brune 77084c577daSBarry Smith .keywords: SNES, SNESQN, type, set, SNESQNType 7711efc8c45SPeter Brune @*/ 7721efc8c45SPeter Brune 7731efc8c45SPeter Brune PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype) 7741efc8c45SPeter Brune { 7751efc8c45SPeter Brune PetscErrorCode ierr; 7761efc8c45SPeter Brune 7771efc8c45SPeter Brune PetscFunctionBegin; 7781efc8c45SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7791efc8c45SPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));CHKERRQ(ierr); 7801efc8c45SPeter Brune PetscFunctionReturn(0); 7811efc8c45SPeter Brune } 7821efc8c45SPeter Brune 7831efc8c45SPeter Brune #undef __FUNCT__ 7841efc8c45SPeter Brune #define __FUNCT__ "SNESQNSetType_QN" 7851efc8c45SPeter Brune PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype) 7861efc8c45SPeter Brune { 7871efc8c45SPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7881efc8c45SPeter Brune 7891efc8c45SPeter Brune PetscFunctionBegin; 7901efc8c45SPeter Brune qn->type = qtype; 7911efc8c45SPeter Brune PetscFunctionReturn(0); 7921efc8c45SPeter Brune } 7931efc8c45SPeter Brune 7944b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 7954b11644fSPeter Brune /*MC 7964b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 7974b11644fSPeter Brune 7986cc8130cSPeter Brune Options Database: 7996cc8130cSPeter Brune 80084c577daSBarry Smith + -snes_qn_m <m> - Number of past states saved for the L-Broyden methods. 80184c577daSBarry Smith + -snes_qn_restart_type <powell,periodic,none> - set the restart type 802*6bdcc836SBarry Smith . -snes_qn_powell_gamma - Angle condition for restart. 8031867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 80484c577daSBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> - QN type 80584c577daSBarry Smith . -snes_qn_scale_type <shanno,none,linesearch,jacobian> - scaling performed on inner Jacobian 80672835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 80784c577daSBarry Smith - -snes_qn_monitor - Monitors the quasi-newton Jacobian. 8086cc8130cSPeter Brune 809b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 810b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 811b8840d6bSPeter Brune updates. 8126cc8130cSPeter Brune 8131867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 8141867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 81584c577daSBarry Smith iteration as the current iteration's values when constructing the approximate Jacobian. The second, composed, 81684c577daSBarry Smith perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 8171867fe5bSPeter Brune 8186cc8130cSPeter Brune References: 81996a0c994SBarry Smith + 1. - Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995. 82096a0c994SBarry Smith . 2. - R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods, 8210a8e8854SPeter Brune Technical Report, Northwestern University, June 1992. 82296a0c994SBarry Smith . 3. - Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE 8230a8e8854SPeter Brune Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 82496a0c994SBarry Smith - 4. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 8254f02bc6aSBarry Smith SIAM Review, 57(4), 2015 8264b11644fSPeter Brune 8274b11644fSPeter Brune Level: beginner 8284b11644fSPeter Brune 82904d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 8306cc8130cSPeter Brune 8314b11644fSPeter Brune M*/ 8324b11644fSPeter Brune #undef __FUNCT__ 8334b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 8348cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 8354b11644fSPeter Brune { 8364b11644fSPeter Brune PetscErrorCode ierr; 8379f83bee8SJed Brown SNES_QN *qn; 8384b11644fSPeter Brune 8394b11644fSPeter Brune PetscFunctionBegin; 8404b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 8414b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 8424b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 8434b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 8445cd83419SPeter Brune snes->ops->view = SNESView_QN; 8454b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 8464b11644fSPeter Brune 84747f26062SPeter Brune snes->pcside = PC_LEFT; 84847f26062SPeter Brune 84942f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 85042f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 85142f4f86dSBarry Smith 85288976e71SPeter Brune if (!snes->tolerancesset) { 8530e444f03SPeter Brune snes->max_funcs = 30000; 8540e444f03SPeter Brune snes->max_its = 10000; 85588976e71SPeter Brune } 8560e444f03SPeter Brune 857b00a9115SJed Brown ierr = PetscNewLog(snes,&qn);CHKERRQ(ierr); 8584b11644fSPeter Brune snes->data = (void*) qn; 8590ecc9e1dSPeter Brune qn->m = 10; 860b21d5a53SPeter Brune qn->scaling = 1.0; 8610298fd71SBarry Smith qn->U = NULL; 8620298fd71SBarry Smith qn->V = NULL; 8630298fd71SBarry Smith qn->dXtdF = NULL; 8640298fd71SBarry Smith qn->dFtdX = NULL; 8650298fd71SBarry Smith qn->dXdFmat = NULL; 8660298fd71SBarry Smith qn->monitor = NULL; 8671c6523dcSPeter Brune qn->singlereduction = PETSC_TRUE; 868b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 86960c08b40SPeter Brune qn->scale_type = SNES_QN_SCALE_DEFAULT; 87060c08b40SPeter Brune qn->restart_type = SNES_QN_RESTART_DEFAULT; 871b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 872ea630c6eSPeter Brune 873bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr); 874bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr); 8751efc8c45SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);CHKERRQ(ierr); 8764b11644fSPeter Brune PetscFunctionReturn(0); 8774b11644fSPeter Brune } 878