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