xref: /petsc/src/snes/impls/qn/qn.c (revision fffbeea80931a16407cd3b4c477c2108c2b8c307)
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 
51efc8c45SPeter Brune const char *const SNESQNScaleTypes[] =        {"DEFAULT","NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0};
61efc8c45SPeter Brune const char *const SNESQNRestartTypes[] =      {"DEFAULT","NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0};
760c08b40SPeter Brune const char *const SNESQNTypes[] =             {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_",0};
8b8840d6bSPeter Brune 
94b11644fSPeter Brune typedef struct {
10b8840d6bSPeter Brune   Vec               *U;                   /* Stored past states (vary from method to method) */
11b8840d6bSPeter Brune   Vec               *V;                   /* Stored past states (vary from method to method) */
128e231d97SPeter Brune   PetscInt          m;                    /* The number of kept previous steps */
135c7a0a03SPeter Brune   PetscReal         *lambda;              /* The line search history of the method */
145c7a0a03SPeter Brune   PetscReal         *norm;                /* norms of the steps */
158e231d97SPeter Brune   PetscScalar       *alpha, *beta;
168e231d97SPeter Brune   PetscScalar       *dXtdF, *dFtdX, *YtdX;
1718aa0c0cSPeter Brune   PetscBool         singlereduction;      /* Aggregated reduction implementation */
188e231d97SPeter Brune   PetscScalar       *dXdFmat;             /* A matrix of values for dX_i dot dF_j */
1944f7e39eSPeter Brune   PetscViewer       monitor;
206bf1b2e5SPeter Brune   PetscReal         powell_gamma;         /* Powell angle restart condition */
216bf1b2e5SPeter Brune   PetscReal         powell_downhill;      /* Powell descent restart condition */
22b21d5a53SPeter Brune   PetscReal         scaling;              /* scaling of H0 */
23b8840d6bSPeter Brune   SNESQNType        type;                 /* the type of quasi-newton method used */
2488f769c5SPeter Brune   SNESQNScaleType   scale_type;           /* the type of scaling used */
250c777b0cSPeter Brune   SNESQNRestartType restart_type;         /* determine the frequency and type of restart conditions */
269f83bee8SJed Brown } SNES_QN;
274b11644fSPeter Brune 
284b11644fSPeter Brune #undef __FUNCT__
29b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_Broyden"
305051edb1SPeter Brune PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D)
312150357eSBarry Smith {
324b11644fSPeter Brune   PetscErrorCode     ierr;
339f83bee8SJed Brown   SNES_QN            *qn = (SNES_QN*)snes->data;
34b8840d6bSPeter Brune   Vec                W   = snes->work[3];
35b8840d6bSPeter Brune   Vec                *U  = qn->U;
36b8840d6bSPeter Brune   KSPConvergedReason kspreason;
37b8840d6bSPeter Brune   PetscInt           m = qn->m;
385c7a0a03SPeter Brune   PetscInt           k,i,j,lits,l = m;
395c7a0a03SPeter Brune   PetscReal          unorm,a,b;
405c7a0a03SPeter Brune   PetscReal          *lambda=qn->lambda;
41b8840d6bSPeter Brune   PetscScalar        gdot;
4259e7931aSPeter Brune   PetscReal          udot;
432150357eSBarry Smith 
44b8840d6bSPeter Brune   PetscFunctionBegin;
45b8840d6bSPeter Brune   if (it < m) l = it;
465c7a0a03SPeter Brune   if (it > 0) {
475c7a0a03SPeter Brune     k = (it-1)%l;
485c7a0a03SPeter Brune     ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr);
49a49fa0d8SPeter Brune     ierr = VecCopy(Xold,U[k]);CHKERRQ(ierr);
505051edb1SPeter Brune     ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr);
515c7a0a03SPeter Brune     if (qn->monitor) {
525c7a0a03SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
535c7a0a03SPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %d of %d by lambda: %14.12e \n",k,l,lambda[k]);CHKERRQ(ierr);
545c7a0a03SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
555c7a0a03SPeter Brune     }
565c7a0a03SPeter Brune   }
57b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
58d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr);
59b8840d6bSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
60b8840d6bSPeter Brune     if (kspreason < 0) {
61b8840d6bSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
62b8840d6bSPeter 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);
63b8840d6bSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
64b8840d6bSPeter Brune         PetscFunctionReturn(0);
65b8840d6bSPeter Brune       }
66b8840d6bSPeter Brune     }
67b8840d6bSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
68b8840d6bSPeter Brune     snes->linear_its += lits;
69b8840d6bSPeter Brune     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
70b8840d6bSPeter Brune   } else {
71b8840d6bSPeter Brune     ierr = VecCopy(D,Y);CHKERRQ(ierr);
72b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
73b8840d6bSPeter Brune   }
74b8840d6bSPeter Brune 
75b8840d6bSPeter Brune   /* inward recursion starting at the first update and working forward */
76b8840d6bSPeter Brune   for (i = 0; i < l-1; i++) {
775c7a0a03SPeter Brune     j = (it+i-l)%l;
785c7a0a03SPeter Brune     k = (it+i-l+1)%l;
795c7a0a03SPeter Brune     ierr = VecNorm(U[j],NORM_2,&unorm);CHKERRQ(ierr);
805c7a0a03SPeter Brune     ierr = VecDot(U[j],Y,&gdot);CHKERRQ(ierr);
8159e7931aSPeter Brune     unorm *= unorm;
8259e7931aSPeter Brune     udot = PetscRealPart(gdot);
835c7a0a03SPeter Brune     a = (lambda[j]/lambda[k]);
845c7a0a03SPeter Brune     b = -(1.-lambda[j]);
8559e7931aSPeter Brune     a *= udot/unorm;
8659e7931aSPeter Brune     b *= udot/unorm;
875c7a0a03SPeter Brune     ierr = VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);CHKERRQ(ierr);
885c7a0a03SPeter Brune 
89b8840d6bSPeter Brune     if (qn->monitor) {
90b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
91e639b251SJed Brown       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d and %d, gdot: %14.12e\n",k,j,PetscRealPart(gdot));CHKERRQ(ierr);
92b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
93b8840d6bSPeter Brune     }
94b8840d6bSPeter Brune   }
95b8840d6bSPeter Brune   if (l > 0) {
96b8840d6bSPeter Brune     k = (it-1)%l;
97b8840d6bSPeter Brune     ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
985c7a0a03SPeter Brune     ierr = VecNorm(U[k],NORM_2,&unorm);CHKERRQ(ierr);
9959e7931aSPeter Brune     unorm *= unorm;
10059e7931aSPeter Brune     udot = PetscRealPart(gdot);
10159e7931aSPeter Brune     a = unorm/(unorm-lambda[k]*udot);
10259e7931aSPeter Brune     b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot);
103b8840d6bSPeter Brune     if (qn->monitor) {
104b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
10559e7931aSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d: a: %14.12e b: %14.12e \n",k,a,b);CHKERRQ(ierr);
106b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
107b8840d6bSPeter Brune     }
1085c7a0a03SPeter Brune     ierr = VecAXPBY(Y,b,a,U[k]);CHKERRQ(ierr);
109b8840d6bSPeter Brune   }
1105c7a0a03SPeter Brune   l = m;
1115c7a0a03SPeter Brune   if (it+1<m)l=it+1;
1125c7a0a03SPeter Brune   k = it%l;
1135c7a0a03SPeter Brune   if (qn->monitor) {
1145c7a0a03SPeter Brune     ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1155c7a0a03SPeter Brune     ierr = PetscViewerASCIIPrintf(qn->monitor, "setting vector %d of %d\n",k,l);CHKERRQ(ierr);
1165c7a0a03SPeter Brune     ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1175c7a0a03SPeter Brune   }
118b8840d6bSPeter Brune   PetscFunctionReturn(0);
119b8840d6bSPeter Brune }
120b8840d6bSPeter Brune 
121b8840d6bSPeter Brune #undef __FUNCT__
122b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden"
1230adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1240adebc6cSBarry Smith {
125b8840d6bSPeter Brune   PetscErrorCode ierr;
126b8840d6bSPeter Brune   SNES_QN        *qn = (SNES_QN*)snes->data;
127b8840d6bSPeter Brune   Vec            W   = snes->work[3];
128b8840d6bSPeter Brune   Vec            *U  = qn->U;
129b8840d6bSPeter Brune   Vec            *T  = qn->V;
130b8840d6bSPeter Brune 
131b8840d6bSPeter Brune   /* ksp thing for jacobian scaling */
132b8840d6bSPeter Brune   KSPConvergedReason kspreason;
133c9e59058SPeter Brune   PetscInt           h,k,j,i,lits;
134b8840d6bSPeter Brune   PetscInt           m = qn->m;
1355051edb1SPeter Brune   PetscScalar        gdot,udot;
136b8840d6bSPeter Brune   PetscInt           l = m;
1370adebc6cSBarry Smith 
138b8840d6bSPeter Brune   PetscFunctionBegin;
139b8840d6bSPeter Brune   if (it < m) l = it;
140b8840d6bSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
141b8840d6bSPeter Brune   if (l > 0) {
142b8840d6bSPeter Brune     k    = (it-1)%l;
143c9e59058SPeter Brune     ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);CHKERRQ(ierr);
144b8840d6bSPeter Brune     ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr);
145b8840d6bSPeter Brune     ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr);
1465051edb1SPeter Brune     ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr);
1475051edb1SPeter Brune     ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr);
148b8840d6bSPeter Brune   }
149b8840d6bSPeter Brune 
150b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
151d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
152b8840d6bSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
153b8840d6bSPeter Brune     if (kspreason < 0) {
154b8840d6bSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
155b8840d6bSPeter 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);
156b8840d6bSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
157b8840d6bSPeter Brune         PetscFunctionReturn(0);
158b8840d6bSPeter Brune       }
159b8840d6bSPeter Brune     }
160b8840d6bSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
161b8840d6bSPeter Brune     snes->linear_its += lits;
162b8840d6bSPeter Brune     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
163b8840d6bSPeter Brune   } else {
164b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
165b8840d6bSPeter Brune   }
1665051edb1SPeter Brune 
1675051edb1SPeter Brune   /* inward recursion starting at the first update and working forward */
1685051edb1SPeter Brune   if (l > 0) {
1695051edb1SPeter Brune     for (i = 0; i < l-1; i++) {
170c9e59058SPeter Brune       j    = (it+i-l)%l;
171c9e59058SPeter Brune       k    = (it+i-l+1)%l;
172c9e59058SPeter Brune       h    = (it-1)%l;
173c9e59058SPeter Brune       ierr = VecDotBegin(U[j],U[h],&gdot);CHKERRQ(ierr);
174c9e59058SPeter Brune       ierr = VecDotBegin(U[j],U[j],&udot);CHKERRQ(ierr);
175c9e59058SPeter Brune       ierr = VecDotEnd(U[j],U[h],&gdot);CHKERRQ(ierr);
176c9e59058SPeter Brune       ierr = VecDotEnd(U[j],U[j],&udot);CHKERRQ(ierr);
1775051edb1SPeter Brune       ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr);
178c9e59058SPeter Brune       ierr = VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);CHKERRQ(ierr);
1795051edb1SPeter Brune       if (qn->monitor) {
1805051edb1SPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1815051edb1SPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
1825051edb1SPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1835051edb1SPeter Brune       }
1845051edb1SPeter Brune     }
1855051edb1SPeter Brune   }
186b8840d6bSPeter Brune   PetscFunctionReturn(0);
187b8840d6bSPeter Brune }
188b8840d6bSPeter Brune 
189b8840d6bSPeter Brune #undef __FUNCT__
190b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS"
1910adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1920adebc6cSBarry Smith {
193b8840d6bSPeter Brune   PetscErrorCode ierr;
194b8840d6bSPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
195b8840d6bSPeter Brune   Vec            W      = snes->work[3];
196b8840d6bSPeter Brune   Vec            *dX    = qn->U;
197b8840d6bSPeter Brune   Vec            *dF    = qn->V;
198bd052dfeSPeter Brune   PetscScalar    *alpha = qn->alpha;
199bd052dfeSPeter Brune   PetscScalar    *beta  = qn->beta;
2008e231d97SPeter Brune   PetscScalar    *dXtdF = qn->dXtdF;
201b8840d6bSPeter Brune   PetscScalar    *dFtdX = qn->dFtdX;
2028e231d97SPeter Brune   PetscScalar    *YtdX  = qn->YtdX;
203bd052dfeSPeter Brune 
2040ecc9e1dSPeter Brune   /* ksp thing for jacobian scaling */
2050ecc9e1dSPeter Brune   KSPConvergedReason kspreason;
2068e231d97SPeter Brune   PetscInt           k,i,j,g,lits;
2074b11644fSPeter Brune   PetscInt           m = qn->m;
2084b11644fSPeter Brune   PetscScalar        t;
2094b11644fSPeter Brune   PetscInt           l = m;
2108e231d97SPeter Brune 
2114b11644fSPeter Brune   PetscFunctionBegin;
2125ba6227bSPeter Brune   if (it < m) l = it;
2131c6523dcSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
214b8840d6bSPeter Brune   if (it > 0) {
215b8840d6bSPeter Brune     k    = (it - 1) % l;
216b8840d6bSPeter Brune     ierr = VecCopy(D,dF[k]);CHKERRQ(ierr);
217b8840d6bSPeter Brune     ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
218b8840d6bSPeter Brune     ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
219b8840d6bSPeter Brune     ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
22018aa0c0cSPeter Brune     if (qn->singlereduction) {
2211c6523dcSPeter Brune       ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
2221c6523dcSPeter Brune       ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
2231c6523dcSPeter Brune       ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr);
2241c6523dcSPeter Brune       ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
2251c6523dcSPeter Brune       ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
2261c6523dcSPeter Brune       ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr);
227b8840d6bSPeter Brune       for (j = 0; j < l; j++) {
228b8840d6bSPeter Brune         H(k, j) = dFtdX[j];
229b8840d6bSPeter Brune         H(j, k) = dXtdF[j];
230b8840d6bSPeter Brune       }
231b8840d6bSPeter Brune       /* copy back over to make the computation of alpha and beta easier */
2321aa26658SKarl Rupp       for (j = 0; j < l; j++) dXtdF[j] = H(j, j);
233b8840d6bSPeter Brune     } else {
234b8840d6bSPeter Brune       ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr);
23523d44fbcSPeter Brune     }
23623d44fbcSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
237b8840d6bSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr);
238b8840d6bSPeter Brune     }
239b8840d6bSPeter Brune   }
240b8840d6bSPeter Brune 
2414b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
2424b11644fSPeter Brune   for (i=0; i<l; i++) {
243b21d5a53SPeter Brune     k = (it-i-1)%l;
24418aa0c0cSPeter Brune     if (qn->singlereduction) {
2458e231d97SPeter Brune       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
2468e231d97SPeter Brune       t = YtdX[k];
2478e231d97SPeter Brune       for (j=0; j<i; j++) {
2488e231d97SPeter Brune         g  = (it-j-1)%l;
249b3cfccb2SPeter Brune         t -= alpha[g]*H(k, g);
2508e231d97SPeter Brune       }
2518e231d97SPeter Brune       alpha[k] = t/H(k,k);
2528e231d97SPeter Brune     } else {
2533af51624SPeter Brune       ierr     = VecDot(dX[k],Y,&t);CHKERRQ(ierr);
2548e231d97SPeter Brune       alpha[k] = t/dXtdF[k];
2558e231d97SPeter Brune     }
25644f7e39eSPeter Brune     if (qn->monitor) {
2573af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2585ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
2593af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
26044f7e39eSPeter Brune     }
2616bf1b2e5SPeter Brune     ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr);
2624b11644fSPeter Brune   }
2634b11644fSPeter Brune 
2640c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
265d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
2660ecc9e1dSPeter Brune     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
2670ecc9e1dSPeter Brune     if (kspreason < 0) {
2680ecc9e1dSPeter Brune       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
2690ecc9e1dSPeter 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);
2700ecc9e1dSPeter Brune         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
2710ecc9e1dSPeter Brune         PetscFunctionReturn(0);
2720ecc9e1dSPeter Brune       }
2730ecc9e1dSPeter Brune     }
2740ecc9e1dSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2750ecc9e1dSPeter Brune     snes->linear_its += lits;
276b8840d6bSPeter Brune     ierr              = VecCopy(W, Y);CHKERRQ(ierr);
2770ecc9e1dSPeter Brune   } else {
278b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
2790ecc9e1dSPeter Brune   }
28018aa0c0cSPeter Brune   if (qn->singlereduction) {
281b8840d6bSPeter Brune     ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr);
2828e231d97SPeter Brune   }
2834b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
284bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
285b21d5a53SPeter Brune     k = (it + i - l) % l;
28618aa0c0cSPeter Brune     if (qn->singlereduction) {
2878e231d97SPeter Brune       t = YtdX[k];
2888e231d97SPeter Brune       for (j = 0; j < i; j++) {
2898e231d97SPeter Brune         g  = (it + j - l) % l;
290b3cfccb2SPeter Brune         t += (alpha[g] - beta[g])*H(g, k);
2918e231d97SPeter Brune       }
2928e231d97SPeter Brune       beta[k] = t / H(k, k);
2938e231d97SPeter Brune     } else {
2946bf1b2e5SPeter Brune       ierr    = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
2958e231d97SPeter Brune       beta[k] = t / dXtdF[k];
2968e231d97SPeter Brune     }
29722d28d08SBarry Smith     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr);
29844f7e39eSPeter Brune     if (qn->monitor) {
2993af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
3005ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
3013af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
30244f7e39eSPeter Brune     }
3034b11644fSPeter Brune   }
3044b11644fSPeter Brune   PetscFunctionReturn(0);
3054b11644fSPeter Brune }
3064b11644fSPeter Brune 
3074b11644fSPeter Brune #undef __FUNCT__
3084b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
3094b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
3104b11644fSPeter Brune {
3114b11644fSPeter Brune   PetscErrorCode      ierr;
3129f83bee8SJed Brown   SNES_QN             *qn = (SNES_QN*) snes->data;
31315f5eeeaSPeter Brune   Vec                 X,Xold;
3140a3905e1SPeter Brune   Vec                 F,W;
31547f26062SPeter Brune   Vec                 Y,D,Dold;
316b8840d6bSPeter Brune   PetscInt            i, i_r;
317335fb975SPeter Brune   PetscReal           fnorm,xnorm,ynorm,gnorm;
3180c777b0cSPeter Brune   PetscBool           lssucceed,powell,periodic;
3191c6523dcSPeter Brune   PetscScalar         DolddotD,DolddotDold;
3200ecc9e1dSPeter Brune   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
321b7281c8aSPeter Brune   SNESConvergedReason reason;
3220ecc9e1dSPeter Brune 
3234b11644fSPeter Brune   /* basically just a regular newton's method except for the application of the jacobian */
3244b11644fSPeter Brune 
3256e111a19SKarl Rupp   PetscFunctionBegin;
326*fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
3279f3a0142SPeter Brune   F    = snes->vec_func;                /* residual vector */
3283af51624SPeter Brune   Y    = snes->vec_sol_update;          /* search direction generated by J^-1D*/
3290a3905e1SPeter Brune   W    = snes->work[3];
330b8840d6bSPeter Brune   X    = snes->vec_sol;                 /* solution vector */
331335fb975SPeter Brune   Xold = snes->work[0];
3323af51624SPeter Brune 
3333af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
334335fb975SPeter Brune   D    = snes->work[1];
335335fb975SPeter Brune   Dold = snes->work[2];
3364b11644fSPeter Brune 
3374b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
3384b11644fSPeter Brune 
339e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3404b11644fSPeter Brune   snes->iter = 0;
3414b11644fSPeter Brune   snes->norm = 0.;
342e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
34347f26062SPeter Brune 
344b28a06ddSPeter Brune   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
34547f26062SPeter Brune     ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr);
346b7281c8aSPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
347b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
348b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
349b7281c8aSPeter Brune       PetscFunctionReturn(0);
350b7281c8aSPeter Brune     }
35147f26062SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
35247f26062SPeter Brune   } else {
353e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
35415f5eeeaSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3554b11644fSPeter Brune       if (snes->domainerror) {
3564b11644fSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
3574b11644fSPeter Brune         PetscFunctionReturn(0);
3584b11644fSPeter Brune       }
3591aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
360c1c75074SPeter Brune 
36147f26062SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
362189a9710SBarry Smith     if (PetscIsInfOrNanReal(fnorm)) {
363189a9710SBarry Smith       snes->reason = SNES_DIVERGED_FNORM_NAN;
364189a9710SBarry Smith       PetscFunctionReturn(0);
365189a9710SBarry Smith     }
36647f26062SPeter Brune   }
367b28a06ddSPeter Brune   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
36847f26062SPeter Brune       ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr);
369b7281c8aSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
370b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
371b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
372b7281c8aSPeter Brune         PetscFunctionReturn(0);
373b7281c8aSPeter Brune       }
37447f26062SPeter Brune   } else {
37547f26062SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
37647f26062SPeter Brune   }
377b28a06ddSPeter Brune 
378e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3794b11644fSPeter Brune   snes->norm = fnorm;
380e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
381a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
3824b11644fSPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
3834b11644fSPeter Brune 
3844b11644fSPeter Brune   /* test convergence */
3854b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3864b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
38770d3b23bSPeter Brune 
3883cf07b75SPeter Brune   if (snes->pc && snes->pcside == PC_RIGHT) {
3893cf07b75SPeter Brune     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
3903cf07b75SPeter Brune     ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
3913cf07b75SPeter Brune     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
392ddd40ce5SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
393ddd40ce5SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
394ddd40ce5SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
395ddd40ce5SPeter Brune       PetscFunctionReturn(0);
396ddd40ce5SPeter Brune     }
397ddd40ce5SPeter Brune     ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
3983cf07b75SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
3993cf07b75SPeter Brune   }
4003cf07b75SPeter Brune 
401f8e15203SPeter Brune   /* scale the initial update */
4020c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
4030ecc9e1dSPeter Brune     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
4048cba9625SJed Brown     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
4050ecc9e1dSPeter Brune   }
4060ecc9e1dSPeter Brune 
4075ba6227bSPeter Brune   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
4080a3905e1SPeter Brune     if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) {
4090a3905e1SPeter Brune       PetscScalar ff,xf;
4100a3905e1SPeter Brune       ierr = VecCopy(Dold,Y);CHKERRQ(ierr);
4110a3905e1SPeter Brune       ierr = VecCopy(Xold,W);CHKERRQ(ierr);
4120a3905e1SPeter Brune       ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr);
4130a3905e1SPeter Brune       ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr);
4140a3905e1SPeter Brune       ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr);
4150a3905e1SPeter Brune       ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr);
4160a3905e1SPeter Brune       ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr);
4170a3905e1SPeter Brune       ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr);
4180a3905e1SPeter Brune       qn->scaling = PetscRealPart(xf)/PetscRealPart(ff);
4190a3905e1SPeter Brune     }
420b8840d6bSPeter Brune     switch (qn->type) {
421b8840d6bSPeter Brune     case SNES_QN_BADBROYDEN:
422b8840d6bSPeter Brune       ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
423b8840d6bSPeter Brune       break;
424b8840d6bSPeter Brune     case SNES_QN_BROYDEN:
4255051edb1SPeter Brune       ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr);
426b8840d6bSPeter Brune       break;
427b8840d6bSPeter Brune     case SNES_QN_LBFGS:
428b8840d6bSPeter Brune       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
429b8840d6bSPeter Brune       break;
430b8840d6bSPeter Brune     }
43170d3b23bSPeter Brune     /* line search for lambda */
43270d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
4330a3905e1SPeter Brune     ierr  = VecCopy(D, Dold);CHKERRQ(ierr);
43415f5eeeaSPeter Brune     ierr  = VecCopy(X, Xold);CHKERRQ(ierr);
435f1c6b773SPeter Brune     ierr  = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
4369f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
4379f3a0142SPeter Brune     if (snes->domainerror) {
4389f3a0142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
4399f3a0142SPeter Brune       PetscFunctionReturn(0);
4409f3a0142SPeter Brune     }
441f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr);
4429f3a0142SPeter Brune     if (!lssucceed) {
4439f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
4449f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
4459f3a0142SPeter Brune         break;
4469f3a0142SPeter Brune       }
4479f3a0142SPeter Brune     }
448f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
4490c777b0cSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
45005ee83feSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
4510ecc9e1dSPeter Brune     }
4523af51624SPeter Brune 
45388d374b2SPeter Brune     /* convergence monitoring */
454b21d5a53SPeter 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);
455b21d5a53SPeter Brune 
456b28a06ddSPeter Brune     if (snes->pc && snes->pcside == PC_RIGHT) {
457b28a06ddSPeter Brune       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
458b28a06ddSPeter Brune       ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
459b28a06ddSPeter Brune       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
460ddd40ce5SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
461ddd40ce5SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
462ddd40ce5SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
463ddd40ce5SPeter Brune         PetscFunctionReturn(0);
464ddd40ce5SPeter Brune       }
465ddd40ce5SPeter Brune       ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
466b28a06ddSPeter Brune     }
467b28a06ddSPeter Brune 
468360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
469360c497dSPeter Brune     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
470360c497dSPeter Brune 
471a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
4728409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
4734b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
474d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
4754b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
476b28a06ddSPeter Brune     if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
47747f26062SPeter Brune       ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr);
478b7281c8aSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
479b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
480b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
481b7281c8aSPeter Brune         PetscFunctionReturn(0);
482b7281c8aSPeter Brune       }
48388d374b2SPeter Brune     } else {
48488d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
48588d374b2SPeter Brune     }
4860c777b0cSPeter Brune     powell = PETSC_FALSE;
4870c777b0cSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_POWELL) {
4886bf1b2e5SPeter Brune       /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
48923c918e7SPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
49023c918e7SPeter Brune         ierr = MatMult(snes->jacobian_pre,Dold,W);CHKERRQ(ierr);
49123c918e7SPeter Brune       } else {
49223c918e7SPeter Brune         ierr = VecCopy(Dold,W);CHKERRQ(ierr);
49323c918e7SPeter Brune       }
49423c918e7SPeter Brune       ierr = VecDotBegin(W, Dold, &DolddotDold);CHKERRQ(ierr);
49523c918e7SPeter Brune       ierr = VecDotBegin(W, D, &DolddotD);CHKERRQ(ierr);
49623c918e7SPeter Brune       ierr = VecDotEnd(W, Dold, &DolddotDold);CHKERRQ(ierr);
49723c918e7SPeter Brune       ierr = VecDotEnd(W, D, &DolddotD);CHKERRQ(ierr);
4980c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
4990c777b0cSPeter Brune     }
5000c777b0cSPeter Brune     periodic = PETSC_FALSE;
501b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
502b8840d6bSPeter Brune       if (i_r>qn->m-1) periodic = PETSC_TRUE;
5030c777b0cSPeter Brune     }
5040c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
5050c777b0cSPeter Brune     if (powell || periodic) {
5065ba6227bSPeter Brune       if (qn->monitor) {
5075ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
508516fe3c3SPeter 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);
5095ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
5105ba6227bSPeter Brune       }
5115ba6227bSPeter Brune       i_r = -1;
5125ba6227bSPeter Brune       /* general purpose update */
5135ba6227bSPeter Brune       if (snes->ops->update) {
5145ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
5155ba6227bSPeter Brune       }
5160c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
5170ecc9e1dSPeter Brune         ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
5188cba9625SJed Brown         ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
5190ecc9e1dSPeter Brune       }
5200ecc9e1dSPeter Brune     }
52170d3b23bSPeter Brune     /* general purpose update */
52270d3b23bSPeter Brune     if (snes->ops->update) {
52370d3b23bSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
52470d3b23bSPeter Brune     }
5255ba6227bSPeter Brune   }
5264b11644fSPeter Brune   if (i == snes->max_its) {
5274b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
5284b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
5294b11644fSPeter Brune   }
5304b11644fSPeter Brune   PetscFunctionReturn(0);
5314b11644fSPeter Brune }
5324b11644fSPeter Brune 
5334b11644fSPeter Brune #undef __FUNCT__
5344b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
5354b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
5364b11644fSPeter Brune {
5379f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
5384b11644fSPeter Brune   PetscErrorCode ierr;
539335fb975SPeter Brune 
5404b11644fSPeter Brune   PetscFunctionBegin;
541b8840d6bSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr);
54259e7931aSPeter Brune   if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr);
543dcca6d9dSJed Brown   ierr = PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);CHKERRQ(ierr);
5448e231d97SPeter Brune 
54518aa0c0cSPeter Brune   if (qn->singlereduction) {
546dcca6d9dSJed Brown     ierr = PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);CHKERRQ(ierr);
5478e231d97SPeter Brune   }
548fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
54960c08b40SPeter Brune   /* set method defaults */
5501efc8c45SPeter Brune   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
55160c08b40SPeter Brune     if (qn->type == SNES_QN_BADBROYDEN) {
55260c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_NONE;
55360c08b40SPeter Brune     } else {
55460c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_SHANNO;
55560c08b40SPeter Brune     }
55660c08b40SPeter Brune   }
5571efc8c45SPeter Brune   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
55860c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
55960c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_POWELL;
56060c08b40SPeter Brune     } else {
56160c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_PERIODIC;
56260c08b40SPeter Brune     }
56360c08b40SPeter Brune   }
56460c08b40SPeter Brune 
5650c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
5668e231d97SPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
5678e231d97SPeter Brune   }
5686c67d002SPeter Brune   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;}
5694b11644fSPeter Brune   PetscFunctionReturn(0);
5704b11644fSPeter Brune }
5714b11644fSPeter Brune 
5724b11644fSPeter Brune #undef __FUNCT__
5734b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
5744b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
5754b11644fSPeter Brune {
5764b11644fSPeter Brune   PetscErrorCode ierr;
5779f83bee8SJed Brown   SNES_QN        *qn;
5780adebc6cSBarry Smith 
5794b11644fSPeter Brune   PetscFunctionBegin;
5804b11644fSPeter Brune   if (snes->data) {
5819f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
582b8840d6bSPeter Brune     if (qn->U) {
583b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr);
5844b11644fSPeter Brune     }
585b8840d6bSPeter Brune     if (qn->V) {
586b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr);
5874b11644fSPeter Brune     }
58818aa0c0cSPeter Brune     if (qn->singlereduction) {
5898e231d97SPeter Brune       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
5908e231d97SPeter Brune     }
5915c7a0a03SPeter Brune     ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr);
5924b11644fSPeter Brune   }
5934b11644fSPeter Brune   PetscFunctionReturn(0);
5944b11644fSPeter Brune }
5954b11644fSPeter Brune 
5964b11644fSPeter Brune #undef __FUNCT__
5974b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
5984b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
5994b11644fSPeter Brune {
6004b11644fSPeter Brune   PetscErrorCode ierr;
6016e111a19SKarl Rupp 
6024b11644fSPeter Brune   PetscFunctionBegin;
6034b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
6044b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
605bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
6064b11644fSPeter Brune   PetscFunctionReturn(0);
6074b11644fSPeter Brune }
6084b11644fSPeter Brune 
6094b11644fSPeter Brune #undef __FUNCT__
6104b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
6114b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
6124b11644fSPeter Brune {
6134b11644fSPeter Brune 
6144b11644fSPeter Brune   PetscErrorCode    ierr;
6152150357eSBarry Smith   SNES_QN           *qn    = (SNES_QN*)snes->data;
61688f769c5SPeter Brune   PetscBool         monflg = PETSC_FALSE,flg;
617f1c6b773SPeter Brune   SNESLineSearch    linesearch;
6182150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
6192150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
6201efc8c45SPeter Brune   SNESQNType        qtype = qn->type;
6212150357eSBarry Smith 
6224b11644fSPeter Brune   PetscFunctionBegin;
6234b11644fSPeter Brune   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
6240298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr);
6250298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr);
6260298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr);
6270298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr);
6280298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr);
62988f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr);
63088f769c5SPeter Brune   if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr);
63188f769c5SPeter Brune 
63288f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr);
63388f769c5SPeter Brune   if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr);
63488f769c5SPeter Brune 
635b8840d6bSPeter Brune   ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,
6361efc8c45SPeter Brune                           (PetscEnum)qtype,(PetscEnum*)&qtype,&flg);CHKERRQ(ierr);
6371efc8c45SPeter Brune   if (flg) {ierr = SNESQNSetType(snes,qtype);CHKERRQ(ierr);}
6384b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
6399e764e56SPeter Brune   if (!snes->linesearch) {
6407601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
64160c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
6421a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
64360c08b40SPeter Brune     } else if (qn->type == SNES_QN_BROYDEN) {
64460c08b40SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);CHKERRQ(ierr);
64560c08b40SPeter Brune     } else {
64660c08b40SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
64760c08b40SPeter Brune     }
6489e764e56SPeter Brune   }
64944f7e39eSPeter Brune   if (monflg) {
650ce94432eSBarry Smith     qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
65144f7e39eSPeter Brune   }
6524b11644fSPeter Brune   PetscFunctionReturn(0);
6534b11644fSPeter Brune }
6544b11644fSPeter Brune 
6555cd83419SPeter Brune #undef __FUNCT__
6565cd83419SPeter Brune #define __FUNCT__ "SNESView_QN"
6575cd83419SPeter Brune static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
6585cd83419SPeter Brune {
6595cd83419SPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
6605cd83419SPeter Brune   PetscBool      iascii;
6615cd83419SPeter Brune   PetscErrorCode ierr;
6625cd83419SPeter Brune 
6635cd83419SPeter Brune   PetscFunctionBegin;
6645cd83419SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
6655cd83419SPeter Brune   if (iascii) {
6665cd83419SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  QN type is %s, restart type is %s, scale type is %s\n",SNESQNTypes[qn->type],SNESQNRestartTypes[qn->restart_type],SNESQNScaleTypes[qn->scale_type]);CHKERRQ(ierr);
6675cd83419SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  Stored subspace size: %d\n", qn->m);CHKERRQ(ierr);
6685cd83419SPeter Brune     if (qn->singlereduction) {
6695cd83419SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  Using the single reduction variant.\n");CHKERRQ(ierr);
6705cd83419SPeter Brune     }
6715cd83419SPeter Brune   }
6725cd83419SPeter Brune   PetscFunctionReturn(0);
6735cd83419SPeter Brune }
67492c02d66SPeter Brune 
6750c777b0cSPeter Brune #undef __FUNCT__
6760c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType"
6770c777b0cSPeter Brune /*@
6780c777b0cSPeter Brune     SNESQNSetRestartType - Sets the restart type for SNESQN.
6790c777b0cSPeter Brune 
6800c777b0cSPeter Brune     Logically Collective on SNES
6810c777b0cSPeter Brune 
6820c777b0cSPeter Brune     Input Parameters:
6830c777b0cSPeter Brune +   snes - the iterative context
6840c777b0cSPeter Brune -   rtype - restart type
6850c777b0cSPeter Brune 
6860c777b0cSPeter Brune     Options Database:
6870c777b0cSPeter Brune +   -snes_qn_restart_type<powell,periodic,none> - set the restart type
688b8840d6bSPeter Brune -   -snes_qn_restart[10] - sets the number of iterations before restart for periodic
6890c777b0cSPeter Brune 
6900c777b0cSPeter Brune     Level: intermediate
6910c777b0cSPeter Brune 
6920c777b0cSPeter Brune     SNESQNRestartTypes:
6930c777b0cSPeter Brune +   SNES_QN_RESTART_NONE - never restart
6940c777b0cSPeter Brune .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
6950c777b0cSPeter Brune -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
6960c777b0cSPeter Brune 
6970c777b0cSPeter Brune     Notes:
6980c777b0cSPeter Brune     The default line search used is the L2 line search and it requires two additional function evaluations.
6990c777b0cSPeter Brune 
7000c777b0cSPeter Brune .keywords: SNES, SNESQN, restart, type, set SNESLineSearch
7010c777b0cSPeter Brune @*/
7022150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
7032150357eSBarry Smith {
7040c777b0cSPeter Brune   PetscErrorCode ierr;
7056e111a19SKarl Rupp 
7060c777b0cSPeter Brune   PetscFunctionBegin;
7070c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7080c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
7090c777b0cSPeter Brune   PetscFunctionReturn(0);
7100c777b0cSPeter Brune }
7110c777b0cSPeter Brune 
7120c777b0cSPeter Brune #undef __FUNCT__
7130c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType"
7140c777b0cSPeter Brune /*@
7150c777b0cSPeter Brune     SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN.
7160c777b0cSPeter Brune 
7170c777b0cSPeter Brune     Logically Collective on SNES
7180c777b0cSPeter Brune 
7190c777b0cSPeter Brune     Input Parameters:
7200c777b0cSPeter Brune +   snes - the iterative context
7210c777b0cSPeter Brune -   stype - scale type
7220c777b0cSPeter Brune 
7230c777b0cSPeter Brune     Options Database:
7240c777b0cSPeter Brune .   -snes_qn_scale_type<shanno,none,linesearch,jacobian>
7250c777b0cSPeter Brune 
7260c777b0cSPeter Brune     Level: intermediate
7270c777b0cSPeter Brune 
7280c777b0cSPeter Brune     SNESQNSelectTypes:
7290c777b0cSPeter Brune +   SNES_QN_SCALE_NONE - don't scale the problem
7300c777b0cSPeter Brune .   SNES_QN_SCALE_SHANNO - use shanno scaling
7310c777b0cSPeter Brune .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
7320c777b0cSPeter Brune -   SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian.
7330c777b0cSPeter Brune 
7340c777b0cSPeter Brune .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch
7350c777b0cSPeter Brune @*/
7360c777b0cSPeter Brune 
7372150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
7382150357eSBarry Smith {
7390c777b0cSPeter Brune   PetscErrorCode ierr;
7406e111a19SKarl Rupp 
7410c777b0cSPeter Brune   PetscFunctionBegin;
7420c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7430c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
7440c777b0cSPeter Brune   PetscFunctionReturn(0);
7450c777b0cSPeter Brune }
7460c777b0cSPeter Brune 
7470c777b0cSPeter Brune #undef __FUNCT__
7480c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN"
7490adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
7500adebc6cSBarry Smith {
7510c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7526e111a19SKarl Rupp 
7530c777b0cSPeter Brune   PetscFunctionBegin;
7540c777b0cSPeter Brune   qn->scale_type = stype;
7550c777b0cSPeter Brune   PetscFunctionReturn(0);
7560c777b0cSPeter Brune }
7570c777b0cSPeter Brune 
7580c777b0cSPeter Brune #undef __FUNCT__
7590c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN"
7600adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
7610adebc6cSBarry Smith {
7620c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7636e111a19SKarl Rupp 
7640c777b0cSPeter Brune   PetscFunctionBegin;
7650c777b0cSPeter Brune   qn->restart_type = rtype;
7660c777b0cSPeter Brune   PetscFunctionReturn(0);
7670c777b0cSPeter Brune }
7680c777b0cSPeter Brune 
7691efc8c45SPeter Brune #undef __FUNCT__
7701efc8c45SPeter Brune #define __FUNCT__ "SNESQNSetType"
7711efc8c45SPeter Brune /*@
7721efc8c45SPeter Brune     SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN.
7731efc8c45SPeter Brune 
7741efc8c45SPeter Brune     Logically Collective on SNES
7751efc8c45SPeter Brune 
7761efc8c45SPeter Brune     Input Parameters:
7771efc8c45SPeter Brune +   snes - the iterative context
7781efc8c45SPeter Brune -   qtype - variant type
7791efc8c45SPeter Brune 
7801efc8c45SPeter Brune     Options Database:
7811efc8c45SPeter Brune .   -snes_qn_scale_type<lbfgs,broyden,badbroyden>
7821efc8c45SPeter Brune 
7831efc8c45SPeter Brune     Level: beginner
7841efc8c45SPeter Brune 
7851efc8c45SPeter Brune     SNESQNTypes:
7861efc8c45SPeter Brune +   SNES_QN_LBFGS - LBFGS variant
7871efc8c45SPeter Brune .   SNES_QN_BROYDEN - Broyden variant
7881efc8c45SPeter Brune -   SNES_QN_BADBROYDEN - Bad Broyden variant
7891efc8c45SPeter Brune 
7901efc8c45SPeter Brune .keywords: SNES, SNESQN, type, set
7911efc8c45SPeter Brune @*/
7921efc8c45SPeter Brune 
7931efc8c45SPeter Brune PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
7941efc8c45SPeter Brune {
7951efc8c45SPeter Brune   PetscErrorCode ierr;
7961efc8c45SPeter Brune 
7971efc8c45SPeter Brune   PetscFunctionBegin;
7981efc8c45SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7991efc8c45SPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));CHKERRQ(ierr);
8001efc8c45SPeter Brune   PetscFunctionReturn(0);
8011efc8c45SPeter Brune }
8021efc8c45SPeter Brune 
8031efc8c45SPeter Brune #undef __FUNCT__
8041efc8c45SPeter Brune #define __FUNCT__ "SNESQNSetType_QN"
8051efc8c45SPeter Brune PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
8061efc8c45SPeter Brune {
8071efc8c45SPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
8081efc8c45SPeter Brune 
8091efc8c45SPeter Brune   PetscFunctionBegin;
8101efc8c45SPeter Brune   qn->type = qtype;
8111efc8c45SPeter Brune   PetscFunctionReturn(0);
8121efc8c45SPeter Brune }
8131efc8c45SPeter Brune 
8144b11644fSPeter Brune /* -------------------------------------------------------------------------- */
8154b11644fSPeter Brune /*MC
8164b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
8174b11644fSPeter Brune 
8186cc8130cSPeter Brune       Options Database:
8196cc8130cSPeter Brune 
8206cc8130cSPeter Brune +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
8211867fe5bSPeter Brune .     -snes_qn_powell_angle - Angle condition for restart.
8221867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
82372835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
82444f7e39eSPeter Brune -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
8256cc8130cSPeter Brune 
826b8840d6bSPeter Brune       Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
827b8840d6bSPeter Brune       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
828b8840d6bSPeter Brune       updates.
8296cc8130cSPeter Brune 
8301867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
8311867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
8321867fe5bSPeter Brune       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
8331867fe5bSPeter Brune       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
8341867fe5bSPeter Brune 
8356cc8130cSPeter Brune       References:
8366cc8130cSPeter Brune 
8370a8e8854SPeter Brune       Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
8386cc8130cSPeter Brune 
8390a8e8854SPeter Brune       R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi-Newton Matrices and their use in Limited Memory Methods,
8400a8e8854SPeter Brune       Technical Report, Northwestern University, June 1992.
8410a8e8854SPeter Brune 
8420a8e8854SPeter Brune       Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
8430a8e8854SPeter Brune       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
844b8840d6bSPeter Brune 
8454b11644fSPeter Brune 
8464b11644fSPeter Brune       Level: beginner
8474b11644fSPeter Brune 
84804d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
8496cc8130cSPeter Brune 
8504b11644fSPeter Brune M*/
8514b11644fSPeter Brune #undef __FUNCT__
8524b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
8538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
8544b11644fSPeter Brune {
8554b11644fSPeter Brune   PetscErrorCode ierr;
8569f83bee8SJed Brown   SNES_QN        *qn;
8574b11644fSPeter Brune 
8584b11644fSPeter Brune   PetscFunctionBegin;
8594b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
8604b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
8614b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
8624b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
8635cd83419SPeter Brune   snes->ops->view           = SNESView_QN;
8644b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
8654b11644fSPeter Brune 
86647f26062SPeter Brune   snes->pcside = PC_LEFT;
86747f26062SPeter Brune 
86842f4f86dSBarry Smith   snes->usespc  = PETSC_TRUE;
86942f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
87042f4f86dSBarry Smith 
87188976e71SPeter Brune   if (!snes->tolerancesset) {
8720e444f03SPeter Brune     snes->max_funcs = 30000;
8730e444f03SPeter Brune     snes->max_its   = 10000;
87488976e71SPeter Brune   }
8750e444f03SPeter Brune 
876b00a9115SJed Brown   ierr                = PetscNewLog(snes,&qn);CHKERRQ(ierr);
8774b11644fSPeter Brune   snes->data          = (void*) qn;
8780ecc9e1dSPeter Brune   qn->m               = 10;
879b21d5a53SPeter Brune   qn->scaling         = 1.0;
8800298fd71SBarry Smith   qn->U               = NULL;
8810298fd71SBarry Smith   qn->V               = NULL;
8820298fd71SBarry Smith   qn->dXtdF           = NULL;
8830298fd71SBarry Smith   qn->dFtdX           = NULL;
8840298fd71SBarry Smith   qn->dXdFmat         = NULL;
8850298fd71SBarry Smith   qn->monitor         = NULL;
8861c6523dcSPeter Brune   qn->singlereduction = PETSC_TRUE;
887b8840d6bSPeter Brune   qn->powell_gamma    = 0.9999;
88860c08b40SPeter Brune   qn->scale_type      = SNES_QN_SCALE_DEFAULT;
88960c08b40SPeter Brune   qn->restart_type    = SNES_QN_RESTART_DEFAULT;
890b8840d6bSPeter Brune   qn->type            = SNES_QN_LBFGS;
891ea630c6eSPeter Brune 
892bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr);
893bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr);
8941efc8c45SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);CHKERRQ(ierr);
8954b11644fSPeter Brune   PetscFunctionReturn(0);
8964b11644fSPeter Brune }
897