Actual source code: taolinesearch.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1:  #include <petsctaolinesearch.h>
  2:  #include <petsc/private/taolinesearchimpl.h>

  4: PetscFunctionList TaoLineSearchList = NULL;

  6: PetscClassId TAOLINESEARCH_CLASSID=0;

  8: PetscLogEvent TAOLINESEARCH_Apply;
  9: PetscLogEvent TAOLINESEARCH_Eval;

 11: /*@C
 12:    TaoLineSearchViewFromOptions - View from Options

 14:    Collective on TaoLineSearch

 16:    Input Parameters:
 17: +  A - the Tao context
 18: .  obj - Optional object
 19: -  name - command line option

 21:    Level: intermediate
 22: .seealso:  TaoLineSearch, TaoLineSearchView, PetscObjectViewFromOptions(), TaoLineSearchCreate()
 23: @*/
 24: PetscErrorCode  TaoLineSearchViewFromOptions(TaoLineSearch A,PetscObject obj,const char name[])
 25: {

 30:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
 31:   return(0);
 32: }

 34: /*@C
 35:   TaoLineSearchView - Prints information about the TaoLineSearch

 37:   Collective on TaoLineSearch

 39:   InputParameters:
 40: + ls - the Tao context
 41: - viewer - visualization context

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

 46:   Notes:
 47:   The available visualization contexts include
 48: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
 49: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
 50:          output where only the first processor opens
 51:          the file.  All other processors send their
 52:          data to the first processor to print.

 54:   Level: beginner

 56: .seealso: PetscViewerASCIIOpen()
 57: @*/

 59: PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
 60: {
 61:   PetscErrorCode          ierr;
 62:   PetscBool               isascii, isstring;
 63:   TaoLineSearchType       type;

 67:   if (!viewer) {
 68:     PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer);
 69:   }

 73:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 74:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
 75:   if (isascii) {
 76:     PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer);
 77:     if (ls->ops->view) {
 78:       PetscViewerASCIIPushTab(viewer);
 79:       (*ls->ops->view)(ls,viewer);
 80:       PetscViewerASCIIPopTab(viewer);
 81:     }
 82:     PetscViewerASCIIPushTab(viewer);
 83:     PetscViewerASCIIPrintf(viewer,"maximum function evaluations=%D\n",ls->max_funcs);
 84:     PetscViewerASCIIPrintf(viewer,"tolerances: ftol=%g, rtol=%g, gtol=%g\n",(double)ls->ftol,(double)ls->rtol,(double)ls->gtol);
 85:     PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D\n",ls->nfeval);
 86:     PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D\n",ls->ngeval);
 87:     PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D\n",ls->nfgeval);

 89:     if (ls->bounded) {
 90:       PetscViewerASCIIPrintf(viewer,"using variable bounds\n");
 91:     }
 92:     PetscViewerASCIIPrintf(viewer,"Termination reason: %d\n",(int)ls->reason);
 93:     PetscViewerASCIIPopTab(viewer);
 94:   } else if (isstring) {
 95:     TaoLineSearchGetType(ls,&type);
 96:     PetscViewerStringSPrintf(viewer," %-3.3s",type);
 97:   }
 98:   return(0);
 99: }

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

105:   Collective

107:   Input Parameter:
108: . comm - MPI communicator

110:   Output Parameter:
111: . newls - the new TaoLineSearch context

113:   Available methods include:
114: + more-thuente
115: . gpcg
116: - unit - Do not perform any line search


119:    Options Database Keys:
120: .   -tao_ls_type - select which method TAO should use

122:    Level: beginner

124: .seealso: TaoLineSearchSetType(), TaoLineSearchApply(), TaoLineSearchDestroy()
125: @*/

