Actual source code: ex5.c

petsc-3.7.3 2016-07-24
Report Typos and Errors
  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: #include <petscts.h>

 27: #define freq 50
 28: #define ws (2*PETSC_PI*freq)
 29: #define MVAbase 100

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

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

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

 58:   Mat      Sol;   /* Solution matrix */
 59:   PetscInt stepnum;   /* Column number of solution matrix */
 60: } AppCtx;

 62: /* Initial values computed by Power flow and initialization */
 63: PetscScalar s = -0.00011577790353;
 64: /*Pw = 0.011064344110238; %Te*wm */
 65: PetscScalar       vwa  = 22.317142184449754;
 66: PetscReal         tmax = 20.0;

 68: /* Saves the solution at each time to a matrix */
 71: PetscErrorCode SaveSolution(TS ts)
 72: {
 73:   PetscErrorCode    ierr;
 74:   AppCtx            *user;
 75:   Vec               X;
 76:   PetscScalar       *mat;
 77:   const PetscScalar *x;
 78:   PetscInt          idx;
 79:   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:   PetscMemcpy(mat+idx+1,x,2*sizeof(PetscScalar));
 90:   MatDenseRestoreArray(user->Sol,&mat);
 91:   VecRestoreArrayRead(X,&x);
 92:   user->stepnum++;
 93:   return(0);
 94: }


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

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

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

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

145: /* Sets the parameters for wind turbine */
146: PetscErrorCode SetWindTurbineParams(AppCtx *user)
147: {
149:   user->Rt  = 35;
150:   user->Ar  = PETSC_PI*user->Rt*user->Rt;
151:   user->nGB = 1.0/89.0;
152:   user->rho = 1.225;
153:   user->Ht  = 1.5;
154:   return(0);
155: }

159: /* Sets the parameters for induction generator */
160: PetscErrorCode SetInductionGeneratorParams(AppCtx *user)
161: {
163:   user->np = 4;
164:   user->Xm = 3.0;
165:   user->Xs = 0.1;
166:   user->Xr = 0.08;
167:   user->Rs = 0.01;
168:   user->Rr = 0.01;
169:   user->Xp = user->Xs + user->Xm*user->Xr/(user->Xm + user->Xr);
170:   user->Hm = 1.0;
171:   user->Te = 0.011063063063251968;
172:   return(0);
173: }

177: /* Computes the power extracted from wind */
178: PetscErrorCode GetWindPower(PetscScalar wm,PetscScalar vw,PetscScalar *Pw,AppCtx *user)
179: {
180:   PetscScalar temp,lambda,lambda_i,cp;

183:   temp     = user->nGB*2*user->Rt*ws/user->np;
184:   lambda   = temp*wm/vw;
185:   lambda_i = 1/(1/lambda + 0.002);
186:   cp       = 0.44*(125/lambda_i - 6.94)*PetscExpScalar(-16.5/lambda_i);
187:   *Pw      = 0.5*user->rho*cp*user->Ar*vw*vw*vw/(MVAbase*1e6);
188:   return(0);
189: }

193: /*
194:      Defines the ODE passed to the ODE solver
195: */
196: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *user)
197: {
198:   PetscErrorCode    ierr;
199:   PetscScalar       *f,wm,Pw,*wd;
200:   const PetscScalar *u,*udot;
201:   PetscInt          stepnum;

204:   TSGetTimeStepNumber(ts,&stepnum);
205:   /*  The next three lines allow us to access the entries of the vectors directly */
206:   VecGetArrayRead(U,&u);
207:   VecGetArrayRead(Udot,&udot);
208:   VecGetArray(F,&f);
209:   VecGetArray(user->wind_data,&wd);

211:   f[0] = user->Tw*udot[0] - wd[stepnum] + u[0];
212:   wm   = 1-u[1];
213:   GetWindPower(wm,u[0],&Pw,user);
214:   f[1] = 2.0*(user->Ht+user->Hm)*udot[1] - Pw/wm + user->Te;

216:   VecRestoreArray(user->wind_data,&wd);
217:   VecRestoreArrayRead(U,&u);
218:   VecRestoreArrayRead(Udot,&udot);
219:   VecRestoreArray(F,&f);
220:   return(0);
221: }

