xref: /petsc/src/ts/tutorials/autodiff/ex16adj.cxx (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 */
339371c9d4SSatish Balay static PetscErrorCode RHSFunctionPassive(TS ts, PetscReal t, Vec X, Vec F, void *ctx) {
34c4762a1bSJed Brown   User               user = (User)ctx;
35c4762a1bSJed Brown   PetscScalar       *f;
36c4762a1bSJed Brown   const PetscScalar *x;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   PetscFunctionBeginUser;
399566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
409566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
41c4762a1bSJed Brown   f[0] = x[1];
42c4762a1bSJed Brown   f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0];
439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
45c4762a1bSJed Brown   PetscFunctionReturn(0);
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /*
49c4762a1bSJed Brown   Trace RHS to mark on tape 1 the dependence of f upon x. This tape is used in generating the
50c4762a1bSJed Brown   Jacobian transform.
51c4762a1bSJed Brown */
529371c9d4SSatish Balay static PetscErrorCode RHSFunctionActive(TS ts, PetscReal t, Vec X, Vec F, void *ctx) {
53c4762a1bSJed Brown   User               user = (User)ctx;
54c4762a1bSJed Brown   PetscScalar       *f;
55c4762a1bSJed Brown   const PetscScalar *x;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   adouble f_a[2]; /* 'active' double for dependent variables */
58c4762a1bSJed Brown   adouble x_a[2]; /* 'active' double for independent variables */
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   PetscFunctionBeginUser;
619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Start of active section */
65c4762a1bSJed Brown   trace_on(1);
669371c9d4SSatish Balay   x_a[0] <<= x[0];
679371c9d4SSatish Balay   x_a[1] <<= x[1]; /* Mark independence */
68c4762a1bSJed Brown   f_a[0] = x_a[1];
69c4762a1bSJed Brown   f_a[1] = user->mu * (1. - x_a[0] * x_a[0]) * x_a[1] - x_a[0];
709371c9d4SSatish Balay   f_a[0] >>= f[0];
719371c9d4SSatish Balay   f_a[1] >>= f[1]; /* Mark dependence */
72c4762a1bSJed Brown   trace_off();
73c4762a1bSJed Brown   /* End of active section */
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
77c4762a1bSJed Brown   PetscFunctionReturn(0);
78c4762a1bSJed Brown }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown /*
81c4762a1bSJed Brown   Trace RHS again to mark on tape 2 the dependence of f upon the parameter mu. This tape is used in
82c4762a1bSJed Brown   generating JacobianP.
83c4762a1bSJed Brown */
849371c9d4SSatish Balay static PetscErrorCode RHSFunctionActiveP(TS ts, PetscReal t, Vec X, Vec F, void *ctx) {
85c4762a1bSJed Brown   User               user = (User)ctx;
86c4762a1bSJed Brown   PetscScalar       *f;
87c4762a1bSJed Brown   const PetscScalar *x;
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   adouble f_a[2];       /* 'active' double for dependent variables */
90c4762a1bSJed Brown   adouble x_a[2], mu_a; /* 'active' double for independent variables */
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   PetscFunctionBeginUser;
939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
949566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Start of active section */
97c4762a1bSJed Brown   trace_on(3);
989371c9d4SSatish Balay   x_a[0] <<= x[0];
999371c9d4SSatish Balay   x_a[1] <<= x[1];
1009371c9d4SSatish Balay   mu_a <<= user->mu; /* Mark independence */
101c4762a1bSJed Brown   f_a[0] = x_a[1];
102c4762a1bSJed Brown   f_a[1] = mu_a * (1. - x_a[0] * x_a[0]) * x_a[1] - x_a[0];
1039371c9d4SSatish Balay   f_a[0] >>= f[0];
1049371c9d4SSatish Balay   f_a[1] >>= f[1]; /* Mark dependence */
105c4762a1bSJed Brown   trace_off();
106c4762a1bSJed Brown   /* End of active section */
107c4762a1bSJed Brown 
1089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
110c4762a1bSJed Brown   PetscFunctionReturn(0);
111c4762a1bSJed Brown }
112c4762a1bSJed Brown 
113c4762a1bSJed Brown /*
114c4762a1bSJed Brown   Compute the Jacobian w.r.t. x using PETSc-ADOL-C driver for explicit TS.
115c4762a1bSJed Brown */
1169371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, void *ctx) {
117c4762a1bSJed Brown   User               user = (User)ctx;
118a8c08197SHong Zhang   const PetscScalar *x;
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   PetscFunctionBeginUser;
1219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1229566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeRHSJacobian(1, A, x, user->adctx));
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
124c4762a1bSJed Brown   PetscFunctionReturn(0);
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown /*
128c4762a1bSJed Brown   Compute the Jacobian w.r.t. mu using PETSc-ADOL-C driver for explicit TS.
129c4762a1bSJed Brown */
1309371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx) {
131c4762a1bSJed Brown   User               user = (User)ctx;
132410585f6SHong Zhang   const PetscScalar *x;
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   PetscFunctionBeginUser;
1359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1369566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeRHSJacobianP(3, A, x, &user->mu, user->adctx));
1379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
138c4762a1bSJed Brown   PetscFunctionReturn(0);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*
142c4762a1bSJed Brown   Monitor timesteps and use interpolation to output at integer multiples of 0.1
143c4762a1bSJed Brown */
1449371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) {
145c4762a1bSJed Brown   const PetscScalar *x;
146c4762a1bSJed Brown   PetscReal          tfinal, dt, tprev;
147c4762a1bSJed Brown   User               user = (User)ctx;
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   PetscFunctionBeginUser;
1509566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1519566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
1529566063dSJacob Faibussowitsch   PetscCall(TSGetPrevTime(ts, &tprev));
1539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
15463a3b9bcSJacob 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])));
1559566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev));
1569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
157c4762a1bSJed Brown   PetscFunctionReturn(0);
158c4762a1bSJed Brown }
159c4762a1bSJed Brown 
1609371c9d4SSatish Balay int main(int argc, char **argv) {
161c4762a1bSJed Brown   TS             ts;   /* nonlinear solver */
162c4762a1bSJed Brown   Vec            x;    /* solution, residual vectors */
163c4762a1bSJed Brown   Mat            A;    /* Jacobian matrix */
164c4762a1bSJed Brown   Mat            Jacp; /* JacobianP matrix */
165c4762a1bSJed Brown   PetscInt       steps;
166c4762a1bSJed Brown   PetscReal      ftime   = 0.5;
167c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
168c4762a1bSJed Brown   PetscScalar   *x_ptr;
169c4762a1bSJed Brown   PetscMPIInt    size;
170c4762a1bSJed Brown   struct _n_User user;
171c4762a1bSJed Brown   AdolcCtx      *adctx;
172c4762a1bSJed Brown   Vec            lambda[2], mu[2], r;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175c4762a1bSJed Brown      Initialize program
176c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177327415f7SBarry Smith   PetscFunctionBeginUser;
1789566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
18008401ef6SPierre Jolivet   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183c4762a1bSJed Brown     Set runtime options and create AdolcCtx
184c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1859566063dSJacob Faibussowitsch   PetscCall(PetscNew(&adctx));
186c4762a1bSJed Brown   user.mu           = 1;
187c4762a1bSJed Brown   user.next_output  = 0.0;
1889371c9d4SSatish Balay   adctx->m          = 2;
1899371c9d4SSatish Balay   adctx->n          = 2;
1909371c9d4SSatish Balay   adctx->p          = 2;
1919371c9d4SSatish Balay   adctx->num_params = 1;
192c4762a1bSJed Brown   user.adctx        = adctx;
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
1959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
199c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2009566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2019566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
2029566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2039566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
2049566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, NULL));
205c4762a1bSJed Brown 
2069566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
2079566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
2089566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jacp));
2099566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Jacp));
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212c4762a1bSJed Brown      Create timestepping solver context
213c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2149566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2159566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSRK));
2169566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, &user));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219c4762a1bSJed Brown      Set initial conditions
220c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2219566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &x_ptr));
2229371c9d4SSatish Balay   x_ptr[0] = 2;
2239371c9d4SSatish Balay   x_ptr[1] = 0.66666654321;
2249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &x_ptr));
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227c4762a1bSJed Brown      Trace just once on each tape and put zeros on Jacobian diagonal
228c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2299566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
2309566063dSJacob Faibussowitsch   PetscCall(RHSFunctionActive(ts, 0., x, r, &user));
2319566063dSJacob Faibussowitsch   PetscCall(RHSFunctionActiveP(ts, 0., x, r, &user));
2329566063dSJacob Faibussowitsch   PetscCall(VecSet(r, 0));
2339566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(A, r, INSERT_VALUES));
2349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237c4762a1bSJed Brown      Set RHS Jacobian for the adjoint integration
238c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2399566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &user));
2409566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
2419566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
242*48a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
2439566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .001));
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /*
246c4762a1bSJed Brown     Have the TS save its trajectory so that TSAdjointSolve() may be used
247c4762a1bSJed Brown   */
2489566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
251c4762a1bSJed Brown      Set runtime options
252c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2539566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
256c4762a1bSJed Brown      Solve nonlinear system
257c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2589566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
2599566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2609566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
26163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime));
2629566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
265c4762a1bSJed Brown      Start the Adjoint model
266c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2679566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
2689566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lambda[1], NULL));
269c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
2709566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &x_ptr));
2719371c9d4SSatish Balay   x_ptr[0] = 1.0;
2729371c9d4SSatish Balay   x_ptr[1] = 0.0;
2739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &x_ptr));
2749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[1], &x_ptr));
2759371c9d4SSatish Balay   x_ptr[0] = 0.0;
2769371c9d4SSatish Balay   x_ptr[1] = 1.0;
2779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[1], &x_ptr));
278c4762a1bSJed Brown 
2799566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp, &mu[0], NULL));
2809566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp, &mu[1], NULL));
2819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &x_ptr));
282c4762a1bSJed Brown   x_ptr[0] = 0.0;
2839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &x_ptr));
2849566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[1], &x_ptr));
285c4762a1bSJed Brown   x_ptr[0] = 0.0;
2869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[1], &x_ptr));
2879566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 2, lambda, mu));
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   /*   Set RHS JacobianP */
2909566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &user));
291c4762a1bSJed Brown 
2929566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
293c4762a1bSJed Brown 
2949566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD));
2959566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[1], PETSC_VIEWER_STDOUT_WORLD));
2969566063dSJacob Faibussowitsch   PetscCall(VecView(mu[0], PETSC_VIEWER_STDOUT_WORLD));
2979566063dSJacob Faibussowitsch   PetscCall(VecView(mu[1], PETSC_VIEWER_STDOUT_WORLD));
298c4762a1bSJed Brown 
299c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
300c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
301c4762a1bSJed Brown      are no longer needed.
302c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jacp));
3059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
3079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[1]));
3089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
3099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[1]));
3109566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3119566063dSJacob Faibussowitsch   PetscCall(PetscFree(adctx));
3129566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
313b122ec5aSJacob Faibussowitsch   return 0;
314c4762a1bSJed Brown }
315c4762a1bSJed Brown 
316c4762a1bSJed Brown /*TEST
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   build:
319c4762a1bSJed Brown     requires: double !complex adolc
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   test:
322c4762a1bSJed Brown     suffix: 1
323c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
324c4762a1bSJed Brown     output_file: output/ex16adj_1.out
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   test:
327c4762a1bSJed Brown     suffix: 2
328c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5
329c4762a1bSJed Brown     output_file: output/ex16adj_2.out
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   test:
332c4762a1bSJed Brown     suffix: 3
333c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor
334c4762a1bSJed Brown     output_file: output/ex16adj_3.out
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   test:
337c4762a1bSJed Brown     suffix: 4
338c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor -mu 5
339c4762a1bSJed Brown     output_file: output/ex16adj_4.out
340c4762a1bSJed Brown 
341c4762a1bSJed Brown TEST*/
342