xref: /petsc/src/tao/linesearch/impls/unit/unit.c (revision f8d70eaab67b53d355abe83957ffb5a3120bcf52)
1af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h>
2a7e14dcfSSatish Balay 
3d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchDestroy_Unit(TaoLineSearch ls)
4d71ae5a4SJacob Faibussowitsch {
5a7e14dcfSSatish Balay   PetscFunctionBegin;
63ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7a7e14dcfSSatish Balay }
8a7e14dcfSSatish Balay 
9d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchSetFromOptions_Unit(TaoLineSearch ls, PetscOptionItems *PetscOptionsObject)
10d71ae5a4SJacob Faibussowitsch {
11a7e14dcfSSatish Balay   PetscFunctionBegin;
123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13a7e14dcfSSatish Balay }
14a7e14dcfSSatish Balay 
15d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchView_Unit(TaoLineSearch ls, PetscViewer viewer)
16d71ae5a4SJacob Faibussowitsch {
17a7e14dcfSSatish Balay   PetscBool isascii;
18a7e14dcfSSatish Balay 
19050fc7a3SBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2163a3b9bcSJacob Faibussowitsch   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Line Search: Unit Step %g.\n", (double)ls->initstep));
223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23a7e14dcfSSatish Balay }
24a7e14dcfSSatish Balay 
25a39c8e28SStefano Zampini /* Take unit step (newx = startx + initstep*step_direction) */
26d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchApply_Unit(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec step_direction)
27d71ae5a4SJacob Faibussowitsch {
28a7e14dcfSSatish Balay   PetscFunctionBegin;
299566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
30a39c8e28SStefano Zampini   ls->step = ls->initstep;
31a39c8e28SStefano Zampini   PetscCall(VecAXPY(x, ls->step, step_direction));
32a39c8e28SStefano Zampini   PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, x, f, g));
33a39c8e28SStefano Zampini   PetscCall(TaoLineSearchMonitor(ls, 1, *f, ls->step));
34a7e14dcfSSatish Balay   ls->reason = TAOLINESEARCH_SUCCESS;
353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36a7e14dcfSSatish Balay }
37a7e14dcfSSatish Balay 
3890b6438dSAlp Dener /*MC
3990b6438dSAlp Dener    TAOLINESEARCHUNIT - Line-search type that disables line search and accepts the unit step length every time
40a7e14dcfSSatish Balay 
41a39c8e28SStefano Zampini   Options Database Keys:
42a39c8e28SStefano Zampini . -tao_ls_stepinit <step> - steplength
43a39c8e28SStefano Zampini 
4490b6438dSAlp Dener    Level: developer
45a7e14dcfSSatish Balay 
46*f8d70eaaSPierre Jolivet .seealso: `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`
4790b6438dSAlp Dener M*/
48d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Unit(TaoLineSearch ls)
49d71ae5a4SJacob Faibussowitsch {
50a7e14dcfSSatish Balay   PetscFunctionBegin;
5183c8fe1dSLisandro Dalcin   ls->ops->setup          = NULL;
5283c8fe1dSLisandro Dalcin   ls->ops->reset          = NULL;
534664e3ffSStefano Zampini   ls->ops->monitor        = NULL;
54a7e14dcfSSatish Balay   ls->ops->apply          = TaoLineSearchApply_Unit;
55a7e14dcfSSatish Balay   ls->ops->view           = TaoLineSearchView_Unit;
56a7e14dcfSSatish Balay   ls->ops->destroy        = TaoLineSearchDestroy_Unit;
57a7e14dcfSSatish Balay   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Unit;
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59a7e14dcfSSatish Balay }
60