Actual source code: taolinesearch.c

  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 a `TaoLineSearch` object based on values in the options database

 14:   Collective

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

 21:   Level: intermediate

 23:   Note:
 24:   See `PetscObjectViewFromOptions()` for available viewer options

 26: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchView()`, `PetscObjectViewFromOptions()`, `TaoLineSearchCreate()`
 27: @*/
 28: PetscErrorCode TaoLineSearchViewFromOptions(TaoLineSearch A, PetscObject obj, const char name[])
 29: {
 30:   PetscFunctionBegin;
 32:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
 33:   PetscFunctionReturn(PETSC_SUCCESS);
 34: }

 36: /*@C
 37:   TaoLineSearchView - Prints information about the `TaoLineSearch`

 39:   Collective

 41:   Input Parameters:
 42: + ls     - the `TaoLineSearch` context
 43: - viewer - visualization context

 45:   Options Database Key:
 46: . -tao_ls_view - Calls `TaoLineSearchView()` at the end of each line search

 48:   Level: beginner

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

 58: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `PetscViewerASCIIOpen()`, `TaoLineSearchViewFromOptions()`
 59: @*/
 60: PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
 61: {
 62:   PetscBool         isascii, isstring;
 63:   TaoLineSearchType type;

 65:   PetscFunctionBegin;
 67:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer));
 69:   PetscCheckSameComm(ls, 1, viewer, 2);

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

 85:     if (ls->bounded) PetscCall(PetscViewerASCIIPrintf(viewer, "using variable bounds\n"));
 86:     PetscCall(PetscViewerASCIIPrintf(viewer, "Termination reason: %d\n", (int)ls->reason));
 87:     PetscCall(PetscViewerASCIIPopTab(viewer));
 88:   } else if (isstring) {
 89:     PetscCall(TaoLineSearchGetType(ls, &type));
 90:     PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type));
 91:   }
 92:   PetscFunctionReturn(PETSC_SUCCESS);
 93: }

 95: /*@C
 96:   TaoLineSearchCreate - Creates a `TaoLineSearch` object.  Algorithms in `Tao` that use
 97:   line-searches will automatically create one so this all is rarely needed

 99:   Collective

101:   Input Parameter:
102: . comm - MPI communicator

104:   Output Parameter:
105: . newls - the new `TaoLineSearch` context

107:   Options Database Key:
108: . -tao_ls_type - select which method `Tao` should use

110:   Level: developer

112: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`, `TaoLineSearchDestroy()`
113: @*/
114: PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
115: {
116:   TaoLineSearch ls;

118:   PetscFunctionBegin;
119:   PetscAssertPointer(newls, 2);
120:   PetscCall(TaoLineSearchInitializePackage());

122:   PetscCall(PetscHeaderCreate(ls, TAOLINESEARCH_CLASSID, "TaoLineSearch", "Linesearch", "Tao", comm, TaoLineSearchDestroy, TaoLineSearchView));
123:   ls->max_funcs = 30;
124:   ls->ftol      = 0.0001;
125:   ls->gtol      = 0.9;
126: #if defined(PETSC_USE_REAL_SINGLE)
127:   ls->rtol = 1.0e-5;
128: #else
129:   ls->rtol = 1.0e-10;
130: #endif
131:   ls->stepmin  = 1.0e-20;
132:   ls->stepmax  = 1.0e+20;
133:   ls->step     = 1.0;
134:   ls->initstep = 1.0;
135:   *newls       = ls;
136:   PetscFunctionReturn(PETSC_SUCCESS);
137: }

