Actual source code: ex3adj.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2: static char help[] = "Adjoint sensitivity analysis of the basic equation for generator stability analysis.\n";


\begin{eqnarray}
\frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
\frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
\end{eqnarray}

 13: /*
 14:   This code demonstrate the TSAdjoint interface to a system of ordinary differential equations with discontinuities.
 15:   It computes the sensitivities of an integral cost function
 16:   \int c*max(0,\theta(t)-u_s)^beta dt
 17:   w.r.t. initial conditions and the parameter P_m.
 18:   Backward Euler method is used for time integration.
 19:   The discontinuities are dealt with TSEvent, which is compatible with TSAdjoint.
 20:  */
 21: #include <petscts.h>

 23: typedef struct {
 24:   PetscScalar H,D,omega_b,omega_s,Pmax,Pmax_ini,Pm,E,V,X,u_s,c;
 25:   PetscInt    beta;
 26:   PetscReal   tf,tcl;
 27: } AppCtx;

 29: /* Event check */
 30: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
 31: {
 32:   AppCtx        *user=(AppCtx*)ctx;

 35:   /* Event for fault-on time */
 36:   fvalue[0] = t - user->tf;
 37:   /* Event for fault-off time */
 38:   fvalue[1] = t - user->tcl;

 40:   return(0);
 41: }

 43: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx)
 44: {
 45:   AppCtx *user=(AppCtx*)ctx;


 49:   if (event_list[0] == 0) {
 50:     if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */
 51:     else user->Pmax = user->Pmax_ini; /* Going backward, reversal of event */
 52:   } else if(event_list[0] == 1) {
 53:     if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault  - this is done by setting Pmax = Pmax_ini */
 54:     else user->Pmax = 0.0; /* Going backward, reversal of event */
 55:   }
 56:   return(0);
 57: }

 59: PetscErrorCode PostStepFunction(TS ts)
 60: {
 61:   PetscErrorCode    ierr;
 62:   Vec               U;
 63:   PetscReal         t;
 64:   const PetscScalar *u;

 67:   TSGetTime(ts,&t);
 68:   TSGetSolution(ts,&U);
 69:   VecGetArrayRead(U,&u);
 70:   PetscPrintf(PETSC_COMM_SELF,"delta(%3.2f) = %8.7f\n",(double)t,(double)u[0]);
 71:   VecRestoreArrayRead(U,&u);

 73:   return(0);
 74: }

 76: /*
 77:      Defines the ODE passed to the ODE solver
 78: */
 79: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 80: {
 81:   PetscErrorCode    ierr;
 82:   PetscScalar       *f,Pmax;
 83:   const PetscScalar *u,*udot;

 86:   /*  The next three lines allow us to access the entries of the vectors directly */
 87:   VecGetArrayRead(U,&u);
 88:   VecGetArrayRead(Udot,&udot);
 89:   VecGetArray(F,&f);
 90:   Pmax = ctx->Pmax;
 91:   f[0] = udot[0] - ctx->omega_b*(u[1] - ctx->omega_s);
 92:   f[1] = 2.0*ctx->H/ctx->omega_s*udot[1] +  Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - ctx->omega_s)- ctx->Pm;

 94:   VecRestoreArrayRead(U,&u);
 95:   VecRestoreArrayRead(Udot,&udot);
 96:   VecRestoreArray(F,&f);
 97:   return(0);
 98: }

100: /*
101:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
102: */
103: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
104: {
105:   PetscErrorCode    ierr;
106:   PetscInt          rowcol[] = {0,1};
107:   PetscScalar       J[2][2],Pmax;
108:   const PetscScalar *u,*udot;

111:   VecGetArrayRead(U,&u);
112:   VecGetArrayRead(Udot,&udot);
113:   Pmax = ctx->Pmax;

115:   J[0][0] = a;                       J[0][1] = -ctx->omega_b;
116:   J[1][1] = 2.0*ctx->H/ctx->omega_s*a + ctx->D;   J[1][0] = Pmax*PetscCosScalar(u[0]);

118:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
119:   VecRestoreArrayRead(U,&u);
120:   VecRestoreArrayRead(Udot,&udot);

122:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
123:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
124:   if (A != B) {
125:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
126:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
127:   }
128:   return(0);
129: }

