xref: /petsc/src/snes/impls/qn/qn.c (revision c688d0420b4e513ff34944d1e1ad7d4e50aafa8d)
1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2 #include <petscdm.h>
3 
4 #define H(i,j)  qn->dXdFmat[i*qn->m + j]
5 
6 const char *const SNESQNScaleTypes[] =        {"DEFAULT","NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0};
7 const char *const SNESQNRestartTypes[] =      {"DEFAULT","NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0};
8 const char *const SNESQNTypes[] =             {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_",0};
9 
10 typedef struct {
11   Vec               *U;                   /* Stored past states (vary from method to method) */
12   Vec               *V;                   /* Stored past states (vary from method to method) */
13   PetscInt          m;                    /* The number of kept previous steps */
14   PetscReal         *lambda;              /* The line search history of the method */
15   PetscReal         *norm;                /* norms of the steps */
16   PetscScalar       *alpha, *beta;
17   PetscScalar       *dXtdF, *dFtdX, *YtdX;
18   PetscBool         singlereduction;      /* Aggregated reduction implementation */
19   PetscScalar       *dXdFmat;             /* A matrix of values for dX_i dot dF_j */
20   PetscViewer       monitor;
21   PetscReal         powell_gamma;         /* Powell angle restart condition */
22   PetscReal         scaling;              /* scaling of H0 */
23   SNESQNType        type;                 /* the type of quasi-newton method used */
24   SNESQNScaleType   scale_type;           /* the type of scaling used */
25   SNESQNRestartType restart_type;         /* determine the frequency and type of restart conditions */
26 } SNES_QN;
27 
28 #undef __FUNCT__
29 #define __FUNCT__ "SNESQNApply_Broyden"
30 PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D)
31 {
32   PetscErrorCode     ierr;
33   SNES_QN            *qn = (SNES_QN*)snes->data;
34   Vec                W   = snes->work[3];
35   Vec                *U  = qn->U;
36   PetscInt           m = qn->m;
37   PetscInt           k,i,j,lits,l = m;
38   PetscReal          unorm,a,b;
39   PetscReal          *lambda=qn->lambda;
40   PetscScalar        gdot;
41   PetscReal          udot;
42 
43   PetscFunctionBegin;
44   if (it < m) l = it;
45   if (it > 0) {
46     k = (it-1)%l;
47     ierr = SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);CHKERRQ(ierr);
48     ierr = VecCopy(Xold,U[k]);CHKERRQ(ierr);
49     ierr = VecAXPY(U[k],-1.0,X);CHKERRQ(ierr);
50     if (qn->monitor) {
51       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
52       ierr = PetscViewerASCIIPrintf(qn->monitor, "scaling vector %D of %D by lambda: %14.12e \n",k,l,(double)lambda[k]);CHKERRQ(ierr);
53       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
54     }
55   }
56   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
57     ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr);
58     SNESCheckKSPSolve(snes);
59     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
60     snes->linear_its += lits;
61     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
62   } else {
63     ierr = VecCopy(D,Y);CHKERRQ(ierr);
64     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
65   }
66 
67   /* inward recursion starting at the first update and working forward */
68   for (i = 0; i < l-1; i++) {
69     j = (it+i-l)%l;
70     k = (it+i-l+1)%l;
71     ierr = VecNorm(U[j],NORM_2,&unorm);CHKERRQ(ierr);
72     ierr = VecDot(U[j],Y,&gdot);CHKERRQ(ierr);
73     unorm *= unorm;
74     udot = PetscRealPart(gdot);
75     a = (lambda[j]/lambda[k]);
76     b = -(1.-lambda[j]);
77     a *= udot/unorm;
78     b *= udot/unorm;
79     ierr = VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);CHKERRQ(ierr);
80 
81     if (qn->monitor) {
82       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
83       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %D and %D, gdot: %14.12e\n",k,j,(double)PetscRealPart(gdot));CHKERRQ(ierr);
84       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
85     }
86   }
87   if (l > 0) {
88     k = (it-1)%l;
89     ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr);
90     ierr = VecNorm(U[k],NORM_2,&unorm);CHKERRQ(ierr);
91     unorm *= unorm;
92     udot = PetscRealPart(gdot);
93     a = unorm/(unorm-lambda[k]*udot);
94     b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot);
95     if (qn->monitor) {
96       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
97       ierr = PetscViewerASCIIPrintf(qn->monitor, "using vector %D: a: %14.12e b: %14.12e \n",k,(double)a,(double)b);CHKERRQ(ierr);
98       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
99     }
100     ierr = VecAXPBY(Y,b,a,U[k]);CHKERRQ(ierr);
101   }
102   l = m;
103   if (it+1<m)l=it+1;
104   k = it%l;
105   if (qn->monitor) {
106     ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
107     ierr = PetscViewerASCIIPrintf(qn->monitor, "setting vector %D of %D\n",k,l);CHKERRQ(ierr);
108     ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
109   }
110   PetscFunctionReturn(0);
111 }
112 
113 #undef __FUNCT__
114 #define __FUNCT__ "SNESQNApply_BadBroyden"
115 PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
116 {
117   PetscErrorCode ierr;
118   SNES_QN        *qn = (SNES_QN*)snes->data;
119   Vec            W   = snes->work[3];
120   Vec            *U  = qn->U;
121   Vec            *T  = qn->V;
122 
123   /* ksp thing for Jacobian scaling */
124   PetscInt           h,k,j,i,lits;
125   PetscInt           m = qn->m;
126   PetscScalar        gdot,udot;
127   PetscInt           l = m;
128 
129   PetscFunctionBegin;
130   if (it < m) l = it;
131   ierr = VecCopy(D,Y);CHKERRQ(ierr);
132   if (l > 0) {
133     k    = (it-1)%l;
134     ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);CHKERRQ(ierr);
135     ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr);
136     ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr);
137     ierr = VecCopy(Xold,T[k]);CHKERRQ(ierr);
138     ierr = VecAXPY(T[k],-1.0,X);CHKERRQ(ierr);
139   }
140 
141   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
142     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
143     SNESCheckKSPSolve(snes);
144     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
145     snes->linear_its += lits;
146     ierr              = VecCopy(W,Y);CHKERRQ(ierr);
147   } else {
148     ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr);
149   }
150 
151   /* inward recursion starting at the first update and working forward */
152   if (l > 0) {
153     for (i = 0; i < l-1; i++) {
154       j    = (it+i-l)%l;
155       k    = (it+i-l+1)%l;
156       h    = (it-1)%l;
157       ierr = VecDotBegin(U[j],U[h],&gdot);CHKERRQ(ierr);
158       ierr = VecDotBegin(U[j],U[j],&udot);CHKERRQ(ierr);
159       ierr = VecDotEnd(U[j],U[h],&gdot);CHKERRQ(ierr);
160       ierr = VecDotEnd(U[j],U[j],&udot);CHKERRQ(ierr);
161       ierr = VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);CHKERRQ(ierr);
162       ierr = VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);CHKERRQ(ierr);
163       if (qn->monitor) {
164         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
165         ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D gdot: %14.12e\n", it, k, (double)PetscRealPart(gdot));CHKERRQ(ierr);
166         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
167       }
168     }
169   }
170   PetscFunctionReturn(0);
171 }
172 
173 #undef __FUNCT__
174 #define __FUNCT__ "SNESQNApply_LBFGS"
175 PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
176 {
177   PetscErrorCode ierr;
178   SNES_QN        *qn    = (SNES_QN*)snes->data;
179   Vec            W      = snes->work[3];
180   Vec            *dX    = qn->U;
181   Vec            *dF    = qn->V;
182   PetscScalar    *alpha = qn->alpha;
183   PetscScalar    *beta  = qn->beta;
184   PetscScalar    *dXtdF = qn->dXtdF;
185   PetscScalar    *dFtdX = qn->dFtdX;
186   PetscScalar    *YtdX  = qn->YtdX;
187 
188   /* ksp thing for Jacobian scaling */
189   PetscInt           k,i,j,g,lits;
190   PetscInt           m = qn->m;
191   PetscScalar        t;
192   PetscInt           l = m;
193 
194   PetscFunctionBegin;
195   if (it < m) l = it;
196   ierr = VecCopy(D,Y);CHKERRQ(ierr);
197   if (it > 0) {
198     k    = (it - 1) % l;
199     ierr = VecCopy(D,dF[k]);CHKERRQ(ierr);
200     ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr);
201     ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
202     ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
203     if (qn->singlereduction) {
204       ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
205       ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
206       ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr);
207       ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr);
208       ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr);
209       ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr);
210       for (j = 0; j < l; j++) {
211         H(k, j) = dFtdX[j];
212         H(j, k) = dXtdF[j];
213       }
214       /* copy back over to make the computation of alpha and beta easier */
215       for (j = 0; j < l; j++) dXtdF[j] = H(j, j);
216     } else {
217       ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr);
218     }
219     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
220       ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr);
221     }
222   }
223 
224   /* outward recursion starting at iteration k's update and working back */
225   for (i=0; i<l; i++) {
226     k = (it-i-1)%l;
227     if (qn->singlereduction) {
228       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
229       t = YtdX[k];
230       for (j=0; j<i; j++) {
231         g  = (it-j-1)%l;
232         t -= alpha[g]*H(k, g);
233       }
234       alpha[k] = t/H(k,k);
235     } else {
236       ierr     = VecDot(dX[k],Y,&t);CHKERRQ(ierr);
237       alpha[k] = t/dXtdF[k];
238     }
239     if (qn->monitor) {
240       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
241       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D alpha:        %14.12e\n", it, k, (double)PetscRealPart(alpha[k]));CHKERRQ(ierr);
242       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
243     }
244     ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr);
245   }
246 
247   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
248     ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr);
249     SNESCheckKSPSolve(snes);
250     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
251     snes->linear_its += lits;
252     ierr              = VecCopy(W, Y);CHKERRQ(ierr);
253   } else {
254     ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
255   }
256   if (qn->singlereduction) {
257     ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr);
258   }
259   /* inward recursion starting at the first update and working forward */
260   for (i = 0; i < l; i++) {
261     k = (it + i - l) % l;
262     if (qn->singlereduction) {
263       t = YtdX[k];
264       for (j = 0; j < i; j++) {
265         g  = (it + j - l) % l;
266         t += (alpha[g] - beta[g])*H(g, k);
267       }
268       beta[k] = t / H(k, k);
269     } else {
270       ierr    = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
271       beta[k] = t / dXtdF[k];
272     }
273     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr);
274     if (qn->monitor) {
275       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
276       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %D k: %D alpha - beta: %14.12e\n", it, k, (double)PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
277       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
278     }
279   }
280   PetscFunctionReturn(0);
281 }
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "SNESSolve_QN"
285 static PetscErrorCode SNESSolve_QN(SNES snes)
286 {
287   PetscErrorCode       ierr;
288   SNES_QN              *qn = (SNES_QN*) snes->data;
289   Vec                  X,Xold;
290   Vec                  F,W;
291   Vec                  Y,D,Dold;
292   PetscInt             i, i_r;
293   PetscReal            fnorm,xnorm,ynorm,gnorm;
294   SNESLineSearchReason lssucceed;
295   PetscBool            powell,periodic;
296   PetscScalar          DolddotD,DolddotDold;
297   SNESConvergedReason  reason;
298 
299   /* basically just a regular newton's method except for the application of the Jacobian */
300 
301   PetscFunctionBegin;
302   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
303 
304   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
305   F    = snes->vec_func;                /* residual vector */
306   Y    = snes->vec_sol_update;          /* search direction generated by J^-1D*/
307   W    = snes->work[3];
308   X    = snes->vec_sol;                 /* solution vector */
309   Xold = snes->work[0];
310 
311   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
312   D    = snes->work[1];
313   Dold = snes->work[2];
314 
315   snes->reason = SNES_CONVERGED_ITERATING;
316 
317   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
318   snes->iter = 0;
319   snes->norm = 0.;
320   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
321 
322   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
323     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
324     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
325     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
326       snes->reason = SNES_DIVERGED_INNER;
327       PetscFunctionReturn(0);
328     }
329     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
330   } else {
331     if (!snes->vec_func_init_set) {
332       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
333     } else snes->vec_func_init_set = PETSC_FALSE;
334 
335     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
336     SNESCheckFunctionNorm(snes,fnorm);
337   }
338   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
339       ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr);
340       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
341       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
342         snes->reason = SNES_DIVERGED_INNER;
343         PetscFunctionReturn(0);
344       }
345   } else {
346     ierr = VecCopy(F,D);CHKERRQ(ierr);
347   }
348 
349   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
350   snes->norm = fnorm;
351   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
352   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
353   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
354 
355   /* test convergence */
356   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
357   if (snes->reason) PetscFunctionReturn(0);
358 
359   if (snes->pc && snes->pcside == PC_RIGHT) {
360     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
361     ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
362     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
363     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
364     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
365       snes->reason = SNES_DIVERGED_INNER;
366       PetscFunctionReturn(0);
367     }
368     ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
369     ierr = VecCopy(F,D);CHKERRQ(ierr);
370   }
371 
372   /* scale the initial update */
373   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
374     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
375     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
376   }
377 
378   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
379     if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) {
380       PetscScalar ff,xf;
381       ierr = VecCopy(Dold,Y);CHKERRQ(ierr);
382       ierr = VecCopy(Xold,W);CHKERRQ(ierr);
383       ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr);
384       ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr);
385       ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr);
386       ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr);
387       ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr);
388       ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr);
389       qn->scaling = PetscRealPart(xf)/PetscRealPart(ff);
390       if (qn->monitor) {
391         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
392         ierr = PetscViewerASCIIPrintf(qn->monitor, "Shanno scaling %D %g\n", i,(double)qn->scaling);CHKERRQ(ierr);
393         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
394       }
395     }
396     switch (qn->type) {
397     case SNES_QN_BADBROYDEN:
398       ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
399       break;
400     case SNES_QN_BROYDEN:
401       ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr);
402       break;
403     case SNES_QN_LBFGS:
404       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
405       break;
406     }
407     /* line search for lambda */
408     ynorm = 1; gnorm = fnorm;
409     ierr  = VecCopy(D, Dold);CHKERRQ(ierr);
410     ierr  = VecCopy(X, Xold);CHKERRQ(ierr);
411     ierr  = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
412     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
413     ierr = SNESLineSearchGetReason(snes->linesearch, &lssucceed);CHKERRQ(ierr);
414     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
415     if (lssucceed) {
416       if (++snes->numFailures >= snes->maxFailures) {
417         snes->reason = SNES_DIVERGED_LINE_SEARCH;
418         break;
419       }
420     }
421     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
422       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
423     }
424 
425     /* convergence monitoring */
426     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);
427 
428     if (snes->pc && snes->pcside == PC_RIGHT) {
429       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
430       ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
431       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
432       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
433       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
434         snes->reason = SNES_DIVERGED_INNER;
435         PetscFunctionReturn(0);
436       }
437       ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
438     }
439 
440     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
441     snes->norm = fnorm;
442 
443     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
444     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
445     /* set parameter for default relative tolerance convergence test */
446     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
447     if (snes->reason) PetscFunctionReturn(0);
448     if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
449       ierr = SNESApplyNPC(snes,X,F,D);CHKERRQ(ierr);
450       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
451       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
452         snes->reason = SNES_DIVERGED_INNER;
453         PetscFunctionReturn(0);
454       }
455     } else {
456       ierr = VecCopy(F, D);CHKERRQ(ierr);
457     }
458     powell = PETSC_FALSE;
459     if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) {
460       /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */
461       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
462         ierr = MatMult(snes->jacobian_pre,Dold,W);CHKERRQ(ierr);
463       } else {
464         ierr = VecCopy(Dold,W);CHKERRQ(ierr);
465       }
466       ierr = VecDotBegin(W, Dold, &DolddotDold);CHKERRQ(ierr);
467       ierr = VecDotBegin(W, D, &DolddotD);CHKERRQ(ierr);
468       ierr = VecDotEnd(W, Dold, &DolddotDold);CHKERRQ(ierr);
469       ierr = VecDotEnd(W, D, &DolddotD);CHKERRQ(ierr);
470       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
471     }
472     periodic = PETSC_FALSE;
473     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
474       if (i_r>qn->m-1) periodic = PETSC_TRUE;
475     }
476     /* restart if either powell or periodic restart is satisfied. */
477     if (powell || periodic) {
478       if (qn->monitor) {
479         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
480         if (powell) {
481           ierr = PetscViewerASCIIPrintf(qn->monitor, "Powell restart! |%14.12e| > %6.4f*|%14.12e| i_r = %D\n", (double)PetscRealPart(DolddotD), (double)qn->powell_gamma, (double)PetscRealPart(DolddotDold),i_r);CHKERRQ(ierr);
482         } else {
483           ierr = PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %D\n", i_r);CHKERRQ(ierr);
484         }
485         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
486       }
487       i_r = -1;
488       /* general purpose update */
489       if (snes->ops->update) {
490         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
491       }
492       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
493         ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
494       }
495     }
496     /* general purpose update */
497     if (snes->ops->update) {
498       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
499     }
500   }
501   if (i == snes->max_its) {
502     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
503     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
504   }
505   PetscFunctionReturn(0);
506 }
507 
508 #undef __FUNCT__
509 #define __FUNCT__ "SNESSetUp_QN"
510 static PetscErrorCode SNESSetUp_QN(SNES snes)
511 {
512   SNES_QN        *qn = (SNES_QN*)snes->data;
513   PetscErrorCode ierr;
514   DM             dm;
515 
516   PetscFunctionBegin;
517 
518   if (!snes->vec_sol) {
519     ierr             = SNESGetDM(snes,&dm);CHKERRQ(ierr);
520     ierr             = DMCreateGlobalVector(dm,&snes->vec_sol);CHKERRQ(ierr);
521   }
522 
523   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr);
524   if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr);
525   ierr = PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);CHKERRQ(ierr);
526 
527   if (qn->singlereduction) {
528     ierr = PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);CHKERRQ(ierr);
529   }
530   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
531   /* set method defaults */
532   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
533     if (qn->type == SNES_QN_BADBROYDEN) {
534       qn->scale_type = SNES_QN_SCALE_NONE;
535     } else {
536       qn->scale_type = SNES_QN_SCALE_SHANNO;
537     }
538   }
539   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
540     if (qn->type == SNES_QN_LBFGS) {
541       qn->restart_type = SNES_QN_RESTART_POWELL;
542     } else {
543       qn->restart_type = SNES_QN_RESTART_PERIODIC;
544     }
545   }
546 
547   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
548     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
549   }
550   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;}
551   PetscFunctionReturn(0);
552 }
553 
554 #undef __FUNCT__
555 #define __FUNCT__ "SNESReset_QN"
556 static PetscErrorCode SNESReset_QN(SNES snes)
557 {
558   PetscErrorCode ierr;
559   SNES_QN        *qn;
560 
561   PetscFunctionBegin;
562   if (snes->data) {
563     qn = (SNES_QN*)snes->data;
564     if (qn->U) {
565       ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr);
566     }
567     if (qn->V) {
568       ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr);
569     }
570     if (qn->singlereduction) {
571       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
572     }
573     ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr);
574   }
575   PetscFunctionReturn(0);
576 }
577 
578 #undef __FUNCT__
579 #define __FUNCT__ "SNESDestroy_QN"
580 static PetscErrorCode SNESDestroy_QN(SNES snes)
581 {
582   PetscErrorCode ierr;
583 
584   PetscFunctionBegin;
585   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
586   ierr = PetscFree(snes->data);CHKERRQ(ierr);
587   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
588   PetscFunctionReturn(0);
589 }
590 
591 #undef __FUNCT__
592 #define __FUNCT__ "SNESSetFromOptions_QN"
593 static PetscErrorCode SNESSetFromOptions_QN(PetscOptionItems *PetscOptionsObject,SNES snes)
594 {
595 
596   PetscErrorCode    ierr;
597   SNES_QN           *qn    = (SNES_QN*)snes->data;
598   PetscBool         monflg = PETSC_FALSE,flg;
599   SNESLineSearch    linesearch;
600   SNESQNRestartType rtype = qn->restart_type;
601   SNESQNScaleType   stype = qn->scale_type;
602   SNESQNType        qtype = qn->type;
603 
604   PetscFunctionBegin;
605   ierr = PetscOptionsHead(PetscOptionsObject,"SNES QN options");CHKERRQ(ierr);
606   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr);
607   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr);
608   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr);
609   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr);
610   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr);
611   if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr);
612 
613   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr);
614   if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr);
615 
616   ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg);CHKERRQ(ierr);
617   if (flg) {ierr = SNESQNSetType(snes,qtype);CHKERRQ(ierr);}
618   ierr = PetscOptionsTail();CHKERRQ(ierr);
619   if (!snes->linesearch) {
620     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
621     if (qn->type == SNES_QN_LBFGS) {
622       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
623     } else if (qn->type == SNES_QN_BROYDEN) {
624       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);CHKERRQ(ierr);
625     } else {
626       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
627     }
628   }
629   if (monflg) {
630     qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
631   }
632   PetscFunctionReturn(0);
633 }
634 
635 #undef __FUNCT__
636 #define __FUNCT__ "SNESView_QN"
637 static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
638 {
639   SNES_QN        *qn    = (SNES_QN*)snes->data;
640   PetscBool      iascii;
641   PetscErrorCode ierr;
642 
643   PetscFunctionBegin;
644   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
645   if (iascii) {
646     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);
647     ierr = PetscViewerASCIIPrintf(viewer,"  Stored subspace size: %D\n", qn->m);CHKERRQ(ierr);
648     if (qn->singlereduction) {
649       ierr = PetscViewerASCIIPrintf(viewer,"  Using the single reduction variant.\n");CHKERRQ(ierr);
650     }
651   }
652   PetscFunctionReturn(0);
653 }
654 
655 #undef __FUNCT__
656 #define __FUNCT__ "SNESQNSetRestartType"
657 /*@
658     SNESQNSetRestartType - Sets the restart type for SNESQN.
659 
660     Logically Collective on SNES
661 
662     Input Parameters:
663 +   snes - the iterative context
664 -   rtype - restart type
665 
666     Options Database:
667 +   -snes_qn_restart_type <powell,periodic,none> - set the restart type
668 -   -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
669 
670     Level: intermediate
671 
672     SNESQNRestartTypes:
673 +   SNES_QN_RESTART_NONE - never restart
674 .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
675 -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
676 
677 .keywords: SNES, SNESQN, restart, type, set SNESLineSearch, SNESQNRestartType
678 @*/
679 PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
680 {
681   PetscErrorCode ierr;
682 
683   PetscFunctionBegin;
684   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
685   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
686   PetscFunctionReturn(0);
687 }
688 
689 #undef __FUNCT__
690 #define __FUNCT__ "SNESQNSetScaleType"
691 /*@
692     SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN.
693 
694     Logically Collective on SNES
695 
696     Input Parameters:
697 +   snes - the iterative context
698 -   stype - scale type
699 
700     Options Database:
701 .   -snes_qn_scale_type <shanno,none,linesearch,jacobian>
702 
703     Level: intermediate
704 
705     SNESQNScaleTypes:
706 +   SNES_QN_SCALE_NONE - don't scale the problem
707 .   SNES_QN_SCALE_SHANNO - use shanno scaling
708 .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
709 -   SNES_QN_SCALE_JACOBIAN - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration
710                              of QN and at ever restart.
711 
712 .keywords: scaling type
713 
714 .seealso: SNES, SNESQN, SNESLineSearch, SNESQNScaleType, SNESetJacobian()
715 @*/
716 
717 PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
718 {
719   PetscErrorCode ierr;
720 
721   PetscFunctionBegin;
722   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
723   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
724   PetscFunctionReturn(0);
725 }
726 
727 #undef __FUNCT__
728 #define __FUNCT__ "SNESQNSetScaleType_QN"
729 PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
730 {
731   SNES_QN *qn = (SNES_QN*)snes->data;
732 
733   PetscFunctionBegin;
734   qn->scale_type = stype;
735   PetscFunctionReturn(0);
736 }
737 
738 #undef __FUNCT__
739 #define __FUNCT__ "SNESQNSetRestartType_QN"
740 PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
741 {
742   SNES_QN *qn = (SNES_QN*)snes->data;
743 
744   PetscFunctionBegin;
745   qn->restart_type = rtype;
746   PetscFunctionReturn(0);
747 }
748 
749 #undef __FUNCT__
750 #define __FUNCT__ "SNESQNSetType"
751 /*@
752     SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN.
753 
754     Logically Collective on SNES
755 
756     Input Parameters:
757 +   snes - the iterative context
758 -   qtype - variant type
759 
760     Options Database:
761 .   -snes_qn_type <lbfgs,broyden,badbroyden>
762 
763     Level: beginner
764 
765     SNESQNTypes:
766 +   SNES_QN_LBFGS - LBFGS variant
767 .   SNES_QN_BROYDEN - Broyden variant
768 -   SNES_QN_BADBROYDEN - Bad Broyden variant
769 
770 .keywords: SNES, SNESQN, type, set, SNESQNType
771 @*/
772 
773 PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
774 {
775   PetscErrorCode ierr;
776 
777   PetscFunctionBegin;
778   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
779   ierr = PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));CHKERRQ(ierr);
780   PetscFunctionReturn(0);
781 }
782 
783 #undef __FUNCT__
784 #define __FUNCT__ "SNESQNSetType_QN"
785 PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
786 {
787   SNES_QN *qn = (SNES_QN*)snes->data;
788 
789   PetscFunctionBegin;
790   qn->type = qtype;
791   PetscFunctionReturn(0);
792 }
793 
794 /* -------------------------------------------------------------------------- */
795 /*MC
796       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
797 
798       Options Database:
799 
800 +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
801 +     -snes_qn_restart_type <powell,periodic,none> - set the restart type
802 .     -snes_qn_powell_gamma - Angle condition for restart.
803 .     -snes_qn_powell_descent - Descent condition for restart.
804 .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
805 .     -snes_qn_scale_type <shanno,none,linesearch,jacobian> - scaling performed on inner Jacobian
806 .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
807 -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.
808 
809       Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
810       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
811       updates.
812 
813       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
814       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
815       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
816       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
817 
818       Uses left nonlinear preconditioning by default.
819 
820       References:
821 +   1. -   Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
822 .   2. -   R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods,
823       Technical Report, Northwestern University, June 1992.
824 .   3. -   Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
825       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
826 -   4. -   Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
827        SIAM Review, 57(4), 2015
828 
829       Level: beginner
830 
831 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
832 
833 M*/
834 #undef __FUNCT__
835 #define __FUNCT__ "SNESCreate_QN"
836 PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
837 {
838   PetscErrorCode ierr;
839   SNES_QN        *qn;
840 
841   PetscFunctionBegin;
842   snes->ops->setup          = SNESSetUp_QN;
843   snes->ops->solve          = SNESSolve_QN;
844   snes->ops->destroy        = SNESDestroy_QN;
845   snes->ops->setfromoptions = SNESSetFromOptions_QN;
846   snes->ops->view           = SNESView_QN;
847   snes->ops->reset          = SNESReset_QN;
848 
849   snes->pcside = PC_LEFT;
850 
851   snes->usespc  = PETSC_TRUE;
852   snes->usesksp = PETSC_FALSE;
853 
854   snes->alwayscomputesfinalresidual = PETSC_TRUE;
855 
856   if (!snes->tolerancesset) {
857     snes->max_funcs = 30000;
858     snes->max_its   = 10000;
859   }
860 
861   ierr                = PetscNewLog(snes,&qn);CHKERRQ(ierr);
862   snes->data          = (void*) qn;
863   qn->m               = 10;
864   qn->scaling         = 1.0;
865   qn->U               = NULL;
866   qn->V               = NULL;
867   qn->dXtdF           = NULL;
868   qn->dFtdX           = NULL;
869   qn->dXdFmat         = NULL;
870   qn->monitor         = NULL;
871   qn->singlereduction = PETSC_TRUE;
872   qn->powell_gamma    = 0.9999;
873   qn->scale_type      = SNES_QN_SCALE_DEFAULT;
874   qn->restart_type    = SNES_QN_RESTART_DEFAULT;
875   qn->type            = SNES_QN_LBFGS;
876 
877   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr);
878   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr);
879   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);CHKERRQ(ierr);
880   PetscFunctionReturn(0);
881 }
882