139: /*@
140:   TaoLineSearchSetUp - Sets up the internal data structures for the later use
141:   of a `TaoLineSearch`

143:   Collective

145:   Input Parameter:
146: . ls - the `TaoLineSearch` context

148:   Level: developer

150:   Note:
151:   The user will not need to explicitly call `TaoLineSearchSetUp()`, as it will
152:   automatically be called in `TaoLineSearchSolve()`.  However, if the user
153:   desires to call it explicitly, it should come after `TaoLineSearchCreate()`
154:   but before `TaoLineSearchApply()`.

156: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchApply()`
157: @*/
158: PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
159: {
160:   const char *default_type = TAOLINESEARCHMT;
161:   PetscBool   flg;

163:   PetscFunctionBegin;
165:   if (ls->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
166:   if (!((PetscObject)ls)->type_name) PetscCall(TaoLineSearchSetType(ls, default_type));
167:   PetscTryTypeMethod(ls, setup);
168:   if (ls->usetaoroutines) {
169:     PetscCall(TaoIsObjectiveDefined(ls->tao, &flg));
170:     ls->hasobjective = flg;
171:     PetscCall(TaoIsGradientDefined(ls->tao, &flg));
172:     ls->hasgradient = flg;
173:     PetscCall(TaoIsObjectiveAndGradientDefined(ls->tao, &flg));
174:     ls->hasobjectiveandgradient = flg;
175:   } else {
176:     if (ls->ops->computeobjective) {
177:       ls->hasobjective = PETSC_TRUE;
178:     } else {
179:       ls->hasobjective = PETSC_FALSE;
180:     }
181:     if (ls->ops->computegradient) {
182:       ls->hasgradient = PETSC_TRUE;
183:     } else {
184:       ls->hasgradient = PETSC_FALSE;
185:     }
186:     if (ls->ops->computeobjectiveandgradient) {
187:       ls->hasobjectiveandgradient = PETSC_TRUE;
188:     } else {
189:       ls->hasobjectiveandgradient = PETSC_FALSE;
190:     }
191:   }
192:   ls->setupcalled = PETSC_TRUE;
193:   PetscFunctionReturn(PETSC_SUCCESS);
194: }

196: /*@
197:   TaoLineSearchReset - Some line searches may carry state information
198:   from one `TaoLineSearchApply()` to the next.  This function resets this
199:   state information.

201:   Collective

203:   Input Parameter:
204: . ls - the `TaoLineSearch` context

206:   Level: developer

208: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchApply()`
209: @*/
210: PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
211: {
212:   PetscFunctionBegin;
214:   PetscTryTypeMethod(ls, reset);
215:   PetscFunctionReturn(PETSC_SUCCESS);
216: }

218: /*@
219:   TaoLineSearchDestroy - Destroys the `TaoLineSearch` context that was created with
220:   `TaoLineSearchCreate()`

222:   Collective

224:   Input Parameter:
225: . ls - the `TaoLineSearch` context

227:   Level: developer

229: .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSolve()`
230: @*/
231: PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
232: {
233:   PetscFunctionBegin;
234:   if (!*ls) PetscFunctionReturn(PETSC_SUCCESS);
236:   if (--((PetscObject)*ls)->refct > 0) {
237:     *ls = NULL;
238:     PetscFunctionReturn(PETSC_SUCCESS);
239:   }
240:   PetscCall(VecDestroy(&(*ls)->stepdirection));
241:   PetscCall(VecDestroy(&(*ls)->start_x));
242:   PetscCall(VecDestroy(&(*ls)->upper));
243:   PetscCall(VecDestroy(&(*ls)->lower));
244:   PetscTryTypeMethod(*ls, destroy);
245:   if ((*ls)->usemonitor) PetscCall(PetscViewerDestroy(&(*ls)->viewer));
246:   PetscCall(PetscHeaderDestroy(ls));
247:   PetscFunctionReturn(PETSC_SUCCESS);
248: }

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

254:   Collective

256:   Input Parameters:
257: + ls - the `TaoLineSearch` context
258: - s  - search direction

260:   Output Parameters:
261: + x          - On input the current solution, on output `x` contains the new solution determined by the line search
262: . f          - On input the objective function value at current solution, on output contains the objective function value at new solution
263: . g          - On input the gradient evaluated at `x`, on output contains the gradient at new solution
264: . steplength - scalar multiplier of s used ( x = x0 + steplength * x)
265: - reason     - `TaoLineSearchConvergedReason` reason why the line-search stopped

267:   Level: advanced

269:   Notes:
270:   The algorithm developer must set up the `TaoLineSearch` with calls to
271:   `TaoLineSearchSetObjectiveRoutine()` and `TaoLineSearchSetGradientRoutine()`,
272:   `TaoLineSearchSetObjectiveAndGradientRoutine()`, or `TaoLineSearchUseTaoRoutines()`.
273:   The latter is done automatically by default and thus requires no user input.

275:   You may or may not need to follow this with a call to
276:   `TaoAddLineSearchCounts()`, depending on whether you want these
277:   evaluations to count toward the total function/gradient evaluations.

279: .seealso: [](ch_tao), `Tao`, `TaoLineSearchConvergedReason`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetType()`,
280:           `TaoLineSearchSetInitialStepLength()`, `TaoAddLineSearchCounts()`
281: @*/
282: PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
283: {
284:   PetscInt low1, low2, low3, high1, high2, high3;

286:   PetscFunctionBegin;
289:   PetscAssertPointer(f, 3);
292:   PetscAssertPointer(reason, 7);
293:   PetscCheckSameComm(ls, 1, x, 2);
294:   PetscCheckSameTypeAndComm(x, 2, g, 4);
295:   PetscCheckSameTypeAndComm(x, 2, s, 5);
296:   PetscCall(VecGetOwnershipRange(x, &low1, &high1));
297:   PetscCall(VecGetOwnershipRange(g, &low2, &high2));
298:   PetscCall(VecGetOwnershipRange(s, &low3, &high3));
299:   PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible vector local lengths");

301:   *reason = TAOLINESEARCH_CONTINUE_ITERATING;
302:   PetscCall(PetscObjectReference((PetscObject)s));
303:   PetscCall(VecDestroy(&ls->stepdirection));
304:   ls->stepdirection = s;

306:   PetscCall(TaoLineSearchSetUp(ls));
307:   ls->nfeval  = 0;
308:   ls->ngeval  = 0;
309:   ls->nfgeval = 0;
310:   /* Check parameter values */
311:   if (ls->ftol < 0.0) {
312:     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: ftol (%g) < 0\n", (double)ls->ftol));
313:     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
314:   }
315:   if (ls->rtol < 0.0) {
316:     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: rtol (%g) < 0\n", (double)ls->rtol));
317:     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
318:   }
319:   if (ls->gtol < 0.0) {
320:     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: gtol (%g) < 0\n", (double)ls->gtol));
321:     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
322:   }
323:   if (ls->stepmin < 0.0) {
324:     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) < 0\n", (double)ls->stepmin));
325:     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
326:   }
327:   if (ls->stepmax < ls->stepmin) {
328:     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n", (double)ls->stepmin, (double)ls->stepmax));
329:     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
330:   }
331:   if (ls->max_funcs < 0) {
332:     PetscCall(PetscInfo(ls, "Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n", ls->max_funcs));
333:     *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
334:   }
335:   if (PetscIsInfOrNanReal(*f)) {
336:     PetscCall(PetscInfo(ls, "Initial Line Search Function Value is Inf or Nan (%g)\n", (double)*f));
337:     *reason = TAOLINESEARCH_FAILED_INFORNAN;
338:   }

340:   PetscCall(PetscObjectReference((PetscObject)x));
341:   PetscCall(VecDestroy(&ls->start_x));
342:   ls->start_x = x;

344:   PetscCall(PetscLogEventBegin(TAOLINESEARCH_Apply, ls, 0, 0, 0));
345:   PetscUseTypeMethod(ls, apply, x, f, g, s);
346:   PetscCall(PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0, 0, 0));
347:   *reason   = ls->reason;
348:   ls->new_f = *f;

350:   if (steplength) *steplength = ls->step;

352:   PetscCall(TaoLineSearchViewFromOptions(ls, NULL, "-tao_ls_view"));
353:   PetscFunctionReturn(PETSC_SUCCESS);
354: }

356: /*@C
357:   TaoLineSearchSetType - Sets the algorithm used in a line search

359:   Collective

361:   Input Parameters:
362: + ls   - the `TaoLineSearch` context
363: - type - the `TaoLineSearchType` selection

365:   Options Database Key:
366: . -tao_ls_type <type> - select which method Tao should use at runtime

368:   Level: beginner

370: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`, `TaoLineSearchCreate()`, `TaoLineSearchGetType()`,
371:           `TaoLineSearchApply()`
372: @*/
373: PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type)
374: {
375:   PetscErrorCode (*r)(TaoLineSearch);
376:   PetscBool flg;

378:   PetscFunctionBegin;
380:   PetscAssertPointer(type, 2);
381:   PetscCall(PetscObjectTypeCompare((PetscObject)ls, type, &flg));
382:   if (flg) PetscFunctionReturn(PETSC_SUCCESS);

384:   PetscCall(PetscFunctionListFind(TaoLineSearchList, type, (void (**)(void)) & r));
385:   PetscCheck(r, PetscObjectComm((PetscObject)ls), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested TaoLineSearch type %s", type);
386:   PetscTryTypeMethod(ls, destroy);
387:   ls->max_funcs = 30;
388:   ls->ftol      = 0.0001;
389:   ls->gtol      = 0.9;
390: #if defined(PETSC_USE_REAL_SINGLE)
391:   ls->rtol = 1.0e-5;
392: #else
393:   ls->rtol = 1.0e-10;
394: #endif
395:   ls->stepmin = 1.0e-20;
396:   ls->stepmax = 1.0e+20;

398:   ls->nfeval              = 0;
399:   ls->ngeval              = 0;
400:   ls->nfgeval             = 0;
401:   ls->ops->setup          = NULL;
402:   ls->ops->apply          = NULL;
403:   ls->ops->view           = NULL;
404:   ls->ops->setfromoptions = NULL;
405:   ls->ops->destroy        = NULL;
406:   ls->setupcalled         = PETSC_FALSE;
407:   PetscCall((*r)(ls));
408:   PetscCall(PetscObjectChangeTypeName((PetscObject)ls, type));
409:   PetscFunctionReturn(PETSC_SUCCESS);
410: }

412: /*@C
413:   TaoLineSearchMonitor - Monitor the line search steps. This routine will otuput the
414:   iteration number, step length, and function value before calling the implementation
415:   specific monitor.

417:   Input Parameters:
418: + ls   - the `TaoLineSearch` context
419: . its  - the current iterate number (>=0)
420: . f    - the current objective function value
421: - step - the step length

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

426:   Level: developer

428: .seealso: `TaoLineSearch`
429: @*/
430: PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step)
431: {
432:   PetscInt tabs;

434:   PetscFunctionBegin;
436:   if (ls->usemonitor) {
437:     PetscCall(PetscViewerASCIIGetTab(ls->viewer, &tabs));
438:     PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel));
439:     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "%3" PetscInt_FMT " LS", its));
440:     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "  Function value: %g,", (double)f));
441:     PetscCall(PetscViewerASCIIPrintf(ls->viewer, "  Step length: %g\n", (double)step));
442:     if (ls->ops->monitor && its > 0) {
443:       PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3));
444:       PetscUseTypeMethod(ls, monitor);
445:     }
446:     PetscCall(PetscViewerASCIISetTab(ls->viewer, tabs));
447:   }
448:   PetscFunctionReturn(PETSC_SUCCESS);
449: }

