Actual source code: ex3sa.c

petsc-3.9.0 2018-04-07
Report Typos and Errors

  2: static char help[] = "Adjoint and tangent linear 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 sensitivity analysis 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 detected with TSEvent.
 20:  */

 22: #include <petscts.h>

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

 30: typedef enum {SA_ADJ, SA_TLM} SAMethod;
 31: static const char *const SAMethods[] = {"ADJ","TLM","SAMethod","SA_",0};

 33: /* Event check */
 34: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
 35: {
 36:   AppCtx        *user=(AppCtx*)ctx;

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

 44:   return(0);
 45: }

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

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

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

 77:   return(0);
 78: }

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

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

 98:   VecRestoreArrayRead(U,&u);
 99:   VecRestoreArrayRead(Udot,&udot);
100:   VecRestoreArray(F,&f);
101:   return(0);
102: }

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

115:   VecGetArrayRead(U,&u);
116:   VecGetArrayRead(Udot,&udot);
117:   Pmax = ctx->Pmax;

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

122:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
123:   VecRestoreArrayRead(U,&u);
124:   VecRestoreArrayRead(Udot,&udot);

126:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
127:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
128:   if (A != B) {
129:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
130:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
131:   }
132:   return(0);
133: }

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

142:   VecGetArray(X,&x);

144:   J[0][0] = 0;
145:   J[1][0] = 1.;
146:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

148:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
149:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
150:   return(0);
151: }

153: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
154: {
155:   PetscErrorCode    ierr;
156:   PetscScalar       *r;
157:   const PetscScalar *u;

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

168: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
169: {
170:   PetscErrorCode    ierr;
171:   PetscScalar       *ry;
172:   const PetscScalar *u;

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

183: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
184: {
185:   PetscErrorCode    ierr;
186:   PetscScalar       *rp;
187:   const PetscScalar *u;

190:   VecGetArrayRead(U,&u);
191:   VecGetArray(drdp[0],&rp);
192:   rp[0] = 0.;
193:   VecRestoreArray(drdp[0],&rp);
194:   VecGetArrayRead(U,&u);
195:   return(0);
196: }

198: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,PetscScalar *val,AppCtx *ctx)
199: {
200:   PetscErrorCode    ierr;
201:   const PetscScalar *x,*y;

204:   VecGetArrayRead(lambda,&x);
205:   VecGetArrayRead(mu,&y);
206:   val[0] = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
207:   VecRestoreArrayRead(lambda,&x);
208:   VecRestoreArrayRead(mu,&y);
209:   return(0);
210: }

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

235:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236:      Initialize program
237:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
239:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
240:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

242:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
243:     Create necessary matrix and vectors
244:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
245:   MatCreate(PETSC_COMM_WORLD,&A);
246:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
247:   MatSetType(A,MATDENSE);
248:   MatSetFromOptions(A);
249:   MatSetUp(A);

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

253:   MatCreate(PETSC_COMM_WORLD,&Jacp);
254:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
255:   MatSetFromOptions(Jacp);
256:   MatSetUp(Jacp);

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

290:     VecGetArray(U,&u);
291:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
292:     u[1] = 1.0;
293:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
294:     n    = 2;
295:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
296:     u[0] += du[0];
297:     u[1] += du[1];
298:     VecRestoreArray(U,&u);
299:     if (flg1 || flg2) {
300:       ctx.tf      = -1;
301:       ctx.tcl     = -1;
302:     }
303:     sa = SA_ADJ;
304:     PetscOptionsEnum("-sa_method","Sensitivity analysis method (adj or tlm)","",SAMethods,(PetscEnum)sa,(PetscEnum*)&sa,NULL);
305:   }
306:   PetscOptionsEnd();

308:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
309:      Create timestepping solver context
310:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
311:   TSCreate(PETSC_COMM_WORLD,&ts);
312:   TSSetProblemType(ts,TS_NONLINEAR);
313:   TSSetType(ts,TSBEULER);
314:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
315:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

317:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
318:      Set initial conditions
319:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
320:   TSSetSolution(ts,U);

