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