451: /*@
452:   TaoLineSearchSetFromOptions - Sets various `TaoLineSearch` parameters from user
453:   options.

455:   Collective

457:   Input Parameter:
458: . ls - the `TaoLineSearch` context

460:   Options Database Keys:
461: + -tao_ls_type <type>     - The algorithm that `TaoLineSearch` uses (more-thuente, gpcg, unit)
462: . -tao_ls_ftol <tol>      - tolerance for sufficient decrease
463: . -tao_ls_gtol <tol>      - tolerance for curvature condition
464: . -tao_ls_rtol <tol>      - relative tolerance for acceptable step
465: . -tao_ls_stepinit <step> - initial steplength allowed
466: . -tao_ls_stepmin <step>  - minimum steplength allowed
467: . -tao_ls_stepmax <step>  - maximum steplength allowed
468: . -tao_ls_max_funcs <n>   - maximum number of function evaluations allowed
469: - -tao_ls_view            - display line-search results to standard output

471:   Level: beginner

473: .seealso: `TaoLineSearch`
474: @*/
475: PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
476: {
477:   const char *default_type = TAOLINESEARCHMT;
478:   char        type[256], monfilename[PETSC_MAX_PATH_LEN];
479:   PetscViewer monviewer;
480:   PetscBool   flg;

482:   PetscFunctionBegin;
484:   PetscObjectOptionsBegin((PetscObject)ls);
485:   if (((PetscObject)ls)->type_name) default_type = ((PetscObject)ls)->type_name;
486:   /* Check for type from options */
487:   PetscCall(PetscOptionsFList("-tao_ls_type", "Tao Line Search type", "TaoLineSearchSetType", TaoLineSearchList, default_type, type, 256, &flg));
488:   if (flg) {
489:     PetscCall(TaoLineSearchSetType(ls, type));
490:   } else if (!((PetscObject)ls)->type_name) {
491:     PetscCall(TaoLineSearchSetType(ls, default_type));
492:   }

494:   PetscCall(PetscOptionsInt("-tao_ls_max_funcs", "max function evals in line search", "", ls->max_funcs, &ls->max_funcs, NULL));
495:   PetscCall(PetscOptionsReal("-tao_ls_ftol", "tol for sufficient decrease", "", ls->ftol, &ls->ftol, NULL));
496:   PetscCall(PetscOptionsReal("-tao_ls_gtol", "tol for curvature condition", "", ls->gtol, &ls->gtol, NULL));
497:   PetscCall(PetscOptionsReal("-tao_ls_rtol", "relative tol for acceptable step", "", ls->rtol, &ls->rtol, NULL));
498:   PetscCall(PetscOptionsReal("-tao_ls_stepmin", "lower bound for step", "", ls->stepmin, &ls->stepmin, NULL));
499:   PetscCall(PetscOptionsReal("-tao_ls_stepmax", "upper bound for step", "", ls->stepmax, &ls->stepmax, NULL));
500:   PetscCall(PetscOptionsReal("-tao_ls_stepinit", "initial step", "", ls->initstep, &ls->initstep, NULL));
501:   PetscCall(PetscOptionsString("-tao_ls_monitor", "enable the basic monitor", "TaoLineSearchSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg));
502:   if (flg) {
503:     PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls), monfilename, &monviewer));
504:     ls->viewer     = monviewer;
505:     ls->usemonitor = PETSC_TRUE;
506:   }
507:   PetscTryTypeMethod(ls, setfromoptions, PetscOptionsObject);
508:   PetscOptionsEnd();
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }

