xref: /petsc/src/ts/tutorials/autodiff/ex16adj_tl.cxx (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Demonstrates tapeless automatic Jacobian generation using ADOL-C for an adjoint sensitivity analysis of the van der Pol equation.\n\
2*c4762a1bSJed Brown Input parameters include:\n\
3*c4762a1bSJed Brown       -mu : stiffness parameter\n\n";
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown /*
6*c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
7*c4762a1bSJed Brown    Concepts: TS^van der Pol equation
8*c4762a1bSJed Brown    Concepts: TS^adjoint sensitivity analysis
9*c4762a1bSJed Brown    Concepts: Automatic differentation using ADOL-C
10*c4762a1bSJed Brown    Concepts: Tapeless automatic differentiation using ADOL-C
11*c4762a1bSJed Brown    Concepts: Automatic differentation w.r.t. a parameter using ADOL-C
12*c4762a1bSJed Brown    Processors: 1
13*c4762a1bSJed Brown */
14*c4762a1bSJed Brown /*
15*c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown    For documentation on ADOL-C, see
18*c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
19*c4762a1bSJed Brown */
20*c4762a1bSJed Brown /* ------------------------------------------------------------------------
21*c4762a1bSJed Brown    See ex16adj for a description of the problem being solved.
22*c4762a1bSJed Brown   ------------------------------------------------------------------------- */
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown #include <petscts.h>
25*c4762a1bSJed Brown #include <petscmat.h>
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown #define ADOLC_TAPELESS
28*c4762a1bSJed Brown #define NUMBER_DIRECTIONS 3
29*c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
30*c4762a1bSJed Brown #include <adolc/adtl.h>
31*c4762a1bSJed Brown using namespace adtl;
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown typedef struct _n_User *User;
34*c4762a1bSJed Brown struct _n_User {
35*c4762a1bSJed Brown   PetscReal mu;
36*c4762a1bSJed Brown   PetscReal next_output;
37*c4762a1bSJed Brown   PetscReal tprev;
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown   /* Automatic differentiation support */
40*c4762a1bSJed Brown   AdolcCtx  *adctx;
41*c4762a1bSJed Brown   Vec       F;
42*c4762a1bSJed Brown };
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown /*
45*c4762a1bSJed Brown   Residual evaluation templated, so as to allow for PetscScalar or adouble
46*c4762a1bSJed Brown   arguments.
47*c4762a1bSJed Brown */
48*c4762a1bSJed Brown template <class T> PetscErrorCode EvaluateResidual(T *x,T mu,T *f)
49*c4762a1bSJed Brown {
50*c4762a1bSJed Brown   PetscFunctionBegin;
51*c4762a1bSJed Brown   f[0] = x[1];
52*c4762a1bSJed Brown   f[1] = mu*(1.-x[0]*x[0])*x[1]-x[0];
53*c4762a1bSJed Brown   PetscFunctionReturn(0);
54*c4762a1bSJed Brown }
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown /*
57*c4762a1bSJed Brown   'Passive' RHS function, used in residual evaluations during the time integration.
58*c4762a1bSJed Brown */
59*c4762a1bSJed Brown static PetscErrorCode RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
60*c4762a1bSJed Brown {
61*c4762a1bSJed Brown   PetscErrorCode    ierr;
62*c4762a1bSJed Brown   User              user = (User)ctx;
63*c4762a1bSJed Brown   PetscScalar       *f,*x;
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown   PetscFunctionBeginUser;
66*c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = EvaluateResidual(x,user->mu,f);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
71*c4762a1bSJed Brown   PetscFunctionReturn(0);
72*c4762a1bSJed Brown }
73*c4762a1bSJed Brown 
74*c4762a1bSJed Brown /*
75*c4762a1bSJed Brown   Compute the Jacobian w.r.t. x using tapeless mode of ADOL-C.
76*c4762a1bSJed Brown */
77*c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
78*c4762a1bSJed Brown {
79*c4762a1bSJed Brown   PetscErrorCode ierr;
80*c4762a1bSJed Brown   User           user = (User)ctx;
81*c4762a1bSJed Brown   PetscScalar    *x,**J;
82*c4762a1bSJed Brown   adouble        f_a[2];      /* 'active' double for dependent variables */
83*c4762a1bSJed Brown   adouble        x_a[2],mu_a; /* 'active' doubles for independent variables */
84*c4762a1bSJed Brown   PetscInt       i,j;
85*c4762a1bSJed Brown 
86*c4762a1bSJed Brown   PetscFunctionBeginUser;
87*c4762a1bSJed Brown 
88*c4762a1bSJed Brown   /* Set values for independent variables and parameters */
89*c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
90*c4762a1bSJed Brown   x_a[0].setValue(x[0]);
91*c4762a1bSJed Brown   x_a[1].setValue(x[1]);
92*c4762a1bSJed Brown   mu_a.setValue(user->mu);
93*c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   /* Set seed matrix as 3x3 identity matrix */
96*c4762a1bSJed Brown   x_a[0].setADValue(0,1.);x_a[0].setADValue(1,0.);x_a[0].setADValue(2,0.);
97*c4762a1bSJed Brown   x_a[1].setADValue(0,0.);x_a[1].setADValue(1,1.);x_a[1].setADValue(2,0.);
98*c4762a1bSJed Brown   mu_a.setADValue(0,0.);mu_a.setADValue(1,0.);mu_a.setADValue(2,1.);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   /* Evaluate residual (on active variables) */
101*c4762a1bSJed Brown   ierr = EvaluateResidual(x_a,mu_a,f_a);CHKERRQ(ierr);
102*c4762a1bSJed Brown 
103*c4762a1bSJed Brown   /* Extract derivatives */
104*c4762a1bSJed Brown   ierr = PetscMalloc1(user->adctx->n,&J);
105*c4762a1bSJed Brown   J[0] = (PetscScalar*) f_a[0].getADValue();
106*c4762a1bSJed Brown   J[1] = (PetscScalar*) f_a[1].getADValue();
107*c4762a1bSJed Brown 
108*c4762a1bSJed Brown   /* Set matrix values */
109*c4762a1bSJed Brown   for (i=0; i<user->adctx->m; i++) {
110*c4762a1bSJed Brown     for (j=0; j<user->adctx->n; j++) {
111*c4762a1bSJed Brown       ierr = MatSetValues(A,1,&i,1,&j,&J[i][j],INSERT_VALUES);CHKERRQ(ierr);
112*c4762a1bSJed Brown     }
113*c4762a1bSJed Brown   }
114*c4762a1bSJed Brown   ierr = PetscFree(J);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
117*c4762a1bSJed Brown   if (A != B) {
118*c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
119*c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
120*c4762a1bSJed Brown   }
121*c4762a1bSJed Brown   PetscFunctionReturn(0);
122*c4762a1bSJed Brown }
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown /*
125*c4762a1bSJed Brown   Compute the Jacobian w.r.t. mu using tapeless mode of ADOL-C.
126*c4762a1bSJed Brown */
127*c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
128*c4762a1bSJed Brown {
129*c4762a1bSJed Brown   PetscErrorCode ierr;
130*c4762a1bSJed Brown   User           user = (User)ctx;
131*c4762a1bSJed Brown   PetscScalar    **J;
132*c4762a1bSJed Brown   PetscScalar    *x;
133*c4762a1bSJed Brown   adouble        f_a[2];      /* 'active' double for dependent variables */
134*c4762a1bSJed Brown   adouble        x_a[2],mu_a; /* 'active' doubles for independent variables */
135*c4762a1bSJed Brown   PetscInt       i,j = 0;
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown   PetscFunctionBeginUser;
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown   /* Set values for independent variables and parameters */
140*c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
141*c4762a1bSJed Brown   x_a[0].setValue(x[0]);
142*c4762a1bSJed Brown   x_a[1].setValue(x[1]);
143*c4762a1bSJed Brown   mu_a.setValue(user->mu);
144*c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
145*c4762a1bSJed Brown 
146*c4762a1bSJed Brown   /* Set seed matrix as 3x3 identity matrix */
147*c4762a1bSJed Brown   x_a[0].setADValue(0,1.);x_a[0].setADValue(1,0.);x_a[0].setADValue(2,0.);
148*c4762a1bSJed Brown   x_a[1].setADValue(0,0.);x_a[1].setADValue(1,1.);x_a[1].setADValue(2,0.);
149*c4762a1bSJed Brown   mu_a.setADValue(0,0.);mu_a.setADValue(1,0.);mu_a.setADValue(2,1.);
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown   /* Evaluate residual (on active variables) */
152*c4762a1bSJed Brown   ierr = EvaluateResidual(x_a,mu_a,f_a);CHKERRQ(ierr);
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown   /* Extract derivatives */
155*c4762a1bSJed Brown   ierr = PetscMalloc1(2,&J);CHKERRQ(ierr);
156*c4762a1bSJed Brown   J[0] = (PetscScalar*) f_a[0].getADValue();
157*c4762a1bSJed Brown   J[1] = (PetscScalar*) f_a[1].getADValue();
158*c4762a1bSJed Brown 
159*c4762a1bSJed Brown   /* Set matrix values */
160*c4762a1bSJed Brown   for (i=0; i<user->adctx->m; i++) {
161*c4762a1bSJed Brown     ierr = MatSetValues(A,1,&i,1,&j,&J[i][user->adctx->n],INSERT_VALUES);CHKERRQ(ierr);
162*c4762a1bSJed Brown   }
163*c4762a1bSJed Brown   ierr = PetscFree(J);CHKERRQ(ierr);
164*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166*c4762a1bSJed Brown   PetscFunctionReturn(0);
167*c4762a1bSJed Brown }
168*c4762a1bSJed Brown 
169*c4762a1bSJed Brown /*
170*c4762a1bSJed Brown   Monitor timesteps and use interpolation to output at integer multiples of 0.1
171*c4762a1bSJed Brown */
172*c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
173*c4762a1bSJed Brown {
174*c4762a1bSJed Brown   PetscErrorCode    ierr;
175*c4762a1bSJed Brown   const PetscScalar *x;
176*c4762a1bSJed Brown   PetscReal         tfinal, dt, tprev;
177*c4762a1bSJed Brown   User              user = (User)ctx;
178*c4762a1bSJed Brown 
179*c4762a1bSJed Brown   PetscFunctionBeginUser;
180*c4762a1bSJed Brown   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
181*c4762a1bSJed Brown   ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr);
182*c4762a1bSJed Brown   ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = 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]));CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
187*c4762a1bSJed Brown   PetscFunctionReturn(0);
188*c4762a1bSJed Brown }
189*c4762a1bSJed Brown 
190*c4762a1bSJed Brown int main(int argc,char **argv)
191*c4762a1bSJed Brown {
192*c4762a1bSJed Brown   TS             ts;            /* nonlinear solver */
193*c4762a1bSJed Brown   Vec            x;             /* solution, residual vectors */
194*c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
195*c4762a1bSJed Brown   Mat            Jacp;          /* JacobianP matrix */
196*c4762a1bSJed Brown   PetscInt       steps;
197*c4762a1bSJed Brown   PetscReal      ftime   = 0.5;
198*c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
199*c4762a1bSJed Brown   PetscScalar    *x_ptr;
200*c4762a1bSJed Brown   PetscMPIInt    size;
201*c4762a1bSJed Brown   struct _n_User user;
202*c4762a1bSJed Brown   AdolcCtx       *adctx;
203*c4762a1bSJed Brown   PetscErrorCode ierr;
204*c4762a1bSJed Brown   Vec            lambda[2],mu[2];
205*c4762a1bSJed Brown 
206*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207*c4762a1bSJed Brown      Initialize program
208*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
210*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
211*c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_MPI_WRONG_SIZE,"This is a uniprocessor example only!");
212*c4762a1bSJed Brown 
213*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214*c4762a1bSJed Brown     Set runtime options and create AdolcCtx
215*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216*c4762a1bSJed Brown   ierr = PetscNew(&adctx);CHKERRQ(ierr);
217*c4762a1bSJed Brown   user.mu          = 1;
218*c4762a1bSJed Brown   user.next_output = 0.0;
219*c4762a1bSJed Brown   adctx->m = 2;adctx->n = 2;adctx->p = 2;
220*c4762a1bSJed Brown   user.adctx = adctx;
221*c4762a1bSJed Brown   adtl::setNumDir(adctx->n+1); /* #indep. variables, plus parameters */
222*c4762a1bSJed Brown 
223*c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr);
224*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr);
225*c4762a1bSJed Brown 
226*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227*c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
228*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
230*c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
231*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
232*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
233*c4762a1bSJed Brown   ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr);
234*c4762a1bSJed Brown 
235*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr);
237*c4762a1bSJed Brown   ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr);
238*c4762a1bSJed Brown   ierr = MatSetUp(Jacp);CHKERRQ(ierr);
239*c4762a1bSJed Brown 
240*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241*c4762a1bSJed Brown      Create timestepping solver context
242*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
244*c4762a1bSJed Brown   ierr = TSSetType(ts,TSRK);CHKERRQ(ierr);
245*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user);CHKERRQ(ierr);
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248*c4762a1bSJed Brown      Set initial conditions
249*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
250*c4762a1bSJed Brown   ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr);
251*c4762a1bSJed Brown   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
252*c4762a1bSJed Brown   ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr);
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
255*c4762a1bSJed Brown      Set RHS Jacobian for the adjoint integration
256*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
257*c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);CHKERRQ(ierr);
258*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
259*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
260*c4762a1bSJed Brown   if (monitor) {
261*c4762a1bSJed Brown     ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr);
262*c4762a1bSJed Brown   }
263*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr);
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown   /*
266*c4762a1bSJed Brown     Have the TS save its trajectory so that TSAdjointSolve() may be used
267*c4762a1bSJed Brown   */
268*c4762a1bSJed Brown   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
271*c4762a1bSJed Brown      Set runtime options
272*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
273*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
274*c4762a1bSJed Brown 
275*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276*c4762a1bSJed Brown      Solve nonlinear system
277*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
278*c4762a1bSJed Brown   ierr = TSSolve(ts,x);CHKERRQ(ierr);
279*c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
280*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
281*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);CHKERRQ(ierr);
282*c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
283*c4762a1bSJed Brown 
284*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
285*c4762a1bSJed Brown      Start the Adjoint model
286*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
287*c4762a1bSJed Brown   ierr = MatCreateVecs(A,&lambda[0],NULL);CHKERRQ(ierr);
288*c4762a1bSJed Brown   ierr = MatCreateVecs(A,&lambda[1],NULL);CHKERRQ(ierr);
289*c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
290*c4762a1bSJed Brown   ierr = VecGetArray(lambda[0],&x_ptr);CHKERRQ(ierr);
291*c4762a1bSJed Brown   x_ptr[0] = 1.0;   x_ptr[1] = 0.0;
292*c4762a1bSJed Brown   ierr = VecRestoreArray(lambda[0],&x_ptr);CHKERRQ(ierr);
293*c4762a1bSJed Brown   ierr = VecGetArray(lambda[1],&x_ptr);CHKERRQ(ierr);
294*c4762a1bSJed Brown   x_ptr[0] = 0.0;   x_ptr[1] = 1.0;
295*c4762a1bSJed Brown   ierr = VecRestoreArray(lambda[1],&x_ptr);CHKERRQ(ierr);
296*c4762a1bSJed Brown 
297*c4762a1bSJed Brown   ierr = MatCreateVecs(Jacp,&mu[0],NULL);CHKERRQ(ierr);
298*c4762a1bSJed Brown   ierr = MatCreateVecs(Jacp,&mu[1],NULL);CHKERRQ(ierr);
299*c4762a1bSJed Brown   ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr);
300*c4762a1bSJed Brown   x_ptr[0] = 0.0;
301*c4762a1bSJed Brown   ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr);
302*c4762a1bSJed Brown   ierr = VecGetArray(mu[1],&x_ptr);CHKERRQ(ierr);
303*c4762a1bSJed Brown   x_ptr[0] = 0.0;
304*c4762a1bSJed Brown   ierr = VecRestoreArray(mu[1],&x_ptr);CHKERRQ(ierr);
305*c4762a1bSJed Brown   ierr = TSSetCostGradients(ts,2,lambda,mu);CHKERRQ(ierr);
306*c4762a1bSJed Brown 
307*c4762a1bSJed Brown 
308*c4762a1bSJed Brown   /*   Set RHS JacobianP */
309*c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user);CHKERRQ(ierr);
310*c4762a1bSJed Brown 
311*c4762a1bSJed Brown   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
312*c4762a1bSJed Brown 
313*c4762a1bSJed Brown   ierr = VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
314*c4762a1bSJed Brown   ierr = VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
315*c4762a1bSJed Brown   ierr = VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
316*c4762a1bSJed Brown   ierr = VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
317*c4762a1bSJed Brown 
318*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
319*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
320*c4762a1bSJed Brown      are no longer needed.
321*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
322*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
323*c4762a1bSJed Brown   ierr = MatDestroy(&Jacp);CHKERRQ(ierr);
324*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
325*c4762a1bSJed Brown   ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
326*c4762a1bSJed Brown   ierr = VecDestroy(&lambda[1]);CHKERRQ(ierr);
327*c4762a1bSJed Brown   ierr = VecDestroy(&mu[0]);CHKERRQ(ierr);
328*c4762a1bSJed Brown   ierr = VecDestroy(&mu[1]);CHKERRQ(ierr);
329*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
330*c4762a1bSJed Brown   ierr = PetscFree(adctx);CHKERRQ(ierr);
331*c4762a1bSJed Brown   ierr = PetscFinalize();
332*c4762a1bSJed Brown   return ierr;
333*c4762a1bSJed Brown }
334*c4762a1bSJed Brown 
335*c4762a1bSJed Brown /*TEST
336*c4762a1bSJed Brown 
337*c4762a1bSJed Brown   build:
338*c4762a1bSJed Brown     requires: double !complex adolc
339*c4762a1bSJed Brown 
340*c4762a1bSJed Brown   test:
341*c4762a1bSJed Brown     suffix: 1
342*c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
343*c4762a1bSJed Brown     output_file: output/ex16adj_tl_1.out
344*c4762a1bSJed Brown 
345*c4762a1bSJed Brown   test:
346*c4762a1bSJed Brown     suffix: 2
347*c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5
348*c4762a1bSJed Brown     output_file: output/ex16adj_tl_2.out
349*c4762a1bSJed Brown 
350*c4762a1bSJed Brown   test:
351*c4762a1bSJed Brown     suffix: 3
352*c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor
353*c4762a1bSJed Brown     output_file: output/ex16adj_tl_3.out
354*c4762a1bSJed Brown 
355*c4762a1bSJed Brown   test:
356*c4762a1bSJed Brown     suffix: 4
357*c4762a1bSJed Brown     args: -ts_max_steps 10 -monitor -mu 5
358*c4762a1bSJed Brown     output_file: output/ex16adj_tl_4.out
359*c4762a1bSJed Brown 
360*c4762a1bSJed Brown TEST*/
361