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