xref: /petsc/src/ts/tutorials/power_grid/ex9.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
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{d \theta}{dt} = \omega_b (\omega - \omega_s)
7c4762a1bSJed Brown                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\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_b, 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 */
RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx * ctx)41d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
42d71ae5a4SJacob Faibussowitsch {
43c4762a1bSJed Brown   const PetscScalar *u;
44c4762a1bSJed Brown   PetscScalar       *f, Pmax;
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(VecGetArray(F, &f));
50c4762a1bSJed 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 */
51c4762a1bSJed Brown   else Pmax = ctx->Pmax;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
54c4762a1bSJed Brown   f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H);
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown /*
62c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
63c4762a1bSJed Brown */
RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx * ctx)64d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx)
65d71ae5a4SJacob Faibussowitsch {
66c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
67c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
68c4762a1bSJed Brown   const PetscScalar *u;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBegin;
719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
72c4762a1bSJed 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 */
73c4762a1bSJed Brown   else Pmax = ctx->Pmax;
74c4762a1bSJed Brown 
759371c9d4SSatish Balay   J[0][0] = 0;
769371c9d4SSatish Balay   J[0][1] = ctx->omega_b;
779371c9d4SSatish Balay   J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H);
789371c9d4SSatish Balay   J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H);
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
85c4762a1bSJed Brown   if (A != B) {
869566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
879566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
88c4762a1bSJed Brown   }
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90c4762a1bSJed Brown }
91c4762a1bSJed Brown 
main(int argc,char ** argv)92d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
93d71ae5a4SJacob Faibussowitsch {
94c4762a1bSJed Brown   TS           ts; /* ODE integrator */
95c4762a1bSJed Brown   Vec          U;  /* solution will be stored here */
96c4762a1bSJed Brown   Mat          A;  /* Jacobian matrix */
97c4762a1bSJed Brown   PetscMPIInt  size;
98c4762a1bSJed Brown   PetscInt     n = 2;
99c4762a1bSJed Brown   AppCtx       ctx;
100c4762a1bSJed Brown   PetscScalar *u;
101c4762a1bSJed Brown   PetscReal    du[2]    = {0.0, 0.0};
102c4762a1bSJed Brown   PetscBool    ensemble = PETSC_FALSE, flg1, flg2;
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105c4762a1bSJed Brown      Initialize program
106c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107327415f7SBarry Smith   PetscFunctionBeginUser;
108*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1103c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113c4762a1bSJed Brown     Create necessary matrix and vectors
114c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1159566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1169566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
1179566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATDENSE));
1189566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1199566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124c4762a1bSJed Brown     Set runtime options
125c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
127c4762a1bSJed Brown   {
128c4762a1bSJed Brown     ctx.omega_b = 1.0;
129c4762a1bSJed Brown     ctx.omega_s = 2.0 * PETSC_PI * 60.0;
130c4762a1bSJed Brown     ctx.H       = 5.0;
1319566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
132c4762a1bSJed Brown     ctx.D = 5.0;
1339566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
134c4762a1bSJed Brown     ctx.E    = 1.1378;
135c4762a1bSJed Brown     ctx.V    = 1.0;
136c4762a1bSJed Brown     ctx.X    = 0.545;
137c4762a1bSJed Brown     ctx.Pmax = ctx.E * ctx.V / ctx.X;
1389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
139c4762a1bSJed Brown     ctx.Pm = 0.9;
1409566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
141c4762a1bSJed Brown     ctx.tf  = 1.0;
142c4762a1bSJed Brown     ctx.tcl = 1.05;
1439566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
1449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
1459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
146c4762a1bSJed Brown     if (ensemble) {
147c4762a1bSJed Brown       ctx.tf  = -1;
148c4762a1bSJed Brown       ctx.tcl = -1;
149c4762a1bSJed Brown     }
150c4762a1bSJed Brown 
1519566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
152c4762a1bSJed Brown     u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
153c4762a1bSJed Brown     u[1] = 1.0;
1549566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
155c4762a1bSJed Brown     n = 2;
1569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
157c4762a1bSJed Brown     u[0] += du[0];
158c4762a1bSJed Brown     u[1] += du[1];
1599566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
160c4762a1bSJed Brown     if (flg1 || flg2) {
161c4762a1bSJed Brown       ctx.tf  = -1;
162c4762a1bSJed Brown       ctx.tcl = -1;
163c4762a1bSJed Brown     }
164c4762a1bSJed Brown   }
165d0609cedSBarry Smith   PetscOptionsEnd();
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168c4762a1bSJed Brown      Create timestepping solver context
169c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1709566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1719566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1729566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSTHETA));
1738434afd1SBarry Smith   PetscCall(TSSetRHSFunction(ts, NULL, (TSRHSFunctionFn *)RHSFunction, &ctx));
1748434afd1SBarry Smith   PetscCall(TSSetRHSJacobian(ts, A, A, (TSRHSJacobianFn *)RHSJacobian, &ctx));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177c4762a1bSJed Brown      Set initial conditions
178c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1799566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182c4762a1bSJed Brown      Set solver options
183c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1849566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 35.0));
1859566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1869566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .01));
1879566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190c4762a1bSJed Brown      Solve nonlinear system
191c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192c4762a1bSJed Brown   if (ensemble) {
193c4762a1bSJed Brown     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
1949566063dSJacob Faibussowitsch       PetscCall(VecGetArray(U, &u));
195c4762a1bSJed Brown       u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
196c4762a1bSJed Brown       u[1] = ctx.omega_s;
197c4762a1bSJed Brown       u[0] += du[0];
198c4762a1bSJed Brown       u[1] += du[1];
1999566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(U, &u));
2009566063dSJacob Faibussowitsch       PetscCall(TSSetTimeStep(ts, .01));
2019566063dSJacob Faibussowitsch       PetscCall(TSSolve(ts, U));
202c4762a1bSJed Brown     }
203c4762a1bSJed Brown   } else {
2049566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, U));
205c4762a1bSJed Brown   }
2069566063dSJacob Faibussowitsch   PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
207c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
209c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2109566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2129566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2139566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
214b122ec5aSJacob Faibussowitsch   return 0;
215c4762a1bSJed Brown }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown /*TEST
218c4762a1bSJed Brown 
219c4762a1bSJed Brown    build:
220c4762a1bSJed Brown      requires: !complex
221c4762a1bSJed Brown 
222c4762a1bSJed Brown    test:
223c4762a1bSJed Brown 
224c4762a1bSJed Brown TEST*/
225