xref: /petsc/src/snes/impls/qn/qn.c (revision aeb49b86d700c2d4a06d6f7a22069c109820b654)
1af0996ceSBarry Smith #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   PetscInt           m = qn->m;
375c7a0a03SPeter Brune   PetscInt           k,i,j,lits,l = m;
385c7a0a03SPeter Brune   PetscReal          unorm,a,b;
395c7a0a03SPeter Brune   PetscReal          *lambda=qn->lambda;
40b8840d6bSPeter Brune   PetscScalar        gdot;
4159e7931aSPeter Brune   PetscReal          udot;
422150357eSBarry Smith 
43b8840d6bSPeter Brune   PetscFunctionBegin;
44b8840d6bSPeter Brune   if (it < m) l = it;
455c7a0a03SPeter Brune   if (it > 0) {
465c7a0a03SPeter Brune     k = (it-1)%l;
475c7a0a03SPeter Brune     ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr);
48a49fa0d8SPeter Brune     ierr = VecCopy(Xold,U[k]);CHKERRQ(ierr);
495051edb1SPeter Brune     ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr);
505c7a0a03SPeter Brune     if (qn->monitor) {
515c7a0a03SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
525c7a0a03SPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %d of %d by lambda: %14.12e \n",k,l,lambda[k]);CHKERRQ(ierr);
535c7a0a03SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
545c7a0a03SPeter Brune     }
555c7a0a03SPeter Brune   }
56b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
57d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr);
58422a814eSBarry Smith     SNESCheckKSPSolve(snes);
59b8840d6bSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
60b8840d6bSPeter Brune     snes->linear_its += lits;
61b8840d6bSPeter Brune     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
62b8840d6bSPeter Brune   } else {
63b8840d6bSPeter Brune     ierr = VecCopy(D,Y);CHKERRQ(ierr);
64b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
65b8840d6bSPeter Brune   }
66b8840d6bSPeter Brune 
67b8840d6bSPeter Brune   /* inward recursion starting at the first update and working forward */
68b8840d6bSPeter Brune   for (i = 0; i < l-1; i++) {
695c7a0a03SPeter Brune     j = (it+i-l)%l;
705c7a0a03SPeter Brune     k = (it+i-l+1)%l;
715c7a0a03SPeter Brune     ierr = VecNorm(U[j],NORM_2,&unorm);CHKERRQ(ierr);
725c7a0a03SPeter Brune     ierr = VecDot(U[j],Y,&gdot);CHKERRQ(ierr);
7359e7931aSPeter Brune     unorm *= unorm;
7459e7931aSPeter Brune     udot = PetscRealPart(gdot);
755c7a0a03SPeter Brune     a = (lambda[j]/lambda[k]);
765c7a0a03SPeter Brune     b = -(1.-lambda[j]);
7759e7931aSPeter Brune     a *= udot/unorm;
7859e7931aSPeter Brune     b *= udot/unorm;
795c7a0a03SPeter Brune     ierr = VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);CHKERRQ(ierr);
805c7a0a03SPeter Brune 
81b8840d6bSPeter Brune     if (qn->monitor) {
82b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
83e639b251SJed Brown       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d and %d, gdot: %14.12e\n",k,j,PetscRealPart(gdot));CHKERRQ(ierr);
84b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
85b8840d6bSPeter Brune     }
86b8840d6bSPeter Brune   }
87b8840d6bSPeter Brune   if (l > 0) {
88b8840d6bSPeter Brune     k = (it-1)%l;
89b8840d6bSPeter Brune     ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
905c7a0a03SPeter Brune     ierr = VecNorm(U[k],NORM_2,&unorm);CHKERRQ(ierr);
9159e7931aSPeter Brune     unorm *= unorm;
9259e7931aSPeter Brune     udot = PetscRealPart(gdot);
9359e7931aSPeter Brune     a = unorm/(unorm-lambda[k]*udot);
9459e7931aSPeter Brune     b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot);
95b8840d6bSPeter Brune     if (qn->monitor) {
96b8840d6bSPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
9759e7931aSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %d: a: %14.12e b: %14.12e \n",k,a,b);CHKERRQ(ierr);
98b8840d6bSPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
99b8840d6bSPeter Brune     }
1005c7a0a03SPeter Brune     ierr = VecAXPBY(Y,b,a,U[k]);CHKERRQ(ierr);
101b8840d6bSPeter Brune   }
1025c7a0a03SPeter Brune   l = m;
1035c7a0a03SPeter Brune   if (it+1<m)l=it+1;
1045c7a0a03SPeter Brune   k = it%l;
1055c7a0a03SPeter Brune   if (qn->monitor) {
1065c7a0a03SPeter Brune     ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1075c7a0a03SPeter Brune     ierr = PetscViewerASCIIPrintf(qn->monitor, "setting vector %d of %d\n",k,l);CHKERRQ(ierr);
1085c7a0a03SPeter Brune     ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1095c7a0a03SPeter Brune   }
110b8840d6bSPeter Brune   PetscFunctionReturn(0);
111b8840d6bSPeter Brune }
112b8840d6bSPeter Brune 
113b8840d6bSPeter Brune #undef __FUNCT__
114b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_BadBroyden"
1150adebc6cSBarry Smith PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1160adebc6cSBarry Smith {
117b8840d6bSPeter Brune   PetscErrorCode ierr;
118b8840d6bSPeter Brune   SNES_QN        *qn = (SNES_QN*)snes->data;
119b8840d6bSPeter Brune   Vec            W   = snes->work[3];
120b8840d6bSPeter Brune   Vec            *U  = qn->U;
121b8840d6bSPeter Brune   Vec            *T  = qn->V;
122b8840d6bSPeter Brune 
12384c577daSBarry Smith   /* ksp thing for Jacobian scaling */
124c9e59058SPeter Brune   PetscInt           h,k,j,i,lits;
125b8840d6bSPeter Brune   PetscInt           m = qn->m;
1265051edb1SPeter Brune   PetscScalar        gdot,udot;
127b8840d6bSPeter Brune   PetscInt           l = m;
1280adebc6cSBarry Smith 
129b8840d6bSPeter Brune   PetscFunctionBegin;
130b8840d6bSPeter Brune   if (it < m) l = it;
131b8840d6bSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
132b8840d6bSPeter Brune   if (l > 0) {
133b8840d6bSPeter Brune     k    = (it-1)%l;
134c9e59058SPeter Brune     ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);CHKERRQ(ierr);
135b8840d6bSPeter Brune     ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr);
136b8840d6bSPeter Brune     ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr);
1375051edb1SPeter Brune     ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr);
1385051edb1SPeter Brune     ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr);
139b8840d6bSPeter Brune   }
140b8840d6bSPeter Brune 
141b8840d6bSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
142d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
143422a814eSBarry Smith     SNESCheckKSPSolve(snes);
144b8840d6bSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
145b8840d6bSPeter Brune     snes->linear_its += lits;
146b8840d6bSPeter Brune     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
147b8840d6bSPeter Brune   } else {
148b8840d6bSPeter Brune     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
149b8840d6bSPeter Brune   }
1505051edb1SPeter Brune 
1515051edb1SPeter Brune   /* inward recursion starting at the first update and working forward */
1525051edb1SPeter Brune   if (l > 0) {
1535051edb1SPeter Brune     for (i = 0; i < l-1; i++) {
154c9e59058SPeter Brune       j    = (it+i-l)%l;
155c9e59058SPeter Brune       k    = (it+i-l+1)%l;
156c9e59058SPeter Brune       h    = (it-1)%l;
157c9e59058SPeter Brune       ierr = VecDotBegin(U[j],U[h],&gdot);CHKERRQ(ierr);
158c9e59058SPeter Brune       ierr = VecDotBegin(U[j],U[j],&udot);CHKERRQ(ierr);
159c9e59058SPeter Brune       ierr = VecDotEnd(U[j],U[h],&gdot);CHKERRQ(ierr);
160c9e59058SPeter Brune       ierr = VecDotEnd(U[j],U[j],&udot);CHKERRQ(ierr);
1615051edb1SPeter Brune       ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr);
162c9e59058SPeter Brune       ierr = VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);CHKERRQ(ierr);
1635051edb1SPeter Brune       if (qn->monitor) {
1645051edb1SPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1655051edb1SPeter Brune         ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr);
1665051edb1SPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
1675051edb1SPeter Brune       }
1685051edb1SPeter Brune     }
1695051edb1SPeter Brune   }
170b8840d6bSPeter Brune   PetscFunctionReturn(0);
171b8840d6bSPeter Brune }
172b8840d6bSPeter Brune 
173b8840d6bSPeter Brune #undef __FUNCT__
174b8840d6bSPeter Brune #define __FUNCT__ "SNESQNApply_LBFGS"
1750adebc6cSBarry Smith PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
1760adebc6cSBarry Smith {
177b8840d6bSPeter Brune   PetscErrorCode ierr;
178b8840d6bSPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
179b8840d6bSPeter Brune   Vec            W      = snes->work[3];
180b8840d6bSPeter Brune   Vec            *dX    = qn->U;
181b8840d6bSPeter Brune   Vec            *dF    = qn->V;
182bd052dfeSPeter Brune   PetscScalar    *alpha = qn->alpha;
183bd052dfeSPeter Brune   PetscScalar    *beta  = qn->beta;
1848e231d97SPeter Brune   PetscScalar    *dXtdF = qn->dXtdF;
185b8840d6bSPeter Brune   PetscScalar    *dFtdX = qn->dFtdX;
1868e231d97SPeter Brune   PetscScalar    *YtdX  = qn->YtdX;
187bd052dfeSPeter Brune 
18884c577daSBarry Smith   /* ksp thing for Jacobian scaling */
1898e231d97SPeter Brune   PetscInt           k,i,j,g,lits;
1904b11644fSPeter Brune   PetscInt           m = qn->m;
1914b11644fSPeter Brune   PetscScalar        t;
1924b11644fSPeter Brune   PetscInt           l = m;
1938e231d97SPeter Brune 
1944b11644fSPeter Brune   PetscFunctionBegin;
1955ba6227bSPeter Brune   if (it < m) l = it;
1961c6523dcSPeter Brune   ierr = VecCopy(D,Y);CHKERRQ(ierr);
197b8840d6bSPeter Brune   if (it > 0) {
198b8840d6bSPeter Brune     k    = (it - 1) % l;
199b8840d6bSPeter Brune     ierr = VecCopy(D,dF[k]);CHKERRQ(ierr);
200b8840d6bSPeter Brune     ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
201b8840d6bSPeter Brune     ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
202b8840d6bSPeter Brune     ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
20318aa0c0cSPeter Brune     if (qn->singlereduction) {
2041c6523dcSPeter Brune       ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
2051c6523dcSPeter Brune       ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
2061c6523dcSPeter Brune       ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr);
2071c6523dcSPeter Brune       ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
2081c6523dcSPeter Brune       ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
2091c6523dcSPeter Brune       ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr);
210b8840d6bSPeter Brune       for (j = 0; j < l; j++) {
211b8840d6bSPeter Brune         H(k, j) = dFtdX[j];
212b8840d6bSPeter Brune         H(j, k) = dXtdF[j];
213b8840d6bSPeter Brune       }
214b8840d6bSPeter Brune       /* copy back over to make the computation of alpha and beta easier */
2151aa26658SKarl Rupp       for (j = 0; j < l; j++) dXtdF[j] = H(j, j);
216b8840d6bSPeter Brune     } else {
217b8840d6bSPeter Brune       ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr);
21823d44fbcSPeter Brune     }
21923d44fbcSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
220b8840d6bSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr);
221b8840d6bSPeter Brune     }
222b8840d6bSPeter Brune   }
223b8840d6bSPeter Brune 
2244b11644fSPeter Brune   /* outward recursion starting at iteration k's update and working back */
2254b11644fSPeter Brune   for (i=0; i<l; i++) {
226b21d5a53SPeter Brune     k = (it-i-1)%l;
22718aa0c0cSPeter Brune     if (qn->singlereduction) {
2288e231d97SPeter Brune       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
2298e231d97SPeter Brune       t = YtdX[k];
2308e231d97SPeter Brune       for (j=0; j<i; j++) {
2318e231d97SPeter Brune         g  = (it-j-1)%l;
232b3cfccb2SPeter Brune         t -= alpha[g]*H(k, g);
2338e231d97SPeter Brune       }
2348e231d97SPeter Brune       alpha[k] = t/H(k,k);
2358e231d97SPeter Brune     } else {
2363af51624SPeter Brune       ierr     = VecDot(dX[k],Y,&t);CHKERRQ(ierr);
2378e231d97SPeter Brune       alpha[k] = t/dXtdF[k];
2388e231d97SPeter Brune     }
23944f7e39eSPeter Brune     if (qn->monitor) {
2403af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2415ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
2423af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
24344f7e39eSPeter Brune     }
2446bf1b2e5SPeter Brune     ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr);
2454b11644fSPeter Brune   }
2464b11644fSPeter Brune 
2470c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
248d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
249422a814eSBarry Smith     SNESCheckKSPSolve(snes);
2500ecc9e1dSPeter Brune     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2510ecc9e1dSPeter Brune     snes->linear_its += lits;
252b8840d6bSPeter Brune     ierr              = VecCopy(W, Y);CHKERRQ(ierr);
2530ecc9e1dSPeter Brune   } else {
254b21d5a53SPeter Brune     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
2550ecc9e1dSPeter Brune   }
25618aa0c0cSPeter Brune   if (qn->singlereduction) {
257b8840d6bSPeter Brune     ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr);
2588e231d97SPeter Brune   }
2594b11644fSPeter Brune   /* inward recursion starting at the first update and working forward */
260bd052dfeSPeter Brune   for (i = 0; i < l; i++) {
261b21d5a53SPeter Brune     k = (it + i - l) % l;
26218aa0c0cSPeter Brune     if (qn->singlereduction) {
2638e231d97SPeter Brune       t = YtdX[k];
2648e231d97SPeter Brune       for (j = 0; j < i; j++) {
2658e231d97SPeter Brune         g  = (it + j - l) % l;
266b3cfccb2SPeter Brune         t += (alpha[g] - beta[g])*H(g, k);
2678e231d97SPeter Brune       }
2688e231d97SPeter Brune       beta[k] = t / H(k, k);
2698e231d97SPeter Brune     } else {
2706bf1b2e5SPeter Brune       ierr    = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
2718e231d97SPeter Brune       beta[k] = t / dXtdF[k];
2728e231d97SPeter Brune     }
27322d28d08SBarry Smith     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr);
27444f7e39eSPeter Brune     if (qn->monitor) {
2753af51624SPeter Brune       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
2765ba6227bSPeter Brune       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
2773af51624SPeter Brune       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
27844f7e39eSPeter Brune     }
2794b11644fSPeter Brune   }
2804b11644fSPeter Brune   PetscFunctionReturn(0);
2814b11644fSPeter Brune }
2824b11644fSPeter Brune 
2834b11644fSPeter Brune #undef __FUNCT__
2844b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN"
2854b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes)
2864b11644fSPeter Brune {
2874b11644fSPeter Brune   PetscErrorCode       ierr;
2889f83bee8SJed Brown   SNES_QN              *qn = (SNES_QN*) snes->data;
28915f5eeeaSPeter Brune   Vec                  X,Xold;
2900a3905e1SPeter Brune   Vec                  F,W;
29147f26062SPeter Brune   Vec                  Y,D,Dold;
292b8840d6bSPeter Brune   PetscInt             i, i_r;
293335fb975SPeter Brune   PetscReal            fnorm,xnorm,ynorm,gnorm;
294422a814eSBarry Smith   SNESLineSearchReason lssucceed;
295422a814eSBarry Smith   PetscBool            powell,periodic;
2961c6523dcSPeter Brune   PetscScalar          DolddotD,DolddotDold;
297b7281c8aSPeter Brune   SNESConvergedReason  reason;
2980ecc9e1dSPeter Brune 
29984c577daSBarry Smith   /* basically just a regular newton's method except for the application of the Jacobian */
3004b11644fSPeter Brune 
3016e111a19SKarl Rupp   PetscFunctionBegin;
302c579b300SPatrick Farrell 
303c579b300SPatrick Farrell   if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
304c579b300SPatrick Farrell     SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
305c579b300SPatrick Farrell   }
306c579b300SPatrick Farrell 
307fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
3089f3a0142SPeter Brune   F    = snes->vec_func;                /* residual vector */
3093af51624SPeter Brune   Y    = snes->vec_sol_update;          /* search direction generated by J^-1D*/
3100a3905e1SPeter Brune   W    = snes->work[3];
311b8840d6bSPeter Brune   X    = snes->vec_sol;                 /* solution vector */
312335fb975SPeter Brune   Xold = snes->work[0];
3133af51624SPeter Brune 
3143af51624SPeter Brune   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
315335fb975SPeter Brune   D    = snes->work[1];
316335fb975SPeter Brune   Dold = snes->work[2];
3174b11644fSPeter Brune 
3184b11644fSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
3194b11644fSPeter Brune 
320e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3214b11644fSPeter Brune   snes->iter = 0;
3224b11644fSPeter Brune   snes->norm = 0.;
323e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
32447f26062SPeter Brune 
325b28a06ddSPeter Brune   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
326be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
327b7281c8aSPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
328b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
329b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
330b7281c8aSPeter Brune       PetscFunctionReturn(0);
331b7281c8aSPeter Brune     }
33247f26062SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
33347f26062SPeter Brune   } else {
334e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
33515f5eeeaSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3361aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
337c1c75074SPeter Brune 
33847f26062SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
339422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
34047f26062SPeter Brune   }
341b28a06ddSPeter Brune   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
342be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr);
343b7281c8aSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
344b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
345b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
346b7281c8aSPeter Brune         PetscFunctionReturn(0);
347b7281c8aSPeter Brune       }
34847f26062SPeter Brune   } else {
34947f26062SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
35047f26062SPeter Brune   }
351b28a06ddSPeter Brune 
352e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3534b11644fSPeter Brune   snes->norm = fnorm;
354e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
355a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
3564b11644fSPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
3574b11644fSPeter Brune 
3584b11644fSPeter Brune   /* test convergence */
3594b11644fSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3604b11644fSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
36170d3b23bSPeter Brune 
3623cf07b75SPeter Brune   if (snes->pc && snes->pcside == PC_RIGHT) {
3633cf07b75SPeter Brune     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
3643cf07b75SPeter Brune     ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
3653cf07b75SPeter Brune     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
366ddd40ce5SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
367ddd40ce5SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
368ddd40ce5SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
369ddd40ce5SPeter Brune       PetscFunctionReturn(0);
370ddd40ce5SPeter Brune     }
371be95d8f1SBarry Smith     ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
3723cf07b75SPeter Brune     ierr = VecCopy(F,D);CHKERRQ(ierr);
3733cf07b75SPeter Brune   }
3743cf07b75SPeter Brune 
375f8e15203SPeter Brune   /* scale the initial update */
3760c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
377d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
378*aeb49b86SAsbjørn Nilsen Riseth     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
3790ecc9e1dSPeter Brune   }
3800ecc9e1dSPeter Brune 
3815ba6227bSPeter Brune   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
3820a3905e1SPeter Brune     if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) {
3830a3905e1SPeter Brune       PetscScalar ff,xf;
3840a3905e1SPeter Brune       ierr = VecCopy(Dold,Y);CHKERRQ(ierr);
3850a3905e1SPeter Brune       ierr = VecCopy(Xold,W);CHKERRQ(ierr);
3860a3905e1SPeter Brune       ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr);
3870a3905e1SPeter Brune       ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr);
3880a3905e1SPeter Brune       ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr);
3890a3905e1SPeter Brune       ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr);
3900a3905e1SPeter Brune       ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr);
3910a3905e1SPeter Brune       ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr);
3920a3905e1SPeter Brune       qn->scaling = PetscRealPart(xf)/PetscRealPart(ff);
3930a3905e1SPeter Brune     }
394b8840d6bSPeter Brune     switch (qn->type) {
395b8840d6bSPeter Brune     case SNES_QN_BADBROYDEN:
396b8840d6bSPeter Brune       ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
397b8840d6bSPeter Brune       break;
398b8840d6bSPeter Brune     case SNES_QN_BROYDEN:
3995051edb1SPeter Brune       ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr);
400b8840d6bSPeter Brune       break;
401b8840d6bSPeter Brune     case SNES_QN_LBFGS:
402b8840d6bSPeter Brune       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
403b8840d6bSPeter Brune       break;
404b8840d6bSPeter Brune     }
40570d3b23bSPeter Brune     /* line search for lambda */
40670d3b23bSPeter Brune     ynorm = 1; gnorm = fnorm;
4070a3905e1SPeter Brune     ierr  = VecCopy(D, Dold);CHKERRQ(ierr);
40815f5eeeaSPeter Brune     ierr  = VecCopy(X, Xold);CHKERRQ(ierr);
409f1c6b773SPeter Brune     ierr  = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
4109f3a0142SPeter Brune     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
411422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lssucceed);CHKERRQ(ierr);
412422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
413422a814eSBarry Smith     if (lssucceed) {
4149f3a0142SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
4159f3a0142SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
4169f3a0142SPeter Brune         break;
4179f3a0142SPeter Brune       }
4189f3a0142SPeter Brune     }
4190c777b0cSPeter Brune     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
42005ee83feSPeter Brune       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
4210ecc9e1dSPeter Brune     }
4223af51624SPeter Brune 
42388d374b2SPeter Brune     /* convergence monitoring */
424b21d5a53SPeter 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);
425b21d5a53SPeter Brune 
426b28a06ddSPeter Brune     if (snes->pc && snes->pcside == PC_RIGHT) {
427b28a06ddSPeter Brune       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
428b28a06ddSPeter Brune       ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
429b28a06ddSPeter Brune       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
430ddd40ce5SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
431ddd40ce5SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
432ddd40ce5SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
433ddd40ce5SPeter Brune         PetscFunctionReturn(0);
434ddd40ce5SPeter Brune       }
435be95d8f1SBarry Smith       ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
436b28a06ddSPeter Brune     }
437b28a06ddSPeter Brune 
438360c497dSPeter Brune     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
43971dbe336SPeter Brune     snes->norm = fnorm;
440360c497dSPeter Brune 
441a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
4428409ca45SMatthew G Knepley     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
4434b11644fSPeter Brune     /* set parameter for default relative tolerance convergence test */
444d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
4454b11644fSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
446b28a06ddSPeter Brune     if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
447be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr);
448b7281c8aSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
449b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
450b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
451b7281c8aSPeter Brune         PetscFunctionReturn(0);
452b7281c8aSPeter Brune       }
45388d374b2SPeter Brune     } else {
45488d374b2SPeter Brune       ierr = VecCopy(F, D);CHKERRQ(ierr);
45588d374b2SPeter Brune     }
4560c777b0cSPeter Brune     powell = PETSC_FALSE;
4570c777b0cSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_POWELL) {
4586bf1b2e5SPeter Brune       /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
45923c918e7SPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
46023c918e7SPeter Brune         ierr = MatMult(snes->jacobian_pre,Dold,W);CHKERRQ(ierr);
46123c918e7SPeter Brune       } else {
46223c918e7SPeter Brune         ierr = VecCopy(Dold,W);CHKERRQ(ierr);
46323c918e7SPeter Brune       }
46423c918e7SPeter Brune       ierr = VecDotBegin(W, Dold, &DolddotDold);CHKERRQ(ierr);
46523c918e7SPeter Brune       ierr = VecDotBegin(W, D, &DolddotD);CHKERRQ(ierr);
46623c918e7SPeter Brune       ierr = VecDotEnd(W, Dold, &DolddotDold);CHKERRQ(ierr);
46723c918e7SPeter Brune       ierr = VecDotEnd(W, D, &DolddotD);CHKERRQ(ierr);
4680c777b0cSPeter Brune       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
4690c777b0cSPeter Brune     }
4700c777b0cSPeter Brune     periodic = PETSC_FALSE;
471b8840d6bSPeter Brune     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
472b8840d6bSPeter Brune       if (i_r>qn->m-1) periodic = PETSC_TRUE;
4730c777b0cSPeter Brune     }
4740c777b0cSPeter Brune     /* restart if either powell or periodic restart is satisfied. */
4750c777b0cSPeter Brune     if (powell || periodic) {
4765ba6227bSPeter Brune       if (qn->monitor) {
4775ba6227bSPeter Brune         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
478516fe3c3SPeter 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);
4795ba6227bSPeter Brune         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
4805ba6227bSPeter Brune       }
4815ba6227bSPeter Brune       i_r = -1;
4825ba6227bSPeter Brune       /* general purpose update */
4835ba6227bSPeter Brune       if (snes->ops->update) {
4845ba6227bSPeter Brune         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
4855ba6227bSPeter Brune       }
4860c777b0cSPeter Brune       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
487d1e9a80fSBarry Smith         ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
4880ecc9e1dSPeter Brune       }
4890ecc9e1dSPeter Brune     }
49070d3b23bSPeter Brune     /* general purpose update */
49170d3b23bSPeter Brune     if (snes->ops->update) {
49270d3b23bSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
49370d3b23bSPeter Brune     }
4945ba6227bSPeter Brune   }
4954b11644fSPeter Brune   if (i == snes->max_its) {
4964b11644fSPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
4974b11644fSPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
4984b11644fSPeter Brune   }
4994b11644fSPeter Brune   PetscFunctionReturn(0);
5004b11644fSPeter Brune }
5014b11644fSPeter Brune 
5024b11644fSPeter Brune #undef __FUNCT__
5034b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN"
5044b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes)
5054b11644fSPeter Brune {
5069f83bee8SJed Brown   SNES_QN        *qn = (SNES_QN*)snes->data;
5074b11644fSPeter Brune   PetscErrorCode ierr;
508335fb975SPeter Brune 
5094b11644fSPeter Brune   PetscFunctionBegin;
510b8840d6bSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr);
51159e7931aSPeter Brune   if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr);
512dcca6d9dSJed Brown   ierr = PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);CHKERRQ(ierr);
5138e231d97SPeter Brune 
51418aa0c0cSPeter Brune   if (qn->singlereduction) {
515dcca6d9dSJed Brown     ierr = PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);CHKERRQ(ierr);
5168e231d97SPeter Brune   }
517fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
51860c08b40SPeter Brune   /* set method defaults */
5191efc8c45SPeter Brune   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
52060c08b40SPeter Brune     if (qn->type == SNES_QN_BADBROYDEN) {
52160c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_NONE;
52260c08b40SPeter Brune     } else {
52360c08b40SPeter Brune       qn->scale_type = SNES_QN_SCALE_SHANNO;
52460c08b40SPeter Brune     }
52560c08b40SPeter Brune   }
5261efc8c45SPeter Brune   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
52760c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
52860c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_POWELL;
52960c08b40SPeter Brune     } else {
53060c08b40SPeter Brune       qn->restart_type = SNES_QN_RESTART_PERIODIC;
53160c08b40SPeter Brune     }
53260c08b40SPeter Brune   }
53360c08b40SPeter Brune 
5340c777b0cSPeter Brune   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
5358e231d97SPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
5368e231d97SPeter Brune   }
5376c67d002SPeter Brune   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;}
5384b11644fSPeter Brune   PetscFunctionReturn(0);
5394b11644fSPeter Brune }
5404b11644fSPeter Brune 
5414b11644fSPeter Brune #undef __FUNCT__
5424b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN"
5434b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes)
5444b11644fSPeter Brune {
5454b11644fSPeter Brune   PetscErrorCode ierr;
5469f83bee8SJed Brown   SNES_QN        *qn;
5470adebc6cSBarry Smith 
5484b11644fSPeter Brune   PetscFunctionBegin;
5494b11644fSPeter Brune   if (snes->data) {
5509f83bee8SJed Brown     qn = (SNES_QN*)snes->data;
551b8840d6bSPeter Brune     if (qn->U) {
552b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr);
5534b11644fSPeter Brune     }
554b8840d6bSPeter Brune     if (qn->V) {
555b8840d6bSPeter Brune       ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr);
5564b11644fSPeter Brune     }
55718aa0c0cSPeter Brune     if (qn->singlereduction) {
5588e231d97SPeter Brune       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
5598e231d97SPeter Brune     }
5605c7a0a03SPeter Brune     ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr);
5614b11644fSPeter Brune   }
5624b11644fSPeter Brune   PetscFunctionReturn(0);
5634b11644fSPeter Brune }
5644b11644fSPeter Brune 
5654b11644fSPeter Brune #undef __FUNCT__
5664b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN"
5674b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes)
5684b11644fSPeter Brune {
5694b11644fSPeter Brune   PetscErrorCode ierr;
5706e111a19SKarl Rupp 
5714b11644fSPeter Brune   PetscFunctionBegin;
5724b11644fSPeter Brune   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
5734b11644fSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
574bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
5754b11644fSPeter Brune   PetscFunctionReturn(0);
5764b11644fSPeter Brune }
5774b11644fSPeter Brune 
5784b11644fSPeter Brune #undef __FUNCT__
5794b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN"
5808c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_QN(PetscOptions *PetscOptionsObject,SNES snes)
5814b11644fSPeter Brune {
5824b11644fSPeter Brune 
5834b11644fSPeter Brune   PetscErrorCode    ierr;
5842150357eSBarry Smith   SNES_QN           *qn    = (SNES_QN*)snes->data;
58588f769c5SPeter Brune   PetscBool         monflg = PETSC_FALSE,flg;
586f1c6b773SPeter Brune   SNESLineSearch    linesearch;
5872150357eSBarry Smith   SNESQNRestartType rtype = qn->restart_type;
5882150357eSBarry Smith   SNESQNScaleType   stype = qn->scale_type;
5891efc8c45SPeter Brune   SNESQNType        qtype = qn->type;
5902150357eSBarry Smith 
5914b11644fSPeter Brune   PetscFunctionBegin;
592e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES QN options");CHKERRQ(ierr);
5930298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr);
5940298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr);
5950298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr);
5960298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr);
5970298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr);
59888f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr);
59988f769c5SPeter Brune   if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr);
60088f769c5SPeter Brune 
60188f769c5SPeter Brune   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr);
60288f769c5SPeter Brune   if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr);
60388f769c5SPeter Brune 
6040fdccdaeSBarry Smith   ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg);CHKERRQ(ierr);
6051efc8c45SPeter Brune   if (flg) {ierr = SNESQNSetType(snes,qtype);CHKERRQ(ierr);}
6064b11644fSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
6079e764e56SPeter Brune   if (!snes->linesearch) {
6087601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
60960c08b40SPeter Brune     if (qn->type == SNES_QN_LBFGS) {
6101a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
61160c08b40SPeter Brune     } else if (qn->type == SNES_QN_BROYDEN) {
61260c08b40SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);CHKERRQ(ierr);
61360c08b40SPeter Brune     } else {
61460c08b40SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
61560c08b40SPeter Brune     }
6169e764e56SPeter Brune   }
61744f7e39eSPeter Brune   if (monflg) {
618ce94432eSBarry Smith     qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
61944f7e39eSPeter Brune   }
6204b11644fSPeter Brune   PetscFunctionReturn(0);
6214b11644fSPeter Brune }
6224b11644fSPeter Brune 
6235cd83419SPeter Brune #undef __FUNCT__
6245cd83419SPeter Brune #define __FUNCT__ "SNESView_QN"
6255cd83419SPeter Brune static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
6265cd83419SPeter Brune {
6275cd83419SPeter Brune   SNES_QN        *qn    = (SNES_QN*)snes->data;
6285cd83419SPeter Brune   PetscBool      iascii;
6295cd83419SPeter Brune   PetscErrorCode ierr;
6305cd83419SPeter Brune 
6315cd83419SPeter Brune   PetscFunctionBegin;
6325cd83419SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
6335cd83419SPeter Brune   if (iascii) {
6345cd83419SPeter 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);
6355cd83419SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  Stored subspace size: %d\n", qn->m);CHKERRQ(ierr);
6365cd83419SPeter Brune     if (qn->singlereduction) {
6375cd83419SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  Using the single reduction variant.\n");CHKERRQ(ierr);
6385cd83419SPeter Brune     }
6395cd83419SPeter Brune   }
6405cd83419SPeter Brune   PetscFunctionReturn(0);
6415cd83419SPeter Brune }
64292c02d66SPeter Brune 
6430c777b0cSPeter Brune #undef __FUNCT__
6440c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType"
6450c777b0cSPeter Brune /*@
6460c777b0cSPeter Brune     SNESQNSetRestartType - Sets the restart type for SNESQN.
6470c777b0cSPeter Brune 
6480c777b0cSPeter Brune     Logically Collective on SNES
6490c777b0cSPeter Brune 
6500c777b0cSPeter Brune     Input Parameters:
6510c777b0cSPeter Brune +   snes - the iterative context
6520c777b0cSPeter Brune -   rtype - restart type
6530c777b0cSPeter Brune 
6540c777b0cSPeter Brune     Options Database:
6550c777b0cSPeter Brune +   -snes_qn_restart_type <powell,periodic,none> - set the restart type
65684c577daSBarry Smith -   -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
6570c777b0cSPeter Brune 
6580c777b0cSPeter Brune     Level: intermediate
6590c777b0cSPeter Brune 
6600c777b0cSPeter Brune     SNESQNRestartTypes:
6610c777b0cSPeter Brune +   SNES_QN_RESTART_NONE - never restart
6620c777b0cSPeter Brune .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
6630c777b0cSPeter Brune -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
6640c777b0cSPeter Brune 
66584c577daSBarry Smith .keywords: SNES, SNESQN, restart, type, set SNESLineSearch, SNESQNRestartType
6660c777b0cSPeter Brune @*/
6672150357eSBarry Smith PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
6682150357eSBarry Smith {
6690c777b0cSPeter Brune   PetscErrorCode ierr;
6706e111a19SKarl Rupp 
6710c777b0cSPeter Brune   PetscFunctionBegin;
6720c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6730c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
6740c777b0cSPeter Brune   PetscFunctionReturn(0);
6750c777b0cSPeter Brune }
6760c777b0cSPeter Brune 
6770c777b0cSPeter Brune #undef __FUNCT__
6780c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType"
6790c777b0cSPeter Brune /*@
68084c577daSBarry Smith     SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN.
6810c777b0cSPeter Brune 
6820c777b0cSPeter Brune     Logically Collective on SNES
6830c777b0cSPeter Brune 
6840c777b0cSPeter Brune     Input Parameters:
6850c777b0cSPeter Brune +   snes - the iterative context
6860c777b0cSPeter Brune -   stype - scale type
6870c777b0cSPeter Brune 
6880c777b0cSPeter Brune     Options Database:
6890c777b0cSPeter Brune .   -snes_qn_scale_type <shanno,none,linesearch,jacobian>
6900c777b0cSPeter Brune 
6910c777b0cSPeter Brune     Level: intermediate
6920c777b0cSPeter Brune 
69384c577daSBarry Smith     SNESQNScaleTypes:
6940c777b0cSPeter Brune +   SNES_QN_SCALE_NONE - don't scale the problem
6950c777b0cSPeter Brune .   SNES_QN_SCALE_SHANNO - use shanno scaling
6960c777b0cSPeter Brune .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
6970c777b0cSPeter Brune -   SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian.
6980c777b0cSPeter Brune 
69984c577daSBarry Smith .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch, SNESQNScaleType
7000c777b0cSPeter Brune @*/
7010c777b0cSPeter Brune 
7022150357eSBarry Smith PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
7032150357eSBarry Smith {
7040c777b0cSPeter Brune   PetscErrorCode ierr;
7056e111a19SKarl Rupp 
7060c777b0cSPeter Brune   PetscFunctionBegin;
7070c777b0cSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7080c777b0cSPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
7090c777b0cSPeter Brune   PetscFunctionReturn(0);
7100c777b0cSPeter Brune }
7110c777b0cSPeter Brune 
7120c777b0cSPeter Brune #undef __FUNCT__
7130c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetScaleType_QN"
7140adebc6cSBarry Smith PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
7150adebc6cSBarry Smith {
7160c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7176e111a19SKarl Rupp 
7180c777b0cSPeter Brune   PetscFunctionBegin;
7190c777b0cSPeter Brune   qn->scale_type = stype;
7200c777b0cSPeter Brune   PetscFunctionReturn(0);
7210c777b0cSPeter Brune }
7220c777b0cSPeter Brune 
7230c777b0cSPeter Brune #undef __FUNCT__
7240c777b0cSPeter Brune #define __FUNCT__ "SNESQNSetRestartType_QN"
7250adebc6cSBarry Smith PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
7260adebc6cSBarry Smith {
7270c777b0cSPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7286e111a19SKarl Rupp 
7290c777b0cSPeter Brune   PetscFunctionBegin;
7300c777b0cSPeter Brune   qn->restart_type = rtype;
7310c777b0cSPeter Brune   PetscFunctionReturn(0);
7320c777b0cSPeter Brune }
7330c777b0cSPeter Brune 
7341efc8c45SPeter Brune #undef __FUNCT__
7351efc8c45SPeter Brune #define __FUNCT__ "SNESQNSetType"
7361efc8c45SPeter Brune /*@
7371efc8c45SPeter Brune     SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN.
7381efc8c45SPeter Brune 
7391efc8c45SPeter Brune     Logically Collective on SNES
7401efc8c45SPeter Brune 
7411efc8c45SPeter Brune     Input Parameters:
7421efc8c45SPeter Brune +   snes - the iterative context
7431efc8c45SPeter Brune -   qtype - variant type
7441efc8c45SPeter Brune 
7451efc8c45SPeter Brune     Options Database:
74684c577daSBarry Smith .   -snes_qn_type <lbfgs,broyden,badbroyden>
7471efc8c45SPeter Brune 
7481efc8c45SPeter Brune     Level: beginner
7491efc8c45SPeter Brune 
7501efc8c45SPeter Brune     SNESQNTypes:
7511efc8c45SPeter Brune +   SNES_QN_LBFGS - LBFGS variant
7521efc8c45SPeter Brune .   SNES_QN_BROYDEN - Broyden variant
7531efc8c45SPeter Brune -   SNES_QN_BADBROYDEN - Bad Broyden variant
7541efc8c45SPeter Brune 
75584c577daSBarry Smith .keywords: SNES, SNESQN, type, set, SNESQNType
7561efc8c45SPeter Brune @*/
7571efc8c45SPeter Brune 
7581efc8c45SPeter Brune PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
7591efc8c45SPeter Brune {
7601efc8c45SPeter Brune   PetscErrorCode ierr;
7611efc8c45SPeter Brune 
7621efc8c45SPeter Brune   PetscFunctionBegin;
7631efc8c45SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7641efc8c45SPeter Brune   ierr = PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));CHKERRQ(ierr);
7651efc8c45SPeter Brune   PetscFunctionReturn(0);
7661efc8c45SPeter Brune }
7671efc8c45SPeter Brune 
7681efc8c45SPeter Brune #undef __FUNCT__
7691efc8c45SPeter Brune #define __FUNCT__ "SNESQNSetType_QN"
7701efc8c45SPeter Brune PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
7711efc8c45SPeter Brune {
7721efc8c45SPeter Brune   SNES_QN *qn = (SNES_QN*)snes->data;
7731efc8c45SPeter Brune 
7741efc8c45SPeter Brune   PetscFunctionBegin;
7751efc8c45SPeter Brune   qn->type = qtype;
7761efc8c45SPeter Brune   PetscFunctionReturn(0);
7771efc8c45SPeter Brune }
7781efc8c45SPeter Brune 
7794b11644fSPeter Brune /* -------------------------------------------------------------------------- */
7804b11644fSPeter Brune /*MC
7814b11644fSPeter Brune       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
7824b11644fSPeter Brune 
7836cc8130cSPeter Brune       Options Database:
7846cc8130cSPeter Brune 
78584c577daSBarry Smith +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
78684c577daSBarry Smith +     -snes_qn_restart_type <powell,periodic,none> - set the restart type
7871867fe5bSPeter Brune .     -snes_qn_powell_angle - Angle condition for restart.
7881867fe5bSPeter Brune .     -snes_qn_powell_descent - Descent condition for restart.
78984c577daSBarry Smith .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
79084c577daSBarry Smith .     -snes_qn_scale_type <shanno,none,linesearch,jacobian> - scaling performed on inner Jacobian
79172835e02SPeter Brune .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
79284c577daSBarry Smith -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.
7936cc8130cSPeter Brune 
794b8840d6bSPeter Brune       Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
795b8840d6bSPeter Brune       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
796b8840d6bSPeter Brune       updates.
7976cc8130cSPeter Brune 
7981867fe5bSPeter Brune       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
7991867fe5bSPeter Brune       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
80084c577daSBarry Smith       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
80184c577daSBarry Smith       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
8021867fe5bSPeter Brune 
8036cc8130cSPeter Brune       References:
8046cc8130cSPeter Brune 
8050a8e8854SPeter Brune       Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
8066cc8130cSPeter Brune 
8070a8e8854SPeter Brune       R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi-Newton Matrices and their use in Limited Memory Methods,
8080a8e8854SPeter Brune       Technical Report, Northwestern University, June 1992.
8090a8e8854SPeter Brune 
8100a8e8854SPeter Brune       Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
8110a8e8854SPeter Brune       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
812b8840d6bSPeter Brune 
8134b11644fSPeter Brune 
8144b11644fSPeter Brune       Level: beginner
8154b11644fSPeter Brune 
81604d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
8176cc8130cSPeter Brune 
8184b11644fSPeter Brune M*/
8194b11644fSPeter Brune #undef __FUNCT__
8204b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN"
8218cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
8224b11644fSPeter Brune {
8234b11644fSPeter Brune   PetscErrorCode ierr;
8249f83bee8SJed Brown   SNES_QN        *qn;
8254b11644fSPeter Brune 
8264b11644fSPeter Brune   PetscFunctionBegin;
8274b11644fSPeter Brune   snes->ops->setup          = SNESSetUp_QN;
8284b11644fSPeter Brune   snes->ops->solve          = SNESSolve_QN;
8294b11644fSPeter Brune   snes->ops->destroy        = SNESDestroy_QN;
8304b11644fSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_QN;
8315cd83419SPeter Brune   snes->ops->view           = SNESView_QN;
8324b11644fSPeter Brune   snes->ops->reset          = SNESReset_QN;
8334b11644fSPeter Brune 
83447f26062SPeter Brune   snes->pcside = PC_LEFT;
83547f26062SPeter Brune 
83642f4f86dSBarry Smith   snes->usespc  = PETSC_TRUE;
83742f4f86dSBarry Smith   snes->usesksp = PETSC_FALSE;
83842f4f86dSBarry Smith 
83988976e71SPeter Brune   if (!snes->tolerancesset) {
8400e444f03SPeter Brune     snes->max_funcs = 30000;
8410e444f03SPeter Brune     snes->max_its   = 10000;
84288976e71SPeter Brune   }
8430e444f03SPeter Brune 
844b00a9115SJed Brown   ierr                = PetscNewLog(snes,&qn);CHKERRQ(ierr);
8454b11644fSPeter Brune   snes->data          = (void*) qn;
8460ecc9e1dSPeter Brune   qn->m               = 10;
847b21d5a53SPeter Brune   qn->scaling         = 1.0;
8480298fd71SBarry Smith   qn->U               = NULL;
8490298fd71SBarry Smith   qn->V               = NULL;
8500298fd71SBarry Smith   qn->dXtdF           = NULL;
8510298fd71SBarry Smith   qn->dFtdX           = NULL;
8520298fd71SBarry Smith   qn->dXdFmat         = NULL;
8530298fd71SBarry Smith   qn->monitor         = NULL;
8541c6523dcSPeter Brune   qn->singlereduction = PETSC_TRUE;
855b8840d6bSPeter Brune   qn->powell_gamma    = 0.9999;
85660c08b40SPeter Brune   qn->scale_type      = SNES_QN_SCALE_DEFAULT;
85760c08b40SPeter Brune   qn->restart_type    = SNES_QN_RESTART_DEFAULT;
858b8840d6bSPeter Brune   qn->type            = SNES_QN_LBFGS;
859ea630c6eSPeter Brune 
860bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr);
861bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr);
8621efc8c45SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);CHKERRQ(ierr);
8634b11644fSPeter Brune   PetscFunctionReturn(0);
8644b11644fSPeter Brune }
865