xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 40cbb1a031ea8f2be4fe2b92dc842b003ad37be3)
1 #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/
2 const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",NULL};
3 
4 static PetscErrorCode SNESReset_NCG(SNES snes)
5 {
6   PetscFunctionBegin;
7   PetscFunctionReturn(0);
8 }
9 
10 /*
11   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
12 
13   Input Parameter:
14 . snes - the SNES context
15 
16   Application Interface Routine: SNESDestroy()
17 */
18 static PetscErrorCode SNESDestroy_NCG(SNES snes)
19 {
20   PetscFunctionBegin;
21   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C",NULL));
22   PetscCall(PetscFree(snes->data));
23   PetscFunctionReturn(0);
24 }
25 
26 /*
27    SNESSetUp_NCG - Sets up the internal data structures for the later use
28    of the SNESNCG nonlinear solver.
29 
30    Input Parameters:
31 +  snes - the SNES context
32 -  x - the solution vector
33 
34    Application Interface Routine: SNESSetUp()
35  */
36 
37 static PetscErrorCode SNESSetUp_NCG(SNES snes)
38 {
39   PetscFunctionBegin;
40   PetscCall(SNESSetWorkVecs(snes,2));
41   PetscCheck(snes->npcside!= PC_RIGHT,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
42   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
43   PetscFunctionReturn(0);
44 }
45 
46 static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
47 {
48   PetscScalar    alpha, ptAp;
49   Vec            X, Y, F, W;
50   SNES           snes;
51   PetscReal      *fnorm, *xnorm, *ynorm;
52 
53   PetscFunctionBegin;
54   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
55   X     = linesearch->vec_sol;
56   W     = linesearch->vec_sol_new;
57   F     = linesearch->vec_func;
58   Y     = linesearch->vec_update;
59   fnorm = &linesearch->fnorm;
60   xnorm = &linesearch->xnorm;
61   ynorm = &linesearch->ynorm;
62 
63   if (!snes->jacobian) {
64     PetscCall(SNESSetUpMatrices(snes));
65   }
66 
67   /*
68 
69    The exact step size for unpreconditioned linear CG is just:
70    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
71    */
72   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
73   PetscCall(VecDot(F, F, &alpha));
74   PetscCall(MatMult(snes->jacobian, Y, W));
75   PetscCall(VecDot(Y, W, &ptAp));
76   alpha = alpha / ptAp;
77   PetscCall(VecAXPY(X,-alpha,Y));
78   PetscCall(SNESComputeFunction(snes,X,F));
79 
80   PetscCall(VecNorm(F,NORM_2,fnorm));
81   PetscCall(VecNorm(X,NORM_2,xnorm));
82   PetscCall(VecNorm(Y,NORM_2,ynorm));
83   PetscFunctionReturn(0);
84 }
85 
86 /*MC
87    SNESLINESEARCHNCGLINEAR - Special line search only for SNESNCG
88 
89    This line search uses the length "as if" the problem is linear (that is what is computed by the linear CG method) using the Jacobian of the function.
90    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) where r (f) is the current residual (function value), p (y) is the current search direction.
91 
92    Notes: This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
93 
94    This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
95 
96    Level: advanced
97 
98 .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
99 M*/
100 
101 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
102 {
103   PetscFunctionBegin;
104   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
105   linesearch->ops->destroy        = NULL;
106   linesearch->ops->setfromoptions = NULL;
107   linesearch->ops->reset          = NULL;
108   linesearch->ops->view           = NULL;
109   linesearch->ops->setup          = NULL;
110   PetscFunctionReturn(0);
111 }
112 
113 /*
114   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
115 
116   Input Parameter:
117 . snes - the SNES context
118 
119   Application Interface Routine: SNESSetFromOptions()
120 */
121 static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes)
122 {
123   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
124   PetscBool      debug = PETSC_FALSE;
125   SNESNCGType    ncgtype=ncg->type;
126   SNESLineSearch linesearch;
127 
128   PetscFunctionBegin;
129   PetscOptionsHeadBegin(PetscOptionsObject,"SNES NCG options");
130   PetscCall(PetscOptionsBool("-snes_ncg_monitor","Monitor the beta values used in the NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
131   if (debug) {
132     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
133   }
134   PetscCall(PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL));
135   PetscCall(SNESNCGSetType(snes, ncgtype));
136   PetscOptionsHeadEnd();
137   if (!snes->linesearch) {
138     PetscCall(SNESGetLineSearch(snes, &linesearch));
139     if (!((PetscObject)linesearch)->type_name) {
140       if (!snes->npc) {
141         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
142       } else {
143         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
144       }
145     }
146   }
147   PetscFunctionReturn(0);
148 }
149 
150 /*
151   SNESView_NCG - Prints info from the SNESNCG data structure.
152 
153   Input Parameters:
154 + SNES - the SNES context
155 - viewer - visualization context
156 
157   Application Interface Routine: SNESView()
158 */
159 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
160 {
161   SNES_NCG      *ncg = (SNES_NCG *) snes->data;
162   PetscBool      iascii;
163 
164   PetscFunctionBegin;
165   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
166   if (iascii) {
167     PetscCall(PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]));
168   }
169   PetscFunctionReturn(0);
170 }
171 
172 /*
173 
174  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
175 
176  This routine is not currently used. I don't know what its intended purpose is.
177 
178  Note it has a hardwired differencing parameter of 1e-5
179 
180  */
181 PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
182 {
183   PetscScalar    ftf, ftg, fty, h;
184 
185   PetscFunctionBegin;
186   PetscCall(VecDot(F, F, &ftf));
187   PetscCall(VecDot(F, Y, &fty));
188   h      = 1e-5*fty / fty;
189   PetscCall(VecCopy(X, W));
190   PetscCall(VecAXPY(W, -h, Y));          /* this is arbitrary */
191   PetscCall(SNESComputeFunction(snes, W, G));
192   PetscCall(VecDot(G, F, &ftg));
193   *ytJtf = (ftg - ftf) / h;
194   PetscFunctionReturn(0);
195 }
196 
197 /*@
198     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
199 
200     Logically Collective on SNES
201 
202     Input Parameters:
203 +   snes - the iterative context
204 -   btype - update type
205 
206     Options Database:
207 .   -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta
208 
209     Level: intermediate
210 
211     SNESNCGSelectTypes:
212 +   SNES_NCG_FR - Fletcher-Reeves update
213 .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
214 .   SNES_NCG_HS - Hestenes-Steifel update
215 .   SNES_NCG_DY - Dai-Yuan update
216 -   SNES_NCG_CD - Conjugate Descent update
217 
218    Notes:
219    SNES_NCG_PRP is the default, and the only one that tolerates generalized search directions.
220 
221    It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
222    that is using -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC()?
223 
224 @*/
225 PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
226 {
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
229   PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));
230   PetscFunctionReturn(0);
231 }
232 
233 static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
234 {
235   SNES_NCG *ncg = (SNES_NCG*)snes->data;
236 
237   PetscFunctionBegin;
238   ncg->type = btype;
239   PetscFunctionReturn(0);
240 }
241 
242 /*
243   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
244 
245   Input Parameters:
246 . snes - the SNES context
247 
248   Output Parameter:
249 . outits - number of iterations until termination
250 
251   Application Interface Routine: SNESSolve()
252 */
253 static PetscErrorCode SNESSolve_NCG(SNES snes)
254 {
255   SNES_NCG             *ncg = (SNES_NCG*)snes->data;
256   Vec                  X,dX,lX,F,dXold;
257   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
258   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
259   PetscInt             maxits, i;
260   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
261   SNESLineSearch       linesearch;
262   SNESConvergedReason  reason;
263 
264   PetscFunctionBegin;
265   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
266 
267   PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite));
268   snes->reason = SNES_CONVERGED_ITERATING;
269 
270   maxits = snes->max_its;            /* maximum number of iterations */
271   X      = snes->vec_sol;            /* X^n */
272   dXold  = snes->work[0];            /* The previous iterate of X */
273   dX     = snes->work[1];            /* the preconditioned direction */
274   lX     = snes->vec_sol_update;     /* the conjugate direction */
275   F      = snes->vec_func;           /* residual vector */
276 
277   PetscCall(SNESGetLineSearch(snes, &linesearch));
278 
279   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
280   snes->iter = 0;
281   snes->norm = 0.;
282   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
283 
284   /* compute the initial function and preconditioned update dX */
285 
286   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
287     PetscCall(SNESApplyNPC(snes,X,NULL,dX));
288     PetscCall(SNESGetConvergedReason(snes->npc,&reason));
289     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
290       snes->reason = SNES_DIVERGED_INNER;
291       PetscFunctionReturn(0);
292     }
293     PetscCall(VecCopy(dX,F));
294     PetscCall(VecNorm(F,NORM_2,&fnorm));
295   } else {
296     if (!snes->vec_func_init_set) {
297       PetscCall(SNESComputeFunction(snes,X,F));
298     } else snes->vec_func_init_set = PETSC_FALSE;
299 
300     /* convergence test */
301     PetscCall(VecNorm(F,NORM_2,&fnorm));
302     SNESCheckFunctionNorm(snes,fnorm);
303     PetscCall(VecCopy(F,dX));
304   }
305   if (snes->npc) {
306     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
307       PetscCall(SNESApplyNPC(snes,X,F,dX));
308       PetscCall(SNESGetConvergedReason(snes->npc,&reason));
309       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
310         snes->reason = SNES_DIVERGED_INNER;
311         PetscFunctionReturn(0);
312       }
313     }
314   }
315   PetscCall(VecCopy(dX,lX));
316   PetscCall(VecDot(dX, dX, &dXdotdX));
317 
318   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
319   snes->norm = fnorm;
320   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
321   PetscCall(SNESLogConvergenceHistory(snes,fnorm,0));
322   PetscCall(SNESMonitor(snes,0,fnorm));
323 
324   /* test convergence */
325   PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP));
326   if (snes->reason) PetscFunctionReturn(0);
327 
328   /* Call general purpose update function */
329   if (snes->ops->update) {
330     PetscCall((*snes->ops->update)(snes, snes->iter));
331   }
332 
333   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
334 
335   for (i = 1; i < maxits + 1; i++) {
336     /* some update types require the old update direction or conjugate direction */
337     if (ncg->type != SNES_NCG_FR) {
338       PetscCall(VecCopy(dX, dXold));
339     }
340     PetscCall(SNESLineSearchApply(linesearch,X,F,&fnorm,lX));
341     PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
342     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
343     if (lsresult) {
344       if (++snes->numFailures >= snes->maxFailures) {
345         snes->reason = SNES_DIVERGED_LINE_SEARCH;
346         PetscFunctionReturn(0);
347       }
348     }
349     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
350       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
351       PetscFunctionReturn(0);
352     }
353     /* Monitor convergence */
354     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
355     snes->iter = i;
356     snes->norm = fnorm;
357     snes->xnorm = xnorm;
358     snes->ynorm = ynorm;
359     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
360     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
361     PetscCall(SNESMonitor(snes,snes->iter,snes->norm));
362 
363     /* Test for convergence */
364     PetscCall((*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP));
365     if (snes->reason) PetscFunctionReturn(0);
366 
367     /* Call general purpose update function */
368     if (snes->ops->update) {
369       PetscCall((*snes->ops->update)(snes, snes->iter));
370     }
371     if (snes->npc) {
372       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
373         PetscCall(SNESApplyNPC(snes,X,NULL,dX));
374         PetscCall(SNESGetConvergedReason(snes->npc,&reason));
375         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
376           snes->reason = SNES_DIVERGED_INNER;
377           PetscFunctionReturn(0);
378         }
379         PetscCall(VecCopy(dX,F));
380       } else {
381         PetscCall(SNESApplyNPC(snes,X,F,dX));
382         PetscCall(SNESGetConvergedReason(snes->npc,&reason));
383         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
384           snes->reason = SNES_DIVERGED_INNER;
385           PetscFunctionReturn(0);
386         }
387       }
388     } else {
389       PetscCall(VecCopy(F,dX));
390     }
391 
392     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
393     switch (ncg->type) {
394     case SNES_NCG_FR: /* Fletcher-Reeves */
395       dXolddotdXold= dXdotdX;
396       PetscCall(VecDot(dX, dX, &dXdotdX));
397       beta         = PetscRealPart(dXdotdX / dXolddotdXold);
398       break;
399     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
400       dXolddotdXold= dXdotdX;
401       PetscCall(VecDotBegin(dX, dX, &dXdotdXold));
402       PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
403       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
404       PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
405       beta         = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
406       if (beta < 0.0) beta = 0.0; /* restart */
407       break;
408     case SNES_NCG_HS: /* Hestenes-Stiefel */
409       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
410       PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
411       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
412       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
413       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
414       PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
415       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
416       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
417       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
418       break;
419     case SNES_NCG_DY: /* Dai-Yuan */
420       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
421       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
422       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
423       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
424       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
425       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
426       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
427       break;
428     case SNES_NCG_CD: /* Conjugate Descent */
429       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
430       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
431       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
432       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
433       beta = PetscRealPart(dXdotdX / lXdotdXold);
434       break;
435     }
436     if (ncg->monitor) {
437       PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
438     }
439     PetscCall(VecAYPX(lX, beta, dX));
440   }
441   PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
442   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
443   PetscFunctionReturn(0);
444 }
445 
446 /*MC
447   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
448 
449   Level: beginner
450 
451   Options Database:
452 +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
453 .   -snes_linesearch_type <cp,l2,basic> - Line search type.
454 -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
455 
456    Notes:
457     This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
458           gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
459           chooses the initial search direction as F(x) for the initial guess x.
460 
461           Only supports left non-linear preconditioning.
462 
463     Default line search is SNESLINESEARCHCP, unless a nonlinear preconditioner is used with -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC() then
464     SNESLINESEARCHL2 is used. Also supports the special purpose line search SNESLINESEARCHNCGLINEAR
465 
466    References:
467 .  * -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
468    SIAM Review, 57(4), 2015
469 
470 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESNCGGetType()`, `SNESLineSearchSetType()`
471 M*/
472 PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
473 {
474   SNES_NCG       * neP;
475 
476   PetscFunctionBegin;
477   snes->ops->destroy        = SNESDestroy_NCG;
478   snes->ops->setup          = SNESSetUp_NCG;
479   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
480   snes->ops->view           = SNESView_NCG;
481   snes->ops->solve          = SNESSolve_NCG;
482   snes->ops->reset          = SNESReset_NCG;
483 
484   snes->usesksp = PETSC_FALSE;
485   snes->usesnpc = PETSC_TRUE;
486   snes->npcside = PC_LEFT;
487 
488   snes->alwayscomputesfinalresidual = PETSC_TRUE;
489 
490   if (!snes->tolerancesset) {
491     snes->max_funcs = 30000;
492     snes->max_its   = 10000;
493     snes->stol      = 1e-20;
494   }
495 
496   PetscCall(PetscNewLog(snes,&neP));
497   snes->data   = (void*) neP;
498   neP->monitor = NULL;
499   neP->type    = SNES_NCG_PRP;
500   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG));
501   PetscFunctionReturn(0);
502 }
503