1c4762a1bSJed Brown static char help[] = "An example of hybrid system using TS event.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown The dynamics is described by the ODE 5c4762a1bSJed Brown u_t = A_i u 6c4762a1bSJed Brown 7c4762a1bSJed Brown where A_1 = [ 1 -100 8c4762a1bSJed Brown 10 1 ], 9c4762a1bSJed Brown A_2 = [ 1 10 10c4762a1bSJed Brown -100 1 ]. 11c4762a1bSJed Brown The index i changes from 1 to 2 when u[1]=2.75u[0] and from 2 to 1 when u[1]=0.36u[0]. 12c4762a1bSJed Brown Initially u=[0 1]^T and i=1. 13c4762a1bSJed Brown 14c4762a1bSJed Brown Reference: 15c4762a1bSJed Brown I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000 16c4762a1bSJed Brown */ 17c4762a1bSJed Brown 18c4762a1bSJed Brown #include <petscts.h> 19c4762a1bSJed Brown 20c4762a1bSJed Brown typedef struct { 21c4762a1bSJed Brown PetscReal lambda1; 22c4762a1bSJed Brown PetscReal lambda2; 23c4762a1bSJed Brown PetscInt mode; /* mode flag*/ 24c4762a1bSJed Brown } AppCtx; 25c4762a1bSJed Brown 26c4762a1bSJed Brown PetscErrorCode FWDRun(TS, Vec, void *); 27c4762a1bSJed Brown 289371c9d4SSatish Balay PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx) { 29c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 30c4762a1bSJed Brown const PetscScalar *u; 31c4762a1bSJed Brown 32c4762a1bSJed Brown PetscFunctionBegin; 339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 34c4762a1bSJed Brown if (actx->mode == 1) { 35c4762a1bSJed Brown fvalue[0] = u[1] - actx->lambda1 * u[0]; 36c4762a1bSJed Brown } else if (actx->mode == 2) { 37c4762a1bSJed Brown fvalue[0] = u[1] - actx->lambda2 * u[0]; 38c4762a1bSJed Brown } 399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 40c4762a1bSJed Brown PetscFunctionReturn(0); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 439371c9d4SSatish Balay PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx) { 44c4762a1bSJed Brown Vec *lambda, *mu; 45c4762a1bSJed Brown PetscScalar *x, *y; 46c4762a1bSJed Brown const PetscScalar *u; 47c4762a1bSJed Brown PetscScalar tmp[2], A1[2][2], A2[2], denorm1, denorm2; 48c4762a1bSJed Brown PetscInt numcost; 49c4762a1bSJed Brown 50c4762a1bSJed Brown PetscFunctionBegin; 519566063dSJacob Faibussowitsch PetscCall(TSGetCostGradients(ts, &numcost, &lambda, &mu)); 529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown if (actx->mode == 2) { 55c4762a1bSJed Brown denorm1 = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]); 56c4762a1bSJed Brown denorm2 = -actx->lambda1 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]); 57c4762a1bSJed Brown A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm1 + 1.; 58c4762a1bSJed Brown A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm1; 59c4762a1bSJed Brown A1[1][0] = 110. * u[1] * 1. / denorm1; 60c4762a1bSJed Brown A1[1][1] = -110. * u[0] * 1. / denorm1 + 1.; 61c4762a1bSJed Brown 62c4762a1bSJed Brown A2[0] = 110. * u[1] * (-u[0]) / denorm2; 63c4762a1bSJed Brown A2[1] = -110. * u[0] * (-u[0]) / denorm2; 64c4762a1bSJed Brown } else { 65c4762a1bSJed Brown denorm2 = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]); 66c4762a1bSJed Brown A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm2 + 1; 67c4762a1bSJed Brown A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm2; 68c4762a1bSJed Brown A1[1][0] = 110. * u[1] * 1. / denorm2; 69c4762a1bSJed Brown A1[1][1] = -110. * u[0] * 1. / denorm2 + 1.; 70c4762a1bSJed Brown 71c4762a1bSJed Brown A2[0] = 0; 72c4762a1bSJed Brown A2[1] = 0; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 76c4762a1bSJed Brown 779566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &x)); 789566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &y)); 79c4762a1bSJed Brown tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1]; 80c4762a1bSJed Brown tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1]; 81c4762a1bSJed Brown y[0] = y[0] + A2[0] * x[0] + A2[1] * x[1]; 82c4762a1bSJed Brown x[0] = tmp[0]; 83c4762a1bSJed Brown x[1] = tmp[1]; 849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &y)); 859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &x)); 86c4762a1bSJed Brown 879566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[1], &x)); 889566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[1], &y)); 89c4762a1bSJed Brown tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1]; 90c4762a1bSJed Brown tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1]; 91c4762a1bSJed Brown y[0] = y[0] + A2[0] * x[0] + A2[1] * x[1]; 92c4762a1bSJed Brown x[0] = tmp[0]; 93c4762a1bSJed Brown x[1] = tmp[1]; 949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[1], &y)); 959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[1], &x)); 96c4762a1bSJed Brown PetscFunctionReturn(0); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 999371c9d4SSatish Balay PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) { 100c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 101c4762a1bSJed Brown 102c4762a1bSJed Brown PetscFunctionBegin; 103*48a46eb9SPierre Jolivet if (!forwardsolve) PetscCall(ShiftGradients(ts, U, actx)); 104c4762a1bSJed Brown if (actx->mode == 1) { 105c4762a1bSJed Brown actx->mode = 2; 1069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Change from mode 1 to 2 at t = %f \n", (double)t)); 107c4762a1bSJed Brown } else if (actx->mode == 2) { 108c4762a1bSJed Brown actx->mode = 1; 1099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Change from mode 2 to 1 at t = %f \n", (double)t)); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown PetscFunctionReturn(0); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* 115c4762a1bSJed Brown Defines the ODE passed to the ODE solver 116c4762a1bSJed Brown */ 1179371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) { 118c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 119c4762a1bSJed Brown PetscScalar *f; 120c4762a1bSJed Brown const PetscScalar *u, *udot; 121c4762a1bSJed Brown 122c4762a1bSJed Brown PetscFunctionBegin; 123c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 1249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 1269566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown if (actx->mode == 1) { 129c4762a1bSJed Brown f[0] = udot[0] - u[0] + 100 * u[1]; 130c4762a1bSJed Brown f[1] = udot[1] - 10 * u[0] - u[1]; 131c4762a1bSJed Brown } else if (actx->mode == 2) { 132c4762a1bSJed Brown f[0] = udot[0] - u[0] - 10 * u[1]; 133c4762a1bSJed Brown f[1] = udot[1] + 100 * u[0] - u[1]; 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 1389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 139c4762a1bSJed Brown PetscFunctionReturn(0); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* 143c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 144c4762a1bSJed Brown */ 1459371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) { 146c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 147c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 148c4762a1bSJed Brown PetscScalar J[2][2]; 149c4762a1bSJed Brown const PetscScalar *u, *udot; 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscFunctionBegin; 1529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown if (actx->mode == 1) { 1569371c9d4SSatish Balay J[0][0] = a - 1; 1579371c9d4SSatish Balay J[0][1] = 100; 1589371c9d4SSatish Balay J[1][0] = -10; 1599371c9d4SSatish Balay J[1][1] = a - 1; 160c4762a1bSJed Brown } else if (actx->mode == 2) { 1619371c9d4SSatish Balay J[0][0] = a - 1; 1629371c9d4SSatish Balay J[0][1] = -10; 1639371c9d4SSatish Balay J[1][0] = 100; 1649371c9d4SSatish Balay J[1][1] = a - 1; 165c4762a1bSJed Brown } 1669566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 167c4762a1bSJed Brown 1689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1729566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 173c4762a1bSJed Brown if (A != B) { 1749566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1759566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown PetscFunctionReturn(0); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 1809371c9d4SSatish Balay int main(int argc, char **argv) { 181c4762a1bSJed Brown TS ts; /* ODE integrator */ 182c4762a1bSJed Brown Vec U; /* solution will be stored here */ 183c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 184c4762a1bSJed Brown Mat Ap; /* dfdp */ 185c4762a1bSJed Brown PetscMPIInt size; 186c4762a1bSJed Brown PetscInt n = 2; 187c4762a1bSJed Brown PetscScalar *u; 188c4762a1bSJed Brown AppCtx app; 189c4762a1bSJed Brown PetscInt direction[1]; 190c4762a1bSJed Brown PetscBool terminate[1]; 191c4762a1bSJed Brown PetscReal delta, tmp[2], sensi[2]; 192c4762a1bSJed Brown 193c4762a1bSJed Brown delta = 1e-8; 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 196c4762a1bSJed Brown Initialize program 197c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 198327415f7SBarry Smith PetscFunctionBeginUser; 1999566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 2009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2013c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 202c4762a1bSJed Brown app.mode = 1; 203c4762a1bSJed Brown app.lambda1 = 2.75; 204c4762a1bSJed Brown app.lambda2 = 0.36; 205d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1 options", ""); 206c4762a1bSJed Brown { 2079566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL)); 2089566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL)); 209c4762a1bSJed Brown } 210d0609cedSBarry Smith PetscOptionsEnd(); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 213c4762a1bSJed Brown Create necessary matrix and vectors 214c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2159566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 2169566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 2179566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 2189566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2199566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 220c4762a1bSJed Brown 2219566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 222c4762a1bSJed Brown 2239566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap)); 2249566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Ap, n, 1, PETSC_DETERMINE, PETSC_DETERMINE)); 2259566063dSJacob Faibussowitsch PetscCall(MatSetType(Ap, MATDENSE)); 2269566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Ap)); 2279566063dSJacob Faibussowitsch PetscCall(MatSetUp(Ap)); 2289566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */ 229c4762a1bSJed Brown 2309566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 231c4762a1bSJed Brown u[0] = 0; 232c4762a1bSJed Brown u[1] = 1; 2339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 234c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 235c4762a1bSJed Brown Create timestepping solver context 236c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2379566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 2389566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 2399566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 2409566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app)); 2419566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app)); 242c4762a1bSJed Brown 243c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 244c4762a1bSJed Brown Set initial conditions 245c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2469566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 249c4762a1bSJed Brown Set solver options 250c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2519566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 0.125)); 2529566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 2539566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1. / 256.)); 2549566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 255c4762a1bSJed Brown 256c4762a1bSJed Brown /* Set directions and terminate flags for the two events */ 257c4762a1bSJed Brown direction[0] = 0; 258c4762a1bSJed Brown terminate[0] = PETSC_FALSE; 2599566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app)); 260c4762a1bSJed Brown 261c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 262c4762a1bSJed Brown Run timestepping solver 263c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2649566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 265c4762a1bSJed Brown 2669566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 267c4762a1bSJed Brown tmp[0] = u[0]; 268c4762a1bSJed Brown tmp[1] = u[1]; 269c4762a1bSJed Brown 270c4762a1bSJed Brown u[0] = 0 + delta; 271c4762a1bSJed Brown u[1] = 1; 2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 273c4762a1bSJed Brown 2749566063dSJacob Faibussowitsch PetscCall(FWDRun(ts, U, (void *)&app)); 275c4762a1bSJed Brown 2769566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 277c4762a1bSJed Brown sensi[0] = (u[0] - tmp[0]) / delta; 278c4762a1bSJed Brown sensi[1] = (u[1] - tmp[1]) / delta; 27963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "d x1(tf) /d x1(t0) = %f d x2(tf) / d x1(t0) = %f \n", (double)sensi[0], (double)sensi[1])); 280c4762a1bSJed Brown u[0] = 0; 281c4762a1bSJed Brown u[1] = 1 + delta; 2829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 283c4762a1bSJed Brown 2849566063dSJacob Faibussowitsch PetscCall(FWDRun(ts, U, (void *)&app)); 285c4762a1bSJed Brown 2869566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 287c4762a1bSJed Brown sensi[0] = (u[0] - tmp[0]) / delta; 288c4762a1bSJed Brown sensi[1] = (u[1] - tmp[1]) / delta; 28963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "d x1(tf) /d x2(t0) = %f d x2(tf) / d x2(t0) = %f \n", (double)sensi[0], (double)sensi[1])); 290c4762a1bSJed Brown u[0] = 0; 291c4762a1bSJed Brown u[1] = 1; 292c4762a1bSJed Brown app.lambda1 = app.lambda1 + delta; 2939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 294c4762a1bSJed Brown 2959566063dSJacob Faibussowitsch PetscCall(FWDRun(ts, U, (void *)&app)); 296c4762a1bSJed Brown 2979566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 298c4762a1bSJed Brown sensi[0] = (u[0] - tmp[0]) / delta; 299c4762a1bSJed Brown sensi[1] = (u[1] - tmp[1]) / delta; 30063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Final gradients: d x1(tf) /d p = %f d x2(tf) / d p = %f \n", (double)sensi[0], (double)sensi[1])); 3019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 302c4762a1bSJed Brown 303c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 304c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 305c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3089566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 309c4762a1bSJed Brown 3109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ap)); 3119566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 312b122ec5aSJacob Faibussowitsch return 0; 313c4762a1bSJed Brown } 314c4762a1bSJed Brown 3159371c9d4SSatish Balay PetscErrorCode FWDRun(TS ts, Vec U0, void *ctx0) { 316c4762a1bSJed Brown Vec U; /* solution will be stored here */ 317c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 318c4762a1bSJed Brown 319c4762a1bSJed Brown PetscFunctionBeginUser; 3209566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 3219566063dSJacob Faibussowitsch PetscCall(VecCopy(U0, U)); 322c4762a1bSJed Brown 323c4762a1bSJed Brown ctx->mode = 1; 324c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 325c4762a1bSJed Brown Run timestepping solver 326c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3279566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 328c4762a1bSJed Brown 3299566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscFunctionReturn(0); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown 334c4762a1bSJed Brown /*TEST 335c4762a1bSJed Brown 336c4762a1bSJed Brown build: 337dfd57a17SPierre Jolivet requires: !defined(PETSC_USE_CXXCOMPLEX) 338c4762a1bSJed Brown 339c4762a1bSJed Brown test: 340c4762a1bSJed Brown args: -ts_event_tol 1e-9 341c4762a1bSJed Brown timeoutfactor: 18 342c4762a1bSJed Brown requires: !single 343c4762a1bSJed Brown 344c4762a1bSJed Brown TEST*/ 345