1af0996ceSBarry Smith #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
2a7e14dcfSSatish Balay
3a7e14dcfSSatish Balay /*@C
4a82e8c82SStefano Zampini TaoSetHessian - Sets the function to compute the Hessian as well as the location to store the matrix.
5a7e14dcfSSatish Balay
620f4b53cSBarry Smith Logically Collective
7a7e14dcfSSatish Balay
8a7e14dcfSSatish Balay Input Parameters:
947450a7bSBarry Smith + tao - the `Tao` context
10a7e14dcfSSatish Balay . H - Matrix used for the hessian
1147450a7bSBarry Smith . Hpre - Matrix that will be used to construct the preconditioner, can be same as `H`
12f4c1ad5cSStefano Zampini . func - Hessian evaluation routine
13a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the
1447450a7bSBarry Smith Hessian evaluation routine (may be `NULL`)
15a7e14dcfSSatish Balay
1620f4b53cSBarry Smith Calling sequence of `func`:
178847d985SBarry Smith + tao - the `Tao` context
18a7e14dcfSSatish Balay . x - input vector
19a7e14dcfSSatish Balay . H - Hessian matrix
2047450a7bSBarry Smith . Hpre - matrix used to construct the preconditioner, usually the same as `H`
21a7e14dcfSSatish Balay - ctx - [optional] user-defined Hessian context
22a7e14dcfSSatish Balay
23a7e14dcfSSatish Balay Level: beginner
24a82e8c82SStefano Zampini
251cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoTypes`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetObjectiveAndGradient()`, `TaoGetHessian()`
26a7e14dcfSSatish Balay @*/
TaoSetHessian(Tao tao,Mat H,Mat Hpre,PetscErrorCode (* func)(Tao tao,Vec x,Mat H,Mat Hpre,PetscCtx ctx),PetscCtx ctx)27*2a8381b2SBarry Smith PetscErrorCode TaoSetHessian(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao tao, Vec x, Mat H, Mat Hpre, PetscCtx ctx), PetscCtx ctx)
28d71ae5a4SJacob Faibussowitsch {
29a7e14dcfSSatish Balay PetscFunctionBegin;
30441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
31a7e14dcfSSatish Balay if (H) {
32a7e14dcfSSatish Balay PetscValidHeaderSpecific(H, MAT_CLASSID, 2);
33a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, H, 2);
34a7e14dcfSSatish Balay }
35a7e14dcfSSatish Balay if (Hpre) {
36a7e14dcfSSatish Balay PetscValidHeaderSpecific(Hpre, MAT_CLASSID, 3);
37a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, Hpre, 3);
38a7e14dcfSSatish Balay }
39a82e8c82SStefano Zampini if (ctx) tao->user_hessP = ctx;
40a82e8c82SStefano Zampini if (func) tao->ops->computehessian = func;
41a7e14dcfSSatish Balay if (H) {
429566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H));
439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->hessian));
44a7e14dcfSSatish Balay tao->hessian = H;
45a7e14dcfSSatish Balay }
46a7e14dcfSSatish Balay if (Hpre) {
479566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Hpre));
489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->hessian_pre));
49a7e14dcfSSatish Balay tao->hessian_pre = Hpre;
50a7e14dcfSSatish Balay }
513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
52a7e14dcfSSatish Balay }
53a7e14dcfSSatish Balay
54a82e8c82SStefano Zampini /*@C
55a82e8c82SStefano Zampini TaoGetHessian - Gets the function to compute the Hessian as well as the location to store the matrix.
56a82e8c82SStefano Zampini
5720f4b53cSBarry Smith Not Collective
58a82e8c82SStefano Zampini
59a82e8c82SStefano Zampini Input Parameter:
6047450a7bSBarry Smith . tao - the `Tao` context
61a82e8c82SStefano Zampini
62a82e8c82SStefano Zampini Output Parameters:
63a82e8c82SStefano Zampini + H - Matrix used for the hessian
6447450a7bSBarry Smith . Hpre - Matrix that will be used to construct the preconditioner, can be the same as `H`
65a82e8c82SStefano Zampini . func - Hessian evaluation routine
66a82e8c82SStefano Zampini - ctx - user-defined context for private data for the Hessian evaluation routine
67a82e8c82SStefano Zampini
6820f4b53cSBarry Smith Calling sequence of `func`:
698847d985SBarry Smith + tao - the `Tao` context
70a82e8c82SStefano Zampini . x - input vector
71a82e8c82SStefano Zampini . H - Hessian matrix
7247450a7bSBarry Smith . Hpre - matrix used to construct the preconditioner, usually the same as `H`
73a82e8c82SStefano Zampini - ctx - [optional] user-defined Hessian context
74a82e8c82SStefano Zampini
75a82e8c82SStefano Zampini Level: beginner
76a82e8c82SStefano Zampini
77e056e8ceSJacob Faibussowitsch .seealso: [](ch_tao), `Tao`, `TaoType`, `TaoGetObjective()`, `TaoGetGradient()`, `TaoGetObjectiveAndGradient()`, `TaoSetHessian()`
78a82e8c82SStefano Zampini @*/
TaoGetHessian(Tao tao,Mat * H,Mat * Hpre,PetscErrorCode (** func)(Tao tao,Vec x,Mat H,Mat Hpre,PetscCtx ctx),PetscCtxRt ctx)79*2a8381b2SBarry Smith PetscErrorCode TaoGetHessian(Tao tao, Mat *H, Mat *Hpre, PetscErrorCode (**func)(Tao tao, Vec x, Mat H, Mat Hpre, PetscCtx ctx), PetscCtxRt ctx)
80d71ae5a4SJacob Faibussowitsch {
81a82e8c82SStefano Zampini PetscFunctionBegin;
82a82e8c82SStefano Zampini PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
83a82e8c82SStefano Zampini if (H) *H = tao->hessian;
84a82e8c82SStefano Zampini if (Hpre) *Hpre = tao->hessian_pre;
85*2a8381b2SBarry Smith if (ctx) *(void **)ctx = tao->user_hessP;
86a82e8c82SStefano Zampini if (func) *func = tao->ops->computehessian;
873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
88a82e8c82SStefano Zampini }
89a82e8c82SStefano Zampini
TaoTestHessian(Tao tao)90d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoTestHessian(Tao tao)
91d71ae5a4SJacob Faibussowitsch {
9209baa881SHong Zhang Mat A, B, C, D, hessian;
9309baa881SHong Zhang Vec x = tao->solution;
9409baa881SHong Zhang PetscReal nrm, gnorm;
9509baa881SHong Zhang PetscReal threshold = 1.e-5;
9609baa881SHong Zhang PetscInt m, n, M, N;
9709baa881SHong Zhang PetscBool complete_print = PETSC_FALSE, test = PETSC_FALSE, flg;
98f49d1c87SHong Zhang PetscViewer viewer, mviewer;
9909baa881SHong Zhang MPI_Comm comm;
10009baa881SHong Zhang PetscInt tabs;
10109baa881SHong Zhang static PetscBool directionsprinted = PETSC_FALSE;
102f49d1c87SHong Zhang PetscViewerFormat format;
10309baa881SHong Zhang
10409baa881SHong Zhang PetscFunctionBegin;
105d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)tao);
1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-tao_test_hessian", "Compare hand-coded and finite difference Hessians", "None", &test));
1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_test_hessian", "Threshold for element difference between hand-coded and finite difference being meaningful", "None", threshold, &threshold, NULL));
1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsViewer("-tao_test_hessian_view", "View difference between hand-coded and finite difference Hessians element entries", "None", &mviewer, &format, &complete_print));
109d0609cedSBarry Smith PetscOptionsEnd();
1103ba16761SJacob Faibussowitsch if (!test) PetscFunctionReturn(PETSC_SUCCESS);
11109baa881SHong Zhang
1129566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
1139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(comm, &viewer));
1149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------- Testing Hessian -------------\n"));
11709baa881SHong Zhang if (!complete_print && !directionsprinted) {
1189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Run with -tao_test_hessian_view and optionally -tao_test_hessian <threshold> to show difference\n"));
1199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " of hand-coded and finite difference Hessian entries greater than <threshold>.\n"));
12009baa881SHong Zhang }
12109baa881SHong Zhang if (!directionsprinted) {
1229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Testing hand-coded Hessian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n"));
1239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " O(1.e-8), the hand-coded Hessian is probably correct.\n"));
12409baa881SHong Zhang directionsprinted = PETSC_TRUE;
12509baa881SHong Zhang }
1261baa6e33SBarry Smith if (complete_print) PetscCall(PetscViewerPushFormat(mviewer, format));
12709baa881SHong Zhang
1289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian, MATMFFD, &flg));
12909baa881SHong Zhang if (!flg) hessian = tao->hessian;
13009baa881SHong Zhang else hessian = tao->hessian_pre;
13109baa881SHong Zhang
13209baa881SHong Zhang while (hessian) {
1332b2f8cc6SPierre Jolivet PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)hessian, &flg, MATSEQAIJ, MATMPIAIJ, MATSEQDENSE, MATMPIDENSE, MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ, ""));
13409baa881SHong Zhang if (flg) {
13509baa881SHong Zhang A = hessian;
1369566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A));
13709baa881SHong Zhang } else {
1389566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(hessian, MATAIJ, &A));
13909baa881SHong Zhang }
14009baa881SHong Zhang
1419566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
1429566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N));
1439566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n));
1449566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, n, M, N));
1459566063dSJacob Faibussowitsch PetscCall(MatSetType(B, ((PetscObject)A)->type_name));
1469566063dSJacob Faibussowitsch PetscCall(MatSetUp(B));
1479566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
14809baa881SHong Zhang
1499566063dSJacob Faibussowitsch PetscCall(TaoDefaultComputeHessian(tao, x, B, B, NULL));
15009baa881SHong Zhang
1519566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &D));
1529566063dSJacob Faibussowitsch PetscCall(MatAYPX(D, -1.0, A, DIFFERENT_NONZERO_PATTERN));
1539566063dSJacob Faibussowitsch PetscCall(MatNorm(D, NORM_FROBENIUS, &nrm));
1549566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_FROBENIUS, &gnorm));
1559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D));
15609baa881SHong Zhang if (!gnorm) gnorm = 1; /* just in case */
1579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ||H - Hfd||_F/||H||_F = %g, ||H - Hfd||_F = %g\n", (double)(nrm / gnorm), (double)nrm));
15809baa881SHong Zhang
15909baa881SHong Zhang if (complete_print) {
1609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Hand-coded Hessian ----------\n"));
1619566063dSJacob Faibussowitsch PetscCall(MatView(A, mviewer));
1629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Finite difference Hessian ----------\n"));
1639566063dSJacob Faibussowitsch PetscCall(MatView(B, mviewer));
16409baa881SHong Zhang }
16509baa881SHong Zhang
16609baa881SHong Zhang if (complete_print) {
16709baa881SHong Zhang PetscInt Istart, Iend, *ccols, bncols, cncols, j, row;
16809baa881SHong Zhang PetscScalar *cvals;
16909baa881SHong Zhang const PetscInt *bcols;
17009baa881SHong Zhang const PetscScalar *bvals;
17109baa881SHong Zhang
1729566063dSJacob Faibussowitsch PetscCall(MatAYPX(B, -1.0, A, DIFFERENT_NONZERO_PATTERN));
1739566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
1749566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m, n, M, N));
1759566063dSJacob Faibussowitsch PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
1769566063dSJacob Faibussowitsch PetscCall(MatSetUp(C));
1779566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
1789566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B, &Istart, &Iend));
17909baa881SHong Zhang
18009baa881SHong Zhang for (row = Istart; row < Iend; row++) {
1819566063dSJacob Faibussowitsch PetscCall(MatGetRow(B, row, &bncols, &bcols, &bvals));
1829566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(bncols, &ccols, bncols, &cvals));
18309baa881SHong Zhang for (j = 0, cncols = 0; j < bncols; j++) {
18409baa881SHong Zhang if (PetscAbsScalar(bvals[j]) > threshold) {
18509baa881SHong Zhang ccols[cncols] = bcols[j];
18609baa881SHong Zhang cvals[cncols] = bvals[j];
18709baa881SHong Zhang cncols += 1;
18809baa881SHong Zhang }
18909baa881SHong Zhang }
19048a46eb9SPierre Jolivet if (cncols) PetscCall(MatSetValues(C, 1, &row, cncols, ccols, cvals, INSERT_VALUES));
1919566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(B, row, &bncols, &bcols, &bvals));
1929566063dSJacob Faibussowitsch PetscCall(PetscFree2(ccols, cvals));
19309baa881SHong Zhang }
1949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
1969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Finite-difference minus hand-coded Hessian with tolerance %g ----------\n", (double)threshold));
1979566063dSJacob Faibussowitsch PetscCall(MatView(C, mviewer));
1989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
19909baa881SHong Zhang }
2009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
2019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
20209baa881SHong Zhang
20309baa881SHong Zhang if (hessian != tao->hessian_pre) {
20409baa881SHong Zhang hessian = tao->hessian_pre;
2059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ---------- Testing Hessian for preconditioner -------------\n"));
20609baa881SHong Zhang } else hessian = NULL;
20709baa881SHong Zhang }
208f49d1c87SHong Zhang if (complete_print) {
2099566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(mviewer));
2109566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&mviewer));
211f49d1c87SHong Zhang }
2129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISetTab(viewer, tabs));
2133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
21409baa881SHong Zhang }
21509baa881SHong Zhang
216cc4c1da9SBarry Smith /*@
217a7e14dcfSSatish Balay TaoComputeHessian - Computes the Hessian matrix that has been
21865ba42b6SBarry Smith set with `TaoSetHessian()`.
219a7e14dcfSSatish Balay
220c3339decSBarry Smith Collective
221a7e14dcfSSatish Balay
222a7e14dcfSSatish Balay Input Parameters:
223f4c1ad5cSStefano Zampini + tao - the Tao solver context
224f4c1ad5cSStefano Zampini - X - input vector
225a7e14dcfSSatish Balay
226a7e14dcfSSatish Balay Output Parameters:
227a7e14dcfSSatish Balay + H - Hessian matrix
2287addb90fSBarry Smith - Hpre - matrix used to construct the preconditioner, usually the same as `H`
229a7e14dcfSSatish Balay
23009baa881SHong Zhang Options Database Keys:
23109baa881SHong Zhang + -tao_test_hessian - compare the user provided Hessian with one compute via finite differences to check for errors
23209baa881SHong Zhang . -tao_test_hessian <numerical value> - display entries in the difference between the user provided Hessian and finite difference Hessian that are greater than a certain value to help users detect errors
233dfe02fe6SHong Zhang - -tao_test_hessian_view - display the user provided Hessian, the finite difference Hessian and the difference between them to help users detect the location of errors in the user provided Hessian
23409baa881SHong Zhang
23547450a7bSBarry Smith Level: developer
23647450a7bSBarry Smith
237a7e14dcfSSatish Balay Notes:
238a7e14dcfSSatish Balay Most users should not need to explicitly call this routine, as it
239a7e14dcfSSatish Balay is used internally within the minimization solvers.
240a7e14dcfSSatish Balay
24165ba42b6SBarry Smith `TaoComputeHessian()` is typically used within optimization algorithms,
24265ba42b6SBarry Smith so most users would not generally call this routine
243a7e14dcfSSatish Balay themselves.
244a7e14dcfSSatish Balay
245e056e8ceSJacob Faibussowitsch Developer Notes:
24665ba42b6SBarry Smith The Hessian test mechanism follows `SNESTestJacobian()`.
24709baa881SHong Zhang
2481cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetHessian()`
249a7e14dcfSSatish Balay @*/
TaoComputeHessian(Tao tao,Vec X,Mat H,Mat Hpre)250d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre)
251d71ae5a4SJacob Faibussowitsch {
252a7e14dcfSSatish Balay PetscFunctionBegin;
253441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
254a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
255a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
256794dad8cSStefano Zampini PetscCheck(tao->ops->computehessian, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TaoSetHessian() not called");
257794dad8cSStefano Zampini
258a7e14dcfSSatish Balay ++tao->nhess;
2599566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
2609566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_HessianEval, tao, X, H, Hpre));
261792fecdfSBarry Smith PetscCallBack("Tao callback Hessian", (*tao->ops->computehessian)(tao, X, H, Hpre, tao->user_hessP));
2629566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_HessianEval, tao, X, H, Hpre));
2639566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
26409baa881SHong Zhang
2659566063dSJacob Faibussowitsch PetscCall(TaoTestHessian(tao));
2663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
267a7e14dcfSSatish Balay }
268a7e14dcfSSatish Balay
269cc4c1da9SBarry Smith /*@
270a7e14dcfSSatish Balay TaoComputeJacobian - Computes the Jacobian matrix that has been
271a7e14dcfSSatish Balay set with TaoSetJacobianRoutine().
272a7e14dcfSSatish Balay
273c3339decSBarry Smith Collective
274a7e14dcfSSatish Balay
275a7e14dcfSSatish Balay Input Parameters:
276f4c1ad5cSStefano Zampini + tao - the Tao solver context
277f4c1ad5cSStefano Zampini - X - input vector
278a7e14dcfSSatish Balay
279a7e14dcfSSatish Balay Output Parameters:
280f4c1ad5cSStefano Zampini + J - Jacobian matrix
2817addb90fSBarry Smith - Jpre - matrix used to compute the preconditioner, often the same as `J`
282a7e14dcfSSatish Balay
28347450a7bSBarry Smith Level: developer
28447450a7bSBarry Smith
285a7e14dcfSSatish Balay Notes:
286a7e14dcfSSatish Balay Most users should not need to explicitly call this routine, as it
287a7e14dcfSSatish Balay is used internally within the minimization solvers.
288a7e14dcfSSatish Balay
28965ba42b6SBarry Smith `TaoComputeJacobian()` is typically used within minimization
290a7e14dcfSSatish Balay implementations, so most users would not generally call this routine
291a7e14dcfSSatish Balay themselves.
292a7e14dcfSSatish Balay
2931cc06b55SBarry Smith .seealso: [](ch_tao), `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianRoutine()`
294a7e14dcfSSatish Balay @*/
TaoComputeJacobian(Tao tao,Vec X,Mat J,Mat Jpre)295d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeJacobian(Tao tao, Vec X, Mat J, Mat Jpre)
296d71ae5a4SJacob Faibussowitsch {
297a7e14dcfSSatish Balay PetscFunctionBegin;
298441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
299a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
300a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
301a7e14dcfSSatish Balay ++tao->njac;
3029566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
3039566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_JacobianEval, tao, X, J, Jpre));
304792fecdfSBarry Smith PetscCallBack("Tao callback Jacobian", (*tao->ops->computejacobian)(tao, X, J, Jpre, tao->user_jacP));
3059566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_JacobianEval, tao, X, J, Jpre));
3069566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
3073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
308a7e14dcfSSatish Balay }
309a7e14dcfSSatish Balay
310cc4c1da9SBarry Smith /*@
3114a48860cSAlp Dener TaoComputeResidualJacobian - Computes the least-squares residual Jacobian matrix that has been
31265ba42b6SBarry Smith set with `TaoSetJacobianResidual()`.
3134a48860cSAlp Dener
314c3339decSBarry Smith Collective
3154a48860cSAlp Dener
3164a48860cSAlp Dener Input Parameters:
3174a48860cSAlp Dener + tao - the Tao solver context
3184a48860cSAlp Dener - X - input vector
3194a48860cSAlp Dener
3204a48860cSAlp Dener Output Parameters:
3214a48860cSAlp Dener + J - Jacobian matrix
3227addb90fSBarry Smith - Jpre - matrix used to compute the preconditioner, often the same as `J`
3234a48860cSAlp Dener
32447450a7bSBarry Smith Level: developer
32547450a7bSBarry Smith
3264a48860cSAlp Dener Notes:
3274a48860cSAlp Dener Most users should not need to explicitly call this routine, as it
3284a48860cSAlp Dener is used internally within the minimization solvers.
3294a48860cSAlp Dener
33065ba42b6SBarry Smith `TaoComputeResidualJacobian()` is typically used within least-squares
3314a48860cSAlp Dener implementations, so most users would not generally call this routine
3324a48860cSAlp Dener themselves.
3334a48860cSAlp Dener
3341cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoComputeResidual()`, `TaoSetJacobianResidual()`
3354a48860cSAlp Dener @*/
TaoComputeResidualJacobian(Tao tao,Vec X,Mat J,Mat Jpre)336d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeResidualJacobian(Tao tao, Vec X, Mat J, Mat Jpre)
337d71ae5a4SJacob Faibussowitsch {
3384a48860cSAlp Dener PetscFunctionBegin;
3394a48860cSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
3404a48860cSAlp Dener PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
3414a48860cSAlp Dener PetscCheckSameComm(tao, 1, X, 2);
3424a48860cSAlp Dener ++tao->njac;
3439566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
3449566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_JacobianEval, tao, X, J, Jpre));
345792fecdfSBarry Smith PetscCallBack("Tao callback least-squares residual Jacobian", (*tao->ops->computeresidualjacobian)(tao, X, J, Jpre, tao->user_lsjacP));
3469566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_JacobianEval, tao, X, J, Jpre));
3479566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3494a48860cSAlp Dener }
3504a48860cSAlp Dener
351cc4c1da9SBarry Smith /*@
352a7e14dcfSSatish Balay TaoComputeJacobianState - Computes the Jacobian matrix that has been
35365ba42b6SBarry Smith set with `TaoSetJacobianStateRoutine()`.
354a7e14dcfSSatish Balay
355c3339decSBarry Smith Collective
356a7e14dcfSSatish Balay
357a7e14dcfSSatish Balay Input Parameters:
35847450a7bSBarry Smith + tao - the `Tao` solver context
359f4c1ad5cSStefano Zampini - X - input vector
360a7e14dcfSSatish Balay
361a7e14dcfSSatish Balay Output Parameters:
3626b867d5aSJose E. Roman + J - Jacobian matrix
36347450a7bSBarry Smith . Jpre - matrix used to construct the preconditioner, often the same as `J`
36465ba42b6SBarry Smith - Jinv - unknown
365a7e14dcfSSatish Balay
366a7e14dcfSSatish Balay Level: developer
367a7e14dcfSSatish Balay
36847450a7bSBarry Smith Note:
36947450a7bSBarry Smith Most users should not need to explicitly call this routine, as it
37047450a7bSBarry Smith is used internally within the optimization algorithms.
37147450a7bSBarry Smith
3721cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
373a7e14dcfSSatish Balay @*/
TaoComputeJacobianState(Tao tao,Vec X,Mat J,Mat Jpre,Mat Jinv)374d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv)
375d71ae5a4SJacob Faibussowitsch {
376a7e14dcfSSatish Balay PetscFunctionBegin;
377441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
378a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
379a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
380a7e14dcfSSatish Balay ++tao->njac_state;
3819566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
3829566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_JacobianEval, tao, X, J, Jpre));
383792fecdfSBarry Smith PetscCallBack("Tao callback Jacobian(state)", (*tao->ops->computejacobianstate)(tao, X, J, Jpre, Jinv, tao->user_jac_stateP));
3849566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_JacobianEval, tao, X, J, Jpre));
3859566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
3863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
387a7e14dcfSSatish Balay }
388a7e14dcfSSatish Balay
389cc4c1da9SBarry Smith /*@
390a7e14dcfSSatish Balay TaoComputeJacobianDesign - Computes the Jacobian matrix that has been
39165ba42b6SBarry Smith set with `TaoSetJacobianDesignRoutine()`.
392a7e14dcfSSatish Balay
393c3339decSBarry Smith Collective
394a7e14dcfSSatish Balay
395a7e14dcfSSatish Balay Input Parameters:
396f4c1ad5cSStefano Zampini + tao - the Tao solver context
397f4c1ad5cSStefano Zampini - X - input vector
398a7e14dcfSSatish Balay
3992fe279fdSBarry Smith Output Parameter:
400f4c1ad5cSStefano Zampini . J - Jacobian matrix
401a7e14dcfSSatish Balay
40247450a7bSBarry Smith Level: developer
40347450a7bSBarry Smith
40447450a7bSBarry Smith Note:
405a7e14dcfSSatish Balay Most users should not need to explicitly call this routine, as it
40665ba42b6SBarry Smith is used internally within the optimization algorithms.
407a7e14dcfSSatish Balay
40842747ad1SJacob Faibussowitsch .seealso: [](ch_tao), `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianDesignRoutine()`, `TaoSetStateDesignIS()`
409a7e14dcfSSatish Balay @*/
TaoComputeJacobianDesign(Tao tao,Vec X,Mat J)410d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J)
411d71ae5a4SJacob Faibussowitsch {
412a7e14dcfSSatish Balay PetscFunctionBegin;
413441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
414a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
415a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
416a7e14dcfSSatish Balay ++tao->njac_design;
4179566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
4189566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_JacobianEval, tao, X, J, NULL));
419792fecdfSBarry Smith PetscCallBack("Tao callback Jacobian(design)", (*tao->ops->computejacobiandesign)(tao, X, J, tao->user_jac_designP));
4209566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_JacobianEval, tao, X, J, NULL));
4219566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
4223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
423a7e14dcfSSatish Balay }
424a7e14dcfSSatish Balay
425a7e14dcfSSatish Balay /*@C
426a7e14dcfSSatish Balay TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix.
427a7e14dcfSSatish Balay
42820f4b53cSBarry Smith Logically Collective
429a7e14dcfSSatish Balay
430a7e14dcfSSatish Balay Input Parameters:
43147450a7bSBarry Smith + tao - the `Tao` context
43247450a7bSBarry Smith . J - Matrix used for the Jacobian
43347450a7bSBarry Smith . Jpre - Matrix that will be used to construct the preconditioner, can be same as `J`
434f4c1ad5cSStefano Zampini . func - Jacobian evaluation routine
435a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the
43647450a7bSBarry Smith Jacobian evaluation routine (may be `NULL`)
437a7e14dcfSSatish Balay
43820f4b53cSBarry Smith Calling sequence of `func`:
43947450a7bSBarry Smith + tao - the `Tao` context
440a7e14dcfSSatish Balay . x - input vector
441a7e14dcfSSatish Balay . J - Jacobian matrix
44247450a7bSBarry Smith . Jpre - matrix used to construct the preconditioner, usually the same as `J`
443a7e14dcfSSatish Balay - ctx - [optional] user-defined Jacobian context
444a7e14dcfSSatish Balay
445a7e14dcfSSatish Balay Level: intermediate
44665ba42b6SBarry Smith
4471cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetGradient()`, `TaoSetObjective()`
448a7e14dcfSSatish Balay @*/
TaoSetJacobianRoutine(Tao tao,Mat J,Mat Jpre,PetscErrorCode (* func)(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx),PetscCtx ctx)449*2a8381b2SBarry Smith PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx), PetscCtx ctx)
450d71ae5a4SJacob Faibussowitsch {
451a7e14dcfSSatish Balay PetscFunctionBegin;
452441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
453a7e14dcfSSatish Balay if (J) {
454a7e14dcfSSatish Balay PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
455a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, J, 2);
456a7e14dcfSSatish Balay }
457a7e14dcfSSatish Balay if (Jpre) {
458a7e14dcfSSatish Balay PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
459a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, Jpre, 3);
460a7e14dcfSSatish Balay }
461ad540459SPierre Jolivet if (ctx) tao->user_jacP = ctx;
462ad540459SPierre Jolivet if (func) tao->ops->computejacobian = func;
463a7e14dcfSSatish Balay if (J) {
4649566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J));
4659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian));
466a7e14dcfSSatish Balay tao->jacobian = J;
467a7e14dcfSSatish Balay }
468a7e14dcfSSatish Balay if (Jpre) {
4699566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jpre));
4709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian_pre));
471a7e14dcfSSatish Balay tao->jacobian_pre = Jpre;
472a7e14dcfSSatish Balay }
4733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
474a7e14dcfSSatish Balay }
475a7e14dcfSSatish Balay
476a7e14dcfSSatish Balay /*@C
4774ffbe8acSAlp Dener TaoSetJacobianResidualRoutine - Sets the function to compute the least-squares residual Jacobian as well as the
4784a48860cSAlp Dener location to store the matrix.
4794a48860cSAlp Dener
48020f4b53cSBarry Smith Logically Collective
4814a48860cSAlp Dener
4824a48860cSAlp Dener Input Parameters:
48347450a7bSBarry Smith + tao - the `Tao` context
4844a48860cSAlp Dener . J - Matrix used for the jacobian
48547450a7bSBarry Smith . Jpre - Matrix that will be used to construct the preconditioner, can be same as `J`
4864a48860cSAlp Dener . func - Jacobian evaluation routine
4874a48860cSAlp Dener - ctx - [optional] user-defined context for private data for the
48847450a7bSBarry Smith Jacobian evaluation routine (may be `NULL`)
4894a48860cSAlp Dener
49020f4b53cSBarry Smith Calling sequence of `func`:
49147450a7bSBarry Smith + tao - the `Tao` context
4924a48860cSAlp Dener . x - input vector
4934a48860cSAlp Dener . J - Jacobian matrix
49447450a7bSBarry Smith . Jpre - matrix used to construct the preconditioner, usually the same as `J`
4954a48860cSAlp Dener - ctx - [optional] user-defined Jacobian context
4964a48860cSAlp Dener
4974a48860cSAlp Dener Level: intermediate
49865ba42b6SBarry Smith
4991cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetGradient()`, `TaoSetObjective()`
5004a48860cSAlp Dener @*/
TaoSetJacobianResidualRoutine(Tao tao,Mat J,Mat Jpre,PetscErrorCode (* func)(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx),PetscCtx ctx)501*2a8381b2SBarry Smith PetscErrorCode TaoSetJacobianResidualRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx), PetscCtx ctx)
502d71ae5a4SJacob Faibussowitsch {
5034a48860cSAlp Dener PetscFunctionBegin;
5044a48860cSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
5054a48860cSAlp Dener if (J) {
5064a48860cSAlp Dener PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
5074a48860cSAlp Dener PetscCheckSameComm(tao, 1, J, 2);
5084a48860cSAlp Dener }
5094a48860cSAlp Dener if (Jpre) {
5104a48860cSAlp Dener PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
5114a48860cSAlp Dener PetscCheckSameComm(tao, 1, Jpre, 3);
5124a48860cSAlp Dener }
513ad540459SPierre Jolivet if (ctx) tao->user_lsjacP = ctx;
514ad540459SPierre Jolivet if (func) tao->ops->computeresidualjacobian = func;
5154a48860cSAlp Dener if (J) {
5169566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J));
5179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->ls_jac));
5184a48860cSAlp Dener tao->ls_jac = J;
5194a48860cSAlp Dener }
5204a48860cSAlp Dener if (Jpre) {
5219566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jpre));
5229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->ls_jac_pre));
5234a48860cSAlp Dener tao->ls_jac_pre = Jpre;
5244a48860cSAlp Dener }
5253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5264a48860cSAlp Dener }
5274a48860cSAlp Dener
5284a48860cSAlp Dener /*@C
529a7e14dcfSSatish Balay TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian
530a7e14dcfSSatish Balay (and its inverse) of the constraint function with respect to the state variables.
53165ba42b6SBarry Smith Used only for PDE-constrained optimization.
532a7e14dcfSSatish Balay
53320f4b53cSBarry Smith Logically Collective
534a7e14dcfSSatish Balay
535a7e14dcfSSatish Balay Input Parameters:
53647450a7bSBarry Smith + tao - the `Tao` context
53747450a7bSBarry Smith . J - Matrix used for the Jacobian
53847450a7bSBarry Smith . Jpre - Matrix that will be used to construct the preconditioner, can be same as `J`. Only used if `Jinv` is `NULL`
53947450a7bSBarry Smith . Jinv - [optional] Matrix used to apply the inverse of the state Jacobian. Use `NULL` to default to PETSc `KSP` solvers to apply the inverse.
540f4c1ad5cSStefano Zampini . func - Jacobian evaluation routine
541a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the
54247450a7bSBarry Smith Jacobian evaluation routine (may be `NULL`)
543a7e14dcfSSatish Balay
54420f4b53cSBarry Smith Calling sequence of `func`:
54547450a7bSBarry Smith + tao - the `Tao` context
546a7e14dcfSSatish Balay . x - input vector
547a7e14dcfSSatish Balay . J - Jacobian matrix
54847450a7bSBarry Smith . Jpre - matrix used to construct the preconditioner, usually the same as `J`
54920f4b53cSBarry Smith . Jinv - inverse of `J`
550a7e14dcfSSatish Balay - ctx - [optional] user-defined Jacobian context
551a7e14dcfSSatish Balay
552a7e14dcfSSatish Balay Level: intermediate
55365ba42b6SBarry Smith
5541cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoComputeJacobianState()`, `TaoSetJacobianDesignRoutine()`, `TaoSetStateDesignIS()`
555a7e14dcfSSatish Balay @*/
TaoSetJacobianStateRoutine(Tao tao,Mat J,Mat Jpre,Mat Jinv,PetscErrorCode (* func)(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,PetscCtx ctx),PetscCtx ctx)556*2a8381b2SBarry Smith PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao tao, Vec x, Mat J, Mat Jpre, Mat Jinv, PetscCtx ctx), PetscCtx ctx)
557d71ae5a4SJacob Faibussowitsch {
558a7e14dcfSSatish Balay PetscFunctionBegin;
559441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
560a7e14dcfSSatish Balay if (J) {
561a7e14dcfSSatish Balay PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
562a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, J, 2);
563a7e14dcfSSatish Balay }
564a7e14dcfSSatish Balay if (Jpre) {
565a7e14dcfSSatish Balay PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
566a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, Jpre, 3);
567a7e14dcfSSatish Balay }
568a7e14dcfSSatish Balay if (Jinv) {
569a7e14dcfSSatish Balay PetscValidHeaderSpecific(Jinv, MAT_CLASSID, 4);
570a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, Jinv, 4);
571a7e14dcfSSatish Balay }
572ad540459SPierre Jolivet if (ctx) tao->user_jac_stateP = ctx;
573ad540459SPierre Jolivet if (func) tao->ops->computejacobianstate = func;
574a7e14dcfSSatish Balay if (J) {
5759566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J));
5769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian_state));
577a7e14dcfSSatish Balay tao->jacobian_state = J;
578a7e14dcfSSatish Balay }
579a7e14dcfSSatish Balay if (Jpre) {
5809566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jpre));
5819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian_state_pre));
582a7e14dcfSSatish Balay tao->jacobian_state_pre = Jpre;
583a7e14dcfSSatish Balay }
584a7e14dcfSSatish Balay if (Jinv) {
5859566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jinv));
5869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian_state_inv));
587a7e14dcfSSatish Balay tao->jacobian_state_inv = Jinv;
588a7e14dcfSSatish Balay }
5893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
590a7e14dcfSSatish Balay }
591a7e14dcfSSatish Balay
592a7e14dcfSSatish Balay /*@C
593a7e14dcfSSatish Balay TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of
594a7e14dcfSSatish Balay the constraint function with respect to the design variables. Used only for
59565ba42b6SBarry Smith PDE-constrained optimization.
596a7e14dcfSSatish Balay
59720f4b53cSBarry Smith Logically Collective
598a7e14dcfSSatish Balay
599a7e14dcfSSatish Balay Input Parameters:
60047450a7bSBarry Smith + tao - the `Tao` context
60147450a7bSBarry Smith . J - Matrix used for the Jacobian
602f4c1ad5cSStefano Zampini . func - Jacobian evaluation routine
603a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the
60447450a7bSBarry Smith Jacobian evaluation routine (may be `NULL`)
605a7e14dcfSSatish Balay
60620f4b53cSBarry Smith Calling sequence of `func`:
60747450a7bSBarry Smith + tao - the `Tao` context
608a7e14dcfSSatish Balay . x - input vector
609a7e14dcfSSatish Balay . J - Jacobian matrix
610a7e14dcfSSatish Balay - ctx - [optional] user-defined Jacobian context
611a7e14dcfSSatish Balay
612a7e14dcfSSatish Balay Level: intermediate
613f4c1ad5cSStefano Zampini
6141cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoComputeJacobianDesign()`, `TaoSetJacobianStateRoutine()`, `TaoSetStateDesignIS()`
615a7e14dcfSSatish Balay @*/
TaoSetJacobianDesignRoutine(Tao tao,Mat J,PetscErrorCode (* func)(Tao tao,Vec x,Mat J,PetscCtx ctx),PetscCtx ctx)616*2a8381b2SBarry Smith PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao tao, Vec x, Mat J, PetscCtx ctx), PetscCtx ctx)
617d71ae5a4SJacob Faibussowitsch {
618a7e14dcfSSatish Balay PetscFunctionBegin;
619441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
620a7e14dcfSSatish Balay if (J) {
621a7e14dcfSSatish Balay PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
622a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, J, 2);
623a7e14dcfSSatish Balay }
624ad540459SPierre Jolivet if (ctx) tao->user_jac_designP = ctx;
625ad540459SPierre Jolivet if (func) tao->ops->computejacobiandesign = func;
626a7e14dcfSSatish Balay if (J) {
6279566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J));
6289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian_design));
629a7e14dcfSSatish Balay tao->jacobian_design = J;
630a7e14dcfSSatish Balay }
6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
632a7e14dcfSSatish Balay }
633a7e14dcfSSatish Balay
634a7e14dcfSSatish Balay /*@
63547450a7bSBarry Smith TaoSetStateDesignIS - Indicate to the `Tao` object which variables in the
636a7e14dcfSSatish Balay solution vector are state variables and which are design. Only applies to
63765ba42b6SBarry Smith PDE-constrained optimization.
638a7e14dcfSSatish Balay
639c3339decSBarry Smith Logically Collective
640a7e14dcfSSatish Balay
641a7e14dcfSSatish Balay Input Parameters:
64247450a7bSBarry Smith + tao - The `Tao` context
643a7e14dcfSSatish Balay . s_is - the index set corresponding to the state variables
644a7e14dcfSSatish Balay - d_is - the index set corresponding to the design variables
645a7e14dcfSSatish Balay
646a7e14dcfSSatish Balay Level: intermediate
647a7e14dcfSSatish Balay
6481cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetJacobianStateRoutine()`, `TaoSetJacobianDesignRoutine()`
649a7e14dcfSSatish Balay @*/
TaoSetStateDesignIS(Tao tao,IS s_is,IS d_is)650d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is)
651d71ae5a4SJacob Faibussowitsch {
65245cf516eSBarry Smith PetscFunctionBegin;
6539566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)s_is));
6549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tao->state_is));
655a7e14dcfSSatish Balay tao->state_is = s_is;
656f4f49eeaSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)d_is));
6579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tao->design_is));
658a7e14dcfSSatish Balay tao->design_is = d_is;
6593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
660a7e14dcfSSatish Balay }
661a7e14dcfSSatish Balay
662cc4c1da9SBarry Smith /*@
663a7e14dcfSSatish Balay TaoComputeJacobianEquality - Computes the Jacobian matrix that has been
66465ba42b6SBarry Smith set with `TaoSetJacobianEqualityRoutine()`.
665a7e14dcfSSatish Balay
666c3339decSBarry Smith Collective
667a7e14dcfSSatish Balay
668a7e14dcfSSatish Balay Input Parameters:
66947450a7bSBarry Smith + tao - the `Tao` solver context
670f4c1ad5cSStefano Zampini - X - input vector
671a7e14dcfSSatish Balay
672a7e14dcfSSatish Balay Output Parameters:
673f4c1ad5cSStefano Zampini + J - Jacobian matrix
67447450a7bSBarry Smith - Jpre - matrix used to construct the preconditioner, often the same as `J`
67547450a7bSBarry Smith
67647450a7bSBarry Smith Level: developer
677a7e14dcfSSatish Balay
678a7e14dcfSSatish Balay Notes:
679a7e14dcfSSatish Balay Most users should not need to explicitly call this routine, as it
68065ba42b6SBarry Smith is used internally within the optimization algorithms.
681a7e14dcfSSatish Balay
6821cc06b55SBarry Smith .seealso: [](ch_tao), `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
683a7e14dcfSSatish Balay @*/
TaoComputeJacobianEquality(Tao tao,Vec X,Mat J,Mat Jpre)684d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre)
685d71ae5a4SJacob Faibussowitsch {
686a7e14dcfSSatish Balay PetscFunctionBegin;
687441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
688a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
689a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
690a7e14dcfSSatish Balay ++tao->njac_equality;
6919566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
6929566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_JacobianEval, tao, X, J, Jpre));
693792fecdfSBarry Smith PetscCallBack("Tao callback Jacobian(equality)", (*tao->ops->computejacobianequality)(tao, X, J, Jpre, tao->user_jac_equalityP));
6949566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_JacobianEval, tao, X, J, Jpre));
6959566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
6963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
697a7e14dcfSSatish Balay }
698a7e14dcfSSatish Balay
699cc4c1da9SBarry Smith /*@
700a7e14dcfSSatish Balay TaoComputeJacobianInequality - Computes the Jacobian matrix that has been
70165ba42b6SBarry Smith set with `TaoSetJacobianInequalityRoutine()`.
702a7e14dcfSSatish Balay
703c3339decSBarry Smith Collective
704a7e14dcfSSatish Balay
705a7e14dcfSSatish Balay Input Parameters:
70647450a7bSBarry Smith + tao - the `Tao` solver context
707f4c1ad5cSStefano Zampini - X - input vector
708a7e14dcfSSatish Balay
709a7e14dcfSSatish Balay Output Parameters:
710f4c1ad5cSStefano Zampini + J - Jacobian matrix
71147450a7bSBarry Smith - Jpre - matrix used to construct the preconditioner
712a7e14dcfSSatish Balay
713a7e14dcfSSatish Balay Level: developer
714a7e14dcfSSatish Balay
71547450a7bSBarry Smith Note:
71647450a7bSBarry Smith Most users should not need to explicitly call this routine, as it
71747450a7bSBarry Smith is used internally within the minimization solvers.
71847450a7bSBarry Smith
7191cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
720a7e14dcfSSatish Balay @*/
TaoComputeJacobianInequality(Tao tao,Vec X,Mat J,Mat Jpre)721d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre)
722d71ae5a4SJacob Faibussowitsch {
723a7e14dcfSSatish Balay PetscFunctionBegin;
724441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
725a7e14dcfSSatish Balay PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
726a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, X, 2);
727a7e14dcfSSatish Balay ++tao->njac_inequality;
7289566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X));
7299566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TAO_JacobianEval, tao, X, J, Jpre));
730792fecdfSBarry Smith PetscCallBack("Tao callback Jacobian (inequality)", (*tao->ops->computejacobianinequality)(tao, X, J, Jpre, tao->user_jac_inequalityP));
7319566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TAO_JacobianEval, tao, X, J, Jpre));
7329566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X));
7333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
734a7e14dcfSSatish Balay }
735a7e14dcfSSatish Balay
736a7e14dcfSSatish Balay /*@C
737a7e14dcfSSatish Balay TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian
738a7e14dcfSSatish Balay (and its inverse) of the constraint function with respect to the equality variables.
73965ba42b6SBarry Smith Used only for PDE-constrained optimization.
740a7e14dcfSSatish Balay
74120f4b53cSBarry Smith Logically Collective
742a7e14dcfSSatish Balay
743a7e14dcfSSatish Balay Input Parameters:
74447450a7bSBarry Smith + tao - the `Tao` context
74547450a7bSBarry Smith . J - Matrix used for the Jacobian
74647450a7bSBarry Smith . Jpre - Matrix that will be used to construct the preconditioner, can be same as `J`.
747f4c1ad5cSStefano Zampini . func - Jacobian evaluation routine
748a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the
74947450a7bSBarry Smith Jacobian evaluation routine (may be `NULL`)
750a7e14dcfSSatish Balay
75120f4b53cSBarry Smith Calling sequence of `func`:
75247450a7bSBarry Smith + tao - the `Tao` context
753a7e14dcfSSatish Balay . x - input vector
754a7e14dcfSSatish Balay . J - Jacobian matrix
75547450a7bSBarry Smith . Jpre - matrix used to construct the preconditioner, usually the same as `J`
756a7e14dcfSSatish Balay - ctx - [optional] user-defined Jacobian context
757a7e14dcfSSatish Balay
758a7e14dcfSSatish Balay Level: intermediate
759f4c1ad5cSStefano Zampini
7601cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoComputeJacobianEquality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetEqualityDesignIS()`
761a7e14dcfSSatish Balay @*/
TaoSetJacobianEqualityRoutine(Tao tao,Mat J,Mat Jpre,PetscErrorCode (* func)(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx),PetscCtx ctx)762*2a8381b2SBarry Smith PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx), PetscCtx ctx)
763d71ae5a4SJacob Faibussowitsch {
764a7e14dcfSSatish Balay PetscFunctionBegin;
765441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
766a7e14dcfSSatish Balay if (J) {
767a7e14dcfSSatish Balay PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
768a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, J, 2);
769a7e14dcfSSatish Balay }
770a7e14dcfSSatish Balay if (Jpre) {
771a7e14dcfSSatish Balay PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
772a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, Jpre, 3);
773a7e14dcfSSatish Balay }
774ad540459SPierre Jolivet if (ctx) tao->user_jac_equalityP = ctx;
775ad540459SPierre Jolivet if (func) tao->ops->computejacobianequality = func;
776a7e14dcfSSatish Balay if (J) {
7779566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J));
7789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian_equality));
779a7e14dcfSSatish Balay tao->jacobian_equality = J;
780a7e14dcfSSatish Balay }
781a7e14dcfSSatish Balay if (Jpre) {
7829566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jpre));
7839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian_equality_pre));
784a7e14dcfSSatish Balay tao->jacobian_equality_pre = Jpre;
785a7e14dcfSSatish Balay }
7863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
787a7e14dcfSSatish Balay }
788a7e14dcfSSatish Balay
789a7e14dcfSSatish Balay /*@C
790c70ced77Spaul.kuehner TaoGetJacobianEqualityRoutine - Gets the function used to compute equality constraint Jacobian.
791c70ced77Spaul.kuehner
792c70ced77Spaul.kuehner Not Collective
793c70ced77Spaul.kuehner
794c70ced77Spaul.kuehner Input Parameter:
795c70ced77Spaul.kuehner . tao - the `Tao` context
796c70ced77Spaul.kuehner
797c70ced77Spaul.kuehner Output Parameters:
798c70ced77Spaul.kuehner + J - the matrix to internally hold the constraint computation
799c70ced77Spaul.kuehner . Jpre - the matrix used to construct the preconditioner
800c70ced77Spaul.kuehner . func - Jacobian evaluation routine
801c70ced77Spaul.kuehner - ctx - the (optional) user-defined context
802c70ced77Spaul.kuehner
803c70ced77Spaul.kuehner Calling sequence of `func`:
804c70ced77Spaul.kuehner + tao - the `Tao` context
805c70ced77Spaul.kuehner . x - input vector
806c70ced77Spaul.kuehner . J - Jacobian matrix
807c70ced77Spaul.kuehner . Jpre - matrix used to construct the preconditioner, usually the same as `J`
808c70ced77Spaul.kuehner - ctx - [optional] user-defined Jacobian context
809c70ced77Spaul.kuehner
810c70ced77Spaul.kuehner Level: intermediate
811c70ced77Spaul.kuehner
812c70ced77Spaul.kuehner .seealso: [](ch_tao), `Tao`, `TaoComputeJacobianEquality()`, `TaoSetJacobianEqualityRoutine()`
813c70ced77Spaul.kuehner @*/
TaoGetJacobianEqualityRoutine(Tao tao,Mat * J,Mat * Jpre,PetscErrorCode (** func)(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx),PetscCtxRt ctx)814*2a8381b2SBarry Smith PetscErrorCode TaoGetJacobianEqualityRoutine(Tao tao, Mat *J, Mat *Jpre, PetscErrorCode (**func)(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx), PetscCtxRt ctx)
815c70ced77Spaul.kuehner {
816c70ced77Spaul.kuehner PetscFunctionBegin;
817c70ced77Spaul.kuehner PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
818c70ced77Spaul.kuehner if (J) *J = tao->jacobian_equality;
819c70ced77Spaul.kuehner if (Jpre) *Jpre = tao->jacobian_equality_pre;
820c70ced77Spaul.kuehner if (func) *func = tao->ops->computejacobianequality;
821*2a8381b2SBarry Smith if (ctx) *(void **)ctx = tao->user_jac_equalityP;
822c70ced77Spaul.kuehner PetscFunctionReturn(PETSC_SUCCESS);
823c70ced77Spaul.kuehner }
824c70ced77Spaul.kuehner
825c70ced77Spaul.kuehner /*@C
826a7e14dcfSSatish Balay TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian
827a7e14dcfSSatish Balay (and its inverse) of the constraint function with respect to the inequality variables.
82865ba42b6SBarry Smith Used only for PDE-constrained optimization.
829a7e14dcfSSatish Balay
83020f4b53cSBarry Smith Logically Collective
831a7e14dcfSSatish Balay
832a7e14dcfSSatish Balay Input Parameters:
83347450a7bSBarry Smith + tao - the `Tao` context
83447450a7bSBarry Smith . J - Matrix used for the Jacobian
83547450a7bSBarry Smith . Jpre - Matrix that will be used to construct the preconditioner, can be same as `J`.
836f4c1ad5cSStefano Zampini . func - Jacobian evaluation routine
837a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the
83847450a7bSBarry Smith Jacobian evaluation routine (may be `NULL`)
839a7e14dcfSSatish Balay
84020f4b53cSBarry Smith Calling sequence of `func`:
84147450a7bSBarry Smith + tao - the `Tao` context
842a7e14dcfSSatish Balay . x - input vector
843a7e14dcfSSatish Balay . J - Jacobian matrix
84447450a7bSBarry Smith . Jpre - matrix used to construct the preconditioner, usually the same as `J`
845a7e14dcfSSatish Balay - ctx - [optional] user-defined Jacobian context
846a7e14dcfSSatish Balay
847a7e14dcfSSatish Balay Level: intermediate
848f4c1ad5cSStefano Zampini
8491cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoComputeJacobianInequality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetInequalityDesignIS()`
850a7e14dcfSSatish Balay @*/
TaoSetJacobianInequalityRoutine(Tao tao,Mat J,Mat Jpre,PetscErrorCode (* func)(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx),PetscCtx ctx)851*2a8381b2SBarry Smith PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx), PetscCtx ctx)
852d71ae5a4SJacob Faibussowitsch {
853a7e14dcfSSatish Balay PetscFunctionBegin;
854441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
855a7e14dcfSSatish Balay if (J) {
856a7e14dcfSSatish Balay PetscValidHeaderSpecific(J, MAT_CLASSID, 2);
857a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, J, 2);
858a7e14dcfSSatish Balay }
859a7e14dcfSSatish Balay if (Jpre) {
860a7e14dcfSSatish Balay PetscValidHeaderSpecific(Jpre, MAT_CLASSID, 3);
861a7e14dcfSSatish Balay PetscCheckSameComm(tao, 1, Jpre, 3);
862a7e14dcfSSatish Balay }
863ad540459SPierre Jolivet if (ctx) tao->user_jac_inequalityP = ctx;
864ad540459SPierre Jolivet if (func) tao->ops->computejacobianinequality = func;
865a7e14dcfSSatish Balay if (J) {
8669566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)J));
8679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian_inequality));
868a7e14dcfSSatish Balay tao->jacobian_inequality = J;
869a7e14dcfSSatish Balay }
870a7e14dcfSSatish Balay if (Jpre) {
8719566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Jpre));
8729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tao->jacobian_inequality_pre));
873a7e14dcfSSatish Balay tao->jacobian_inequality_pre = Jpre;
874a7e14dcfSSatish Balay }
8753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
876a7e14dcfSSatish Balay }
877d625fd55Spaul.kuehner
878d625fd55Spaul.kuehner /*@C
879d625fd55Spaul.kuehner TaoGetJacobianInequalityRoutine - Gets the function used to compute inequality constraint Jacobian.
880d625fd55Spaul.kuehner
881d625fd55Spaul.kuehner Not Collective
882d625fd55Spaul.kuehner
883d625fd55Spaul.kuehner Input Parameter:
884d625fd55Spaul.kuehner . tao - the `Tao` context
885d625fd55Spaul.kuehner
886d625fd55Spaul.kuehner Output Parameters:
887d625fd55Spaul.kuehner + J - the matrix to internally hold the constraint computation
888d625fd55Spaul.kuehner . Jpre - the matrix used to construct the preconditioner
889d625fd55Spaul.kuehner . func - Jacobian evaluation routine
890d625fd55Spaul.kuehner - ctx - the (optional) user-defined context
891d625fd55Spaul.kuehner
892d625fd55Spaul.kuehner Calling sequence of `func`:
893d625fd55Spaul.kuehner + tao - the `Tao` context
894d625fd55Spaul.kuehner . x - input vector
895d625fd55Spaul.kuehner . J - Jacobian matrix
896d625fd55Spaul.kuehner . Jpre - matrix used to construct the preconditioner, usually the same as `J`
897d625fd55Spaul.kuehner - ctx - [optional] user-defined Jacobian context
898d625fd55Spaul.kuehner
899d625fd55Spaul.kuehner Level: intermediate
900d625fd55Spaul.kuehner
901d625fd55Spaul.kuehner .seealso: [](ch_tao), `Tao`, `TaoComputeJacobianInequality()`, `TaoSetJacobianInequalityRoutine()`
902d625fd55Spaul.kuehner @*/
TaoGetJacobianInequalityRoutine(Tao tao,Mat * J,Mat * Jpre,PetscErrorCode (** func)(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx),PetscCtxRt ctx)903*2a8381b2SBarry Smith PetscErrorCode TaoGetJacobianInequalityRoutine(Tao tao, Mat *J, Mat *Jpre, PetscErrorCode (**func)(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx), PetscCtxRt ctx)
904d625fd55Spaul.kuehner {
905d625fd55Spaul.kuehner PetscFunctionBegin;
906d625fd55Spaul.kuehner PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
907d625fd55Spaul.kuehner if (J) *J = tao->jacobian_inequality;
908d625fd55Spaul.kuehner if (Jpre) *Jpre = tao->jacobian_inequality_pre;
909d625fd55Spaul.kuehner if (func) *func = tao->ops->computejacobianinequality;
910*2a8381b2SBarry Smith if (ctx) *(void **)ctx = tao->user_jac_inequalityP;
911d625fd55Spaul.kuehner PetscFunctionReturn(PETSC_SUCCESS);
912d625fd55Spaul.kuehner }
913