Actual source code: taolinesearch.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: #include <petsctaolinesearch.h> /*I "petsctaolinesearch.h" I*/
  2: #include <petsc-private/taolinesearchimpl.h>

  4: PetscBool TaoLineSearchInitialized = PETSC_FALSE;
  5: PetscFunctionList TaoLineSearchList = NULL;

  7: PetscClassId TAOLINESEARCH_CLASSID=0;
  8: PetscLogEvent TaoLineSearch_ApplyEvent = 0, TaoLineSearch_EvalEvent=0;

 12: /*@C
 13:   TaoLineSearchView - Prints information about the TaoLineSearch

 15:   Collective on TaoLineSearch

 17:   InputParameters:
 18: + ls - the Tao context
 19: - viewer - visualization context

 21:   Options Database Key:
 22: . -tao_ls_view - Calls TaoLineSearchView() at the end of each line search

 24:   Notes:
 25:   The available visualization contexts include
 26: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
 27: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
 28:          output where only the first processor opens
 29:          the file.  All other processors send their
 30:          data to the first processor to print.

 32:   Level: beginner

 34: .seealso: PetscViewerASCIIOpen()
 35: @*/

 37: PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
 38: {
 39:   PetscErrorCode          ierr;
 40:   PetscBool               isascii, isstring;
 41:   const TaoLineSearchType type;
 42:   PetscBool               destroyviewer = PETSC_FALSE;

 46:   if (!viewer) {
 47:     destroyviewer = PETSC_TRUE;
 48:     PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer);
 49:   }

 53:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 54:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
 55:   if (isascii) {
 56:     if (((PetscObject)ls)->prefix) {
 57:       PetscViewerASCIIPrintf(viewer,"TaoLineSearch Object:(%s)\n",((PetscObject)ls)->prefix);
 58:     } else {
 59:       PetscViewerASCIIPrintf(viewer,"TaoLineSearch Object:\n");
 60:     }
 61:     PetscViewerASCIIPushTab(viewer);
 62:     TaoLineSearchGetType(ls,&type);
 63:     if (type) {
 64:       PetscViewerASCIIPrintf(viewer,"type: %s\n",type);
 65:     } else {
 66:       PetscViewerASCIIPrintf(viewer,"type: not set yet\n");
 67:     }
 68:     if (ls->ops->view) {
 69:       PetscViewerASCIIPushTab(viewer);
 70:       (*ls->ops->view)(ls,viewer);
 71:       PetscViewerASCIIPopTab(viewer);
 72:     }
 73:     PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);
 74:     PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);
 75:     PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);
 76:     PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);
 77:     PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);

 79:     if (ls->bounded) {
 80:       PetscViewerASCIIPrintf(viewer,"using variable bounds\n");
 81:     }
 82:     PetscViewerASCIIPrintf(viewer,"Termination reason: %D\n",(int)ls->reason);
 83:     PetscViewerASCIIPopTab(viewer);

 85:   } else if (isstring) {
 86:     TaoLineSearchGetType(ls,&type);
 87:     PetscViewerStringSPrintf(viewer," %-3.3s",type);
 88:   }
 89:   if (destroyviewer) {
 90:     PetscViewerDestroy(&viewer);
 91:   }
 92:   return(0);
 93: }

 97: /*@C
 98:   TaoLineSearchCreate - Creates a TAO Line Search object.  Algorithms in TAO that use
 99:   line-searches will automatically create one.

101:   Collective on MPI_Comm

103:   Input Parameter:
104: . comm - MPI communicator

106:   Output Parameter:
107: . newls - the new TaoLineSearch context

109:   Available methods include:
110: + more-thuente
111: . gpcg
112: - unit - Do not perform any line search


115:    Options Database Keys:
116: .   -tao_ls_type - select which method TAO should use

118:    Level: beginner

120: .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy()
121: @*/

