Actual source code: ex3adj.c

petsc-3.7.2 2016-06-05
Report Typos and Errors
  2: static char help[] = "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 */
 32: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
 33: {
 34:   AppCtx        *user=(AppCtx*)ctx;

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

 42:   return(0);
 43: }

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


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

 65: PetscErrorCode PostStepFunction(TS ts)
 66: {
 67:   PetscErrorCode    ierr;
 68:   Vec               U;
 69:   PetscReal         t;
 70:   const PetscScalar *u;

 73:   TSGetTime(ts,&t);
 74:   TSGetSolution(ts,&U);
 75:   VecGetArrayRead(U,&u);
 76:   PetscPrintf(PETSC_COMM_SELF,"delta(%3.2f) = %8.7f\n",(double)t,(double)u[0]);
 77:   VecRestoreArrayRead(U,&u);

 79:   return(0);
 80: }

 84: /*
 85:      Defines the ODE passed to the ODE solver
 86: */
 87: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 88: {
 89:   PetscErrorCode    ierr;
 90:   PetscScalar       *f,Pmax;
 91:   const PetscScalar *u,*udot;

 94:   /*  The next three lines allow us to access the entries of the vectors directly */
 95:   VecGetArrayRead(U,&u);
 96:   VecGetArrayRead(Udot,&udot);
 97:   VecGetArray(F,&f);
 98:   Pmax = ctx->Pmax;
 99:   f[0] = udot[0] - ctx->omega_b*(u[1] - ctx->omega_s);
100:   f[1] = 2.0*ctx->H/ctx->omega_s*udot[1] +  Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - ctx->omega_s)- ctx->Pm;

102:   VecRestoreArrayRead(U,&u);
103:   VecRestoreArrayRead(Udot,&udot);
104:   VecRestoreArray(F,&f);
105:   return(0);
106: }

110: /*
111:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
112: */
113: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
114: {
115:   PetscErrorCode    ierr;
116:   PetscInt          rowcol[] = {0,1};
117:   PetscScalar       J[2][2],Pmax;
118:   const PetscScalar *u,*udot;

121:   VecGetArrayRead(U,&u);
122:   VecGetArrayRead(Udot,&udot);
123:   Pmax = ctx->Pmax;

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

128:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
129:   VecRestoreArrayRead(U,&u);
130:   VecRestoreArrayRead(Udot,&udot);

132:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
133:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
134:   if (A != B) {
135:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
136:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
137:   }
138:   return(0);
139: }

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

150:   VecGetArray(X,&x);

152:   J[0][0] = 0;
153:   J[1][0] = 1.;
154:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

156:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
157:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
158:   return(0);
159: }

163: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
164: {
165:   PetscErrorCode    ierr;
166:   PetscScalar       *r;
167:   const PetscScalar *u;

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

180: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
181: {
182:   PetscErrorCode    ierr;
183:   PetscScalar       *ry;
184:   const PetscScalar *u;

187:   VecGetArrayRead(U,&u);
188:   VecGetArray(drdy[0],&ry);
189:   ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
190:   VecRestoreArray(drdy[0],&ry);
191:   VecRestoreArrayRead(U,&u);
192:   return(0);
193: }

197: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
198: {
199:   PetscErrorCode    ierr;
200:   PetscScalar       *rp;
201:   const PetscScalar *u;

204:   VecGetArrayRead(U,&u);
205:   VecGetArray(drdp[0],&rp);
206:   rp[0] = 0.;
207:   VecRestoreArray(drdp[0],&rp);
208:   VecGetArrayRead(U,&u);
209:   return(0);
210: }

214: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
215: {
216:   PetscErrorCode    ierr;
217:   PetscScalar       sensip;
218:   const PetscScalar *x,*y;

221:   VecGetArrayRead(lambda,&x);
222:   VecGetArrayRead(mu,&y);
223:   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
224:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip);
225:   VecRestoreArrayRead(lambda,&x);
226:   VecRestoreArrayRead(mu,&y);
227:   return(0);
228: }

