xref: /petsc/src/ts/tutorials/ex44.c (revision ca4445c7a2f5ca546b532f08b853c371604af17c)
1c4762a1bSJed Brown static char help[] = "Parallel bouncing ball example formulated as a second-order system to test TS event feature.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   The dynamics of the bouncing ball with drag coefficient Cd is described by the ODE
5c4762a1bSJed Brown 
6c4762a1bSJed Brown       u_tt = -9.8 - 1/2 Cd (u_t)^2 sign(u_t)
7c4762a1bSJed Brown 
8c4762a1bSJed Brown   There are two events set in this example. The first one checks for the ball hitting the
9c4762a1bSJed Brown   ground (u = 0). Every time the ball hits the ground, its velocity u_t is attenuated by
10c4762a1bSJed Brown   a restitution coefficient Cr. 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   PetscReal Cd; /* drag coefficient */
17c4762a1bSJed Brown   PetscReal Cr; /* restitution coefficient */
18c4762a1bSJed Brown   PetscInt  bounces;
19c4762a1bSJed Brown   PetscInt  maxbounces;
20c4762a1bSJed Brown } AppCtx;
21c4762a1bSJed Brown 
22d71ae5a4SJacob Faibussowitsch static PetscErrorCode Event(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx)
23d71ae5a4SJacob Faibussowitsch {
24c4762a1bSJed Brown   AppCtx            *app = (AppCtx *)ctx;
25c4762a1bSJed Brown   Vec                V;
26c4762a1bSJed Brown   const PetscScalar *u, *v;
27c4762a1bSJed Brown 
287510d9b0SBarry Smith   PetscFunctionBeginUser;
29c4762a1bSJed Brown   /* Event for ball height */
309566063dSJacob Faibussowitsch   PetscCall(TS2GetSolution(ts, &U, &V));
319566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V, &v));
33c4762a1bSJed Brown   fvalue[0] = u[0];
34c4762a1bSJed Brown   /* Event for number of bounces */
35c4762a1bSJed Brown   fvalue[1] = app->maxbounces - app->bounces;
369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V, &v));
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39c4762a1bSJed Brown }
40c4762a1bSJed Brown 
41d71ae5a4SJacob Faibussowitsch static PetscErrorCode PostEvent(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx)
42d71ae5a4SJacob Faibussowitsch {
43c4762a1bSJed Brown   AppCtx      *app = (AppCtx *)ctx;
44c4762a1bSJed Brown   Vec          V;
45c4762a1bSJed Brown   PetscScalar *u, *v;
46c4762a1bSJed Brown   PetscMPIInt  rank;
47c4762a1bSJed Brown 
487510d9b0SBarry Smith   PetscFunctionBeginUser;
493ba16761SJacob Faibussowitsch   if (!nevents) PetscFunctionReturn(PETSC_SUCCESS);
509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
51c4762a1bSJed Brown   if (event_list[0] == 0) {
529566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor [%d]: Ball hit the ground at t = %5.2f seconds\n", rank, (double)t));
53c4762a1bSJed Brown     /* Set new initial conditions with .9 attenuation */
549566063dSJacob Faibussowitsch     PetscCall(TS2GetSolution(ts, &U, &V));
559566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
569566063dSJacob Faibussowitsch     PetscCall(VecGetArray(V, &v));
579371c9d4SSatish Balay     u[0] = 0.0;
589371c9d4SSatish Balay     v[0] = -app->Cr * v[0];
599566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
609566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(V, &v));
61c4762a1bSJed Brown     app->bounces++;
62c4762a1bSJed Brown   } else if (event_list[0] == 1) {
6363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor [%d]: Ball bounced %" PetscInt_FMT " times\n", rank, app->bounces));
64c4762a1bSJed Brown   }
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
68d71ae5a4SJacob Faibussowitsch static PetscErrorCode I2Function(TS ts, PetscReal t, Vec U, Vec V, Vec A, Vec F, void *ctx)
69d71ae5a4SJacob Faibussowitsch {
70c4762a1bSJed Brown   AppCtx            *app = (AppCtx *)ctx;
71c4762a1bSJed Brown   const PetscScalar *u, *v, *a;
72c4762a1bSJed Brown   PetscScalar        Res, *f;
73c4762a1bSJed Brown 
747510d9b0SBarry Smith   PetscFunctionBeginUser;
759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V, &v));
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(A, &a));
78c4762a1bSJed Brown   Res = a[0] + 9.8 + 0.5 * app->Cd * v[0] * v[0] * PetscSignReal(PetscRealPart(v[0]));
799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V, &v));
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(A, &a));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
84c4762a1bSJed Brown   f[0] = Res;
859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89d71ae5a4SJacob Faibussowitsch static PetscErrorCode I2Jacobian(TS ts, PetscReal t, Vec U, Vec V, Vec A, PetscReal shiftV, PetscReal shiftA, Mat J, Mat P, void *ctx)
90d71ae5a4SJacob Faibussowitsch {
91c4762a1bSJed Brown   AppCtx            *app = (AppCtx *)ctx;
92c4762a1bSJed Brown   const PetscScalar *u, *v, *a;
93c4762a1bSJed Brown   PetscInt           i;
94c4762a1bSJed Brown   PetscScalar        Jac;
95c4762a1bSJed Brown 
967510d9b0SBarry Smith   PetscFunctionBeginUser;
979566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
989566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V, &v));
999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(A, &a));
100c4762a1bSJed Brown   Jac = shiftA + shiftV * app->Cd * v[0];
1019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V, &v));
1039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(A, &a));
104c4762a1bSJed Brown 
1059566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(P, &i, NULL));
1069566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, i, i, Jac, INSERT_VALUES));
1079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
109c4762a1bSJed Brown   if (J != P) {
1109566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
1119566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
112c4762a1bSJed Brown   }
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
116d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
117d71ae5a4SJacob Faibussowitsch {
118c4762a1bSJed Brown   TS           ts;   /* ODE integrator */
119c4762a1bSJed Brown   Vec          U, V; /* solution will be stored here */
120c4762a1bSJed Brown   Vec          F;    /* residual vector */
121c4762a1bSJed Brown   Mat          J;    /* Jacobian matrix */
122c4762a1bSJed Brown   PetscMPIInt  rank;
123c4762a1bSJed Brown   PetscScalar *u, *v;
124c4762a1bSJed Brown   AppCtx       app;
125c4762a1bSJed Brown   PetscInt     direction[2];
126c4762a1bSJed Brown   PetscBool    terminate[2];
127c4762a1bSJed Brown   TSAdapt      adapt;
128c4762a1bSJed Brown 
129327415f7SBarry Smith   PetscFunctionBeginUser;
1309566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   app.Cd         = 0.0;
134c4762a1bSJed Brown   app.Cr         = 0.9;
135c4762a1bSJed Brown   app.bounces    = 0;
136c4762a1bSJed Brown   app.maxbounces = 10;
137d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex44 options", "");
1389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-Cd", "Drag coefficient", "", app.Cd, &app.Cd, NULL));
1399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-Cr", "Restitution coefficient", "", app.Cr, &app.Cr, NULL));
1409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-maxbounces", "Maximum number of bounces", "", app.maxbounces, &app.maxbounces, NULL));
141d0609cedSBarry Smith   PetscOptionsEnd();
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1449566063dSJacob Faibussowitsch   /*PetscCall(TSSetSaveTrajectory(ts));*/
1459566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1469566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSALPHA2));
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, PETSC_INFINITY));
1499566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.1));
1509566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1519566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
1529566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5));
153c4762a1bSJed Brown 
1549371c9d4SSatish Balay   direction[0] = -1;
1559371c9d4SSatish Balay   terminate[0] = PETSC_FALSE;
1569371c9d4SSatish Balay   direction[1] = -1;
1579371c9d4SSatish Balay   terminate[1] = PETSC_TRUE;
1589566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts, 2, direction, terminate, Event, PostEvent, &app));
159c4762a1bSJed Brown 
1609566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, 1, 1, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, 0, NULL, &J));
1619566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
1629566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
1639566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(J, NULL, &F));
1649566063dSJacob Faibussowitsch   PetscCall(TSSetI2Function(ts, F, I2Function, &app));
1659566063dSJacob Faibussowitsch   PetscCall(TSSetI2Jacobian(ts, J, J, I2Jacobian, &app));
1669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&F));
1679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
168c4762a1bSJed Brown 
1699566063dSJacob Faibussowitsch   PetscCall(TSGetI2Jacobian(ts, &J, NULL, NULL, NULL));
1709566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(J, &U, NULL));
1719566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(J, &V, NULL));
1729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
1739566063dSJacob Faibussowitsch   PetscCall(VecGetArray(V, &v));
1749371c9d4SSatish Balay   u[0] = 5.0 * rank;
1759371c9d4SSatish Balay   v[0] = 20.0;
1769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
1779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(V, &v));
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch   PetscCall(TS2SetSolution(ts, U, V));
1809566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1819566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, NULL));
182c4762a1bSJed Brown 
1839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&V));
1859566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
186c4762a1bSJed Brown 
1879566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
188b122ec5aSJacob Faibussowitsch   return 0;
189c4762a1bSJed Brown }
190c4762a1bSJed Brown 
191c4762a1bSJed Brown /*TEST
192c4762a1bSJed Brown 
193c4762a1bSJed Brown     test:
194c4762a1bSJed Brown       suffix: a
195*ca4445c7SIlya Fursov       args: -ts_alpha_radius {{1.0 0.5}} -ts_max_time 50
196c4762a1bSJed Brown       output_file: output/ex44.out
197c4762a1bSJed Brown 
198c4762a1bSJed Brown     test:
199c4762a1bSJed Brown       suffix: b
200*ca4445c7SIlya Fursov       args: -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic -ts_max_time 50
201*ca4445c7SIlya Fursov       output_file: output/ex44.out
202*ca4445c7SIlya Fursov 
203*ca4445c7SIlya Fursov     test:
204*ca4445c7SIlya Fursov       suffix: bmf
205*ca4445c7SIlya Fursov       args: -snes_mf_operator -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic -ts_max_time 50
206c4762a1bSJed Brown       output_file: output/ex44.out
207c4762a1bSJed Brown 
208c4762a1bSJed Brown     test:
209c4762a1bSJed Brown       suffix: 2
210c4762a1bSJed Brown       nsize: 2
211*ca4445c7SIlya Fursov       args: -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic -ts_max_time 50
212c4762a1bSJed Brown       output_file: output/ex44_2.out
213c4762a1bSJed Brown       filter: sort -b
214c4762a1bSJed Brown       filter_output: sort -b
215c4762a1bSJed Brown 
216e6b32627SLisandro Dalcin     test:
217e6b32627SLisandro Dalcin       requires: !single
218e6b32627SLisandro Dalcin       args: -ts_dt 0.25 -ts_adapt_type basic -ts_adapt_wnormtype INFINITY -ts_adapt_monitor
219e6b32627SLisandro Dalcin       args: -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false
220e6b32627SLisandro Dalcin 
221c4762a1bSJed Brown TEST*/
222