123: PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
124: {
126:   TaoLineSearch  ls;

130:   *newls = NULL;

132:  #ifndef PETSC_USE_DYNAMIC_LIBRARIES
133:   TaoLineSearchInitializePackage();
134:  #endif

136:   PetscHeaderCreate(ls,_p_TaoLineSearch,struct _TaoLineSearchOps,TAOLINESEARCH_CLASSID,"TaoLineSearch",0,0,comm,TaoLineSearchDestroy,TaoLineSearchView);
137:   ls->bounded = 0;
138:   ls->max_funcs=30;
139:   ls->ftol = 0.0001;
140:   ls->gtol = 0.9;
141: #if defined(PETSC_USE_REAL_SINGLE)
142:   ls->rtol = 1.0e-5;
143: #else
144:   ls->rtol = 1.0e-10;
145: #endif
146:   ls->stepmin=1.0e-20;
147:   ls->stepmax=1.0e+20;
148:   ls->step=1.0;
149:   ls->nfeval=0;
150:   ls->ngeval=0;
151:   ls->nfgeval=0;

153:   ls->ops->computeobjective=0;
154:   ls->ops->computegradient=0;
155:   ls->ops->computeobjectiveandgradient=0;
156:   ls->ops->computeobjectiveandgts=0;
157:   ls->ops->setup=0;
158:   ls->ops->apply=0;
159:   ls->ops->view=0;
160:   ls->ops->setfromoptions=0;
161:   ls->ops->reset=0;
162:   ls->ops->destroy=0;
163:   ls->setupcalled=PETSC_FALSE;
164:   ls->usetaoroutines=PETSC_FALSE;
165:   *newls = ls;
166:   return(0);
167: }

171: /*@
172:   TaoLineSearchSetUp - Sets up the internal data structures for the later use
173:   of a Tao solver

175:   Collective on ls

177:   Input Parameters:
178: . ls - the TaoLineSearch context

180:   Notes:
181:   The user will not need to explicitly call TaoLineSearchSetUp(), as it will
182:   automatically be called in TaoLineSearchSolve().  However, if the user
183:   desires to call it explicitly, it should come after TaoLineSearchCreate()
184:   but before TaoLineSearchApply().

186:   Level: developer

188: .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
189: @*/

191: PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
192: {
194:   const char     *default_type=TAOLINESEARCHMT;
195:   PetscBool      flg;

199:   if (ls->setupcalled) return(0);
200:   if (!((PetscObject)ls)->type_name) {
201:     TaoLineSearchSetType(ls,default_type);
202:   }
203:   if (ls->ops->setup) {
204:     (*ls->ops->setup)(ls);
205:   }
206:   if (ls->usetaoroutines) {
207:     TaoIsObjectiveDefined(ls->tao,&flg);
208:     ls->hasobjective = flg;
209:     TaoIsGradientDefined(ls->tao,&flg);
210:     ls->hasgradient = flg;
211:     TaoIsObjectiveAndGradientDefined(ls->tao,&flg);
212:     ls->hasobjectiveandgradient = flg;
213:   } else {
214:     if (ls->ops->computeobjective) {
215:       ls->hasobjective = PETSC_TRUE;
216:     } else {
217:       ls->hasobjective = PETSC_FALSE;
218:     }
219:     if (ls->ops->computegradient) {
220:       ls->hasgradient = PETSC_TRUE;
221:     } else {
222:       ls->hasgradient = PETSC_FALSE;
223:     }
224:     if (ls->ops->computeobjectiveandgradient) {
225:       ls->hasobjectiveandgradient = PETSC_TRUE;
226:     } else {
227:       ls->hasobjectiveandgradient = PETSC_FALSE;
228:     }
229:   }
230:   ls->setupcalled = PETSC_TRUE;
231:   return(0);
232: }

236: /*@
237:   TaoLineSearchReset - Some line searches may carry state information
238:   from one TaoLineSearchApply() to the next.  This function resets this
239:   state information.

241:   Collective on TaoLineSearch

243:   Input Parameter:
244: . ls - the TaoLineSearch context

246:   Level: developer

248: .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
249: @*/
250: PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
251: {

256:   if (ls->ops->reset) {
257:     (*ls->ops->reset)(ls);
258:   }
259:   return(0);
260: }

264: /*@
265:   TaoLineSearchDestroy - Destroys the TAO context that was created with
266:   TaoLineSearchCreate()

268:   Collective on TaoLineSearch

270:   Input Parameter
271: . ls - the TaoLineSearch context

273:   Level: beginner

275: .seealse: TaoLineSearchCreate(), TaoLineSearchSolve()
276: @*/
277: PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
278: {

282:   if (!*ls) return(0);
284:   if (--((PetscObject)*ls)->refct > 0) {*ls=0; return(0);}
285:   VecDestroy(&(*ls)->stepdirection);
286:   VecDestroy(&(*ls)->start_x);
287:   if ((*ls)->ops->destroy) {
288:     (*(*ls)->ops->destroy)(*ls);
289:   }
290:   PetscHeaderDestroy(ls);
291:   return(0);
292: }

