Actual source code: linesearchbasic.c
1: #include <petsc/private/linesearchimpl.h>
2: #include <petsc/private/snesimpl.h>
4: static PetscErrorCode SNESLineSearchApply_Basic(SNESLineSearch linesearch)
5: {
6: PetscBool changed_y, changed_w;
7: Vec X, F, Y, W;
8: SNES snes;
9: PetscReal gnorm, xnorm, ynorm, lambda, fnorm;
11: PetscFunctionBegin;
12: PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL));
13: PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
14: PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
15: PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
16: PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
18: /* precheck */
19: PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
21: /* update */
22: PetscCall(VecWAXPY(W, -lambda, Y, X));
23: if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
25: /* postcheck */
26: PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
27: if (changed_y) {
28: if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X));
29: if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
30: }
31: if (linesearch->norms || snes->iter < snes->max_its - 1) {
32: PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
33: PetscCall(VecNorm(F, NORM_2, &fnorm));
34: if (PetscIsInfOrNanReal(fnorm)) {
35: PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_DOMAIN));
36: PetscFunctionReturn(PETSC_SUCCESS);
37: }
38: }
39: if (linesearch->norms) {
40: PetscCall(VecNormBegin(Y, NORM_2, &linesearch->ynorm));
41: PetscCall(VecNormBegin(W, NORM_2, &linesearch->xnorm));
42: PetscCall(VecNormEnd(Y, NORM_2, &linesearch->ynorm));
43: PetscCall(VecNormEnd(W, NORM_2, &linesearch->xnorm));
45: if (linesearch->ops->vinorm) {
46: linesearch->fnorm = gnorm;
48: PetscCall((*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm));
49: } else linesearch->fnorm = fnorm;
50: }
52: /* copy the solution over */
53: PetscCall(VecCopy(W, X));
54: PetscFunctionReturn(PETSC_SUCCESS);
55: }
57: /*MC
58: SNESLINESEARCHBASIC - This line search implementation is not a line
59: search at all; it simply uses the full step $x_{k+1} = x_k - \lambda Y_k$ with $\lambda=1$.
60: Alternatively, $\lambda$ can be configured to be a constant damping factor by setting `snes_linesearch_damping`.
61: Thus, this routine is intended for methods with well-scaled updates; i.e. Newton's method (`SNESNEWTONLS`), on
62: well-behaved problems. Also named `SNESLINESEARCHNONE`.
64: Options Database Keys:
65: + -snes_linesearch_damping <1.0> - step length is scaled by this factor
66: - -snes_linesearch_norms <true> - whether to compute norms or not (`SNESLineSearchSetComputeNorms()`)
68: Note:
69: For methods with ill-scaled updates (`SNESNRICHARDSON`, `SNESNCG`), a small
70: damping parameter may yield satisfactory, but slow convergence, despite
71: the lack of the line search.
73: Level: advanced
75: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLineSearchSetDamping()`, `SNESLineSearchSetComputeNorms()`
76: M*/
77: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch)
78: {
79: PetscFunctionBegin;
80: linesearch->ops->apply = SNESLineSearchApply_Basic;
81: linesearch->ops->destroy = NULL;
82: linesearch->ops->setfromoptions = NULL;
83: linesearch->ops->reset = NULL;
84: linesearch->ops->view = NULL;
85: linesearch->ops->setup = NULL;
86: PetscFunctionReturn(PETSC_SUCCESS);
87: }