Actual source code: taolinesearch.c

petsc-3.11.4 2019-09-28
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:   TaoLineSearchView - Prints information about the TaoLineSearch

 14:   Collective on TaoLineSearch

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

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

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

 31:   Level: beginner

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

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

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

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

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

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

 82:   Collective on MPI_Comm

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

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

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


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

 99:    Level: beginner

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

104: PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
105: {
107:   TaoLineSearch  ls;

111:   *newls = NULL;

113:   TaoLineSearchInitializePackage();

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

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

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

154:   Collective on ls

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

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

165:   Level: developer

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

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

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

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

218:   Collective on TaoLineSearch

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

223:   Level: developer

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

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

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

243:   Collective on TaoLineSearch

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

248:   Level: beginner

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

257:   if (!*ls) return(0);
259:   if (--((PetscObject)*ls)->refct > 0) {*ls=0; return(0);}
260:   VecDestroy(&(*ls)->stepdirection);
261:   VecDestroy(&(*ls)->start_x);
262:   if ((*ls)->ops->destroy) {
263:     (*(*ls)->ops->destroy)(*ls);
264:   }
265:   if ((*ls)->usemonitor) {
266:     PetscViewerDestroy(&(*ls)->viewer);
267:   }
268:   PetscHeaderDestroy(ls);
269:   return(0);
270: }

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

275:   Collective on TaoLineSearch

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

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

291:   reason will be set to one of:

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

304:   Note:
305:   The algorithm developer must set up the TaoLineSearch with calls to
306:   TaoLineSearchSetObjectiveRoutine() and TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), or TaoLineSearchUseTaoRoutines()

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

313:   Level: beginner

315:   .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchSetInitialStepLength(), TaoAddLineSearchCounts()
316:  @*/

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

324:   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
334:   VecGetOwnershipRange(x, &low1, &high1);
335:   VecGetOwnershipRange(g, &low2, &high2);
336:   VecGetOwnershipRange(s, &low3, &high3);
337:   if ( low1!= low2 || low1!= low3 || high1!= high2 || high1!= high3) SETERRQ(PETSC_COMM_SELF,1,"InCompatible vector local lengths");

339:   PetscObjectReference((PetscObject)s);
340:   VecDestroy(&ls->stepdirection);
341:   ls->stepdirection = s;

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

378:   PetscObjectReference((PetscObject)x);
379:   VecDestroy(&ls->start_x);
380:   ls->start_x = x;

382:   PetscLogEventBegin(TAOLINESEARCH_Apply,ls,0,0,0);
383:   (*ls->ops->apply)(ls,x,f,g,s);
384:   PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0,0,0);
385:   *reason=ls->reason;
386:   ls->new_f = *f;

388:   if (steplength) {
389:     *steplength=ls->step;
390:   }

392:   TaoLineSearchViewFromOptions(ls,NULL,"-tao_ls_view");
393:   return(0);
394: }

396: /*@C
397:    TaoLineSearchSetType - Sets the algorithm used in a line search

399:    Collective on TaoLineSearch

401:    Input Parameters:
402: +  ls - the TaoLineSearch context
403: -  type - a known method

405:   Available methods include:
406: + more-thuente
407: . gpcg
408: - unit - Do not perform any line search


411:   Options Database Keys:
412: .   -tao_ls_type - select which method TAO should use

414:   Level: beginner


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

419: @*/

