Actual source code: unit.c
petsc-3.9.4 2018-09-11
2: #include <petsc/private/taolinesearchimpl.h>
4: static PetscErrorCode TaoLineSearchDestroy_Unit(TaoLineSearch ls)
5: {
8: PetscFree(ls->data);
9: return(0);
10: }
12: static PetscErrorCode TaoLineSearchSetFromOptions_Unit(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
13: {
16: PetscOptionsHead(PetscOptionsObject,"No Unit line search options");
17: PetscOptionsTail();
18: return(0);
19: }
21: static PetscErrorCode TaoLineSearchView_Unit(TaoLineSearch ls,PetscViewer viewer)
22: {
24: PetscBool isascii;
27: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
28: if (isascii) {
29: ierr=PetscViewerASCIIPrintf(viewer," Line Search: Unit Step.\n");
30: }
31: return(0);
32: }
34: static PetscErrorCode TaoLineSearchApply_Unit(TaoLineSearch ls,Vec x,PetscReal *f,Vec g,Vec step_direction)
35: {
37: PetscReal ftry;
38: PetscReal startf = *f;
41: /* Take unit step (newx = startx + 1.0*step_direction) */
42: VecAXPY(x,1.0,step_direction);
43: TaoLineSearchComputeObjectiveAndGradient(ls,x,&ftry,g);
44: PetscInfo1(ls,"Tao Apply Unit Step: %4.4e\n",1.0);
45: if (startf < ftry){
46: PetscInfo2(ls,"Tao Apply Unit Step, FINCREASE: F old:= %12.10e, F new: %12.10e\n",(double)startf,(double)ftry);
47: }
48: *f = ftry;
49: ls->step = 1.0;
50: ls->reason=TAOLINESEARCH_SUCCESS;
51: return(0);
52: }
54: /*@C
55: TaoCreateUnitLineSearch - Always use step length of 1.0
57: Input Parameters:
58: . tao - Tao context
60: Level: advanced
62: .keywords: Tao, linesearch
63: @*/
64: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Unit(TaoLineSearch ls)
65: {
67: ls->ops->setup = 0;
68: ls->ops->reset = 0;
69: ls->ops->apply = TaoLineSearchApply_Unit;
70: ls->ops->view = TaoLineSearchView_Unit;
71: ls->ops->destroy = TaoLineSearchDestroy_Unit;
72: ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Unit;
73: return(0);
74: }