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