xref: /petsc/src/ts/tutorials/autodiff/ex16adj_tl.cxx (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
1c4762a1bSJed Brown static char help[] = "Demonstrates tapeless 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: Tapeless automatic differentiation using ADOL-C
11c4762a1bSJed Brown    Concepts: Automatic differentation w.r.t. a parameter using ADOL-C
12c4762a1bSJed Brown    Processors: 1
13c4762a1bSJed Brown */
14c4762a1bSJed Brown /*
15c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
16c4762a1bSJed Brown 
17c4762a1bSJed Brown    For documentation on ADOL-C, see
18c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
19c4762a1bSJed Brown */
20c4762a1bSJed Brown /* ------------------------------------------------------------------------
21c4762a1bSJed Brown    See ex16adj for a description of the problem being solved.
22c4762a1bSJed Brown   ------------------------------------------------------------------------- */
23c4762a1bSJed Brown 
24c4762a1bSJed Brown #include <petscts.h>
25c4762a1bSJed Brown #include <petscmat.h>
26c4762a1bSJed Brown 
27c4762a1bSJed Brown #define ADOLC_TAPELESS
28c4762a1bSJed Brown #define NUMBER_DIRECTIONS 3
29c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
30c4762a1bSJed Brown #include <adolc/adtl.h>
31c4762a1bSJed Brown using namespace adtl;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown typedef struct _n_User *User;
34c4762a1bSJed Brown struct _n_User {
35c4762a1bSJed Brown   PetscReal mu;
36c4762a1bSJed Brown   PetscReal next_output;
37c4762a1bSJed Brown   PetscReal tprev;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   /* Automatic differentiation support */
40c4762a1bSJed Brown   AdolcCtx  *adctx;
41c4762a1bSJed Brown   Vec       F;
42c4762a1bSJed Brown };
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*
45c4762a1bSJed Brown   Residual evaluation templated, so as to allow for PetscScalar or adouble
46c4762a1bSJed Brown   arguments.
47c4762a1bSJed Brown */
48a8c08197SHong Zhang template <class T> PetscErrorCode EvaluateResidual(const T *x,T mu,T *f)
49c4762a1bSJed Brown {
50c4762a1bSJed Brown   PetscFunctionBegin;
51c4762a1bSJed Brown   f[0] = x[1];
52c4762a1bSJed Brown   f[1] = mu*(1.-x[0]*x[0])*x[1]-x[0];
53c4762a1bSJed Brown   PetscFunctionReturn(0);
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*
57c4762a1bSJed Brown   'Passive' RHS function, used in residual evaluations during the time integration.
58c4762a1bSJed Brown */
59c4762a1bSJed Brown static PetscErrorCode RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   User              user = (User)ctx;
62a8c08197SHong Zhang   PetscScalar       *f;
63a8c08197SHong Zhang   const PetscScalar *x;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   PetscFunctionBeginUser;
669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
689566063dSJacob Faibussowitsch   PetscCall(EvaluateResidual(x,user->mu,f));
699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
71c4762a1bSJed Brown   PetscFunctionReturn(0);
72c4762a1bSJed Brown }
73c4762a1bSJed Brown 
74c4762a1bSJed Brown /*
75c4762a1bSJed Brown   Compute the Jacobian w.r.t. x using tapeless mode of ADOL-C.
76c4762a1bSJed Brown */
77c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
78c4762a1bSJed Brown {
79c4762a1bSJed Brown   User              user = (User)ctx;
80a8c08197SHong Zhang   PetscScalar       **J;
81a8c08197SHong Zhang   const PetscScalar *x;
82c4762a1bSJed Brown   adouble           f_a[2];      /* 'active' double for dependent variables */
83c4762a1bSJed Brown   adouble           x_a[2],mu_a; /* 'active' doubles for independent variables */
84c4762a1bSJed Brown   PetscInt          i,j;
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   PetscFunctionBeginUser;
87c4762a1bSJed Brown   /* Set values for independent variables and parameters */
889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
89c4762a1bSJed Brown   x_a[0].setValue(x[0]);
90c4762a1bSJed Brown   x_a[1].setValue(x[1]);
91c4762a1bSJed Brown   mu_a.setValue(user->mu);
929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* Set seed matrix as 3x3 identity matrix */
95c4762a1bSJed Brown   x_a[0].setADValue(0,1.);x_a[0].setADValue(1,0.);x_a[0].setADValue(2,0.);
96c4762a1bSJed Brown   x_a[1].setADValue(0,0.);x_a[1].setADValue(1,1.);x_a[1].setADValue(2,0.);
97c4762a1bSJed Brown   mu_a.setADValue(0,0.);mu_a.setADValue(1,0.);mu_a.setADValue(2,1.);
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* Evaluate residual (on active variables) */
1009566063dSJacob Faibussowitsch   PetscCall(EvaluateResidual(x_a,mu_a,f_a));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* Extract derivatives */
1039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->adctx->n,&J));
104c4762a1bSJed Brown   J[0] = (PetscScalar*) f_a[0].getADValue();
105c4762a1bSJed Brown   J[1] = (PetscScalar*) f_a[1].getADValue();
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* Set matrix values */
108c4762a1bSJed Brown   for (i=0; i<user->adctx->m; i++) {
109c4762a1bSJed Brown     for (j=0; j<user->adctx->n; j++) {
1109566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,1,&j,&J[i][j],INSERT_VALUES));
111c4762a1bSJed Brown     }
112c4762a1bSJed Brown   }
1139566063dSJacob Faibussowitsch   PetscCall(PetscFree(J));
1149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
116c4762a1bSJed Brown   if (A != B) {
1179566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1189566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
119c4762a1bSJed Brown   }
120c4762a1bSJed Brown   PetscFunctionReturn(0);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown /*
124c4762a1bSJed Brown   Compute the Jacobian w.r.t. mu using tapeless mode of ADOL-C.
125c4762a1bSJed Brown */
126c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
127c4762a1bSJed Brown {
128c4762a1bSJed Brown   User           user = (User)ctx;
129c4762a1bSJed Brown   PetscScalar    **J;
130c4762a1bSJed Brown   PetscScalar    *x;
131c4762a1bSJed Brown   adouble        f_a[2];      /* 'active' double for dependent variables */
132c4762a1bSJed Brown   adouble        x_a[2],mu_a; /* 'active' doubles for independent variables */
133c4762a1bSJed Brown   PetscInt       i,j = 0;
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   PetscFunctionBeginUser;
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* Set values for independent variables and parameters */
1389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,&x));
139c4762a1bSJed Brown   x_a[0].setValue(x[0]);
140c4762a1bSJed Brown   x_a[1].setValue(x[1]);
141c4762a1bSJed Brown   mu_a.setValue(user->mu);
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,&x));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Set seed matrix as 3x3 identity matrix */
145c4762a1bSJed Brown   x_a[0].setADValue(0,1.);x_a[0].setADValue(1,0.);x_a[0].setADValue(2,0.);
146c4762a1bSJed Brown   x_a[1].setADValue(0,0.);x_a[1].setADValue(1,1.);x_a[1].setADValue(2,0.);
147c4762a1bSJed Brown   mu_a.setADValue(0,0.);mu_a.setADValue(1,0.);mu_a.setADValue(2,1.);
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* Evaluate residual (on active variables) */
1509566063dSJacob Faibussowitsch   PetscCall(EvaluateResidual(x_a,mu_a,f_a));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* Extract derivatives */
1539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2,&J));
154c4762a1bSJed Brown   J[0] = (PetscScalar*) f_a[0].getADValue();
155c4762a1bSJed Brown   J[1] = (PetscScalar*) f_a[1].getADValue();
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* Set matrix values */
158c4762a1bSJed Brown   for (i=0; i<user->adctx->m; i++) {
1599566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&i,1,&j,&J[i][user->adctx->n],INSERT_VALUES));
160c4762a1bSJed Brown   }
1619566063dSJacob Faibussowitsch   PetscCall(PetscFree(J));
1629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
164c4762a1bSJed Brown   PetscFunctionReturn(0);
165c4762a1bSJed Brown }
166c4762a1bSJed Brown 
167c4762a1bSJed Brown /*
168c4762a1bSJed Brown   Monitor timesteps and use interpolation to output at integer multiples of 0.1
169c4762a1bSJed Brown */
170c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
171c4762a1bSJed Brown {
172c4762a1bSJed Brown   const PetscScalar *x;
173c4762a1bSJed Brown   PetscReal         tfinal, dt, tprev;
174c4762a1bSJed Brown   User              user = (User)ctx;
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   PetscFunctionBeginUser;
1779566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts,&dt));
1789566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts,&tfinal));
1799566063dSJacob Faibussowitsch   PetscCall(TSGetPrevTime(ts,&tprev));
1809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
1819566063dSJacob Faibussowitsch   PetscCall(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])));
1829566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev));
1839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
184c4762a1bSJed Brown   PetscFunctionReturn(0);
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown int main(int argc,char **argv)
188c4762a1bSJed Brown {
189c4762a1bSJed Brown   TS             ts;            /* nonlinear solver */
190c4762a1bSJed Brown   Vec            x;             /* solution, residual vectors */
191c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
192c4762a1bSJed Brown   Mat            Jacp;          /* JacobianP matrix */
193c4762a1bSJed Brown   PetscInt       steps;
194c4762a1bSJed Brown   PetscReal      ftime   = 0.5;
195c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
196c4762a1bSJed Brown   PetscScalar    *x_ptr;
197c4762a1bSJed Brown   PetscMPIInt    size;
198c4762a1bSJed Brown   struct _n_User user;
199c4762a1bSJed Brown   AdolcCtx       *adctx;
200c4762a1bSJed Brown   Vec            lambda[2],mu[2];
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203c4762a1bSJed Brown      Initialize program
204c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2059566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
2069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
207*08401ef6SPierre Jolivet   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown     Set runtime options and create AdolcCtx
211c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2129566063dSJacob Faibussowitsch   PetscCall(PetscNew(&adctx));
213c4762a1bSJed Brown   user.mu          = 1;
214c4762a1bSJed Brown   user.next_output = 0.0;
215c4762a1bSJed Brown   adctx->m = 2;adctx->n = 2;adctx->p = 2;
216c4762a1bSJed Brown   user.adctx = adctx;
217c4762a1bSJed Brown   adtl::setNumDir(adctx->n+1); /* #indep. variables, plus parameters */
218c4762a1bSJed Brown 
2199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
2209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
224c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2259566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
2269566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2));
2279566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2289566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
2299566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&x,NULL));
230c4762a1bSJed Brown 
2319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&Jacp));
2329566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
2339566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jacp));
2349566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Jacp));
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237c4762a1bSJed Brown      Create timestepping solver context
238c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2399566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
2409566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSRK));
2419566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown      Set initial conditions
245c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&x_ptr));
247c4762a1bSJed Brown   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
2489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&x_ptr));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
251c4762a1bSJed Brown      Set RHS Jacobian for the adjoint integration
252c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2539566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,A,A,RHSJacobian,&user));
2549566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
2559566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
256c4762a1bSJed Brown   if (monitor) {
2579566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,Monitor,&user,NULL));
258c4762a1bSJed Brown   }
2599566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,.001));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   /*
262c4762a1bSJed Brown     Have the TS save its trajectory so that TSAdjointSolve() may be used
263c4762a1bSJed Brown   */
2649566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267c4762a1bSJed Brown      Set runtime options
268c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2699566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272c4762a1bSJed Brown      Solve nonlinear system
273c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2749566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,x));
2759566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
2769566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
2779566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime));
2789566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
281c4762a1bSJed Brown      Start the Adjoint model
282c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2839566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&lambda[0],NULL));
2849566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&lambda[1],NULL));
285c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
2869566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0],&x_ptr));
287c4762a1bSJed Brown   x_ptr[0] = 1.0;   x_ptr[1] = 0.0;
2889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0],&x_ptr));
2899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[1],&x_ptr));
290c4762a1bSJed Brown   x_ptr[0] = 0.0;   x_ptr[1] = 1.0;
2919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[1],&x_ptr));
292c4762a1bSJed Brown 
2939566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp,&mu[0],NULL));
2949566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp,&mu[1],NULL));
2959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0],&x_ptr));
296c4762a1bSJed Brown   x_ptr[0] = 0.0;
2979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0],&x_ptr));
2989566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[1],&x_ptr));
299c4762a1bSJed Brown   x_ptr[0] = 0.0;
3009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[1],&x_ptr));
3019566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,2,lambda,mu));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   /*   Set RHS JacobianP */
3049566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user));
305c4762a1bSJed Brown 
3069566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
307c4762a1bSJed Brown 
3089566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD));
3099566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD));
3109566063dSJacob Faibussowitsch   PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD));
3119566063dSJacob Faibussowitsch   PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD));
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
314c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
315c4762a1bSJed Brown      are no longer needed.
316c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jacp));
3199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
3219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[1]));
3229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
3239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[1]));
3249566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3259566063dSJacob Faibussowitsch   PetscCall(PetscFree(adctx));
3269566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
327b122ec5aSJacob Faibussowitsch   return 0;
328c4762a1bSJed Brown }
329c4762a1bSJed Brown 
330c4762a1bSJed Brown /*TEST
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   build:
333c4762a1bSJed Brown     requires: double !complex adolc
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   test:
336c4762a1bSJed Brown     suffix: 1
337c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
338c4762a1bSJed Brown     output_file: output/ex16adj_tl_1.out
339c4762a1bSJed Brown 
340c4762a1bSJed Brown   test:
341c4762a1bSJed Brown     suffix: 2
342c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5
343c4762a1bSJed Brown     output_file: output/ex16adj_tl_2.out
344c4762a1bSJed Brown 
345c4762a1bSJed Brown   test:
346c4762a1bSJed Brown     suffix: 3
347c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor
348c4762a1bSJed Brown     output_file: output/ex16adj_tl_3.out
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   test:
351c4762a1bSJed Brown     suffix: 4
352c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor -mu 5
353c4762a1bSJed Brown     output_file: output/ex16adj_tl_4.out
354c4762a1bSJed Brown 
355c4762a1bSJed Brown TEST*/
356