xref: /petsc/src/ts/tutorials/autodiff/adr_ex1.cxx (revision d0e6bf2ad94dcc89b258ce16c7987200a4714786)
1c4762a1bSJed Brown static char help[] = "Demonstrates automatic Jacobian generation using ADOL-C for a nonlinear reaction problem from chemistry.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
5c4762a1bSJed Brown 
6c4762a1bSJed Brown    For documentation on ADOL-C, see
7c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
8c4762a1bSJed Brown */
9c4762a1bSJed Brown /* ------------------------------------------------------------------------
10c4762a1bSJed Brown   See ../advection-diffusion-reaction/ex1 for a description of the problem
11c4762a1bSJed Brown   ------------------------------------------------------------------------- */
12c4762a1bSJed Brown #include <petscts.h>
13c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
14c4762a1bSJed Brown #include <adolc/adolc.h>
15c4762a1bSJed Brown 
16c4762a1bSJed Brown typedef struct {
17c4762a1bSJed Brown   PetscScalar k;
18c4762a1bSJed Brown   Vec         initialsolution;
19c4762a1bSJed Brown   AdolcCtx   *adctx; /* Automatic differentiation support */
20c4762a1bSJed Brown } AppCtx;
21c4762a1bSJed Brown 
IFunctionView(AppCtx * ctx,PetscViewer v)22d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionView(AppCtx *ctx, PetscViewer v)
23d71ae5a4SJacob Faibussowitsch {
24c4762a1bSJed Brown   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(v, &ctx->k, 1, PETSC_SCALAR));
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27c4762a1bSJed Brown }
28c4762a1bSJed Brown 
IFunctionLoad(AppCtx ** ctx,PetscViewer v)29d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionLoad(AppCtx **ctx, PetscViewer v)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
339566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(v, &(*ctx)->k, 1, NULL, PETSC_SCALAR));
343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown /*
38c4762a1bSJed Brown   Defines the ODE passed to the ODE solver
39c4762a1bSJed Brown */
IFunctionPassive(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx * ctx)40d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionPassive(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
41d71ae5a4SJacob Faibussowitsch {
42c4762a1bSJed Brown   PetscScalar       *f;
43c4762a1bSJed Brown   const PetscScalar *u, *udot;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   PetscFunctionBegin;
46c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
50c4762a1bSJed Brown   f[0] = udot[0] + ctx->k * u[0] * u[1];
51c4762a1bSJed Brown   f[1] = udot[1] + ctx->k * u[0] * u[1];
52c4762a1bSJed Brown   f[2] = udot[2] - ctx->k * u[0] * u[1];
539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*
60c4762a1bSJed Brown   'Active' ADOL-C annotated version, marking dependence upon u.
61c4762a1bSJed Brown */
IFunctionActive1(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx * ctx)62d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionActive1(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
63d71ae5a4SJacob Faibussowitsch {
64c4762a1bSJed Brown   PetscScalar       *f;
65c4762a1bSJed Brown   const PetscScalar *u, *udot;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   adouble f_a[3]; /* 'active' double for dependent variables */
68c4762a1bSJed Brown   adouble u_a[3]; /* 'active' double for independent variables */
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBegin;
71c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
729566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* Start of active section */
77c4762a1bSJed Brown   trace_on(1);
789371c9d4SSatish Balay   u_a[0] <<= u[0];
799371c9d4SSatish Balay   u_a[1] <<= u[1];
809371c9d4SSatish Balay   u_a[2] <<= u[2]; /* Mark independence */
81c4762a1bSJed Brown   f_a[0] = udot[0] + ctx->k * u_a[0] * u_a[1];
82c4762a1bSJed Brown   f_a[1] = udot[1] + ctx->k * u_a[0] * u_a[1];
83c4762a1bSJed Brown   f_a[2] = udot[2] - ctx->k * u_a[0] * u_a[1];
849371c9d4SSatish Balay   f_a[0] >>= f[0];
859371c9d4SSatish Balay   f_a[1] >>= f[1];
869371c9d4SSatish Balay   f_a[2] >>= f[2]; /* Mark dependence */
87c4762a1bSJed Brown   trace_off();
88c4762a1bSJed Brown   /* End of active section */
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown /*
97c4762a1bSJed Brown   'Active' ADOL-C annotated version, marking dependence upon udot.
98c4762a1bSJed Brown */
IFunctionActive2(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx * ctx)99d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionActive2(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
100d71ae5a4SJacob Faibussowitsch {
101c4762a1bSJed Brown   PetscScalar       *f;
102c4762a1bSJed Brown   const PetscScalar *u, *udot;
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   adouble f_a[3];    /* 'active' double for dependent variables */
105c4762a1bSJed Brown   adouble udot_a[3]; /* 'active' double for independent variables */
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   PetscFunctionBegin;
108c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /* Start of active section */
114c4762a1bSJed Brown   trace_on(2);
1159371c9d4SSatish Balay   udot_a[0] <<= udot[0];
1169371c9d4SSatish Balay   udot_a[1] <<= udot[1];
1179371c9d4SSatish Balay   udot_a[2] <<= udot[2]; /* Mark independence */
118c4762a1bSJed Brown   f_a[0] = udot_a[0] + ctx->k * u[0] * u[1];
119c4762a1bSJed Brown   f_a[1] = udot_a[1] + ctx->k * u[0] * u[1];
120c4762a1bSJed Brown   f_a[2] = udot_a[2] - ctx->k * u[0] * u[1];
1219371c9d4SSatish Balay   f_a[0] >>= f[0];
1229371c9d4SSatish Balay   f_a[1] >>= f[1];
1239371c9d4SSatish Balay   f_a[2] >>= f[2]; /* Mark dependence */
124c4762a1bSJed Brown   trace_off();
125c4762a1bSJed Brown   /* End of active section */
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131c4762a1bSJed Brown }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown /*
134c4762a1bSJed Brown  Defines the Jacobian of the ODE passed to the ODE solver, using the PETSc-ADOL-C driver for
135c4762a1bSJed Brown  implicit TS.
136c4762a1bSJed Brown */
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx * ctx)137d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
138d71ae5a4SJacob Faibussowitsch {
139c4762a1bSJed Brown   AppCtx            *appctx = (AppCtx *)ctx;
140410585f6SHong Zhang   const PetscScalar *u;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   PetscFunctionBegin;
1439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1449566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeIJacobian(1, 2, A, u, a, appctx->adctx));
1459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
147c4762a1bSJed Brown }
148c4762a1bSJed Brown 
149c4762a1bSJed Brown /*
150c4762a1bSJed Brown      Defines the exact (analytic) solution to the ODE
151c4762a1bSJed Brown */
Solution(TS ts,PetscReal t,Vec U,AppCtx * ctx)152d71ae5a4SJacob Faibussowitsch static PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *ctx)
153d71ae5a4SJacob Faibussowitsch {
154c4762a1bSJed Brown   const PetscScalar *uinit;
155c4762a1bSJed Brown   PetscScalar       *u, d0, q;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   PetscFunctionBegin;
1589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->initialsolution, &uinit));
1599566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
160c4762a1bSJed Brown   d0 = uinit[0] - uinit[1];
161c4762a1bSJed Brown   if (d0 == 0.0) q = ctx->k * t;
162c4762a1bSJed Brown   else q = (1.0 - PetscExpScalar(-ctx->k * t * d0)) / d0;
163c4762a1bSJed Brown   u[0] = uinit[0] / (1.0 + uinit[1] * q);
164c4762a1bSJed Brown   u[1] = u[0] - d0;
165c4762a1bSJed Brown   u[2] = uinit[1] + uinit[2] - u[1];
1669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
1679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->initialsolution, &uinit));
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
169c4762a1bSJed Brown }
170c4762a1bSJed Brown 
main(int argc,char ** argv)171d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
172d71ae5a4SJacob Faibussowitsch {
173c4762a1bSJed Brown   TS                ts;         /* ODE integrator */
174c4762a1bSJed Brown   Vec               U, Udot, R; /* solution, derivative, residual */
175c4762a1bSJed Brown   Mat               A;          /* Jacobian matrix */
176c4762a1bSJed Brown   PetscMPIInt       size;
177c4762a1bSJed Brown   PetscInt          n = 3;
178c4762a1bSJed Brown   AppCtx            ctx;
179c4762a1bSJed Brown   AdolcCtx         *adctx;
180c4762a1bSJed Brown   PetscScalar      *u;
181c4762a1bSJed Brown   const char *const names[] = {"U1", "U2", "U3", NULL};
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184c4762a1bSJed Brown      Initialize program
185c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186327415f7SBarry Smith   PetscFunctionBeginUser;
1879566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
18908401ef6SPierre Jolivet   PetscCheck(size <= 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Only for sequential runs");
1909566063dSJacob Faibussowitsch   PetscCall(PetscNew(&adctx));
1919371c9d4SSatish Balay   adctx->m  = n;
1929371c9d4SSatish Balay   adctx->n  = n;
1939371c9d4SSatish Balay   adctx->p  = n;
194c4762a1bSJed Brown   ctx.adctx = adctx;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197c4762a1bSJed Brown     Create necessary matrix and vectors
198c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1999566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2009566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
2019566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2029566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
203c4762a1bSJed Brown 
2049566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207c4762a1bSJed Brown     Set runtime options
208c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209c4762a1bSJed Brown   ctx.k = .9;
2109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-k", &ctx.k, NULL));
2119566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &ctx.initialsolution));
2129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx.initialsolution, &u));
213c4762a1bSJed Brown   u[0] = 1;
214c4762a1bSJed Brown   u[1] = .7;
215c4762a1bSJed Brown   u[2] = 0;
2169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx.initialsolution, &u));
2179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetVec(NULL, NULL, "-initial", ctx.initialsolution, NULL));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220c4762a1bSJed Brown      Create timestepping solver context
221c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2229566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2239566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2249566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW));
225*8434afd1SBarry Smith   PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunctionPassive, &ctx));
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228c4762a1bSJed Brown      Set initial conditions
229c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2309566063dSJacob Faibussowitsch   PetscCall(Solution(ts, 0, U, &ctx));
2319566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234c4762a1bSJed Brown      Trace just once for each tape
235c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &Udot));
2379566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &R));
2389566063dSJacob Faibussowitsch   PetscCall(IFunctionActive1(ts, 0., U, Udot, R, &ctx));
2399566063dSJacob Faibussowitsch   PetscCall(IFunctionActive2(ts, 0., U, Udot, R, &ctx));
2409566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
2419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Udot));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown      Set Jacobian
245c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246*8434afd1SBarry Smith   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobianFn *)IJacobian, &ctx));
247*8434afd1SBarry Smith   PetscCall(TSSetSolutionFunction(ts, (TSSolutionFn *)Solution, &ctx));
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   {
250c4762a1bSJed Brown     DM    dm;
251c4762a1bSJed Brown     void *ptr;
2529566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
2539566063dSJacob Faibussowitsch     PetscCall(PetscDLSym(NULL, "IFunctionView", &ptr));
2549566063dSJacob Faibussowitsch     PetscCall(PetscDLSym(NULL, "IFunctionLoad", &ptr));
2559566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunctionSerialize(dm, (PetscErrorCode (*)(void *, PetscViewer))IFunctionView, (PetscErrorCode (*)(void **, PetscViewer))IFunctionLoad));
2569566063dSJacob Faibussowitsch     PetscCall(DMTSSetIJacobianSerialize(dm, (PetscErrorCode (*)(void *, PetscViewer))IFunctionView, (PetscErrorCode (*)(void **, PetscViewer))IFunctionLoad));
257c4762a1bSJed Brown   }
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260c4762a1bSJed Brown      Set solver options
261c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2629566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .001));
2639566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 1000));
2649566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 20.0));
2659566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2669566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
2679566063dSJacob Faibussowitsch   PetscCall(TSMonitorLGSetVariableNames(ts, names));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
270c4762a1bSJed Brown      Solve nonlinear system
271c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2729566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
273c4762a1bSJed Brown 
2749566063dSJacob Faibussowitsch   PetscCall(TSView(ts, PETSC_VIEWER_BINARY_WORLD));
275c4762a1bSJed Brown 
276c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
277c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
278c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.initialsolution));
2809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2829566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2839566063dSJacob Faibussowitsch   PetscCall(PetscFree(adctx));
2849566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
285b122ec5aSJacob Faibussowitsch   return 0;
286c4762a1bSJed Brown }
287c4762a1bSJed Brown 
288c4762a1bSJed Brown /*TEST
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   build:
291c4762a1bSJed Brown     requires: double !complex adolc
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   test:
294c4762a1bSJed Brown     suffix: 1
295c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
296c4762a1bSJed Brown     output_file: output/adr_ex1_1.out
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   test:
299c4762a1bSJed Brown     suffix: 2
300c4762a1bSJed Brown     args: -ts_max_steps 1 -snes_test_jacobian
301c4762a1bSJed Brown     output_file: output/adr_ex1_2.out
302c4762a1bSJed Brown 
303c4762a1bSJed Brown TEST*/
304