1 static char help[] = "Adjoint sensitivity of a hybrid system with state-dependent switchings.\n"; 2 3 /* 4 The dynamics is described by the ODE 5 u_t = A_i u 6 7 where A_1 = [ 1 -100 8 10 1 ], 9 A_2 = [ 1 10 10 -100 1 ]. 11 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]. 12 Initially u=[0 1]^T and i=1. 13 14 References: 15 + * - 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 16 - * - I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000 17 */ 18 19 #include <petscts.h> 20 21 typedef struct { 22 PetscScalar lambda1; 23 PetscScalar lambda2; 24 PetscInt mode; /* mode flag*/ 25 } AppCtx; 26 27 PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx) { 28 AppCtx *actx = (AppCtx *)ctx; 29 const PetscScalar *u; 30 31 PetscFunctionBegin; 32 PetscCall(VecGetArrayRead(U, &u)); 33 if (actx->mode == 1) { 34 fvalue[0] = u[1] - actx->lambda1 * u[0]; 35 } else if (actx->mode == 2) { 36 fvalue[0] = u[1] - actx->lambda2 * u[0]; 37 } 38 PetscCall(VecRestoreArrayRead(U, &u)); 39 PetscFunctionReturn(0); 40 } 41 42 PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx) { 43 Vec *lambda, *mu; 44 PetscScalar *x, *y; 45 const PetscScalar *u; 46 PetscScalar tmp[2], A1[2][2], A2[2], denorm; 47 PetscInt numcost; 48 49 PetscFunctionBegin; 50 PetscCall(TSGetCostGradients(ts, &numcost, &lambda, &mu)); 51 PetscCall(VecGetArrayRead(U, &u)); 52 53 if (actx->mode == 2) { 54 denorm = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]); 55 A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm + 1.; 56 A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm; 57 A1[1][0] = 110. * u[1] * 1. / denorm; 58 A1[1][1] = -110. * u[0] * 1. / denorm + 1.; 59 60 A2[0] = 110. * u[1] * (-u[0]) / denorm; 61 A2[1] = -110. * u[0] * (-u[0]) / denorm; 62 } else { 63 denorm = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]); 64 A1[0][0] = 110. * u[1] * (actx->lambda2) / denorm + 1; 65 A1[0][1] = -110. * u[0] * (actx->lambda2) / denorm; 66 A1[1][0] = -110. * u[1] * 1. / denorm; 67 A1[1][1] = 110. * u[0] * 1. / denorm + 1.; 68 69 A2[0] = 0; 70 A2[1] = 0; 71 } 72 73 PetscCall(VecRestoreArrayRead(U, &u)); 74 75 PetscCall(VecGetArray(lambda[0], &x)); 76 PetscCall(VecGetArray(mu[0], &y)); 77 tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1]; 78 tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1]; 79 y[0] = y[0] + A2[0] * x[0] + A2[1] * x[1]; 80 x[0] = tmp[0]; 81 x[1] = tmp[1]; 82 PetscCall(VecRestoreArray(mu[0], &y)); 83 PetscCall(VecRestoreArray(lambda[0], &x)); 84 85 PetscCall(VecGetArray(lambda[1], &x)); 86 PetscCall(VecGetArray(mu[1], &y)); 87 tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1]; 88 tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1]; 89 y[0] = y[0] + A2[0] * x[0] + A2[1] * x[1]; 90 x[0] = tmp[0]; 91 x[1] = tmp[1]; 92 PetscCall(VecRestoreArray(mu[1], &y)); 93 PetscCall(VecRestoreArray(lambda[1], &x)); 94 PetscFunctionReturn(0); 95 } 96 97 PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) { 98 AppCtx *actx = (AppCtx *)ctx; 99 100 PetscFunctionBegin; 101 /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */ 102 if (!forwardsolve) { PetscCall(ShiftGradients(ts, U, actx)); } 103 if (actx->mode == 1) { 104 actx->mode = 2; 105 /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t)); */ 106 } else if (actx->mode == 2) { 107 actx->mode = 1; 108 /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t)); */ 109 } 110 PetscFunctionReturn(0); 111 } 112 113 /* 114 Defines the ODE passed to the ODE solver 115 */ 116 static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) { 117 AppCtx *actx = (AppCtx *)ctx; 118 PetscScalar *f; 119 const PetscScalar *u, *udot; 120 121 PetscFunctionBegin; 122 /* The next three lines allow us to access the entries of the vectors directly */ 123 PetscCall(VecGetArrayRead(U, &u)); 124 PetscCall(VecGetArrayRead(Udot, &udot)); 125 PetscCall(VecGetArray(F, &f)); 126 127 if (actx->mode == 1) { 128 f[0] = udot[0] - u[0] + 100 * u[1]; 129 f[1] = udot[1] - 10 * u[0] - u[1]; 130 } else if (actx->mode == 2) { 131 f[0] = udot[0] - u[0] - 10 * u[1]; 132 f[1] = udot[1] + 100 * u[0] - u[1]; 133 } 134 135 PetscCall(VecRestoreArrayRead(U, &u)); 136 PetscCall(VecRestoreArrayRead(Udot, &udot)); 137 PetscCall(VecRestoreArray(F, &f)); 138 PetscFunctionReturn(0); 139 } 140 141 /* 142 Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 143 */ 144 static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) { 145 AppCtx *actx = (AppCtx *)ctx; 146 PetscInt rowcol[] = {0, 1}; 147 PetscScalar J[2][2]; 148 const PetscScalar *u, *udot; 149 150 PetscFunctionBegin; 151 PetscCall(VecGetArrayRead(U, &u)); 152 PetscCall(VecGetArrayRead(Udot, &udot)); 153 154 if (actx->mode == 1) { 155 J[0][0] = a - 1; 156 J[0][1] = 100; 157 J[1][0] = -10; 158 J[1][1] = a - 1; 159 } else if (actx->mode == 2) { 160 J[0][0] = a - 1; 161 J[0][1] = -10; 162 J[1][0] = 100; 163 J[1][1] = a - 1; 164 } 165 PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 166 167 PetscCall(VecRestoreArrayRead(U, &u)); 168 PetscCall(VecRestoreArrayRead(Udot, &udot)); 169 170 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 171 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 172 if (A != B) { 173 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 174 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 175 } 176 PetscFunctionReturn(0); 177 } 178 179 /* Matrix JacobianP is constant so that it only needs to be evaluated once */ 180 static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx) { 181 PetscFunctionBeginUser; 182 PetscFunctionReturn(0); 183 } 184 185 int main(int argc, char **argv) { 186 TS ts; /* ODE integrator */ 187 Vec U; /* solution will be stored here */ 188 Mat A; /* Jacobian matrix */ 189 Mat Ap; /* dfdp */ 190 PetscMPIInt size; 191 PetscInt n = 2; 192 PetscScalar *u, *v; 193 AppCtx app; 194 PetscInt direction[1]; 195 PetscBool terminate[1]; 196 Vec lambda[2], mu[2]; 197 PetscReal tend; 198 199 FILE *f; 200 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 201 Initialize program 202 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 203 PetscFunctionBeginUser; 204 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 205 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 206 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 207 app.mode = 1; 208 app.lambda1 = 2.75; 209 app.lambda2 = 0.36; 210 tend = 0.125; 211 PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1adj options", ""); 212 { 213 PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL)); 214 PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL)); 215 PetscCall(PetscOptionsReal("-tend", "", "", tend, &tend, NULL)); 216 } 217 PetscOptionsEnd(); 218 219 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 220 Create necessary matrix and vectors 221 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 222 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 223 PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 224 PetscCall(MatSetType(A, MATDENSE)); 225 PetscCall(MatSetFromOptions(A)); 226 PetscCall(MatSetUp(A)); 227 228 PetscCall(MatCreateVecs(A, &U, NULL)); 229 230 PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap)); 231 PetscCall(MatSetSizes(Ap, n, 1, PETSC_DETERMINE, PETSC_DETERMINE)); 232 PetscCall(MatSetType(Ap, MATDENSE)); 233 PetscCall(MatSetFromOptions(Ap)); 234 PetscCall(MatSetUp(Ap)); 235 PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */ 236 237 PetscCall(VecGetArray(U, &u)); 238 u[0] = 0; 239 u[1] = 1; 240 PetscCall(VecRestoreArray(U, &u)); 241 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 242 Create timestepping solver context 243 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 244 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 245 PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 246 PetscCall(TSSetType(ts, TSCN)); 247 PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app)); 248 PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app)); 249 PetscCall(TSSetRHSJacobianP(ts, Ap, RHSJacobianP, &app)); 250 251 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 252 Set initial conditions 253 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 254 PetscCall(TSSetSolution(ts, U)); 255 256 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 257 Save trajectory of solution so that TSAdjointSolve() may be used 258 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 259 PetscCall(TSSetSaveTrajectory(ts)); 260 261 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 262 Set solver options 263 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 264 PetscCall(TSSetMaxTime(ts, tend)); 265 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 266 PetscCall(TSSetTimeStep(ts, 1. / 256.)); 267 PetscCall(TSSetFromOptions(ts)); 268 269 /* Set directions and terminate flags for the two events */ 270 direction[0] = 0; 271 terminate[0] = PETSC_FALSE; 272 PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app)); 273 274 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 275 Run timestepping solver 276 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 277 PetscCall(TSSolve(ts, U)); 278 279 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 280 Adjoint model starts here 281 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 282 PetscCall(MatCreateVecs(A, &lambda[0], NULL)); 283 PetscCall(MatCreateVecs(A, &lambda[1], NULL)); 284 /* Set initial conditions for the adjoint integration */ 285 PetscCall(VecZeroEntries(lambda[0])); 286 PetscCall(VecZeroEntries(lambda[1])); 287 PetscCall(VecGetArray(lambda[0], &u)); 288 u[0] = 1.; 289 PetscCall(VecRestoreArray(lambda[0], &u)); 290 PetscCall(VecGetArray(lambda[1], &u)); 291 u[1] = 1.; 292 PetscCall(VecRestoreArray(lambda[1], &u)); 293 294 PetscCall(MatCreateVecs(Ap, &mu[0], NULL)); 295 PetscCall(MatCreateVecs(Ap, &mu[1], NULL)); 296 PetscCall(VecZeroEntries(mu[0])); 297 PetscCall(VecZeroEntries(mu[1])); 298 PetscCall(TSSetCostGradients(ts, 2, lambda, mu)); 299 300 PetscCall(TSAdjointSolve(ts)); 301 302 /* 303 PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD)); 304 PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD)); 305 PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD)); 306 PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD)); 307 */ 308 PetscCall(VecGetArray(mu[0], &u)); 309 PetscCall(VecGetArray(mu[1], &v)); 310 f = fopen("adj_mu.out", "a"); 311 PetscCall(PetscFPrintf(PETSC_COMM_WORLD, f, "%20.15lf %20.15lf %20.15lf\n", (double)tend, (double)PetscRealPart(u[0]), (double)PetscRealPart(v[0]))); 312 PetscCall(VecRestoreArray(mu[0], &u)); 313 PetscCall(VecRestoreArray(mu[1], &v)); 314 fclose(f); 315 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 316 Free work space. All PETSc objects should be destroyed when they are no longer needed. 317 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 318 PetscCall(MatDestroy(&A)); 319 PetscCall(VecDestroy(&U)); 320 PetscCall(TSDestroy(&ts)); 321 322 PetscCall(MatDestroy(&Ap)); 323 PetscCall(VecDestroy(&lambda[0])); 324 PetscCall(VecDestroy(&lambda[1])); 325 PetscCall(VecDestroy(&mu[0])); 326 PetscCall(VecDestroy(&mu[1])); 327 PetscCall(PetscFinalize()); 328 return 0; 329 } 330 331 /*TEST 332 333 build: 334 requires: !complex 335 336 test: 337 args: -ts_monitor -ts_adjoint_monitor 338 339 TEST*/ 340