131: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
132: {
134:   PetscInt       row[] = {0,1},col[]={0};
135:   PetscScalar    *x,J[2][1];

138:   VecGetArray(X,&x);

140:   J[0][0] = 0;
141:   J[1][0] = 1.;
142:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

144:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
145:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
146:   return(0);
147: }

149: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
150: {
151:   PetscErrorCode    ierr;
152:   PetscScalar       *r;
153:   const PetscScalar *u;

156:   VecGetArrayRead(U,&u);
157:   VecGetArray(R,&r);
158:   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);
159:   VecRestoreArray(R,&r);
160:   VecRestoreArrayRead(U,&u);
161:   return(0);
162: }

164: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
165: {
166:   PetscErrorCode    ierr;
167:   PetscScalar       *ry;
168:   const PetscScalar *u;

171:   VecGetArrayRead(U,&u);
172:   VecGetArray(drdy[0],&ry);
173:   ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
174:   VecRestoreArray(drdy[0],&ry);
175:   VecRestoreArrayRead(U,&u);
176:   return(0);
177: }

179: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
180: {
181:   PetscErrorCode    ierr;
182:   PetscScalar       *rp;
183:   const PetscScalar *u;

186:   VecGetArrayRead(U,&u);
187:   VecGetArray(drdp[0],&rp);
188:   rp[0] = 0.;
189:   VecRestoreArray(drdp[0],&rp);
190:   VecGetArrayRead(U,&u);
191:   return(0);
192: }

194: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
195: {
196:   PetscErrorCode    ierr;
197:   PetscScalar       sensip;
198:   const PetscScalar *x,*y;

201:   VecGetArrayRead(lambda,&x);
202:   VecGetArrayRead(mu,&y);
203:   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
204:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip);
205:   VecRestoreArrayRead(lambda,&x);
206:   VecRestoreArrayRead(mu,&y);
207:   return(0);
208: }

