Actual source code: ex1adj.c

  1: static char help[] = "Adjoint sensitivity of a hybrid system with state-dependent switchings.\n";

  3: /*
  4:   The dynamics is described by the ODE
  5:                   u_t = A_i u

  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.

 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: */

 20: #include <petscts.h>

 22: typedef struct {
 23:   PetscScalar lambda1;
 24:   PetscScalar lambda2;
 25:   PetscInt    mode; /* mode flag*/
 26: } AppCtx;

 28: PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal *fvalue, void *ctx)
 29: {
 30:   AppCtx            *actx = (AppCtx *)ctx;
 31:   const PetscScalar *u;

 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: }

 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], denorm;
 50:   PetscInt           numcost;

 52:   PetscFunctionBegin;
 53:   PetscCall(TSGetCostGradients(ts, &numcost, &lambda, &mu));
 54:   PetscCall(VecGetArrayRead(U, &u));

 56:   if (actx->mode == 2) {
 57:     denorm   = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]);
 58:     A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm + 1.;
 59:     A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm;
 60:     A1[1][0] = 110. * u[1] * 1. / denorm;
 61:     A1[1][1] = -110. * u[0] * 1. / denorm + 1.;

 63:     A2[0] = 110. * u[1] * (-u[0]) / denorm;
 64:     A2[1] = -110. * u[0] * (-u[0]) / denorm;
 65:   } else {
 66:     denorm   = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
 67:     A1[0][0] = 110. * u[1] * (actx->lambda2) / denorm + 1;
 68:     A1[0][1] = -110. * u[0] * (actx->lambda2) / denorm;
 69:     A1[1][0] = -110. * u[1] * 1. / denorm;
 70:     A1[1][1] = 110. * u[0] * 1. / denorm + 1.;

 72:     A2[0] = 0;
 73:     A2[1] = 0;
 74:   }

 76:   PetscCall(VecRestoreArrayRead(U, &u));

 78:   PetscCall(VecGetArray(lambda[0], &x));
 79:   PetscCall(VecGetArray(mu[0], &y));
 80:   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
 81:   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
 82:   y[0]   = y[0] + A2[0] * x[0] + A2[1] * x[1];
 83:   x[0]   = tmp[0];
 84:   x[1]   = tmp[1];
 85:   PetscCall(VecRestoreArray(mu[0], &y));
 86:   PetscCall(VecRestoreArray(lambda[0], &x));

 88:   PetscCall(VecGetArray(lambda[1], &x));
 89:   PetscCall(VecGetArray(mu[1], &y));
 90:   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
 91:   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
 92:   y[0]   = y[0] + A2[0] * x[0] + A2[1] * x[1];
 93:   x[0]   = tmp[0];
 94:   x[1]   = tmp[1];
 95:   PetscCall(VecRestoreArray(mu[1], &y));
 96:   PetscCall(VecRestoreArray(lambda[1], &x));
 97:   PetscFunctionReturn(PETSC_SUCCESS);
 98: }

100: PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx)
101: {
102:   AppCtx *actx = (AppCtx *)ctx;

104:   PetscFunctionBegin;
105:   /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */
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: }

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;

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));

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:   }

140:   PetscCall(VecRestoreArrayRead(U, &u));
141:   PetscCall(VecRestoreArrayRead(Udot, &udot));
142:   PetscCall(VecRestoreArray(F, &f));
143:   PetscFunctionReturn(PETSC_SUCCESS);
144: }

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;

156:   PetscFunctionBegin;
157:   PetscCall(VecGetArrayRead(U, &u));
158:   PetscCall(VecGetArrayRead(Udot, &udot));

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));

173:   PetscCall(VecRestoreArrayRead(U, &u));
174:   PetscCall(VecRestoreArrayRead(Udot, &udot));

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: }

185: /* Matrix JacobianP is constant so that it only needs to be evaluated once */
186: static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx)
187: {
188:   PetscFunctionBeginUser;
189:   PetscFunctionReturn(PETSC_SUCCESS);
190: }

192: int main(int argc, char **argv)
193: {
194:   TS           ts; /* ODE integrator */
195:   Vec          U;  /* solution will be stored here */
196:   Mat          A;  /* Jacobian matrix */
197:   Mat          Ap; /* dfdp */
198:   PetscMPIInt  size;
199:   PetscInt     n = 2;
200:   PetscScalar *u, *v;
201:   AppCtx       app;
202:   PetscInt     direction[1];
203:   PetscBool    terminate[1];
204:   Vec          lambda[2], mu[2];
205:   PetscReal    tend;

207:   FILE *f;
208:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209:      Initialize program
210:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211:   PetscFunctionBeginUser;
212:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
213:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
214:   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
215:   app.mode    = 1;
216:   app.lambda1 = 2.75;
217:   app.lambda2 = 0.36;
218:   tend        = 0.125;
219:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1adj options", "");
220:   {
221:     PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL));
222:     PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL));
223:     PetscCall(PetscOptionsReal("-tend", "", "", tend, &tend, NULL));
224:   }
225:   PetscOptionsEnd();

227:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228:     Create necessary matrix and vectors
229:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
230:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
231:   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
232:   PetscCall(MatSetType(A, MATDENSE));
233:   PetscCall(MatSetFromOptions(A));
234:   PetscCall(MatSetUp(A));

236:   PetscCall(MatCreateVecs(A, &U, NULL));

238:   PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap));
239:   PetscCall(MatSetSizes(Ap, n, 1, PETSC_DETERMINE, PETSC_DETERMINE));
240:   PetscCall(MatSetType(Ap, MATDENSE));
241:   PetscCall(MatSetFromOptions(Ap));
242:   PetscCall(MatSetUp(Ap));
243:   PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */

245:   PetscCall(VecGetArray(U, &u));
246:   u[0] = 0;
247:   u[1] = 1;
248:   PetscCall(VecRestoreArray(U, &u));
249:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
250:      Create timestepping solver context
251:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
252:   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
253:   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
254:   PetscCall(TSSetType(ts, TSCN));
255:   PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, &app));
256:   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobianFn *)IJacobian, &app));
257:   PetscCall(TSSetRHSJacobianP(ts, Ap, RHSJacobianP, &app));

259:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260:      Set initial conditions
261:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
262:   PetscCall(TSSetSolution(ts, U));

264:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
265:     Save trajectory of solution so that TSAdjointSolve() may be used
266:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
267:   PetscCall(TSSetSaveTrajectory(ts));

269:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
270:      Set solver options
271:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
272:   PetscCall(TSSetMaxTime(ts, tend));
273:   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
274:   PetscCall(TSSetTimeStep(ts, 1. / 256.));
275:   PetscCall(TSSetFromOptions(ts));

277:   /* Set directions and terminate flags for the two events */
278:   direction[0] = 0;
279:   terminate[0] = PETSC_FALSE;
280:   PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));

282:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283:      Run timestepping solver
284:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
285:   PetscCall(TSSolve(ts, U));

287:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
288:      Adjoint model starts here
289:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
290:   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
291:   PetscCall(MatCreateVecs(A, &lambda[1], NULL));
292:   /*   Set initial conditions for the adjoint integration */
293:   PetscCall(VecZeroEntries(lambda[0]));
294:   PetscCall(VecZeroEntries(lambda[1]));
295:   PetscCall(VecGetArray(lambda[0], &u));
296:   u[0] = 1.;
297:   PetscCall(VecRestoreArray(lambda[0], &u));
298:   PetscCall(VecGetArray(lambda[1], &u));
299:   u[1] = 1.;
300:   PetscCall(VecRestoreArray(lambda[1], &u));

302:   PetscCall(MatCreateVecs(Ap, &mu[0], NULL));
303:   PetscCall(MatCreateVecs(Ap, &mu[1], NULL));
304:   PetscCall(VecZeroEntries(mu[0]));
305:   PetscCall(VecZeroEntries(mu[1]));
306:   PetscCall(TSSetCostGradients(ts, 2, lambda, mu));

308:   PetscCall(TSAdjointSolve(ts));

310:   /*
311:   PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD));
312:   PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD));
313:   PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD));
314:   PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD));
315:   */
316:   PetscCall(VecGetArray(mu[0], &u));
317:   PetscCall(VecGetArray(mu[1], &v));
318:   f = fopen("adj_mu.out", "a");
319:   PetscCall(PetscFPrintf(PETSC_COMM_WORLD, f, "%20.15lf %20.15lf %20.15lf\n", (double)tend, (double)PetscRealPart(u[0]), (double)PetscRealPart(v[0])));
320:   PetscCall(VecRestoreArray(mu[0], &u));
321:   PetscCall(VecRestoreArray(mu[1], &v));
322:   fclose(f);
323:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
324:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
325:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
326:   PetscCall(MatDestroy(&A));
327:   PetscCall(VecDestroy(&U));
328:   PetscCall(TSDestroy(&ts));

330:   PetscCall(MatDestroy(&Ap));
331:   PetscCall(VecDestroy(&lambda[0]));
332:   PetscCall(VecDestroy(&lambda[1]));
333:   PetscCall(VecDestroy(&mu[0]));
334:   PetscCall(VecDestroy(&mu[1]));
335:   PetscCall(PetscFinalize());
336:   return 0;
337: }

339: /*TEST

341:    build:
342:       requires: !complex

344:    test:
345:       args: -ts_monitor -ts_adjoint_monitor

347: TEST*/