xref: /petsc/src/ts/tests/ex25.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
46*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
47*d71ae5a4SJacob Faibussowitsch {
48c4762a1bSJed Brown   PetscInt       cycle;
49c4762a1bSJed Brown   PetscErrorCode ierr;
50c4762a1bSJed Brown 
519371c9d4SSatish Balay   ierr = MPI_Init(&argc, &argv);
529371c9d4SSatish Balay   if (ierr) return ierr;
53c4762a1bSJed Brown   for (cycle = 0; cycle < 4; cycle++) {
54c4762a1bSJed Brown     ierr = Brusselator(argc, argv, cycle);
55c4762a1bSJed Brown     if (ierr) return 1;
56c4762a1bSJed Brown   }
57c4762a1bSJed Brown   ierr = MPI_Finalize();
58c4762a1bSJed Brown   return ierr;
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61*d71ae5a4SJacob Faibussowitsch PetscErrorCode Brusselator(int argc, char **argv, PetscInt cycle)
62*d71ae5a4SJacob Faibussowitsch {
63c4762a1bSJed Brown   TS                ts; /* nonlinear solver */
64c4762a1bSJed Brown   Vec               X;  /* solution, residual vectors */
65c4762a1bSJed Brown   Mat               J;  /* Jacobian matrix */
66c4762a1bSJed Brown   PetscInt          steps, mx;
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 
72327415f7SBarry Smith   PetscFunctionBeginUser;
739566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
77c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
789566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 11, 2, 2, NULL, &da));
799566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
809566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83c4762a1bSJed Brown      Extract global vectors from DMDA;
84c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
859566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &X));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Initialize user application context */
88d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Advection-reaction options", "");
89c4762a1bSJed Brown   {
90c4762a1bSJed Brown     user.A      = 1;
91c4762a1bSJed Brown     user.B      = 3;
92c4762a1bSJed Brown     user.alpha  = 0.1;
93c4762a1bSJed Brown     user.uleft  = 1;
94c4762a1bSJed Brown     user.uright = 1;
95c4762a1bSJed Brown     user.vleft  = 3;
96c4762a1bSJed Brown     user.vright = 3;
979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-A", "Reaction rate", "", user.A, &user.A, NULL));
989566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-B", "Reaction rate", "", user.B, &user.B, NULL));
999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-alpha", "Diffusion coefficient", "", user.alpha, &user.alpha, NULL));
1009566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-uleft", "Dirichlet boundary condition", "", user.uleft, &user.uleft, NULL));
1019566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-uright", "Dirichlet boundary condition", "", user.uright, &user.uright, NULL));
1029566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-vleft", "Dirichlet boundary condition", "", user.vleft, &user.vleft, NULL));
1039566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-vright", "Dirichlet boundary condition", "", user.vright, &user.vright, NULL));
104c4762a1bSJed Brown   }
105d0609cedSBarry Smith   PetscOptionsEnd();
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown      Create timestepping solver context
109c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1109566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1119566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
1129566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSARKIMEX));
1139566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, &user));
1149566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &user));
1159566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da, MATAIJ));
1169566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da, &J));
1179566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &user));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   ftime = 1.0;
1209566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
1219566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124c4762a1bSJed Brown      Set initial conditions
125c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1269566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts, X, &user));
1279566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
1289566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &mx));
129c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(mx / 2 - 1);
130c4762a1bSJed Brown   dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
131c4762a1bSJed Brown   dt *= PetscPowRealInt(0.2, cycle);    /* Shrink the time step in convergence study. */
1329566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
1339566063dSJacob Faibussowitsch   PetscCall(TSSetTolerances(ts, 1e-3 * PetscPowRealInt(0.5, cycle), NULL, 1e-3 * PetscPowRealInt(0.5, cycle), NULL));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136c4762a1bSJed Brown      Set runtime options
137c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1389566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown      Solve nonlinear system
142c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1439566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
1449566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
1459566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
1469566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
1479566063dSJacob Faibussowitsch   PetscCall(VecMin(X, NULL, &xmin));
1489566063dSJacob Faibussowitsch   PetscCall(VecMax(X, NULL, &xmax));
14963a3b9bcSJacob 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));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152c4762a1bSJed Brown      Free work space.
153c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1569566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1579566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
159b122ec5aSJacob Faibussowitsch   return 0;
160c4762a1bSJed Brown }
161c4762a1bSJed Brown 
162*d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ptr)
163*d71ae5a4SJacob Faibussowitsch {
164c4762a1bSJed Brown   User          user = (User)ptr;
165c4762a1bSJed Brown   DM            da;
166c4762a1bSJed Brown   DMDALocalInfo info;
167c4762a1bSJed Brown   PetscInt      i;
168c4762a1bSJed Brown   Field        *x, *xdot, *f;
169c4762a1bSJed Brown   PetscReal     hx;
170c4762a1bSJed Brown   Vec           Xloc;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   PetscFunctionBeginUser;
1739566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1749566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
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   */
1839566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
1849566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
1859566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* Get pointers to vector data */
1889566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, Xloc, &x));
1899566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, Xdot, &xdot));
1909566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
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 */
2079566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, Xloc, &x));
2089566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, Xdot, &xdot));
2099566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
2109566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
211c4762a1bSJed Brown   PetscFunctionReturn(0);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214*d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr)
215*d71ae5a4SJacob Faibussowitsch {
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 
223c4762a1bSJed Brown   PetscFunctionBeginUser;
2249566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2259566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
226c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(info.mx - 1);
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* Get pointers to vector data */
2299566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, X, &x));
2309566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
233c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
234c4762a1bSJed Brown     PetscScalar u = x[i].u, v = x[i].v;
235c4762a1bSJed Brown     f[i].u = hx * (user->A + u * u * v - (user->B + 1) * u);
236c4762a1bSJed Brown     f[i].v = hx * (user->B * u - u * u * v);
237c4762a1bSJed Brown   }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   /* Restore vectors */
2409566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, X, &x));
2419566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
242c4762a1bSJed Brown   PetscFunctionReturn(0);
243c4762a1bSJed Brown }
244c4762a1bSJed Brown 
245c4762a1bSJed Brown /* --------------------------------------------------------------------- */
246c4762a1bSJed Brown /*
247c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
248c4762a1bSJed Brown */
249*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ptr)
250*d71ae5a4SJacob Faibussowitsch {
251c4762a1bSJed Brown   User          user = (User)ptr;
252c4762a1bSJed Brown   DMDALocalInfo info;
253c4762a1bSJed Brown   PetscInt      i;
254c4762a1bSJed Brown   PetscReal     hx;
255c4762a1bSJed Brown   DM            da;
256c4762a1bSJed Brown   Field        *x, *xdot;
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   PetscFunctionBeginUser;
2599566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2609566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
261c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(info.mx - 1);
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /* Get pointers to vector data */
2649566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, X, &x));
2659566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, Xdot, &xdot));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
268c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
269c4762a1bSJed Brown     if (i == 0 || i == info.mx - 1) {
270c4762a1bSJed Brown       const PetscInt    row = i, col = i;
2719371c9d4SSatish Balay       const PetscScalar vals[2][2] = {
2729371c9d4SSatish Balay         {hx, 0 },
2739371c9d4SSatish Balay         {0,  hx}
2749371c9d4SSatish Balay       };
2759566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 1, &col, &vals[0][0], INSERT_VALUES));
276c4762a1bSJed Brown     } else {
277c4762a1bSJed Brown       const PetscInt    row = i, col[] = {i - 1, i, i + 1};
278c4762a1bSJed Brown       const PetscScalar dxxL = -user->alpha / hx, dxx0 = 2. * user->alpha / hx, dxxR = -user->alpha / hx;
2799371c9d4SSatish Balay       const PetscScalar vals[2][3][2] = {
2809371c9d4SSatish Balay         {{dxxL, 0}, {a * hx + dxx0, 0}, {dxxR, 0}},
2819371c9d4SSatish Balay         {{0, dxxL}, {0, a * hx + dxx0}, {0, dxxR}}
2829371c9d4SSatish Balay       };
2839566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 3, col, &vals[0][0][0], INSERT_VALUES));
284c4762a1bSJed Brown     }
285c4762a1bSJed Brown   }
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /* Restore vectors */
2889566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, X, &x));
2899566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, Xdot, &xdot));
290c4762a1bSJed Brown 
2919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
2929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
293c4762a1bSJed Brown   if (J != Jpre) {
2949566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
2959566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
296c4762a1bSJed Brown   }
297c4762a1bSJed Brown   PetscFunctionReturn(0);
298c4762a1bSJed Brown }
299c4762a1bSJed Brown 
300*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(TS ts, Vec X, void *ctx)
301*d71ae5a4SJacob Faibussowitsch {
302c4762a1bSJed Brown   User          user = (User)ctx;
303c4762a1bSJed Brown   DM            da;
304c4762a1bSJed Brown   PetscInt      i;
305c4762a1bSJed Brown   DMDALocalInfo info;
306c4762a1bSJed Brown   Field        *x;
307c4762a1bSJed Brown   PetscReal     hx;
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   PetscFunctionBeginUser;
3109566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
3119566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
312c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(info.mx - 1);
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   /* Get pointers to vector data */
3159566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
318c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
319c4762a1bSJed Brown     PetscReal xi = i * hx;
320c4762a1bSJed Brown     x[i].u       = user->uleft * (1. - xi) + user->uright * xi + PetscSinReal(2. * PETSC_PI * xi);
321c4762a1bSJed Brown     x[i].v       = user->vleft * (1. - xi) + user->vright * xi;
322c4762a1bSJed Brown   }
3239566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
324c4762a1bSJed Brown   PetscFunctionReturn(0);
325c4762a1bSJed Brown }
326c4762a1bSJed Brown 
327c4762a1bSJed Brown /*TEST
328c4762a1bSJed Brown 
329c4762a1bSJed Brown     test:
330c4762a1bSJed Brown       args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3
331c4762a1bSJed Brown 
332c4762a1bSJed Brown     test:
333c4762a1bSJed Brown       suffix: 2
334c4762a1bSJed Brown       args:   -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3
335c4762a1bSJed Brown 
336c4762a1bSJed Brown TEST*/
337