Actual source code: linesearchbasic.c
petsc-3.5.4 2015-05-23
1: #include <petsc-private/linesearchimpl.h>
2: #include <petsc-private/snesimpl.h>
6: static PetscErrorCode SNESLineSearchApply_Basic(SNESLineSearch linesearch)
7: {
8: PetscBool changed_y, changed_w;
10: Vec X, F, Y, W;
11: SNES snes;
12: PetscReal gnorm, xnorm, ynorm, lambda;
13: PetscBool domainerror;
16: SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL);
17: SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);
18: SNESLineSearchGetLambda(linesearch, &lambda);
19: SNESLineSearchGetSNES(linesearch, &snes);
20: SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);
22: /* precheck */
23: SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);
25: /* update */
26: VecWAXPY(W,-lambda,Y,X);
27: if (linesearch->ops->viproject) {
28: (*linesearch->ops->viproject)(snes, W);
29: }
31: /* postcheck */
32: SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);
33: if (changed_y) {
34: VecWAXPY(W,-lambda,Y,X);
35: if (linesearch->ops->viproject) {
36: (*linesearch->ops->viproject)(snes, W);
37: }
38: }
39: if (linesearch->norms || snes->iter < snes->max_its-1) {
40: (*linesearch->ops->snesfunc)(snes,W,F);
41: SNESGetFunctionDomainError(snes, &domainerror);
42: if (domainerror) {
43: SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
44: return(0);
45: }
46: }
48: if (linesearch->norms) {
49: if (!linesearch->ops->vinorm) VecNormBegin(F, NORM_2, &linesearch->fnorm);
50: VecNormBegin(Y, NORM_2, &linesearch->ynorm);
51: VecNormBegin(W, NORM_2, &linesearch->xnorm);
52: if (!linesearch->ops->vinorm) VecNormEnd(F, NORM_2, &linesearch->fnorm);
53: VecNormEnd(Y, NORM_2, &linesearch->ynorm);
54: VecNormEnd(W, NORM_2, &linesearch->xnorm);
56: if (linesearch->ops->vinorm) {
57: linesearch->fnorm = gnorm;
59: (*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm);
60: } else {
61: VecNorm(F,NORM_2,&linesearch->fnorm);
62: }
63: }
65: /* copy the solution over */
66: VecCopy(W, X);
67: return(0);
68: }
72: /*MC
73: SNESLINESEARCHBASIC - This line search implementation is not a line
74: search at all; it simply uses the full step. Thus, this routine is intended
75: for methods with well-scaled updates; i.e. Newton's method (SNESNEWTONLS), on
76: well-behaved problems.
78: Options Database Keys:
79: + -snes_linesearch_damping <damping> search vector is scaled by this amount, default is 1.0
80: - -snes_linesearch_norms <flag> whether to compute norms or not, default is true
82: Notes:
83: For methods with ill-scaled updates (SNESNRICHARDSON, SNESNCG), a small
84: damping parameter may yield satisfactory but slow convergence despite
85: the simplicity of the line search.
87: Level: advanced
89: .keywords: SNES, SNESLineSearch, damping
91: .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
92: M*/
93: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch)
94: {
96: linesearch->ops->apply = SNESLineSearchApply_Basic;
97: linesearch->ops->destroy = NULL;
98: linesearch->ops->setfromoptions = NULL;
99: linesearch->ops->reset = NULL;
100: linesearch->ops->view = NULL;
101: linesearch->ops->setup = NULL;
102: return(0);
103: }