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