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