296: /*@
297:   TaoLineSearchApply - Performs a line-search in a given step direction.  Criteria for acceptable step length depends on the line-search algorithm chosen

299:   Collective on TaoLineSearch

301:   Input Parameters:
302: + ls - the Tao context
303: . x - The current solution (on output x contains the new solution determined by the line search)
304: . f - objective function value at current solution (on output contains the objective function value at new solution)
305: . g - gradient evaluated at x (on output contains the gradient at new solution)
306: - s - search direction

308:   Output Parameters:
309: + x - new solution
310: . f - objective function value at x
311: . g - gradient vector at x
312: . steplength - scalar multiplier of s used ( x = x0 + steplength * x )
313: - reason - reason why the line-search stopped

315:   reason will be set to one of:

317: + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
318: . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
319: . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
320: . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
321: . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
322: . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
323: . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
324: . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
325: . TAOLINESEARCH_HALTED_OTHER - any other reason
326: - TAOLINESEARCH_SUCCESS - successful line search

328:   Note:
329:   The algorithm developer must set up the TaoLineSearch with calls to
330:   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines()

332:   Note:
333:   You may or may not need to follow this with a call to
334:   TaoAddLineSearchCounts(), depending on whether you want these
335:   evaluations to count toward the total function/gradient evaluations.

337:   Level: beginner

339:   .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts()
340:  @*/

342: PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
343: {
345:   PetscInt       low1,low2,low3,high1,high2,high3;

348:   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
358:   VecGetOwnershipRange(x, &low1, &high1);
359:   VecGetOwnershipRange(g, &low2, &high2);
360:   VecGetOwnershipRange(s, &low3, &high3);
361:   if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths");

363:   PetscObjectReference((PetscObject)s);
364:   VecDestroy(&ls->stepdirection);
365:   ls->stepdirection = s;

367:   TaoLineSearchSetUp(ls);
368:   if (!ls->ops->apply) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine");
369:   ls->nfeval=0;
370:   ls->ngeval=0;
371:   ls->nfgeval=0;
372:   /* Check parameter values */
373:   if (ls->ftol < 0.0) {
374:     PetscInfo1(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol);
375:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
376:   }
377:   if (ls->rtol < 0.0) {
378:     PetscInfo1(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol);
379:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
380:   }
381:   if (ls->gtol < 0.0) {
382:     PetscInfo1(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol);
383:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
384:   }
385:   if (ls->stepmin < 0.0) {
386:     PetscInfo1(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin);
387:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
388:   }
389:   if (ls->stepmax < ls->stepmin) {
390:     PetscInfo2(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax);
391:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
392:   }
393:   if (ls->max_funcs < 0) {
394:     PetscInfo1(ls,"Bad Line Search Parameter: max_funcs (%D) < 0\n",ls->max_funcs);
395:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
396:   }
397:   if (PetscIsInfOrNanReal(*f)) {
398:     PetscInfo1(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f);
399:     *reason=TAOLINESEARCH_FAILED_INFORNAN;
400:   }

402:   PetscObjectReference((PetscObject)x);
403:   VecDestroy(&ls->start_x);
404:   ls->start_x = x;

406:   PetscLogEventBegin(TaoLineSearch_ApplyEvent,ls,0,0,0);
407:   (*ls->ops->apply)(ls,x,f,g,s);
408:   PetscLogEventEnd(TaoLineSearch_ApplyEvent, ls, 0,0,0);
409:   *reason=ls->reason;
410:   ls->new_f = *f;

412:   if (steplength) {
413:     *steplength=ls->step;
414:   }

416:   TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");
417:   return(0);
418: }

422: /*@C
423:    TaoLineSearchSetType - Sets the algorithm used in a line search

425:    Collective on TaoLineSearch

427:    Input Parameters:
428: +  ls - the TaoLineSearch context
429: -  type - a known method

431:   Available methods include:
432: + more-thuente
433: . gpcg
434: - unit - Do not perform any line search


437:   Options Database Keys:
438: .   -tao_ls_type - select which method TAO should use

440:   Level: beginner


443: .seealso: TaoLineSearchCreate(), TaoLineSearchGetType(), TaoLineSearchApply()

445: @*/