232: int main(int argc,char **argv)
233: {
234:   TS             ts;            /* ODE integrator */
235:   Vec            U;             /* solution will be stored here */
236:   Mat            A;             /* Jacobian matrix */
237:   Mat            Jacp;          /* Jacobian matrix */
239:   PetscMPIInt    size;
240:   PetscInt       n = 2;
241:   AppCtx         ctx;
242:   PetscScalar    *u;
243:   PetscReal      du[2] = {0.0,0.0};
244:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
245:   PetscReal      ftime;
246:   PetscInt       steps;
247:   PetscScalar    *x_ptr,*y_ptr;
248:   Vec            lambda[1],q,mu[1];
249:   PetscInt       direction[2];
250:   PetscBool      terminate[2];

252:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
253:      Initialize program
254:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
255:   PetscInitialize(&argc,&argv,(char*)0,help);
256:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
257:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

259:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260:     Create necessary matrix and vectors
261:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
262:   MatCreate(PETSC_COMM_WORLD,&A);
263:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
264:   MatSetType(A,MATDENSE);
265:   MatSetFromOptions(A);
266:   MatSetUp(A);

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

270:   MatCreate(PETSC_COMM_WORLD,&Jacp);
271:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
272:   MatSetFromOptions(Jacp);
273:   MatSetUp(Jacp);

275:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276:     Set runtime options
277:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
278:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
279:   {
280:     ctx.beta    = 2;
281:     ctx.c       = 10000.0;
282:     ctx.u_s     = 1.0;
283:     ctx.omega_s = 1.0;
284:     ctx.omega_b = 120.0*PETSC_PI;
285:     ctx.H       = 5.0;
286:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
287:     ctx.D       = 5.0;
288:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
289:     ctx.E       = 1.1378;
290:     ctx.V       = 1.0;
291:     ctx.X       = 0.545;
292:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
293:     ctx.Pmax_ini = ctx.Pmax;
294:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
295:     ctx.Pm      = 1.1;
296:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
297:     ctx.tf      = 0.1;
298:     ctx.tcl     = 0.2;
299:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
300:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
301:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
302:     if (ensemble) {
303:       ctx.tf      = -1;
304:       ctx.tcl     = -1;
305:     }

307:     VecGetArray(U,&u);
308:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
309:     u[1] = 1.0;
310:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
311:     n    = 2;
312:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
313:     u[0] += du[0];
314:     u[1] += du[1];
315:     VecRestoreArray(U,&u);
316:     if (flg1 || flg2) {
317:       ctx.tf      = -1;
318:       ctx.tcl     = -1;
319:     }
320:   }
321:   PetscOptionsEnd();

323:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
324:      Create timestepping solver context
325:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
326:   TSCreate(PETSC_COMM_WORLD,&ts);
327:   TSSetProblemType(ts,TS_NONLINEAR);
328:   TSSetType(ts,TSBEULER);
329:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
330:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);
331:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

333:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
334:      Set initial conditions
335:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
336:   TSSetSolution(ts,U);

338:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
339:     Save trajectory of solution so that TSAdjointSolve() may be used
340:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
341:   TSSetSaveTrajectory(ts);

343:   MatCreateVecs(A,&lambda[0],NULL);
344:   /*   Set initial conditions for the adjoint integration */
345:   VecGetArray(lambda[0],&y_ptr);
346:   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
347:   VecRestoreArray(lambda[0],&y_ptr);

349:   MatCreateVecs(Jacp,&mu[0],NULL);
350:   VecGetArray(mu[0],&x_ptr);
351:   x_ptr[0] = -1.0;
352:   VecRestoreArray(mu[0],&x_ptr);
353:   TSSetCostGradients(ts,1,lambda,mu);
354:   TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
355:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
356:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_TRUE,&ctx);

358:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
359:      Set solver options
360:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
361:   TSSetDuration(ts,PETSC_DEFAULT,1.0);
362:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
363:   TSSetInitialTimeStep(ts,0.0,.005);
364:   TSSetFromOptions(ts);

366:   direction[0] = direction[1] = 1;
367:   terminate[0] = terminate[1] = PETSC_FALSE;

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

371:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
372:      Solve nonlinear system
373:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
374:   if (ensemble) {
375:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
376:       VecGetArray(U,&u);
377:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
378:       u[1] = ctx.omega_s;
379:       u[0] += du[0];
380:       u[1] += du[1];
381:       VecRestoreArray(U,&u);
382:       TSSetInitialTimeStep(ts,0.0,.005);
383:       TSSolve(ts,U);
384:     }
385:   } else {
386:     TSSolve(ts,U);
387:   }
388:   VecView(U,PETSC_VIEWER_STDOUT_WORLD);
389:   TSGetSolveTime(ts,&ftime);
390:   TSGetTimeStepNumber(ts,&steps);

392:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
393:      Adjoint model starts here
394:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
395:   /*   Set initial conditions for the adjoint integration */
396:   VecGetArray(lambda[0],&y_ptr);
397:   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
398:   VecRestoreArray(lambda[0],&y_ptr);

400:   VecGetArray(mu[0],&x_ptr);
401:   x_ptr[0] = -1.0;
402:   VecRestoreArray(mu[0],&x_ptr);

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

407:   TSAdjointSolve(ts);

409:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n");
410:   VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
411:   VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
412:   TSGetCostIntegral(ts,&q);
413:   VecView(q,PETSC_VIEWER_STDOUT_WORLD);
414:   VecGetArray(q,&x_ptr);
415:   PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));
416:   VecRestoreArray(q,&x_ptr);

418:   ComputeSensiP(lambda[0],mu[0],&ctx);

420:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
421:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
422:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
423:   MatDestroy(&A);
424:   MatDestroy(&Jacp);
425:   VecDestroy(&U);
426:   VecDestroy(&lambda[0]);
427:   VecDestroy(&mu[0]);
428:   TSDestroy(&ts);

430:   PetscFinalize();
431:   return(0);
432: }