xref: /petsc/src/ts/tests/ex25.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
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   PetscErrorCode    ierr;
67c4762a1bSJed Brown   DM                da;
68c4762a1bSJed Brown   PetscReal         ftime,hx,dt,xmax,xmin;
69c4762a1bSJed Brown   struct _User      user;       /* user-defined work context */
70c4762a1bSJed Brown   TSConvergedReason reason;
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
76c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da);CHKERRQ(ierr);
78c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
79c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82c4762a1bSJed Brown      Extract global vectors from DMDA;
83c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr);
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Initialize user application context */
87*1e1ea65dSPierre Jolivet   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");CHKERRQ(ierr);
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;
96c4762a1bSJed Brown     ierr        = PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL);CHKERRQ(ierr);
97c4762a1bSJed Brown     ierr        = PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL);CHKERRQ(ierr);
98c4762a1bSJed Brown     ierr        = PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL);CHKERRQ(ierr);
99c4762a1bSJed Brown     ierr        = PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL);CHKERRQ(ierr);
100c4762a1bSJed Brown     ierr        = PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL);CHKERRQ(ierr);
101c4762a1bSJed Brown     ierr        = PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL);CHKERRQ(ierr);
102c4762a1bSJed Brown     ierr        = PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL);CHKERRQ(ierr);
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107c4762a1bSJed Brown      Create timestepping solver context
108c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
110c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr);
112c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,FormIFunction,&user);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr);
116c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&user);CHKERRQ(ierr);
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   ftime = 1.0;
119c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123c4762a1bSJed Brown      Set initial conditions
124c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125c4762a1bSJed Brown   ierr = FormInitialSolution(ts,X,&user);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = VecGetSize(X,&mx);CHKERRQ(ierr);
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. */
131c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
132c4762a1bSJed Brown   ierr = TSSetTolerances(ts,1e-3*PetscPowRealInt(0.5,cycle),NULL,1e-3*PetscPowRealInt(0.5,cycle),NULL);CHKERRQ(ierr);
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135c4762a1bSJed Brown      Set runtime options
136c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140c4762a1bSJed Brown      Solve nonlinear system
141c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142c4762a1bSJed Brown   ierr = TSSolve(ts,X);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr);
146c4762a1bSJed Brown   ierr = VecMin(X,NULL,&xmin);CHKERRQ(ierr);
147c4762a1bSJed Brown   ierr = VecMax(X,NULL,&xmax);CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after % 3D steps. Range [%6.4f,%6.4f]\n",TSConvergedReasons[reason],(double)ftime,steps,(double)xmin,(double)xmax);CHKERRQ(ierr);
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151c4762a1bSJed Brown      Free work space.
152c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
154c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
155c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
156c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
157c4762a1bSJed Brown   ierr = PetscFinalize();
158c4762a1bSJed Brown   return ierr;
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   PetscErrorCode ierr;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   PetscFunctionBeginUser;
173c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
174c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
175c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /*
178c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
179c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
180c4762a1bSJed Brown      By placing code between these two statements, computations can be
181c4762a1bSJed Brown      done while messages are in transition.
182c4762a1bSJed Brown   */
183c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
184c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
185c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* Get pointers to vector data */
188c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Xloc,&x);CHKERRQ(ierr);
189c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Xdot,&xdot);CHKERRQ(ierr);
190c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
193c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
194c4762a1bSJed Brown     if (i == 0) {
195c4762a1bSJed Brown       f[i].u = hx * (x[i].u - user->uleft);
196c4762a1bSJed Brown       f[i].v = hx * (x[i].v - user->vleft);
197c4762a1bSJed Brown     } else if (i == info.mx-1) {
198c4762a1bSJed Brown       f[i].u = hx * (x[i].u - user->uright);
199c4762a1bSJed Brown       f[i].v = hx * (x[i].v - user->vright);
200c4762a1bSJed Brown     } else {
201c4762a1bSJed Brown       f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx;
202c4762a1bSJed Brown       f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx;
203c4762a1bSJed Brown     }
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* Restore vectors */
207c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Xloc,&x);CHKERRQ(ierr);
208c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Xdot,&xdot);CHKERRQ(ierr);
209c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
210c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
211c4762a1bSJed Brown   PetscFunctionReturn(0);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
215c4762a1bSJed Brown {
216c4762a1bSJed Brown   User           user = (User)ptr;
217c4762a1bSJed Brown   DM             da;
218c4762a1bSJed Brown   DMDALocalInfo  info;
219c4762a1bSJed Brown   PetscInt       i;
220c4762a1bSJed Brown   PetscReal      hx;
221c4762a1bSJed Brown   Field          *x,*f;
222c4762a1bSJed Brown   PetscErrorCode ierr;
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   PetscFunctionBeginUser;
225c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
226c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
227c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* Get pointers to vector data */
230c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,X,&x);CHKERRQ(ierr);
231c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
234c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
235c4762a1bSJed Brown     PetscScalar u = x[i].u,v = x[i].v;
236c4762a1bSJed Brown     f[i].u = hx*(user->A + u*u*v - (user->B+1)*u);
237c4762a1bSJed Brown     f[i].v = hx*(user->B*u - u*u*v);
238c4762a1bSJed Brown   }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   /* Restore vectors */
241c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,X,&x);CHKERRQ(ierr);
242c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
243c4762a1bSJed Brown   PetscFunctionReturn(0);
244c4762a1bSJed Brown }
245c4762a1bSJed Brown 
246c4762a1bSJed Brown /* --------------------------------------------------------------------- */
247c4762a1bSJed Brown /*
248c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
249c4762a1bSJed Brown */
250c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
251c4762a1bSJed Brown {
252c4762a1bSJed Brown   User           user = (User)ptr;
253c4762a1bSJed Brown   PetscErrorCode ierr;
254c4762a1bSJed Brown   DMDALocalInfo  info;
255c4762a1bSJed Brown   PetscInt       i;
256c4762a1bSJed Brown   PetscReal      hx;
257c4762a1bSJed Brown   DM             da;
258c4762a1bSJed Brown   Field          *x,*xdot;
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   PetscFunctionBeginUser;
261c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
262c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
263c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /* Get pointers to vector data */
266c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,X,&x);CHKERRQ(ierr);
267c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Xdot,&xdot);CHKERRQ(ierr);
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
270c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
271c4762a1bSJed Brown     if (i == 0 || i == info.mx-1) {
272c4762a1bSJed Brown       const PetscInt    row        = i,col = i;
273c4762a1bSJed Brown       const PetscScalar vals[2][2] = {{hx,0},{0,hx}};
274c4762a1bSJed Brown       ierr = MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES);CHKERRQ(ierr);
275c4762a1bSJed Brown     } else {
276c4762a1bSJed Brown       const PetscInt    row           = i,col[] = {i-1,i,i+1};
277c4762a1bSJed Brown       const PetscScalar dxxL          = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
278c4762a1bSJed Brown       const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
279c4762a1bSJed Brown                                          {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
280c4762a1bSJed Brown       ierr = MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);CHKERRQ(ierr);
281c4762a1bSJed Brown     }
282c4762a1bSJed Brown   }
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   /* Restore vectors */
285c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,X,&x);CHKERRQ(ierr);
286c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Xdot,&xdot);CHKERRQ(ierr);
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
289c4762a1bSJed Brown   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
290c4762a1bSJed Brown   if (J != Jpre) {
291c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
292c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown   PetscFunctionReturn(0);
295c4762a1bSJed Brown }
296c4762a1bSJed Brown 
297c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
298c4762a1bSJed Brown {
299c4762a1bSJed Brown   User           user = (User)ctx;
300c4762a1bSJed Brown   DM             da;
301c4762a1bSJed Brown   PetscInt       i;
302c4762a1bSJed Brown   DMDALocalInfo  info;
303c4762a1bSJed Brown   Field          *x;
304c4762a1bSJed Brown   PetscReal      hx;
305c4762a1bSJed Brown   PetscErrorCode ierr;
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   PetscFunctionBeginUser;
308c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
309c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
310c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(info.mx-1);
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   /* Get pointers to vector data */
313c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
316c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
317c4762a1bSJed Brown     PetscReal xi = i*hx;
318c4762a1bSJed Brown     x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi);
319c4762a1bSJed Brown     x[i].v = user->vleft*(1.-xi) + user->vright*xi;
320c4762a1bSJed Brown   }
321c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
322c4762a1bSJed Brown   PetscFunctionReturn(0);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown 
325c4762a1bSJed Brown /*TEST
326c4762a1bSJed Brown 
327c4762a1bSJed Brown     test:
328c4762a1bSJed Brown       args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3
329c4762a1bSJed Brown 
330c4762a1bSJed Brown     test:
331c4762a1bSJed Brown       suffix: 2
332c4762a1bSJed Brown       args:   -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3
333c4762a1bSJed Brown 
334c4762a1bSJed Brown TEST*/
335c4762a1bSJed Brown 
336