Actual source code: ex9adj.c

petsc-3.8.3 2017-12-09
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;

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

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

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

 69:   /*  The next three lines allow us to access the entries of the vectors directly */
 70:   VecGetArrayRead(U,&u);
 71:   VecGetArray(F,&f);
 72:   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 */
 73:   else Pmax = ctx->Pmax;
 74: 
 75:   f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
 76:   f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H);

 78:   VecRestoreArrayRead(U,&u);
 79:   VecRestoreArray(F,&f);
 80:   return(0);
 81: }

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

 94:   VecGetArrayRead(U,&u);
 95:   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 */
 96:   else Pmax = ctx->Pmax;

 98:   J[0][0] = 0;                                  J[0][1] = ctx->omega_b;
 99:   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);

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

104:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
105:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
106:   if (A != B) {
107:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
108:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
109:   }
110:   return(0);
111: }

113: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
114: {
116:   PetscInt       row[] = {0,1},col[]={0};
117:   PetscScalar    J[2][1];
118:   AppCtx         *ctx=(AppCtx*)ctx0;

121:   J[0][0] = 0;
122:   J[1][0] = ctx->omega_s/(2.0*ctx->H);
123:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
124:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
125:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
126:   return(0);
127: }

129: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
130: {
131:   PetscErrorCode    ierr;
132:   PetscScalar       *r;
133:   const PetscScalar *u;

136:   VecGetArrayRead(U,&u);
137:   VecGetArray(R,&r);
138:   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);
139:   VecRestoreArray(R,&r);
140:   VecRestoreArrayRead(U,&u);
141:   return(0);
142: }

144: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
145: {
146:   PetscErrorCode    ierr;
147:   PetscScalar       *ry;
148:   const PetscScalar *u;

151:   VecGetArrayRead(U,&u);
152:   VecGetArray(drdy[0],&ry);
153:   ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
154:   VecRestoreArray(drdy[0],&ry);
155:   VecRestoreArrayRead(U,&u);
156:   return(0);
157: }

159: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
160: {
161:   PetscErrorCode    ierr;
162:   PetscScalar       *rp;
163:   const PetscScalar *u;

166:   VecGetArrayRead(U,&u);
167:   VecGetArray(drdp[0],&rp);
168:   rp[0] = 0.;
169:   VecRestoreArray(drdp[0],&rp);
170:   VecGetArrayRead(U,&u);
171:   return(0);
172: }

174: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
175: {
176:   PetscErrorCode    ierr;
177:   PetscScalar       sensip;
178:   const PetscScalar *x,*y;
179: 
181:   VecGetArrayRead(lambda,&x);
182:   VecGetArrayRead(mu,&y);
183:   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
184:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip);
185:   VecRestoreArrayRead(lambda,&x);
186:   VecRestoreArrayRead(mu,&y);
187:   return(0);
188: }