447: PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, const TaoLineSearchType type)
448: {
450:   PetscErrorCode (*r)(TaoLineSearch);
451:   PetscBool      flg;

456:   PetscObjectTypeCompare((PetscObject)ls, type, &flg);
457:   if (flg) return(0);

459:   PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r);
460:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type);
461:   if (ls->ops->destroy) {
462:     (*(ls)->ops->destroy)(ls);
463:   }
464:   ls->max_funcs=30;
465:   ls->ftol = 0.0001;
466:   ls->gtol = 0.9;
467: #if defined(PETSC_USE_REAL_SINGLE)
468:   ls->rtol = 1.0e-5;
469: #else
470:   ls->rtol = 1.0e-10;
471: #endif
472:   ls->stepmin=1.0e-20;
473:   ls->stepmax=1.0e+20;

475:   ls->nfeval=0;
476:   ls->ngeval=0;
477:   ls->nfgeval=0;
478:   ls->ops->setup=0;
479:   ls->ops->apply=0;
480:   ls->ops->view=0;
481:   ls->ops->setfromoptions=0;
482:   ls->ops->destroy=0;
483:   ls->setupcalled = PETSC_FALSE;
484:   (*r)(ls);
485:   PetscObjectChangeTypeName((PetscObject)ls, type);
486:   return(0);
487: }

491: /*@
492:   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
493:   options.

495:   Collective on TaoLineSearch

497:   Input Paremeter:
498: . ls - the TaoLineSearch context

500:   Options Database Keys:
501: + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
502: . -tao_ls_ftol <tol> - tolerance for sufficient decrease
503: . -tao_ls_gtol <tol> - tolerance for curvature condition
504: . -tao_ls_rtol <tol> - relative tolerance for acceptable step
505: . -tao_ls_stepmin <step> - minimum steplength allowed
506: . -tao_ls_stepmax <step> - maximum steplength allowed
507: . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
508: - -tao_ls_view - display line-search results to standard output

510:   Level: beginner
511: @*/
512: PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
513: {
515:   const char     *default_type=TAOLINESEARCHMT;
516:   char           type[256];
517:   PetscBool      flg;

521:   PetscObjectOptionsBegin((PetscObject)ls);
522:   if (!TaoLineSearchInitialized) {
523:     TaoLineSearchInitializePackage();
524:   }
525:   if (((PetscObject)ls)->type_name) {
526:     default_type = ((PetscObject)ls)->type_name;
527:   }
528:   /* Check for type from options */
529:   PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);
530:   if (flg) {
531:     TaoLineSearchSetType(ls,type);
532:   } else if (!((PetscObject)ls)->type_name) {
533:     TaoLineSearchSetType(ls,default_type);
534:   }

536:   PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,0);
537:   PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,0);
538:   PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,0);
539:   PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,0);
540:   PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,0);
541:   PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,0);
542:   if (ls->ops->setfromoptions) {
543:     (*ls->ops->setfromoptions)(ls);
544:   }
545:   PetscOptionsEnd();
546:   return(0);
547: }

551: /*@C
552:   TaoLineSearchGetType - Gets the current line search algorithm

554:   Not Collective

556:   Input Parameter:
557: . ls - the TaoLineSearch context

559:   Output Paramter:
560: . type - the line search algorithm in effect

562:   Level: developer

564: @*/
565: PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, const TaoLineSearchType *type)
566: {
570:   *type = ((PetscObject)ls)->type_name;
571:   return(0);
572: }

576: /*@
577:   TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation
578:   routines used by the line search in last application (not cumulative).

580:   Not Collective

582:   Input Parameter:
583: . ls - the TaoLineSearch context

585:   Output Parameters:
586: + nfeval   - number of function evaluations
587: . ngeval   - number of gradient evaluations
588: - nfgeval  - number of function/gradient evaluations

590:   Level: intermediate

592:   Note:
593:   If the line search is using the Tao objective and gradient
594:   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
595:   is already counting the number of evaluations.

597: @*/
598: PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
599: {
602:   *nfeval = ls->nfeval;
603:   *ngeval = ls->ngeval;
604:   *nfgeval = ls->nfgeval;
605:   return(0);
606: }

