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