Actual source code: ex9adj.c
petsc-3.6.4 2016-04-12
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: {
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" ,t,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 *x,J[2][1];
126: AppCtx *ctx=(AppCtx*)ctx0;
129: VecGetArray (X,&x);
131: J[0][0] = 0;
132: J[1][0] = ctx->omega_s/(2.0*ctx->H);
133: MatSetValues (A,2,row,1,col,&J[0][0],INSERT_VALUES );
135: MatAssemblyBegin (A,MAT_FINAL_ASSEMBLY);
136: MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);
137: return (0);
138: }
142: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
143: {
144: PetscErrorCode ierr;
145: PetscScalar *r;
146: const PetscScalar *u;
149: VecGetArrayRead (U,&u);
150: VecGetArray (R,&r);
151: r[0] = ctx->c*PetscPowScalarInt(PetscMax (0., u[0]-ctx->u_s),ctx->beta);
152: VecRestoreArray (R,&r);
153: VecRestoreArrayRead (U,&u);
154: return (0);
155: }
159: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
160: {
161: PetscErrorCode ierr;
162: PetscScalar *ry;
163: const PetscScalar *u;
166: VecGetArrayRead (U,&u);
167: VecGetArray (drdy[0],&ry);
168: ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax (0., u[0]-ctx->u_s),ctx->beta-1);
169: VecRestoreArray (drdy[0],&ry);
170: VecRestoreArrayRead (U,&u);
171: return (0);
172: }
176: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
177: {
178: PetscErrorCode ierr;
179: PetscScalar *rp;
180: const PetscScalar *u;
183: VecGetArrayRead (U,&u);
184: VecGetArray (drdp[0],&rp);
185: rp[0] = 0.;
186: VecRestoreArray (drdp[0],&rp);
187: VecGetArrayRead (U,&u);
188: return (0);
189: }
193: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
194: {
195: PetscErrorCode ierr;
196: PetscScalar sensip;
197: const PetscScalar *x,*y;
198:
200: VecGetArrayRead (lambda,&x);
201: VecGetArrayRead (mu,&y);
202: sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
203: PetscPrintf (PETSC_COMM_WORLD ,"\n sensitivity wrt parameter pm: %.7f \n" ,(double)sensip);
204: VecRestoreArrayRead (lambda,&x);
205: VecRestoreArrayRead (mu,&y);
206: return (0);
207: }
211: int main(int argc,char **argv)
212: {
213: TS ts; /* ODE integrator */
214: Vec U; /* solution will be stored here */
215: Mat A; /* Jacobian matrix */
216: Mat Jacp; /* Jacobian matrix */
218: PetscMPIInt size;
219: PetscInt n = 2;
220: AppCtx ctx;
221: PetscScalar *u;
222: PetscReal du[2] = {0.0,0.0};
223: PetscBool ensemble = PETSC_FALSE ,flg1,flg2;
224: PetscReal ftime;
225: PetscInt steps;
226: PetscScalar *x_ptr,*y_ptr;
227: Vec lambda[1],q,mu[1];
229: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230: Initialize program
231: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
232: PetscInitialize (&argc,&argv,(char*)0,help);
233: MPI_Comm_size (PETSC_COMM_WORLD ,&size);
234: if (size > 1) SETERRQ (PETSC_COMM_WORLD ,PETSC_ERR_SUP,"Only for sequential runs" );
236: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237: Create necessary matrix and vectors
238: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
239: MatCreate (PETSC_COMM_WORLD ,&A);
240: MatSetSizes (A,n,n,PETSC_DETERMINE ,PETSC_DETERMINE );
241: MatSetType (A,MATDENSE );
242: MatSetFromOptions (A);
243: MatSetUp (A);
245: MatCreateVecs (A,&U,NULL);
247: MatCreate (PETSC_COMM_WORLD ,&Jacp);
248: MatSetSizes (Jacp,PETSC_DECIDE ,PETSC_DECIDE ,2,1);
249: MatSetFromOptions (Jacp);
250: MatSetUp (Jacp);
252: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
253: Set runtime options
254: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
255: PetscOptionsBegin (PETSC_COMM_WORLD ,NULL,"Swing equation options" ,"" );
256: {
257: ctx.beta = 2;
258: ctx.c = 10000.0;
259: ctx.u_s = 1.0;
260: ctx.omega_s = 1.0;
261: ctx.omega_b = 120.0*PETSC_PI;
262: ctx.H = 5.0;
263: PetscOptionsScalar ("-Inertia" ,"" ,"" ,ctx.H,&ctx.H,NULL);
264: ctx.D = 5.0;
265: PetscOptionsScalar ("-D" ,"" ,"" ,ctx.D,&ctx.D,NULL);
266: ctx.E = 1.1378;
267: ctx.V = 1.0;
268: ctx.X = 0.545;
269: ctx.Pmax = ctx.E*ctx.V/ctx.X;;
270: PetscOptionsScalar ("-Pmax" ,"" ,"" ,ctx.Pmax,&ctx.Pmax,NULL);
271: ctx.Pm = 1.1;
272: PetscOptionsScalar ("-Pm" ,"" ,"" ,ctx.Pm,&ctx.Pm,NULL);
273: ctx.tf = 0.1;
274: ctx.tcl = 0.2;
275: PetscOptionsReal ("-tf" ,"Time to start fault" ,"" ,ctx.tf,&ctx.tf,NULL);
276: PetscOptionsReal ("-tcl" ,"Time to end fault" ,"" ,ctx.tcl,&ctx.tcl,NULL);
277: PetscOptionsBool ("-ensemble" ,"Run ensemble of different initial conditions" ,"" ,ensemble,&ensemble,NULL);
278: if (ensemble) {
279: ctx.tf = -1;
280: ctx.tcl = -1;
281: }
283: VecGetArray (U,&u);
284: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
285: u[1] = 1.0;
286: PetscOptionsRealArray ("-u" ,"Initial solution" ,"" ,u,&n,&flg1);
287: n = 2;
288: PetscOptionsRealArray ("-du" ,"Perturbation in initial solution" ,"" ,du,&n,&flg2);
289: u[0] += du[0];
290: u[1] += du[1];
291: VecRestoreArray (U,&u);
292: if (flg1 || flg2) {
293: ctx.tf = -1;
294: ctx.tcl = -1;
295: }
296: }
297: PetscOptionsEnd ();
299: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
300: Create timestepping solver context
301: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
302: TSCreate (PETSC_COMM_WORLD ,&ts);
303: TSSetProblemType (ts,TS_NONLINEAR);
304: TSSetType (ts,TSRK );
305: TSSetRHSFunction (ts,NULL,(TSRHSFunction)RHSFunction,&ctx);
306: TSSetRHSJacobian (ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);
307: TSSetExactFinalTime (ts,TS_EXACTFINALTIME_MATCHSTEP);
309: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
310: Set initial conditions
311: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
312: TSSetSolution (ts,U);
314: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
315: Save trajectory of solution so that TSAdjointSolve () may be used
316: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
317: TSSetSaveTrajectory (ts);
319: MatCreateVecs (A,&lambda[0],NULL);
320: /* Set initial conditions for the adjoint integration */
321: VecGetArray (lambda[0],&y_ptr);
322: y_ptr[0] = 0.0; y_ptr[1] = 0.0;
323: VecRestoreArray (lambda[0],&y_ptr);
325: MatCreateVecs (Jacp,&mu[0],NULL);
326: VecGetArray (mu[0],&x_ptr);
327: x_ptr[0] = -1.0;
328: VecRestoreArray (mu[0],&x_ptr);
329: TSSetCostGradients (ts,1,lambda,mu);
330: TSSetCostIntegrand (ts,1,(PetscErrorCode (*)(TS ,PetscReal ,Vec ,Vec ,void*))CostIntegrand,
331: (PetscErrorCode (*)(TS ,PetscReal ,Vec ,Vec *,void*))DRDYFunction,
332: (PetscErrorCode (*)(TS ,PetscReal ,Vec ,Vec *,void*))DRDPFunction,&ctx);
334: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
335: Set solver options
336: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
337: TSSetDuration (ts,PETSC_DEFAULT ,10.0);
338: TSSetInitialTimeStep (ts,0.0,.01);
339: TSSetFromOptions (ts);
341: #if 0
342: TSSetPostStep (ts,PostStepFunction);
343: #endif
344: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
345: Solve nonlinear system
346: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
347: if (ensemble) {
348: for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
349: VecGetArray (U,&u);
350: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
351: u[1] = ctx.omega_s;
352: u[0] += du[0];
353: u[1] += du[1];
354: VecRestoreArray (U,&u);
355: TSSetInitialTimeStep (ts,0.0,.01);
356: TSSolve (ts,U);
357: }
358: } else {
359: TSSolve (ts,U);
360: }
361: VecView (U,PETSC_VIEWER_STDOUT_WORLD );
362: TSGetSolveTime (ts,&ftime);
363: TSGetTimeStepNumber (ts,&steps);
365: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
366: Adjoint model starts here
367: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
368: /* Set initial conditions for the adjoint integration */
369: VecGetArray (lambda[0],&y_ptr);
370: y_ptr[0] = 0.0; y_ptr[1] = 0.0;
371: VecRestoreArray (lambda[0],&y_ptr);
373: VecGetArray (mu[0],&x_ptr);
374: x_ptr[0] = -1.0;
375: VecRestoreArray (mu[0],&x_ptr);
377: /* Set RHS JacobianP */
378: TSAdjointSetRHSJacobian (ts,Jacp,RHSJacobianP,&ctx);
380: TSAdjointSolve (ts);
382: PetscPrintf (PETSC_COMM_WORLD ,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0] d[Psi(tf)]/d[omega0]\n" );
383: VecView (lambda[0],PETSC_VIEWER_STDOUT_WORLD );
384: VecView (mu[0],PETSC_VIEWER_STDOUT_WORLD );
385: TSGetCostIntegral (ts,&q);
386: VecView (q,PETSC_VIEWER_STDOUT_WORLD );
387: VecGetArray (q,&x_ptr);
388: PetscPrintf (PETSC_COMM_WORLD ,"\n cost function=%g\n" ,(double)(x_ptr[0]-ctx.Pm));
389:
390: ComputeSensiP(lambda[0],mu[0],&ctx);
392: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
393: Free work space. All PETSc objects should be destroyed when they are no longer needed.
394: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
395: MatDestroy (&A);
396: MatDestroy (&Jacp);
397: VecDestroy (&U);
398: VecDestroy (&lambda[0]);
399: VecDestroy (&mu[0]);
400: TSDestroy (&ts);
401: PetscFinalize ();
402: return (0);
403: }