xref: /petsc/src/ts/tutorials/hybrid/ex1fd.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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   const PetscScalar *u;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   PetscFunctionBegin;
349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
35c4762a1bSJed Brown   if (actx->mode == 1) {
36c4762a1bSJed Brown     fvalue[0] = u[1]-actx->lambda1*u[0];
37c4762a1bSJed Brown   } else if (actx->mode == 2) {
38c4762a1bSJed Brown     fvalue[0] = u[1]-actx->lambda2*u[0];
39c4762a1bSJed Brown   }
409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
41c4762a1bSJed Brown   PetscFunctionReturn(0);
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown PetscErrorCode ShiftGradients(TS ts,Vec U,AppCtx *actx)
45c4762a1bSJed Brown {
46c4762a1bSJed Brown   Vec               *lambda,*mu;
47c4762a1bSJed Brown   PetscScalar       *x,*y;
48c4762a1bSJed Brown   const PetscScalar *u;
49c4762a1bSJed Brown   PetscScalar       tmp[2],A1[2][2],A2[2],denorm1,denorm2;
50c4762a1bSJed Brown   PetscInt          numcost;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(TSGetCostGradients(ts,&numcost,&lambda,&mu));
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   if (actx->mode==2) {
57c4762a1bSJed Brown     denorm1 = -actx->lambda1*(u[0]-100.*u[1])+1.*(10.*u[0]+u[1]);
58c4762a1bSJed Brown     denorm2 = -actx->lambda1*(u[0]+10.*u[1])+1.*(-100.*u[0]+u[1]);
59c4762a1bSJed Brown     A1[0][0] = 110.*u[1]*(-actx->lambda1)/denorm1+1.;
60c4762a1bSJed Brown     A1[0][1] = -110.*u[0]*(-actx->lambda1)/denorm1;
61c4762a1bSJed Brown     A1[1][0] = 110.*u[1]*1./denorm1;
62c4762a1bSJed Brown     A1[1][1] = -110.*u[0]*1./denorm1+1.;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown     A2[0] = 110.*u[1]*(-u[0])/denorm2;
65c4762a1bSJed Brown     A2[1] = -110.*u[0]*(-u[0])/denorm2;
66c4762a1bSJed Brown   } else {
67c4762a1bSJed Brown     denorm2 = -actx->lambda2*(u[0]+10.*u[1])+1.*(-100.*u[0]+u[1]);
68c4762a1bSJed Brown     A1[0][0] = 110.*u[1]*(-actx->lambda1)/denorm2+1;
69c4762a1bSJed Brown     A1[0][1] = -110.*u[0]*(-actx->lambda1)/denorm2;
70c4762a1bSJed Brown     A1[1][0] = 110.*u[1]*1./denorm2;
71c4762a1bSJed Brown     A1[1][1] = -110.*u[0]*1./denorm2+1.;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown     A2[0] = 0;
74c4762a1bSJed Brown     A2[1] = 0;
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0],&x));
809566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0],&y));
81c4762a1bSJed Brown   tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1];
82c4762a1bSJed Brown   tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1];
83c4762a1bSJed Brown   y[0]   = y[0] + A2[0]*x[0]+A2[1]*x[1];
84c4762a1bSJed Brown   x[0]   = tmp[0];
85c4762a1bSJed Brown   x[1]   = tmp[1];
869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0],&y));
879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0],&x));
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[1],&x));
909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[1],&y));
91c4762a1bSJed Brown   tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1];
92c4762a1bSJed Brown   tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1];
93c4762a1bSJed Brown   y[0]   = y[0] + A2[0]*x[0]+A2[1]*x[1];
94c4762a1bSJed Brown   x[0]   = tmp[0];
95c4762a1bSJed Brown   x[1]   = tmp[1];
969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[1],&y));
979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[1],&x));
98c4762a1bSJed Brown   PetscFunctionReturn(0);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
102c4762a1bSJed Brown {
103c4762a1bSJed Brown   AppCtx         *actx=(AppCtx*)ctx;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   PetscFunctionBegin;
106c4762a1bSJed Brown   if (!forwardsolve) {
1079566063dSJacob Faibussowitsch     PetscCall(ShiftGradients(ts,U,actx));
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown   if (actx->mode == 1) {
110c4762a1bSJed Brown     actx->mode = 2;
1119566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t));
112c4762a1bSJed Brown   } else if (actx->mode == 2) {
113c4762a1bSJed Brown     actx->mode = 1;
1149566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t));
115c4762a1bSJed Brown   }
116c4762a1bSJed Brown   PetscFunctionReturn(0);
117c4762a1bSJed Brown }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown /*
120c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
121c4762a1bSJed Brown */
122c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
123c4762a1bSJed Brown {
124c4762a1bSJed Brown   AppCtx            *actx=(AppCtx*)ctx;
125c4762a1bSJed Brown   PetscScalar       *f;
126c4762a1bSJed Brown   const PetscScalar *u,*udot;
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   PetscFunctionBegin;
129c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1319566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
1329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   if (actx->mode == 1) {
135c4762a1bSJed Brown     f[0] = udot[0]-u[0]+100*u[1];
136c4762a1bSJed Brown     f[1] = udot[1]-10*u[0]-u[1];
137c4762a1bSJed Brown   } else if (actx->mode == 2) {
138c4762a1bSJed Brown     f[0] = udot[0]-u[0]-10*u[1];
139c4762a1bSJed Brown     f[1] = udot[1]+100*u[0]-u[1];
140c4762a1bSJed Brown   }
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
1449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
145c4762a1bSJed Brown   PetscFunctionReturn(0);
146c4762a1bSJed Brown }
147c4762a1bSJed Brown 
148c4762a1bSJed Brown /*
149c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
150c4762a1bSJed Brown */
151c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
152c4762a1bSJed Brown {
153c4762a1bSJed Brown   AppCtx            *actx=(AppCtx*)ctx;
154c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
155c4762a1bSJed Brown   PetscScalar       J[2][2];
156c4762a1bSJed Brown   const PetscScalar *u,*udot;
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   PetscFunctionBegin;
1599566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   if (actx->mode == 1) {
163c4762a1bSJed Brown     J[0][0] = a-1;                       J[0][1] = 100;
164c4762a1bSJed Brown     J[1][0] = -10;                       J[1][1] = a-1;
165c4762a1bSJed Brown   } else if (actx->mode == 2) {
166c4762a1bSJed Brown     J[0][0] = a-1;                       J[0][1] = -10;
167c4762a1bSJed Brown     J[1][0] = 100;                       J[1][1] = a-1;
168c4762a1bSJed Brown   }
1699566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
170c4762a1bSJed Brown 
1719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
176c4762a1bSJed Brown   if (A != B) {
1779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1789566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
179c4762a1bSJed Brown   }
180c4762a1bSJed Brown   PetscFunctionReturn(0);
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown int main(int argc,char **argv)
184c4762a1bSJed Brown {
185c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
186c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
187c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
188c4762a1bSJed Brown   Mat            Ap;            /* dfdp */
189c4762a1bSJed Brown   PetscMPIInt    size;
190c4762a1bSJed Brown   PetscInt       n = 2;
191c4762a1bSJed Brown   PetscScalar    *u;
192c4762a1bSJed Brown   AppCtx         app;
193c4762a1bSJed Brown   PetscInt       direction[1];
194c4762a1bSJed Brown   PetscBool      terminate[1];
195c4762a1bSJed Brown   PetscReal      delta,tmp[2],sensi[2];
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   delta = 1e-8;
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200c4762a1bSJed Brown      Initialize program
201c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2029566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
2039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
2043c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
205c4762a1bSJed Brown   app.mode = 1;
206c4762a1bSJed Brown   app.lambda1 = 2.75;
207c4762a1bSJed Brown   app.lambda2 = 0.36;
208d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex1 options","");
209c4762a1bSJed Brown   {
2109566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda1","","",app.lambda1,&app.lambda1,NULL));
2119566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda2","","",app.lambda2,&app.lambda2,NULL));
212c4762a1bSJed Brown   }
213d0609cedSBarry Smith   PetscOptionsEnd();
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216c4762a1bSJed Brown     Create necessary matrix and vectors
217c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2189566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
2199566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
2209566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATDENSE));
2219566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2229566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
223c4762a1bSJed Brown 
2249566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&U,NULL));
225c4762a1bSJed Brown 
2269566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&Ap));
2279566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Ap,n,1,PETSC_DETERMINE,PETSC_DETERMINE));
2289566063dSJacob Faibussowitsch   PetscCall(MatSetType(Ap,MATDENSE));
2299566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Ap));
2309566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Ap));
2319566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */
232c4762a1bSJed Brown 
2339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u));
234c4762a1bSJed Brown   u[0] = 0;
235c4762a1bSJed Brown   u[1] = 1;
2369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u));
237c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238c4762a1bSJed Brown      Create timestepping solver context
239c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2409566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
2419566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
2429566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSCN));
2439566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,(TSIFunction)IFunction,&app));
2449566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&app));
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
247c4762a1bSJed Brown      Set initial conditions
248c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2499566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,U));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252c4762a1bSJed Brown      Set solver options
253c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2549566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,0.125));
2559566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
2569566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,1./256.));
2579566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Set directions and terminate flags for the two events */
260c4762a1bSJed Brown   direction[0] = 0;
261c4762a1bSJed Brown   terminate[0] = PETSC_FALSE;
2629566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts,1,direction,terminate,EventFunction,PostEventFunction,(void*)&app));
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
265c4762a1bSJed Brown      Run timestepping solver
266c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2679566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,U));
268c4762a1bSJed Brown 
2699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u));
270c4762a1bSJed Brown   tmp[0] = u[0];
271c4762a1bSJed Brown   tmp[1] = u[1];
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   u[0] = 0+delta;
274c4762a1bSJed Brown   u[1] = 1;
2759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u));
276c4762a1bSJed Brown 
2779566063dSJacob Faibussowitsch   PetscCall(FWDRun(ts,U,(void*)&app));
278c4762a1bSJed Brown 
2799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u));
280c4762a1bSJed Brown   sensi[0] = (u[0]-tmp[0])/delta;
281c4762a1bSJed Brown   sensi[1] = (u[1]-tmp[1])/delta;
282*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"d x1(tf) /d x1(t0) = %f d x2(tf) / d x1(t0) = %f \n",(double)sensi[0],(double)sensi[1]));
283c4762a1bSJed Brown   u[0] = 0;
284c4762a1bSJed Brown   u[1] = 1+delta;
2859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u));
286c4762a1bSJed Brown 
2879566063dSJacob Faibussowitsch   PetscCall(FWDRun(ts,U,(void*)&app));
288c4762a1bSJed Brown 
2899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u));
290c4762a1bSJed Brown   sensi[0] = (u[0]-tmp[0])/delta;
291c4762a1bSJed Brown   sensi[1] = (u[1]-tmp[1])/delta;
292*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"d x1(tf) /d x2(t0) = %f d x2(tf) / d x2(t0) = %f \n",(double)sensi[0],(double)sensi[1]));
293c4762a1bSJed Brown   u[0] = 0;
294c4762a1bSJed Brown   u[1] = 1;
295c4762a1bSJed Brown   app.lambda1 = app.lambda1+delta;
2969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u));
297c4762a1bSJed Brown 
2989566063dSJacob Faibussowitsch   PetscCall(FWDRun(ts,U,(void*)&app));
299c4762a1bSJed Brown 
3009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u));
301c4762a1bSJed Brown   sensi[0] = (u[0]-tmp[0])/delta;
302c4762a1bSJed Brown   sensi[1] = (u[1]-tmp[1])/delta;
303*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"Final gradients: d x1(tf) /d p = %f d x2(tf) / d p = %f \n",(double)sensi[0],(double)sensi[1]));
3049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u));
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
307c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
308c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3099566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3119566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
312c4762a1bSJed Brown 
3139566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ap));
3149566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
315b122ec5aSJacob Faibussowitsch   return 0;
316c4762a1bSJed Brown }
317c4762a1bSJed Brown 
318c4762a1bSJed Brown PetscErrorCode FWDRun(TS ts, Vec U0, void *ctx0)
319c4762a1bSJed Brown {
320c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
321c4762a1bSJed Brown   AppCtx         *ctx=(AppCtx*)ctx0;
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   PetscFunctionBeginUser;
3249566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts,&U));
3259566063dSJacob Faibussowitsch   PetscCall(VecCopy(U0,U));
326c4762a1bSJed Brown 
327c4762a1bSJed Brown   ctx->mode = 1;
328c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
329c4762a1bSJed Brown      Run timestepping solver
330c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3319566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
332c4762a1bSJed Brown 
3339566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,U));
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   PetscFunctionReturn(0);
336c4762a1bSJed Brown }
337c4762a1bSJed Brown 
338c4762a1bSJed Brown /*TEST
339c4762a1bSJed Brown 
340c4762a1bSJed Brown    build:
341dfd57a17SPierre Jolivet       requires: !defined(PETSC_USE_CXXCOMPLEX)
342c4762a1bSJed Brown 
343c4762a1bSJed Brown    test:
344c4762a1bSJed Brown       args: -ts_event_tol 1e-9
345c4762a1bSJed Brown       timeoutfactor: 18
346c4762a1bSJed Brown       requires: !single
347c4762a1bSJed Brown 
348c4762a1bSJed Brown TEST*/
349