Actual source code: taolinesearch.c

petsc-3.7.5 2017-01-01
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;

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

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

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

 83:   } else if (isstring) {
 84:     TaoLineSearchGetType(ls,&type);
 85:     PetscViewerStringSPrintf(viewer," %-3.3s",type);
 86:   }
 87:   return(0);
 88: }

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

 96:   Collective on MPI_Comm

 98:   Input Parameter:
 99: . comm - MPI communicator

101:   Output Parameter:
102: . newls - the new TaoLineSearch context

104:   Available methods include:
105: + more-thuente
106: . gpcg
107: - unit - Do not perform any line search


110:    Options Database Keys:
111: .   -tao_ls_type - select which method TAO should use

113:    Level: beginner

115: .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy()
116: @*/

118: PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
119: {
121:   TaoLineSearch  ls;

125:   *newls = NULL;

127:  #ifndef PETSC_USE_DYNAMIC_LIBRARIES
128:   TaoLineSearchInitializePackage();
129:  #endif

131:   PetscHeaderCreate(ls,TAOLINESEARCH_CLASSID,"TaoLineSearch","Linesearch","Tao",comm,TaoLineSearchDestroy,TaoLineSearchView);
132:   ls->bounded = 0;
133:   ls->max_funcs=30;
134:   ls->ftol = 0.0001;
135:   ls->gtol = 0.9;
136: #if defined(PETSC_USE_REAL_SINGLE)
137:   ls->rtol = 1.0e-5;
138: #else
139:   ls->rtol = 1.0e-10;
140: #endif
141:   ls->stepmin=1.0e-20;
142:   ls->stepmax=1.0e+20;
143:   ls->step=1.0;
144:   ls->nfeval=0;
145:   ls->ngeval=0;
146:   ls->nfgeval=0;

148:   ls->ops->computeobjective=0;
149:   ls->ops->computegradient=0;
150:   ls->ops->computeobjectiveandgradient=0;
151:   ls->ops->computeobjectiveandgts=0;
152:   ls->ops->setup=0;
153:   ls->ops->apply=0;
154:   ls->ops->view=0;
155:   ls->ops->setfromoptions=0;
156:   ls->ops->reset=0;
157:   ls->ops->destroy=0;
158:   ls->setupcalled=PETSC_FALSE;
159:   ls->usetaoroutines=PETSC_FALSE;
160:   *newls = ls;
161:   return(0);
162: }

166: /*@
167:   TaoLineSearchSetUp - Sets up the internal data structures for the later use
168:   of a Tao solver

170:   Collective on ls

172:   Input Parameters:
173: . ls - the TaoLineSearch context

175:   Notes:
176:   The user will not need to explicitly call TaoLineSearchSetUp(), as it will
177:   automatically be called in TaoLineSearchSolve().  However, if the user
178:   desires to call it explicitly, it should come after TaoLineSearchCreate()
179:   but before TaoLineSearchApply().

181:   Level: developer

183: .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
184: @*/