127: PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
128: {
130:   TaoLineSearch  ls;

134:   *newls = NULL;

136:   TaoLineSearchInitializePackage();

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

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

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

177:   Collective on ls

179:   Input Parameters:
180: . ls - the TaoLineSearch context

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

188:   Level: developer

190: .seealso: TaoLineSearchCreate(), TaoLineSearchApply()
191: @*/

193: PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
194: {
196:   const char     *default_type=TAOLINESEARCHMT;
197:   PetscBool      flg;

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

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

241:   Collective on TaoLineSearch

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

246:   Level: developer

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

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

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

266:   Collective on TaoLineSearch

268:   Input Parameter:
269: . ls - the TaoLineSearch context

271:   Level: beginner

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

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

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

298:   Collective on TaoLineSearch

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

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

314:   reason will be set to one of:

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

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

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

336:   Level: beginner

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

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

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

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

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

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

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

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

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

419: /*@C
420:    TaoLineSearchSetType - Sets the algorithm used in a line search

422:    Collective on TaoLineSearch

424:    Input Parameters:
425: +  ls - the TaoLineSearch context
426: -  type - the TaoLineSearchType selection

428:   Available methods include:
429: +  more-thuente - line search with a cubic model enforcing the strong Wolfe/curvature condition
430: .  armijo - simple backtracking line search enforcing only the sufficient decrease condition
431: -  unit - do not perform a line search and always accept unit step length

433:   Options Database Keys:
434: .  -tao_ls_type <more-thuente, armijo, unit> - select which method TAO should use at runtime

436:   Level: beginner

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

440: @*/

442: PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, 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(PetscObjectComm((PetscObject)ls),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: }

484: /*@C
485:   TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the
486:   iteration number, step length, and function value before calling the implementation 
487:   specific monitor.

489:    Input Parameters:
490: +  ls - the TaoLineSearch context
491: .  its - the current iterate number (>=0)
492: .  f - the current objective function value
493: -  step - the step length

495:    Options Database Key:
496: .  -tao_ls_monitor - Use the default monitor, which prints statistics to standard output

498:    Level: developer

500: @*/
501: PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step)
502: {
504:   PetscInt       tabs;

508:   if (ls->usemonitor) {
509:     PetscViewerASCIIGetTab(ls->viewer, &tabs);
510:     PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel);
511:     PetscViewerASCIIPrintf(ls->viewer, "%3D LS", its);
512:     PetscViewerASCIIPrintf(ls->viewer, "  Function value: %g,", (double)f);
513:     PetscViewerASCIIPrintf(ls->viewer, "  Step length: %g\n", (double)step);
514:     if (ls->ops->monitor && its > 0) {
515:       PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3);
516:       (*ls->ops->monitor)(ls);
517:     }
518:     PetscViewerASCIISetTab(ls->viewer, tabs);
519:   }
520:   return(0);
521: }

523: /*@
524:   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
525:   options.

527:   Collective on TaoLineSearch

529:   Input Paremeter:
530: . ls - the TaoLineSearch context

532:   Options Database Keys:
533: + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
534: . -tao_ls_ftol <tol> - tolerance for sufficient decrease
535: . -tao_ls_gtol <tol> - tolerance for curvature condition
536: . -tao_ls_rtol <tol> - relative tolerance for acceptable step
537: . -tao_ls_stepmin <step> - minimum steplength allowed
538: . -tao_ls_stepmax <step> - maximum steplength allowed
539: . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
540: - -tao_ls_view - display line-search results to standard output

542:   Level: beginner
543: @*/
544: PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
545: {
547:   const char     *default_type=TAOLINESEARCHMT;
548:   char           type[256],monfilename[PETSC_MAX_PATH_LEN];
549:   PetscViewer    monviewer;
550:   PetscBool      flg;

554:   PetscObjectOptionsBegin((PetscObject)ls);
555:   if (((PetscObject)ls)->type_name) {
556:     default_type = ((PetscObject)ls)->type_name;
557:   }
558:   /* Check for type from options */
559:   PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);
560:   if (flg) {
561:     TaoLineSearchSetType(ls,type);
562:   } else if (!((PetscObject)ls)->type_name) {
563:     TaoLineSearchSetType(ls,default_type);
564:   }

566:   PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL);
567:   PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL);
568:   PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL);
569:   PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL);
570:   PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL);
571:   PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL);
572:   PetscOptionsString("-tao_ls_monitor","enable the basic monitor","TaoLineSearchSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
573:   if (flg) {
574:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls),monfilename,&monviewer);
575:     ls->viewer = monviewer;
576:     ls->usemonitor = PETSC_TRUE;
577:   }
578:   if (ls->ops->setfromoptions) {
579:     (*ls->ops->setfromoptions)(PetscOptionsObject,ls);
580:   }
581:   PetscOptionsEnd();
582:   return(0);
583: }

