xref: /petsc/src/snes/impls/qn/qn.c (revision 07b62357f41014a58be1c7cccdd6f30e8aca3e32)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2fced5a79SAsbjørn Nilsen Riseth #include <petscdm.h>
34b11644fSPeter Brune 
48e231d97SPeter Brune #define H(i,j)  qn->dXdFmat[i*qn->m + j]
58e231d97SPeter Brune 
61efc8c45SPeter Brune const char *const SNESQNScaleTypes[] =        {"DEFAULT","NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0};
71efc8c45SPeter Brune const char *const SNESQNRestartTypes[] =      {"DEFAULT","NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0};
860c08b40SPeter Brune const char *const SNESQNTypes[] =             {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_",0};
9b8840d6bSPeter Brune 
104b11644fSPeter Brune typedef struct {
11b8840d6bSPeter Brune   Vec               *U;                   /* Stored past states (vary from method to method) */
12b8840d6bSPeter Brune   Vec               *V;                   /* Stored past states (vary from method to method) */
138e231d97SPeter Brune   PetscInt          m;                    /* The number of kept previous steps */
145c7a0a03SPeter Brune   PetscReal         *lambda;              /* The line search history of the method */
155c7a0a03SPeter Brune   PetscReal         *norm;                /* norms of the steps */
168e231d97SPeter Brune   PetscScalar       *alpha, *beta;
178e231d97SPeter Brune   PetscScalar       *dXtdF, *dFtdX, *YtdX;
1818aa0c0cSPeter Brune   PetscBool         singlereduction;      /* Aggregated reduction implementation */
198e231d97SPeter Brune   PetscScalar       *dXdFmat;             /* A matrix of values for dX_i dot dF_j */
2044f7e39eSPeter Brune   PetscViewer       monitor;
216bf1b2e5SPeter Brune   PetscReal         powell_gamma;         /* Powell angle restart condition */
22b21d5a53SPeter Brune   PetscReal         scaling;              /* scaling of H0 */
23b8840d6bSPeter Brune   SNESQNType        type;                 /* the type of quasi-newton method used */
2488f769c5SPeter Brune   SNESQNScaleType   scale_type;           /* the type of scaling used */
250c777b0cSPeter Brune   SNESQNRestartType restart_type;         /* determine the frequency and type of restart conditions */
269f83bee8SJed Brown } SNES_QN;
274b11644fSPeter Brune 
285051edb1SPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D)
292150357eSBarry Smith {
304b11644fSPeter Brune   PetscErrorCode     ierr;
319f83bee8SJed Brown   SNES_QN            *qn = (SNES_QN*)snes->data;
32b8840d6bSPeter Brune   Vec                W   = snes->work[3];
33b8840d6bSPeter Brune   Vec                *U  = qn->U;
34b8840d6bSPeter Brune   PetscInt           m = qn->m;
357f41f16fSJed Brown   PetscInt           k,i,j,l = m;
365c7a0a03SPeter Brune   PetscReal          unorm,a,b;
375c7a0a03SPeter Brune   PetscReal          *lambda=qn->lambda;
38b8840d6bSPeter Brune   PetscScalar        gdot;
3959e7931aSPeter Brune   PetscReal          udot;
402150357eSBarry Smith 
41b8840d6bSPeter Brune   PetscFunctionBegin;
42b8840d6bSPeter Brune   if (it < m) l = it;
435c7a0a03SPeter Brune   if (it > 0) {
445c7a0a03SPeter Brune     k = (it-1)%l;
455c7a0a03SPeter Brune     ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr);
46a49fa0d8SPeter Brune     ierr = VecCopy(Xold,U[k]);CHKERRQ(ierr);
475051edb1SPeter Brune     ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr);
485c7a0a03SPeter Brune     if (qn->monitor) {
495c7a0a03SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
506bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %D of %D by lambda: %14.12e \n",k,l,(double)lambda[k]);CHKERRQ(ierr);
515c7a0a03SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
525c7a0a03SPeter Brune     }
535c7a0a03SPeter Brune   }
54b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
55d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr);
56422a814eSBarry Smith     SNESCheckKSPSolve(snes);
57b8840d6bSPeter Brune     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
58b8840d6bSPeter Brune   } else {
59b8840d6bSPeter Brune     ierr = VecCopy(D,Y);CHKERRQ(ierr);
60b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
61b8840d6bSPeter Brune   }
62b8840d6bSPeter Brune 
63b8840d6bSPeter Brune   /* inward recursion starting at the first update and working forward */
64b8840d6bSPeter Brune   for (i = 0; i < l-1; i++) {
655c7a0a03SPeter Brune     j = (it+i-l)%l;
665c7a0a03SPeter Brune     k = (it+i-l+1)%l;
675c7a0a03SPeter Brune     ierr = VecNorm(U[j],NORM_2,&unorm);CHKERRQ(ierr);
685c7a0a03SPeter Brune     ierr = VecDot(U[j],Y,&gdot);CHKERRQ(ierr);
6959e7931aSPeter Brune     unorm *= unorm;
7059e7931aSPeter Brune     udot = PetscRealPart(gdot);
715c7a0a03SPeter Brune     a = (lambda[j]/lambda[k]);
725c7a0a03SPeter Brune     b = -(1.-lambda[j]);
7359e7931aSPeter Brune     a *= udot/unorm;
7459e7931aSPeter Brune     b *= udot/unorm;
755c7a0a03SPeter Brune     ierr = VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);CHKERRQ(ierr);
765c7a0a03SPeter Brune 
77b8840d6bSPeter Brune     if (qn->monitor) {
78b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
796bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %D and %D, gdot: %14.12e\n",k,j,(double)PetscRealPart(gdot));CHKERRQ(ierr);
80b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
81b8840d6bSPeter Brune     }
82b8840d6bSPeter Brune   }
83b8840d6bSPeter Brune   if (l > 0) {
84b8840d6bSPeter Brune     k = (it-1)%l;
85b8840d6bSPeter Brune     ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
865c7a0a03SPeter Brune     ierr = VecNorm(U[k],NORM_2,&unorm);CHKERRQ(ierr);
8759e7931aSPeter Brune     unorm *= unorm;
8859e7931aSPeter Brune     udot = PetscRealPart(gdot);
8959e7931aSPeter Brune     a = unorm/(unorm-lambda[k]*udot);
9059e7931aSPeter Brune     b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot);
91b8840d6bSPeter Brune     if (qn->monitor) {
92b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
936bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %D: a: %14.12e b: %14.12e \n",k,(double)a,(double)b);CHKERRQ(ierr);
94b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
95b8840d6bSPeter Brune     }
965c7a0a03SPeter Brune     ierr = VecAXPBY(Y,b,a,U[k]);CHKERRQ(ierr);
97b8840d6bSPeter Brune   }
985c7a0a03SPeter Brune   l = m;
995c7a0a03SPeter Brune   if (it+1<m)l=it+1;
1005c7a0a03SPeter Brune   k = it%l;
1015c7a0a03SPeter Brune   if (qn->monitor) {
1025c7a0a03SPeter Brune     ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1036bdcc836SBarry Smith     ierr = PetscViewerASCIIPrintf(qn->monitor, "setting vector %D of %D\n",k,l);CHKERRQ(ierr);
1045c7a0a03SPeter Brune     ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1055c7a0a03SPeter Brune   }
106b8840d6bSPeter Brune   PetscFunctionReturn(0);
107b8840d6bSPeter Brune }
108b8840d6bSPeter Brune 
1090adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1100adebc6cSBarry Smith {
111b8840d6bSPeter Brune   PetscErrorCode ierr;
112b8840d6bSPeter Brune   SNES_QN        *qn = (SNES_QN*)snes->data;
113b8840d6bSPeter Brune   Vec            W   = snes->work[3];
114b8840d6bSPeter Brune   Vec            *U  = qn->U;
115b8840d6bSPeter Brune   Vec            *T  = qn->V;
116b8840d6bSPeter Brune 
11784c577daSBarry Smith   /* ksp thing for Jacobian scaling */
1187f41f16fSJed Brown   PetscInt           h,k,j,i;
119b8840d6bSPeter Brune   PetscInt           m = qn->m;
1205051edb1SPeter Brune   PetscScalar        gdot,udot;
121b8840d6bSPeter Brune   PetscInt           l = m;
1220adebc6cSBarry Smith 
123b8840d6bSPeter Brune   PetscFunctionBegin;
124b8840d6bSPeter Brune   if (it < m) l = it;
125b8840d6bSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
126b8840d6bSPeter Brune   if (l > 0) {
127b8840d6bSPeter Brune     k    = (it-1)%l;
128c9e59058SPeter Brune     ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);CHKERRQ(ierr);
129b8840d6bSPeter Brune     ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr);
130b8840d6bSPeter Brune     ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr);
1315051edb1SPeter Brune     ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr);
1325051edb1SPeter Brune     ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr);
133b8840d6bSPeter Brune   }
134b8840d6bSPeter Brune 
135b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
136d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
137422a814eSBarry Smith     SNESCheckKSPSolve(snes);
138b8840d6bSPeter Brune     ierr = VecCopy(W,Y);CHKERRQ(ierr);
139b8840d6bSPeter Brune   } else {
140b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
141b8840d6bSPeter Brune   }
1425051edb1SPeter Brune 
1435051edb1SPeter Brune   /* inward recursion starting at the first update and working forward */
1445051edb1SPeter Brune   if (l > 0) {
1455051edb1SPeter Brune     for (i = 0; i < l-1; i++) {
146c9e59058SPeter Brune       j    = (it+i-l)%l;
147c9e59058SPeter Brune       k    = (it+i-l+1)%l;
148c9e59058SPeter Brune       h    = (it-1)%l;
149c9e59058SPeter Brune       ierr = VecDotBegin(U[j],U[h],&gdot);CHKERRQ(ierr);
150c9e59058SPeter Brune       ierr = VecDotBegin(U[j],U[j],&udot);CHKERRQ(ierr);
151c9e59058SPeter Brune       ierr = VecDotEnd(U[j],U[h],&gdot);CHKERRQ(ierr);
152c9e59058SPeter Brune       ierr = VecDotEnd(U[j],U[j],&udot);CHKERRQ(ierr);
1535051edb1SPeter Brune       ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr);
154c9e59058SPeter Brune       ierr = VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);CHKERRQ(ierr);
1555051edb1SPeter Brune       if (qn->monitor) {
1565051edb1SPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1576bdcc836SBarry Smith         ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D gdot: %14.12e\n", it, k, (double)PetscRealPart(gdot));CHKERRQ(ierr);
1585051edb1SPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1595051edb1SPeter Brune       }
1605051edb1SPeter Brune     }
1615051edb1SPeter Brune   }
162b8840d6bSPeter Brune   PetscFunctionReturn(0);
163b8840d6bSPeter Brune }
164b8840d6bSPeter Brune 
1650adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1660adebc6cSBarry Smith {
167b8840d6bSPeter Brune   PetscErrorCode ierr;
168b8840d6bSPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
169b8840d6bSPeter Brune   Vec            W      = snes->work[3];
170b8840d6bSPeter Brune   Vec            *dX    = qn->U;
171b8840d6bSPeter Brune   Vec            *dF    = qn->V;
172bd052dfeSPeter Brune   PetscScalar    *alpha = qn->alpha;
173bd052dfeSPeter Brune   PetscScalar    *beta  = qn->beta;
1748e231d97SPeter Brune   PetscScalar    *dXtdF = qn->dXtdF;
175b8840d6bSPeter Brune   PetscScalar    *dFtdX = qn->dFtdX;
1768e231d97SPeter Brune   PetscScalar    *YtdX  = qn->YtdX;
177bd052dfeSPeter Brune 
17884c577daSBarry Smith   /* ksp thing for Jacobian scaling */
1797f41f16fSJed Brown   PetscInt           k,i,j,g;
1804b11644fSPeter Brune   PetscInt           m = qn->m;
1814b11644fSPeter Brune   PetscScalar        t;
1824b11644fSPeter Brune   PetscInt           l = m;
1838e231d97SPeter Brune 
1844b11644fSPeter Brune   PetscFunctionBegin;
1855ba6227bSPeter Brune   if (it < m) l = it;
1861c6523dcSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
187b8840d6bSPeter Brune   if (it > 0) {
188b8840d6bSPeter Brune     k    = (it - 1) % l;
189b8840d6bSPeter Brune     ierr = VecCopy(D,dF[k]);CHKERRQ(ierr);
190b8840d6bSPeter Brune     ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
191b8840d6bSPeter Brune     ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
192b8840d6bSPeter Brune     ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
19318aa0c0cSPeter Brune     if (qn->singlereduction) {
1941c6523dcSPeter Brune       ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
1951c6523dcSPeter Brune       ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
1961c6523dcSPeter Brune       ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr);
1971c6523dcSPeter Brune       ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
1981c6523dcSPeter Brune       ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
1991c6523dcSPeter Brune       ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr);
200b8840d6bSPeter Brune       for (j = 0; j < l; j++) {
201b8840d6bSPeter Brune         H(k, j) = dFtdX[j];
202b8840d6bSPeter Brune         H(j, k) = dXtdF[j];
203b8840d6bSPeter Brune       }
204b8840d6bSPeter Brune       /* copy back over to make the computation of alpha and beta easier */
2051aa26658SKarl Rupp       for (j = 0; j < l; j++) dXtdF[j] = H(j, j);
206b8840d6bSPeter Brune     } else {
207b8840d6bSPeter Brune       ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr);
20823d44fbcSPeter Brune     }
20923d44fbcSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
210b8840d6bSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr);
211b8840d6bSPeter Brune     }
212b8840d6bSPeter Brune   }
213b8840d6bSPeter Brune 
2144b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
2154b11644fSPeter Brune   for (i=0; i<l; i++) {
216b21d5a53SPeter Brune     k = (it-i-1)%l;
21718aa0c0cSPeter Brune     if (qn->singlereduction) {
2188e231d97SPeter Brune       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
2198e231d97SPeter Brune       t = YtdX[k];
2208e231d97SPeter Brune       for (j=0; j<i; j++) {
2218e231d97SPeter Brune         g  = (it-j-1)%l;
222b3cfccb2SPeter Brune         t -= alpha[g]*H(k, g);
2238e231d97SPeter Brune       }
2248e231d97SPeter Brune       alpha[k] = t/H(k,k);
2258e231d97SPeter Brune     } else {
2263af51624SPeter Brune       ierr     = VecDot(dX[k],Y,&t);CHKERRQ(ierr);
2278e231d97SPeter Brune       alpha[k] = t/dXtdF[k];
2288e231d97SPeter Brune     }
22944f7e39eSPeter Brune     if (qn->monitor) {
2303af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2316bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D alpha:        %14.12e\n", it, k, (double)PetscRealPart(alpha[k]));CHKERRQ(ierr);
2323af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
23344f7e39eSPeter Brune     }
2346bf1b2e5SPeter Brune     ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr);
2354b11644fSPeter Brune   }
2364b11644fSPeter Brune 
2370c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
238d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
239422a814eSBarry Smith     SNESCheckKSPSolve(snes);
240b8840d6bSPeter Brune     ierr              = VecCopy(W, Y);CHKERRQ(ierr);
2410ecc9e1dSPeter Brune   } else {
242b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
2430ecc9e1dSPeter Brune   }
24418aa0c0cSPeter Brune   if (qn->singlereduction) {
245b8840d6bSPeter Brune     ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr);
2468e231d97SPeter Brune   }
2474b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
248bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
249b21d5a53SPeter Brune     k = (it + i - l) % l;
25018aa0c0cSPeter Brune     if (qn->singlereduction) {
2518e231d97SPeter Brune       t = YtdX[k];
2528e231d97SPeter Brune       for (j = 0; j < i; j++) {
2538e231d97SPeter Brune         g  = (it + j - l) % l;
254b3cfccb2SPeter Brune         t += (alpha[g] - beta[g])*H(g, k);
2558e231d97SPeter Brune       }
2568e231d97SPeter Brune       beta[k] = t / H(k, k);
2578e231d97SPeter Brune     } else {
2586bf1b2e5SPeter Brune       ierr    = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
2598e231d97SPeter Brune       beta[k] = t / dXtdF[k];
2608e231d97SPeter Brune     }
26122d28d08SBarry Smith     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr);
26244f7e39eSPeter Brune     if (qn->monitor) {
2633af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2646bdcc836SBarry Smith       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D alpha - beta: %14.12e\n", it, k, (double)PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
2653af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
26644f7e39eSPeter Brune     }
2674b11644fSPeter Brune   }
2684b11644fSPeter Brune   PetscFunctionReturn(0);
2694b11644fSPeter Brune }
2704b11644fSPeter Brune 
2714b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
2724b11644fSPeter Brune {
2734b11644fSPeter Brune   PetscErrorCode       ierr;
2749f83bee8SJed Brown   SNES_QN              *qn = (SNES_QN*) snes->data;
27515f5eeeaSPeter Brune   Vec                  X,Xold;
2760a3905e1SPeter Brune   Vec                  F,W;
27747f26062SPeter Brune   Vec                  Y,D,Dold;
278b8840d6bSPeter Brune   PetscInt             i, i_r;
279335fb975SPeter Brune   PetscReal            fnorm,xnorm,ynorm,gnorm;
280422a814eSBarry Smith   SNESLineSearchReason lssucceed;
281422a814eSBarry Smith   PetscBool            powell,periodic;
2821c6523dcSPeter Brune   PetscScalar          DolddotD,DolddotDold;
283b7281c8aSPeter Brune   SNESConvergedReason  reason;
2840ecc9e1dSPeter Brune 
28584c577daSBarry Smith   /* basically just a regular newton's method except for the application of the Jacobian */
2864b11644fSPeter Brune 
2876e111a19SKarl Rupp   PetscFunctionBegin;
2886c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
289c579b300SPatrick Farrell 
290fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
2919f3a0142SPeter Brune   F    = snes->vec_func;                /* residual vector */
2923af51624SPeter Brune   Y    = snes->vec_sol_update;          /* search direction generated by J^-1D*/
2930a3905e1SPeter Brune   W    = snes->work[3];
294b8840d6bSPeter Brune   X    = snes->vec_sol;                 /* solution vector */
295335fb975SPeter Brune   Xold = snes->work[0];
2963af51624SPeter Brune 
2973af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
298335fb975SPeter Brune   D    = snes->work[1];
299335fb975SPeter Brune   Dold = snes->work[2];
3004b11644fSPeter Brune 
3014b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
3024b11644fSPeter Brune 
303e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3044b11644fSPeter Brune   snes->iter = 0;
3054b11644fSPeter Brune   snes->norm = 0.;
306e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
30747f26062SPeter Brune 
308efd4aadfSBarry Smith   if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
309be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
310efd4aadfSBarry Smith     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
311b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
312b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
313b7281c8aSPeter Brune       PetscFunctionReturn(0);
314b7281c8aSPeter Brune     }
31547f26062SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
31647f26062SPeter Brune   } else {
317e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
31815f5eeeaSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3191aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
320c1c75074SPeter Brune 
32147f26062SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
322422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
32347f26062SPeter Brune   }
324efd4aadfSBarry Smith   if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
325be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr);
326efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
327b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
328b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
329b7281c8aSPeter Brune         PetscFunctionReturn(0);
330b7281c8aSPeter Brune       }
33147f26062SPeter Brune   } else {
33247f26062SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
33347f26062SPeter Brune   }
334b28a06ddSPeter Brune 
335e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3364b11644fSPeter Brune   snes->norm = fnorm;
337e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
338a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
3394b11644fSPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
3404b11644fSPeter Brune 
3414b11644fSPeter Brune   /* test convergence */
3424b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3434b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
34470d3b23bSPeter Brune 
345efd4aadfSBarry Smith   if (snes->npc && snes->npcside== PC_RIGHT) {
346efd4aadfSBarry Smith     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,X,0,0);CHKERRQ(ierr);
347efd4aadfSBarry Smith     ierr = SNESSolve(snes->npc,snes->vec_rhs,X);CHKERRQ(ierr);
348efd4aadfSBarry Smith     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,X,0,0);CHKERRQ(ierr);
349efd4aadfSBarry Smith     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
350ddd40ce5SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
351ddd40ce5SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
352ddd40ce5SPeter Brune       PetscFunctionReturn(0);
353ddd40ce5SPeter Brune     }
354be95d8f1SBarry Smith     ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
3553cf07b75SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
3563cf07b75SPeter Brune   }
3573cf07b75SPeter Brune 
358f8e15203SPeter Brune   /* scale the initial update */
3590c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
360d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
361*07b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
362aeb49b86SAsbjørn Nilsen Riseth     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
3630ecc9e1dSPeter Brune   }
3640ecc9e1dSPeter Brune 
3655ba6227bSPeter Brune   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
3660a3905e1SPeter Brune     if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) {
3670a3905e1SPeter Brune       PetscScalar ff,xf;
3680a3905e1SPeter Brune       ierr = VecCopy(Dold,Y);CHKERRQ(ierr);
3690a3905e1SPeter Brune       ierr = VecCopy(Xold,W);CHKERRQ(ierr);
3700a3905e1SPeter Brune       ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr);
3710a3905e1SPeter Brune       ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr);
3720a3905e1SPeter Brune       ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr);
3730a3905e1SPeter Brune       ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr);
3740a3905e1SPeter Brune       ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr);
3750a3905e1SPeter Brune       ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr);
3760a3905e1SPeter Brune       qn->scaling = PetscRealPart(xf)/PetscRealPart(ff);
3776bdcc836SBarry Smith       if (qn->monitor) {
3786bdcc836SBarry Smith         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
3796bdcc836SBarry Smith         ierr = PetscViewerASCIIPrintf(qn->monitor, "Shanno scaling %D %g\n", i,(double)qn->scaling);CHKERRQ(ierr);
3806bdcc836SBarry Smith         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
3816bdcc836SBarry Smith       }
3820a3905e1SPeter Brune     }
383b8840d6bSPeter Brune     switch (qn->type) {
384b8840d6bSPeter Brune     case SNES_QN_BADBROYDEN:
385b8840d6bSPeter Brune       ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
386b8840d6bSPeter Brune       break;
387b8840d6bSPeter Brune     case SNES_QN_BROYDEN:
3885051edb1SPeter Brune       ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr);
389b8840d6bSPeter Brune       break;
390b8840d6bSPeter Brune     case SNES_QN_LBFGS:
391b8840d6bSPeter Brune       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
392b8840d6bSPeter Brune       break;
393b8840d6bSPeter Brune     }
39470d3b23bSPeter Brune     /* line search for lambda */
39570d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
3960a3905e1SPeter Brune     ierr  = VecCopy(D, Dold);CHKERRQ(ierr);
39715f5eeeaSPeter Brune     ierr  = VecCopy(X, Xold);CHKERRQ(ierr);
398f1c6b773SPeter Brune     ierr  = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
3999f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
400422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lssucceed);CHKERRQ(ierr);
401422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
402422a814eSBarry Smith     if (lssucceed) {
4039f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
4049f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
4059f3a0142SPeter Brune         break;
4069f3a0142SPeter Brune       }
4079f3a0142SPeter Brune     }
4080c777b0cSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
40905ee83feSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
4100ecc9e1dSPeter Brune     }
4113af51624SPeter Brune 
41288d374b2SPeter Brune     /* convergence monitoring */
413b21d5a53SPeter 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);
414b21d5a53SPeter Brune 
415efd4aadfSBarry Smith     if (snes->npc && snes->npcside== PC_RIGHT) {
416efd4aadfSBarry Smith       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,X,0,0);CHKERRQ(ierr);
417efd4aadfSBarry Smith       ierr = SNESSolve(snes->npc,snes->vec_rhs,X);CHKERRQ(ierr);
418efd4aadfSBarry Smith       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,X,0,0);CHKERRQ(ierr);
419efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
420ddd40ce5SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
421ddd40ce5SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
422ddd40ce5SPeter Brune         PetscFunctionReturn(0);
423ddd40ce5SPeter Brune       }
424be95d8f1SBarry Smith       ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
425b28a06ddSPeter Brune     }
426b28a06ddSPeter Brune 
427360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
42871dbe336SPeter Brune     snes->norm = fnorm;
429c1e67a49SFande Kong     snes->xnorm = xnorm;
430c1e67a49SFande Kong     snes->ynorm = ynorm;
431360c497dSPeter Brune 
432a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
4338409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
4344b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
435d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
4364b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
437efd4aadfSBarry Smith     if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
438be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr);
439efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
440b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
441b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
442b7281c8aSPeter Brune         PetscFunctionReturn(0);
443b7281c8aSPeter Brune       }
44488d374b2SPeter Brune     } else {
44588d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
44688d374b2SPeter Brune     }
4470c777b0cSPeter Brune     powell = PETSC_FALSE;
4486bdcc836SBarry Smith     if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) {
4496bdcc836SBarry Smith       /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */
45023c918e7SPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
45123c918e7SPeter Brune         ierr = MatMult(snes->jacobian_pre,Dold,W);CHKERRQ(ierr);
45223c918e7SPeter Brune       } else {
45323c918e7SPeter Brune         ierr = VecCopy(Dold,W);CHKERRQ(ierr);
45423c918e7SPeter Brune       }
45523c918e7SPeter Brune       ierr = VecDotBegin(W, Dold, &DolddotDold);CHKERRQ(ierr);
45623c918e7SPeter Brune       ierr = VecDotBegin(W, D, &DolddotD);CHKERRQ(ierr);
45723c918e7SPeter Brune       ierr = VecDotEnd(W, Dold, &DolddotDold);CHKERRQ(ierr);
45823c918e7SPeter Brune       ierr = VecDotEnd(W, D, &DolddotD);CHKERRQ(ierr);
4590c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
4600c777b0cSPeter Brune     }
4610c777b0cSPeter Brune     periodic = PETSC_FALSE;
462b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
463b8840d6bSPeter Brune       if (i_r>qn->m-1) periodic = PETSC_TRUE;
4640c777b0cSPeter Brune     }
4650c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
4660c777b0cSPeter Brune     if (powell || periodic) {
4675ba6227bSPeter Brune       if (qn->monitor) {
4685ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
4696bdcc836SBarry Smith         if (powell) {
4706bdcc836SBarry Smith           ierr = PetscViewerASCIIPrintf(qn->monitor, "Powell restart! |%14.12e| > %6.4f*|%14.12e| i_r = %D\n", (double)PetscRealPart(DolddotD), (double)qn->powell_gamma, (double)PetscRealPart(DolddotDold),i_r);CHKERRQ(ierr);
4716bdcc836SBarry Smith         } else {
4726bdcc836SBarry Smith           ierr = PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %D\n", i_r);CHKERRQ(ierr);
4736bdcc836SBarry Smith         }
4745ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
4755ba6227bSPeter Brune       }
4765ba6227bSPeter Brune       i_r = -1;
4775ba6227bSPeter Brune       /* general purpose update */
4785ba6227bSPeter Brune       if (snes->ops->update) {
4795ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
4805ba6227bSPeter Brune       }
4810c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
482d1e9a80fSBarry Smith         ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
483*07b62357SFande Kong         SNESCheckJacobianDomainerror(snes);
4840ecc9e1dSPeter Brune       }
4850ecc9e1dSPeter Brune     }
48670d3b23bSPeter Brune     /* general purpose update */
48770d3b23bSPeter Brune     if (snes->ops->update) {
48870d3b23bSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
48970d3b23bSPeter Brune     }
4905ba6227bSPeter Brune   }
4914b11644fSPeter Brune   if (i == snes->max_its) {
4924b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
4934b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
4944b11644fSPeter Brune   }
4954b11644fSPeter Brune   PetscFunctionReturn(0);
4964b11644fSPeter Brune }
4974b11644fSPeter Brune 
4984b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
4994b11644fSPeter Brune {
5009f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
5014b11644fSPeter Brune   PetscErrorCode ierr;
502fced5a79SAsbjørn Nilsen Riseth   DM             dm;
503335fb975SPeter Brune 
5044b11644fSPeter Brune   PetscFunctionBegin;
505fced5a79SAsbjørn Nilsen Riseth 
506fced5a79SAsbjørn Nilsen Riseth   if (!snes->vec_sol) {
507fced5a79SAsbjørn Nilsen Riseth     ierr             = SNESGetDM(snes,&dm);CHKERRQ(ierr);
508fced5a79SAsbjørn Nilsen Riseth     ierr             = DMCreateGlobalVector(dm,&snes->vec_sol);CHKERRQ(ierr);
509fced5a79SAsbjørn Nilsen Riseth   }
510fced5a79SAsbjørn Nilsen Riseth 
511b8840d6bSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr);
51259e7931aSPeter Brune   if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr);
513dcca6d9dSJed Brown   ierr = PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);CHKERRQ(ierr);
5148e231d97SPeter Brune 
51518aa0c0cSPeter Brune   if (qn->singlereduction) {
516dcca6d9dSJed Brown     ierr = PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);CHKERRQ(ierr);
5178e231d97SPeter Brune   }
518fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
51960c08b40SPeter Brune   /* set method defaults */
5201efc8c45SPeter Brune   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
52160c08b40SPeter Brune     if (qn->type == SNES_QN_BADBROYDEN) {
52260c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_NONE;
52360c08b40SPeter Brune     } else {
52460c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_SHANNO;
52560c08b40SPeter Brune     }
52660c08b40SPeter Brune   }
5271efc8c45SPeter Brune   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
52860c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
52960c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_POWELL;
53060c08b40SPeter Brune     } else {
53160c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_PERIODIC;
53260c08b40SPeter Brune     }
53360c08b40SPeter Brune   }
53460c08b40SPeter Brune 
5350c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
5368e231d97SPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
5378e231d97SPeter Brune   }
538efd4aadfSBarry Smith   if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;}
5394b11644fSPeter Brune   PetscFunctionReturn(0);
5404b11644fSPeter Brune }
5414b11644fSPeter Brune 
5424b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
5434b11644fSPeter Brune {
5444b11644fSPeter Brune   PetscErrorCode ierr;
5459f83bee8SJed Brown   SNES_QN        *qn;
5460adebc6cSBarry Smith 
5474b11644fSPeter Brune   PetscFunctionBegin;
5484b11644fSPeter Brune   if (snes->data) {
5499f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
550b8840d6bSPeter Brune     if (qn->U) {
551b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr);
5524b11644fSPeter Brune     }
553b8840d6bSPeter Brune     if (qn->V) {
554b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr);
5554b11644fSPeter Brune     }
55618aa0c0cSPeter Brune     if (qn->singlereduction) {
5578e231d97SPeter Brune       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
5588e231d97SPeter Brune     }
5595c7a0a03SPeter Brune     ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr);
5604b11644fSPeter Brune   }
5614b11644fSPeter Brune   PetscFunctionReturn(0);
5624b11644fSPeter Brune }
5634b11644fSPeter Brune 
5644b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
5654b11644fSPeter Brune {
5664b11644fSPeter Brune   PetscErrorCode ierr;
5676e111a19SKarl Rupp 
5684b11644fSPeter Brune   PetscFunctionBegin;
5694b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
5704b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
571bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
5724b11644fSPeter Brune   PetscFunctionReturn(0);
5734b11644fSPeter Brune }
5744b11644fSPeter Brune 
5754416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_QN(PetscOptionItems *PetscOptionsObject,SNES snes)
5764b11644fSPeter Brune {
5774b11644fSPeter Brune 
5784b11644fSPeter Brune   PetscErrorCode    ierr;
5792150357eSBarry Smith   SNES_QN           *qn    = (SNES_QN*)snes->data;
58088f769c5SPeter Brune   PetscBool         monflg = PETSC_FALSE,flg;
581f1c6b773SPeter Brune   SNESLineSearch    linesearch;
5822150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
5832150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
5841efc8c45SPeter Brune   SNESQNType        qtype = qn->type;
5852150357eSBarry Smith 
5864b11644fSPeter Brune   PetscFunctionBegin;
587e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES QN options");CHKERRQ(ierr);
5880298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr);
5890298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr);
5900298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr);
5910298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr);
59288f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr);
59388f769c5SPeter Brune   if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr);
59488f769c5SPeter Brune 
59588f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr);
59688f769c5SPeter Brune   if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr);
59788f769c5SPeter Brune 
5980fdccdaeSBarry Smith   ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg);CHKERRQ(ierr);
5991efc8c45SPeter Brune   if (flg) {ierr = SNESQNSetType(snes,qtype);CHKERRQ(ierr);}
6004b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
6019e764e56SPeter Brune   if (!snes->linesearch) {
6027601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
60360c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
6041a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
60560c08b40SPeter Brune     } else if (qn->type == SNES_QN_BROYDEN) {
60660c08b40SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);CHKERRQ(ierr);
60760c08b40SPeter Brune     } else {
60860c08b40SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
60960c08b40SPeter Brune     }
6109e764e56SPeter Brune   }
61144f7e39eSPeter Brune   if (monflg) {
612ce94432eSBarry Smith     qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
61344f7e39eSPeter Brune   }
6144b11644fSPeter Brune   PetscFunctionReturn(0);
6154b11644fSPeter Brune }
6164b11644fSPeter Brune 
6175cd83419SPeter Brune static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
6185cd83419SPeter Brune {
6195cd83419SPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
6205cd83419SPeter Brune   PetscBool      iascii;
6215cd83419SPeter Brune   PetscErrorCode ierr;
6225cd83419SPeter Brune 
6235cd83419SPeter Brune   PetscFunctionBegin;
6245cd83419SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
6255cd83419SPeter Brune   if (iascii) {
626efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, restart type is %s, scale type is %s\n",SNESQNTypes[qn->type],SNESQNRestartTypes[qn->restart_type],SNESQNScaleTypes[qn->scale_type]);CHKERRQ(ierr);
6276bdcc836SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Stored subspace size: %D\n", qn->m);CHKERRQ(ierr);
6285cd83419SPeter Brune     if (qn->singlereduction) {
6295cd83419SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  Using the single reduction variant.\n");CHKERRQ(ierr);
6305cd83419SPeter Brune     }
6315cd83419SPeter Brune   }
6325cd83419SPeter Brune   PetscFunctionReturn(0);
6335cd83419SPeter Brune }
63492c02d66SPeter Brune 
6350c777b0cSPeter Brune /*@
6360c777b0cSPeter Brune     SNESQNSetRestartType - Sets the restart type for SNESQN.
6370c777b0cSPeter Brune 
6380c777b0cSPeter Brune     Logically Collective on SNES
6390c777b0cSPeter Brune 
6400c777b0cSPeter Brune     Input Parameters:
6410c777b0cSPeter Brune +   snes - the iterative context
6420c777b0cSPeter Brune -   rtype - restart type
6430c777b0cSPeter Brune 
6440c777b0cSPeter Brune     Options Database:
6450c777b0cSPeter Brune +   -snes_qn_restart_type <powell,periodic,none> - set the restart type
64684c577daSBarry Smith -   -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
6470c777b0cSPeter Brune 
6480c777b0cSPeter Brune     Level: intermediate
6490c777b0cSPeter Brune 
6500c777b0cSPeter Brune     SNESQNRestartTypes:
6510c777b0cSPeter Brune +   SNES_QN_RESTART_NONE - never restart
6520c777b0cSPeter Brune .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
6530c777b0cSPeter Brune -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
6540c777b0cSPeter Brune 
65584c577daSBarry Smith .keywords: SNES, SNESQN, restart, type, set SNESLineSearch, SNESQNRestartType
6560c777b0cSPeter Brune @*/
6572150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
6582150357eSBarry Smith {
6590c777b0cSPeter Brune   PetscErrorCode ierr;
6606e111a19SKarl Rupp 
6610c777b0cSPeter Brune   PetscFunctionBegin;
6620c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6630c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
6640c777b0cSPeter Brune   PetscFunctionReturn(0);
6650c777b0cSPeter Brune }
6660c777b0cSPeter Brune 
6670c777b0cSPeter Brune /*@
66884c577daSBarry Smith     SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN.
6690c777b0cSPeter Brune 
6700c777b0cSPeter Brune     Logically Collective on SNES
6710c777b0cSPeter Brune 
6720c777b0cSPeter Brune     Input Parameters:
6730c777b0cSPeter Brune +   snes - the iterative context
6740c777b0cSPeter Brune -   stype - scale type
6750c777b0cSPeter Brune 
6760c777b0cSPeter Brune     Options Database:
6770c777b0cSPeter Brune .   -snes_qn_scale_type <shanno,none,linesearch,jacobian>
6780c777b0cSPeter Brune 
6790c777b0cSPeter Brune     Level: intermediate
6800c777b0cSPeter Brune 
68184c577daSBarry Smith     SNESQNScaleTypes:
6820c777b0cSPeter Brune +   SNES_QN_SCALE_NONE - don't scale the problem
6830c777b0cSPeter Brune .   SNES_QN_SCALE_SHANNO - use shanno scaling
6840c777b0cSPeter Brune .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
685a01a0525SBarry Smith -   SNES_QN_SCALE_JACOBIAN - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration
686a01a0525SBarry Smith                              of QN and at ever restart.
6870c777b0cSPeter Brune 
688a01a0525SBarry Smith .keywords: scaling type
689a01a0525SBarry Smith 
690a01a0525SBarry Smith .seealso: SNES, SNESQN, SNESLineSearch, SNESQNScaleType, SNESetJacobian()
6910c777b0cSPeter Brune @*/
6920c777b0cSPeter Brune 
6932150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
6942150357eSBarry Smith {
6950c777b0cSPeter Brune   PetscErrorCode ierr;
6966e111a19SKarl Rupp 
6970c777b0cSPeter Brune   PetscFunctionBegin;
6980c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6990c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
7000c777b0cSPeter Brune   PetscFunctionReturn(0);
7010c777b0cSPeter Brune }
7020c777b0cSPeter Brune 
7030adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
7040adebc6cSBarry Smith {
7050c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7066e111a19SKarl Rupp 
7070c777b0cSPeter Brune   PetscFunctionBegin;
7080c777b0cSPeter Brune   qn->scale_type = stype;
7090c777b0cSPeter Brune   PetscFunctionReturn(0);
7100c777b0cSPeter Brune }
7110c777b0cSPeter Brune 
7120adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
7130adebc6cSBarry Smith {
7140c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7156e111a19SKarl Rupp 
7160c777b0cSPeter Brune   PetscFunctionBegin;
7170c777b0cSPeter Brune   qn->restart_type = rtype;
7180c777b0cSPeter Brune   PetscFunctionReturn(0);
7190c777b0cSPeter Brune }
7200c777b0cSPeter Brune 
7211efc8c45SPeter Brune /*@
7221efc8c45SPeter Brune     SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN.
7231efc8c45SPeter Brune 
7241efc8c45SPeter Brune     Logically Collective on SNES
7251efc8c45SPeter Brune 
7261efc8c45SPeter Brune     Input Parameters:
7271efc8c45SPeter Brune +   snes - the iterative context
7281efc8c45SPeter Brune -   qtype - variant type
7291efc8c45SPeter Brune 
7301efc8c45SPeter Brune     Options Database:
73184c577daSBarry Smith .   -snes_qn_type <lbfgs,broyden,badbroyden>
7321efc8c45SPeter Brune 
7331efc8c45SPeter Brune     Level: beginner
7341efc8c45SPeter Brune 
7351efc8c45SPeter Brune     SNESQNTypes:
7361efc8c45SPeter Brune +   SNES_QN_LBFGS - LBFGS variant
7371efc8c45SPeter Brune .   SNES_QN_BROYDEN - Broyden variant
7381efc8c45SPeter Brune -   SNES_QN_BADBROYDEN - Bad Broyden variant
7391efc8c45SPeter Brune 
74084c577daSBarry Smith .keywords: SNES, SNESQN, type, set, SNESQNType
7411efc8c45SPeter Brune @*/
7421efc8c45SPeter Brune 
7431efc8c45SPeter Brune PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
7441efc8c45SPeter Brune {
7451efc8c45SPeter Brune   PetscErrorCode ierr;
7461efc8c45SPeter Brune 
7471efc8c45SPeter Brune   PetscFunctionBegin;
7481efc8c45SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7491efc8c45SPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));CHKERRQ(ierr);
7501efc8c45SPeter Brune   PetscFunctionReturn(0);
7511efc8c45SPeter Brune }
7521efc8c45SPeter Brune 
7531efc8c45SPeter Brune PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
7541efc8c45SPeter Brune {
7551efc8c45SPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7561efc8c45SPeter Brune 
7571efc8c45SPeter Brune   PetscFunctionBegin;
7581efc8c45SPeter Brune   qn->type = qtype;
7591efc8c45SPeter Brune   PetscFunctionReturn(0);
7601efc8c45SPeter Brune }
7611efc8c45SPeter Brune 
7624b11644fSPeter Brune /* -------------------------------------------------------------------------- */
7634b11644fSPeter Brune /*MC
7644b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
7654b11644fSPeter Brune 
7666cc8130cSPeter Brune       Options Database:
7676cc8130cSPeter Brune 
76884c577daSBarry Smith +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
76984c577daSBarry Smith +     -snes_qn_restart_type <powell,periodic,none> - set the restart type
7706bdcc836SBarry Smith .     -snes_qn_powell_gamma - Angle condition for restart.
7711867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
77284c577daSBarry Smith .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
77384c577daSBarry Smith .     -snes_qn_scale_type <shanno,none,linesearch,jacobian> - scaling performed on inner Jacobian
77472835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
77584c577daSBarry Smith -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.
7766cc8130cSPeter Brune 
77795452b02SPatrick Sanan       Notes:
77895452b02SPatrick Sanan     This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
779b8840d6bSPeter Brune       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
780b8840d6bSPeter Brune       updates.
7816cc8130cSPeter Brune 
7821867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
7831867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
78484c577daSBarry Smith       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
78584c577daSBarry Smith       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
7861867fe5bSPeter Brune 
7872d547940SBarry Smith       Uses left nonlinear preconditioning by default.
7882d547940SBarry Smith 
7896cc8130cSPeter Brune       References:
79096a0c994SBarry Smith +   1. -   Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
79196a0c994SBarry Smith .   2. -   R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods,
7920a8e8854SPeter Brune       Technical Report, Northwestern University, June 1992.
79396a0c994SBarry Smith .   3. -   Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
7940a8e8854SPeter Brune       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
79596a0c994SBarry Smith -   4. -   Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
7964f02bc6aSBarry Smith        SIAM Review, 57(4), 2015
7974b11644fSPeter Brune 
7984b11644fSPeter Brune       Level: beginner
7994b11644fSPeter Brune 
80004d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
8016cc8130cSPeter Brune 
8024b11644fSPeter Brune M*/
8038cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
8044b11644fSPeter Brune {
8054b11644fSPeter Brune   PetscErrorCode ierr;
8069f83bee8SJed Brown   SNES_QN        *qn;
8074b11644fSPeter Brune 
8084b11644fSPeter Brune   PetscFunctionBegin;
8094b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
8104b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
8114b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
8124b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
8135cd83419SPeter Brune   snes->ops->view           = SNESView_QN;
8144b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
8154b11644fSPeter Brune 
816efd4aadfSBarry Smith   snes->npcside= PC_LEFT;
81747f26062SPeter Brune 
818efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
81942f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
82042f4f86dSBarry Smith 
8214fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
8224fc747eaSLawrence Mitchell 
82388976e71SPeter Brune   if (!snes->tolerancesset) {
8240e444f03SPeter Brune     snes->max_funcs = 30000;
8250e444f03SPeter Brune     snes->max_its   = 10000;
82688976e71SPeter Brune   }
8270e444f03SPeter Brune 
828b00a9115SJed Brown   ierr                = PetscNewLog(snes,&qn);CHKERRQ(ierr);
8294b11644fSPeter Brune   snes->data          = (void*) qn;
8300ecc9e1dSPeter Brune   qn->m               = 10;
831b21d5a53SPeter Brune   qn->scaling         = 1.0;
8320298fd71SBarry Smith   qn->U               = NULL;
8330298fd71SBarry Smith   qn->V               = NULL;
8340298fd71SBarry Smith   qn->dXtdF           = NULL;
8350298fd71SBarry Smith   qn->dFtdX           = NULL;
8360298fd71SBarry Smith   qn->dXdFmat         = NULL;
8370298fd71SBarry Smith   qn->monitor         = NULL;
8381c6523dcSPeter Brune   qn->singlereduction = PETSC_TRUE;
839b8840d6bSPeter Brune   qn->powell_gamma    = 0.9999;
84060c08b40SPeter Brune   qn->scale_type      = SNES_QN_SCALE_DEFAULT;
84160c08b40SPeter Brune   qn->restart_type    = SNES_QN_RESTART_DEFAULT;
842b8840d6bSPeter Brune   qn->type            = SNES_QN_LBFGS;
843ea630c6eSPeter Brune 
844bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr);
845bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr);
8461efc8c45SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);CHKERRQ(ierr);
8474b11644fSPeter Brune   PetscFunctionReturn(0);
8484b11644fSPeter Brune }
849