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