xref: /petsc/src/snes/linesearch/impls/basic/linesearchbasic.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 #include <petsc/private/linesearchimpl.h>
2 #include <petsc/private/snesimpl.h>
3 
4 static PetscErrorCode SNESLineSearchApply_Basic(SNESLineSearch linesearch) {
5   PetscBool changed_y, changed_w;
6   Vec       X, F, Y, W;
7   SNES      snes;
8   PetscReal gnorm, xnorm, ynorm, lambda;
9   PetscBool domainerror;
10 
11   PetscFunctionBegin;
12   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL));
13   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
14   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
15   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
16   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
17 
18   /* precheck */
19   PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
20 
21   /* update */
22   PetscCall(VecWAXPY(W, -lambda, Y, X));
23   if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
24 
25   /* postcheck */
26   PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
27   if (changed_y) {
28     PetscCall(VecWAXPY(W, -lambda, Y, X));
29     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
30   }
31   if (linesearch->norms || snes->iter < snes->max_its - 1) {
32     PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
33     PetscCall(SNESGetFunctionDomainError(snes, &domainerror));
34     if (domainerror) {
35       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_DOMAIN));
36       PetscFunctionReturn(0);
37     }
38   }
39 
40   if (linesearch->norms) {
41     if (!linesearch->ops->vinorm) PetscCall(VecNormBegin(F, NORM_2, &linesearch->fnorm));
42     PetscCall(VecNormBegin(Y, NORM_2, &linesearch->ynorm));
43     PetscCall(VecNormBegin(W, NORM_2, &linesearch->xnorm));
44     if (!linesearch->ops->vinorm) PetscCall(VecNormEnd(F, NORM_2, &linesearch->fnorm));
45     PetscCall(VecNormEnd(Y, NORM_2, &linesearch->ynorm));
46     PetscCall(VecNormEnd(W, NORM_2, &linesearch->xnorm));
47 
48     if (linesearch->ops->vinorm) {
49       linesearch->fnorm = gnorm;
50 
51       PetscCall((*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm));
52     } else {
53       PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
54     }
55   }
56 
57   /* copy the solution over */
58   PetscCall(VecCopy(W, X));
59   PetscFunctionReturn(0);
60 }
61 
62 /*MC
63    SNESLINESEARCHBASIC - This line search implementation is not a line
64    search at all; it simply uses the full step.  Thus, this routine is intended
65    for methods with well-scaled updates; i.e. Newton's method (SNESNEWTONLS), on
66    well-behaved problems. Also named as SNESLINESEARCHNONE
67 
68    Options Database Keys:
69 +   -snes_linesearch_damping <damping> - search vector is scaled by this amount, default is 1.0
70 -   -snes_linesearch_norms <flag> - whether to compute norms or not, default is true (SNESLineSearchSetComputeNorms())
71 
72    Notes:
73    For methods with ill-scaled updates (SNESNRICHARDSON, SNESNCG), a small
74    damping parameter may yield satisfactory but slow convergence despite
75    the simplicity of the line search.
76 
77    Level: advanced
78 
79 .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLineSearchSetDamping()`, `SNESLineSearchSetComputeNorms()`
80 M*/
81 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch) {
82   PetscFunctionBegin;
83   linesearch->ops->apply          = SNESLineSearchApply_Basic;
84   linesearch->ops->destroy        = NULL;
85   linesearch->ops->setfromoptions = NULL;
86   linesearch->ops->reset          = NULL;
87   linesearch->ops->view           = NULL;
88   linesearch->ops->setup          = NULL;
89   PetscFunctionReturn(0);
90 }
91