xref: /petsc/src/ts/tutorials/power_grid/ex5.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 
2 static char help[] = "Basic equation for an induction generator driven by a wind turbine.\n";
3 
4 /*F
5 \begin{eqnarray}
6           T_w\frac{dv_w}{dt} & = & v_w - v_we \\
7           2(H_t+H_m)\frac{ds}{dt} & = & P_w - P_e
8 \end{eqnarray}
9 F*/
10 /*
11  - Pw is the power extracted from the wind turbine given by
12            Pw = 0.5*\rho*cp*Ar*vw^3
13 
14  - The wind speed time series is modeled using a Weibull distribution and then
15    passed through a low pass filter (with time constant T_w).
16  - v_we is the wind speed data calculated using Weibull distribution while v_w is
17    the output of the filter.
18  - P_e is assumed as constant electrical torque
19 
20  - This example does not work with adaptive time stepping!
21 
22 Reference:
23 Power System Modeling and Scripting - F. Milano
24 */
25 
26 #include <petscts.h>
27 
28 #define freq    50
29 #define ws      (2 * PETSC_PI * freq)
30 #define MVAbase 100
31 
32 typedef struct {
33   /* Parameters for wind speed model */
34   PetscInt  nsamples;  /* Number of wind samples */
35   PetscReal cw;        /* Scale factor for Weibull distribution */
36   PetscReal kw;        /* Shape factor for Weibull distribution */
37   Vec       wind_data; /* Vector to hold wind speeds */
38   Vec       t_wind;    /* Vector to hold wind speed times */
39   PetscReal Tw;        /* Filter time constant */
40 
41   /* Wind turbine parameters */
42   PetscScalar Rt;  /* Rotor radius */
43   PetscScalar Ar;  /* Area swept by rotor (pi*R*R) */
44   PetscReal   nGB; /* Gear box ratio */
45   PetscReal   Ht;  /* Turbine inertia constant */
46   PetscReal   rho; /* Atmospheric pressure */
47 
48   /* Induction generator parameters */
49   PetscInt    np; /* Number of poles */
50   PetscReal   Xm; /* Magnetizing reactance */
51   PetscReal   Xs; /* Stator Reactance */
52   PetscReal   Xr; /* Rotor reactance */
53   PetscReal   Rs; /* Stator resistance */
54   PetscReal   Rr; /* Rotor resistance */
55   PetscReal   Hm; /* Motor inertia constant */
56   PetscReal   Xp; /* Xs + Xm*Xr/(Xm + Xr) */
57   PetscScalar Te; /* Electrical Torque */
58 
59   Mat      Sol;     /* Solution matrix */
60   PetscInt stepnum; /* Column number of solution matrix */
61 } AppCtx;
62 
63 /* Initial values computed by Power flow and initialization */
64 PetscScalar s    = -0.00011577790353;
65 /*Pw = 0.011064344110238; %Te*wm */
66 PetscScalar vwa  = 22.317142184449754;
67 PetscReal   tmax = 20.0;
68 
69 /* Saves the solution at each time to a matrix */
70 PetscErrorCode SaveSolution(TS ts) {
71   AppCtx            *user;
72   Vec                X;
73   PetscScalar       *mat;
74   const PetscScalar *x;
75   PetscInt           idx;
76   PetscReal          t;
77 
78   PetscFunctionBegin;
79   PetscCall(TSGetApplicationContext(ts, &user));
80   PetscCall(TSGetTime(ts, &t));
81   PetscCall(TSGetSolution(ts, &X));
82   idx = 3 * user->stepnum;
83   PetscCall(MatDenseGetArray(user->Sol, &mat));
84   PetscCall(VecGetArrayRead(X, &x));
85   mat[idx] = t;
86   PetscCall(PetscArraycpy(mat + idx + 1, x, 2));
87   PetscCall(MatDenseRestoreArray(user->Sol, &mat));
88   PetscCall(VecRestoreArrayRead(X, &x));
89   user->stepnum++;
90   PetscFunctionReturn(0);
91 }
92 
93 /* Computes the wind speed using Weibull distribution */
94 PetscErrorCode WindSpeeds(AppCtx *user) {
95   PetscScalar *x, *t, avg_dev, sum;
96   PetscInt     i;
97 
98   PetscFunctionBegin;
99   user->cw       = 5;
100   user->kw       = 2; /* Rayleigh distribution */
101   user->nsamples = 2000;
102   user->Tw       = 0.2;
103   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Wind Speed Options", "");
104   {
105     PetscCall(PetscOptionsReal("-cw", "", "", user->cw, &user->cw, NULL));
106     PetscCall(PetscOptionsReal("-kw", "", "", user->kw, &user->kw, NULL));
107     PetscCall(PetscOptionsInt("-nsamples", "", "", user->nsamples, &user->nsamples, NULL));
108     PetscCall(PetscOptionsReal("-Tw", "", "", user->Tw, &user->Tw, NULL));
109   }
110   PetscOptionsEnd();
111   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->wind_data));
112   PetscCall(VecSetSizes(user->wind_data, PETSC_DECIDE, user->nsamples));
113   PetscCall(VecSetFromOptions(user->wind_data));
114   PetscCall(VecDuplicate(user->wind_data, &user->t_wind));
115 
116   PetscCall(VecGetArray(user->t_wind, &t));
117   for (i = 0; i < user->nsamples; i++) t[i] = (i + 1) * tmax / user->nsamples;
118   PetscCall(VecRestoreArray(user->t_wind, &t));
119 
120   /* Wind speed deviation = (-log(rand)/cw)^(1/kw) */
121   PetscCall(VecSetRandom(user->wind_data, NULL));
122   PetscCall(VecLog(user->wind_data));
123   PetscCall(VecScale(user->wind_data, -1 / user->cw));
124   PetscCall(VecGetArray(user->wind_data, &x));
125   for (i = 0; i < user->nsamples; i++) x[i] = PetscPowScalar(x[i], (1 / user->kw));
126   PetscCall(VecRestoreArray(user->wind_data, &x));
127   PetscCall(VecSum(user->wind_data, &sum));
128   avg_dev = sum / user->nsamples;
129   /* Wind speed (t) = (1 + wind speed deviation(t) - avg_dev)*average wind speed */
130   PetscCall(VecShift(user->wind_data, (1 - avg_dev)));
131   PetscCall(VecScale(user->wind_data, vwa));
132   PetscFunctionReturn(0);
133 }
134 
135 /* Sets the parameters for wind turbine */
136 PetscErrorCode SetWindTurbineParams(AppCtx *user) {
137   PetscFunctionBegin;
138   user->Rt  = 35;
139   user->Ar  = PETSC_PI * user->Rt * user->Rt;
140   user->nGB = 1.0 / 89.0;
141   user->rho = 1.225;
142   user->Ht  = 1.5;
143   PetscFunctionReturn(0);
144 }
145 
146 /* Sets the parameters for induction generator */
147 PetscErrorCode SetInductionGeneratorParams(AppCtx *user) {
148   PetscFunctionBegin;
149   user->np = 4;
150   user->Xm = 3.0;
151   user->Xs = 0.1;
152   user->Xr = 0.08;
153   user->Rs = 0.01;
154   user->Rr = 0.01;
155   user->Xp = user->Xs + user->Xm * user->Xr / (user->Xm + user->Xr);
156   user->Hm = 1.0;
157   user->Te = 0.011063063063251968;
158   PetscFunctionReturn(0);
159 }
160 
161 /* Computes the power extracted from wind */
162 PetscErrorCode GetWindPower(PetscScalar wm, PetscScalar vw, PetscScalar *Pw, AppCtx *user) {
163   PetscScalar temp, lambda, lambda_i, cp;
164 
165   PetscFunctionBegin;
166   temp     = user->nGB * 2 * user->Rt * ws / user->np;
167   lambda   = temp * wm / vw;
168   lambda_i = 1 / (1 / lambda + 0.002);
169   cp       = 0.44 * (125 / lambda_i - 6.94) * PetscExpScalar(-16.5 / lambda_i);
170   *Pw      = 0.5 * user->rho * cp * user->Ar * vw * vw * vw / (MVAbase * 1e6);
171   PetscFunctionReturn(0);
172 }
173 
174 /*
175      Defines the ODE passed to the ODE solver
176 */
177 static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *user) {
178   PetscScalar       *f, wm, Pw, *wd;
179   const PetscScalar *u, *udot;
180   PetscInt           stepnum;
181 
182   PetscFunctionBegin;
183   PetscCall(TSGetStepNumber(ts, &stepnum));
184   /*  The next three lines allow us to access the entries of the vectors directly */
185   PetscCall(VecGetArrayRead(U, &u));
186   PetscCall(VecGetArrayRead(Udot, &udot));
187   PetscCall(VecGetArray(F, &f));
188   PetscCall(VecGetArray(user->wind_data, &wd));
189 
190   f[0] = user->Tw * udot[0] - wd[stepnum] + u[0];
191   wm   = 1 - u[1];
192   PetscCall(GetWindPower(wm, u[0], &Pw, user));
193   f[1] = 2.0 * (user->Ht + user->Hm) * udot[1] - Pw / wm + user->Te;
194 
195   PetscCall(VecRestoreArray(user->wind_data, &wd));
196   PetscCall(VecRestoreArrayRead(U, &u));
197   PetscCall(VecRestoreArrayRead(Udot, &udot));
198   PetscCall(VecRestoreArray(F, &f));
199   PetscFunctionReturn(0);
200 }
201 
202 int main(int argc, char **argv) {
203   TS                 ts; /* ODE integrator */
204   Vec                U;  /* solution will be stored here */
205   Mat                A;  /* Jacobian matrix */
206   PetscMPIInt        size;
207   PetscInt           n = 2, idx;
208   AppCtx             user;
209   PetscScalar       *u;
210   SNES               snes;
211   PetscScalar       *mat;
212   const PetscScalar *x, *rmat;
213   Mat                B;
214   PetscScalar       *amat;
215   PetscViewer        viewer;
216 
217   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218      Initialize program
219      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
220   PetscFunctionBeginUser;
221   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
222   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
223   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
224 
225   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226     Create necessary matrix and vectors
227     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
228   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
229   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
230   PetscCall(MatSetFromOptions(A));
231   PetscCall(MatSetUp(A));
232 
233   PetscCall(MatCreateVecs(A, &U, NULL));
234 
235   /* Create wind speed data using Weibull distribution */
236   PetscCall(WindSpeeds(&user));
237   /* Set parameters for wind turbine and induction generator */
238   PetscCall(SetWindTurbineParams(&user));
239   PetscCall(SetInductionGeneratorParams(&user));
240 
241   PetscCall(VecGetArray(U, &u));
242   u[0] = vwa;
243   u[1] = s;
244   PetscCall(VecRestoreArray(U, &u));
245 
246   /* Create matrix to save solutions at each time step */
247   user.stepnum = 0;
248 
249   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, 3, 2010, NULL, &user.Sol));
250 
251   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252      Create timestepping solver context
253      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
254   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
255   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
256   PetscCall(TSSetType(ts, TSBEULER));
257   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &user));
258 
259   PetscCall(TSGetSNES(ts, &snes));
260   PetscCall(SNESSetJacobian(snes, A, A, SNESComputeJacobianDefault, NULL));
261   /*  PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&user)); */
262   PetscCall(TSSetApplicationContext(ts, &user));
263 
264   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
265      Set initial conditions
266    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
267   PetscCall(TSSetSolution(ts, U));
268 
269   /* Save initial solution */
270   idx = 3 * user.stepnum;
271 
272   PetscCall(MatDenseGetArray(user.Sol, &mat));
273   PetscCall(VecGetArrayRead(U, &x));
274 
275   mat[idx] = 0.0;
276 
277   PetscCall(PetscArraycpy(mat + idx + 1, x, 2));
278   PetscCall(MatDenseRestoreArray(user.Sol, &mat));
279   PetscCall(VecRestoreArrayRead(U, &x));
280   user.stepnum++;
281 
282   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283      Set solver options
284    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
285   PetscCall(TSSetMaxTime(ts, 20.0));
286   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
287   PetscCall(TSSetTimeStep(ts, .01));
288   PetscCall(TSSetFromOptions(ts));
289   PetscCall(TSSetPostStep(ts, SaveSolution));
290   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
291      Solve nonlinear system
292      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
293   PetscCall(TSSolve(ts, U));
294 
295   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, 3, user.stepnum, NULL, &B));
296   PetscCall(MatDenseGetArrayRead(user.Sol, &rmat));
297   PetscCall(MatDenseGetArray(B, &amat));
298   PetscCall(PetscArraycpy(amat, rmat, user.stepnum * 3));
299   PetscCall(MatDenseRestoreArray(B, &amat));
300   PetscCall(MatDenseRestoreArrayRead(user.Sol, &rmat));
301 
302   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "out.bin", FILE_MODE_WRITE, &viewer));
303   PetscCall(MatView(B, viewer));
304   PetscCall(PetscViewerDestroy(&viewer));
305   PetscCall(MatDestroy(&user.Sol));
306   PetscCall(MatDestroy(&B));
307   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
308      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
309    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
310   PetscCall(VecDestroy(&user.wind_data));
311   PetscCall(VecDestroy(&user.t_wind));
312   PetscCall(MatDestroy(&A));
313   PetscCall(VecDestroy(&U));
314   PetscCall(TSDestroy(&ts));
315 
316   PetscCall(PetscFinalize());
317   return 0;
318 }
319 
320 /*TEST
321 
322    build:
323       requires: !complex
324 
325    test:
326 
327 TEST*/
328