xref: /petsc/src/snes/impls/qn/qn.c (revision 6950c3364dbdc26b53fdae8a802a5c652fef0a3b)
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;
302   Vec                 Y,D,Dold;
303   PetscInt            i, i_r;
304   PetscReal           fnorm,xnorm,ynorm,gnorm;
305   PetscBool           lssucceed,powell,periodic;
306   PetscScalar         DolddotD,DolddotDold,DdotD,YdotD;
307   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
308   SNESConvergedReason reason;
309 
310   /* basically just a regular newton's method except for the application of the jacobian */
311 
312   PetscFunctionBegin;
313   F = snes->vec_func;                   /* residual vector */
314   Y = snes->vec_sol_update;             /* search direction generated by J^-1D*/
315 
316   X    = snes->vec_sol;                 /* solution vector */
317   Xold = snes->work[0];
318 
319   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
320   D    = snes->work[1];
321   Dold = snes->work[2];
322 
323   snes->reason = SNES_CONVERGED_ITERATING;
324 
325   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
326   snes->iter = 0;
327   snes->norm = 0.;
328   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
329 
330   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
331     ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr);
332     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
333     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
334       snes->reason = SNES_DIVERGED_INNER;
335       PetscFunctionReturn(0);
336     }
337     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
338   } else {
339     if (!snes->vec_func_init_set) {
340       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
341       if (snes->domainerror) {
342         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
343         PetscFunctionReturn(0);
344       }
345     } else snes->vec_func_init_set = PETSC_FALSE;
346     if (!snes->norm_init_set) {
347       /* convergence test */
348       ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
349       if (PetscIsInfOrNanReal(fnorm)) {
350         snes->reason = SNES_DIVERGED_FNORM_NAN;
351         PetscFunctionReturn(0);
352       }
353     } else {
354       fnorm               = snes->norm_init;
355       snes->norm_init_set = PETSC_FALSE;
356     }
357   }
358   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
359       ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr);
360       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
361       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
362         snes->reason = SNES_DIVERGED_INNER;
363         PetscFunctionReturn(0);
364       }
365   } else {
366     ierr = VecCopy(F,D);CHKERRQ(ierr);
367   }
368 
369   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
370   snes->norm = fnorm;
371   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
372   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
373   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
374 
375   /* set parameter for default relative tolerance convergence test */
376   snes->ttol = fnorm*snes->rtol;
377   /* test convergence */
378   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
379   if (snes->reason) PetscFunctionReturn(0);
380 
381   if (snes->pc && snes->pcside == PC_RIGHT) {
382     Vec FPC;
383     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
384     ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
385     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
386     ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr);
387     ierr = SNESGetFunctionNorm(snes->pc,&fnorm);CHKERRQ(ierr);
388     ierr = VecCopy(FPC,F);CHKERRQ(ierr);
389     ierr = VecCopy(F,D);CHKERRQ(ierr);
390   }
391 
392   /* scale the initial update */
393   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
394     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
395     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
396   }
397 
398   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
399     switch (qn->type) {
400     case SNES_QN_BADBROYDEN:
401       ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
402       break;
403     case SNES_QN_BROYDEN:
404       ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
405       break;
406     case SNES_QN_LBFGS:
407       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr);
408       break;
409     }
410     /* line search for lambda */
411     ynorm = 1; gnorm = fnorm;
412     ierr  = VecCopy(D, Dold);CHKERRQ(ierr);
413     ierr  = VecCopy(X, Xold);CHKERRQ(ierr);
414     ierr  = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
415     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
416     if (snes->domainerror) {
417       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
418       PetscFunctionReturn(0);
419     }
420     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr);
421     if (!lssucceed) {
422       if (++snes->numFailures >= snes->maxFailures) {
423         snes->reason = SNES_DIVERGED_LINE_SEARCH;
424         break;
425       }
426     }
427     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
428     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
429       ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr);
430     }
431 
432     /* convergence monitoring */
433     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);
434 
435     if (snes->pc && snes->pcside == PC_RIGHT) {
436       Vec FPC;
437       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
438       ierr = SNESSolve(snes->pc,snes->vec_rhs,X);CHKERRQ(ierr);
439       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);CHKERRQ(ierr);
440       ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr);
441       ierr = SNESGetFunctionNorm(snes->pc,&fnorm);CHKERRQ(ierr);
442       ierr = VecCopy(FPC,F);CHKERRQ(ierr);
443     }
444 
445     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
446     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
447 
448     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
449     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
450     /* set parameter for default relative tolerance convergence test */
451     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
452     if (snes->reason) PetscFunctionReturn(0);
453     if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
454       ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr);
455       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
456       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
457         snes->reason = SNES_DIVERGED_INNER;
458         PetscFunctionReturn(0);
459       }
460     } else {
461       ierr = VecCopy(F, D);CHKERRQ(ierr);
462     }
463 
464     powell = PETSC_FALSE;
465     if (qn->restart_type == SNES_QN_RESTART_POWELL) {
466       /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
467       ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
468       ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr);
469       ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr);
470       ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr);
471       ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr);
472       ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr);
473       ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr);
474       ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr);
475       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
476     }
477     periodic = PETSC_FALSE;
478     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
479       if (i_r>qn->m-1) periodic = PETSC_TRUE;
480     }
481     /* restart if either powell or periodic restart is satisfied. */
482     if (powell || periodic) {
483       if (qn->monitor) {
484         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
485         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);
486         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
487       }
488       i_r = -1;
489       /* general purpose update */
490       if (snes->ops->update) {
491         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
492       }
493       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
494         ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
495         ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
496       }
497     }
498     /* general purpose update */
499     if (snes->ops->update) {
500       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
501     }
502   }
503   if (i == snes->max_its) {
504     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
505     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
506   }
507   PetscFunctionReturn(0);
508 }
509 
510 #undef __FUNCT__
511 #define __FUNCT__ "SNESSetUp_QN"
512 static PetscErrorCode SNESSetUp_QN(SNES snes)
513 {
514   SNES_QN        *qn = (SNES_QN*)snes->data;
515   PetscErrorCode ierr;
516 
517   PetscFunctionBegin;
518   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr);
519   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr);
520   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha,
521                       qn->m, PetscScalar, &qn->beta,
522                       qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr);
523 
524   if (qn->singlereduction) {
525     ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat,
526                         qn->m, PetscScalar, &qn->dFtdX,
527                         qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr);
528   }
529   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
530 
531   /* set up the line search */
532   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
533     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
534   }
535   PetscFunctionReturn(0);
536 }
537 
538 #undef __FUNCT__
539 #define __FUNCT__ "SNESReset_QN"
540 static PetscErrorCode SNESReset_QN(SNES snes)
541 {
542   PetscErrorCode ierr;
543   SNES_QN        *qn;
544 
545   PetscFunctionBegin;
546   if (snes->data) {
547     qn = (SNES_QN*)snes->data;
548     if (qn->U) {
549       ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr);
550     }
551     if (qn->V) {
552       ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr);
553     }
554     if (qn->singlereduction) {
555       ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr);
556     }
557     ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr);
558   }
559   PetscFunctionReturn(0);
560 }
561 
562 #undef __FUNCT__
563 #define __FUNCT__ "SNESDestroy_QN"
564 static PetscErrorCode SNESDestroy_QN(SNES snes)
565 {
566   PetscErrorCode ierr;
567 
568   PetscFunctionBegin;
569   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
570   ierr = PetscFree(snes->data);CHKERRQ(ierr);
571   ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr);
572   PetscFunctionReturn(0);
573 }
574 
575 #undef __FUNCT__
576 #define __FUNCT__ "SNESSetFromOptions_QN"
577 static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
578 {
579 
580   PetscErrorCode    ierr;
581   SNES_QN           *qn    = (SNES_QN*)snes->data;
582   PetscBool         monflg = PETSC_FALSE,flg;
583   SNESLineSearch    linesearch;
584   SNESQNRestartType rtype = qn->restart_type;
585   SNESQNScaleType   stype = qn->scale_type;
586 
587   PetscFunctionBegin;
588   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
589   ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr);
590   ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr);
591   ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr);
592   ierr = PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr);
593   ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr);
594   ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr);
595   if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr);
596 
597   ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr);
598   if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr);
599 
600   ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,
601                           (PetscEnum)qn->type,(PetscEnum*)&qn->type,NULL);CHKERRQ(ierr);
602   ierr = PetscOptionsTail();CHKERRQ(ierr);
603   if (!snes->linesearch) {
604     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
605     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
606   }
607   if (monflg) {
608     qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
609   }
610   PetscFunctionReturn(0);
611 }
612 
613 
614 #undef __FUNCT__
615 #define __FUNCT__ "SNESQNSetRestartType"
616 /*@
617     SNESQNSetRestartType - Sets the restart type for SNESQN.
618 
619     Logically Collective on SNES
620 
621     Input Parameters:
622 +   snes - the iterative context
623 -   rtype - restart type
624 
625     Options Database:
626 +   -snes_qn_restart_type<powell,periodic,none> - set the restart type
627 -   -snes_qn_restart[10] - sets the number of iterations before restart for periodic
628 
629     Level: intermediate
630 
631     SNESQNRestartTypes:
632 +   SNES_QN_RESTART_NONE - never restart
633 .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
634 -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
635 
636     Notes:
637     The default line search used is the L2 line search and it requires two additional function evaluations.
638 
639 .keywords: SNES, SNESQN, restart, type, set SNESLineSearch
640 @*/
641 PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
642 {
643   PetscErrorCode ierr;
644 
645   PetscFunctionBegin;
646   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
647   ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr);
648   PetscFunctionReturn(0);
649 }
650 
651 #undef __FUNCT__
652 #define __FUNCT__ "SNESQNSetScaleType"
653 /*@
654     SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN.
655 
656     Logically Collective on SNES
657 
658     Input Parameters:
659 +   snes - the iterative context
660 -   stype - scale type
661 
662     Options Database:
663 .   -snes_qn_scale_type<shanno,none,linesearch,jacobian>
664 
665     Level: intermediate
666 
667     SNESQNSelectTypes:
668 +   SNES_QN_SCALE_NONE - don't scale the problem
669 .   SNES_QN_SCALE_SHANNO - use shanno scaling
670 .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
671 -   SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian.
672 
673 .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch
674 @*/
675 
676 PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
677 {
678   PetscErrorCode ierr;
679 
680   PetscFunctionBegin;
681   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
682   ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr);
683   PetscFunctionReturn(0);
684 }
685 
686 #undef __FUNCT__
687 #define __FUNCT__ "SNESQNSetScaleType_QN"
688 PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
689 {
690   SNES_QN *qn = (SNES_QN*)snes->data;
691 
692   PetscFunctionBegin;
693   qn->scale_type = stype;
694   PetscFunctionReturn(0);
695 }
696 
697 #undef __FUNCT__
698 #define __FUNCT__ "SNESQNSetRestartType_QN"
699 PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
700 {
701   SNES_QN *qn = (SNES_QN*)snes->data;
702 
703   PetscFunctionBegin;
704   qn->restart_type = rtype;
705   PetscFunctionReturn(0);
706 }
707 
708 /* -------------------------------------------------------------------------- */
709 /*MC
710       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
711 
712       Options Database:
713 
714 +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
715 .     -snes_qn_powell_angle - Angle condition for restart.
716 .     -snes_qn_powell_descent - Descent condition for restart.
717 .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
718 -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
719 
720       Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
721       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
722       updates.
723 
724       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
725       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
726       iteration as the current iteration's values when constructing the approximate jacobian.  The second, composed,
727       perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
728 
729       References:
730 
731       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
732       International Journal of Computer Mathematics, vol. 86, 2009.
733 
734       Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker
735       SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
736 
737 
738       Level: beginner
739 
740 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
741 
742 M*/
743 #undef __FUNCT__
744 #define __FUNCT__ "SNESCreate_QN"
745 PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
746 {
747   PetscErrorCode ierr;
748   SNES_QN        *qn;
749 
750   PetscFunctionBegin;
751   snes->ops->setup          = SNESSetUp_QN;
752   snes->ops->solve          = SNESSolve_QN;
753   snes->ops->destroy        = SNESDestroy_QN;
754   snes->ops->setfromoptions = SNESSetFromOptions_QN;
755   snes->ops->view           = 0;
756   snes->ops->reset          = SNESReset_QN;
757 
758   snes->pcside = PC_LEFT;
759   snes->functype = SNES_FUNCTION_PRECONDITIONED;
760 
761   snes->usespc  = PETSC_TRUE;
762   snes->usesksp = PETSC_FALSE;
763 
764   if (!snes->tolerancesset) {
765     snes->max_funcs = 30000;
766     snes->max_its   = 10000;
767   }
768 
769   ierr                = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
770   snes->data          = (void*) qn;
771   qn->m               = 10;
772   qn->scaling         = 1.0;
773   qn->U               = NULL;
774   qn->V               = NULL;
775   qn->dXtdF           = NULL;
776   qn->dFtdX           = NULL;
777   qn->dXdFmat         = NULL;
778   qn->monitor         = NULL;
779   qn->singlereduction = PETSC_FALSE;
780   qn->powell_gamma    = 0.9999;
781   qn->scale_type      = SNES_QN_SCALE_SHANNO;
782   qn->restart_type    = SNES_QN_RESTART_POWELL;
783   qn->type            = SNES_QN_LBFGS;
784 
785   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr);
786   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr);
787   PetscFunctionReturn(0);
788 }
789