225: int main(int argc,char **argv)
226: {
227:   TS             ts;            /* ODE integrator */
228:   Vec            U;             /* solution will be stored here */
229:   Mat            A;             /* Jacobian matrix */
231:   PetscMPIInt    size;
232:   PetscInt       n = 2,idx;
233:   AppCtx         user;
234:   PetscScalar    *u;
235:   SNES           snes;
236:   PetscScalar       *mat;
237:   const PetscScalar *x;
238:   Mat         B;
239:   PetscScalar *amat;
240:   PetscViewer viewer;



244:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
245:      Initialize program
246:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
247:   PetscInitialize(&argc,&argv,(char*)0,help);
248:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
249:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

251:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252:     Create necessary matrix and vectors
253:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
254:   MatCreate(PETSC_COMM_WORLD,&A);
255:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
256:   MatSetFromOptions(A);
257:   MatSetUp(A);

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

261:   /* Create wind speed data using Weibull distribution */
262:   WindSpeeds(&user);
263:   /* Set parameters for wind turbine and induction generator */
264:   SetWindTurbineParams(&user);
265:   SetInductionGeneratorParams(&user);

267:   VecGetArray(U,&u);
268:   u[0] = vwa;
269:   u[1] = s;
270:   VecRestoreArray(U,&u);

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

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

277:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
278:      Create timestepping solver context
279:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
280:   TSCreate(PETSC_COMM_WORLD,&ts);
281:   TSSetProblemType(ts,TS_NONLINEAR);
282:   TSSetType(ts,TSBEULER);
283:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user);

285:   TSGetSNES(ts,&snes);
286:   SNESSetJacobian(snes,A,A,SNESComputeJacobianDefault,NULL);
287:   /*  TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&user); */
288:   TSSetApplicationContext(ts,&user);

290:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
291:      Set initial conditions
292:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
293:   TSSetSolution(ts,U);

295:   /* Save initial solution */
296:   idx=3*user.stepnum;

298:   MatDenseGetArray(user.Sol,&mat);
299:   VecGetArrayRead(U,&x);

301:   mat[idx] = 0.0;

303:   PetscMemcpy(mat+idx+1,x,2*sizeof(PetscScalar));
304:   MatDenseRestoreArray(user.Sol,&mat);
305:   VecRestoreArrayRead(U,&x);
306:   user.stepnum++;


309:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
310:      Set solver options
311:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
312:   TSSetDuration(ts,2000,20.0);
313:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
314:   TSSetInitialTimeStep(ts,0.0,.01);
315:   TSSetFromOptions(ts);
316:   TSSetPostStep(ts,SaveSolution);
317:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
318:      Solve nonlinear system
319:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
320:   TSSolve(ts,U);

322:   MatCreateSeqDense(PETSC_COMM_SELF,3,user.stepnum,NULL,&B);
323:   MatDenseGetArray(user.Sol,&mat);
324:   MatDenseGetArray(B,&amat);
325:   PetscMemcpy(amat,mat,user.stepnum*3*sizeof(PetscScalar));
326:   MatDenseRestoreArray(B,&amat);
327:   MatDenseRestoreArray(user.Sol,&mat);

329:   PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer);
330:   MatView(B,viewer);
331:   PetscViewerDestroy(&viewer);
332:   MatDestroy(&user.Sol);
333:   MatDestroy(&B);
334:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
335:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
336:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
337:   VecDestroy(&user.wind_data);
338:   VecDestroy(&user.t_wind);
339:   MatDestroy(&A);
340:   VecDestroy(&U);
341:   TSDestroy(&ts);

343:   PetscFinalize();
344:   return(0);
345: }