610: /*@
611:   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
612:   Tao evaluation routines.

614:   Not Collective

616:   Input Parameter:
617: . ls - the TaoLineSearch context

619:   Output Parameter:
620: . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
621:         otherwise PETSC_FALSE

623:   Level: developer
624: @*/
625: PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
626: {
629:   *flg = ls->usetaoroutines;
630:   return(0);
631: }

635: /*@C
636:   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search

638:   Logically Collective on TaoLineSearch

640:   Input Parameter:
641: + ls - the TaoLineSearch context
642: . func - the objective function evaluation routine
643: - ctx - the (optional) user-defined context for private data

645:   Calling sequence of func:
646: $      func (TaoLinesearch ls, Vec x, PetscReal *f, void *ctx);

648: + x - input vector
649: . f - function value
650: - ctx (optional) user-defined context

652:   Level: beginner

654:   Note:
655:   Use this routine only if you want the line search objective
656:   evaluation routine to be different from the Tao's objective
657:   evaluation routine. If you use this routine you must also set
658:   the line search gradient and/or function/gradient routine.

660:   Note:
661:   Some algorithms (lcl, gpcg) set their own objective routine for the
662:   line search, application programmers should be wary of overriding the
663:   default objective routine.

665: .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
666: @*/
667: PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
668: {

672:   ls->ops->computeobjective=func;
673:   if (ctx) ls->userctx_func=ctx;
674:   ls->usetaoroutines=PETSC_FALSE;
675:   return(0);
676: }

680: /*@C
681:   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search

683:   Logically Collective on TaoLineSearch

685:   Input Parameter:
686: + ls - the TaoLineSearch context
687: . func - the gradient evaluation routine
688: - ctx - the (optional) user-defined context for private data

690:   Calling sequence of func:
691: $      func (TaoLinesearch ls, Vec x, Vec g, void *ctx);

693: + x - input vector
694: . g - gradient vector
695: - ctx (optional) user-defined context

697:   Level: beginner

699:   Note:
700:   Use this routine only if you want the line search gradient
701:   evaluation routine to be different from the Tao's gradient
702:   evaluation routine. If you use this routine you must also set
703:   the line search function and/or function/gradient routine.

705:   Note:
706:   Some algorithms (lcl, gpcg) set their own gradient routine for the
707:   line search, application programmers should be wary of overriding the
708:   default gradient routine.

710: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
711: @*/
712: PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
713: {
716:   ls->ops->computegradient=func;
717:   if (ctx) ls->userctx_grad=ctx;
718:   ls->usetaoroutines=PETSC_FALSE;
719:   return(0);
720: }

724: /*@C
725:   TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search

727:   Logically Collective on TaoLineSearch

729:   Input Parameter:
730: + ls - the TaoLineSearch context
731: . func - the objective and gradient evaluation routine
732: - ctx - the (optional) user-defined context for private data

734:   Calling sequence of func:
735: $      func (TaoLinesearch ls, Vec x, PetscReal *f, Vec g, void *ctx);

737: + x - input vector
738: . f - function value
739: . g - gradient vector
740: - ctx (optional) user-defined context

742:   Level: beginner

744:   Note:
745:   Use this routine only if you want the line search objective and gradient
746:   evaluation routines to be different from the Tao's objective
747:   and gradient evaluation routines.

749:   Note:
750:   Some algorithms (lcl, gpcg) set their own objective routine for the
751:   line search, application programmers should be wary of overriding the
752:   default objective routine.

754: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines()
755: @*/
756: PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx)
757: {
760:   ls->ops->computeobjectiveandgradient=func;
761:   if (ctx) ls->userctx_funcgrad=ctx;
762:   ls->usetaoroutines = PETSC_FALSE;
763:   return(0);
764: }

