xref: /petsc/src/snes/impls/qn/qn.c (revision 7cf18bfc0af3bd20f57ca1ab9dbfc3d774cdee5d)
1 #include <private/snesimpl.h>
2 
3 typedef struct {
4   Vec * dX;           /* The change in X */
5   Vec * dF;           /* The change in F */
6   PetscInt m;         /* the number of kept previous steps */
7   PetscScalar * alpha;
8   PetscScalar * beta;
9   PetscScalar * rho;
10 } QNContext;
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "LBGFSApplyJinv_Private"
14 PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec g, Vec z) {
15 
16   PetscErrorCode ierr;
17 
18   QNContext * qn = (QNContext *)snes->data;
19 
20   Vec * dX = qn->dX;
21   Vec * dF = qn->dF;
22 
23   PetscScalar * alpha = qn->alpha;
24   PetscScalar * beta = qn->beta;
25   PetscScalar * rho = qn->rho;
26 
27   PetscInt k, i;
28   PetscInt m = qn->m;
29   PetscScalar t;
30   PetscInt l = m;
31 
32   PetscFunctionBegin;
33 
34   if (it < m) l = it;
35 
36   ierr = VecCopy(g, z);CHKERRQ(ierr);
37 
38   /* outward recursion starting at iteration k's update and working back */
39   for (i = 0; i < l; i++) {
40     k = (it - i - 1) % l;
41     ierr = VecDot(dX[k], z, &t);CHKERRQ(ierr);
42     alpha[k] = t*rho[k];
43     /*
44     ierr = PetscPrintf(PETSC_COMM_WORLD, "%d, alpha %g\n", k, alpha[k]);CHKERRQ(ierr);
45      */
46     ierr = VecAXPY(z, -alpha[k], dF[k]);CHKERRQ(ierr);
47   }
48 
49   /* inner application of the initial inverse jacobian approximation */
50   /* right now it's just the identity. Nothing needs to go here. */
51 
52   /* inward recursion starting at the first update and working forward*/
53   for (i = 0; i < l; i++) {
54     k = (it + i - l) % l;
55     ierr = VecDot(dF[k], z, &t);CHKERRQ(ierr);
56     beta[k] = rho[k]*t;
57     ierr = VecAXPY(z, (alpha[k] - beta[k]), dX[k]);
58     /*
59     ierr = PetscPrintf(PETSC_COMM_WORLD, "%d, alpha - beta %g\n", k, alpha[k] - beta[k]);CHKERRQ(ierr);
60      */
61   }
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "SNESLineSearchQuadratic_QN"
67 PetscErrorCode  SNESLineSearchQuadratic_QN(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal xnorm,Vec G,Vec W,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
68 {
69   PetscInt         i;
70   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
71   PetscReal        norms[3];
72   PetscReal        alpha,a,b;
73   PetscErrorCode   ierr;
74 
75   PetscFunctionBegin;
76   norms[0]  = fnorm;
77   for(i=1; i < 3; ++i) {
78     ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
79     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
80     ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr);
81   }
82   for(i = 0; i < 3; ++i) {
83     norms[i] = PetscSqr(norms[i]);
84   }
85   /* Fit a quadratic:
86        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
87        a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1)
88        b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1)
89        c = y_0
90        x_min = -b/2a
91 
92        If we let x_{0,1,2} = 0, 0.5, 1.0
93        a = 2 y_2 - 4 y_1 + 2 y_0
94        b =  -y_2 + 4 y_1 - 3 y_0
95        c =   y_0
96   */
97   a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1]));
98   b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]);
99   /* Check for positive a (concave up) */
100   if (a >= 0.0) {
101     alpha = -b/(2.0*a);
102     alpha = PetscMin(alpha, alphas[2]);
103     alpha = PetscMax(alpha, alphas[0]);
104   } else {
105     alpha = 1.0;
106   }
107   if (snes->ls_monitor) {
108     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
109     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr);
110     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
111   }
112   ierr = VecCopy(X, W);CHKERRQ(ierr);
113   ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr);
114   if (alpha != 1.0) {
115     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
116     ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr);
117   } else {
118     *gnorm = PetscSqrtReal(norms[2]);
119   }
120   if (alpha == 0.0) *flag = PETSC_FALSE;
121   else              *flag = PETSC_TRUE;
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "SNESSolve_QN"
127 static PetscErrorCode SNESSolve_QN(SNES snes)
128 {
129 
130   PetscErrorCode ierr;
131   QNContext * qn = (QNContext*) snes->data;
132 
133   Vec X, Xold;
134   Vec F, Fold, G;
135   Vec W, Y;
136 
137   PetscInt i, k;
138 
139   PetscReal fnorm, xnorm = 0, ynorm, gnorm;
140   PetscInt m = qn->m;
141   PetscBool lssucceed;
142 
143   PetscScalar rhosc;
144 
145   Vec * dX = qn->dX;
146   Vec * dF = qn->dF;
147   PetscScalar * rho = qn->rho;
148 
149   /* basically just a regular newton's method except for the application of the jacobian */
150   PetscFunctionBegin;
151 
152   X             = snes->vec_sol;        /* solution vector */
153   F             = snes->vec_func;       /* residual vector */
154   Y             = snes->vec_sol_update; /* search direction */
155   G             = snes->work[0];
156   W             = snes->work[1];
157   Xold          = snes->work[2];
158   Fold          = snes->work[3];
159 
160   snes->reason = SNES_CONVERGED_ITERATING;
161 
162   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
163   snes->iter = 0;
164   snes->norm = 0.;
165   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
166   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
167   if (snes->domainerror) {
168     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
169     PetscFunctionReturn(0);
170   }
171   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
172   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
173   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
174   snes->norm = fnorm;
175   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
176   SNESLogConvHistory(snes,fnorm,0);
177   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
178 
179   /* set parameter for default relative tolerance convergence test */
180    snes->ttol = fnorm*snes->rtol;
181   /* test convergence */
182   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
183   if (snes->reason) PetscFunctionReturn(0);
184 
185   /* initialize the search direction as steepest descent */
186   ierr = VecCopy(F, Y);CHKERRQ(ierr);
187  for(i = 0; i < snes->max_its; i++) {
188    /* line search for lambda */
189    ynorm = 1; gnorm = fnorm;
190    ierr = VecCopy(F, Fold);CHKERRQ(ierr);
191    ierr = VecCopy(X, Xold);CHKERRQ(ierr);
192    ierr = VecScale(Y, -1.0);CHKERRQ(ierr);
193    ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
194    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);
195    if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
196    if (snes->domainerror) {
197      snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
198      PetscFunctionReturn(0);
199    }
200    if (!lssucceed) {
201      if (++snes->numFailures >= snes->maxFailures) {
202        snes->reason = SNES_DIVERGED_LINE_SEARCH;
203        break;
204      }
205    }
206    /* Update function and solution vectors */
207    fnorm = gnorm;
208    ierr = VecCopy(G,F);CHKERRQ(ierr);
209    ierr = VecCopy(W,X);CHKERRQ(ierr);
210 
211    ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
212    snes->iter = i+1;
213    snes->norm = fnorm;
214    ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
215    SNESLogConvHistory(snes,snes->norm,snes->iter);
216    ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
217    /* set parameter for default relative tolerance convergence test */
218    ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
219    if (snes->reason) PetscFunctionReturn(0);
220 
221    /* set the differences */
222    k = i % m;
223    ierr = VecCopy(F, dF[k]);CHKERRQ(ierr);
224    ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr);
225    ierr = VecCopy(X, dX[k]);CHKERRQ(ierr);
226    ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr);
227    ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr);
228    rho[k] = 1. / rhosc;
229 
230    /* general purpose update */
231    if (snes->ops->update) {
232      ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
233    }
234 
235    /* apply the current iteration of the approximate jacobian in order to get the next search direction*/
236    ierr = LBGFSApplyJinv_Private(snes, i, F, Y);CHKERRQ(ierr);
237  }
238   if (i == snes->max_its) {
239     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr);
240     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
241   }
242   PetscFunctionReturn(0);
243 }
244 
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "SNESSetUp_QN"
248 static PetscErrorCode SNESSetUp_QN(SNES snes)
249 {
250   QNContext * qn = (QNContext *)snes->data;
251   PetscErrorCode ierr;
252   PetscFunctionBegin;
253   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr);
254   ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr);
255   ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr);
256   ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr);
257   PetscFunctionReturn(0);
258 }
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "SNESReset_QN"
262 static PetscErrorCode SNESReset_QN(SNES snes)
263 {
264   PetscErrorCode ierr;
265   QNContext * qn;
266   PetscFunctionBegin;
267   if (snes->data) {
268     qn = (QNContext *)snes->data;
269     if (qn->dX) {
270       ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr);
271     }
272     if (qn->dF) {
273       ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr);
274     }
275     ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr);
276   }
277   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
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   PetscFunctionReturn(0);
290 }
291 
292 #undef __FUNCT__
293 #define __FUNCT__ "SNESSetFromOptions_QN"
294 static PetscErrorCode SNESSetFromOptions_QN(SNES snes)
295 {
296 
297   PetscErrorCode ierr;
298   QNContext * qn;
299 
300   PetscFunctionBegin;
301 
302   qn = (QNContext *)snes->data;
303 
304   ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr);
305   ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr);
306   ierr = PetscOptionsTail();CHKERRQ(ierr);
307   PetscFunctionReturn(0);
308 }
309 
310 EXTERN_C_BEGIN
311 #undef __FUNCT__
312 #define __FUNCT__ "SNESLineSearchSetType_QN"
313 PetscErrorCode  SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type)
314 {
315   PetscErrorCode ierr;
316   PetscFunctionBegin;
317 
318   switch (type) {
319   case SNES_LS_BASIC:
320     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
321     break;
322   case SNES_LS_BASIC_NONORMS:
323     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
324     break;
325   case SNES_LS_QUADRATIC:
326     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_QN,PETSC_NULL);CHKERRQ(ierr);
327     break;
328   default:
329     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
330     break;
331   }
332   snes->ls_type = type;
333   PetscFunctionReturn(0);
334 }
335 EXTERN_C_END
336 
337 
338 /* -------------------------------------------------------------------------- */
339 /*MC
340       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
341 
342       Options Database:
343 
344 +     -snes_qn_m - Number of past states saved for the L-Broyden methods.
345 
346       Notes: This implements the L-BFGS algorithm for the solution of F(x) = 0 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->dF = PETSC_NULL;
386 
387   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr);
388   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
389 
390   PetscFunctionReturn(0);
391 }
392 EXTERN_C_END
393