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