Actual source code: taolinesearch.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1:  #include <petsctaolinesearch.h>
  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;

 10: /*@C
 11:   TaoLineSearchView - Prints information about the TaoLineSearch

 13:   Collective on TaoLineSearch

 15:   InputParameters:
 16: + ls - the Tao context
 17: - viewer - visualization context

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

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

 30:   Level: beginner

 32: .seealso: PetscViewerASCIIOpen()
 33: @*/

 35: PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
 36: {
 37:   PetscErrorCode          ierr;
 38:   PetscBool               isascii, isstring;
 39:   const TaoLineSearchType type;

 43:   if (!viewer) {
 44:     PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer);
 45:   }

 49:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 50:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
 51:   if (isascii) {
 52:     PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer);
 53:     if (ls->ops->view) {
 54:       PetscViewerASCIIPushTab(viewer);
 55:       (*ls->ops->view)(ls,viewer);
 56:       PetscViewerASCIIPopTab(viewer);
 57:     }
 58:     PetscViewerASCIIPushTab(viewer);
 59:     PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);
 60:     PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);
 61:     PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);
 62:     PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);
 63:     PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);

 65:     if (ls->bounded) {
 66:       PetscViewerASCIIPrintf(viewer,"using variable bounds\n");
 67:     }
 68:     PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason);
 69:     PetscViewerASCIIPopTab(viewer);
 70:   } else if (isstring) {
 71:     TaoLineSearchGetType(ls,&type);
 72:     PetscViewerStringSPrintf(viewer," %-3.3s",type);
 73:   }
 74:   return(0);
 75: }

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

 81:   Collective on MPI_Comm

 83:   Input Parameter:
 84: . comm - MPI communicator

 86:   Output Parameter:
 87: . newls - the new TaoLineSearch context

 89:   Available methods include:
 90: + more-thuente
 91: . gpcg
 92: - unit - Do not perform any line search


 95:    Options Database Keys:
 96: .   -tao_ls_type - select which method TAO should use

 98:    Level: beginner

100: .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy()
101: @*/

103: PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
104: {
106:   TaoLineSearch  ls;

110:   *newls = NULL;

112:  #ifndef PETSC_USE_DYNAMIC_LIBRARIES
113:   TaoLineSearchInitializePackage();
114:  #endif

116:   PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView);
117:   ls->bounded = 0;
118:   ls->max_funcs=30;
119:   ls->ftol = 0.0001;
120:   ls->gtol = 0.9;
121: #if defined(PETSC_USE_REAL_SINGLE)
122:   ls->rtol = 1.0e-5;
123: #else
124:   ls->rtol = 1.0e-10;
125: #endif
126:   ls->stepmin=1.0e-20;
127:   ls->stepmax=1.0e+20;
128:   ls->step=1.0;
129:   ls->nfeval=0;
130:   ls->ngeval=0;
131:   ls->nfgeval=0;

133:   ls->ops->computeobjective=0;
134:   ls->ops->computegradient=0;
135:   ls->ops->computeobjectiveandgradient=0;
136:   ls->ops->computeobjectiveandgts=0;
137:   ls->ops->setup=0;
138:   ls->ops->apply=0;
139:   ls->ops->view=0;
140:   ls->ops->setfromoptions=0;
141:   ls->ops->reset=0;
142:   ls->ops->destroy=0;
143:   ls->setupcalled=PETSC_FALSE;
144:   ls->usetaoroutines=PETSC_FALSE;
145:   *newls = ls;
146:   return(0);
147: }

149: /*@
150:   TaoLineSearchSetUp - Sets up the internal data structures for the later use
151:   of a Tao solver

153:   Collective on ls

155:   Input Parameters:
156: . ls - the TaoLineSearch context

158:   Notes:
159:   The user will not need to explicitly call TaoLineSearchSetUp(), as it will
160:   automatically be called in TaoLineSearchSolve().  However, if the user
161:   desires to call it explicitly, it should come after TaoLineSearchCreate()
162:   but before TaoLineSearchApply().

164:   Level: developer

166: .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
167: @*/

