1c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*F
4c4762a1bSJed Brown
5c4762a1bSJed Brown \begin{eqnarray}
6c4762a1bSJed Brown \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - \frac{EV}{X} \sin(\theta) -D(\omega - \omega_s)\\
7c4762a1bSJed Brown \frac{d \theta}{dt} = \omega - \omega_s
8c4762a1bSJed Brown \end{eqnarray}
9c4762a1bSJed Brown
10c4762a1bSJed Brown Ensemble of initial conditions
11c4762a1bSJed Brown ./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
12c4762a1bSJed Brown
13c4762a1bSJed Brown Fault at .1 seconds
14c4762a1bSJed Brown ./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
15c4762a1bSJed Brown
16c4762a1bSJed Brown Initial conditions same as when fault is ended
17c4762a1bSJed Brown ./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
18c4762a1bSJed Brown
19c4762a1bSJed Brown F*/
20c4762a1bSJed Brown
21c4762a1bSJed Brown /*
22c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this
23c4762a1bSJed Brown file automatically includes:
24c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors
25c4762a1bSJed Brown petscmat.h - matrices
26c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods
27c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners
28c4762a1bSJed Brown petscksp.h - linear solvers
29c4762a1bSJed Brown */
30c4762a1bSJed Brown
31c4762a1bSJed Brown #include <petscts.h>
32c4762a1bSJed Brown
33c4762a1bSJed Brown typedef struct {
34c4762a1bSJed Brown PetscScalar H, D, omega_s, Pmax, Pm, E, V, X;
35c4762a1bSJed Brown PetscReal tf, tcl;
36c4762a1bSJed Brown } AppCtx;
37c4762a1bSJed Brown
38c4762a1bSJed Brown /*
39c4762a1bSJed Brown Defines the ODE passed to the ODE solver
40c4762a1bSJed Brown */
IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx * ctx)41d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
42d71ae5a4SJacob Faibussowitsch {
43c4762a1bSJed Brown PetscScalar *f, Pmax;
44c4762a1bSJed Brown const PetscScalar *u, *udot;
45c4762a1bSJed Brown
46c4762a1bSJed Brown PetscFunctionBegin;
47c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */
489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot));
509566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
51c4762a1bSJed Brown 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 */
52c4762a1bSJed Brown else if (t >= ctx->tcl) Pmax = ctx->E / 0.745;
53c4762a1bSJed Brown else Pmax = ctx->Pmax;
54c4762a1bSJed Brown f[0] = udot[0] - ctx->omega_s * (u[1] - 1.0);
55c4762a1bSJed Brown f[1] = 2.0 * ctx->H * udot[1] + Pmax * PetscSinScalar(u[0]) + ctx->D * (u[1] - 1.0) - ctx->Pm;
56c4762a1bSJed Brown
579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot));
599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
61c4762a1bSJed Brown }
62c4762a1bSJed Brown
63c4762a1bSJed Brown /*
64c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
65c4762a1bSJed Brown */
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx * ctx)66d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
67d71ae5a4SJacob Faibussowitsch {
68c4762a1bSJed Brown PetscInt rowcol[] = {0, 1};
69c4762a1bSJed Brown PetscScalar J[2][2], Pmax;
70c4762a1bSJed Brown const PetscScalar *u, *udot;
71c4762a1bSJed Brown
72c4762a1bSJed Brown PetscFunctionBegin;
739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot));
75c4762a1bSJed Brown 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 */
76c4762a1bSJed Brown else if (t >= ctx->tcl) Pmax = ctx->E / 0.745;
77c4762a1bSJed Brown else Pmax = ctx->Pmax;
78c4762a1bSJed Brown
799371c9d4SSatish Balay J[0][0] = a;
809371c9d4SSatish Balay J[0][1] = -ctx->omega_s;
819371c9d4SSatish Balay J[1][1] = 2.0 * ctx->H * a + ctx->D;
829371c9d4SSatish Balay J[1][0] = Pmax * PetscCosScalar(u[0]);
83c4762a1bSJed Brown
849566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot));
87c4762a1bSJed Brown
889566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
899566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
90c4762a1bSJed Brown if (A != B) {
919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
93c4762a1bSJed Brown }
943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown
PostStep(TS ts)97d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStep(TS ts)
98d71ae5a4SJacob Faibussowitsch {
99c4762a1bSJed Brown Vec X;
100c4762a1bSJed Brown PetscReal t;
101c4762a1bSJed Brown
102c4762a1bSJed Brown PetscFunctionBegin;
1039566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t));
104c4762a1bSJed Brown if (t >= .2) {
1059566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X));
1069566063dSJacob Faibussowitsch PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
107c4762a1bSJed Brown exit(0);
108c4762a1bSJed Brown /* results in initial conditions after fault of -u 0.496792,1.00932 */
109c4762a1bSJed Brown }
1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
111c4762a1bSJed Brown }
112c4762a1bSJed Brown
main(int argc,char ** argv)113d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
114d71ae5a4SJacob Faibussowitsch {
115c4762a1bSJed Brown TS ts; /* ODE integrator */
116c4762a1bSJed Brown Vec U; /* solution will be stored here */
117c4762a1bSJed Brown Mat A; /* Jacobian matrix */
118c4762a1bSJed Brown PetscMPIInt size;
119c4762a1bSJed Brown PetscInt n = 2;
120c4762a1bSJed Brown AppCtx ctx;
121c4762a1bSJed Brown PetscScalar *u;
122c4762a1bSJed Brown PetscReal du[2] = {0.0, 0.0};
123c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE, flg1, flg2;
124c4762a1bSJed Brown
125c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown Initialize program
127c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128327415f7SBarry Smith PetscFunctionBeginUser;
129c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1313c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
132c4762a1bSJed Brown
133c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134c4762a1bSJed Brown Create necessary matrix and vectors
135c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1369566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1379566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
1389566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE));
1399566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
1409566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
141c4762a1bSJed Brown
1429566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL));
143c4762a1bSJed Brown
144c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown Set runtime options
146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
148c4762a1bSJed Brown {
149c4762a1bSJed Brown ctx.omega_s = 2.0 * PETSC_PI * 60.0;
150c4762a1bSJed Brown ctx.H = 5.0;
1519566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
152c4762a1bSJed Brown ctx.D = 5.0;
1539566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
154c4762a1bSJed Brown ctx.E = 1.1378;
155c4762a1bSJed Brown ctx.V = 1.0;
156c4762a1bSJed Brown ctx.X = 0.545;
157c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X;
1589566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
159c4762a1bSJed Brown ctx.Pm = 0.9;
1609566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
161c4762a1bSJed Brown ctx.tf = 1.0;
162c4762a1bSJed Brown ctx.tcl = 1.05;
1639566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
1649566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
1659566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
166c4762a1bSJed Brown if (ensemble) {
167c4762a1bSJed Brown ctx.tf = -1;
168c4762a1bSJed Brown ctx.tcl = -1;
169c4762a1bSJed Brown }
170c4762a1bSJed Brown
1719566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u));
172c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
173c4762a1bSJed Brown u[1] = 1.0;
1749566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
175c4762a1bSJed Brown n = 2;
1769566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
177c4762a1bSJed Brown u[0] += du[0];
178c4762a1bSJed Brown u[1] += du[1];
1799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u));
180c4762a1bSJed Brown if (flg1 || flg2) {
181c4762a1bSJed Brown ctx.tf = -1;
182c4762a1bSJed Brown ctx.tcl = -1;
183c4762a1bSJed Brown }
184c4762a1bSJed Brown }
185d0609cedSBarry Smith PetscOptionsEnd();
186c4762a1bSJed Brown
187c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188c4762a1bSJed Brown Create timestepping solver context
189c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1909566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1919566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1929566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSROSW));
1938434afd1SBarry Smith PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, &ctx));
1948434afd1SBarry Smith PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobianFn *)IJacobian, &ctx));
195c4762a1bSJed Brown
196c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197c4762a1bSJed Brown Set initial conditions
198c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1999566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U));
200c4762a1bSJed Brown
201c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202c4762a1bSJed Brown Set solver options
203c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2049566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 35.0));
2059566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2069566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01));
2079566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
2089566063dSJacob Faibussowitsch /* PetscCall(TSSetPostStep(ts,PostStep)); */
209c4762a1bSJed Brown
210c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211c4762a1bSJed Brown Solve nonlinear system
212c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213c4762a1bSJed Brown if (ensemble) {
214c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
2159566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u));
216c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
217c4762a1bSJed Brown u[1] = ctx.omega_s;
218c4762a1bSJed Brown u[0] += du[0];
219c4762a1bSJed Brown u[1] += du[1];
2209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u));
2219566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01));
2229566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
223c4762a1bSJed Brown }
224c4762a1bSJed Brown } else {
2259566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
226c4762a1bSJed Brown }
227c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed.
229c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
2319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U));
2329566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
2339566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
234b122ec5aSJacob Faibussowitsch return 0;
235c4762a1bSJed Brown }
236c4762a1bSJed Brown
237c4762a1bSJed Brown /*TEST
238c4762a1bSJed Brown
239c4762a1bSJed Brown build:
240c4762a1bSJed Brown requires: !complex
241c4762a1bSJed Brown
242c4762a1bSJed Brown test:
243*188af4bfSBarry Smith args: -nox -ts_time_step 10
2443886731fSPierre Jolivet output_file: output/empty.out
245c4762a1bSJed Brown
246c4762a1bSJed Brown TEST*/
247