Actual source code: ex5.c


  2: static char help[] = "Basic equation for an induction generator driven by a wind turbine.\n";

\begin{eqnarray}
T_w\frac{dv_w}{dt} & = & v_w - v_we \\
2(H_t+H_m)\frac{ds}{dt} & = & P_w - P_e
\end{eqnarray}
 10: /*
 11:  - Pw is the power extracted from the wind turbine given by
 12:            Pw = 0.5*\rho*cp*Ar*vw^3

 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

 20:  - This example does not work with adaptive time stepping!

 22: Reference:
 23: Power System Modeling and Scripting - F. Milano
 24: */
 25: /*T

 27: T*/

 29: #include <petscts.h>

 31: #define freq 50
 32: #define ws (2*PETSC_PI*freq)
 33: #define MVAbase 100

 35: typedef struct {
 36:   /* Parameters for wind speed model */
 37:   PetscInt  nsamples; /* Number of wind samples */
 38:   PetscReal cw;   /* Scale factor for Weibull distribution */
 39:   PetscReal kw;   /* Shape factor for Weibull distribution */
 40:   Vec       wind_data; /* Vector to hold wind speeds */
 41:   Vec       t_wind; /* Vector to hold wind speed times */
 42:   PetscReal Tw;     /* Filter time constant */

 44:   /* Wind turbine parameters */
 45:   PetscScalar Rt; /* Rotor radius */
 46:   PetscScalar Ar; /* Area swept by rotor (pi*R*R) */
 47:   PetscReal   nGB; /* Gear box ratio */
 48:   PetscReal   Ht;  /* Turbine inertia constant */
 49:   PetscReal   rho; /* Atmospheric pressure */

 51:   /* Induction generator parameters */
 52:   PetscInt    np; /* Number of poles */
 53:   PetscReal   Xm; /* Magnetizing reactance */
 54:   PetscReal   Xs; /* Stator Reactance */
 55:   PetscReal   Xr; /* Rotor reactance */
 56:   PetscReal   Rs; /* Stator resistance */
 57:   PetscReal   Rr; /* Rotor resistance */
 58:   PetscReal   Hm; /* Motor inertia constant */
 59:   PetscReal   Xp; /* Xs + Xm*Xr/(Xm + Xr) */
 60:   PetscScalar Te; /* Electrical Torque */

 62:   Mat      Sol;   /* Solution matrix */
 63:   PetscInt stepnum;   /* Column number of solution matrix */
 64: } AppCtx;

 66: /* Initial values computed by Power flow and initialization */
 67: PetscScalar s = -0.00011577790353;
 68: /*Pw = 0.011064344110238; %Te*wm */
 69: PetscScalar       vwa  = 22.317142184449754;
 70: PetscReal         tmax = 20.0;

 72: /* Saves the solution at each time to a matrix */
 73: PetscErrorCode SaveSolution(TS ts)
 74: {
 75:   AppCtx            *user;
 76:   Vec               X;
 77:   PetscScalar       *mat;
 78:   const PetscScalar *x;
 79:   PetscInt          idx;
 80:   PetscReal         t;

 82:   TSGetApplicationContext(ts,&user);
 83:   TSGetTime(ts,&t);
 84:   TSGetSolution(ts,&X);
 85:   idx      =  3*user->stepnum;
 86:   MatDenseGetArray(user->Sol,&mat);
 87:   VecGetArrayRead(X,&x);
 88:   mat[idx] = t;
 89:   PetscArraycpy(mat+idx+1,x,2);
 90:   MatDenseRestoreArray(user->Sol,&mat);
 91:   VecRestoreArrayRead(X,&x);
 92:   user->stepnum++;
 93:   return 0;
 94: }

 96: /* Computes the wind speed using Weibull distribution */
 97: PetscErrorCode WindSpeeds(AppCtx *user)
 98: {
100:   PetscScalar    *x,*t,avg_dev,sum;
101:   PetscInt       i;

103:   user->cw       = 5;
104:   user->kw       = 2; /* Rayleigh distribution */
105:   user->nsamples = 2000;
106:   user->Tw       = 0.2;
107:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Wind Speed Options","");
108:   {
109:     PetscOptionsReal("-cw","","",user->cw,&user->cw,NULL);
110:     PetscOptionsReal("-kw","","",user->kw,&user->kw,NULL);
111:     PetscOptionsInt("-nsamples","","",user->nsamples,&user->nsamples,NULL);
112:     PetscOptionsReal("-Tw","","",user->Tw,&user->Tw,NULL);
113:   }
114:   PetscOptionsEnd();
115:   VecCreate(PETSC_COMM_WORLD,&user->wind_data);
116:   VecSetSizes(user->wind_data,PETSC_DECIDE,user->nsamples);
117:   VecSetFromOptions(user->wind_data);
118:   VecDuplicate(user->wind_data,&user->t_wind);

120:   VecGetArray(user->t_wind,&t);
121:   for (i=0; i < user->nsamples; i++) t[i] = (i+1)*tmax/user->nsamples;
122:   VecRestoreArray(user->t_wind,&t);

124:   /* Wind speed deviation = (-log(rand)/cw)^(1/kw) */
125:   VecSetRandom(user->wind_data,NULL);
126:   VecLog(user->wind_data);
127:   VecScale(user->wind_data,-1/user->cw);
128:   VecGetArray(user->wind_data,&x);
129:   for (i=0;i < user->nsamples;i++) x[i] = PetscPowScalar(x[i],(1/user->kw));
130:   VecRestoreArray(user->wind_data,&x);
131:   VecSum(user->wind_data,&sum);
132:   avg_dev = sum/user->nsamples;
133:   /* Wind speed (t) = (1 + wind speed deviation(t) - avg_dev)*average wind speed */
134:   VecShift(user->wind_data,(1-avg_dev));
135:   VecScale(user->wind_data,vwa);
136:   return 0;
137: }

