xref: /petsc/src/ts/tests/ex25.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static const char help[] = "Call PetscInitialize multiple times.\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    This example is based on the Brusselator tutorial of the same name, but tests multiple calls to PetscInitialize().
4c4762a1bSJed Brown    This is a bad "convergence study" because it only compares min and max values of the solution rather than comparing
5c4762a1bSJed Brown    norms of the errors.  For convergence studies, we recommend invoking PetscInitialize() only once and comparing norms
6c4762a1bSJed Brown    of errors (perhaps estimated using an accurate reference solution).
7c4762a1bSJed Brown 
8c4762a1bSJed Brown    Time-dependent Brusselator reaction-diffusion PDE in 1d. Demonstrates IMEX methods and multiple solves.
9c4762a1bSJed Brown 
10c4762a1bSJed Brown    u_t - alpha u_xx = A + u^2 v - (B+1) u
11c4762a1bSJed Brown    v_t - alpha v_xx = B u - u^2 v
12c4762a1bSJed Brown    0 < x < 1;
13c4762a1bSJed Brown    A = 1, B = 3, alpha = 1/10
14c4762a1bSJed Brown 
15c4762a1bSJed Brown    Initial conditions:
16c4762a1bSJed Brown    u(x,0) = 1 + sin(2 pi x)
17c4762a1bSJed Brown    v(x,0) = 3
18c4762a1bSJed Brown 
19c4762a1bSJed Brown    Boundary conditions:
20c4762a1bSJed Brown    u(0,t) = u(1,t) = 1
21c4762a1bSJed Brown    v(0,t) = v(1,t) = 3
22c4762a1bSJed Brown */
23c4762a1bSJed Brown 
24c4762a1bSJed Brown #include <petscdm.h>
25c4762a1bSJed Brown #include <petscdmda.h>
26c4762a1bSJed Brown #include <petscts.h>
27c4762a1bSJed Brown 
28c4762a1bSJed Brown typedef struct {
29c4762a1bSJed Brown   PetscScalar u,v;
30c4762a1bSJed Brown } Field;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown typedef struct _User *User;
33c4762a1bSJed Brown struct _User {
34c4762a1bSJed Brown   PetscReal A,B;                /* Reaction coefficients */
35c4762a1bSJed Brown   PetscReal alpha;              /* Diffusion coefficient */
36c4762a1bSJed Brown   PetscReal uleft,uright;       /* Dirichlet boundary conditions */
37c4762a1bSJed Brown   PetscReal vleft,vright;       /* Dirichlet boundary conditions */
38c4762a1bSJed Brown };
39c4762a1bSJed Brown 
40c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
41c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
42c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
43c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
44c4762a1bSJed Brown static int Brusselator(int,char**,PetscInt);
45c4762a1bSJed Brown 
46c4762a1bSJed Brown int main(int argc,char **argv)
47c4762a1bSJed Brown {
48c4762a1bSJed Brown   PetscInt       cycle;
49c4762a1bSJed Brown   PetscErrorCode ierr;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   ierr = MPI_Init(&argc,&argv);if (ierr) return ierr;
52c4762a1bSJed Brown   for (cycle=0; cycle<4; cycle++) {
53c4762a1bSJed Brown     ierr = Brusselator(argc,argv,cycle);
54c4762a1bSJed Brown     if (ierr) return 1;
55c4762a1bSJed Brown   }
56c4762a1bSJed Brown   ierr = MPI_Finalize();
57c4762a1bSJed Brown   return ierr;
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown PetscErrorCode Brusselator(int argc,char **argv,PetscInt cycle)
61c4762a1bSJed Brown {
62c4762a1bSJed Brown   TS                ts;         /* nonlinear solver */
63c4762a1bSJed Brown   Vec               X;          /* solution, residual vectors */
64c4762a1bSJed Brown   Mat               J;          /* Jacobian matrix */
65c4762a1bSJed Brown   PetscInt          steps,mx;
66c4762a1bSJed Brown   DM                da;
67c4762a1bSJed Brown   PetscReal         ftime,hx,dt,xmax,xmin;
68c4762a1bSJed Brown   struct _User      user;       /* user-defined work context */
69c4762a1bSJed Brown   TSConvergedReason reason;
70c4762a1bSJed Brown 
71*327415f7SBarry Smith   PetscFunctionBeginUser;
729566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
76c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
779566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da));
789566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
799566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82c4762a1bSJed Brown      Extract global vectors from DMDA;
83c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
849566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&X));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Initialize user application context */
87d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
88c4762a1bSJed Brown   {
89c4762a1bSJed Brown     user.A      = 1;
90c4762a1bSJed Brown     user.B      = 3;
91c4762a1bSJed Brown     user.alpha  = 0.1;
92c4762a1bSJed Brown     user.uleft  = 1;
93c4762a1bSJed Brown     user.uright = 1;
94c4762a1bSJed Brown     user.vleft  = 3;
95c4762a1bSJed Brown     user.vright = 3;
969566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL));
979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL));
989566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL));
999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL));
1009566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL));
1019566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL));
1029566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL));
103c4762a1bSJed Brown   }
104d0609cedSBarry Smith   PetscOptionsEnd();
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107c4762a1bSJed Brown      Create timestepping solver context
108c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1099566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1109566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts,da));
1119566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSARKIMEX));
1129566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user));
1139566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,FormIFunction,&user));
1149566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da,MATAIJ));
1159566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da,&J));
1169566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,J,J,FormIJacobian,&user));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   ftime = 1.0;
1199566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
1209566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123c4762a1bSJed Brown      Set initial conditions
124c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1259566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts,X,&user));
1269566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,X));
1279566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X,&mx));
128c4762a1bSJed Brown   hx = 1.0/(PetscReal)(mx/2-1);
129c4762a1bSJed Brown   dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
130c4762a1bSJed Brown   dt *= PetscPowRealInt(0.2,cycle);     /* Shrink the time step in convergence study. */
1319566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
1329566063dSJacob Faibussowitsch   PetscCall(TSSetTolerances(ts,1e-3*PetscPowRealInt(0.5,cycle),NULL,1e-3*PetscPowRealInt(0.5,cycle),NULL));
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135c4762a1bSJed Brown      Set runtime options
136c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1379566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140c4762a1bSJed Brown      Solve nonlinear system
141c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1429566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,X));
1439566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
1449566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
1459566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts,&reason));
1469566063dSJacob Faibussowitsch   PetscCall(VecMin(X,NULL,&xmin));
1479566063dSJacob Faibussowitsch   PetscCall(VecMax(X,NULL,&xmax));
14863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after % 3" PetscInt_FMT " steps. Range [%6.4f,%6.4f]\n",TSConvergedReasons[reason],(double)ftime,steps,(double)xmin,(double)xmax));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151c4762a1bSJed Brown      Free work space.
152c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1559566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1569566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
158b122ec5aSJacob Faibussowitsch   return 0;
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
162c4762a1bSJed Brown {
163c4762a1bSJed Brown   User           user = (User)ptr;
164c4762a1bSJed Brown   DM             da;
165c4762a1bSJed Brown   DMDALocalInfo  info;
166c4762a1bSJed Brown   PetscInt       i;
167c4762a1bSJed Brown   Field          *x,*xdot,*f;
168c4762a1bSJed Brown   PetscReal      hx;
169c4762a1bSJed Brown   Vec            Xloc;
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   PetscFunctionBeginUser;
1729566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
1739566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
174c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /*
177c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
178c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
179c4762a1bSJed Brown      By placing code between these two statements, computations can be
180c4762a1bSJed Brown      done while messages are in transition.
181c4762a1bSJed Brown   */
1829566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&Xloc));
1839566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
1849566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* Get pointers to vector data */
1879566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Xloc,&x));
1889566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Xdot,&xdot));
1899566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
192c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
193c4762a1bSJed Brown     if (i == 0) {
194c4762a1bSJed Brown       f[i].u = hx * (x[i].u - user->uleft);
195c4762a1bSJed Brown       f[i].v = hx * (x[i].v - user->vleft);
196c4762a1bSJed Brown     } else if (i == info.mx-1) {
197c4762a1bSJed Brown       f[i].u = hx * (x[i].u - user->uright);
198c4762a1bSJed Brown       f[i].v = hx * (x[i].v - user->vright);
199c4762a1bSJed Brown     } else {
200c4762a1bSJed Brown       f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx;
201c4762a1bSJed Brown       f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx;
202c4762a1bSJed Brown     }
203c4762a1bSJed Brown   }
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /* Restore vectors */
2069566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Xloc,&x));
2079566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Xdot,&xdot));
2089566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
2099566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&Xloc));
210c4762a1bSJed Brown   PetscFunctionReturn(0);
211c4762a1bSJed Brown }
212c4762a1bSJed Brown 
213c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
214c4762a1bSJed Brown {
215c4762a1bSJed Brown   User           user = (User)ptr;
216c4762a1bSJed Brown   DM             da;
217c4762a1bSJed Brown   DMDALocalInfo  info;
218c4762a1bSJed Brown   PetscInt       i;
219c4762a1bSJed Brown   PetscReal      hx;
220c4762a1bSJed Brown   Field          *x,*f;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   PetscFunctionBeginUser;
2239566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
2249566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
225c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /* Get pointers to vector data */
2289566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,X,&x));
2299566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
232c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
233c4762a1bSJed Brown     PetscScalar u = x[i].u,v = x[i].v;
234c4762a1bSJed Brown     f[i].u = hx*(user->A + u*u*v - (user->B+1)*u);
235c4762a1bSJed Brown     f[i].v = hx*(user->B*u - u*u*v);
236c4762a1bSJed Brown   }
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   /* Restore vectors */
2399566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,X,&x));
2409566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
241c4762a1bSJed Brown   PetscFunctionReturn(0);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown /* --------------------------------------------------------------------- */
245c4762a1bSJed Brown /*
246c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
247c4762a1bSJed Brown */
248c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
249c4762a1bSJed Brown {
250c4762a1bSJed Brown   User           user = (User)ptr;
251c4762a1bSJed Brown   DMDALocalInfo  info;
252c4762a1bSJed Brown   PetscInt       i;
253c4762a1bSJed Brown   PetscReal      hx;
254c4762a1bSJed Brown   DM             da;
255c4762a1bSJed Brown   Field          *x,*xdot;
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   PetscFunctionBeginUser;
2589566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
2599566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
260c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /* Get pointers to vector data */
2639566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,X,&x));
2649566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Xdot,&xdot));
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
267c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
268c4762a1bSJed Brown     if (i == 0 || i == info.mx-1) {
269c4762a1bSJed Brown       const PetscInt    row        = i,col = i;
270c4762a1bSJed Brown       const PetscScalar vals[2][2] = {{hx,0},{0,hx}};
2719566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES));
272c4762a1bSJed Brown     } else {
273c4762a1bSJed Brown       const PetscInt    row           = i,col[] = {i-1,i,i+1};
274c4762a1bSJed Brown       const PetscScalar dxxL          = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
275c4762a1bSJed Brown       const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
276c4762a1bSJed Brown                                          {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
2779566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES));
278c4762a1bSJed Brown     }
279c4762a1bSJed Brown   }
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   /* Restore vectors */
2829566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,X,&x));
2839566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Xdot,&xdot));
284c4762a1bSJed Brown 
2859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY));
2869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY));
287c4762a1bSJed Brown   if (J != Jpre) {
2889566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
2899566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
290c4762a1bSJed Brown   }
291c4762a1bSJed Brown   PetscFunctionReturn(0);
292c4762a1bSJed Brown }
293c4762a1bSJed Brown 
294c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
295c4762a1bSJed Brown {
296c4762a1bSJed Brown   User           user = (User)ctx;
297c4762a1bSJed Brown   DM             da;
298c4762a1bSJed Brown   PetscInt       i;
299c4762a1bSJed Brown   DMDALocalInfo  info;
300c4762a1bSJed Brown   Field          *x;
301c4762a1bSJed Brown   PetscReal      hx;
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   PetscFunctionBeginUser;
3049566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
3059566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
306c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /* Get pointers to vector data */
3099566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,X,&x));
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
312c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
313c4762a1bSJed Brown     PetscReal xi = i*hx;
314c4762a1bSJed Brown     x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi);
315c4762a1bSJed Brown     x[i].v = user->vleft*(1.-xi) + user->vright*xi;
316c4762a1bSJed Brown   }
3179566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,X,&x));
318c4762a1bSJed Brown   PetscFunctionReturn(0);
319c4762a1bSJed Brown }
320c4762a1bSJed Brown 
321c4762a1bSJed Brown /*TEST
322c4762a1bSJed Brown 
323c4762a1bSJed Brown     test:
324c4762a1bSJed Brown       args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3
325c4762a1bSJed Brown 
326c4762a1bSJed Brown     test:
327c4762a1bSJed Brown       suffix: 2
328c4762a1bSJed Brown       args:   -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3
329c4762a1bSJed Brown 
330c4762a1bSJed Brown TEST*/
331