186: PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
187: {
189:   const char     *default_type=TAOLINESEARCHMT;
190:   PetscBool      flg;

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

231: /*@
232:   TaoLineSearchReset - Some line searches may carry state information
233:   from one TaoLineSearchApply() to the next.  This function resets this
234:   state information.

236:   Collective on TaoLineSearch

238:   Input Parameter:
239: . ls - the TaoLineSearch context

241:   Level: developer

243: .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
244: @*/
245: PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
246: {

251:   if (ls->ops->reset) {
252:     (*ls->ops->reset)(ls);
253:   }
254:   return(0);
255: }

259: /*@
260:   TaoLineSearchDestroy - Destroys the TAO context that was created with
261:   TaoLineSearchCreate()

263:   Collective on TaoLineSearch

265:   Input Parameter
266: . ls - the TaoLineSearch context

268:   Level: beginner

270: .seealse: TaoLineSearchCreate(), TaoLineSearchSolve()
271: @*/
272: PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
273: {

277:   if (!*ls) return(0);
279:   if (--((PetscObject)*ls)->refct > 0) {*ls=0; return(0);}
280:   VecDestroy(&(*ls)->stepdirection);
281:   VecDestroy(&(*ls)->start_x);
282:   if ((*ls)->ops->destroy) {
283:     (*(*ls)->ops->destroy)(*ls);
284:   }
285:   PetscHeaderDestroy(ls);
286:   return(0);
287: }

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

294:   Collective on TaoLineSearch

296:   Input Parameters:
297: + ls - the Tao context
298: . x - The current solution (on output x contains the new solution determined by the line search)
299: . f - objective function value at current solution (on output contains the objective function value at new solution)
300: . g - gradient evaluated at x (on output contains the gradient at new solution)
301: - s - search direction

303:   Output Parameters:
304: + x - new solution
305: . f - objective function value at x
306: . g - gradient vector at x
307: . steplength - scalar multiplier of s used ( x = x0 + steplength * x )
308: - reason - reason why the line-search stopped

310:   reason will be set to one of:

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

323:   Note:
324:   The algorithm developer must set up the TaoLineSearch with calls to
325:   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines()

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

332:   Level: beginner

334:   .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts()
335:  @*/

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

343:   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
353:   VecGetOwnershipRange(x, &low1, &high1);
354:   VecGetOwnershipRange(g, &low2, &high2);
355:   VecGetOwnershipRange(s, &low3, &high3);
356:   if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths");

358:   PetscObjectReference((PetscObject)s);
359:   VecDestroy(&ls->stepdirection);
360:   ls->stepdirection = s;

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

397:   PetscObjectReference((PetscObject)x);
398:   VecDestroy(&ls->start_x);
399:   ls->start_x = x;

401:   PetscLogEventBegin(TaoLineSearch_ApplyEvent,ls,0,0,0);
402:   (*ls->ops->apply)(ls,x,f,g,s);
403:   PetscLogEventEnd(TaoLineSearch_ApplyEvent, ls, 0,0,0);
404:   *reason=ls->reason;
405:   ls->new_f = *f;

407:   if (steplength) {
408:     *steplength=ls->step;
409:   }

411:   TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");
412:   return(0);
413: }

417: /*@C
418:    TaoLineSearchSetType - Sets the algorithm used in a line search

420:    Collective on TaoLineSearch

422:    Input Parameters:
423: +  ls - the TaoLineSearch context
424: -  type - a known method

426:   Available methods include:
427: + more-thuente
428: . gpcg
429: - unit - Do not perform any line search


432:   Options Database Keys:
433: .   -tao_ls_type - select which method TAO should use

435:   Level: beginner


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

440: @*/

442: PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, const TaoLineSearchType type)
443: {
445:   PetscErrorCode (*r)(TaoLineSearch);
446:   PetscBool      flg;

451:   PetscObjectTypeCompare((PetscObject)ls, type, &flg);
452:   if (flg) return(0);

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

470:   ls->nfeval=0;
471:   ls->ngeval=0;
472:   ls->nfgeval=0;
473:   ls->ops->setup=0;
474:   ls->ops->apply=0;
475:   ls->ops->view=0;
476:   ls->ops->setfromoptions=0;
477:   ls->ops->destroy=0;
478:   ls->setupcalled = PETSC_FALSE;
479:   (*r)(ls);
480:   PetscObjectChangeTypeName((PetscObject)ls, type);
481:   return(0);
482: }

486: /*@
487:   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
488:   options.

490:   Collective on TaoLineSearch

492:   Input Paremeter:
493: . ls - the TaoLineSearch context

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

505:   Level: beginner
506: @*/
507: PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
508: {
510:   const char     *default_type=TAOLINESEARCHMT;
511:   char           type[256];
512:   PetscBool      flg;

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

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

546: /*@C
547:   TaoLineSearchGetType - Gets the current line search algorithm

549:   Not Collective

551:   Input Parameter:
552: . ls - the TaoLineSearch context

554:   Output Paramter:
555: . type - the line search algorithm in effect

557:   Level: developer

559: @*/
560: PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, const TaoLineSearchType *type)
561: {
565:   *type = ((PetscObject)ls)->type_name;
566:   return(0);
567: }

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

575:   Not Collective

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

580:   Output Parameters:
581: + nfeval   - number of function evaluations
582: . ngeval   - number of gradient evaluations
583: - nfgeval  - number of function/gradient evaluations

585:   Level: intermediate

587:   Note:
588:   If the line search is using the Tao objective and gradient
589:   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
590:   is already counting the number of evaluations.

592: @*/
593: PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
594: {
597:   *nfeval = ls->nfeval;
598:   *ngeval = ls->ngeval;
599:   *nfgeval = ls->nfgeval;
600:   return(0);
601: }

605: /*@
606:   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
607:   Tao evaluation routines.

609:   Not Collective

611:   Input Parameter:
612: . ls - the TaoLineSearch context

614:   Output Parameter:
615: . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
616:         otherwise PETSC_FALSE

618:   Level: developer
619: @*/
620: PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
621: {
624:   *flg = ls->usetaoroutines;
625:   return(0);
626: }

