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