169: PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
170: {
172:   const char     *default_type=TAOLINESEARCHMT;
173:   PetscBool      flg;

177:   if (ls->setupcalled) return(0);
178:   if (!((PetscObject)ls)->type_name) {
179:     TaoLineSearchSetType(ls,default_type);
180:   }
181:   if (ls->ops->setup) {
182:     (*ls->ops->setup)(ls);
183:   }
184:   if (ls->usetaoroutines) {
185:     TaoIsObjectiveDefined(ls->tao,&flg);
186:     ls->hasobjective = flg;
187:     TaoIsGradientDefined(ls->tao,&flg);
188:     ls->hasgradient = flg;
189:     TaoIsObjectiveAndGradientDefined(ls->tao,&flg);
190:     ls->hasobjectiveandgradient = flg;
191:   } else {
192:     if (ls->ops->computeobjective) {
193:       ls->hasobjective = PETSC_TRUE;
194:     } else {
195:       ls->hasobjective = PETSC_FALSE;
196:     }
197:     if (ls->ops->computegradient) {
198:       ls->hasgradient = PETSC_TRUE;
199:     } else {
200:       ls->hasgradient = PETSC_FALSE;
201:     }
202:     if (ls->ops->computeobjectiveandgradient) {
203:       ls->hasobjectiveandgradient = PETSC_TRUE;
204:     } else {
205:       ls->hasobjectiveandgradient = PETSC_FALSE;
206:     }
207:   }
208:   ls->setupcalled = PETSC_TRUE;
209:   return(0);
210: }

212: /*@
213:   TaoLineSearchReset - Some line searches may carry state information
214:   from one TaoLineSearchApply() to the next.  This function resets this
215:   state information.

217:   Collective on TaoLineSearch

219:   Input Parameter:
220: . ls - the TaoLineSearch context

222:   Level: developer

224: .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
225: @*/
226: PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
227: {

232:   if (ls->ops->reset) {
233:     (*ls->ops->reset)(ls);
234:   }
235:   return(0);
236: }

238: /*@
239:   TaoLineSearchDestroy - Destroys the TAO context that was created with
240:   TaoLineSearchCreate()

242:   Collective on TaoLineSearch

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

247:   Level: beginner

249: .seealse: TaoLineSearchCreate(), TaoLineSearchSolve()
250: @*/
251: PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
252: {

256:   if (!*ls) return(0);
258:   if (--((PetscObject)*ls)->refct > 0) {*ls=0; return(0);}
259:   VecDestroy(&(*ls)->stepdirection);
260:   VecDestroy(&(*ls)->start_x);
261:   if ((*ls)->ops->destroy) {
262:     (*(*ls)->ops->destroy)(*ls);
263:   }
264:   PetscHeaderDestroy(ls);
265:   return(0);
266: }

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

271:   Collective on TaoLineSearch

273:   Input Parameters:
274: + ls - the Tao context
275: . x - The current solution (on output x contains the new solution determined by the line search)
276: . f - objective function value at current solution (on output contains the objective function value at new solution)
277: . g - gradient evaluated at x (on output contains the gradient at new solution)
278: - s - search direction

280:   Output Parameters:
281: + x - new solution
282: . f - objective function value at x
283: . g - gradient vector at x
284: . steplength - scalar multiplier of s used ( x = x0 + steplength * x )
285: - reason - reason why the line-search stopped

287:   reason will be set to one of:

289: + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
290: . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
291: . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
292: . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
293: . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
294: . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
295: . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
296: . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
297: . TAOLINESEARCH_HALTED_OTHER - any other reason
298: - TAOLINESEARCH_SUCCESS - successful line search

300:   Note:
301:   The algorithm developer must set up the TaoLineSearch with calls to
302:   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines()

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

309:   Level: beginner

311:   .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts()
312:  @*/

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

320:   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
330:   VecGetOwnershipRange(x, &low1, &high1);
331:   VecGetOwnershipRange(g, &low2, &high2);
332:   VecGetOwnershipRange(s, &low3, &high3);
333:   if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths");

335:   PetscObjectReference((PetscObject)s);
336:   VecDestroy(&ls->stepdirection);
337:   ls->stepdirection = s;

