xref: /petsc/src/ts/tutorials/ex25.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static const char help[] = "Time-dependent Brusselator reaction-diffusion PDE in 1d. Demonstrates IMEX methods.\n";
2*c4762a1bSJed Brown /*
3*c4762a1bSJed Brown    u_t - alpha u_xx = A + u^2 v - (B+1) u
4*c4762a1bSJed Brown    v_t - alpha v_xx = B u - u^2 v
5*c4762a1bSJed Brown    0 < x < 1;
6*c4762a1bSJed Brown    A = 1, B = 3, alpha = 1/50
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown    Initial conditions:
9*c4762a1bSJed Brown    u(x,0) = 1 + sin(2 pi x)
10*c4762a1bSJed Brown    v(x,0) = 3
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown    Boundary conditions:
13*c4762a1bSJed Brown    u(0,t) = u(1,t) = 1
14*c4762a1bSJed Brown    v(0,t) = v(1,t) = 3
15*c4762a1bSJed Brown */
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown #include <petscdm.h>
18*c4762a1bSJed Brown #include <petscdmda.h>
19*c4762a1bSJed Brown #include <petscts.h>
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown typedef struct {
22*c4762a1bSJed Brown   PetscScalar u,v;
23*c4762a1bSJed Brown } Field;
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown typedef struct _User *User;
26*c4762a1bSJed Brown struct _User {
27*c4762a1bSJed Brown   PetscReal A,B;                /* Reaction coefficients */
28*c4762a1bSJed Brown   PetscReal alpha;              /* Diffusion coefficient */
29*c4762a1bSJed Brown   PetscReal uleft,uright;       /* Dirichlet boundary conditions */
30*c4762a1bSJed Brown   PetscReal vleft,vright;       /* Dirichlet boundary conditions */
31*c4762a1bSJed Brown };
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
34*c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
35*c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
36*c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown int main(int argc,char **argv)
39*c4762a1bSJed Brown {
40*c4762a1bSJed Brown   TS                ts;         /* nonlinear solver */
41*c4762a1bSJed Brown   Vec               X;          /* solution, residual vectors */
42*c4762a1bSJed Brown   Mat               J;          /* Jacobian matrix */
43*c4762a1bSJed Brown   PetscInt          steps,mx;
44*c4762a1bSJed Brown   PetscErrorCode    ierr;
45*c4762a1bSJed Brown   DM                da;
46*c4762a1bSJed Brown   PetscReal         ftime,hx,dt;
47*c4762a1bSJed Brown   struct _User      user;       /* user-defined work context */
48*c4762a1bSJed Brown   TSConvergedReason reason;
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
51*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
53*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59*c4762a1bSJed Brown      Extract global vectors from DMDA;
60*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr);
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown   /* Initialize user application context */
64*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
65*c4762a1bSJed Brown   {
66*c4762a1bSJed Brown     user.A      = 1;
67*c4762a1bSJed Brown     user.B      = 3;
68*c4762a1bSJed Brown     user.alpha  = 0.02;
69*c4762a1bSJed Brown     user.uleft  = 1;
70*c4762a1bSJed Brown     user.uright = 1;
71*c4762a1bSJed Brown     user.vleft  = 3;
72*c4762a1bSJed Brown     user.vright = 3;
73*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL);CHKERRQ(ierr);
74*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL);CHKERRQ(ierr);
75*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL);CHKERRQ(ierr);
76*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL);CHKERRQ(ierr);
77*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL);CHKERRQ(ierr);
78*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL);CHKERRQ(ierr);
79*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL);CHKERRQ(ierr);
80*c4762a1bSJed Brown   }
81*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84*c4762a1bSJed Brown      Create timestepping solver context
85*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,FormIFunction,&user);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&user);CHKERRQ(ierr);
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown   ftime = 10.0;
97*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101*c4762a1bSJed Brown      Set initial conditions
102*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103*c4762a1bSJed Brown   ierr = FormInitialSolution(ts,X,&user);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = VecGetSize(X,&mx);CHKERRQ(ierr);
106*c4762a1bSJed Brown   hx = 1.0/(PetscReal)(mx/2-1);
107*c4762a1bSJed Brown   dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
108*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111*c4762a1bSJed Brown      Set runtime options
112*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116*c4762a1bSJed Brown      Solve nonlinear system
117*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118*c4762a1bSJed Brown   ierr = TSSolve(ts,X);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);CHKERRQ(ierr);
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125*c4762a1bSJed Brown      Free work space.
126*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
130*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = PetscFinalize();
132*c4762a1bSJed Brown   return ierr;
133*c4762a1bSJed Brown }
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
136*c4762a1bSJed Brown {
137*c4762a1bSJed Brown   User           user = (User)ptr;
138*c4762a1bSJed Brown   DM             da;
139*c4762a1bSJed Brown   DMDALocalInfo  info;
140*c4762a1bSJed Brown   PetscInt       i;
141*c4762a1bSJed Brown   Field          *x,*xdot,*f;
142*c4762a1bSJed Brown   PetscReal      hx;
143*c4762a1bSJed Brown   Vec            Xloc;
144*c4762a1bSJed Brown   PetscErrorCode ierr;
145*c4762a1bSJed Brown 
146*c4762a1bSJed Brown   PetscFunctionBeginUser;
147*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
149*c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown   /*
152*c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
153*c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
154*c4762a1bSJed Brown      By placing code between these two statements, computations can be
155*c4762a1bSJed Brown      done while messages are in transition.
156*c4762a1bSJed Brown   */
157*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown   /* Get pointers to vector data */
162*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Xloc,&x);CHKERRQ(ierr);
163*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Xdot,&xdot);CHKERRQ(ierr);
164*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
167*c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
168*c4762a1bSJed Brown     if (i == 0) {
169*c4762a1bSJed Brown       f[i].u = hx * (x[i].u - user->uleft);
170*c4762a1bSJed Brown       f[i].v = hx * (x[i].v - user->vleft);
171*c4762a1bSJed Brown     } else if (i == info.mx-1) {
172*c4762a1bSJed Brown       f[i].u = hx * (x[i].u - user->uright);
173*c4762a1bSJed Brown       f[i].v = hx * (x[i].v - user->vright);
174*c4762a1bSJed Brown     } else {
175*c4762a1bSJed Brown       f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx;
176*c4762a1bSJed Brown       f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx;
177*c4762a1bSJed Brown     }
178*c4762a1bSJed Brown   }
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   /* Restore vectors */
181*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Xloc,&x);CHKERRQ(ierr);
182*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Xdot,&xdot);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
185*c4762a1bSJed Brown   PetscFunctionReturn(0);
186*c4762a1bSJed Brown }
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
189*c4762a1bSJed Brown {
190*c4762a1bSJed Brown   User           user = (User)ptr;
191*c4762a1bSJed Brown   DM             da;
192*c4762a1bSJed Brown   DMDALocalInfo  info;
193*c4762a1bSJed Brown   PetscInt       i;
194*c4762a1bSJed Brown   PetscReal      hx;
195*c4762a1bSJed Brown   Field          *x,*f;
196*c4762a1bSJed Brown   PetscErrorCode ierr;
197*c4762a1bSJed Brown 
198*c4762a1bSJed Brown   PetscFunctionBeginUser;
199*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
200*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
201*c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
202*c4762a1bSJed Brown 
203*c4762a1bSJed Brown   /* Get pointers to vector data */
204*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,X,&x);CHKERRQ(ierr);
205*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
206*c4762a1bSJed Brown 
207*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
208*c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
209*c4762a1bSJed Brown     PetscScalar u = x[i].u,v = x[i].v;
210*c4762a1bSJed Brown     f[i].u = hx*(user->A + u*u*v - (user->B+1)*u);
211*c4762a1bSJed Brown     f[i].v = hx*(user->B*u - u*u*v);
212*c4762a1bSJed Brown   }
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown   /* Restore vectors */
215*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,X,&x);CHKERRQ(ierr);
216*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
217*c4762a1bSJed Brown   PetscFunctionReturn(0);
218*c4762a1bSJed Brown }
219*c4762a1bSJed Brown 
220*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
221*c4762a1bSJed Brown /*
222*c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
223*c4762a1bSJed Brown */
224*c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
225*c4762a1bSJed Brown {
226*c4762a1bSJed Brown   User           user = (User)ptr;
227*c4762a1bSJed Brown   PetscErrorCode ierr;
228*c4762a1bSJed Brown   DMDALocalInfo  info;
229*c4762a1bSJed Brown   PetscInt       i;
230*c4762a1bSJed Brown   PetscReal      hx;
231*c4762a1bSJed Brown   DM             da;
232*c4762a1bSJed Brown   Field          *x,*xdot;
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown   PetscFunctionBeginUser;
235*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
237*c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
238*c4762a1bSJed Brown 
239*c4762a1bSJed Brown   /* Get pointers to vector data */
240*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,X,&x);CHKERRQ(ierr);
241*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Xdot,&xdot);CHKERRQ(ierr);
242*c4762a1bSJed Brown 
243*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
244*c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
245*c4762a1bSJed Brown     if (i == 0 || i == info.mx-1) {
246*c4762a1bSJed Brown       const PetscInt    row        = i,col = i;
247*c4762a1bSJed Brown       const PetscScalar vals[2][2] = {{hx,0},{0,hx}};
248*c4762a1bSJed Brown       ierr = MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES);CHKERRQ(ierr);
249*c4762a1bSJed Brown     } else {
250*c4762a1bSJed Brown       const PetscInt    row           = i,col[] = {i-1,i,i+1};
251*c4762a1bSJed Brown       const PetscScalar dxxL          = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
252*c4762a1bSJed Brown       const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
253*c4762a1bSJed Brown                                          {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
254*c4762a1bSJed Brown       ierr = MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);CHKERRQ(ierr);
255*c4762a1bSJed Brown     }
256*c4762a1bSJed Brown   }
257*c4762a1bSJed Brown 
258*c4762a1bSJed Brown   /* Restore vectors */
259*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,X,&x);CHKERRQ(ierr);
260*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Xdot,&xdot);CHKERRQ(ierr);
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown   ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
263*c4762a1bSJed Brown   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
264*c4762a1bSJed Brown   if (J != Jpre) {
265*c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
266*c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
267*c4762a1bSJed Brown   }
268*c4762a1bSJed Brown   PetscFunctionReturn(0);
269*c4762a1bSJed Brown }
270*c4762a1bSJed Brown 
271*c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
272*c4762a1bSJed Brown {
273*c4762a1bSJed Brown   User           user = (User)ctx;
274*c4762a1bSJed Brown   DM             da;
275*c4762a1bSJed Brown   PetscInt       i;
276*c4762a1bSJed Brown   DMDALocalInfo  info;
277*c4762a1bSJed Brown   Field          *x;
278*c4762a1bSJed Brown   PetscReal      hx;
279*c4762a1bSJed Brown   PetscErrorCode ierr;
280*c4762a1bSJed Brown 
281*c4762a1bSJed Brown   PetscFunctionBeginUser;
282*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
283*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
284*c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
285*c4762a1bSJed Brown 
286*c4762a1bSJed Brown   /* Get pointers to vector data */
287*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
288*c4762a1bSJed Brown 
289*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
290*c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
291*c4762a1bSJed Brown     PetscReal xi = i*hx;
292*c4762a1bSJed Brown     x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi);
293*c4762a1bSJed Brown     x[i].v = user->vleft*(1.-xi) + user->vright*xi;
294*c4762a1bSJed Brown   }
295*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
296*c4762a1bSJed Brown   PetscFunctionReturn(0);
297*c4762a1bSJed Brown }
298*c4762a1bSJed Brown 
299*c4762a1bSJed Brown /*TEST
300*c4762a1bSJed Brown 
301*c4762a1bSJed Brown     build:
302*c4762a1bSJed Brown       requires: c99
303*c4762a1bSJed Brown 
304*c4762a1bSJed Brown     test:
305*c4762a1bSJed Brown       args: -nox -da_grid_x 20 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type 2p -ts_dt 5e-2 -ts_adapt_type none
306*c4762a1bSJed Brown 
307*c4762a1bSJed Brown TEST*/
308