xref: /petsc/src/snes/linesearch/impls/basic/linesearchbasic.c (revision 6a388c361ff2893f4363aecae9e02aa0b4e91a18)
1bf7f4e0aSPeter Brune #include <private/linesearchimpl.h>
2bf7f4e0aSPeter Brune #include <private/snesimpl.h>
3bf7f4e0aSPeter Brune 
4bf7f4e0aSPeter Brune #undef __FUNCT__
5bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchApply_Basic"
6bf7f4e0aSPeter Brune 
7bf7f4e0aSPeter Brune /*@C
8bf7f4e0aSPeter Brune    SNESLineSearchBasic - This routine is not a line search at all;
9bf7f4e0aSPeter Brune    it simply uses the full step.  Thus, this routine is intended
10bf7f4e0aSPeter Brune    to serve as a template and is not recommended for general use.
11bf7f4e0aSPeter Brune 
12bf7f4e0aSPeter Brune    Logically Collective on SNES and Vec
13bf7f4e0aSPeter Brune 
14bf7f4e0aSPeter Brune    Input Parameters:
15bf7f4e0aSPeter Brune +  snes - nonlinear context
16bf7f4e0aSPeter Brune .  lsctx - optional context for line search (not used here)
17bf7f4e0aSPeter Brune .  x - current iterate
18bf7f4e0aSPeter Brune .  f - residual evaluated at x
19bf7f4e0aSPeter Brune .  y - search direction
20bf7f4e0aSPeter Brune .  fnorm - 2-norm of f
21bf7f4e0aSPeter Brune -  xnorm - norm of x if known, otherwise 0
22bf7f4e0aSPeter Brune 
23bf7f4e0aSPeter Brune    Output Parameters:
24bf7f4e0aSPeter Brune +  g - residual evaluated at new iterate y
25bf7f4e0aSPeter Brune .  w - new iterate
26bf7f4e0aSPeter Brune .  gnorm - 2-norm of g
27bf7f4e0aSPeter Brune .  ynorm - 2-norm of search length
28bf7f4e0aSPeter Brune -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
29bf7f4e0aSPeter Brune 
30bf7f4e0aSPeter Brune    Options Database Key:
31bf7f4e0aSPeter Brune .  -snes_ls basic - Activates SNESLineSearchNo()
32bf7f4e0aSPeter Brune 
33bf7f4e0aSPeter Brune    Level: advanced
34bf7f4e0aSPeter Brune 
35bf7f4e0aSPeter Brune .keywords: SNES, nonlinear, line search, cubic
36bf7f4e0aSPeter Brune 
37bf7f4e0aSPeter Brune .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
38bf7f4e0aSPeter Brune           SNESLineSearchSet(), SNESLineSearchNoNorms()
39bf7f4e0aSPeter Brune @*/
40bf7f4e0aSPeter Brune PetscErrorCode  LineSearchApply_Basic(LineSearch linesearch)
41bf7f4e0aSPeter Brune {
42bf7f4e0aSPeter Brune   PetscBool      changed_y, changed_w;
43bf7f4e0aSPeter Brune   PetscErrorCode ierr;
44*6a388c36SPeter Brune   Vec            X, F, Y, W;
45*6a388c36SPeter Brune   SNES           snes;
46*6a388c36SPeter Brune   PetscReal      gnorm, xnorm, ynorm, lambda;
47*6a388c36SPeter Brune   PetscBool      domainerror;
48bf7f4e0aSPeter Brune 
49bf7f4e0aSPeter Brune   PetscFunctionBegin;
50*6a388c36SPeter Brune 
51*6a388c36SPeter Brune   ierr = LineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr);
52*6a388c36SPeter Brune   ierr = LineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
53*6a388c36SPeter Brune   ierr = LineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
54*6a388c36SPeter Brune   ierr = LineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
55*6a388c36SPeter Brune   ierr = LineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
56bf7f4e0aSPeter Brune 
57bf7f4e0aSPeter Brune   /* precheck */
58bf7f4e0aSPeter Brune   ierr = LineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr);
59bf7f4e0aSPeter Brune 
60bf7f4e0aSPeter Brune   /* update */
61*6a388c36SPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
62bf7f4e0aSPeter Brune 
63bf7f4e0aSPeter Brune   /* postcheck */
64bf7f4e0aSPeter Brune   ierr = LineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr);
65bf7f4e0aSPeter Brune   if (changed_y) {
66*6a388c36SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
67bf7f4e0aSPeter Brune   }
68bf7f4e0aSPeter Brune   ierr = SNESComputeFunction(snes,W,F);CHKERRQ(ierr);
69*6a388c36SPeter Brune   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
70*6a388c36SPeter Brune   if (domainerror) {
71*6a388c36SPeter Brune     ierr = LineSearchSetSuccess(linesearch, PETSC_FALSE);
72bf7f4e0aSPeter Brune     PetscFunctionReturn(0);
73bf7f4e0aSPeter Brune   }
74*6a388c36SPeter Brune 
75*6a388c36SPeter Brune   ierr = LineSearchComputeNorms(linesearch);CHKERRQ(ierr);
76*6a388c36SPeter Brune 
77bf7f4e0aSPeter Brune   /* copy the solution over */
78bf7f4e0aSPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
79bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
80bf7f4e0aSPeter Brune }
81bf7f4e0aSPeter Brune 
82bf7f4e0aSPeter Brune EXTERN_C_BEGIN
83bf7f4e0aSPeter Brune #undef __FUNCT__
84bf7f4e0aSPeter Brune #define __FUNCT__ "LineSearchCreate_Basic"
85bf7f4e0aSPeter Brune PetscErrorCode LineSearchCreate_Basic(LineSearch linesearch)
86bf7f4e0aSPeter Brune {
87bf7f4e0aSPeter Brune   PetscFunctionBegin;
88bf7f4e0aSPeter Brune   linesearch->ops->apply          = LineSearchApply_Basic;
89bf7f4e0aSPeter Brune   linesearch->ops->destroy        = PETSC_NULL;
90bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = PETSC_NULL;
91bf7f4e0aSPeter Brune   linesearch->ops->reset          = PETSC_NULL;
92bf7f4e0aSPeter Brune   linesearch->ops->view           = PETSC_NULL;
93bf7f4e0aSPeter Brune   linesearch->ops->setup          = PETSC_NULL;
94bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
95bf7f4e0aSPeter Brune }
96bf7f4e0aSPeter Brune EXTERN_C_END
97