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