139: /* Sets the parameters for wind turbine */
140: PetscErrorCode SetWindTurbineParams(AppCtx *user)
141: {
142:   user->Rt  = 35;
143:   user->Ar  = PETSC_PI*user->Rt*user->Rt;
144:   user->nGB = 1.0/89.0;
145:   user->rho = 1.225;
146:   user->Ht  = 1.5;
147:   return 0;
148: }

150: /* Sets the parameters for induction generator */
151: PetscErrorCode SetInductionGeneratorParams(AppCtx *user)
152: {
153:   user->np = 4;
154:   user->Xm = 3.0;
155:   user->Xs = 0.1;
156:   user->Xr = 0.08;
157:   user->Rs = 0.01;
158:   user->Rr = 0.01;
159:   user->Xp = user->Xs + user->Xm*user->Xr/(user->Xm + user->Xr);
160:   user->Hm = 1.0;
161:   user->Te = 0.011063063063251968;
162:   return 0;
163: }

165: /* Computes the power extracted from wind */
166: PetscErrorCode GetWindPower(PetscScalar wm,PetscScalar vw,PetscScalar *Pw,AppCtx *user)
167: {
168:   PetscScalar temp,lambda,lambda_i,cp;

170:   temp     = user->nGB*2*user->Rt*ws/user->np;
171:   lambda   = temp*wm/vw;
172:   lambda_i = 1/(1/lambda + 0.002);
173:   cp       = 0.44*(125/lambda_i - 6.94)*PetscExpScalar(-16.5/lambda_i);
174:   *Pw      = 0.5*user->rho*cp*user->Ar*vw*vw*vw/(MVAbase*1e6);
175:   return 0;
176: }

178: /*
179:      Defines the ODE passed to the ODE solver
180: */
181: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *user)
182: {
183:   PetscScalar       *f,wm,Pw,*wd;
184:   const PetscScalar *u,*udot;
185:   PetscInt          stepnum;

187:   TSGetStepNumber(ts,&stepnum);
188:   /*  The next three lines allow us to access the entries of the vectors directly */
189:   VecGetArrayRead(U,&u);
190:   VecGetArrayRead(Udot,&udot);
191:   VecGetArray(F,&f);
192:   VecGetArray(user->wind_data,&wd);

194:   f[0] = user->Tw*udot[0] - wd[stepnum] + u[0];
195:   wm   = 1-u[1];
196:   GetWindPower(wm,u[0],&Pw,user);
197:   f[1] = 2.0*(user->Ht+user->Hm)*udot[1] - Pw/wm + user->Te;

199:   VecRestoreArray(user->wind_data,&wd);
200:   VecRestoreArrayRead(U,&u);
201:   VecRestoreArrayRead(Udot,&udot);
202:   VecRestoreArray(F,&f);
203:   return 0;
204: }

