xref: /petsc/src/ts/tutorials/power_grid/ex5.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
835f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetApplicationContext(ts,&user));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts,&t));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ts,&X));
86c4762a1bSJed Brown   idx      =  3*user->stepnum;
875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArray(user->Sol,&mat));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
89c4762a1bSJed Brown   mat[idx] = t;
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(mat+idx+1,x,2));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArray(user->Sol,&mat));
925f80ce2aSJacob 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   {
1115f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-cw","","",user->cw,&user->cw,NULL));
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-kw","","",user->kw,&user->kw,NULL));
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsInt("-nsamples","","",user->nsamples,&user->nsamples,NULL));
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-Tw","","",user->Tw,&user->Tw,NULL));
115c4762a1bSJed Brown   }
116c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user->wind_data));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(user->wind_data,PETSC_DECIDE,user->nsamples));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->wind_data));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->wind_data,&user->t_wind));
121c4762a1bSJed Brown 
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->t_wind,&t));
123c4762a1bSJed Brown   for (i=0; i < user->nsamples; i++) t[i] = (i+1)*tmax/user->nsamples;
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->t_wind,&t));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* Wind speed deviation = (-log(rand)/cw)^(1/kw) */
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(user->wind_data,NULL));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecLog(user->wind_data));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(user->wind_data,-1/user->cw));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->wind_data,&x));
131c4762a1bSJed Brown   for (i=0;i < user->nsamples;i++) x[i] = PetscPowScalar(x[i],(1/user->kw));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->wind_data,&x));
1335f80ce2aSJacob 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 */
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecShift(user->wind_data,(1-avg_dev)));
1375f80ce2aSJacob 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;
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&stepnum));
194c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Udot,&udot));
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
1985f80ce2aSJacob 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];
2025f80ce2aSJacob 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 
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->wind_data,&wd));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Udot,&udot));
2085f80ce2aSJacob 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   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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
2325f80ce2aSJacob Faibussowitsch   CHKERRMPI(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     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
2405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
2415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
242c4762a1bSJed Brown 
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&U,NULL));
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /* Create wind speed data using Weibull distribution */
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(WindSpeeds(&user));
247c4762a1bSJed Brown   /* Set parameters for wind turbine and induction generator */
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(SetWindTurbineParams(&user));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(SetInductionGeneratorParams(&user));
250c4762a1bSJed Brown 
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(U,&u));
252c4762a1bSJed Brown   u[0] = vwa;
253c4762a1bSJed Brown   u[1] = s;
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(U,&u));
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* Create matrix to save solutions at each time step */
257c4762a1bSJed Brown   user.stepnum = 0;
258c4762a1bSJed Brown 
2595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,3,2010,NULL,&user.Sol));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262c4762a1bSJed Brown      Create timestepping solver context
263c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSBEULER));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user));
268c4762a1bSJed Brown 
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSNES(ts,&snes));
2705f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,A,A,SNESComputeJacobianDefault,NULL));
2715f80ce2aSJacob Faibussowitsch   /*  CHKERRQ(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&user)); */
2725f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetApplicationContext(ts,&user));
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
275c4762a1bSJed Brown      Set initial conditions
276c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2775f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,U));
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   /* Save initial solution */
280c4762a1bSJed Brown   idx=3*user.stepnum;
281c4762a1bSJed Brown 
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArray(user.Sol,&mat));
2835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&x));
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   mat[idx] = 0.0;
286c4762a1bSJed Brown 
2875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(mat+idx+1,x,2));
2885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArray(user.Sol,&mat));
2895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&x));
290c4762a1bSJed Brown   user.stepnum++;
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
293c4762a1bSJed Brown      Set solver options
294c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2955f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,20.0));
2965f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
2975f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.01));
2985f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetPostStep(ts,SaveSolution));
300c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
301c4762a1bSJed Brown      Solve nonlinear system
302c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3035f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,U));
304c4762a1bSJed Brown 
3055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,3,user.stepnum,NULL,&B));
3065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArrayRead(user.Sol,&rmat));
3075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArray(B,&amat));
3085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(amat,rmat,user.stepnum*3));
3095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArray(B,&amat));
3105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArrayRead(user.Sol,&rmat));
311c4762a1bSJed Brown 
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer));
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,viewer));
3145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
3155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.Sol));
3165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
317c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
318c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
319c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.wind_data));
3215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.t_wind));
3225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
3235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
3245f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
325c4762a1bSJed Brown 
326*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
327*b122ec5aSJacob 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