1 #include <petsctaolinesearch.h> /*I "petsctaolinesearch.h" I*/
2 #include <petsc/private/taolinesearchimpl.h>
3
4 PetscFunctionList TaoLineSearchList = NULL;
5
6 PetscClassId TAOLINESEARCH_CLASSID = 0;
7
8 PetscLogEvent TAOLINESEARCH_Apply;
9 PetscLogEvent TAOLINESEARCH_Eval;
10
11 const char *const TaoLineSearchConvergedReasons_Shifted[] = {"FAILED_ASCENT", "FAILED_BADPARAMETER", "FAILED_INFORNAN", "CONTINUE_ITERATING", "SUCCESS", "SUCCESS_USER", "HALTED_OTHER", "HALTED_MAXFCN", "HALTED_UPPERBOUND", "HALTED_LOWERBOUND", "HALTED_RTOL", "HALTED_USER", "TaoLineSearchConvergedReason", "TAOLINESEARCH_", NULL};
12 const char *const *TaoLineSearchConvergedReasons = TaoLineSearchConvergedReasons_Shifted + 3;
13
14 /*@
15 TaoLineSearchViewFromOptions - View a `TaoLineSearch` object based on values in the options database
16
17 Collective
18
19 Input Parameters:
20 + A - the `Tao` context
21 . obj - Optional object
22 - name - command line option
23
24 Level: intermediate
25
26 Note:
27 See `PetscObjectViewFromOptions()` for available viewer options
28
29 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchView()`, `PetscObjectViewFromOptions()`, `TaoLineSearchCreate()`
30 @*/
TaoLineSearchViewFromOptions(TaoLineSearch A,PetscObject obj,const char name[])31 PetscErrorCode TaoLineSearchViewFromOptions(TaoLineSearch A, PetscObject obj, const char name[])
32 {
33 PetscFunctionBegin;
34 PetscValidHeaderSpecific(A, TAOLINESEARCH_CLASSID, 1);
35 PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
36 PetscFunctionReturn(PETSC_SUCCESS);
37 }
38
39 /*@
40 TaoLineSearchView - Prints information about the `TaoLineSearch`
41
42 Collective
43
44 Input Parameters:
45 + ls - the `TaoLineSearch` context
46 - viewer - visualization context
47
48 Options Database Key:
49 . -tao_ls_view - Calls `TaoLineSearchView()` at the end of each line search
50
51 Level: beginner
52
53 Notes:
54 The available visualization contexts include
55 + `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
56 - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
57 output where only the first processor opens
58 the file. All other processors send their
59 data to the first processor to print.
60
61 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `PetscViewerASCIIOpen()`, `TaoLineSearchViewFromOptions()`
62 @*/
TaoLineSearchView(TaoLineSearch ls,PetscViewer viewer)63 PetscErrorCode TaoLineSearchView(TaoLineSearch ls, PetscViewer viewer)
64 {
65 PetscBool isascii, isstring;
66 TaoLineSearchType type;
67
68 PetscFunctionBegin;
69 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
70 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)ls)->comm, &viewer));
71 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
72 PetscCheckSameComm(ls, 1, viewer, 2);
73
74 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
75 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
76 if (isascii) {
77 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ls, viewer));
78 PetscCall(PetscViewerASCIIPushTab(viewer));
79 PetscTryTypeMethod(ls, view, viewer);
80 PetscCall(PetscViewerASCIIPopTab(viewer));
81 PetscCall(PetscViewerASCIIPushTab(viewer));
82 PetscCall(PetscViewerASCIIPrintf(viewer, "maximum function evaluations=%" PetscInt_FMT "\n", ls->max_funcs));
83 PetscCall(PetscViewerASCIIPrintf(viewer, "tolerances: ftol=%g, rtol=%g, gtol=%g\n", (double)ls->ftol, (double)ls->rtol, (double)ls->gtol));
84 PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT "\n", ls->nfeval));
85 PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT "\n", ls->ngeval));
86 PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT "\n", ls->nfgeval));
87
88 if (ls->bounded) PetscCall(PetscViewerASCIIPrintf(viewer, "using variable bounds\n"));
89 PetscCall(PetscViewerASCIIPrintf(viewer, "Termination reason: %s\n", TaoLineSearchConvergedReasons[ls->reason]));
90 PetscCall(PetscViewerASCIIPopTab(viewer));
91 } else if (isstring) {
92 PetscCall(TaoLineSearchGetType(ls, &type));
93 PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type));
94 }
95 PetscFunctionReturn(PETSC_SUCCESS);
96 }
97
98 /*@C
99 TaoLineSearchCreate - Creates a `TaoLineSearch` object. Algorithms in `Tao` that use
100 line-searches will automatically create one so this all is rarely needed
101
102 Collective
103
104 Input Parameter:
105 . comm - MPI communicator
106
107 Output Parameter:
108 . newls - the new `TaoLineSearch` context
109
110 Options Database Key:
111 . -tao_ls_type - select which method `Tao` should use
112
113 Level: developer
114
115 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`, `TaoLineSearchDestroy()`
116 @*/
TaoLineSearchCreate(MPI_Comm comm,TaoLineSearch * newls)117 PetscErrorCode TaoLineSearchCreate(MPI_Comm comm, TaoLineSearch *newls)
118 {
119 TaoLineSearch ls;
120
121 PetscFunctionBegin;
122 PetscAssertPointer(newls, 2);
123 PetscCall(TaoLineSearchInitializePackage());
124
125 PetscCall(PetscHeaderCreate(ls, TAOLINESEARCH_CLASSID, "TaoLineSearch", "Linesearch", "Tao", comm, TaoLineSearchDestroy, TaoLineSearchView));
126 ls->max_funcs = 30;
127 ls->ftol = 0.0001;
128 ls->gtol = 0.9;
129 #if defined(PETSC_USE_REAL_SINGLE)
130 ls->rtol = 1.0e-5;
131 #else
132 ls->rtol = 1.0e-10;
133 #endif
134 ls->stepmin = 1.0e-20;
135 ls->stepmax = 1.0e+20;
136 ls->step = 1.0;
137 ls->initstep = 1.0;
138 *newls = ls;
139 PetscFunctionReturn(PETSC_SUCCESS);
140 }
141
142 /*@
143 TaoLineSearchSetUp - Sets up the internal data structures for the later use
144 of a `TaoLineSearch`
145
146 Collective
147
148 Input Parameter:
149 . ls - the `TaoLineSearch` context
150
151 Level: developer
152
153 Note:
154 The user will not need to explicitly call `TaoLineSearchSetUp()`, as it will
155 automatically be called in `TaoLineSearchSolve()`. However, if the user
156 desires to call it explicitly, it should come after `TaoLineSearchCreate()`
157 but before `TaoLineSearchApply()`.
158
159 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchApply()`
160 @*/
TaoLineSearchSetUp(TaoLineSearch ls)161 PetscErrorCode TaoLineSearchSetUp(TaoLineSearch ls)
162 {
163 const char *default_type = TAOLINESEARCHMT;
164 PetscBool flg;
165
166 PetscFunctionBegin;
167 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
168 if (ls->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
169 if (!((PetscObject)ls)->type_name) PetscCall(TaoLineSearchSetType(ls, default_type));
170 PetscTryTypeMethod(ls, setup);
171 if (ls->usetaoroutines) {
172 PetscCall(TaoIsObjectiveDefined(ls->tao, &flg));
173 ls->hasobjective = flg;
174 PetscCall(TaoIsGradientDefined(ls->tao, &flg));
175 ls->hasgradient = flg;
176 PetscCall(TaoIsObjectiveAndGradientDefined(ls->tao, &flg));
177 ls->hasobjectiveandgradient = flg;
178 } else {
179 if (ls->ops->computeobjective) {
180 ls->hasobjective = PETSC_TRUE;
181 } else {
182 ls->hasobjective = PETSC_FALSE;
183 }
184 if (ls->ops->computegradient) {
185 ls->hasgradient = PETSC_TRUE;
186 } else {
187 ls->hasgradient = PETSC_FALSE;
188 }
189 if (ls->ops->computeobjectiveandgradient) {
190 ls->hasobjectiveandgradient = PETSC_TRUE;
191 } else {
192 ls->hasobjectiveandgradient = PETSC_FALSE;
193 }
194 }
195 ls->setupcalled = PETSC_TRUE;
196 PetscFunctionReturn(PETSC_SUCCESS);
197 }
198
199 /*@
200 TaoLineSearchReset - Some line searches may carry state information
201 from one `TaoLineSearchApply()` to the next. This function resets this
202 state information.
203
204 Collective
205
206 Input Parameter:
207 . ls - the `TaoLineSearch` context
208
209 Level: developer
210
211 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchApply()`
212 @*/
TaoLineSearchReset(TaoLineSearch ls)213 PetscErrorCode TaoLineSearchReset(TaoLineSearch ls)
214 {
215 PetscFunctionBegin;
216 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
217 PetscTryTypeMethod(ls, reset);
218 PetscFunctionReturn(PETSC_SUCCESS);
219 }
220
221 /*@
222 TaoLineSearchDestroy - Destroys the `TaoLineSearch` context that was created with
223 `TaoLineSearchCreate()`
224
225 Collective
226
227 Input Parameter:
228 . ls - the `TaoLineSearch` context
229
230 Level: developer
231
232 .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSolve()`
233 @*/
TaoLineSearchDestroy(TaoLineSearch * ls)234 PetscErrorCode TaoLineSearchDestroy(TaoLineSearch *ls)
235 {
236 PetscFunctionBegin;
237 if (!*ls) PetscFunctionReturn(PETSC_SUCCESS);
238 PetscValidHeaderSpecific(*ls, TAOLINESEARCH_CLASSID, 1);
239 if (--((PetscObject)*ls)->refct > 0) {
240 *ls = NULL;
241 PetscFunctionReturn(PETSC_SUCCESS);
242 }
243 PetscCall(VecDestroy(&(*ls)->stepdirection));
244 PetscCall(VecDestroy(&(*ls)->start_x));
245 PetscCall(VecDestroy(&(*ls)->upper));
246 PetscCall(VecDestroy(&(*ls)->lower));
247 PetscTryTypeMethod(*ls, destroy);
248 if ((*ls)->usemonitor) PetscCall(PetscViewerDestroy(&(*ls)->viewer));
249 PetscCall(PetscHeaderDestroy(ls));
250 PetscFunctionReturn(PETSC_SUCCESS);
251 }
252
253 /*@
254 TaoLineSearchApply - Performs a line-search in a given step direction.
255 Criteria for acceptable step length depends on the line-search algorithm chosen
256
257 Collective
258
259 Input Parameters:
260 + ls - the `TaoLineSearch` context
261 - s - search direction
262
263 Output Parameters:
264 + x - On input the current solution, on output `x` contains the new solution determined by the line search
265 . f - On input the objective function value at current solution, on output contains the objective function value at new solution
266 . g - On input the gradient evaluated at `x`, on output contains the gradient at new solution
267 . steplength - scalar multiplier of s used ( x = x0 + steplength * x)
268 - reason - `TaoLineSearchConvergedReason` reason why the line-search stopped
269
270 Level: advanced
271
272 Notes:
273 The algorithm developer must set up the `TaoLineSearch` with calls to
274 `TaoLineSearchSetObjectiveRoutine()` and `TaoLineSearchSetGradientRoutine()`,
275 `TaoLineSearchSetObjectiveAndGradientRoutine()`, or `TaoLineSearchUseTaoRoutines()`.
276 The latter is done automatically by default and thus requires no user input.
277
278 You may or may not need to follow this with a call to
279 `TaoAddLineSearchCounts()`, depending on whether you want these
280 evaluations to count toward the total function/gradient evaluations.
281
282 .seealso: [](ch_tao), `Tao`, `TaoLineSearchConvergedReason`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetType()`,
283 `TaoLineSearchSetInitialStepLength()`, `TaoAddLineSearchCounts()`
284 @*/
TaoLineSearchApply(TaoLineSearch ls,Vec x,PetscReal * f,Vec g,Vec s,PetscReal * steplength,TaoLineSearchConvergedReason * reason)285 PetscErrorCode TaoLineSearchApply(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
286 {
287 PetscInt low1, low2, low3, high1, high2, high3;
288
289 PetscFunctionBegin;
290 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
291 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
292 PetscAssertPointer(f, 3);
293 PetscValidHeaderSpecific(g, VEC_CLASSID, 4);
294 PetscValidHeaderSpecific(s, VEC_CLASSID, 5);
295 PetscAssertPointer(reason, 7);
296 PetscCheckSameComm(ls, 1, x, 2);
297 PetscCheckSameTypeAndComm(x, 2, g, 4);
298 PetscCheckSameTypeAndComm(x, 2, s, 5);
299 PetscCall(VecGetOwnershipRange(x, &low1, &high1));
300 PetscCall(VecGetOwnershipRange(g, &low2, &high2));
301 PetscCall(VecGetOwnershipRange(s, &low3, &high3));
302 PetscCheck(low1 == low2 && low1 == low3 && high1 == high2 && high1 == high3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible vector local lengths");
303
304 *reason = TAOLINESEARCH_CONTINUE_ITERATING;
305 PetscCall(PetscObjectReference((PetscObject)s));
306 PetscCall(VecDestroy(&ls->stepdirection));
307 ls->stepdirection = s;
308
309 PetscCall(TaoLineSearchSetUp(ls));
310 ls->nfeval = 0;
311 ls->ngeval = 0;
312 ls->nfgeval = 0;
313 /* Check parameter values */
314 if (ls->ftol < 0.0) {
315 PetscCall(PetscInfo(ls, "Bad Line Search Parameter: ftol (%g) < 0\n", (double)ls->ftol));
316 *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
317 }
318 if (ls->rtol < 0.0) {
319 PetscCall(PetscInfo(ls, "Bad Line Search Parameter: rtol (%g) < 0\n", (double)ls->rtol));
320 *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
321 }
322 if (ls->gtol < 0.0) {
323 PetscCall(PetscInfo(ls, "Bad Line Search Parameter: gtol (%g) < 0\n", (double)ls->gtol));
324 *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
325 }
326 if (ls->stepmin < 0.0) {
327 PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) < 0\n", (double)ls->stepmin));
328 *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
329 }
330 if (ls->stepmax < ls->stepmin) {
331 PetscCall(PetscInfo(ls, "Bad Line Search Parameter: stepmin (%g) > stepmax (%g)\n", (double)ls->stepmin, (double)ls->stepmax));
332 *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
333 }
334 if (ls->max_funcs < 0) {
335 PetscCall(PetscInfo(ls, "Bad Line Search Parameter: max_funcs (%" PetscInt_FMT ") < 0\n", ls->max_funcs));
336 *reason = TAOLINESEARCH_FAILED_BADPARAMETER;
337 }
338 if (PetscIsInfOrNanReal(*f)) {
339 PetscCall(PetscInfo(ls, "Initial Line Search Function Value is infinity or NaN (%g)\n", (double)*f));
340 *reason = TAOLINESEARCH_FAILED_INFORNAN;
341 }
342
343 PetscCall(PetscObjectReference((PetscObject)x));
344 PetscCall(VecDestroy(&ls->start_x));
345 ls->start_x = x;
346
347 PetscCall(PetscLogEventBegin(TAOLINESEARCH_Apply, ls, 0, 0, 0));
348 PetscUseTypeMethod(ls, apply, x, f, g, s);
349 PetscCall(PetscLogEventEnd(TAOLINESEARCH_Apply, ls, 0, 0, 0));
350 *reason = ls->reason;
351 ls->new_f = *f;
352
353 if (steplength) *steplength = ls->step;
354
355 PetscCall(TaoLineSearchViewFromOptions(ls, NULL, "-tao_ls_view"));
356 PetscFunctionReturn(PETSC_SUCCESS);
357 }
358
359 /*@
360 TaoLineSearchSetType - Sets the algorithm used in a line search
361
362 Collective
363
364 Input Parameters:
365 + ls - the `TaoLineSearch` context
366 - type - the `TaoLineSearchType` selection
367
368 Options Database Key:
369 . -tao_ls_type <type> - select which method Tao should use at runtime
370
371 Level: beginner
372
373 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`, `TaoLineSearchCreate()`, `TaoLineSearchGetType()`,
374 `TaoLineSearchApply()`
375 @*/
TaoLineSearchSetType(TaoLineSearch ls,TaoLineSearchType type)376 PetscErrorCode TaoLineSearchSetType(TaoLineSearch ls, TaoLineSearchType type)
377 {
378 PetscErrorCode (*r)(TaoLineSearch);
379 PetscBool flg;
380
381 PetscFunctionBegin;
382 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
383 PetscAssertPointer(type, 2);
384 PetscCall(PetscObjectTypeCompare((PetscObject)ls, type, &flg));
385 if (flg) PetscFunctionReturn(PETSC_SUCCESS);
386
387 PetscCall(PetscFunctionListFind(TaoLineSearchList, type, &r));
388 PetscCheck(r, PetscObjectComm((PetscObject)ls), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested TaoLineSearch type %s", type);
389 PetscTryTypeMethod(ls, destroy);
390 ls->max_funcs = 30;
391 ls->ftol = 0.0001;
392 ls->gtol = 0.9;
393 #if defined(PETSC_USE_REAL_SINGLE)
394 ls->rtol = 1.0e-5;
395 #else
396 ls->rtol = 1.0e-10;
397 #endif
398 ls->stepmin = 1.0e-20;
399 ls->stepmax = 1.0e+20;
400
401 ls->nfeval = 0;
402 ls->ngeval = 0;
403 ls->nfgeval = 0;
404 ls->ops->setup = NULL;
405 ls->ops->apply = NULL;
406 ls->ops->view = NULL;
407 ls->ops->setfromoptions = NULL;
408 ls->ops->destroy = NULL;
409 ls->setupcalled = PETSC_FALSE;
410 PetscCall((*r)(ls));
411 PetscCall(PetscObjectChangeTypeName((PetscObject)ls, type));
412 PetscFunctionReturn(PETSC_SUCCESS);
413 }
414
415 /*@C
416 TaoLineSearchMonitor - Monitor the line search steps. This routine will output the
417 iteration number, step length, and function value before calling the implementation
418 specific monitor.
419
420 Input Parameters:
421 + ls - the `TaoLineSearch` context
422 . its - the current iterate number (>=0)
423 . f - the current objective function value
424 - step - the step length
425
426 Options Database Key:
427 . -tao_ls_monitor - Use the default monitor, which prints statistics to standard output
428
429 Level: developer
430
431 .seealso: `TaoLineSearch`
432 @*/
TaoLineSearchMonitor(TaoLineSearch ls,PetscInt its,PetscReal f,PetscReal step)433 PetscErrorCode TaoLineSearchMonitor(TaoLineSearch ls, PetscInt its, PetscReal f, PetscReal step)
434 {
435 PetscInt tabs;
436
437 PetscFunctionBegin;
438 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
439 if (ls->usemonitor) {
440 PetscCall(PetscViewerASCIIGetTab(ls->viewer, &tabs));
441 PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel));
442 PetscCall(PetscViewerASCIIPrintf(ls->viewer, "%3" PetscInt_FMT " LS", its));
443 PetscCall(PetscViewerASCIIPrintf(ls->viewer, " Function value: %g,", (double)f));
444 PetscCall(PetscViewerASCIIPrintf(ls->viewer, " Step length: %g\n", (double)step));
445 if (ls->ops->monitor && its > 0) {
446 PetscCall(PetscViewerASCIISetTab(ls->viewer, ((PetscObject)ls)->tablevel + 3));
447 PetscUseTypeMethod(ls, monitor);
448 }
449 PetscCall(PetscViewerASCIISetTab(ls->viewer, tabs));
450 }
451 PetscFunctionReturn(PETSC_SUCCESS);
452 }
453
454 /*@
455 TaoLineSearchSetFromOptions - Sets various `TaoLineSearch` parameters from user
456 options.
457
458 Collective
459
460 Input Parameter:
461 . ls - the `TaoLineSearch` context
462
463 Options Database Keys:
464 + -tao_ls_type <type> - The algorithm that `TaoLineSearch` uses (more-thuente, gpcg, unit)
465 . -tao_ls_ftol <tol> - tolerance for sufficient decrease
466 . -tao_ls_gtol <tol> - tolerance for curvature condition
467 . -tao_ls_rtol <tol> - relative tolerance for acceptable step
468 . -tao_ls_stepinit <step> - initial steplength allowed
469 . -tao_ls_stepmin <step> - minimum steplength allowed
470 . -tao_ls_stepmax <step> - maximum steplength allowed
471 . -tao_ls_max_funcs <n> - maximum number of function evaluations allowed
472 - -tao_ls_view - display line-search results to standard output
473
474 Level: beginner
475
476 .seealso: `TaoLineSearch`
477 @*/
TaoLineSearchSetFromOptions(TaoLineSearch ls)478 PetscErrorCode TaoLineSearchSetFromOptions(TaoLineSearch ls)
479 {
480 const char *default_type = TAOLINESEARCHMT;
481 char type[256], monfilename[PETSC_MAX_PATH_LEN];
482 PetscViewer monviewer;
483 PetscBool flg;
484
485 PetscFunctionBegin;
486 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
487 PetscObjectOptionsBegin((PetscObject)ls);
488 if (((PetscObject)ls)->type_name) default_type = ((PetscObject)ls)->type_name;
489 /* Check for type from options */
490 PetscCall(PetscOptionsFList("-tao_ls_type", "Tao Line Search type", "TaoLineSearchSetType", TaoLineSearchList, default_type, type, 256, &flg));
491 if (flg) {
492 PetscCall(TaoLineSearchSetType(ls, type));
493 } else if (!((PetscObject)ls)->type_name) {
494 PetscCall(TaoLineSearchSetType(ls, default_type));
495 }
496
497 PetscCall(PetscOptionsInt("-tao_ls_max_funcs", "max function evals in line search", "", ls->max_funcs, &ls->max_funcs, NULL));
498 PetscCall(PetscOptionsReal("-tao_ls_ftol", "tol for sufficient decrease", "", ls->ftol, &ls->ftol, NULL));
499 PetscCall(PetscOptionsReal("-tao_ls_gtol", "tol for curvature condition", "", ls->gtol, &ls->gtol, NULL));
500 PetscCall(PetscOptionsReal("-tao_ls_rtol", "relative tol for acceptable step", "", ls->rtol, &ls->rtol, NULL));
501 PetscCall(PetscOptionsReal("-tao_ls_stepmin", "lower bound for step", "", ls->stepmin, &ls->stepmin, NULL));
502 PetscCall(PetscOptionsReal("-tao_ls_stepmax", "upper bound for step", "", ls->stepmax, &ls->stepmax, NULL));
503 PetscCall(PetscOptionsReal("-tao_ls_stepinit", "initial step", "", ls->initstep, &ls->initstep, NULL));
504 PetscCall(PetscOptionsString("-tao_ls_monitor", "enable the basic monitor", "TaoLineSearchSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg));
505 if (flg) {
506 PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ls), monfilename, &monviewer));
507 ls->viewer = monviewer;
508 ls->usemonitor = PETSC_TRUE;
509 }
510 PetscTryTypeMethod(ls, setfromoptions, PetscOptionsObject);
511 PetscOptionsEnd();
512 PetscFunctionReturn(PETSC_SUCCESS);
513 }
514
515 /*@
516 TaoLineSearchGetType - Gets the current line search algorithm
517
518 Not Collective
519
520 Input Parameter:
521 . ls - the `TaoLineSearch` context
522
523 Output Parameter:
524 . type - the line search algorithm in effect
525
526 Level: developer
527
528 .seealso: `TaoLineSearch`
529 @*/
TaoLineSearchGetType(TaoLineSearch ls,TaoLineSearchType * type)530 PetscErrorCode TaoLineSearchGetType(TaoLineSearch ls, TaoLineSearchType *type)
531 {
532 PetscFunctionBegin;
533 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
534 PetscAssertPointer(type, 2);
535 *type = ((PetscObject)ls)->type_name;
536 PetscFunctionReturn(PETSC_SUCCESS);
537 }
538
539 /*@
540 TaoLineSearchGetNumberFunctionEvaluations - Gets the number of function and gradient evaluation
541 routines used by the line search in last application (not cumulative).
542
543 Not Collective
544
545 Input Parameter:
546 . ls - the `TaoLineSearch` context
547
548 Output Parameters:
549 + nfeval - number of function evaluations
550 . ngeval - number of gradient evaluations
551 - nfgeval - number of function/gradient evaluations
552
553 Level: intermediate
554
555 Note:
556 If the line search is using the `Tao` objective and gradient
557 routines directly (see `TaoLineSearchUseTaoRoutines()`), then the `Tao`
558 is already counting the number of evaluations.
559
560 .seealso: `TaoLineSearch`
561 @*/
TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls,PetscInt * nfeval,PetscInt * ngeval,PetscInt * nfgeval)562 PetscErrorCode TaoLineSearchGetNumberFunctionEvaluations(TaoLineSearch ls, PetscInt *nfeval, PetscInt *ngeval, PetscInt *nfgeval)
563 {
564 PetscFunctionBegin;
565 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
566 *nfeval = ls->nfeval;
567 *ngeval = ls->ngeval;
568 *nfgeval = ls->nfgeval;
569 PetscFunctionReturn(PETSC_SUCCESS);
570 }
571
572 /*@
573 TaoLineSearchIsUsingTaoRoutines - Checks whether the line search is using
574 the standard `Tao` evaluation routines.
575
576 Not Collective
577
578 Input Parameter:
579 . ls - the `TaoLineSearch` context
580
581 Output Parameter:
582 . flg - `PETSC_TRUE` if the line search is using `Tao` evaluation routines,
583 otherwise `PETSC_FALSE`
584
585 Level: developer
586
587 .seealso: `TaoLineSearch`
588 @*/
TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls,PetscBool * flg)589 PetscErrorCode TaoLineSearchIsUsingTaoRoutines(TaoLineSearch ls, PetscBool *flg)
590 {
591 PetscFunctionBegin;
592 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
593 *flg = ls->usetaoroutines;
594 PetscFunctionReturn(PETSC_SUCCESS);
595 }
596
597 /*@C
598 TaoLineSearchSetObjectiveRoutine - Sets the function evaluation routine for the line search
599
600 Logically Collective
601
602 Input Parameters:
603 + ls - the `TaoLineSearch` context
604 . func - the objective function evaluation routine
605 - ctx - the (optional) user-defined context for private data
606
607 Calling sequence of `func`:
608 + ls - the line search context
609 . x - input vector
610 . f - function value
611 - ctx - (optional) user-defined context
612
613 Level: advanced
614
615 Notes:
616 Use this routine only if you want the line search objective
617 evaluation routine to be different from the `Tao`'s objective
618 evaluation routine. If you use this routine you must also set
619 the line search gradient and/or function/gradient routine.
620
621 Some algorithms (lcl, gpcg) set their own objective routine for the
622 line search, application programmers should be wary of overriding the
623 default objective routine.
624
625 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
626 @*/
TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls,PetscErrorCode (* func)(TaoLineSearch ls,Vec x,PetscReal * f,PetscCtx ctx),PetscCtx ctx)627 PetscErrorCode TaoLineSearchSetObjectiveRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *f, PetscCtx ctx), PetscCtx ctx)
628 {
629 PetscFunctionBegin;
630 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
631
632 ls->ops->computeobjective = func;
633 if (ctx) ls->userctx_func = ctx;
634 ls->usetaoroutines = PETSC_FALSE;
635 PetscFunctionReturn(PETSC_SUCCESS);
636 }
637
638 /*@C
639 TaoLineSearchSetGradientRoutine - Sets the gradient evaluation routine for the line search
640
641 Logically Collective
642
643 Input Parameters:
644 + ls - the `TaoLineSearch` context
645 . func - the gradient evaluation routine
646 - ctx - the (optional) user-defined context for private data
647
648 Calling sequence of `func`:
649 + ls - the linesearch object
650 . x - input vector
651 . g - gradient vector
652 - ctx - (optional) user-defined context
653
654 Level: beginner
655
656 Note:
657 Use this routine only if you want the line search gradient
658 evaluation routine to be different from the `Tao`'s gradient
659 evaluation routine. If you use this routine you must also set
660 the line search function and/or function/gradient routine.
661
662 Some algorithms (lcl, gpcg) set their own gradient routine for the
663 line search, application programmers should be wary of overriding the
664 default gradient routine.
665
666 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetObjectiveAndGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
667 @*/
TaoLineSearchSetGradientRoutine(TaoLineSearch ls,PetscErrorCode (* func)(TaoLineSearch ls,Vec x,Vec g,PetscCtx ctx),PetscCtx ctx)668 PetscErrorCode TaoLineSearchSetGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec g, PetscCtx ctx), PetscCtx ctx)
669 {
670 PetscFunctionBegin;
671 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
672 ls->ops->computegradient = func;
673 if (ctx) ls->userctx_grad = ctx;
674 ls->usetaoroutines = PETSC_FALSE;
675 PetscFunctionReturn(PETSC_SUCCESS);
676 }
677
678 /*@C
679 TaoLineSearchSetObjectiveAndGradientRoutine - Sets the objective/gradient evaluation routine for the line search
680
681 Logically Collective
682
683 Input Parameters:
684 + ls - the `TaoLineSearch` context
685 . func - the objective and gradient evaluation routine
686 - ctx - the (optional) user-defined context for private data
687
688 Calling sequence of `func`:
689 + ls - the linesearch object
690 . x - input vector
691 . f - function value
692 . g - gradient vector
693 - ctx - (optional) user-defined context
694
695 Level: beginner
696
697 Note:
698 Use this routine only if you want the line search objective and gradient
699 evaluation routines to be different from the `Tao`'s objective
700 and gradient evaluation routines.
701
702 Some algorithms (lcl, gpcg) set their own objective routine for the
703 line search, application programmers should be wary of overriding the
704 default objective routine.
705
706 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjectiveRoutine()`, `TaoLineSearchSetGradientRoutine()`, `TaoLineSearchUseTaoRoutines()`
707 @*/
TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls,PetscErrorCode (* func)(TaoLineSearch ls,Vec x,PetscReal * f,Vec g,PetscCtx ctx),PetscCtx ctx)708 PetscErrorCode TaoLineSearchSetObjectiveAndGradientRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscCtx ctx), PetscCtx ctx)
709 {
710 PetscFunctionBegin;
711 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
712 ls->ops->computeobjectiveandgradient = func;
713 if (ctx) ls->userctx_funcgrad = ctx;
714 ls->usetaoroutines = PETSC_FALSE;
715 PetscFunctionReturn(PETSC_SUCCESS);
716 }
717
718 /*@C
719 TaoLineSearchSetObjectiveAndGTSRoutine - Sets the objective and
720 (gradient'*stepdirection) evaluation routine for the line search.
721
722 Logically Collective
723
724 Input Parameters:
725 + ls - the `TaoLineSearch` context
726 . func - the objective and gradient evaluation routine
727 - ctx - the (optional) user-defined context for private data
728
729 Calling sequence of `func`:
730 + ls - the linesearch context
731 . x - input vector
732 . s - step direction
733 . f - function value
734 . gts - inner product of gradient and step direction vectors
735 - ctx - (optional) user-defined context
736
737 Level: advanced
738
739 Notes:
740 Sometimes it is more efficient to compute the inner product of the gradient and the step
741 direction than it is to compute the gradient, and this is all the line search typically needs
742 of the gradient.
743
744 The gradient will still need to be computed at the end of the line
745 search, so you will still need to set a line search gradient evaluation
746 routine
747
748 Bounded line searches (those used in bounded optimization algorithms)
749 don't use g's directly, but rather (g'x - g'x0)/steplength. You can get the
750 x0 and steplength with `TaoLineSearchGetStartingVector()` and `TaoLineSearchGetStepLength()`
751
752 Some algorithms (lcl, gpcg) set their own objective routine for the
753 line search, application programmers should be wary of overriding the
754 default objective routine.
755
756 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`, `TaoLineSearchSetObjective()`, `TaoLineSearchSetGradient()`, `TaoLineSearchUseTaoRoutines()`
757 @*/
TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls,PetscErrorCode (* func)(TaoLineSearch ls,Vec x,Vec s,PetscReal * f,PetscReal * gts,PetscCtx ctx),PetscCtx ctx)758 PetscErrorCode TaoLineSearchSetObjectiveAndGTSRoutine(TaoLineSearch ls, PetscErrorCode (*func)(TaoLineSearch ls, Vec x, Vec s, PetscReal *f, PetscReal *gts, PetscCtx ctx), PetscCtx ctx)
759 {
760 PetscFunctionBegin;
761 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
762 ls->ops->computeobjectiveandgts = func;
763 if (ctx) ls->userctx_funcgts = ctx;
764 ls->usegts = PETSC_TRUE;
765 ls->usetaoroutines = PETSC_FALSE;
766 PetscFunctionReturn(PETSC_SUCCESS);
767 }
768
769 /*@
770 TaoLineSearchUseTaoRoutines - Informs the `TaoLineSearch` to use the
771 objective and gradient evaluation routines from the given `Tao` object. The default.
772
773 Logically Collective
774
775 Input Parameters:
776 + ls - the `TaoLineSearch` context
777 - ts - the `Tao` context with defined objective/gradient evaluation routines
778
779 Level: developer
780
781 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchCreate()`
782 @*/
TaoLineSearchUseTaoRoutines(TaoLineSearch ls,Tao ts)783 PetscErrorCode TaoLineSearchUseTaoRoutines(TaoLineSearch ls, Tao ts)
784 {
785 PetscFunctionBegin;
786 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
787 PetscValidHeaderSpecific(ts, TAO_CLASSID, 2);
788 ls->tao = ts;
789 ls->usetaoroutines = PETSC_TRUE;
790 PetscFunctionReturn(PETSC_SUCCESS);
791 }
792
793 /*@
794 TaoLineSearchComputeObjective - Computes the objective function value at a given point
795
796 Collective
797
798 Input Parameters:
799 + ls - the `TaoLineSearch` context
800 - x - input vector
801
802 Output Parameter:
803 . f - Objective value at `x`
804
805 Level: developer
806
807 Note:
808 `TaoLineSearchComputeObjective()` is typically used within line searches
809 so most users would not generally call this routine themselves.
810
811 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
812 @*/
TaoLineSearchComputeObjective(TaoLineSearch ls,Vec x,PetscReal * f)813 PetscErrorCode TaoLineSearchComputeObjective(TaoLineSearch ls, Vec x, PetscReal *f)
814 {
815 Vec gdummy;
816 PetscReal gts;
817
818 PetscFunctionBegin;
819 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
820 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
821 PetscAssertPointer(f, 3);
822 PetscCheckSameComm(ls, 1, x, 2);
823 if (ls->usetaoroutines) {
824 PetscCall(TaoComputeObjective(ls->tao, x, f));
825 } else {
826 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");
827 PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
828 if (ls->ops->computeobjective) PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func));
829 else if (ls->ops->computeobjectiveandgradient) {
830 PetscCall(VecDuplicate(x, &gdummy));
831 PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgradient)(ls, x, f, gdummy, ls->userctx_funcgrad));
832 PetscCall(VecDestroy(&gdummy));
833 } else PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, >s, ls->userctx_funcgts));
834 PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
835 }
836 ls->nfeval++;
837 PetscFunctionReturn(PETSC_SUCCESS);
838 }
839
840 /*@
841 TaoLineSearchComputeObjectiveAndGradient - Computes the objective function value at a given point
842
843 Collective
844
845 Input Parameters:
846 + ls - the `TaoLineSearch` context
847 - x - input vector
848
849 Output Parameters:
850 + f - Objective value at `x`
851 - g - Gradient vector at `x`
852
853 Level: developer
854
855 Note:
856 `TaoLineSearchComputeObjectiveAndGradient()` is typically used within line searches
857 so most users would not generally call this routine themselves.
858
859 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchSetObjectiveRoutine()`
860 @*/
TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls,Vec x,PetscReal * f,Vec g)861 PetscErrorCode TaoLineSearchComputeObjectiveAndGradient(TaoLineSearch ls, Vec x, PetscReal *f, Vec g)
862 {
863 PetscFunctionBegin;
864 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
865 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
866 PetscAssertPointer(f, 3);
867 PetscValidHeaderSpecific(g, VEC_CLASSID, 4);
868 PetscCheckSameComm(ls, 1, x, 2);
869 PetscCheckSameComm(ls, 1, g, 4);
870 if (ls->usetaoroutines) {
871 PetscCall(TaoComputeObjectiveAndGradient(ls->tao, x, f, g));
872 } else {
873 PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
874 if (ls->ops->computeobjectiveandgradient) PetscCallBack("TaoLineSearch callback objective/gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, f, g, ls->userctx_funcgrad));
875 else {
876 PetscCallBack("TaoLineSearch callback objective", (*ls->ops->computeobjective)(ls, x, f, ls->userctx_func));
877 PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad));
878 }
879 PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
880 PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f)));
881 }
882 ls->nfgeval++;
883 PetscFunctionReturn(PETSC_SUCCESS);
884 }
885
886 /*@
887 TaoLineSearchComputeGradient - Computes the gradient of the objective function
888
889 Collective
890
891 Input Parameters:
892 + ls - the `TaoLineSearch` context
893 - x - input vector
894
895 Output Parameter:
896 . g - gradient vector
897
898 Level: developer
899
900 Note:
901 `TaoComputeGradient()` is typically used within line searches
902 so most users would not generally call this routine themselves.
903
904 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeObjective()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetGradient()`
905 @*/
TaoLineSearchComputeGradient(TaoLineSearch ls,Vec x,Vec g)906 PetscErrorCode TaoLineSearchComputeGradient(TaoLineSearch ls, Vec x, Vec g)
907 {
908 PetscReal fdummy;
909
910 PetscFunctionBegin;
911 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
912 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
913 PetscValidHeaderSpecific(g, VEC_CLASSID, 3);
914 PetscCheckSameComm(ls, 1, x, 2);
915 PetscCheckSameComm(ls, 1, g, 3);
916 if (ls->usetaoroutines) {
917 PetscCall(TaoComputeGradient(ls->tao, x, g));
918 } else {
919 PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
920 if (ls->ops->computegradient) PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computegradient)(ls, x, g, ls->userctx_grad));
921 else PetscCallBack("TaoLineSearch callback gradient", (*ls->ops->computeobjectiveandgradient)(ls, x, &fdummy, g, ls->userctx_funcgrad));
922 PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
923 }
924 ls->ngeval++;
925 PetscFunctionReturn(PETSC_SUCCESS);
926 }
927
928 /*@
929 TaoLineSearchComputeObjectiveAndGTS - Computes the objective function value and inner product of gradient and
930 step direction at a given point
931
932 Collective
933
934 Input Parameters:
935 + ls - the `TaoLineSearch` context
936 - x - input vector
937
938 Output Parameters:
939 + f - Objective value at `x`
940 - gts - inner product of gradient and step direction at `x`
941
942 Level: developer
943
944 Note:
945 `TaoLineSearchComputeObjectiveAndGTS()` is typically used within line searches
946 so most users would not generally call this routine themselves.
947
948 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchComputeGradient()`, `TaoLineSearchComputeObjectiveAndGradient()`, `TaoLineSearchSetObjectiveRoutine()`
949 @*/
TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls,Vec x,PetscReal * f,PetscReal * gts)950 PetscErrorCode TaoLineSearchComputeObjectiveAndGTS(TaoLineSearch ls, Vec x, PetscReal *f, PetscReal *gts)
951 {
952 PetscFunctionBegin;
953 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
954 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
955 PetscAssertPointer(f, 3);
956 PetscAssertPointer(gts, 4);
957 PetscCheckSameComm(ls, 1, x, 2);
958 PetscCall(PetscLogEventBegin(TAOLINESEARCH_Eval, ls, 0, 0, 0));
959 PetscCallBack("TaoLineSearch callback objective/gts", (*ls->ops->computeobjectiveandgts)(ls, x, ls->stepdirection, f, gts, ls->userctx_funcgts));
960 PetscCall(PetscLogEventEnd(TAOLINESEARCH_Eval, ls, 0, 0, 0));
961 PetscCall(PetscInfo(ls, "TaoLineSearch Function evaluation: %14.12e\n", (double)(*f)));
962 ls->nfeval++;
963 PetscFunctionReturn(PETSC_SUCCESS);
964 }
965
966 /*@
967 TaoLineSearchGetSolution - Returns the solution to the line search
968
969 Collective
970
971 Input Parameter:
972 . ls - the `TaoLineSearch` context
973
974 Output Parameters:
975 + x - the new solution
976 . f - the objective function value at `x`
977 . g - the gradient at `x`
978 . steplength - the multiple of the step direction taken by the line search
979 - reason - the reason why the line search terminated
980
981 Level: developer
982
983 .seealso: `TaoLineSearchGetStartingVector()`, `TaoLineSearchGetStepDirection()`
984 @*/
TaoLineSearchGetSolution(TaoLineSearch ls,Vec x,PetscReal * f,Vec g,PetscReal * steplength,TaoLineSearchConvergedReason * reason)985 PetscErrorCode TaoLineSearchGetSolution(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, PetscReal *steplength, TaoLineSearchConvergedReason *reason)
986 {
987 PetscFunctionBegin;
988 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
989 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
990 PetscAssertPointer(f, 3);
991 PetscValidHeaderSpecific(g, VEC_CLASSID, 4);
992 PetscAssertPointer(reason, 6);
993 if (ls->new_x) PetscCall(VecCopy(ls->new_x, x));
994 *f = ls->new_f;
995 if (ls->new_g) PetscCall(VecCopy(ls->new_g, g));
996 if (steplength) *steplength = ls->step;
997 *reason = ls->reason;
998 PetscFunctionReturn(PETSC_SUCCESS);
999 }
1000
1001 /*@
1002 TaoLineSearchGetStartingVector - Gets a the initial point of the line
1003 search.
1004
1005 Not Collective
1006
1007 Input Parameter:
1008 . ls - the `TaoLineSearch` context
1009
1010 Output Parameter:
1011 . x - The initial point of the line search
1012
1013 Level: advanced
1014
1015 .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStepDirection()`
1016 @*/
TaoLineSearchGetStartingVector(TaoLineSearch ls,Vec * x)1017 PetscErrorCode TaoLineSearchGetStartingVector(TaoLineSearch ls, Vec *x)
1018 {
1019 PetscFunctionBegin;
1020 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1021 if (x) *x = ls->start_x;
1022 PetscFunctionReturn(PETSC_SUCCESS);
1023 }
1024
1025 /*@
1026 TaoLineSearchGetStepDirection - Gets the step direction of the line
1027 search.
1028
1029 Not Collective
1030
1031 Input Parameter:
1032 . ls - the `TaoLineSearch` context
1033
1034 Output Parameter:
1035 . s - the step direction of the line search
1036
1037 Level: advanced
1038
1039 .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStartingVector()`
1040 @*/
TaoLineSearchGetStepDirection(TaoLineSearch ls,Vec * s)1041 PetscErrorCode TaoLineSearchGetStepDirection(TaoLineSearch ls, Vec *s)
1042 {
1043 PetscFunctionBegin;
1044 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1045 if (s) *s = ls->stepdirection;
1046 PetscFunctionReturn(PETSC_SUCCESS);
1047 }
1048
1049 /*@
1050 TaoLineSearchGetFullStepObjective - Returns the objective function value at the full step. Useful for some minimization algorithms.
1051
1052 Not Collective
1053
1054 Input Parameter:
1055 . ls - the `TaoLineSearch` context
1056
1057 Output Parameter:
1058 . f_fullstep - the objective value at the full step length
1059
1060 Level: developer
1061
1062 .seealso: `TaoLineSearchGetSolution()`, `TaoLineSearchGetStartingVector()`, `TaoLineSearchGetStepDirection()`
1063 @*/
TaoLineSearchGetFullStepObjective(TaoLineSearch ls,PetscReal * f_fullstep)1064 PetscErrorCode TaoLineSearchGetFullStepObjective(TaoLineSearch ls, PetscReal *f_fullstep)
1065 {
1066 PetscFunctionBegin;
1067 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1068 *f_fullstep = ls->f_fullstep;
1069 PetscFunctionReturn(PETSC_SUCCESS);
1070 }
1071
1072 /*@
1073 TaoLineSearchSetVariableBounds - Sets the upper and lower bounds for a bounded line search
1074
1075 Logically Collective
1076
1077 Input Parameters:
1078 + ls - the `TaoLineSearch` context
1079 . xl - vector of lower bounds
1080 - xu - vector of upper bounds
1081
1082 Level: beginner
1083
1084 Note:
1085 If the variable bounds are not set with this routine, then
1086 `PETSC_NINFINITY` and `PETSC_INFINITY` are assumed
1087
1088 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoSetVariableBounds()`, `TaoLineSearchCreate()`
1089 @*/
TaoLineSearchSetVariableBounds(TaoLineSearch ls,Vec xl,Vec xu)1090 PetscErrorCode TaoLineSearchSetVariableBounds(TaoLineSearch ls, Vec xl, Vec xu)
1091 {
1092 PetscFunctionBegin;
1093 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1094 if (xl) PetscValidHeaderSpecific(xl, VEC_CLASSID, 2);
1095 if (xu) PetscValidHeaderSpecific(xu, VEC_CLASSID, 3);
1096 PetscCall(PetscObjectReference((PetscObject)xl));
1097 PetscCall(PetscObjectReference((PetscObject)xu));
1098 PetscCall(VecDestroy(&ls->lower));
1099 PetscCall(VecDestroy(&ls->upper));
1100 ls->lower = xl;
1101 ls->upper = xu;
1102 ls->bounded = (PetscBool)(xl || xu);
1103 PetscFunctionReturn(PETSC_SUCCESS);
1104 }
1105
1106 /*@
1107 TaoLineSearchSetInitialStepLength - Sets the initial step length of a line
1108 search. If this value is not set then 1.0 is assumed.
1109
1110 Logically Collective
1111
1112 Input Parameters:
1113 + ls - the `TaoLineSearch` context
1114 - s - the initial step size
1115
1116 Level: intermediate
1117
1118 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchGetStepLength()`, `TaoLineSearchApply()`
1119 @*/
TaoLineSearchSetInitialStepLength(TaoLineSearch ls,PetscReal s)1120 PetscErrorCode TaoLineSearchSetInitialStepLength(TaoLineSearch ls, PetscReal s)
1121 {
1122 PetscFunctionBegin;
1123 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1124 PetscValidLogicalCollectiveReal(ls, s, 2);
1125 ls->initstep = s;
1126 PetscFunctionReturn(PETSC_SUCCESS);
1127 }
1128
1129 /*@
1130 TaoLineSearchGetStepLength - Get the current step length
1131
1132 Not Collective
1133
1134 Input Parameter:
1135 . ls - the `TaoLineSearch` context
1136
1137 Output Parameter:
1138 . s - the current step length
1139
1140 Level: intermediate
1141
1142 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetInitialStepLength()`, `TaoLineSearchApply()`
1143 @*/
TaoLineSearchGetStepLength(TaoLineSearch ls,PetscReal * s)1144 PetscErrorCode TaoLineSearchGetStepLength(TaoLineSearch ls, PetscReal *s)
1145 {
1146 PetscFunctionBegin;
1147 PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
1148 *s = ls->step;
1149 PetscFunctionReturn(PETSC_SUCCESS);
1150 }
1151
1152 /*@C
1153 TaoLineSearchRegister - Adds a line-search algorithm to the registry
1154
1155 Not Collective, No Fortran Support
1156
1157 Input Parameters:
1158 + sname - name of a new user-defined solver
1159 - func - routine to Create method context
1160
1161 Example Usage:
1162 .vb
1163 TaoLineSearchRegister("my_linesearch", MyLinesearchCreate);
1164 .ve
1165
1166 Then, your solver can be chosen with the procedural interface via
1167 .vb
1168 TaoLineSearchSetType(ls, "my_linesearch")
1169 .ve
1170 or at runtime via the option
1171 .vb
1172 -tao_ls_type my_linesearch
1173 .ve
1174
1175 Level: developer
1176
1177 Note:
1178 `TaoLineSearchRegister()` may be called multiple times to add several user-defined solvers.
1179
1180 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`
1181 @*/
TaoLineSearchRegister(const char sname[],PetscErrorCode (* func)(TaoLineSearch))1182 PetscErrorCode TaoLineSearchRegister(const char sname[], PetscErrorCode (*func)(TaoLineSearch))
1183 {
1184 PetscFunctionBegin;
1185 PetscCall(TaoLineSearchInitializePackage());
1186 PetscCall(PetscFunctionListAdd(&TaoLineSearchList, sname, func));
1187 PetscFunctionReturn(PETSC_SUCCESS);
1188 }
1189
1190 /*@
1191 TaoLineSearchAppendOptionsPrefix - Appends to the prefix used for searching
1192 for all `TaoLineSearch` options in the database.
1193
1194 Collective
1195
1196 Input Parameters:
1197 + ls - the `TaoLineSearch` solver context
1198 - p - the prefix string to prepend to all line search requests
1199
1200 Level: advanced
1201
1202 Notes:
1203 A hyphen (-) must NOT be given at the beginning of the prefix name.
1204 The first character of all runtime options is AUTOMATICALLY the hyphen.
1205
1206 This is inherited from the `Tao` object so rarely needs to be set
1207
1208 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()`
1209 @*/
TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls,const char p[])1210 PetscErrorCode TaoLineSearchAppendOptionsPrefix(TaoLineSearch ls, const char p[])
1211 {
1212 return PetscObjectAppendOptionsPrefix((PetscObject)ls, p);
1213 }
1214
1215 /*@
1216 TaoLineSearchGetOptionsPrefix - Gets the prefix used for searching for all
1217 `TaoLineSearch` options in the database
1218
1219 Not Collective
1220
1221 Input Parameter:
1222 . ls - the `TaoLineSearch` context
1223
1224 Output Parameter:
1225 . p - pointer to the prefix string used is returned
1226
1227 Level: advanced
1228
1229 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchSetOptionsPrefix()`, `TaoLineSearchAppendOptionsPrefix()`
1230 @*/
TaoLineSearchGetOptionsPrefix(TaoLineSearch ls,const char * p[])1231 PetscErrorCode TaoLineSearchGetOptionsPrefix(TaoLineSearch ls, const char *p[])
1232 {
1233 return PetscObjectGetOptionsPrefix((PetscObject)ls, p);
1234 }
1235
1236 /*@
1237 TaoLineSearchSetOptionsPrefix - Sets the prefix used for searching for all
1238 `TaoLineSearch` options in the database.
1239
1240 Logically Collective
1241
1242 Input Parameters:
1243 + ls - the `TaoLineSearch` context
1244 - p - the prefix string to prepend to all `ls` option requests
1245
1246 Level: advanced
1247
1248 Notes:
1249 A hyphen (-) must NOT be given at the beginning of the prefix name.
1250 The first character of all runtime options is AUTOMATICALLY the hyphen.
1251
1252 This is inherited from the `Tao` object so rarely needs to be set
1253
1254 For example, to distinguish between the runtime options for two
1255 different line searches, one could call
1256 .vb
1257 TaoLineSearchSetOptionsPrefix(ls1,"sys1_")
1258 TaoLineSearchSetOptionsPrefix(ls2,"sys2_")
1259 .ve
1260
1261 This would enable use of different options for each system, such as
1262 .vb
1263 -sys1_tao_ls_type mt
1264 -sys2_tao_ls_type armijo
1265 .ve
1266
1267 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchAppendOptionsPrefix()`, `TaoLineSearchGetOptionsPrefix()`
1268 @*/
TaoLineSearchSetOptionsPrefix(TaoLineSearch ls,const char p[])1269 PetscErrorCode TaoLineSearchSetOptionsPrefix(TaoLineSearch ls, const char p[])
1270 {
1271 return PetscObjectSetOptionsPrefix((PetscObject)ls, p);
1272 }
1273