206: int main(int argc,char **argv)
207: {
208:   TS                ts;            /* ODE integrator */
209:   Vec               U;             /* solution will be stored here */
210:   Mat               A;             /* Jacobian matrix */
211:   PetscMPIInt       size;
212:   PetscInt          n = 2,idx;
213:   AppCtx            user;
214:   PetscScalar       *u;
215:   SNES              snes;
216:   PetscScalar       *mat;
217:   const PetscScalar *x,*rmat;
218:   Mat               B;
219:   PetscScalar       *amat;
220:   PetscViewer       viewer;

222:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223:      Initialize program
224:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
225:   PetscInitialize(&argc,&argv,(char*)0,help);
226:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

229:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230:     Create necessary matrix and vectors
231:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
232:   MatCreate(PETSC_COMM_WORLD,&A);
233:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
234:   MatSetFromOptions(A);
235:   MatSetUp(A);

237:   MatCreateVecs(A,&U,NULL);

239:   /* Create wind speed data using Weibull distribution */
240:   WindSpeeds(&user);
241:   /* Set parameters for wind turbine and induction generator */
242:   SetWindTurbineParams(&user);
243:   SetInductionGeneratorParams(&user);

245:   VecGetArray(U,&u);
246:   u[0] = vwa;
247:   u[1] = s;
248:   VecRestoreArray(U,&u);

250:   /* Create matrix to save solutions at each time step */
251:   user.stepnum = 0;

253:   MatCreateSeqDense(PETSC_COMM_SELF,3,2010,NULL,&user.Sol);

255:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
256:      Create timestepping solver context
257:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
258:   TSCreate(PETSC_COMM_WORLD,&ts);
259:   TSSetProblemType(ts,TS_NONLINEAR);
260:   TSSetType(ts,TSBEULER);
261:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user);

263:   TSGetSNES(ts,&snes);
264:   SNESSetJacobian(snes,A,A,SNESComputeJacobianDefault,NULL);
265:   /*  TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&user); */
266:   TSSetApplicationContext(ts,&user);

268:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269:      Set initial conditions
270:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
271:   TSSetSolution(ts,U);

273:   /* Save initial solution */
274:   idx=3*user.stepnum;

276:   MatDenseGetArray(user.Sol,&mat);
277:   VecGetArrayRead(U,&x);

279:   mat[idx] = 0.0;

281:   PetscArraycpy(mat+idx+1,x,2);
282:   MatDenseRestoreArray(user.Sol,&mat);
283:   VecRestoreArrayRead(U,&x);
284:   user.stepnum++;

286:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
287:      Set solver options
288:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
289:   TSSetMaxTime(ts,20.0);
290:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
291:   TSSetTimeStep(ts,.01);
292:   TSSetFromOptions(ts);
293:   TSSetPostStep(ts,SaveSolution);
294:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
295:      Solve nonlinear system
296:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
297:   TSSolve(ts,U);

299:   MatCreateSeqDense(PETSC_COMM_SELF,3,user.stepnum,NULL,&B);
300:   MatDenseGetArrayRead(user.Sol,&rmat);
301:   MatDenseGetArray(B,&amat);
302:   PetscArraycpy(amat,rmat,user.stepnum*3);
303:   MatDenseRestoreArray(B,&amat);
304:   MatDenseRestoreArrayRead(user.Sol,&rmat);

306:   PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer);
307:   MatView(B,viewer);
308:   PetscViewerDestroy(&viewer);
309:   MatDestroy(&user.Sol);
310:   MatDestroy(&B);
311:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
312:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
313:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
314:   VecDestroy(&user.wind_data);
315:   VecDestroy(&user.t_wind);
316:   MatDestroy(&A);
317:   VecDestroy(&U);
318:   TSDestroy(&ts);

320:   PetscFinalize();
321:   return 0;
322: }

324: /*TEST

326:    build:
327:       requires: !complex

329:    test:

331: TEST*/