Actual source code: ex2.c
2: static char help[] = "Basic equation for generator stability analysis.\n";
4: /*F
6: \begin{eqnarray}
7: \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - \frac{EV}{X} \sin(\theta) -D(\omega - \omega_s)\\
8: \frac{d \theta}{dt} = \omega - \omega_s
9: \end{eqnarray}
11: Ensemble of initial conditions
12: ./ex2 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly
14: Fault at .1 seconds
15: ./ex2 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly
17: Initial conditions same as when fault is ended
18: ./ex2 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -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_s, Pmax, Pm, E, V, X;
36: PetscReal tf, tcl;
37: } AppCtx;
39: /*
40: Defines the ODE passed to the ODE solver
41: */
42: static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
43: {
44: PetscScalar *f, Pmax;
45: const PetscScalar *u, *udot;
47: PetscFunctionBegin;
48: /* The next three lines allow us to access the entries of the vectors directly */
49: PetscCall(VecGetArrayRead(U, &u));
50: PetscCall(VecGetArrayRead(Udot, &udot));
51: PetscCall(VecGetArray(F, &f));
52: 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 */
53: else if (t >= ctx->tcl) Pmax = ctx->E / 0.745;
54: else Pmax = ctx->Pmax;
55: f[0] = udot[0] - ctx->omega_s * (u[1] - 1.0);
56: f[1] = 2.0 * ctx->H * udot[1] + Pmax * PetscSinScalar(u[0]) + ctx->D * (u[1] - 1.0) - ctx->Pm;
58: PetscCall(VecRestoreArrayRead(U, &u));
59: PetscCall(VecRestoreArrayRead(Udot, &udot));
60: PetscCall(VecRestoreArray(F, &f));
61: PetscFunctionReturn(PETSC_SUCCESS);
62: }
64: /*
65: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
66: */
67: static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
68: {
69: PetscInt rowcol[] = {0, 1};
70: PetscScalar J[2][2], Pmax;
71: const PetscScalar *u, *udot;
73: PetscFunctionBegin;
74: PetscCall(VecGetArrayRead(U, &u));
75: PetscCall(VecGetArrayRead(Udot, &udot));
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 if (t >= ctx->tcl) Pmax = ctx->E / 0.745;
78: else Pmax = ctx->Pmax;
80: J[0][0] = a;
81: J[0][1] = -ctx->omega_s;
82: J[1][1] = 2.0 * ctx->H * a + ctx->D;
83: J[1][0] = Pmax * PetscCosScalar(u[0]);
85: PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
86: PetscCall(VecRestoreArrayRead(U, &u));
87: PetscCall(VecRestoreArrayRead(Udot, &udot));
89: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
90: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
91: if (A != B) {
92: PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
93: PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
94: }
95: PetscFunctionReturn(PETSC_SUCCESS);
96: }
98: PetscErrorCode PostStep(TS ts)
99: {
100: Vec X;
101: PetscReal t;
103: PetscFunctionBegin;
104: PetscCall(TSGetTime(ts, &t));
105: if (t >= .2) {
106: PetscCall(TSGetSolution(ts, &X));
107: PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
108: exit(0);
109: /* results in initial conditions after fault of -u 0.496792,1.00932 */
110: }
111: PetscFunctionReturn(PETSC_SUCCESS);
112: }
114: int main(int argc, char **argv)
115: {
116: TS ts; /* ODE integrator */
117: Vec U; /* solution will be stored here */
118: Mat A; /* Jacobian matrix */
119: PetscMPIInt size;
120: PetscInt n = 2;
121: AppCtx ctx;
122: PetscScalar *u;
123: PetscReal du[2] = {0.0, 0.0};
124: PetscBool ensemble = PETSC_FALSE, flg1, flg2;
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Initialize program
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: PetscFunctionBeginUser;
130: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
131: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
132: PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Create necessary matrix and vectors
136: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
138: PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
139: PetscCall(MatSetType(A, MATDENSE));
140: PetscCall(MatSetFromOptions(A));
141: PetscCall(MatSetUp(A));
143: PetscCall(MatCreateVecs(A, &U, NULL));
145: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: Set runtime options
147: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148: PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
149: {
150: ctx.omega_s = 2.0 * PETSC_PI * 60.0;
151: ctx.H = 5.0;
152: PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
153: ctx.D = 5.0;
154: PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
155: ctx.E = 1.1378;
156: ctx.V = 1.0;
157: ctx.X = 0.545;
158: ctx.Pmax = ctx.E * ctx.V / ctx.X;
159: PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
160: ctx.Pm = 0.9;
161: PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
162: ctx.tf = 1.0;
163: ctx.tcl = 1.05;
164: PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
165: PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
166: PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
167: if (ensemble) {
168: ctx.tf = -1;
169: ctx.tcl = -1;
170: }
172: PetscCall(VecGetArray(U, &u));
173: u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
174: u[1] = 1.0;
175: PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
176: n = 2;
177: PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
178: u[0] += du[0];
179: u[1] += du[1];
180: PetscCall(VecRestoreArray(U, &u));
181: if (flg1 || flg2) {
182: ctx.tf = -1;
183: ctx.tcl = -1;
184: }
185: }
186: PetscOptionsEnd();
188: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189: Create timestepping solver context
190: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191: PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
192: PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
193: PetscCall(TSSetType(ts, TSROSW));
194: PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &ctx));
195: PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &ctx));
197: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198: Set initial conditions
199: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200: PetscCall(TSSetSolution(ts, U));
202: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: Set solver options
204: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205: PetscCall(TSSetMaxTime(ts, 35.0));
206: PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
207: PetscCall(TSSetTimeStep(ts, .01));
208: PetscCall(TSSetFromOptions(ts));
209: /* PetscCall(TSSetPostStep(ts,PostStep)); */
211: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212: Solve nonlinear system
213: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
214: if (ensemble) {
215: for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
216: PetscCall(VecGetArray(U, &u));
217: u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
218: u[1] = ctx.omega_s;
219: u[0] += du[0];
220: u[1] += du[1];
221: PetscCall(VecRestoreArray(U, &u));
222: PetscCall(TSSetTimeStep(ts, .01));
223: PetscCall(TSSolve(ts, U));
224: }
225: } else {
226: PetscCall(TSSolve(ts, U));
227: }
228: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229: Free work space. All PETSc objects should be destroyed when they are no longer needed.
230: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231: PetscCall(MatDestroy(&A));
232: PetscCall(VecDestroy(&U));
233: PetscCall(TSDestroy(&ts));
234: PetscCall(PetscFinalize());
235: return 0;
236: }
238: /*TEST
240: build:
241: requires: !complex
243: test:
244: args: -nox -ts_dt 10
246: TEST*/