xref: /petsc/src/ts/tests/ex4.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown        The Problem:
3c4762a1bSJed Brown            Solve the convection-diffusion equation:
4c4762a1bSJed Brown 
5c4762a1bSJed Brown              u_t+a*(u_x+u_y)=epsilon*(u_xx+u_yy)
6c4762a1bSJed Brown              u=0   at x=0, y=0
7c4762a1bSJed Brown              u_x=0 at x=1
8c4762a1bSJed Brown              u_y=0 at y=1
9c4762a1bSJed Brown              u = exp(-20.0*(pow(x-0.5,2.0)+pow(y-0.5,2.0))) at t=0
10c4762a1bSJed Brown 
11c4762a1bSJed Brown        This program tests the routine of computing the Jacobian by the
12c4762a1bSJed Brown        finite difference method as well as PETSc.
13c4762a1bSJed Brown 
14c4762a1bSJed Brown */
15c4762a1bSJed Brown 
16c4762a1bSJed Brown static char help[] = "Solve the convection-diffusion equation. \n\n";
17c4762a1bSJed Brown 
18c4762a1bSJed Brown #include <petscts.h>
19c4762a1bSJed Brown 
209371c9d4SSatish Balay typedef struct {
21c4762a1bSJed Brown   PetscInt  m;       /* the number of mesh points in x-direction */
22c4762a1bSJed Brown   PetscInt  n;       /* the number of mesh points in y-direction */
23c4762a1bSJed Brown   PetscReal dx;      /* the grid space in x-direction */
24c4762a1bSJed Brown   PetscReal dy;      /* the grid space in y-direction */
25c4762a1bSJed Brown   PetscReal a;       /* the convection coefficient    */
26c4762a1bSJed Brown   PetscReal epsilon; /* the diffusion coefficient     */
27c4762a1bSJed Brown   PetscReal tfinal;
28c4762a1bSJed Brown } Data;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
31c4762a1bSJed Brown extern PetscErrorCode Initial(Vec, void *);
32c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
33c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
34c4762a1bSJed Brown extern PetscErrorCode PostStep(TS);
35c4762a1bSJed Brown 
369371c9d4SSatish Balay int main(int argc, char **argv) {
37c4762a1bSJed Brown   PetscInt      time_steps = 100, iout, NOUT = 1;
38c4762a1bSJed Brown   Vec           global;
39c4762a1bSJed Brown   PetscReal     dt, ftime, ftime_original;
40c4762a1bSJed Brown   TS            ts;
41c4762a1bSJed Brown   PetscViewer   viewfile;
42c4762a1bSJed Brown   Mat           J = 0;
43c4762a1bSJed Brown   Vec           x;
44c4762a1bSJed Brown   Data          data;
45c4762a1bSJed Brown   PetscInt      mn;
46c4762a1bSJed Brown   PetscBool     flg;
47c4762a1bSJed Brown   MatColoring   mc;
48c4762a1bSJed Brown   ISColoring    iscoloring;
49c4762a1bSJed Brown   MatFDColoring matfdcoloring        = 0;
50c4762a1bSJed Brown   PetscBool     fd_jacobian_coloring = PETSC_FALSE;
51c4762a1bSJed Brown   SNES          snes;
52c4762a1bSJed Brown   KSP           ksp;
53c4762a1bSJed Brown   PC            pc;
54c4762a1bSJed Brown 
55327415f7SBarry Smith   PetscFunctionBeginUser;
569566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* set data */
59c4762a1bSJed Brown   data.m       = 9;
60c4762a1bSJed Brown   data.n       = 9;
61c4762a1bSJed Brown   data.a       = 1.0;
62c4762a1bSJed Brown   data.epsilon = 0.1;
63c4762a1bSJed Brown   data.dx      = 1.0 / (data.m + 1.0);
64c4762a1bSJed Brown   data.dy      = 1.0 / (data.n + 1.0);
65c4762a1bSJed Brown   mn           = (data.m) * (data.n);
669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-time", &time_steps, NULL));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* set initial conditions */
699566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &global));
709566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(global, PETSC_DECIDE, mn));
719566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(global));
729566063dSJacob Faibussowitsch   PetscCall(Initial(global, &data));
739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(global, &x));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* create timestep context */
769566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
779566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, Monitor, &data, NULL));
789566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSEULER));
79c4762a1bSJed Brown   dt             = 0.1;
80c4762a1bSJed Brown   ftime_original = data.tfinal = 1.0;
81c4762a1bSJed Brown 
829566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
839566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps));
849566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime_original));
859566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
869566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, global));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* set user provided RHSFunction and RHSJacobian */
899566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &data));
909566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
919566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, mn, mn));
929566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
939566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 5, NULL));
949566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(J, 5, NULL, 5, NULL));
95c4762a1bSJed Brown 
969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-ts_fd", &flg));
97c4762a1bSJed Brown   if (!flg) {
989566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, &data));
99c4762a1bSJed Brown   } else {
1009566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
1019566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-fd_color", &fd_jacobian_coloring));
102c4762a1bSJed Brown     if (fd_jacobian_coloring) { /* Use finite differences with coloring */
103c4762a1bSJed Brown       /* Get data structure of J */
104c4762a1bSJed Brown       PetscBool pc_diagonal;
1059566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(NULL, NULL, "-pc_diagonal", &pc_diagonal));
106c4762a1bSJed Brown       if (pc_diagonal) { /* the preconditioner of J is a diagonal matrix */
107c4762a1bSJed Brown         PetscInt    rstart, rend, i;
108c4762a1bSJed Brown         PetscScalar zero = 0.0;
1099566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(J, &rstart, &rend));
110*48a46eb9SPierre Jolivet         for (i = rstart; i < rend; i++) PetscCall(MatSetValues(J, 1, &i, 1, &i, &zero, INSERT_VALUES));
1119566063dSJacob Faibussowitsch         PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1129566063dSJacob Faibussowitsch         PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
113c4762a1bSJed Brown       } else {
114c4762a1bSJed Brown         /* Fill the structure using the expensive SNESComputeJacobianDefault. Temporarily set up the TS so we can call this function */
1159566063dSJacob Faibussowitsch         PetscCall(TSSetType(ts, TSBEULER));
1169566063dSJacob Faibussowitsch         PetscCall(TSSetUp(ts));
1179566063dSJacob Faibussowitsch         PetscCall(SNESComputeJacobianDefault(snes, x, J, J, ts));
118c4762a1bSJed Brown       }
119c4762a1bSJed Brown 
120c4762a1bSJed Brown       /* create coloring context */
1219566063dSJacob Faibussowitsch       PetscCall(MatColoringCreate(J, &mc));
1229566063dSJacob Faibussowitsch       PetscCall(MatColoringSetType(mc, MATCOLORINGSL));
1239566063dSJacob Faibussowitsch       PetscCall(MatColoringSetFromOptions(mc));
1249566063dSJacob Faibussowitsch       PetscCall(MatColoringApply(mc, &iscoloring));
1259566063dSJacob Faibussowitsch       PetscCall(MatColoringDestroy(&mc));
1269566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring));
1279566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode(*)(void))SNESTSFormFunction, ts));
1289566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
1299566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(J, iscoloring, matfdcoloring));
1309566063dSJacob Faibussowitsch       PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, matfdcoloring));
1319566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&iscoloring));
132c4762a1bSJed Brown     } else { /* Use finite differences (slow) */
1339566063dSJacob Faibussowitsch       PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefault, NULL));
134c4762a1bSJed Brown     }
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* Pick up a Petsc preconditioner */
138c4762a1bSJed Brown   /* one can always set method or preconditioner during the run time */
1399566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
1409566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
1419566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
1429566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCJACOBI));
1439566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1469566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Test TSSetPostStep() */
1499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_PostStep", &flg));
1501baa6e33SBarry Smith   if (flg) PetscCall(TSSetPostStep(ts, PostStep));
151c4762a1bSJed Brown 
1529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-NOUT", &NOUT, NULL));
153c4762a1bSJed Brown   for (iout = 1; iout <= NOUT; iout++) {
1549566063dSJacob Faibussowitsch     PetscCall(TSSetMaxSteps(ts, time_steps));
1559566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts, iout * ftime_original / NOUT));
1569566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, global));
1579566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &ftime));
1589566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, ftime));
1599566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
160c4762a1bSJed Brown   }
161c4762a1bSJed Brown   /* Interpolate solution at tfinal */
1629566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &global));
1639566063dSJacob Faibussowitsch   PetscCall(TSInterpolate(ts, ftime_original, global));
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-matlab_view", &flg));
166c4762a1bSJed Brown   if (flg) { /* print solution into a MATLAB file */
1679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "out.m", &viewfile));
1689566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewfile, PETSC_VIEWER_ASCII_MATLAB));
1699566063dSJacob Faibussowitsch     PetscCall(VecView(global, viewfile));
1709566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewfile));
1719566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewfile));
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* free the memories */
1759566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
1779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1799566063dSJacob Faibussowitsch   if (fd_jacobian_coloring) PetscCall(MatFDColoringDestroy(&matfdcoloring));
1809566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
181b122ec5aSJacob Faibussowitsch   return 0;
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown /* -------------------------------------------------------------------*/
185c4762a1bSJed Brown /* the initial function */
1869371c9d4SSatish Balay PetscReal f_ini(PetscReal x, PetscReal y) {
187c4762a1bSJed Brown   PetscReal f;
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   f = PetscExpReal(-20.0 * (PetscPowRealInt(x - 0.5, 2) + PetscPowRealInt(y - 0.5, 2)));
190c4762a1bSJed Brown   return f;
191c4762a1bSJed Brown }
192c4762a1bSJed Brown 
1939371c9d4SSatish Balay PetscErrorCode Initial(Vec global, void *ctx) {
194c4762a1bSJed Brown   Data        *data = (Data *)ctx;
195c4762a1bSJed Brown   PetscInt     m, row, col;
196c4762a1bSJed Brown   PetscReal    x, y, dx, dy;
197c4762a1bSJed Brown   PetscScalar *localptr;
198c4762a1bSJed Brown   PetscInt     i, mybase, myend, locsize;
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   PetscFunctionBeginUser;
201c4762a1bSJed Brown   /* make the local  copies of parameters */
202c4762a1bSJed Brown   m  = data->m;
203c4762a1bSJed Brown   dx = data->dx;
204c4762a1bSJed Brown   dy = data->dy;
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* determine starting point of each processor */
2079566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(global, &mybase, &myend));
2089566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global, &locsize));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /* Initialize the array */
2119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(global, &localptr));
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   for (i = 0; i < locsize; i++) {
214c4762a1bSJed Brown     row         = 1 + (mybase + i) - ((mybase + i) / m) * m;
215c4762a1bSJed Brown     col         = (mybase + i) / m + 1;
216c4762a1bSJed Brown     x           = dx * row;
217c4762a1bSJed Brown     y           = dy * col;
218c4762a1bSJed Brown     localptr[i] = f_ini(x, y);
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown 
2219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(global, &localptr));
222c4762a1bSJed Brown   PetscFunctionReturn(0);
223c4762a1bSJed Brown }
224c4762a1bSJed Brown 
2259371c9d4SSatish Balay PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec global, void *ctx) {
226c4762a1bSJed Brown   VecScatter         scatter;
227c4762a1bSJed Brown   IS                 from, to;
228c4762a1bSJed Brown   PetscInt           i, n, *idx, nsteps, maxsteps;
229c4762a1bSJed Brown   Vec                tmp_vec;
230c4762a1bSJed Brown   const PetscScalar *tmp;
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   PetscFunctionBeginUser;
2339566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &nsteps));
234c4762a1bSJed Brown   /* display output at selected time steps */
2359566063dSJacob Faibussowitsch   PetscCall(TSGetMaxSteps(ts, &maxsteps));
236c4762a1bSJed Brown   if (nsteps % 10 != 0) PetscFunctionReturn(0);
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   /* Get the size of the vector */
2399566063dSJacob Faibussowitsch   PetscCall(VecGetSize(global, &n));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   /* Set the index sets */
2429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idx));
243c4762a1bSJed Brown   for (i = 0; i < n; i++) idx[i] = i;
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /* Create local sequential vectors */
2469566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &tmp_vec));
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   /* Create scatter context */
2499566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &from));
2509566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &to));
2519566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global, from, tmp_vec, to, &scatter));
2529566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
2539566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
254c4762a1bSJed Brown 
2559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_vec, &tmp));
25663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "At t[%" PetscInt_FMT "] =%14.2e u= %14.2e at the center \n", nsteps, (double)time, (double)PetscRealPart(tmp[n / 2])));
2579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_vec, &tmp));
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
2609566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
2619566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
2629566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
2639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_vec));
264c4762a1bSJed Brown   PetscFunctionReturn(0);
265c4762a1bSJed Brown }
266c4762a1bSJed Brown 
2679371c9d4SSatish Balay PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat BB, void *ptr) {
268c4762a1bSJed Brown   Data       *data = (Data *)ptr;
269c4762a1bSJed Brown   PetscScalar v[5];
270c4762a1bSJed Brown   PetscInt    idx[5], i, j, row;
271c4762a1bSJed Brown   PetscInt    m, n, mn;
272c4762a1bSJed Brown   PetscReal   dx, dy, a, epsilon, xc, xl, xr, yl, yr;
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   PetscFunctionBeginUser;
275c4762a1bSJed Brown   m       = data->m;
276c4762a1bSJed Brown   n       = data->n;
277c4762a1bSJed Brown   mn      = m * n;
278c4762a1bSJed Brown   dx      = data->dx;
279c4762a1bSJed Brown   dy      = data->dy;
280c4762a1bSJed Brown   a       = data->a;
281c4762a1bSJed Brown   epsilon = data->epsilon;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   xc = -2.0 * epsilon * (1.0 / (dx * dx) + 1.0 / (dy * dy));
284c4762a1bSJed Brown   xl = 0.5 * a / dx + epsilon / (dx * dx);
285c4762a1bSJed Brown   xr = -0.5 * a / dx + epsilon / (dx * dx);
286c4762a1bSJed Brown   yl = 0.5 * a / dy + epsilon / (dy * dy);
287c4762a1bSJed Brown   yr = -0.5 * a / dy + epsilon / (dy * dy);
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   row    = 0;
2909371c9d4SSatish Balay   v[0]   = xc;
2919371c9d4SSatish Balay   v[1]   = xr;
2929371c9d4SSatish Balay   v[2]   = yr;
2939371c9d4SSatish Balay   idx[0] = 0;
2949371c9d4SSatish Balay   idx[1] = 2;
2959371c9d4SSatish Balay   idx[2] = m;
2969566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 3, idx, v, INSERT_VALUES));
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   row    = m - 1;
2999371c9d4SSatish Balay   v[0]   = 2.0 * xl;
3009371c9d4SSatish Balay   v[1]   = xc;
3019371c9d4SSatish Balay   v[2]   = yr;
3029371c9d4SSatish Balay   idx[0] = m - 2;
3039371c9d4SSatish Balay   idx[1] = m - 1;
3049371c9d4SSatish Balay   idx[2] = m - 1 + m;
3059566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 3, idx, v, INSERT_VALUES));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   for (i = 1; i < m - 1; i++) {
308c4762a1bSJed Brown     row    = i;
3099371c9d4SSatish Balay     v[0]   = xl;
3109371c9d4SSatish Balay     v[1]   = xc;
3119371c9d4SSatish Balay     v[2]   = xr;
3129371c9d4SSatish Balay     v[3]   = yr;
3139371c9d4SSatish Balay     idx[0] = i - 1;
3149371c9d4SSatish Balay     idx[1] = i;
3159371c9d4SSatish Balay     idx[2] = i + 1;
3169371c9d4SSatish Balay     idx[3] = i + m;
3179566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
318c4762a1bSJed Brown   }
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   for (j = 1; j < n - 1; j++) {
321c4762a1bSJed Brown     row    = j * m;
3229371c9d4SSatish Balay     v[0]   = xc;
3239371c9d4SSatish Balay     v[1]   = xr;
3249371c9d4SSatish Balay     v[2]   = yl;
3259371c9d4SSatish Balay     v[3]   = yr;
3269371c9d4SSatish Balay     idx[0] = j * m;
3279371c9d4SSatish Balay     idx[1] = j * m;
3289371c9d4SSatish Balay     idx[2] = j * m - m;
3299371c9d4SSatish Balay     idx[3] = j * m + m;
3309566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
331c4762a1bSJed Brown 
332c4762a1bSJed Brown     row    = j * m + m - 1;
3339371c9d4SSatish Balay     v[0]   = xc;
3349371c9d4SSatish Balay     v[1]   = 2.0 * xl;
3359371c9d4SSatish Balay     v[2]   = yl;
3369371c9d4SSatish Balay     v[3]   = yr;
3379371c9d4SSatish Balay     idx[0] = j * m + m - 1;
3389371c9d4SSatish Balay     idx[1] = j * m + m - 1 - 1;
3399371c9d4SSatish Balay     idx[2] = j * m + m - 1 - m;
3409371c9d4SSatish Balay     idx[3] = j * m + m - 1 + m;
3419566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
342c4762a1bSJed Brown 
343c4762a1bSJed Brown     for (i = 1; i < m - 1; i++) {
344c4762a1bSJed Brown       row    = j * m + i;
3459371c9d4SSatish Balay       v[0]   = xc;
3469371c9d4SSatish Balay       v[1]   = xl;
3479371c9d4SSatish Balay       v[2]   = xr;
3489371c9d4SSatish Balay       v[3]   = yl;
3499371c9d4SSatish Balay       v[4]   = yr;
3509371c9d4SSatish Balay       idx[0] = j * m + i;
3519371c9d4SSatish Balay       idx[1] = j * m + i - 1;
3529371c9d4SSatish Balay       idx[2] = j * m + i + 1;
3539371c9d4SSatish Balay       idx[3] = j * m + i - m;
354c4762a1bSJed Brown       idx[4] = j * m + i + m;
3559566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &row, 5, idx, v, INSERT_VALUES));
356c4762a1bSJed Brown     }
357c4762a1bSJed Brown   }
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   row    = mn - m;
3609371c9d4SSatish Balay   v[0]   = xc;
3619371c9d4SSatish Balay   v[1]   = xr;
3629371c9d4SSatish Balay   v[2]   = 2.0 * yl;
3639371c9d4SSatish Balay   idx[0] = mn - m;
3649371c9d4SSatish Balay   idx[1] = mn - m + 1;
3659371c9d4SSatish Balay   idx[2] = mn - m - m;
3669566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 3, idx, v, INSERT_VALUES));
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   row    = mn - 1;
3699371c9d4SSatish Balay   v[0]   = xc;
3709371c9d4SSatish Balay   v[1]   = 2.0 * xl;
3719371c9d4SSatish Balay   v[2]   = 2.0 * yl;
3729371c9d4SSatish Balay   idx[0] = mn - 1;
3739371c9d4SSatish Balay   idx[1] = mn - 2;
3749371c9d4SSatish Balay   idx[2] = mn - 1 - m;
3759566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   for (i = 1; i < m - 1; i++) {
378c4762a1bSJed Brown     row    = mn - m + i;
3799371c9d4SSatish Balay     v[0]   = xl;
3809371c9d4SSatish Balay     v[1]   = xc;
3819371c9d4SSatish Balay     v[2]   = xr;
3829371c9d4SSatish Balay     v[3]   = 2.0 * yl;
3839371c9d4SSatish Balay     idx[0] = mn - m + i - 1;
3849371c9d4SSatish Balay     idx[1] = mn - m + i;
3859371c9d4SSatish Balay     idx[2] = mn - m + i + 1;
3869371c9d4SSatish Balay     idx[3] = mn - m + i - m;
3879566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
388c4762a1bSJed Brown   }
389c4762a1bSJed Brown 
3909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   PetscFunctionReturn(0);
394c4762a1bSJed Brown }
395c4762a1bSJed Brown 
396c4762a1bSJed Brown /* globalout = -a*(u_x+u_y) + epsilon*(u_xx+u_yy) */
3979371c9d4SSatish Balay PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx) {
398c4762a1bSJed Brown   Data              *data = (Data *)ctx;
399c4762a1bSJed Brown   PetscInt           m, n, mn;
400c4762a1bSJed Brown   PetscReal          dx, dy;
401c4762a1bSJed Brown   PetscReal          xc, xl, xr, yl, yr;
402c4762a1bSJed Brown   PetscReal          a, epsilon;
403c4762a1bSJed Brown   PetscScalar       *outptr;
404c4762a1bSJed Brown   const PetscScalar *inptr;
405c4762a1bSJed Brown   PetscInt           i, j, len;
406c4762a1bSJed Brown   IS                 from, to;
407c4762a1bSJed Brown   PetscInt          *idx;
408c4762a1bSJed Brown   VecScatter         scatter;
409c4762a1bSJed Brown   Vec                tmp_in, tmp_out;
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   PetscFunctionBeginUser;
412c4762a1bSJed Brown   m       = data->m;
413c4762a1bSJed Brown   n       = data->n;
414c4762a1bSJed Brown   mn      = m * n;
415c4762a1bSJed Brown   dx      = data->dx;
416c4762a1bSJed Brown   dy      = data->dy;
417c4762a1bSJed Brown   a       = data->a;
418c4762a1bSJed Brown   epsilon = data->epsilon;
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   xc = -2.0 * epsilon * (1.0 / (dx * dx) + 1.0 / (dy * dy));
421c4762a1bSJed Brown   xl = 0.5 * a / dx + epsilon / (dx * dx);
422c4762a1bSJed Brown   xr = -0.5 * a / dx + epsilon / (dx * dx);
423c4762a1bSJed Brown   yl = 0.5 * a / dy + epsilon / (dy * dy);
424c4762a1bSJed Brown   yr = -0.5 * a / dy + epsilon / (dy * dy);
425c4762a1bSJed Brown 
426c4762a1bSJed Brown   /* Get the length of parallel vector */
4279566063dSJacob Faibussowitsch   PetscCall(VecGetSize(globalin, &len));
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   /* Set the index sets */
4309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len, &idx));
431c4762a1bSJed Brown   for (i = 0; i < len; i++) idx[i] = i;
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   /* Create local sequential vectors */
4349566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, len, &tmp_in));
4359566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tmp_in, &tmp_out));
436c4762a1bSJed Brown 
437c4762a1bSJed Brown   /* Create scatter context */
4389566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, len, idx, PETSC_COPY_VALUES, &from));
4399566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, len, idx, PETSC_COPY_VALUES, &to));
4409566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(globalin, from, tmp_in, to, &scatter));
4419566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
4429566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
4439566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   /*Extract income array - include ghost points */
4469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_in, &inptr));
447c4762a1bSJed Brown 
448c4762a1bSJed Brown   /* Extract outcome array*/
4499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(tmp_out, &outptr));
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   outptr[0]     = xc * inptr[0] + xr * inptr[1] + yr * inptr[m];
452c4762a1bSJed Brown   outptr[m - 1] = 2.0 * xl * inptr[m - 2] + xc * inptr[m - 1] + yr * inptr[m - 1 + m];
4539371c9d4SSatish Balay   for (i = 1; i < m - 1; i++) { outptr[i] = xc * inptr[i] + xl * inptr[i - 1] + xr * inptr[i + 1] + yr * inptr[i + m]; }
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   for (j = 1; j < n - 1; j++) {
456c4762a1bSJed Brown     outptr[j * m]         = xc * inptr[j * m] + xr * inptr[j * m + 1] + yl * inptr[j * m - m] + yr * inptr[j * m + m];
457c4762a1bSJed Brown     outptr[j * m + m - 1] = xc * inptr[j * m + m - 1] + 2.0 * xl * inptr[j * m + m - 1 - 1] + yl * inptr[j * m + m - 1 - m] + yr * inptr[j * m + m - 1 + m];
4589371c9d4SSatish Balay     for (i = 1; i < m - 1; i++) { outptr[j * m + i] = xc * inptr[j * m + i] + xl * inptr[j * m + i - 1] + xr * inptr[j * m + i + 1] + yl * inptr[j * m + i - m] + yr * inptr[j * m + i + m]; }
459c4762a1bSJed Brown   }
460c4762a1bSJed Brown 
461c4762a1bSJed Brown   outptr[mn - m] = xc * inptr[mn - m] + xr * inptr[mn - m + 1] + 2.0 * yl * inptr[mn - m - m];
462c4762a1bSJed Brown   outptr[mn - 1] = 2.0 * xl * inptr[mn - 2] + xc * inptr[mn - 1] + 2.0 * yl * inptr[mn - 1 - m];
4639371c9d4SSatish Balay   for (i = 1; i < m - 1; i++) { outptr[mn - m + i] = xc * inptr[mn - m + i] + xl * inptr[mn - m + i - 1] + xr * inptr[mn - m + i + 1] + 2 * yl * inptr[mn - m + i - m]; }
464c4762a1bSJed Brown 
4659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_in, &inptr));
4669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(tmp_out, &outptr));
467c4762a1bSJed Brown 
4689566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(tmp_out, from, globalout, to, &scatter));
4699566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
4709566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
471c4762a1bSJed Brown 
472c4762a1bSJed Brown   /* Destroy idx aand scatter */
4739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_in));
4749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_out));
4759566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
4769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
4779566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
478c4762a1bSJed Brown 
4799566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
480c4762a1bSJed Brown   PetscFunctionReturn(0);
481c4762a1bSJed Brown }
482c4762a1bSJed Brown 
4839371c9d4SSatish Balay PetscErrorCode PostStep(TS ts) {
484c4762a1bSJed Brown   PetscReal t;
485c4762a1bSJed Brown 
486c4762a1bSJed Brown   PetscFunctionBeginUser;
4879566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
4889566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  PostStep, t: %g\n", (double)t));
489c4762a1bSJed Brown   PetscFunctionReturn(0);
490c4762a1bSJed Brown }
491c4762a1bSJed Brown 
492c4762a1bSJed Brown /*TEST
493c4762a1bSJed Brown 
494c4762a1bSJed Brown     test:
495c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
496c4762a1bSJed Brown       output_file: output/ex4.out
497c4762a1bSJed Brown 
498c4762a1bSJed Brown     test:
499c4762a1bSJed Brown       suffix: 2
500c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
501c4762a1bSJed Brown       nsize: 2
502c4762a1bSJed Brown       output_file: output/ex4.out
503c4762a1bSJed Brown 
504c4762a1bSJed Brown     test:
505c4762a1bSJed Brown       suffix: 3
506c4762a1bSJed Brown       args: -ts_fd -ts_type cn
507c4762a1bSJed Brown 
508c4762a1bSJed Brown     test:
509c4762a1bSJed Brown       suffix: 4
510c4762a1bSJed Brown       args: -ts_fd -ts_type cn
511c4762a1bSJed Brown       output_file: output/ex4_3.out
512c4762a1bSJed Brown       nsize: 2
513c4762a1bSJed Brown 
514c4762a1bSJed Brown     test:
515c4762a1bSJed Brown       suffix: 5
516c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
517c4762a1bSJed Brown       output_file: output/ex4.out
518c4762a1bSJed Brown 
519c4762a1bSJed Brown     test:
520c4762a1bSJed Brown       suffix: 6
521c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
522c4762a1bSJed Brown       output_file: output/ex4.out
523c4762a1bSJed Brown       nsize: 2
524c4762a1bSJed Brown 
525c4762a1bSJed Brown     test:
526c4762a1bSJed Brown       suffix: 7
527c4762a1bSJed Brown       requires: !single
528c4762a1bSJed Brown       args: -ts_fd -ts_type beuler -test_PostStep -ts_dt .1
529c4762a1bSJed Brown 
530c4762a1bSJed Brown     test:
531c4762a1bSJed Brown       suffix: 8
532c4762a1bSJed Brown       requires: !single
533c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 5dp -ts_dt .01 -ts_adapt_type none -ts_view
534c4762a1bSJed Brown 
535c4762a1bSJed Brown TEST*/
536