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