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