xref: /petsc/src/snes/impls/qn/qn.c (revision b5c23020607e2274b96d6554e5ed7d43c908fac4)
1 #include <private/snesimpl.h>
2 
3 
4 
5 typedef struct {
6   Vec         *dX;              /* The change in X */
7   Vec         *dF;              /* The change in F */
8   PetscInt    m;                /* the number of kept previous steps */
9   PetscScalar *alpha;
10   PetscScalar *beta;
11   PetscScalar *rho;
12   PetscViewer monitor;
13   PetscReal   powell_gamma;     /* Powell angle restart condition */
14   PetscReal   powell_downhill;  /* Powell descent restart condition */
15   PetscReal   scaling;          /* scaling of H0 */
16 } SNES_QN;
17 
18 #undef __FUNCT__
19 #define __FUNCT__ "LBGFSApplyJinv_Private"
20 PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) {
21 
22   PetscErrorCode ierr;
23 
24   SNES_QN *qn = (SNES_QN*)snes->data;
25 
26   Vec *dX = qn->dX;
27   Vec *dF = qn->dF;
28 
29   PetscScalar *alpha = qn->alpha;
30   PetscScalar *beta = qn->beta;
31   PetscScalar *rho = qn->rho;
32 
33   PetscInt k, i;
34   PetscInt m = qn->m;
35   PetscScalar t;
36   PetscInt l = m;
37 
38   PetscFunctionBegin;
39 
40   ierr = VecCopy(D, Y);CHKERRQ(ierr);
41 
42   if (it < m) l = it;
43 
44   /* outward recursion starting at iteration k's update and working back */
45   for (i = 0; i < l; i++) {
46     k = (it - i - 1) % l;
47     ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr);
48     alpha[k] = t*rho[k];
49     if (qn->monitor) {
50       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
51       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr);
52       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
53     }
54     ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr);
55   }
56 
57   ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr);
58 
59   /* inward recursion starting at the first update and working forward */
60   for (i = 0; i < l; i++) {
61     k = (it + i - l) % l;
62     ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr);
63     beta[k] = rho[k]*t;
64     ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
65     if (qn->monitor) {
66       ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
67       ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr);
68       ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
69     }
70   }
71   PetscFunctionReturn(0);
72 }
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "SNESSolve_QN"
76 static PetscErrorCode SNESSolve_QN(SNES snes)
77 {
78 
79   PetscErrorCode ierr;
80   SNES_QN *qn = (SNES_QN*) snes->data;
81 
82   Vec X, Xold;
83   Vec F, G, B;
84   Vec W, Y, FPC, Fold;
85   SNESConvergedReason reason;
86   PetscInt i, i_r, k;
87 
88   PetscReal fnorm, xnorm = 0, ynorm, gnorm;
89   PetscInt m = qn->m;
90   PetscBool lssucceed, changed;
91 
92   PetscScalar rhosc;
93 
94   Vec *dX = qn->dX;
95   Vec *dF = qn->dF;
96   PetscScalar *rho = qn->rho;
97   PetscScalar FolddotF, FolddotFold, FdotF, YdotF, a;
98 
99   /* basically just a regular newton's method except for the application of the jacobian */
100   PetscFunctionBegin;
101 
102   X             = snes->vec_sol;        /* solution vector */
103   F             = snes->vec_func;       /* residual vector */
104   Y             = snes->vec_sol_update; /* search direction generated by J^-1D*/
105   B             = snes->vec_rhs;
106   G             = snes->work[0];
107   W             = snes->work[1];
108   Xold          = snes->work[2];
109 
110   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
111   Fold          = snes->work[3];
112 
113   snes->reason = SNES_CONVERGED_ITERATING;
114 
115   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
116   snes->iter = 0;
117   snes->norm = 0.;
118   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
119   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
120   if (snes->domainerror) {
121     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
122     PetscFunctionReturn(0);
123   }
124   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
125   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
126   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
127   snes->norm = fnorm;
128   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
129   SNESLogConvHistory(snes,fnorm,0);
130   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
131 
132   /* set parameter for default relative tolerance convergence test */
133    snes->ttol = fnorm*snes->rtol;
134   /* test convergence */
135   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
136   if (snes->reason) PetscFunctionReturn(0);
137 
138   /* initialize the search direction as steepest descent */
139   ierr = VecCopy(F, Y);CHKERRQ(ierr);
140 
141   for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
142     /* line search for lambda */
143     ynorm = 1; gnorm = fnorm;
144     ierr = VecCopy(F, Fold);CHKERRQ(ierr);
145     ierr = VecCopy(X, Xold);CHKERRQ(ierr);
146     ierr = SNESLineSearchPreCheckApply(snes,X,Y,&changed);CHKERRQ(ierr);
147     ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
148     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
149     if (snes->domainerror) {
150       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
151       PetscFunctionReturn(0);
152       }
153     if (!lssucceed) {
154       if (++snes->numFailures >= snes->maxFailures) {
155         snes->reason = SNES_DIVERGED_LINE_SEARCH;
156         break;
157       }
158     }
159 
160     /* Update function and solution vectors */
161     fnorm = gnorm;
162     ierr = VecCopy(G,F);CHKERRQ(ierr);
163     ierr = VecCopy(W,X);CHKERRQ(ierr);
164 
165     if (snes->pc) {
166       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
167       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
168       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
169         snes->reason = SNES_DIVERGED_INNER;
170         PetscFunctionReturn(0);
171       }
172       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
173       ierr = VecCopy(FPC, F);CHKERRQ(ierr);
174       ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
175     }
176 
177     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);
178 
179     ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr);
180     ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr);
181 
182     SNESLogConvHistory(snes,snes->norm,snes->iter);
183     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
184     /* set parameter for default relative tolerance convergence test */
185     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
186     if (snes->reason) PetscFunctionReturn(0);
187 
188     /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
189     ierr = VecDot(Fold, Fold, &FolddotFold);CHKERRQ(ierr);
190     ierr = VecDot(Fold, F, &FolddotF);CHKERRQ(ierr);
191     ierr = VecDot(F, F, &FdotF);CHKERRQ(ierr);
192     ierr = VecDot(Y, F, &YdotF);CHKERRQ(ierr);
193     if (PetscAbs(PetscRealPart(FolddotF)) > qn->powell_gamma*PetscAbs(PetscRealPart(FolddotFold))) {
194       if (qn->monitor) {
195         ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
196         ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e|\n", k, PetscRealPart(FolddotF), qn->powell_gamma, PetscRealPart(FolddotFold));CHKERRQ(ierr);
197         ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr);
198       }
199       i_r = -1;
200       /* general purpose update */
201       if (snes->ops->update) {
202         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
203       }
204       ierr = VecCopy(F, Y);CHKERRQ(ierr);
205     } else {
206       /* set the differences */
207       k = i_r % m;
208       ierr = VecCopy(F, dF[k]);CHKERRQ(ierr);
209       ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr);
210       ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
211       ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
212       ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr);
213 
214       /* set scaling to be shanno scaling */
215       rho[k] = 1. / rhosc;
216       ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr);
217       qn->scaling = PetscRealPart(rhosc) / PetscRealPart(a);
218 
219       /* general purpose update */
220       if (snes->ops->update) {
221         ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
222       }
223       /* apply the current iteration of the approximate jacobian in order to get the next search direction*/
224       ierr = LBGFSApplyJinv_Private(snes, i_r+1, F, Y);CHKERRQ(ierr);
225     }
226 }
227   if (i == snes->max_its) {
228     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
229     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
230   }
231   PetscFunctionReturn(0);
232 }
233 
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "SNESSetUp_QN"
237 static PetscErrorCode SNESSetUp_QN(SNES snes)
238 {
239   SNES_QN *qn = (SNES_QN*)snes->data;
240   PetscErrorCode ierr;
241   PetscFunctionBegin;
242   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
243   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr);
244   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr);
245   ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr);
246   PetscFunctionReturn(0);
247 }
248 
249 #undef __FUNCT__
250 #define __FUNCT__ "SNESReset_QN"
251 static PetscErrorCode SNESReset_QN(SNES snes)
252 {
253   PetscErrorCode ierr;
254   SNES_QN *qn;
255   PetscFunctionBegin;
256   if (snes->data) {
257     qn = (SNES_QN*)snes->data;
258     if (qn->dX) {
259       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
260     }
261     if (qn->dF) {
262       ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr);
263     }
264     ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr);
265   }
266   PetscFunctionReturn(0);
267 }
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "SNESDestroy_QN"
271 static PetscErrorCode SNESDestroy_QN(SNES snes)
272 {
273   PetscErrorCode ierr;
274   PetscFunctionBegin;
275   ierr = SNESReset_QN(snes);CHKERRQ(ierr);
276   ierr = PetscFree(snes->data);CHKERRQ(ierr);
277   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr);
278   PetscFunctionReturn(0);
279 }
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "SNESSetFromOptions_QN"
283 static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
284 {
285 
286   PetscErrorCode ierr;
287   SNES_QN *qn;
288   PetscBool monflg = PETSC_FALSE;
289   PetscFunctionBegin;
290 
291   qn = (SNES_QN*)snes->data;
292 
293   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
294   ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
295   ierr = PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance",             "SNES", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr);
296   ierr = PetscOptionsReal("-snes_qn_powell_downhill", "Powell descent tolerance",        "SNES", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr);
297   ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods",              "SNES", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr);
298   ierr = PetscOptionsTail();CHKERRQ(ierr);
299   if (monflg) {
300     qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
301   }
302   PetscFunctionReturn(0);
303 }
304 
305 EXTERN_C_BEGIN
306 #undef __FUNCT__
307 #define __FUNCT__ "SNESLineSearchSetType_QN"
308 PetscErrorCode  SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type)
309 {
310   PetscErrorCode ierr;
311   PetscFunctionBegin;
312 
313   switch (type) {
314   case SNES_LS_BASIC:
315     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
316     break;
317   case SNES_LS_BASIC_NONORMS:
318     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
319     break;
320   case SNES_LS_QUADRATIC:
321     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
322     break;
323   case SNES_LS_SECANT:
324     ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr);
325     break;
326   default:
327     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
328     break;
329   }
330   snes->ls_type = type;
331   PetscFunctionReturn(0);
332 }
333 EXTERN_C_END
334 
335 
336 /* -------------------------------------------------------------------------- */
337 /*MC
338       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
339 
340       Options Database:
341 
342 +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
343 -     -snes_qn_monitor - Monitors the quasi-newton jacobian.
344 
345       Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to
346       form the approximate inverse Jacobian using a series of multiplicative rank-one updates.  This will eventually be
347       generalized to implement several limited-memory Broyden methods.
348 
349       References:
350 
351       L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed,
352       International Journal of Computer Mathematics, vol. 86, 2009.
353 
354 
355       Level: beginner
356 
357 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
358 
359 M*/
360 EXTERN_C_BEGIN
361 #undef __FUNCT__
362 #define __FUNCT__ "SNESCreate_QN"
363 PetscErrorCode  SNESCreate_QN(SNES snes)
364 {
365 
366   PetscErrorCode ierr;
367   SNES_QN *qn;
368 
369   PetscFunctionBegin;
370   snes->ops->setup           = SNESSetUp_QN;
371   snes->ops->solve           = SNESSolve_QN;
372   snes->ops->destroy         = SNESDestroy_QN;
373   snes->ops->setfromoptions  = SNESSetFromOptions_QN;
374   snes->ops->view            = 0;
375   snes->ops->reset           = SNESReset_QN;
376 
377   snes->usespc          = PETSC_TRUE;
378   snes->usesksp         = PETSC_FALSE;
379 
380   snes->max_funcs = 30000;
381   snes->max_its   = 10000;
382 
383   ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr);
384   snes->data = (void *) qn;
385   qn->m       = 10;
386   qn->scaling = 1.0;
387   qn->dX      = PETSC_NULL;
388   qn->dF      = PETSC_NULL;
389   qn->monitor = PETSC_NULL;
390   qn->powell_gamma   = 0.8;
391   qn->powell_downhill = 0.4;
392 
393   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr);
394   ierr = SNESLineSearchSetType(snes, SNES_LS_SECANT);CHKERRQ(ierr);
395 
396   PetscFunctionReturn(0);
397 }
398 EXTERN_C_END
399