xref: /petsc/src/tao/linesearch/impls/unit/unit.c (revision 6c23d07575f07c9a033c9cb13c1a90726fd2593a)
1 #include "petscvec.h"
2 #include "taosolver.h"
3 #include "tao-private/taolinesearch_impl.h"
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "TaoLineSearchDestroy_Unit"
7 static PetscErrorCode TaoLineSearchDestroy_Unit(TaoLineSearch ls)
8 {
9   PetscErrorCode ierr;
10   PetscFunctionBegin;
11   ierr = PetscFree(ls->data); CHKERRQ(ierr);
12   PetscFunctionReturn(0);
13 }
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "TaoLineSearchSetFromOptions_Unit"
17 static PetscErrorCode TaoLineSearchSetFromOptions_Unit(TaoLineSearch ls)
18 {
19   PetscErrorCode ierr;
20   PetscFunctionBegin;
21   ierr = PetscOptionsHead("No Unit line search options");CHKERRQ(ierr);
22   ierr = PetscOptionsTail();CHKERRQ(ierr);
23   PetscFunctionReturn(0);
24 }
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "TaoLineSearchView_Unit"
28 static PetscErrorCode TaoLineSearchView_Unit(TaoLineSearch ls,PetscViewer viewer)
29 {
30 
31   PetscErrorCode ierr;
32   PetscBool isascii;
33   PetscFunctionBegin;
34 
35   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii); CHKERRQ(ierr);
36   if (isascii) {
37       ierr=PetscViewerASCIIPrintf(viewer,"  Line Search: Unit Step.\n");CHKERRQ(ierr);
38   }
39   PetscFunctionReturn(0);
40 }
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "TaoLineSearchApply_Unit"
44 static PetscErrorCode TaoLineSearchApply_Unit(TaoLineSearch ls,Vec x,PetscReal *f,Vec g,Vec step_direction)
45 {
46   PetscErrorCode   ierr;
47   PetscReal ftry;
48   PetscReal startf = *f;
49 
50   PetscFunctionBegin;
51 
52   /* Take unit step (newx = startx + 1.0*step_direction) */
53   ierr = VecAXPY(x,1.0,step_direction);CHKERRQ(ierr);
54 
55   ierr = TaoLineSearchComputeObjectiveAndGradient(ls,x,&ftry,g); CHKERRQ(ierr);
56   ierr = PetscInfo1(ls,"Tao Apply Unit Step: %4.4e\n",1.0);
57          CHKERRQ(ierr);
58   if (startf < ftry){
59     ierr = PetscInfo2(ls,"Tao Apply Unit Step, FINCREASE: F old:= %12.10e, F new: %12.10e\n",startf,ftry); CHKERRQ(ierr);
60   }
61   *f = ftry;
62   ls->step = 1.0;
63   ls->reason=TAOLINESEARCH_SUCCESS;
64   PetscFunctionReturn(0);
65 }
66 
67 EXTERN_C_BEGIN
68 #undef __FUNCT__
69 #define __FUNCT__ "TaoLineSearchCreate_Unit"
70 /*@C
71    TaoCreateUnitLineSearch - Always use step length of 1.0
72 
73    Input Parameters:
74 .  tao - TaoSolver context
75 
76 
77    Level: advanced
78 
79 .keywords: TaoSolver, linesearch
80 @*/
81 PetscErrorCode TaoLineSearchCreate_Unit(TaoLineSearch ls)
82 {
83 
84   PetscFunctionBegin;
85   ls->ops->setup = 0;
86   ls->ops->apply = TaoLineSearchApply_Unit;
87   ls->ops->view = TaoLineSearchView_Unit;
88   ls->ops->destroy = TaoLineSearchDestroy_Unit;
89   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Unit;
90 
91   PetscFunctionReturn(0);
92 }
93 EXTERN_C_END
94 
95