xref: /petsc/src/tao/linesearch/impls/unit/unit.c (revision 050fc7a34c6d42d9b367a0bc7bbb4624be070c15)
1 #include <petscvec.h>
2 #include <taosolver.h>
3 #include <petsc-private/taolinesearchimpl.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   PetscErrorCode ierr;
31   PetscBool      isascii;
32 
33   PetscFunctionBegin;
34   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
35   if (isascii) {
36     ierr=PetscViewerASCIIPrintf(viewer,"  Line Search: Unit Step.\n");CHKERRQ(ierr);
37   }
38   PetscFunctionReturn(0);
39 }
40 
41 #undef __FUNCT__
42 #define __FUNCT__ "TaoLineSearchApply_Unit"
43 static PetscErrorCode TaoLineSearchApply_Unit(TaoLineSearch ls,Vec x,PetscReal *f,Vec g,Vec step_direction)
44 {
45   PetscErrorCode ierr;
46   PetscReal      ftry;
47   PetscReal      startf = *f;
48 
49   PetscFunctionBegin;
50   /* Take unit step (newx = startx + 1.0*step_direction) */
51   ierr = VecAXPY(x,1.0,step_direction);CHKERRQ(ierr);
52   ierr = TaoLineSearchComputeObjectiveAndGradient(ls,x,&ftry,g);CHKERRQ(ierr);
53   ierr = PetscInfo1(ls,"Tao Apply Unit Step: %4.4e\n",1.0);CHKERRQ(ierr);
54   if (startf < ftry){
55     ierr = PetscInfo2(ls,"Tao Apply Unit Step, FINCREASE: F old:= %12.10e, F new: %12.10e\n",startf,ftry);CHKERRQ(ierr);
56   }
57   *f = ftry;
58   ls->step = 1.0;
59   ls->reason=TAOLINESEARCH_SUCCESS;
60   PetscFunctionReturn(0);
61 }
62 
63 EXTERN_C_BEGIN
64 #undef __FUNCT__
65 #define __FUNCT__ "TaoLineSearchCreate_Unit"
66 /*@C
67    TaoCreateUnitLineSearch - Always use step length of 1.0
68 
69    Input Parameters:
70 .  tao - TaoSolver context
71 
72    Level: advanced
73 
74 .keywords: TaoSolver, linesearch
75 @*/
76 PetscErrorCode TaoLineSearchCreate_Unit(TaoLineSearch ls)
77 {
78   PetscFunctionBegin;
79   ls->ops->setup = 0;
80   ls->ops->apply = TaoLineSearchApply_Unit;
81   ls->ops->view = TaoLineSearchView_Unit;
82   ls->ops->destroy = TaoLineSearchDestroy_Unit;
83   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Unit;
84   PetscFunctionReturn(0);
85 }
86 EXTERN_C_END
87 
88