Actual source code: taolinesearch.c

petsc-3.8.4 2018-03-24
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:     if (((PetscObject)ls)->prefix) {
 53:       PetscViewerASCIIPrintf(viewer,"TaoLineSearch Object:(%s)\n",((PetscObject)ls)->prefix);
 54:     } else {
 55:       PetscViewerASCIIPrintf(viewer,"TaoLineSearch Object:\n");
 56:     }
 57:     PetscViewerASCIIPushTab(viewer);
 58:     TaoLineSearchGetType(ls,&type);
 59:     if (type) {
 60:       PetscViewerASCIIPrintf(viewer,"type: %s\n",type);
 61:     } else {
 62:       PetscViewerASCIIPrintf(viewer,"type: not set yet\n");
 63:     }
 64:     if (ls->ops->view) {
 65:       PetscViewerASCIIPushTab(viewer);
 66:       (*ls->ops->view)(ls,viewer);
 67:       PetscViewerASCIIPopTab(viewer);
 68:     }
 69:     PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);
 70:     PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);
 71:     PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);
 72:     PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);
 73:     PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);

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

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

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

 92:   Collective on MPI_Comm

 94:   Input Parameter:
 95: . comm - MPI communicator

 97:   Output Parameter:
 98: . newls - the new TaoLineSearch context

100:   Available methods include:
101: + more-thuente
102: . gpcg
103: - unit - Do not perform any line search


106:    Options Database Keys:
107: .   -tao_ls_type - select which method TAO should use

109:    Level: beginner

111: .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy()
112: @*/

