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