xref: /petsc/src/snes/linesearch/impls/basic/linesearchbasic.c (revision e3feb28dfd2e346befa5c90158d87b9d63f5a38c)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h>
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>
3bf7f4e0aSPeter Brune 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_Basic(SNESLineSearch linesearch)
5d71ae5a4SJacob Faibussowitsch {
6bf7f4e0aSPeter Brune   PetscBool changed_y, changed_w;
76a388c36SPeter Brune   Vec       X, F, Y, W;
86a388c36SPeter Brune   SNES      snes;
9*e3feb28dSStefano Zampini   PetscReal gnorm, xnorm, ynorm, lambda, fnorm = 0.0;
10bf7f4e0aSPeter Brune 
11bf7f4e0aSPeter Brune   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL));
139566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
149566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
159566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
169566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
17bf7f4e0aSPeter Brune 
18bf7f4e0aSPeter Brune   /* precheck */
199566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
20bf7f4e0aSPeter Brune 
21bf7f4e0aSPeter Brune   /* update */
229566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(W, -lambda, Y, X));
231baa6e33SBarry Smith   if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
24bf7f4e0aSPeter Brune 
25bf7f4e0aSPeter Brune   /* postcheck */
269566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
276b095a85SStefano Zampini   if (changed_y) {
286b095a85SStefano Zampini     if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X));
291baa6e33SBarry Smith     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
30bf7f4e0aSPeter Brune   }
315fc09b23SPeter Brune   if (linesearch->norms || snes->iter < snes->max_its - 1) {
329566063dSJacob Faibussowitsch     PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
33a286dac1SStefano Zampini     PetscCall(VecNorm(F, NORM_2, &fnorm));
345fc09b23SPeter Brune   }
355fc09b23SPeter Brune   if (linesearch->norms) {
369566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(Y, NORM_2, &linesearch->ynorm));
379566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(W, NORM_2, &linesearch->xnorm));
389566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(Y, NORM_2, &linesearch->ynorm));
399566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(W, NORM_2, &linesearch->xnorm));
405fc09b23SPeter Brune 
419bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
429bd66eb0SPeter Brune       linesearch->fnorm = gnorm;
43f5af7f23SKarl Rupp 
449566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm));
45a286dac1SStefano Zampini     } else linesearch->fnorm = fnorm;
465fc09b23SPeter Brune   }
47*e3feb28dSStefano Zampini   if (PetscIsInfOrNanReal(fnorm)) {
48*e3feb28dSStefano Zampini     PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_DOMAIN));
49*e3feb28dSStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
50*e3feb28dSStefano Zampini   }
516a388c36SPeter Brune 
52bf7f4e0aSPeter Brune   /* copy the solution over */
539566063dSJacob Faibussowitsch   PetscCall(VecCopy(W, X));
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55bf7f4e0aSPeter Brune }
56bf7f4e0aSPeter Brune 
57954494b2SJed Brown /*MC
581a4f838cSPeter Brune    SNESLINESEARCHBASIC - This line search implementation is not a line
59a99ef635SJonas Heinzmann    search at all; it simply uses the full step $x_{k+1} = x_k - \lambda Y_k$ with $\lambda=1$.
60a99ef635SJonas Heinzmann    Alternatively, $\lambda$ can be configured to be a constant damping factor by setting `snes_linesearch_damping`.
61a99ef635SJonas Heinzmann    Thus, this routine is intended for methods with well-scaled updates; i.e. Newton's method (`SNESNEWTONLS`), on
62a99ef635SJonas Heinzmann    well-behaved problems. Also named `SNESLINESEARCHNONE`.
63cd7522eaSPeter Brune 
64cd7522eaSPeter Brune    Options Database Keys:
65a99ef635SJonas Heinzmann +  -snes_linesearch_damping <1.0>     - step length is scaled by this factor
66a99ef635SJonas Heinzmann -  -snes_linesearch_norms <true>      - whether to compute norms or not (`SNESLineSearchSetComputeNorms()`)
67cd7522eaSPeter Brune 
68f6dfbefdSBarry Smith    Note:
69f6dfbefdSBarry Smith    For methods with ill-scaled updates (`SNESNRICHARDSON`, `SNESNCG`), a small
7089dfbbd5SBarry Smith    damping parameter may yield satisfactory, but slow convergence, despite
71f6dfbefdSBarry Smith    the lack of the line search.
72954494b2SJed Brown 
73954494b2SJed Brown    Level: advanced
74954494b2SJed Brown 
75420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLineSearchSetDamping()`, `SNESLineSearchSetComputeNorms()`
76954494b2SJed Brown M*/
77d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch)
78d71ae5a4SJacob Faibussowitsch {
79bf7f4e0aSPeter Brune   PetscFunctionBegin;
80f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_Basic;
810298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
820298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
830298fd71SBarry Smith   linesearch->ops->reset          = NULL;
840298fd71SBarry Smith   linesearch->ops->view           = NULL;
850298fd71SBarry Smith   linesearch->ops->setup          = NULL;
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87bf7f4e0aSPeter Brune }
88