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