585: /*@C
586:   TaoLineSearchGetType - Gets the current line search algorithm

588:   Not Collective

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

593:   Output Parameter:
594: . type - the line search algorithm in effect

596:   Level: developer

598: @*/
599: PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type)
600: {
604:   *type = ((PetscObject)ls)->type_name;
605:   return(0);
606: }

608: /*@
609:   TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation
610:   routines used by the line search in last Section 1.5 Writing Application Codes with PETSc (not cumulative).

612:   Not Collective

614:   Input Parameter:
615: . ls - the TaoLineSearch context

617:   Output Parameters:
618: + nfeval   - number of function evaluations
619: . ngeval   - number of gradient evaluations
620: - nfgeval  - number of function/gradient evaluations

622:   Level: intermediate

624:   Note:
625:   If the line search is using the Tao objective and gradient
626:   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
627:   is already counting the number of evaluations.

629: @*/
630: PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
631: {
634:   *nfeval = ls->nfeval;
635:   *ngeval = ls->ngeval;
636:   *nfgeval = ls->nfgeval;
637:   return(0);
638: }

640: /*@
641:   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
642:   Tao evaluation routines.

644:   Not Collective

646:   Input Parameter:
647: . ls - the TaoLineSearch context

649:   Output Parameter:
650: . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
651:         otherwise PETSC_FALSE

653:   Level: developer
654: @*/
655: PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
656: {
659:   *flg = ls->usetaoroutines;
660:   return(0);
661: }

663: /*@C
664:   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search

666:   Logically Collective on TaoLineSearch

668:   Input Parameter:
669: + ls - the TaoLineSearch context
670: . func - the objective function evaluation routine
671: - ctx - the (optional) user-defined context for private data

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

676: + x - input vector
677: . f - function value
678: - ctx (optional) user-defined context

680:   Level: beginner

682:   Note:
683:   Use this routine only if you want the line search objective
684:   evaluation routine to be different from the Tao's objective
685:   evaluation routine. If you use this routine you must also set
686:   the line search gradient and/or function/gradient routine.

688:   Note:
689:   Some algorithms (lcl, gpcg) set their own objective routine for the
690:   line search, Section 1.5 Writing Application Codes with PETSc programmers should be wary of overriding the
691:   default objective routine.

693: .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
694: @*/
695: PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
696: {

700:   ls->ops->computeobjective=func;
701:   if (ctx) ls->userctx_func=ctx;
702:   ls->usetaoroutines=PETSC_FALSE;
703:   return(0);
704: }

706: /*@C
707:   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search

709:   Logically Collective on TaoLineSearch

711:   Input Parameter:
712: + ls - the TaoLineSearch context
713: . func - the gradient evaluation routine
714: - ctx - the (optional) user-defined context for private data

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

719: + x - input vector
720: . g - gradient vector
721: - ctx (optional) user-defined context

723:   Level: beginner

725:   Note:
726:   Use this routine only if you want the line search gradient
727:   evaluation routine to be different from the Tao's gradient
728:   evaluation routine. If you use this routine you must also set
729:   the line search function and/or function/gradient routine.

731:   Note:
732:   Some algorithms (lcl, gpcg) set their own gradient routine for the
733:   line search, Section 1.5 Writing Application Codes with PETSc programmers should be wary of overriding the
734:   default gradient routine.

736: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
737: @*/
738: PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
739: {
742:   ls->ops->computegradient=func;
743:   if (ctx) ls->userctx_grad=ctx;
744:   ls->usetaoroutines=PETSC_FALSE;
745:   return(0);
746: }

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

751:   Logically Collective on TaoLineSearch

753:   Input Parameter:
754: + ls - the TaoLineSearch context
755: . func - the objective and gradient evaluation routine
756: - ctx - the (optional) user-defined context for private data

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

761: + x - input vector
762: . f - function value
763: . g - gradient vector
764: - ctx (optional) user-defined context

766:   Level: beginner

768:   Note:
769:   Use this routine only if you want the line search objective and gradient
770:   evaluation routines to be different from the Tao's objective
771:   and gradient evaluation routines.

