1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Basic equation for an induction generator driven by a wind turbine.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*F 5c4762a1bSJed Brown \begin{eqnarray} 6c4762a1bSJed Brown T_w\frac{dv_w}{dt} & = & v_w - v_we \\ 7c4762a1bSJed Brown 2(H_t+H_m)\frac{ds}{dt} & = & P_w - P_e 8c4762a1bSJed Brown \end{eqnarray} 9c4762a1bSJed Brown F*/ 10c4762a1bSJed Brown /* 11c4762a1bSJed Brown - Pw is the power extracted from the wind turbine given by 12c4762a1bSJed Brown Pw = 0.5*\rho*cp*Ar*vw^3 13c4762a1bSJed Brown 14c4762a1bSJed Brown - The wind speed time series is modeled using a Weibull distribution and then 15c4762a1bSJed Brown passed through a low pass filter (with time constant T_w). 16c4762a1bSJed Brown - v_we is the wind speed data calculated using Weibull distribution while v_w is 17c4762a1bSJed Brown the output of the filter. 18c4762a1bSJed Brown - P_e is assumed as constant electrical torque 19c4762a1bSJed Brown 20c4762a1bSJed Brown - This example does not work with adaptive time stepping! 21c4762a1bSJed Brown 22c4762a1bSJed Brown Reference: 23c4762a1bSJed Brown Power System Modeling and Scripting - F. Milano 24c4762a1bSJed Brown */ 25c4762a1bSJed Brown 26c4762a1bSJed Brown #include <petscts.h> 27c4762a1bSJed Brown 28c4762a1bSJed Brown #define freq 50 29c4762a1bSJed Brown #define ws (2*PETSC_PI*freq) 30c4762a1bSJed Brown #define MVAbase 100 31c4762a1bSJed Brown 32c4762a1bSJed Brown typedef struct { 33c4762a1bSJed Brown /* Parameters for wind speed model */ 34c4762a1bSJed Brown PetscInt nsamples; /* Number of wind samples */ 35c4762a1bSJed Brown PetscReal cw; /* Scale factor for Weibull distribution */ 36c4762a1bSJed Brown PetscReal kw; /* Shape factor for Weibull distribution */ 37c4762a1bSJed Brown Vec wind_data; /* Vector to hold wind speeds */ 38c4762a1bSJed Brown Vec t_wind; /* Vector to hold wind speed times */ 39c4762a1bSJed Brown PetscReal Tw; /* Filter time constant */ 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Wind turbine parameters */ 42c4762a1bSJed Brown PetscScalar Rt; /* Rotor radius */ 43c4762a1bSJed Brown PetscScalar Ar; /* Area swept by rotor (pi*R*R) */ 44c4762a1bSJed Brown PetscReal nGB; /* Gear box ratio */ 45c4762a1bSJed Brown PetscReal Ht; /* Turbine inertia constant */ 46c4762a1bSJed Brown PetscReal rho; /* Atmospheric pressure */ 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Induction generator parameters */ 49c4762a1bSJed Brown PetscInt np; /* Number of poles */ 50c4762a1bSJed Brown PetscReal Xm; /* Magnetizing reactance */ 51c4762a1bSJed Brown PetscReal Xs; /* Stator Reactance */ 52c4762a1bSJed Brown PetscReal Xr; /* Rotor reactance */ 53c4762a1bSJed Brown PetscReal Rs; /* Stator resistance */ 54c4762a1bSJed Brown PetscReal Rr; /* Rotor resistance */ 55c4762a1bSJed Brown PetscReal Hm; /* Motor inertia constant */ 56c4762a1bSJed Brown PetscReal Xp; /* Xs + Xm*Xr/(Xm + Xr) */ 57c4762a1bSJed Brown PetscScalar Te; /* Electrical Torque */ 58c4762a1bSJed Brown 59c4762a1bSJed Brown Mat Sol; /* Solution matrix */ 60c4762a1bSJed Brown PetscInt stepnum; /* Column number of solution matrix */ 61c4762a1bSJed Brown } AppCtx; 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Initial values computed by Power flow and initialization */ 64c4762a1bSJed Brown PetscScalar s = -0.00011577790353; 65c4762a1bSJed Brown /*Pw = 0.011064344110238; %Te*wm */ 66c4762a1bSJed Brown PetscScalar vwa = 22.317142184449754; 67c4762a1bSJed Brown PetscReal tmax = 20.0; 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Saves the solution at each time to a matrix */ 70c4762a1bSJed Brown PetscErrorCode SaveSolution(TS ts) 71c4762a1bSJed Brown { 72c4762a1bSJed Brown AppCtx *user; 73c4762a1bSJed Brown Vec X; 74c4762a1bSJed Brown PetscScalar *mat; 75c4762a1bSJed Brown const PetscScalar *x; 76c4762a1bSJed Brown PetscInt idx; 77c4762a1bSJed Brown PetscReal t; 78c4762a1bSJed Brown 79c4762a1bSJed Brown PetscFunctionBegin; 809566063dSJacob Faibussowitsch PetscCall(TSGetApplicationContext(ts,&user)); 819566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&t)); 829566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts,&X)); 83c4762a1bSJed Brown idx = 3*user->stepnum; 849566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(user->Sol,&mat)); 859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 86c4762a1bSJed Brown mat[idx] = t; 879566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat+idx+1,x,2)); 889566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(user->Sol,&mat)); 899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 90c4762a1bSJed Brown user->stepnum++; 91c4762a1bSJed Brown PetscFunctionReturn(0); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Computes the wind speed using Weibull distribution */ 95c4762a1bSJed Brown PetscErrorCode WindSpeeds(AppCtx *user) 96c4762a1bSJed Brown { 97c4762a1bSJed Brown PetscScalar *x,*t,avg_dev,sum; 98c4762a1bSJed Brown PetscInt i; 99c4762a1bSJed Brown 100c4762a1bSJed Brown PetscFunctionBegin; 101c4762a1bSJed Brown user->cw = 5; 102c4762a1bSJed Brown user->kw = 2; /* Rayleigh distribution */ 103c4762a1bSJed Brown user->nsamples = 2000; 104c4762a1bSJed Brown user->Tw = 0.2; 105d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Wind Speed Options",""); 106c4762a1bSJed Brown { 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cw","","",user->cw,&user->cw,NULL)); 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-kw","","",user->kw,&user->kw,NULL)); 1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-nsamples","","",user->nsamples,&user->nsamples,NULL)); 1109566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-Tw","","",user->Tw,&user->Tw,NULL)); 111c4762a1bSJed Brown } 112d0609cedSBarry Smith PetscOptionsEnd(); 1139566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->wind_data)); 1149566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->wind_data,PETSC_DECIDE,user->nsamples)); 1159566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->wind_data)); 1169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->wind_data,&user->t_wind)); 117c4762a1bSJed Brown 1189566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->t_wind,&t)); 119c4762a1bSJed Brown for (i=0; i < user->nsamples; i++) t[i] = (i+1)*tmax/user->nsamples; 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->t_wind,&t)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* Wind speed deviation = (-log(rand)/cw)^(1/kw) */ 1239566063dSJacob Faibussowitsch PetscCall(VecSetRandom(user->wind_data,NULL)); 1249566063dSJacob Faibussowitsch PetscCall(VecLog(user->wind_data)); 1259566063dSJacob Faibussowitsch PetscCall(VecScale(user->wind_data,-1/user->cw)); 1269566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->wind_data,&x)); 127c4762a1bSJed Brown for (i=0;i < user->nsamples;i++) x[i] = PetscPowScalar(x[i],(1/user->kw)); 1289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->wind_data,&x)); 1299566063dSJacob Faibussowitsch PetscCall(VecSum(user->wind_data,&sum)); 130c4762a1bSJed Brown avg_dev = sum/user->nsamples; 131c4762a1bSJed Brown /* Wind speed (t) = (1 + wind speed deviation(t) - avg_dev)*average wind speed */ 1329566063dSJacob Faibussowitsch PetscCall(VecShift(user->wind_data,(1-avg_dev))); 1339566063dSJacob Faibussowitsch PetscCall(VecScale(user->wind_data,vwa)); 134c4762a1bSJed Brown PetscFunctionReturn(0); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Sets the parameters for wind turbine */ 138c4762a1bSJed Brown PetscErrorCode SetWindTurbineParams(AppCtx *user) 139c4762a1bSJed Brown { 140c4762a1bSJed Brown PetscFunctionBegin; 141c4762a1bSJed Brown user->Rt = 35; 142c4762a1bSJed Brown user->Ar = PETSC_PI*user->Rt*user->Rt; 143c4762a1bSJed Brown user->nGB = 1.0/89.0; 144c4762a1bSJed Brown user->rho = 1.225; 145c4762a1bSJed Brown user->Ht = 1.5; 146c4762a1bSJed Brown PetscFunctionReturn(0); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Sets the parameters for induction generator */ 150c4762a1bSJed Brown PetscErrorCode SetInductionGeneratorParams(AppCtx *user) 151c4762a1bSJed Brown { 152c4762a1bSJed Brown PetscFunctionBegin; 153c4762a1bSJed Brown user->np = 4; 154c4762a1bSJed Brown user->Xm = 3.0; 155c4762a1bSJed Brown user->Xs = 0.1; 156c4762a1bSJed Brown user->Xr = 0.08; 157c4762a1bSJed Brown user->Rs = 0.01; 158c4762a1bSJed Brown user->Rr = 0.01; 159c4762a1bSJed Brown user->Xp = user->Xs + user->Xm*user->Xr/(user->Xm + user->Xr); 160c4762a1bSJed Brown user->Hm = 1.0; 161c4762a1bSJed Brown user->Te = 0.011063063063251968; 162c4762a1bSJed Brown PetscFunctionReturn(0); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* Computes the power extracted from wind */ 166c4762a1bSJed Brown PetscErrorCode GetWindPower(PetscScalar wm,PetscScalar vw,PetscScalar *Pw,AppCtx *user) 167c4762a1bSJed Brown { 168c4762a1bSJed Brown PetscScalar temp,lambda,lambda_i,cp; 169c4762a1bSJed Brown 170c4762a1bSJed Brown PetscFunctionBegin; 171c4762a1bSJed Brown temp = user->nGB*2*user->Rt*ws/user->np; 172c4762a1bSJed Brown lambda = temp*wm/vw; 173c4762a1bSJed Brown lambda_i = 1/(1/lambda + 0.002); 174c4762a1bSJed Brown cp = 0.44*(125/lambda_i - 6.94)*PetscExpScalar(-16.5/lambda_i); 175c4762a1bSJed Brown *Pw = 0.5*user->rho*cp*user->Ar*vw*vw*vw/(MVAbase*1e6); 176c4762a1bSJed Brown PetscFunctionReturn(0); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* 180c4762a1bSJed Brown Defines the ODE passed to the ODE solver 181c4762a1bSJed Brown */ 182c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *user) 183c4762a1bSJed Brown { 184c4762a1bSJed Brown PetscScalar *f,wm,Pw,*wd; 185c4762a1bSJed Brown const PetscScalar *u,*udot; 186c4762a1bSJed Brown PetscInt stepnum; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBegin; 1899566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&stepnum)); 190c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 1919566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 1929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot,&udot)); 1939566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 1949566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->wind_data,&wd)); 195c4762a1bSJed Brown 196c4762a1bSJed Brown f[0] = user->Tw*udot[0] - wd[stepnum] + u[0]; 197c4762a1bSJed Brown wm = 1-u[1]; 1989566063dSJacob Faibussowitsch PetscCall(GetWindPower(wm,u[0],&Pw,user)); 199c4762a1bSJed Brown f[1] = 2.0*(user->Ht+user->Hm)*udot[1] - Pw/wm + user->Te; 200c4762a1bSJed Brown 2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->wind_data,&wd)); 2029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 2039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot,&udot)); 2049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 205c4762a1bSJed Brown PetscFunctionReturn(0); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208c4762a1bSJed Brown int main(int argc,char **argv) 209c4762a1bSJed Brown { 210c4762a1bSJed Brown TS ts; /* ODE integrator */ 211c4762a1bSJed Brown Vec U; /* solution will be stored here */ 212c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 213c4762a1bSJed Brown PetscMPIInt size; 214c4762a1bSJed Brown PetscInt n = 2,idx; 215c4762a1bSJed Brown AppCtx user; 216c4762a1bSJed Brown PetscScalar *u; 217c4762a1bSJed Brown SNES snes; 218c4762a1bSJed Brown PetscScalar *mat; 219c4762a1bSJed Brown const PetscScalar *x,*rmat; 220c4762a1bSJed Brown Mat B; 221c4762a1bSJed Brown PetscScalar *amat; 222c4762a1bSJed Brown PetscViewer viewer; 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225c4762a1bSJed Brown Initialize program 226c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 227*327415f7SBarry Smith PetscFunctionBeginUser; 2289566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 2299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 2303c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 233c4762a1bSJed Brown Create necessary matrix and vectors 234c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2359566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 2369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 2379566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2389566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 239c4762a1bSJed Brown 2409566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&U,NULL)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* Create wind speed data using Weibull distribution */ 2439566063dSJacob Faibussowitsch PetscCall(WindSpeeds(&user)); 244c4762a1bSJed Brown /* Set parameters for wind turbine and induction generator */ 2459566063dSJacob Faibussowitsch PetscCall(SetWindTurbineParams(&user)); 2469566063dSJacob Faibussowitsch PetscCall(SetInductionGeneratorParams(&user)); 247c4762a1bSJed Brown 2489566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u)); 249c4762a1bSJed Brown u[0] = vwa; 250c4762a1bSJed Brown u[1] = s; 2519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u)); 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* Create matrix to save solutions at each time step */ 254c4762a1bSJed Brown user.stepnum = 0; 255c4762a1bSJed Brown 2569566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,3,2010,NULL,&user.Sol)); 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 259c4762a1bSJed Brown Create timestepping solver context 260c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2619566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 2629566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 2639566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSBEULER)); 2649566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user)); 265c4762a1bSJed Brown 2669566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 2679566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,A,A,SNESComputeJacobianDefault,NULL)); 2689566063dSJacob Faibussowitsch /* PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&user)); */ 2699566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts,&user)); 270c4762a1bSJed Brown 271c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 272c4762a1bSJed Brown Set initial conditions 273c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2749566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,U)); 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* Save initial solution */ 277c4762a1bSJed Brown idx=3*user.stepnum; 278c4762a1bSJed Brown 2799566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(user.Sol,&mat)); 2809566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&x)); 281c4762a1bSJed Brown 282c4762a1bSJed Brown mat[idx] = 0.0; 283c4762a1bSJed Brown 2849566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat+idx+1,x,2)); 2859566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(user.Sol,&mat)); 2869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&x)); 287c4762a1bSJed Brown user.stepnum++; 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 290c4762a1bSJed Brown Set solver options 291c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2929566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,20.0)); 2939566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 2949566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,.01)); 2959566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 2969566063dSJacob Faibussowitsch PetscCall(TSSetPostStep(ts,SaveSolution)); 297c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 298c4762a1bSJed Brown Solve nonlinear system 299c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3009566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,U)); 301c4762a1bSJed Brown 3029566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,3,user.stepnum,NULL,&B)); 3039566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(user.Sol,&rmat)); 3049566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(B,&amat)); 3059566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(amat,rmat,user.stepnum*3)); 3069566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(B,&amat)); 3079566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(user.Sol,&rmat)); 308c4762a1bSJed Brown 3099566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer)); 3109566063dSJacob Faibussowitsch PetscCall(MatView(B,viewer)); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 3129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.Sol)); 3139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 314c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 315c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 316c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.wind_data)); 3189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.t_wind)); 3199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3219566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 322c4762a1bSJed Brown 3239566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 324b122ec5aSJacob Faibussowitsch return 0; 325c4762a1bSJed Brown } 326c4762a1bSJed Brown 327c4762a1bSJed Brown /*TEST 328c4762a1bSJed Brown 329c4762a1bSJed Brown build: 330c4762a1bSJed Brown requires: !complex 331c4762a1bSJed Brown 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown 334c4762a1bSJed Brown TEST*/ 335