1c4762a1bSJed Brown static char help[] = "Demonstrates tapeless 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
18c4762a1bSJed Brown #define ADOLC_TAPELESS
19c4762a1bSJed Brown #define NUMBER_DIRECTIONS 3
20c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
21c4762a1bSJed Brown #include <adolc/adtl.h>
22c4762a1bSJed Brown using namespace adtl;
23c4762a1bSJed Brown
24c4762a1bSJed Brown typedef struct _n_User *User;
25c4762a1bSJed Brown struct _n_User {
26c4762a1bSJed Brown PetscReal mu;
27c4762a1bSJed Brown PetscReal next_output;
28c4762a1bSJed Brown PetscReal tprev;
29c4762a1bSJed Brown
30c4762a1bSJed Brown /* Automatic differentiation support */
31c4762a1bSJed Brown AdolcCtx *adctx;
32c4762a1bSJed Brown Vec F;
33c4762a1bSJed Brown };
34c4762a1bSJed Brown
35c4762a1bSJed Brown /*
36c4762a1bSJed Brown Residual evaluation templated, so as to allow for PetscScalar or adouble
37c4762a1bSJed Brown arguments.
38c4762a1bSJed Brown */
399371c9d4SSatish Balay template <class T>
EvaluateResidual(const T * x,T mu,T * f)40d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateResidual(const T *x, T mu, T *f)
41d71ae5a4SJacob Faibussowitsch {
42c4762a1bSJed Brown PetscFunctionBegin;
43c4762a1bSJed Brown f[0] = x[1];
44c4762a1bSJed Brown f[1] = mu * (1. - x[0] * x[0]) * x[1] - x[0];
453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
46c4762a1bSJed Brown }
47c4762a1bSJed Brown
48c4762a1bSJed Brown /*
49c4762a1bSJed Brown 'Passive' RHS function, used in residual evaluations during the time integration.
50c4762a1bSJed Brown */
RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,PetscCtx ctx)51*2a8381b2SBarry Smith static PetscErrorCode RHSFunctionPassive(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx)
52d71ae5a4SJacob Faibussowitsch {
53c4762a1bSJed Brown User user = (User)ctx;
54a8c08197SHong Zhang PetscScalar *f;
55a8c08197SHong Zhang const PetscScalar *x;
56c4762a1bSJed Brown
57c4762a1bSJed Brown PetscFunctionBeginUser;
589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
599566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
609566063dSJacob Faibussowitsch PetscCall(EvaluateResidual(x, user->mu, f));
619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
629566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown
66c4762a1bSJed Brown /*
67c4762a1bSJed Brown Compute the Jacobian w.r.t. x using tapeless mode of ADOL-C.
68c4762a1bSJed Brown */
RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,PetscCtx ctx)69*2a8381b2SBarry Smith static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, PetscCtx ctx)
70d71ae5a4SJacob Faibussowitsch {
71c4762a1bSJed Brown User user = (User)ctx;
72a8c08197SHong Zhang PetscScalar **J;
73a8c08197SHong Zhang const PetscScalar *x;
74c4762a1bSJed Brown adouble f_a[2]; /* 'active' double for dependent variables */
75c4762a1bSJed Brown adouble x_a[2], mu_a; /* 'active' doubles for independent variables */
76c4762a1bSJed Brown PetscInt i, j;
77c4762a1bSJed Brown
78c4762a1bSJed Brown PetscFunctionBeginUser;
79c4762a1bSJed Brown /* Set values for independent variables and parameters */
809566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
81c4762a1bSJed Brown x_a[0].setValue(x[0]);
82c4762a1bSJed Brown x_a[1].setValue(x[1]);
83c4762a1bSJed Brown mu_a.setValue(user->mu);
849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
85c4762a1bSJed Brown
86c4762a1bSJed Brown /* Set seed matrix as 3x3 identity matrix */
879371c9d4SSatish Balay x_a[0].setADValue(0, 1.);
889371c9d4SSatish Balay x_a[0].setADValue(1, 0.);
899371c9d4SSatish Balay x_a[0].setADValue(2, 0.);
909371c9d4SSatish Balay x_a[1].setADValue(0, 0.);
919371c9d4SSatish Balay x_a[1].setADValue(1, 1.);
929371c9d4SSatish Balay x_a[1].setADValue(2, 0.);
939371c9d4SSatish Balay mu_a.setADValue(0, 0.);
949371c9d4SSatish Balay mu_a.setADValue(1, 0.);
959371c9d4SSatish Balay mu_a.setADValue(2, 1.);
96c4762a1bSJed Brown
97c4762a1bSJed Brown /* Evaluate residual (on active variables) */
989566063dSJacob Faibussowitsch PetscCall(EvaluateResidual(x_a, mu_a, f_a));
99c4762a1bSJed Brown
100c4762a1bSJed Brown /* Extract derivatives */
1019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->adctx->n, &J));
102c4762a1bSJed Brown J[0] = (PetscScalar *)f_a[0].getADValue();
103c4762a1bSJed Brown J[1] = (PetscScalar *)f_a[1].getADValue();
104c4762a1bSJed Brown
105c4762a1bSJed Brown /* Set matrix values */
106c4762a1bSJed Brown for (i = 0; i < user->adctx->m; i++) {
10748a46eb9SPierre Jolivet for (j = 0; j < user->adctx->n; j++) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES));
108c4762a1bSJed Brown }
1099566063dSJacob Faibussowitsch PetscCall(PetscFree(J));
1109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
112c4762a1bSJed Brown if (A != B) {
1139566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1149566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
115c4762a1bSJed Brown }
1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
117c4762a1bSJed Brown }
118c4762a1bSJed Brown
119c4762a1bSJed Brown /*
120c4762a1bSJed Brown Compute the Jacobian w.r.t. mu using tapeless mode of ADOL-C.
121c4762a1bSJed Brown */
RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,PetscCtx ctx)122*2a8381b2SBarry Smith static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, PetscCtx ctx)
123d71ae5a4SJacob Faibussowitsch {
124c4762a1bSJed Brown User user = (User)ctx;
125c4762a1bSJed Brown PetscScalar **J;
126c4762a1bSJed Brown PetscScalar *x;
127c4762a1bSJed Brown adouble f_a[2]; /* 'active' double for dependent variables */
128c4762a1bSJed Brown adouble x_a[2], mu_a; /* 'active' doubles for independent variables */
129c4762a1bSJed Brown PetscInt i, j = 0;
130c4762a1bSJed Brown
131c4762a1bSJed Brown PetscFunctionBeginUser;
132c4762a1bSJed Brown /* Set values for independent variables and parameters */
1339566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x));
134c4762a1bSJed Brown x_a[0].setValue(x[0]);
135c4762a1bSJed Brown x_a[1].setValue(x[1]);
136c4762a1bSJed Brown mu_a.setValue(user->mu);
1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x));
138c4762a1bSJed Brown
139c4762a1bSJed Brown /* Set seed matrix as 3x3 identity matrix */
1409371c9d4SSatish Balay x_a[0].setADValue(0, 1.);
1419371c9d4SSatish Balay x_a[0].setADValue(1, 0.);
1429371c9d4SSatish Balay x_a[0].setADValue(2, 0.);
1439371c9d4SSatish Balay x_a[1].setADValue(0, 0.);
1449371c9d4SSatish Balay x_a[1].setADValue(1, 1.);
1459371c9d4SSatish Balay x_a[1].setADValue(2, 0.);
1469371c9d4SSatish Balay mu_a.setADValue(0, 0.);
1479371c9d4SSatish Balay mu_a.setADValue(1, 0.);
1489371c9d4SSatish Balay mu_a.setADValue(2, 1.);
149c4762a1bSJed Brown
150c4762a1bSJed Brown /* Evaluate residual (on active variables) */
1519566063dSJacob Faibussowitsch PetscCall(EvaluateResidual(x_a, mu_a, f_a));
152c4762a1bSJed Brown
153c4762a1bSJed Brown /* Extract derivatives */
1549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &J));
155c4762a1bSJed Brown J[0] = (PetscScalar *)f_a[0].getADValue();
156c4762a1bSJed Brown J[1] = (PetscScalar *)f_a[1].getADValue();
157c4762a1bSJed Brown
158c4762a1bSJed Brown /* Set matrix values */
15948a46eb9SPierre Jolivet for (i = 0; i < user->adctx->m; i++) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][user->adctx->n], INSERT_VALUES));
1609566063dSJacob Faibussowitsch PetscCall(PetscFree(J));
1619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown
166c4762a1bSJed Brown /*
167c4762a1bSJed Brown Monitor timesteps and use interpolation to output at integer multiples of 0.1
168c4762a1bSJed Brown */
Monitor(TS ts,PetscInt step,PetscReal t,Vec X,PetscCtx ctx)169*2a8381b2SBarry Smith static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, PetscCtx ctx)
170d71ae5a4SJacob Faibussowitsch {
171c4762a1bSJed Brown const PetscScalar *x;
172c4762a1bSJed Brown PetscReal tfinal, dt, tprev;
173c4762a1bSJed Brown User user = (User)ctx;
174c4762a1bSJed Brown
175c4762a1bSJed Brown PetscFunctionBeginUser;
1769566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt));
1779566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal));
1789566063dSJacob Faibussowitsch PetscCall(TSGetPrevTime(ts, &tprev));
1799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
18063a3b9bcSJacob 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])));
1819566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev));
1829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
184c4762a1bSJed Brown }
185c4762a1bSJed Brown
main(int argc,char ** argv)186d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
187d71ae5a4SJacob Faibussowitsch {
188c4762a1bSJed Brown TS ts; /* nonlinear solver */
189c4762a1bSJed Brown Vec x; /* solution, residual vectors */
190c4762a1bSJed Brown Mat A; /* Jacobian matrix */
191c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix */
192c4762a1bSJed Brown PetscInt steps;
193c4762a1bSJed Brown PetscReal ftime = 0.5;
194c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE;
195c4762a1bSJed Brown PetscScalar *x_ptr;
196c4762a1bSJed Brown PetscMPIInt size;
197c4762a1bSJed Brown struct _n_User user;
198c4762a1bSJed Brown AdolcCtx *adctx;
199c4762a1bSJed Brown Vec lambda[2], mu[2];
200c4762a1bSJed Brown
201c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202c4762a1bSJed Brown Initialize program
203c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204327415f7SBarry Smith PetscFunctionBeginUser;
2059566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
20708401ef6SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
208c4762a1bSJed Brown
209c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown Set runtime options and create AdolcCtx
211c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2129566063dSJacob Faibussowitsch PetscCall(PetscNew(&adctx));
213c4762a1bSJed Brown user.mu = 1;
214c4762a1bSJed Brown user.next_output = 0.0;
2159371c9d4SSatish Balay adctx->m = 2;
2169371c9d4SSatish Balay adctx->n = 2;
2179371c9d4SSatish Balay adctx->p = 2;
218c4762a1bSJed Brown user.adctx = adctx;
219c4762a1bSJed Brown adtl::setNumDir(adctx->n + 1); /* #indep. variables, plus parameters */
220c4762a1bSJed Brown
2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
223c4762a1bSJed Brown
224c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process
226c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2289566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
2299566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
2309566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
2319566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL));
232c4762a1bSJed Brown
2339566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
2349566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
2359566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jacp));
2369566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jacp));
237c4762a1bSJed Brown
238c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239c4762a1bSJed Brown Create timestepping solver context
240c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2419566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2429566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK));
2439566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, &user));
244c4762a1bSJed Brown
245c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
246c4762a1bSJed Brown Set initial conditions
247c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2489566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr));
2499371c9d4SSatish Balay x_ptr[0] = 2;
2509371c9d4SSatish Balay x_ptr[1] = 0.66666654321;
2519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr));
252c4762a1bSJed Brown
253c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
254c4762a1bSJed Brown Set RHS Jacobian for the adjoint integration
255c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2569566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &user));
2579566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime));
2589566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
25948a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
2609566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .001));
261c4762a1bSJed Brown
262c4762a1bSJed Brown /*
263c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used
264c4762a1bSJed Brown */
2659566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts));
266c4762a1bSJed Brown
267c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268c4762a1bSJed Brown Set runtime options
269c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2709566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
271c4762a1bSJed Brown
272c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
273c4762a1bSJed Brown Solve nonlinear system
274c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2759566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x));
2769566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime));
2779566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
27863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime));
2799566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
280c4762a1bSJed Brown
281c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
282c4762a1bSJed Brown Start the Adjoint model
283c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2849566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[0], NULL));
2859566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[1], NULL));
286c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */
2879566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &x_ptr));
2889371c9d4SSatish Balay x_ptr[0] = 1.0;
2899371c9d4SSatish Balay x_ptr[1] = 0.0;
2909566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &x_ptr));
2919566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[1], &x_ptr));
2929371c9d4SSatish Balay x_ptr[0] = 0.0;
2939371c9d4SSatish Balay x_ptr[1] = 1.0;
2949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[1], &x_ptr));
295c4762a1bSJed Brown
2969566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp, &mu[0], NULL));
2979566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp, &mu[1], NULL));
2989566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &x_ptr));
299c4762a1bSJed Brown x_ptr[0] = 0.0;
3009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &x_ptr));
3019566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[1], &x_ptr));
302c4762a1bSJed Brown x_ptr[0] = 0.0;
3039566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[1], &x_ptr));
3049566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 2, lambda, mu));
305c4762a1bSJed Brown
306c4762a1bSJed Brown /* Set RHS JacobianP */
3079566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &user));
308c4762a1bSJed Brown
3099566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts));
310c4762a1bSJed Brown
3119566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD));
3129566063dSJacob Faibussowitsch PetscCall(VecView(lambda[1], PETSC_VIEWER_STDOUT_WORLD));
3139566063dSJacob Faibussowitsch PetscCall(VecView(mu[0], PETSC_VIEWER_STDOUT_WORLD));
3149566063dSJacob Faibussowitsch PetscCall(VecView(mu[1], PETSC_VIEWER_STDOUT_WORLD));
315c4762a1bSJed Brown
316c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
317c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
318c4762a1bSJed Brown are no longer needed.
319c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
3219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jacp));
3229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
3239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0]));
3249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[1]));
3259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0]));
3269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[1]));
3279566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
3289566063dSJacob Faibussowitsch PetscCall(PetscFree(adctx));
3299566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
330b122ec5aSJacob Faibussowitsch return 0;
331c4762a1bSJed Brown }
332c4762a1bSJed Brown
333c4762a1bSJed Brown /*TEST
334c4762a1bSJed Brown
335c4762a1bSJed Brown build:
336c4762a1bSJed Brown requires: double !complex adolc
337c4762a1bSJed Brown
338c4762a1bSJed Brown test:
339c4762a1bSJed Brown suffix: 1
340c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
341c4762a1bSJed Brown output_file: output/ex16adj_tl_1.out
342c4762a1bSJed Brown
343c4762a1bSJed Brown test:
344c4762a1bSJed Brown suffix: 2
345c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5
346c4762a1bSJed Brown output_file: output/ex16adj_tl_2.out
347c4762a1bSJed Brown
348c4762a1bSJed Brown test:
349c4762a1bSJed Brown suffix: 3
350c4762a1bSJed Brown args: -ts_max_steps 10 -monitor
351c4762a1bSJed Brown output_file: output/ex16adj_tl_3.out
352c4762a1bSJed Brown
353c4762a1bSJed Brown test:
354c4762a1bSJed Brown suffix: 4
355c4762a1bSJed Brown args: -ts_max_steps 10 -monitor -mu 5
356c4762a1bSJed Brown output_file: output/ex16adj_tl_4.out
357c4762a1bSJed Brown
358c4762a1bSJed Brown TEST*/
359