10c777b0cSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 24b11644fSPeter Brune 38e231d97SPeter Brune #define H(i,j) qn->dXdFmat[i*qn->m + j] 48e231d97SPeter Brune 56a6fc655SJed Brown const char *const SNESQNScaleTypes[] = {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0}; 66a6fc655SJed Brown const char *const SNESQNRestartTypes[] = {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0}; 7b8840d6bSPeter Brune const char *const SNESQNTypes[] = {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_RESTART_",0}; 8b8840d6bSPeter Brune 9b8840d6bSPeter Brune typedef enum {SNES_QN_LBFGS = 0, 10b8840d6bSPeter Brune SNES_QN_BROYDEN = 1, 110ad7597dSKarl Rupp SNES_QN_BADBROYDEN = 2 120ad7597dSKarl Rupp } SNESQNType; 136bf1b2e5SPeter Brune 144b11644fSPeter Brune typedef struct { 15b8840d6bSPeter Brune Vec *U; /* Stored past states (vary from method to method) */ 16b8840d6bSPeter Brune Vec *V; /* Stored past states (vary from method to method) */ 178e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 185c7a0a03SPeter Brune PetscReal *lambda; /* The line search history of the method */ 195c7a0a03SPeter Brune PetscReal *norm; /* norms of the steps */ 208e231d97SPeter Brune PetscScalar *alpha, *beta; 218e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 2218aa0c0cSPeter Brune PetscBool singlereduction; /* Aggregated reduction implementation */ 238e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 2444f7e39eSPeter Brune PetscViewer monitor; 256bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 266bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 27b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 28b8840d6bSPeter Brune SNESQNType type; /* the type of quasi-newton method used */ 2988f769c5SPeter Brune SNESQNScaleType scale_type; /* the type of scaling used */ 300c777b0cSPeter Brune SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 319f83bee8SJed Brown } SNES_QN; 324b11644fSPeter Brune 334b11644fSPeter Brune #undef __FUNCT__ 34b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden" 35*5051edb1SPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D) 362150357eSBarry Smith { 374b11644fSPeter Brune PetscErrorCode ierr; 389f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 39b8840d6bSPeter Brune Vec W = snes->work[3]; 40b8840d6bSPeter Brune Vec *U = qn->U; 41b8840d6bSPeter Brune KSPConvergedReason kspreason; 42b8840d6bSPeter Brune PetscInt m = qn->m; 435c7a0a03SPeter Brune PetscInt k,i,j,lits,l = m; 445c7a0a03SPeter Brune PetscReal unorm,a,b; 455c7a0a03SPeter Brune PetscReal *lambda=qn->lambda; 46b8840d6bSPeter Brune PetscScalar gdot; 4759e7931aSPeter Brune PetscReal udot; 482150357eSBarry Smith 49b8840d6bSPeter Brune PetscFunctionBegin; 50b8840d6bSPeter Brune if (it < m) l = it; 515c7a0a03SPeter Brune if (it > 0) { 525c7a0a03SPeter Brune k = (it-1)%l; 535c7a0a03SPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr); 54*5051edb1SPeter Brune ierr = VecCopy(U[k],Xold);CHKERRQ(ierr); 55*5051edb1SPeter Brune ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr); 565c7a0a03SPeter Brune if (qn->monitor) { 575c7a0a03SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 585c7a0a03SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %d of %d by lambda: %14.12e \n",k,l,lambda[k]);CHKERRQ(ierr); 595c7a0a03SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 605c7a0a03SPeter Brune } 615c7a0a03SPeter Brune } 62b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 63d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr); 64b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 65b8840d6bSPeter Brune if (kspreason < 0) { 66b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 67b8840d6bSPeter Brune ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 68b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 69b8840d6bSPeter Brune PetscFunctionReturn(0); 70b8840d6bSPeter Brune } 71b8840d6bSPeter Brune } 72b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 73b8840d6bSPeter Brune snes->linear_its += lits; 74b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 75b8840d6bSPeter Brune } else { 76b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 77b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 78b8840d6bSPeter Brune } 79b8840d6bSPeter Brune 80b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 81b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 825c7a0a03SPeter Brune j = (it+i-l)%l; 835c7a0a03SPeter Brune k = (it+i-l+1)%l; 845c7a0a03SPeter Brune 855c7a0a03SPeter Brune ierr = VecNorm(U[j],NORM_2,&unorm);CHKERRQ(ierr); 865c7a0a03SPeter Brune ierr = VecDot(U[j],Y,&gdot);CHKERRQ(ierr); 8759e7931aSPeter Brune unorm *= unorm; 8859e7931aSPeter Brune udot = PetscRealPart(gdot); 895c7a0a03SPeter Brune a = (lambda[j]/lambda[k]); 905c7a0a03SPeter Brune b = -(1.-lambda[j]); 9159e7931aSPeter Brune a *= udot/unorm; 9259e7931aSPeter Brune b *= udot/unorm; 935c7a0a03SPeter Brune ierr = VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);CHKERRQ(ierr); 945c7a0a03SPeter Brune 95b8840d6bSPeter Brune if (qn->monitor) { 96b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 975c7a0a03SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d and %d, gdot: %14.12e\n",k,j,gdot);CHKERRQ(ierr); 98b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 99b8840d6bSPeter Brune } 100b8840d6bSPeter Brune } 101b8840d6bSPeter Brune if (l > 0) { 102b8840d6bSPeter Brune k = (it-1)%l; 103b8840d6bSPeter Brune ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 1045c7a0a03SPeter Brune ierr = VecNorm(U[k],NORM_2,&unorm);CHKERRQ(ierr); 10559e7931aSPeter Brune unorm *= unorm; 10659e7931aSPeter Brune udot = PetscRealPart(gdot); 10759e7931aSPeter Brune a = unorm/(unorm-lambda[k]*udot); 10859e7931aSPeter Brune b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot); 109b8840d6bSPeter Brune if (qn->monitor) { 110b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 11159e7931aSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d: a: %14.12e b: %14.12e \n",k,a,b);CHKERRQ(ierr); 112b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 113b8840d6bSPeter Brune } 1145c7a0a03SPeter Brune ierr = VecAXPBY(Y,b,a,U[k]);CHKERRQ(ierr); 115b8840d6bSPeter Brune } 1165c7a0a03SPeter Brune l = m; 1175c7a0a03SPeter Brune if (it+1<m)l=it+1; 1185c7a0a03SPeter Brune k = it%l; 1195c7a0a03SPeter Brune if (qn->monitor) { 1205c7a0a03SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1215c7a0a03SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "setting vector %d of %d\n",k,l);CHKERRQ(ierr); 1225c7a0a03SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1235c7a0a03SPeter Brune } 1245c7a0a03SPeter Brune ierr = VecCopy(Y,U[k]);CHKERRQ(ierr); 125b8840d6bSPeter Brune PetscFunctionReturn(0); 126b8840d6bSPeter Brune } 127b8840d6bSPeter Brune 128b8840d6bSPeter Brune #undef __FUNCT__ 129b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden" 1300adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1310adebc6cSBarry Smith { 132b8840d6bSPeter Brune PetscErrorCode ierr; 133b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 134b8840d6bSPeter Brune Vec W = snes->work[3]; 135b8840d6bSPeter Brune Vec *U = qn->U; 136b8840d6bSPeter Brune Vec *T = qn->V; 137b8840d6bSPeter Brune 138b8840d6bSPeter Brune /* ksp thing for jacobian scaling */ 139b8840d6bSPeter Brune KSPConvergedReason kspreason; 140*5051edb1SPeter Brune PetscInt k,j,i,lits; 141b8840d6bSPeter Brune PetscInt m = qn->m; 142*5051edb1SPeter Brune PetscScalar gdot,udot; 143b8840d6bSPeter Brune PetscInt l = m; 1440adebc6cSBarry Smith 145b8840d6bSPeter Brune PetscFunctionBegin; 146b8840d6bSPeter Brune if (it < m) l = it; 147b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 148b8840d6bSPeter Brune if (l > 0) { 149b8840d6bSPeter Brune k = (it-1)%l; 150b8840d6bSPeter Brune ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 151b8840d6bSPeter Brune ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 152*5051edb1SPeter Brune ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr); 153*5051edb1SPeter Brune ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr); 154b8840d6bSPeter Brune if (qn->monitor) { 155b8840d6bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 156b8840d6bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 157b8840d6bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 158b8840d6bSPeter Brune } 159b8840d6bSPeter Brune } 160b8840d6bSPeter Brune 161b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 162d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 163b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 164b8840d6bSPeter Brune if (kspreason < 0) { 165b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 166b8840d6bSPeter Brune ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 167b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 168b8840d6bSPeter Brune PetscFunctionReturn(0); 169b8840d6bSPeter Brune } 170b8840d6bSPeter Brune } 171b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 172b8840d6bSPeter Brune snes->linear_its += lits; 173b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 174b8840d6bSPeter Brune } else { 175b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 176b8840d6bSPeter Brune } 177*5051edb1SPeter Brune 178*5051edb1SPeter Brune /* inward recursion starting at the first update and working forward */ 179*5051edb1SPeter Brune if (l > 0) { 180*5051edb1SPeter Brune for (i = 0; i < l-1; i++) { 181*5051edb1SPeter Brune k = (it+i-l)%l; 182*5051edb1SPeter Brune j = (it-1)%l; 183*5051edb1SPeter Brune ierr = VecDotBegin(U[k],U[j],&gdot);CHKERRQ(ierr); 184*5051edb1SPeter Brune ierr = VecDotBegin(U[k],U[k],&udot);CHKERRQ(ierr); 185*5051edb1SPeter Brune ierr = VecDotEnd(U[k],U[j],&gdot);CHKERRQ(ierr); 186*5051edb1SPeter Brune ierr = VecDotEnd(U[k],U[k],&udot);CHKERRQ(ierr); 187*5051edb1SPeter Brune ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr); 188*5051edb1SPeter Brune if (qn->monitor) { 189*5051edb1SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 190*5051edb1SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 191*5051edb1SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 192*5051edb1SPeter Brune } 193*5051edb1SPeter Brune } 194*5051edb1SPeter Brune } 195b8840d6bSPeter Brune PetscFunctionReturn(0); 196b8840d6bSPeter Brune } 197b8840d6bSPeter Brune 198b8840d6bSPeter Brune #undef __FUNCT__ 199b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 2000adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 2010adebc6cSBarry Smith { 202b8840d6bSPeter Brune PetscErrorCode ierr; 203b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 204b8840d6bSPeter Brune Vec W = snes->work[3]; 205b8840d6bSPeter Brune Vec *dX = qn->U; 206b8840d6bSPeter Brune Vec *dF = qn->V; 207bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 208bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 2098e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 210b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 2118e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 212bd052dfeSPeter Brune 2130ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 2140ecc9e1dSPeter Brune KSPConvergedReason kspreason; 2158e231d97SPeter Brune PetscInt k,i,j,g,lits; 2164b11644fSPeter Brune PetscInt m = qn->m; 2174b11644fSPeter Brune PetscScalar t; 2184b11644fSPeter Brune PetscInt l = m; 2198e231d97SPeter Brune 2204b11644fSPeter Brune PetscFunctionBegin; 2215ba6227bSPeter Brune if (it < m) l = it; 2221c6523dcSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 223b8840d6bSPeter Brune if (it > 0) { 224b8840d6bSPeter Brune k = (it - 1) % l; 225b8840d6bSPeter Brune ierr = VecCopy(D,dF[k]);CHKERRQ(ierr); 226b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 227b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 228b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 22918aa0c0cSPeter Brune if (qn->singlereduction) { 23023d44fbcSPeter Brune PetscScalar dFtdF; 2311c6523dcSPeter Brune ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2321c6523dcSPeter Brune ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2331c6523dcSPeter Brune ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr); 23423d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) {ierr = VecDotBegin(dF[k],dF[k],&dFtdF);CHKERRQ(ierr);} 2351c6523dcSPeter Brune ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2361c6523dcSPeter Brune ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2371c6523dcSPeter Brune ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr); 23823d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 23923d44fbcSPeter Brune ierr = VecDotEnd(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 24023d44fbcSPeter Brune qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(dFtdF); 24123d44fbcSPeter Brune } 242b8840d6bSPeter Brune for (j = 0; j < l; j++) { 243b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 244b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 245b8840d6bSPeter Brune } 246b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 2471aa26658SKarl Rupp for (j = 0; j < l; j++) dXtdF[j] = H(j, j); 248b8840d6bSPeter Brune } else { 249b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 250b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 25167392de3SBarry Smith PetscReal dFtdF; 25267392de3SBarry Smith ierr = VecDotRealPart(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 25367392de3SBarry Smith qn->scaling = PetscRealPart(dXtdF[k])/dFtdF; 25423d44fbcSPeter Brune } 25523d44fbcSPeter Brune } 25623d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 257b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 258b8840d6bSPeter Brune } 259b8840d6bSPeter Brune } 260b8840d6bSPeter Brune 2614b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2624b11644fSPeter Brune for (i=0; i<l; i++) { 263b21d5a53SPeter Brune k = (it-i-1)%l; 26418aa0c0cSPeter Brune if (qn->singlereduction) { 2658e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2668e231d97SPeter Brune t = YtdX[k]; 2678e231d97SPeter Brune for (j=0; j<i; j++) { 2688e231d97SPeter Brune g = (it-j-1)%l; 269b3cfccb2SPeter Brune t -= alpha[g]*H(k, g); 2708e231d97SPeter Brune } 2718e231d97SPeter Brune alpha[k] = t/H(k,k); 2728e231d97SPeter Brune } else { 2733af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2748e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2758e231d97SPeter Brune } 27644f7e39eSPeter Brune if (qn->monitor) { 2773af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2785ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2793af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 28044f7e39eSPeter Brune } 2816bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2824b11644fSPeter Brune } 2834b11644fSPeter Brune 2840c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 285d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 2860ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2870ecc9e1dSPeter Brune if (kspreason < 0) { 2880ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 2890ecc9e1dSPeter Brune ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 2900ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 2910ecc9e1dSPeter Brune PetscFunctionReturn(0); 2920ecc9e1dSPeter Brune } 2930ecc9e1dSPeter Brune } 2940ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2950ecc9e1dSPeter Brune snes->linear_its += lits; 296b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2970ecc9e1dSPeter Brune } else { 298b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2990ecc9e1dSPeter Brune } 30018aa0c0cSPeter Brune if (qn->singlereduction) { 301b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 3028e231d97SPeter Brune } 3034b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 304bd052dfeSPeter Brune for (i = 0; i < l; i++) { 305b21d5a53SPeter Brune k = (it + i - l) % l; 30618aa0c0cSPeter Brune if (qn->singlereduction) { 3078e231d97SPeter Brune t = YtdX[k]; 3088e231d97SPeter Brune for (j = 0; j < i; j++) { 3098e231d97SPeter Brune g = (it + j - l) % l; 310b3cfccb2SPeter Brune t += (alpha[g] - beta[g])*H(g, k); 3118e231d97SPeter Brune } 3128e231d97SPeter Brune beta[k] = t / H(k, k); 3138e231d97SPeter Brune } else { 3146bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 3158e231d97SPeter Brune beta[k] = t / dXtdF[k]; 3168e231d97SPeter Brune } 31722d28d08SBarry Smith ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 31844f7e39eSPeter Brune if (qn->monitor) { 3193af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3205ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 3213af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 32244f7e39eSPeter Brune } 3234b11644fSPeter Brune } 3244b11644fSPeter Brune PetscFunctionReturn(0); 3254b11644fSPeter Brune } 3264b11644fSPeter Brune 3274b11644fSPeter Brune #undef __FUNCT__ 3284b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 3294b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 3304b11644fSPeter Brune { 3314b11644fSPeter Brune PetscErrorCode ierr; 3329f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 33315f5eeeaSPeter Brune Vec X,Xold; 33447f26062SPeter Brune Vec F; 33547f26062SPeter Brune Vec Y,D,Dold; 336b8840d6bSPeter Brune PetscInt i, i_r; 337335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 3380c777b0cSPeter Brune PetscBool lssucceed,powell,periodic; 3391c6523dcSPeter Brune PetscScalar DolddotD,DolddotDold; 3400ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 341b7281c8aSPeter Brune SNESConvergedReason reason; 3420ecc9e1dSPeter Brune 3434b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 3444b11644fSPeter Brune 3456e111a19SKarl Rupp PetscFunctionBegin; 3469f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3473af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 348b8840d6bSPeter Brune 349b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 350335fb975SPeter Brune Xold = snes->work[0]; 3513af51624SPeter Brune 3523af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 353335fb975SPeter Brune D = snes->work[1]; 354335fb975SPeter Brune Dold = snes->work[2]; 3554b11644fSPeter Brune 3564b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3574b11644fSPeter Brune 358ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3594b11644fSPeter Brune snes->iter = 0; 3604b11644fSPeter Brune snes->norm = 0.; 361ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 36247f26062SPeter Brune 363b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 36447f26062SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 365b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 366b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 367b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 368b7281c8aSPeter Brune PetscFunctionReturn(0); 369b7281c8aSPeter Brune } 37047f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 37147f26062SPeter Brune } else { 372e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 37315f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3744b11644fSPeter Brune if (snes->domainerror) { 3754b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3764b11644fSPeter Brune PetscFunctionReturn(0); 3774b11644fSPeter Brune } 3781aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 379c1c75074SPeter Brune 38047f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 381189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 382189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 383189a9710SBarry Smith PetscFunctionReturn(0); 384189a9710SBarry Smith } 38547f26062SPeter Brune } 386b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 38747f26062SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 388b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 389b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 390b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 391b7281c8aSPeter Brune PetscFunctionReturn(0); 392b7281c8aSPeter Brune } 39347f26062SPeter Brune } else { 39447f26062SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 39547f26062SPeter Brune } 396b28a06ddSPeter Brune 397ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3984b11644fSPeter Brune snes->norm = fnorm; 399ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 400a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 4014b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 4024b11644fSPeter Brune 4034b11644fSPeter Brune /* test convergence */ 4044b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4054b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 40670d3b23bSPeter Brune 4073cf07b75SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 4083cf07b75SPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 4093cf07b75SPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 4103cf07b75SPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 411ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 412ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 413ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 414ddd40ce5SPeter Brune PetscFunctionReturn(0); 415ddd40ce5SPeter Brune } 416ddd40ce5SPeter Brune ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 4173cf07b75SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 4183cf07b75SPeter Brune } 4193cf07b75SPeter Brune 420f8e15203SPeter Brune /* scale the initial update */ 4210c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4220ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4238cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 4240ecc9e1dSPeter Brune } 4250ecc9e1dSPeter Brune 4265ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 427b8840d6bSPeter Brune switch (qn->type) { 428b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 429b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 430*5051edb1SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 431b8840d6bSPeter Brune break; 432b8840d6bSPeter Brune case SNES_QN_BROYDEN: 433*5051edb1SPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr); 434b8840d6bSPeter Brune break; 435b8840d6bSPeter Brune case SNES_QN_LBFGS: 436b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 437*5051edb1SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 438b8840d6bSPeter Brune break; 439b8840d6bSPeter Brune } 44070d3b23bSPeter Brune /* line search for lambda */ 44170d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 44215f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 443f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4449f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 4459f3a0142SPeter Brune if (snes->domainerror) { 4469f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 4479f3a0142SPeter Brune PetscFunctionReturn(0); 4489f3a0142SPeter Brune } 449f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 4509f3a0142SPeter Brune if (!lssucceed) { 4519f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4529f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4539f3a0142SPeter Brune break; 4549f3a0142SPeter Brune } 4559f3a0142SPeter Brune } 456f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 4570c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 45805ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4590ecc9e1dSPeter Brune } 4603af51624SPeter Brune 46188d374b2SPeter Brune /* convergence monitoring */ 462b21d5a53SPeter 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); 463b21d5a53SPeter Brune 464b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 465b28a06ddSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 466b28a06ddSPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 467b28a06ddSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 468ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 469ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 470ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 471ddd40ce5SPeter Brune PetscFunctionReturn(0); 472ddd40ce5SPeter Brune } 473ddd40ce5SPeter Brune ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 474b28a06ddSPeter Brune } 475b28a06ddSPeter Brune 476360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 477360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 478360c497dSPeter Brune 479a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 4808409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4814b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 482d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4834b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 484b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 48547f26062SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 486b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 487b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 488b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 489b7281c8aSPeter Brune PetscFunctionReturn(0); 490b7281c8aSPeter Brune } 49188d374b2SPeter Brune } else { 49288d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 49388d374b2SPeter Brune } 49488d374b2SPeter Brune 4950c777b0cSPeter Brune powell = PETSC_FALSE; 4960c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 4976bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 4988e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4998e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 5008e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 5018e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 5020c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 5030c777b0cSPeter Brune } 5040c777b0cSPeter Brune periodic = PETSC_FALSE; 505b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 506b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 5070c777b0cSPeter Brune } 5080c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 5090c777b0cSPeter Brune if (powell || periodic) { 5105ba6227bSPeter Brune if (qn->monitor) { 5115ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 512516fe3c3SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);CHKERRQ(ierr); 5135ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 5145ba6227bSPeter Brune } 5155ba6227bSPeter Brune i_r = -1; 5165ba6227bSPeter Brune /* general purpose update */ 5175ba6227bSPeter Brune if (snes->ops->update) { 5185ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 5195ba6227bSPeter Brune } 5200c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5210ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 5228cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 5230ecc9e1dSPeter Brune } 5240ecc9e1dSPeter Brune } 52570d3b23bSPeter Brune /* general purpose update */ 52670d3b23bSPeter Brune if (snes->ops->update) { 52770d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 52870d3b23bSPeter Brune } 5295ba6227bSPeter Brune } 5304b11644fSPeter Brune if (i == snes->max_its) { 5314b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 5324b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 5334b11644fSPeter Brune } 5344b11644fSPeter Brune PetscFunctionReturn(0); 5354b11644fSPeter Brune } 5364b11644fSPeter Brune 5374b11644fSPeter Brune #undef __FUNCT__ 5384b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5394b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5404b11644fSPeter Brune { 5419f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5424b11644fSPeter Brune PetscErrorCode ierr; 543335fb975SPeter Brune 5444b11644fSPeter Brune PetscFunctionBegin; 545b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 54659e7931aSPeter Brune if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 5475c7a0a03SPeter Brune ierr = PetscMalloc4(qn->m, PetscScalar, &qn->alpha, 5488e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 5495c7a0a03SPeter Brune qn->m, PetscScalar, &qn->dXtdF, 5505c7a0a03SPeter Brune qn->m, PetscReal, &qn->lambda);CHKERRQ(ierr); 5518e231d97SPeter Brune 55218aa0c0cSPeter Brune if (qn->singlereduction) { 5538e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 5548e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 5558e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 5568e231d97SPeter Brune } 557fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 558335fb975SPeter Brune 559335fb975SPeter Brune /* set up the line search */ 5600c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5618e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5628e231d97SPeter Brune } 5636c67d002SPeter Brune 5646c67d002SPeter Brune if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;} 5656c67d002SPeter Brune 5664b11644fSPeter Brune PetscFunctionReturn(0); 5674b11644fSPeter Brune } 5684b11644fSPeter Brune 5694b11644fSPeter Brune #undef __FUNCT__ 5704b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5714b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5724b11644fSPeter Brune { 5734b11644fSPeter Brune PetscErrorCode ierr; 5749f83bee8SJed Brown SNES_QN *qn; 5750adebc6cSBarry Smith 5764b11644fSPeter Brune PetscFunctionBegin; 5774b11644fSPeter Brune if (snes->data) { 5789f83bee8SJed Brown qn = (SNES_QN*)snes->data; 579b8840d6bSPeter Brune if (qn->U) { 580b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5814b11644fSPeter Brune } 582b8840d6bSPeter Brune if (qn->V) { 583b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5844b11644fSPeter Brune } 58518aa0c0cSPeter Brune if (qn->singlereduction) { 5868e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5878e231d97SPeter Brune } 5885c7a0a03SPeter Brune ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr); 5894b11644fSPeter Brune } 5904b11644fSPeter Brune PetscFunctionReturn(0); 5914b11644fSPeter Brune } 5924b11644fSPeter Brune 5934b11644fSPeter Brune #undef __FUNCT__ 5944b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5954b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5964b11644fSPeter Brune { 5974b11644fSPeter Brune PetscErrorCode ierr; 5986e111a19SKarl Rupp 5994b11644fSPeter Brune PetscFunctionBegin; 6004b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 6014b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 602bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr); 6034b11644fSPeter Brune PetscFunctionReturn(0); 6044b11644fSPeter Brune } 6054b11644fSPeter Brune 6064b11644fSPeter Brune #undef __FUNCT__ 6074b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 6084b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 6094b11644fSPeter Brune { 6104b11644fSPeter Brune 6114b11644fSPeter Brune PetscErrorCode ierr; 6122150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 61388f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 614f1c6b773SPeter Brune SNESLineSearch linesearch; 6152150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 6162150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 6172150357eSBarry Smith 6184b11644fSPeter Brune PetscFunctionBegin; 6194b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 6200298fd71SBarry Smith ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr); 6210298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr); 6220298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr); 6230298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr); 6240298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr); 62588f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 62688f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 62788f769c5SPeter Brune 62888f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 62988f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 63088f769c5SPeter Brune 631b8840d6bSPeter Brune ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 6320298fd71SBarry Smith (PetscEnum)qn->type,(PetscEnum*)&qn->type,NULL);CHKERRQ(ierr); 6334b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6349e764e56SPeter Brune if (!snes->linesearch) { 6357601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 6361a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 6379e764e56SPeter Brune } 63844f7e39eSPeter Brune if (monflg) { 639ce94432eSBarry Smith qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 64044f7e39eSPeter Brune } 6414b11644fSPeter Brune PetscFunctionReturn(0); 6424b11644fSPeter Brune } 6434b11644fSPeter Brune 64492c02d66SPeter Brune 6450c777b0cSPeter Brune #undef __FUNCT__ 6460c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6470c777b0cSPeter Brune /*@ 6480c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6490c777b0cSPeter Brune 6500c777b0cSPeter Brune Logically Collective on SNES 6510c777b0cSPeter Brune 6520c777b0cSPeter Brune Input Parameters: 6530c777b0cSPeter Brune + snes - the iterative context 6540c777b0cSPeter Brune - rtype - restart type 6550c777b0cSPeter Brune 6560c777b0cSPeter Brune Options Database: 6570c777b0cSPeter Brune + -snes_qn_restart_type<powell,periodic,none> - set the restart type 658b8840d6bSPeter Brune - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 6590c777b0cSPeter Brune 6600c777b0cSPeter Brune Level: intermediate 6610c777b0cSPeter Brune 6620c777b0cSPeter Brune SNESQNRestartTypes: 6630c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6640c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6650c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6660c777b0cSPeter Brune 6670c777b0cSPeter Brune Notes: 6680c777b0cSPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 6690c777b0cSPeter Brune 6700c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 6710c777b0cSPeter Brune @*/ 6722150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6732150357eSBarry Smith { 6740c777b0cSPeter Brune PetscErrorCode ierr; 6756e111a19SKarl Rupp 6760c777b0cSPeter Brune PetscFunctionBegin; 6770c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6780c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6790c777b0cSPeter Brune PetscFunctionReturn(0); 6800c777b0cSPeter Brune } 6810c777b0cSPeter Brune 6820c777b0cSPeter Brune #undef __FUNCT__ 6830c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6840c777b0cSPeter Brune /*@ 6850c777b0cSPeter Brune SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 6860c777b0cSPeter Brune 6870c777b0cSPeter Brune Logically Collective on SNES 6880c777b0cSPeter Brune 6890c777b0cSPeter Brune Input Parameters: 6900c777b0cSPeter Brune + snes - the iterative context 6910c777b0cSPeter Brune - stype - scale type 6920c777b0cSPeter Brune 6930c777b0cSPeter Brune Options Database: 6940c777b0cSPeter Brune . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 6950c777b0cSPeter Brune 6960c777b0cSPeter Brune Level: intermediate 6970c777b0cSPeter Brune 6980c777b0cSPeter Brune SNESQNSelectTypes: 6990c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 7000c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 7010c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 7020c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 7030c777b0cSPeter Brune 7040c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 7050c777b0cSPeter Brune @*/ 7060c777b0cSPeter Brune 7072150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 7082150357eSBarry Smith { 7090c777b0cSPeter Brune PetscErrorCode ierr; 7106e111a19SKarl Rupp 7110c777b0cSPeter Brune PetscFunctionBegin; 7120c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7130c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 7140c777b0cSPeter Brune PetscFunctionReturn(0); 7150c777b0cSPeter Brune } 7160c777b0cSPeter Brune 7170c777b0cSPeter Brune #undef __FUNCT__ 7180c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 7190adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 7200adebc6cSBarry Smith { 7210c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7226e111a19SKarl Rupp 7230c777b0cSPeter Brune PetscFunctionBegin; 7240c777b0cSPeter Brune qn->scale_type = stype; 7250c777b0cSPeter Brune PetscFunctionReturn(0); 7260c777b0cSPeter Brune } 7270c777b0cSPeter Brune 7280c777b0cSPeter Brune #undef __FUNCT__ 7290c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7300adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 7310adebc6cSBarry Smith { 7320c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7336e111a19SKarl Rupp 7340c777b0cSPeter Brune PetscFunctionBegin; 7350c777b0cSPeter Brune qn->restart_type = rtype; 7360c777b0cSPeter Brune PetscFunctionReturn(0); 7370c777b0cSPeter Brune } 7380c777b0cSPeter Brune 7394b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 7404b11644fSPeter Brune /*MC 7414b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 7424b11644fSPeter Brune 7436cc8130cSPeter Brune Options Database: 7446cc8130cSPeter Brune 7456cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 7461867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 7471867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 74872835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 74944f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 7506cc8130cSPeter Brune 751b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 752b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 753b8840d6bSPeter Brune updates. 7546cc8130cSPeter Brune 7551867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 7561867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 7571867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 7581867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 7591867fe5bSPeter Brune 7606cc8130cSPeter Brune References: 7616cc8130cSPeter Brune 7620a8e8854SPeter Brune Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995. 7636cc8130cSPeter Brune 7640a8e8854SPeter Brune R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi-Newton Matrices and their use in Limited Memory Methods, 7650a8e8854SPeter Brune Technical Report, Northwestern University, June 1992. 7660a8e8854SPeter Brune 7670a8e8854SPeter Brune Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE 7680a8e8854SPeter Brune Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 769b8840d6bSPeter Brune 7704b11644fSPeter Brune 7714b11644fSPeter Brune Level: beginner 7724b11644fSPeter Brune 77304d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 7746cc8130cSPeter Brune 7754b11644fSPeter Brune M*/ 7764b11644fSPeter Brune #undef __FUNCT__ 7774b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 7788cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 7794b11644fSPeter Brune { 7804b11644fSPeter Brune PetscErrorCode ierr; 7819f83bee8SJed Brown SNES_QN *qn; 7824b11644fSPeter Brune 7834b11644fSPeter Brune PetscFunctionBegin; 7844b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 7854b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 7864b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 7874b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 7884b11644fSPeter Brune snes->ops->view = 0; 7894b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 7904b11644fSPeter Brune 79147f26062SPeter Brune snes->pcside = PC_LEFT; 79247f26062SPeter Brune 79342f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 79442f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 79542f4f86dSBarry Smith 79688976e71SPeter Brune if (!snes->tolerancesset) { 7970e444f03SPeter Brune snes->max_funcs = 30000; 7980e444f03SPeter Brune snes->max_its = 10000; 79988976e71SPeter Brune } 8000e444f03SPeter Brune 8019f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 8024b11644fSPeter Brune snes->data = (void*) qn; 8030ecc9e1dSPeter Brune qn->m = 10; 804b21d5a53SPeter Brune qn->scaling = 1.0; 8050298fd71SBarry Smith qn->U = NULL; 8060298fd71SBarry Smith qn->V = NULL; 8070298fd71SBarry Smith qn->dXtdF = NULL; 8080298fd71SBarry Smith qn->dFtdX = NULL; 8090298fd71SBarry Smith qn->dXdFmat = NULL; 8100298fd71SBarry Smith qn->monitor = NULL; 8111c6523dcSPeter Brune qn->singlereduction = PETSC_TRUE; 812b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 8130c777b0cSPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 8140c777b0cSPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 815b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 816ea630c6eSPeter Brune 817bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr); 818bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr); 8194b11644fSPeter Brune PetscFunctionReturn(0); 8204b11644fSPeter Brune } 821