1*c4762a1bSJed Brown static char help[] = "Serial bouncing ball example to test TS event feature.\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown /* 4*c4762a1bSJed Brown The dynamics of the bouncing ball is described by the ODE 5*c4762a1bSJed Brown u1_t = u2 6*c4762a1bSJed Brown u2_t = -9.8 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 (u1 = 0). Every time the ball hits the ground, its velocity u2 is attenuated by 10*c4762a1bSJed Brown a factor of 0.9. 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 PetscInt maxbounces; 17*c4762a1bSJed Brown PetscInt nbounces; 18*c4762a1bSJed Brown } AppCtx; 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx) 21*c4762a1bSJed Brown { 22*c4762a1bSJed Brown AppCtx *app=(AppCtx*)ctx; 23*c4762a1bSJed Brown PetscErrorCode ierr; 24*c4762a1bSJed Brown const PetscScalar *u; 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown PetscFunctionBegin; 27*c4762a1bSJed Brown /* Event for ball height */ 28*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 29*c4762a1bSJed Brown fvalue[0] = u[0]; 30*c4762a1bSJed Brown /* Event for number of bounces */ 31*c4762a1bSJed Brown fvalue[1] = app->maxbounces - app->nbounces; 32*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 33*c4762a1bSJed Brown PetscFunctionReturn(0); 34*c4762a1bSJed Brown } 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 37*c4762a1bSJed Brown { 38*c4762a1bSJed Brown AppCtx *app=(AppCtx*)ctx; 39*c4762a1bSJed Brown PetscErrorCode ierr; 40*c4762a1bSJed Brown PetscScalar *u; 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown PetscFunctionBegin; 43*c4762a1bSJed Brown if (event_list[0] == 0) { 44*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Ball hit the ground at t = %5.2f seconds\n",(double)t);CHKERRQ(ierr); 45*c4762a1bSJed Brown /* Set new initial conditions with .9 attenuation */ 46*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 47*c4762a1bSJed Brown u[0] = 0.0; 48*c4762a1bSJed Brown u[1] = -0.9*u[1]; 49*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 50*c4762a1bSJed Brown } else if (event_list[0] == 1) { 51*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Ball bounced %D times\n",app->nbounces);CHKERRQ(ierr); 52*c4762a1bSJed Brown } 53*c4762a1bSJed Brown app->nbounces++; 54*c4762a1bSJed Brown PetscFunctionReturn(0); 55*c4762a1bSJed Brown } 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown /* 58*c4762a1bSJed Brown Defines the ODE passed to the ODE solver in explicit form: U_t = F(U) 59*c4762a1bSJed Brown */ 60*c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx) 61*c4762a1bSJed Brown { 62*c4762a1bSJed Brown PetscErrorCode ierr; 63*c4762a1bSJed Brown PetscScalar *f; 64*c4762a1bSJed Brown const PetscScalar *u; 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown PetscFunctionBegin; 67*c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 68*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown f[0] = u[1]; 72*c4762a1bSJed Brown f[1] = - 9.8; 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 76*c4762a1bSJed Brown PetscFunctionReturn(0); 77*c4762a1bSJed Brown } 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown /* 80*c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of the Jacobian. 81*c4762a1bSJed Brown */ 82*c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx) 83*c4762a1bSJed Brown { 84*c4762a1bSJed Brown PetscErrorCode ierr; 85*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 86*c4762a1bSJed Brown PetscScalar J[2][2]; 87*c4762a1bSJed Brown const PetscScalar *u; 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown PetscFunctionBegin; 90*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown J[0][0] = 0.0; J[0][1] = 1.0; 93*c4762a1bSJed Brown J[1][0] = 0.0; J[1][1] = 0.0; 94*c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 100*c4762a1bSJed Brown if (A != B) { 101*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 103*c4762a1bSJed Brown } 104*c4762a1bSJed Brown PetscFunctionReturn(0); 105*c4762a1bSJed Brown } 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown /* 108*c4762a1bSJed Brown Defines the ODE passed to the ODE solver in implicit form: F(U_t,U) = 0 109*c4762a1bSJed Brown */ 110*c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 111*c4762a1bSJed Brown { 112*c4762a1bSJed Brown PetscErrorCode ierr; 113*c4762a1bSJed Brown PetscScalar *f; 114*c4762a1bSJed Brown const PetscScalar *u,*udot; 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown PetscFunctionBegin; 117*c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 118*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown f[0] = udot[0] - u[1]; 123*c4762a1bSJed Brown f[1] = udot[1] + 9.8; 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 128*c4762a1bSJed Brown PetscFunctionReturn(0); 129*c4762a1bSJed Brown } 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown /* 132*c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 133*c4762a1bSJed Brown */ 134*c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 135*c4762a1bSJed Brown { 136*c4762a1bSJed Brown PetscErrorCode ierr; 137*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 138*c4762a1bSJed Brown PetscScalar J[2][2]; 139*c4762a1bSJed Brown const PetscScalar *u,*udot; 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown PetscFunctionBegin; 142*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown J[0][0] = a; J[0][1] = -1.0; 146*c4762a1bSJed Brown J[1][0] = 0.0; J[1][1] = a; 147*c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 150*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 154*c4762a1bSJed Brown if (A != B) { 155*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157*c4762a1bSJed Brown } 158*c4762a1bSJed Brown PetscFunctionReturn(0); 159*c4762a1bSJed Brown } 160*c4762a1bSJed Brown 161*c4762a1bSJed Brown int main(int argc,char **argv) 162*c4762a1bSJed Brown { 163*c4762a1bSJed Brown TS ts; /* ODE integrator */ 164*c4762a1bSJed Brown Vec U; /* solution will be stored here */ 165*c4762a1bSJed Brown PetscErrorCode ierr; 166*c4762a1bSJed Brown PetscMPIInt size; 167*c4762a1bSJed Brown PetscInt n = 2; 168*c4762a1bSJed Brown PetscScalar *u; 169*c4762a1bSJed Brown AppCtx app; 170*c4762a1bSJed Brown PetscInt direction[2]; 171*c4762a1bSJed Brown PetscBool terminate[2]; 172*c4762a1bSJed Brown PetscBool rhs_form=PETSC_FALSE,hist=PETSC_TRUE; 173*c4762a1bSJed Brown TSAdapt adapt; 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 176*c4762a1bSJed Brown Initialize program 177*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 178*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 179*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 180*c4762a1bSJed Brown if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown app.nbounces = 0; 183*c4762a1bSJed Brown app.maxbounces = 10; 184*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex40 options","");CHKERRQ(ierr); 185*c4762a1bSJed Brown ierr = PetscOptionsInt("-maxbounces","","",app.maxbounces,&app.maxbounces,NULL);CHKERRQ(ierr); 186*c4762a1bSJed Brown ierr = PetscOptionsBool("-test_adapthistory","","",hist,&hist,NULL);CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 190*c4762a1bSJed Brown Create timestepping solver context 191*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 192*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 193*c4762a1bSJed Brown ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr); 194*c4762a1bSJed Brown 195*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 196*c4762a1bSJed Brown Set ODE routines 197*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 198*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 199*c4762a1bSJed Brown /* Users are advised against the following branching and code duplication. 200*c4762a1bSJed Brown For problems without a mass matrix like the one at hand, the RHSFunction 201*c4762a1bSJed Brown (and companion RHSJacobian) interface is enough to support both explicit 202*c4762a1bSJed Brown and implicit timesteppers. This tutorial example also deals with the 203*c4762a1bSJed Brown IFunction/IJacobian interface for demonstration and testing purposes. */ 204*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-rhs-form",&rhs_form,NULL);CHKERRQ(ierr); 205*c4762a1bSJed Brown if (rhs_form) { 206*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,NULL);CHKERRQ(ierr); 207*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,NULL);CHKERRQ(ierr); 208*c4762a1bSJed Brown } else { 209*c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 210*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 211*c4762a1bSJed Brown ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 212*c4762a1bSJed Brown ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 213*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 214*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 215*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,NULL);CHKERRQ(ierr); 216*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,IJacobian,NULL);CHKERRQ(ierr); 217*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 218*c4762a1bSJed Brown } 219*c4762a1bSJed Brown 220*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 221*c4762a1bSJed Brown Set initial conditions 222*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 223*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&U);CHKERRQ(ierr); 224*c4762a1bSJed Brown ierr = VecSetSizes(U,n,PETSC_DETERMINE);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = VecSetUp(U);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 227*c4762a1bSJed Brown u[0] = 0.0; 228*c4762a1bSJed Brown u[1] = 20.0; 229*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 230*c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 231*c4762a1bSJed Brown 232*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 233*c4762a1bSJed Brown Set solver options 234*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 235*c4762a1bSJed Brown ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,30.0);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.1);CHKERRQ(ierr); 239*c4762a1bSJed Brown /* The adapative time step controller could take very large timesteps resulting in 240*c4762a1bSJed Brown the same event occuring multiple times in the same interval. A maximum step size 241*c4762a1bSJed Brown limit is enforced here to avoid this issue. */ 242*c4762a1bSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 243*c4762a1bSJed Brown ierr = TSAdaptSetStepLimits(adapt,0.0,0.5);CHKERRQ(ierr); 244*c4762a1bSJed Brown 245*c4762a1bSJed Brown /* Set directions and terminate flags for the two events */ 246*c4762a1bSJed Brown direction[0] = -1; direction[1] = -1; 247*c4762a1bSJed Brown terminate[0] = PETSC_FALSE; terminate[1] = PETSC_TRUE; 248*c4762a1bSJed Brown ierr = TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&app);CHKERRQ(ierr); 249*c4762a1bSJed Brown 250*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 251*c4762a1bSJed Brown 252*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 253*c4762a1bSJed Brown Run timestepping solver 254*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 255*c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown if (hist) { /* replay following history */ 258*c4762a1bSJed Brown TSTrajectory tj; 259*c4762a1bSJed Brown PetscReal tf,t0,dt; 260*c4762a1bSJed Brown 261*c4762a1bSJed Brown app.nbounces = 0; 262*c4762a1bSJed Brown ierr = TSGetTime(ts,&tf);CHKERRQ(ierr); 263*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,tf);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr); 265*c4762a1bSJed Brown ierr = TSRestartStep(ts);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 267*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 268*c4762a1bSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 269*c4762a1bSJed Brown ierr = TSAdaptSetType(adapt,TSADAPTHISTORY);CHKERRQ(ierr); 270*c4762a1bSJed Brown ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr); 271*c4762a1bSJed Brown ierr = TSAdaptHistorySetTrajectory(adapt,tj,PETSC_FALSE);CHKERRQ(ierr); 272*c4762a1bSJed Brown ierr = TSAdaptHistoryGetStep(adapt,0,&t0,&dt);CHKERRQ(ierr); 273*c4762a1bSJed Brown /* this example fails with single (or smaller) precision */ 274*c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL__FP16) 275*c4762a1bSJed Brown ierr = TSAdaptSetType(adapt,TSADAPTBASIC);CHKERRQ(ierr); 276*c4762a1bSJed Brown ierr = TSAdaptSetStepLimits(adapt,0.0,0.5);CHKERRQ(ierr); 277*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 278*c4762a1bSJed Brown #endif 279*c4762a1bSJed Brown ierr = TSSetTime(ts,t0);CHKERRQ(ierr); 280*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 281*c4762a1bSJed Brown ierr = TSResetTrajectory(ts);CHKERRQ(ierr); 282*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 283*c4762a1bSJed Brown u[0] = 0.0; 284*c4762a1bSJed Brown u[1] = 20.0; 285*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 286*c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 287*c4762a1bSJed Brown } 288*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 289*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 290*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 291*c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 292*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 293*c4762a1bSJed Brown 294*c4762a1bSJed Brown ierr = PetscFinalize(); 295*c4762a1bSJed Brown return ierr; 296*c4762a1bSJed Brown } 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown /*TEST 299*c4762a1bSJed Brown 300*c4762a1bSJed Brown test: 301*c4762a1bSJed Brown suffix: a 302*c4762a1bSJed Brown args: -snes_stol 1e-4 -ts_trajectory_dirname ex40_a_dir 303*c4762a1bSJed Brown output_file: output/ex40.out 304*c4762a1bSJed Brown 305*c4762a1bSJed Brown test: 306*c4762a1bSJed Brown suffix: b 307*c4762a1bSJed Brown args: -ts_type arkimex -snes_stol 1e-4 -ts_trajectory_dirname ex40_b_dir 308*c4762a1bSJed Brown output_file: output/ex40.out 309*c4762a1bSJed Brown 310*c4762a1bSJed Brown test: 311*c4762a1bSJed Brown suffix: c 312*c4762a1bSJed Brown args: -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_c_dir 313*c4762a1bSJed Brown output_file: output/ex40.out 314*c4762a1bSJed Brown 315*c4762a1bSJed Brown test: 316*c4762a1bSJed Brown suffix: d 317*c4762a1bSJed Brown args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_d_dir 318*c4762a1bSJed Brown output_file: output/ex40.out 319*c4762a1bSJed Brown 320*c4762a1bSJed Brown test: 321*c4762a1bSJed Brown suffix: e 322*c4762a1bSJed Brown args: -ts_type bdf -ts_adapt_dt_max 0.025 -ts_max_steps 1500 -ts_trajectory_dirname ex40_e_dir 323*c4762a1bSJed Brown output_file: output/ex40.out 324*c4762a1bSJed Brown 325*c4762a1bSJed Brown test: 326*c4762a1bSJed Brown suffix: f 327*c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 3bs -ts_trajectory_dirname ex40_f_dir 328*c4762a1bSJed Brown output_file: output/ex40.out 329*c4762a1bSJed Brown 330*c4762a1bSJed Brown test: 331*c4762a1bSJed Brown suffix: g 332*c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 5bs -ts_trajectory_dirname ex40_g_dir 333*c4762a1bSJed Brown output_file: output/ex40.out 334*c4762a1bSJed Brown 335*c4762a1bSJed Brown test: 336*c4762a1bSJed Brown suffix: h 337*c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 6vr -ts_trajectory_dirname ex40_h_dir 338*c4762a1bSJed Brown output_file: output/ex40.out 339*c4762a1bSJed Brown 340*c4762a1bSJed Brown test: 341*c4762a1bSJed Brown suffix: i 342*c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 7vr -ts_trajectory_dirname ex40_i_dir 343*c4762a1bSJed Brown output_file: output/ex40.out 344*c4762a1bSJed Brown 345*c4762a1bSJed Brown test: 346*c4762a1bSJed Brown suffix: j 347*c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 8vr -ts_trajectory_dirname ex40_j_dir 348*c4762a1bSJed Brown output_file: output/ex40.out 349*c4762a1bSJed Brown 350*c4762a1bSJed Brown TEST*/ 351