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