421: PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type)
422: {
424:   PetscErrorCode (*r)(TaoLineSearch);
425:   PetscBool      flg;

430:   PetscObjectTypeCompare((PetscObject)ls, type, &flg);
431:   if (flg) return(0);

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

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

463: /*@C
464:   TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the
465:   iteration number, step length, and function value before calling the implementation 
466:   specific monitor.

468:    Input Parameters:
469: +  ls - the TaoLineSearch context
470: .  its - the current iterate number (>=0)
471: .  f - the current objective function value
472: -  step - the step length

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

477:    Level: developer

479: @*/
480: PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step)
481: {
483:   PetscInt       tabs;

487:   if (ls->usemonitor) {
488:     PetscViewerASCIIGetTab(ls->viewer, &tabs);
489:     PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel);
490:     PetscViewerASCIIPrintf(ls->viewer, "%3D LS", its);
491:     PetscViewerASCIIPrintf(ls->viewer, "  Function value: %g,", (double)f);
492:     PetscViewerASCIIPrintf(ls->viewer, "  Step length: %g\n", (double)step);
493:     if (ls->ops->monitor && its > 0) {
494:       PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3);
495:       (*ls->ops->monitor)(ls);
496:     }
497:     PetscViewerASCIISetTab(ls->viewer, tabs);
498:   }
499:   return(0);
500: }

502: /*@
503:   TaoLineSearchSetFromOptions - Sets various TaoLineSearch parameters from user
504:   options.

506:   Collective on TaoLineSearch

508:   Input Paremeter:
509: . ls - the TaoLineSearch context

511:   Options Database Keys:
512: + -tao_ls_type <type> - The algorithm that TAO uses (more-thuente, gpcg, unit)
513: . -tao_ls_ftol <tol> - tolerance for sufficient decrease
514: . -tao_ls_gtol <tol> - tolerance for curvature condition
515: . -tao_ls_rtol <tol> - relative tolerance for acceptable step
516: . -tao_ls_stepmin <step> - minimum steplength allowed
517: . -tao_ls_stepmax <step> - maximum steplength allowed
518: . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
519: - -tao_ls_view - display line-search results to standard output

521:   Level: beginner
522: @*/
523: PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
524: {
526:   const char     *default_type=TAOLINESEARCHMT;
527:   char           type[256],monfilename[PETSC_MAX_PATH_LEN];
528:   PetscViewer    monviewer;
529:   PetscBool      flg;

533:   PetscObjectOptionsBegin((PetscObject)ls);
534:   if (((PetscObject)ls)->type_name) {
535:     default_type = ((PetscObject)ls)->type_name;
536:   }
537:   /* Check for type from options */
538:   PetscOptionsFList("-tao_ls_type","Tao Line Search type","TaoLineSearchSetType",TaoLineSearchList,default_type,type,256,&flg);
539:   if (flg) {
540:     TaoLineSearchSetType(ls,type);
541:   } else if (!((PetscObject)ls)->type_name) {
542:     TaoLineSearchSetType(ls,default_type);
543:   }

545:   PetscOptionsInt("-tao_ls_max_funcs","max function evals in line search","",ls->max_funcs,&ls->max_funcs,NULL);
546:   PetscOptionsReal("-tao_ls_ftol","tol for sufficient decrease","",ls->ftol,&ls->ftol,NULL);
547:   PetscOptionsReal("-tao_ls_gtol","tol for curvature condition","",ls->gtol,&ls->gtol,NULL);
548:   PetscOptionsReal("-tao_ls_rtol","relative tol for acceptable step","",ls->rtol,&ls->rtol,NULL);
549:   PetscOptionsReal("-tao_ls_stepmin","lower bound for step","",ls->stepmin,&ls->stepmin,NULL);
550:   PetscOptionsReal("-tao_ls_stepmax","upper bound for step","",ls->stepmax,&ls->stepmax,NULL);
551:   PetscOptionsString("-tao_ls_monitor","enable the basic monitor","TaoLineSearchSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
552:   if (flg) {
553:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls),monfilename,&monviewer);
554:     ls->viewer = monviewer;
555:     ls->usemonitor = PETSC_TRUE;
556:   }
557:   if (ls->ops->setfromoptions) {
558:     (*ls->ops->setfromoptions)(PetscOptionsObject,ls);
559:   }
560:   PetscOptionsEnd();
561:   return(0);
562: }

564: /*@C
565:   TaoLineSearchGetType - Gets the current line search algorithm

567:   Not Collective

569:   Input Parameter:
570: . ls - the TaoLineSearch context

572:   Output Paramter:
573: . type - the line search algorithm in effect

575:   Level: developer

577: @*/
578: PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type)
579: {
583:   *type = ((PetscObject)ls)->type_name;
584:   return(0);
585: }

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