339:   TaoLineSearchSetUp(ls);
340:   if (!ls->ops->apply) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search Object does not have 'apply' routine");
341:   ls->nfeval=0;
342:   ls->ngeval=0;
343:   ls->nfgeval=0;
344:   /* Check parameter values */
345:   if (ls->ftol < 0.0) {
346:     PetscInfo1(ls,"Bad Line Search Parameter: ftol (%g) < 0\n",(double)ls->ftol);
347:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
348:   }
349:   if (ls->rtol < 0.0) {
350:     PetscInfo1(ls,"Bad Line Search Parameter: rtol (%g) < 0\n",(double)ls->rtol);
351:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
352:   }
353:   if (ls->gtol < 0.0) {
354:     PetscInfo1(ls,"Bad Line Search Parameter: gtol (%g) < 0\n",(double)ls->gtol);
355:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
356:   }
357:   if (ls->stepmin < 0.0) {
358:     PetscInfo1(ls,"Bad Line Search Parameter: stepmin (%g) < 0\n",(double)ls->stepmin);
359:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
360:   }
361:   if (ls->stepmax < ls->stepmin) {
362:     PetscInfo2(ls,"Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n",(double)ls->stepmin,(double)ls->stepmax);
363:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
364:   }
365:   if (ls->max_funcs < 0) {
366:     PetscInfo1(ls,"Bad Line Search Parameter: max_funcs (%D) < 0\n",ls->max_funcs);
367:     *reason=TAOLINESEARCH_FAILED_BADPARAMETER;
368:   }
369:   if (PetscIsInfOrNanReal(*f)) {
370:     PetscInfo1(ls,"Initial Line Search Function Value is Inf or Nan (%g)\n",(double)*f);
371:     *reason=TAOLINESEARCH_FAILED_INFORNAN;
372:   }

374:   PetscObjectReference((PetscObject)x);
375:   VecDestroy(&ls->start_x);
376:   ls->start_x = x;

378:   PetscLogEventBegin(TaoLineSearch_ApplyEvent,ls,0,0,0);
379:   (*ls->ops->apply)(ls,x,f,g,s);
380:   PetscLogEventEnd(TaoLineSearch_ApplyEvent, ls, 0,0,0);
381:   *reason=ls->reason;
382:   ls->new_f = *f;

384:   if (steplength) {
385:     *steplength=ls->step;
386:   }

388:   TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");
389:   return(0);
390: }

392: /*@C
393:    TaoLineSearchSetType - Sets the algorithm used in a line search

395:    Collective on TaoLineSearch

397:    Input Parameters:
398: +  ls - the TaoLineSearch context
399: -  type - a known method

401:   Available methods include:
402: + more-thuente
403: . gpcg
404: - unit - Do not perform any line search


407:   Options Database Keys:
408: .   -tao_ls_type - select which method TAO should use

410:   Level: beginner


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

415: @*/

417: PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, const TaoLineSearchType type)
418: {
420:   PetscErrorCode (*r)(TaoLineSearch);
421:   PetscBool      flg;

426:   PetscObjectTypeCompare((PetscObject)ls, type, &flg);
427:   if (flg) return(0);

429:   PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r);
430:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type);
431:   if (ls->ops->destroy) {
432:     (*(ls)->ops->destroy)(ls);
433:   }
434:   ls->max_funcs=30;
435:   ls->ftol = 0.0001;
436:   ls->gtol = 0.9;
437: #if defined(PETSC_USE_REAL_SINGLE)
438:   ls->rtol = 1.0e-5;
439: #else
440:   ls->rtol = 1.0e-10;
441: #endif
442:   ls->stepmin=1.0e-20;
443:   ls->stepmax=1.0e+20;

445:   ls->nfeval=0;
446:   ls->ngeval=0;
447:   ls->nfgeval=0;
448:   ls->ops->setup=0;
449:   ls->ops->apply=0;
450:   ls->ops->view=0;
451:   ls->ops->setfromoptions=0;
452:   ls->ops->destroy=0;
453:   ls->setupcalled = PETSC_FALSE;
454:   (*r)(ls);
455:   PetscObjectChangeTypeName((PetscObject)ls, type);
456:   return(0);
457: }