773:   Note:
774:   Some algorithms (lcl, gpcg) set their own objective routine for the
775:   line search, Section 1.5 Writing Application Codes with PETSc programmers should be wary of overriding the
776:   default objective routine.

778: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetGradientRoutine(), TaoLineSearchUseTaoRoutines()
779: @*/
780: PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal *, Vec g, void*), void *ctx)
781: {
784:   ls->ops->computeobjectiveandgradient=func;
785:   if (ctx) ls->userctx_funcgrad=ctx;
786:   ls->usetaoroutines = PETSC_FALSE;
787:   return(0);
788: }

790: /*@C
791:   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
792:   (gradient'*stepdirection) evaluation routine for the line search.
793:   Sometimes it is more efficient to compute the inner product of the gradient
794:   and the step direction than it is to compute the gradient, and this is all
795:   the line search typically needs of the gradient.

797:   Logically Collective on TaoLineSearch

799:   Input Parameter:
800: + ls - the TaoLineSearch context
801: . func - the objective and gradient evaluation routine
802: - ctx - the (optional) user-defined context for private data

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

807: + x - input vector
808: . s - step direction
809: . f - function value
810: . gts - inner product of gradient and step direction vectors
811: - ctx (optional) user-defined context

813:   Note: The gradient will still need to be computed at the end of the line
814:   search, so you will still need to set a line search gradient evaluation
815:   routine

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

821:   Level: advanced

823:   Note:
824:   Some algorithms (lcl, gpcg) set their own objective routine for the
825:   line search, Section 1.5 Writing Application Codes with PETSc programmers should be wary of overriding the
826:   default objective routine.

828: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjective(), TaoLineSearchSetGradient(), TaoLineSearchUseTaoRoutines()
829: @*/
830: PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *, PetscReal *, void*), void *ctx)
831: {
834:   ls->ops->computeobjectiveandgts=func;
835:   if (ctx) ls->userctx_funcgts=ctx;
836:   ls->usegts = PETSC_TRUE;
837:   ls->usetaoroutines=PETSC_FALSE;
838:   return(0);
839: }

841: /*@C
842:   TaoLineSearchUseTaoRoutines - Informs the TaoLineSearch to use the
843:   objective and gradient evaluation routines from the given Tao object.

845:   Logically Collective on TaoLineSearch

847:   Input Parameter:
848: + ls - the TaoLineSearch context
849: - ts - the Tao context with defined objective/gradient evaluation routines

851:   Level: developer

853: .seealso: TaoLineSearchCreate()
854: @*/
855: PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
856: {
860:   ls->tao = ts;
861:   ls->usetaoroutines=PETSC_TRUE;
862:   return(0);
863: }

865: /*@
866:   TaoLineSearchComputeObjective - Computes the objective function value at a given point

868:   Collective on TaoLineSearch

870:   Input Parameters:
871: + ls - the TaoLineSearch context
872: - x - input vector

874:   Output Parameter:
875: . f - Objective value at X

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

881:   Level: developer

883: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
884: @*/
885: PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
886: {
888:   Vec            gdummy;
889:   PetscReal      gts;

896:   if (ls->usetaoroutines) {
897:     TaoComputeObjective(ls->tao,x,f);
898:   } else {
899:     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient && !ls->ops->computeobjectiveandgts) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
900:     PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);
901:     PetscStackPush("TaoLineSearch user objective routine");
902:     if (ls->ops->computeobjective) {
903:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
904:     } else if (ls->ops->computeobjectiveandgradient) {
905:       VecDuplicate(x,&gdummy);
906:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,gdummy,ls->userctx_funcgrad);
907:       VecDestroy(&gdummy);
908:     } else {
909:       (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,&gts,ls->userctx_funcgts);
910:     }
911:     PetscStackPop;
912:     PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);
913:   }
914:   ls->nfeval++;
915:   return(0);
916: }

