xref: /petsc/src/ts/tutorials/ex20adj.c (revision 5ab1ac2bae2c53bf5666daef2efaf6cb34bd2dd4)
1c4762a1bSJed Brown static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
5c4762a1bSJed Brown    Concepts: TS^van der Pol equation DAE equivalent
6c4762a1bSJed Brown    Concepts: TS^adjoint sensitivity analysis
7c4762a1bSJed Brown    Processors: 1
8c4762a1bSJed Brown */
9c4762a1bSJed Brown /* ------------------------------------------------------------------------
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    This program solves the van der Pol DAE ODE equivalent
12c4762a1bSJed Brown       [ u_1' ] = [          u_2                ]  (2)
13c4762a1bSJed Brown       [ u_2' ]   [ \mu ((1 - u_1^2) u_2 - u_1) ]
14c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
15c4762a1bSJed Brown        u_1(0) = 2, u_2(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2),
16c4762a1bSJed Brown    and
17c4762a1bSJed Brown        \mu = 10^6 ( y'(0) ~ -0.6666665432100101).,
18c4762a1bSJed Brown    and computes the sensitivities of the final solution w.r.t. initial conditions and parameter \mu with the implicit theta method and its discrete adjoint.
19c4762a1bSJed Brown 
20c4762a1bSJed Brown    Notes:
21c4762a1bSJed Brown    This code demonstrates the TSAdjoint interface to a DAE system.
22c4762a1bSJed Brown 
23c4762a1bSJed Brown    The user provides the implicit right-hand-side function
24*5ab1ac2bSHong Zhang    [ F(u',u,t) ] = [u' - f(u,t)] = [ u_1'] - [        u_2             ]
25c4762a1bSJed Brown                                    [ u_2']   [ \mu ((1-u_1^2)u_2-u_1) ]
26c4762a1bSJed Brown 
27*5ab1ac2bSHong Zhang    and the Jacobian of F (from the PETSc user manual)
28c4762a1bSJed Brown 
29*5ab1ac2bSHong Zhang               dF   dF
30*5ab1ac2bSHong Zhang    J(F) = a * -- + --
31c4762a1bSJed Brown               du'  du
32c4762a1bSJed Brown 
33*5ab1ac2bSHong Zhang    and the JacobianP of the explicit right-hand side of (2) f(u,t) ( which is equivalent to -F(0,u,t) ).
34c4762a1bSJed Brown    df   [       0               ]
35c4762a1bSJed Brown    -- = [                       ]
36c4762a1bSJed Brown    dp   [ (1 - u_1^2) u_2 - u_1 ].
37c4762a1bSJed Brown 
38c4762a1bSJed Brown    See ex20.c for more details on the Jacobian.
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   ------------------------------------------------------------------------- */
41c4762a1bSJed Brown #include <petscts.h>
42c4762a1bSJed Brown #include <petsctao.h>
43c4762a1bSJed Brown 
44c4762a1bSJed Brown typedef struct _n_User *User;
45c4762a1bSJed Brown struct _n_User {
46c4762a1bSJed Brown   PetscReal mu;
47c4762a1bSJed Brown   PetscReal next_output;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /* Sensitivity analysis support */
50c4762a1bSJed Brown   PetscInt  steps;
51c4762a1bSJed Brown   PetscReal ftime;
52c4762a1bSJed Brown   Mat       A;                   /* Jacobian matrix */
53c4762a1bSJed Brown   Mat       Jacp;                /* JacobianP matrix */
54c4762a1bSJed Brown   Vec       U,lambda[2],mup[2];  /* adjoint variables */
55c4762a1bSJed Brown };
56c4762a1bSJed Brown 
57c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE  -------------------- */
58c4762a1bSJed Brown 
59c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   PetscErrorCode    ierr;
62c4762a1bSJed Brown   User              user = (User)ctx;
63c4762a1bSJed Brown   PetscScalar       *f;
64c4762a1bSJed Brown   const PetscScalar *u;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   PetscFunctionBeginUser;
67c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
69c4762a1bSJed Brown   f[0] = u[1];
70c4762a1bSJed Brown   f[1] = user->mu*((1.-u[0]*u[0])*u[1]-u[0]);
71c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
73c4762a1bSJed Brown   PetscFunctionReturn(0);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
77c4762a1bSJed Brown {
78c4762a1bSJed Brown   PetscErrorCode    ierr;
79c4762a1bSJed Brown   User              user = (User)ctx;
80c4762a1bSJed Brown   PetscReal         mu   = user->mu;
81c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
82c4762a1bSJed Brown   PetscScalar       J[2][2];
83c4762a1bSJed Brown   const PetscScalar *u;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   PetscFunctionBeginUser;
86c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
87c4762a1bSJed Brown   J[0][0] = 0;
88c4762a1bSJed Brown   J[1][0] = -mu*(2.0*u[1]*u[0]+1.);
89c4762a1bSJed Brown   J[0][1] = 1.0;
90c4762a1bSJed Brown   J[1][1] = mu*(1.0-u[0]*u[0]);
91c4762a1bSJed Brown   ierr    = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94c4762a1bSJed Brown   if (A != B) {
95c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
96c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
97c4762a1bSJed Brown   }
98c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
99c4762a1bSJed Brown   PetscFunctionReturn(0);
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE  -------------------- */
103c4762a1bSJed Brown 
104c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
105c4762a1bSJed Brown {
106c4762a1bSJed Brown   PetscErrorCode    ierr;
107c4762a1bSJed Brown   User              user = (User)ctx;
108c4762a1bSJed Brown   const PetscScalar *u,*udot;
109c4762a1bSJed Brown   PetscScalar       *f;
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   PetscFunctionBeginUser;
112c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
115c4762a1bSJed Brown   f[0] = udot[0] - u[1];
116c4762a1bSJed Brown   f[1] = udot[1] - user->mu*((1.0-u[0]*u[0])*u[1] - u[0]);
117c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
120c4762a1bSJed Brown   PetscFunctionReturn(0);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
124c4762a1bSJed Brown {
125c4762a1bSJed Brown   PetscErrorCode    ierr;
126c4762a1bSJed Brown   User              user     = (User)ctx;
127c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
128c4762a1bSJed Brown   PetscScalar       J[2][2];
129c4762a1bSJed Brown   const PetscScalar *u;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   PetscFunctionBeginUser;
132c4762a1bSJed Brown   ierr    = VecGetArrayRead(U,&u);CHKERRQ(ierr);
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   J[0][0] = a;     J[0][1] =  -1.0;
135c4762a1bSJed Brown   J[1][0] = user->mu*(2.0*u[0]*u[1] + 1.0);   J[1][1] = a - user->mu*(1.0-u[0]*u[0]);
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
138c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
142c4762a1bSJed Brown   if (B && A != B) {
143c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
144c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown   PetscFunctionReturn(0);
147c4762a1bSJed Brown }
148c4762a1bSJed Brown 
149c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
150c4762a1bSJed Brown {
151c4762a1bSJed Brown   PetscErrorCode    ierr;
152c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[]={0};
153c4762a1bSJed Brown   PetscScalar       J[2][1];
154c4762a1bSJed Brown   const PetscScalar *u;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   PetscFunctionBeginUser;
157c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
158c4762a1bSJed Brown   J[0][0] = 0;
159c4762a1bSJed Brown   J[1][0] = (1.-u[0]*u[0])*u[1]-u[0];
160c4762a1bSJed Brown   ierr    = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
163c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
164c4762a1bSJed Brown   PetscFunctionReturn(0);
165c4762a1bSJed Brown }
166c4762a1bSJed Brown 
167c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
168c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx)
169c4762a1bSJed Brown {
170c4762a1bSJed Brown   PetscErrorCode    ierr;
171c4762a1bSJed Brown   const PetscScalar *u;
172c4762a1bSJed Brown   PetscReal         tfinal, dt;
173c4762a1bSJed Brown   User              user = (User)ctx;
174c4762a1bSJed Brown   Vec               interpolatedU;
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   PetscFunctionBeginUser;
177c4762a1bSJed Brown   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
178c4762a1bSJed Brown   ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr);
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
181c4762a1bSJed Brown     ierr = VecDuplicate(U,&interpolatedU);CHKERRQ(ierr);
182c4762a1bSJed Brown     ierr = TSInterpolate(ts,user->next_output,interpolatedU);CHKERRQ(ierr);
183c4762a1bSJed Brown     ierr = VecGetArrayRead(interpolatedU,&u);CHKERRQ(ierr);
184c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"[%g] %D TS %g (dt = %g) X %g %g\n",
185c4762a1bSJed Brown                        (double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(u[0]),
186c4762a1bSJed Brown                        (double)PetscRealPart(u[1]));CHKERRQ(ierr);
187c4762a1bSJed Brown     ierr = VecRestoreArrayRead(interpolatedU,&u);CHKERRQ(ierr);
188c4762a1bSJed Brown     ierr = VecDestroy(&interpolatedU);CHKERRQ(ierr);
189c4762a1bSJed Brown     user->next_output += 0.1;
190c4762a1bSJed Brown   }
191c4762a1bSJed Brown   PetscFunctionReturn(0);
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown int main(int argc,char **argv)
195c4762a1bSJed Brown {
196c4762a1bSJed Brown   TS             ts;
197c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE,implicitform = PETSC_TRUE;
198c4762a1bSJed Brown   PetscScalar    *x_ptr,*y_ptr,derp;
199c4762a1bSJed Brown   PetscMPIInt    size;
200c4762a1bSJed Brown   struct _n_User user;
201c4762a1bSJed Brown   PetscErrorCode ierr;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204c4762a1bSJed Brown      Initialize program
205c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
207c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
208c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211c4762a1bSJed Brown     Set runtime options
212c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213c4762a1bSJed Brown   user.next_output = 0.0;
214c4762a1bSJed Brown   user.mu          = 1.0e3;
215c4762a1bSJed Brown   user.steps       = 0;
216c4762a1bSJed Brown   user.ftime       = 0.5;
217c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr);
218c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr);
219c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);CHKERRQ(ierr);
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
223c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user.A);CHKERRQ(ierr);
225c4762a1bSJed Brown   ierr = MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
226c4762a1bSJed Brown   ierr = MatSetFromOptions(user.A);CHKERRQ(ierr);
227c4762a1bSJed Brown   ierr = MatSetUp(user.A);CHKERRQ(ierr);
228c4762a1bSJed Brown   ierr = MatCreateVecs(user.A,&user.U,NULL);CHKERRQ(ierr);
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr);
231c4762a1bSJed Brown   ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = MatSetUp(user.Jacp);CHKERRQ(ierr);
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236c4762a1bSJed Brown      Create timestepping solver context
237c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
239c4762a1bSJed Brown   ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
240c4762a1bSJed Brown   if (implicitform) {
241c4762a1bSJed Brown     ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr);
242c4762a1bSJed Brown     ierr = TSSetIJacobian(ts,user.A,user.A,IJacobian,&user);CHKERRQ(ierr);
243c4762a1bSJed Brown     ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
244c4762a1bSJed Brown   } else {
245c4762a1bSJed Brown     ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr);
246c4762a1bSJed Brown     ierr = TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user);CHKERRQ(ierr);
247c4762a1bSJed Brown     ierr = TSSetType(ts,TSRK);CHKERRQ(ierr);
248c4762a1bSJed Brown   }
249c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user);CHKERRQ(ierr);
250c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,user.ftime);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr);
252c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
253c4762a1bSJed Brown   if (monitor) {
254c4762a1bSJed Brown     ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr);
255c4762a1bSJed Brown   }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258c4762a1bSJed Brown      Set initial conditions
259c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
260c4762a1bSJed Brown   ierr = VecGetArray(user.U,&x_ptr);CHKERRQ(ierr);
261c4762a1bSJed Brown   x_ptr[0] = 2.0;
262c4762a1bSJed Brown   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
263c4762a1bSJed Brown   ierr = VecRestoreArray(user.U,&x_ptr);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr);
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
268c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
269c4762a1bSJed Brown   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272c4762a1bSJed Brown      Set runtime options
273c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
274c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
275c4762a1bSJed Brown 
276c4762a1bSJed Brown   ierr = TSSolve(ts,user.U);CHKERRQ(ierr);
277c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&user.ftime);CHKERRQ(ierr);
278c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&user.steps);CHKERRQ(ierr);
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
281c4762a1bSJed Brown      Adjoint model starts here
282c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
283c4762a1bSJed Brown   ierr = MatCreateVecs(user.A,&user.lambda[0],NULL);CHKERRQ(ierr);
284c4762a1bSJed Brown   /* Set initial conditions for the adjoint integration */
285c4762a1bSJed Brown   ierr = VecGetArray(user.lambda[0],&y_ptr);CHKERRQ(ierr);
286c4762a1bSJed Brown   y_ptr[0] = 1.0; y_ptr[1] = 0.0;
287c4762a1bSJed Brown   ierr = VecRestoreArray(user.lambda[0],&y_ptr);CHKERRQ(ierr);
288c4762a1bSJed Brown   ierr = MatCreateVecs(user.A,&user.lambda[1],NULL);CHKERRQ(ierr);
289c4762a1bSJed Brown   ierr = VecGetArray(user.lambda[1],&y_ptr);CHKERRQ(ierr);
290c4762a1bSJed Brown   y_ptr[0] = 0.0; y_ptr[1] = 1.0;
291c4762a1bSJed Brown   ierr = VecRestoreArray(user.lambda[1],&y_ptr);CHKERRQ(ierr);
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   ierr = MatCreateVecs(user.Jacp,&user.mup[0],NULL);CHKERRQ(ierr);
294c4762a1bSJed Brown   ierr = VecGetArray(user.mup[0],&x_ptr);CHKERRQ(ierr);
295c4762a1bSJed Brown   x_ptr[0] = 0.0;
296c4762a1bSJed Brown   ierr = VecRestoreArray(user.mup[0],&x_ptr);CHKERRQ(ierr);
297c4762a1bSJed Brown   ierr = MatCreateVecs(user.Jacp,&user.mup[1],NULL);CHKERRQ(ierr);
298c4762a1bSJed Brown   ierr = VecGetArray(user.mup[1],&x_ptr);CHKERRQ(ierr);
299c4762a1bSJed Brown   x_ptr[0] = 0.0;
300c4762a1bSJed Brown   ierr = VecRestoreArray(user.mup[1],&x_ptr);CHKERRQ(ierr);
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   ierr = TSSetCostGradients(ts,2,user.lambda,user.mup);CHKERRQ(ierr);
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[y(tf)]/d[y0]  d[y(tf)]/d[z0]\n");CHKERRQ(ierr);
307c4762a1bSJed Brown   ierr = VecView(user.lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
308c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[z(tf)]/d[y0]  d[z(tf)]/d[z0]\n");CHKERRQ(ierr);
309c4762a1bSJed Brown   ierr = VecView(user.lambda[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   ierr = VecGetArray(user.mup[0],&x_ptr);CHKERRQ(ierr);
312c4762a1bSJed Brown   ierr = VecGetArray(user.lambda[0],&y_ptr);CHKERRQ(ierr);
313c4762a1bSJed Brown   derp = y_ptr[1]*(-10.0/(81.0*user.mu*user.mu)+2.0*292.0/(2187.0*user.mu*user.mu*user.mu))+x_ptr[0];
314c4762a1bSJed Brown   ierr = VecRestoreArray(user.mup[0],&x_ptr);CHKERRQ(ierr);
315c4762a1bSJed Brown   ierr = VecRestoreArray(user.lambda[0],&y_ptr);CHKERRQ(ierr);
316c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameters: d[y(tf)]/d[mu]\n%g\n",(double)PetscRealPart(derp));CHKERRQ(ierr);
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   ierr = VecGetArray(user.mup[1],&x_ptr);CHKERRQ(ierr);
319c4762a1bSJed Brown   ierr = VecGetArray(user.lambda[1],&y_ptr);CHKERRQ(ierr);
320c4762a1bSJed Brown   derp = y_ptr[1]*(-10.0/(81.0*user.mu*user.mu)+2.0*292.0/(2187.0*user.mu*user.mu*user.mu))+x_ptr[0];
321c4762a1bSJed Brown   ierr = VecRestoreArray(user.mup[1],&x_ptr);CHKERRQ(ierr);
322c4762a1bSJed Brown   ierr = VecRestoreArray(user.lambda[1],&y_ptr);CHKERRQ(ierr);
323c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensivitity wrt parameters: d[z(tf)]/d[mu]\n%g\n",(double)PetscRealPart(derp));CHKERRQ(ierr);
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
326c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
327c4762a1bSJed Brown      are no longer needed.
328c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
329c4762a1bSJed Brown   ierr = MatDestroy(&user.A);CHKERRQ(ierr);
330c4762a1bSJed Brown   ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr);
331c4762a1bSJed Brown   ierr = VecDestroy(&user.U);CHKERRQ(ierr);
332c4762a1bSJed Brown   ierr = VecDestroy(&user.lambda[0]);CHKERRQ(ierr);
333c4762a1bSJed Brown   ierr = VecDestroy(&user.lambda[1]);CHKERRQ(ierr);
334c4762a1bSJed Brown   ierr = VecDestroy(&user.mup[0]);CHKERRQ(ierr);
335c4762a1bSJed Brown   ierr = VecDestroy(&user.mup[1]);CHKERRQ(ierr);
336c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   ierr = PetscFinalize();
339c4762a1bSJed Brown   return(ierr);
340c4762a1bSJed Brown }
341c4762a1bSJed Brown 
342c4762a1bSJed Brown /*TEST
343c4762a1bSJed Brown 
344c4762a1bSJed Brown     test:
345c4762a1bSJed Brown       requires: revolve
346c4762a1bSJed Brown       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -viewer_binary_skip_info -ts_dt 0.001 -mu 100000
347c4762a1bSJed Brown 
348c4762a1bSJed Brown     test:
349c4762a1bSJed Brown       suffix: 2
350c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only
351c4762a1bSJed Brown 
352c4762a1bSJed Brown     test:
353c4762a1bSJed Brown       suffix: 3
354c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 0
355c4762a1bSJed Brown       output_file: output/ex20adj_2.out
356c4762a1bSJed Brown 
357c4762a1bSJed Brown     test:
358c4762a1bSJed Brown       suffix: 4
359c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack
360c4762a1bSJed Brown       output_file: output/ex20adj_2.out
361c4762a1bSJed Brown 
362c4762a1bSJed Brown     test:
363c4762a1bSJed Brown       suffix: 5
364c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack
365c4762a1bSJed Brown       output_file: output/ex20adj_2.out
366c4762a1bSJed Brown 
367c4762a1bSJed Brown     test:
368c4762a1bSJed Brown       suffix: 6
369c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 0
370c4762a1bSJed Brown       output_file: output/ex20adj_2.out
371c4762a1bSJed Brown 
372c4762a1bSJed Brown     test:
373c4762a1bSJed Brown       suffix: 7
374c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0
375c4762a1bSJed Brown       output_file: output/ex20adj_2.out
376c4762a1bSJed Brown 
377c4762a1bSJed Brown     test:
378c4762a1bSJed Brown       suffix: 8
379c4762a1bSJed Brown       requires: revolve
380c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_solution_only -ts_trajectory_monitor
381c4762a1bSJed Brown       output_file: output/ex20adj_3.out
382c4762a1bSJed Brown 
383c4762a1bSJed Brown     test:
384c4762a1bSJed Brown       suffix: 9
385c4762a1bSJed Brown       requires: revolve
386c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_solution_only 0 -ts_trajectory_monitor
387c4762a1bSJed Brown       output_file: output/ex20adj_4.out
388c4762a1bSJed Brown 
389c4762a1bSJed Brown     test:
390c4762a1bSJed Brown       requires: revolve
391c4762a1bSJed Brown       suffix: 10
392c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_revolve_online -ts_trajectory_solution_only
393c4762a1bSJed Brown       output_file: output/ex20adj_2.out
394c4762a1bSJed Brown 
395c4762a1bSJed Brown     test:
396c4762a1bSJed Brown       requires: revolve
397c4762a1bSJed Brown       suffix: 11
398c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_revolve_online -ts_trajectory_solution_only 0
399c4762a1bSJed Brown       output_file: output/ex20adj_2.out
400c4762a1bSJed Brown 
401c4762a1bSJed Brown     test:
402c4762a1bSJed Brown       suffix: 12
403c4762a1bSJed Brown       requires: revolve
404c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only
405c4762a1bSJed Brown       output_file: output/ex20adj_2.out
406c4762a1bSJed Brown 
407c4762a1bSJed Brown     test:
408c4762a1bSJed Brown       suffix: 13
409c4762a1bSJed Brown       requires: revolve
410c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only 0
411c4762a1bSJed Brown       output_file: output/ex20adj_2.out
412c4762a1bSJed Brown 
413c4762a1bSJed Brown     test:
414c4762a1bSJed Brown       suffix: 14
415c4762a1bSJed Brown       requires: revolve
416c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack
417c4762a1bSJed Brown       output_file: output/ex20adj_2.out
418c4762a1bSJed Brown 
419c4762a1bSJed Brown     test:
420c4762a1bSJed Brown       suffix: 15
421c4762a1bSJed Brown       requires: revolve
422c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 0
423c4762a1bSJed Brown       output_file: output/ex20adj_2.out
424c4762a1bSJed Brown 
425c4762a1bSJed Brown     test:
426c4762a1bSJed Brown       suffix: 16
427c4762a1bSJed Brown       requires: revolve
428c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack
429c4762a1bSJed Brown       output_file: output/ex20adj_2.out
430c4762a1bSJed Brown 
431c4762a1bSJed Brown     test:
432c4762a1bSJed Brown       suffix: 17
433c4762a1bSJed Brown       requires: revolve
434c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0
435c4762a1bSJed Brown       output_file: output/ex20adj_2.out
436c4762a1bSJed Brown 
437c4762a1bSJed Brown     test:
438c4762a1bSJed Brown       suffix: 18
439c4762a1bSJed Brown       requires: revolve
440c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack
441c4762a1bSJed Brown       output_file: output/ex20adj_2.out
442c4762a1bSJed Brown 
443c4762a1bSJed Brown     test:
444c4762a1bSJed Brown       suffix: 19
445c4762a1bSJed Brown       requires: revolve
446c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack
447c4762a1bSJed Brown       output_file: output/ex20adj_2.out
448c4762a1bSJed Brown 
449c4762a1bSJed Brown     test:
450c4762a1bSJed Brown       suffix: 20
451c4762a1bSJed Brown       requires: revolve
452c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only 0
453c4762a1bSJed Brown       output_file: output/ex20adj_2.out
454c4762a1bSJed Brown 
455c4762a1bSJed Brown     test:
456c4762a1bSJed Brown       suffix: 21
457c4762a1bSJed Brown       requires: revolve
458c4762a1bSJed Brown       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0
459c4762a1bSJed Brown       output_file: output/ex20adj_2.out
460c4762a1bSJed Brown 
461c4762a1bSJed Brown     test:
462c4762a1bSJed Brown       suffix: 22
463c4762a1bSJed Brown       args: -ts_type beuler -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only
464c4762a1bSJed Brown       output_file: output/ex20adj_2.out
465c4762a1bSJed Brown TEST*/
466