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