512: /*@C
513:   TaoLineSearchGetType - Gets the current line search algorithm

515:   Not Collective

517:   Input Parameter:
518: . ls - the `TaoLineSearch` context

520:   Output Parameter:
521: . type - the line search algorithm in effect

523:   Level: developer

525: .seealso: `TaoLineSearch`
526: @*/
527: PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type)
528: {
529:   PetscFunctionBegin;
531:   PetscAssertPointer(type, 2);
532:   *type = ((PetscObject)ls)->type_name;
533:   PetscFunctionReturn(PETSC_SUCCESS);
534: }

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

540:   Not Collective

542:   Input Parameter:
543: . ls - the `TaoLineSearch` context

545:   Output Parameters:
546: + nfeval  - number of function evaluations
547: . ngeval  - number of gradient evaluations
548: - nfgeval - number of function/gradient evaluations

550:   Level: intermediate

552:   Note:
553:   If the line search is using the `Tao` objective and gradient
554:   routines directly (see `TaoLineSearchUseTaoRoutines()`), then the `Tao`
555:   is already counting the number of evaluations.

557: .seealso: `TaoLineSearch`
558: @*/
559: PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
560: {
561:   PetscFunctionBegin;
563:   *nfeval  = ls->nfeval;
564:   *ngeval  = ls->ngeval;
565:   *nfgeval = ls->nfgeval;
566:   PetscFunctionReturn(PETSC_SUCCESS);
567: }