114: PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
115: {
117:   TaoLineSearch  ls;

121:   *newls = NULL;

123:  #ifndef PETSC_USE_DYNAMIC_LIBRARIES
124:   TaoLineSearchInitializePackage();
125:  #endif

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

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

160: /*@
161:   TaoLineSearchSetUp - Sets up the internal data structures for the later use
162:   of a Tao solver

164:   Collective on ls

166:   Input Parameters:
167: . ls - the TaoLineSearch context

169:   Notes:
170:   The user will not need to explicitly call TaoLineSearchSetUp(), as it will
171:   automatically be called in TaoLineSearchSolve().  However, if the user
172:   desires to call it explicitly, it should come after TaoLineSearchCreate()
173:   but before TaoLineSearchApply().

175:   Level: developer

177: .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
178: @*/

180: PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
181: {
183:   const char     *default_type=TAOLINESEARCHMT;
184:   PetscBool      flg;

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

223: /*@
224:   TaoLineSearchReset - Some line searches may carry state information
225:   from one TaoLineSearchApply() to the next.  This function resets this
226:   state information.

228:   Collective on TaoLineSearch

230:   Input Parameter:
231: . ls - the TaoLineSearch context

233:   Level: developer

235: .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
236: @*/
237: PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
238: {

243:   if (ls->ops->reset) {
244:     (*ls->ops->reset)(ls);
245:   }
246:   return(0);
247: }

249: /*@
250:   TaoLineSearchDestroy - Destroys the TAO context that was created with
251:   TaoLineSearchCreate()

253:   Collective on TaoLineSearch

255:   Input Parameter
256: . ls - the TaoLineSearch context

258:   Level: beginner

260: .seealse: TaoLineSearchCreate(), TaoLineSearchSolve()
261: @*/
262: PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
263: {

267:   if (!*ls) return(0);
269:   if (--((PetscObject)*ls)->refct > 0) {*ls=0; return(0);}
270:   VecDestroy(&(*ls)->stepdirection);
271:   VecDestroy(&(*ls)->start_x);
272:   if ((*ls)->ops->destroy) {
273:     (*(*ls)->ops->destroy)(*ls);
274:   }
275:   PetscHeaderDestroy(ls);
276:   return(0);
277: }

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

282:   Collective on TaoLineSearch

284:   Input Parameters:
285: + ls - the Tao context
286: . x - The current solution (on output x contains the new solution determined by the line search)
287: . f - objective function value at current solution (on output contains the objective function value at new solution)
288: . g - gradient evaluated at x (on output contains the gradient at new solution)
289: - s - search direction

291:   Output Parameters:
292: + x - new solution
293: . f - objective function value at x
294: . g - gradient vector at x
295: . steplength - scalar multiplier of s used ( x = x0 + steplength * x )
296: - reason - reason why the line-search stopped

298:   reason will be set to one of:

300: + TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
301: . TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
302: . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
303: . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
304: . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
305: . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
306: . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance
307: . TAOLINESEARCH_HALTED_USER - user can set this reason to stop line search
308: . TAOLINESEARCH_HALTED_OTHER - any other reason
309: - TAOLINESEARCH_SUCCESS - successful line search

311:   Note:
312:   The algorithm developer must set up the TaoLineSearch with calls to
313:   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines()

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

320:   Level: beginner

322:   .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts()
323:  @*/

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

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

346:   PetscObjectReference((PetscObject)s);
347:   VecDestroy(&ls->stepdirection);
348:   ls->stepdirection = s;

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

385:   PetscObjectReference((PetscObject)x);
386:   VecDestroy(&ls->start_x);
387:   ls->start_x = x;

389:   PetscLogEventBegin(TaoLineSearch_ApplyEvent,ls,0,0,0);
390:   (*ls->ops->apply)(ls,x,f,g,s);
391:   PetscLogEventEnd(TaoLineSearch_ApplyEvent, ls, 0,0,0);
392:   *reason=ls->reason;
393:   ls->new_f = *f;

395:   if (steplength) {
396:     *steplength=ls->step;
397:   }

399:   TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");
400:   return(0);
401: }

403: /*@C
404:    TaoLineSearchSetType - Sets the algorithm used in a line search

406:    Collective on TaoLineSearch

408:    Input Parameters:
409: +  ls - the TaoLineSearch context
410: -  type - a known method

412:   Available methods include:
413: + more-thuente
414: . gpcg
415: - unit - Do not perform any line search


418:   Options Database Keys:
419: .   -tao_ls_type - select which method TAO should use

421:   Level: beginner


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

426: @*/

428: PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, const TaoLineSearchType type)
429: {
431:   PetscErrorCode (*r)(TaoLineSearch);
432:   PetscBool      flg;

437:   PetscObjectTypeCompare((PetscObject)ls, type, &flg);
438:   if (flg) return(0);

440:   PetscFunctionListFind(TaoLineSearchList,type, (void (**)(void)) &r);
441:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested TaoLineSearch type %s",type);
442:   if (ls->ops->destroy) {
443:     (*(ls)->ops->destroy)(ls);
444:   }
445:   ls->max_funcs=30;
446:   ls->ftol = 0.0001;
447:   ls->gtol = 0.9;
448: #if defined(PETSC_USE_REAL_SINGLE)
449:   ls->rtol = 1.0e-5;
450: #else
451:   ls->rtol = 1.0e-10;
452: #endif
453:   ls->stepmin=1.0e-20;
454:   ls->stepmax=1.0e+20;

456:   ls->nfeval=0;
457:   ls->ngeval=0;
458:   ls->nfgeval=0;
459:   ls->ops->setup=0;
460:   ls->ops->apply=0;
461:   ls->ops->view=0;
462:   ls->ops->setfromoptions=0;
463:   ls->ops->destroy=0;
464:   ls->setupcalled = PETSC_FALSE;
465:   (*r)(ls);
466:   PetscObjectChangeTypeName((PetscObject)ls, type);
467:   return(0);
468: }

470: /*@
471:   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
472:   options.

474:   Collective on TaoLineSearch

476:   Input Paremeter:
477: . ls - the TaoLineSearch context

479:   Options Database Keys:
480: + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
481: . -tao_ls_ftol <tol> - tolerance for sufficient decrease
482: . -tao_ls_gtol <tol> - tolerance for curvature condition
483: . -tao_ls_rtol <tol> - relative tolerance for acceptable step
484: . -tao_ls_stepmin <step> - minimum steplength allowed
485: . -tao_ls_stepmax <step> - maximum steplength allowed
486: . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
487: - -tao_ls_view - display line-search results to standard output

489:   Level: beginner
490: @*/
491: PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
492: {
494:   const char     *default_type=TAOLINESEARCHMT;
495:   char           type[256];
496:   PetscBool      flg;

500:   PetscObjectOptionsBegin((PetscObject)ls);
501:   if (!TaoLineSearchInitialized) {
502:     TaoLineSearchInitializePackage();
503:   }
504:   if (((PetscObject)ls)->type_name) {
505:     default_type = ((PetscObject)ls)->type_name;
506:   }
507:   /* Check for type from options */
508:   PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);
509:   if (flg) {
510:     TaoLineSearchSetType(ls,type);
511:   } else if (!((PetscObject)ls)->type_name) {
512:     TaoLineSearchSetType(ls,default_type);
513:   }