768: /*@C
769:   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
770:   (gradient'*stepdirection) evaluation routine for the line search.
771:   Sometimes it is more efficient to compute the inner product of the gradient
772:   and the step direction than it is to compute the gradient, and this is all
773:   the line search typically needs of the gradient.

775:   Logically Collective on TaoLineSearch

777:   Input Parameter:
778: + ls - the TaoLineSearch context
779: . func - the objective and gradient evaluation routine
780: - ctx - the (optional) user-defined context for private data

782:   Calling sequence of func:
783: $      func (TaoLinesearch ls, Vec x, PetscReal *f, PetscReal *gts, void *ctx);

785: + x - input vector
786: . s - step direction
787: . f - function value
788: . gts - inner product of gradient and step direction vectors
789: - ctx (optional) user-defined context

791:   Note: The gradient will still need to be computed at the end of the line
792:   search, so you will still need to set a line search gradient evaluation
793:   routine

795:   Note: Bounded line searches (those used in bounded optimization algorithms)
796:   don't use g's directly, but rather (g'x - g'x0)/steplength.  You can get the
797:   x0 and steplength with TaoLineSearchGetStartingVector() and TaoLineSearchGetStepLength()

799:   Level: advanced

801:   Note:
802:   Some algorithms (lcl, gpcg) set their own objective routine for the
803:   line search, application programmers should be wary of overriding the
804:   default objective routine.

806: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines()
807: @*/
808: PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx)
809: {
812:   ls->ops->computeobjectiveandgts=func;
813:   if (ctx) ls->userctx_funcgts=ctx;
814:   ls->usegts = PETSC_TRUE;
815:   ls->usetaoroutines=PETSC_FALSE;
816:   return(0);
817: }

821: /*@C
822:   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
823:   objective and gradient evaluation routines from the given Tao object.

825:   Logically Collective on TaoLineSearch

827:   Input Parameter:
828: + ls - the TaoLineSearch context
829: - ts - the Tao context with defined objective/gradient evaluation routines

831:   Level: developer

833: .seealso: TaoLineSearchCreate()
834: @*/
835: PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
836: {
840:   ls->tao = ts;
841:   ls->usetaoroutines=PETSC_TRUE;
842:   return(0);
843: }

847: /*@
848:   TaoLineSearchComputeObjective - Computes the objective function value at a given point

850:   Collective on TaoLineSearch

852:   Input Parameters:
853: + ls - the TaoLineSearch context
854: - x - input vector

856:   Output Parameter:
857: . f - Objective value at X

859:   Notes: TaoLineSearchComputeObjective() is typically used within line searches
860:   so most users would not generally call this routine themselves.

862:   Level: developer

864: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
865: @*/
866: PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
867: {
869:   Vec            gdummy;
870:   PetscReal      gts;

877:   if (ls->usetaoroutines) {
878:     TaoComputeObjective(ls->tao,x,f);
879:   } else {
880:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
881:     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient && !ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
882:     PetscStackPush("TaoLineSearch user objective routine");
883:     if (ls->ops->computeobjective) {
884:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
885:     } else if (ls->ops->computeobjectiveandgradient) {
886:       VecDuplicate(x,&gdummy);
887:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);
888:       VecDestroy(&gdummy);
889:     } else {
890:       (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,&gts,ls->userctx_funcgts);
891:     }
892:     PetscStackPop;
893:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
894:   }
895:   ls->nfeval++;
896:   return(0);
897: }

901: /*@
902:   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point

904:   Collective on Tao

906:   Input Parameters:
907: + ls - the TaoLineSearch context
908: - x - input vector

910:   Output Parameter:
911: + f - Objective value at X
912: - g - Gradient vector at X

914:   Notes: TaoLineSearchComputeObjectiveAndGradient() is typically used within line searches
915:   so most users would not generally call this routine themselves.

917:   Level: developer

919: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
920: @*/
921: PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
922: {

932:   if (ls->usetaoroutines) {
933:       TaoComputeObjectiveAndGradient(ls->tao,x,f,g);
934:   } else {
935:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
936:     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
937:     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set");

939:     PetscStackPush("TaoLineSearch user objective/gradient routine");
940:     if (ls->ops->computeobjectiveandgradient) {
941:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);
942:     } else {
943:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
944:       (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);
945:     }
946:     PetscStackPop;
947:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
948:     PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
949:   }
950:   ls->nfgeval++;
951:   return(0);
952: }

956: /*@
957:   TaoLineSearchComputeGradient - Computes the gradient of the objective function

959:   Collective on TaoLineSearch

961:   Input Parameters:
962: + ls - the TaoLineSearch context
963: - x - input vector

965:   Output Parameter:
966: . g - gradient vector

968:   Notes: TaoComputeGradient() is typically used within line searches
969:   so most users would not generally call this routine themselves.

971:   Level: developer

973: .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient()
974: @*/
975: PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
976: {
978:   PetscReal      fdummy;

986:   if (ls->usetaoroutines) {
987:     TaoComputeGradient(ls->tao,x,g);
988:   } else {
989:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
990:     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set");
991:     PetscStackPush("TaoLineSearch user gradient routine");
992:     if (ls->ops->computegradient) {
993:       (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);
994:     } else {
995:       (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);
996:     }
997:     PetscStackPop;
998:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
999:   }
1000:   ls->ngeval++;
1001:   return(0);
1002: }