591:   Not Collective

593:   Input Parameter:
594: . ls - the TaoLineSearch context

596:   Output Parameters:
597: + nfeval   - number of function evaluations
598: . ngeval   - number of gradient evaluations
599: - nfgeval  - number of function/gradient evaluations

601:   Level: intermediate

603:   Note:
604:   If the line search is using the Tao objective and gradient
605:   routines directly (see TaoLineSearchUseTaoRoutines()), then TAO
606:   is already counting the number of evaluations.

608: @*/
609: PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
610: {
613:   *nfeval = ls->nfeval;
614:   *ngeval = ls->ngeval;
615:   *nfgeval = ls->nfgeval;
616:   return(0);
617: }

619: /*@
620:   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
621:   Tao evaluation routines.

623:   Not Collective

625:   Input Parameter:
626: . ls - the TaoLineSearch context

628:   Output Parameter:
629: . flg - PETSC_TRUE if the line search is using Tao evaluation routines,
630:         otherwise PETSC_FALSE

632:   Level: developer
633: @*/
634: PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
635: {
638:   *flg = ls->usetaoroutines;
639:   return(0);
640: }

642: /*@C
643:   TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search

645:   Logically Collective on TaoLineSearch

647:   Input Parameter:
648: + ls - the TaoLineSearch context
649: . func - the objective function evaluation routine
650: - ctx - the (optional) user-defined context for private data

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

655: + x - input vector
656: . f - function value
657: - ctx (optional) user-defined context

659:   Level: beginner

661:   Note:
662:   Use this routine only if you want the line search objective
663:   evaluation routine to be different from the Tao's objective
664:   evaluation routine. If you use this routine you must also set
665:   the line search gradient and/or function/gradient routine.

667:   Note:
668:   Some algorithms (lcl, gpcg) set their own objective routine for the
669:   line search, application programmers should be wary of overriding the
670:   default objective routine.

672: .seealso: TaoLineSearchCreate(), TaoLineSearchSetGradientRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
673: @*/
674: PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, PetscReal*, void*), void *ctx)
675: {

679:   ls->ops->computeobjective=func;
680:   if (ctx) ls->userctx_func=ctx;
681:   ls->usetaoroutines=PETSC_FALSE;
682:   return(0);
683: }

685: /*@C
686:   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search

688:   Logically Collective on TaoLineSearch

690:   Input Parameter:
691: + ls - the TaoLineSearch context
692: . func - the gradient evaluation routine
693: - ctx - the (optional) user-defined context for private data

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

698: + x - input vector
699: . g - gradient vector
700: - ctx (optional) user-defined context

702:   Level: beginner

704:   Note:
705:   Use this routine only if you want the line search gradient
706:   evaluation routine to be different from the Tao's gradient
707:   evaluation routine. If you use this routine you must also set
708:   the line search function and/or function/gradient routine.

710:   Note:
711:   Some algorithms (lcl, gpcg) set their own gradient routine for the
712:   line search, application programmers should be wary of overriding the
713:   default gradient routine.

715: .seealso: TaoLineSearchCreate(), TaoLineSearchSetObjectiveRoutine(), TaoLineSearchSetObjectiveAndGradientRoutine(), TaoLineSearchUseTaoRoutines()
716: @*/
717: PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode(*func)(TaoLineSearch ls, Vec x, Vec g, void*), void *ctx)
718: {
721:   ls->ops->computegradient=func;
722:   if (ctx) ls->userctx_grad=ctx;
723:   ls->usetaoroutines=PETSC_FALSE;
724:   return(0);
725: }

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

730:   Logically Collective on TaoLineSearch

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

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

740: + x - input vector
741: . f - function value
742: . g - gradient vector
743: - ctx (optional) user-defined context

745:   Level: beginner

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

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

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

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

776:   Logically Collective on TaoLineSearch

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

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

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

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

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

800:   Level: advanced

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

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

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

