xref: /petsc/src/snes/linesearch/impls/cp/linesearchcp.c (revision bbd5d0b317af04eb83c39e72581bf86c64bd741a)
1*bbd5d0b3SPeter Brune #include <private/linesearchimpl.h>
2*bbd5d0b3SPeter Brune #include <private/snesimpl.h>
3*bbd5d0b3SPeter Brune 
4*bbd5d0b3SPeter Brune #undef __FUNCT__
5*bbd5d0b3SPeter Brune #define __FUNCT__ "LineSearchApply_CP"
6*bbd5d0b3SPeter Brune 
7*bbd5d0b3SPeter Brune /*@C
8*bbd5d0b3SPeter Brune    LineSearchCP - This routine is not a line search at all;
9*bbd5d0b3SPeter Brune    it simply uses the full step.  Thus, this routine is intended
10*bbd5d0b3SPeter Brune    to serve as a template and is not recommended for general use.
11*bbd5d0b3SPeter Brune 
12*bbd5d0b3SPeter Brune    Logically Collective on SNES and Vec
13*bbd5d0b3SPeter Brune 
14*bbd5d0b3SPeter Brune    Input Parameters:
15*bbd5d0b3SPeter Brune +  snes - nonlinear context
16*bbd5d0b3SPeter Brune .  lsctx - optional context for line search (not used here)
17*bbd5d0b3SPeter Brune .  x - current iterate
18*bbd5d0b3SPeter Brune .  f - residual evaluated at x
19*bbd5d0b3SPeter Brune .  y - search direction
20*bbd5d0b3SPeter Brune .  fnorm - 2-norm of f
21*bbd5d0b3SPeter Brune -  xnorm - norm of x if known, otherwise 0
22*bbd5d0b3SPeter Brune 
23*bbd5d0b3SPeter Brune    Output Parameters:
24*bbd5d0b3SPeter Brune +  g - residual evaluated at new iterate y
25*bbd5d0b3SPeter Brune .  w - new iterate
26*bbd5d0b3SPeter Brune .  gnorm - 2-norm of g
27*bbd5d0b3SPeter Brune .  ynorm - 2-norm of search length
28*bbd5d0b3SPeter Brune -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
29*bbd5d0b3SPeter Brune 
30*bbd5d0b3SPeter Brune    Options Database Key:
31*bbd5d0b3SPeter Brune .  -snes_ls basic - Activates SNESLineSearchNo()
32*bbd5d0b3SPeter Brune 
33*bbd5d0b3SPeter Brune    Level: advanced
34*bbd5d0b3SPeter Brune 
35*bbd5d0b3SPeter Brune .keywords: SNES, nonlinear, line search, cubic
36*bbd5d0b3SPeter Brune 
37*bbd5d0b3SPeter Brune .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
38*bbd5d0b3SPeter Brune           SNESLineSearchSet(), SNESLineSearchNoNorms()
39*bbd5d0b3SPeter Brune @*/
40*bbd5d0b3SPeter Brune PetscErrorCode  LineSearchApply_CP(LineSearch linesearch)
41*bbd5d0b3SPeter Brune {
42*bbd5d0b3SPeter Brune 
43*bbd5d0b3SPeter Brune   PetscBool      changed_y, changed_w;
44*bbd5d0b3SPeter Brune   PetscErrorCode ierr;
45*bbd5d0b3SPeter Brune   Vec             X = linesearch->vec_sol;
46*bbd5d0b3SPeter Brune   Vec             F  = linesearch->vec_func;
47*bbd5d0b3SPeter Brune   Vec             Y  = linesearch->vec_update;
48*bbd5d0b3SPeter Brune   Vec             W =  linesearch->vec_sol_new;
49*bbd5d0b3SPeter Brune   SNES            snes  = linesearch->snes;
50*bbd5d0b3SPeter Brune   PetscReal       *gnorm  = &linesearch->fnorm;
51*bbd5d0b3SPeter Brune   PetscReal       *ynorm  = &linesearch->ynorm;
52*bbd5d0b3SPeter Brune   PetscReal       *xnorm =  &linesearch->xnorm;
53*bbd5d0b3SPeter Brune 
54*bbd5d0b3SPeter Brune   PetscReal       lambda, lambda_old, lambda_update, delLambda;
55*bbd5d0b3SPeter Brune   PetscScalar     fty, fty_old;
56*bbd5d0b3SPeter Brune   PetscInt        i;
57*bbd5d0b3SPeter Brune 
58*bbd5d0b3SPeter Brune   PetscFunctionBegin;
59*bbd5d0b3SPeter Brune   /* precheck */
60*bbd5d0b3SPeter Brune   ierr = LineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr);
61*bbd5d0b3SPeter Brune   lambda = linesearch->lambda;
62*bbd5d0b3SPeter Brune   lambda_old = 0.0;
63*bbd5d0b3SPeter Brune   ierr = VecDot(F, Y, &fty_old);CHKERRQ(ierr);
64*bbd5d0b3SPeter Brune   linesearch->success = PETSC_TRUE;
65*bbd5d0b3SPeter Brune   for (i = 0; i < linesearch->max_its; i++) {
66*bbd5d0b3SPeter Brune 
67*bbd5d0b3SPeter Brune     /* compute the norm at lambda */
68*bbd5d0b3SPeter Brune     ierr = VecCopy(X, W);CHKERRQ(ierr);
69*bbd5d0b3SPeter Brune     ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr);
70*bbd5d0b3SPeter Brune     ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr);
71*bbd5d0b3SPeter Brune 
72*bbd5d0b3SPeter Brune     ierr = VecDot(F, Y, &fty);CHKERRQ(ierr);
73*bbd5d0b3SPeter Brune 
74*bbd5d0b3SPeter Brune     delLambda    = lambda - lambda_old;
75*bbd5d0b3SPeter Brune     if (PetscAbsReal(delLambda) < linesearch->steptol) break;
76*bbd5d0b3SPeter Brune     if (linesearch->monitor) {
77*bbd5d0b3SPeter Brune       ierr = PetscViewerASCIIAddTab(linesearch->monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
78*bbd5d0b3SPeter Brune       ierr = PetscViewerASCIIPrintf(linesearch->monitor,"    Line search: lambdas = [%g, %g], ftys = [%g, %g]\n",
79*bbd5d0b3SPeter Brune                                     lambda, lambda_old, PetscRealPart(fty), PetscRealPart(fty_old));CHKERRQ(ierr);
80*bbd5d0b3SPeter Brune       ierr = PetscViewerASCIISubtractTab(linesearch->monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
81*bbd5d0b3SPeter Brune     }
82*bbd5d0b3SPeter Brune 
83*bbd5d0b3SPeter Brune     /* compute the search direction */
84*bbd5d0b3SPeter Brune     lambda_update =  PetscRealPart((fty*lambda_old - fty_old*lambda) / (fty - fty_old));
85*bbd5d0b3SPeter Brune     if (PetscIsInfOrNanScalar(lambda_update)) break;
86*bbd5d0b3SPeter Brune     if (lambda_update > linesearch->maxstep) {
87*bbd5d0b3SPeter Brune       break;
88*bbd5d0b3SPeter Brune     }
89*bbd5d0b3SPeter Brune 
90*bbd5d0b3SPeter Brune     /* compute the new state of the line search */
91*bbd5d0b3SPeter Brune     lambda_old = lambda;
92*bbd5d0b3SPeter Brune     lambda = lambda_update;
93*bbd5d0b3SPeter Brune     fty_old = fty;
94*bbd5d0b3SPeter Brune   }
95*bbd5d0b3SPeter Brune   /* construct the solution */
96*bbd5d0b3SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
97*bbd5d0b3SPeter Brune   ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr);
98*bbd5d0b3SPeter Brune 
99*bbd5d0b3SPeter Brune   /* postcheck */
100*bbd5d0b3SPeter Brune   ierr = LineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr);
101*bbd5d0b3SPeter Brune   if (changed_y) {
102*bbd5d0b3SPeter Brune     ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr);
103*bbd5d0b3SPeter Brune   } else {
104*bbd5d0b3SPeter Brune     ierr = VecCopy(W, X);CHKERRQ(ierr);
105*bbd5d0b3SPeter Brune   }
106*bbd5d0b3SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
107*bbd5d0b3SPeter Brune   if (snes->domainerror) {
108*bbd5d0b3SPeter Brune     linesearch->success = PETSC_FALSE;
109*bbd5d0b3SPeter Brune     PetscFunctionReturn(0);
110*bbd5d0b3SPeter Brune   }
111*bbd5d0b3SPeter Brune 
112*bbd5d0b3SPeter Brune   ierr = VecNormBegin(F, NORM_2, gnorm);CHKERRQ(ierr);
113*bbd5d0b3SPeter Brune   ierr = VecNormBegin(X, NORM_2, xnorm);CHKERRQ(ierr);
114*bbd5d0b3SPeter Brune   ierr = VecNormBegin(Y, NORM_2, ynorm);CHKERRQ(ierr);
115*bbd5d0b3SPeter Brune   ierr = VecNormEnd(F, NORM_2, gnorm);CHKERRQ(ierr);
116*bbd5d0b3SPeter Brune   ierr = VecNormEnd(X, NORM_2, xnorm);CHKERRQ(ierr);
117*bbd5d0b3SPeter Brune   ierr = VecNormEnd(Y, NORM_2, ynorm);CHKERRQ(ierr);
118*bbd5d0b3SPeter Brune 
119*bbd5d0b3SPeter Brune   linesearch->lambda = lambda;
120*bbd5d0b3SPeter Brune   if (linesearch->monitor) {
121*bbd5d0b3SPeter Brune     ierr = PetscViewerASCIIAddTab(linesearch->monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
122*bbd5d0b3SPeter Brune     ierr = PetscViewerASCIIPrintf(linesearch->monitor,"    Line search terminated: lambda = %g, fnorms = %g\n", lambda, linesearch->fnorm);CHKERRQ(ierr);
123*bbd5d0b3SPeter Brune     ierr = PetscViewerASCIISubtractTab(linesearch->monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
124*bbd5d0b3SPeter Brune   }
125*bbd5d0b3SPeter Brune   if (lambda <= linesearch->steptol) {
126*bbd5d0b3SPeter Brune     linesearch->success = PETSC_FALSE;
127*bbd5d0b3SPeter Brune   }
128*bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
129*bbd5d0b3SPeter Brune }
130*bbd5d0b3SPeter Brune 
131*bbd5d0b3SPeter Brune EXTERN_C_BEGIN
132*bbd5d0b3SPeter Brune #undef __FUNCT__
133*bbd5d0b3SPeter Brune #define __FUNCT__ "LineSearchCreate_CP"
134*bbd5d0b3SPeter Brune PetscErrorCode LineSearchCreate_CP(LineSearch linesearch)
135*bbd5d0b3SPeter Brune {
136*bbd5d0b3SPeter Brune   PetscFunctionBegin;
137*bbd5d0b3SPeter Brune   linesearch->ops->apply          = LineSearchApply_CP;
138*bbd5d0b3SPeter Brune   linesearch->ops->destroy        = PETSC_NULL;
139*bbd5d0b3SPeter Brune   linesearch->ops->setfromoptions = PETSC_NULL;
140*bbd5d0b3SPeter Brune   linesearch->ops->reset          = PETSC_NULL;
141*bbd5d0b3SPeter Brune   linesearch->ops->view           = PETSC_NULL;
142*bbd5d0b3SPeter Brune   linesearch->ops->setup          = PETSC_NULL;
143*bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
144*bbd5d0b3SPeter Brune }
145*bbd5d0b3SPeter Brune EXTERN_C_END
146