1006: /*@
1007:   TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and step direction at a given point

1009:   Collective on Tao

1011:   Input Parameters:
1012: + ls - the TaoLineSearch context
1013: - x - input vector

1015:   Output Parameter:
1016: + f - Objective value at X
1017: - gts - inner product of gradient and step direction at X

1019:   Notes: TaoLineSearchComputeObjectiveAndGTS() is typically used within line searches
1020:   so most users would not generally call this routine themselves.

1022:   Level: developer

1024: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
1025: @*/
1026: PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
1027: {
1035:   PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
1036:   if (!ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set");
1037:   PetscStackPush("TaoLineSearch user objective/gts routine");
1038:   (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);
1039:   PetscStackPop;
1040:   PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
1041:   PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
1042:   ls->nfeval++;
1043:   return(0);
1044: }

1048: /*@
1049:   TaoLineSearchGetSolution - Returns the solution to the line search

1051:   Collective on TaoLineSearch

1053:   Input Parameter:
1054: . ls - the TaoLineSearch context

1056:   Output Parameter:
1057: + x - the new solution
1058: . f - the objective function value at x
1059: . g - the gradient at x
1060: . steplength - the multiple of the step direction taken by the line search
1061: - reason - the reason why the line search terminated

1063:   reason will be set to one of:

1065: + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
1066: . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
1067: . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
1068: . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
1069: . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
1070: . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
1071: . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance

1073: . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
1074: . TAOLINESEARCH_HALTED_OTHER - any other reason

1076: + TAOLINESEARCH_SUCCESS - successful line search

1078:   Level: developer

1080: @*/
1081: PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
1082: {

1091:   if (ls->new_x) {
1092:     VecCopy(ls->new_x,x);
1093:   }
1094:   *f = ls->new_f;
1095:   if (ls->new_g) {
1096:     VecCopy(ls->new_g,g);
1097:   }
1098:   if (steplength) {
1099:     *steplength=ls->step;
1100:   }
1101:   *reason = ls->reason;
1102:   return(0);
1103: }

1107: /*@
1108:   TaoLineSearchGetStartingVector - Gets a the initial point of the line
1109:   search.

1111:   Not Collective

1113:   Input Parameter:
1114: . ls - the TaoLineSearch context

1116:   Output Parameter:
1117: . x - The initial point of the line search

1119:   Level: intermediate
1120: @*/
1121: PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1122: {
1125:   if (x) {
1126:     *x = ls->start_x;
1127:   }
1128:   return(0);
1129: }

1133: /*@
1134:   TaoLineSearchGetStepDirection - Gets the step direction of the line
1135:   search.

1137:   Not Collective

1139:   Input Parameter:
1140: . ls - the TaoLineSearch context

1142:   Output Parameter:
1143: . s - the step direction of the line search

1145:   Level: advanced
1146: @*/
1147: PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1148: {
1151:   if (s) {
1152:     *s = ls->stepdirection;
1153:   }
1154:   return(0);

1156: }

1160: /*@
1161:   TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step.  Useful for some minimization algorithms.

1163:   Not Collective

1165:   Input Parameter:
1166: . ls - the TaoLineSearch context

1168:   Output Parameter:
1169: . f - the objective value at the full step length

1171:   Level: developer
1172: @*/

1174: PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1175: {
1178:   *f_fullstep = ls->f_fullstep;
1179:   return(0);
1180: }

1184: /*@
1185:   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.

1187:   Logically Collective on Tao

1189:   Input Parameters:
1190: + ls - the TaoLineSearch context
1191: . xl  - vector of lower bounds
1192: - xu  - vector of upper bounds

1194:   Note: If the variable bounds are not set with this routine, then
1195:   PETSC_NINFINITY and PETSC_INFINITY are assumed

1197:   Level: beginner

1199: .seealso: TaoSetVariableBounds(), TaoLineSearchCreate()
1200: @*/
1201: PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
1202: {
1207:   ls->lower = xl;
1208:   ls->upper = xu;
1209:   ls->bounded = 1;
1210:   return(0);
1211: }

1215: /*@
1216:   TaoLineSearchSetInitialStepLength - Sets the initial step length of a line
1217:   search.  If this value is not set then 1.0 is assumed.

1219:   Logically Collective on TaoLineSearch

1221:   Input Parameters:
1222: + ls - the TaoLineSearch context
1223: - s - the initial step size

1225:   Level: intermediate

1227: .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply()
1228: @*/
1229: PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
1230: {
1233:   ls->initstep = s;
1234:   return(0);
1235: }

1239: /*@
1240:   TaoLineSearchGetStepLength - Get the current step length

1242:   Not Collective

1244:   Input Parameters:
1245: . ls - the TaoLineSearch context

1247:   Output Parameters
1248: . s - the current step length

1250:   Level: beginner

1252: .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply()
1253: @*/
1254: PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s)
1255: {
1258:   *s = ls->step;
1259:   return(0);
1260: }