630: /*@C
631:   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search

633:   Logically Collective on TaoLineSearch

635:   Input Parameter:
636: + ls - the TaoLineSearch context
637: . func - the objective function evaluation routine
638: - ctx - the (optional) user-defined context for private data

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

643: + x - input vector
644: . f - function value
645: - ctx (optional) user-defined context

647:   Level: beginner

649:   Note:
650:   Use this routine only if you want the line search objective
651:   evaluation routine to be different from the Tao's objective
652:   evaluation routine. If you use this routine you must also set
653:   the line search gradient and/or function/gradient routine.

655:   Note:
656:   Some algorithms (lcl, gpcg) set their own objective routine for the
657:   line search, application programmers should be wary of overriding the
658:   default objective routine.

660: .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
661: @*/
662: PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
663: {

667:   ls->ops->computeobjective=func;
668:   if (ctx) ls->userctx_func=ctx;
669:   ls->usetaoroutines=PETSC_FALSE;
670:   return(0);
671: }

675: /*@C
676:   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search

678:   Logically Collective on TaoLineSearch

680:   Input Parameter:
681: + ls - the TaoLineSearch context
682: . func - the gradient evaluation routine
683: - ctx - the (optional) user-defined context for private data

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

688: + x - input vector
689: . g - gradient vector
690: - ctx (optional) user-defined context

692:   Level: beginner

694:   Note:
695:   Use this routine only if you want the line search gradient
696:   evaluation routine to be different from the Tao's gradient
697:   evaluation routine. If you use this routine you must also set
698:   the line search function and/or function/gradient routine.

700:   Note:
701:   Some algorithms (lcl, gpcg) set their own gradient routine for the
702:   line search, application programmers should be wary of overriding the
703:   default gradient routine.

705: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
706: @*/
707: PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
708: {
711:   ls->ops->computegradient=func;
712:   if (ctx) ls->userctx_grad=ctx;
713:   ls->usetaoroutines=PETSC_FALSE;
714:   return(0);
715: }

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

722:   Logically Collective on TaoLineSearch

724:   Input Parameter:
725: + ls - the TaoLineSearch context
726: . func - the objective and gradient evaluation routine
727: - ctx - the (optional) user-defined context for private data

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

732: + x - input vector
733: . f - function value
734: . g - gradient vector
735: - ctx (optional) user-defined context

737:   Level: beginner

739:   Note:
740:   Use this routine only if you want the line search objective and gradient
741:   evaluation routines to be different from the Tao's objective
742:   and gradient evaluation routines.

744:   Note:
745:   Some algorithms (lcl, gpcg) set their own objective routine for the
746:   line search, application programmers should be wary of overriding the
747:   default objective routine.

749: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines()
750: @*/
751: PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx)
752: {
755:   ls->ops->computeobjectiveandgradient=func;
756:   if (ctx) ls->userctx_funcgrad=ctx;
757:   ls->usetaoroutines = PETSC_FALSE;
758:   return(0);
759: }

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

770:   Logically Collective on TaoLineSearch

772:   Input Parameter:
773: + ls - the TaoLineSearch context
774: . func - the objective and gradient evaluation routine
775: - ctx - the (optional) user-defined context for private data

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

780: + x - input vector
781: . s - step direction
782: . f - function value
783: . gts - inner product of gradient and step direction vectors
784: - ctx (optional) user-defined context

786:   Note: The gradient will still need to be computed at the end of the line
787:   search, so you will still need to set a line search gradient evaluation
788:   routine

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

794:   Level: advanced

796:   Note:
797:   Some algorithms (lcl, gpcg) set their own objective routine for the
798:   line search, application programmers should be wary of overriding the
799:   default objective routine.

801: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines()
802: @*/
803: PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx)
804: {
807:   ls->ops->computeobjectiveandgts=func;
808:   if (ctx) ls->userctx_funcgts=ctx;
809:   ls->usegts = PETSC_TRUE;
810:   ls->usetaoroutines=PETSC_FALSE;
811:   return(0);
812: }

816: /*@C
817:   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
818:   objective and gradient evaluation routines from the given Tao object.

820:   Logically Collective on TaoLineSearch

822:   Input Parameter:
823: + ls - the TaoLineSearch context
824: - ts - the Tao context with defined objective/gradient evaluation routines

826:   Level: developer

828: .seealso: TaoLineSearchCreate()
829: @*/
830: PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
831: {
835:   ls->tao = ts;
836:   ls->usetaoroutines=PETSC_TRUE;
837:   return(0);
838: }

