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 22c4762a1bSJed Brown static PetscErrorCode Event(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx) 23c4762a1bSJed Brown { 24c4762a1bSJed Brown AppCtx *app = (AppCtx*)ctx; 25c4762a1bSJed Brown Vec V; 26c4762a1bSJed Brown const PetscScalar *u,*v; 27c4762a1bSJed Brown PetscErrorCode ierr; 28c4762a1bSJed Brown 29c4762a1bSJed Brown PetscFunctionBegin; 30c4762a1bSJed Brown /* Event for ball height */ 31c4762a1bSJed Brown ierr = TS2GetSolution(ts,&U,&V);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = VecGetArrayRead(V,&v);CHKERRQ(ierr); 34c4762a1bSJed Brown fvalue[0] = u[0]; 35c4762a1bSJed Brown /* Event for number of bounces */ 36c4762a1bSJed Brown fvalue[1] = app->maxbounces - app->bounces; 37c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = VecRestoreArrayRead(V,&v);CHKERRQ(ierr); 39c4762a1bSJed Brown PetscFunctionReturn(0); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42c4762a1bSJed Brown static PetscErrorCode PostEvent(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 43c4762a1bSJed Brown { 44c4762a1bSJed Brown AppCtx *app = (AppCtx*)ctx; 45c4762a1bSJed Brown Vec V; 46c4762a1bSJed Brown PetscScalar *u,*v; 47c4762a1bSJed Brown PetscMPIInt rank; 48c4762a1bSJed Brown PetscErrorCode ierr; 49c4762a1bSJed Brown 50c4762a1bSJed Brown PetscFunctionBegin; 51c4762a1bSJed Brown if (!nevents) PetscFunctionReturn(0); 52c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 53c4762a1bSJed Brown if (event_list[0] == 0) { 54c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Processor [%d]: Ball hit the ground at t = %5.2f seconds\n",rank,(double)t);CHKERRQ(ierr); 55c4762a1bSJed Brown /* Set new initial conditions with .9 attenuation */ 56c4762a1bSJed Brown ierr = TS2GetSolution(ts,&U,&V);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = VecGetArray(V,&v);CHKERRQ(ierr); 59c4762a1bSJed Brown u[0] = 0.0; v[0] = -app->Cr*v[0]; 60c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = VecRestoreArray(V,&v);CHKERRQ(ierr); 62c4762a1bSJed Brown app->bounces++; 63c4762a1bSJed Brown } else if (event_list[0] == 1) { 64c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Processor [%d]: Ball bounced %D times\n",rank,app->bounces);CHKERRQ(ierr); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown PetscFunctionReturn(0); 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown static PetscErrorCode I2Function(TS ts,PetscReal t,Vec U,Vec V,Vec A,Vec F,void *ctx) 70c4762a1bSJed Brown { 71c4762a1bSJed Brown AppCtx *app = (AppCtx*)ctx; 72c4762a1bSJed Brown const PetscScalar *u,*v,*a; 73c4762a1bSJed Brown PetscScalar Res,*f; 74c4762a1bSJed Brown PetscErrorCode ierr; 75c4762a1bSJed Brown 76c4762a1bSJed Brown PetscFunctionBegin; 77c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 78c4762a1bSJed Brown ierr = VecGetArrayRead(V,&v);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = VecGetArrayRead(A,&a);CHKERRQ(ierr); 80c4762a1bSJed Brown Res = a[0] + 9.8 + 0.5 * app->Cd * v[0]*v[0] * PetscSignReal(PetscRealPart(v[0])); 81c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = VecRestoreArrayRead(V,&v);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = VecRestoreArrayRead(A,&a);CHKERRQ(ierr); 84c4762a1bSJed Brown 85c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 86c4762a1bSJed Brown f[0] = Res; 87c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 88c4762a1bSJed Brown PetscFunctionReturn(0); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown static PetscErrorCode I2Jacobian(TS ts,PetscReal t,Vec U,Vec V,Vec A,PetscReal shiftV,PetscReal shiftA,Mat J,Mat P,void *ctx) 92c4762a1bSJed Brown { 93c4762a1bSJed Brown AppCtx *app = (AppCtx*)ctx; 94c4762a1bSJed Brown const PetscScalar *u,*v,*a; 95c4762a1bSJed Brown PetscInt i; 96c4762a1bSJed Brown PetscScalar Jac; 97c4762a1bSJed Brown PetscErrorCode ierr; 98c4762a1bSJed Brown 99c4762a1bSJed Brown PetscFunctionBegin; 100c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = VecGetArrayRead(V,&v);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = VecGetArrayRead(A,&a);CHKERRQ(ierr); 103c4762a1bSJed Brown Jac = shiftA + shiftV * app->Cd * v[0]; 104c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = VecRestoreArrayRead(V,&v);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = VecRestoreArrayRead(A,&a);CHKERRQ(ierr); 107c4762a1bSJed Brown 108c4762a1bSJed Brown ierr = MatGetOwnershipRange(P,&i,NULL);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = MatSetValue(P,i,i,Jac,INSERT_VALUES);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 112c4762a1bSJed Brown if (J != P) { 113c4762a1bSJed Brown ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown PetscFunctionReturn(0); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown int main(int argc,char **argv) 120c4762a1bSJed Brown { 121c4762a1bSJed Brown TS ts; /* ODE integrator */ 122c4762a1bSJed Brown Vec U,V; /* solution will be stored here */ 123c4762a1bSJed Brown Vec F; /* residual vector */ 124c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 125c4762a1bSJed Brown PetscMPIInt rank; 126c4762a1bSJed Brown PetscScalar *u,*v; 127c4762a1bSJed Brown AppCtx app; 128c4762a1bSJed Brown PetscInt direction[2]; 129c4762a1bSJed Brown PetscBool terminate[2]; 130c4762a1bSJed Brown TSAdapt adapt; 131c4762a1bSJed Brown PetscErrorCode ierr; 132c4762a1bSJed Brown 133c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 134c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 135c4762a1bSJed Brown 136c4762a1bSJed Brown app.Cd = 0.0; 137c4762a1bSJed Brown app.Cr = 0.9; 138c4762a1bSJed Brown app.bounces = 0; 139c4762a1bSJed Brown app.maxbounces = 10; 140c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex44 options","");CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = PetscOptionsReal("-Cd","Drag coefficient","",app.Cd,&app.Cd,NULL);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = PetscOptionsReal("-Cr","Restitution coefficient","",app.Cr,&app.Cr,NULL);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = PetscOptionsInt("-maxbounces","Maximum number of bounces","",app.maxbounces,&app.maxbounces,NULL);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 145c4762a1bSJed Brown 146c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 147c4762a1bSJed Brown /*ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);*/ 148c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = TSSetType(ts,TSALPHA2);CHKERRQ(ierr); 150c4762a1bSJed Brown 151c4762a1bSJed Brown ierr = TSSetMaxTime(ts,PETSC_INFINITY);CHKERRQ(ierr); 152c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.1);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = TSAdaptSetStepLimits(adapt,0.0,0.5);CHKERRQ(ierr); 156c4762a1bSJed Brown 157c4762a1bSJed Brown direction[0] = -1; terminate[0] = PETSC_FALSE; 158c4762a1bSJed Brown direction[1] = -1; terminate[1] = PETSC_TRUE; 159c4762a1bSJed Brown ierr = TSSetEventHandler(ts,2,direction,terminate,Event,PostEvent,&app);CHKERRQ(ierr); 160c4762a1bSJed Brown 161c4762a1bSJed Brown ierr = MatCreateAIJ(PETSC_COMM_WORLD,1,1,PETSC_DECIDE,PETSC_DECIDE,1,NULL,0,NULL,&J);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = MatSetFromOptions(J);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = MatCreateVecs(J,NULL,&F);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = TSSetI2Function(ts,F,I2Function,&app);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = TSSetI2Jacobian(ts,J,J,I2Jacobian,&app);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = VecDestroy(&F);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 169c4762a1bSJed Brown 170c4762a1bSJed Brown ierr = TSGetI2Jacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = MatCreateVecs(J,&U,NULL);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = MatCreateVecs(J,&V,NULL);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 174c4762a1bSJed Brown ierr = VecGetArray(V,&v);CHKERRQ(ierr); 175c4762a1bSJed Brown u[0] = 5.0*rank; v[0] = 20.0; 176c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = VecRestoreArray(V,&v);CHKERRQ(ierr); 178c4762a1bSJed Brown 179c4762a1bSJed Brown ierr = TS2SetSolution(ts,U,V);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = TSSolve(ts,NULL);CHKERRQ(ierr); 182c4762a1bSJed Brown 183c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = VecDestroy(&V);CHKERRQ(ierr); 185c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 186c4762a1bSJed Brown 187c4762a1bSJed Brown ierr = PetscFinalize(); 188c4762a1bSJed Brown return ierr; 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown /*TEST 192c4762a1bSJed Brown 193c4762a1bSJed Brown test: 194c4762a1bSJed Brown suffix: a 195*e6b32627SLisandro Dalcin args: -ts_alpha_radius {{1.0 0.5}} 196c4762a1bSJed Brown output_file: output/ex44.out 197c4762a1bSJed Brown 198c4762a1bSJed Brown test: 199c4762a1bSJed Brown suffix: b 200c4762a1bSJed Brown args: -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic 201c4762a1bSJed Brown output_file: output/ex44.out 202c4762a1bSJed Brown 203c4762a1bSJed Brown test: 204c4762a1bSJed Brown suffix: 2 205c4762a1bSJed Brown nsize: 2 206c4762a1bSJed Brown args: -ts_rtol 0 -ts_atol 1e-1 -ts_adapt_type basic 207c4762a1bSJed Brown output_file: output/ex44_2.out 208c4762a1bSJed Brown filter: sort -b 209c4762a1bSJed Brown filter_output: sort -b 210c4762a1bSJed Brown 211*e6b32627SLisandro Dalcin test: 212*e6b32627SLisandro Dalcin requires: !single 213*e6b32627SLisandro Dalcin args: -ts_dt 0.25 -ts_adapt_type basic -ts_adapt_wnormtype INFINITY -ts_adapt_monitor 214*e6b32627SLisandro Dalcin args: -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false 215*e6b32627SLisandro Dalcin 216c4762a1bSJed Brown TEST*/ 217