569: /*@
570:   TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
571:   the standard `Tao` evaluation routines.

573:   Not Collective

575:   Input Parameter:
576: . ls - the `TaoLineSearch` context

578:   Output Parameter:
579: . flg - `PETSC_TRUE` if the line search is using `Tao` evaluation routines,
580:         otherwise `PETSC_FALSE`

582:   Level: developer

584: .seealso: `TaoLineSearch`
585: @*/
586: PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
587: {
588:   PetscFunctionBegin;
590:   *flg = ls->usetaoroutines;
591:   PetscFunctionReturn(PETSC_SUCCESS);
592: }

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

597:   Logically Collective

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

604:   Calling sequence of `func`:
605: + ls  - the line search context
606: . x   - input vector
607: . f   - function value
608: - ctx - (optional) user-defined context

610:   Level: advanced

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

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

622: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
623: @*/
624: PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *f, void *ctx), void *ctx)
625: {
626:   PetscFunctionBegin;

629:   ls->ops->computeobjective = func;
630:   if (ctx) ls->userctx_func = ctx;
631:   ls->usetaoroutines = PETSC_FALSE;
632:   PetscFunctionReturn(PETSC_SUCCESS);
633: }

635: /*@C
636:   TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search

638:   Logically Collective

640:   Input Parameters:
641: + ls   - the `TaoLineSearch` context
642: . func - the gradient evaluation routine
643: - ctx  - the (optional) user-defined context for private data

645:   Calling sequence of `func`:
646: + ls  - the linesearch object
647: . x   - input vector
648: . g   - gradient vector
649: - ctx - (optional) user-defined context

651:   Level: beginner

653:   Note:
654:   Use this routine only if you want the line search gradient
655:   evaluation routine to be different from the `Tao`'s gradient
656:   evaluation routine. If you use this routine you must also set
657:   the line search function and/or function/gradient routine.

659:   Some algorithms (lcl, gpcg) set their own gradient routine for the
660:   line search, application programmers should be wary of overriding the
661:   default gradient routine.

663: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
664: @*/
665: PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec g, void *ctx), void *ctx)
666: {
667:   PetscFunctionBegin;
669:   ls->ops->computegradient = func;
670:   if (ctx) ls->userctx_grad = ctx;
671:   ls->usetaoroutines = PETSC_FALSE;
672:   PetscFunctionReturn(PETSC_SUCCESS);
673: }

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

678:   Logically Collective

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

685:   Calling sequence of `func`:
686: + ls  - the linesearch object
687: . x   - input vector
688: . f   - function value
689: . g   - gradient vector
690: - ctx - (optional) user-defined context

692:   Level: beginner

694:   Note:
695:   Use this routine only if you want the line search objective and gradient
696:   evaluation routines to be different from the `Tao`'s objective
697:   and gradient evaluation routines.

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

703: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
704: @*/
705: PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, void *ctx), void *ctx)
706: {
707:   PetscFunctionBegin;
709:   ls->ops->computeobjectiveandgradient = func;
710:   if (ctx) ls->userctx_funcgrad = ctx;
711:   ls->usetaoroutines = PETSC_FALSE;
712:   PetscFunctionReturn(PETSC_SUCCESS);
713: }

