xref: /petsc/src/snes/linesearch/impls/basic/linesearchbasic.c (revision 04d7464b71f14a73adf7cb3f664404cc6936f6d1)
1b45d2f2cSJed Brown #include <petsc-private/linesearchimpl.h>
2b45d2f2cSJed 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 */
247b1df9c1SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch,X,Y,&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 */
337b1df9c1SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&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   }
405fc09b23SPeter Brune   if (linesearch->norms || snes->iter < snes->max_its-1) {
41bf7f4e0aSPeter Brune     ierr = SNESComputeFunction(snes,W,F);CHKERRQ(ierr);
426a388c36SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
436a388c36SPeter Brune     if (domainerror) {
44f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
45bf7f4e0aSPeter Brune       PetscFunctionReturn(0);
46bf7f4e0aSPeter Brune     }
475fc09b23SPeter Brune   }
486a388c36SPeter Brune 
495fc09b23SPeter Brune   if (linesearch->norms) {
505fc09b23SPeter Brune     if (!linesearch->ops->vinorm) VecNormBegin(F, NORM_2, &linesearch->fnorm);
515fc09b23SPeter Brune     ierr = VecNormBegin(Y, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
525fc09b23SPeter Brune     ierr = VecNormBegin(W, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
535fc09b23SPeter Brune     if (!linesearch->ops->vinorm) VecNormEnd(F, NORM_2, &linesearch->fnorm);
545fc09b23SPeter Brune     ierr = VecNormEnd(Y, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
555fc09b23SPeter Brune     ierr = VecNormEnd(W, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
565fc09b23SPeter Brune 
579bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
589bd66eb0SPeter Brune       linesearch->fnorm = gnorm;
599bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm);CHKERRQ(ierr);
609bd66eb0SPeter Brune     } else {
619bd66eb0SPeter Brune       ierr = VecNorm(F,NORM_2,&linesearch->fnorm);CHKERRQ(ierr);
629bd66eb0SPeter Brune     }
635fc09b23SPeter Brune   }
646a388c36SPeter Brune 
65bf7f4e0aSPeter Brune   /* copy the solution over */
66bf7f4e0aSPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
67bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
68bf7f4e0aSPeter Brune }
69bf7f4e0aSPeter Brune 
70bf7f4e0aSPeter Brune #undef __FUNCT__
71f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_Basic"
72954494b2SJed Brown /*MC
731a4f838cSPeter Brune    SNESLINESEARCHBASIC - This line search implementation is not a line
74cd7522eaSPeter Brune    search at all; it simply uses the full step.  Thus, this routine is intended
75*04d7464bSBarry Smith    for methods with well-scaled updates; i.e. Newton's method (SNESNEWTONLS), on
76cd7522eaSPeter Brune    well-behaved problems.
77cd7522eaSPeter Brune 
78cd7522eaSPeter Brune    Options Database Keys:
795fc09b23SPeter Brune +   -snes_linesearch_damping (1.0) damping parameter.
805fc09b23SPeter Brune -   -snes_linesearch_norms (true) whether to compute norms or not.
81cd7522eaSPeter Brune 
82cd7522eaSPeter Brune    Notes:
83cd7522eaSPeter Brune    For methods with ill-scaled updates (SNESNRICHARDSON, SNESNCG), a small
84cd7522eaSPeter Brune    damping parameter may yield satisfactory but slow convergence despite
85cd7522eaSPeter Brune    the simplicity of the line search.
86954494b2SJed Brown 
87954494b2SJed Brown    Level: advanced
88954494b2SJed Brown 
89f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
90954494b2SJed Brown 
91f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
92954494b2SJed Brown M*/
93f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch)
94bf7f4e0aSPeter Brune {
95bf7f4e0aSPeter Brune   PetscFunctionBegin;
96f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_Basic;
97bf7f4e0aSPeter Brune   linesearch->ops->destroy        = PETSC_NULL;
98bf7f4e0aSPeter Brune   linesearch->ops->setfromoptions = PETSC_NULL;
99bf7f4e0aSPeter Brune   linesearch->ops->reset          = PETSC_NULL;
100bf7f4e0aSPeter Brune   linesearch->ops->view           = PETSC_NULL;
101bf7f4e0aSPeter Brune   linesearch->ops->setup          = PETSC_NULL;
102bf7f4e0aSPeter Brune   PetscFunctionReturn(0);
103bf7f4e0aSPeter Brune }
104