xref: /petsc/src/snes/linesearch/impls/basic/linesearchbasic.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
1*b45d2f2cSJed Brown #include <petsc-private/linesearchimpl.h>
2*b45d2f2cSJed Brown #include <petsc-private/snesimpl.h>
3bf7f4e0aSPeter Brune 
4f40b411bSPeter Brune #undef __FUNCT__
5f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_Basic"
6f1c6b773SPeter Brune static PetscErrorCode  SNESLineSearchApply_Basic(SNESLineSearch linesearch)
7bf7f4e0aSPeter Brune {
8bf7f4e0aSPeter Brune   PetscBool      changed_y, changed_w;
9bf7f4e0aSPeter Brune   PetscErrorCode ierr;
106a388c36SPeter Brune   Vec            X, F, Y, W;
116a388c36SPeter Brune   SNES           snes;
126a388c36SPeter Brune   PetscReal      gnorm, xnorm, ynorm, lambda;
136a388c36SPeter Brune   PetscBool      domainerror;
14bf7f4e0aSPeter Brune 
15bf7f4e0aSPeter Brune   PetscFunctionBegin;
166a388c36SPeter Brune 
17f1c6b773SPeter Brune   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, PETSC_NULL);CHKERRQ(ierr);
18f1c6b773SPeter Brune   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
19f1c6b773SPeter Brune   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
20f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
21f1c6b773SPeter Brune   ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
22bf7f4e0aSPeter Brune 
23bf7f4e0aSPeter Brune   /* precheck */
24f1c6b773SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr);
25bf7f4e0aSPeter Brune 
26bf7f4e0aSPeter Brune   /* update */
276a388c36SPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
289bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
299bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
309bd66eb0SPeter Brune   }
31bf7f4e0aSPeter Brune 
32bf7f4e0aSPeter Brune   /* postcheck */
33f1c6b773SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr);
34bf7f4e0aSPeter Brune   if (changed_y) {
356a388c36SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
369bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
379bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
389bd66eb0SPeter Brune     }
39bf7f4e0aSPeter Brune   }
40bf7f4e0aSPeter Brune   ierr = SNESComputeFunction(snes,W,F);CHKERRQ(ierr);
416a388c36SPeter Brune   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
426a388c36SPeter Brune   if (domainerror) {
43f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
44bf7f4e0aSPeter Brune     PetscFunctionReturn(0);
45bf7f4e0aSPeter Brune   }
466a388c36SPeter Brune 
479bd66eb0SPeter Brune   ierr = VecNorm(Y, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
489bd66eb0SPeter Brune   ierr = VecNorm(W, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
499bd66eb0SPeter Brune   if (linesearch->ops->vinorm) {
509bd66eb0SPeter Brune     linesearch->fnorm = gnorm;
519bd66eb0SPeter Brune     ierr = (*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm);CHKERRQ(ierr);
529bd66eb0SPeter Brune   } else {
539bd66eb0SPeter Brune     ierr = VecNorm(F,NORM_2,&linesearch->fnorm);CHKERRQ(ierr);
549bd66eb0SPeter Brune   }
559bd66eb0SPeter Brune 
569bd66eb0SPeter Brune   /*
57f1c6b773SPeter Brune   ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr);
589bd66eb0SPeter Brune    */
596a388c36SPeter Brune 
60bf7f4e0aSPeter Brune   /* copy the solution over */
61bf7f4e0aSPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
62bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
63bf7f4e0aSPeter Brune }
64bf7f4e0aSPeter Brune 
65bf7f4e0aSPeter Brune #undef __FUNCT__
66f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_Basic"
67954494b2SJed Brown /*MC
68cd7522eaSPeter Brune    SNES_LINESEARCH_BASIC - This line search implementation is not a line
69cd7522eaSPeter Brune    search at all; it simply uses the full step.  Thus, this routine is intended
70cd7522eaSPeter Brune    for methods with well-scaled updates; i.e. Newton's method (SNESLS), in on
71cd7522eaSPeter Brune    well-behaved problems.
72cd7522eaSPeter Brune 
73cd7522eaSPeter Brune    Options Database Keys:
74cd7522eaSPeter Brune    -snes_linesearch_damping (1.0) damping parameter.
75cd7522eaSPeter Brune 
76cd7522eaSPeter Brune    Notes:
77cd7522eaSPeter Brune    For methods with ill-scaled updates (SNESNRICHARDSON, SNESNCG), a small
78cd7522eaSPeter Brune    damping parameter may yield satisfactory but slow convergence despite
79cd7522eaSPeter Brune    the simplicity of the line search.
80954494b2SJed Brown 
81954494b2SJed Brown    Level: advanced
82954494b2SJed Brown 
83f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
84954494b2SJed Brown 
85f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
86954494b2SJed Brown M*/
87f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch)
88bf7f4e0aSPeter Brune {
89bf7f4e0aSPeter Brune   PetscFunctionBegin;
90f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_Basic;
91bf7f4e0aSPeter Brune   linesearch->ops->destroy        = PETSC_NULL;
92bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = PETSC_NULL;
93bf7f4e0aSPeter Brune   linesearch->ops->reset          = PETSC_NULL;
94bf7f4e0aSPeter Brune   linesearch->ops->view           = PETSC_NULL;
95bf7f4e0aSPeter Brune   linesearch->ops->setup          = PETSC_NULL;
96bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
97bf7f4e0aSPeter Brune }
98