Actual source code: ex9adj.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: static char help[] = "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}



Ensemble of initial conditions
./ex9 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly

Fault at .1 seconds
./ex9 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly

Initial conditions same as when fault is ended
./ex9 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly


 25: /*
 26:    Include "petscts.h" so that we can use TS solvers.  Note that this
 27:    file automatically includes:
 28:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 29:      petscmat.h - matrices
 30:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 31:      petscviewer.h - viewers               petscpc.h  - preconditioners
 32:      petscksp.h   - linear solvers
 33: */
 34: #include <petscts.h>

 36: typedef struct {
 37:   PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c;
 38:   PetscInt    beta;
 39:   PetscReal   tf,tcl;
 40: } AppCtx;

 44: PetscErrorCode PostStepFunction(TS ts)
 45: {
 46:   PetscErrorCode    ierr;
 47:   Vec               U;
 48:   PetscReal         t;
 49:   const PetscScalar *u;

 52:   TSGetTime(ts,&t);
 53:   TSGetSolution(ts,&U);
 54:   VecGetArrayRead(U,&u);
 55:   PetscPrintf(PETSC_COMM_SELF,"delta(%3.2f) = %8.7f\n",(double)t,(double)u[0]);
 56:   VecRestoreArrayRead(U,&u);
 57: 
 58:   return(0);
 59: }

 63: /*
 64:      Defines the ODE passed to the ODE solver
 65: */
 66: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
 67: {
 68:   PetscErrorCode    ierr;
 69:   PetscScalar       *f,Pmax;
 70:   const PetscScalar *u;

 73:   /*  The next three lines allow us to access the entries of the vectors directly */
 74:   VecGetArrayRead(U,&u);
 75:   VecGetArray(F,&f);
 76:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
 77:   else Pmax = ctx->Pmax;
 78: 
 79:   f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
 80:   f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H);

 82:   VecRestoreArrayRead(U,&u);
 83:   VecRestoreArray(F,&f);
 84:   return(0);
 85: }

 89: /*
 90:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 91: */
 92: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
 93: {
 94:   PetscErrorCode    ierr;
 95:   PetscInt          rowcol[] = {0,1};
 96:   PetscScalar       J[2][2],Pmax;
 97:   const PetscScalar *u;

100:   VecGetArrayRead(U,&u);
101:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
102:   else Pmax = ctx->Pmax;

104:   J[0][0] = 0;                                  J[0][1] = ctx->omega_b;
105:   J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H);  J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H);

107:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
108:   VecRestoreArrayRead(U,&u);

110:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
111:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
112:   if (A != B) {
113:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
114:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
115:   }
116:   return(0);
117: }

121: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
122: {
124:   PetscInt       row[] = {0,1},col[]={0};
125:   PetscScalar    J[2][1];
126:   AppCtx         *ctx=(AppCtx*)ctx0;

129:   J[0][0] = 0;
130:   J[1][0] = ctx->omega_s/(2.0*ctx->H);
131:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
132:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
133:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
134:   return(0);
135: }

139: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
140: {
141:   PetscErrorCode    ierr;
142:   PetscScalar       *r;
143:   const PetscScalar *u;

146:   VecGetArrayRead(U,&u);
147:   VecGetArray(R,&r);
148:   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);
149:   VecRestoreArray(R,&r);
150:   VecRestoreArrayRead(U,&u);
151:   return(0);
152: }

156: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
157: {
158:   PetscErrorCode    ierr;
159:   PetscScalar       *ry;
160:   const PetscScalar *u;

163:   VecGetArrayRead(U,&u);
164:   VecGetArray(drdy[0],&ry);
165:   ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
166:   VecRestoreArray(drdy[0],&ry);
167:   VecRestoreArrayRead(U,&u);
168:   return(0);
169: }

173: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
174: {
175:   PetscErrorCode    ierr;
176:   PetscScalar       *rp;
177:   const PetscScalar *u;

180:   VecGetArrayRead(U,&u);
181:   VecGetArray(drdp[0],&rp);
182:   rp[0] = 0.;
183:   VecRestoreArray(drdp[0],&rp);
184:   VecGetArrayRead(U,&u);
185:   return(0);
186: }

190: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
191: {
192:   PetscErrorCode    ierr;
193:   PetscScalar       sensip;
194:   const PetscScalar *x,*y;
195: 
197:   VecGetArrayRead(lambda,&x);
198:   VecGetArrayRead(mu,&y);
199:   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
200:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip);
201:   VecRestoreArrayRead(lambda,&x);
202:   VecRestoreArrayRead(mu,&y);
203:   return(0);
204: }

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

226:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227:      Initialize program
228:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229:   PetscInitialize(&argc,&argv,(char*)0,help);
230:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
231:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

233:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234:     Create necessary matrix and vectors
235:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236:   MatCreate(PETSC_COMM_WORLD,&A);
237:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
238:   MatSetType(A,MATDENSE);
239:   MatSetFromOptions(A);
240:   MatSetUp(A);

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

244:   MatCreate(PETSC_COMM_WORLD,&Jacp);
245:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
246:   MatSetFromOptions(Jacp);
247:   MatSetUp(Jacp);

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

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

296:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
297:      Create timestepping solver context
298:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
299:   TSCreate(PETSC_COMM_WORLD,&ts);
300:   TSSetProblemType(ts,TS_NONLINEAR);
301:   TSSetType(ts,TSRK);
302:   TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);
303:   TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);
304:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

306:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
307:      Set initial conditions
308:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
309:   TSSetSolution(ts,U);

311:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
312:     Save trajectory of solution so that TSAdjointSolve() may be used
313:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
314:   TSSetSaveTrajectory(ts);

316:   MatCreateVecs(A,&lambda[0],NULL);
317:   /*   Set initial conditions for the adjoint integration */
318:   VecGetArray(lambda[0],&y_ptr);
319:   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
320:   VecRestoreArray(lambda[0],&y_ptr);

322:   MatCreateVecs(Jacp,&mu[0],NULL);
323:   VecGetArray(mu[0],&x_ptr);
324:   x_ptr[0] = -1.0;
325:   VecRestoreArray(mu[0],&x_ptr);
326:   TSSetCostGradients(ts,1,lambda,mu);
327:   TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
328:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
329:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_TRUE,&ctx);

331:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
332:      Set solver options
333:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
334:   TSSetDuration(ts,PETSC_DEFAULT,10.0);
335:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
336:   TSSetInitialTimeStep(ts,0.0,.01);
337:   TSSetFromOptions(ts);

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

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

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

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

375:   TSAdjointSolve(ts);

377:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n");
378:   VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
379:   VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
380:   TSGetCostIntegral(ts,&q);
381:   VecView(q,PETSC_VIEWER_STDOUT_WORLD);
382:   VecGetArray(q,&x_ptr);
383:   PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));
384:   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(0);
399: }