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