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