210: int main(int argc,char **argv)
211: {
212:   TS             ts;            /* ODE integrator */
213:   Vec            U;             /* solution will be stored here */
214:   Mat            A;             /* Jacobian matrix */
215:   Mat            Jacp;          /* Jacobian matrix */
217:   PetscMPIInt    size;
218:   PetscInt       n = 2;
219:   AppCtx         ctx;
220:   PetscScalar    *u;
221:   PetscReal      du[2] = {0.0,0.0};
222:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
223:   PetscReal      ftime;
224:   PetscInt       steps;
225:   PetscScalar    *x_ptr,*y_ptr;
226:   Vec            lambda[1],q,mu[1];
227:   PetscInt       direction[2];
228:   PetscBool      terminate[2];

230:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231:      Initialize program
232:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
233:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
234:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
235:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

237:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238:     Create necessary matrix and vectors
239:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
240:   MatCreate(PETSC_COMM_WORLD,&A);
241:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
242:   MatSetType(A,MATDENSE);
243:   MatSetFromOptions(A);
244:   MatSetUp(A);

246:   MatCreateVecs(A,&U,NULL);

248:   MatCreate(PETSC_COMM_WORLD,&Jacp);
249:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
250:   MatSetFromOptions(Jacp);
251:   MatSetUp(Jacp);

253:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
254:     Set runtime options
255:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
256:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
257:   {
258:     ctx.beta    = 2;
259:     ctx.c       = 10000.0;
260:     ctx.u_s     = 1.0;
261:     ctx.omega_s = 1.0;
262:     ctx.omega_b = 120.0*PETSC_PI;
263:     ctx.H       = 5.0;
264:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
265:     ctx.D       = 5.0;
266:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
267:     ctx.E       = 1.1378;
268:     ctx.V       = 1.0;
269:     ctx.X       = 0.545;
270:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
271:     ctx.Pmax_ini = ctx.Pmax;
272:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
273:     ctx.Pm      = 1.1;
274:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
275:     ctx.tf      = 0.1;
276:     ctx.tcl     = 0.2;
277:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
278:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
279:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
280:     if (ensemble) {
281:       ctx.tf      = -1;
282:       ctx.tcl     = -1;
283:     }

285:     VecGetArray(U,&u);
286:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
287:     u[1] = 1.0;
288:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
289:     n    = 2;
290:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
291:     u[0] += du[0];
292:     u[1] += du[1];
293:     VecRestoreArray(U,&u);
294:     if (flg1 || flg2) {
295:       ctx.tf      = -1;
296:       ctx.tcl     = -1;
297:     }
298:   }
299:   PetscOptionsEnd();

301:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
302:      Create timestepping solver context
303:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
304:   TSCreate(PETSC_COMM_WORLD,&ts);
305:   TSSetProblemType(ts,TS_NONLINEAR);
306:   TSSetType(ts,TSBEULER);
307:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
308:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

310:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
311:      Set initial conditions
312:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
313:   TSSetSolution(ts,U);

315:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
316:     Save trajectory of solution so that TSAdjointSolve() may be used
317:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
318:   TSSetSaveTrajectory(ts);

320:   MatCreateVecs(A,&lambda[0],NULL);
321:   MatCreateVecs(Jacp,&mu[0],NULL);
322:   TSSetCostGradients(ts,1,lambda,mu);
323:   TSSetCostIntegrand(ts,1,NULL,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
324:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
325:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_TRUE,&ctx);

327:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
328:      Set solver options
329:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
330:   TSSetMaxTime(ts,1.0);
331:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
332:   TSSetTimeStep(ts,0.03125);
333:   TSSetFromOptions(ts);

335:   direction[0] = direction[1] = 1;
336:   terminate[0] = terminate[1] = PETSC_FALSE;

338:   TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&ctx);

340:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
341:      Solve nonlinear system
342:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
343:   if (ensemble) {
344:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
345:       VecGetArray(U,&u);
346:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
347:       u[1] = ctx.omega_s;
348:       u[0] += du[0];
349:       u[1] += du[1];
350:       VecRestoreArray(U,&u);
351:       TSSetTimeStep(ts,0.03125);
352:       TSSolve(ts,U);
353:     }
354:   } else {
355:     TSSolve(ts,U);
356:   }
357:   VecView(U,PETSC_VIEWER_STDOUT_WORLD);
358:   TSGetSolveTime(ts,&ftime);
359:   TSGetStepNumber(ts,&steps);

361:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
362:      Adjoint model starts here
363:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
364:   /*   Set initial conditions for the adjoint integration */
365:   VecGetArray(lambda[0],&y_ptr);
366:   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
367:   VecRestoreArray(lambda[0],&y_ptr);

369:   VecGetArray(mu[0],&x_ptr);
370:   x_ptr[0] = -1.0;
371:   VecRestoreArray(mu[0],&x_ptr);

373:   /*   Set RHS JacobianP */
374:   TSAdjointSetRHSJacobian(ts,Jacp,RHSJacobianP,&ctx);

376:   TSAdjointSolve(ts);

378:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n");
379:   VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
380:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameters (partial derivative): d[Psi(tf)]/d[pm]  d[Psi(tf)]/d[pm]\n");
381:   VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
382:   TSGetCostIntegral(ts,&q);
383:   VecGetArray(q,&x_ptr);
384:   PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));
385:   VecRestoreArray(q,&x_ptr);
386:   ComputeSensiP(lambda[0],mu[0],&ctx);

388:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
389:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
390:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
391:   MatDestroy(&A);
392:   MatDestroy(&Jacp);
393:   VecDestroy(&U);
394:   VecDestroy(&lambda[0]);
395:   VecDestroy(&mu[0]);
396:   TSDestroy(&ts);
397:   PetscFinalize();
398:   return ierr;
399: }