xref: /petsc/src/snes/linesearch/impls/cp/linesearchcp.c (revision cd7522eaf66b5dfe973598b9eb081fb87c780851)
1bbd5d0b3SPeter Brune #include <private/linesearchimpl.h>
26a388c36SPeter Brune #include <petscsnes.h>
3bbd5d0b3SPeter Brune 
4bbd5d0b3SPeter Brune #undef __FUNCT__
5f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_CP"
6f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchApply_CP(SNESLineSearch linesearch)
7bbd5d0b3SPeter Brune {
86a388c36SPeter Brune   PetscBool      changed_y, changed_w, domainerror;
9bbd5d0b3SPeter Brune   PetscErrorCode ierr;
106a388c36SPeter Brune   Vec             X, Y, F, W;
116a388c36SPeter Brune   SNES            snes;
12516fe3c3SPeter Brune   PetscReal       xnorm, ynorm, gnorm, steptol, atol, rtol, ltol, maxstep;
13bbd5d0b3SPeter Brune 
14bbd5d0b3SPeter Brune   PetscReal       lambda, lambda_old, lambda_update, delLambda;
15bbd5d0b3SPeter Brune   PetscScalar     fty, fty_old;
166a388c36SPeter Brune   PetscInt        i, max_its;
176a388c36SPeter Brune 
186a388c36SPeter Brune   PetscViewer     monitor;
19bbd5d0b3SPeter Brune 
20bbd5d0b3SPeter Brune   PetscFunctionBegin;
216a388c36SPeter Brune 
22f1c6b773SPeter Brune   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr);
23f1c6b773SPeter Brune   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
24f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
25f1c6b773SPeter Brune   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
26f1c6b773SPeter Brune   ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, &ltol, &max_its);CHKERRQ(ierr);
27f1c6b773SPeter Brune   ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
28f1c6b773SPeter Brune   ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr);
296a388c36SPeter Brune 
30bbd5d0b3SPeter Brune   /* precheck */
31f1c6b773SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr);
32bbd5d0b3SPeter Brune   lambda_old = 0.0;
33bbd5d0b3SPeter Brune   ierr = VecDot(F, Y, &fty_old);CHKERRQ(ierr);
346a388c36SPeter Brune   for (i = 0; i < max_its; i++) {
35bbd5d0b3SPeter Brune 
36bbd5d0b3SPeter Brune     /* compute the norm at lambda */
37bbd5d0b3SPeter Brune     ierr = VecCopy(X, W);CHKERRQ(ierr);
38bbd5d0b3SPeter Brune     ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr);
399bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
409bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
419bd66eb0SPeter Brune     }
42bbd5d0b3SPeter Brune     ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr);
43bbd5d0b3SPeter Brune 
44bbd5d0b3SPeter Brune     ierr = VecDot(F, Y, &fty);CHKERRQ(ierr);
45bbd5d0b3SPeter Brune 
46bbd5d0b3SPeter Brune     delLambda    = lambda - lambda_old;
476a388c36SPeter Brune     if (PetscAbsReal(delLambda) < steptol) break;
486a388c36SPeter Brune     if (monitor) {
496a388c36SPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
506a388c36SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: lambdas = [%g, %g], ftys = [%g, %g]\n",
51bbd5d0b3SPeter Brune                                     lambda, lambda_old, PetscRealPart(fty), PetscRealPart(fty_old));CHKERRQ(ierr);
526a388c36SPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
53bbd5d0b3SPeter Brune     }
54bbd5d0b3SPeter Brune 
55bbd5d0b3SPeter Brune     /* compute the search direction */
56bbd5d0b3SPeter Brune     lambda_update =  PetscRealPart((fty*lambda_old - fty_old*lambda) / (fty - fty_old));
57bbd5d0b3SPeter Brune     if (PetscIsInfOrNanScalar(lambda_update)) break;
586a388c36SPeter Brune     if (lambda_update > maxstep) {
59bbd5d0b3SPeter Brune       break;
60bbd5d0b3SPeter Brune     }
61bbd5d0b3SPeter Brune 
62bbd5d0b3SPeter Brune     /* compute the new state of the line search */
63bbd5d0b3SPeter Brune     lambda_old = lambda;
64bbd5d0b3SPeter Brune     lambda = lambda_update;
65bbd5d0b3SPeter Brune     fty_old = fty;
66bbd5d0b3SPeter Brune   }
67bbd5d0b3SPeter Brune   /* construct the solution */
68bbd5d0b3SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
69bbd5d0b3SPeter Brune   ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr);
709bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
719bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
729bd66eb0SPeter Brune   }
73bbd5d0b3SPeter Brune   /* postcheck */
74f1c6b773SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr);
75bbd5d0b3SPeter Brune   if (changed_y) {
76bbd5d0b3SPeter Brune     ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr);
779bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
789bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, X);CHKERRQ(ierr);
799bd66eb0SPeter Brune     }
80bbd5d0b3SPeter Brune   } else {
81bbd5d0b3SPeter Brune     ierr = VecCopy(W, X);CHKERRQ(ierr);
82bbd5d0b3SPeter Brune   }
83bbd5d0b3SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
846a388c36SPeter Brune   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
856a388c36SPeter Brune   if (domainerror) {
86f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
87bbd5d0b3SPeter Brune     PetscFunctionReturn(0);
88bbd5d0b3SPeter Brune   }
89bbd5d0b3SPeter Brune 
90f1c6b773SPeter Brune   ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr);
91f1c6b773SPeter Brune   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
92bbd5d0b3SPeter Brune 
93f1c6b773SPeter Brune   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
946a388c36SPeter Brune 
956a388c36SPeter Brune   if (monitor) {
966a388c36SPeter Brune     ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
976a388c36SPeter Brune     ierr = PetscViewerASCIIPrintf(monitor,"    Line search terminated: lambda = %g, fnorms = %g\n", lambda, gnorm);CHKERRQ(ierr);
986a388c36SPeter Brune     ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
99bbd5d0b3SPeter Brune   }
1006a388c36SPeter Brune   if (lambda <= steptol) {
101f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
102bbd5d0b3SPeter Brune   }
103bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
104bbd5d0b3SPeter Brune }
105bbd5d0b3SPeter Brune 
106bbd5d0b3SPeter Brune #undef __FUNCT__
107f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_CP"
108954494b2SJed Brown /*MC
109*cd7522eaSPeter Brune    SNES_LINESEARCH_CP - Critical point line search. This line search assumes that there exists some
110*cd7522eaSPeter Brune    artificial G(x) for which the SNESFunction F(x) = grad G(x).  Therefore, this line search seeks
111*cd7522eaSPeter Brune    to find roots of f^ty via a secant method.
112*cd7522eaSPeter Brune 
113*cd7522eaSPeter Brune    Options Database Keys:
114*cd7522eaSPeter Brune .  -snes_linesearch_damping - initial trial step length
115*cd7522eaSPeter Brune 
116*cd7522eaSPeter Brune    Notes:
117*cd7522eaSPeter Brune    This method is the preferred line search for SNESQN and SNESNCG.
118954494b2SJed Brown 
119954494b2SJed Brown    Level: advanced
120954494b2SJed Brown 
121f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
122954494b2SJed Brown 
123f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
124954494b2SJed Brown M*/
125f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch)
126bbd5d0b3SPeter Brune {
127bbd5d0b3SPeter Brune   PetscFunctionBegin;
128f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_CP;
129bbd5d0b3SPeter Brune   linesearch->ops->destroy        = PETSC_NULL;
130bbd5d0b3SPeter Brune   linesearch->ops->setfromoptions = PETSC_NULL;
131bbd5d0b3SPeter Brune   linesearch->ops->reset          = PETSC_NULL;
132bbd5d0b3SPeter Brune   linesearch->ops->view           = PETSC_NULL;
133bbd5d0b3SPeter Brune   linesearch->ops->setup          = PETSC_NULL;
134bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
135bbd5d0b3SPeter Brune }
136