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