515:   PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL);
516:   PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL);
517:   PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL);
518:   PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL);
519:   PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL);
520:   PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL);
521:   if (ls->ops->setfromoptions) {
522:     (*ls->ops->setfromoptions)(PetscOptionsObject,ls);
523:   }
524:   PetscOptionsEnd();
525:   return(0);
526: }

528: /*@C
529:   TaoLineSearchGetType - Gets the current line search algorithm

531:   Not Collective

533:   Input Parameter:
534: . ls - the TaoLineSearch context

536:   Output Paramter:
537: . type - the line search algorithm in effect

539:   Level: developer

541: @*/
542: PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, const TaoLineSearchType *type)
543: {
547:   *type = ((PetscObject)ls)->type_name;
548:   return(0);
549: }

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

555:   Not Collective

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

560:   Output Parameters:
561: + nfeval   - number of function evaluations
562: . ngeval   - number of gradient evaluations
563: - nfgeval  - number of function/gradient evaluations

565:   Level: intermediate

567:   Note:
568:   If the line search is using the Tao objective and gradient
569:   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
570:   is already counting the number of evaluations.

572: @*/
573: PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
574: {
577:   *nfeval = ls->nfeval;
578:   *ngeval = ls->ngeval;
579:   *nfgeval = ls->nfgeval;
580:   return(0);
581: }

583: /*@
584:   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
585:   Tao evaluation routines.

587:   Not Collective

589:   Input Parameter:
590: . ls - the TaoLineSearch context

592:   Output Parameter:
593: . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
594:         otherwise PETSC_FALSE

596:   Level: developer
597: @*/
598: PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
599: {
602:   *flg = ls->usetaoroutines;
603:   return(0);
604: }

606: /*@C
607:   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search

609:   Logically Collective on TaoLineSearch

611:   Input Parameter:
612: + ls - the TaoLineSearch context
613: . func - the objective function evaluation routine
614: - ctx - the (optional) user-defined context for private data

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

619: + x - input vector
620: . f - function value
621: - ctx (optional) user-defined context

623:   Level: beginner

625:   Note:
626:   Use this routine only if you want the line search objective
627:   evaluation routine to be different from the Tao's objective
628:   evaluation routine. If you use this routine you must also set
629:   the line search gradient and/or function/gradient routine.

631:   Note:
632:   Some algorithms (lcl, gpcg) set their own objective routine for the
633:   line search, application programmers should be wary of overriding the
634:   default objective routine.

636: .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
637: @*/
638: PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
639: {

643:   ls->ops->computeobjective=func;
644:   if (ctx) ls->userctx_func=ctx;
645:   ls->usetaoroutines=PETSC_FALSE;
646:   return(0);
647: }

649: /*@C
650:   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search

652:   Logically Collective on TaoLineSearch

654:   Input Parameter:
655: + ls - the TaoLineSearch context
656: . func - the gradient evaluation routine
657: - ctx - the (optional) user-defined context for private data

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

662: + x - input vector
663: . g - gradient vector
664: - ctx (optional) user-defined context

666:   Level: beginner

668:   Note:
669:   Use this routine only if you want the line search gradient
670:   evaluation routine to be different from the Tao's gradient
671:   evaluation routine. If you use this routine you must also set
672:   the line search function and/or function/gradient routine.

674:   Note:
675:   Some algorithms (lcl, gpcg) set their own gradient routine for the
676:   line search, application programmers should be wary of overriding the
677:   default gradient routine.

679: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
680: @*/
681: PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
682: {
685:   ls->ops->computegradient=func;
686:   if (ctx) ls->userctx_grad=ctx;
687:   ls->usetaoroutines=PETSC_FALSE;
688:   return(0);
689: }

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

694:   Logically Collective on TaoLineSearch

