xref: /petsc/src/ts/tutorials/ex22.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
1c4762a1bSJed Brown static const char help[] = "Time-dependent advection-reaction PDE in 1d, demonstrates IMEX methods.\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    u_t + a1*u_x = -k1*u + k2*v + s1
4c4762a1bSJed Brown    v_t + a2*v_x = k1*u - k2*v + s2
5c4762a1bSJed Brown    0 < x < 1;
6c4762a1bSJed Brown    a1 = 1, k1 = 10^6, s1 = 0,
7c4762a1bSJed Brown    a2 = 0, k2 = 2*k1, s2 = 1
8c4762a1bSJed Brown 
9c4762a1bSJed Brown    Initial conditions:
10c4762a1bSJed Brown    u(x,0) = 1 + s2*x
11c4762a1bSJed Brown    v(x,0) = k0/k1*u(x,0) + s1/k1
12c4762a1bSJed Brown 
13c4762a1bSJed Brown    Upstream boundary conditions:
14c4762a1bSJed Brown    u(0,t) = 1-sin(12*t)^4
15c4762a1bSJed Brown 
16c4762a1bSJed Brown    Useful command-line parameters:
17c4762a1bSJed Brown    -ts_arkimex_fully_implicit - treats advection implicitly, only relevant with -ts_type arkimex (default)
18c4762a1bSJed Brown    -ts_type rosw - use Rosenbrock-W method (linearly implicit IMEX, amortizes assembly and preconditioner setup)
19c4762a1bSJed Brown */
20c4762a1bSJed Brown 
21c4762a1bSJed Brown #include <petscdm.h>
22c4762a1bSJed Brown #include <petscdmda.h>
23c4762a1bSJed Brown #include <petscts.h>
24c4762a1bSJed Brown 
25c4762a1bSJed Brown typedef PetscScalar Field[2];
26c4762a1bSJed Brown 
27c4762a1bSJed Brown typedef struct _User *User;
28c4762a1bSJed Brown struct _User {
29c4762a1bSJed Brown   PetscReal a[2];              /* Advection speeds */
30c4762a1bSJed Brown   PetscReal k[2];              /* Reaction coefficients */
31c4762a1bSJed Brown   PetscReal s[2];              /* Source terms */
32c4762a1bSJed Brown };
33c4762a1bSJed Brown 
34c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
35c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
36c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
37c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
38c4762a1bSJed Brown 
39c4762a1bSJed Brown int main(int argc,char **argv)
40c4762a1bSJed Brown {
41c4762a1bSJed Brown   TS                ts;         /* time integrator */
42c4762a1bSJed Brown   SNES              snes;       /* nonlinear solver */
43c4762a1bSJed Brown   SNESLineSearch    linesearch; /* line search */
44c4762a1bSJed Brown   Vec               X;          /* solution, residual vectors */
45c4762a1bSJed Brown   Mat               J;          /* Jacobian matrix */
46c4762a1bSJed Brown   PetscInt          steps,mx;
47c4762a1bSJed Brown   PetscErrorCode    ierr;
48c4762a1bSJed Brown   DM                da;
49c4762a1bSJed Brown   PetscReal         ftime,dt;
50c4762a1bSJed Brown   struct _User      user;       /* user-defined work context */
51c4762a1bSJed Brown   TSConvergedReason reason;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
57c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63c4762a1bSJed Brown      Extract global vectors from DMDA;
64c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* Initialize user application context */
68*1e1ea65dSPierre Jolivet   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");CHKERRQ(ierr);
69c4762a1bSJed Brown   {
70c4762a1bSJed Brown     user.a[0] = 1;           ierr = PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],NULL);CHKERRQ(ierr);
71c4762a1bSJed Brown     user.a[1] = 0;           ierr = PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],NULL);CHKERRQ(ierr);
72c4762a1bSJed Brown     user.k[0] = 1e6;         ierr = PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],NULL);CHKERRQ(ierr);
73c4762a1bSJed Brown     user.k[1] = 2*user.k[0]; ierr = PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],NULL);CHKERRQ(ierr);
74c4762a1bSJed Brown     user.s[0] = 0;           ierr = PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],NULL);CHKERRQ(ierr);
75c4762a1bSJed Brown     user.s[1] = 1;           ierr = PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],NULL);CHKERRQ(ierr);
76c4762a1bSJed Brown   }
77c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80c4762a1bSJed Brown      Create timestepping solver context
81c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
83c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
84c4762a1bSJed Brown   ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr);
85c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);CHKERRQ(ierr);
86c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,FormIFunction,&user);CHKERRQ(ierr);
87c4762a1bSJed Brown   ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr);
88c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&user);CHKERRQ(ierr);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* A line search in the nonlinear solve can fail due to ill-conditioning unless an absolute tolerance is set. Since
92c4762a1bSJed Brown    * this problem is linear, we deactivate the line search. For a linear problem, it is usually recommended to also use
93c4762a1bSJed Brown    * SNESSetType(snes,SNESKSPONLY). */
94c4762a1bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
95c4762a1bSJed Brown   ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
96c4762a1bSJed Brown   ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   ftime = .1;
99c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
100c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103c4762a1bSJed Brown      Set initial conditions
104c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105c4762a1bSJed Brown   ierr = FormInitialSolution(ts,X,&user);CHKERRQ(ierr);
106c4762a1bSJed Brown   ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
107c4762a1bSJed Brown   ierr = VecGetSize(X,&mx);CHKERRQ(ierr);
108c4762a1bSJed Brown   dt   = .1 * PetscMax(user.a[0],user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */
109c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112c4762a1bSJed Brown      Set runtime options
113c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117c4762a1bSJed Brown      Solve nonlinear system
118c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119c4762a1bSJed Brown   ierr = TSSolve(ts,X);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
121c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
122c4762a1bSJed Brown   ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr);
123c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);CHKERRQ(ierr);
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown      Free work space.
127c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
129c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
130c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
131c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
132c4762a1bSJed Brown   ierr = PetscFinalize();
133c4762a1bSJed Brown   return ierr;
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
137c4762a1bSJed Brown {
138c4762a1bSJed Brown   User           user = (User)ptr;
139c4762a1bSJed Brown   DM             da;
140c4762a1bSJed Brown   DMDALocalInfo  info;
141c4762a1bSJed Brown   PetscInt       i;
142c4762a1bSJed Brown   Field          *f;
143c4762a1bSJed Brown   const Field    *x,*xdot;
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 
150c4762a1bSJed Brown   /* Get pointers to vector data */
151c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,X,(void*)&x);CHKERRQ(ierr);
152c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Xdot,(void*)&xdot);CHKERRQ(ierr);
153c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
156c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
157c4762a1bSJed Brown     f[i][0] = xdot[i][0] + user->k[0]*x[i][0] - user->k[1]*x[i][1] - user->s[0];
158c4762a1bSJed Brown     f[i][1] = xdot[i][1] - user->k[0]*x[i][0] + user->k[1]*x[i][1] - user->s[1];
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* Restore vectors */
162c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,X,(void*)&x);CHKERRQ(ierr);
163c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
165c4762a1bSJed Brown   PetscFunctionReturn(0);
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
169c4762a1bSJed Brown {
170c4762a1bSJed Brown   User           user = (User)ptr;
171c4762a1bSJed Brown   DM             da;
172c4762a1bSJed Brown   Vec            Xloc;
173c4762a1bSJed Brown   DMDALocalInfo  info;
174c4762a1bSJed Brown   PetscInt       i,j;
175c4762a1bSJed Brown   PetscReal      hx;
176c4762a1bSJed Brown   Field          *f;
177c4762a1bSJed Brown   const Field    *x;
178c4762a1bSJed Brown   PetscErrorCode ierr;
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   PetscFunctionBeginUser;
181c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
182c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
183c4762a1bSJed Brown   hx   = 1.0/(PetscReal)info.mx;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /*
186c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
187c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
188c4762a1bSJed Brown      By placing code between these two statements, computations can be
189c4762a1bSJed Brown      done while messages are in transition.
190c4762a1bSJed Brown   */
191c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
192c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
193c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* Get pointers to vector data */
196c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Xloc,(void*)&x);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
200c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
201c4762a1bSJed Brown     const PetscReal *a  = user->a;
202c4762a1bSJed Brown     PetscReal       u0t[2];
203c4762a1bSJed Brown     u0t[0] = 1.0 - PetscPowRealInt(PetscSinReal(12*t),4);
204c4762a1bSJed Brown     u0t[1] = 0.0;
205c4762a1bSJed Brown     for (j=0; j<2; j++) {
206c4762a1bSJed Brown       if (i == 0)              f[i][j] = a[j]/hx*(1./3*u0t[j] + 0.5*x[i][j] - x[i+1][j] + 1./6*x[i+2][j]);
207c4762a1bSJed Brown       else if (i == 1)         f[i][j] = a[j]/hx*(-1./12*u0t[j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
208c4762a1bSJed Brown       else if (i == info.mx-2) f[i][j] = a[j]/hx*(-1./6*x[i-2][j] + x[i-1][j] - 0.5*x[i][j] - 1./3*x[i+1][j]);
209c4762a1bSJed Brown       else if (i == info.mx-1) f[i][j] = a[j]/hx*(-x[i][j] + x[i-1][j]);
210c4762a1bSJed Brown       else                     f[i][j] = a[j]/hx*(-1./12*x[i-2][j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
211c4762a1bSJed Brown     }
212c4762a1bSJed Brown   }
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* Restore vectors */
215c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Xloc,(void*)&x);CHKERRQ(ierr);
216c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
217c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
218c4762a1bSJed Brown   PetscFunctionReturn(0);
219c4762a1bSJed Brown }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown /* --------------------------------------------------------------------- */
222c4762a1bSJed Brown /*
223c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
224c4762a1bSJed Brown */
225c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
226c4762a1bSJed Brown {
227c4762a1bSJed Brown   User           user = (User)ptr;
228c4762a1bSJed Brown   PetscErrorCode ierr;
229c4762a1bSJed Brown   DMDALocalInfo  info;
230c4762a1bSJed Brown   PetscInt       i;
231c4762a1bSJed Brown   DM             da;
232c4762a1bSJed Brown   const 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 
238c4762a1bSJed Brown   /* Get pointers to vector data */
239c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,X,(void*)&x);CHKERRQ(ierr);
240c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Xdot,(void*)&xdot);CHKERRQ(ierr);
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
243c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
244c4762a1bSJed Brown     const PetscReal *k = user->k;
245c4762a1bSJed Brown     PetscScalar     v[2][2];
246c4762a1bSJed Brown 
247c4762a1bSJed Brown     v[0][0] = a + k[0]; v[0][1] =  -k[1];
248c4762a1bSJed Brown     v[1][0] =    -k[0]; v[1][1] = a+k[1];
249c4762a1bSJed Brown     ierr    = MatSetValuesBlocked(Jpre,1,&i,1,&i,&v[0][0],INSERT_VALUES);CHKERRQ(ierr);
250c4762a1bSJed Brown   }
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /* Restore vectors */
253c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,X,(void*)&x);CHKERRQ(ierr);
254c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot);CHKERRQ(ierr);
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257c4762a1bSJed Brown   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258c4762a1bSJed Brown   if (J != Jpre) {
259c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
261c4762a1bSJed Brown   }
262c4762a1bSJed Brown   PetscFunctionReturn(0);
263c4762a1bSJed Brown }
264c4762a1bSJed Brown 
265c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
266c4762a1bSJed Brown {
267c4762a1bSJed Brown   User           user = (User)ctx;
268c4762a1bSJed Brown   DM             da;
269c4762a1bSJed Brown   PetscInt       i;
270c4762a1bSJed Brown   DMDALocalInfo  info;
271c4762a1bSJed Brown   Field          *x;
272c4762a1bSJed Brown   PetscReal      hx;
273c4762a1bSJed Brown   PetscErrorCode ierr;
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   PetscFunctionBeginUser;
276c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
277c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
278c4762a1bSJed Brown   hx   = 1.0/(PetscReal)info.mx;
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   /* Get pointers to vector data */
281c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
284c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
285c4762a1bSJed Brown     PetscReal r = (i+1)*hx,ik = user->k[1] != 0.0 ? 1.0/user->k[1] : 1.0;
286c4762a1bSJed Brown     x[i][0] = 1 + user->s[1]*r;
287c4762a1bSJed Brown     x[i][1] = user->k[0]*ik*x[i][0] + user->s[1]*ik;
288c4762a1bSJed Brown   }
289c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
290c4762a1bSJed Brown   PetscFunctionReturn(0);
291c4762a1bSJed Brown }
292c4762a1bSJed Brown 
293c4762a1bSJed Brown /*TEST
294c4762a1bSJed Brown 
295c4762a1bSJed Brown     test:
296c4762a1bSJed Brown       args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_arkimex_type 4 -ts_adapt_type none -ts_dt .005 -ts_max_time .1
297c4762a1bSJed Brown       requires: !single
298c4762a1bSJed Brown 
299c4762a1bSJed Brown     test:
300c4762a1bSJed Brown       suffix: 2
301c4762a1bSJed Brown       args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_dt 1e-3 -ts_adapt_type none -ts_dt .005 -ts_max_time .1
302c4762a1bSJed Brown       nsize: 2
303c4762a1bSJed Brown 
304c4762a1bSJed Brown     test:
305c4762a1bSJed Brown       suffix: 3
306c4762a1bSJed Brown       args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type ra34pw2 -ts_dt 5e-3 -ts_max_time .1  -ts_adapt_type none
307c4762a1bSJed Brown       nsize: 2
308c4762a1bSJed Brown 
309c4762a1bSJed Brown     test:
310c4762a1bSJed Brown       suffix: 4
311c4762a1bSJed Brown       args: -ts_type eimex -da_grid_x 200 -ts_eimex_order_adapt true -ts_dt 0.001 -ts_monitor -ts_max_steps 100
312c4762a1bSJed Brown       filter: sed "s/ITS/TIME/g"
313c4762a1bSJed Brown       nsize: 2
314c4762a1bSJed Brown 
315c4762a1bSJed Brown TEST*/
316