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 22*9371c9d4SSatish Balay static PetscErrorCode Event(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx) { 23c4762a1bSJed Brown AppCtx *app = (AppCtx *)ctx; 24c4762a1bSJed Brown Vec V; 25c4762a1bSJed Brown const PetscScalar *u, *v; 26c4762a1bSJed Brown 277510d9b0SBarry Smith PetscFunctionBeginUser; 28c4762a1bSJed Brown /* Event for ball height */ 299566063dSJacob Faibussowitsch PetscCall(TS2GetSolution(ts, &U, &V)); 309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(V, &v)); 32c4762a1bSJed Brown fvalue[0] = u[0]; 33c4762a1bSJed Brown /* Event for number of bounces */ 34c4762a1bSJed Brown fvalue[1] = app->maxbounces - app->bounces; 359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(V, &v)); 37c4762a1bSJed Brown PetscFunctionReturn(0); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40*9371c9d4SSatish Balay static PetscErrorCode PostEvent(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) { 41c4762a1bSJed Brown AppCtx *app = (AppCtx *)ctx; 42c4762a1bSJed Brown Vec V; 43c4762a1bSJed Brown PetscScalar *u, *v; 44c4762a1bSJed Brown PetscMPIInt rank; 45c4762a1bSJed Brown 467510d9b0SBarry Smith PetscFunctionBeginUser; 47c4762a1bSJed Brown if (!nevents) PetscFunctionReturn(0); 489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 49c4762a1bSJed Brown if (event_list[0] == 0) { 509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor [%d]: Ball hit the ground at t = %5.2f seconds\n", rank, (double)t)); 51c4762a1bSJed Brown /* Set new initial conditions with .9 attenuation */ 529566063dSJacob Faibussowitsch PetscCall(TS2GetSolution(ts, &U, &V)); 539566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 549566063dSJacob Faibussowitsch PetscCall(VecGetArray(V, &v)); 55*9371c9d4SSatish Balay u[0] = 0.0; 56*9371c9d4SSatish Balay v[0] = -app->Cr * v[0]; 579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 589566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(V, &v)); 59c4762a1bSJed Brown app->bounces++; 60c4762a1bSJed Brown } else if (event_list[0] == 1) { 6163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor [%d]: Ball bounced %" PetscInt_FMT " times\n", rank, app->bounces)); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown PetscFunctionReturn(0); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66*9371c9d4SSatish Balay static PetscErrorCode I2Function(TS ts, PetscReal t, Vec U, Vec V, Vec A, Vec F, void *ctx) { 67c4762a1bSJed Brown AppCtx *app = (AppCtx *)ctx; 68c4762a1bSJed Brown const PetscScalar *u, *v, *a; 69c4762a1bSJed Brown PetscScalar Res, *f; 70c4762a1bSJed Brown 717510d9b0SBarry Smith PetscFunctionBeginUser; 729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(V, &v)); 749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(A, &a)); 75c4762a1bSJed Brown Res = a[0] + 9.8 + 0.5 * app->Cd * v[0] * v[0] * PetscSignReal(PetscRealPart(v[0])); 769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(V, &v)); 789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(A, &a)); 79c4762a1bSJed Brown 809566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 81c4762a1bSJed Brown f[0] = Res; 829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 83c4762a1bSJed Brown PetscFunctionReturn(0); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86*9371c9d4SSatish Balay static PetscErrorCode I2Jacobian(TS ts, PetscReal t, Vec U, Vec V, Vec A, PetscReal shiftV, PetscReal shiftA, Mat J, Mat P, void *ctx) { 87c4762a1bSJed Brown AppCtx *app = (AppCtx *)ctx; 88c4762a1bSJed Brown const PetscScalar *u, *v, *a; 89c4762a1bSJed Brown PetscInt i; 90c4762a1bSJed Brown PetscScalar Jac; 91c4762a1bSJed Brown 927510d9b0SBarry Smith PetscFunctionBeginUser; 939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(V, &v)); 959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(A, &a)); 96c4762a1bSJed Brown Jac = shiftA + shiftV * app->Cd * v[0]; 979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(V, &v)); 999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(A, &a)); 100c4762a1bSJed Brown 1019566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(P, &i, NULL)); 1029566063dSJacob Faibussowitsch PetscCall(MatSetValue(P, i, i, Jac, INSERT_VALUES)); 1039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 105c4762a1bSJed Brown if (J != P) { 1069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY)); 1079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY)); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown PetscFunctionReturn(0); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112*9371c9d4SSatish Balay int main(int argc, char **argv) { 113c4762a1bSJed Brown TS ts; /* ODE integrator */ 114c4762a1bSJed Brown Vec U, V; /* solution will be stored here */ 115c4762a1bSJed Brown Vec F; /* residual vector */ 116c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 117c4762a1bSJed Brown PetscMPIInt rank; 118c4762a1bSJed Brown PetscScalar *u, *v; 119c4762a1bSJed Brown AppCtx app; 120c4762a1bSJed Brown PetscInt direction[2]; 121c4762a1bSJed Brown PetscBool terminate[2]; 122c4762a1bSJed Brown TSAdapt adapt; 123c4762a1bSJed Brown 124327415f7SBarry Smith PetscFunctionBeginUser; 1259566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown app.Cd = 0.0; 129c4762a1bSJed Brown app.Cr = 0.9; 130c4762a1bSJed Brown app.bounces = 0; 131c4762a1bSJed Brown app.maxbounces = 10; 132d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex44 options", ""); 1339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-Cd", "Drag coefficient", "", app.Cd, &app.Cd, NULL)); 1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-Cr", "Restitution coefficient", "", app.Cr, &app.Cr, NULL)); 1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-maxbounces", "Maximum number of bounces", "", app.maxbounces, &app.maxbounces, NULL)); 136d0609cedSBarry Smith PetscOptionsEnd(); 137c4762a1bSJed Brown 1389566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1399566063dSJacob Faibussowitsch /*PetscCall(TSSetSaveTrajectory(ts));*/ 1409566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1419566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSALPHA2)); 142c4762a1bSJed Brown 1439566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, PETSC_INFINITY)); 1449566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.1)); 1459566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1469566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 1479566063dSJacob Faibussowitsch PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5)); 148c4762a1bSJed Brown 149*9371c9d4SSatish Balay direction[0] = -1; 150*9371c9d4SSatish Balay terminate[0] = PETSC_FALSE; 151*9371c9d4SSatish Balay direction[1] = -1; 152*9371c9d4SSatish Balay terminate[1] = PETSC_TRUE; 1539566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts, 2, direction, terminate, Event, PostEvent, &app)); 154c4762a1bSJed Brown 1559566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, 1, 1, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, 0, NULL, &J)); 1569566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 1579566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 1589566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J, NULL, &F)); 1599566063dSJacob Faibussowitsch PetscCall(TSSetI2Function(ts, F, I2Function, &app)); 1609566063dSJacob Faibussowitsch PetscCall(TSSetI2Jacobian(ts, J, J, I2Jacobian, &app)); 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F)); 1629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 163c4762a1bSJed Brown 1649566063dSJacob Faibussowitsch PetscCall(TSGetI2Jacobian(ts, &J, NULL, NULL, NULL)); 1659566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J, &U, NULL)); 1669566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J, &V, NULL)); 1679566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 1689566063dSJacob Faibussowitsch PetscCall(VecGetArray(V, &v)); 169*9371c9d4SSatish Balay u[0] = 5.0 * rank; 170*9371c9d4SSatish Balay v[0] = 20.0; 1719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 1729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(V, &v)); 173c4762a1bSJed Brown 1749566063dSJacob Faibussowitsch PetscCall(TS2SetSolution(ts, U, V)); 1759566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1769566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, NULL)); 177c4762a1bSJed Brown 1789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 1799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V)); 1809566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 181c4762a1bSJed Brown 1829566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 183b122ec5aSJacob Faibussowitsch return 0; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown /*TEST 187c4762a1bSJed Brown 188c4762a1bSJed Brown test: 189c4762a1bSJed Brown suffix: a 190e6b32627SLisandro Dalcin args: -ts_alpha_radius {{1.0 0.5}} 191c4762a1bSJed Brown output_file: output/ex44.out 192c4762a1bSJed Brown 193c4762a1bSJed Brown test: 194c4762a1bSJed Brown suffix: b 195c4762a1bSJed Brown args: -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic 196c4762a1bSJed Brown output_file: output/ex44.out 197c4762a1bSJed Brown 198c4762a1bSJed Brown test: 199c4762a1bSJed Brown suffix: 2 200c4762a1bSJed Brown nsize: 2 201c4762a1bSJed Brown args: -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic 202c4762a1bSJed Brown output_file: output/ex44_2.out 203c4762a1bSJed Brown filter: sort -b 204c4762a1bSJed Brown filter_output: sort -b 205c4762a1bSJed Brown 206e6b32627SLisandro Dalcin test: 207e6b32627SLisandro Dalcin requires: !single 208e6b32627SLisandro Dalcin args: -ts_dt 0.25 -ts_adapt_type basic -ts_adapt_wnormtype INFINITY -ts_adapt_monitor 209e6b32627SLisandro Dalcin args: -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false 210e6b32627SLisandro Dalcin 211c4762a1bSJed Brown TEST*/ 212