696:   Input Parameter:
697: + ls - the TaoLineSearch context
698: . func - the objective and gradient evaluation routine
699: - ctx - the (optional) user-defined context for private data

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

704: + x - input vector
705: . f - function value
706: . g - gradient vector
707: - ctx (optional) user-defined context

709:   Level: beginner

711:   Note:
712:   Use this routine only if you want the line search objective and gradient
713:   evaluation routines to be different from the Tao's objective
714:   and gradient evaluation routines.

716:   Note:
717:   Some algorithms (lcl, gpcg) set their own objective routine for the
718:   line search, application programmers should be wary of overriding the
719:   default objective routine.

721: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines()
722: @*/
723: PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx)
724: {
727:   ls->ops->computeobjectiveandgradient=func;
728:   if (ctx) ls->userctx_funcgrad=ctx;
729:   ls->usetaoroutines = PETSC_FALSE;
730:   return(0);
731: }

733: /*@C
734:   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
735:   (gradient'*stepdirection) evaluation routine for the line search.
736:   Sometimes it is more efficient to compute the inner product of the gradient
737:   and the step direction than it is to compute the gradient, and this is all
738:   the line search typically needs of the gradient.

740:   Logically Collective on TaoLineSearch

742:   Input Parameter:
743: + ls - the TaoLineSearch context
744: . func - the objective and gradient evaluation routine
745: - ctx - the (optional) user-defined context for private data

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

750: + x - input vector
751: . s - step direction
752: . f - function value
753: . gts - inner product of gradient and step direction vectors
754: - ctx (optional) user-defined context

756:   Note: The gradient will still need to be computed at the end of the line
757:   search, so you will still need to set a line search gradient evaluation
758:   routine

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

764:   Level: advanced

766:   Note:
767:   Some algorithms (lcl, gpcg) set their own objective routine for the
768:   line search, application programmers should be wary of overriding the
769:   default objective routine.

771: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines()
772: @*/
773: PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx)
774: {
777:   ls->ops->computeobjectiveandgts=func;
778:   if (ctx) ls->userctx_funcgts=ctx;
779:   ls->usegts = PETSC_TRUE;
780:   ls->usetaoroutines=PETSC_FALSE;
781:   return(0);
782: }

784: /*@C
785:   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
786:   objective and gradient evaluation routines from the given Tao object.

788:   Logically Collective on TaoLineSearch

790:   Input Parameter:
791: + ls - the TaoLineSearch context
792: - ts - the Tao context with defined objective/gradient evaluation routines

794:   Level: developer

796: .seealso: TaoLineSearchCreate()
797: @*/
798: PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
799: {
803:   ls->tao = ts;
804:   ls->usetaoroutines=PETSC_TRUE;
805:   return(0);
806: }

808: /*@
809:   TaoLineSearchComputeObjective - Computes the objective function value at a given point

811:   Collective on TaoLineSearch

813:   Input Parameters:
814: + ls - the TaoLineSearch context
815: - x - input vector

817:   Output Parameter:
818: . f - Objective value at X

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

823:   Level: developer

825: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
826: @*/
827: PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
828: {
830:   Vec            gdummy;
831:   PetscReal      gts;

838:   if (ls->usetaoroutines) {
839:     TaoComputeObjective(ls->tao,x,f);
840:   } else {
841:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
842:     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");
843:     PetscStackPush("TaoLineSearch user objective routine");
844:     if (ls->ops->computeobjective) {
845:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
846:     } else if (ls->ops->computeobjectiveandgradient) {
847:       VecDuplicate(x,&gdummy);
848:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);
849:       VecDestroy(&gdummy);
850:     } else {
851:       (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,&gts,ls->userctx_funcgts);
852:     }
853:     PetscStackPop;
854:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
855:   }
856:   ls->nfeval++;
857:   return(0);
858: }

860: /*@
861:   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point

863:   Collective on Tao

865:   Input Parameters:
866: + ls - the TaoLineSearch context
867: - x - input vector

869:   Output Parameter:
870: + f - Objective value at X
871: - g - Gradient vector at X

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

876:   Level: developer

878: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
879: @*/
880: PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
881: {

891:   if (ls->usetaoroutines) {
892:       TaoComputeObjectiveAndGradient(ls->tao,x,f,g);
893:   } else {
894:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
895:     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
896:     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set");

898:     PetscStackPush("TaoLineSearch user objective/gradient routine");
899:     if (ls->ops->computeobjectiveandgradient) {
900:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);
901:     } else {
902:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
903:       (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);
904:     }
905:     PetscStackPop;
906:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
907:     PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
908:   }
909:   ls->nfgeval++;
910:   return(0);
911: }

