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
8fa9584fbSIlya Fursov There is one event set in this example, which checks for the ball hitting the
9c4762a1bSJed Brown ground (u1 = 0). Every time the ball hits the ground, its velocity u2 is attenuated by
10fa9584fbSIlya Fursov a factor of 0.9. On reaching the limit on the number of ball bounces,
11fa9584fbSIlya Fursov the TS run is requested to terminate from the PostEvent() callback.
12c4762a1bSJed Brown */
13c4762a1bSJed Brown
14c4762a1bSJed Brown #include <petscts.h>
15c4762a1bSJed Brown
16c4762a1bSJed Brown typedef struct {
17c4762a1bSJed Brown PetscInt maxbounces;
18c4762a1bSJed Brown PetscInt nbounces;
19c4762a1bSJed Brown } AppCtx;
20c4762a1bSJed Brown
EventFunction(TS ts,PetscReal t,Vec U,PetscReal * fvalue,PetscCtx ctx)21*2a8381b2SBarry Smith PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal *fvalue, PetscCtx ctx)
22d71ae5a4SJacob Faibussowitsch {
23c4762a1bSJed Brown const PetscScalar *u;
24c4762a1bSJed Brown
257510d9b0SBarry Smith PetscFunctionBeginUser;
26c4762a1bSJed Brown /* Event for ball height */
279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
28fa9584fbSIlya Fursov fvalue[0] = PetscRealPart(u[0]);
299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
31c4762a1bSJed Brown }
32c4762a1bSJed Brown
PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,PetscCtx ctx)33*2a8381b2SBarry Smith PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, PetscCtx ctx)
34d71ae5a4SJacob Faibussowitsch {
35c4762a1bSJed Brown AppCtx *app = (AppCtx *)ctx;
36c4762a1bSJed Brown PetscScalar *u;
37c4762a1bSJed Brown
387510d9b0SBarry Smith PetscFunctionBeginUser;
399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Ball hit the ground at t = %5.2f seconds\n", (double)t));
40c4762a1bSJed Brown /* Set new initial conditions with .9 attenuation */
419566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u));
42c4762a1bSJed Brown u[0] = 0.0;
43c4762a1bSJed Brown u[1] = -0.9 * u[1];
449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u));
45c4762a1bSJed Brown app->nbounces++;
46fa9584fbSIlya Fursov if (app->nbounces >= app->maxbounces) {
47fa9584fbSIlya Fursov PetscCall(PetscPrintf(PETSC_COMM_SELF, "Ball bounced %" PetscInt_FMT " times\n", app->nbounces));
48fa9584fbSIlya Fursov PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); // request TS to terminate; since the program is serial, no need to sync this call
49fa9584fbSIlya Fursov }
503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
51c4762a1bSJed Brown }
52c4762a1bSJed Brown
53c4762a1bSJed Brown /*
54c4762a1bSJed Brown Defines the ODE passed to the ODE solver in explicit form: U_t = F(U)
55c4762a1bSJed Brown */
RHSFunction(TS ts,PetscReal t,Vec U,Vec F,PetscCtx ctx)56*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, PetscCtx ctx)
57d71ae5a4SJacob Faibussowitsch {
58c4762a1bSJed Brown PetscScalar *f;
59c4762a1bSJed Brown const PetscScalar *u;
60c4762a1bSJed Brown
617510d9b0SBarry Smith PetscFunctionBeginUser;
62fa9584fbSIlya Fursov /* The following lines allow us to access the entries of the vectors directly */
639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
649566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
65c4762a1bSJed Brown
66c4762a1bSJed Brown f[0] = u[1];
67c4762a1bSJed Brown f[1] = -9.8;
68c4762a1bSJed Brown
699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
72c4762a1bSJed Brown }
73c4762a1bSJed Brown
74c4762a1bSJed Brown /*
75c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of the Jacobian.
76c4762a1bSJed Brown */
RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,PetscCtx ctx)77*2a8381b2SBarry Smith static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, PetscCtx ctx)
78d71ae5a4SJacob Faibussowitsch {
79c4762a1bSJed Brown PetscInt rowcol[] = {0, 1};
80c4762a1bSJed Brown PetscScalar J[2][2];
81c4762a1bSJed Brown const PetscScalar *u;
82c4762a1bSJed Brown
837510d9b0SBarry Smith PetscFunctionBeginUser;
849566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
85c4762a1bSJed Brown
869371c9d4SSatish Balay J[0][0] = 0.0;
879371c9d4SSatish Balay J[0][1] = 1.0;
889371c9d4SSatish Balay J[1][0] = 0.0;
899371c9d4SSatish Balay J[1][1] = 0.0;
909566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
91c4762a1bSJed Brown
929566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
93c4762a1bSJed Brown
949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
96fa9584fbSIlya Fursov if (A != B) {
97fa9584fbSIlya Fursov PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
98fa9584fbSIlya Fursov PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
99c4762a1bSJed Brown }
1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
101c4762a1bSJed Brown }
102c4762a1bSJed Brown
103c4762a1bSJed Brown /*
104c4762a1bSJed Brown Defines the ODE passed to the ODE solver in implicit form: F(U_t,U) = 0
105c4762a1bSJed Brown */
IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,PetscCtx ctx)106*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx)
107d71ae5a4SJacob Faibussowitsch {
108c4762a1bSJed Brown PetscScalar *f;
109c4762a1bSJed Brown const PetscScalar *u, *udot;
110c4762a1bSJed Brown
1117510d9b0SBarry Smith PetscFunctionBeginUser;
112c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */
1139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
1149566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot));
1159566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
116c4762a1bSJed Brown
117c4762a1bSJed Brown f[0] = udot[0] - u[1];
118c4762a1bSJed Brown f[1] = udot[1] + 9.8;
119c4762a1bSJed Brown
1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
1219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot));
1229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
124c4762a1bSJed Brown }
125c4762a1bSJed Brown
126c4762a1bSJed Brown /*
127c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
128c4762a1bSJed Brown */
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,PetscCtx ctx)129*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
130d71ae5a4SJacob Faibussowitsch {
131c4762a1bSJed Brown PetscInt rowcol[] = {0, 1};
132c4762a1bSJed Brown PetscScalar J[2][2];
133c4762a1bSJed Brown const PetscScalar *u, *udot;
134c4762a1bSJed Brown
1357510d9b0SBarry Smith PetscFunctionBeginUser;
1369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
1379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot));
138c4762a1bSJed Brown
1399371c9d4SSatish Balay J[0][0] = a;
1409371c9d4SSatish Balay J[0][1] = -1.0;
1419371c9d4SSatish Balay J[1][0] = 0.0;
1429371c9d4SSatish Balay J[1][1] = a;
1439566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
144c4762a1bSJed Brown
1459566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
1469566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot));
147c4762a1bSJed Brown
1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
150fa9584fbSIlya Fursov if (A != B) {
151fa9584fbSIlya Fursov PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
152fa9584fbSIlya Fursov PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
153c4762a1bSJed Brown }
1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
155c4762a1bSJed Brown }
156c4762a1bSJed Brown
main(int argc,char ** argv)157d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
158d71ae5a4SJacob Faibussowitsch {
159c4762a1bSJed Brown TS ts; /* ODE integrator */
160c4762a1bSJed Brown Vec U; /* solution will be stored here */
161c4762a1bSJed Brown PetscMPIInt size;
162c4762a1bSJed Brown PetscInt n = 2;
163c4762a1bSJed Brown PetscScalar *u;
164c4762a1bSJed Brown AppCtx app;
165fa9584fbSIlya Fursov PetscInt direction[1];
166fa9584fbSIlya Fursov PetscBool terminate[1];
167c4762a1bSJed Brown PetscBool rhs_form = PETSC_FALSE, hist = PETSC_TRUE;
168c4762a1bSJed Brown TSAdapt adapt;
169c4762a1bSJed Brown
170c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171c4762a1bSJed Brown Initialize program
172c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173327415f7SBarry Smith PetscFunctionBeginUser;
174c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1763c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
177c4762a1bSJed Brown
178c4762a1bSJed Brown app.nbounces = 0;
179c4762a1bSJed Brown app.maxbounces = 10;
180d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex40 options", "");
1819566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-maxbounces", "", "", app.maxbounces, &app.maxbounces, NULL));
1829566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_adapthistory", "", "", hist, &hist, NULL));
183d0609cedSBarry Smith PetscOptionsEnd();
184c4762a1bSJed Brown
185ca4445c7SIlya Fursov Mat A; /* Jacobian matrix */
186ca4445c7SIlya Fursov PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
187ca4445c7SIlya Fursov PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
188ca4445c7SIlya Fursov PetscCall(MatSetType(A, MATDENSE));
189ca4445c7SIlya Fursov PetscCall(MatSetFromOptions(A));
190ca4445c7SIlya Fursov PetscCall(MatSetUp(A));
191c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192c4762a1bSJed Brown Create timestepping solver context
193c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1949566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1959566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSROSW));
196c4762a1bSJed Brown
197c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown Set ODE routines
199c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2009566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
201c4762a1bSJed Brown /* Users are advised against the following branching and code duplication.
202c4762a1bSJed Brown For problems without a mass matrix like the one at hand, the RHSFunction
203c4762a1bSJed Brown (and companion RHSJacobian) interface is enough to support both explicit
204c4762a1bSJed Brown and implicit timesteppers. This tutorial example also deals with the
205c4762a1bSJed Brown IFunction/IJacobian interface for demonstration and testing purposes. */
2069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-rhs-form", &rhs_form, NULL));
207c4762a1bSJed Brown if (rhs_form) {
2089566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL));
209ca4445c7SIlya Fursov PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
210c4762a1bSJed Brown } else {
2119566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, NULL));
2129566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, NULL));
213c4762a1bSJed Brown }
214c4762a1bSJed Brown
215c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216c4762a1bSJed Brown Set initial conditions
217c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2189566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &U));
2199566063dSJacob Faibussowitsch PetscCall(VecSetSizes(U, n, PETSC_DETERMINE));
2209566063dSJacob Faibussowitsch PetscCall(VecSetUp(U));
2219566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u));
222c4762a1bSJed Brown u[0] = 0.0;
223c4762a1bSJed Brown u[1] = 20.0;
2249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u));
2259566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U));
226c4762a1bSJed Brown
227c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228c4762a1bSJed Brown Set solver options
229c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2309566063dSJacob Faibussowitsch if (hist) PetscCall(TSSetSaveTrajectory(ts));
2319566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 30.0));
2329566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2339566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.1));
234fa9584fbSIlya Fursov /* The adaptive time step controller could take very large timesteps
235fa9584fbSIlya Fursov jumping over the next event zero-crossing point. A maximum step size
236c4762a1bSJed Brown limit is enforced here to avoid this issue. */
2379566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt));
2389566063dSJacob Faibussowitsch PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5));
239c4762a1bSJed Brown
240fa9584fbSIlya Fursov /* Set direction and terminate flag for the event */
2419371c9d4SSatish Balay direction[0] = -1;
2429371c9d4SSatish Balay terminate[0] = PETSC_FALSE;
243fa9584fbSIlya Fursov PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
244c4762a1bSJed Brown
2459566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
246c4762a1bSJed Brown
247c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248c4762a1bSJed Brown Run timestepping solver
249c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2509566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
251c4762a1bSJed Brown
252c4762a1bSJed Brown if (hist) { /* replay following history */
253c4762a1bSJed Brown TSTrajectory tj;
254c4762a1bSJed Brown PetscReal tf, t0, dt;
255c4762a1bSJed Brown
256c4762a1bSJed Brown app.nbounces = 0;
2579566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &tf));
2589566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, tf));
2599566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0));
2609566063dSJacob Faibussowitsch PetscCall(TSRestartStep(ts));
2619566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2629566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
2639566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt));
2649566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt, TSADAPTHISTORY));
2659566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts, &tj));
2669566063dSJacob Faibussowitsch PetscCall(TSAdaptHistorySetTrajectory(adapt, tj, PETSC_FALSE));
2679566063dSJacob Faibussowitsch PetscCall(TSAdaptHistoryGetStep(adapt, 0, &t0, &dt));
268c4762a1bSJed Brown /* this example fails with single (or smaller) precision */
269cc404788SPierre Jolivet #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
270ca4445c7SIlya Fursov /*
271ca4445c7SIlya Fursov In the first TSSolve() the final time 'tf' is the event location found after a few event handler iterations.
272ca4445c7SIlya Fursov If 'tf' is set as the max time for the second run, the TS solver may approach this point by
273ca4445c7SIlya Fursov slightly different steps, resulting in a slightly different solution and fvalue[] at 'tf',
274ca4445c7SIlya Fursov so that the event may not be triggered at 'tf' anymore. Fix: apply safety factor 1.05
275ca4445c7SIlya Fursov */
276ca4445c7SIlya Fursov PetscCall(TSSetMaxTime(ts, tf * 1.05));
2779566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt, TSADAPTBASIC));
2789566063dSJacob Faibussowitsch PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5));
2799566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
280c4762a1bSJed Brown #endif
2819566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, t0));
2829566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt));
2839566063dSJacob Faibussowitsch PetscCall(TSResetTrajectory(ts));
2849566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u));
285c4762a1bSJed Brown u[0] = 0.0;
286c4762a1bSJed Brown u[1] = 20.0;
2879566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u));
2889566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
289c4762a1bSJed Brown }
290c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
291c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed.
292c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
293ca4445c7SIlya Fursov PetscCall(MatDestroy(&A));
2949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U));
2959566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
296c4762a1bSJed Brown
2979566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
298b122ec5aSJacob Faibussowitsch return 0;
299c4762a1bSJed Brown }
300c4762a1bSJed Brown
301c4762a1bSJed Brown /*TEST
302c4762a1bSJed Brown
303c4762a1bSJed Brown test:
304c4762a1bSJed Brown suffix: a
305c4762a1bSJed Brown args: -snes_stol 1e-4 -ts_trajectory_dirname ex40_a_dir
306c4762a1bSJed Brown output_file: output/ex40.out
307c4762a1bSJed Brown
308c4762a1bSJed Brown test:
309c4762a1bSJed Brown suffix: b
310c4762a1bSJed Brown args: -ts_type arkimex -snes_stol 1e-4 -ts_trajectory_dirname ex40_b_dir
311c4762a1bSJed Brown output_file: output/ex40.out
312c4762a1bSJed Brown
313c4762a1bSJed Brown test:
314c4762a1bSJed Brown suffix: c
315ca4445c7SIlya 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
316ca4445c7SIlya Fursov output_file: output/ex40.out
317ca4445c7SIlya Fursov
318ca4445c7SIlya Fursov test:
319ca4445c7SIlya Fursov suffix: cr
320ca4445c7SIlya 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
321ca4445c7SIlya Fursov output_file: output/ex40.out
322ca4445c7SIlya Fursov
323ca4445c7SIlya Fursov test:
324ca4445c7SIlya Fursov suffix: crmf
325ca4445c7SIlya 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
326c4762a1bSJed Brown output_file: output/ex40.out
327c4762a1bSJed Brown
328c4762a1bSJed Brown test:
329c4762a1bSJed Brown suffix: d
330c4762a1bSJed Brown args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_d_dir
331c4762a1bSJed Brown output_file: output/ex40.out
332c4762a1bSJed Brown
333c4762a1bSJed Brown test:
334c4762a1bSJed Brown suffix: e
335c4762a1bSJed Brown args: -ts_type bdf -ts_adapt_dt_max 0.025 -ts_max_steps 1500 -ts_trajectory_dirname ex40_e_dir
336c4762a1bSJed Brown output_file: output/ex40.out
337c4762a1bSJed Brown
338c4762a1bSJed Brown test:
339c4762a1bSJed Brown suffix: f
340c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 3bs -ts_trajectory_dirname ex40_f_dir
341c4762a1bSJed Brown output_file: output/ex40.out
342c4762a1bSJed Brown
343c4762a1bSJed Brown test:
344c4762a1bSJed Brown suffix: g
345c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 5bs -ts_trajectory_dirname ex40_g_dir
346c4762a1bSJed Brown output_file: output/ex40.out
347c4762a1bSJed Brown
348c4762a1bSJed Brown test:
349c4762a1bSJed Brown suffix: h
350c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 6vr -ts_trajectory_dirname ex40_h_dir
351c4762a1bSJed Brown output_file: output/ex40.out
352c4762a1bSJed Brown
353c4762a1bSJed Brown test:
354c4762a1bSJed Brown suffix: i
355c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 7vr -ts_trajectory_dirname ex40_i_dir
356c4762a1bSJed Brown output_file: output/ex40.out
357c4762a1bSJed Brown
358c4762a1bSJed Brown test:
359c4762a1bSJed Brown suffix: j
360c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 8vr -ts_trajectory_dirname ex40_j_dir
3615385558fSLisandro Dalcin output_file: output/ex40.out
3625385558fSLisandro Dalcin
3635385558fSLisandro Dalcin test:
3645385558fSLisandro Dalcin suffix: k
3655385558fSLisandro Dalcin args: -ts_type theta -ts_adapt_type dsp -ts_trajectory_dirname ex40_k_dir
3665385558fSLisandro Dalcin output_file: output/ex40.out
3675385558fSLisandro Dalcin
3685385558fSLisandro Dalcin test:
3695385558fSLisandro Dalcin suffix: l
3705385558fSLisandro Dalcin args: -rhs-form -ts_type rk -ts_rk_type 2a -ts_trajectory_dirname ex40_l_dir
3715385558fSLisandro Dalcin args: -ts_adapt_type dsp -ts_adapt_always_accept {{false true}} -ts_adapt_dt_min 0.01
372c4762a1bSJed Brown output_file: output/ex40.out
373c4762a1bSJed Brown
37408c4bdbbSLisandro Dalcin test:
37508c4bdbbSLisandro Dalcin suffix: m
37608c4bdbbSLisandro Dalcin args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -test_adapthistory false
37708c4bdbbSLisandro Dalcin args: -ts_max_time 10 -ts_exact_final_time {{STEPOVER MATCHSTEP INTERPOLATE}}
37808c4bdbbSLisandro Dalcin
37908c4bdbbSLisandro Dalcin test:
38008c4bdbbSLisandro Dalcin requires: !single
38108c4bdbbSLisandro Dalcin suffix: n
38208c4bdbbSLisandro Dalcin args: -test_adapthistory false
38308c4bdbbSLisandro Dalcin args: -ts_type alpha -ts_alpha_radius 1.0 -ts_view
384188af4bfSBarry Smith args: -ts_time_step 0.25 -ts_adapt_type basic -ts_adapt_wnormtype INFINITY -ts_adapt_monitor
385188af4bfSBarry Smith args: -ts_max_steps 1 -ts_max_step_rejections {{0 1 2}separate_output} -ts_error_if_step_fails false
38608c4bdbbSLisandro Dalcin
387e7685601SHong Zhang test:
388e7685601SHong Zhang requires: !single
389e7685601SHong Zhang suffix: o
390e7685601SHong Zhang args: -rhs-form -ts_type rk -ts_rk_type 2b -ts_trajectory_dirname ex40_o_dir
391e7685601SHong Zhang output_file: output/ex40.out
392c4762a1bSJed Brown TEST*/
393