824:   Logically Collective on TaoLineSearch

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

830:   Level: developer

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

844: /*@
845:   TaoLineSearchComputeObjective - Computes the objective function value at a given point

847:   Collective on TaoLineSearch

849:   Input Parameters:
850: + ls - the TaoLineSearch context
851: - x - input vector

853:   Output Parameter:
854: . f - Objective value at X

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

860:   Level: developer

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

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

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

900:   Collective on Tao

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

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

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

914:   Level: developer

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

929:   if (ls->usetaoroutines) {
930:       TaoComputeObjectiveAndGradient(ls->tao,x,f,g);
931:   } else {
932:     if (!ls->ops->computeobjective && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have objective function set");
933:     if (!ls->ops->computegradient  && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient function set");
934:     PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);
935:     PetscStackPush("TaoLineSearch user objective/gradient routine");
936:     if (ls->ops->computeobjectiveandgradient) {
937:       (*ls->ops->computeobjectiveandgradient)(ls,x,f,g,ls->userctx_funcgrad);
938:     } else {
939:       (*ls->ops->computeobjective)(ls,x,f,ls->userctx_func);
940:       (*ls->ops->computegradient)(ls,x,g,ls->userctx_grad);
941:     }
942:     PetscStackPop;
943:     PetscLogEventEnd(TAOLINESEARCH_Eval,ls,0,0,0);
944:     PetscInfo1(ls,"TaoLineSearch Function evaluation: %14.12e\n",(double)(*f));
945:   }
946:   ls->nfgeval++;
947:   return(0);
948: }

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

953:   Collective on TaoLineSearch

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

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

962:   Notes:
963:     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:     if (!ls->ops->computegradient && !ls->ops->computeobjectiveandgradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Line Search does not have gradient functions set");
985:     PetscLogEventBegin(TAOLINESEARCH_Eval,ls,0,0,0);
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_Eval,ls,0,0,0);
994:   }
995:   ls->ngeval++;
996:   return(0);
997: }

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

1002:   Collective on Tao

1004:   Input Parameters:
1005: + ls - the TaoLineSearch context
1006: - x - input vector

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

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

1016:   Level: developer

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

1040: /*@
1041:   TaoLineSearchGetSolution - Returns the solution to the line search

1043:   Collective on TaoLineSearch

1045:   Input Parameter:
1046: . ls - the TaoLineSearch context

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

1055:   reason will be set to one of:

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

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

1068: + TAOLINESEARCH_SUCCESS - successful line search

1070:   Level: developer

1072: @*/
1073: PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
1074: {

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

1097: /*@
1098:   TaoLineSearchGetStartingVector - Gets a the initial point of the line
1099:   search.

1101:   Not Collective

1103:   Input Parameter:
1104: . ls - the TaoLineSearch context

1106:   Output Parameter:
1107: . x - The initial point of the line search

1109:   Level: intermediate
1110: @*/
1111: PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1112: {
1115:   if (x) {
1116:     *x = ls->start_x;
1117:   }
1118:   return(0);
1119: }

1121: /*@
1122:   TaoLineSearchGetStepDirection - Gets the step direction of the line
1123:   search.

1125:   Not Collective

1127:   Input Parameter:
1128: . ls - the TaoLineSearch context

1130:   Output Parameter:
1131: . s - the step direction of the line search

1133:   Level: advanced
1134: @*/
1135: PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1136: {
1139:   if (s) {
1140:     *s = ls->stepdirection;
1141:   }
1142:   return(0);

1144: }

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

1149:   Not Collective

1151:   Input Parameter:
1152: . ls - the TaoLineSearch context

1154:   Output Parameter:
1155: . f - the objective value at the full step length

1157:   Level: developer
1158: @*/

1160: PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1161: {
1164:   *f_fullstep = ls->f_fullstep;
1165:   return(0);
1166: }

1168: /*@
1169:   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds.

1171:   Logically Collective on Tao

1173:   Input Parameters:
1174: + ls - the TaoLineSearch context
1175: . xl  - vector of lower bounds
1176: - xu  - vector of upper bounds

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

1181:   Level: beginner

1183: .seealso: TaoSetVariableBounds(), TaoLineSearchCreate()
1184: @*/
1185: PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl, Vec xu)
1186: {
1191:   ls->lower = xl;
1192:   ls->upper = xu;
1193:   ls->bounded = 1;
1194:   return(0);
1195: }

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

