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 { 72 AppCtx *user; 73 Vec X; 74 PetscScalar *mat; 75 const PetscScalar *x; 76 PetscInt idx; 77 PetscReal t; 78 79 PetscFunctionBegin; 80 PetscCall(TSGetApplicationContext(ts,&user)); 81 PetscCall(TSGetTime(ts,&t)); 82 PetscCall(TSGetSolution(ts,&X)); 83 idx = 3*user->stepnum; 84 PetscCall(MatDenseGetArray(user->Sol,&mat)); 85 PetscCall(VecGetArrayRead(X,&x)); 86 mat[idx] = t; 87 PetscCall(PetscArraycpy(mat+idx+1,x,2)); 88 PetscCall(MatDenseRestoreArray(user->Sol,&mat)); 89 PetscCall(VecRestoreArrayRead(X,&x)); 90 user->stepnum++; 91 PetscFunctionReturn(0); 92 } 93 94 /* Computes the wind speed using Weibull distribution */ 95 PetscErrorCode WindSpeeds(AppCtx *user) 96 { 97 PetscScalar *x,*t,avg_dev,sum; 98 PetscInt i; 99 100 PetscFunctionBegin; 101 user->cw = 5; 102 user->kw = 2; /* Rayleigh distribution */ 103 user->nsamples = 2000; 104 user->Tw = 0.2; 105 PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Wind Speed Options",""); 106 { 107 PetscCall(PetscOptionsReal("-cw","","",user->cw,&user->cw,NULL)); 108 PetscCall(PetscOptionsReal("-kw","","",user->kw,&user->kw,NULL)); 109 PetscCall(PetscOptionsInt("-nsamples","","",user->nsamples,&user->nsamples,NULL)); 110 PetscCall(PetscOptionsReal("-Tw","","",user->Tw,&user->Tw,NULL)); 111 } 112 PetscOptionsEnd(); 113 PetscCall(VecCreate(PETSC_COMM_WORLD,&user->wind_data)); 114 PetscCall(VecSetSizes(user->wind_data,PETSC_DECIDE,user->nsamples)); 115 PetscCall(VecSetFromOptions(user->wind_data)); 116 PetscCall(VecDuplicate(user->wind_data,&user->t_wind)); 117 118 PetscCall(VecGetArray(user->t_wind,&t)); 119 for (i=0; i < user->nsamples; i++) t[i] = (i+1)*tmax/user->nsamples; 120 PetscCall(VecRestoreArray(user->t_wind,&t)); 121 122 /* Wind speed deviation = (-log(rand)/cw)^(1/kw) */ 123 PetscCall(VecSetRandom(user->wind_data,NULL)); 124 PetscCall(VecLog(user->wind_data)); 125 PetscCall(VecScale(user->wind_data,-1/user->cw)); 126 PetscCall(VecGetArray(user->wind_data,&x)); 127 for (i=0;i < user->nsamples;i++) x[i] = PetscPowScalar(x[i],(1/user->kw)); 128 PetscCall(VecRestoreArray(user->wind_data,&x)); 129 PetscCall(VecSum(user->wind_data,&sum)); 130 avg_dev = sum/user->nsamples; 131 /* Wind speed (t) = (1 + wind speed deviation(t) - avg_dev)*average wind speed */ 132 PetscCall(VecShift(user->wind_data,(1-avg_dev))); 133 PetscCall(VecScale(user->wind_data,vwa)); 134 PetscFunctionReturn(0); 135 } 136 137 /* Sets the parameters for wind turbine */ 138 PetscErrorCode SetWindTurbineParams(AppCtx *user) 139 { 140 PetscFunctionBegin; 141 user->Rt = 35; 142 user->Ar = PETSC_PI*user->Rt*user->Rt; 143 user->nGB = 1.0/89.0; 144 user->rho = 1.225; 145 user->Ht = 1.5; 146 PetscFunctionReturn(0); 147 } 148 149 /* Sets the parameters for induction generator */ 150 PetscErrorCode SetInductionGeneratorParams(AppCtx *user) 151 { 152 PetscFunctionBegin; 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 PetscFunctionReturn(0); 163 } 164 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; 169 170 PetscFunctionBegin; 171 temp = user->nGB*2*user->Rt*ws/user->np; 172 lambda = temp*wm/vw; 173 lambda_i = 1/(1/lambda + 0.002); 174 cp = 0.44*(125/lambda_i - 6.94)*PetscExpScalar(-16.5/lambda_i); 175 *Pw = 0.5*user->rho*cp*user->Ar*vw*vw*vw/(MVAbase*1e6); 176 PetscFunctionReturn(0); 177 } 178 179 /* 180 Defines the ODE passed to the ODE solver 181 */ 182 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *user) 183 { 184 PetscScalar *f,wm,Pw,*wd; 185 const PetscScalar *u,*udot; 186 PetscInt stepnum; 187 188 PetscFunctionBegin; 189 PetscCall(TSGetStepNumber(ts,&stepnum)); 190 /* The next three lines allow us to access the entries of the vectors directly */ 191 PetscCall(VecGetArrayRead(U,&u)); 192 PetscCall(VecGetArrayRead(Udot,&udot)); 193 PetscCall(VecGetArray(F,&f)); 194 PetscCall(VecGetArray(user->wind_data,&wd)); 195 196 f[0] = user->Tw*udot[0] - wd[stepnum] + u[0]; 197 wm = 1-u[1]; 198 PetscCall(GetWindPower(wm,u[0],&Pw,user)); 199 f[1] = 2.0*(user->Ht+user->Hm)*udot[1] - Pw/wm + user->Te; 200 201 PetscCall(VecRestoreArray(user->wind_data,&wd)); 202 PetscCall(VecRestoreArrayRead(U,&u)); 203 PetscCall(VecRestoreArrayRead(Udot,&udot)); 204 PetscCall(VecRestoreArray(F,&f)); 205 PetscFunctionReturn(0); 206 } 207 208 int main(int argc,char **argv) 209 { 210 TS ts; /* ODE integrator */ 211 Vec U; /* solution will be stored here */ 212 Mat A; /* Jacobian matrix */ 213 PetscMPIInt size; 214 PetscInt n = 2,idx; 215 AppCtx user; 216 PetscScalar *u; 217 SNES snes; 218 PetscScalar *mat; 219 const PetscScalar *x,*rmat; 220 Mat B; 221 PetscScalar *amat; 222 PetscViewer viewer; 223 224 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225 Initialize program 226 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 227 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 228 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 229 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 230 231 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 232 Create necessary matrix and vectors 233 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 234 PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 235 PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 236 PetscCall(MatSetFromOptions(A)); 237 PetscCall(MatSetUp(A)); 238 239 PetscCall(MatCreateVecs(A,&U,NULL)); 240 241 /* Create wind speed data using Weibull distribution */ 242 PetscCall(WindSpeeds(&user)); 243 /* Set parameters for wind turbine and induction generator */ 244 PetscCall(SetWindTurbineParams(&user)); 245 PetscCall(SetInductionGeneratorParams(&user)); 246 247 PetscCall(VecGetArray(U,&u)); 248 u[0] = vwa; 249 u[1] = s; 250 PetscCall(VecRestoreArray(U,&u)); 251 252 /* Create matrix to save solutions at each time step */ 253 user.stepnum = 0; 254 255 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,3,2010,NULL,&user.Sol)); 256 257 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 258 Create timestepping solver context 259 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 260 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 261 PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 262 PetscCall(TSSetType(ts,TSBEULER)); 263 PetscCall(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user)); 264 265 PetscCall(TSGetSNES(ts,&snes)); 266 PetscCall(SNESSetJacobian(snes,A,A,SNESComputeJacobianDefault,NULL)); 267 /* PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&user)); */ 268 PetscCall(TSSetApplicationContext(ts,&user)); 269 270 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 271 Set initial conditions 272 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 273 PetscCall(TSSetSolution(ts,U)); 274 275 /* Save initial solution */ 276 idx=3*user.stepnum; 277 278 PetscCall(MatDenseGetArray(user.Sol,&mat)); 279 PetscCall(VecGetArrayRead(U,&x)); 280 281 mat[idx] = 0.0; 282 283 PetscCall(PetscArraycpy(mat+idx+1,x,2)); 284 PetscCall(MatDenseRestoreArray(user.Sol,&mat)); 285 PetscCall(VecRestoreArrayRead(U,&x)); 286 user.stepnum++; 287 288 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 289 Set solver options 290 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 291 PetscCall(TSSetMaxTime(ts,20.0)); 292 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 293 PetscCall(TSSetTimeStep(ts,.01)); 294 PetscCall(TSSetFromOptions(ts)); 295 PetscCall(TSSetPostStep(ts,SaveSolution)); 296 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 297 Solve nonlinear system 298 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 299 PetscCall(TSSolve(ts,U)); 300 301 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,3,user.stepnum,NULL,&B)); 302 PetscCall(MatDenseGetArrayRead(user.Sol,&rmat)); 303 PetscCall(MatDenseGetArray(B,&amat)); 304 PetscCall(PetscArraycpy(amat,rmat,user.stepnum*3)); 305 PetscCall(MatDenseRestoreArray(B,&amat)); 306 PetscCall(MatDenseRestoreArrayRead(user.Sol,&rmat)); 307 308 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer)); 309 PetscCall(MatView(B,viewer)); 310 PetscCall(PetscViewerDestroy(&viewer)); 311 PetscCall(MatDestroy(&user.Sol)); 312 PetscCall(MatDestroy(&B)); 313 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 314 Free work space. All PETSc objects should be destroyed when they are no longer needed. 315 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 316 PetscCall(VecDestroy(&user.wind_data)); 317 PetscCall(VecDestroy(&user.t_wind)); 318 PetscCall(MatDestroy(&A)); 319 PetscCall(VecDestroy(&U)); 320 PetscCall(TSDestroy(&ts)); 321 322 PetscCall(PetscFinalize()); 323 return 0; 324 } 325 326 /*TEST 327 328 build: 329 requires: !complex 330 331 test: 332 333 TEST*/ 334