xref: /petsc/src/tao/interface/taosolver_fg.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
2 
3 /*@
4   TaoSetSolution - Sets the vector holding the initial guess for the solve
5 
6   Logically collective on tao
7 
8   Input Parameters:
9 + tao - the Tao context
10 - x0  - the initial guess
11 
12   Level: beginner
13 .seealso: `Tao`, `TaoCreate()`, `TaoSolve()`, `TaoGetSolution()`
14 @*/
15 PetscErrorCode TaoSetSolution(Tao tao, Vec x0)
16 {
17   PetscFunctionBegin;
18   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
19   if (x0) PetscValidHeaderSpecific(x0,VEC_CLASSID,2);
20   PetscCall(PetscObjectReference((PetscObject)x0));
21   PetscCall(VecDestroy(&tao->solution));
22   tao->solution = x0;
23   PetscFunctionReturn(0);
24 }
25 
26 PetscErrorCode TaoTestGradient(Tao tao,Vec x,Vec g1)
27 {
28   Vec               g2,g3;
29   PetscBool         complete_print = PETSC_FALSE,test = PETSC_FALSE;
30   PetscReal         hcnorm,fdnorm,hcmax,fdmax,diffmax,diffnorm;
31   PetscScalar       dot;
32   MPI_Comm          comm;
33   PetscViewer       viewer,mviewer;
34   PetscViewerFormat format;
35   PetscInt          tabs;
36   static PetscBool  directionsprinted = PETSC_FALSE;
37 
38   PetscFunctionBegin;
39   PetscObjectOptionsBegin((PetscObject)tao);
40   PetscCall(PetscOptionsName("-tao_test_gradient","Compare hand-coded and finite difference Gradients","None",&test));
41   PetscCall(PetscOptionsViewer("-tao_test_gradient_view","View difference between hand-coded and finite difference Gradients element entries","None",&mviewer,&format,&complete_print));
42   PetscOptionsEnd();
43   if (!test) {
44     if (complete_print) {
45       PetscCall(PetscViewerDestroy(&mviewer));
46     }
47     PetscFunctionReturn(0);
48   }
49 
50   PetscCall(PetscObjectGetComm((PetscObject)tao,&comm));
51   PetscCall(PetscViewerASCIIGetStdout(comm,&viewer));
52   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
53   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
54   PetscCall(PetscViewerASCIIPrintf(viewer,"  ---------- Testing Gradient -------------\n"));
55   if (!complete_print && !directionsprinted) {
56     PetscCall(PetscViewerASCIIPrintf(viewer,"  Run with -tao_test_gradient_view and optionally -tao_test_gradient <threshold> to show difference\n"));
57     PetscCall(PetscViewerASCIIPrintf(viewer,"    of hand-coded and finite difference gradient entries greater than <threshold>.\n"));
58   }
59   if (!directionsprinted) {
60     PetscCall(PetscViewerASCIIPrintf(viewer,"  Testing hand-coded Gradient, if (for double precision runs) ||G - Gfd||/||G|| is\n"));
61     PetscCall(PetscViewerASCIIPrintf(viewer,"    O(1.e-8), the hand-coded Gradient is probably correct.\n"));
62     directionsprinted = PETSC_TRUE;
63   }
64   if (complete_print) PetscCall(PetscViewerPushFormat(mviewer,format));
65 
66   PetscCall(VecDuplicate(x,&g2));
67   PetscCall(VecDuplicate(x,&g3));
68 
69   /* Compute finite difference gradient, assume the gradient is already computed by TaoComputeGradient() and put into g1 */
70   PetscCall(TaoDefaultComputeGradient(tao,x,g2,NULL));
71 
72   PetscCall(VecNorm(g2,NORM_2,&fdnorm));
73   PetscCall(VecNorm(g1,NORM_2,&hcnorm));
74   PetscCall(VecNorm(g2,NORM_INFINITY,&fdmax));
75   PetscCall(VecNorm(g1,NORM_INFINITY,&hcmax));
76   PetscCall(VecDot(g1,g2,&dot));
77   PetscCall(VecCopy(g1,g3));
78   PetscCall(VecAXPY(g3,-1.0,g2));
79   PetscCall(VecNorm(g3,NORM_2,&diffnorm));
80   PetscCall(VecNorm(g3,NORM_INFINITY,&diffmax));
81   PetscCall(PetscViewerASCIIPrintf(viewer,"  ||Gfd|| %g, ||G|| = %g, angle cosine = (Gfd'G)/||Gfd||||G|| = %g\n", (double)fdnorm, (double)hcnorm, (double)(PetscRealPart(dot)/(fdnorm*hcnorm))));
82   PetscCall(PetscViewerASCIIPrintf(viewer,"  2-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n",(double)(diffnorm/PetscMax(hcnorm,fdnorm)),(double)diffnorm));
83   PetscCall(PetscViewerASCIIPrintf(viewer,"  max-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n",(double)(diffmax/PetscMax(hcmax,fdmax)),(double)diffmax));
84 
85   if (complete_print) {
86     PetscCall(PetscViewerASCIIPrintf(viewer,"  Hand-coded gradient ----------\n"));
87     PetscCall(VecView(g1,mviewer));
88     PetscCall(PetscViewerASCIIPrintf(viewer,"  Finite difference gradient ----------\n"));
89     PetscCall(VecView(g2,mviewer));
90     PetscCall(PetscViewerASCIIPrintf(viewer,"  Hand-coded minus finite-difference gradient ----------\n"));
91     PetscCall(VecView(g3,mviewer));
92   }
93   PetscCall(VecDestroy(&g2));
94   PetscCall(VecDestroy(&g3));
95 
96   if (complete_print) {
97     PetscCall(PetscViewerPopFormat(mviewer));
98     PetscCall(PetscViewerDestroy(&mviewer));
99   }
100   PetscCall(PetscViewerASCIISetTab(viewer,tabs));
101   PetscFunctionReturn(0);
102 }
103 
104 /*@
105   TaoComputeGradient - Computes the gradient of the objective function
106 
107   Collective on tao
108 
109   Input Parameters:
110 + tao - the Tao context
111 - X - input vector
112 
113   Output Parameter:
114 . G - gradient vector
115 
116   Options Database Keys:
117 +    -tao_test_gradient - compare the user provided gradient with one compute via finite differences to check for errors
118 -    -tao_test_gradient_view - display the user provided gradient, the finite difference gradient and the difference between them to help users detect the location of errors in the user provided gradient
119 
120   Note:
121     `TaoComputeGradient()` is typically used within the implementation of the optimization method,
122   so most users would not generally call this routine themselves.
123 
124   Level: developer
125 
126 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetGradient()`
127 @*/
128 PetscErrorCode TaoComputeGradient(Tao tao, Vec X, Vec G)
129 {
130   PetscReal      dummy;
131 
132   PetscFunctionBegin;
133   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
134   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
135   PetscValidHeaderSpecific(G,VEC_CLASSID,3);
136   PetscCheckSameComm(tao,1,X,2);
137   PetscCheckSameComm(tao,1,G,3);
138   PetscCall(VecLockReadPush(X));
139   if (tao->ops->computegradient) {
140     PetscCall(PetscLogEventBegin(TAO_GradientEval,tao,X,G,NULL));
141     PetscCallBack("Tao callback gradient",(*tao->ops->computegradient)(tao,X,G,tao->user_gradP));
142     PetscCall(PetscLogEventEnd(TAO_GradientEval,tao,X,G,NULL));
143     tao->ngrads++;
144   } else if (tao->ops->computeobjectiveandgradient) {
145     PetscCall(PetscLogEventBegin(TAO_ObjGradEval,tao,X,G,NULL));
146     PetscCallBack("Tao callback objective/gradient",(*tao->ops->computeobjectiveandgradient)(tao,X,&dummy,G,tao->user_objgradP));
147     PetscCall(PetscLogEventEnd(TAO_ObjGradEval,tao,X,G,NULL));
148     tao->nfuncgrads++;
149   } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetGradient() has not been called");
150   PetscCall(VecLockReadPop(X));
151 
152   PetscCall(TaoTestGradient(tao,X,G));
153   PetscFunctionReturn(0);
154 }
155 
156 /*@
157   TaoComputeObjective - Computes the objective function value at a given point
158 
159   Collective on tao
160 
161   Input Parameters:
162 + tao - the Tao context
163 - X - input vector
164 
165   Output Parameter:
166 . f - Objective value at X
167 
168   Note:
169     `TaoComputeObjective()` is typically used within the implementation of the optimization algorithm
170   so most users would not generally call this routine themselves.
171 
172   Level: developer
173 
174 .seealso: `Tao`, `TaoComputeGradient()`, `TaoComputeObjectiveAndGradient()`, `TaoSetObjective()`
175 @*/
176 PetscErrorCode TaoComputeObjective(Tao tao, Vec X, PetscReal *f)
177 {
178   Vec            temp;
179 
180   PetscFunctionBegin;
181   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
182   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
183   PetscCheckSameComm(tao,1,X,2);
184   PetscCall(VecLockReadPush(X));
185   if (tao->ops->computeobjective) {
186     PetscCall(PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL));
187     PetscCallBack("Tao callback objective",(*tao->ops->computeobjective)(tao,X,f,tao->user_objP));
188     PetscCall(PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL));
189     tao->nfuncs++;
190   } else if (tao->ops->computeobjectiveandgradient) {
191     PetscCall(PetscInfo(tao,"Duplicating variable vector in order to call func/grad routine\n"));
192     PetscCall(VecDuplicate(X,&temp));
193     PetscCall(PetscLogEventBegin(TAO_ObjGradEval,tao,X,NULL,NULL));
194     PetscCallBack("Tao callback objective/gradient",(*tao->ops->computeobjectiveandgradient)(tao,X,f,temp,tao->user_objgradP));
195     PetscCall(PetscLogEventEnd(TAO_ObjGradEval,tao,X,NULL,NULL));
196     PetscCall(VecDestroy(&temp));
197     tao->nfuncgrads++;
198   } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjective() has not been called");
199   PetscCall(PetscInfo(tao,"TAO Function evaluation: %20.19e\n",(double)(*f)));
200   PetscCall(VecLockReadPop(X));
201   PetscFunctionReturn(0);
202 }
203 
204 /*@
205   TaoComputeObjectiveAndGradient - Computes the objective function value at a given point
206 
207   Collective on tao
208 
209   Input Parameters:
210 + tao - the Tao context
211 - X - input vector
212 
213   Output Parameters:
214 + f - Objective value at X
215 - g - Gradient vector at X
216 
217   Note:
218     `TaoComputeObjectiveAndGradient()` is typically used within the implementation of the optimization algorithm,
219   so most users would not generally call this routine themselves.
220 
221   Level: developer
222 
223 .seealso: `TaoComputeGradient()`, `TaoComputeObjectiveAndGradient()`, `TaoSetObjective()`
224 @*/
225 PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G)
226 {
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
229   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
230   PetscValidHeaderSpecific(G,VEC_CLASSID,4);
231   PetscCheckSameComm(tao,1,X,2);
232   PetscCheckSameComm(tao,1,G,4);
233   PetscCall(VecLockReadPush(X));
234   if (tao->ops->computeobjectiveandgradient) {
235     PetscCall(PetscLogEventBegin(TAO_ObjGradEval,tao,X,G,NULL));
236     if (tao->ops->computegradient == TaoDefaultComputeGradient) {
237       PetscCall(TaoComputeObjective(tao,X,f));
238       PetscCall(TaoDefaultComputeGradient(tao,X,G,NULL));
239     } else {
240       PetscCallBack("Tao callback objective/gradient",(*tao->ops->computeobjectiveandgradient)(tao,X,f,G,tao->user_objgradP));
241     }
242     PetscCall(PetscLogEventEnd(TAO_ObjGradEval,tao,X,G,NULL));
243     tao->nfuncgrads++;
244   } else if (tao->ops->computeobjective && tao->ops->computegradient) {
245     PetscCall(PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL));
246     PetscCallBack("Tao callback objective",(*tao->ops->computeobjective)(tao,X,f,tao->user_objP));
247     PetscCall(PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL));
248     tao->nfuncs++;
249     PetscCall(PetscLogEventBegin(TAO_GradientEval,tao,X,G,NULL));
250     PetscCallBack("Tao callback gradient",(*tao->ops->computegradient)(tao,X,G,tao->user_gradP));
251     PetscCall(PetscLogEventEnd(TAO_GradientEval,tao,X,G,NULL));
252     tao->ngrads++;
253   } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjective() or TaoSetGradient() not set");
254   PetscCall(PetscInfo(tao,"TAO Function evaluation: %20.19e\n",(double)(*f)));
255   PetscCall(VecLockReadPop(X));
256 
257   PetscCall(TaoTestGradient(tao,X,G));
258   PetscFunctionReturn(0);
259 }
260 
261 /*@C
262   TaoSetObjective - Sets the function evaluation routine for minimization
263 
264   Logically collective on tao
265 
266   Input Parameters:
267 + tao - the Tao context
268 . func - the objective function
269 - ctx - [optional] user-defined context for private data for the function evaluation
270         routine (may be NULL)
271 
272   Calling sequence of func:
273 $      func (Tao tao, Vec x, PetscReal *f, void *ctx);
274 
275 + x - input vector
276 . f - function value
277 - ctx - [optional] user-defined function context
278 
279   Level: beginner
280 
281 .seealso: `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoGetObjective()`
282 @*/
283 PetscErrorCode TaoSetObjective(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*,void*),void *ctx)
284 {
285   PetscFunctionBegin;
286   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
287   if (ctx) tao->user_objP = ctx;
288   if (func) tao->ops->computeobjective = func;
289   PetscFunctionReturn(0);
290 }
291 
292 /*@C
293   TaoGetObjective - Gets the function evaluation routine for the function to be minimized
294 
295   Not collective
296 
297   Input Parameter:
298 . tao - the Tao context
299 
300   Output Parameters
301 + func - the objective function
302 - ctx - the user-defined context for private data for the function evaluation
303 
304   Calling sequence of func:
305 $      func (Tao tao, Vec x, PetscReal *f, void *ctx);
306 
307 + x - input vector
308 . f - function value
309 - ctx - [optional] user-defined function context
310 
311   Level: beginner
312 
313 .seealso: `Tao`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjective()`
314 @*/
315 PetscErrorCode TaoGetObjective(Tao tao, PetscErrorCode (**func)(Tao, Vec, PetscReal*,void*),void **ctx)
316 {
317   PetscFunctionBegin;
318   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
319   if (func) *func = tao->ops->computeobjective;
320   if (ctx) *ctx = tao->user_objP;
321   PetscFunctionReturn(0);
322 }
323 
324 /*@C
325   TaoSetResidualRoutine - Sets the residual evaluation routine for least-square applications
326 
327   Logically collective on tao
328 
329   Input Parameters:
330 + tao - the Tao context
331 . func - the residual evaluation routine
332 - ctx - [optional] user-defined context for private data for the function evaluation
333         routine (may be NULL)
334 
335   Calling sequence of func:
336 $      func (Tao tao, Vec x, Vec f, void *ctx);
337 
338 + x - input vector
339 . f - function value vector
340 - ctx - [optional] user-defined function context
341 
342   Level: beginner
343 
344 .seealso: `Tao`, `TaoSetObjective()`, `TaoSetJacobianRoutine()`
345 @*/
346 PetscErrorCode TaoSetResidualRoutine(Tao tao, Vec res, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
347 {
348   PetscFunctionBegin;
349   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
350   PetscValidHeaderSpecific(res,VEC_CLASSID,2);
351   PetscCall(PetscObjectReference((PetscObject)res));
352   if (tao->ls_res) {
353     PetscCall(VecDestroy(&tao->ls_res));
354   }
355   tao->ls_res = res;
356   tao->user_lsresP = ctx;
357   tao->ops->computeresidual = func;
358 
359   PetscFunctionReturn(0);
360 }
361 
362 /*@
363   TaoSetResidualWeights - Give weights for the residual values. A vector can be used if only diagonal terms are used, otherwise a matrix can be give.
364    If this function is not provided, or if sigma_v and sigma_w are both NULL, then the identity matrix will be used for weights.
365 
366   Collective on tao
367 
368   Input Parameters:
369 + tao - the Tao context
370 . sigma_v - vector of weights (diagonal terms only)
371 . n       - the number of weights (if using off-diagonal)
372 . rows    - index list of rows for sigma_w
373 . cols    - index list of columns for sigma_w
374 - vals - array of weights
375 
376   Note: Either sigma_v or sigma_w (or both) should be NULL
377 
378   Level: intermediate
379 
380 .seealso: `Tao`, `TaoSetResidualRoutine()`
381 @*/
382 PetscErrorCode TaoSetResidualWeights(Tao tao, Vec sigma_v, PetscInt n, PetscInt *rows, PetscInt *cols, PetscReal *vals)
383 {
384   PetscInt       i;
385 
386   PetscFunctionBegin;
387   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
388   if (sigma_v) PetscValidHeaderSpecific(sigma_v,VEC_CLASSID,2);
389   PetscCall(PetscObjectReference((PetscObject)sigma_v));
390   PetscCall(VecDestroy(&tao->res_weights_v));
391   tao->res_weights_v = sigma_v;
392   if (vals) {
393     PetscCall(PetscFree(tao->res_weights_rows));
394     PetscCall(PetscFree(tao->res_weights_cols));
395     PetscCall(PetscFree(tao->res_weights_w));
396     PetscCall(PetscMalloc1(n,&tao->res_weights_rows));
397     PetscCall(PetscMalloc1(n,&tao->res_weights_cols));
398     PetscCall(PetscMalloc1(n,&tao->res_weights_w));
399     tao->res_weights_n = n;
400     for (i=0;i<n;i++) {
401       tao->res_weights_rows[i] = rows[i];
402       tao->res_weights_cols[i] = cols[i];
403       tao->res_weights_w[i] = vals[i];
404     }
405   } else {
406     tao->res_weights_n = 0;
407     tao->res_weights_rows = NULL;
408     tao->res_weights_cols = NULL;
409   }
410   PetscFunctionReturn(0);
411 }
412 
413 /*@
414   TaoComputeResidual - Computes a least-squares residual vector at a given point
415 
416   Collective on tao
417 
418   Input Parameters:
419 + tao - the Tao context
420 - X - input vector
421 
422   Output Parameter:
423 . f - Objective vector at X
424 
425   Notes:
426     `TaoComputeResidual()` is typically used within the implementation of the optimization algorithm,
427   so most users would not generally call this routine themselves.
428 
429   Level: advanced
430 
431 .seealso: `Tao`, `TaoSetResidualRoutine()`
432 @*/
433 PetscErrorCode TaoComputeResidual(Tao tao, Vec X, Vec F)
434 {
435   PetscFunctionBegin;
436   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
437   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
438   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
439   PetscCheckSameComm(tao,1,X,2);
440   PetscCheckSameComm(tao,1,F,3);
441   if (tao->ops->computeresidual) {
442     PetscCall(PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL));
443     PetscCallBack("Tao callback least-squares residual",(*tao->ops->computeresidual)(tao,X,F,tao->user_lsresP));
444     PetscCall(PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL));
445     tao->nfuncs++;
446   } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetResidualRoutine() has not been called");
447   PetscCall(PetscInfo(tao,"TAO least-squares residual evaluation.\n"));
448   PetscFunctionReturn(0);
449 }
450 
451 /*@C
452   TaoSetGradient - Sets the gradient evaluation routine for the function to be optimized
453 
454   Logically collective on tao
455 
456   Input Parameters:
457 + tao - the Tao context
458 . g - [optional] the vector to internally hold the gradient computation
459 . func - the gradient function
460 - ctx - [optional] user-defined context for private data for the gradient evaluation
461         routine (may be NULL)
462 
463   Calling sequence of func:
464 $      func (Tao tao, Vec x, Vec g, void *ctx);
465 
466 + x - input vector
467 . g - gradient value (output)
468 - ctx - [optional] user-defined function context
469 
470   Level: beginner
471 
472 .seealso: `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoGetGradient()`
473 @*/
474 PetscErrorCode TaoSetGradient(Tao tao, Vec g, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
475 {
476   PetscFunctionBegin;
477   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
478   if (g) {
479     PetscValidHeaderSpecific(g,VEC_CLASSID,2);
480     PetscCheckSameComm(tao,1,g,2);
481     PetscCall(PetscObjectReference((PetscObject)g));
482     PetscCall(VecDestroy(&tao->gradient));
483     tao->gradient = g;
484   }
485   if (func) tao->ops->computegradient = func;
486   if (ctx) tao->user_gradP = ctx;
487   PetscFunctionReturn(0);
488 }
489 
490 /*@C
491   TaoGetGradient - Gets the gradient evaluation routine for the function being optimized
492 
493   Not collective
494 
495   Input Parameter:
496 . tao - the Tao context
497 
498   Output Parameters:
499 + g - the vector to internally hold the gradient computation
500 . func - the gradient function
501 - ctx - user-defined context for private data for the gradient evaluation routine
502 
503   Calling sequence of func:
504 $      func (Tao tao, Vec x, Vec g, void *ctx);
505 
506 + x - input vector
507 . g - gradient value (output)
508 - ctx - [optional] user-defined function context
509 
510   Level: beginner
511 
512 .seealso: `Tao`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetGradient()`
513 @*/
514 PetscErrorCode TaoGetGradient(Tao tao, Vec *g, PetscErrorCode (**func)(Tao, Vec, Vec, void*),void **ctx)
515 {
516   PetscFunctionBegin;
517   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
518   if (g) *g = tao->gradient;
519   if (func) *func = tao->ops->computegradient;
520   if (ctx) *ctx = tao->user_gradP;
521   PetscFunctionReturn(0);
522 }
523 
524 /*@C
525   TaoSetObjectiveAndGradient - Sets a combined objective function and gradient evaluation routine for the function to be optimized
526 
527   Logically collective on tao
528 
529   Input Parameters:
530 + tao - the Tao context
531 . g - [optional] the vector to internally hold the gradient computation
532 . func - the gradient function
533 - ctx - [optional] user-defined context for private data for the gradient evaluation
534         routine (may be NULL)
535 
536   Calling sequence of func:
537 $      func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx);
538 
539 + x - input vector
540 . f - objective value (output)
541 . g - gradient value (output)
542 - ctx - [optional] user-defined function context
543 
544   Level: beginner
545 
546   Note:
547   For some optimization methods using a combined function can be more eifficient.
548 
549 .seealso: `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetGradient()`, `TaoGetObjectiveAndGradient()`
550 @*/
551 PetscErrorCode TaoSetObjectiveAndGradient(Tao tao, Vec g, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx)
552 {
553   PetscFunctionBegin;
554   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
555   if (g) {
556     PetscValidHeaderSpecific(g,VEC_CLASSID,2);
557     PetscCheckSameComm(tao,1,g,2);
558     PetscCall(PetscObjectReference((PetscObject)g));
559     PetscCall(VecDestroy(&tao->gradient));
560     tao->gradient = g;
561   }
562   if (ctx) tao->user_objgradP = ctx;
563   if (func) tao->ops->computeobjectiveandgradient = func;
564   PetscFunctionReturn(0);
565 }
566 
567 /*@C
568   TaoGetObjectiveAndGradient - Gets the combined objective function and gradient evaluation routine for the function to be optimized
569 
570   Not collective
571 
572   Input Parameter:
573 . tao - the Tao context
574 
575   Output Parameters:
576 + g - the vector to internally hold the gradient computation
577 . func - the gradient function
578 - ctx - user-defined context for private data for the gradient evaluation routine
579 
580   Calling sequence of func:
581 $      func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx);
582 
583 + x - input vector
584 . f - objective value (output)
585 . g - gradient value (output)
586 - ctx - [optional] user-defined function context
587 
588   Level: beginner
589 
590 .seealso: `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`
591 @*/
592 PetscErrorCode TaoGetObjectiveAndGradient(Tao tao, Vec *g, PetscErrorCode (**func)(Tao, Vec, PetscReal*, Vec, void*), void **ctx)
593 {
594   PetscFunctionBegin;
595   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
596   if (g) *g = tao->gradient;
597   if (func) *func = tao->ops->computeobjectiveandgradient;
598   if (ctx) *ctx = tao->user_objgradP;
599   PetscFunctionReturn(0);
600 }
601 
602 /*@
603   TaoIsObjectiveDefined - Checks to see if the user has
604   declared an objective-only routine.  Useful for determining when
605   it is appropriate to call `TaoComputeObjective()` or
606   `TaoComputeObjectiveAndGradient()`
607 
608   Not collective
609 
610   Input Parameter:
611 . tao - the Tao context
612 
613   Output Parameter:
614 . flg - `PETSC_TRUE` if function routine is set by user, `PETSC_FALSE` otherwise
615 
616   Level: developer
617 
618 .seealso: `Tao`, `TaoSetObjective()`, `TaoIsGradientDefined()`, `TaoIsObjectiveAndGradientDefined()`
619 @*/
620 PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg)
621 {
622   PetscFunctionBegin;
623   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
624   if (tao->ops->computeobjective == NULL) *flg = PETSC_FALSE;
625   else *flg = PETSC_TRUE;
626   PetscFunctionReturn(0);
627 }
628 
629 /*@
630   TaoIsGradientDefined - Checks to see if the user has
631   declared an objective-only routine.  Useful for determining when
632   it is appropriate to call `TaoComputeGradient()` or
633   `TaoComputeGradientAndGradient()`
634 
635   Not Collective
636 
637   Input Parameter:
638 . tao - the Tao context
639 
640   Output Parameter:
641 . flg - `PETSC_TRUE` if function routine is set by user, `PETSC_FALSE` otherwise
642 
643   Level: developer
644 
645 .seealso: `TaoSetGradient()`, `TaoIsObjectiveDefined()`, `TaoIsObjectiveAndGradientDefined()`
646 @*/
647 PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg)
648 {
649   PetscFunctionBegin;
650   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
651   if (tao->ops->computegradient == NULL) *flg = PETSC_FALSE;
652   else *flg = PETSC_TRUE;
653   PetscFunctionReturn(0);
654 }
655 
656 /*@
657   TaoIsObjectiveAndGradientDefined - Checks to see if the user has
658   declared a joint objective/gradient routine.  Useful for determining when
659   it is appropriate to call `TaoComputeObjective()` or
660   `TaoComputeObjectiveAndGradient()`
661 
662   Not Collective
663 
664   Input Parameter:
665 . tao - the Tao context
666 
667   Output Parameter:
668 . flg - `PETSC_TRUE` if function routine is set by user, `PETSC_FALSE` otherwise
669 
670   Level: developer
671 
672 .seealso: `TaoSetObjectiveAndGradient()`, `TaoIsObjectiveDefined()`, `TaoIsGradientDefined()`
673 @*/
674 PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg)
675 {
676   PetscFunctionBegin;
677   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
678   if (tao->ops->computeobjectiveandgradient == NULL) *flg = PETSC_FALSE;
679   else *flg = PETSC_TRUE;
680   PetscFunctionReturn(0);
681 }
682