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