1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h>
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>
3bf7f4e0aSPeter Brune
SNESLineSearchApply_Basic(SNESLineSearch linesearch)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;
9e3feb28dSStefano 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));
34*76c63389SBarry Smith SNESLineSearchCheckFunctionDomainError(snes, linesearch, fnorm);
355fc09b23SPeter Brune }
365fc09b23SPeter Brune if (linesearch->norms) {
379566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &linesearch->ynorm));
389566063dSJacob Faibussowitsch PetscCall(VecNormBegin(W, NORM_2, &linesearch->xnorm));
399566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &linesearch->ynorm));
409566063dSJacob Faibussowitsch PetscCall(VecNormEnd(W, NORM_2, &linesearch->xnorm));
415fc09b23SPeter Brune
429bd66eb0SPeter Brune if (linesearch->ops->vinorm) {
439bd66eb0SPeter Brune linesearch->fnorm = gnorm;
44f5af7f23SKarl Rupp
459566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm));
46a286dac1SStefano Zampini } else linesearch->fnorm = fnorm;
475fc09b23SPeter Brune }
486a388c36SPeter Brune
49bf7f4e0aSPeter Brune /* copy the solution over */
509566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X));
513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
52bf7f4e0aSPeter Brune }
53bf7f4e0aSPeter Brune
54954494b2SJed Brown /*MC
551a4f838cSPeter Brune SNESLINESEARCHBASIC - This line search implementation is not a line
56a99ef635SJonas Heinzmann search at all; it simply uses the full step $x_{k+1} = x_k - \lambda Y_k$ with $\lambda=1$.
57a99ef635SJonas Heinzmann Alternatively, $\lambda$ can be configured to be a constant damping factor by setting `snes_linesearch_damping`.
58a99ef635SJonas Heinzmann Thus, this routine is intended for methods with well-scaled updates; i.e. Newton's method (`SNESNEWTONLS`), on
59a99ef635SJonas Heinzmann well-behaved problems. Also named `SNESLINESEARCHNONE`.
60cd7522eaSPeter Brune
61cd7522eaSPeter Brune Options Database Keys:
62a99ef635SJonas Heinzmann + -snes_linesearch_damping <1.0> - step length is scaled by this factor
63a99ef635SJonas Heinzmann - -snes_linesearch_norms <true> - whether to compute norms or not (`SNESLineSearchSetComputeNorms()`)
64cd7522eaSPeter Brune
65f6dfbefdSBarry Smith Note:
66f6dfbefdSBarry Smith For methods with ill-scaled updates (`SNESNRICHARDSON`, `SNESNCG`), a small
6789dfbbd5SBarry Smith damping parameter may yield satisfactory, but slow convergence, despite
68f6dfbefdSBarry Smith the lack of the line search.
69954494b2SJed Brown
70954494b2SJed Brown Level: advanced
71954494b2SJed Brown
72420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLineSearchSetDamping()`, `SNESLineSearchSetComputeNorms()`
73954494b2SJed Brown M*/
SNESLineSearchCreate_Basic(SNESLineSearch linesearch)74d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch)
75d71ae5a4SJacob Faibussowitsch {
76bf7f4e0aSPeter Brune PetscFunctionBegin;
77f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_Basic;
780298fd71SBarry Smith linesearch->ops->destroy = NULL;
790298fd71SBarry Smith linesearch->ops->setfromoptions = NULL;
800298fd71SBarry Smith linesearch->ops->reset = NULL;
810298fd71SBarry Smith linesearch->ops->view = NULL;
820298fd71SBarry Smith linesearch->ops->setup = NULL;
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
84bf7f4e0aSPeter Brune }
85