xref: /petsc/src/ts/tutorials/autodiff/ex16adj.cxx (revision 4e8208cbcbc709572b8abe32f33c78b69c819375) !
1c4762a1bSJed Brown static char help[] = "Demonstrates automatic Jacobian generation using ADOL-C for an adjoint sensitivity analysis of the van der Pol equation.\n\
2c4762a1bSJed Brown Input parameters include:\n\
3c4762a1bSJed Brown       -mu : stiffness parameter\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
7c4762a1bSJed Brown 
8c4762a1bSJed Brown    For documentation on ADOL-C, see
9c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
10c4762a1bSJed Brown */
11c4762a1bSJed Brown /* ------------------------------------------------------------------------
12c4762a1bSJed Brown    See ex16adj for a description of the problem being solved.
13c4762a1bSJed Brown   ------------------------------------------------------------------------- */
14c4762a1bSJed Brown 
15c4762a1bSJed Brown #include <petscts.h>
16c4762a1bSJed Brown #include <petscmat.h>
17c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
18c4762a1bSJed Brown #include <adolc/adolc.h>
19c4762a1bSJed Brown 
20c4762a1bSJed Brown typedef struct _n_User *User;
21c4762a1bSJed Brown struct _n_User {
22c4762a1bSJed Brown   PetscReal mu;
23c4762a1bSJed Brown   PetscReal next_output;
24c4762a1bSJed Brown   PetscReal tprev;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Automatic differentiation support */
27c4762a1bSJed Brown   AdolcCtx *adctx;
28c4762a1bSJed Brown };
29c4762a1bSJed Brown 
30c4762a1bSJed Brown /*
31c4762a1bSJed Brown   'Passive' RHS function, used in residual evaluations during the time integration.
32c4762a1bSJed Brown */
RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,PetscCtx ctx)33*2a8381b2SBarry Smith static PetscErrorCode RHSFunctionPassive(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx)
34d71ae5a4SJacob Faibussowitsch {
35c4762a1bSJed Brown   User               user = (User)ctx;
36c4762a1bSJed Brown   PetscScalar       *f;
37c4762a1bSJed Brown   const PetscScalar *x;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   PetscFunctionBeginUser;
409566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
42c4762a1bSJed Brown   f[0] = x[1];
43c4762a1bSJed Brown   f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0];
449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /*
50c4762a1bSJed Brown   Trace RHS to mark on tape 1 the dependence of f upon x. This tape is used in generating the
51c4762a1bSJed Brown   Jacobian transform.
52c4762a1bSJed Brown */
RHSFunctionActive(TS ts,PetscReal t,Vec X,Vec F,PetscCtx ctx)53*2a8381b2SBarry Smith static PetscErrorCode RHSFunctionActive(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx)
54d71ae5a4SJacob Faibussowitsch {
55c4762a1bSJed Brown   User               user = (User)ctx;
56c4762a1bSJed Brown   PetscScalar       *f;
57c4762a1bSJed Brown   const PetscScalar *x;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   adouble f_a[2]; /* 'active' double for dependent variables */
60c4762a1bSJed Brown   adouble x_a[2]; /* 'active' double for independent variables */
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   PetscFunctionBeginUser;
639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* Start of active section */
67c4762a1bSJed Brown   trace_on(1);
689371c9d4SSatish Balay   x_a[0] <<= x[0];
699371c9d4SSatish Balay   x_a[1] <<= x[1]; /* Mark independence */
70c4762a1bSJed Brown   f_a[0] = x_a[1];
71c4762a1bSJed Brown   f_a[1] = user->mu * (1. - x_a[0] * x_a[0]) * x_a[1] - x_a[0];
729371c9d4SSatish Balay   f_a[0] >>= f[0];
739371c9d4SSatish Balay   f_a[1] >>= f[1]; /* Mark dependence */
74c4762a1bSJed Brown   trace_off();
75c4762a1bSJed Brown   /* End of active section */
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80c4762a1bSJed Brown }
81c4762a1bSJed Brown 
82c4762a1bSJed Brown /*
83c4762a1bSJed Brown   Trace RHS again to mark on tape 2 the dependence of f upon the parameter mu. This tape is used in
84c4762a1bSJed Brown   generating JacobianP.
85c4762a1bSJed Brown */
RHSFunctionActiveP(TS ts,PetscReal t,Vec X,Vec F,PetscCtx ctx)86*2a8381b2SBarry Smith static PetscErrorCode RHSFunctionActiveP(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx)
87d71ae5a4SJacob Faibussowitsch {
88c4762a1bSJed Brown   User               user = (User)ctx;
89c4762a1bSJed Brown   PetscScalar       *f;
90c4762a1bSJed Brown   const PetscScalar *x;
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   adouble f_a[2];       /* 'active' double for dependent variables */
93c4762a1bSJed Brown   adouble x_a[2], mu_a; /* 'active' double for independent variables */
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   PetscFunctionBeginUser;
969566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
979566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* Start of active section */
100c4762a1bSJed Brown   trace_on(3);
1019371c9d4SSatish Balay   x_a[0] <<= x[0];
1029371c9d4SSatish Balay   x_a[1] <<= x[1];
1039371c9d4SSatish Balay   mu_a <<= user->mu; /* Mark independence */
104c4762a1bSJed Brown   f_a[0] = x_a[1];
105c4762a1bSJed Brown   f_a[1] = mu_a * (1. - x_a[0] * x_a[0]) * x_a[1] - x_a[0];
1069371c9d4SSatish Balay   f_a[0] >>= f[0];
1079371c9d4SSatish Balay   f_a[1] >>= f[1]; /* Mark dependence */
108c4762a1bSJed Brown   trace_off();
109c4762a1bSJed Brown   /* End of active section */
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown /*
117c4762a1bSJed Brown   Compute the Jacobian w.r.t. x using PETSc-ADOL-C driver for explicit TS.
118c4762a1bSJed Brown */
RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,PetscCtx ctx)119*2a8381b2SBarry Smith static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, PetscCtx ctx)
120d71ae5a4SJacob Faibussowitsch {
121c4762a1bSJed Brown   User               user = (User)ctx;
122a8c08197SHong Zhang   const PetscScalar *x;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   PetscFunctionBeginUser;
1259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1269566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeRHSJacobian(1, A, x, user->adctx));
1279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown /*
132c4762a1bSJed Brown   Compute the Jacobian w.r.t. mu using PETSc-ADOL-C driver for explicit TS.
133c4762a1bSJed Brown */
RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,PetscCtx ctx)134*2a8381b2SBarry Smith static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, PetscCtx ctx)
135d71ae5a4SJacob Faibussowitsch {
136c4762a1bSJed Brown   User               user = (User)ctx;
137410585f6SHong Zhang   const PetscScalar *x;
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   PetscFunctionBeginUser;
1409566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1419566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeRHSJacobianP(3, A, x, &user->mu, user->adctx));
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown /*
147c4762a1bSJed Brown   Monitor timesteps and use interpolation to output at integer multiples of 0.1
148c4762a1bSJed Brown */
Monitor(TS ts,PetscInt step,PetscReal t,Vec X,PetscCtx ctx)149*2a8381b2SBarry Smith static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, PetscCtx ctx)
150d71ae5a4SJacob Faibussowitsch {
151c4762a1bSJed Brown   const PetscScalar *x;
152c4762a1bSJed Brown   PetscReal          tfinal, dt, tprev;
153c4762a1bSJed Brown   User               user = (User)ctx;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   PetscFunctionBeginUser;
1569566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1579566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
1589566063dSJacob Faibussowitsch   PetscCall(TSGetPrevTime(ts, &tprev));
1599566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
16063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%.1f] %" PetscInt_FMT " TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n", (double)user->next_output, step, (double)t, (double)dt, (double)PetscRealPart(x[0]), (double)PetscRealPart(x[1])));
1619566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev));
1629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
main(int argc,char ** argv)166d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
167d71ae5a4SJacob Faibussowitsch {
168c4762a1bSJed Brown   TS             ts;   /* nonlinear solver */
169c4762a1bSJed Brown   Vec            x;    /* solution, residual vectors */
170c4762a1bSJed Brown   Mat            A;    /* Jacobian matrix */
171c4762a1bSJed Brown   Mat            Jacp; /* JacobianP matrix */
172c4762a1bSJed Brown   PetscInt       steps;
173c4762a1bSJed Brown   PetscReal      ftime   = 0.5;
174c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
175c4762a1bSJed Brown   PetscScalar   *x_ptr;
176c4762a1bSJed Brown   PetscMPIInt    size;
177c4762a1bSJed Brown   struct _n_User user;
178c4762a1bSJed Brown   AdolcCtx      *adctx;
179c4762a1bSJed Brown   Vec            lambda[2], mu[2], r;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182c4762a1bSJed Brown      Initialize program
183c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
184327415f7SBarry Smith   PetscFunctionBeginUser;
1859566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
18708401ef6SPierre Jolivet   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190c4762a1bSJed Brown     Set runtime options and create AdolcCtx
191c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1929566063dSJacob Faibussowitsch   PetscCall(PetscNew(&adctx));
193c4762a1bSJed Brown   user.mu           = 1;
194c4762a1bSJed Brown   user.next_output  = 0.0;
1959371c9d4SSatish Balay   adctx->m          = 2;
1969371c9d4SSatish Balay   adctx->n          = 2;
1979371c9d4SSatish Balay   adctx->p          = 2;
1989371c9d4SSatish Balay   adctx->num_params = 1;
199c4762a1bSJed Brown   user.adctx        = adctx;
200c4762a1bSJed Brown 
2019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
2029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
206c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2079566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2089566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
2099566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2109566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
2119566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, NULL));
212c4762a1bSJed Brown 
2139566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
2149566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
2159566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jacp));
2169566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Jacp));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219c4762a1bSJed Brown      Create timestepping solver context
220c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2219566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2229566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSRK));
2239566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, &user));
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226c4762a1bSJed Brown      Set initial conditions
227c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2289566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &x_ptr));
2299371c9d4SSatish Balay   x_ptr[0] = 2;
2309371c9d4SSatish Balay   x_ptr[1] = 0.66666654321;
2319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &x_ptr));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234c4762a1bSJed Brown      Trace just once on each tape and put zeros on Jacobian diagonal
235c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
2379566063dSJacob Faibussowitsch   PetscCall(RHSFunctionActive(ts, 0., x, r, &user));
2389566063dSJacob Faibussowitsch   PetscCall(RHSFunctionActiveP(ts, 0., x, r, &user));
2399566063dSJacob Faibussowitsch   PetscCall(VecSet(r, 0));
2409566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(A, r, INSERT_VALUES));
2419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown      Set RHS Jacobian for the adjoint integration
245c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2469566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &user));
2479566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
2489566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
24948a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
2509566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .001));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /*
253c4762a1bSJed Brown     Have the TS save its trajectory so that TSAdjointSolve() may be used
254c4762a1bSJed Brown   */
2559566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258c4762a1bSJed Brown      Set runtime options
259c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2609566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
263c4762a1bSJed Brown      Solve nonlinear system
264c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2659566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
2669566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2679566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
26863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime));
2699566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272c4762a1bSJed Brown      Start the Adjoint model
273c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2749566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
2759566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lambda[1], NULL));
276c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
2779566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &x_ptr));
2789371c9d4SSatish Balay   x_ptr[0] = 1.0;
2799371c9d4SSatish Balay   x_ptr[1] = 0.0;
2809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &x_ptr));
2819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[1], &x_ptr));
2829371c9d4SSatish Balay   x_ptr[0] = 0.0;
2839371c9d4SSatish Balay   x_ptr[1] = 1.0;
2849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[1], &x_ptr));
285c4762a1bSJed Brown 
2869566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp, &mu[0], NULL));
2879566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp, &mu[1], NULL));
2889566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &x_ptr));
289c4762a1bSJed Brown   x_ptr[0] = 0.0;
2909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &x_ptr));
2919566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[1], &x_ptr));
292c4762a1bSJed Brown   x_ptr[0] = 0.0;
2939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[1], &x_ptr));
2949566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 2, lambda, mu));
295c4762a1bSJed Brown 
296c4762a1bSJed Brown   /*   Set RHS JacobianP */
2979566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &user));
298c4762a1bSJed Brown 
2999566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
300c4762a1bSJed Brown 
3019566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD));
3029566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[1], PETSC_VIEWER_STDOUT_WORLD));
3039566063dSJacob Faibussowitsch   PetscCall(VecView(mu[0], PETSC_VIEWER_STDOUT_WORLD));
3049566063dSJacob Faibussowitsch   PetscCall(VecView(mu[1], PETSC_VIEWER_STDOUT_WORLD));
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
307c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
308c4762a1bSJed Brown      are no longer needed.
309c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3109566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3119566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jacp));
3129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
3149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[1]));
3159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
3169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[1]));
3179566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3189566063dSJacob Faibussowitsch   PetscCall(PetscFree(adctx));
3199566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
320b122ec5aSJacob Faibussowitsch   return 0;
321c4762a1bSJed Brown }
322c4762a1bSJed Brown 
323c4762a1bSJed Brown /*TEST
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   build:
326c4762a1bSJed Brown     requires: double !complex adolc
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   test:
329c4762a1bSJed Brown     suffix: 1
330c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
331c4762a1bSJed Brown     output_file: output/ex16adj_1.out
332c4762a1bSJed Brown 
333c4762a1bSJed Brown   test:
334c4762a1bSJed Brown     suffix: 2
335c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5
336c4762a1bSJed Brown     output_file: output/ex16adj_2.out
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   test:
339c4762a1bSJed Brown     suffix: 3
340c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor
341c4762a1bSJed Brown     output_file: output/ex16adj_3.out
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   test:
344c4762a1bSJed Brown     suffix: 4
345c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor -mu 5
346c4762a1bSJed Brown     output_file: output/ex16adj_4.out
347c4762a1bSJed Brown 
348c4762a1bSJed Brown TEST*/
349