1c4762a1bSJed Brown static char help[] = "Trajectory sensitivity of a hybrid system with state-dependent switchings.\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 References: 15606c0280SSatish Balay + * - H. Zhang, S. Abhyankar, E. Constantinescu, M. Mihai, Discrete Adjoint Sensitivity Analysis of Hybrid Dynamical Systems With Switching, IEEE Transactions on Circuits and Systems I: Regular Papers, 64(5), May 2017 16606c0280SSatish Balay - * - I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000 17c4762a1bSJed Brown */ 18c4762a1bSJed Brown 19c4762a1bSJed Brown #include <petscts.h> 20c4762a1bSJed Brown 21c4762a1bSJed Brown typedef struct { 22c4762a1bSJed Brown PetscScalar lambda1; 23c4762a1bSJed Brown PetscScalar lambda2; 24c4762a1bSJed Brown PetscInt mode; /* mode flag*/ 25c4762a1bSJed Brown PetscReal print_time; 26c4762a1bSJed Brown } AppCtx; 27c4762a1bSJed Brown 28*9371c9d4SSatish Balay PetscErrorCode MyMonitor(TS ts, PetscInt stepnum, PetscReal time, Vec U, void *ctx) { 29c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 30c4762a1bSJed Brown Mat sp; 31c4762a1bSJed Brown PetscScalar *u; 32c4762a1bSJed Brown PetscInt nump; 33c4762a1bSJed Brown FILE *f; 34c4762a1bSJed Brown 35c4762a1bSJed Brown PetscFunctionBegin; 36c4762a1bSJed Brown if (time >= actx->print_time) { 37c4762a1bSJed Brown actx->print_time += 1. / 256.; 389566063dSJacob Faibussowitsch PetscCall(TSForwardGetSensitivities(ts, &nump, &sp)); 399566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(sp, 2, &u)); 40c4762a1bSJed Brown f = fopen("fwd_sp.out", "a"); 4163a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, f, "%20.15lf %20.15lf %20.15lf\n", (double)time, (double)PetscRealPart(u[0]), (double)PetscRealPart(u[1]))); 429566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(sp, &u)); 43c4762a1bSJed Brown fclose(f); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown PetscFunctionReturn(0); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48*9371c9d4SSatish Balay PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx) { 49c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 50c4762a1bSJed Brown const PetscScalar *u; 51c4762a1bSJed Brown 52c4762a1bSJed Brown PetscFunctionBegin; 539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 54c4762a1bSJed Brown if (actx->mode == 1) { 55c4762a1bSJed Brown fvalue[0] = u[1] - actx->lambda1 * u[0]; 56c4762a1bSJed Brown } else if (actx->mode == 2) { 57c4762a1bSJed Brown fvalue[0] = u[1] - actx->lambda2 * u[0]; 58c4762a1bSJed Brown } 599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 60c4762a1bSJed Brown PetscFunctionReturn(0); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63*9371c9d4SSatish Balay PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx) { 64c4762a1bSJed Brown Mat sp; 65c4762a1bSJed Brown PetscScalar *x; 66c4762a1bSJed Brown PetscScalar *u; 67c4762a1bSJed Brown PetscScalar tmp[2], A1[2][2], A2[2], denorm; 68c4762a1bSJed Brown PetscInt nump; 69c4762a1bSJed Brown 70c4762a1bSJed Brown PetscFunctionBegin; 719566063dSJacob Faibussowitsch PetscCall(TSForwardGetSensitivities(ts, &nump, &sp)); 729566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown if (actx->mode == 1) { 75c4762a1bSJed Brown denorm = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]); 76c4762a1bSJed Brown A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm + 1.; 77c4762a1bSJed Brown A1[1][0] = -110. * u[0] * (-actx->lambda1) / denorm; 78c4762a1bSJed Brown A1[0][1] = 110. * u[1] * 1. / denorm; 79c4762a1bSJed Brown A1[1][1] = -110. * u[0] * 1. / denorm + 1.; 80c4762a1bSJed Brown 81c4762a1bSJed Brown A2[0] = 110. * u[1] * (-u[0]) / denorm; 82c4762a1bSJed Brown A2[1] = -110. * u[0] * (-u[0]) / denorm; 83c4762a1bSJed Brown } else { 84c4762a1bSJed Brown denorm = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]); 85c4762a1bSJed Brown A1[0][0] = 110. * u[1] * (actx->lambda2) / denorm + 1; 86c4762a1bSJed Brown A1[1][0] = -110. * u[0] * (actx->lambda2) / denorm; 87c4762a1bSJed Brown A1[0][1] = -110. * u[1] * 1. / denorm; 88c4762a1bSJed Brown A1[1][1] = 110. * u[0] * 1. / denorm + 1.; 89c4762a1bSJed Brown 90c4762a1bSJed Brown A2[0] = 0; 91c4762a1bSJed Brown A2[1] = 0; 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 95c4762a1bSJed Brown 969566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(sp, 0, &x)); 97c4762a1bSJed Brown tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1]; 98c4762a1bSJed Brown tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1]; 99c4762a1bSJed Brown x[0] = tmp[0]; 100c4762a1bSJed Brown x[1] = tmp[1]; 1019566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(sp, &x)); 102c4762a1bSJed Brown 1039566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(sp, 1, &x)); 104c4762a1bSJed Brown tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1]; 105c4762a1bSJed Brown tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1]; 106c4762a1bSJed Brown x[0] = tmp[0]; 107c4762a1bSJed Brown x[1] = tmp[1]; 1089566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(sp, &x)); 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(sp, 2, &x)); 111c4762a1bSJed Brown tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1]; 112c4762a1bSJed Brown tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1]; 113c4762a1bSJed Brown x[0] = tmp[0] + A2[0]; 114c4762a1bSJed Brown x[1] = tmp[1] + A2[1]; 1159566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(sp, &x)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscFunctionReturn(0); 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120*9371c9d4SSatish Balay PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) { 121c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 122c4762a1bSJed Brown 123c4762a1bSJed Brown PetscFunctionBegin; 1249566063dSJacob Faibussowitsch /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */ 1259566063dSJacob Faibussowitsch PetscCall(ShiftGradients(ts, U, actx)); 126c4762a1bSJed Brown if (actx->mode == 1) { 127c4762a1bSJed Brown actx->mode = 2; 1289566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t)); */ 129c4762a1bSJed Brown } else if (actx->mode == 2) { 130c4762a1bSJed Brown actx->mode = 1; 1319566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t)); */ 132c4762a1bSJed Brown } 133c4762a1bSJed Brown PetscFunctionReturn(0); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* 137c4762a1bSJed Brown Defines the ODE passed to the ODE solver 138c4762a1bSJed Brown */ 139*9371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) { 140c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 141c4762a1bSJed Brown PetscScalar *f; 142c4762a1bSJed Brown const PetscScalar *u, *udot; 143c4762a1bSJed Brown 144c4762a1bSJed Brown PetscFunctionBegin; 145c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 1469566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 1489566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown if (actx->mode == 1) { 151c4762a1bSJed Brown f[0] = udot[0] - u[0] + 100 * u[1]; 152c4762a1bSJed Brown f[1] = udot[1] - 10 * u[0] - u[1]; 153c4762a1bSJed Brown } else if (actx->mode == 2) { 154c4762a1bSJed Brown f[0] = udot[0] - u[0] - 10 * u[1]; 155c4762a1bSJed Brown f[1] = udot[1] + 100 * u[0] - u[1]; 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 1589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 1609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 161c4762a1bSJed Brown PetscFunctionReturn(0); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* 165c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 166c4762a1bSJed Brown */ 167*9371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) { 168c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 169c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 170c4762a1bSJed Brown PetscScalar J[2][2]; 171c4762a1bSJed Brown const PetscScalar *u, *udot; 172c4762a1bSJed Brown 173c4762a1bSJed Brown PetscFunctionBegin; 1749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown if (actx->mode == 1) { 178*9371c9d4SSatish Balay J[0][0] = a - 1; 179*9371c9d4SSatish Balay J[0][1] = 100; 180*9371c9d4SSatish Balay J[1][0] = -10; 181*9371c9d4SSatish Balay J[1][1] = a - 1; 182c4762a1bSJed Brown } else if (actx->mode == 2) { 183*9371c9d4SSatish Balay J[0][0] = a - 1; 184*9371c9d4SSatish Balay J[0][1] = -10; 185*9371c9d4SSatish Balay J[1][0] = 100; 186*9371c9d4SSatish Balay J[1][1] = a - 1; 187c4762a1bSJed Brown } 1889566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 189c4762a1bSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 192c4762a1bSJed Brown 1939566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1949566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 195c4762a1bSJed Brown if (A != B) { 1969566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1979566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 198c4762a1bSJed Brown } 199c4762a1bSJed Brown PetscFunctionReturn(0); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* Matrix JacobianP is constant so that it only needs to be evaluated once */ 203*9371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat Ap, void *ctx) { 204c4762a1bSJed Brown PetscFunctionBeginUser; 205c4762a1bSJed Brown PetscFunctionReturn(0); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208*9371c9d4SSatish Balay int main(int argc, char **argv) { 209c4762a1bSJed Brown TS ts; /* ODE integrator */ 210c4762a1bSJed Brown Vec U; /* solution will be stored here */ 211c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 212c4762a1bSJed Brown Mat Ap; /* Jacobian dfdp */ 213c4762a1bSJed Brown PetscMPIInt size; 214c4762a1bSJed Brown PetscInt n = 2; 215c4762a1bSJed Brown PetscScalar *u; 216c4762a1bSJed Brown AppCtx app; 217c4762a1bSJed Brown PetscInt direction[1]; 218c4762a1bSJed Brown PetscBool terminate[1]; 219c4762a1bSJed Brown Mat sp; 220c4762a1bSJed Brown PetscReal tend; 221c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 222c4762a1bSJed Brown Initialize program 223c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 224327415f7SBarry Smith PetscFunctionBeginUser; 2259566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 2269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2273c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 228c4762a1bSJed Brown app.mode = 1; 229c4762a1bSJed Brown app.lambda1 = 2.75; 230c4762a1bSJed Brown app.lambda2 = 0.36; 231c4762a1bSJed Brown app.print_time = 1. / 256.; 232c4762a1bSJed Brown tend = 0.125; 233d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1fwd options", ""); 234c4762a1bSJed Brown { 2359566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL)); 2369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL)); 2379566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tend", "", "", tend, &tend, NULL)); 238c4762a1bSJed Brown } 239d0609cedSBarry Smith PetscOptionsEnd(); 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 242c4762a1bSJed Brown Create necessary matrix and vectors 243c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2449566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 2459566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 2469566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 2479566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2489566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 249c4762a1bSJed Brown 2509566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 251c4762a1bSJed Brown 2529566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap)); 2539566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Ap, n, 3, PETSC_DETERMINE, PETSC_DETERMINE)); 2549566063dSJacob Faibussowitsch PetscCall(MatSetType(Ap, MATDENSE)); 2559566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Ap)); 2569566063dSJacob Faibussowitsch PetscCall(MatSetUp(Ap)); 2579566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Ap)); 258c4762a1bSJed Brown 2599566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, n, 3, NULL, &sp)); 2609566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(sp)); 2619566063dSJacob Faibussowitsch PetscCall(MatShift(sp, 1.0)); 262c4762a1bSJed Brown 2639566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 264c4762a1bSJed Brown u[0] = 0; 265c4762a1bSJed Brown u[1] = 1; 2669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 269c4762a1bSJed Brown Create timestepping solver context 270c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2719566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 2729566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 2739566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 2749566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app)); 2759566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app)); 276c4762a1bSJed Brown 277c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 278c4762a1bSJed Brown Set initial conditions 279c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2809566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 281c4762a1bSJed Brown 2829566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(ts, 3, sp)); 283c4762a1bSJed Brown /* Set RHS JacobianP */ 2849566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts, Ap, RHSJacobianP, &app)); 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 287c4762a1bSJed Brown Set solver options 288c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2899566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, tend)); 2909566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 2919566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1. / 256.)); 2929566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MyMonitor, &app, PETSC_NULL)); 2939566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 294c4762a1bSJed Brown 295c4762a1bSJed Brown /* Set directions and terminate flags for the two events */ 296c4762a1bSJed Brown direction[0] = 0; 297c4762a1bSJed Brown terminate[0] = PETSC_FALSE; 2989566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app)); 299c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 300c4762a1bSJed Brown Run timestepping solver 301c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3029566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 303c4762a1bSJed Brown 304c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 305c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 306c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3099566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 310c4762a1bSJed Brown 3119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ap)); 3129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sp)); 3139566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 314d0609cedSBarry Smith return (0); 315c4762a1bSJed Brown } 316c4762a1bSJed Brown 317c4762a1bSJed Brown /*TEST 318c4762a1bSJed Brown 319c4762a1bSJed Brown build: 320c4762a1bSJed Brown requires: !complex 321c4762a1bSJed Brown 322c4762a1bSJed Brown test: 323c4762a1bSJed Brown args: -ts_monitor 324c4762a1bSJed Brown 325c4762a1bSJed Brown TEST*/ 326