xref: /petsc/src/ts/tutorials/ex40.c (revision 08c4bdbb682ef9cc0b4e244fe978db2ed35e9be3)
1c4762a1bSJed Brown static char help[] = "Serial bouncing ball example to test TS event feature.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   The dynamics of the bouncing ball is described by the ODE
5c4762a1bSJed Brown                   u1_t = u2
6c4762a1bSJed Brown                   u2_t = -9.8
7c4762a1bSJed Brown 
8c4762a1bSJed Brown   There are two events set in this example. The first one checks for the ball hitting the
9c4762a1bSJed Brown   ground (u1 = 0). Every time the ball hits the ground, its velocity u2 is attenuated by
10c4762a1bSJed Brown   a factor of 0.9. The second event sets a limit on the number of ball bounces.
11c4762a1bSJed Brown */
12c4762a1bSJed Brown 
13c4762a1bSJed Brown #include <petscts.h>
14c4762a1bSJed Brown 
15c4762a1bSJed Brown typedef struct {
16c4762a1bSJed Brown   PetscInt maxbounces;
17c4762a1bSJed Brown   PetscInt nbounces;
18c4762a1bSJed Brown } AppCtx;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx)
21c4762a1bSJed Brown {
22c4762a1bSJed Brown   AppCtx            *app=(AppCtx*)ctx;
23c4762a1bSJed Brown   PetscErrorCode    ierr;
24c4762a1bSJed Brown   const PetscScalar *u;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   PetscFunctionBegin;
27c4762a1bSJed Brown   /* Event for ball height */
28c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
29c4762a1bSJed Brown   fvalue[0] = u[0];
30c4762a1bSJed Brown   /* Event for number of bounces */
31c4762a1bSJed Brown   fvalue[1] = app->maxbounces - app->nbounces;
32c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
33c4762a1bSJed Brown   PetscFunctionReturn(0);
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
36c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
37c4762a1bSJed Brown {
38c4762a1bSJed Brown   AppCtx         *app=(AppCtx*)ctx;
39c4762a1bSJed Brown   PetscErrorCode ierr;
40c4762a1bSJed Brown   PetscScalar    *u;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   PetscFunctionBegin;
43c4762a1bSJed Brown   if (event_list[0] == 0) {
44c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Ball hit the ground at t = %5.2f seconds\n",(double)t);CHKERRQ(ierr);
45c4762a1bSJed Brown     /* Set new initial conditions with .9 attenuation */
46c4762a1bSJed Brown     ierr = VecGetArray(U,&u);CHKERRQ(ierr);
47c4762a1bSJed Brown     u[0] =  0.0;
48c4762a1bSJed Brown     u[1] = -0.9*u[1];
49c4762a1bSJed Brown     ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
50c4762a1bSJed Brown   } else if (event_list[0] == 1) {
51c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Ball bounced %D times\n",app->nbounces);CHKERRQ(ierr);
52c4762a1bSJed Brown   }
53c4762a1bSJed Brown   app->nbounces++;
54c4762a1bSJed Brown   PetscFunctionReturn(0);
55c4762a1bSJed Brown }
56c4762a1bSJed Brown 
57c4762a1bSJed Brown /*
58c4762a1bSJed Brown      Defines the ODE passed to the ODE solver in explicit form: U_t = F(U)
59c4762a1bSJed Brown */
60c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
61c4762a1bSJed Brown {
62c4762a1bSJed Brown   PetscErrorCode    ierr;
63c4762a1bSJed Brown   PetscScalar       *f;
64c4762a1bSJed Brown   const PetscScalar *u;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   PetscFunctionBegin;
67c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
68c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   f[0] = u[1];
72c4762a1bSJed Brown   f[1] = - 9.8;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
76c4762a1bSJed Brown   PetscFunctionReturn(0);
77c4762a1bSJed Brown }
78c4762a1bSJed Brown 
79c4762a1bSJed Brown /*
80c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of the Jacobian.
81c4762a1bSJed Brown */
82c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
83c4762a1bSJed Brown {
84c4762a1bSJed Brown   PetscErrorCode    ierr;
85c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
86c4762a1bSJed Brown   PetscScalar       J[2][2];
87c4762a1bSJed Brown   const PetscScalar *u;
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   PetscFunctionBegin;
90c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   J[0][0] = 0.0;     J[0][1] = 1.0;
93c4762a1bSJed Brown   J[1][0] = 0.0;     J[1][1] = 0.0;
94c4762a1bSJed Brown   ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100c4762a1bSJed Brown   if (A != B) {
101c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
102c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown   PetscFunctionReturn(0);
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /*
108c4762a1bSJed Brown      Defines the ODE passed to the ODE solver in implicit form: F(U_t,U) = 0
109c4762a1bSJed Brown */
110c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
111c4762a1bSJed Brown {
112c4762a1bSJed Brown   PetscErrorCode    ierr;
113c4762a1bSJed Brown   PetscScalar       *f;
114c4762a1bSJed Brown   const PetscScalar *u,*udot;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   PetscFunctionBegin;
117c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
118c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   f[0] = udot[0] - u[1];
123c4762a1bSJed Brown   f[1] = udot[1] + 9.8;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
128c4762a1bSJed Brown   PetscFunctionReturn(0);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown /*
132c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
133c4762a1bSJed Brown */
134c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
135c4762a1bSJed Brown {
136c4762a1bSJed Brown   PetscErrorCode    ierr;
137c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
138c4762a1bSJed Brown   PetscScalar       J[2][2];
139c4762a1bSJed Brown   const PetscScalar *u,*udot;
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   PetscFunctionBegin;
142c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   J[0][0] = a;      J[0][1] = -1.0;
146c4762a1bSJed Brown   J[1][0] = 0.0;    J[1][1] = a;
147c4762a1bSJed Brown   ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
153c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154c4762a1bSJed Brown   if (A != B) {
155c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157c4762a1bSJed Brown   }
158c4762a1bSJed Brown   PetscFunctionReturn(0);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown int main(int argc,char **argv)
162c4762a1bSJed Brown {
163c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
164c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
165c4762a1bSJed Brown   PetscErrorCode ierr;
166c4762a1bSJed Brown   PetscMPIInt    size;
167c4762a1bSJed Brown   PetscInt       n = 2;
168c4762a1bSJed Brown   PetscScalar    *u;
169c4762a1bSJed Brown   AppCtx         app;
170c4762a1bSJed Brown   PetscInt       direction[2];
171c4762a1bSJed Brown   PetscBool      terminate[2];
172c4762a1bSJed Brown   PetscBool      rhs_form=PETSC_FALSE,hist=PETSC_TRUE;
173c4762a1bSJed Brown   TSAdapt        adapt;
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176c4762a1bSJed Brown      Initialize program
177c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
179c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
180c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   app.nbounces = 0;
183c4762a1bSJed Brown   app.maxbounces = 10;
184c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex40 options","");CHKERRQ(ierr);
185c4762a1bSJed Brown   ierr = PetscOptionsInt("-maxbounces","","",app.maxbounces,&app.maxbounces,NULL);CHKERRQ(ierr);
186c4762a1bSJed Brown   ierr = PetscOptionsBool("-test_adapthistory","","",hist,&hist,NULL);CHKERRQ(ierr);
187c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190c4762a1bSJed Brown      Create timestepping solver context
191c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
193c4762a1bSJed Brown   ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr);
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196c4762a1bSJed Brown      Set ODE routines
197c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
199c4762a1bSJed Brown   /* Users are advised against the following branching and code duplication.
200c4762a1bSJed Brown      For problems without a mass matrix like the one at hand, the RHSFunction
201c4762a1bSJed Brown      (and companion RHSJacobian) interface is enough to support both explicit
202c4762a1bSJed Brown      and implicit timesteppers. This tutorial example also deals with the
203c4762a1bSJed Brown      IFunction/IJacobian interface for demonstration and testing purposes. */
204c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-rhs-form",&rhs_form,NULL);CHKERRQ(ierr);
205c4762a1bSJed Brown   if (rhs_form) {
206c4762a1bSJed Brown     ierr = TSSetRHSFunction(ts,NULL,RHSFunction,NULL);CHKERRQ(ierr);
207c4762a1bSJed Brown     ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,NULL);CHKERRQ(ierr);
208c4762a1bSJed Brown   } else {
209c4762a1bSJed Brown     Mat A; /* Jacobian matrix */
210c4762a1bSJed Brown     ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
211c4762a1bSJed Brown     ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
212c4762a1bSJed Brown     ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
213c4762a1bSJed Brown     ierr = MatSetFromOptions(A);CHKERRQ(ierr);
214c4762a1bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
215c4762a1bSJed Brown     ierr = TSSetIFunction(ts,NULL,IFunction,NULL);CHKERRQ(ierr);
216c4762a1bSJed Brown     ierr = TSSetIJacobian(ts,A,A,IJacobian,NULL);CHKERRQ(ierr);
217c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
218c4762a1bSJed Brown   }
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221c4762a1bSJed Brown      Set initial conditions
222c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
223c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&U);CHKERRQ(ierr);
224c4762a1bSJed Brown   ierr = VecSetSizes(U,n,PETSC_DETERMINE);CHKERRQ(ierr);
225c4762a1bSJed Brown   ierr = VecSetUp(U);CHKERRQ(ierr);
226c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
227c4762a1bSJed Brown   u[0] = 0.0;
228c4762a1bSJed Brown   u[1] = 20.0;
229c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
230c4762a1bSJed Brown   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233c4762a1bSJed Brown      Set solver options
234c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
235*08c4bdbbSLisandro Dalcin   if (hist) {ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);}
236c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,30.0);CHKERRQ(ierr);
237c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
238c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.1);CHKERRQ(ierr);
239c4762a1bSJed Brown   /* The adapative time step controller could take very large timesteps resulting in
240c4762a1bSJed Brown      the same event occuring multiple times in the same interval. A maximum step size
241c4762a1bSJed Brown      limit is enforced here to avoid this issue. */
242c4762a1bSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
243c4762a1bSJed Brown   ierr = TSAdaptSetStepLimits(adapt,0.0,0.5);CHKERRQ(ierr);
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /* Set directions and terminate flags for the two events */
246c4762a1bSJed Brown   direction[0] = -1;            direction[1] = -1;
247c4762a1bSJed Brown   terminate[0] = PETSC_FALSE;   terminate[1] = PETSC_TRUE;
248c4762a1bSJed Brown   ierr = TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&app);CHKERRQ(ierr);
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
253c4762a1bSJed Brown      Run timestepping solver
254c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
255c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   if (hist) { /* replay following history */
258c4762a1bSJed Brown     TSTrajectory tj;
259c4762a1bSJed Brown     PetscReal    tf,t0,dt;
260c4762a1bSJed Brown 
261c4762a1bSJed Brown     app.nbounces = 0;
262c4762a1bSJed Brown     ierr = TSGetTime(ts,&tf);CHKERRQ(ierr);
263c4762a1bSJed Brown     ierr = TSSetMaxTime(ts,tf);CHKERRQ(ierr);
264c4762a1bSJed Brown     ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
265c4762a1bSJed Brown     ierr = TSRestartStep(ts);CHKERRQ(ierr);
266c4762a1bSJed Brown     ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
267c4762a1bSJed Brown     ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
268c4762a1bSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
269c4762a1bSJed Brown     ierr = TSAdaptSetType(adapt,TSADAPTHISTORY);CHKERRQ(ierr);
270c4762a1bSJed Brown     ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr);
271c4762a1bSJed Brown     ierr = TSAdaptHistorySetTrajectory(adapt,tj,PETSC_FALSE);CHKERRQ(ierr);
272c4762a1bSJed Brown     ierr = TSAdaptHistoryGetStep(adapt,0,&t0,&dt);CHKERRQ(ierr);
273c4762a1bSJed Brown     /* this example fails with single (or smaller) precision */
274c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL__FP16)
275c4762a1bSJed Brown     ierr = TSAdaptSetType(adapt,TSADAPTBASIC);CHKERRQ(ierr);
276c4762a1bSJed Brown     ierr = TSAdaptSetStepLimits(adapt,0.0,0.5);CHKERRQ(ierr);
277c4762a1bSJed Brown     ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
278c4762a1bSJed Brown #endif
279c4762a1bSJed Brown     ierr = TSSetTime(ts,t0);CHKERRQ(ierr);
280c4762a1bSJed Brown     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
281c4762a1bSJed Brown     ierr = TSResetTrajectory(ts);CHKERRQ(ierr);
282c4762a1bSJed Brown     ierr = VecGetArray(U,&u);CHKERRQ(ierr);
283c4762a1bSJed Brown     u[0] = 0.0;
284c4762a1bSJed Brown     u[1] = 20.0;
285c4762a1bSJed Brown     ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
286c4762a1bSJed Brown     ierr = TSSolve(ts,U);CHKERRQ(ierr);
287c4762a1bSJed Brown   }
288c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
289c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
290c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
291c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
292c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   ierr = PetscFinalize();
295c4762a1bSJed Brown   return ierr;
296c4762a1bSJed Brown }
297c4762a1bSJed Brown 
298c4762a1bSJed Brown /*TEST
299c4762a1bSJed Brown 
300c4762a1bSJed Brown     test:
301c4762a1bSJed Brown       suffix: a
302c4762a1bSJed Brown       args: -snes_stol 1e-4 -ts_trajectory_dirname ex40_a_dir
303c4762a1bSJed Brown       output_file: output/ex40.out
304c4762a1bSJed Brown 
305c4762a1bSJed Brown     test:
306c4762a1bSJed Brown       suffix: b
307c4762a1bSJed Brown       args: -ts_type arkimex -snes_stol 1e-4 -ts_trajectory_dirname ex40_b_dir
308c4762a1bSJed Brown       output_file: output/ex40.out
309c4762a1bSJed Brown 
310c4762a1bSJed Brown     test:
311c4762a1bSJed Brown       suffix: c
312c4762a1bSJed Brown       args: -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_c_dir
313c4762a1bSJed Brown       output_file: output/ex40.out
314c4762a1bSJed Brown 
315c4762a1bSJed Brown     test:
316c4762a1bSJed Brown       suffix: d
317c4762a1bSJed Brown       args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_d_dir
318c4762a1bSJed Brown       output_file: output/ex40.out
319c4762a1bSJed Brown 
320c4762a1bSJed Brown     test:
321c4762a1bSJed Brown       suffix: e
322c4762a1bSJed Brown       args: -ts_type bdf -ts_adapt_dt_max 0.025 -ts_max_steps 1500 -ts_trajectory_dirname ex40_e_dir
323c4762a1bSJed Brown       output_file: output/ex40.out
324c4762a1bSJed Brown 
325c4762a1bSJed Brown     test:
326c4762a1bSJed Brown       suffix: f
327c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 3bs -ts_trajectory_dirname ex40_f_dir
328c4762a1bSJed Brown       output_file: output/ex40.out
329c4762a1bSJed Brown 
330c4762a1bSJed Brown     test:
331c4762a1bSJed Brown       suffix: g
332c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 5bs -ts_trajectory_dirname ex40_g_dir
333c4762a1bSJed Brown       output_file: output/ex40.out
334c4762a1bSJed Brown 
335c4762a1bSJed Brown     test:
336c4762a1bSJed Brown       suffix: h
337c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 6vr -ts_trajectory_dirname ex40_h_dir
338c4762a1bSJed Brown       output_file: output/ex40.out
339c4762a1bSJed Brown 
340c4762a1bSJed Brown     test:
341c4762a1bSJed Brown       suffix: i
342c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 7vr -ts_trajectory_dirname ex40_i_dir
343c4762a1bSJed Brown       output_file: output/ex40.out
344c4762a1bSJed Brown 
345c4762a1bSJed Brown     test:
346c4762a1bSJed Brown       suffix: j
347c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 8vr -ts_trajectory_dirname ex40_j_dir
3485385558fSLisandro Dalcin       output_file: output/ex40.out
3495385558fSLisandro Dalcin 
3505385558fSLisandro Dalcin     test:
3515385558fSLisandro Dalcin       suffix: k
3525385558fSLisandro Dalcin       args: -ts_type theta -ts_adapt_type dsp -ts_trajectory_dirname ex40_k_dir
3535385558fSLisandro Dalcin       output_file: output/ex40.out
3545385558fSLisandro Dalcin 
3555385558fSLisandro Dalcin     test:
3565385558fSLisandro Dalcin       suffix: l
3575385558fSLisandro Dalcin       args: -rhs-form -ts_type rk -ts_rk_type 2a -ts_trajectory_dirname ex40_l_dir
3585385558fSLisandro Dalcin       args: -ts_adapt_type dsp -ts_adapt_always_accept {{false true}} -ts_adapt_dt_min 0.01
359c4762a1bSJed Brown       output_file: output/ex40.out
360c4762a1bSJed Brown 
361*08c4bdbbSLisandro Dalcin     test:
362*08c4bdbbSLisandro Dalcin       suffix: m
363*08c4bdbbSLisandro Dalcin       args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -test_adapthistory false
364*08c4bdbbSLisandro Dalcin       args: -ts_max_time 10 -ts_exact_final_time {{STEPOVER MATCHSTEP INTERPOLATE}}
365*08c4bdbbSLisandro Dalcin 
366*08c4bdbbSLisandro Dalcin     test:
367*08c4bdbbSLisandro Dalcin       requires: !single
368*08c4bdbbSLisandro Dalcin       suffix: n
369*08c4bdbbSLisandro Dalcin       args: -test_adapthistory false
370*08c4bdbbSLisandro Dalcin       args: -ts_type alpha -ts_alpha_radius 1.0 -ts_view
371*08c4bdbbSLisandro Dalcin       args: -ts_dt 0.25 -ts_adapt_type basic -ts_adapt_wnormtype INFINITY -ts_adapt_monitor
372*08c4bdbbSLisandro Dalcin       args: -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false
373*08c4bdbbSLisandro Dalcin 
374c4762a1bSJed Brown TEST*/
375