xref: /petsc/src/snes/impls/qn/qn.c (revision e4cb33bb7dbdbae9285fba102465ca0f1dcb3977)
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   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
327   F    = snes->vec_func;                /* residual vector */
328   Y    = snes->vec_sol_update;          /* search direction generated by J^-1D*/
329   W    = snes->work[3];
330   X    = snes->vec_sol;                 /* solution vector */
331   Xold = snes->work[0];
332 
333   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
334   D    = snes->work[1];
335   Dold = snes->work[2];
336 
337   snes->reason = SNES_CONVERGED_ITERATING;
338 
339   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
340   snes->iter = 0;
341   snes->norm = 0.;
342   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
343 
344   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
345     ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr);
346     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
347     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
348       snes->reason = SNES_DIVERGED_INNER;
349       PetscFunctionReturn(0);
350     }
351     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
352   } else {
353     if (!snes->vec_func_init_set) {
354       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
355       if (snes->domainerror) {
356         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
357         PetscFunctionReturn(0);
358       }
359     } else snes->vec_func_init_set = PETSC_FALSE;
360 
361     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
362     if (PetscIsInfOrNanReal(fnorm)) {
363       snes->reason = SNES_DIVERGED_FNORM_NAN;
364       PetscFunctionReturn(0);
365     }
366   }
367   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
368       ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr);
369       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
370       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
371         snes->reason = SNES_DIVERGED_INNER;
372         PetscFunctionReturn(0);
373       }
374   } else {
375     ierr = VecCopy(F,D);CHKERRQ(ierr);
376   }
377 
378   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
379   snes->norm = fnorm;
380   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
381   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
382   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
383 
384   /* test convergence */
385   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
386   if (snes->reason) PetscFunctionReturn(0);
387 
388   if (snes->pc && snes->pcside == PC_RIGHT) {
389     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
390     ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
391     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
392     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
393     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
394       snes->reason = SNES_DIVERGED_INNER;
395       PetscFunctionReturn(0);
396     }
397     ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
398     ierr = VecCopy(F,D);CHKERRQ(ierr);
399   }
400 
401   /* scale the initial update */
402   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
403     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
404     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
405   }
406 
407   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
408     if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) {
409       PetscScalar ff,xf;
410       ierr = VecCopy(Dold,Y);CHKERRQ(ierr);
411       ierr = VecCopy(Xold,W);CHKERRQ(ierr);
412       ierr = VecAXPY(Y,-1.0,D);CHKERRQ(ierr);
413       ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr);
414       ierr = VecDotBegin(Y,Y,&ff);CHKERRQ(ierr);
415       ierr = VecDotBegin(W,Y,&xf);CHKERRQ(ierr);
416       ierr = VecDotEnd(Y,Y,&ff);CHKERRQ(ierr);
417       ierr = VecDotEnd(W,Y,&xf);CHKERRQ(ierr);
418       qn->scaling = PetscRealPart(xf)/PetscRealPart(ff);
419     }
420     switch (qn->type) {
421     case SNES_QN_BADBROYDEN:
422       ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
423       break;
424     case SNES_QN_BROYDEN:
425       ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);CHKERRQ(ierr);
426       break;
427     case SNES_QN_LBFGS:
428       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
429       break;
430     }
431     /* line search for lambda */
432     ynorm = 1; gnorm = fnorm;
433     ierr  = VecCopy(D, Dold);CHKERRQ(ierr);
434     ierr  = VecCopy(X, Xold);CHKERRQ(ierr);
435     ierr  = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
436     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
437     if (snes->domainerror) {
438       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
439       PetscFunctionReturn(0);
440     }
441     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr);
442     if (!lssucceed) {
443       if (++snes->numFailures >= snes->maxFailures) {
444         snes->reason = SNES_DIVERGED_LINE_SEARCH;
445         break;
446       }
447     }
448     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
449     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
450       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
451     }
452 
453     /* convergence monitoring */
454     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);
455 
456     if (snes->pc && snes->pcside == PC_RIGHT) {
457       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
458       ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
459       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
460       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
461       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
462         snes->reason = SNES_DIVERGED_INNER;
463         PetscFunctionReturn(0);
464       }
465       ierr = SNESGetPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
466     }
467 
468     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
469     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
470 
471     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
472     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
473     /* set parameter for default relative tolerance convergence test */
474     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
475     if (snes->reason) PetscFunctionReturn(0);
476     if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
477       ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr);
478       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
479       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
480         snes->reason = SNES_DIVERGED_INNER;
481         PetscFunctionReturn(0);
482       }
483     } else {
484       ierr = VecCopy(F, D);CHKERRQ(ierr);
485     }
486     powell = PETSC_FALSE;
487     if (qn->restart_type == SNES_QN_RESTART_POWELL) {
488       /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
489       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
490         ierr = MatMult(snes->jacobian_pre,Dold,W);CHKERRQ(ierr);
491       } else {
492         ierr = VecCopy(Dold,W);CHKERRQ(ierr);
493       }
494       ierr = VecDotBegin(W, Dold, &DolddotDold);CHKERRQ(ierr);
495       ierr = VecDotBegin(W, D, &DolddotD);CHKERRQ(ierr);
496       ierr = VecDotEnd(W, Dold, &DolddotDold);CHKERRQ(ierr);
497       ierr = VecDotEnd(W, D, &DolddotD);CHKERRQ(ierr);
498       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
499     }
500     periodic = PETSC_FALSE;
501     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
502       if (i_r>qn->m-1) periodic = PETSC_TRUE;
503     }
504     /* restart if either powell or periodic restart is satisfied. */
505     if (powell || periodic) {
506       if (qn->monitor) {
507         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
508         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);
509         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
510       }
511       i_r = -1;
512       /* general purpose update */
513       if (snes->ops->update) {
514         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
515       }
516       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
517         ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
518         ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
519       }
520     }
521     /* general purpose update */
522     if (snes->ops->update) {
523       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
524     }
525   }
526   if (i == snes->max_its) {
527     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
528     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
529   }
530   PetscFunctionReturn(0);
531 }
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "SNESSetUp_QN"
535 static PetscErrorCode SNESSetUp_QN(SNES snes)
536 {
537   SNES_QN        *qn = (SNES_QN*)snes->data;
538   PetscErrorCode ierr;
539 
540   PetscFunctionBegin;
541   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr);
542   if (qn->type != SNES_QN_BROYDEN) ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr);
543   ierr = PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);CHKERRQ(ierr);
544 
545   if (qn->singlereduction) {
546     ierr = PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);CHKERRQ(ierr);
547   }
548   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
549   /* set method defaults */
550   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
551     if (qn->type == SNES_QN_BADBROYDEN) {
552       qn->scale_type = SNES_QN_SCALE_NONE;
553     } else {
554       qn->scale_type = SNES_QN_SCALE_SHANNO;
555     }
556   }
557   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
558     if (qn->type == SNES_QN_LBFGS) {
559       qn->restart_type = SNES_QN_RESTART_POWELL;
560     } else {
561       qn->restart_type = SNES_QN_RESTART_PERIODIC;
562     }
563   }
564 
565   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
566     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
567   }
568   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;}
569   PetscFunctionReturn(0);
570 }
571 
572 #undef __FUNCT__
573 #define __FUNCT__ "SNESReset_QN"
574 static PetscErrorCode SNESReset_QN(SNES snes)
575 {
576   PetscErrorCode ierr;
577   SNES_QN        *qn;
578 
579   PetscFunctionBegin;
580   if (snes->data) {
581     qn = (SNES_QN*)snes->data;
582     if (qn->U) {
583       ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr);
584     }
585     if (qn->V) {
586       ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr);
587     }
588     if (qn->singlereduction) {
589       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
590     }
591     ierr = PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);CHKERRQ(ierr);
592   }
593   PetscFunctionReturn(0);
594 }
595 
596 #undef __FUNCT__
597 #define __FUNCT__ "SNESDestroy_QN"
598 static PetscErrorCode SNESDestroy_QN(SNES snes)
599 {
600   PetscErrorCode ierr;
601 
602   PetscFunctionBegin;
603   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
604   ierr = PetscFree(snes->data);CHKERRQ(ierr);
605   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
606   PetscFunctionReturn(0);
607 }
608 
609 #undef __FUNCT__
610 #define __FUNCT__ "SNESSetFromOptions_QN"
611 static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
612 {
613 
614   PetscErrorCode    ierr;
615   SNES_QN           *qn    = (SNES_QN*)snes->data;
616   PetscBool         monflg = PETSC_FALSE,flg;
617   SNESLineSearch    linesearch;
618   SNESQNRestartType rtype = qn->restart_type;
619   SNESQNScaleType   stype = qn->scale_type;
620   SNESQNType        qtype = qn->type;
621 
622   PetscFunctionBegin;
623   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
624   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr);
625   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr);
626   ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr);
627   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr);
628   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr);
629   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr);
630   if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr);
631 
632   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr);
633   if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr);
634 
635   ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,
636                           (PetscEnum)qtype,(PetscEnum*)&qtype,&flg);CHKERRQ(ierr);
637   if (flg) {ierr = SNESQNSetType(snes,qtype);CHKERRQ(ierr);}
638   ierr = PetscOptionsTail();CHKERRQ(ierr);
639   if (!snes->linesearch) {
640     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
641     if (qn->type == SNES_QN_LBFGS) {
642       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
643     } else if (qn->type == SNES_QN_BROYDEN) {
644       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);CHKERRQ(ierr);
645     } else {
646       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
647     }
648   }
649   if (monflg) {
650     qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
651   }
652   PetscFunctionReturn(0);
653 }
654 
655 #undef __FUNCT__
656 #define __FUNCT__ "SNESView_QN"
657 static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
658 {
659   SNES_QN        *qn    = (SNES_QN*)snes->data;
660   PetscBool      iascii;
661   PetscErrorCode ierr;
662 
663   PetscFunctionBegin;
664   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
665   if (iascii) {
666     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);
667     ierr = PetscViewerASCIIPrintf(viewer,"  Stored subspace size: %d\n", qn->m);CHKERRQ(ierr);
668     if (qn->singlereduction) {
669       ierr = PetscViewerASCIIPrintf(viewer,"  Using the single reduction variant.\n");CHKERRQ(ierr);
670     }
671   }
672   PetscFunctionReturn(0);
673 }
674 
675 #undef __FUNCT__
676 #define __FUNCT__ "SNESQNSetRestartType"
677 /*@
678     SNESQNSetRestartType - Sets the restart type for SNESQN.
679 
680     Logically Collective on SNES
681 
682     Input Parameters:
683 +   snes - the iterative context
684 -   rtype - restart type
685 
686     Options Database:
687 +   -snes_qn_restart_type<powell,periodic,none> - set the restart type
688 -   -snes_qn_restart[10] - sets the number of iterations before restart for periodic
689 
690     Level: intermediate
691 
692     SNESQNRestartTypes:
693 +   SNES_QN_RESTART_NONE - never restart
694 .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
695 -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
696 
697     Notes:
698     The default line search used is the L2 line search and it requires two additional function evaluations.
699 
700 .keywords: SNES, SNESQN, restart, type, set SNESLineSearch
701 @*/
702 PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
703 {
704   PetscErrorCode ierr;
705 
706   PetscFunctionBegin;
707   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
708   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
709   PetscFunctionReturn(0);
710 }
711 
712 #undef __FUNCT__
713 #define __FUNCT__ "SNESQNSetScaleType"
714 /*@
715     SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN.
716 
717     Logically Collective on SNES
718 
719     Input Parameters:
720 +   snes - the iterative context
721 -   stype - scale type
722 
723     Options Database:
724 .   -snes_qn_scale_type<shanno,none,linesearch,jacobian>
725 
726     Level: intermediate
727 
728     SNESQNSelectTypes:
729 +   SNES_QN_SCALE_NONE - don't scale the problem
730 .   SNES_QN_SCALE_SHANNO - use shanno scaling
731 .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
732 -   SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian.
733 
734 .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch
735 @*/
736 
737 PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
738 {
739   PetscErrorCode ierr;
740 
741   PetscFunctionBegin;
742   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
743   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
744   PetscFunctionReturn(0);
745 }
746 
747 #undef __FUNCT__
748 #define __FUNCT__ "SNESQNSetScaleType_QN"
749 PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
750 {
751   SNES_QN *qn = (SNES_QN*)snes->data;
752 
753   PetscFunctionBegin;
754   qn->scale_type = stype;
755   PetscFunctionReturn(0);
756 }
757 
758 #undef __FUNCT__
759 #define __FUNCT__ "SNESQNSetRestartType_QN"
760 PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
761 {
762   SNES_QN *qn = (SNES_QN*)snes->data;
763 
764   PetscFunctionBegin;
765   qn->restart_type = rtype;
766   PetscFunctionReturn(0);
767 }
768 
769 #undef __FUNCT__
770 #define __FUNCT__ "SNESQNSetType"
771 /*@
772     SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN.
773 
774     Logically Collective on SNES
775 
776     Input Parameters:
777 +   snes - the iterative context
778 -   qtype - variant type
779 
780     Options Database:
781 .   -snes_qn_scale_type<lbfgs,broyden,badbroyden>
782 
783     Level: beginner
784 
785     SNESQNTypes:
786 +   SNES_QN_LBFGS - LBFGS variant
787 .   SNES_QN_BROYDEN - Broyden variant
788 -   SNES_QN_BADBROYDEN - Bad Broyden variant
789 
790 .keywords: SNES, SNESQN, type, set
791 @*/
792 
793 PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
794 {
795   PetscErrorCode ierr;
796 
797   PetscFunctionBegin;
798   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
799   ierr = PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));CHKERRQ(ierr);
800   PetscFunctionReturn(0);
801 }
802 
803 #undef __FUNCT__
804 #define __FUNCT__ "SNESQNSetType_QN"
805 PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
806 {
807   SNES_QN *qn = (SNES_QN*)snes->data;
808 
809   PetscFunctionBegin;
810   qn->type = qtype;
811   PetscFunctionReturn(0);
812 }
813 
814 /* -------------------------------------------------------------------------- */
815 /*MC
816       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
817 
818       Options Database:
819 
820 +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
821 .     -snes_qn_powell_angle - Angle condition for restart.
822 .     -snes_qn_powell_descent - Descent condition for restart.
823 .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
824 -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
825 
826       Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
827       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
828       updates.
829 
830       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
831       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
832       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
833       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
834 
835       References:
836 
837       Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
838 
839       R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi-Newton Matrices and their use in Limited Memory Methods,
840       Technical Report, Northwestern University, June 1992.
841 
842       Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
843       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
844 
845 
846       Level: beginner
847 
848 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
849 
850 M*/
851 #undef __FUNCT__
852 #define __FUNCT__ "SNESCreate_QN"
853 PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
854 {
855   PetscErrorCode ierr;
856   SNES_QN        *qn;
857 
858   PetscFunctionBegin;
859   snes->ops->setup          = SNESSetUp_QN;
860   snes->ops->solve          = SNESSolve_QN;
861   snes->ops->destroy        = SNESDestroy_QN;
862   snes->ops->setfromoptions = SNESSetFromOptions_QN;
863   snes->ops->view           = SNESView_QN;
864   snes->ops->reset          = SNESReset_QN;
865 
866   snes->pcside = PC_LEFT;
867 
868   snes->usespc  = PETSC_TRUE;
869   snes->usesksp = PETSC_FALSE;
870 
871   if (!snes->tolerancesset) {
872     snes->max_funcs = 30000;
873     snes->max_its   = 10000;
874   }
875 
876   ierr                = PetscNewLog(snes,&qn);CHKERRQ(ierr);
877   snes->data          = (void*) qn;
878   qn->m               = 10;
879   qn->scaling         = 1.0;
880   qn->U               = NULL;
881   qn->V               = NULL;
882   qn->dXtdF           = NULL;
883   qn->dFtdX           = NULL;
884   qn->dXdFmat         = NULL;
885   qn->monitor         = NULL;
886   qn->singlereduction = PETSC_TRUE;
887   qn->powell_gamma    = 0.9999;
888   qn->scale_type      = SNES_QN_SCALE_DEFAULT;
889   qn->restart_type    = SNES_QN_RESTART_DEFAULT;
890   qn->type            = SNES_QN_LBFGS;
891 
892   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr);
893   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr);
894   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);CHKERRQ(ierr);
895   PetscFunctionReturn(0);
896 }
897