1264: /*MC
1265:    TaoLineSearchRegister - Adds a line-search algorithm to the registry

1267:    Not collective

1269:    Input Parameters:
1270: +  sname - name of a new user-defined solver
1271: -  func - routine to Create method context

1273:    Notes:
1274:    TaoLineSearchRegister() may be called multiple times to add several user-defined solvers.

1276:    Sample usage:
1277: .vb
1278:    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
1279: .ve

1281:    Then, your solver can be chosen with the procedural interface via
1282: $     TaoLineSearchSetType(ls,"my_linesearch")
1283:    or at runtime via the option
1284: $     -tao_ls_type my_linesearch

1286:    Level: developer

1288: .seealso: TaoLineSearchRegisterDestroy()
1289: M*/
1290: PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1291: {
1294:   PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);
1295:   return(0);
1296: }

1300: /*@C
1301:    TaoLineSearchRegisterDestroy - Frees the list of line-search algorithms that were
1302:    registered by TaoLineSearchRegister().

1304:    Not Collective

1306:    Level: developer

1308: .seealso: TaoLineSearchRegister()
1309: @*/
1310: PetscErrorCode TaoLineSearchRegisterDestroy(void)
1311: {
1314:   PetscFunctionListDestroy(&TaoLineSearchList);
1315:   TaoLineSearchInitialized = PETSC_FALSE;
1316:   return(0);
1317: }

1321: /*@C
1322:    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1323:    for all TaoLineSearch options in the database.


1326:    Collective on TaoLineSearch

1328:    Input Parameters:
1329: +  ls - the TaoLineSearch solver context
1330: -  prefix - the prefix string to prepend to all line search requests

1332:    Notes:
1333:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1334:    The first character of all runtime options is AUTOMATICALLY the hyphen.


1337:    Level: advanced

1339: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1340: @*/
1341: PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1342: {
1343:   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
1344: }

1348: /*@C
1349:   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1350:   TaoLineSearch options in the database

1352:   Not Collective

1354:   Input Parameters:
1355: . ls - the TaoLineSearch context

1357:   Output Parameters:
1358: . prefix - pointer to the prefix string used is returned

1360:   Notes: On the fortran side, the user should pass in a string 'prefix' of
1361:   sufficient length to hold the prefix.

1363:   Level: advanced

1365: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix()
1366: @*/
1367: PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1368: {
1369:   return PetscObjectGetOptionsPrefix((PetscObject)ls,p);
1370: }

1374: /*@C
1375:    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1376:    TaoLineSearch options in the database.


1379:    Logically Collective on TaoLineSearch

1381:    Input Parameters:
1382: +  ls - the TaoLineSearch context
1383: -  prefix - the prefix string to prepend to all TAO option requests

1385:    Notes:
1386:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1387:    The first character of all runtime options is AUTOMATICALLY the hyphen.

1389:    For example, to distinguish between the runtime options for two
1390:    different line searches, one could call
1391: .vb
1392:       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1393:       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1394: .ve

1396:    This would enable use of different options for each system, such as
1397: .vb
1398:       -sys1_tao_ls_type mt
1399:       -sys2_tao_ls_type armijo
1400: .ve

1402:    Level: advanced

1404: .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1405: @*/

1407: PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1408: {
1409:   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
1410: }