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