190: int main(int argc,char **argv)
191: {
192:   TS             ts;            /* ODE integrator */
193:   Vec            U;             /* solution will be stored here */
194:   Mat            A;             /* Jacobian matrix */
195:   Mat            Jacp;          /* Jacobian matrix */
197:   PetscMPIInt    size;
198:   PetscInt       n = 2;
199:   AppCtx         ctx;
200:   PetscScalar    *u;
201:   PetscReal      du[2] = {0.0,0.0};
202:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
203:   PetscReal      ftime;
204:   PetscInt       steps;
205:   PetscScalar    *x_ptr,*y_ptr;
206:   Vec            lambda[1],q,mu[1];

208:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209:      Initialize program
210:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
212:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
213:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

215:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216:     Create necessary matrix and vectors
217:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
218:   MatCreate(PETSC_COMM_WORLD,&A);
219:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
220:   MatSetType(A,MATDENSE);
221:   MatSetFromOptions(A);
222:   MatSetUp(A);

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

226:   MatCreate(PETSC_COMM_WORLD,&Jacp);
227:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
228:   MatSetFromOptions(Jacp);
229:   MatSetUp(Jacp);

231:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232:     Set runtime options
233:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
234:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
235:   {
236:     ctx.beta    = 2;
237:     ctx.c       = 10000.0;
238:     ctx.u_s     = 1.0;
239:     ctx.omega_s = 1.0;
240:     ctx.omega_b = 120.0*PETSC_PI;
241:     ctx.H       = 5.0;
242:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
243:     ctx.D       = 5.0;
244:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
245:     ctx.E       = 1.1378;
246:     ctx.V       = 1.0;
247:     ctx.X       = 0.545;
248:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;;
249:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
250:     ctx.Pm      = 1.1;
251:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
252:     ctx.tf      = 0.1;
253:     ctx.tcl     = 0.2;
254:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
255:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
256:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
257:     if (ensemble) {
258:       ctx.tf      = -1;
259:       ctx.tcl     = -1;
260:     }

262:     VecGetArray(U,&u);
263:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
264:     u[1] = 1.0;
265:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
266:     n    = 2;
267:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
268:     u[0] += du[0];
269:     u[1] += du[1];
270:     VecRestoreArray(U,&u);
271:     if (flg1 || flg2) {
272:       ctx.tf      = -1;
273:       ctx.tcl     = -1;
274:     }
275:   }
276:   PetscOptionsEnd();

278:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
279:      Create timestepping solver context
280:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
281:   TSCreate(PETSC_COMM_WORLD,&ts);
282:   TSSetProblemType(ts,TS_NONLINEAR);
283:   TSSetType(ts,TSRK);
284:   TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);
285:   TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);

287:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
288:      Set initial conditions
289:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
290:   TSSetSolution(ts,U);

292:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
293:     Save trajectory of solution so that TSAdjointSolve() may be used
294:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
295:   TSSetSaveTrajectory(ts);

297:   MatCreateVecs(A,&lambda[0],NULL);
298:   /*   Set initial conditions for the adjoint integration */
299:   VecGetArray(lambda[0],&y_ptr);
300:   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
301:   VecRestoreArray(lambda[0],&y_ptr);

303:   MatCreateVecs(Jacp,&mu[0],NULL);
304:   VecGetArray(mu[0],&x_ptr);
305:   x_ptr[0] = -1.0;
306:   VecRestoreArray(mu[0],&x_ptr);
307:   TSSetCostGradients(ts,1,lambda,mu);
308:   TSSetCostIntegrand(ts,1,NULL,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
309:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
310:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_TRUE,&ctx);

312:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
313:      Set solver options
314:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
315:   TSSetMaxTime(ts,10.0);
316:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
317:   TSSetTimeStep(ts,.01);
318:   TSSetFromOptions(ts);

320:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
321:      Solve nonlinear system
322:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
323:   if (ensemble) {
324:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
325:       VecGetArray(U,&u);
326:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
327:       u[1] = ctx.omega_s;
328:       u[0] += du[0];
329:       u[1] += du[1];
330:       VecRestoreArray(U,&u);
331:       TSSetTimeStep(ts,.01);
332:       TSSolve(ts,U);
333:     }
334:   } else {
335:     TSSolve(ts,U);
336:   }
337:   VecView(U,PETSC_VIEWER_STDOUT_WORLD);
338:   TSGetSolveTime(ts,&ftime);
339:   TSGetStepNumber(ts,&steps);

341:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
342:      Adjoint model starts here
343:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
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:   VecGetArray(mu[0],&x_ptr);
350:   x_ptr[0] = -1.0;
351:   VecRestoreArray(mu[0],&x_ptr);

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

356:   TSAdjointSolve(ts);

358:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n");
359:   VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
360:   VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
361:   TSGetCostIntegral(ts,&q);
362:   VecView(q,PETSC_VIEWER_STDOUT_WORLD);
363:   VecGetArray(q,&x_ptr);
364:   PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));
365:   VecRestoreArray(q,&x_ptr);

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

369:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
370:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
371:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
372:   MatDestroy(&A);
373:   MatDestroy(&Jacp);
374:   VecDestroy(&U);
375:   VecDestroy(&lambda[0]);
376:   VecDestroy(&mu[0]);
377:   TSDestroy(&ts);
378:   PetscFinalize();
379:   return ierr;
380: }