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