Actual source code: ex9.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: ./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_b, 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 RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
43: {
44: const PetscScalar *u;
45: PetscScalar *f, Pmax;
47: PetscFunctionBegin;
48: /* The next three lines allow us to access the entries of the vectors directly */
49: PetscCall(VecGetArrayRead(U, &u));
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 Pmax = ctx->Pmax;
54: f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
55: f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H);
57: PetscCall(VecRestoreArrayRead(U, &u));
58: PetscCall(VecRestoreArray(F, &f));
59: PetscFunctionReturn(PETSC_SUCCESS);
60: }
62: /*
63: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
64: */
65: static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx)
66: {
67: PetscInt rowcol[] = {0, 1};
68: PetscScalar J[2][2], Pmax;
69: const PetscScalar *u;
71: PetscFunctionBegin;
72: PetscCall(VecGetArrayRead(U, &u));
73: 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 */
74: else Pmax = ctx->Pmax;
76: J[0][0] = 0;
77: J[0][1] = ctx->omega_b;
78: J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H);
79: J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H);
81: PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
82: PetscCall(VecRestoreArrayRead(U, &u));
84: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
85: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
86: if (A != B) {
87: PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
88: PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
89: }
90: PetscFunctionReturn(PETSC_SUCCESS);
91: }
93: int main(int argc, char **argv)
94: {
95: TS ts; /* ODE integrator */
96: Vec U; /* solution will be stored here */
97: Mat A; /* Jacobian matrix */
98: PetscMPIInt size;
99: PetscInt n = 2;
100: AppCtx ctx;
101: PetscScalar *u;
102: PetscReal du[2] = {0.0, 0.0};
103: PetscBool ensemble = PETSC_FALSE, flg1, flg2;
105: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: Initialize program
107: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: PetscFunctionBeginUser;
109: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
110: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
111: PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
113: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114: Create necessary matrix and vectors
115: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
117: PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
118: PetscCall(MatSetType(A, MATDENSE));
119: PetscCall(MatSetFromOptions(A));
120: PetscCall(MatSetUp(A));
122: PetscCall(MatCreateVecs(A, &U, NULL));
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Set runtime options
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127: PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
128: {
129: ctx.omega_b = 1.0;
130: ctx.omega_s = 2.0 * PETSC_PI * 60.0;
131: ctx.H = 5.0;
132: PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
133: ctx.D = 5.0;
134: PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
135: ctx.E = 1.1378;
136: ctx.V = 1.0;
137: ctx.X = 0.545;
138: ctx.Pmax = ctx.E * ctx.V / ctx.X;
139: PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
140: ctx.Pm = 0.9;
141: PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
142: ctx.tf = 1.0;
143: ctx.tcl = 1.05;
144: PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
145: PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
146: PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
147: if (ensemble) {
148: ctx.tf = -1;
149: ctx.tcl = -1;
150: }
152: PetscCall(VecGetArray(U, &u));
153: u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
154: u[1] = 1.0;
155: PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
156: n = 2;
157: PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
158: u[0] += du[0];
159: u[1] += du[1];
160: PetscCall(VecRestoreArray(U, &u));
161: if (flg1 || flg2) {
162: ctx.tf = -1;
163: ctx.tcl = -1;
164: }
165: }
166: PetscOptionsEnd();
168: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169: Create timestepping solver context
170: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
171: PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
172: PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
173: PetscCall(TSSetType(ts, TSTHETA));
174: PetscCall(TSSetRHSFunction(ts, NULL, (TSRHSFunction)RHSFunction, &ctx));
175: PetscCall(TSSetRHSJacobian(ts, A, A, (TSRHSJacobian)RHSJacobian, &ctx));
177: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178: Set initial conditions
179: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180: PetscCall(TSSetSolution(ts, U));
182: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183: Set solver options
184: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185: PetscCall(TSSetMaxTime(ts, 35.0));
186: PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
187: PetscCall(TSSetTimeStep(ts, .01));
188: PetscCall(TSSetFromOptions(ts));
190: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191: Solve nonlinear system
192: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193: if (ensemble) {
194: for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
195: PetscCall(VecGetArray(U, &u));
196: u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
197: u[1] = ctx.omega_s;
198: u[0] += du[0];
199: u[1] += du[1];
200: PetscCall(VecRestoreArray(U, &u));
201: PetscCall(TSSetTimeStep(ts, .01));
202: PetscCall(TSSolve(ts, U));
203: }
204: } else {
205: PetscCall(TSSolve(ts, U));
206: }
207: PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
208: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209: Free work space. All PETSc objects should be destroyed when they are no longer needed.
210: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211: PetscCall(MatDestroy(&A));
212: PetscCall(VecDestroy(&U));
213: PetscCall(TSDestroy(&ts));
214: PetscCall(PetscFinalize());
215: return 0;
216: }
218: /*TEST
220: build:
221: requires: !complex
223: test:
225: TEST*/