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 8fa9584fbSIlya Fursov There is one event set in this example, which 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 10fa9584fbSIlya Fursov a restitution coefficient Cr. 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 PetscReal Cd; /* drag coefficient */ 18c4762a1bSJed Brown PetscReal Cr; /* restitution coefficient */ 19c4762a1bSJed Brown PetscInt bounces; 20c4762a1bSJed Brown PetscInt maxbounces; 21c4762a1bSJed Brown } AppCtx; 22c4762a1bSJed Brown 23fa9584fbSIlya Fursov static PetscErrorCode Event(TS ts, PetscReal t, Vec U, PetscReal *fvalue, void *ctx) 24d71ae5a4SJacob Faibussowitsch { 25c4762a1bSJed Brown Vec V; 26fa9584fbSIlya Fursov const PetscScalar *u; 27c4762a1bSJed Brown 287510d9b0SBarry Smith PetscFunctionBeginUser; 29c4762a1bSJed Brown /* Event for ball height */ 309566063dSJacob Faibussowitsch PetscCall(TS2GetSolution(ts, &U, &V)); 319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 32fa9584fbSIlya Fursov fvalue[0] = PetscRealPart(u[0]); 339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown 37d71ae5a4SJacob Faibussowitsch static PetscErrorCode PostEvent(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) 38d71ae5a4SJacob Faibussowitsch { 39c4762a1bSJed Brown AppCtx *app = (AppCtx *)ctx; 40c4762a1bSJed Brown Vec V; 41c4762a1bSJed Brown PetscScalar *u, *v; 42c4762a1bSJed Brown PetscMPIInt rank; 43fa9584fbSIlya Fursov PetscBool inflag = PETSC_FALSE, outflag; 44c4762a1bSJed Brown 457510d9b0SBarry Smith PetscFunctionBeginUser; 469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 47fa9584fbSIlya Fursov if (nevents > 0) { 489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor [%d]: Ball hit the ground at t = %5.2f seconds\n", rank, (double)t)); 49fa9584fbSIlya Fursov /* Set new initial conditions with attenuation Cr */ 509566063dSJacob Faibussowitsch PetscCall(TS2GetSolution(ts, &U, &V)); 519566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 529566063dSJacob Faibussowitsch PetscCall(VecGetArray(V, &v)); 539371c9d4SSatish Balay u[0] = 0.0; 549371c9d4SSatish Balay v[0] = -app->Cr * v[0]; 559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(V, &v)); 57c4762a1bSJed Brown app->bounces++; 58c4762a1bSJed Brown } 59fa9584fbSIlya Fursov if (app->bounces >= app->maxbounces) { // 'app->bounces' may be different on different processes 60fa9584fbSIlya Fursov PetscCall(PetscPrintf(PETSC_COMM_SELF, "Processor [%d]: Ball bounced %" PetscInt_FMT " times\n", rank, app->bounces)); 61fa9584fbSIlya Fursov inflag = PETSC_TRUE; // current process requested to terminate 62fa9584fbSIlya Fursov } 63*5440e5dcSBarry Smith PetscCallMPI(MPIU_Allreduce(&inflag, &outflag, 1, MPI_C_BOOL, MPI_LOR, PetscObjectComm((PetscObject)ts))); 64fa9584fbSIlya Fursov if (outflag) PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); // request TS to terminate, sync on all ranks 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(P, MAT_FINAL_ASSEMBLY)); 1089566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY)); 109fa9584fbSIlya Fursov if (J != P) { 110fa9584fbSIlya Fursov PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 111fa9584fbSIlya Fursov PetscCall(MatAssemblyEnd(J, 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; 125fa9584fbSIlya Fursov PetscInt direction[1]; 126fa9584fbSIlya Fursov PetscBool terminate[1]; 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; 156fa9584fbSIlya Fursov PetscCall(TSSetEventHandler(ts, 1, direction, terminate, Event, PostEvent, &app)); // each process has one event-function defined 157c4762a1bSJed Brown 1589566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, 1, 1, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, 0, NULL, &J)); 1599566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 1609566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 1619566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J, NULL, &F)); 1629566063dSJacob Faibussowitsch PetscCall(TSSetI2Function(ts, F, I2Function, &app)); 1639566063dSJacob Faibussowitsch PetscCall(TSSetI2Jacobian(ts, J, J, I2Jacobian, &app)); 1649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F)); 1659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(TSGetI2Jacobian(ts, &J, NULL, NULL, NULL)); 1689566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J, &U, NULL)); 1699566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J, &V, NULL)); 1709566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 1719566063dSJacob Faibussowitsch PetscCall(VecGetArray(V, &v)); 1729371c9d4SSatish Balay u[0] = 5.0 * rank; 1739371c9d4SSatish Balay v[0] = 20.0; 1749566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 1759566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(V, &v)); 176c4762a1bSJed Brown 1779566063dSJacob Faibussowitsch PetscCall(TS2SetSolution(ts, U, V)); 1789566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1799566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, NULL)); 180c4762a1bSJed Brown 1819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 1829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V)); 1839566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 184c4762a1bSJed Brown 1859566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 186b122ec5aSJacob Faibussowitsch return 0; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown /*TEST 190c4762a1bSJed Brown 191c4762a1bSJed Brown test: 192c4762a1bSJed Brown suffix: a 193ca4445c7SIlya Fursov args: -ts_alpha_radius {{1.0 0.5}} -ts_max_time 50 194c4762a1bSJed Brown output_file: output/ex44.out 195c4762a1bSJed Brown 196c4762a1bSJed Brown test: 197c4762a1bSJed Brown suffix: b 198ca4445c7SIlya Fursov args: -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic -ts_max_time 50 199ca4445c7SIlya Fursov output_file: output/ex44.out 200ca4445c7SIlya Fursov 201ca4445c7SIlya Fursov test: 202ca4445c7SIlya Fursov suffix: bmf 203ca4445c7SIlya Fursov args: -snes_mf_operator -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic -ts_max_time 50 204c4762a1bSJed Brown output_file: output/ex44.out 205c4762a1bSJed Brown 206c4762a1bSJed Brown test: 207c4762a1bSJed Brown suffix: 2 208c4762a1bSJed Brown nsize: 2 209ca4445c7SIlya Fursov args: -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic -ts_max_time 50 210c4762a1bSJed Brown output_file: output/ex44_2.out 211c4762a1bSJed Brown filter: sort -b 212c4762a1bSJed Brown filter_output: sort -b 213c4762a1bSJed Brown 214e6b32627SLisandro Dalcin test: 215e6b32627SLisandro Dalcin requires: !single 216e6b32627SLisandro Dalcin args: -ts_dt 0.25 -ts_adapt_type basic -ts_adapt_wnormtype INFINITY -ts_adapt_monitor 217e6b32627SLisandro Dalcin args: -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false 218e6b32627SLisandro Dalcin 219c4762a1bSJed Brown TEST*/ 220