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