459: /*@
460:   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
461:   options.

463:   Collective on TaoLineSearch

465:   Input Paremeter:
466: . ls - the TaoLineSearch context

468:   Options Database Keys:
469: + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
470: . -tao_ls_ftol <tol> - tolerance for sufficient decrease
471: . -tao_ls_gtol <tol> - tolerance for curvature condition
472: . -tao_ls_rtol <tol> - relative tolerance for acceptable step
473: . -tao_ls_stepmin <step> - minimum steplength allowed
474: . -tao_ls_stepmax <step> - maximum steplength allowed
475: . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
476: - -tao_ls_view - display line-search results to standard output

478:   Level: beginner
479: @*/
480: PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
481: {
483:   const char     *default_type=TAOLINESEARCHMT;
484:   char           type[256];
485:   PetscBool      flg;

489:   PetscObjectOptionsBegin((PetscObject)ls);
490:   if (!TaoLineSearchInitialized) {
491:     TaoLineSearchInitializePackage();
492:   }
493:   if (((PetscObject)ls)->type_name) {
494:     default_type = ((PetscObject)ls)->type_name;
495:   }
496:   /* Check for type from options */
497:   PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);
498:   if (flg) {
499:     TaoLineSearchSetType(ls,type);
500:   } else if (!((PetscObject)ls)->type_name) {
501:     TaoLineSearchSetType(ls,default_type);
502:   }

504:   PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL);
505:   PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL);
506:   PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL);
507:   PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL);
508:   PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL);
509:   PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL);
510:   if (ls->ops->setfromoptions) {
511:     (*ls->ops->setfromoptions)(PetscOptionsObject,ls);
512:   }
513:   PetscOptionsEnd();
514:   return(0);
515: }

517: /*@C
518:   TaoLineSearchGetType - Gets the current line search algorithm

520:   Not Collective

522:   Input Parameter:
523: . ls - the TaoLineSearch context

525:   Output Paramter:
526: . type - the line search algorithm in effect

528:   Level: developer

530: @*/
531: PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, const TaoLineSearchType *type)
532: {
536:   *type = ((PetscObject)ls)->type_name;
537:   return(0);
538: }

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

544:   Not Collective

546:   Input Parameter:
547: . ls - the TaoLineSearch context

549:   Output Parameters:
550: + nfeval   - number of function evaluations
551: . ngeval   - number of gradient evaluations
552: - nfgeval  - number of function/gradient evaluations

554:   Level: intermediate

556:   Note:
557:   If the line search is using the Tao objective and gradient
558:   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
559:   is already counting the number of evaluations.

561: @*/
562: PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
563: {
566:   *nfeval = ls->nfeval;
567:   *ngeval = ls->ngeval;
568:   *nfgeval = ls->nfgeval;
569:   return(0);
570: }

572: /*@
573:   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
574:   Tao evaluation routines.

576:   Not Collective

578:   Input Parameter:
579: . ls - the TaoLineSearch context

581:   Output Parameter:
582: . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
583:         otherwise PETSC_FALSE

585:   Level: developer
586: @*/
587: PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
588: {
591:   *flg = ls->usetaoroutines;
592:   return(0);
593: }

595: /*@C
596:   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search

598:   Logically Collective on TaoLineSearch

600:   Input Parameter:
601: + ls - the TaoLineSearch context
602: . func - the objective function evaluation routine
603: - ctx - the (optional) user-defined context for private data

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

608: + x - input vector
609: . f - function value
610: - ctx (optional) user-defined context

612:   Level: beginner

614:   Note:
615:   Use this routine only if you want the line search objective
616:   evaluation routine to be different from the Tao's objective
617:   evaluation routine. If you use this routine you must also set
618:   the line search gradient and/or function/gradient routine.

620:   Note:
621:   Some algorithms (lcl, gpcg) set their own objective routine for the
622:   line search, application programmers should be wary of overriding the
623:   default objective routine.

625: .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
626: @*/
627: PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
628: {

632:   ls->ops->computeobjective=func;
633:   if (ctx) ls->userctx_func=ctx;
634:   ls->usetaoroutines=PETSC_FALSE;
635:   return(0);
636: }

