xref: /petsc/src/snes/linesearch/impls/basic/linesearchbasic.c (revision 9bd66eb0dc31efbcccaa04c7c4d5de114a77b297)
1 #include <private/linesearchimpl.h>
2 #include <private/snesimpl.h>
3 
4 /*MC
5    PetscLineSearchBasic - This routine is not a line search at all;
6    it simply uses the full step.  Thus, this routine is intended
7    to serve as a template and is not recommended for general use.
8 
9    Level: advanced
10 
11 .keywords: SNES, PetscLineSearch, damping
12 
13 .seealso: PetscLineSearchCreate(), PetscLineSearchSetType()
14 M*/
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "PetscLineSearchApply_Basic"
18 
19 PetscErrorCode  PetscLineSearchApply_Basic(PetscLineSearch linesearch)
20 {
21   PetscBool      changed_y, changed_w;
22   PetscErrorCode ierr;
23   Vec            X, F, Y, W;
24   SNES           snes;
25   PetscReal      gnorm, xnorm, ynorm, lambda;
26   PetscBool      domainerror;
27 
28   PetscFunctionBegin;
29 
30   ierr = PetscLineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr);
31   ierr = PetscLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
32   ierr = PetscLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
33   ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
34   ierr = PetscLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
35 
36   /* precheck */
37   ierr = PetscLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr);
38 
39   /* update */
40   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
41   if (linesearch->ops->viproject) {
42     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
43   }
44 
45   /* postcheck */
46   ierr = PetscLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr);
47   if (changed_y) {
48     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
49     if (linesearch->ops->viproject) {
50       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
51     }
52   }
53   ierr = SNESComputeFunction(snes,W,F);CHKERRQ(ierr);
54   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
55   if (domainerror) {
56     ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);
57     PetscFunctionReturn(0);
58   }
59 
60   ierr = VecNorm(Y, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
61   ierr = VecNorm(W, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
62   if (linesearch->ops->vinorm) {
63     linesearch->fnorm = gnorm;
64     ierr = (*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm);CHKERRQ(ierr);
65   } else {
66     ierr = VecNorm(F,NORM_2,&linesearch->fnorm);CHKERRQ(ierr);
67   }
68 
69   /*
70   ierr = PetscLineSearchComputeNorms(linesearch);CHKERRQ(ierr);
71    */
72 
73   /* copy the solution over */
74   ierr = VecCopy(W, X);CHKERRQ(ierr);
75   PetscFunctionReturn(0);
76 }
77 
78 EXTERN_C_BEGIN
79 #undef __FUNCT__
80 #define __FUNCT__ "PetscLineSearchCreate_Basic"
81 PetscErrorCode PetscLineSearchCreate_Basic(PetscLineSearch linesearch)
82 {
83   PetscFunctionBegin;
84   linesearch->ops->apply          = PetscLineSearchApply_Basic;
85   linesearch->ops->destroy        = PETSC_NULL;
86   linesearch->ops->setfromoptions = PETSC_NULL;
87   linesearch->ops->reset          = PETSC_NULL;
88   linesearch->ops->view           = PETSC_NULL;
89   linesearch->ops->setup          = PETSC_NULL;
90   PetscFunctionReturn(0);
91 }
92 EXTERN_C_END
93