xref: /petsc/src/tao/linesearch/impls/unit/unit.c (revision 9566063d113dddea24716c546802770db7481bc0)
1f273b952SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchDestroy_Unit(TaoLineSearch ls)
5a7e14dcfSSatish Balay {
6a7e14dcfSSatish Balay   PetscFunctionBegin;
7*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(ls->data));
8a7e14dcfSSatish Balay   PetscFunctionReturn(0);
9a7e14dcfSSatish Balay }
10a7e14dcfSSatish Balay 
114416b707SBarry Smith static PetscErrorCode TaoLineSearchSetFromOptions_Unit(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
12a7e14dcfSSatish Balay {
13a7e14dcfSSatish Balay   PetscFunctionBegin;
14*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"No Unit line search options"));
15*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
16a7e14dcfSSatish Balay   PetscFunctionReturn(0);
17a7e14dcfSSatish Balay }
18a7e14dcfSSatish Balay 
19a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchView_Unit(TaoLineSearch ls,PetscViewer viewer)
20a7e14dcfSSatish Balay {
21a7e14dcfSSatish Balay   PetscBool      isascii;
22a7e14dcfSSatish Balay 
23050fc7a3SBarry Smith   PetscFunctionBegin;
24*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
25a7e14dcfSSatish Balay   if (isascii) {
26*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Line Search: Unit Step.\n"));
27a7e14dcfSSatish Balay   }
28a7e14dcfSSatish Balay   PetscFunctionReturn(0);
29a7e14dcfSSatish Balay }
30a7e14dcfSSatish Balay 
31a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchApply_Unit(TaoLineSearch ls,Vec x,PetscReal *f,Vec g,Vec step_direction)
32a7e14dcfSSatish Balay {
33a7e14dcfSSatish Balay   PetscReal      ftry;
34a7e14dcfSSatish Balay   PetscReal      startf = *f;
35a7e14dcfSSatish Balay 
36a7e14dcfSSatish Balay   PetscFunctionBegin;
37a7e14dcfSSatish Balay   /* Take unit step (newx = startx + 1.0*step_direction) */
38*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
39*9566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x,1.0,step_direction));
40*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls,x,&ftry,g));
41*9566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchMonitor(ls, 1, *f, 1.0));
42*9566063dSJacob Faibussowitsch   PetscCall(PetscInfo(ls,"Tao Apply Unit Step: %4.4e\n",1.0));
43a7e14dcfSSatish Balay   if (startf < ftry) {
44*9566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Tao Apply Unit Step, FINCREASE: F old:= %12.10e, F new: %12.10e\n",(double)startf,(double)ftry));
45a7e14dcfSSatish Balay   }
46a7e14dcfSSatish Balay   *f = ftry;
47a7e14dcfSSatish Balay   ls->step = 1.0;
48a7e14dcfSSatish Balay   ls->reason=TAOLINESEARCH_SUCCESS;
49a7e14dcfSSatish Balay   PetscFunctionReturn(0);
50a7e14dcfSSatish Balay }
51a7e14dcfSSatish Balay 
5290b6438dSAlp Dener /*MC
5390b6438dSAlp Dener    TAOLINESEARCHUNIT - Line-search type that disables line search and accepts the unit step length every time
54a7e14dcfSSatish Balay 
5590b6438dSAlp Dener    Level: developer
56a7e14dcfSSatish Balay 
5790b6438dSAlp Dener .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchApply()
58a7e14dcfSSatish Balay 
59441846f8SBarry Smith .keywords: Tao, linesearch
6090b6438dSAlp Dener M*/
61728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Unit(TaoLineSearch ls)
62a7e14dcfSSatish Balay {
63a7e14dcfSSatish Balay   PetscFunctionBegin;
6483c8fe1dSLisandro Dalcin   ls->ops->setup = NULL;
6583c8fe1dSLisandro Dalcin   ls->ops->reset = NULL;
664664e3ffSStefano Zampini   ls->ops->monitor = NULL;
67a7e14dcfSSatish Balay   ls->ops->apply = TaoLineSearchApply_Unit;
68a7e14dcfSSatish Balay   ls->ops->view = TaoLineSearchView_Unit;
69a7e14dcfSSatish Balay   ls->ops->destroy = TaoLineSearchDestroy_Unit;
70a7e14dcfSSatish Balay   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Unit;
71a7e14dcfSSatish Balay   PetscFunctionReturn(0);
72a7e14dcfSSatish Balay }
73