xref: /petsc/src/tao/interface/taosolver_hj.c (revision f49d1c87faf37a381abb3ef099353a5ededcc2fa)
1af0996ceSBarry Smith #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
2a7e14dcfSSatish Balay 
3a7e14dcfSSatish Balay /*@C
4a7e14dcfSSatish Balay    TaoSetHessianRoutine - Sets the function to compute the Hessian as well as the location to store the matrix.
5a7e14dcfSSatish Balay 
6441846f8SBarry Smith    Logically collective on Tao
7a7e14dcfSSatish Balay 
8a7e14dcfSSatish Balay    Input Parameters:
9441846f8SBarry Smith +  tao  - the Tao context
10a7e14dcfSSatish Balay .  H    - Matrix used for the hessian
11a7e14dcfSSatish Balay .  Hpre - Matrix that will be used operated on by preconditioner, can be same as H
12f4c1ad5cSStefano Zampini .  func - Hessian evaluation routine
13a7e14dcfSSatish Balay -  ctx  - [optional] user-defined context for private data for the
146c23d075SBarry Smith          Hessian evaluation routine (may be NULL)
15a7e14dcfSSatish Balay 
16f4c1ad5cSStefano Zampini    Calling sequence of func:
17f4c1ad5cSStefano Zampini $    func(Tao tao,Vec x,Mat H,Mat Hpre,void *ctx);
18a7e14dcfSSatish Balay 
19441846f8SBarry Smith +  tao  - the Tao  context
20a7e14dcfSSatish Balay .  x    - input vector
21a7e14dcfSSatish Balay .  H    - Hessian matrix
22a7e14dcfSSatish Balay .  Hpre - preconditioner matrix, usually the same as H
23a7e14dcfSSatish Balay -  ctx  - [optional] user-defined Hessian context
24a7e14dcfSSatish Balay 
25a7e14dcfSSatish Balay    Level: beginner
26a7e14dcfSSatish Balay @*/
27ffad9901SBarry Smith PetscErrorCode TaoSetHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
28a7e14dcfSSatish Balay {
29a7e14dcfSSatish Balay   PetscErrorCode ierr;
3045cf516eSBarry Smith 
31a7e14dcfSSatish Balay   PetscFunctionBegin;
32441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
33a7e14dcfSSatish Balay   if (H) {
34a7e14dcfSSatish Balay     PetscValidHeaderSpecific(H,MAT_CLASSID,2);
35a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,H,2);
36a7e14dcfSSatish Balay   }
37a7e14dcfSSatish Balay   if (Hpre) {
38a7e14dcfSSatish Balay     PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3);
39a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,Hpre,3);
40a7e14dcfSSatish Balay   }
41a7e14dcfSSatish Balay   if (ctx) {
42a7e14dcfSSatish Balay     tao->user_hessP = ctx;
43a7e14dcfSSatish Balay   }
44a7e14dcfSSatish Balay   if (func) {
45a7e14dcfSSatish Balay     tao->ops->computehessian = func;
46a7e14dcfSSatish Balay   }
47a7e14dcfSSatish Balay   if (H) {
48a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)H);CHKERRQ(ierr);
4945cf516eSBarry Smith     ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr);
50a7e14dcfSSatish Balay     tao->hessian = H;
51a7e14dcfSSatish Balay   }
52a7e14dcfSSatish Balay   if (Hpre) {
53a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)Hpre);CHKERRQ(ierr);
5445cf516eSBarry Smith     ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr);
55a7e14dcfSSatish Balay     tao->hessian_pre = Hpre;
56a7e14dcfSSatish Balay   }
57a7e14dcfSSatish Balay   PetscFunctionReturn(0);
58a7e14dcfSSatish Balay }
59a7e14dcfSSatish Balay 
6009baa881SHong Zhang PetscErrorCode TaoTestHessian(Tao tao)
6109baa881SHong Zhang {
6209baa881SHong Zhang   Mat               A,B,C,D,hessian;
6309baa881SHong Zhang   Vec               x = tao->solution;
6409baa881SHong Zhang   PetscErrorCode    ierr;
6509baa881SHong Zhang   PetscReal         nrm,gnorm;
6609baa881SHong Zhang   PetscReal         threshold = 1.e-5;
6709baa881SHong Zhang   PetscInt          m,n,M,N;
6809baa881SHong Zhang   PetscBool         complete_print = PETSC_FALSE,test = PETSC_FALSE,flg;
69*f49d1c87SHong Zhang   PetscViewer       viewer,mviewer;
7009baa881SHong Zhang   MPI_Comm          comm;
7109baa881SHong Zhang   PetscInt          tabs;
7209baa881SHong Zhang   static PetscBool  directionsprinted = PETSC_FALSE;
73*f49d1c87SHong Zhang   PetscViewerFormat format;
7409baa881SHong Zhang 
7509baa881SHong Zhang   PetscFunctionBegin;
7609baa881SHong Zhang   ierr = PetscObjectOptionsBegin((PetscObject)tao);CHKERRQ(ierr);
7709baa881SHong Zhang   ierr = PetscOptionsName("-tao_test_hessian","Compare hand-coded and finite difference Hessians","None",&test);CHKERRQ(ierr);
7809baa881SHong Zhang   ierr = PetscOptionsReal("-tao_test_hessian", "Threshold for element difference between hand-coded and finite difference being meaningful","None",threshold,&threshold,NULL);CHKERRQ(ierr);
79*f49d1c87SHong Zhang   ierr = PetscOptionsViewer("-tao_test_hessian_view","View difference between hand-coded and finite difference Hessians element entries","None",&mviewer,&format,&complete_print);CHKERRQ(ierr);
8009baa881SHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
8109baa881SHong Zhang   if (!test) PetscFunctionReturn(0);
8209baa881SHong Zhang 
8309baa881SHong Zhang   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
8409baa881SHong Zhang   ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr);
8509baa881SHong Zhang   ierr = PetscViewerASCIIGetTab(viewer, &tabs);CHKERRQ(ierr);
8609baa881SHong Zhang   ierr = PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel);CHKERRQ(ierr);
8709baa881SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ---------- Testing Hessian -------------\n");CHKERRQ(ierr);
8809baa881SHong Zhang   if (!complete_print && !directionsprinted) {
8909baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Run with -tao_test_hessian_view and optionally -tao_test_hessian <threshold> to show difference\n");CHKERRQ(ierr);
9009baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    of hand-coded and finite difference Hessian entries greater than <threshold>.\n");CHKERRQ(ierr);
9109baa881SHong Zhang   }
9209baa881SHong Zhang   if (!directionsprinted) {
9309baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Testing hand-coded Hessian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");CHKERRQ(ierr);
9409baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    O(1.e-8), the hand-coded Hessian is probably correct.\n");CHKERRQ(ierr);
9509baa881SHong Zhang     directionsprinted = PETSC_TRUE;
9609baa881SHong Zhang   }
97*f49d1c87SHong Zhang   if (complete_print) {
98*f49d1c87SHong Zhang     ierr = PetscViewerPushFormat(mviewer,format);CHKERRQ(ierr);
99*f49d1c87SHong Zhang   }
10009baa881SHong Zhang 
10109baa881SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)tao->hessian,MATMFFD,&flg);CHKERRQ(ierr);
10209baa881SHong Zhang   if (!flg) hessian = tao->hessian;
10309baa881SHong Zhang   else hessian = tao->hessian_pre;
10409baa881SHong Zhang 
10509baa881SHong Zhang   while (hessian) {
10609baa881SHong Zhang     ierr = PetscObjectTypeCompareAny((PetscObject)hessian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");CHKERRQ(ierr);
10709baa881SHong Zhang     if (flg) {
10809baa881SHong Zhang       A    = hessian;
10909baa881SHong Zhang       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
11009baa881SHong Zhang     } else {
11109baa881SHong Zhang       ierr = MatComputeExplicitOperator(hessian,&A);CHKERRQ(ierr);
11209baa881SHong Zhang     }
11309baa881SHong Zhang 
11409baa881SHong Zhang     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
11509baa881SHong Zhang     ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
11609baa881SHong Zhang     ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
11709baa881SHong Zhang     ierr = MatSetSizes(B,m,n,M,N);CHKERRQ(ierr);
11809baa881SHong Zhang     ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
11909baa881SHong Zhang     ierr = MatSetUp(B);CHKERRQ(ierr);
12009baa881SHong Zhang     ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
12109baa881SHong Zhang 
12209baa881SHong Zhang     ierr = TaoDefaultComputeHessian(tao,x,B,B,NULL);CHKERRQ(ierr);
12309baa881SHong Zhang 
12409baa881SHong Zhang     ierr = MatDuplicate(B,MAT_COPY_VALUES,&D);CHKERRQ(ierr);
12509baa881SHong Zhang     ierr = MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
12609baa881SHong Zhang     ierr = MatNorm(D,NORM_FROBENIUS,&nrm);CHKERRQ(ierr);
12709baa881SHong Zhang     ierr = MatNorm(A,NORM_FROBENIUS,&gnorm);CHKERRQ(ierr);
12809baa881SHong Zhang     ierr = MatDestroy(&D);CHKERRQ(ierr);
12909baa881SHong Zhang     if (!gnorm) gnorm = 1; /* just in case */
13009baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ||H - Hfd||_F/||H||_F = %g, ||H - Hfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);CHKERRQ(ierr);
13109baa881SHong Zhang 
13209baa881SHong Zhang     if (complete_print) {
13309baa881SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Hand-coded Hessian ----------\n");CHKERRQ(ierr);
134*f49d1c87SHong Zhang       ierr = MatView(hessian,mviewer);CHKERRQ(ierr);
13509baa881SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Finite difference Hessian ----------\n");CHKERRQ(ierr);
136*f49d1c87SHong Zhang       ierr = MatView(B,mviewer);CHKERRQ(ierr);
13709baa881SHong Zhang     }
13809baa881SHong Zhang 
13909baa881SHong Zhang     if (complete_print) {
14009baa881SHong Zhang       PetscInt          Istart, Iend, *ccols, bncols, cncols, j, row;
14109baa881SHong Zhang       PetscScalar       *cvals;
14209baa881SHong Zhang       const PetscInt    *bcols;
14309baa881SHong Zhang       const PetscScalar *bvals;
14409baa881SHong Zhang 
14509baa881SHong Zhang       ierr = MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
14609baa881SHong Zhang       ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
14709baa881SHong Zhang       ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
14809baa881SHong Zhang       ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
14909baa881SHong Zhang       ierr = MatSetUp(C);CHKERRQ(ierr);
15009baa881SHong Zhang       ierr = MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
15109baa881SHong Zhang       ierr = MatGetOwnershipRange(B,&Istart,&Iend);CHKERRQ(ierr);
15209baa881SHong Zhang 
15309baa881SHong Zhang       for (row = Istart; row < Iend; row++) {
15409baa881SHong Zhang         ierr = MatGetRow(B,row,&bncols,&bcols,&bvals);CHKERRQ(ierr);
15509baa881SHong Zhang         ierr = PetscMalloc2(bncols,&ccols,bncols,&cvals);CHKERRQ(ierr);
15609baa881SHong Zhang         for (j = 0, cncols = 0; j < bncols; j++) {
15709baa881SHong Zhang           if (PetscAbsScalar(bvals[j]) > threshold) {
15809baa881SHong Zhang             ccols[cncols] = bcols[j];
15909baa881SHong Zhang             cvals[cncols] = bvals[j];
16009baa881SHong Zhang             cncols += 1;
16109baa881SHong Zhang           }
16209baa881SHong Zhang         }
16309baa881SHong Zhang         if (cncols) {
16409baa881SHong Zhang           ierr = MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);CHKERRQ(ierr);
16509baa881SHong Zhang         }
16609baa881SHong Zhang         ierr = MatRestoreRow(B,row,&bncols,&bcols,&bvals);CHKERRQ(ierr);
16709baa881SHong Zhang         ierr = PetscFree2(ccols,cvals);CHKERRQ(ierr);
16809baa881SHong Zhang       }
16909baa881SHong Zhang       ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17009baa881SHong Zhang       ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17109baa881SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Finite-difference minus hand-coded Hessian with tolerance %g ----------\n",(double)threshold);CHKERRQ(ierr);
172*f49d1c87SHong Zhang       ierr = MatView(C,mviewer);CHKERRQ(ierr);
17309baa881SHong Zhang       ierr = MatDestroy(&C);CHKERRQ(ierr);
17409baa881SHong Zhang     }
17509baa881SHong Zhang     ierr = MatDestroy(&A);CHKERRQ(ierr);
17609baa881SHong Zhang     ierr = MatDestroy(&B);CHKERRQ(ierr);
17709baa881SHong Zhang 
17809baa881SHong Zhang     if (hessian != tao->hessian_pre) {
17909baa881SHong Zhang       hessian = tao->hessian_pre;
18009baa881SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ---------- Testing Hessian for preconditioner -------------\n");CHKERRQ(ierr);
18109baa881SHong Zhang     } else hessian = NULL;
18209baa881SHong Zhang   }
183*f49d1c87SHong Zhang   if (complete_print) {
184*f49d1c87SHong Zhang     ierr = PetscViewerPopFormat(mviewer);CHKERRQ(ierr);
185*f49d1c87SHong Zhang   }
18609baa881SHong Zhang   ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
18709baa881SHong Zhang   PetscFunctionReturn(0);
18809baa881SHong Zhang }
18909baa881SHong Zhang 
190a7e14dcfSSatish Balay /*@C
191a7e14dcfSSatish Balay    TaoComputeHessian - Computes the Hessian matrix that has been
192a7e14dcfSSatish Balay    set with TaoSetHessianRoutine().
193a7e14dcfSSatish Balay 
194441846f8SBarry Smith    Collective on Tao
195a7e14dcfSSatish Balay 
196a7e14dcfSSatish Balay    Input Parameters:
197f4c1ad5cSStefano Zampini +  tao - the Tao solver context
198f4c1ad5cSStefano Zampini -  X   - input vector
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay    Output Parameters:
201a7e14dcfSSatish Balay +  H    - Hessian matrix
202aa6c7ce3SBarry Smith -  Hpre - Preconditioning matrix
203a7e14dcfSSatish Balay 
20409baa881SHong Zhang    Options Database Keys:
20509baa881SHong Zhang +     -tao_test_hessian - compare the user provided Hessian with one compute via finite differences to check for errors
20609baa881SHong 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
207dfe02fe6SHong 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
20809baa881SHong Zhang 
209a7e14dcfSSatish Balay    Notes:
210a7e14dcfSSatish Balay    Most users should not need to explicitly call this routine, as it
211a7e14dcfSSatish Balay    is used internally within the minimization solvers.
212a7e14dcfSSatish Balay 
213a7e14dcfSSatish Balay    TaoComputeHessian() is typically used within minimization
214a7e14dcfSSatish Balay    implementations, so most users would not generally call this routine
215a7e14dcfSSatish Balay    themselves.
216a7e14dcfSSatish Balay 
21709baa881SHong Zhang    Developer Notes:
21809baa881SHong Zhang    The Hessian test mechanism follows SNESTestJacobian().
21909baa881SHong Zhang 
220a7e14dcfSSatish Balay    Level: developer
221a7e14dcfSSatish Balay 
222f4c1ad5cSStefano Zampini .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetHessianRoutine()
223a7e14dcfSSatish Balay @*/
224ffad9901SBarry Smith PetscErrorCode TaoComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre)
225a7e14dcfSSatish Balay {
226a7e14dcfSSatish Balay   PetscErrorCode ierr;
22787f595a5SBarry Smith 
228a7e14dcfSSatish Balay   PetscFunctionBegin;
229441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
230a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
231a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
232f4c1ad5cSStefano Zampini   if (!tao->ops->computehessian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetHessianRoutine() first");
233a7e14dcfSSatish Balay   ++tao->nhess;
234f4c1ad5cSStefano Zampini   ierr = VecLockPush(X);CHKERRQ(ierr);
23594ab13aaSBarry Smith   ierr = PetscLogEventBegin(Tao_HessianEval,tao,X,H,Hpre);CHKERRQ(ierr);
236441846f8SBarry Smith   PetscStackPush("Tao user Hessian function");
237ffad9901SBarry Smith   ierr = (*tao->ops->computehessian)(tao,X,H,Hpre,tao->user_hessP);CHKERRQ(ierr);
238a7e14dcfSSatish Balay   PetscStackPop;
23994ab13aaSBarry Smith   ierr = PetscLogEventEnd(Tao_HessianEval,tao,X,H,Hpre);CHKERRQ(ierr);
240f4c1ad5cSStefano Zampini   ierr = VecLockPop(X);CHKERRQ(ierr);
24109baa881SHong Zhang 
24209baa881SHong Zhang   ierr = TaoTestHessian(tao);CHKERRQ(ierr);
243a7e14dcfSSatish Balay   PetscFunctionReturn(0);
244a7e14dcfSSatish Balay }
245a7e14dcfSSatish Balay 
246a7e14dcfSSatish Balay /*@C
247a7e14dcfSSatish Balay    TaoComputeJacobian - Computes the Jacobian matrix that has been
248a7e14dcfSSatish Balay    set with TaoSetJacobianRoutine().
249a7e14dcfSSatish Balay 
250441846f8SBarry Smith    Collective on Tao
251a7e14dcfSSatish Balay 
252a7e14dcfSSatish Balay    Input Parameters:
253f4c1ad5cSStefano Zampini +  tao - the Tao solver context
254f4c1ad5cSStefano Zampini -  X   - input vector
255a7e14dcfSSatish Balay 
256a7e14dcfSSatish Balay    Output Parameters:
257f4c1ad5cSStefano Zampini +  J    - Jacobian matrix
258f4c1ad5cSStefano Zampini -  Jpre - Preconditioning matrix
259a7e14dcfSSatish Balay 
260a7e14dcfSSatish Balay    Notes:
261a7e14dcfSSatish Balay    Most users should not need to explicitly call this routine, as it
262a7e14dcfSSatish Balay    is used internally within the minimization solvers.
263a7e14dcfSSatish Balay 
264a7e14dcfSSatish Balay    TaoComputeJacobian() is typically used within minimization
265a7e14dcfSSatish Balay    implementations, so most users would not generally call this routine
266a7e14dcfSSatish Balay    themselves.
267a7e14dcfSSatish Balay 
268a7e14dcfSSatish Balay    Level: developer
269a7e14dcfSSatish Balay 
270f4c1ad5cSStefano Zampini .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianRoutine()
271a7e14dcfSSatish Balay @*/
272ffad9901SBarry Smith PetscErrorCode TaoComputeJacobian(Tao tao, Vec X, Mat J, Mat Jpre)
273a7e14dcfSSatish Balay {
274a7e14dcfSSatish Balay   PetscErrorCode ierr;
27587f595a5SBarry Smith 
276a7e14dcfSSatish Balay   PetscFunctionBegin;
277441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
278a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
279a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
28087f595a5SBarry Smith   if (!tao->ops->computejacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobian() first");
281a7e14dcfSSatish Balay   ++tao->njac;
282f4c1ad5cSStefano Zampini   ierr = VecLockPush(X);CHKERRQ(ierr);
28394ab13aaSBarry Smith   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
284441846f8SBarry Smith   PetscStackPush("Tao user Jacobian function");
285ffad9901SBarry Smith   ierr = (*tao->ops->computejacobian)(tao,X,J,Jpre,tao->user_jacP);CHKERRQ(ierr);
286a7e14dcfSSatish Balay   PetscStackPop;
28794ab13aaSBarry Smith   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
288f4c1ad5cSStefano Zampini   ierr = VecLockPop(X);CHKERRQ(ierr);
289a7e14dcfSSatish Balay   PetscFunctionReturn(0);
290a7e14dcfSSatish Balay }
291a7e14dcfSSatish Balay 
292a7e14dcfSSatish Balay /*@C
293a7e14dcfSSatish Balay    TaoComputeJacobianState - Computes the Jacobian matrix that has been
294a7e14dcfSSatish Balay    set with TaoSetJacobianStateRoutine().
295a7e14dcfSSatish Balay 
296441846f8SBarry Smith    Collective on Tao
297a7e14dcfSSatish Balay 
298a7e14dcfSSatish Balay    Input Parameters:
299f4c1ad5cSStefano Zampini +  tao - the Tao solver context
300f4c1ad5cSStefano Zampini -  X   - input vector
301a7e14dcfSSatish Balay 
302a7e14dcfSSatish Balay    Output Parameters:
303f4c1ad5cSStefano Zampini +  Jpre - Jacobian matrix
304f4c1ad5cSStefano Zampini -  Jinv - Preconditioning matrix
305a7e14dcfSSatish Balay 
306a7e14dcfSSatish Balay    Notes:
307a7e14dcfSSatish Balay    Most users should not need to explicitly call this routine, as it
308a7e14dcfSSatish Balay    is used internally within the minimization solvers.
309a7e14dcfSSatish Balay 
310a7e14dcfSSatish Balay    TaoComputeJacobianState() is typically used within minimization
311a7e14dcfSSatish Balay    implementations, so most users would not generally call this routine
312a7e14dcfSSatish Balay    themselves.
313a7e14dcfSSatish Balay 
314a7e14dcfSSatish Balay    Level: developer
315a7e14dcfSSatish Balay 
316a7e14dcfSSatish Balay .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
317a7e14dcfSSatish Balay @*/
318ffad9901SBarry Smith PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv)
319a7e14dcfSSatish Balay {
320a7e14dcfSSatish Balay   PetscErrorCode ierr;
32145cf516eSBarry Smith 
322a7e14dcfSSatish Balay   PetscFunctionBegin;
323441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
324a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
325a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
32687f595a5SBarry Smith   if (!tao->ops->computejacobianstate) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first");
327a7e14dcfSSatish Balay   ++tao->njac_state;
328f4c1ad5cSStefano Zampini   ierr = VecLockPush(X);CHKERRQ(ierr);
32994ab13aaSBarry Smith   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
330441846f8SBarry Smith   PetscStackPush("Tao user Jacobian(state) function");
331ffad9901SBarry Smith   ierr = (*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP);CHKERRQ(ierr);
332a7e14dcfSSatish Balay   PetscStackPop;
33394ab13aaSBarry Smith   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
334f4c1ad5cSStefano Zampini   ierr = VecLockPop(X);CHKERRQ(ierr);
335a7e14dcfSSatish Balay   PetscFunctionReturn(0);
336a7e14dcfSSatish Balay }
337a7e14dcfSSatish Balay 
338a7e14dcfSSatish Balay /*@C
339a7e14dcfSSatish Balay    TaoComputeJacobianDesign - Computes the Jacobian matrix that has been
340a7e14dcfSSatish Balay    set with TaoSetJacobianDesignRoutine().
341a7e14dcfSSatish Balay 
342441846f8SBarry Smith    Collective on Tao
343a7e14dcfSSatish Balay 
344a7e14dcfSSatish Balay    Input Parameters:
345f4c1ad5cSStefano Zampini +  tao - the Tao solver context
346f4c1ad5cSStefano Zampini -  X   - input vector
347a7e14dcfSSatish Balay 
348a7e14dcfSSatish Balay    Output Parameters:
349f4c1ad5cSStefano Zampini .  J - Jacobian matrix
350a7e14dcfSSatish Balay 
351a7e14dcfSSatish Balay    Notes:
352a7e14dcfSSatish Balay    Most users should not need to explicitly call this routine, as it
353a7e14dcfSSatish Balay    is used internally within the minimization solvers.
354a7e14dcfSSatish Balay 
355a7e14dcfSSatish Balay    TaoComputeJacobianDesign() is typically used within minimization
356a7e14dcfSSatish Balay    implementations, so most users would not generally call this routine
357a7e14dcfSSatish Balay    themselves.
358a7e14dcfSSatish Balay 
359a7e14dcfSSatish Balay    Level: developer
360a7e14dcfSSatish Balay 
361a7e14dcfSSatish Balay .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianDesignRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
362a7e14dcfSSatish Balay @*/
36394ab13aaSBarry Smith PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J)
364a7e14dcfSSatish Balay {
365a7e14dcfSSatish Balay   PetscErrorCode ierr;
36687f595a5SBarry Smith 
367a7e14dcfSSatish Balay   PetscFunctionBegin;
368441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
369a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
370a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
37187f595a5SBarry Smith   if (!tao->ops->computejacobiandesign) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first");
372a7e14dcfSSatish Balay   ++tao->njac_design;
373f4c1ad5cSStefano Zampini   ierr = VecLockPush(X);CHKERRQ(ierr);
37494ab13aaSBarry Smith   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr);
375441846f8SBarry Smith   PetscStackPush("Tao user Jacobian(design) function");
376a7e14dcfSSatish Balay   ierr = (*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP);CHKERRQ(ierr);
377a7e14dcfSSatish Balay   PetscStackPop;
37894ab13aaSBarry Smith   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr);
379f4c1ad5cSStefano Zampini   ierr = VecLockPop(X);CHKERRQ(ierr);
380a7e14dcfSSatish Balay   PetscFunctionReturn(0);
381a7e14dcfSSatish Balay }
382a7e14dcfSSatish Balay 
383a7e14dcfSSatish Balay /*@C
384a7e14dcfSSatish Balay    TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix.
385a7e14dcfSSatish Balay 
386441846f8SBarry Smith    Logically collective on Tao
387a7e14dcfSSatish Balay 
388a7e14dcfSSatish Balay    Input Parameters:
389441846f8SBarry Smith +  tao  - the Tao context
390a7e14dcfSSatish Balay .  J    - Matrix used for the jacobian
391a7e14dcfSSatish Balay .  Jpre - Matrix that will be used operated on by preconditioner, can be same as J
392f4c1ad5cSStefano Zampini .  func - Jacobian evaluation routine
393a7e14dcfSSatish Balay -  ctx  - [optional] user-defined context for private data for the
3946c23d075SBarry Smith           Jacobian evaluation routine (may be NULL)
395a7e14dcfSSatish Balay 
396f4c1ad5cSStefano Zampini    Calling sequence of func:
397f4c1ad5cSStefano Zampini $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
398a7e14dcfSSatish Balay 
399441846f8SBarry Smith +  tao  - the Tao  context
400a7e14dcfSSatish Balay .  x    - input vector
401a7e14dcfSSatish Balay .  J    - Jacobian matrix
402f4c1ad5cSStefano Zampini .  Jpre - preconditioning matrix, usually the same as J
403a7e14dcfSSatish Balay -  ctx  - [optional] user-defined Jacobian context
404a7e14dcfSSatish Balay 
405a7e14dcfSSatish Balay    Level: intermediate
406a7e14dcfSSatish Balay @*/
407ffad9901SBarry Smith PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
408a7e14dcfSSatish Balay {
409a7e14dcfSSatish Balay   PetscErrorCode ierr;
410f4c1ad5cSStefano Zampini 
411a7e14dcfSSatish Balay   PetscFunctionBegin;
412441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
413a7e14dcfSSatish Balay   if (J) {
414a7e14dcfSSatish Balay     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
415a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,J,2);
416a7e14dcfSSatish Balay   }
417a7e14dcfSSatish Balay   if (Jpre) {
418a7e14dcfSSatish Balay     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
419a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,Jpre,3);
420a7e14dcfSSatish Balay   }
421a7e14dcfSSatish Balay   if (ctx) {
422a7e14dcfSSatish Balay     tao->user_jacP = ctx;
423a7e14dcfSSatish Balay   }
424a7e14dcfSSatish Balay   if (func) {
425a7e14dcfSSatish Balay     tao->ops->computejacobian = func;
426a7e14dcfSSatish Balay   }
427a7e14dcfSSatish Balay   if (J) {
428a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
42945cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian);CHKERRQ(ierr);
430a7e14dcfSSatish Balay     tao->jacobian = J;
431a7e14dcfSSatish Balay   }
432a7e14dcfSSatish Balay   if (Jpre) {
433a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
43445cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian_pre);CHKERRQ(ierr);
435a7e14dcfSSatish Balay     tao->jacobian_pre=Jpre;
436a7e14dcfSSatish Balay   }
437a7e14dcfSSatish Balay   PetscFunctionReturn(0);
438a7e14dcfSSatish Balay }
439a7e14dcfSSatish Balay 
440a7e14dcfSSatish Balay /*@C
441a7e14dcfSSatish Balay    TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian
442a7e14dcfSSatish Balay    (and its inverse) of the constraint function with respect to the state variables.
443a7e14dcfSSatish Balay    Used only for pde-constrained optimization.
444a7e14dcfSSatish Balay 
445441846f8SBarry Smith    Logically collective on Tao
446a7e14dcfSSatish Balay 
447a7e14dcfSSatish Balay    Input Parameters:
448441846f8SBarry Smith +  tao  - the Tao context
449a7e14dcfSSatish Balay .  J    - Matrix used for the jacobian
4506c23d075SBarry Smith .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.  Only used if Jinv is NULL
4516c23d075SBarry 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.
452f4c1ad5cSStefano Zampini .  func - Jacobian evaluation routine
453a7e14dcfSSatish Balay -  ctx  - [optional] user-defined context for private data for the
4546c23d075SBarry Smith           Jacobian evaluation routine (may be NULL)
455a7e14dcfSSatish Balay 
456f4c1ad5cSStefano Zampini    Calling sequence of func:
457f4c1ad5cSStefano Zampini $    func(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,void *ctx);
458a7e14dcfSSatish Balay 
459441846f8SBarry Smith +  tao  - the Tao  context
460a7e14dcfSSatish Balay .  x    - input vector
461a7e14dcfSSatish Balay .  J    - Jacobian matrix
462a7e14dcfSSatish Balay .  Jpre - preconditioner matrix, usually the same as J
463a7e14dcfSSatish Balay .  Jinv - inverse of J
464a7e14dcfSSatish Balay -  ctx  - [optional] user-defined Jacobian context
465a7e14dcfSSatish Balay 
466a7e14dcfSSatish Balay    Level: intermediate
467f4c1ad5cSStefano Zampini .seealso: TaoComputeJacobianState(), TaoSetJacobianDesignRoutine(), TaoSetStateDesignIS()
468a7e14dcfSSatish Balay @*/
469ffad9901SBarry Smith PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat, void*), void *ctx)
470a7e14dcfSSatish Balay {
471a7e14dcfSSatish Balay   PetscErrorCode ierr;
472f4c1ad5cSStefano Zampini 
473a7e14dcfSSatish Balay   PetscFunctionBegin;
474441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
475a7e14dcfSSatish Balay   if (J) {
476a7e14dcfSSatish Balay     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
477a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,J,2);
478a7e14dcfSSatish Balay   }
479a7e14dcfSSatish Balay   if (Jpre) {
480a7e14dcfSSatish Balay     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
481a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,Jpre,3);
482a7e14dcfSSatish Balay   }
483a7e14dcfSSatish Balay   if (Jinv) {
484a7e14dcfSSatish Balay     PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4);
485a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,Jinv,4);
486a7e14dcfSSatish Balay   }
487a7e14dcfSSatish Balay   if (ctx) {
488a7e14dcfSSatish Balay     tao->user_jac_stateP = ctx;
489a7e14dcfSSatish Balay   }
490a7e14dcfSSatish Balay   if (func) {
491a7e14dcfSSatish Balay     tao->ops->computejacobianstate = func;
492a7e14dcfSSatish Balay   }
493a7e14dcfSSatish Balay   if (J) {
494a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
49545cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian_state);CHKERRQ(ierr);
496a7e14dcfSSatish Balay     tao->jacobian_state = J;
497a7e14dcfSSatish Balay   }
498a7e14dcfSSatish Balay   if (Jpre) {
499a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
50045cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian_state_pre);CHKERRQ(ierr);
501a7e14dcfSSatish Balay     tao->jacobian_state_pre=Jpre;
502a7e14dcfSSatish Balay   }
503a7e14dcfSSatish Balay   if (Jinv) {
504a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)Jinv);CHKERRQ(ierr);
50545cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian_state_inv);CHKERRQ(ierr);
506a7e14dcfSSatish Balay     tao->jacobian_state_inv=Jinv;
507a7e14dcfSSatish Balay   }
508a7e14dcfSSatish Balay   PetscFunctionReturn(0);
509a7e14dcfSSatish Balay }
510a7e14dcfSSatish Balay 
511a7e14dcfSSatish Balay /*@C
512a7e14dcfSSatish Balay    TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of
513a7e14dcfSSatish Balay    the constraint function with respect to the design variables.  Used only for
514a7e14dcfSSatish Balay    pde-constrained optimization.
515a7e14dcfSSatish Balay 
516441846f8SBarry Smith    Logically collective on Tao
517a7e14dcfSSatish Balay 
518a7e14dcfSSatish Balay    Input Parameters:
519441846f8SBarry Smith +  tao  - the Tao context
520a7e14dcfSSatish Balay .  J    - Matrix used for the jacobian
521f4c1ad5cSStefano Zampini .  func - Jacobian evaluation routine
522a7e14dcfSSatish Balay -  ctx  - [optional] user-defined context for private data for the
5236c23d075SBarry Smith           Jacobian evaluation routine (may be NULL)
524a7e14dcfSSatish Balay 
525f4c1ad5cSStefano Zampini    Calling sequence of func:
526f4c1ad5cSStefano Zampini $    func(Tao tao,Vec x,Mat J,void *ctx);
527a7e14dcfSSatish Balay 
528441846f8SBarry Smith +  tao - the Tao  context
529a7e14dcfSSatish Balay .  x   - input vector
530a7e14dcfSSatish Balay .  J   - Jacobian matrix
531a7e14dcfSSatish Balay -  ctx - [optional] user-defined Jacobian context
532a7e14dcfSSatish Balay 
533a7e14dcfSSatish Balay    Level: intermediate
534f4c1ad5cSStefano Zampini 
535a7e14dcfSSatish Balay .seealso: TaoComputeJacobianDesign(), TaoSetJacobianStateRoutine(), TaoSetStateDesignIS()
536a7e14dcfSSatish Balay @*/
53794ab13aaSBarry Smith PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx)
538a7e14dcfSSatish Balay {
539a7e14dcfSSatish Balay   PetscErrorCode ierr;
54045cf516eSBarry Smith 
541a7e14dcfSSatish Balay   PetscFunctionBegin;
542441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
543a7e14dcfSSatish Balay   if (J) {
544a7e14dcfSSatish Balay     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
545a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,J,2);
546a7e14dcfSSatish Balay   }
547a7e14dcfSSatish Balay   if (ctx) {
548a7e14dcfSSatish Balay     tao->user_jac_designP = ctx;
549a7e14dcfSSatish Balay   }
550a7e14dcfSSatish Balay   if (func) {
551a7e14dcfSSatish Balay     tao->ops->computejacobiandesign = func;
552a7e14dcfSSatish Balay   }
553a7e14dcfSSatish Balay   if (J) {
554a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
55545cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian_design);CHKERRQ(ierr);
556a7e14dcfSSatish Balay     tao->jacobian_design = J;
557a7e14dcfSSatish Balay   }
558a7e14dcfSSatish Balay   PetscFunctionReturn(0);
559a7e14dcfSSatish Balay }
560a7e14dcfSSatish Balay 
561a7e14dcfSSatish Balay /*@
562441846f8SBarry Smith    TaoSetStateDesignIS - Indicate to the Tao which variables in the
563a7e14dcfSSatish Balay    solution vector are state variables and which are design.  Only applies to
564a7e14dcfSSatish Balay    pde-constrained optimization.
565a7e14dcfSSatish Balay 
566441846f8SBarry Smith    Logically Collective on Tao
567a7e14dcfSSatish Balay 
568a7e14dcfSSatish Balay    Input Parameters:
569441846f8SBarry Smith +  tao  - The Tao context
570a7e14dcfSSatish Balay .  s_is - the index set corresponding to the state variables
571a7e14dcfSSatish Balay -  d_is - the index set corresponding to the design variables
572a7e14dcfSSatish Balay 
573a7e14dcfSSatish Balay    Level: intermediate
574a7e14dcfSSatish Balay 
575a7e14dcfSSatish Balay .seealso: TaoSetJacobianStateRoutine(), TaoSetJacobianDesignRoutine()
576a7e14dcfSSatish Balay @*/
577441846f8SBarry Smith PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is)
578a7e14dcfSSatish Balay {
579a7e14dcfSSatish Balay   PetscErrorCode ierr;
58045cf516eSBarry Smith 
58145cf516eSBarry Smith   PetscFunctionBegin;
58245cf516eSBarry Smith   ierr = PetscObjectReference((PetscObject)s_is);CHKERRQ(ierr);
58345cf516eSBarry Smith   ierr = ISDestroy(&tao->state_is);CHKERRQ(ierr);
584a7e14dcfSSatish Balay   tao->state_is = s_is;
58545cf516eSBarry Smith   ierr = PetscObjectReference((PetscObject)(d_is));CHKERRQ(ierr);
58645cf516eSBarry Smith   ierr = ISDestroy(&tao->design_is);CHKERRQ(ierr);
587a7e14dcfSSatish Balay   tao->design_is = d_is;
588a7e14dcfSSatish Balay   PetscFunctionReturn(0);
589a7e14dcfSSatish Balay }
590a7e14dcfSSatish Balay 
591a7e14dcfSSatish Balay /*@C
592a7e14dcfSSatish Balay    TaoComputeJacobianEquality - Computes the Jacobian matrix that has been
593a7e14dcfSSatish Balay    set with TaoSetJacobianEqualityRoutine().
594a7e14dcfSSatish Balay 
595441846f8SBarry Smith    Collective on Tao
596a7e14dcfSSatish Balay 
597a7e14dcfSSatish Balay    Input Parameters:
598f4c1ad5cSStefano Zampini +  tao - the Tao solver context
599f4c1ad5cSStefano Zampini -  X   - input vector
600a7e14dcfSSatish Balay 
601a7e14dcfSSatish Balay    Output Parameters:
602f4c1ad5cSStefano Zampini +  J    - Jacobian matrix
603f4c1ad5cSStefano Zampini -  Jpre - Preconditioning matrix
604a7e14dcfSSatish Balay 
605a7e14dcfSSatish Balay    Notes:
606a7e14dcfSSatish Balay    Most users should not need to explicitly call this routine, as it
607a7e14dcfSSatish Balay    is used internally within the minimization solvers.
608a7e14dcfSSatish Balay 
609a7e14dcfSSatish Balay    Level: developer
610a7e14dcfSSatish Balay 
611a7e14dcfSSatish Balay .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
612a7e14dcfSSatish Balay @*/
613ffad9901SBarry Smith PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre)
614a7e14dcfSSatish Balay {
615a7e14dcfSSatish Balay   PetscErrorCode ierr;
61645cf516eSBarry Smith 
617a7e14dcfSSatish Balay   PetscFunctionBegin;
618441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
619a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
620a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
62187f595a5SBarry Smith   if (!tao->ops->computejacobianequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first");
622a7e14dcfSSatish Balay   ++tao->njac_equality;
623f4c1ad5cSStefano Zampini   ierr = VecLockPush(X);CHKERRQ(ierr);
62494ab13aaSBarry Smith   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
625441846f8SBarry Smith   PetscStackPush("Tao user Jacobian(equality) function");
626ffad9901SBarry Smith   ierr = (*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP);CHKERRQ(ierr);
627a7e14dcfSSatish Balay   PetscStackPop;
62894ab13aaSBarry Smith   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
629f4c1ad5cSStefano Zampini   ierr = VecLockPop(X);CHKERRQ(ierr);
630a7e14dcfSSatish Balay   PetscFunctionReturn(0);
631a7e14dcfSSatish Balay }
632a7e14dcfSSatish Balay 
633a7e14dcfSSatish Balay /*@C
634a7e14dcfSSatish Balay    TaoComputeJacobianInequality - Computes the Jacobian matrix that has been
635a7e14dcfSSatish Balay    set with TaoSetJacobianInequalityRoutine().
636a7e14dcfSSatish Balay 
637441846f8SBarry Smith    Collective on Tao
638a7e14dcfSSatish Balay 
639a7e14dcfSSatish Balay    Input Parameters:
640f4c1ad5cSStefano Zampini +  tao - the Tao solver context
641f4c1ad5cSStefano Zampini -  X   - input vector
642a7e14dcfSSatish Balay 
643a7e14dcfSSatish Balay    Output Parameters:
644f4c1ad5cSStefano Zampini +  J    - Jacobian matrix
645f4c1ad5cSStefano Zampini -  Jpre - Preconditioning matrix
646a7e14dcfSSatish Balay 
647a7e14dcfSSatish Balay    Notes:
648a7e14dcfSSatish Balay    Most users should not need to explicitly call this routine, as it
649a7e14dcfSSatish Balay    is used internally within the minimization solvers.
650a7e14dcfSSatish Balay 
651a7e14dcfSSatish Balay    Level: developer
652a7e14dcfSSatish Balay 
653a7e14dcfSSatish Balay .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
654a7e14dcfSSatish Balay @*/
655ffad9901SBarry Smith PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre)
656a7e14dcfSSatish Balay {
657a7e14dcfSSatish Balay   PetscErrorCode ierr;
65887f595a5SBarry Smith 
659a7e14dcfSSatish Balay   PetscFunctionBegin;
660441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
661a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
662a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
66387f595a5SBarry Smith   if (!tao->ops->computejacobianinequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first");
664a7e14dcfSSatish Balay   ++tao->njac_inequality;
665f4c1ad5cSStefano Zampini   ierr = VecLockPush(X);CHKERRQ(ierr);
66694ab13aaSBarry Smith   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
667441846f8SBarry Smith   PetscStackPush("Tao user Jacobian(inequality) function");
668ffad9901SBarry Smith   ierr = (*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP);CHKERRQ(ierr);
669a7e14dcfSSatish Balay   PetscStackPop;
67094ab13aaSBarry Smith   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
671f4c1ad5cSStefano Zampini   ierr = VecLockPop(X);CHKERRQ(ierr);
672a7e14dcfSSatish Balay   PetscFunctionReturn(0);
673a7e14dcfSSatish Balay }
674a7e14dcfSSatish Balay 
675a7e14dcfSSatish Balay /*@C
676a7e14dcfSSatish Balay    TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian
677a7e14dcfSSatish Balay    (and its inverse) of the constraint function with respect to the equality variables.
678a7e14dcfSSatish Balay    Used only for pde-constrained optimization.
679a7e14dcfSSatish Balay 
680441846f8SBarry Smith    Logically collective on Tao
681a7e14dcfSSatish Balay 
682a7e14dcfSSatish Balay    Input Parameters:
683441846f8SBarry Smith +  tao  - the Tao context
684a7e14dcfSSatish Balay .  J    - Matrix used for the jacobian
685a7e14dcfSSatish Balay .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
686f4c1ad5cSStefano Zampini .  func - Jacobian evaluation routine
687a7e14dcfSSatish Balay -  ctx  - [optional] user-defined context for private data for the
6886c23d075SBarry Smith           Jacobian evaluation routine (may be NULL)
689a7e14dcfSSatish Balay 
690f4c1ad5cSStefano Zampini    Calling sequence of func:
691f4c1ad5cSStefano Zampini $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
692a7e14dcfSSatish Balay 
693441846f8SBarry Smith +  tao  - the Tao  context
694a7e14dcfSSatish Balay .  x    - input vector
695a7e14dcfSSatish Balay .  J    - Jacobian matrix
696a7e14dcfSSatish Balay .  Jpre - preconditioner matrix, usually the same as J
697a7e14dcfSSatish Balay -  ctx  - [optional] user-defined Jacobian context
698a7e14dcfSSatish Balay 
699a7e14dcfSSatish Balay    Level: intermediate
700f4c1ad5cSStefano Zampini 
701f4c1ad5cSStefano Zampini .seealso: TaoComputeJacobianEquality(), TaoSetJacobianDesignRoutine(), TaoSetEqualityDesignIS()
702a7e14dcfSSatish Balay @*/
703ffad9901SBarry Smith PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
704a7e14dcfSSatish Balay {
705a7e14dcfSSatish Balay   PetscErrorCode ierr;
70645cf516eSBarry Smith 
707a7e14dcfSSatish Balay   PetscFunctionBegin;
708441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
709a7e14dcfSSatish Balay   if (J) {
710a7e14dcfSSatish Balay     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
711a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,J,2);
712a7e14dcfSSatish Balay   }
713a7e14dcfSSatish Balay   if (Jpre) {
714a7e14dcfSSatish Balay     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
715a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,Jpre,3);
716a7e14dcfSSatish Balay   }
717a7e14dcfSSatish Balay   if (ctx) {
718a7e14dcfSSatish Balay     tao->user_jac_equalityP = ctx;
719a7e14dcfSSatish Balay   }
720a7e14dcfSSatish Balay   if (func) {
721a7e14dcfSSatish Balay     tao->ops->computejacobianequality = func;
722a7e14dcfSSatish Balay   }
723a7e14dcfSSatish Balay   if (J) {
724a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
72545cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian_equality);CHKERRQ(ierr);
726a7e14dcfSSatish Balay     tao->jacobian_equality = J;
727a7e14dcfSSatish Balay   }
728a7e14dcfSSatish Balay   if (Jpre) {
729a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
73045cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian_equality_pre);CHKERRQ(ierr);
731a7e14dcfSSatish Balay     tao->jacobian_equality_pre=Jpre;
732a7e14dcfSSatish Balay   }
733a7e14dcfSSatish Balay   PetscFunctionReturn(0);
734a7e14dcfSSatish Balay }
735a7e14dcfSSatish Balay 
736a7e14dcfSSatish Balay /*@C
737a7e14dcfSSatish Balay    TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian
738a7e14dcfSSatish Balay    (and its inverse) of the constraint function with respect to the inequality variables.
739a7e14dcfSSatish Balay    Used only for pde-constrained optimization.
740a7e14dcfSSatish Balay 
741441846f8SBarry Smith    Logically collective on Tao
742a7e14dcfSSatish Balay 
743a7e14dcfSSatish Balay    Input Parameters:
744441846f8SBarry Smith +  tao  - the Tao context
745a7e14dcfSSatish Balay .  J    - Matrix used for the jacobian
746a7e14dcfSSatish Balay .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
747f4c1ad5cSStefano Zampini .  func - Jacobian evaluation routine
748a7e14dcfSSatish Balay -  ctx  - [optional] user-defined context for private data for the
7496c23d075SBarry Smith           Jacobian evaluation routine (may be NULL)
750a7e14dcfSSatish Balay 
751f4c1ad5cSStefano Zampini    Calling sequence of func:
752f4c1ad5cSStefano Zampini $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
753a7e14dcfSSatish Balay 
754441846f8SBarry Smith +  tao  - the Tao  context
755a7e14dcfSSatish Balay .  x    - input vector
756a7e14dcfSSatish Balay .  J    - Jacobian matrix
757a7e14dcfSSatish Balay .  Jpre - preconditioner matrix, usually the same as J
758a7e14dcfSSatish Balay -  ctx  - [optional] user-defined Jacobian context
759a7e14dcfSSatish Balay 
760a7e14dcfSSatish Balay    Level: intermediate
761f4c1ad5cSStefano Zampini 
762f4c1ad5cSStefano Zampini .seealso: TaoComputeJacobianInequality(), TaoSetJacobianDesignRoutine(), TaoSetInequalityDesignIS()
763a7e14dcfSSatish Balay @*/
764ffad9901SBarry Smith PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
765a7e14dcfSSatish Balay {
766a7e14dcfSSatish Balay   PetscErrorCode ierr;
767f4c1ad5cSStefano Zampini 
768a7e14dcfSSatish Balay   PetscFunctionBegin;
769441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
770a7e14dcfSSatish Balay   if (J) {
771a7e14dcfSSatish Balay     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
772a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,J,2);
773a7e14dcfSSatish Balay   }
774a7e14dcfSSatish Balay   if (Jpre) {
775a7e14dcfSSatish Balay     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
776a7e14dcfSSatish Balay     PetscCheckSameComm(tao,1,Jpre,3);
777a7e14dcfSSatish Balay   }
778a7e14dcfSSatish Balay   if (ctx) {
779a7e14dcfSSatish Balay     tao->user_jac_inequalityP = ctx;
780a7e14dcfSSatish Balay   }
781a7e14dcfSSatish Balay   if (func) {
782a7e14dcfSSatish Balay     tao->ops->computejacobianinequality = func;
783a7e14dcfSSatish Balay   }
784a7e14dcfSSatish Balay   if (J) {
785a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
78645cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian_inequality);CHKERRQ(ierr);
787a7e14dcfSSatish Balay     tao->jacobian_inequality = J;
788a7e14dcfSSatish Balay   }
789a7e14dcfSSatish Balay   if (Jpre) {
790a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
79145cf516eSBarry Smith     ierr = MatDestroy(&tao->jacobian_inequality_pre);CHKERRQ(ierr);
792a7e14dcfSSatish Balay     tao->jacobian_inequality_pre=Jpre;
793a7e14dcfSSatish Balay   }
794a7e14dcfSSatish Balay   PetscFunctionReturn(0);
795a7e14dcfSSatish Balay }
796