xref: /petsc/src/snes/linesearch/impls/basic/linesearchbasic.c (revision f1c6b773534bc9d3f9e058625aa590d993f2f137)
1 #include <private/linesearchimpl.h>
2 #include <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, &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, &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   ierr = SNESComputeFunction(snes,W,F);CHKERRQ(ierr);
41   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
42   if (domainerror) {
43     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
44     PetscFunctionReturn(0);
45   }
46 
47   ierr = VecNorm(Y, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
48   ierr = VecNorm(W, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
49   if (linesearch->ops->vinorm) {
50     linesearch->fnorm = gnorm;
51     ierr = (*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm);CHKERRQ(ierr);
52   } else {
53     ierr = VecNorm(F,NORM_2,&linesearch->fnorm);CHKERRQ(ierr);
54   }
55 
56   /*
57   ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr);
58    */
59 
60   /* copy the solution over */
61   ierr = VecCopy(W, X);CHKERRQ(ierr);
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "SNESLineSearchCreate_Basic"
67 /*MC
68    SNES_LINESEARCH_BASIC - This routine is not a line search at all;
69    it simply uses the full step.  Thus, this routine is intended
70    to serve as a template and is not recommended for general use.
71 
72    Level: advanced
73 
74 .keywords: SNES, SNESLineSearch, damping
75 
76 .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
77 M*/
78 PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch)
79 {
80   PetscFunctionBegin;
81   linesearch->ops->apply          = SNESLineSearchApply_Basic;
82   linesearch->ops->destroy        = PETSC_NULL;
83   linesearch->ops->setfromoptions = PETSC_NULL;
84   linesearch->ops->reset          = PETSC_NULL;
85   linesearch->ops->view           = PETSC_NULL;
86   linesearch->ops->setup          = PETSC_NULL;
87   PetscFunctionReturn(0);
88 }
89