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