Actual source code: unit.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  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:   TaoLineSearchMonitor(ls, 0, *f, 0.0);
 43:   VecAXPY(x,1.0,step_direction);
 44:   TaoLineSearchComputeObjectiveAndGradient(ls,x,&ftry,g);
 45:   TaoLineSearchMonitor(ls, 1, *f, 1.0);
 46:   PetscInfo1(ls,"Tao Apply Unit Step: %4.4e\n",1.0);
 47:   if (startf < ftry){
 48:     PetscInfo2(ls,"Tao Apply Unit Step, FINCREASE: F old:= %12.10e, F new: %12.10e\n",(double)startf,(double)ftry);
 49:   }
 50:   *f = ftry;
 51:   ls->step = 1.0;
 52:   ls->reason=TAOLINESEARCH_SUCCESS;
 53:   return(0);
 54: }

 56: /*@C
 57:    TaoCreateUnitLineSearch - Always use step length of 1.0

 59:    Input Parameters:
 60: .  tao - Tao context

 62:    Level: advanced

 64: .keywords: Tao, linesearch
 65: @*/
 66: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_Unit(TaoLineSearch ls)
 67: {
 69:   ls->ops->setup = 0;
 70:   ls->ops->reset = 0;
 71:   ls->ops->apply = TaoLineSearchApply_Unit;
 72:   ls->ops->view = TaoLineSearchView_Unit;
 73:   ls->ops->destroy = TaoLineSearchDestroy_Unit;
 74:   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_Unit;
 75:   return(0);
 76: }