1af0996ceSBarry Smith #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
2a7e14dcfSSatish Balay
3a7e14dcfSSatish Balay /*@
4a82e8c82SStefano Zampini TaoSetSolution - Sets the vector holding the initial guess for the solve
5a7e14dcfSSatish Balay
620f4b53cSBarry Smith Logically Collective
7a7e14dcfSSatish Balay
8a7e14dcfSSatish Balay Input Parameters:
947450a7bSBarry Smith + tao - the `Tao` context
10a7e14dcfSSatish Balay - x0 - the initial guess
11a7e14dcfSSatish Balay
12a7e14dcfSSatish Balay Level: beginner
1347450a7bSBarry Smith
141cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`, `TaoGetSolution()`
15a7e14dcfSSatish Balay @*/
TaoSetSolution(Tao tao,Vec x0)16d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetSolution(Tao tao, Vec x0)
17d71ae5a4SJacob Faibussowitsch {
18a7e14dcfSSatish Balay PetscFunctionBegin;
19441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
20a82e8c82SStefano Zampini if (x0) PetscValidHeaderSpecific(x0, VEC_CLASSID, 2);
219566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)x0));
229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->solution));
23a7e14dcfSSatish Balay tao->solution = x0;
243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
25a7e14dcfSSatish Balay }
26a7e14dcfSSatish Balay
TaoTestGradient(Tao tao,Vec x,Vec g1)27d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoTestGradient(Tao tao, Vec x, Vec g1)
28d71ae5a4SJacob Faibussowitsch {
29412cdd55SHong Zhang Vec g2, g3;
3009baa881SHong Zhang PetscBool complete_print = PETSC_FALSE, test = PETSC_FALSE;
3109baa881SHong Zhang PetscReal hcnorm, fdnorm, hcmax, fdmax, diffmax, diffnorm;
3209baa881SHong Zhang PetscScalar dot;
3309baa881SHong Zhang MPI_Comm comm;
34913eda9aSHong Zhang PetscViewer viewer, mviewer;
35913eda9aSHong Zhang PetscViewerFormat format;
3609baa881SHong Zhang PetscInt tabs;
3709baa881SHong Zhang static PetscBool directionsprinted = PETSC_FALSE;
3809baa881SHong Zhang
3909baa881SHong Zhang PetscFunctionBegin;
40d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)tao);
419566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-tao_test_gradient", "Compare hand-coded and finite difference Gradients", "None", &test));
429566063dSJacob Faibussowitsch PetscCall(PetscOptionsViewer("-tao_test_gradient_view", "View difference between hand-coded and finite difference Gradients element entries", "None", &mviewer, &format, &complete_print));
43d0609cedSBarry Smith PetscOptionsEnd();
442f4b6201SAlp Dener if (!test) {
4548a46eb9SPierre Jolivet if (complete_print) PetscCall(PetscViewerDestroy(&mviewer));
463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
472f4b6201SAlp Dener }
4809baa881SHong Zhang
499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(comm, &viewer));
519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------- Testing Gradient -------------\n"));
5409baa881SHong Zhang if (!complete_print && !directionsprinted) {
559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Run with -tao_test_gradient_view and optionally -tao_test_gradient <threshold> to show difference\n"));
569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " of hand-coded and finite difference gradient entries greater than <threshold>.\n"));
5709baa881SHong Zhang }
5809baa881SHong Zhang if (!directionsprinted) {
599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Testing hand-coded Gradient, if (for double precision runs) ||G - Gfd||/||G|| is\n"));
609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " O(1.e-8), the hand-coded Gradient is probably correct.\n"));
6109baa881SHong Zhang directionsprinted = PETSC_TRUE;
6209baa881SHong Zhang }
631baa6e33SBarry Smith if (complete_print) PetscCall(PetscViewerPushFormat(mviewer, format));
6409baa881SHong Zhang
659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &g2));
669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &g3));
6709baa881SHong Zhang
6809baa881SHong Zhang /* Compute finite difference gradient, assume the gradient is already computed by TaoComputeGradient() and put into g1 */
699566063dSJacob Faibussowitsch PetscCall(TaoDefaultComputeGradient(tao, x, g2, NULL));
7009baa881SHong Zhang
719566063dSJacob Faibussowitsch PetscCall(VecNorm(g2, NORM_2, &fdnorm));
729566063dSJacob Faibussowitsch PetscCall(VecNorm(g1, NORM_2, &hcnorm));
739566063dSJacob Faibussowitsch PetscCall(VecNorm(g2, NORM_INFINITY, &fdmax));
749566063dSJacob Faibussowitsch PetscCall(VecNorm(g1, NORM_INFINITY, &hcmax));
759566063dSJacob Faibussowitsch PetscCall(VecDot(g1, g2, &dot));
769566063dSJacob Faibussowitsch PetscCall(VecCopy(g1, g3));
779566063dSJacob Faibussowitsch PetscCall(VecAXPY(g3, -1.0, g2));
789566063dSJacob Faibussowitsch PetscCall(VecNorm(g3, NORM_2, &diffnorm));
799566063dSJacob Faibussowitsch PetscCall(VecNorm(g3, NORM_INFINITY, &diffmax));
809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ||Gfd|| %g, ||G|| = %g, angle cosine = (Gfd'G)/||Gfd||||G|| = %g\n", (double)fdnorm, (double)hcnorm, (double)(PetscRealPart(dot) / (fdnorm * hcnorm))));
819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " 2-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n", (double)(diffnorm / PetscMax(hcnorm, fdnorm)), (double)diffnorm));
829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " max-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n", (double)(diffmax / PetscMax(hcmax, fdmax)), (double)diffmax));
8309baa881SHong Zhang
8409baa881SHong Zhang if (complete_print) {
859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Hand-coded gradient ----------\n"));
869566063dSJacob Faibussowitsch PetscCall(VecView(g1, mviewer));
879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Finite difference gradient ----------\n"));
889566063dSJacob Faibussowitsch PetscCall(VecView(g2, mviewer));
899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Hand-coded minus finite-difference gradient ----------\n"));
909566063dSJacob Faibussowitsch PetscCall(VecView(g3, mviewer));
9109baa881SHong Zhang }
929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g2));
939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g3));
94913eda9aSHong Zhang
95913eda9aSHong Zhang if (complete_print) {
969566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(mviewer));
979566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&mviewer));
98913eda9aSHong Zhang }
999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10109baa881SHong Zhang }
10209baa881SHong Zhang
103a7e14dcfSSatish Balay /*@
104a7e14dcfSSatish Balay TaoComputeGradient - Computes the gradient of the objective function
105a7e14dcfSSatish Balay
106c3339decSBarry Smith Collective
107a7e14dcfSSatish Balay
108a7e14dcfSSatish Balay Input Parameters:
10947450a7bSBarry Smith + tao - the `Tao` context
110a7e14dcfSSatish Balay - X - input vector
111a7e14dcfSSatish Balay
112a7e14dcfSSatish Balay Output Parameter:
113a7e14dcfSSatish Balay . G - gradient vector
114a7e14dcfSSatish Balay
11509baa881SHong Zhang Options Database Keys:
11609baa881SHong Zhang + -tao_test_gradient - compare the user provided gradient with one compute via finite differences to check for errors
117dfe02fe6SHong Zhang - -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
11809baa881SHong Zhang
11947450a7bSBarry Smith Level: developer
12047450a7bSBarry Smith
12165ba42b6SBarry Smith Note:
12265ba42b6SBarry Smith `TaoComputeGradient()` is typically used within the implementation of the optimization method,
123a7e14dcfSSatish Balay so most users would not generally call this routine themselves.
124a7e14dcfSSatish Balay
1251cc06b55SBarry Smith .seealso: [](ch_tao), `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetGradient()`
126a7e14dcfSSatish Balay @*/
TaoComputeGradient(Tao tao,Vec X,Vec G)127d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeGradient(Tao tao, Vec X, Vec G)
128d71ae5a4SJacob Faibussowitsch {
129a7e14dcfSSatish Balay PetscReal dummy;
13087f595a5SBarry Smith
131a7e14dcfSSatish Balay PetscFunctionBegin;
132441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
133a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
134064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(G, VEC_CLASSID, 3);
135a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
136a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, G, 3);
1379566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
138a7e14dcfSSatish Balay if (tao->ops->computegradient) {
1399566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_GradientEval, tao, X, G, NULL));
140792fecdfSBarry Smith PetscCallBack("Tao callback gradient", (*tao->ops->computegradient)(tao, X, G, tao->user_gradP));
1419566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_GradientEval, tao, X, G, NULL));
142a7e14dcfSSatish Balay tao->ngrads++;
143a7e14dcfSSatish Balay } else if (tao->ops->computeobjectiveandgradient) {
1449566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_ObjGradEval, tao, X, G, NULL));
145792fecdfSBarry Smith PetscCallBack("Tao callback objective/gradient", (*tao->ops->computeobjectiveandgradient)(tao, X, &dummy, G, tao->user_objgradP));
1469566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_ObjGradEval, tao, X, G, NULL));
147a7e14dcfSSatish Balay tao->nfuncgrads++;
148a82e8c82SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TaoSetGradient() has not been called");
1499566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
15009baa881SHong Zhang
1519566063dSJacob Faibussowitsch PetscCall(TaoTestGradient(tao, X, G));
1523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
153a7e14dcfSSatish Balay }
154a7e14dcfSSatish Balay
155a7e14dcfSSatish Balay /*@
156a7e14dcfSSatish Balay TaoComputeObjective - Computes the objective function value at a given point
157a7e14dcfSSatish Balay
158c3339decSBarry Smith Collective
159a7e14dcfSSatish Balay
160a7e14dcfSSatish Balay Input Parameters:
16147450a7bSBarry Smith + tao - the `Tao` context
162a7e14dcfSSatish Balay - X - input vector
163a7e14dcfSSatish Balay
164a7e14dcfSSatish Balay Output Parameter:
165a7e14dcfSSatish Balay . f - Objective value at X
166a7e14dcfSSatish Balay
16747450a7bSBarry Smith Level: developer
16847450a7bSBarry Smith
16965ba42b6SBarry Smith Note:
17065ba42b6SBarry Smith `TaoComputeObjective()` is typically used within the implementation of the optimization algorithm
171a7e14dcfSSatish Balay so most users would not generally call this routine themselves.
172a7e14dcfSSatish Balay
1731cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoComputeGradient()`, `TaoComputeObjectiveAndGradient()`, `TaoSetObjective()`
174a7e14dcfSSatish Balay @*/
TaoComputeObjective(Tao tao,Vec X,PetscReal * f)175d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeObjective(Tao tao, Vec X, PetscReal *f)
176d71ae5a4SJacob Faibussowitsch {
177a7e14dcfSSatish Balay Vec temp;
17887f595a5SBarry Smith
179a7e14dcfSSatish Balay PetscFunctionBegin;
180441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
181a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
182a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
1839566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
184a7e14dcfSSatish Balay if (tao->ops->computeobjective) {
1859566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_ObjectiveEval, tao, X, NULL, NULL));
186792fecdfSBarry Smith PetscCallBack("Tao callback objective", (*tao->ops->computeobjective)(tao, X, f, tao->user_objP));
1879566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_ObjectiveEval, tao, X, NULL, NULL));
188a7e14dcfSSatish Balay tao->nfuncs++;
189a7e14dcfSSatish Balay } else if (tao->ops->computeobjectiveandgradient) {
1909566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Duplicating variable vector in order to call func/grad routine\n"));
1919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &temp));
1929566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_ObjGradEval, tao, X, NULL, NULL));
193792fecdfSBarry Smith PetscCallBack("Tao callback objective/gradient", (*tao->ops->computeobjectiveandgradient)(tao, X, f, temp, tao->user_objgradP));
1949566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_ObjGradEval, tao, X, NULL, NULL));
1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&temp));
196a7e14dcfSSatish Balay tao->nfuncgrads++;
197a82e8c82SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TaoSetObjective() has not been called");
1989566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "TAO Function evaluation: %20.19e\n", (double)(*f)));
1999566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
201a7e14dcfSSatish Balay }
202a7e14dcfSSatish Balay
203a7e14dcfSSatish Balay /*@
204a7e14dcfSSatish Balay TaoComputeObjectiveAndGradient - Computes the objective function value at a given point
205a7e14dcfSSatish Balay
206c3339decSBarry Smith Collective
207a7e14dcfSSatish Balay
208a7e14dcfSSatish Balay Input Parameters:
20947450a7bSBarry Smith + tao - the `Tao` context
210a7e14dcfSSatish Balay - X - input vector
211a7e14dcfSSatish Balay
212d8d19677SJose E. Roman Output Parameters:
21320f4b53cSBarry Smith + f - Objective value at `X`
214e056e8ceSJacob Faibussowitsch - G - Gradient vector at `X`
215a7e14dcfSSatish Balay
21647450a7bSBarry Smith Level: developer
21747450a7bSBarry Smith
21865ba42b6SBarry Smith Note:
21965ba42b6SBarry Smith `TaoComputeObjectiveAndGradient()` is typically used within the implementation of the optimization algorithm,
220a7e14dcfSSatish Balay so most users would not generally call this routine themselves.
221a7e14dcfSSatish Balay
22242747ad1SJacob Faibussowitsch .seealso: [](ch_tao), `TaoComputeGradient()`, `TaoSetObjective()`
223a7e14dcfSSatish Balay @*/
TaoComputeObjectiveAndGradient(Tao tao,Vec X,PetscReal * f,Vec G)224d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G)
225d71ae5a4SJacob Faibussowitsch {
226a7e14dcfSSatish Balay PetscFunctionBegin;
227441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
228a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
229a7e14dcfSSatish Balay PetscValidHeaderSpecific(G, VEC_CLASSID, 4);
230a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
231a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, G, 4);
2329566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
233a7e14dcfSSatish Balay if (tao->ops->computeobjectiveandgradient) {
2349566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_ObjGradEval, tao, X, G, NULL));
235f4c1ad5cSStefano Zampini if (tao->ops->computegradient == TaoDefaultComputeGradient) {
2369566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, X, f));
2379566063dSJacob Faibussowitsch PetscCall(TaoDefaultComputeGradient(tao, X, G, NULL));
238794dad8cSStefano Zampini } else PetscCallBack("Tao callback objective/gradient", (*tao->ops->computeobjectiveandgradient)(tao, X, f, G, tao->user_objgradP));
2399566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_ObjGradEval, tao, X, G, NULL));
240a7e14dcfSSatish Balay tao->nfuncgrads++;
241a7e14dcfSSatish Balay } else if (tao->ops->computeobjective && tao->ops->computegradient) {
2429566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_ObjectiveEval, tao, X, NULL, NULL));
243792fecdfSBarry Smith PetscCallBack("Tao callback objective", (*tao->ops->computeobjective)(tao, X, f, tao->user_objP));
2449566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_ObjectiveEval, tao, X, NULL, NULL));
245a7e14dcfSSatish Balay tao->nfuncs++;
2469566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_GradientEval, tao, X, G, NULL));
247792fecdfSBarry Smith PetscCallBack("Tao callback gradient", (*tao->ops->computegradient)(tao, X, G, tao->user_gradP));
2489566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_GradientEval, tao, X, G, NULL));
249a7e14dcfSSatish Balay tao->ngrads++;
250a82e8c82SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TaoSetObjective() or TaoSetGradient() not set");
2519566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "TAO Function evaluation: %20.19e\n", (double)(*f)));
2529566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
2531657496cSHong Zhang
2549566063dSJacob Faibussowitsch PetscCall(TaoTestGradient(tao, X, G));
2553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
256a7e14dcfSSatish Balay }
257a7e14dcfSSatish Balay
258a7e14dcfSSatish Balay /*@C
259a82e8c82SStefano Zampini TaoSetObjective - Sets the function evaluation routine for minimization
260a7e14dcfSSatish Balay
26120f4b53cSBarry Smith Logically Collective
262a7e14dcfSSatish Balay
263d8d19677SJose E. Roman Input Parameters:
26447450a7bSBarry Smith + tao - the `Tao` context
265a7e14dcfSSatish Balay . func - the objective function
266a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the function evaluation
26747450a7bSBarry Smith routine (may be `NULL`)
268a7e14dcfSSatish Balay
26920f4b53cSBarry Smith Calling sequence of `func`:
27020f4b53cSBarry Smith + tao - the optimizer
27120f4b53cSBarry Smith . x - input vector
272a7e14dcfSSatish Balay . f - function value
273a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
274a7e14dcfSSatish Balay
275a7e14dcfSSatish Balay Level: beginner
276a7e14dcfSSatish Balay
2771cc06b55SBarry Smith .seealso: [](ch_tao), `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoGetObjective()`
278a7e14dcfSSatish Balay @*/
TaoSetObjective(Tao tao,PetscErrorCode (* func)(Tao tao,Vec x,PetscReal * f,PetscCtx ctx),PetscCtx ctx)279*2a8381b2SBarry Smith PetscErrorCode TaoSetObjective(Tao tao, PetscErrorCode (*func)(Tao tao, Vec x, PetscReal *f, PetscCtx ctx), PetscCtx ctx)
280d71ae5a4SJacob Faibussowitsch {
281a7e14dcfSSatish Balay PetscFunctionBegin;
282441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
283a82e8c82SStefano Zampini if (ctx) tao->user_objP = ctx;
284a82e8c82SStefano Zampini if (func) tao->ops->computeobjective = func;
2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
286a82e8c82SStefano Zampini }
287a82e8c82SStefano Zampini
288a82e8c82SStefano Zampini /*@C
28965ba42b6SBarry Smith TaoGetObjective - Gets the function evaluation routine for the function to be minimized
290a82e8c82SStefano Zampini
29120f4b53cSBarry Smith Not Collective
292a82e8c82SStefano Zampini
293a82e8c82SStefano Zampini Input Parameter:
29447450a7bSBarry Smith . tao - the `Tao` context
295a82e8c82SStefano Zampini
2962fe279fdSBarry Smith Output Parameters:
297a82e8c82SStefano Zampini + func - the objective function
298a82e8c82SStefano Zampini - ctx - the user-defined context for private data for the function evaluation
299a82e8c82SStefano Zampini
30020f4b53cSBarry Smith Calling sequence of `func`:
30120f4b53cSBarry Smith + tao - the optimizer
30220f4b53cSBarry Smith . x - input vector
303a82e8c82SStefano Zampini . f - function value
304a82e8c82SStefano Zampini - ctx - [optional] user-defined function context
305a82e8c82SStefano Zampini
306a82e8c82SStefano Zampini Level: beginner
307a82e8c82SStefano Zampini
3081cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjective()`
309a82e8c82SStefano Zampini @*/
TaoGetObjective(Tao tao,PetscErrorCode (** func)(Tao tao,Vec x,PetscReal * f,PetscCtx ctx),PetscCtxRt ctx)310*2a8381b2SBarry Smith PetscErrorCode TaoGetObjective(Tao tao, PetscErrorCode (**func)(Tao tao, Vec x, PetscReal *f, PetscCtx ctx), PetscCtxRt ctx)
311d71ae5a4SJacob Faibussowitsch {
312a82e8c82SStefano Zampini PetscFunctionBegin;
313a82e8c82SStefano Zampini PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
314a82e8c82SStefano Zampini if (func) *func = tao->ops->computeobjective;
315*2a8381b2SBarry Smith if (ctx) *(void **)ctx = tao->user_objP;
3163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
317a7e14dcfSSatish Balay }
318a7e14dcfSSatish Balay
319a7e14dcfSSatish Balay /*@C
3204a48860cSAlp Dener TaoSetResidualRoutine - Sets the residual evaluation routine for least-square applications
321a7e14dcfSSatish Balay
32220f4b53cSBarry Smith Logically Collective
323a7e14dcfSSatish Balay
324d8d19677SJose E. Roman Input Parameters:
32547450a7bSBarry Smith + tao - the `Tao` context
3263b242c63SJacob Faibussowitsch . res - the residual vector
3274a48860cSAlp Dener . func - the residual evaluation routine
328a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the function evaluation
32947450a7bSBarry Smith routine (may be `NULL`)
330a7e14dcfSSatish Balay
33120f4b53cSBarry Smith Calling sequence of `func`:
33220f4b53cSBarry Smith + tao - the optimizer
33320f4b53cSBarry Smith . x - input vector
334e056e8ceSJacob Faibussowitsch . res - function value vector
335a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
336a7e14dcfSSatish Balay
337a7e14dcfSSatish Balay Level: beginner
338a7e14dcfSSatish Balay
3391cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetObjective()`, `TaoSetJacobianRoutine()`
340a7e14dcfSSatish Balay @*/
TaoSetResidualRoutine(Tao tao,Vec res,PetscErrorCode (* func)(Tao tao,Vec x,Vec res,PetscCtx ctx),PetscCtx ctx)341*2a8381b2SBarry Smith PetscErrorCode TaoSetResidualRoutine(Tao tao, Vec res, PetscErrorCode (*func)(Tao tao, Vec x, Vec res, PetscCtx ctx), PetscCtx ctx)
342d71ae5a4SJacob Faibussowitsch {
343a7e14dcfSSatish Balay PetscFunctionBegin;
344441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
3454a48860cSAlp Dener PetscValidHeaderSpecific(res, VEC_CLASSID, 2);
3469566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)res));
34748a46eb9SPierre Jolivet if (tao->ls_res) PetscCall(VecDestroy(&tao->ls_res));
3484a48860cSAlp Dener tao->ls_res = res;
3494ffbe8acSAlp Dener tao->user_lsresP = ctx;
3504a48860cSAlp Dener tao->ops->computeresidual = func;
3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
352737f463aSAlp Dener }
353737f463aSAlp Dener
354737f463aSAlp Dener /*@
35565ba42b6SBarry Smith TaoSetResidualWeights - Give weights for the residual values. A vector can be used if only diagonal terms are used, otherwise a matrix can be give.
356737f463aSAlp Dener
357c3339decSBarry Smith Collective
358737f463aSAlp Dener
359737f463aSAlp Dener Input Parameters:
36047450a7bSBarry Smith + tao - the `Tao` context
361737f463aSAlp Dener . sigma_v - vector of weights (diagonal terms only)
362737f463aSAlp Dener . n - the number of weights (if using off-diagonal)
36347450a7bSBarry Smith . rows - index list of rows for `sigma_v`
36447450a7bSBarry Smith . cols - index list of columns for `sigma_v`
365737f463aSAlp Dener - vals - array of weights
366737f463aSAlp Dener
367737f463aSAlp Dener Level: intermediate
368737f463aSAlp Dener
3693b242c63SJacob Faibussowitsch Notes:
3703b242c63SJacob Faibussowitsch If this function is not provided, or if `sigma_v` and `vals` are both `NULL`, then the
3713b242c63SJacob Faibussowitsch identity matrix will be used for weights.
3723b242c63SJacob Faibussowitsch
37347450a7bSBarry Smith Either `sigma_v` or `vals` should be `NULL`
37447450a7bSBarry Smith
3751cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetResidualRoutine()`
376737f463aSAlp Dener @*/
TaoSetResidualWeights(Tao tao,Vec sigma_v,PetscInt n,PetscInt * rows,PetscInt * cols,PetscReal * vals)377d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetResidualWeights(Tao tao, Vec sigma_v, PetscInt n, PetscInt *rows, PetscInt *cols, PetscReal *vals)
378d71ae5a4SJacob Faibussowitsch {
379737f463aSAlp Dener PetscInt i;
380a82e8c82SStefano Zampini
381737f463aSAlp Dener PetscFunctionBegin;
382737f463aSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
383a82e8c82SStefano Zampini if (sigma_v) PetscValidHeaderSpecific(sigma_v, VEC_CLASSID, 2);
3849566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sigma_v));
3859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->res_weights_v));
3864ffbe8acSAlp Dener tao->res_weights_v = sigma_v;
387737f463aSAlp Dener if (vals) {
3889566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->res_weights_rows));
3899566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->res_weights_cols));
3909566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->res_weights_w));
3919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tao->res_weights_rows));
3929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tao->res_weights_cols));
3939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tao->res_weights_w));
394737f463aSAlp Dener tao->res_weights_n = n;
395737f463aSAlp Dener for (i = 0; i < n; i++) {
396737f463aSAlp Dener tao->res_weights_rows[i] = rows[i];
397737f463aSAlp Dener tao->res_weights_cols[i] = cols[i];
398737f463aSAlp Dener tao->res_weights_w[i] = vals[i];
399737f463aSAlp Dener }
400737f463aSAlp Dener } else {
401737f463aSAlp Dener tao->res_weights_n = 0;
40283c8fe1dSLisandro Dalcin tao->res_weights_rows = NULL;
40383c8fe1dSLisandro Dalcin tao->res_weights_cols = NULL;
404737f463aSAlp Dener }
4053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
406a7e14dcfSSatish Balay }
407a7e14dcfSSatish Balay
4088b7a9b22SJason Sarich /*@
4094a48860cSAlp Dener TaoComputeResidual - Computes a least-squares residual vector at a given point
410a7e14dcfSSatish Balay
411c3339decSBarry Smith Collective
412a7e14dcfSSatish Balay
413a7e14dcfSSatish Balay Input Parameters:
41447450a7bSBarry Smith + tao - the `Tao` context
415a7e14dcfSSatish Balay - X - input vector
416a7e14dcfSSatish Balay
417a7e14dcfSSatish Balay Output Parameter:
418e056e8ceSJacob Faibussowitsch . F - Objective vector at `X`
419a7e14dcfSSatish Balay
42047450a7bSBarry Smith Level: advanced
42147450a7bSBarry Smith
42295452b02SPatrick Sanan Notes:
42365ba42b6SBarry Smith `TaoComputeResidual()` is typically used within the implementation of the optimization algorithm,
424a7e14dcfSSatish Balay so most users would not generally call this routine themselves.
425a7e14dcfSSatish Balay
4261cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetResidualRoutine()`
427a7e14dcfSSatish Balay @*/
TaoComputeResidual(Tao tao,Vec X,Vec F)428d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeResidual(Tao tao, Vec X, Vec F)
429d71ae5a4SJacob Faibussowitsch {
430a7e14dcfSSatish Balay PetscFunctionBegin;
431441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
432a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
433a7e14dcfSSatish Balay PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
434a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
435a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, F, 3);
436966bd95aSPierre Jolivet PetscCheck(tao->ops->computeresidual, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TaoSetResidualRoutine() has not been called");
4379566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_ObjectiveEval, tao, X, NULL, NULL));
438792fecdfSBarry Smith PetscCallBack("Tao callback least-squares residual", (*tao->ops->computeresidual)(tao, X, F, tao->user_lsresP));
4399566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_ObjectiveEval, tao, X, NULL, NULL));
440a7e14dcfSSatish Balay tao->nfuncs++;
4419566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "TAO least-squares residual evaluation.\n"));
4423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
443a7e14dcfSSatish Balay }
444a7e14dcfSSatish Balay
445a7e14dcfSSatish Balay /*@C
44665ba42b6SBarry Smith TaoSetGradient - Sets the gradient evaluation routine for the function to be optimized
447a7e14dcfSSatish Balay
44820f4b53cSBarry Smith Logically Collective
449a7e14dcfSSatish Balay
450d8d19677SJose E. Roman Input Parameters:
45147450a7bSBarry Smith + tao - the `Tao` context
452a82e8c82SStefano Zampini . g - [optional] the vector to internally hold the gradient computation
453a7e14dcfSSatish Balay . func - the gradient function
454a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the gradient evaluation
45547450a7bSBarry Smith routine (may be `NULL`)
456a7e14dcfSSatish Balay
45720f4b53cSBarry Smith Calling sequence of `func`:
45820f4b53cSBarry Smith + tao - the optimization solver
45920f4b53cSBarry Smith . x - input vector
460a7e14dcfSSatish Balay . g - gradient value (output)
461a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
462a7e14dcfSSatish Balay
463a7e14dcfSSatish Balay Level: beginner
464a7e14dcfSSatish Balay
4651cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoGetGradient()`
466a7e14dcfSSatish Balay @*/
TaoSetGradient(Tao tao,Vec g,PetscErrorCode (* func)(Tao tao,Vec x,Vec g,PetscCtx ctx),PetscCtx ctx)467*2a8381b2SBarry Smith PetscErrorCode TaoSetGradient(Tao tao, Vec g, PetscErrorCode (*func)(Tao tao, Vec x, Vec g, PetscCtx ctx), PetscCtx ctx)
468d71ae5a4SJacob Faibussowitsch {
469a7e14dcfSSatish Balay PetscFunctionBegin;
470441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
471a82e8c82SStefano Zampini if (g) {
472a82e8c82SStefano Zampini PetscValidHeaderSpecific(g, VEC_CLASSID, 2);
473a82e8c82SStefano Zampini PetscCheckSameComm(tao, 1, g, 2);
4749566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)g));
4759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->gradient));
476a82e8c82SStefano Zampini tao->gradient = g;
477a82e8c82SStefano Zampini }
478a82e8c82SStefano Zampini if (func) tao->ops->computegradient = func;
479a82e8c82SStefano Zampini if (ctx) tao->user_gradP = ctx;
4803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
481a7e14dcfSSatish Balay }
482a7e14dcfSSatish Balay
483a7e14dcfSSatish Balay /*@C
48465ba42b6SBarry Smith TaoGetGradient - Gets the gradient evaluation routine for the function being optimized
485a82e8c82SStefano Zampini
48620f4b53cSBarry Smith Not Collective
487a82e8c82SStefano Zampini
488a82e8c82SStefano Zampini Input Parameter:
48947450a7bSBarry Smith . tao - the `Tao` context
490a82e8c82SStefano Zampini
491a82e8c82SStefano Zampini Output Parameters:
492a82e8c82SStefano Zampini + g - the vector to internally hold the gradient computation
493a82e8c82SStefano Zampini . func - the gradient function
494a82e8c82SStefano Zampini - ctx - user-defined context for private data for the gradient evaluation routine
495a82e8c82SStefano Zampini
49620f4b53cSBarry Smith Calling sequence of `func`:
49720f4b53cSBarry Smith + tao - the optimizer
49820f4b53cSBarry Smith . x - input vector
499a82e8c82SStefano Zampini . g - gradient value (output)
500a82e8c82SStefano Zampini - ctx - [optional] user-defined function context
501a82e8c82SStefano Zampini
502a82e8c82SStefano Zampini Level: beginner
503a82e8c82SStefano Zampini
5041cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetGradient()`
505a82e8c82SStefano Zampini @*/
TaoGetGradient(Tao tao,Vec * g,PetscErrorCode (** func)(Tao tao,Vec x,Vec g,PetscCtx ctx),PetscCtxRt ctx)506*2a8381b2SBarry Smith PetscErrorCode TaoGetGradient(Tao tao, Vec *g, PetscErrorCode (**func)(Tao tao, Vec x, Vec g, PetscCtx ctx), PetscCtxRt ctx)
507d71ae5a4SJacob Faibussowitsch {
508a82e8c82SStefano Zampini PetscFunctionBegin;
509a82e8c82SStefano Zampini PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
510a82e8c82SStefano Zampini if (g) *g = tao->gradient;
511a82e8c82SStefano Zampini if (func) *func = tao->ops->computegradient;
512*2a8381b2SBarry Smith if (ctx) *(void **)ctx = tao->user_gradP;
5133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
514a82e8c82SStefano Zampini }
515a82e8c82SStefano Zampini
516a82e8c82SStefano Zampini /*@C
51765ba42b6SBarry Smith TaoSetObjectiveAndGradient - Sets a combined objective function and gradient evaluation routine for the function to be optimized
518a7e14dcfSSatish Balay
51920f4b53cSBarry Smith Logically Collective
520a7e14dcfSSatish Balay
521d8d19677SJose E. Roman Input Parameters:
52247450a7bSBarry Smith + tao - the `Tao` context
523a82e8c82SStefano Zampini . g - [optional] the vector to internally hold the gradient computation
524a7e14dcfSSatish Balay . func - the gradient function
525a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the gradient evaluation
52647450a7bSBarry Smith routine (may be `NULL`)
527a7e14dcfSSatish Balay
52820f4b53cSBarry Smith Calling sequence of `func`:
52920f4b53cSBarry Smith + tao - the optimization object
53020f4b53cSBarry Smith . x - input vector
53117477c02SJason Sarich . f - objective value (output)
532a7e14dcfSSatish Balay . g - gradient value (output)
533a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
534a7e14dcfSSatish Balay
535a7e14dcfSSatish Balay Level: beginner
536a7e14dcfSSatish Balay
53765ba42b6SBarry Smith Note:
53865ba42b6SBarry Smith For some optimization methods using a combined function can be more eifficient.
53965ba42b6SBarry Smith
5401cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetGradient()`, `TaoGetObjectiveAndGradient()`
541a7e14dcfSSatish Balay @*/
TaoSetObjectiveAndGradient(Tao tao,Vec g,PetscErrorCode (* func)(Tao tao,Vec x,PetscReal * f,Vec g,PetscCtx ctx),PetscCtx ctx)542*2a8381b2SBarry Smith PetscErrorCode TaoSetObjectiveAndGradient(Tao tao, Vec g, PetscErrorCode (*func)(Tao tao, Vec x, PetscReal *f, Vec g, PetscCtx ctx), PetscCtx ctx)
543d71ae5a4SJacob Faibussowitsch {
544a82e8c82SStefano Zampini PetscFunctionBegin;
545a82e8c82SStefano Zampini PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
546a82e8c82SStefano Zampini if (g) {
547a82e8c82SStefano Zampini PetscValidHeaderSpecific(g, VEC_CLASSID, 2);
548a82e8c82SStefano Zampini PetscCheckSameComm(tao, 1, g, 2);
5499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)g));
5509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->gradient));
551a82e8c82SStefano Zampini tao->gradient = g;
552a82e8c82SStefano Zampini }
553a82e8c82SStefano Zampini if (ctx) tao->user_objgradP = ctx;
554a82e8c82SStefano Zampini if (func) tao->ops->computeobjectiveandgradient = func;
5553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
556a82e8c82SStefano Zampini }
557a82e8c82SStefano Zampini
558a82e8c82SStefano Zampini /*@C
55965ba42b6SBarry Smith TaoGetObjectiveAndGradient - Gets the combined objective function and gradient evaluation routine for the function to be optimized
560a82e8c82SStefano Zampini
56120f4b53cSBarry Smith Not Collective
562a82e8c82SStefano Zampini
563a82e8c82SStefano Zampini Input Parameter:
56447450a7bSBarry Smith . tao - the `Tao` context
565a82e8c82SStefano Zampini
566a82e8c82SStefano Zampini Output Parameters:
567817da375SSatish Balay + g - the vector to internally hold the gradient computation
568a82e8c82SStefano Zampini . func - the gradient function
569a82e8c82SStefano Zampini - ctx - user-defined context for private data for the gradient evaluation routine
570a82e8c82SStefano Zampini
57120f4b53cSBarry Smith Calling sequence of `func`:
57220f4b53cSBarry Smith + tao - the optimizer
57320f4b53cSBarry Smith . x - input vector
574a82e8c82SStefano Zampini . f - objective value (output)
575a82e8c82SStefano Zampini . g - gradient value (output)
576a82e8c82SStefano Zampini - ctx - [optional] user-defined function context
577a82e8c82SStefano Zampini
578a82e8c82SStefano Zampini Level: beginner
579a82e8c82SStefano Zampini
5801cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`
581a82e8c82SStefano Zampini @*/
TaoGetObjectiveAndGradient(Tao tao,Vec * g,PetscErrorCode (** func)(Tao tao,Vec x,PetscReal * f,Vec g,PetscCtx ctx),PetscCtxRt ctx)582*2a8381b2SBarry Smith PetscErrorCode TaoGetObjectiveAndGradient(Tao tao, Vec *g, PetscErrorCode (**func)(Tao tao, Vec x, PetscReal *f, Vec g, PetscCtx ctx), PetscCtxRt ctx)
583d71ae5a4SJacob Faibussowitsch {
584a7e14dcfSSatish Balay PetscFunctionBegin;
585441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
586a82e8c82SStefano Zampini if (g) *g = tao->gradient;
587a82e8c82SStefano Zampini if (func) *func = tao->ops->computeobjectiveandgradient;
588*2a8381b2SBarry Smith if (ctx) *(void **)ctx = tao->user_objgradP;
5893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
590a7e14dcfSSatish Balay }
591a7e14dcfSSatish Balay
592a7e14dcfSSatish Balay /*@
593a82e8c82SStefano Zampini TaoIsObjectiveDefined - Checks to see if the user has
594a7e14dcfSSatish Balay declared an objective-only routine. Useful for determining when
59565ba42b6SBarry Smith it is appropriate to call `TaoComputeObjective()` or
59665ba42b6SBarry Smith `TaoComputeObjectiveAndGradient()`
597a7e14dcfSSatish Balay
59820f4b53cSBarry Smith Not Collective
599a7e14dcfSSatish Balay
600a82e8c82SStefano Zampini Input Parameter:
60147450a7bSBarry Smith . tao - the `Tao` context
602a82e8c82SStefano Zampini
603a82e8c82SStefano Zampini Output Parameter:
60465ba42b6SBarry Smith . flg - `PETSC_TRUE` if function routine is set by user, `PETSC_FALSE` otherwise
605a82e8c82SStefano Zampini
606a7e14dcfSSatish Balay Level: developer
607a7e14dcfSSatish Balay
6081cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetObjective()`, `TaoIsGradientDefined()`, `TaoIsObjectiveAndGradientDefined()`
609a7e14dcfSSatish Balay @*/
TaoIsObjectiveDefined(Tao tao,PetscBool * flg)610d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg)
611d71ae5a4SJacob Faibussowitsch {
612a7e14dcfSSatish Balay PetscFunctionBegin;
613441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
61483c8fe1dSLisandro Dalcin if (tao->ops->computeobjective == NULL) *flg = PETSC_FALSE;
61545cf516eSBarry Smith else *flg = PETSC_TRUE;
6163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
617a7e14dcfSSatish Balay }
618a7e14dcfSSatish Balay
619a7e14dcfSSatish Balay /*@
620a82e8c82SStefano Zampini TaoIsGradientDefined - Checks to see if the user has
621a7e14dcfSSatish Balay declared an objective-only routine. Useful for determining when
62265ba42b6SBarry Smith it is appropriate to call `TaoComputeGradient()` or
62365ba42b6SBarry Smith `TaoComputeGradientAndGradient()`
624a7e14dcfSSatish Balay
625a7e14dcfSSatish Balay Not Collective
626a7e14dcfSSatish Balay
627a82e8c82SStefano Zampini Input Parameter:
62847450a7bSBarry Smith . tao - the `Tao` context
629a82e8c82SStefano Zampini
630a82e8c82SStefano Zampini Output Parameter:
63165ba42b6SBarry Smith . flg - `PETSC_TRUE` if function routine is set by user, `PETSC_FALSE` otherwise
632a82e8c82SStefano Zampini
633a7e14dcfSSatish Balay Level: developer
634a7e14dcfSSatish Balay
6351cc06b55SBarry Smith .seealso: [](ch_tao), `TaoSetGradient()`, `TaoIsObjectiveDefined()`, `TaoIsObjectiveAndGradientDefined()`
636a7e14dcfSSatish Balay @*/
TaoIsGradientDefined(Tao tao,PetscBool * flg)637d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg)
638d71ae5a4SJacob Faibussowitsch {
639a7e14dcfSSatish Balay PetscFunctionBegin;
640441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
64183c8fe1dSLisandro Dalcin if (tao->ops->computegradient == NULL) *flg = PETSC_FALSE;
64245cf516eSBarry Smith else *flg = PETSC_TRUE;
6433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
644a7e14dcfSSatish Balay }
645a7e14dcfSSatish Balay
646a7e14dcfSSatish Balay /*@
647a82e8c82SStefano Zampini TaoIsObjectiveAndGradientDefined - Checks to see if the user has
648a7e14dcfSSatish Balay declared a joint objective/gradient routine. Useful for determining when
64965ba42b6SBarry Smith it is appropriate to call `TaoComputeObjective()` or
65065ba42b6SBarry Smith `TaoComputeObjectiveAndGradient()`
651a7e14dcfSSatish Balay
652a7e14dcfSSatish Balay Not Collective
653a7e14dcfSSatish Balay
654a82e8c82SStefano Zampini Input Parameter:
65547450a7bSBarry Smith . tao - the `Tao` context
656a82e8c82SStefano Zampini
657a82e8c82SStefano Zampini Output Parameter:
65865ba42b6SBarry Smith . flg - `PETSC_TRUE` if function routine is set by user, `PETSC_FALSE` otherwise
659a82e8c82SStefano Zampini
660a7e14dcfSSatish Balay Level: developer
661a7e14dcfSSatish Balay
6621cc06b55SBarry Smith .seealso: [](ch_tao), `TaoSetObjectiveAndGradient()`, `TaoIsObjectiveDefined()`, `TaoIsGradientDefined()`
663a7e14dcfSSatish Balay @*/
TaoIsObjectiveAndGradientDefined(Tao tao,PetscBool * flg)664d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg)
665d71ae5a4SJacob Faibussowitsch {
666a7e14dcfSSatish Balay PetscFunctionBegin;
667441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
66883c8fe1dSLisandro Dalcin if (tao->ops->computeobjectiveandgradient == NULL) *flg = PETSC_FALSE;
66945cf516eSBarry Smith else *flg = PETSC_TRUE;
6703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
671a7e14dcfSSatish Balay }
672