1201:   Logically Collective on TaoLineSearch

1203:   Input Parameters:
1204: + ls - the TaoLineSearch context
1205: - s - the initial step size

1207:   Level: intermediate

1209: .seealso: TaoLineSearchGetStepLength(), TaoLineSearchApply()
1210: @*/
1211: PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)
1212: {
1215:   ls->initstep = s;
1216:   return(0);
1217: }

1219: /*@
1220:   TaoLineSearchGetStepLength - Get the current step length

1222:   Not Collective

1224:   Input Parameters:
1225: . ls - the TaoLineSearch context

1227:   Output Parameters
1228: . s - the current step length

1230:   Level: beginner

1232: .seealso: TaoLineSearchSetInitialStepLength(), TaoLineSearchApply()
1233: @*/
1234: PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal *s)
1235: {
1238:   *s = ls->step;
1239:   return(0);
1240: }

1242: /*@C
1243:    TaoLineSearchRegister - Adds a line-search algorithm to the registry

1245:    Not collective

1247:    Input Parameters:
1248: +  sname - name of a new user-defined solver
1249: -  func - routine to Create method context

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

1254:    Sample usage:
1255: .vb
1256:    TaoLineSearchRegister("my_linesearch",MyLinesearchCreate);
1257: .ve

1259:    Then, your solver can be chosen with the procedural interface via
1260: $     TaoLineSearchSetType(ls,"my_linesearch")
1261:    or at runtime via the option
1262: $     -tao_ls_type my_linesearch

1264:    Level: developer

1266: @*/
1267: PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1268: {
1271:   TaoLineSearchInitializePackage();
1272:   PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func);
1273:   return(0);
1274: }

1276: /*@C
1277:    TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1278:    for all TaoLineSearch options in the database.


1281:    Collective on TaoLineSearch

1283:    Input Parameters:
1284: +  ls - the TaoLineSearch solver context
1285: -  prefix - the prefix string to prepend to all line search requests

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


1292:    Level: advanced

1294: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1295: @*/
1296: PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1297: {
1298:   return PetscObjectAppendOptionsPrefix((PetscObject)ls,p);
1299: }

1301: /*@C
1302:   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1303:   TaoLineSearch options in the database

1305:   Not Collective

1307:   Input Parameters:
1308: . ls - the TaoLineSearch context

1310:   Output Parameters:
1311: . prefix - pointer to the prefix string used is returned

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

1317:   Level: advanced

1319: .seealso: TaoLineSearchSetOptionsPrefix(), TaoLineSearchAppendOptionsPrefix()
1320: @*/
1321: PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1322: {
1323:   return PetscObjectGetOptionsPrefix((PetscObject)ls,p);
1324: }

1326: /*@C
1327:    TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1328:    TaoLineSearch options in the database.


1331:    Logically Collective on TaoLineSearch

1333:    Input Parameters:
1334: +  ls - the TaoLineSearch context
1335: -  prefix - the prefix string to prepend to all TAO option requests

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

1341:    For example, to distinguish between the runtime options for two
1342:    different line searches, one could call
1343: .vb
1344:       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1345:       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1346: .ve

1348:    This would enable use of different options for each system, such as
1349: .vb
1350:       -sys1_tao_ls_type mt
1351:       -sys2_tao_ls_type armijo
1352: .ve

1354:    Level: advanced

1356: .seealso: TaoLineSearchAppendOptionsPrefix(), TaoLineSearchGetOptionsPrefix()
1357: @*/

1359: PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1360: {
1361:   return PetscObjectSetOptionsPrefix((PetscObject)ls,p);
1362: }