xref: /petsc/src/ts/tutorials/phasefield/biharmonic2.c (revision 2a1887a77e7b2c6e00dd0ba96d1387c839460237)
1c4762a1bSJed Brown static char help[] = "Solves biharmonic equation in 1d.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   Solves the equation biharmonic equation in split form
5c4762a1bSJed Brown 
6c4762a1bSJed Brown     w = -kappa \Delta u
7c4762a1bSJed Brown     u_t =  \Delta w
8c4762a1bSJed Brown     -1  <= u <= 1
9c4762a1bSJed Brown     Periodic boundary conditions
10c4762a1bSJed Brown 
11c4762a1bSJed Brown Evolve the biharmonic heat equation with bounds:  (same as biharmonic)
12c4762a1bSJed Brown ---------------
13*188af4bfSBarry Smith ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution  -pc_type lu  -draw_pause .1 -snes_converged_reason  -ts_type beuler  -da_refine 5 -draw_fields 1 -ts_time_step 9.53674e-9
14c4762a1bSJed Brown 
15c4762a1bSJed Brown     w = -kappa \Delta u  + u^3  - u
16c4762a1bSJed Brown     u_t =  \Delta w
17c4762a1bSJed Brown     -1  <= u <= 1
18c4762a1bSJed Brown     Periodic boundary conditions
19c4762a1bSJed Brown 
20c4762a1bSJed Brown Evolve the Cahn-Hillard equations: (this fails after a few timesteps 12/17/2017)
21c4762a1bSJed Brown ---------------
22*188af4bfSBarry Smith ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution  -pc_type lu  -draw_pause .1 -snes_converged_reason   -ts_type beuler    -da_refine 6  -draw_fields 1  -kappa .00001 -ts_time_step 5.96046e-06 -cahn-hillard
23c4762a1bSJed Brown 
24c4762a1bSJed Brown */
25c4762a1bSJed Brown #include <petscdm.h>
26c4762a1bSJed Brown #include <petscdmda.h>
27c4762a1bSJed Brown #include <petscts.h>
28c4762a1bSJed Brown #include <petscdraw.h>
29c4762a1bSJed Brown 
30c4762a1bSJed Brown /*
31c4762a1bSJed Brown    User-defined routines
32c4762a1bSJed Brown */
33c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS, PetscReal, Vec, Vec, Vec, void *), FormInitialSolution(DM, Vec, PetscReal);
349371c9d4SSatish Balay typedef struct {
359371c9d4SSatish Balay   PetscBool cahnhillard;
369371c9d4SSatish Balay   PetscReal kappa;
379371c9d4SSatish Balay   PetscInt  energy;
389371c9d4SSatish Balay   PetscReal tol;
399371c9d4SSatish Balay   PetscReal theta;
409371c9d4SSatish Balay   PetscReal theta_c;
419371c9d4SSatish Balay } UserCtx;
42c4762a1bSJed Brown 
main(int argc,char ** argv)43d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
44d71ae5a4SJacob Faibussowitsch {
45c4762a1bSJed Brown   TS            ts;   /* nonlinear solver */
46c4762a1bSJed Brown   Vec           x, r; /* solution, residual vectors */
47c4762a1bSJed Brown   Mat           J;    /* Jacobian matrix */
48c4762a1bSJed Brown   PetscInt      steps, Mx;
49c4762a1bSJed Brown   DM            da;
50c4762a1bSJed Brown   MatFDColoring matfdcoloring;
51c4762a1bSJed Brown   ISColoring    iscoloring;
52c4762a1bSJed Brown   PetscReal     dt;
53c4762a1bSJed Brown   PetscReal     vbounds[] = {-100000, 100000, -1.1, 1.1};
54c4762a1bSJed Brown   SNES          snes;
55c4762a1bSJed Brown   UserCtx       ctx;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58c4762a1bSJed Brown      Initialize program
59c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60327415f7SBarry Smith   PetscFunctionBeginUser;
61c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
62c4762a1bSJed Brown   ctx.kappa = 1.0;
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", &ctx.kappa, NULL));
64c4762a1bSJed Brown   ctx.cahnhillard = PETSC_FALSE;
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-cahn-hillard", &ctx.cahnhillard, NULL));
679566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 2, vbounds));
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 600, 600));
69c4762a1bSJed Brown   ctx.energy = 1;
709566063dSJacob Faibussowitsch   /*PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL));*/
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-energy", &ctx.energy, NULL));
72c4762a1bSJed Brown   ctx.tol = 1.0e-8;
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &ctx.tol, NULL));
74c4762a1bSJed Brown   ctx.theta   = .001;
75c4762a1bSJed Brown   ctx.theta_c = 1.0;
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta", &ctx.theta, NULL));
779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta_c", &ctx.theta_c, NULL));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
81c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
829566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10, 2, 2, NULL, &da));
839566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
849566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
859566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "Biharmonic heat equation: w = -kappa*u_xx"));
869566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "Biharmonic heat equation: u"));
879566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
88c4762a1bSJed Brown   dt = 1.0 / (10. * ctx.kappa * Mx * Mx * Mx * Mx);
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
92c4762a1bSJed Brown      vectors that are the same types
93c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
949566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
959566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98c4762a1bSJed Brown      Create timestepping solver context
99c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1009566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1019566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
1029566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1039566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, FormFunction, &ctx));
1049566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, .02));
1059566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_INTERPOLATE));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
109c4762a1bSJed Brown 
110c4762a1bSJed Brown <     Set Jacobian matrix data structure and default Jacobian evaluation
111c4762a1bSJed Brown      routine. User can override with:
112c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
113c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
1147addb90fSBarry Smith      -snes_mf_operator : form matrix used to construct the preconditioner as set by the user,
115c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
116c4762a1bSJed Brown                          products within Newton-Krylov method
117c4762a1bSJed Brown 
118c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1199566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
1209566063dSJacob Faibussowitsch   PetscCall(DMCreateColoring(da, IS_COLORING_GLOBAL, &iscoloring));
1219566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da, MATAIJ));
1229566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da, &J));
1239566063dSJacob Faibussowitsch   PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring));
1242ba42892SBarry Smith   PetscCall(MatFDColoringSetFunction(matfdcoloring, (MatFDColoringFn *)SNESTSFormFunction, ts));
1259566063dSJacob Faibussowitsch   PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
1269566063dSJacob Faibussowitsch   PetscCall(MatFDColoringSetUp(J, iscoloring, matfdcoloring));
1279566063dSJacob Faibussowitsch   PetscCall(ISColoringDestroy(&iscoloring));
1289566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, matfdcoloring));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131c4762a1bSJed Brown      Customize nonlinear solver
132c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1339566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136c4762a1bSJed Brown      Set initial conditions
137c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1389566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(da, x, ctx.kappa));
1399566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
1409566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143c4762a1bSJed Brown      Set runtime options
144c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1459566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148c4762a1bSJed Brown      Solve nonlinear system
149c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1509566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
1519566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
155c4762a1bSJed Brown      are no longer needed.
156c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1589566063dSJacob Faibussowitsch   PetscCall(MatFDColoringDestroy(&matfdcoloring));
1599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
1619566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1629566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
165b122ec5aSJacob Faibussowitsch   return 0;
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
1689371c9d4SSatish Balay typedef struct {
1699371c9d4SSatish Balay   PetscScalar w, u;
1709371c9d4SSatish Balay } Field;
171c4762a1bSJed Brown /* ------------------------------------------------------------------- */
172c4762a1bSJed Brown /*
173c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
174c4762a1bSJed Brown 
175c4762a1bSJed Brown    Input Parameters:
176c4762a1bSJed Brown .  ts - the TS context
177c4762a1bSJed Brown .  X - input vector
178c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
179c4762a1bSJed Brown 
180c4762a1bSJed Brown    Output Parameter:
181c4762a1bSJed Brown .  F - function vector
182c4762a1bSJed Brown  */
FormFunction(TS ts,PetscReal ftime,Vec X,Vec Xdot,Vec F,void * ptr)183d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec Xdot, Vec F, void *ptr)
184d71ae5a4SJacob Faibussowitsch {
185c4762a1bSJed Brown   DM        da;
186c4762a1bSJed Brown   PetscInt  i, Mx, xs, xm;
187c4762a1bSJed Brown   PetscReal hx, sx;
188c4762a1bSJed Brown   Field    *x, *xdot, *f;
189c4762a1bSJed Brown   Vec       localX, localXdot;
190c4762a1bSJed Brown   UserCtx  *ctx = (UserCtx *)ptr;
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   PetscFunctionBegin;
1939566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1949566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
1959566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localXdot));
1969566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
197c4762a1bSJed Brown 
1989371c9d4SSatish Balay   hx = 1.0 / (PetscReal)Mx;
1999371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /*
202c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
203c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
204c4762a1bSJed Brown      By placing code between these two statements, computations can be
205c4762a1bSJed Brown      done while messages are in transition.
206c4762a1bSJed Brown   */
2079566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
2089566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
2099566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, Xdot, INSERT_VALUES, localXdot));
2109566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, Xdot, INSERT_VALUES, localXdot));
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /*
213c4762a1bSJed Brown      Get pointers to vector data
214c4762a1bSJed Brown   */
2159566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localX, &x));
2169566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localXdot, &xdot));
2179566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /*
220c4762a1bSJed Brown      Get local grid boundaries
221c4762a1bSJed Brown   */
2229566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   /*
225c4762a1bSJed Brown      Compute function over the locally owned part of the grid
226c4762a1bSJed Brown   */
227c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
228c4762a1bSJed Brown     f[i].w = x[i].w + ctx->kappa * (x[i - 1].u + x[i + 1].u - 2.0 * x[i].u) * sx;
229c4762a1bSJed Brown     if (ctx->cahnhillard) {
230c4762a1bSJed Brown       switch (ctx->energy) {
231d71ae5a4SJacob Faibussowitsch       case 1: /* double well */
232d71ae5a4SJacob Faibussowitsch         f[i].w += -x[i].u * x[i].u * x[i].u + x[i].u;
233d71ae5a4SJacob Faibussowitsch         break;
234d71ae5a4SJacob Faibussowitsch       case 2: /* double obstacle */
235d71ae5a4SJacob Faibussowitsch         f[i].w += x[i].u;
236d71ae5a4SJacob Faibussowitsch         break;
237c4762a1bSJed Brown       case 3: /* logarithmic */
238c4762a1bSJed Brown         if (PetscRealPart(x[i].u) < -1.0 + 2.0 * ctx->tol) f[i].w += .5 * ctx->theta * (-PetscLogReal(ctx->tol) + PetscLogScalar((1.0 - x[i].u) / 2.0)) + ctx->theta_c * x[i].u;
239c4762a1bSJed Brown         else if (PetscRealPart(x[i].u) > 1.0 - 2.0 * ctx->tol) f[i].w += .5 * ctx->theta * (-PetscLogScalar((1.0 + x[i].u) / 2.0) + PetscLogReal(ctx->tol)) + ctx->theta_c * x[i].u;
240c4762a1bSJed Brown         else f[i].w += .5 * ctx->theta * (-PetscLogScalar((1.0 + x[i].u) / 2.0) + PetscLogScalar((1.0 - x[i].u) / 2.0)) + ctx->theta_c * x[i].u;
241c4762a1bSJed Brown         break;
242c4762a1bSJed Brown       }
243c4762a1bSJed Brown     }
244c4762a1bSJed Brown     f[i].u = xdot[i].u - (x[i - 1].w + x[i + 1].w - 2.0 * x[i].w) * sx;
245c4762a1bSJed Brown   }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /*
248c4762a1bSJed Brown      Restore vectors
249c4762a1bSJed Brown   */
2509566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localXdot, &xdot));
2519566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
2529566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
2539566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
2549566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localXdot));
2553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown /* ------------------------------------------------------------------- */
FormInitialSolution(DM da,Vec X,PetscReal kappa)259d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(DM da, Vec X, PetscReal kappa)
260d71ae5a4SJacob Faibussowitsch {
261c4762a1bSJed Brown   PetscInt  i, xs, xm, Mx, xgs, xgm;
262c4762a1bSJed Brown   Field    *x;
263c4762a1bSJed Brown   PetscReal hx, xx, r, sx;
264c4762a1bSJed Brown   Vec       Xg;
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   PetscFunctionBegin;
2679566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   hx = 1.0 / (PetscReal)Mx;
270c4762a1bSJed Brown   sx = 1.0 / (hx * hx);
271c4762a1bSJed Brown 
272c4762a1bSJed Brown   /*
273c4762a1bSJed Brown      Get pointers to vector data
274c4762a1bSJed Brown   */
2759566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da, &Xg));
2769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xg, &x));
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   /*
279c4762a1bSJed Brown      Get local grid boundaries
280c4762a1bSJed Brown   */
2819566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
2829566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &xgs, NULL, NULL, &xgm, NULL, NULL));
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   /*
285c4762a1bSJed Brown      Compute u function over the locally owned part of the grid including ghost points
286c4762a1bSJed Brown   */
287c4762a1bSJed Brown   for (i = xgs; i < xgs + xgm; i++) {
288c4762a1bSJed Brown     xx = i * hx;
289c4762a1bSJed Brown     r  = PetscSqrtReal((xx - .5) * (xx - .5));
290c4762a1bSJed Brown     if (r < .125) x[i].u = 1.0;
291c4762a1bSJed Brown     else x[i].u = -.50;
292c4762a1bSJed Brown     /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */
293c4762a1bSJed Brown     x[i].w = 0;
294c4762a1bSJed Brown   }
295c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) x[i].w = -kappa * (x[i - 1].u + x[i + 1].u - 2.0 * x[i].u) * sx;
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   /*
298c4762a1bSJed Brown      Restore vectors
299c4762a1bSJed Brown   */
3009566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xg, &x));
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   /* Grab only the global part of the vector */
3039566063dSJacob Faibussowitsch   PetscCall(VecSet(X, 0));
3049566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da, Xg, ADD_VALUES, X));
3059566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da, Xg, ADD_VALUES, X));
3069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Xg));
3073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
308c4762a1bSJed Brown }
309c4762a1bSJed Brown 
310c4762a1bSJed Brown /*TEST
311c4762a1bSJed Brown 
312c4762a1bSJed Brown    build:
313c4762a1bSJed Brown      requires: !complex !single
314c4762a1bSJed Brown 
315c4762a1bSJed Brown    test:
316*188af4bfSBarry Smith      args: -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type beuler -da_refine 5 -ts_time_step 9.53674e-9 -ts_max_steps 50
317c4762a1bSJed Brown      requires: x
318c4762a1bSJed Brown 
319c4762a1bSJed Brown TEST*/
320