918: /*@
919:   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point

921:   Collective on Tao

923:   Input Parameters:
924: + ls - the TaoLineSearch context
925: - x - input vector

927:   Output Parameter:
928: + f - Objective value at X
929: - g - Gradient vector at X

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

935:   Level: developer

937: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
938: @*/
939: PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
940: {

950:   if (ls->usetaoroutines) {
951:       TaoComputeObjectiveAndGradient(ls->tao,x,f,g);
952:   } else {
953:     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
954:     if (!ls->ops->computegradient  && !ls->ops->computeobjectiveandgradient) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set");
955:     PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);
956:     PetscStackPush("TaoLineSearch user objective/gradient routine");
957:     if (ls->ops->computeobjectiveandgradient) {
958:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);
959:     } else {
960:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
961:       (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);
962:     }
963:     PetscStackPop;
964:     PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);
965:     PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
966:   }
967:   ls->nfgeval++;
968:   return(0);
969: }

971: /*@
972:   TaoLineSearchComputeGradient - Computes the gradient of the objective function

974:   Collective on TaoLineSearch

976:   Input Parameters:
977: + ls - the TaoLineSearch context
978: - x - input vector

980:   Output Parameter:
981: . g - gradient vector

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

987:   Level: developer

989: .seealso: TaoLineSearchComputeObjective(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetGradient()
990: @*/
991: PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
992: {
994:   PetscReal      fdummy;

1002:   if (ls->usetaoroutines) {
1003:     TaoComputeGradient(ls->tao,x,g);
1004:   } else {
1005:     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set");
1006:     PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);
1007:     PetscStackPush("TaoLineSearch user gradient routine");
1008:     if (ls->ops->computegradient) {
1009:       (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);
1010:     } else {
1011:       (*ls->ops->computeobjectiveandgradient)(ls,x,&fdummy,g,ls->userctx_funcgrad);
1012:     }
1013:     PetscStackPop;
1014:     PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);
1015:   }
1016:   ls->ngeval++;
1017:   return(0);
1018: }

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

1023:   Collective on Tao

1025:   Input Parameters:
1026: + ls - the TaoLineSearch context
1027: - x - input vector

1029:   Output Parameter:
1030: + f - Objective value at X
1031: - gts - inner product of gradient and step direction at X

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

1037:   Level: developer

1039: .seealso: TaoLineSearchComputeGradient(), TaoLineSearchComputeObjectiveAndGradient(), TaoLineSearchSetObjectiveRoutine()
1040: @*/
1041: PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
1042: {
1050:   if (!ls->ops->computeobjectiveandgts) SETERRQ(PetscObjectComm((PetscObject)ls),PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective and gts function set");
1051:   PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);
1052:   PetscStackPush("TaoLineSearch user objective/gts routine");
1053:   (*ls->ops->computeobjectiveandgts)(ls,x,ls->stepdirection,f,gts,ls->userctx_funcgts);
1054:   PetscStackPop;
1055:   PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);
1056:   PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
1057:   ls->nfeval++;
1058:   return(0);
1059: }

1061: /*@
1062:   TaoLineSearchGetSolution - Returns the solution to the line search

1064:   Collective on TaoLineSearch

1066:   Input Parameter:
1067: . ls - the TaoLineSearch context

1069:   Output Parameter:
1070: + x - the new solution
1071: . f - the objective function value at x
1072: . g - the gradient at x
1073: . steplength - the multiple of the step direction taken by the line search
1074: - reason - the reason why the line search terminated

1076:   reason will be set to one of:

1078: + TAOLINESEARCH_FAILED_INFORNAN - function evaluation gives Inf or Nan value
1079: . TAOLINESEARCH_FAILED_BADPARAMETER - negative value set as parameter
1080: . TAOLINESEARCH_FAILED_ASCENT - initial line search step * g is not descent direction
1081: . TAOLINESEARCH_HALTED_MAXFCN - maximum number of function evaluation reached
1082: . TAOLINESEARCH_HALTED_UPPERBOUND - step is at upper bound
1083: . TAOLINESEARCH_HALTED_LOWERBOUND - step is at lower bound
1084: . TAOLINESEARCH_HALTED_RTOL - range of uncertainty is smaller than given tolerance

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

1089: - TAOLINESEARCH_SUCCESS - successful line search

1091:   Level: developer

1093: @*/
1094: PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
1095: {

1104:   if (ls->new_x) {
1105:     VecCopy(ls->new_x,x);
1106:   }
1107:   *f = ls->new_f;
1108:   if (ls->new_g) {
1109:     VecCopy(ls->new_g,g);
1110:   }
1111:   if (steplength) {
1112:     *steplength=ls->step;
1113:   }
1114:   *reason = ls->reason;
1115:   return(0);
1116: }