842: /*@
843:   TaoLineSearchComputeObjective - Computes the objective function value at a given point

845:   Collective on TaoLineSearch

847:   Input Parameters:
848: + ls - the TaoLineSearch context
849: - x - input vector

851:   Output Parameter:
852: . f - Objective value at X

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

857:   Level: developer

859: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
860: @*/
861: PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
862: {
864:   Vec            gdummy;
865:   PetscReal      gts;

872:   if (ls->usetaoroutines) {
873:     TaoComputeObjective(ls->tao,x,f);
874:   } else {
875:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
876:     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");
877:     PetscStackPush("TaoLineSearch user objective routine");
878:     if (ls->ops->computeobjective) {
879:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
880:     } else if (ls->ops->computeobjectiveandgradient) {
881:       VecDuplicate(x,&gdummy);
882:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);
883:       VecDestroy(&gdummy);
884:     } else {
885:       (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,&gts,ls->userctx_funcgts);
886:     }
887:     PetscStackPop;
888:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
889:   }
890:   ls->nfeval++;
891:   return(0);
892: }

896: /*@
897:   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point

899:   Collective on Tao

901:   Input Parameters:
902: + ls - the TaoLineSearch context
903: - x - input vector

905:   Output Parameter:
906: + f - Objective value at X
907: - g - Gradient vector at X

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

912:   Level: developer

914: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
915: @*/
916: PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
917: {

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

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

951: /*@
952:   TaoLineSearchComputeGradient - Computes the gradient of the objective function

954:   Collective on TaoLineSearch

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

960:   Output Parameter:
961: . g - gradient vector

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

966:   Level: developer

968: .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient()
969: @*/
970: PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
971: {
973:   PetscReal      fdummy;

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

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

1004:   Collective on Tao

1006:   Input Parameters:
1007: + ls - the TaoLineSearch context
1008: - x - input vector

1010:   Output Parameter:
1011: + f - Objective value at X
1012: - gts - inner product of gradient and step direction at X

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

1017:   Level: developer

1019: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
1020: @*/
1021: PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
1022: {
1030:   PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
1031:   if (!ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set");
1032:   PetscStackPush("TaoLineSearch user objective/gts routine");
1033:   (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);
1034:   PetscStackPop;
1035:   PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
1036:   PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
1037:   ls->nfeval++;
1038:   return(0);
1039: }

1043: /*@
1044:   TaoLineSearchGetSolution - Returns the solution to the line search

1046:   Collective on TaoLineSearch

1048:   Input Parameter:
1049: . ls - the TaoLineSearch context

1051:   Output Parameter:
1052: + x - the new solution
1053: . f - the objective function value at x
1054: . g - the gradient at x
1055: . steplength - the multiple of the step direction taken by the line search
1056: - reason - the reason why the line search terminated

1058:   reason will be set to one of:

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

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

1071: + TAOLINESEARCH_SUCCESS - successful line search

1073:   Level: developer

1075: @*/
1076: PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
1077: {

1086:   if (ls->new_x) {
1087:     VecCopy(ls->new_x,x);
1088:   }
1089:   *f = ls->new_f;
1090:   if (ls->new_g) {
1091:     VecCopy(ls->new_g,g);
1092:   }
1093:   if (steplength) {
1094:     *steplength=ls->step;
1095:   }
1096:   *reason = ls->reason;
1097:   return(0);
1098: }

1102: /*@
1103:   TaoLineSearchGetStartingVector - Gets a the initial point of the line
1104:   search.

1106:   Not Collective

1108:   Input Parameter:
1109: . ls - the TaoLineSearch context

1111:   Output Parameter:
1112: . x - The initial point of the line search

1114:   Level: intermediate
1115: @*/
1116: PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1117: {
1120:   if (x) {
1121:     *x = ls->start_x;
1122:   }
1123:   return(0);
1124: }

1128: /*@
1129:   TaoLineSearchGetStepDirection - Gets the step direction of the line
1130:   search.

1132:   Not Collective

1134:   Input Parameter:
1135: . ls - the TaoLineSearch context

1137:   Output Parameter:
1138: . s - the step direction of the line search

1140:   Level: advanced
1141: @*/
1142: PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1143: {
1146:   if (s) {
1147:     *s = ls->stepdirection;
1148:   }
1149:   return(0);

1151: }

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

1158:   Not Collective

1160:   Input Parameter:
1161: . ls - the TaoLineSearch context

1163:   Output Parameter:
1164: . f - the objective value at the full step length

1166:   Level: developer
1167: @*/

1169: PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1170: {
1173:   *f_fullstep = ls->f_fullstep;
1174:   return(0);
1175: }

1179: /*@
1180:   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.

1182:   Logically Collective on Tao

1184:   Input Parameters:
1185: + ls - the TaoLineSearch context
1186: . xl  - vector of lower bounds
1187: - xu  - vector of upper bounds

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

1192:   Level: beginner

1194: .seealso: TaoSetVariableBounds(), TaoLineSearchCreate()
1195: @*/
1196: PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
1197: {
1202:   ls->lower = xl;
1203:   ls->upper = xu;
1204:   ls->bounded = 1;
1205:   return(0);
1206: }

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

1214:   Logically Collective on TaoLineSearch

1216:   Input Parameters:
1217: + ls - the TaoLineSearch context
1218: - s - the initial step size

1220:   Level: intermediate

1222: .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply()
1223: @*/
1224: PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
1225: {
1228:   ls->initstep = s;
1229:   return(0);
1230: }

1234: /*@
1235:   TaoLineSearchGetStepLength - Get the current step length

1237:   Not Collective

1239:   Input Parameters:
1240: . ls - the TaoLineSearch context

1242:   Output Parameters
1243: . s - the current step length

1245:   Level: beginner

1247: .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply()
1248: @*/
1249: PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s)
1250: {
1253:   *s = ls->step;
1254:   return(0);
1255: }

1259: /*MC
1260:    TaoLineSearchRegister - Adds a line-search algorithm to the registry

1262:    Not collective

1264:    Input Parameters:
1265: +  sname - name of a new user-defined solver
1266: -  func - routine to Create method context

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

1271:    Sample usage:
1272: .vb
1273:    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
1274: .ve

1276:    Then, your solver can be chosen with the procedural interface via
1277: $     TaoLineSearchSetType(ls,"my_linesearch")
1278:    or at runtime via the option
1279: $     -tao_ls_type my_linesearch

1281:    Level: developer

1283: .seealso: TaoLineSearchRegisterDestroy()
1284: M*/
1285: PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1286: {
1289:   PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);
1290:   return(0);
1291: }

1295: /*@C
1296:    TaoLineSearchRegisterDestroy - Frees the list of line-search algorithms that were
1297:    registered by TaoLineSearchRegister().

1299:    Not Collective

1301:    Level: developer

1303: .seealso: TaoLineSearchRegister()
1304: @*/
1305: PetscErrorCode TaoLineSearchRegisterDestroy(void)
1306: {
1309:   PetscFunctionListDestroy(&TaoLineSearchList);
1310:   TaoLineSearchInitialized = PETSC_FALSE;
1311:   return(0);
1312: }

1316: /*@C
1317:    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1318:    for all TaoLineSearch options in the database.


1321:    Collective on TaoLineSearch

1323:    Input Parameters:
1324: +  ls - the TaoLineSearch solver context
1325: -  prefix - the prefix string to prepend to all line search requests

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


1332:    Level: advanced

1334: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1335: @*/
1336: PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1337: {
1338:   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
1339: }

1343: /*@C
1344:   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1345:   TaoLineSearch options in the database

1347:   Not Collective

1349:   Input Parameters:
1350: . ls - the TaoLineSearch context

1352:   Output Parameters:
1353: . prefix - pointer to the prefix string used is returned

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

1358:   Level: advanced

1360: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix()
1361: @*/
1362: PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1363: {
1364:   return PetscObjectGetOptionsPrefix((PetscObject)ls,p);
1365: }

1369: /*@C
1370:    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1371:    TaoLineSearch options in the database.


1374:    Logically Collective on TaoLineSearch

1376:    Input Parameters:
1377: +  ls - the TaoLineSearch context
1378: -  prefix - the prefix string to prepend to all TAO option requests

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

1384:    For example, to distinguish between the runtime options for two
1385:    different line searches, one could call
1386: .vb
1387:       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1388:       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1389: .ve

1391:    This would enable use of different options for each system, such as
1392: .vb
1393:       -sys1_tao_ls_type mt
1394:       -sys2_tao_ls_type armijo
1395: .ve

1397:    Level: advanced

1399: .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1400: @*/

1402: PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1403: {
1404:   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
1405: }