1af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h>
2a7e14dcfSSatish Balay
TaoLineSearchView_Unit(TaoLineSearch ls,PetscViewer viewer)3d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchView_Unit(TaoLineSearch ls, PetscViewer viewer)
4d71ae5a4SJacob Faibussowitsch {
5a7e14dcfSSatish Balay PetscBool isascii;
6a7e14dcfSSatish Balay
7050fc7a3SBarry Smith PetscFunctionBegin;
89566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
963a3b9bcSJacob Faibussowitsch if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, " Line Search: Unit Step %g.\n", (double)ls->initstep));
103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11a7e14dcfSSatish Balay }
12a7e14dcfSSatish Balay
13a39c8e28SStefano Zampini /* Take unit step (newx = startx + initstep*step_direction) */
TaoLineSearchApply_Unit(TaoLineSearch ls,Vec x,PetscReal * f,Vec g,Vec step_direction)14d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchApply_Unit(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec step_direction)
15d71ae5a4SJacob Faibussowitsch {
16a7e14dcfSSatish Balay PetscFunctionBegin;
179566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
18a39c8e28SStefano Zampini ls->step = ls->initstep;
19a39c8e28SStefano Zampini PetscCall(VecAXPY(x, ls->step, step_direction));
20a39c8e28SStefano Zampini PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, x, f, g));
21a39c8e28SStefano Zampini PetscCall(TaoLineSearchMonitor(ls, 1, *f, ls->step));
22a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_SUCCESS;
233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
24a7e14dcfSSatish Balay }
25a7e14dcfSSatish Balay
2690b6438dSAlp Dener /*MC
2790b6438dSAlp Dener TAOLINESEARCHUNIT - Line-search type that disables line search and accepts the unit step length every time
28a7e14dcfSSatish Balay
29a39c8e28SStefano Zampini Options Database Keys:
30a39c8e28SStefano Zampini . -tao_ls_stepinit <step> - steplength
31a39c8e28SStefano Zampini
3290b6438dSAlp Dener Level: developer
33a7e14dcfSSatish Balay
34*f8d70eaaSPierre Jolivet .seealso: `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`
3590b6438dSAlp Dener M*/
TaoLineSearchCreate_Unit(TaoLineSearch ls)36d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Unit(TaoLineSearch ls)
37d71ae5a4SJacob Faibussowitsch {
38a7e14dcfSSatish Balay PetscFunctionBegin;
3983c8fe1dSLisandro Dalcin ls->ops->setup = NULL;
4083c8fe1dSLisandro Dalcin ls->ops->reset = NULL;
414664e3ffSStefano Zampini ls->ops->monitor = NULL;
42a7e14dcfSSatish Balay ls->ops->apply = TaoLineSearchApply_Unit;
43a7e14dcfSSatish Balay ls->ops->view = TaoLineSearchView_Unit;
443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
45a7e14dcfSSatish Balay }
46