xref: /petsc/src/ts/tutorials/autodiff/ex16adj.cxx (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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    Concepts: TS^time-dependent nonlinear problems
7c4762a1bSJed Brown    Concepts: TS^van der Pol equation
8c4762a1bSJed Brown    Concepts: TS^adjoint sensitivity analysis
9c4762a1bSJed Brown    Concepts: Automatic differentation using ADOL-C
10c4762a1bSJed Brown    Concepts: Automatic differentation w.r.t. a parameter using ADOL-C
11c4762a1bSJed Brown    Processors: 1
12c4762a1bSJed Brown */
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
15c4762a1bSJed Brown 
16c4762a1bSJed Brown    For documentation on ADOL-C, see
17c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
18c4762a1bSJed Brown */
19c4762a1bSJed Brown /* ------------------------------------------------------------------------
20c4762a1bSJed Brown    See ex16adj for a description of the problem being solved.
21c4762a1bSJed Brown   ------------------------------------------------------------------------- */
22c4762a1bSJed Brown 
23c4762a1bSJed Brown #include <petscts.h>
24c4762a1bSJed Brown #include <petscmat.h>
25c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
26c4762a1bSJed Brown #include <adolc/adolc.h>
27c4762a1bSJed Brown 
28c4762a1bSJed Brown typedef struct _n_User *User;
29c4762a1bSJed Brown struct _n_User {
30c4762a1bSJed Brown   PetscReal mu;
31c4762a1bSJed Brown   PetscReal next_output;
32c4762a1bSJed Brown   PetscReal tprev;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /* Automatic differentiation support */
35c4762a1bSJed Brown   AdolcCtx  *adctx;
36c4762a1bSJed Brown };
37c4762a1bSJed Brown 
38c4762a1bSJed Brown /*
39c4762a1bSJed Brown   'Passive' RHS function, used in residual evaluations during the time integration.
40c4762a1bSJed Brown */
41c4762a1bSJed Brown static PetscErrorCode RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
42c4762a1bSJed Brown {
43c4762a1bSJed Brown   User              user = (User)ctx;
44c4762a1bSJed Brown   PetscScalar       *f;
45c4762a1bSJed Brown   const PetscScalar *x;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   PetscFunctionBeginUser;
485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
50c4762a1bSJed Brown   f[0] = x[1];
51c4762a1bSJed Brown   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
54c4762a1bSJed Brown   PetscFunctionReturn(0);
55c4762a1bSJed Brown }
56c4762a1bSJed Brown 
57c4762a1bSJed Brown /*
58c4762a1bSJed Brown   Trace RHS to mark on tape 1 the dependence of f upon x. This tape is used in generating the
59c4762a1bSJed Brown   Jacobian transform.
60c4762a1bSJed Brown */
61c4762a1bSJed Brown static PetscErrorCode RHSFunctionActive(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
62c4762a1bSJed Brown {
63c4762a1bSJed Brown   User              user = (User)ctx;
64c4762a1bSJed Brown   PetscScalar       *f;
65c4762a1bSJed Brown   const PetscScalar *x;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   adouble           f_a[2]; /* 'active' double for dependent variables */
68c4762a1bSJed Brown   adouble           x_a[2]; /* 'active' double for independent variables */
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBeginUser;
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* Start of active section */
75c4762a1bSJed Brown   trace_on(1);
76c4762a1bSJed Brown   x_a[0] <<= x[0];x_a[1] <<= x[1]; /* Mark independence */
77c4762a1bSJed Brown   f_a[0] = x_a[1];
78c4762a1bSJed Brown   f_a[1] = user->mu*(1.-x_a[0]*x_a[0])*x_a[1]-x_a[0];
79c4762a1bSJed Brown   f_a[0] >>= f[0];f_a[1] >>= f[1]; /* Mark dependence */
80c4762a1bSJed Brown   trace_off();
81c4762a1bSJed Brown   /* End of active section */
82c4762a1bSJed Brown 
835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
85c4762a1bSJed Brown   PetscFunctionReturn(0);
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown /*
89c4762a1bSJed Brown   Trace RHS again to mark on tape 2 the dependence of f upon the parameter mu. This tape is used in
90c4762a1bSJed Brown   generating JacobianP.
91c4762a1bSJed Brown */
92c4762a1bSJed Brown static PetscErrorCode RHSFunctionActiveP(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
93c4762a1bSJed Brown {
94c4762a1bSJed Brown   User              user = (User)ctx;
95c4762a1bSJed Brown   PetscScalar       *f;
96c4762a1bSJed Brown   const PetscScalar *x;
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   adouble           f_a[2];      /* 'active' double for dependent variables */
99c4762a1bSJed Brown   adouble           x_a[2],mu_a; /* 'active' double for independent variables */
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   PetscFunctionBeginUser;
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* Start of active section */
106c4762a1bSJed Brown   trace_on(3);
107c4762a1bSJed Brown   x_a[0] <<= x[0];x_a[1] <<= x[1];mu_a <<= user->mu; /* Mark independence */
108c4762a1bSJed Brown   f_a[0] = x_a[1];
109c4762a1bSJed Brown   f_a[1] = mu_a*(1.-x_a[0]*x_a[0])*x_a[1]-x_a[0];
110c4762a1bSJed Brown   f_a[0] >>= f[0];f_a[1] >>= f[1];                   /* Mark dependence */
111c4762a1bSJed Brown   trace_off();
112c4762a1bSJed Brown   /* End of active section */
113c4762a1bSJed Brown 
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
116c4762a1bSJed Brown   PetscFunctionReturn(0);
117c4762a1bSJed Brown }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown /*
120c4762a1bSJed Brown   Compute the Jacobian w.r.t. x using PETSc-ADOL-C driver for explicit TS.
121c4762a1bSJed Brown */
122c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
123c4762a1bSJed Brown {
124c4762a1bSJed Brown   User              user = (User)ctx;
125a8c08197SHong Zhang   const PetscScalar *x;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   PetscFunctionBeginUser;
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscAdolcComputeRHSJacobian(1,A,x,user->adctx));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
131c4762a1bSJed Brown   PetscFunctionReturn(0);
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown /*
135c4762a1bSJed Brown   Compute the Jacobian w.r.t. mu using PETSc-ADOL-C driver for explicit TS.
136c4762a1bSJed Brown */
137c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
138c4762a1bSJed Brown {
139c4762a1bSJed Brown   User              user = (User)ctx;
140410585f6SHong Zhang   const PetscScalar *x;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   PetscFunctionBeginUser;
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscAdolcComputeRHSJacobianP(3,A,x,&user->mu,user->adctx));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
146c4762a1bSJed Brown   PetscFunctionReturn(0);
147c4762a1bSJed Brown }
148c4762a1bSJed Brown 
149c4762a1bSJed Brown /*
150c4762a1bSJed Brown   Monitor timesteps and use interpolation to output at integer multiples of 0.1
151c4762a1bSJed Brown */
152c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   const PetscScalar *x;
155c4762a1bSJed Brown   PetscReal         tfinal, dt, tprev;
156c4762a1bSJed Brown   User              user = (User)ctx;
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   PetscFunctionBeginUser;
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTimeStep(ts,&dt));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetMaxTime(ts,&tfinal));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetPrevTime(ts,&tprev));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D 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])));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown int main(int argc,char **argv)
170c4762a1bSJed Brown {
171c4762a1bSJed Brown   TS             ts;            /* nonlinear solver */
172c4762a1bSJed Brown   Vec            x;             /* solution, residual vectors */
173c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
174c4762a1bSJed Brown   Mat            Jacp;          /* JacobianP matrix */
175c4762a1bSJed Brown   PetscInt       steps;
176c4762a1bSJed Brown   PetscReal      ftime   = 0.5;
177c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
178c4762a1bSJed Brown   PetscScalar    *x_ptr;
179c4762a1bSJed Brown   PetscMPIInt    size;
180c4762a1bSJed Brown   struct _n_User user;
181c4762a1bSJed Brown   AdolcCtx       *adctx;
182c4762a1bSJed Brown   Vec            lambda[2],mu[2],r;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185c4762a1bSJed Brown      Initialize program
186c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
1885f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1892c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192c4762a1bSJed Brown     Set runtime options and create AdolcCtx
193c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&adctx));
195c4762a1bSJed Brown   user.mu          = 1;
196c4762a1bSJed Brown   user.next_output = 0.0;
197c4762a1bSJed Brown   adctx->m = 2;adctx->n = 2;adctx->p = 2;adctx->num_params = 1;
198c4762a1bSJed Brown   user.adctx = adctx;
199c4762a1bSJed Brown 
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
205c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&x,NULL));
211c4762a1bSJed Brown 
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Jacp));
2135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(Jacp));
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(Jacp));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218c4762a1bSJed Brown      Create timestepping solver context
219c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
2215f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSRK));
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user));
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225c4762a1bSJed Brown      Set initial conditions
226c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(x,&x_ptr));
228c4762a1bSJed Brown   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(x,&x_ptr));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232c4762a1bSJed Brown      Trace just once on each tape and put zeros on Jacobian diagonal
233c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(RHSFunctionActive(ts,0.,x,r,&user));
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(RHSFunctionActiveP(ts,0.,x,r,&user));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(r,0));
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalSet(A,r,INSERT_VALUES));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242c4762a1bSJed Brown      Set RHS Jacobian for the adjoint integration
243c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,&user));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,ftime));
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
247c4762a1bSJed Brown   if (monitor) {
2485f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(ts,Monitor,&user,NULL));
249c4762a1bSJed Brown   }
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.001));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /*
253c4762a1bSJed Brown     Have the TS save its trajectory so that TSAdjointSolve() may be used
254c4762a1bSJed Brown   */
2555f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(ts));
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258c4762a1bSJed Brown      Set runtime options
259c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
263c4762a1bSJed Brown      Solve nonlinear system
264c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&ftime));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&steps));
2685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime));
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272c4762a1bSJed Brown      Start the Adjoint model
273c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&lambda[0],NULL));
2755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&lambda[1],NULL));
276c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
2775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(lambda[0],&x_ptr));
278c4762a1bSJed Brown   x_ptr[0] = 1.0;   x_ptr[1] = 0.0;
2795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(lambda[0],&x_ptr));
2805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(lambda[1],&x_ptr));
281c4762a1bSJed Brown   x_ptr[0] = 0.0;   x_ptr[1] = 1.0;
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(lambda[1],&x_ptr));
283c4762a1bSJed Brown 
2845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(Jacp,&mu[0],NULL));
2855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(Jacp,&mu[1],NULL));
2865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(mu[0],&x_ptr));
287c4762a1bSJed Brown   x_ptr[0] = 0.0;
2885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(mu[0],&x_ptr));
2895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(mu[1],&x_ptr));
290c4762a1bSJed Brown   x_ptr[0] = 0.0;
2915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(mu[1],&x_ptr));
2925f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,2,lambda,mu));
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   /*   Set RHS JacobianP */
2955f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user));
296c4762a1bSJed Brown 
2975f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
298c4762a1bSJed Brown 
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD));
3005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD));
3015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD));
3025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD));
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
305c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
306c4762a1bSJed Brown      are no longer needed.
307c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
3095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Jacp));
3105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
3115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&lambda[0]));
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&lambda[1]));
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&mu[0]));
3145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&mu[1]));
3155f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
3165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(adctx));
317*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
318*b122ec5aSJacob Faibussowitsch   return 0;
319c4762a1bSJed Brown }
320c4762a1bSJed Brown 
321c4762a1bSJed Brown /*TEST
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   build:
324c4762a1bSJed Brown     requires: double !complex adolc
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   test:
327c4762a1bSJed Brown     suffix: 1
328c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
329c4762a1bSJed Brown     output_file: output/ex16adj_1.out
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   test:
332c4762a1bSJed Brown     suffix: 2
333c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5
334c4762a1bSJed Brown     output_file: output/ex16adj_2.out
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   test:
337c4762a1bSJed Brown     suffix: 3
338c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor
339c4762a1bSJed Brown     output_file: output/ex16adj_3.out
340c4762a1bSJed Brown 
341c4762a1bSJed Brown   test:
342c4762a1bSJed Brown     suffix: 4
343c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor -mu 5
344c4762a1bSJed Brown     output_file: output/ex16adj_4.out
345c4762a1bSJed Brown 
346c4762a1bSJed Brown TEST*/
347