638: /*@C
639:   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search

641:   Logically Collective on TaoLineSearch

643:   Input Parameter:
644: + ls - the TaoLineSearch context
645: . func - the gradient evaluation routine
646: - ctx - the (optional) user-defined context for private data

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

651: + x - input vector
652: . g - gradient vector
653: - ctx (optional) user-defined context

655:   Level: beginner

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

663:   Note:
664:   Some algorithms (lcl, gpcg) set their own gradient routine for the
665:   line search, application programmers should be wary of overriding the
666:   default gradient routine.

668: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
669: @*/
670: PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
671: {
674:   ls->ops->computegradient=func;
675:   if (ctx) ls->userctx_grad=ctx;
676:   ls->usetaoroutines=PETSC_FALSE;
677:   return(0);
678: }

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

683:   Logically Collective on TaoLineSearch

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

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

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

698:   Level: beginner

700:   Note:
701:   Use this routine only if you want the line search objective and gradient
702:   evaluation routines to be different from the Tao's objective
703:   and gradient evaluation routines.

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

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

722: /*@C
723:   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
724:   (gradient'*stepdirection) evaluation routine for the line search.
725:   Sometimes it is more efficient to compute the inner product of the gradient
726:   and the step direction than it is to compute the gradient, and this is all
727:   the line search typically needs of the gradient.

729:   Logically Collective on TaoLineSearch

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

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

739: + x - input vector
740: . s - step direction
741: . f - function value
742: . gts - inner product of gradient and step direction vectors
743: - ctx (optional) user-defined context

745:   Note: The gradient will still need to be computed at the end of the line
746:   search, so you will still need to set a line search gradient evaluation
747:   routine

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

753:   Level: advanced

755:   Note:
756:   Some algorithms (lcl, gpcg) set their own objective routine for the
757:   line search, application programmers should be wary of overriding the
758:   default objective routine.

760: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines()
761: @*/
762: PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx)
763: {
766:   ls->ops->computeobjectiveandgts=func;
767:   if (ctx) ls->userctx_funcgts=ctx;
768:   ls->usegts = PETSC_TRUE;
769:   ls->usetaoroutines=PETSC_FALSE;
770:   return(0);
771: }

773: /*@C
774:   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
775:   objective and gradient evaluation routines from the given Tao object.

777:   Logically Collective on TaoLineSearch

779:   Input Parameter:
780: + ls - the TaoLineSearch context
781: - ts - the Tao context with defined objective/gradient evaluation routines

783:   Level: developer

785: .seealso: TaoLineSearchCreate()
786: @*/
787: PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
788: {
792:   ls->tao = ts;
793:   ls->usetaoroutines=PETSC_TRUE;
794:   return(0);
795: }

797: /*@
798:   TaoLineSearchComputeObjective - Computes the objective function value at a given point

800:   Collective on TaoLineSearch

802:   Input Parameters:
803: + ls - the TaoLineSearch context
804: - x - input vector

806:   Output Parameter:
807: . f - Objective value at X

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

812:   Level: developer

814: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
815: @*/
816: PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
817: {
819:   Vec            gdummy;
820:   PetscReal      gts;

827:   if (ls->usetaoroutines) {
828:     TaoComputeObjective(ls->tao,x,f);
829:   } else {
830:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
831:     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");
832:     PetscStackPush("TaoLineSearch user objective routine");
833:     if (ls->ops->computeobjective) {
834:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
835:     } else if (ls->ops->computeobjectiveandgradient) {
836:       VecDuplicate(x,&gdummy);
837:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);
838:       VecDestroy(&gdummy);
839:     } else {
840:       (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,&gts,ls->userctx_funcgts);
841:     }
842:     PetscStackPop;
843:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
844:   }
845:   ls->nfeval++;
846:   return(0);
847: }