1118: /*@
1119:   TaoLineSearchGetStartingVector - Gets a the initial point of the line
1120:   search.

1122:   Not Collective

1124:   Input Parameter:
1125: . ls - the TaoLineSearch context

1127:   Output Parameter:
1128: . x - The initial point of the line search

1130:   Level: intermediate
1131: @*/
1132: PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1133: {
1136:   if (x) {
1137:     *x = ls->start_x;
1138:   }
1139:   return(0);
1140: }

1142: /*@
1143:   TaoLineSearchGetStepDirection - Gets the step direction of the line
1144:   search.

1146:   Not Collective

1148:   Input Parameter:
1149: . ls - the TaoLineSearch context

1151:   Output Parameter:
1152: . s - the step direction of the line search

1154:   Level: advanced
1155: @*/
1156: PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1157: {
1160:   if (s) {
1161:     *s = ls->stepdirection;
1162:   }
1163:   return(0);

1165: }

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

1170:   Not Collective

1172:   Input Parameter:
1173: . ls - the TaoLineSearch context

1175:   Output Parameter:
1176: . f - the objective value at the full step length

1178:   Level: developer
1179: @*/

1181: PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1182: {
1185:   *f_fullstep = ls->f_fullstep;
1186:   return(0);
1187: }

1189: /*@
1190:   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.

1192:   Logically Collective on Tao

1194:   Input Parameters:
1195: + ls - the TaoLineSearch context
1196: . xl  - vector of lower bounds
1197: - xu  - vector of upper bounds

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

1202:   Level: beginner

1204: .seealso: TaoSetVariableBounds(), TaoLineSearchCreate()
1205: @*/
1206: PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
1207: {
1212:   ls->lower = xl;
1213:   ls->upper = xu;
1214:   ls->bounded = 1;
1215:   return(0);
1216: }

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

1222:   Logically Collective on TaoLineSearch

1224:   Input Parameters:
1225: + ls - the TaoLineSearch context
1226: - s - the initial step size

1228:   Level: intermediate

1230: .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply()
1231: @*/
1232: PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
1233: {
1236:   ls->initstep = s;
1237:   return(0);
1238: }

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

1243:   Not Collective

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

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

1251:   Level: beginner

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

1263: /*@C
1264:    TaoLineSearchRegister - Adds a line-search algorithm to the registry

1266:    Not collective

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

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

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

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

1285:    Level: developer

1287: @*/
1288: PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1289: {
1292:   TaoLineSearchInitializePackage();
1293:   PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);
1294:   return(0);
1295: }

1297: /*@C
1298:    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1299:    for all TaoLineSearch options in the database.


1302:    Collective on TaoLineSearch

1304:    Input Parameters:
1305: +  ls - the TaoLineSearch solver context
1306: -  prefix - the prefix string to prepend to all line search requests

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


1313:    Level: advanced

1315: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1316: @*/
1317: PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1318: {
1319:   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
1320: }

1322: /*@C
1323:   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1324:   TaoLineSearch options in the database

1326:   Not Collective

1328:   Input Parameters:
1329: . ls - the TaoLineSearch context

1331:   Output Parameters:
1332: . prefix - pointer to the prefix string used is returned

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

1338:   Level: advanced

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

1347: /*@C
1348:    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1349:    TaoLineSearch options in the database.


1352:    Logically Collective on TaoLineSearch

1354:    Input Parameters:
1355: +  ls - the TaoLineSearch context
1356: -  prefix - the prefix string to prepend to all TAO option requests

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

1362:    For example, to distinguish between the runtime options for two
1363:    different line searches, one could call
1364: .vb
1365:       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1366:       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1367: .ve

1369:    This would enable use of different options for each system, such as
1370: .vb
1371:       -sys1_tao_ls_type mt
1372:       -sys2_tao_ls_type armijo
1373: .ve

1375:    Level: advanced

1377: .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1378: @*/

1380: PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1381: {
1382:   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
1383: }