913: /*@
914:   TaoLineSearchComputeGradient - Computes the gradient of the objective function

916:   Collective on TaoLineSearch

918:   Input Parameters:
919: + ls - the TaoLineSearch context
920: - x - input vector

922:   Output Parameter:
923: . g - gradient vector

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

928:   Level: developer

930: .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient()
931: @*/
932: PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
933: {
935:   PetscReal      fdummy;

943:   if (ls->usetaoroutines) {
944:     TaoComputeGradient(ls->tao,x,g);
945:   } else {
946:     PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
947:     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set");
948:     PetscStackPush("TaoLineSearch user gradient routine");
949:     if (ls->ops->computegradient) {
950:       (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);
951:     } else {
952:       (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);
953:     }
954:     PetscStackPop;
955:     PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
956:   }
957:   ls->ngeval++;
958:   return(0);
959: }

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

964:   Collective on Tao

966:   Input Parameters:
967: + ls - the TaoLineSearch context
968: - x - input vector

970:   Output Parameter:
971: + f - Objective value at X
972: - gts - inner product of gradient and step direction at X

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

977:   Level: developer

979: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
980: @*/
981: PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
982: {
990:   PetscLogEventBegin(TaoLineSearch_EvalEvent,ls,0,0,0);
991:   if (!ls->ops->computeobjectiveandgts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set");
992:   PetscStackPush("TaoLineSearch user objective/gts routine");
993:   (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);
994:   PetscStackPop;
995:   PetscLogEventEnd(TaoLineSearch_EvalEvent,ls,0,0,0);
996:   PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
997:   ls->nfeval++;
998:   return(0);
999: }

1001: /*@
1002:   TaoLineSearchGetSolution - Returns the solution to the line search

1004:   Collective on TaoLineSearch

1006:   Input Parameter:
1007: . ls - the TaoLineSearch context

1009:   Output Parameter:
1010: + x - the new solution
1011: . f - the objective function value at x
1012: . g - the gradient at x
1013: . steplength - the multiple of the step direction taken by the line search
1014: - reason - the reason why the line search terminated

1016:   reason will be set to one of:

1018: + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
1019: . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
1020: . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
1021: . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
1022: . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
1023: . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
1024: . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance

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

1029: + TAOLINESEARCH_SUCCESS - successful line search

1031:   Level: developer

1033: @*/
1034: PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
1035: {

1044:   if (ls->new_x) {
1045:     VecCopy(ls->new_x,x);
1046:   }
1047:   *f = ls->new_f;
1048:   if (ls->new_g) {
1049:     VecCopy(ls->new_g,g);
1050:   }
1051:   if (steplength) {
1052:     *steplength=ls->step;
1053:   }
1054:   *reason = ls->reason;
1055:   return(0);
1056: }

1058: /*@
1059:   TaoLineSearchGetStartingVector - Gets a the initial point of the line
1060:   search.

1062:   Not Collective

1064:   Input Parameter:
1065: . ls - the TaoLineSearch context

1067:   Output Parameter:
1068: . x - The initial point of the line search

1070:   Level: intermediate
1071: @*/
1072: PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1073: {
1076:   if (x) {
1077:     *x = ls->start_x;
1078:   }
1079:   return(0);
1080: }

1082: /*@
1083:   TaoLineSearchGetStepDirection - Gets the step direction of the line
1084:   search.

1086:   Not Collective

1088:   Input Parameter:
1089: . ls - the TaoLineSearch context

1091:   Output Parameter:
1092: . s - the step direction of the line search

1094:   Level: advanced
1095: @*/
1096: PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1097: {
1100:   if (s) {
1101:     *s = ls->stepdirection;
1102:   }
1103:   return(0);

1105: }

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

1110:   Not Collective

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

1115:   Output Parameter:
1116: . f - the objective value at the full step length

1118:   Level: developer
1119: @*/

1121: PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1122: {
1125:   *f_fullstep = ls->f_fullstep;
1126:   return(0);
1127: }

1129: /*@
1130:   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.

1132:   Logically Collective on Tao

1134:   Input Parameters:
1135: + ls - the TaoLineSearch context
1136: . xl  - vector of lower bounds
1137: - xu  - vector of upper bounds

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

1142:   Level: beginner

1144: .seealso: TaoSetVariableBounds(), TaoLineSearchCreate()
1145: @*/
1146: PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
1147: {
1152:   ls->lower = xl;
1153:   ls->upper = xu;
1154:   ls->bounded = 1;
1155:   return(0);
1156: }

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

1162:   Logically Collective on TaoLineSearch

1164:   Input Parameters:
1165: + ls - the TaoLineSearch context
1166: - s - the initial step size

1168:   Level: intermediate

1170: .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply()
1171: @*/
1172: PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
1173: {
1176:   ls->initstep = s;
1177:   return(0);
1178: }

1180: /*@
1181:   TaoLineSearchGetStepLength - Get the current step length

1183:   Not Collective

1185:   Input Parameters:
1186: . ls - the TaoLineSearch context

1188:   Output Parameters
1189: . s - the current step length

1191:   Level: beginner

1193: .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply()
1194: @*/
1195: PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s)
1196: {
1199:   *s = ls->step;
1200:   return(0);
1201: }

