1c4762a1bSJed Brown static char help[] = "Adjoint 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 } AppCtx; 26c4762a1bSJed Brown 27fa9584fbSIlya Fursov PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal *fvalue, void *ctx) 28d71ae5a4SJacob Faibussowitsch { 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) { 35fa9584fbSIlya Fursov fvalue[0] = PetscRealPart(u[1] - actx->lambda1 * u[0]); 36c4762a1bSJed Brown } else if (actx->mode == 2) { 37fa9584fbSIlya Fursov fvalue[0] = PetscRealPart(u[1] - actx->lambda2 * u[0]); 38c4762a1bSJed Brown } 399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43d71ae5a4SJacob Faibussowitsch PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx) 44d71ae5a4SJacob Faibussowitsch { 45c4762a1bSJed Brown Vec *lambda, *mu; 46c4762a1bSJed Brown PetscScalar *x, *y; 47c4762a1bSJed Brown const PetscScalar *u; 48c4762a1bSJed Brown PetscScalar tmp[2], A1[2][2], A2[2], denorm; 49c4762a1bSJed Brown PetscInt numcost; 50c4762a1bSJed Brown 51c4762a1bSJed Brown PetscFunctionBegin; 529566063dSJacob Faibussowitsch PetscCall(TSGetCostGradients(ts, &numcost, &lambda, &mu)); 539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown if (actx->mode == 2) { 56c4762a1bSJed Brown denorm = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]); 57c4762a1bSJed Brown A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm + 1.; 58c4762a1bSJed Brown A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm; 59c4762a1bSJed Brown A1[1][0] = 110. * u[1] * 1. / denorm; 60c4762a1bSJed Brown A1[1][1] = -110. * u[0] * 1. / denorm + 1.; 61c4762a1bSJed Brown 62c4762a1bSJed Brown A2[0] = 110. * u[1] * (-u[0]) / denorm; 63c4762a1bSJed Brown A2[1] = -110. * u[0] * (-u[0]) / denorm; 64c4762a1bSJed Brown } else { 65c4762a1bSJed Brown denorm = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]); 66c4762a1bSJed Brown A1[0][0] = 110. * u[1] * (actx->lambda2) / denorm + 1; 67c4762a1bSJed Brown A1[0][1] = -110. * u[0] * (actx->lambda2) / denorm; 68c4762a1bSJed Brown A1[1][0] = -110. * u[1] * 1. / denorm; 69c4762a1bSJed Brown A1[1][1] = 110. * u[0] * 1. / denorm + 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)); 963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99d71ae5a4SJacob Faibussowitsch PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) 100d71ae5a4SJacob Faibussowitsch { 101c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 102c4762a1bSJed Brown 103c4762a1bSJed Brown PetscFunctionBegin; 1049566063dSJacob Faibussowitsch /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */ 10548a46eb9SPierre Jolivet if (!forwardsolve) PetscCall(ShiftGradients(ts, U, actx)); 106c4762a1bSJed Brown if (actx->mode == 1) { 107c4762a1bSJed Brown actx->mode = 2; 1089566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t)); */ 109c4762a1bSJed Brown } else if (actx->mode == 2) { 110c4762a1bSJed Brown actx->mode = 1; 1119566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t)); */ 112c4762a1bSJed Brown } 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* 117c4762a1bSJed Brown Defines the ODE passed to the ODE solver 118c4762a1bSJed Brown */ 119d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) 120d71ae5a4SJacob Faibussowitsch { 121c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 122c4762a1bSJed Brown PetscScalar *f; 123c4762a1bSJed Brown const PetscScalar *u, *udot; 124c4762a1bSJed Brown 125c4762a1bSJed Brown PetscFunctionBegin; 126c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 1279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 1299566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown if (actx->mode == 1) { 132c4762a1bSJed Brown f[0] = udot[0] - u[0] + 100 * u[1]; 133c4762a1bSJed Brown f[1] = udot[1] - 10 * u[0] - u[1]; 134c4762a1bSJed Brown } else if (actx->mode == 2) { 135c4762a1bSJed Brown f[0] = udot[0] - u[0] - 10 * u[1]; 136c4762a1bSJed Brown f[1] = udot[1] + 100 * u[0] - u[1]; 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 1419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 1423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* 146c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 147c4762a1bSJed Brown */ 148d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) 149d71ae5a4SJacob Faibussowitsch { 150c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx; 151c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 152c4762a1bSJed Brown PetscScalar J[2][2]; 153c4762a1bSJed Brown const PetscScalar *u, *udot; 154c4762a1bSJed Brown 155c4762a1bSJed Brown PetscFunctionBegin; 1569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown if (actx->mode == 1) { 1609371c9d4SSatish Balay J[0][0] = a - 1; 1619371c9d4SSatish Balay J[0][1] = 100; 1629371c9d4SSatish Balay J[1][0] = -10; 1639371c9d4SSatish Balay J[1][1] = a - 1; 164c4762a1bSJed Brown } else if (actx->mode == 2) { 1659371c9d4SSatish Balay J[0][0] = a - 1; 1669371c9d4SSatish Balay J[0][1] = -10; 1679371c9d4SSatish Balay J[1][0] = 100; 1689371c9d4SSatish Balay J[1][1] = a - 1; 169c4762a1bSJed Brown } 1709566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 171c4762a1bSJed Brown 1729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 174c4762a1bSJed Brown 1759566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1769566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 177c4762a1bSJed Brown if (A != B) { 1789566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1799566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 180c4762a1bSJed Brown } 1813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* Matrix JacobianP is constant so that it only needs to be evaluated once */ 185d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx) 186d71ae5a4SJacob Faibussowitsch { 187c4762a1bSJed Brown PetscFunctionBeginUser; 1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 192d71ae5a4SJacob Faibussowitsch { 193c4762a1bSJed Brown TS ts; /* ODE integrator */ 194c4762a1bSJed Brown Vec U; /* solution will be stored here */ 195c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 196c4762a1bSJed Brown Mat Ap; /* dfdp */ 197c4762a1bSJed Brown PetscMPIInt size; 198c4762a1bSJed Brown PetscInt n = 2; 199c4762a1bSJed Brown PetscScalar *u, *v; 200c4762a1bSJed Brown AppCtx app; 201c4762a1bSJed Brown PetscInt direction[1]; 202c4762a1bSJed Brown PetscBool terminate[1]; 203c4762a1bSJed Brown Vec lambda[2], mu[2]; 204c4762a1bSJed Brown PetscReal tend; 205c4762a1bSJed Brown 206c4762a1bSJed Brown FILE *f; 207c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 208c4762a1bSJed Brown Initialize program 209c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 210327415f7SBarry Smith PetscFunctionBeginUser; 2119566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 2129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2133c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 214c4762a1bSJed Brown app.mode = 1; 215c4762a1bSJed Brown app.lambda1 = 2.75; 216c4762a1bSJed Brown app.lambda2 = 0.36; 217c4762a1bSJed Brown tend = 0.125; 218d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1adj options", ""); 219c4762a1bSJed Brown { 2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL)); 2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tend", "", "", tend, &tend, NULL)); 223c4762a1bSJed Brown } 224d0609cedSBarry Smith PetscOptionsEnd(); 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 227c4762a1bSJed Brown Create necessary matrix and vectors 228c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2299566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 2309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 2319566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 2329566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2339566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 234c4762a1bSJed Brown 2359566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 236c4762a1bSJed Brown 2379566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap)); 2389566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Ap, n, 1, PETSC_DETERMINE, PETSC_DETERMINE)); 2399566063dSJacob Faibussowitsch PetscCall(MatSetType(Ap, MATDENSE)); 2409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Ap)); 2419566063dSJacob Faibussowitsch PetscCall(MatSetUp(Ap)); 2429566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */ 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 245c4762a1bSJed Brown u[0] = 0; 246c4762a1bSJed Brown u[1] = 1; 2479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 248c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 249c4762a1bSJed Brown Create timestepping solver context 250c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2519566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 2529566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 2539566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 2549566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app)); 2559566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app)); 2569566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts, Ap, RHSJacobianP, &app)); 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 259c4762a1bSJed Brown Set initial conditions 260c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2619566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 264c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 265c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2669566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 269c4762a1bSJed Brown Set solver options 270c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2719566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, tend)); 2729566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 2739566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1. / 256.)); 2749566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* Set directions and terminate flags for the two events */ 277c4762a1bSJed Brown direction[0] = 0; 278c4762a1bSJed Brown terminate[0] = PETSC_FALSE; 2799566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app)); 280c4762a1bSJed Brown 281c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 282c4762a1bSJed Brown Run timestepping solver 283c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2849566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 287c4762a1bSJed Brown Adjoint model starts here 288c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2899566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[0], NULL)); 2909566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[1], NULL)); 291c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 2929566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(lambda[0])); 2939566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(lambda[1])); 2949566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &u)); 295c4762a1bSJed Brown u[0] = 1.; 2969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &u)); 2979566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[1], &u)); 298c4762a1bSJed Brown u[1] = 1.; 2999566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[1], &u)); 300c4762a1bSJed Brown 3019566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Ap, &mu[0], NULL)); 3029566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Ap, &mu[1], NULL)); 3039566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mu[0])); 3049566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mu[1])); 3059566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 2, lambda, mu)); 306c4762a1bSJed Brown 3079566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 308c4762a1bSJed Brown 309c4762a1bSJed Brown /* 3109566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD)); 3119566063dSJacob Faibussowitsch PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD)); 3129566063dSJacob Faibussowitsch PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD)); 3139566063dSJacob Faibussowitsch PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD)); 314c4762a1bSJed Brown */ 3159566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &u)); 3169566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[1], &v)); 317c4762a1bSJed Brown f = fopen("adj_mu.out", "a"); 31863a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, f, "%20.15lf %20.15lf %20.15lf\n", (double)tend, (double)PetscRealPart(u[0]), (double)PetscRealPart(v[0]))); 3199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &u)); 3209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[1], &v)); 321c4762a1bSJed Brown fclose(f); 322c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 323c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 324c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3279566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 328c4762a1bSJed Brown 3299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ap)); 3309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 3319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[1])); 3329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0])); 3339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[1])); 3349566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 335b122ec5aSJacob Faibussowitsch return 0; 336c4762a1bSJed Brown } 337c4762a1bSJed Brown 338c4762a1bSJed Brown /*TEST 339c4762a1bSJed Brown 340c4762a1bSJed Brown build: 341c4762a1bSJed Brown requires: !complex 342c4762a1bSJed Brown 343c4762a1bSJed Brown test: 344*fe4ad979SIlya Fursov args: -ts_monitor -ts_adjoint_monitor 345c4762a1bSJed Brown 346c4762a1bSJed Brown TEST*/ 347