xref: /petsc/src/ts/tutorials/power_grid/ex2.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 
2 static char help[] = "Basic equation for generator stability analysis.\n";
3 
4 /*F
5 
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}
10 
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
13 
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
16 
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
19 
20 F*/
21 
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 */
31 
32 #include <petscts.h>
33 
34 typedef struct {
35   PetscScalar H, D, omega_s, Pmax, Pm, E, V, X;
36   PetscReal   tf, tcl;
37 } AppCtx;
38 
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   PetscScalar       *f, Pmax;
44   const PetscScalar *u, *udot;
45 
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;
56 
57   PetscCall(VecRestoreArrayRead(U, &u));
58   PetscCall(VecRestoreArrayRead(Udot, &udot));
59   PetscCall(VecRestoreArray(F, &f));
60   PetscFunctionReturn(0);
61 }
62 
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   PetscInt           rowcol[] = {0, 1};
68   PetscScalar        J[2][2], Pmax;
69   const PetscScalar *u, *udot;
70 
71   PetscFunctionBegin;
72   PetscCall(VecGetArrayRead(U, &u));
73   PetscCall(VecGetArrayRead(Udot, &udot));
74   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 */
75   else if (t >= ctx->tcl) Pmax = ctx->E / 0.745;
76   else Pmax = ctx->Pmax;
77 
78   J[0][0] = a;
79   J[0][1] = -ctx->omega_s;
80   J[1][1] = 2.0 * ctx->H * a + ctx->D;
81   J[1][0] = Pmax * PetscCosScalar(u[0]);
82 
83   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
84   PetscCall(VecRestoreArrayRead(U, &u));
85   PetscCall(VecRestoreArrayRead(Udot, &udot));
86 
87   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
88   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
89   if (A != B) {
90     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
91     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
92   }
93   PetscFunctionReturn(0);
94 }
95 
96 PetscErrorCode PostStep(TS ts) {
97   Vec       X;
98   PetscReal t;
99 
100   PetscFunctionBegin;
101   PetscCall(TSGetTime(ts, &t));
102   if (t >= .2) {
103     PetscCall(TSGetSolution(ts, &X));
104     PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
105     exit(0);
106     /* results in initial conditions after fault of -u 0.496792,1.00932 */
107   }
108   PetscFunctionReturn(0);
109 }
110 
111 int main(int argc, char **argv) {
112   TS           ts; /* ODE integrator */
113   Vec          U;  /* solution will be stored here */
114   Mat          A;  /* Jacobian matrix */
115   PetscMPIInt  size;
116   PetscInt     n = 2;
117   AppCtx       ctx;
118   PetscScalar *u;
119   PetscReal    du[2]    = {0.0, 0.0};
120   PetscBool    ensemble = PETSC_FALSE, flg1, flg2;
121 
122   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123      Initialize program
124      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125   PetscFunctionBeginUser;
126   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
127   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
128   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
129 
130   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131     Create necessary matrix and vectors
132     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
134   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
135   PetscCall(MatSetType(A, MATDENSE));
136   PetscCall(MatSetFromOptions(A));
137   PetscCall(MatSetUp(A));
138 
139   PetscCall(MatCreateVecs(A, &U, NULL));
140 
141   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142     Set runtime options
143     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
145   {
146     ctx.omega_s = 2.0 * PETSC_PI * 60.0;
147     ctx.H       = 5.0;
148     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
149     ctx.D = 5.0;
150     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
151     ctx.E    = 1.1378;
152     ctx.V    = 1.0;
153     ctx.X    = 0.545;
154     ctx.Pmax = ctx.E * ctx.V / ctx.X;
155     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
156     ctx.Pm = 0.9;
157     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
158     ctx.tf  = 1.0;
159     ctx.tcl = 1.05;
160     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
161     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
162     PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
163     if (ensemble) {
164       ctx.tf  = -1;
165       ctx.tcl = -1;
166     }
167 
168     PetscCall(VecGetArray(U, &u));
169     u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
170     u[1] = 1.0;
171     PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
172     n = 2;
173     PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
174     u[0] += du[0];
175     u[1] += du[1];
176     PetscCall(VecRestoreArray(U, &u));
177     if (flg1 || flg2) {
178       ctx.tf  = -1;
179       ctx.tcl = -1;
180     }
181   }
182   PetscOptionsEnd();
183 
184   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185      Create timestepping solver context
186      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
188   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
189   PetscCall(TSSetType(ts, TSROSW));
190   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &ctx));
191   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &ctx));
192 
193   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194      Set initial conditions
195    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
196   PetscCall(TSSetSolution(ts, U));
197 
198   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199      Set solver options
200    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201   PetscCall(TSSetMaxTime(ts, 35.0));
202   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
203   PetscCall(TSSetTimeStep(ts, .01));
204   PetscCall(TSSetFromOptions(ts));
205   /* PetscCall(TSSetPostStep(ts,PostStep));  */
206 
207   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208      Solve nonlinear system
209      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210   if (ensemble) {
211     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
212       PetscCall(VecGetArray(U, &u));
213       u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
214       u[1] = ctx.omega_s;
215       u[0] += du[0];
216       u[1] += du[1];
217       PetscCall(VecRestoreArray(U, &u));
218       PetscCall(TSSetTimeStep(ts, .01));
219       PetscCall(TSSolve(ts, U));
220     }
221   } else {
222     PetscCall(TSSolve(ts, U));
223   }
224   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
226    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227   PetscCall(MatDestroy(&A));
228   PetscCall(VecDestroy(&U));
229   PetscCall(TSDestroy(&ts));
230   PetscCall(PetscFinalize());
231   return 0;
232 }
233 
234 /*TEST
235 
236    build:
237       requires: !complex
238 
239    test:
240       args: -nox -ts_dt 10
241 
242 TEST*/
243