715: /*@C
716:   TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
717:   (gradient'*stepdirection) evaluation routine for the line search.

719:   Logically Collective

721:   Input Parameters:
722: + ls   - the `TaoLineSearch` context
723: . func - the objective and gradient evaluation routine
724: - ctx  - the (optional) user-defined context for private data

726:   Calling sequence of `func`:
727: + ls  - the linesearch context
728: . x   - input vector
729: . s   - step direction
730: . f   - function value
731: . gts - inner product of gradient and step direction vectors
732: - ctx - (optional) user-defined context

734:   Level: advanced

736:   Notes:
737:   Sometimes it is more efficient to compute the inner product of the gradient and the step
738:   direction than it is to compute the gradient, and this is all the line search typically needs
739:   of the gradient.

741:   The gradient will still need to be computed at the end of the line
742:   search, so you will still need to set a line search gradient evaluation
743:   routine

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

749:   Some algorithms (lcl, gpcg) set their own objective routine for the
750:   line search, application programmers should be wary of overriding the
751:   default objective routine.

753: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjective()`, `TaoLineSearchSetGradient()`, `TaoLineSearchUseTaoRoutines()`
754: @*/
755: PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *f, PetscReal *gts, void *ctx), void *ctx)
756: {
757:   PetscFunctionBegin;
759:   ls->ops->computeobjectiveandgts = func;
760:   if (ctx) ls->userctx_funcgts = ctx;
761:   ls->usegts         = PETSC_TRUE;
762:   ls->usetaoroutines = PETSC_FALSE;
763:   PetscFunctionReturn(PETSC_SUCCESS);
764: }

766: /*@C
767:   TaoLineSearchUseTaoRoutines - Informs the `TaoLineSearch` to use the
768:   objective and gradient evaluation routines from the given `Tao` object. The default.

770:   Logically Collective

772:   Input Parameters:
773: + ls - the `TaoLineSearch` context
774: - ts - the `Tao` context with defined objective/gradient evaluation routines

776:   Level: developer

778: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`
779: @*/
780: PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
781: {
782:   PetscFunctionBegin;
785:   ls->tao            = ts;
786:   ls->usetaoroutines = PETSC_TRUE;
787:   PetscFunctionReturn(PETSC_SUCCESS);
788: }

790: /*@
791:   TaoLineSearchComputeObjective - Computes the objective function value at a given point

793:   Collective

795:   Input Parameters:
796: + ls - the `TaoLineSearch` context
797: - x  - input vector

799:   Output Parameter:
800: . f - Objective value at `x`

802:   Level: developer

804:   Note:
805:   `TaoLineSearchComputeObjective()` is typically used within line searches
806:   so most users would not generally call this routine themselves.

808: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
809: @*/
810: PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
811: {
812:   Vec       gdummy;
813:   PetscReal gts;

815:   PetscFunctionBegin;
818:   PetscAssertPointer(f, 3);
819:   PetscCheckSameComm(ls, 1, x, 2);
820:   if (ls->usetaoroutines) {
821:     PetscCall(TaoComputeObjective(ls->tao, x, f));
822:   } else {
823:     PetscCheck(ls->ops->computeobjective || ls->ops->computeobjectiveandgradient || ls->ops->computeobjectiveandgts, PetscObjectComm((PetscObject)ls), PETSC_ERR_ARG_WRONGSTATE, "Line Search does not have objective function set");
824:     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
825:     if (ls->ops->computeobjective) PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func));
826:     else if (ls->ops->computeobjectiveandgradient) {
827:       PetscCall(VecDuplicate(x, &gdummy));
828:       PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgradient)(ls, x, f, gdummy, ls->userctx_funcgrad));
829:       PetscCall(VecDestroy(&gdummy));
830:     } else PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, &gts, ls->userctx_funcgts));
831:     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
832:   }
833:   ls->nfeval++;
834:   PetscFunctionReturn(PETSC_SUCCESS);
835: }

837: /*@
838:   TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point

840:   Collective

842:   Input Parameters:
843: + ls - the `TaoLineSearch` context
844: - x  - input vector

846:   Output Parameters:
847: + f - Objective value at `x`
848: - g - Gradient vector at `x`

850:   Level: developer

852:   Note:
853:   `TaoLineSearchComputeObjectiveAndGradient()` is typically used within line searches
854:   so most users would not generally call this routine themselves.

856: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchSetObjectiveRoutine()`
857: @*/
858: PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
859: {
860:   PetscFunctionBegin;
863:   PetscAssertPointer(f, 3);
865:   PetscCheckSameComm(ls, 1, x, 2);
866:   PetscCheckSameComm(ls, 1, g, 4);
867:   if (ls->usetaoroutines) {
868:     PetscCall(TaoComputeObjectiveAndGradient(ls->tao, x, f, g));
869:   } else {
870:     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
871:     if (ls->ops->computeobjectiveandgradient) PetscCallBack("TaoLineSearch callback objective/gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, f, g, ls->userctx_funcgrad));
872:     else {
873:       PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func));
874:       PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad));
875:     }
876:     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
877:     PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f)));
878:   }
879:   ls->nfgeval++;
880:   PetscFunctionReturn(PETSC_SUCCESS);
881: }

