xref: /petsc/src/ts/tutorials/hybrid/ex1fd.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
1c4762a1bSJed Brown static char help[] = "An example of hybrid system using TS event.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   The dynamics is described by the ODE
5c4762a1bSJed Brown                   u_t = A_i u
6c4762a1bSJed Brown 
7c4762a1bSJed Brown   where A_1 = [ 1  -100
8c4762a1bSJed Brown                 10  1  ],
9c4762a1bSJed Brown         A_2 = [ 1    10
10c4762a1bSJed Brown                -100  1 ].
11c4762a1bSJed Brown   The index i changes from 1 to 2 when u[1]=2.75u[0] and from 2 to 1 when u[1]=0.36u[0].
12c4762a1bSJed Brown   Initially u=[0 1]^T and i=1.
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   Reference:
15c4762a1bSJed Brown   I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000
16c4762a1bSJed Brown */
17c4762a1bSJed Brown 
18c4762a1bSJed Brown #include <petscts.h>
19c4762a1bSJed Brown 
20c4762a1bSJed Brown typedef struct {
21c4762a1bSJed Brown   PetscReal lambda1;
22c4762a1bSJed Brown   PetscReal lambda2;
23c4762a1bSJed Brown   PetscInt  mode;  /* mode flag*/
24c4762a1bSJed Brown } AppCtx;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown PetscErrorCode FWDRun(TS, Vec, void *);
27c4762a1bSJed Brown 
28c4762a1bSJed Brown PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx)
29c4762a1bSJed Brown {
30c4762a1bSJed Brown   AppCtx            *actx=(AppCtx*)ctx;
31c4762a1bSJed Brown   PetscErrorCode    ierr;
32c4762a1bSJed Brown   const PetscScalar *u;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   PetscFunctionBegin;
35c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
36c4762a1bSJed Brown   if (actx->mode == 1) {
37c4762a1bSJed Brown     fvalue[0] = u[1]-actx->lambda1*u[0];
38c4762a1bSJed Brown   } else if (actx->mode == 2) {
39c4762a1bSJed Brown     fvalue[0] = u[1]-actx->lambda2*u[0];
40c4762a1bSJed Brown   }
41c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
42c4762a1bSJed Brown   PetscFunctionReturn(0);
43c4762a1bSJed Brown }
44c4762a1bSJed Brown 
45c4762a1bSJed Brown PetscErrorCode ShiftGradients(TS ts,Vec U,AppCtx *actx)
46c4762a1bSJed Brown {
47c4762a1bSJed Brown   Vec               *lambda,*mu;
48c4762a1bSJed Brown   PetscScalar       *x,*y;
49c4762a1bSJed Brown   const PetscScalar *u;
50c4762a1bSJed Brown   PetscErrorCode    ierr;
51c4762a1bSJed Brown   PetscScalar       tmp[2],A1[2][2],A2[2],denorm1,denorm2;
52c4762a1bSJed Brown   PetscInt          numcost;
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   PetscFunctionBegin;
55c4762a1bSJed Brown   ierr = TSGetCostGradients(ts,&numcost,&lambda,&mu);CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   if (actx->mode==2) {
59c4762a1bSJed Brown     denorm1 = -actx->lambda1*(u[0]-100.*u[1])+1.*(10.*u[0]+u[1]);
60c4762a1bSJed Brown     denorm2 = -actx->lambda1*(u[0]+10.*u[1])+1.*(-100.*u[0]+u[1]);
61c4762a1bSJed Brown     A1[0][0] = 110.*u[1]*(-actx->lambda1)/denorm1+1.;
62c4762a1bSJed Brown     A1[0][1] = -110.*u[0]*(-actx->lambda1)/denorm1;
63c4762a1bSJed Brown     A1[1][0] = 110.*u[1]*1./denorm1;
64c4762a1bSJed Brown     A1[1][1] = -110.*u[0]*1./denorm1+1.;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown     A2[0] = 110.*u[1]*(-u[0])/denorm2;
67c4762a1bSJed Brown     A2[1] = -110.*u[0]*(-u[0])/denorm2;
68c4762a1bSJed Brown   } else {
69c4762a1bSJed Brown     denorm2 = -actx->lambda2*(u[0]+10.*u[1])+1.*(-100.*u[0]+u[1]);
70c4762a1bSJed Brown     A1[0][0] = 110.*u[1]*(-actx->lambda1)/denorm2+1;
71c4762a1bSJed Brown     A1[0][1] = -110.*u[0]*(-actx->lambda1)/denorm2;
72c4762a1bSJed Brown     A1[1][0] = 110.*u[1]*1./denorm2;
73c4762a1bSJed Brown     A1[1][1] = -110.*u[0]*1./denorm2+1.;
74c4762a1bSJed Brown 
75c4762a1bSJed Brown     A2[0] = 0;
76c4762a1bSJed Brown     A2[1] = 0;
77c4762a1bSJed Brown   }
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   ierr   = VecGetArray(lambda[0],&x);CHKERRQ(ierr);
82c4762a1bSJed Brown   ierr   = VecGetArray(mu[0],&y);CHKERRQ(ierr);
83c4762a1bSJed Brown   tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1];
84c4762a1bSJed Brown   tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1];
85c4762a1bSJed Brown   y[0]   = y[0] + A2[0]*x[0]+A2[1]*x[1];
86c4762a1bSJed Brown   x[0]   = tmp[0];
87c4762a1bSJed Brown   x[1]   = tmp[1];
88c4762a1bSJed Brown   ierr   = VecRestoreArray(mu[0],&y);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr   = VecRestoreArray(lambda[0],&x);CHKERRQ(ierr);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   ierr   = VecGetArray(lambda[1],&x);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr   = VecGetArray(mu[1],&y);CHKERRQ(ierr);
93c4762a1bSJed Brown   tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1];
94c4762a1bSJed Brown   tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1];
95c4762a1bSJed Brown   y[0]   = y[0] + A2[0]*x[0]+A2[1]*x[1];
96c4762a1bSJed Brown   x[0]   = tmp[0];
97c4762a1bSJed Brown   x[1]   = tmp[1];
98c4762a1bSJed Brown   ierr   = VecRestoreArray(mu[1],&y);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr   = VecRestoreArray(lambda[1],&x);CHKERRQ(ierr);
100c4762a1bSJed Brown   PetscFunctionReturn(0);
101c4762a1bSJed Brown }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
104c4762a1bSJed Brown {
105c4762a1bSJed Brown   AppCtx         *actx=(AppCtx*)ctx;
106c4762a1bSJed Brown   PetscErrorCode ierr;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   PetscFunctionBegin;
109c4762a1bSJed Brown   if (!forwardsolve) {
110c4762a1bSJed Brown     ierr = ShiftGradients(ts,U,actx);CHKERRQ(ierr);
111c4762a1bSJed Brown   }
112c4762a1bSJed Brown   if (actx->mode == 1) {
113c4762a1bSJed Brown     actx->mode = 2;
114c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t);CHKERRQ(ierr);
115c4762a1bSJed Brown   } else if (actx->mode == 2) {
116c4762a1bSJed Brown     actx->mode = 1;
117c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t);CHKERRQ(ierr);
118c4762a1bSJed Brown   }
119c4762a1bSJed Brown   PetscFunctionReturn(0);
120c4762a1bSJed Brown }
121c4762a1bSJed Brown 
122c4762a1bSJed Brown /*
123c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
124c4762a1bSJed Brown */
125c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
126c4762a1bSJed Brown {
127c4762a1bSJed Brown   AppCtx            *actx=(AppCtx*)ctx;
128c4762a1bSJed Brown   PetscErrorCode    ierr;
129c4762a1bSJed Brown   PetscScalar       *f;
130c4762a1bSJed Brown   const PetscScalar *u,*udot;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   PetscFunctionBegin;
133c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
134c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
135c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
136c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   if (actx->mode == 1) {
139c4762a1bSJed Brown     f[0] = udot[0]-u[0]+100*u[1];
140c4762a1bSJed Brown     f[1] = udot[1]-10*u[0]-u[1];
141c4762a1bSJed Brown   } else if (actx->mode == 2) {
142c4762a1bSJed Brown     f[0] = udot[0]-u[0]-10*u[1];
143c4762a1bSJed Brown     f[1] = udot[1]+100*u[0]-u[1];
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
147c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
149c4762a1bSJed Brown   PetscFunctionReturn(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown /*
153c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
154c4762a1bSJed Brown */
155c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
156c4762a1bSJed Brown {
157c4762a1bSJed Brown   AppCtx            *actx=(AppCtx*)ctx;
158c4762a1bSJed Brown   PetscErrorCode    ierr;
159c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
160c4762a1bSJed Brown   PetscScalar       J[2][2];
161c4762a1bSJed Brown   const PetscScalar *u,*udot;
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   PetscFunctionBegin;
164c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   if (actx->mode == 1) {
168c4762a1bSJed Brown     J[0][0] = a-1;                       J[0][1] = 100;
169c4762a1bSJed Brown     J[1][0] = -10;                       J[1][1] = a-1;
170c4762a1bSJed Brown   } else if (actx->mode == 2) {
171c4762a1bSJed Brown     J[0][0] = a-1;                       J[0][1] = -10;
172c4762a1bSJed Brown     J[1][0] = 100;                       J[1][1] = a-1;
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown   ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
177c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
180c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
181c4762a1bSJed Brown   if (A != B) {
182c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
183c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
184c4762a1bSJed Brown   }
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown int main(int argc,char **argv)
189c4762a1bSJed Brown {
190c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
191c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
192c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
193c4762a1bSJed Brown   Mat            Ap;            /* dfdp */
194c4762a1bSJed Brown   PetscErrorCode ierr;
195c4762a1bSJed Brown   PetscMPIInt    size;
196c4762a1bSJed Brown   PetscInt       n = 2;
197c4762a1bSJed Brown   PetscScalar    *u;
198c4762a1bSJed Brown   AppCtx         app;
199c4762a1bSJed Brown   PetscInt       direction[1];
200c4762a1bSJed Brown   PetscBool      terminate[1];
201c4762a1bSJed Brown   PetscReal      delta,tmp[2],sensi[2];
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   delta = 1e-8;
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206c4762a1bSJed Brown      Initialize program
207c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
208c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
209ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
210c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
211c4762a1bSJed Brown   app.mode = 1;
212c4762a1bSJed Brown   app.lambda1 = 2.75;
213c4762a1bSJed Brown   app.lambda2 = 0.36;
214c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex1 options","");CHKERRQ(ierr);
215c4762a1bSJed Brown   {
216c4762a1bSJed Brown     ierr = PetscOptionsReal("-lambda1","","",app.lambda1,&app.lambda1,NULL);CHKERRQ(ierr);
217c4762a1bSJed Brown     ierr = PetscOptionsReal("-lambda2","","",app.lambda2,&app.lambda2,NULL);CHKERRQ(ierr);
218c4762a1bSJed Brown   }
219c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222c4762a1bSJed Brown     Create necessary matrix and vectors
223c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
225c4762a1bSJed Brown   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
226c4762a1bSJed Brown   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
227c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
228c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&Ap);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = MatSetSizes(Ap,n,1,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = MatSetType(Ap,MATDENSE);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = MatSetFromOptions(Ap);CHKERRQ(ierr);
236c4762a1bSJed Brown   ierr = MatSetUp(Ap);CHKERRQ(ierr);
237c4762a1bSJed Brown   ierr = MatZeroEntries(Ap);CHKERRQ(ierr); /* initialize to zeros */
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
240c4762a1bSJed Brown   u[0] = 0;
241c4762a1bSJed Brown   u[1] = 1;
242c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
243c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown      Create timestepping solver context
245c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
247c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
248c4762a1bSJed Brown   ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
249c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,(TSIFunction)IFunction,&app);CHKERRQ(ierr);
250c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&app);CHKERRQ(ierr);
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
253c4762a1bSJed Brown      Set initial conditions
254c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
255c4762a1bSJed Brown   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258c4762a1bSJed Brown      Set solver options
259c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
260c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,0.125);CHKERRQ(ierr);
261c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
262c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,1./256.);CHKERRQ(ierr);
263c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /* Set directions and terminate flags for the two events */
266c4762a1bSJed Brown   direction[0] = 0;
267c4762a1bSJed Brown   terminate[0] = PETSC_FALSE;
268c4762a1bSJed Brown   ierr = TSSetEventHandler(ts,1,direction,terminate,EventFunction,PostEventFunction,(void*)&app);CHKERRQ(ierr);
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
271c4762a1bSJed Brown      Run timestepping solver
272c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
273c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
276c4762a1bSJed Brown   tmp[0] = u[0];
277c4762a1bSJed Brown   tmp[1] = u[1];
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   u[0] = 0+delta;
280c4762a1bSJed Brown   u[1] = 1;
281c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   ierr = FWDRun(ts,U,(void*)&app);CHKERRQ(ierr);
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
286c4762a1bSJed Brown   sensi[0] = (u[0]-tmp[0])/delta;
287c4762a1bSJed Brown   sensi[1] = (u[1]-tmp[1])/delta;
288c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"d x1(tf) /d x1(t0) = %f d x2(tf) / d x1(t0) = %f \n",sensi[0],sensi[1]);CHKERRQ(ierr);
289c4762a1bSJed Brown   u[0] = 0;
290c4762a1bSJed Brown   u[1] = 1+delta;
291c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   ierr = FWDRun(ts,U,(void*)&app);CHKERRQ(ierr);
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
296c4762a1bSJed Brown   sensi[0] = (u[0]-tmp[0])/delta;
297c4762a1bSJed Brown   sensi[1] = (u[1]-tmp[1])/delta;
298c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"d x1(tf) /d x2(t0) = %f d x2(tf) / d x2(t0) = %f \n",sensi[0],sensi[1]);CHKERRQ(ierr);
299c4762a1bSJed Brown   u[0] = 0;
300c4762a1bSJed Brown   u[1] = 1;
301c4762a1bSJed Brown   app.lambda1 = app.lambda1+delta;
302c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   ierr = FWDRun(ts,U,(void*)&app);CHKERRQ(ierr);
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
307c4762a1bSJed Brown   sensi[0] = (u[0]-tmp[0])/delta;
308c4762a1bSJed Brown   sensi[1] = (u[1]-tmp[1])/delta;
309c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"Final gradients: d x1(tf) /d p = %f d x2(tf) / d p = %f \n",sensi[0],sensi[1]);CHKERRQ(ierr);
310c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
313c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
314c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
315c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
316c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
317c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   ierr = MatDestroy(&Ap);CHKERRQ(ierr);
320c4762a1bSJed Brown   ierr = PetscFinalize();
321c4762a1bSJed Brown   return ierr;
322c4762a1bSJed Brown }
323c4762a1bSJed Brown 
324c4762a1bSJed Brown PetscErrorCode FWDRun(TS ts, Vec U0, void *ctx0)
325c4762a1bSJed Brown {
326c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
327c4762a1bSJed Brown   PetscErrorCode ierr;
328c4762a1bSJed Brown   AppCtx         *ctx=(AppCtx*)ctx0;
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   PetscFunctionBeginUser;
331c4762a1bSJed Brown   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
332c4762a1bSJed Brown   ierr = VecCopy(U0,U);CHKERRQ(ierr);
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   ctx->mode = 1;
335c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
336c4762a1bSJed Brown      Run timestepping solver
337c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
338c4762a1bSJed Brown   ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr);
339c4762a1bSJed Brown 
340c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
341c4762a1bSJed Brown 
342c4762a1bSJed Brown   PetscFunctionReturn(0);
343c4762a1bSJed Brown }
344c4762a1bSJed Brown 
345c4762a1bSJed Brown /*TEST
346c4762a1bSJed Brown 
347c4762a1bSJed Brown    build:
348*dfd57a17SPierre Jolivet       requires: !defined(PETSC_USE_CXXCOMPLEX)
349c4762a1bSJed Brown 
350c4762a1bSJed Brown    test:
351c4762a1bSJed Brown       args: -ts_event_tol 1e-9
352c4762a1bSJed Brown       timeoutfactor: 18
353c4762a1bSJed Brown       requires: !single
354c4762a1bSJed Brown 
355c4762a1bSJed Brown TEST*/
356