1203: /*MC
1204:    TaoLineSearchRegister - Adds a line-search algorithm to the registry

1206:    Not collective

1208:    Input Parameters:
1209: +  sname - name of a new user-defined solver
1210: -  func - routine to Create method context

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

1215:    Sample usage:
1216: .vb
1217:    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
1218: .ve

1220:    Then, your solver can be chosen with the procedural interface via
1221: $     TaoLineSearchSetType(ls,"my_linesearch")
1222:    or at runtime via the option
1223: $     -tao_ls_type my_linesearch

1225:    Level: developer

1227: .seealso: TaoLineSearchRegisterDestroy()
1228: M*/
1229: PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1230: {
1233:   PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);
1234:   return(0);
1235: }

1237: /*@C
1238:    TaoLineSearchRegisterDestroy - Frees the list of line-search algorithms that were
1239:    registered by TaoLineSearchRegister().

1241:    Not Collective

1243:    Level: developer

1245: .seealso: TaoLineSearchRegister()
1246: @*/
1247: PetscErrorCode TaoLineSearchRegisterDestroy(void)
1248: {
1251:   PetscFunctionListDestroy(&TaoLineSearchList);
1252:   TaoLineSearchInitialized = PETSC_FALSE;
1253:   return(0);
1254: }

1256: /*@C
1257:    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1258:    for all TaoLineSearch options in the database.


1261:    Collective on TaoLineSearch

1263:    Input Parameters:
1264: +  ls - the TaoLineSearch solver context
1265: -  prefix - the prefix string to prepend to all line search requests

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


1272:    Level: advanced

1274: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1275: @*/
1276: PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1277: {
1278:   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
1279: }

1281: /*@C
1282:   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1283:   TaoLineSearch options in the database

1285:   Not Collective

1287:   Input Parameters:
1288: . ls - the TaoLineSearch context

1290:   Output Parameters:
1291: . prefix - pointer to the prefix string used is returned

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

1296:   Level: advanced

1298: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix()
1299: @*/
1300: PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1301: {
1302:   return PetscObjectGetOptionsPrefix((PetscObject)ls,p);
1303: }

1305: /*@C
1306:    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1307:    TaoLineSearch options in the database.


1310:    Logically Collective on TaoLineSearch

1312:    Input Parameters:
1313: +  ls - the TaoLineSearch context
1314: -  prefix - the prefix string to prepend to all TAO option requests

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

1320:    For example, to distinguish between the runtime options for two
1321:    different line searches, one could call
1322: .vb
1323:       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1324:       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1325: .ve

1327:    This would enable use of different options for each system, such as
1328: .vb
1329:       -sys1_tao_ls_type mt
1330:       -sys2_tao_ls_type armijo
1331: .ve

1333:    Level: advanced

1335: .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1336: @*/

1338: PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1339: {
1340:   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
1341: }