849: /*@
850:   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point

852:   Collective on Tao

854:   Input Parameters:
855: + ls - the TaoLineSearch context
856: - x - input vector

858:   Output Parameter:
859: + f - Objective value at X
860: - g - Gradient vector at X

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

865:   Level: developer

867: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
868: @*/
869: PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
870: {

880:   if (ls->usetaoroutines) {
881:       TaoComputeObjectiveAndGradient(ls->tao,x,f,g);
882:   } else {
883:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
884:     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
885:     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set");

887:     PetscStackPush("TaoLineSearch user objective/gradient routine");
888:     if (ls->ops->computeobjectiveandgradient) {
889:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);
890:     } else {
891:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
892:       (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);
893:     }
894:     PetscStackPop;
895:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
896:     PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
897:   }
898:   ls->nfgeval++;
899:   return(0);
900: }

902: /*@
903:   TaoLineSearchComputeGradient - Computes the gradient of the objective function

905:   Collective on TaoLineSearch

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

911:   Output Parameter:
912: . g - gradient vector

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

917:   Level: developer

919: .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient()
920: @*/
921: PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
922: {
924:   PetscReal      fdummy;

932:   if (ls->usetaoroutines) {
933:     TaoComputeGradient(ls->tao,x,g);
934:   } else {
935:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
936:     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set");
937:     PetscStackPush("TaoLineSearch user gradient routine");
938:     if (ls->ops->computegradient) {
939:       (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);
940:     } else {
941:       (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);
942:     }
943:     PetscStackPop;
944:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
945:   }
946:   ls->ngeval++;
947:   return(0);
948: }

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

953:   Collective on Tao

955:   Input Parameters:
956: + ls - the TaoLineSearch context
957: - x - input vector

959:   Output Parameter:
960: + f - Objective value at X
961: - gts - inner product of gradient and step direction at X

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

966:   Level: developer

968: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
969: @*/
970: PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
971: {
979:   PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
980:   if (!ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set");
981:   PetscStackPush("TaoLineSearch user objective/gts routine");
982:   (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);
983:   PetscStackPop;
984:   PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
985:   PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
986:   ls->nfeval++;
987:   return(0);
988: }

990: /*@
991:   TaoLineSearchGetSolution - Returns the solution to the line search

993:   Collective on TaoLineSearch

995:   Input Parameter:
996: . ls - the TaoLineSearch context

998:   Output Parameter:
999: + x - the new solution
1000: . f - the objective function value at x
1001: . g - the gradient at x
1002: . steplength - the multiple of the step direction taken by the line search
1003: - reason - the reason why the line search terminated

1005:   reason will be set to one of:

1007: + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
1008: . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
1009: . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
1010: . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
1011: . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
1012: . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
1013: . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance

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

1018: + TAOLINESEARCH_SUCCESS - successful line search

1020:   Level: developer

1022: @*/
1023: PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
1024: {

1033:   if (ls->new_x) {
1034:     VecCopy(ls->new_x,x);
1035:   }
1036:   *f = ls->new_f;
1037:   if (ls->new_g) {
1038:     VecCopy(ls->new_g,g);
1039:   }
1040:   if (steplength) {
1041:     *steplength=ls->step;
1042:   }
1043:   *reason = ls->reason;
1044:   return(0);
1045: }

1047: /*@
1048:   TaoLineSearchGetStartingVector - Gets a the initial point of the line
1049:   search.

1051:   Not Collective

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

1056:   Output Parameter:
1057: . x - The initial point of the line search

1059:   Level: intermediate
1060: @*/
1061: PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1062: {
1065:   if (x) {
1066:     *x = ls->start_x;
1067:   }
1068:   return(0);
1069: }

1071: /*@
1072:   TaoLineSearchGetStepDirection - Gets the step direction of the line
1073:   search.

1075:   Not Collective

1077:   Input Parameter:
1078: . ls - the TaoLineSearch context

1080:   Output Parameter:
1081: . s - the step direction of the line search

1083:   Level: advanced
1084: @*/
1085: PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1086: {
1089:   if (s) {
1090:     *s = ls->stepdirection;
1091:   }
1092:   return(0);

1094: }

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

1099:   Not Collective

1101:   Input Parameter:
1102: . ls - the TaoLineSearch context

1104:   Output Parameter:
1105: . f - the objective value at the full step length

1107:   Level: developer
1108: @*/

1110: PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1111: {
1114:   *f_fullstep = ls->f_fullstep;
1115:   return(0);
1116: }

1118: /*@
1119:   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.

1121:   Logically Collective on Tao

1123:   Input Parameters:
1124: + ls - the TaoLineSearch context
1125: . xl  - vector of lower bounds
1126: - xu  - vector of upper bounds

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

1131:   Level: beginner

1133: .seealso: TaoSetVariableBounds(), TaoLineSearchCreate()
1134: @*/
1135: PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
1136: {
1141:   ls->lower = xl;
1142:   ls->upper = xu;
1143:   ls->bounded = 1;
1144:   return(0);
1145: }

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

1151:   Logically Collective on TaoLineSearch

1153:   Input Parameters:
1154: + ls - the TaoLineSearch context
1155: - s - the initial step size

1157:   Level: intermediate

1159: .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply()
1160: @*/
1161: PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
1162: {
1165:   ls->initstep = s;
1166:   return(0);
1167: }

1169: /*@
1170:   TaoLineSearchGetStepLength - Get the current step length

1172:   Not Collective

1174:   Input Parameters:
1175: . ls - the TaoLineSearch context

1177:   Output Parameters
1178: . s - the current step length

1180:   Level: beginner

1182: .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply()
1183: @*/
1184: PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s)
1185: {
1188:   *s = ls->step;
1189:   return(0);
1190: }

