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