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;
10: PetscBool domainerror;
12: PetscFunctionBegin;
13: PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL));
14: PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
15: PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
16: PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
17: PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
19: /* precheck */
20: PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
22: /* update */
23: PetscCall(VecWAXPY(W, -lambda, Y, X));
24: if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
26: /* postcheck */
27: PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
28: if (changed_y && !changed_w) {
29: PetscCall(VecWAXPY(W, -lambda, Y, X));
30: if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
31: }
32: if (linesearch->norms || snes->iter < snes->max_its - 1) {
33: PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
34: PetscCall(SNESGetFunctionDomainError(snes, &domainerror));
35: if (domainerror) {
36: PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_DOMAIN));
37: PetscFunctionReturn(PETSC_SUCCESS);
38: }
39: }
41: if (linesearch->norms) {
42: if (!linesearch->ops->vinorm) PetscCall(VecNormBegin(F, NORM_2, &linesearch->fnorm));
43: PetscCall(VecNormBegin(Y, NORM_2, &linesearch->ynorm));
44: PetscCall(VecNormBegin(W, NORM_2, &linesearch->xnorm));
45: if (!linesearch->ops->vinorm) PetscCall(VecNormEnd(F, NORM_2, &linesearch->fnorm));
46: PetscCall(VecNormEnd(Y, NORM_2, &linesearch->ynorm));
47: PetscCall(VecNormEnd(W, NORM_2, &linesearch->xnorm));
49: if (linesearch->ops->vinorm) {
50: linesearch->fnorm = gnorm;
52: PetscCall((*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm));
53: } else {
54: PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm));
55: }
56: }
58: /* copy the solution over */
59: PetscCall(VecCopy(W, X));
60: PetscFunctionReturn(PETSC_SUCCESS);
61: }
63: /*MC
64: SNESLINESEARCHBASIC - This line search implementation is not a line
65: search at all; it simply uses the full step. Thus, this routine is intended
66: for methods with well-scaled updates; i.e. Newton's method (`SNESNEWTONLS`), on
67: well-behaved problems. Also named as `SNESLINESEARCHNONE`
69: Options Database Keys:
70: + -snes_linesearch_damping <damping> - search vector is scaled by this amount, default is 1.0
71: - -snes_linesearch_norms <flag> - whether to compute norms or not, default is true (SNESLineSearchSetComputeNorms())
73: Note:
74: For methods with ill-scaled updates (`SNESNRICHARDSON`, `SNESNCG`), a small
75: damping parameter may yield satisfactory but slow convergence despite
76: the lack of the line search.
78: Level: advanced
80: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLineSearchSetDamping()`, `SNESLineSearchSetComputeNorms()`
81: M*/
82: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch)
83: {
84: PetscFunctionBegin;
85: linesearch->ops->apply = SNESLineSearchApply_Basic;
86: linesearch->ops->destroy = NULL;
87: linesearch->ops->setfromoptions = NULL;
88: linesearch->ops->reset = NULL;
89: linesearch->ops->view = NULL;
90: linesearch->ops->setup = NULL;
91: PetscFunctionReturn(PETSC_SUCCESS);
92: }