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