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