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