883: /*@
884:   TaoLineSearchComputeGradient - Computes the gradient of the objective function

886:   Collective

888:   Input Parameters:
889: + ls - the `TaoLineSearch` context
890: - x  - input vector

892:   Output Parameter:
893: . g - gradient vector

895:   Level: developer

897:   Note:
898:   `TaoComputeGradient()` is typically used within line searches
899:   so most users would not generally call this routine themselves.

901: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeObjective()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetGradient()`
902: @*/
903: PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
904: {
905:   PetscReal fdummy;

907:   PetscFunctionBegin;
911:   PetscCheckSameComm(ls, 1, x, 2);
912:   PetscCheckSameComm(ls, 1, g, 3);
913:   if (ls->usetaoroutines) {
914:     PetscCall(TaoComputeGradient(ls->tao, x, g));
915:   } else {
916:     PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
917:     if (ls->ops->computegradient) PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad));
918:     else PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, &fdummy, g, ls->userctx_funcgrad));
919:     PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
920:   }
921:   ls->ngeval++;
922:   PetscFunctionReturn(PETSC_SUCCESS);
923: }

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

929:   Collective

931:   Input Parameters:
932: + ls - the `TaoLineSearch` context
933: - x  - input vector

935:   Output Parameters:
936: + f   - Objective value at `x`
937: - gts - inner product of gradient and step direction at `x`

939:   Level: developer

941:   Note:
942:   `TaoLineSearchComputeObjectiveAndGTS()` is typically used within line searches
943:   so most users would not generally call this routine themselves.

945: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
946: @*/
947: PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
948: {
949:   PetscFunctionBegin;
952:   PetscAssertPointer(f, 3);
953:   PetscAssertPointer(gts, 4);
954:   PetscCheckSameComm(ls, 1, x, 2);
955:   PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
956:   PetscCallBack("TaoLineSearch callback objective/gts", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, gts, ls->userctx_funcgts));
957:   PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
958:   PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f)));
959:   ls->nfeval++;
960:   PetscFunctionReturn(PETSC_SUCCESS);
961: }

963: /*@
964:   TaoLineSearchGetSolution - Returns the solution to the line search

966:   Collective

968:   Input Parameter:
969: . ls - the `TaoLineSearch` context

971:   Output Parameters:
972: + x          - the new solution
973: . f          - the objective function value at `x`
974: . g          - the gradient at `x`
975: . steplength - the multiple of the step direction taken by the line search
976: - reason     - the reason why the line search terminated

978:   Level: developer

980: .seealso: `TaoLineSearchGetStartingVector()`, `TaoLineSearchGetStepDirection()`
981: @*/
982: PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
983: {
984:   PetscFunctionBegin;
987:   PetscAssertPointer(f, 3);
989:   PetscAssertPointer(reason, 6);
990:   if (ls->new_x) PetscCall(VecCopy(ls->new_x, x));
991:   *f = ls->new_f;
992:   if (ls->new_g) PetscCall(VecCopy(ls->new_g, g));
993:   if (steplength) *steplength = ls->step;
994:   *reason = ls->reason;
995:   PetscFunctionReturn(PETSC_SUCCESS);
996: }

998: /*@
999:   TaoLineSearchGetStartingVector - Gets a the initial point of the line
1000:   search.

1002:   Not Collective

1004:   Input Parameter:
1005: . ls - the `TaoLineSearch` context

1007:   Output Parameter:
1008: . x - The initial point of the line search

1010:   Level: advanced

1012: .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStepDirection()`
1013: @*/
1014: PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1015: {
1016:   PetscFunctionBegin;
1018:   if (x) *x = ls->start_x;
1019:   PetscFunctionReturn(PETSC_SUCCESS);
1020: }

1022: /*@
1023:   TaoLineSearchGetStepDirection - Gets the step direction of the line
1024:   search.

1026:   Not Collective

1028:   Input Parameter:
1029: . ls - the `TaoLineSearch` context

1031:   Output Parameter:
1032: . s - the step direction of the line search

1034:   Level: advanced

1036: .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStartingVector()`
1037: @*/
1038: PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1039: {
1040:   PetscFunctionBegin;
1042:   if (s) *s = ls->stepdirection;
1043:   PetscFunctionReturn(PETSC_SUCCESS);
1044: }

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

1049:   Not Collective

1051:   Input Parameter:
1052: . ls - the `TaoLineSearch` context

1054:   Output Parameter:
1055: . f_fullstep - the objective value at the full step length

1057:   Level: developer

