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*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetApplicationContext(ts,&user)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts,&t)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolution(ts,&X)); 86c4762a1bSJed Brown idx = 3*user->stepnum; 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(user->Sol,&mat)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 89c4762a1bSJed Brown mat[idx] = t; 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(mat+idx+1,x,2)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(user->Sol,&mat)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 109c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Wind Speed Options","");CHKERRQ(ierr); 110c4762a1bSJed Brown { 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-cw","","",user->cw,&user->cw,NULL)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-kw","","",user->kw,&user->kw,NULL)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-nsamples","","",user->nsamples,&user->nsamples,NULL)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-Tw","","",user->Tw,&user->Tw,NULL)); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user->wind_data)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(user->wind_data,PETSC_DECIDE,user->nsamples)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(user->wind_data)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(user->wind_data,&user->t_wind)); 121c4762a1bSJed Brown 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user->t_wind,&t)); 123c4762a1bSJed Brown for (i=0; i < user->nsamples; i++) t[i] = (i+1)*tmax/user->nsamples; 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user->t_wind,&t)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* Wind speed deviation = (-log(rand)/cw)^(1/kw) */ 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(user->wind_data,NULL)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecLog(user->wind_data)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(user->wind_data,-1/user->cw)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user->wind_data,&x)); 131c4762a1bSJed Brown for (i=0;i < user->nsamples;i++) x[i] = PetscPowScalar(x[i],(1/user->kw)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user->wind_data,&x)); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(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*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(user->wind_data,(1-avg_dev))); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(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*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&stepnum)); 194c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Udot,&udot)); 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(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*5f80ce2aSJacob Faibussowitsch CHKERRQ(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*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user->wind_data,&wd)); 206*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(Udot,&udot)); 208*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 PetscErrorCode ierr; 218c4762a1bSJed Brown PetscMPIInt size; 219c4762a1bSJed Brown PetscInt n = 2,idx; 220c4762a1bSJed Brown AppCtx user; 221c4762a1bSJed Brown PetscScalar *u; 222c4762a1bSJed Brown SNES snes; 223c4762a1bSJed Brown PetscScalar *mat; 224c4762a1bSJed Brown const PetscScalar *x,*rmat; 225c4762a1bSJed Brown Mat B; 226c4762a1bSJed Brown PetscScalar *amat; 227c4762a1bSJed Brown PetscViewer viewer; 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 230c4762a1bSJed Brown Initialize program 231c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 232c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 233*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 2343c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 237c4762a1bSJed Brown Create necessary matrix and vectors 238c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 239*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 242*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 243c4762a1bSJed Brown 244*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&U,NULL)); 245c4762a1bSJed Brown 246c4762a1bSJed Brown /* Create wind speed data using Weibull distribution */ 247*5f80ce2aSJacob Faibussowitsch CHKERRQ(WindSpeeds(&user)); 248c4762a1bSJed Brown /* Set parameters for wind turbine and induction generator */ 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetWindTurbineParams(&user)); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetInductionGeneratorParams(&user)); 251c4762a1bSJed Brown 252*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(U,&u)); 253c4762a1bSJed Brown u[0] = vwa; 254c4762a1bSJed Brown u[1] = s; 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(U,&u)); 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* Create matrix to save solutions at each time step */ 258c4762a1bSJed Brown user.stepnum = 0; 259c4762a1bSJed Brown 260*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,3,2010,NULL,&user.Sol)); 261c4762a1bSJed Brown 262c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 263c4762a1bSJed Brown Create timestepping solver context 264c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 265*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 266*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 267*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSBEULER)); 268*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user)); 269c4762a1bSJed Brown 270*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSNES(ts,&snes)); 271*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,A,A,SNESComputeJacobianDefault,NULL)); 272*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&user)); */ 273*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetApplicationContext(ts,&user)); 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 276c4762a1bSJed Brown Set initial conditions 277c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 278*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,U)); 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* Save initial solution */ 281c4762a1bSJed Brown idx=3*user.stepnum; 282c4762a1bSJed Brown 283*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(user.Sol,&mat)); 284*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&x)); 285c4762a1bSJed Brown 286c4762a1bSJed Brown mat[idx] = 0.0; 287c4762a1bSJed Brown 288*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(mat+idx+1,x,2)); 289*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(user.Sol,&mat)); 290*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&x)); 291c4762a1bSJed Brown user.stepnum++; 292c4762a1bSJed Brown 293c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 294c4762a1bSJed Brown Set solver options 295c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 296*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,20.0)); 297*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 298*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,.01)); 299*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 300*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetPostStep(ts,SaveSolution)); 301c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 302c4762a1bSJed Brown Solve nonlinear system 303c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 304*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,U)); 305c4762a1bSJed Brown 306*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,3,user.stepnum,NULL,&B)); 307*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(user.Sol,&rmat)); 308*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(B,&amat)); 309*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(amat,rmat,user.stepnum*3)); 310*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(B,&amat)); 311*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(user.Sol,&rmat)); 312c4762a1bSJed Brown 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer)); 314*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,viewer)); 315*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 316*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user.Sol)); 317*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 318c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 319c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 320c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 321*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user.wind_data)); 322*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user.t_wind)); 323*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 324*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); 325*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 326c4762a1bSJed Brown 327c4762a1bSJed Brown ierr = PetscFinalize(); 328c4762a1bSJed Brown return ierr; 329c4762a1bSJed Brown } 330c4762a1bSJed Brown 331c4762a1bSJed Brown /*TEST 332c4762a1bSJed Brown 333c4762a1bSJed Brown build: 334c4762a1bSJed Brown requires: !complex 335c4762a1bSJed Brown 336c4762a1bSJed Brown test: 337c4762a1bSJed Brown 338c4762a1bSJed Brown TEST*/ 339