322:   /*   Set RHS JacobianP */
323:   TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&ctx);
324:   if (sa == SA_ADJ) {
325:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
326:       Save trajectory of solution so that TSAdjointSolve() may be used
327:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
328:     TSSetSaveTrajectory(ts);

330:     MatCreateVecs(A,&lambda[0],NULL);
331:     MatCreateVecs(Jacp,&mu[0],NULL);
332:     TSSetCostGradients(ts,1,lambda,mu);
333:     TSSetCostIntegrand(ts,1,NULL,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
334:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
335:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_TRUE,&ctx);
336:   }

338:   if (sa == SA_TLM) {
339:     PetscScalar val[2];
340:     PetscInt    row[]={0,1},col[]={0};

342:     VecCreate(PETSC_COMM_WORLD,&qgrad[0]);
343:     VecSetSizes(qgrad[0],PETSC_DECIDE,1);
344:     VecSetFromOptions(qgrad[0]);

346:     MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&sp);
347:     TSForwardSetSensitivities(ts,1,sp);
348:     TSForwardSetIntegralGradients(ts,1,qgrad);
349:     TSSetCostIntegrand(ts,1,NULL,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
350:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
351:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_TRUE,&ctx);
352:     val[0] = 1./PetscSqrtScalar(1.-(ctx.Pm/ctx.Pmax)*(ctx.Pm/ctx.Pmax))/ctx.Pmax;
353:     val[1] = 0.0;
354:     MatSetValues(sp,2,row,1,col,val,INSERT_VALUES);
355:     MatAssemblyBegin(sp,MAT_FINAL_ASSEMBLY);
356:     MatAssemblyEnd(sp,MAT_FINAL_ASSEMBLY);
357:   }

359:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
360:      Set solver options
361:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
362:   TSSetMaxTime(ts,1.0);
363:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
364:   TSSetTimeStep(ts,0.03125);
365:   TSSetFromOptions(ts);

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

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

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

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

402:     VecGetArray(mu[0],&x_ptr);
403:     x_ptr[0] = 0.0;
404:     VecRestoreArray(mu[0],&x_ptr);

406:     TSAdjointSolve(ts);

408:     PetscPrintf(PETSC_COMM_WORLD,"\n lambda: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n");
409:     VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
410:     PetscPrintf(PETSC_COMM_WORLD,"\n mu: d[Psi(tf)]/d[pm]  d[Psi(tf)]/d[pm]\n");
411:     VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
412:     TSGetCostIntegral(ts,&q);
413:     VecGetArray(q,&x_ptr);
414:     PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));
415:     VecRestoreArray(q,&x_ptr);
416:     ComputeSensiP(lambda[0],mu[0],grad,&ctx);
417:     PetscPrintf(PETSC_COMM_WORLD,"\n gradient=%g\n",(double)grad[0]);
418:     VecDestroy(&lambda[0]);
419:     VecDestroy(&mu[0]);
420:   }

422:   if (sa == SA_TLM) {
423:     PetscPrintf(PETSC_COMM_WORLD,"\n trajectory sensitivity: d[Psi(tf)]/d[pm]  d[Psi(tf)]/d[pm]\n");
424:     MatView(sp,PETSC_VIEWER_STDOUT_WORLD);
425:     TSGetCostIntegral(ts,&q);
426:     VecGetArray(q,&s_ptr);
427:     PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(s_ptr[0]-ctx.Pm));
428:     VecRestoreArray(q,&s_ptr);
429:     VecGetArray(qgrad[0],&s_ptr);
430:     PetscPrintf(PETSC_COMM_WORLD,"\n gradient=%g\n",(double)s_ptr[0]);
431:     VecRestoreArray(qgrad[0],&s_ptr);
432:     VecDestroy(&qgrad[0]);
433:     MatDestroy(&sp);
434:   }
435:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
436:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
437:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
438:   MatDestroy(&A);
439:   MatDestroy(&Jacp);
440:   VecDestroy(&U);
441:   TSDestroy(&ts);
442:   PetscFinalize();
443:   return ierr;
444: }


447: /*TEST

449:    build:
450:       requires: !complex !single

452:    test:
453:       args: -sa_method adj -viewer_binary_skip_info -ts_type cn -pc_type lu

455:    test:
456:       suffix: 2
457:       args: -sa_method tlm -ts_type cn -pc_type lu

459: TEST*/