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