10c777b0cSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 24b11644fSPeter Brune 38e231d97SPeter Brune #define H(i,j) qn->dXdFmat[i*qn->m + j] 48e231d97SPeter Brune 56a6fc655SJed Brown const char *const SNESQNScaleTypes[] = {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0}; 66a6fc655SJed Brown const char *const SNESQNRestartTypes[] = {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0}; 7b8840d6bSPeter Brune const char *const SNESQNTypes[] = {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_RESTART_",0}; 8b8840d6bSPeter Brune 9b8840d6bSPeter Brune typedef enum {SNES_QN_LBFGS = 0, 10b8840d6bSPeter Brune SNES_QN_BROYDEN = 1, 110ad7597dSKarl Rupp SNES_QN_BADBROYDEN = 2 120ad7597dSKarl Rupp } SNESQNType; 136bf1b2e5SPeter Brune 144b11644fSPeter Brune typedef struct { 15b8840d6bSPeter Brune Vec *U; /* Stored past states (vary from method to method) */ 16b8840d6bSPeter Brune Vec *V; /* Stored past states (vary from method to method) */ 178e231d97SPeter Brune PetscInt m; /* The number of kept previous steps */ 185c7a0a03SPeter Brune PetscReal *lambda; /* The line search history of the method */ 195c7a0a03SPeter Brune PetscReal *norm; /* norms of the steps */ 208e231d97SPeter Brune PetscScalar *alpha, *beta; 218e231d97SPeter Brune PetscScalar *dXtdF, *dFtdX, *YtdX; 2218aa0c0cSPeter Brune PetscBool singlereduction; /* Aggregated reduction implementation */ 238e231d97SPeter Brune PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 2444f7e39eSPeter Brune PetscViewer monitor; 256bf1b2e5SPeter Brune PetscReal powell_gamma; /* Powell angle restart condition */ 266bf1b2e5SPeter Brune PetscReal powell_downhill; /* Powell descent restart condition */ 27b21d5a53SPeter Brune PetscReal scaling; /* scaling of H0 */ 28b8840d6bSPeter Brune SNESQNType type; /* the type of quasi-newton method used */ 2988f769c5SPeter Brune SNESQNScaleType scale_type; /* the type of scaling used */ 300c777b0cSPeter Brune SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 319f83bee8SJed Brown } SNES_QN; 324b11644fSPeter Brune 334b11644fSPeter Brune #undef __FUNCT__ 34b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden" 355051edb1SPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D) 362150357eSBarry Smith { 374b11644fSPeter Brune PetscErrorCode ierr; 389f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 39b8840d6bSPeter Brune Vec W = snes->work[3]; 40b8840d6bSPeter Brune Vec *U = qn->U; 41b8840d6bSPeter Brune KSPConvergedReason kspreason; 42b8840d6bSPeter Brune PetscInt m = qn->m; 435c7a0a03SPeter Brune PetscInt k,i,j,lits,l = m; 445c7a0a03SPeter Brune PetscReal unorm,a,b; 455c7a0a03SPeter Brune PetscReal *lambda=qn->lambda; 46b8840d6bSPeter Brune PetscScalar gdot; 4759e7931aSPeter Brune PetscReal udot; 482150357eSBarry Smith 49b8840d6bSPeter Brune PetscFunctionBegin; 50b8840d6bSPeter Brune if (it < m) l = it; 515c7a0a03SPeter Brune if (it > 0) { 525c7a0a03SPeter Brune k = (it-1)%l; 535c7a0a03SPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr); 545051edb1SPeter Brune ierr = VecCopy(U[k],Xold);CHKERRQ(ierr); 555051edb1SPeter Brune ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr); 565c7a0a03SPeter Brune if (qn->monitor) { 575c7a0a03SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 585c7a0a03SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %d of %d by lambda: %14.12e \n",k,l,lambda[k]);CHKERRQ(ierr); 595c7a0a03SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 605c7a0a03SPeter Brune } 615c7a0a03SPeter Brune } 62b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 63d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr); 64b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 65b8840d6bSPeter Brune if (kspreason < 0) { 66b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 67b8840d6bSPeter Brune ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 68b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 69b8840d6bSPeter Brune PetscFunctionReturn(0); 70b8840d6bSPeter Brune } 71b8840d6bSPeter Brune } 72b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 73b8840d6bSPeter Brune snes->linear_its += lits; 74b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 75b8840d6bSPeter Brune } else { 76b8840d6bSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 77b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 78b8840d6bSPeter Brune } 79b8840d6bSPeter Brune 80b8840d6bSPeter Brune /* inward recursion starting at the first update and working forward */ 81b8840d6bSPeter Brune for (i = 0; i < l-1; i++) { 825c7a0a03SPeter Brune j = (it+i-l)%l; 835c7a0a03SPeter Brune k = (it+i-l+1)%l; 845c7a0a03SPeter Brune 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; 1405051edb1SPeter Brune PetscInt k,j,i,lits; 141b8840d6bSPeter Brune PetscInt m = qn->m; 1425051edb1SPeter 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); 1525051edb1SPeter Brune ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr); 1535051edb1SPeter Brune ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr); 154b8840d6bSPeter Brune } 155b8840d6bSPeter Brune 156b8840d6bSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 157d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 158b8840d6bSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 159b8840d6bSPeter Brune if (kspreason < 0) { 160b8840d6bSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 161b8840d6bSPeter 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); 162b8840d6bSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 163b8840d6bSPeter Brune PetscFunctionReturn(0); 164b8840d6bSPeter Brune } 165b8840d6bSPeter Brune } 166b8840d6bSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 167b8840d6bSPeter Brune snes->linear_its += lits; 168b8840d6bSPeter Brune ierr = VecCopy(W,Y);CHKERRQ(ierr); 169b8840d6bSPeter Brune } else { 170b8840d6bSPeter Brune ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 171b8840d6bSPeter Brune } 1725051edb1SPeter Brune 1735051edb1SPeter Brune /* inward recursion starting at the first update and working forward */ 1745051edb1SPeter Brune if (l > 0) { 1755051edb1SPeter Brune for (i = 0; i < l-1; i++) { 1765051edb1SPeter Brune k = (it+i-l)%l; 1775051edb1SPeter Brune j = (it-1)%l; 1785051edb1SPeter Brune ierr = VecDotBegin(U[k],U[j],&gdot);CHKERRQ(ierr); 1795051edb1SPeter Brune ierr = VecDotBegin(U[k],U[k],&udot);CHKERRQ(ierr); 1805051edb1SPeter Brune ierr = VecDotEnd(U[k],U[j],&gdot);CHKERRQ(ierr); 1815051edb1SPeter Brune ierr = VecDotEnd(U[k],U[k],&udot);CHKERRQ(ierr); 1825051edb1SPeter Brune ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr); 1835051edb1SPeter Brune if (qn->monitor) { 1845051edb1SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1855051edb1SPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 1865051edb1SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 1875051edb1SPeter Brune } 1885051edb1SPeter Brune } 1895051edb1SPeter Brune } 190b8840d6bSPeter Brune PetscFunctionReturn(0); 191b8840d6bSPeter Brune } 192b8840d6bSPeter Brune 193b8840d6bSPeter Brune #undef __FUNCT__ 194b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS" 1950adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 1960adebc6cSBarry Smith { 197b8840d6bSPeter Brune PetscErrorCode ierr; 198b8840d6bSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 199b8840d6bSPeter Brune Vec W = snes->work[3]; 200b8840d6bSPeter Brune Vec *dX = qn->U; 201b8840d6bSPeter Brune Vec *dF = qn->V; 202bd052dfeSPeter Brune PetscScalar *alpha = qn->alpha; 203bd052dfeSPeter Brune PetscScalar *beta = qn->beta; 2048e231d97SPeter Brune PetscScalar *dXtdF = qn->dXtdF; 205b8840d6bSPeter Brune PetscScalar *dFtdX = qn->dFtdX; 2068e231d97SPeter Brune PetscScalar *YtdX = qn->YtdX; 207bd052dfeSPeter Brune 2080ecc9e1dSPeter Brune /* ksp thing for jacobian scaling */ 2090ecc9e1dSPeter Brune KSPConvergedReason kspreason; 2108e231d97SPeter Brune PetscInt k,i,j,g,lits; 2114b11644fSPeter Brune PetscInt m = qn->m; 2124b11644fSPeter Brune PetscScalar t; 2134b11644fSPeter Brune PetscInt l = m; 2148e231d97SPeter Brune 2154b11644fSPeter Brune PetscFunctionBegin; 2165ba6227bSPeter Brune if (it < m) l = it; 2171c6523dcSPeter Brune ierr = VecCopy(D,Y);CHKERRQ(ierr); 218b8840d6bSPeter Brune if (it > 0) { 219b8840d6bSPeter Brune k = (it - 1) % l; 220b8840d6bSPeter Brune ierr = VecCopy(D,dF[k]);CHKERRQ(ierr); 221b8840d6bSPeter Brune ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 222b8840d6bSPeter Brune ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 223b8840d6bSPeter Brune ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 22418aa0c0cSPeter Brune if (qn->singlereduction) { 2251c6523dcSPeter Brune ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2261c6523dcSPeter Brune ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2271c6523dcSPeter Brune ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr); 2281c6523dcSPeter Brune ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr); 2291c6523dcSPeter Brune ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr); 2301c6523dcSPeter Brune ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr); 231b8840d6bSPeter Brune for (j = 0; j < l; j++) { 232b8840d6bSPeter Brune H(k, j) = dFtdX[j]; 233b8840d6bSPeter Brune H(j, k) = dXtdF[j]; 234b8840d6bSPeter Brune } 235b8840d6bSPeter Brune /* copy back over to make the computation of alpha and beta easier */ 2361aa26658SKarl Rupp for (j = 0; j < l; j++) dXtdF[j] = H(j, j); 237b8840d6bSPeter Brune } else { 238b8840d6bSPeter Brune ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 23923d44fbcSPeter Brune } 24023d44fbcSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 241b8840d6bSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 242b8840d6bSPeter Brune } 243b8840d6bSPeter Brune } 244b8840d6bSPeter Brune 2454b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 2464b11644fSPeter Brune for (i=0; i<l; i++) { 247b21d5a53SPeter Brune k = (it-i-1)%l; 24818aa0c0cSPeter Brune if (qn->singlereduction) { 2498e231d97SPeter Brune /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 2508e231d97SPeter Brune t = YtdX[k]; 2518e231d97SPeter Brune for (j=0; j<i; j++) { 2528e231d97SPeter Brune g = (it-j-1)%l; 253b3cfccb2SPeter Brune t -= alpha[g]*H(k, g); 2548e231d97SPeter Brune } 2558e231d97SPeter Brune alpha[k] = t/H(k,k); 2568e231d97SPeter Brune } else { 2573af51624SPeter Brune ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 2588e231d97SPeter Brune alpha[k] = t/dXtdF[k]; 2598e231d97SPeter Brune } 26044f7e39eSPeter Brune if (qn->monitor) { 2613af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 2625ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 2633af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 26444f7e39eSPeter Brune } 2656bf1b2e5SPeter Brune ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 2664b11644fSPeter Brune } 2674b11644fSPeter Brune 2680c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 269d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 2700ecc9e1dSPeter Brune ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2710ecc9e1dSPeter Brune if (kspreason < 0) { 2720ecc9e1dSPeter Brune if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 2730ecc9e1dSPeter 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); 2740ecc9e1dSPeter Brune snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 2750ecc9e1dSPeter Brune PetscFunctionReturn(0); 2760ecc9e1dSPeter Brune } 2770ecc9e1dSPeter Brune } 2780ecc9e1dSPeter Brune ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2790ecc9e1dSPeter Brune snes->linear_its += lits; 280b8840d6bSPeter Brune ierr = VecCopy(W, Y);CHKERRQ(ierr); 2810ecc9e1dSPeter Brune } else { 282b21d5a53SPeter Brune ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 2830ecc9e1dSPeter Brune } 28418aa0c0cSPeter Brune if (qn->singlereduction) { 285b8840d6bSPeter Brune ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 2868e231d97SPeter Brune } 2874b11644fSPeter Brune /* inward recursion starting at the first update and working forward */ 288bd052dfeSPeter Brune for (i = 0; i < l; i++) { 289b21d5a53SPeter Brune k = (it + i - l) % l; 29018aa0c0cSPeter Brune if (qn->singlereduction) { 2918e231d97SPeter Brune t = YtdX[k]; 2928e231d97SPeter Brune for (j = 0; j < i; j++) { 2938e231d97SPeter Brune g = (it + j - l) % l; 294b3cfccb2SPeter Brune t += (alpha[g] - beta[g])*H(g, k); 2958e231d97SPeter Brune } 2968e231d97SPeter Brune beta[k] = t / H(k, k); 2978e231d97SPeter Brune } else { 2986bf1b2e5SPeter Brune ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 2998e231d97SPeter Brune beta[k] = t / dXtdF[k]; 3008e231d97SPeter Brune } 30122d28d08SBarry Smith ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 30244f7e39eSPeter Brune if (qn->monitor) { 3033af51624SPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 3045ba6227bSPeter Brune ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 3053af51624SPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 30644f7e39eSPeter Brune } 3074b11644fSPeter Brune } 3084b11644fSPeter Brune PetscFunctionReturn(0); 3094b11644fSPeter Brune } 3104b11644fSPeter Brune 3114b11644fSPeter Brune #undef __FUNCT__ 3124b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 3134b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 3144b11644fSPeter Brune { 3154b11644fSPeter Brune PetscErrorCode ierr; 3169f83bee8SJed Brown SNES_QN *qn = (SNES_QN*) snes->data; 31715f5eeeaSPeter Brune Vec X,Xold; 318*0a3905e1SPeter Brune Vec F,W; 31947f26062SPeter Brune Vec Y,D,Dold; 320b8840d6bSPeter Brune PetscInt i, i_r; 321335fb975SPeter Brune PetscReal fnorm,xnorm,ynorm,gnorm; 3220c777b0cSPeter Brune PetscBool lssucceed,powell,periodic; 3231c6523dcSPeter Brune PetscScalar DolddotD,DolddotDold; 3240ecc9e1dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 325b7281c8aSPeter Brune SNESConvergedReason reason; 3260ecc9e1dSPeter Brune 3274b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 3284b11644fSPeter Brune 3296e111a19SKarl Rupp PetscFunctionBegin; 3309f3a0142SPeter Brune F = snes->vec_func; /* residual vector */ 3313af51624SPeter Brune Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 332*0a3905e1SPeter Brune W = snes->work[3]; 333b8840d6bSPeter Brune X = snes->vec_sol; /* solution vector */ 334335fb975SPeter Brune Xold = snes->work[0]; 3353af51624SPeter Brune 3363af51624SPeter Brune /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 337335fb975SPeter Brune D = snes->work[1]; 338335fb975SPeter Brune Dold = snes->work[2]; 3394b11644fSPeter Brune 3404b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 3414b11644fSPeter Brune 342ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3434b11644fSPeter Brune snes->iter = 0; 3444b11644fSPeter Brune snes->norm = 0.; 345ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 34647f26062SPeter Brune 347b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 34847f26062SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 349b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 350b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 351b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 352b7281c8aSPeter Brune PetscFunctionReturn(0); 353b7281c8aSPeter Brune } 35447f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 35547f26062SPeter Brune } else { 356e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 35715f5eeeaSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3584b11644fSPeter Brune if (snes->domainerror) { 3594b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3604b11644fSPeter Brune PetscFunctionReturn(0); 3614b11644fSPeter Brune } 3621aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 363c1c75074SPeter Brune 36447f26062SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 365189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 366189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 367189a9710SBarry Smith PetscFunctionReturn(0); 368189a9710SBarry Smith } 36947f26062SPeter Brune } 370b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 37147f26062SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 372b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 373b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 374b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 375b7281c8aSPeter Brune PetscFunctionReturn(0); 376b7281c8aSPeter Brune } 37747f26062SPeter Brune } else { 37847f26062SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 37947f26062SPeter Brune } 380b28a06ddSPeter Brune 381ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3824b11644fSPeter Brune snes->norm = fnorm; 383ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 384a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 3854b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 3864b11644fSPeter Brune 3874b11644fSPeter Brune /* test convergence */ 3884b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3894b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 39070d3b23bSPeter Brune 3913cf07b75SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 3923cf07b75SPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 3933cf07b75SPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 3943cf07b75SPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 395ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 396ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 397ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 398ddd40ce5SPeter Brune PetscFunctionReturn(0); 399ddd40ce5SPeter Brune } 400ddd40ce5SPeter Brune ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 4013cf07b75SPeter Brune ierr = VecCopy(F,D);CHKERRQ(ierr); 4023cf07b75SPeter Brune } 4033cf07b75SPeter Brune 404f8e15203SPeter Brune /* scale the initial update */ 4050c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 4060ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 4078cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 4080ecc9e1dSPeter Brune } 4090ecc9e1dSPeter Brune 4105ba6227bSPeter Brune for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 411*0a3905e1SPeter Brune if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) { 412*0a3905e1SPeter Brune PetscScalar ff,xf; 413*0a3905e1SPeter Brune ierr = VecCopy(Dold,Y);CHKERRQ(ierr); 414*0a3905e1SPeter Brune ierr = VecCopy(Xold,W);CHKERRQ(ierr); 415*0a3905e1SPeter Brune ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr); 416*0a3905e1SPeter Brune ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 417*0a3905e1SPeter Brune ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr); 418*0a3905e1SPeter Brune ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr); 419*0a3905e1SPeter Brune ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr); 420*0a3905e1SPeter Brune ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr); 421*0a3905e1SPeter Brune qn->scaling = PetscRealPart(xf)/PetscRealPart(ff); 422*0a3905e1SPeter Brune } 423b8840d6bSPeter Brune switch (qn->type) { 424b8840d6bSPeter Brune case SNES_QN_BADBROYDEN: 425b8840d6bSPeter Brune ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 426b8840d6bSPeter Brune break; 427b8840d6bSPeter Brune case SNES_QN_BROYDEN: 4285051edb1SPeter Brune ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr); 429b8840d6bSPeter Brune break; 430b8840d6bSPeter Brune case SNES_QN_LBFGS: 431b8840d6bSPeter Brune SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 432b8840d6bSPeter Brune break; 433b8840d6bSPeter Brune } 43470d3b23bSPeter Brune /* line search for lambda */ 43570d3b23bSPeter Brune ynorm = 1; gnorm = fnorm; 436*0a3905e1SPeter Brune ierr = VecCopy(D, Dold);CHKERRQ(ierr); 43715f5eeeaSPeter Brune ierr = VecCopy(X, Xold);CHKERRQ(ierr); 438f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 4399f3a0142SPeter Brune if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 4409f3a0142SPeter Brune if (snes->domainerror) { 4419f3a0142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 4429f3a0142SPeter Brune PetscFunctionReturn(0); 4439f3a0142SPeter Brune } 444f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 4459f3a0142SPeter Brune if (!lssucceed) { 4469f3a0142SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 4479f3a0142SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 4489f3a0142SPeter Brune break; 4499f3a0142SPeter Brune } 4509f3a0142SPeter Brune } 451f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 4520c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 45305ee83feSPeter Brune ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 4540ecc9e1dSPeter Brune } 4553af51624SPeter Brune 45688d374b2SPeter Brune /* convergence monitoring */ 457b21d5a53SPeter 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); 458b21d5a53SPeter Brune 459b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 460b28a06ddSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 461b28a06ddSPeter Brune ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr); 462b28a06ddSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr); 463ddd40ce5SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 464ddd40ce5SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 465ddd40ce5SPeter Brune snes->reason = SNES_DIVERGED_INNER; 466ddd40ce5SPeter Brune PetscFunctionReturn(0); 467ddd40ce5SPeter Brune } 468ddd40ce5SPeter Brune ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 469b28a06ddSPeter Brune } 470b28a06ddSPeter Brune 471360c497dSPeter Brune ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 472360c497dSPeter Brune ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 473360c497dSPeter Brune 474a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 4758409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 4764b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 477d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 4784b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 479b28a06ddSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 48047f26062SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 481b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 482b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 483b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 484b7281c8aSPeter Brune PetscFunctionReturn(0); 485b7281c8aSPeter Brune } 48688d374b2SPeter Brune } else { 48788d374b2SPeter Brune ierr = VecCopy(F, D);CHKERRQ(ierr); 48888d374b2SPeter Brune } 4890c777b0cSPeter Brune powell = PETSC_FALSE; 4900c777b0cSPeter Brune if (qn->restart_type == SNES_QN_RESTART_POWELL) { 4916bf1b2e5SPeter Brune /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 4928e231d97SPeter Brune ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4938e231d97SPeter Brune ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 4948e231d97SPeter Brune ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 4958e231d97SPeter Brune ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 4960c777b0cSPeter Brune if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 4970c777b0cSPeter Brune } 4980c777b0cSPeter Brune periodic = PETSC_FALSE; 499b8840d6bSPeter Brune if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 500b8840d6bSPeter Brune if (i_r>qn->m-1) periodic = PETSC_TRUE; 5010c777b0cSPeter Brune } 5020c777b0cSPeter Brune /* restart if either powell or periodic restart is satisfied. */ 5030c777b0cSPeter Brune if (powell || periodic) { 5045ba6227bSPeter Brune if (qn->monitor) { 5055ba6227bSPeter Brune ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 506516fe3c3SPeter 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); 5075ba6227bSPeter Brune ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 5085ba6227bSPeter Brune } 5095ba6227bSPeter Brune i_r = -1; 5105ba6227bSPeter Brune /* general purpose update */ 5115ba6227bSPeter Brune if (snes->ops->update) { 5125ba6227bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 5135ba6227bSPeter Brune } 5140c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5150ecc9e1dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 5168cba9625SJed Brown ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 5170ecc9e1dSPeter Brune } 5180ecc9e1dSPeter Brune } 51970d3b23bSPeter Brune /* general purpose update */ 52070d3b23bSPeter Brune if (snes->ops->update) { 52170d3b23bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 52270d3b23bSPeter Brune } 5235ba6227bSPeter Brune } 5244b11644fSPeter Brune if (i == snes->max_its) { 5254b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 5264b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 5274b11644fSPeter Brune } 5284b11644fSPeter Brune PetscFunctionReturn(0); 5294b11644fSPeter Brune } 5304b11644fSPeter Brune 5314b11644fSPeter Brune #undef __FUNCT__ 5324b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 5334b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 5344b11644fSPeter Brune { 5359f83bee8SJed Brown SNES_QN *qn = (SNES_QN*)snes->data; 5364b11644fSPeter Brune PetscErrorCode ierr; 537335fb975SPeter Brune 5384b11644fSPeter Brune PetscFunctionBegin; 539b8840d6bSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 54059e7931aSPeter Brune if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 5415c7a0a03SPeter Brune ierr = PetscMalloc4(qn->m, PetscScalar, &qn->alpha, 5428e231d97SPeter Brune qn->m, PetscScalar, &qn->beta, 5435c7a0a03SPeter Brune qn->m, PetscScalar, &qn->dXtdF, 5445c7a0a03SPeter Brune qn->m, PetscReal, &qn->lambda);CHKERRQ(ierr); 5458e231d97SPeter Brune 54618aa0c0cSPeter Brune if (qn->singlereduction) { 5478e231d97SPeter Brune ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 5488e231d97SPeter Brune qn->m, PetscScalar, &qn->dFtdX, 5498e231d97SPeter Brune qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 5508e231d97SPeter Brune } 551fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 552335fb975SPeter Brune 553335fb975SPeter Brune /* set up the line search */ 5540c777b0cSPeter Brune if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 5558e231d97SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 5568e231d97SPeter Brune } 5576c67d002SPeter Brune 5586c67d002SPeter Brune if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;} 5596c67d002SPeter Brune 5604b11644fSPeter Brune PetscFunctionReturn(0); 5614b11644fSPeter Brune } 5624b11644fSPeter Brune 5634b11644fSPeter Brune #undef __FUNCT__ 5644b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 5654b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 5664b11644fSPeter Brune { 5674b11644fSPeter Brune PetscErrorCode ierr; 5689f83bee8SJed Brown SNES_QN *qn; 5690adebc6cSBarry Smith 5704b11644fSPeter Brune PetscFunctionBegin; 5714b11644fSPeter Brune if (snes->data) { 5729f83bee8SJed Brown qn = (SNES_QN*)snes->data; 573b8840d6bSPeter Brune if (qn->U) { 574b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 5754b11644fSPeter Brune } 576b8840d6bSPeter Brune if (qn->V) { 577b8840d6bSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 5784b11644fSPeter Brune } 57918aa0c0cSPeter Brune if (qn->singlereduction) { 5808e231d97SPeter Brune ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 5818e231d97SPeter Brune } 5825c7a0a03SPeter Brune ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr); 5834b11644fSPeter Brune } 5844b11644fSPeter Brune PetscFunctionReturn(0); 5854b11644fSPeter Brune } 5864b11644fSPeter Brune 5874b11644fSPeter Brune #undef __FUNCT__ 5884b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 5894b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 5904b11644fSPeter Brune { 5914b11644fSPeter Brune PetscErrorCode ierr; 5926e111a19SKarl Rupp 5934b11644fSPeter Brune PetscFunctionBegin; 5944b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 5954b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 596bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr); 5974b11644fSPeter Brune PetscFunctionReturn(0); 5984b11644fSPeter Brune } 5994b11644fSPeter Brune 6004b11644fSPeter Brune #undef __FUNCT__ 6014b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 6024b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 6034b11644fSPeter Brune { 6044b11644fSPeter Brune 6054b11644fSPeter Brune PetscErrorCode ierr; 6062150357eSBarry Smith SNES_QN *qn = (SNES_QN*)snes->data; 60788f769c5SPeter Brune PetscBool monflg = PETSC_FALSE,flg; 608f1c6b773SPeter Brune SNESLineSearch linesearch; 6092150357eSBarry Smith SNESQNRestartType rtype = qn->restart_type; 6102150357eSBarry Smith SNESQNScaleType stype = qn->scale_type; 6112150357eSBarry Smith 6124b11644fSPeter Brune PetscFunctionBegin; 6134b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 6140298fd71SBarry Smith ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr); 6150298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr); 6160298fd71SBarry Smith ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr); 6170298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr); 6180298fd71SBarry Smith ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr); 61988f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 62088f769c5SPeter Brune if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 62188f769c5SPeter Brune 62288f769c5SPeter Brune ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 62388f769c5SPeter Brune if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 62488f769c5SPeter Brune 625b8840d6bSPeter Brune ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 6260298fd71SBarry Smith (PetscEnum)qn->type,(PetscEnum*)&qn->type,NULL);CHKERRQ(ierr); 6274b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 6289e764e56SPeter Brune if (!snes->linesearch) { 6297601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 6301a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 6319e764e56SPeter Brune } 63244f7e39eSPeter Brune if (monflg) { 633ce94432eSBarry Smith qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 63444f7e39eSPeter Brune } 6354b11644fSPeter Brune PetscFunctionReturn(0); 6364b11644fSPeter Brune } 6374b11644fSPeter Brune 63892c02d66SPeter Brune 6390c777b0cSPeter Brune #undef __FUNCT__ 6400c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType" 6410c777b0cSPeter Brune /*@ 6420c777b0cSPeter Brune SNESQNSetRestartType - Sets the restart type for SNESQN. 6430c777b0cSPeter Brune 6440c777b0cSPeter Brune Logically Collective on SNES 6450c777b0cSPeter Brune 6460c777b0cSPeter Brune Input Parameters: 6470c777b0cSPeter Brune + snes - the iterative context 6480c777b0cSPeter Brune - rtype - restart type 6490c777b0cSPeter Brune 6500c777b0cSPeter Brune Options Database: 6510c777b0cSPeter Brune + -snes_qn_restart_type<powell,periodic,none> - set the restart type 652b8840d6bSPeter Brune - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 6530c777b0cSPeter Brune 6540c777b0cSPeter Brune Level: intermediate 6550c777b0cSPeter Brune 6560c777b0cSPeter Brune SNESQNRestartTypes: 6570c777b0cSPeter Brune + SNES_QN_RESTART_NONE - never restart 6580c777b0cSPeter Brune . SNES_QN_RESTART_POWELL - restart based upon descent criteria 6590c777b0cSPeter Brune - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 6600c777b0cSPeter Brune 6610c777b0cSPeter Brune Notes: 6620c777b0cSPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 6630c777b0cSPeter Brune 6640c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 6650c777b0cSPeter Brune @*/ 6662150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 6672150357eSBarry Smith { 6680c777b0cSPeter Brune PetscErrorCode ierr; 6696e111a19SKarl Rupp 6700c777b0cSPeter Brune PetscFunctionBegin; 6710c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6720c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 6730c777b0cSPeter Brune PetscFunctionReturn(0); 6740c777b0cSPeter Brune } 6750c777b0cSPeter Brune 6760c777b0cSPeter Brune #undef __FUNCT__ 6770c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType" 6780c777b0cSPeter Brune /*@ 6790c777b0cSPeter Brune SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 6800c777b0cSPeter Brune 6810c777b0cSPeter Brune Logically Collective on SNES 6820c777b0cSPeter Brune 6830c777b0cSPeter Brune Input Parameters: 6840c777b0cSPeter Brune + snes - the iterative context 6850c777b0cSPeter Brune - stype - scale type 6860c777b0cSPeter Brune 6870c777b0cSPeter Brune Options Database: 6880c777b0cSPeter Brune . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 6890c777b0cSPeter Brune 6900c777b0cSPeter Brune Level: intermediate 6910c777b0cSPeter Brune 6920c777b0cSPeter Brune SNESQNSelectTypes: 6930c777b0cSPeter Brune + SNES_QN_SCALE_NONE - don't scale the problem 6940c777b0cSPeter Brune . SNES_QN_SCALE_SHANNO - use shanno scaling 6950c777b0cSPeter Brune . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 6960c777b0cSPeter Brune - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 6970c777b0cSPeter Brune 6980c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 6990c777b0cSPeter Brune @*/ 7000c777b0cSPeter Brune 7012150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 7022150357eSBarry Smith { 7030c777b0cSPeter Brune PetscErrorCode ierr; 7046e111a19SKarl Rupp 7050c777b0cSPeter Brune PetscFunctionBegin; 7060c777b0cSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7070c777b0cSPeter Brune ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 7080c777b0cSPeter Brune PetscFunctionReturn(0); 7090c777b0cSPeter Brune } 7100c777b0cSPeter Brune 7110c777b0cSPeter Brune #undef __FUNCT__ 7120c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN" 7130adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 7140adebc6cSBarry Smith { 7150c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7166e111a19SKarl Rupp 7170c777b0cSPeter Brune PetscFunctionBegin; 7180c777b0cSPeter Brune qn->scale_type = stype; 7190c777b0cSPeter Brune PetscFunctionReturn(0); 7200c777b0cSPeter Brune } 7210c777b0cSPeter Brune 7220c777b0cSPeter Brune #undef __FUNCT__ 7230c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN" 7240adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 7250adebc6cSBarry Smith { 7260c777b0cSPeter Brune SNES_QN *qn = (SNES_QN*)snes->data; 7276e111a19SKarl Rupp 7280c777b0cSPeter Brune PetscFunctionBegin; 7290c777b0cSPeter Brune qn->restart_type = rtype; 7300c777b0cSPeter Brune PetscFunctionReturn(0); 7310c777b0cSPeter Brune } 7320c777b0cSPeter Brune 7334b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 7344b11644fSPeter Brune /*MC 7354b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 7364b11644fSPeter Brune 7376cc8130cSPeter Brune Options Database: 7386cc8130cSPeter Brune 7396cc8130cSPeter Brune + -snes_qn_m - Number of past states saved for the L-Broyden methods. 7401867fe5bSPeter Brune . -snes_qn_powell_angle - Angle condition for restart. 7411867fe5bSPeter Brune . -snes_qn_powell_descent - Descent condition for restart. 74272835e02SPeter Brune . -snes_linesearch_type <cp, l2, basic> - Type of line search. 74344f7e39eSPeter Brune - -snes_qn_monitor - Monitors the quasi-newton jacobian. 7446cc8130cSPeter Brune 745b8840d6bSPeter Brune Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 746b8840d6bSPeter Brune previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 747b8840d6bSPeter Brune updates. 7486cc8130cSPeter Brune 7491867fe5bSPeter Brune When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 7501867fe5bSPeter Brune these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 7511867fe5bSPeter Brune iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 7521867fe5bSPeter Brune perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 7531867fe5bSPeter Brune 7546cc8130cSPeter Brune References: 7556cc8130cSPeter Brune 7560a8e8854SPeter Brune Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995. 7576cc8130cSPeter Brune 7580a8e8854SPeter Brune R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi-Newton Matrices and their use in Limited Memory Methods, 7590a8e8854SPeter Brune Technical Report, Northwestern University, June 1992. 7600a8e8854SPeter Brune 7610a8e8854SPeter Brune Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE 7620a8e8854SPeter Brune Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 763b8840d6bSPeter Brune 7644b11644fSPeter Brune 7654b11644fSPeter Brune Level: beginner 7664b11644fSPeter Brune 76704d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 7686cc8130cSPeter Brune 7694b11644fSPeter Brune M*/ 7704b11644fSPeter Brune #undef __FUNCT__ 7714b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 7728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 7734b11644fSPeter Brune { 7744b11644fSPeter Brune PetscErrorCode ierr; 7759f83bee8SJed Brown SNES_QN *qn; 7764b11644fSPeter Brune 7774b11644fSPeter Brune PetscFunctionBegin; 7784b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 7794b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 7804b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 7814b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 7824b11644fSPeter Brune snes->ops->view = 0; 7834b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 7844b11644fSPeter Brune 78547f26062SPeter Brune snes->pcside = PC_LEFT; 78647f26062SPeter Brune 78742f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 78842f4f86dSBarry Smith snes->usesksp = PETSC_FALSE; 78942f4f86dSBarry Smith 79088976e71SPeter Brune if (!snes->tolerancesset) { 7910e444f03SPeter Brune snes->max_funcs = 30000; 7920e444f03SPeter Brune snes->max_its = 10000; 79388976e71SPeter Brune } 7940e444f03SPeter Brune 7959f83bee8SJed Brown ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 7964b11644fSPeter Brune snes->data = (void*) qn; 7970ecc9e1dSPeter Brune qn->m = 10; 798b21d5a53SPeter Brune qn->scaling = 1.0; 7990298fd71SBarry Smith qn->U = NULL; 8000298fd71SBarry Smith qn->V = NULL; 8010298fd71SBarry Smith qn->dXtdF = NULL; 8020298fd71SBarry Smith qn->dFtdX = NULL; 8030298fd71SBarry Smith qn->dXdFmat = NULL; 8040298fd71SBarry Smith qn->monitor = NULL; 8051c6523dcSPeter Brune qn->singlereduction = PETSC_TRUE; 806b8840d6bSPeter Brune qn->powell_gamma = 0.9999; 8070c777b0cSPeter Brune qn->scale_type = SNES_QN_SCALE_SHANNO; 8080c777b0cSPeter Brune qn->restart_type = SNES_QN_RESTART_POWELL; 809b8840d6bSPeter Brune qn->type = SNES_QN_LBFGS; 810ea630c6eSPeter Brune 811bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr); 812bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr); 8134b11644fSPeter Brune PetscFunctionReturn(0); 8144b11644fSPeter Brune } 815