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