xref: /petsc/src/ts/tutorials/autodiff/adr_ex1.cxx (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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    Concepts: TS^time-dependent nonlinear problems
5c4762a1bSJed Brown    Concepts: Automatic differentiation using ADOL-C
6c4762a1bSJed Brown    Processors: 1
7c4762a1bSJed Brown */
8c4762a1bSJed Brown /*
9c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    For documentation on ADOL-C, see
12c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
13c4762a1bSJed Brown */
14c4762a1bSJed Brown /* ------------------------------------------------------------------------
15c4762a1bSJed Brown   See ../advection-diffusion-reaction/ex1 for a description of the problem
16c4762a1bSJed Brown   ------------------------------------------------------------------------- */
17c4762a1bSJed Brown #include <petscts.h>
18c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
19c4762a1bSJed Brown #include <adolc/adolc.h>
20c4762a1bSJed Brown 
21c4762a1bSJed Brown typedef struct {
22c4762a1bSJed Brown   PetscScalar k;
23c4762a1bSJed Brown   Vec         initialsolution;
24c4762a1bSJed Brown   AdolcCtx    *adctx; /* Automatic differentiation support */
25c4762a1bSJed Brown } AppCtx;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown PetscErrorCode IFunctionView(AppCtx *ctx,PetscViewer v)
28c4762a1bSJed Brown {
29c4762a1bSJed Brown   PetscErrorCode ierr;
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   PetscFunctionBegin;
32c4762a1bSJed Brown   ierr = PetscViewerBinaryWrite(v,&ctx->k,1,PETSC_SCALAR);CHKERRQ(ierr);
33c4762a1bSJed Brown   PetscFunctionReturn(0);
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
36c4762a1bSJed Brown PetscErrorCode IFunctionLoad(AppCtx **ctx,PetscViewer v)
37c4762a1bSJed Brown {
38c4762a1bSJed Brown   PetscErrorCode ierr;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   PetscFunctionBegin;
41c4762a1bSJed Brown   ierr = PetscNew(ctx);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = PetscViewerBinaryRead(v,&(*ctx)->k,1,NULL,PETSC_SCALAR);CHKERRQ(ierr);
43c4762a1bSJed Brown   PetscFunctionReturn(0);
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /*
47c4762a1bSJed Brown   Defines the ODE passed to the ODE solver
48c4762a1bSJed Brown */
49c4762a1bSJed Brown PetscErrorCode IFunctionPassive(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
50c4762a1bSJed Brown {
51c4762a1bSJed Brown   PetscErrorCode    ierr;
52c4762a1bSJed Brown   PetscScalar       *f;
53c4762a1bSJed Brown   const PetscScalar *u,*udot;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   PetscFunctionBegin;
56c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
57c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
58c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
60c4762a1bSJed Brown   f[0] = udot[0] + ctx->k*u[0]*u[1];
61c4762a1bSJed Brown   f[1] = udot[1] + ctx->k*u[0]*u[1];
62c4762a1bSJed Brown   f[2] = udot[2] - ctx->k*u[0]*u[1];
63c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
64c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
66c4762a1bSJed Brown   PetscFunctionReturn(0);
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
69c4762a1bSJed Brown /*
70c4762a1bSJed Brown   'Active' ADOL-C annotated version, marking dependence upon u.
71c4762a1bSJed Brown */
72c4762a1bSJed Brown PetscErrorCode IFunctionActive1(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
73c4762a1bSJed Brown {
74c4762a1bSJed Brown   PetscErrorCode    ierr;
75c4762a1bSJed Brown   PetscScalar       *f;
76c4762a1bSJed Brown   const PetscScalar *u,*udot;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   adouble           f_a[3]; /* 'active' double for dependent variables */
79c4762a1bSJed Brown   adouble           u_a[3]; /* 'active' double for independent variables */
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   PetscFunctionBegin;
82c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
83c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
84c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
85c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Start of active section */
88c4762a1bSJed Brown   trace_on(1);
89c4762a1bSJed Brown   u_a[0] <<= u[0]; u_a[1] <<= u[1]; u_a[2] <<= u[2]; /* Mark independence */
90c4762a1bSJed Brown   f_a[0] = udot[0] + ctx->k*u_a[0]*u_a[1];
91c4762a1bSJed Brown   f_a[1] = udot[1] + ctx->k*u_a[0]*u_a[1];
92c4762a1bSJed Brown   f_a[2] = udot[2] - ctx->k*u_a[0]*u_a[1];
93c4762a1bSJed Brown   f_a[0] >>= f[0]; f_a[1] >>= f[1]; f_a[2] >>= f[2]; /* Mark dependence */
94c4762a1bSJed Brown   trace_off();
95c4762a1bSJed Brown   /* End of active section */
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
100c4762a1bSJed Brown   PetscFunctionReturn(0);
101c4762a1bSJed Brown }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown /*
104c4762a1bSJed Brown   'Active' ADOL-C annotated version, marking dependence upon udot.
105c4762a1bSJed Brown */
106c4762a1bSJed Brown PetscErrorCode IFunctionActive2(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
107c4762a1bSJed Brown {
108c4762a1bSJed Brown   PetscErrorCode    ierr;
109c4762a1bSJed Brown   PetscScalar       *f;
110c4762a1bSJed Brown   const PetscScalar *u,*udot;
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   adouble           f_a[3];    /* 'active' double for dependent variables */
113c4762a1bSJed Brown   adouble           udot_a[3]; /* 'active' double for independent variables */
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   PetscFunctionBegin;
116c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
117c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* Start of active section */
122c4762a1bSJed Brown   trace_on(2);
123c4762a1bSJed Brown   udot_a[0] <<= udot[0]; udot_a[1] <<= udot[1]; udot_a[2] <<= udot[2]; /* Mark independence */
124c4762a1bSJed Brown   f_a[0] = udot_a[0] + ctx->k*u[0]*u[1];
125c4762a1bSJed Brown   f_a[1] = udot_a[1] + ctx->k*u[0]*u[1];
126c4762a1bSJed Brown   f_a[2] = udot_a[2] - ctx->k*u[0]*u[1];
127c4762a1bSJed Brown   f_a[0] >>= f[0]; f_a[1] >>= f[1]; f_a[2] >>= f[2];                   /* Mark dependence */
128c4762a1bSJed Brown   trace_off();
129c4762a1bSJed Brown   /* End of active section */
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
132c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
133c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
134c4762a1bSJed Brown   PetscFunctionReturn(0);
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown /*
138c4762a1bSJed Brown  Defines the Jacobian of the ODE passed to the ODE solver, using the PETSc-ADOL-C driver for
139c4762a1bSJed Brown  implicit TS.
140c4762a1bSJed Brown */
141c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
142c4762a1bSJed Brown {
143c4762a1bSJed Brown   PetscErrorCode    ierr;
144c4762a1bSJed Brown   AppCtx            *appctx = (AppCtx*)ctx;
145410585f6SHong Zhang   const PetscScalar *u;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   PetscFunctionBegin;
148a8c08197SHong Zhang   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
149c4762a1bSJed Brown   ierr = PetscAdolcComputeIJacobian(1,2,A,u,a,appctx->adctx);CHKERRQ(ierr);
150a8c08197SHong Zhang   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
151c4762a1bSJed Brown   PetscFunctionReturn(0);
152c4762a1bSJed Brown }
153c4762a1bSJed Brown 
154c4762a1bSJed Brown /*
155c4762a1bSJed Brown      Defines the exact (analytic) solution to the ODE
156c4762a1bSJed Brown */
157c4762a1bSJed Brown static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
158c4762a1bSJed Brown {
159c4762a1bSJed Brown   PetscErrorCode    ierr;
160c4762a1bSJed Brown   const PetscScalar *uinit;
161c4762a1bSJed Brown   PetscScalar       *u,d0,q;
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   PetscFunctionBegin;
164c4762a1bSJed Brown   ierr = VecGetArrayRead(ctx->initialsolution,&uinit);CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
166c4762a1bSJed Brown   d0   = uinit[0] - uinit[1];
167c4762a1bSJed Brown   if (d0 == 0.0) q = ctx->k*t;
168c4762a1bSJed Brown   else q = (1.0 - PetscExpScalar(-ctx->k*t*d0))/d0;
169c4762a1bSJed Brown   u[0] = uinit[0]/(1.0 + uinit[1]*q);
170c4762a1bSJed Brown   u[1] = u[0] - d0;
171c4762a1bSJed Brown   u[2] = uinit[1] + uinit[2] - u[1];
172c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
173c4762a1bSJed Brown   ierr = VecRestoreArrayRead(ctx->initialsolution,&uinit);CHKERRQ(ierr);
174c4762a1bSJed Brown   PetscFunctionReturn(0);
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown int main(int argc,char **argv)
178c4762a1bSJed Brown {
179c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
180c4762a1bSJed Brown   Vec            U,Udot,R;      /* solution, derivative, residual */
181c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
182c4762a1bSJed Brown   PetscErrorCode ierr;
183c4762a1bSJed Brown   PetscMPIInt    size;
184c4762a1bSJed Brown   PetscInt       n = 3;
185c4762a1bSJed Brown   AppCtx         ctx;
186c4762a1bSJed Brown   AdolcCtx       *adctx;
187c4762a1bSJed Brown   PetscScalar    *u;
188c4762a1bSJed Brown   const char     * const names[] = {"U1","U2","U3",NULL};
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191c4762a1bSJed Brown      Initialize program
192c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
194ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
195*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
196c4762a1bSJed Brown   ierr = PetscNew(&adctx);CHKERRQ(ierr);
197c4762a1bSJed Brown   adctx->m = n;adctx->n = n;adctx->p = n;
198c4762a1bSJed Brown   ctx.adctx = adctx;
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201c4762a1bSJed Brown     Create necessary matrix and vectors
202c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
204c4762a1bSJed Brown   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
205c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
206c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211c4762a1bSJed Brown     Set runtime options
212c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213c4762a1bSJed Brown   ctx.k = .9;
214c4762a1bSJed Brown   ierr  = PetscOptionsGetScalar(NULL,NULL,"-k",&ctx.k,NULL);CHKERRQ(ierr);
215c4762a1bSJed Brown   ierr  = VecDuplicate(U,&ctx.initialsolution);CHKERRQ(ierr);
216c4762a1bSJed Brown   ierr  = VecGetArray(ctx.initialsolution,&u);CHKERRQ(ierr);
217c4762a1bSJed Brown   u[0]  = 1;
218c4762a1bSJed Brown   u[1]  = .7;
219c4762a1bSJed Brown   u[2]  = 0;
220c4762a1bSJed Brown   ierr  = VecRestoreArray(ctx.initialsolution,&u);CHKERRQ(ierr);
221c4762a1bSJed Brown   ierr  = PetscOptionsGetVec(NULL,NULL,"-initial",ctx.initialsolution,NULL);CHKERRQ(ierr);
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
224c4762a1bSJed Brown      Create timestepping solver context
225c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
226c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
227c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
228c4762a1bSJed Brown   ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr);
229c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunctionPassive,&ctx);CHKERRQ(ierr);
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232c4762a1bSJed Brown      Set initial conditions
233c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
234c4762a1bSJed Brown   ierr = Solution(ts,0,U,&ctx);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238c4762a1bSJed Brown      Trace just once for each tape
239c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
240c4762a1bSJed Brown   ierr = VecDuplicate(U,&Udot);CHKERRQ(ierr);
241c4762a1bSJed Brown   ierr = VecDuplicate(U,&R);CHKERRQ(ierr);
242c4762a1bSJed Brown   ierr = IFunctionActive1(ts,0.,U,Udot,R,&ctx);CHKERRQ(ierr);
243c4762a1bSJed Brown   ierr = IFunctionActive2(ts,0.,U,Udot,R,&ctx);CHKERRQ(ierr);
244c4762a1bSJed Brown   ierr = VecDestroy(&R);CHKERRQ(ierr);
245c4762a1bSJed Brown   ierr = VecDestroy(&Udot);CHKERRQ(ierr);
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248c4762a1bSJed Brown      Set Jacobian
249c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
250c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = TSSetSolutionFunction(ts,(TSSolutionFunction)Solution,&ctx);CHKERRQ(ierr);
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   {
254c4762a1bSJed Brown     DM   dm;
255c4762a1bSJed Brown     void *ptr;
256c4762a1bSJed Brown     ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
257c4762a1bSJed Brown     ierr = PetscDLSym(NULL,"IFunctionView",&ptr);CHKERRQ(ierr);
258c4762a1bSJed Brown     ierr = PetscDLSym(NULL,"IFunctionLoad",&ptr);CHKERRQ(ierr);
259c4762a1bSJed Brown     ierr = DMTSSetIFunctionSerialize(dm,(PetscErrorCode (*)(void*,PetscViewer))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer))IFunctionLoad);CHKERRQ(ierr);
260c4762a1bSJed Brown     ierr = DMTSSetIJacobianSerialize(dm,(PetscErrorCode (*)(void*,PetscViewer))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer))IFunctionLoad);CHKERRQ(ierr);
261c4762a1bSJed Brown   }
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
264c4762a1bSJed Brown      Set solver options
265c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
266c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr);
267c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,1000);CHKERRQ(ierr);
268c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,20.0);CHKERRQ(ierr);
269c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
270c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
271c4762a1bSJed Brown   ierr = TSMonitorLGSetVariableNames(ts,names);CHKERRQ(ierr);
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
274c4762a1bSJed Brown      Solve nonlinear system
275c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
276c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   ierr = TSView(ts,PETSC_VIEWER_BINARY_WORLD);CHKERRQ(ierr);
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
281c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
282c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
283c4762a1bSJed Brown   ierr = VecDestroy(&ctx.initialsolution);CHKERRQ(ierr);
284c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
285c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
286c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
287c4762a1bSJed Brown   ierr = PetscFree(adctx);CHKERRQ(ierr);
288c4762a1bSJed Brown   ierr = PetscFinalize();
289c4762a1bSJed Brown   return ierr;
290c4762a1bSJed Brown }
291c4762a1bSJed Brown 
292c4762a1bSJed Brown /*TEST
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   build:
295c4762a1bSJed Brown     requires: double !complex adolc
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   test:
298c4762a1bSJed Brown     suffix: 1
299c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
300c4762a1bSJed Brown     output_file: output/adr_ex1_1.out
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   test:
303c4762a1bSJed Brown     suffix: 2
304c4762a1bSJed Brown     args: -ts_max_steps 1 -snes_test_jacobian
305c4762a1bSJed Brown     output_file: output/adr_ex1_2.out
306c4762a1bSJed Brown 
307c4762a1bSJed Brown TEST*/
308