xref: /petsc/src/ts/tutorials/ex40.c (revision ca4445c7a2f5ca546b532f08b853c371604af17c)
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 
20d71ae5a4SJacob Faibussowitsch PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx)
21d71ae5a4SJacob Faibussowitsch {
22c4762a1bSJed Brown   AppCtx            *app = (AppCtx *)ctx;
23c4762a1bSJed Brown   const PetscScalar *u;
24c4762a1bSJed Brown 
257510d9b0SBarry Smith   PetscFunctionBeginUser;
26c4762a1bSJed Brown   /* Event for ball height */
279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
28c4762a1bSJed Brown   fvalue[0] = u[0];
29c4762a1bSJed Brown   /* Event for number of bounces */
30c4762a1bSJed Brown   fvalue[1] = app->maxbounces - app->nbounces;
319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33c4762a1bSJed Brown }
34c4762a1bSJed Brown 
35d71ae5a4SJacob Faibussowitsch PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx)
36d71ae5a4SJacob Faibussowitsch {
37c4762a1bSJed Brown   AppCtx      *app = (AppCtx *)ctx;
38c4762a1bSJed Brown   PetscScalar *u;
39c4762a1bSJed Brown 
407510d9b0SBarry Smith   PetscFunctionBeginUser;
41c4762a1bSJed Brown   if (event_list[0] == 0) {
429566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Ball hit the ground at t = %5.2f seconds\n", (double)t));
43c4762a1bSJed Brown     /* Set new initial conditions with .9 attenuation */
449566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
45c4762a1bSJed Brown     u[0] = 0.0;
46c4762a1bSJed Brown     u[1] = -0.9 * u[1];
479566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
48c4762a1bSJed Brown   } else if (event_list[0] == 1) {
4963a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Ball bounced %" PetscInt_FMT " times\n", app->nbounces));
50c4762a1bSJed Brown   }
51c4762a1bSJed Brown   app->nbounces++;
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown /*
56c4762a1bSJed Brown      Defines the ODE passed to the ODE solver in explicit form: U_t = F(U)
57c4762a1bSJed Brown */
58d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx)
59d71ae5a4SJacob Faibussowitsch {
60c4762a1bSJed Brown   PetscScalar       *f;
61c4762a1bSJed Brown   const PetscScalar *u;
62c4762a1bSJed Brown 
637510d9b0SBarry Smith   PetscFunctionBeginUser;
64c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
669566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   f[0] = u[1];
69c4762a1bSJed Brown   f[1] = -9.8;
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown /*
77c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of the Jacobian.
78c4762a1bSJed Brown */
79d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx)
80d71ae5a4SJacob Faibussowitsch {
81c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
82c4762a1bSJed Brown   PetscScalar        J[2][2];
83c4762a1bSJed Brown   const PetscScalar *u;
84c4762a1bSJed Brown 
857510d9b0SBarry Smith   PetscFunctionBeginUser;
869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
87c4762a1bSJed Brown 
889371c9d4SSatish Balay   J[0][0] = 0.0;
899371c9d4SSatish Balay   J[0][1] = 1.0;
909371c9d4SSatish Balay   J[1][0] = 0.0;
919371c9d4SSatish Balay   J[1][1] = 0.0;
929566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
95c4762a1bSJed Brown 
969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
98c4762a1bSJed Brown   if (A != B) {
999566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1009566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
101c4762a1bSJed Brown   }
1023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
103c4762a1bSJed Brown }
104c4762a1bSJed Brown 
105c4762a1bSJed Brown /*
106c4762a1bSJed Brown      Defines the ODE passed to the ODE solver in implicit form: F(U_t,U) = 0
107c4762a1bSJed Brown */
108d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
109d71ae5a4SJacob Faibussowitsch {
110c4762a1bSJed Brown   PetscScalar       *f;
111c4762a1bSJed Brown   const PetscScalar *u, *udot;
112c4762a1bSJed Brown 
1137510d9b0SBarry Smith   PetscFunctionBeginUser;
114c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1159566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   f[0] = udot[0] - u[1];
120c4762a1bSJed Brown   f[1] = udot[1] + 9.8;
121c4762a1bSJed Brown 
1229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126c4762a1bSJed Brown }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown /*
129c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
130c4762a1bSJed Brown */
131d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx)
132d71ae5a4SJacob Faibussowitsch {
133c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
134c4762a1bSJed Brown   PetscScalar        J[2][2];
135c4762a1bSJed Brown   const PetscScalar *u, *udot;
136c4762a1bSJed Brown 
1377510d9b0SBarry Smith   PetscFunctionBeginUser;
1389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1399566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
140c4762a1bSJed Brown 
1419371c9d4SSatish Balay   J[0][0] = a;
1429371c9d4SSatish Balay   J[0][1] = -1.0;
1439371c9d4SSatish Balay   J[1][0] = 0.0;
1449371c9d4SSatish Balay   J[1][1] = a;
1459566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
146c4762a1bSJed Brown 
1479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
149c4762a1bSJed Brown 
1509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
152c4762a1bSJed Brown   if (A != B) {
1539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1549566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
155c4762a1bSJed Brown   }
1563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
157c4762a1bSJed Brown }
158c4762a1bSJed Brown 
159d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
160d71ae5a4SJacob Faibussowitsch {
161c4762a1bSJed Brown   TS           ts; /* ODE integrator */
162c4762a1bSJed Brown   Vec          U;  /* solution will be stored here */
163c4762a1bSJed Brown   PetscMPIInt  size;
164c4762a1bSJed Brown   PetscInt     n = 2;
165c4762a1bSJed Brown   PetscScalar *u;
166c4762a1bSJed Brown   AppCtx       app;
167c4762a1bSJed Brown   PetscInt     direction[2];
168c4762a1bSJed Brown   PetscBool    terminate[2];
169c4762a1bSJed Brown   PetscBool    rhs_form = PETSC_FALSE, hist = PETSC_TRUE;
170c4762a1bSJed Brown   TSAdapt      adapt;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173c4762a1bSJed Brown      Initialize program
174c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175327415f7SBarry Smith   PetscFunctionBeginUser;
1769566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1783c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   app.nbounces   = 0;
181c4762a1bSJed Brown   app.maxbounces = 10;
182d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex40 options", "");
1839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-maxbounces", "", "", app.maxbounces, &app.maxbounces, NULL));
1849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_adapthistory", "", "", hist, &hist, NULL));
185d0609cedSBarry Smith   PetscOptionsEnd();
186c4762a1bSJed Brown 
187*ca4445c7SIlya Fursov   Mat A; /* Jacobian matrix */
188*ca4445c7SIlya Fursov   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
189*ca4445c7SIlya Fursov   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
190*ca4445c7SIlya Fursov   PetscCall(MatSetType(A, MATDENSE));
191*ca4445c7SIlya Fursov   PetscCall(MatSetFromOptions(A));
192*ca4445c7SIlya Fursov   PetscCall(MatSetUp(A));
193c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194c4762a1bSJed Brown      Create timestepping solver context
195c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1969566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1979566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200c4762a1bSJed Brown      Set ODE routines
201c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2029566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
203c4762a1bSJed Brown   /* Users are advised against the following branching and code duplication.
204c4762a1bSJed Brown      For problems without a mass matrix like the one at hand, the RHSFunction
205c4762a1bSJed Brown      (and companion RHSJacobian) interface is enough to support both explicit
206c4762a1bSJed Brown      and implicit timesteppers. This tutorial example also deals with the
207c4762a1bSJed Brown      IFunction/IJacobian interface for demonstration and testing purposes. */
2089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-rhs-form", &rhs_form, NULL));
209c4762a1bSJed Brown   if (rhs_form) {
2109566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL));
211*ca4445c7SIlya Fursov     PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
212c4762a1bSJed Brown   } else {
2139566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction, NULL));
2149566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, A, A, IJacobian, NULL));
215c4762a1bSJed Brown   }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218c4762a1bSJed Brown      Set initial conditions
219c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2209566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &U));
2219566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(U, n, PETSC_DETERMINE));
2229566063dSJacob Faibussowitsch   PetscCall(VecSetUp(U));
2239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
224c4762a1bSJed Brown   u[0] = 0.0;
225c4762a1bSJed Brown   u[1] = 20.0;
2269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
2279566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230c4762a1bSJed Brown      Set solver options
231c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2329566063dSJacob Faibussowitsch   if (hist) PetscCall(TSSetSaveTrajectory(ts));
2339566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 30.0));
2349566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2359566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.1));
236a5b23f4aSJose E. Roman   /* The adaptive time step controller could take very large timesteps resulting in
237a5b23f4aSJose E. Roman      the same event occurring multiple times in the same interval. A maximum step size
238c4762a1bSJed Brown      limit is enforced here to avoid this issue. */
2399566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
2409566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5));
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* Set directions and terminate flags for the two events */
2439371c9d4SSatish Balay   direction[0] = -1;
2449371c9d4SSatish Balay   direction[1] = -1;
2459371c9d4SSatish Balay   terminate[0] = PETSC_FALSE;
2469371c9d4SSatish Balay   terminate[1] = PETSC_TRUE;
2479566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
248c4762a1bSJed Brown 
2499566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252c4762a1bSJed Brown      Run timestepping solver
253c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2549566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   if (hist) { /* replay following history */
257c4762a1bSJed Brown     TSTrajectory tj;
258c4762a1bSJed Brown     PetscReal    tf, t0, dt;
259c4762a1bSJed Brown 
260c4762a1bSJed Brown     app.nbounces = 0;
2619566063dSJacob Faibussowitsch     PetscCall(TSGetTime(ts, &tf));
2629566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts, tf));
2639566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, 0));
2649566063dSJacob Faibussowitsch     PetscCall(TSRestartStep(ts));
2659566063dSJacob Faibussowitsch     PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2669566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
2679566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
2689566063dSJacob Faibussowitsch     PetscCall(TSAdaptSetType(adapt, TSADAPTHISTORY));
2699566063dSJacob Faibussowitsch     PetscCall(TSGetTrajectory(ts, &tj));
2709566063dSJacob Faibussowitsch     PetscCall(TSAdaptHistorySetTrajectory(adapt, tj, PETSC_FALSE));
2719566063dSJacob Faibussowitsch     PetscCall(TSAdaptHistoryGetStep(adapt, 0, &t0, &dt));
272c4762a1bSJed Brown     /* this example fails with single (or smaller) precision */
273cc404788SPierre Jolivet #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
274*ca4445c7SIlya Fursov     /*
275*ca4445c7SIlya Fursov        In the first TSSolve() the final time 'tf' is the event location found after a few event handler iterations.
276*ca4445c7SIlya Fursov        If 'tf' is set as the max time for the second run, the TS solver may approach this point by
277*ca4445c7SIlya Fursov        slightly different steps, resulting in a slightly different solution and fvalue[] at 'tf',
278*ca4445c7SIlya Fursov        so that the event may not be triggered at 'tf' anymore. Fix: apply safety factor 1.05
279*ca4445c7SIlya Fursov     */
280*ca4445c7SIlya Fursov     PetscCall(TSSetMaxTime(ts, tf * 1.05));
2819566063dSJacob Faibussowitsch     PetscCall(TSAdaptSetType(adapt, TSADAPTBASIC));
2829566063dSJacob Faibussowitsch     PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5));
2839566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
284c4762a1bSJed Brown #endif
2859566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, t0));
2869566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
2879566063dSJacob Faibussowitsch     PetscCall(TSResetTrajectory(ts));
2889566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
289c4762a1bSJed Brown     u[0] = 0.0;
290c4762a1bSJed Brown     u[1] = 20.0;
2919566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
2929566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, U));
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
295c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
296c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
297*ca4445c7SIlya Fursov   PetscCall(MatDestroy(&A));
2989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2999566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
300c4762a1bSJed Brown 
3019566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
302b122ec5aSJacob Faibussowitsch   return 0;
303c4762a1bSJed Brown }
304c4762a1bSJed Brown 
305c4762a1bSJed Brown /*TEST
306c4762a1bSJed Brown 
307c4762a1bSJed Brown     test:
308c4762a1bSJed Brown       suffix: a
309c4762a1bSJed Brown       args: -snes_stol 1e-4 -ts_trajectory_dirname ex40_a_dir
310c4762a1bSJed Brown       output_file: output/ex40.out
311c4762a1bSJed Brown 
312c4762a1bSJed Brown     test:
313c4762a1bSJed Brown       suffix: b
314c4762a1bSJed Brown       args: -ts_type arkimex -snes_stol 1e-4 -ts_trajectory_dirname ex40_b_dir
315c4762a1bSJed Brown       output_file: output/ex40.out
316c4762a1bSJed Brown 
317c4762a1bSJed Brown     test:
318c4762a1bSJed Brown       suffix: c
319*ca4445c7SIlya Fursov       args: -snes_mf_operator -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_c_dir
320*ca4445c7SIlya Fursov       output_file: output/ex40.out
321*ca4445c7SIlya Fursov 
322*ca4445c7SIlya Fursov     test:
323*ca4445c7SIlya Fursov       suffix: cr
324*ca4445c7SIlya Fursov       args: -rhs-form -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_cr_dir
325*ca4445c7SIlya Fursov       output_file: output/ex40.out
326*ca4445c7SIlya Fursov 
327*ca4445c7SIlya Fursov     test:
328*ca4445c7SIlya Fursov       suffix: crmf
329*ca4445c7SIlya Fursov       args: -rhs-form -snes_mf_operator -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_crmf_dir
330c4762a1bSJed Brown       output_file: output/ex40.out
331c4762a1bSJed Brown 
332c4762a1bSJed Brown     test:
333c4762a1bSJed Brown       suffix: d
334c4762a1bSJed Brown       args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_d_dir
335c4762a1bSJed Brown       output_file: output/ex40.out
336c4762a1bSJed Brown 
337c4762a1bSJed Brown     test:
338c4762a1bSJed Brown       suffix: e
339c4762a1bSJed Brown       args: -ts_type bdf -ts_adapt_dt_max 0.025 -ts_max_steps 1500 -ts_trajectory_dirname ex40_e_dir
340c4762a1bSJed Brown       output_file: output/ex40.out
341c4762a1bSJed Brown 
342c4762a1bSJed Brown     test:
343c4762a1bSJed Brown       suffix: f
344c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 3bs -ts_trajectory_dirname ex40_f_dir
345c4762a1bSJed Brown       output_file: output/ex40.out
346c4762a1bSJed Brown 
347c4762a1bSJed Brown     test:
348c4762a1bSJed Brown       suffix: g
349c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 5bs -ts_trajectory_dirname ex40_g_dir
350c4762a1bSJed Brown       output_file: output/ex40.out
351c4762a1bSJed Brown 
352c4762a1bSJed Brown     test:
353c4762a1bSJed Brown       suffix: h
354c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 6vr -ts_trajectory_dirname ex40_h_dir
355c4762a1bSJed Brown       output_file: output/ex40.out
356c4762a1bSJed Brown 
357c4762a1bSJed Brown     test:
358c4762a1bSJed Brown       suffix: i
359c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 7vr -ts_trajectory_dirname ex40_i_dir
360c4762a1bSJed Brown       output_file: output/ex40.out
361c4762a1bSJed Brown 
362c4762a1bSJed Brown     test:
363c4762a1bSJed Brown       suffix: j
364c4762a1bSJed Brown       args: -rhs-form -ts_type rk -ts_rk_type 8vr -ts_trajectory_dirname ex40_j_dir
3655385558fSLisandro Dalcin       output_file: output/ex40.out
3665385558fSLisandro Dalcin 
3675385558fSLisandro Dalcin     test:
3685385558fSLisandro Dalcin       suffix: k
3695385558fSLisandro Dalcin       args: -ts_type theta -ts_adapt_type dsp -ts_trajectory_dirname ex40_k_dir
3705385558fSLisandro Dalcin       output_file: output/ex40.out
3715385558fSLisandro Dalcin 
3725385558fSLisandro Dalcin     test:
3735385558fSLisandro Dalcin       suffix: l
3745385558fSLisandro Dalcin       args: -rhs-form -ts_type rk -ts_rk_type 2a -ts_trajectory_dirname ex40_l_dir
3755385558fSLisandro Dalcin       args: -ts_adapt_type dsp -ts_adapt_always_accept {{false true}} -ts_adapt_dt_min 0.01
376c4762a1bSJed Brown       output_file: output/ex40.out
377c4762a1bSJed Brown 
37808c4bdbbSLisandro Dalcin     test:
37908c4bdbbSLisandro Dalcin       suffix: m
38008c4bdbbSLisandro Dalcin       args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -test_adapthistory false
38108c4bdbbSLisandro Dalcin       args: -ts_max_time 10 -ts_exact_final_time {{STEPOVER MATCHSTEP INTERPOLATE}}
38208c4bdbbSLisandro Dalcin 
38308c4bdbbSLisandro Dalcin     test:
38408c4bdbbSLisandro Dalcin       requires: !single
38508c4bdbbSLisandro Dalcin       suffix: n
38608c4bdbbSLisandro Dalcin       args: -test_adapthistory false
38708c4bdbbSLisandro Dalcin       args: -ts_type alpha -ts_alpha_radius 1.0 -ts_view
38808c4bdbbSLisandro Dalcin       args: -ts_dt 0.25 -ts_adapt_type basic -ts_adapt_wnormtype INFINITY -ts_adapt_monitor
38908c4bdbbSLisandro Dalcin       args: -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false
39008c4bdbbSLisandro Dalcin 
391e7685601SHong Zhang     test:
392e7685601SHong Zhang       requires: !single
393e7685601SHong Zhang       suffix: o
394e7685601SHong Zhang       args: -rhs-form -ts_type rk -ts_rk_type 2b -ts_trajectory_dirname ex40_o_dir
395e7685601SHong Zhang       output_file: output/ex40.out
396c4762a1bSJed Brown TEST*/
397