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 1447450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoCreate()`, `TaoSolve()`, `TaoGetSolution()` 15a7e14dcfSSatish Balay @*/ 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 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 12547450a7bSBarry Smith .seealso: [](chapter_tao), `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetGradient()` 126a7e14dcfSSatish Balay @*/ 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 17347450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoComputeGradient()`, `TaoComputeObjectiveAndGradient()`, `TaoSetObjective()` 174a7e14dcfSSatish Balay @*/ 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` 21420f4b53cSBarry Smith - 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 22247450a7bSBarry Smith .seealso: [](chapter_tao), `TaoComputeGradient()`, `TaoComputeObjectiveAndGradient()`, `TaoSetObjective()` 223a7e14dcfSSatish Balay @*/ 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 $ PetscErrorCode func(Tao tao, Vec x, PetscReal *f, void *ctx); 27120f4b53cSBarry Smith + tao - the optimizer 27220f4b53cSBarry Smith . x - input vector 273a7e14dcfSSatish Balay . f - function value 274a7e14dcfSSatish Balay - ctx - [optional] user-defined function context 275a7e14dcfSSatish Balay 276a7e14dcfSSatish Balay Level: beginner 277a7e14dcfSSatish Balay 27847450a7bSBarry Smith .seealso: [](chapter_tao), `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoGetObjective()` 279a7e14dcfSSatish Balay @*/ 280d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetObjective(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, void *), void *ctx) 281d71ae5a4SJacob Faibussowitsch { 282a7e14dcfSSatish Balay PetscFunctionBegin; 283441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 284a82e8c82SStefano Zampini if (ctx) tao->user_objP = ctx; 285a82e8c82SStefano Zampini if (func) tao->ops->computeobjective = func; 2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 287a82e8c82SStefano Zampini } 288a82e8c82SStefano Zampini 289a82e8c82SStefano Zampini /*@C 29065ba42b6SBarry Smith TaoGetObjective - Gets the function evaluation routine for the function to be minimized 291a82e8c82SStefano Zampini 29220f4b53cSBarry Smith Not Collective 293a82e8c82SStefano Zampini 294a82e8c82SStefano Zampini Input Parameter: 29547450a7bSBarry Smith . tao - the `Tao` context 296a82e8c82SStefano Zampini 297*2fe279fdSBarry Smith Output Parameters: 298a82e8c82SStefano Zampini + func - the objective function 299a82e8c82SStefano Zampini - ctx - the user-defined context for private data for the function evaluation 300a82e8c82SStefano Zampini 30120f4b53cSBarry Smith Calling sequence of `func`: 302*2fe279fdSBarry Smith $ PetscErrorCode func(Tao tao, Vec x, PetscReal *f, void *ctx) 30320f4b53cSBarry Smith + tao - the optimizer 30420f4b53cSBarry Smith . x - input vector 305a82e8c82SStefano Zampini . f - function value 306a82e8c82SStefano Zampini - ctx - [optional] user-defined function context 307a82e8c82SStefano Zampini 308a82e8c82SStefano Zampini Level: beginner 309a82e8c82SStefano Zampini 31047450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjective()` 311a82e8c82SStefano Zampini @*/ 312d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetObjective(Tao tao, PetscErrorCode (**func)(Tao, Vec, PetscReal *, void *), void **ctx) 313d71ae5a4SJacob Faibussowitsch { 314a82e8c82SStefano Zampini PetscFunctionBegin; 315a82e8c82SStefano Zampini PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 316a82e8c82SStefano Zampini if (func) *func = tao->ops->computeobjective; 317a82e8c82SStefano Zampini if (ctx) *ctx = tao->user_objP; 3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 319a7e14dcfSSatish Balay } 320a7e14dcfSSatish Balay 321a7e14dcfSSatish Balay /*@C 3224a48860cSAlp Dener TaoSetResidualRoutine - Sets the residual evaluation routine for least-square applications 323a7e14dcfSSatish Balay 32420f4b53cSBarry Smith Logically Collective 325a7e14dcfSSatish Balay 326d8d19677SJose E. Roman Input Parameters: 32747450a7bSBarry Smith + tao - the `Tao` context 3284a48860cSAlp Dener . func - the residual evaluation routine 329a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the function evaluation 33047450a7bSBarry Smith routine (may be `NULL`) 331a7e14dcfSSatish Balay 33220f4b53cSBarry Smith Calling sequence of `func`: 33320f4b53cSBarry Smith $ PetscErrorCode func(Tao tao, Vec x, Vec f, void *ctx); 33420f4b53cSBarry Smith + tao - the optimizer 33520f4b53cSBarry Smith . x - input vector 336a7e14dcfSSatish Balay . f - function value vector 337a7e14dcfSSatish Balay - ctx - [optional] user-defined function context 338a7e14dcfSSatish Balay 339a7e14dcfSSatish Balay Level: beginner 340a7e14dcfSSatish Balay 34147450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoSetObjective()`, `TaoSetJacobianRoutine()` 342a7e14dcfSSatish Balay @*/ 343d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetResidualRoutine(Tao tao, Vec res, PetscErrorCode (*func)(Tao, Vec, Vec, void *), void *ctx) 344d71ae5a4SJacob Faibussowitsch { 345a7e14dcfSSatish Balay PetscFunctionBegin; 346441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 3474a48860cSAlp Dener PetscValidHeaderSpecific(res, VEC_CLASSID, 2); 3489566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)res)); 34948a46eb9SPierre Jolivet if (tao->ls_res) PetscCall(VecDestroy(&tao->ls_res)); 3504a48860cSAlp Dener tao->ls_res = res; 3514ffbe8acSAlp Dener tao->user_lsresP = ctx; 3524a48860cSAlp Dener tao->ops->computeresidual = func; 353737f463aSAlp Dener 3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 355737f463aSAlp Dener } 356737f463aSAlp Dener 357737f463aSAlp Dener /*@ 35865ba42b6SBarry 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. 35947450a7bSBarry Smith If this function is not provided, or if `sigma_v` and `vals` are both `NULL`, then the identity matrix will be used for weights. 360737f463aSAlp Dener 361c3339decSBarry Smith Collective 362737f463aSAlp Dener 363737f463aSAlp Dener Input Parameters: 36447450a7bSBarry Smith + tao - the `Tao` context 365737f463aSAlp Dener . sigma_v - vector of weights (diagonal terms only) 366737f463aSAlp Dener . n - the number of weights (if using off-diagonal) 36747450a7bSBarry Smith . rows - index list of rows for `sigma_v` 36847450a7bSBarry Smith . cols - index list of columns for `sigma_v` 369737f463aSAlp Dener - vals - array of weights 370737f463aSAlp Dener 371737f463aSAlp Dener Level: intermediate 372737f463aSAlp Dener 37347450a7bSBarry Smith Note: 37447450a7bSBarry Smith Either `sigma_v` or `vals` should be `NULL` 37547450a7bSBarry Smith 37647450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoSetResidualRoutine()` 377737f463aSAlp Dener @*/ 378d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetResidualWeights(Tao tao, Vec sigma_v, PetscInt n, PetscInt *rows, PetscInt *cols, PetscReal *vals) 379d71ae5a4SJacob Faibussowitsch { 380737f463aSAlp Dener PetscInt i; 381a82e8c82SStefano Zampini 382737f463aSAlp Dener PetscFunctionBegin; 383737f463aSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 384a82e8c82SStefano Zampini if (sigma_v) PetscValidHeaderSpecific(sigma_v, VEC_CLASSID, 2); 3859566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sigma_v)); 3869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->res_weights_v)); 3874ffbe8acSAlp Dener tao->res_weights_v = sigma_v; 388737f463aSAlp Dener if (vals) { 3899566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->res_weights_rows)); 3909566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->res_weights_cols)); 3919566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->res_weights_w)); 3929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tao->res_weights_rows)); 3939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tao->res_weights_cols)); 3949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tao->res_weights_w)); 395737f463aSAlp Dener tao->res_weights_n = n; 396737f463aSAlp Dener for (i = 0; i < n; i++) { 397737f463aSAlp Dener tao->res_weights_rows[i] = rows[i]; 398737f463aSAlp Dener tao->res_weights_cols[i] = cols[i]; 399737f463aSAlp Dener tao->res_weights_w[i] = vals[i]; 400737f463aSAlp Dener } 401737f463aSAlp Dener } else { 402737f463aSAlp Dener tao->res_weights_n = 0; 40383c8fe1dSLisandro Dalcin tao->res_weights_rows = NULL; 40483c8fe1dSLisandro Dalcin tao->res_weights_cols = NULL; 405737f463aSAlp Dener } 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 407a7e14dcfSSatish Balay } 408a7e14dcfSSatish Balay 4098b7a9b22SJason Sarich /*@ 4104a48860cSAlp Dener TaoComputeResidual - Computes a least-squares residual vector at a given point 411a7e14dcfSSatish Balay 412c3339decSBarry Smith Collective 413a7e14dcfSSatish Balay 414a7e14dcfSSatish Balay Input Parameters: 41547450a7bSBarry Smith + tao - the `Tao` context 416a7e14dcfSSatish Balay - X - input vector 417a7e14dcfSSatish Balay 418a7e14dcfSSatish Balay Output Parameter: 41920f4b53cSBarry Smith . f - Objective vector at `X` 420a7e14dcfSSatish Balay 42147450a7bSBarry Smith Level: advanced 42247450a7bSBarry Smith 42395452b02SPatrick Sanan Notes: 42465ba42b6SBarry Smith `TaoComputeResidual()` is typically used within the implementation of the optimization algorithm, 425a7e14dcfSSatish Balay so most users would not generally call this routine themselves. 426a7e14dcfSSatish Balay 42747450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoSetResidualRoutine()` 428a7e14dcfSSatish Balay @*/ 429d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeResidual(Tao tao, Vec X, Vec F) 430d71ae5a4SJacob Faibussowitsch { 431a7e14dcfSSatish Balay PetscFunctionBegin; 432441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 433a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 434a7e14dcfSSatish Balay PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 435a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2); 436a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, F, 3); 4374a48860cSAlp Dener if (tao->ops->computeresidual) { 4389566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_ObjectiveEval, tao, X, NULL, NULL)); 439792fecdfSBarry Smith PetscCallBack("Tao callback least-squares residual", (*tao->ops->computeresidual)(tao, X, F, tao->user_lsresP)); 4409566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_ObjectiveEval, tao, X, NULL, NULL)); 441a7e14dcfSSatish Balay tao->nfuncs++; 442691b26d3SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TaoSetResidualRoutine() has not been called"); 4439566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "TAO least-squares residual evaluation.\n")); 4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 445a7e14dcfSSatish Balay } 446a7e14dcfSSatish Balay 447a7e14dcfSSatish Balay /*@C 44865ba42b6SBarry Smith TaoSetGradient - Sets the gradient evaluation routine for the function to be optimized 449a7e14dcfSSatish Balay 45020f4b53cSBarry Smith Logically Collective 451a7e14dcfSSatish Balay 452d8d19677SJose E. Roman Input Parameters: 45347450a7bSBarry Smith + tao - the `Tao` context 454a82e8c82SStefano Zampini . g - [optional] the vector to internally hold the gradient computation 455a7e14dcfSSatish Balay . func - the gradient function 456a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the gradient evaluation 45747450a7bSBarry Smith routine (may be `NULL`) 458a7e14dcfSSatish Balay 45920f4b53cSBarry Smith Calling sequence of `func`: 46020f4b53cSBarry Smith $ PetscErrorCode func(Tao tao, Vec x, Vec g, void *ctx); 46120f4b53cSBarry Smith + tao - the optimization solver 46220f4b53cSBarry Smith . x - input vector 463a7e14dcfSSatish Balay . g - gradient value (output) 464a7e14dcfSSatish Balay - ctx - [optional] user-defined function context 465a7e14dcfSSatish Balay 466a7e14dcfSSatish Balay Level: beginner 467a7e14dcfSSatish Balay 46847450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoGetGradient()` 469a7e14dcfSSatish Balay @*/ 470d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetGradient(Tao tao, Vec g, PetscErrorCode (*func)(Tao, Vec, Vec, void *), void *ctx) 471d71ae5a4SJacob Faibussowitsch { 472a7e14dcfSSatish Balay PetscFunctionBegin; 473441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 474a82e8c82SStefano Zampini if (g) { 475a82e8c82SStefano Zampini PetscValidHeaderSpecific(g, VEC_CLASSID, 2); 476a82e8c82SStefano Zampini PetscCheckSameComm(tao, 1, g, 2); 4779566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)g)); 4789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->gradient)); 479a82e8c82SStefano Zampini tao->gradient = g; 480a82e8c82SStefano Zampini } 481a82e8c82SStefano Zampini if (func) tao->ops->computegradient = func; 482a82e8c82SStefano Zampini if (ctx) tao->user_gradP = ctx; 4833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 484a7e14dcfSSatish Balay } 485a7e14dcfSSatish Balay 486a7e14dcfSSatish Balay /*@C 48765ba42b6SBarry Smith TaoGetGradient - Gets the gradient evaluation routine for the function being optimized 488a82e8c82SStefano Zampini 48920f4b53cSBarry Smith Not Collective 490a82e8c82SStefano Zampini 491a82e8c82SStefano Zampini Input Parameter: 49247450a7bSBarry Smith . tao - the `Tao` context 493a82e8c82SStefano Zampini 494a82e8c82SStefano Zampini Output Parameters: 495a82e8c82SStefano Zampini + g - the vector to internally hold the gradient computation 496a82e8c82SStefano Zampini . func - the gradient function 497a82e8c82SStefano Zampini - ctx - user-defined context for private data for the gradient evaluation routine 498a82e8c82SStefano Zampini 49920f4b53cSBarry Smith Calling sequence of `func`: 50020f4b53cSBarry Smith $ PetscErrorCode func(Tao tao, Vec x, Vec g, void *ctx); 50120f4b53cSBarry Smith + tao - the optimizer 50220f4b53cSBarry Smith . x - input vector 503a82e8c82SStefano Zampini . g - gradient value (output) 504a82e8c82SStefano Zampini - ctx - [optional] user-defined function context 505a82e8c82SStefano Zampini 506a82e8c82SStefano Zampini Level: beginner 507a82e8c82SStefano Zampini 50847450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`, `TaoSetGradient()` 509a82e8c82SStefano Zampini @*/ 510d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetGradient(Tao tao, Vec *g, PetscErrorCode (**func)(Tao, Vec, Vec, void *), void **ctx) 511d71ae5a4SJacob Faibussowitsch { 512a82e8c82SStefano Zampini PetscFunctionBegin; 513a82e8c82SStefano Zampini PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 514a82e8c82SStefano Zampini if (g) *g = tao->gradient; 515a82e8c82SStefano Zampini if (func) *func = tao->ops->computegradient; 516a82e8c82SStefano Zampini if (ctx) *ctx = tao->user_gradP; 5173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 518a82e8c82SStefano Zampini } 519a82e8c82SStefano Zampini 520a82e8c82SStefano Zampini /*@C 52165ba42b6SBarry Smith TaoSetObjectiveAndGradient - Sets a combined objective function and gradient evaluation routine for the function to be optimized 522a7e14dcfSSatish Balay 52320f4b53cSBarry Smith Logically Collective 524a7e14dcfSSatish Balay 525d8d19677SJose E. Roman Input Parameters: 52647450a7bSBarry Smith + tao - the `Tao` context 527a82e8c82SStefano Zampini . g - [optional] the vector to internally hold the gradient computation 528a7e14dcfSSatish Balay . func - the gradient function 529a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the gradient evaluation 53047450a7bSBarry Smith routine (may be `NULL`) 531a7e14dcfSSatish Balay 53220f4b53cSBarry Smith Calling sequence of `func`: 53320f4b53cSBarry Smith $ PetscErrorCode func(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx); 53420f4b53cSBarry Smith + tao - the optimization object 53520f4b53cSBarry Smith . x - input vector 53617477c02SJason Sarich . f - objective value (output) 537a7e14dcfSSatish Balay . g - gradient value (output) 538a7e14dcfSSatish Balay - ctx - [optional] user-defined function context 539a7e14dcfSSatish Balay 540a7e14dcfSSatish Balay Level: beginner 541a7e14dcfSSatish Balay 54265ba42b6SBarry Smith Note: 54365ba42b6SBarry Smith For some optimization methods using a combined function can be more eifficient. 54465ba42b6SBarry Smith 54547450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetHessian()`, `TaoSetGradient()`, `TaoGetObjectiveAndGradient()` 546a7e14dcfSSatish Balay @*/ 547d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetObjectiveAndGradient(Tao tao, Vec g, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void *), void *ctx) 548d71ae5a4SJacob Faibussowitsch { 549a82e8c82SStefano Zampini PetscFunctionBegin; 550a82e8c82SStefano Zampini PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 551a82e8c82SStefano Zampini if (g) { 552a82e8c82SStefano Zampini PetscValidHeaderSpecific(g, VEC_CLASSID, 2); 553a82e8c82SStefano Zampini PetscCheckSameComm(tao, 1, g, 2); 5549566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)g)); 5559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->gradient)); 556a82e8c82SStefano Zampini tao->gradient = g; 557a82e8c82SStefano Zampini } 558a82e8c82SStefano Zampini if (ctx) tao->user_objgradP = ctx; 559a82e8c82SStefano Zampini if (func) tao->ops->computeobjectiveandgradient = func; 5603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 561a82e8c82SStefano Zampini } 562a82e8c82SStefano Zampini 563a82e8c82SStefano Zampini /*@C 56465ba42b6SBarry Smith TaoGetObjectiveAndGradient - Gets the combined objective function and gradient evaluation routine for the function to be optimized 565a82e8c82SStefano Zampini 56620f4b53cSBarry Smith Not Collective 567a82e8c82SStefano Zampini 568a82e8c82SStefano Zampini Input Parameter: 56947450a7bSBarry Smith . tao - the `Tao` context 570a82e8c82SStefano Zampini 571a82e8c82SStefano Zampini Output Parameters: 572817da375SSatish Balay + g - the vector to internally hold the gradient computation 573a82e8c82SStefano Zampini . func - the gradient function 574a82e8c82SStefano Zampini - ctx - user-defined context for private data for the gradient evaluation routine 575a82e8c82SStefano Zampini 57620f4b53cSBarry Smith Calling sequence of `func`: 57720f4b53cSBarry Smith $ PetscErrorCode func(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx); 57820f4b53cSBarry Smith + tao - the optimizer 57920f4b53cSBarry Smith . x - input vector 580a82e8c82SStefano Zampini . f - objective value (output) 581a82e8c82SStefano Zampini . g - gradient value (output) 582a82e8c82SStefano Zampini - ctx - [optional] user-defined function context 583a82e8c82SStefano Zampini 584a82e8c82SStefano Zampini Level: beginner 585a82e8c82SStefano Zampini 58647450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()` 587a82e8c82SStefano Zampini @*/ 588d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetObjectiveAndGradient(Tao tao, Vec *g, PetscErrorCode (**func)(Tao, Vec, PetscReal *, Vec, void *), void **ctx) 589d71ae5a4SJacob Faibussowitsch { 590a7e14dcfSSatish Balay PetscFunctionBegin; 591441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 592a82e8c82SStefano Zampini if (g) *g = tao->gradient; 593a82e8c82SStefano Zampini if (func) *func = tao->ops->computeobjectiveandgradient; 594a82e8c82SStefano Zampini if (ctx) *ctx = tao->user_objgradP; 5953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 596a7e14dcfSSatish Balay } 597a7e14dcfSSatish Balay 598a7e14dcfSSatish Balay /*@ 599a82e8c82SStefano Zampini TaoIsObjectiveDefined - Checks to see if the user has 600a7e14dcfSSatish Balay declared an objective-only routine. Useful for determining when 60165ba42b6SBarry Smith it is appropriate to call `TaoComputeObjective()` or 60265ba42b6SBarry Smith `TaoComputeObjectiveAndGradient()` 603a7e14dcfSSatish Balay 60420f4b53cSBarry Smith Not Collective 605a7e14dcfSSatish Balay 606a82e8c82SStefano Zampini Input Parameter: 60747450a7bSBarry Smith . tao - the `Tao` context 608a82e8c82SStefano Zampini 609a82e8c82SStefano Zampini Output Parameter: 61065ba42b6SBarry Smith . flg - `PETSC_TRUE` if function routine is set by user, `PETSC_FALSE` otherwise 611a82e8c82SStefano Zampini 612a7e14dcfSSatish Balay Level: developer 613a7e14dcfSSatish Balay 61447450a7bSBarry Smith .seealso: [](chapter_tao), `Tao`, `TaoSetObjective()`, `TaoIsGradientDefined()`, `TaoIsObjectiveAndGradientDefined()` 615a7e14dcfSSatish Balay @*/ 616d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg) 617d71ae5a4SJacob Faibussowitsch { 618a7e14dcfSSatish Balay PetscFunctionBegin; 619441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 62083c8fe1dSLisandro Dalcin if (tao->ops->computeobjective == NULL) *flg = PETSC_FALSE; 62145cf516eSBarry Smith else *flg = PETSC_TRUE; 6223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 623a7e14dcfSSatish Balay } 624a7e14dcfSSatish Balay 625a7e14dcfSSatish Balay /*@ 626a82e8c82SStefano Zampini TaoIsGradientDefined - Checks to see if the user has 627a7e14dcfSSatish Balay declared an objective-only routine. Useful for determining when 62865ba42b6SBarry Smith it is appropriate to call `TaoComputeGradient()` or 62965ba42b6SBarry Smith `TaoComputeGradientAndGradient()` 630a7e14dcfSSatish Balay 631a7e14dcfSSatish Balay Not Collective 632a7e14dcfSSatish Balay 633a82e8c82SStefano Zampini Input Parameter: 63447450a7bSBarry Smith . tao - the `Tao` context 635a82e8c82SStefano Zampini 636a82e8c82SStefano Zampini Output Parameter: 63765ba42b6SBarry Smith . flg - `PETSC_TRUE` if function routine is set by user, `PETSC_FALSE` otherwise 638a82e8c82SStefano Zampini 639a7e14dcfSSatish Balay Level: developer 640a7e14dcfSSatish Balay 64147450a7bSBarry Smith .seealso: [](chapter_tao), `TaoSetGradient()`, `TaoIsObjectiveDefined()`, `TaoIsObjectiveAndGradientDefined()` 642a7e14dcfSSatish Balay @*/ 643d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg) 644d71ae5a4SJacob Faibussowitsch { 645a7e14dcfSSatish Balay PetscFunctionBegin; 646441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 64783c8fe1dSLisandro Dalcin if (tao->ops->computegradient == NULL) *flg = PETSC_FALSE; 64845cf516eSBarry Smith else *flg = PETSC_TRUE; 6493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 650a7e14dcfSSatish Balay } 651a7e14dcfSSatish Balay 652a7e14dcfSSatish Balay /*@ 653a82e8c82SStefano Zampini TaoIsObjectiveAndGradientDefined - Checks to see if the user has 654a7e14dcfSSatish Balay declared a joint objective/gradient routine. Useful for determining when 65565ba42b6SBarry Smith it is appropriate to call `TaoComputeObjective()` or 65665ba42b6SBarry Smith `TaoComputeObjectiveAndGradient()` 657a7e14dcfSSatish Balay 658a7e14dcfSSatish Balay Not Collective 659a7e14dcfSSatish Balay 660a82e8c82SStefano Zampini Input Parameter: 66147450a7bSBarry Smith . tao - the `Tao` context 662a82e8c82SStefano Zampini 663a82e8c82SStefano Zampini Output Parameter: 66465ba42b6SBarry Smith . flg - `PETSC_TRUE` if function routine is set by user, `PETSC_FALSE` otherwise 665a82e8c82SStefano Zampini 666a7e14dcfSSatish Balay Level: developer 667a7e14dcfSSatish Balay 66847450a7bSBarry Smith .seealso: [](chapter_tao), `TaoSetObjectiveAndGradient()`, `TaoIsObjectiveDefined()`, `TaoIsGradientDefined()` 669a7e14dcfSSatish Balay @*/ 670d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg) 671d71ae5a4SJacob Faibussowitsch { 672a7e14dcfSSatish Balay PetscFunctionBegin; 673441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 67483c8fe1dSLisandro Dalcin if (tao->ops->computeobjectiveandgradient == NULL) *flg = PETSC_FALSE; 67545cf516eSBarry Smith else *flg = PETSC_TRUE; 6763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 677a7e14dcfSSatish Balay } 678