1059: .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStartingVector()`, `TaoLineSearchGetStepDirection()`
1060: @*/
1061: PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1062: {
1063:   PetscFunctionBegin;
1065:   *f_fullstep = ls->f_fullstep;
1066:   PetscFunctionReturn(PETSC_SUCCESS);
1067: }

1069: /*@
1070:   TaoLineSearchSetVariableBounds - Sets the upper and lower bounds for a bounded line search

1072:   Logically Collective

1074:   Input Parameters:
1075: + ls - the `TaoLineSearch` context
1076: . xl - vector of lower bounds
1077: - xu - vector of upper bounds

1079:   Level: beginner

1081:   Note:
1082:   If the variable bounds are not set with this routine, then
1083:   `PETSC_NINFINITY` and `PETSC_INFINITY` are assumed

1085: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoSetVariableBounds()`, `TaoLineSearchCreate()`
1086: @*/
1087: PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls, Vec xl, Vec xu)
1088: {
1089:   PetscFunctionBegin;
1093:   PetscCall(PetscObjectReference((PetscObject)xl));
1094:   PetscCall(PetscObjectReference((PetscObject)xu));
1095:   PetscCall(VecDestroy(&ls->lower));
1096:   PetscCall(VecDestroy(&ls->upper));
1097:   ls->lower   = xl;
1098:   ls->upper   = xu;
1099:   ls->bounded = (PetscBool)(xl || xu);
1100:   PetscFunctionReturn(PETSC_SUCCESS);
1101: }

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

1107:   Logically Collective

1109:   Input Parameters:
1110: + ls - the `TaoLineSearch` context
1111: - s  - the initial step size

1113:   Level: intermediate

1115: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchGetStepLength()`, `TaoLineSearchApply()`
1116: @*/
1117: PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls, PetscReal s)
1118: {
1119:   PetscFunctionBegin;
1122:   ls->initstep = s;
1123:   PetscFunctionReturn(PETSC_SUCCESS);
1124: }

1126: /*@
1127:   TaoLineSearchGetStepLength - Get the current step length

1129:   Not Collective

1131:   Input Parameter:
1132: . ls - the `TaoLineSearch` context

1134:   Output Parameter:
1135: . s - the current step length

1137:   Level: intermediate

1139: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetInitialStepLength()`, `TaoLineSearchApply()`
1140: @*/
1141: PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls, PetscReal *s)
1142: {
1143:   PetscFunctionBegin;
1145:   *s = ls->step;
1146:   PetscFunctionReturn(PETSC_SUCCESS);
1147: }

1149: /*@C
1150:   TaoLineSearchRegister - Adds a line-search algorithm to the registry

1152:   Not Collective

1154:   Input Parameters:
1155: + sname - name of a new user-defined solver
1156: - func  - routine to Create method context

1158:   Example Usage:
1159: .vb
1160:    TaoLineSearchRegister("my_linesearch", MyLinesearchCreate);
1161: .ve

1163:   Then, your solver can be chosen with the procedural interface via
1164: $     TaoLineSearchSetType(ls, "my_linesearch")
1165:   or at runtime via the option
1166: $     -tao_ls_type my_linesearch

1168:   Level: developer

1170:   Note:
1171:   `TaoLineSearchRegister()` may be called multiple times to add several user-defined solvers.

1173: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`
1174: @*/
1175: PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1176: {
1177:   PetscFunctionBegin;
1178:   PetscCall(TaoLineSearchInitializePackage());
1179:   PetscCall(PetscFunctionListAdd(&TaoLineSearchList, sname, (void (*)(void))func));
1180:   PetscFunctionReturn(PETSC_SUCCESS);
1181: }

1183: /*@C
1184:   TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1185:   for all `TaoLineSearch` options in the database.

1187:   Collective

1189:   Input Parameters:
1190: + ls - the `TaoLineSearch` solver context
1191: - p  - the prefix string to prepend to all line search requests

1193:   Level: advanced

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

1199:   This is inherited from the `Tao` object so rarely needs to be set

1201: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()`
1202: @*/
1203: PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1204: {
1205:   return PetscObjectAppendOptionsPrefix((PetscObject)ls, p);
1206: }

1208: /*@C
1209:   TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1210:   `TaoLineSearch` options in the database

1212:   Not Collective

1214:   Input Parameter:
1215: . ls - the `TaoLineSearch` context

1217:   Output Parameter:
1218: . p - pointer to the prefix string used is returned

1220:   Level: advanced

1222:   Fortran Notes:
1223:   The user should pass in a string 'prefix' of
1224:   sufficient length to hold the prefix.

1226: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchAppendOptionsPrefix()`
1227: @*/
1228: PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1229: {
1230:   return PetscObjectGetOptionsPrefix((PetscObject)ls, p);
1231: }

1233: /*@C
1234:   TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1235:   `TaoLineSearch` options in the database.

1237:   Logically Collective

1239:   Input Parameters:
1240: + ls - the `TaoLineSearch` context
1241: - p  - the prefix string to prepend to all `ls` option requests

1243:   Level: advanced

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

1249:   This is inherited from the `Tao` object so rarely needs to be set

1251:   For example, to distinguish between the runtime options for two
1252:   different line searches, one could call
1253: .vb
1254:       TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1255:       TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1256: .ve

1258:   This would enable use of different options for each system, such as
1259: .vb
1260:       -sys1_tao_ls_type mt
1261:       -sys2_tao_ls_type armijo
1262: .ve

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