Actual source code: ex9adj.c
2: static char help[] = "Basic equation for generator stability analysis.\n";
4: /*F
6: \begin{eqnarray}
7: \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
8: \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
9: \end{eqnarray}
11: Ensemble of initial conditions
12: ./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
14: Fault at .1 seconds
15: ./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
17: Initial conditions same as when fault is ended
18: ./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
20: F*/
22: /*
23: Include "petscts.h" so that we can use TS solvers. Note that this
24: file automatically includes:
25: petscsys.h - base PETSc routines petscvec.h - vectors
26: petscmat.h - matrices
27: petscis.h - index sets petscksp.h - Krylov subspace methods
28: petscviewer.h - viewers petscpc.h - preconditioners
29: petscksp.h - linear solvers
30: */
32: #include <petscts.h>
34: typedef struct {
35: PetscScalar H, D, omega_b, omega_s, Pmax, Pm, E, V, X, u_s, c;
36: PetscInt beta;
37: PetscReal tf, tcl;
38: } AppCtx;
40: PetscErrorCode PostStepFunction(TS ts)
41: {
42: Vec U;
43: PetscReal t;
44: const PetscScalar *u;
46: TSGetTime(ts, &t);
47: TSGetSolution(ts, &U);
48: VecGetArrayRead(U, &u);
49: PetscPrintf(PETSC_COMM_SELF, "delta(%3.2f) = %8.7f\n", (double)t, (double)u[0]);
50: VecRestoreArrayRead(U, &u);
51: return 0;
52: }
54: /*
55: Defines the ODE passed to the ODE solver
56: */
57: static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
58: {
59: PetscScalar *f, Pmax;
60: const PetscScalar *u;
62: /* The next three lines allow us to access the entries of the vectors directly */
63: VecGetArrayRead(U, &u);
64: VecGetArray(F, &f);
65: 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 */
66: else Pmax = ctx->Pmax;
68: f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
69: f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H);
71: VecRestoreArrayRead(U, &u);
72: VecRestoreArray(F, &f);
73: return 0;
74: }
76: /*
77: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
78: */
79: static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx)
80: {
81: PetscInt rowcol[] = {0, 1};
82: PetscScalar J[2][2], Pmax;
83: const PetscScalar *u;
85: VecGetArrayRead(U, &u);
86: 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 */
87: else Pmax = ctx->Pmax;
89: J[0][0] = 0;
90: J[0][1] = ctx->omega_b;
91: J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H);
92: J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H);
94: MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES);
95: VecRestoreArrayRead(U, &u);
97: MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
98: MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
99: if (A != B) {
100: MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
101: MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
102: }
103: return 0;
104: }
106: static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0)
107: {
108: PetscInt row[] = {0, 1}, col[] = {0};
109: PetscScalar J[2][1];
110: AppCtx *ctx = (AppCtx *)ctx0;
113: J[0][0] = 0;
114: J[1][0] = ctx->omega_s / (2.0 * ctx->H);
115: MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES);
116: MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
117: MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
118: return 0;
119: }
121: static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx)
122: {
123: PetscScalar *r;
124: const PetscScalar *u;
126: VecGetArrayRead(U, &u);
127: VecGetArray(R, &r);
128: r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta);
129: VecRestoreArray(R, &r);
130: VecRestoreArrayRead(U, &u);
131: return 0;
132: }
134: static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx)
135: {
136: PetscScalar ru[1];
137: const PetscScalar *u;
138: PetscInt row[] = {0}, col[] = {0};
140: VecGetArrayRead(U, &u);
141: ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1);
142: VecRestoreArrayRead(U, &u);
143: MatSetValues(DRDU, 1, row, 1, col, ru, INSERT_VALUES);
144: MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY);
145: MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY);
146: return 0;
147: }
149: static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, AppCtx *ctx)
150: {
151: MatZeroEntries(DRDP);
152: MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY);
153: MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY);
154: return 0;
155: }
157: PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx)
158: {
159: PetscScalar sensip;
160: const PetscScalar *x, *y;
162: VecGetArrayRead(lambda, &x);
163: VecGetArrayRead(mu, &y);
164: sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0];
165: PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt parameter pm: %.7f \n", (double)sensip);
166: VecRestoreArrayRead(lambda, &x);
167: VecRestoreArrayRead(mu, &y);
168: return 0;
169: }
171: int main(int argc, char **argv)
172: {
173: TS ts, quadts; /* ODE integrator */
174: Vec U; /* solution will be stored here */
175: Mat A; /* Jacobian matrix */
176: Mat Jacp; /* Jacobian matrix */
177: Mat DRDU, DRDP;
178: PetscMPIInt size;
179: PetscInt n = 2;
180: AppCtx ctx;
181: PetscScalar *u;
182: PetscReal du[2] = {0.0, 0.0};
183: PetscBool ensemble = PETSC_FALSE, flg1, flg2;
184: PetscReal ftime;
185: PetscInt steps;
186: PetscScalar *x_ptr, *y_ptr;
187: Vec lambda[1], q, mu[1];
189: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190: Initialize program
191: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193: PetscInitialize(&argc, &argv, (char *)0, help);
194: MPI_Comm_size(PETSC_COMM_WORLD, &size);
197: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198: Create necessary matrix and vectors
199: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200: MatCreate(PETSC_COMM_WORLD, &A);
201: MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE);
202: MatSetType(A, MATDENSE);
203: MatSetFromOptions(A);
204: MatSetUp(A);
206: MatCreateVecs(A, &U, NULL);
208: MatCreate(PETSC_COMM_WORLD, &Jacp);
209: MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1);
210: MatSetFromOptions(Jacp);
211: MatSetUp(Jacp);
213: MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &DRDP);
214: MatSetUp(DRDP);
215: MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, &DRDU);
216: MatSetUp(DRDU);
218: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219: Set runtime options
220: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
221: PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
222: {
223: ctx.beta = 2;
224: ctx.c = 10000.0;
225: ctx.u_s = 1.0;
226: ctx.omega_s = 1.0;
227: ctx.omega_b = 120.0 * PETSC_PI;
228: ctx.H = 5.0;
229: PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL);
230: ctx.D = 5.0;
231: PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL);
232: ctx.E = 1.1378;
233: ctx.V = 1.0;
234: ctx.X = 0.545;
235: ctx.Pmax = ctx.E * ctx.V / ctx.X;
236: PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL);
237: ctx.Pm = 1.1;
238: PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL);
239: ctx.tf = 0.1;
240: ctx.tcl = 0.2;
241: PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL);
242: PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL);
243: PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL);
244: if (ensemble) {
245: ctx.tf = -1;
246: ctx.tcl = -1;
247: }
249: VecGetArray(U, &u);
250: u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
251: u[1] = 1.0;
252: PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1);
253: n = 2;
254: PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2);
255: u[0] += du[0];
256: u[1] += du[1];
257: VecRestoreArray(U, &u);
258: if (flg1 || flg2) {
259: ctx.tf = -1;
260: ctx.tcl = -1;
261: }
262: }
263: PetscOptionsEnd();
265: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
266: Create timestepping solver context
267: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
268: TSCreate(PETSC_COMM_WORLD, &ts);
269: TSSetProblemType(ts, TS_NONLINEAR);
270: TSSetEquationType(ts, TS_EQ_ODE_EXPLICIT); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
271: TSSetType(ts, TSRK);
272: TSSetRHSFunction(ts, NULL, (TSRHSFunction)RHSFunction, &ctx);
273: TSSetRHSJacobian(ts, A, A, (TSRHSJacobian)RHSJacobian, &ctx);
274: TSCreateQuadratureTS(ts, PETSC_TRUE, &quadts);
275: TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, &ctx);
276: TSSetRHSJacobian(quadts, DRDU, DRDU, (TSRHSJacobian)DRDUJacobianTranspose, &ctx);
277: TSSetRHSJacobianP(quadts, DRDP, (TSRHSJacobianP)DRDPJacobianTranspose, &ctx);
279: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
280: Set initial conditions
281: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
282: TSSetSolution(ts, U);
284: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
285: Save trajectory of solution so that TSAdjointSolve() may be used
286: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
287: TSSetSaveTrajectory(ts);
289: MatCreateVecs(A, &lambda[0], NULL);
290: /* Set initial conditions for the adjoint integration */
291: VecGetArray(lambda[0], &y_ptr);
292: y_ptr[0] = 0.0;
293: y_ptr[1] = 0.0;
294: VecRestoreArray(lambda[0], &y_ptr);
296: MatCreateVecs(Jacp, &mu[0], NULL);
297: VecGetArray(mu[0], &x_ptr);
298: x_ptr[0] = -1.0;
299: VecRestoreArray(mu[0], &x_ptr);
300: TSSetCostGradients(ts, 1, lambda, mu);
302: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
303: Set solver options
304: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
305: TSSetMaxTime(ts, 10.0);
306: TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
307: TSSetTimeStep(ts, .01);
308: TSSetFromOptions(ts);
310: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
311: Solve nonlinear system
312: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
313: if (ensemble) {
314: for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
315: VecGetArray(U, &u);
316: u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
317: u[1] = ctx.omega_s;
318: u[0] += du[0];
319: u[1] += du[1];
320: VecRestoreArray(U, &u);
321: TSSetTimeStep(ts, .01);
322: TSSolve(ts, U);
323: }
324: } else {
325: TSSolve(ts, U);
326: }
327: VecView(U, PETSC_VIEWER_STDOUT_WORLD);
328: TSGetSolveTime(ts, &ftime);
329: TSGetStepNumber(ts, &steps);
331: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
332: Adjoint model starts here
333: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
334: /* Set initial conditions for the adjoint integration */
335: VecGetArray(lambda[0], &y_ptr);
336: y_ptr[0] = 0.0;
337: y_ptr[1] = 0.0;
338: VecRestoreArray(lambda[0], &y_ptr);
340: VecGetArray(mu[0], &x_ptr);
341: x_ptr[0] = -1.0;
342: VecRestoreArray(mu[0], &x_ptr);
344: /* Set RHS JacobianP */
345: TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &ctx);
347: TSAdjointSolve(ts);
349: PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0] d[Psi(tf)]/d[omega0]\n");
350: VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD);
351: VecView(mu[0], PETSC_VIEWER_STDOUT_WORLD);
352: TSGetCostIntegral(ts, &q);
353: VecView(q, PETSC_VIEWER_STDOUT_WORLD);
354: VecGetArray(q, &x_ptr);
355: PetscPrintf(PETSC_COMM_WORLD, "\n cost function=%g\n", (double)(x_ptr[0] - ctx.Pm));
356: VecRestoreArray(q, &x_ptr);
358: ComputeSensiP(lambda[0], mu[0], &ctx);
360: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
361: Free work space. All PETSc objects should be destroyed when they are no longer needed.
362: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
363: MatDestroy(&A);
364: MatDestroy(&Jacp);
365: MatDestroy(&DRDU);
366: MatDestroy(&DRDP);
367: VecDestroy(&U);
368: VecDestroy(&lambda[0]);
369: VecDestroy(&mu[0]);
370: TSDestroy(&ts);
371: PetscFinalize();
372: return 0;
373: }
375: /*TEST
377: build:
378: requires: !complex
380: test:
381: args: -viewer_binary_skip_info -ts_adapt_type none
383: TEST*/