Actual source code: ex2.c
1: static char help[] = "Basic equation for generator stability analysis.\n";
3: /*F
5: \begin{eqnarray}
6: \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - \frac{EV}{X} \sin(\theta) -D(\omega - \omega_s)\\
7: \frac{d \theta}{dt} = \omega - \omega_s
8: \end{eqnarray}
10: Ensemble of initial conditions
11: ./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
13: Fault at .1 seconds
14: ./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
16: Initial conditions same as when fault is ended
17: ./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
19: F*/
21: /*
22: Include "petscts.h" so that we can use TS solvers. Note that this
23: file automatically includes:
24: petscsys.h - base PETSc routines petscvec.h - vectors
25: petscmat.h - matrices
26: petscis.h - index sets petscksp.h - Krylov subspace methods
27: petscviewer.h - viewers petscpc.h - preconditioners
28: petscksp.h - linear solvers
29: */
31: #include <petscts.h>
33: typedef struct {
34: PetscScalar H, D, omega_s, Pmax, Pm, E, V, X;
35: PetscReal tf, tcl;
36: } AppCtx;
38: /*
39: Defines the ODE passed to the ODE solver
40: */
41: static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
42: {
43: PetscScalar *f, Pmax;
44: const PetscScalar *u, *udot;
46: PetscFunctionBegin;
47: /* The next three lines allow us to access the entries of the vectors directly */
48: PetscCall(VecGetArrayRead(U, &u));
49: PetscCall(VecGetArrayRead(Udot, &udot));
50: PetscCall(VecGetArray(F, &f));
51: 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 */
52: else if (t >= ctx->tcl) Pmax = ctx->E / 0.745;
53: else Pmax = ctx->Pmax;
54: f[0] = udot[0] - ctx->omega_s * (u[1] - 1.0);
55: f[1] = 2.0 * ctx->H * udot[1] + Pmax * PetscSinScalar(u[0]) + ctx->D * (u[1] - 1.0) - ctx->Pm;
57: PetscCall(VecRestoreArrayRead(U, &u));
58: PetscCall(VecRestoreArrayRead(Udot, &udot));
59: PetscCall(VecRestoreArray(F, &f));
60: PetscFunctionReturn(PETSC_SUCCESS);
61: }
63: /*
64: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
65: */
66: static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
67: {
68: PetscInt rowcol[] = {0, 1};
69: PetscScalar J[2][2], Pmax;
70: const PetscScalar *u, *udot;
72: PetscFunctionBegin;
73: PetscCall(VecGetArrayRead(U, &u));
74: PetscCall(VecGetArrayRead(Udot, &udot));
75: 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 */
76: else if (t >= ctx->tcl) Pmax = ctx->E / 0.745;
77: else Pmax = ctx->Pmax;
79: J[0][0] = a;
80: J[0][1] = -ctx->omega_s;
81: J[1][1] = 2.0 * ctx->H * a + ctx->D;
82: J[1][0] = Pmax * PetscCosScalar(u[0]);
84: PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
85: PetscCall(VecRestoreArrayRead(U, &u));
86: PetscCall(VecRestoreArrayRead(Udot, &udot));
88: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
89: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
90: if (A != B) {
91: PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
92: PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
93: }
94: PetscFunctionReturn(PETSC_SUCCESS);
95: }
97: PetscErrorCode PostStep(TS ts)
98: {
99: Vec X;
100: PetscReal t;
102: PetscFunctionBegin;
103: PetscCall(TSGetTime(ts, &t));
104: if (t >= .2) {
105: PetscCall(TSGetSolution(ts, &X));
106: PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
107: exit(0);
108: /* results in initial conditions after fault of -u 0.496792,1.00932 */
109: }
110: PetscFunctionReturn(PETSC_SUCCESS);
111: }
113: int main(int argc, char **argv)
114: {
115: TS ts; /* ODE integrator */
116: Vec U; /* solution will be stored here */
117: Mat A; /* Jacobian matrix */
118: PetscMPIInt size;
119: PetscInt n = 2;
120: AppCtx ctx;
121: PetscScalar *u;
122: PetscReal du[2] = {0.0, 0.0};
123: PetscBool ensemble = PETSC_FALSE, flg1, flg2;
125: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: Initialize program
127: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128: PetscFunctionBeginUser;
129: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
130: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
131: PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
133: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134: Create necessary matrix and vectors
135: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
137: PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
138: PetscCall(MatSetType(A, MATDENSE));
139: PetscCall(MatSetFromOptions(A));
140: PetscCall(MatSetUp(A));
142: PetscCall(MatCreateVecs(A, &U, NULL));
144: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145: Set runtime options
146: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147: PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
148: {
149: ctx.omega_s = 2.0 * PETSC_PI * 60.0;
150: ctx.H = 5.0;
151: PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
152: ctx.D = 5.0;
153: PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
154: ctx.E = 1.1378;
155: ctx.V = 1.0;
156: ctx.X = 0.545;
157: ctx.Pmax = ctx.E * ctx.V / ctx.X;
158: PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
159: ctx.Pm = 0.9;
160: PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
161: ctx.tf = 1.0;
162: ctx.tcl = 1.05;
163: PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
164: PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
165: PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
166: if (ensemble) {
167: ctx.tf = -1;
168: ctx.tcl = -1;
169: }
171: PetscCall(VecGetArray(U, &u));
172: u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
173: u[1] = 1.0;
174: PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
175: n = 2;
176: PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
177: u[0] += du[0];
178: u[1] += du[1];
179: PetscCall(VecRestoreArray(U, &u));
180: if (flg1 || flg2) {
181: ctx.tf = -1;
182: ctx.tcl = -1;
183: }
184: }
185: PetscOptionsEnd();
187: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188: Create timestepping solver context
189: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
190: PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
191: PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
192: PetscCall(TSSetType(ts, TSROSW));
193: PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, &ctx));
194: PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobianFn *)IJacobian, &ctx));
196: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197: Set initial conditions
198: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
199: PetscCall(TSSetSolution(ts, U));
201: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202: Set solver options
203: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204: PetscCall(TSSetMaxTime(ts, 35.0));
205: PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
206: PetscCall(TSSetTimeStep(ts, .01));
207: PetscCall(TSSetFromOptions(ts));
208: /* PetscCall(TSSetPostStep(ts,PostStep)); */
210: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: Solve nonlinear system
212: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213: if (ensemble) {
214: for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
215: PetscCall(VecGetArray(U, &u));
216: u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
217: u[1] = ctx.omega_s;
218: u[0] += du[0];
219: u[1] += du[1];
220: PetscCall(VecRestoreArray(U, &u));
221: PetscCall(TSSetTimeStep(ts, .01));
222: PetscCall(TSSolve(ts, U));
223: }
224: } else {
225: PetscCall(TSSolve(ts, U));
226: }
227: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228: Free work space. All PETSc objects should be destroyed when they are no longer needed.
229: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
230: PetscCall(MatDestroy(&A));
231: PetscCall(VecDestroy(&U));
232: PetscCall(TSDestroy(&ts));
233: PetscCall(PetscFinalize());
234: return 0;
235: }
237: /*TEST
239: build:
240: requires: !complex
242: test:
243: args: -nox -ts_dt 10
245: TEST*/