Actual source code: unit.c
petsc-3.7.3 2016-08-01
2: #include <petsc/private/taolinesearchimpl.h>
6: static PetscErrorCode TaoLineSearchDestroy_Unit(TaoLineSearch ls)
7: {
10: PetscFree(ls->data);
11: return(0);
12: }
16: static PetscErrorCode TaoLineSearchSetFromOptions_Unit(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
17: {
20: PetscOptionsHead(PetscOptionsObject,"No Unit line search options");
21: PetscOptionsTail();
22: return(0);
23: }
27: static PetscErrorCode TaoLineSearchView_Unit(TaoLineSearch ls,PetscViewer viewer)
28: {
30: PetscBool isascii;
33: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
34: if (isascii) {
35: ierr=PetscViewerASCIIPrintf(viewer," Line Search: Unit Step.\n");
36: }
37: return(0);
38: }
42: static PetscErrorCode TaoLineSearchApply_Unit(TaoLineSearch ls,Vec x,PetscReal *f,Vec g,Vec step_direction)
43: {
45: PetscReal ftry;
46: PetscReal startf = *f;
49: /* Take unit step (newx = startx + 1.0*step_direction) */
50: VecAXPY(x,1.0,step_direction);
51: TaoLineSearchComputeObjectiveAndGradient(ls,x,&ftry,g);
52: PetscInfo1(ls,"Tao Apply Unit Step: %4.4e\n",1.0);
53: if (startf < ftry){
54: PetscInfo2(ls,"Tao Apply Unit Step, FINCREASE: F old:= %12.10e, F new: %12.10e\n",(double)startf,(double)ftry);
55: }
56: *f = ftry;
57: ls->step = 1.0;
58: ls->reason=TAOLINESEARCH_SUCCESS;
59: return(0);
60: }
64: /*@C
65: TaoCreateUnitLineSearch - Always use step length of 1.0
67: Input Parameters:
68: . tao - Tao context
70: Level: advanced
72: .keywords: Tao, linesearch
73: @*/
74: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Unit(TaoLineSearch ls)
75: {
77: ls->ops->setup = 0;
78: ls->ops->reset = 0;
79: ls->ops->apply = TaoLineSearchApply_Unit;
80: ls->ops->view = TaoLineSearchView_Unit;
81: ls->ops->destroy = TaoLineSearchDestroy_Unit;
82: ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Unit;
83: return(0);
84: }