1c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d. Modified from ex13.c for illustrating how to solve DAEs. \n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown u_t = uxx + uyy
4c4762a1bSJed Brown 0 < x < 1, 0 < y < 1;
5c4762a1bSJed Brown At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125
6c4762a1bSJed Brown u(x,y) = 0.0 if r >= .125
7c4762a1bSJed Brown
8c4762a1bSJed Brown Boundary conditions:
9c4762a1bSJed Brown Drichlet BC:
10c4762a1bSJed Brown At x=0, x=1, y=0, y=1: u = 0.0
11c4762a1bSJed Brown
12c4762a1bSJed Brown Neumann BC:
13c4762a1bSJed Brown At x=0, x=1: du(x,y,t)/dx = 0
14c4762a1bSJed Brown At y=0, y=1: du(x,y,t)/dy = 0
15c4762a1bSJed Brown
16c4762a1bSJed Brown mpiexec -n 2 ./ex15 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
17c4762a1bSJed Brown ./ex15 -da_grid_x 40 -da_grid_y 40 -draw_pause .1 -boundary 1 -ts_monitor_draw_solution
18c4762a1bSJed Brown ./ex15 -da_grid_x 40 -da_grid_y 40 -draw_pause .1 -boundary 1 -Jtype 2 -nstencilpts 9
19c4762a1bSJed Brown
20c4762a1bSJed Brown */
21c4762a1bSJed Brown
22c4762a1bSJed Brown #include <petscdm.h>
23c4762a1bSJed Brown #include <petscdmda.h>
24c4762a1bSJed Brown #include <petscts.h>
25c4762a1bSJed Brown
26c4762a1bSJed Brown /*
27c4762a1bSJed Brown User-defined data structures and routines
28c4762a1bSJed Brown */
29c4762a1bSJed Brown
30c4762a1bSJed Brown /* AppCtx: used by FormIFunction() and FormIJacobian() */
31c4762a1bSJed Brown typedef struct {
32c4762a1bSJed Brown DM da;
33c4762a1bSJed Brown PetscInt nstencilpts; /* number of stencil points: 5 or 9 */
34c4762a1bSJed Brown PetscReal c;
35c4762a1bSJed Brown PetscInt boundary; /* Type of boundary condition */
36c4762a1bSJed Brown PetscBool viewJacobian;
37c4762a1bSJed Brown } AppCtx;
38c4762a1bSJed Brown
39c4762a1bSJed Brown extern PetscErrorCode FormIFunction(TS, PetscReal, Vec, Vec, Vec, void *);
40c4762a1bSJed Brown extern PetscErrorCode FormIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
41c4762a1bSJed Brown extern PetscErrorCode FormInitialSolution(Vec, void *);
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 u, r; /* solution, residual vectors */
47c4762a1bSJed Brown Mat J, Jmf = NULL; /* Jacobian matrices */
48c4762a1bSJed Brown DM da;
49c4762a1bSJed Brown PetscReal dt;
50c4762a1bSJed Brown AppCtx user; /* user-defined work context */
51c4762a1bSJed Brown SNES snes;
52c4762a1bSJed Brown PetscInt Jtype; /* Jacobian type
53c4762a1bSJed Brown 0: user provide Jacobian;
54c4762a1bSJed Brown 1: slow finite difference;
55c4762a1bSJed Brown 2: fd with coloring; */
56c4762a1bSJed Brown
57327415f7SBarry Smith PetscFunctionBeginUser;
58c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
59c4762a1bSJed Brown /* Initialize user application context */
60c4762a1bSJed Brown user.da = NULL;
61c4762a1bSJed Brown user.nstencilpts = 5;
62c4762a1bSJed Brown user.c = -30.0;
63c4762a1bSJed Brown user.boundary = 0; /* 0: Drichlet BC; 1: Neumann BC */
64c4762a1bSJed Brown user.viewJacobian = PETSC_FALSE;
65c4762a1bSJed Brown
669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nstencilpts", &user.nstencilpts, NULL));
679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-boundary", &user.boundary, NULL));
689566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-viewJacobian", &user.viewJacobian));
69c4762a1bSJed Brown
70c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
72c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73c4762a1bSJed Brown if (user.nstencilpts == 5) {
749566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 11, 11, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
75c4762a1bSJed Brown } else if (user.nstencilpts == 9) {
769566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 11, 11, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
7763a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "nstencilpts %" PetscInt_FMT " is not supported", user.nstencilpts);
789566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
799566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
80c4762a1bSJed Brown user.da = da;
81c4762a1bSJed Brown
82c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83c4762a1bSJed Brown Extract global vectors from DMDA;
84c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
859566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &u));
869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r));
87c4762a1bSJed Brown
88c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89c4762a1bSJed Brown Create timestepping solver context
90c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
919566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
929566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
939566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER));
949566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
959566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, r, FormIFunction, &user));
969566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.0));
979566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
98c4762a1bSJed Brown
99c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100c4762a1bSJed Brown Set initial conditions
101c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1029566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(u, &user));
1039566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, u));
104c4762a1bSJed Brown dt = .01;
1059566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt));
106c4762a1bSJed Brown
107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown Set Jacobian evaluation routine
109c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1109566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATAIJ));
1119566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &J));
112c4762a1bSJed Brown Jtype = 0;
1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-Jtype", &Jtype, NULL));
114c4762a1bSJed Brown if (Jtype == 0) { /* use user provided Jacobian evaluation routine */
11563a3b9bcSJacob Faibussowitsch PetscCheck(user.nstencilpts == 5, PETSC_COMM_WORLD, PETSC_ERR_SUP, "user Jacobian routine FormIJacobian() does not support nstencilpts=%" PetscInt_FMT, user.nstencilpts);
1169566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &user));
117c4762a1bSJed Brown } else { /* use finite difference Jacobian J as preconditioner and '-snes_mf_operator' for Mat*vec */
1189566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes));
1199566063dSJacob Faibussowitsch PetscCall(MatCreateSNESMF(snes, &Jmf));
120c4762a1bSJed Brown if (Jtype == 1) { /* slow finite difference J; */
1219566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Jmf, J, SNESComputeJacobianDefault, NULL));
122c4762a1bSJed Brown } else if (Jtype == 2) { /* Use coloring to compute finite difference J efficiently */
1239566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Jmf, J, SNESComputeJacobianDefaultColor, 0));
124c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Jtype is not supported");
125c4762a1bSJed Brown }
126c4762a1bSJed Brown
127c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128c4762a1bSJed Brown Sets various TS parameters from user options
129c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1309566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
131c4762a1bSJed Brown
132c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133c4762a1bSJed Brown Solve nonlinear system
134c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1359566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u));
136c4762a1bSJed Brown
137c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown Free work space.
139c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jmf));
1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
1439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
1449566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1459566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
146c4762a1bSJed Brown
1479566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
148b122ec5aSJacob Faibussowitsch return 0;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown
151c4762a1bSJed Brown /* --------------------------------------------------------------------- */
152c4762a1bSJed Brown /*
153c4762a1bSJed Brown FormIFunction = Udot - RHSFunction
154c4762a1bSJed Brown */
FormIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,PetscCtx ctx)155*2a8381b2SBarry Smith PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx)
156d71ae5a4SJacob Faibussowitsch {
157c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx;
158300f1712SStefano Zampini DM da = user->da;
159c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym;
160c4762a1bSJed Brown PetscReal hx, hy, sx, sy;
161c4762a1bSJed Brown PetscScalar u, uxx, uyy, **uarray, **f, **udot;
162c4762a1bSJed Brown Vec localU;
163c4762a1bSJed Brown
164c4762a1bSJed Brown PetscFunctionBeginUser;
1659566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
1669566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
167c4762a1bSJed Brown
1689371c9d4SSatish Balay hx = 1.0 / (PetscReal)(Mx - 1);
1699371c9d4SSatish Balay sx = 1.0 / (hx * hx);
1709371c9d4SSatish Balay hy = 1.0 / (PetscReal)(My - 1);
1719371c9d4SSatish Balay sy = 1.0 / (hy * hy);
1723c633725SBarry Smith PetscCheck(user->nstencilpts != 9 || hx == hy, PETSC_COMM_WORLD, PETSC_ERR_SUP, "hx must equal hy when nstencilpts = 9 for this example");
173c4762a1bSJed Brown
174c4762a1bSJed Brown /*
175c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
176c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
177c4762a1bSJed Brown By placing code between these two statements, computations can be
178c4762a1bSJed Brown done while messages are in transition.
179c4762a1bSJed Brown */
1809566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
1819566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
182c4762a1bSJed Brown
183c4762a1bSJed Brown /* Get pointers to vector data */
1849566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &uarray));
1859566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
1869566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Udot, &udot));
187c4762a1bSJed Brown
188c4762a1bSJed Brown /* Get local grid boundaries */
1899566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
190c4762a1bSJed Brown
191c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
192c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
193c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
194c4762a1bSJed Brown /* Boundary conditions */
195c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
196c4762a1bSJed Brown if (user->boundary == 0) { /* Drichlet BC */
197c4762a1bSJed Brown f[j][i] = uarray[j][i]; /* F = U */
198c4762a1bSJed Brown } else { /* Neumann BC */
199c4762a1bSJed Brown if (i == 0 && j == 0) { /* SW corner */
200c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j + 1][i + 1];
201c4762a1bSJed Brown } else if (i == Mx - 1 && j == 0) { /* SE corner */
202c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j + 1][i - 1];
203c4762a1bSJed Brown } else if (i == 0 && j == My - 1) { /* NW corner */
204c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j - 1][i + 1];
205c4762a1bSJed Brown } else if (i == Mx - 1 && j == My - 1) { /* NE corner */
206c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j - 1][i - 1];
207c4762a1bSJed Brown } else if (i == 0) { /* Left */
208c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j][i + 1];
209c4762a1bSJed Brown } else if (i == Mx - 1) { /* Right */
210c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j][i - 1];
211c4762a1bSJed Brown } else if (j == 0) { /* Bottom */
212c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j + 1][i];
213c4762a1bSJed Brown } else if (j == My - 1) { /* Top */
214c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j - 1][i];
215c4762a1bSJed Brown }
216c4762a1bSJed Brown }
217c4762a1bSJed Brown } else { /* Interior */
218c4762a1bSJed Brown u = uarray[j][i];
219c4762a1bSJed Brown /* 5-point stencil */
220c4762a1bSJed Brown uxx = (-2.0 * u + uarray[j][i - 1] + uarray[j][i + 1]);
221c4762a1bSJed Brown uyy = (-2.0 * u + uarray[j - 1][i] + uarray[j + 1][i]);
222c4762a1bSJed Brown if (user->nstencilpts == 9) {
223c4762a1bSJed Brown /* 9-point stencil: assume hx=hy */
224c4762a1bSJed Brown uxx = 2.0 * uxx / 3.0 + (0.5 * (uarray[j - 1][i - 1] + uarray[j - 1][i + 1] + uarray[j + 1][i - 1] + uarray[j + 1][i + 1]) - 2.0 * u) / 6.0;
225c4762a1bSJed Brown uyy = 2.0 * uyy / 3.0 + (0.5 * (uarray[j - 1][i - 1] + uarray[j - 1][i + 1] + uarray[j + 1][i - 1] + uarray[j + 1][i + 1]) - 2.0 * u) / 6.0;
226c4762a1bSJed Brown }
227c4762a1bSJed Brown f[j][i] = udot[j][i] - (uxx * sx + uyy * sy);
228c4762a1bSJed Brown }
229c4762a1bSJed Brown }
230c4762a1bSJed Brown }
231c4762a1bSJed Brown
232c4762a1bSJed Brown /* Restore vectors */
2339566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &uarray));
2349566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
2359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Udot, &udot));
2369566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
2379566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(11.0 * ym * xm));
2383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
239c4762a1bSJed Brown }
240c4762a1bSJed Brown
241c4762a1bSJed Brown /* --------------------------------------------------------------------- */
242c4762a1bSJed Brown /*
243c4762a1bSJed Brown FormIJacobian() - Compute IJacobian = dF/dU + a dF/dUdot
244c4762a1bSJed Brown This routine is not used with option '-use_coloring'
245c4762a1bSJed Brown */
FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,PetscCtx ctx)246*2a8381b2SBarry Smith PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat J, Mat Jpre, PetscCtx ctx)
247d71ae5a4SJacob Faibussowitsch {
248c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym, nc;
249c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx;
250300f1712SStefano Zampini DM da = user->da;
251c4762a1bSJed Brown MatStencil col[5], row;
252c4762a1bSJed Brown PetscScalar vals[5], hx, hy, sx, sy;
253c4762a1bSJed Brown
254c4762a1bSJed Brown PetscFunctionBeginUser;
2559566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
2569566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
257c4762a1bSJed Brown
2589371c9d4SSatish Balay hx = 1.0 / (PetscReal)(Mx - 1);
2599371c9d4SSatish Balay sx = 1.0 / (hx * hx);
2609371c9d4SSatish Balay hy = 1.0 / (PetscReal)(My - 1);
2619371c9d4SSatish Balay sy = 1.0 / (hy * hy);
262c4762a1bSJed Brown
263c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
264c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
265c4762a1bSJed Brown nc = 0;
2669371c9d4SSatish Balay row.j = j;
2679371c9d4SSatish Balay row.i = i;
268c4762a1bSJed Brown if (user->boundary == 0 && (i == 0 || i == Mx - 1 || j == 0 || j == My - 1)) {
2699371c9d4SSatish Balay col[nc].j = j;
2709371c9d4SSatish Balay col[nc].i = i;
2719371c9d4SSatish Balay vals[nc++] = 1.0;
272c4762a1bSJed Brown
273c4762a1bSJed Brown } else if (user->boundary > 0 && i == 0) { /* Left Neumann */
2749371c9d4SSatish Balay col[nc].j = j;
2759371c9d4SSatish Balay col[nc].i = i;
2769371c9d4SSatish Balay vals[nc++] = 1.0;
2779371c9d4SSatish Balay col[nc].j = j;
2789371c9d4SSatish Balay col[nc].i = i + 1;
2799371c9d4SSatish Balay vals[nc++] = -1.0;
280c4762a1bSJed Brown } else if (user->boundary > 0 && i == Mx - 1) { /* Right Neumann */
2819371c9d4SSatish Balay col[nc].j = j;
2829371c9d4SSatish Balay col[nc].i = i;
2839371c9d4SSatish Balay vals[nc++] = 1.0;
2849371c9d4SSatish Balay col[nc].j = j;
2859371c9d4SSatish Balay col[nc].i = i - 1;
2869371c9d4SSatish Balay vals[nc++] = -1.0;
287c4762a1bSJed Brown } else if (user->boundary > 0 && j == 0) { /* Bottom Neumann */
2889371c9d4SSatish Balay col[nc].j = j;
2899371c9d4SSatish Balay col[nc].i = i;
2909371c9d4SSatish Balay vals[nc++] = 1.0;
2919371c9d4SSatish Balay col[nc].j = j + 1;
2929371c9d4SSatish Balay col[nc].i = i;
2939371c9d4SSatish Balay vals[nc++] = -1.0;
294c4762a1bSJed Brown } else if (user->boundary > 0 && j == My - 1) { /* Top Neumann */
2959371c9d4SSatish Balay col[nc].j = j;
2969371c9d4SSatish Balay col[nc].i = i;
2979371c9d4SSatish Balay vals[nc++] = 1.0;
2989371c9d4SSatish Balay col[nc].j = j - 1;
2999371c9d4SSatish Balay col[nc].i = i;
3009371c9d4SSatish Balay vals[nc++] = -1.0;
301c4762a1bSJed Brown } else { /* Interior */
3029371c9d4SSatish Balay col[nc].j = j - 1;
3039371c9d4SSatish Balay col[nc].i = i;
3049371c9d4SSatish Balay vals[nc++] = -sy;
3059371c9d4SSatish Balay col[nc].j = j;
3069371c9d4SSatish Balay col[nc].i = i - 1;
3079371c9d4SSatish Balay vals[nc++] = -sx;
3089371c9d4SSatish Balay col[nc].j = j;
3099371c9d4SSatish Balay col[nc].i = i;
3109371c9d4SSatish Balay vals[nc++] = 2.0 * (sx + sy) + a;
3119371c9d4SSatish Balay col[nc].j = j;
3129371c9d4SSatish Balay col[nc].i = i + 1;
3139371c9d4SSatish Balay vals[nc++] = -sx;
3149371c9d4SSatish Balay col[nc].j = j + 1;
3159371c9d4SSatish Balay col[nc].i = i;
3169371c9d4SSatish Balay vals[nc++] = -sy;
317c4762a1bSJed Brown }
3189566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, vals, INSERT_VALUES));
319c4762a1bSJed Brown }
320c4762a1bSJed Brown }
3219566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
3229566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
323c4762a1bSJed Brown if (J != Jpre) {
3249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
3259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
326c4762a1bSJed Brown }
327c4762a1bSJed Brown
328c4762a1bSJed Brown if (user->viewJacobian) {
3299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)Jpre), "Jpre:\n"));
3309566063dSJacob Faibussowitsch PetscCall(MatView(Jpre, PETSC_VIEWER_STDOUT_WORLD));
331c4762a1bSJed Brown }
3323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
333c4762a1bSJed Brown }
334c4762a1bSJed Brown
335c4762a1bSJed Brown /* ------------------------------------------------------------------- */
FormInitialSolution(Vec U,void * ptr)336d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(Vec U, void *ptr)
337d71ae5a4SJacob Faibussowitsch {
338c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
339c4762a1bSJed Brown DM da = user->da;
340c4762a1bSJed Brown PetscReal c = user->c;
341c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My;
342c4762a1bSJed Brown PetscScalar **u;
343c4762a1bSJed Brown PetscReal hx, hy, x, y, r;
344c4762a1bSJed Brown
345c4762a1bSJed Brown PetscFunctionBeginUser;
3469566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
347c4762a1bSJed Brown
348c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1);
349c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1);
350c4762a1bSJed Brown
351c4762a1bSJed Brown /* Get pointers to vector data */
3529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
353c4762a1bSJed Brown
354c4762a1bSJed Brown /* Get local grid boundaries */
3559566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
356c4762a1bSJed Brown
357c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
358c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
359c4762a1bSJed Brown y = j * hy;
360c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
361c4762a1bSJed Brown x = i * hx;
362c4762a1bSJed Brown r = PetscSqrtReal((x - .5) * (x - .5) + (y - .5) * (y - .5));
363c4762a1bSJed Brown if (r < .125) u[j][i] = PetscExpReal(c * r * r * r);
364c4762a1bSJed Brown else u[j][i] = 0.0;
365c4762a1bSJed Brown }
366c4762a1bSJed Brown }
367c4762a1bSJed Brown
368c4762a1bSJed Brown /* Restore vectors */
3699566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
3703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
371c4762a1bSJed Brown }
372c4762a1bSJed Brown
373c4762a1bSJed Brown /*TEST
374c4762a1bSJed Brown
375c4762a1bSJed Brown test:
376c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -ts_monitor
377c4762a1bSJed Brown
378c4762a1bSJed Brown test:
379c4762a1bSJed Brown suffix: 2
380c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 2 -ts_monitor
381c4762a1bSJed Brown
382c4762a1bSJed Brown test:
383c4762a1bSJed Brown suffix: 3
384c4762a1bSJed Brown requires: !single
385c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor
386c4762a1bSJed Brown
387c4762a1bSJed Brown test:
388c4762a1bSJed Brown suffix: 4
389c4762a1bSJed Brown requires: !single
390c4762a1bSJed Brown nsize: 2
391c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor
392c4762a1bSJed Brown
393c4762a1bSJed Brown test:
394c4762a1bSJed Brown suffix: 5
395c4762a1bSJed Brown nsize: 1
396c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 1 -ts_monitor
397c4762a1bSJed Brown
398c4762a1bSJed Brown TEST*/
399