1192: /*MC
1193:    TaoLineSearchRegister - Adds a line-search algorithm to the registry

1195:    Not collective

1197:    Input Parameters:
1198: +  sname - name of a new user-defined solver
1199: -  func - routine to Create method context

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

1204:    Sample usage:
1205: .vb
1206:    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
1207: .ve

1209:    Then, your solver can be chosen with the procedural interface via
1210: $     TaoLineSearchSetType(ls,"my_linesearch")
1211:    or at runtime via the option
1212: $     -tao_ls_type my_linesearch

1214:    Level: developer

1216: .seealso: TaoLineSearchRegisterDestroy()
1217: M*/
1218: PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1219: {
1222:   PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);
1223:   return(0);
1224: }

1226: /*@C
1227:    TaoLineSearchRegisterDestroy - Frees the list of line-search algorithms that were
1228:    registered by TaoLineSearchRegister().

1230:    Not Collective

1232:    Level: developer

1234: .seealso: TaoLineSearchRegister()
1235: @*/
1236: PetscErrorCode TaoLineSearchRegisterDestroy(void)
1237: {
1240:   PetscFunctionListDestroy(&TaoLineSearchList);
1241:   TaoLineSearchInitialized = PETSC_FALSE;
1242:   return(0);
1243: }

1245: /*@C
1246:    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1247:    for all TaoLineSearch options in the database.


1250:    Collective on TaoLineSearch

1252:    Input Parameters:
1253: +  ls - the TaoLineSearch solver context
1254: -  prefix - the prefix string to prepend to all line search requests

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


1261:    Level: advanced

1263: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1264: @*/
1265: PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1266: {
1267:   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
1268: }

1270: /*@C
1271:   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1272:   TaoLineSearch options in the database

1274:   Not Collective

1276:   Input Parameters:
1277: . ls - the TaoLineSearch context

1279:   Output Parameters:
1280: . prefix - pointer to the prefix string used is returned

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

1285:   Level: advanced

1287: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix()
1288: @*/
1289: PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1290: {
1291:   return PetscObjectGetOptionsPrefix((PetscObject)ls,p);
1292: }

1294: /*@C
1295:    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1296:    TaoLineSearch options in the database.


1299:    Logically Collective on TaoLineSearch

1301:    Input Parameters:
1302: +  ls - the TaoLineSearch context
1303: -  prefix - the prefix string to prepend to all TAO option requests

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

1309:    For example, to distinguish between the runtime options for two
1310:    different line searches, one could call
1311: .vb
1312:       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1313:       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1314: .ve

1316:    This would enable use of different options for each system, such as
1317: .vb
1318:       -sys1_tao_ls_type mt
1319:       -sys2_tao_ls_type armijo
1320: .ve

1322:    Level: advanced

1324: .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1325: @*/

1327: PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1328: {
1329:   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
1330: }