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