xref: /petsc/src/tao/linesearch/impls/unit/unit.c (revision df4cd43f92eaa320656440c40edb1046daee8f75)
1 
2 #include <petsc/private/taolinesearchimpl.h>
3 
4 static PetscErrorCode TaoLineSearchDestroy_Unit(TaoLineSearch ls)
5 {
6   PetscFunctionBegin;
7   PetscFunctionReturn(PETSC_SUCCESS);
8 }
9 
10 static PetscErrorCode TaoLineSearchSetFromOptions_Unit(TaoLineSearch ls, PetscOptionItems *PetscOptionsObject)
11 {
12   PetscFunctionBegin;
13   PetscFunctionReturn(PETSC_SUCCESS);
14 }
15 
16 static PetscErrorCode TaoLineSearchView_Unit(TaoLineSearch ls, PetscViewer viewer)
17 {
18   PetscBool isascii;
19 
20   PetscFunctionBegin;
21   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
22   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Line Search: Unit Step %g.\n", (double)ls->initstep));
23   PetscFunctionReturn(PETSC_SUCCESS);
24 }
25 
26 /* Take unit step (newx = startx + initstep*step_direction) */
27 static PetscErrorCode TaoLineSearchApply_Unit(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec step_direction)
28 {
29   PetscFunctionBegin;
30   PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
31   ls->step = ls->initstep;
32   PetscCall(VecAXPY(x, ls->step, step_direction));
33   PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, x, f, g));
34   PetscCall(TaoLineSearchMonitor(ls, 1, *f, ls->step));
35   ls->reason = TAOLINESEARCH_SUCCESS;
36   PetscFunctionReturn(PETSC_SUCCESS);
37 }
38 
39 /*MC
40    TAOLINESEARCHUNIT - Line-search type that disables line search and accepts the unit step length every time
41 
42   Options Database Keys:
43 . -tao_ls_stepinit <step> - steplength
44 
45    Level: developer
46 
47 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`
48 
49 .keywords: Tao, linesearch
50 M*/
51 PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Unit(TaoLineSearch ls)
52 {
53   PetscFunctionBegin;
54   ls->ops->setup          = NULL;
55   ls->ops->reset          = NULL;
56   ls->ops->monitor        = NULL;
57   ls->ops->apply          = TaoLineSearchApply_Unit;
58   ls->ops->view           = TaoLineSearchView_Unit;
59   ls->ops->destroy        = TaoLineSearchDestroy_Unit;
60   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Unit;
61   PetscFunctionReturn(PETSC_SUCCESS);
62 }
63