xref: /petsc/src/ts/tests/ex4.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 
main(int argc,char ** argv)36d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
37d71ae5a4SJacob Faibussowitsch {
38c4762a1bSJed Brown   PetscInt      time_steps = 100, iout, NOUT = 1;
39c4762a1bSJed Brown   Vec           global;
40c4762a1bSJed Brown   PetscReal     dt, ftime, ftime_original;
41c4762a1bSJed Brown   TS            ts;
42c4762a1bSJed Brown   PetscViewer   viewfile;
43c4762a1bSJed Brown   Mat           J = 0;
44c4762a1bSJed Brown   Vec           x;
45c4762a1bSJed Brown   Data          data;
46c4762a1bSJed Brown   PetscInt      mn;
47c4762a1bSJed Brown   PetscBool     flg;
48c4762a1bSJed Brown   MatColoring   mc;
49c4762a1bSJed Brown   ISColoring    iscoloring;
50c4762a1bSJed Brown   MatFDColoring matfdcoloring        = 0;
51c4762a1bSJed Brown   PetscBool     fd_jacobian_coloring = PETSC_FALSE;
52c4762a1bSJed Brown   SNES          snes;
53c4762a1bSJed Brown   KSP           ksp;
54c4762a1bSJed Brown   PC            pc;
55c4762a1bSJed Brown 
56327415f7SBarry Smith   PetscFunctionBeginUser;
57c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* set data */
60c4762a1bSJed Brown   data.m       = 9;
61c4762a1bSJed Brown   data.n       = 9;
62c4762a1bSJed Brown   data.a       = 1.0;
63c4762a1bSJed Brown   data.epsilon = 0.1;
64c4762a1bSJed Brown   data.dx      = 1.0 / (data.m + 1.0);
65c4762a1bSJed Brown   data.dy      = 1.0 / (data.n + 1.0);
66c4762a1bSJed Brown   mn           = (data.m) * (data.n);
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-time", &time_steps, NULL));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* set initial conditions */
709566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &global));
719566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(global, PETSC_DECIDE, mn));
729566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(global));
739566063dSJacob Faibussowitsch   PetscCall(Initial(global, &data));
749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(global, &x));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* create timestep context */
779566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
789566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, Monitor, &data, NULL));
799566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSEULER));
80c4762a1bSJed Brown   dt             = 0.1;
81c4762a1bSJed Brown   ftime_original = data.tfinal = 1.0;
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
849566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, time_steps));
859566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime_original));
869566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
879566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, global));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* set user provided RHSFunction and RHSJacobian */
909566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &data));
919566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
929566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, mn, mn));
939566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
949566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 5, NULL));
959566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(J, 5, NULL, 5, NULL));
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-ts_fd", &flg));
98c4762a1bSJed Brown   if (!flg) {
999566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, &data));
100c4762a1bSJed Brown   } else {
1019566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
1029566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-fd_color", &fd_jacobian_coloring));
103c4762a1bSJed Brown     if (fd_jacobian_coloring) { /* Use finite differences with coloring */
104c4762a1bSJed Brown       /* Get data structure of J */
105c4762a1bSJed Brown       PetscBool pc_diagonal;
1069566063dSJacob Faibussowitsch       PetscCall(PetscOptionsHasName(NULL, NULL, "-pc_diagonal", &pc_diagonal));
107c4762a1bSJed Brown       if (pc_diagonal) { /* the preconditioner of J is a diagonal matrix */
108c4762a1bSJed Brown         PetscInt    rstart, rend, i;
109c4762a1bSJed Brown         PetscScalar zero = 0.0;
1109566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(J, &rstart, &rend));
11148a46eb9SPierre Jolivet         for (i = rstart; i < rend; i++) PetscCall(MatSetValues(J, 1, &i, 1, &i, &zero, INSERT_VALUES));
1129566063dSJacob Faibussowitsch         PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1139566063dSJacob Faibussowitsch         PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
114c4762a1bSJed Brown       } else {
115c4762a1bSJed Brown         /* Fill the structure using the expensive SNESComputeJacobianDefault. Temporarily set up the TS so we can call this function */
1169566063dSJacob Faibussowitsch         PetscCall(TSSetType(ts, TSBEULER));
1179566063dSJacob Faibussowitsch         PetscCall(TSSetUp(ts));
1189566063dSJacob Faibussowitsch         PetscCall(SNESComputeJacobianDefault(snes, x, J, J, ts));
119c4762a1bSJed Brown       }
120c4762a1bSJed Brown 
121c4762a1bSJed Brown       /* create coloring context */
1229566063dSJacob Faibussowitsch       PetscCall(MatColoringCreate(J, &mc));
1239566063dSJacob Faibussowitsch       PetscCall(MatColoringSetType(mc, MATCOLORINGSL));
1249566063dSJacob Faibussowitsch       PetscCall(MatColoringSetFromOptions(mc));
1259566063dSJacob Faibussowitsch       PetscCall(MatColoringApply(mc, &iscoloring));
1269566063dSJacob Faibussowitsch       PetscCall(MatColoringDestroy(&mc));
1279566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring));
1282ba42892SBarry Smith       PetscCall(MatFDColoringSetFunction(matfdcoloring, (MatFDColoringFn *)SNESTSFormFunction, ts));
1299566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
1309566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(J, iscoloring, matfdcoloring));
1319566063dSJacob Faibussowitsch       PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, matfdcoloring));
1329566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&iscoloring));
133c4762a1bSJed Brown     } else { /* Use finite differences (slow) */
1349566063dSJacob Faibussowitsch       PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefault, NULL));
135c4762a1bSJed Brown     }
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown 
138f0b74427SPierre Jolivet   /* Pick up a PETSc preconditioner */
139c4762a1bSJed Brown   /* one can always set method or preconditioner during the run time */
1409566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
1419566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
1429566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
1439566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCJACOBI));
1449566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
145c4762a1bSJed Brown 
1469566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1479566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* Test TSSetPostStep() */
1509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_PostStep", &flg));
1511baa6e33SBarry Smith   if (flg) PetscCall(TSSetPostStep(ts, PostStep));
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-NOUT", &NOUT, NULL));
154c4762a1bSJed Brown   for (iout = 1; iout <= NOUT; iout++) {
1559566063dSJacob Faibussowitsch     PetscCall(TSSetMaxSteps(ts, time_steps));
1569566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts, iout * ftime_original / NOUT));
1579566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, global));
1589566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &ftime));
1599566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, ftime));
1609566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
161c4762a1bSJed Brown   }
162c4762a1bSJed Brown   /* Interpolate solution at tfinal */
1639566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &global));
1649566063dSJacob Faibussowitsch   PetscCall(TSInterpolate(ts, ftime_original, global));
165c4762a1bSJed Brown 
1669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-matlab_view", &flg));
167c4762a1bSJed Brown   if (flg) { /* print solution into a MATLAB file */
1689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "out.m", &viewfile));
1699566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewfile, PETSC_VIEWER_ASCII_MATLAB));
1709566063dSJacob Faibussowitsch     PetscCall(VecView(global, viewfile));
1719566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewfile));
1729566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewfile));
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* free the memories */
1769566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
1789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1809566063dSJacob Faibussowitsch   if (fd_jacobian_coloring) PetscCall(MatFDColoringDestroy(&matfdcoloring));
1819566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
182b122ec5aSJacob Faibussowitsch   return 0;
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown /* -------------------------------------------------------------------*/
186c4762a1bSJed Brown /* the initial function */
f_ini(PetscReal x,PetscReal y)187d71ae5a4SJacob Faibussowitsch PetscReal f_ini(PetscReal x, PetscReal y)
188d71ae5a4SJacob Faibussowitsch {
189c4762a1bSJed Brown   PetscReal f;
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   f = PetscExpReal(-20.0 * (PetscPowRealInt(x - 0.5, 2) + PetscPowRealInt(y - 0.5, 2)));
192c4762a1bSJed Brown   return f;
193c4762a1bSJed Brown }
194c4762a1bSJed Brown 
Initial(Vec global,PetscCtx ctx)195*2a8381b2SBarry Smith PetscErrorCode Initial(Vec global, PetscCtx ctx)
196d71ae5a4SJacob Faibussowitsch {
197c4762a1bSJed Brown   Data        *data = (Data *)ctx;
198c4762a1bSJed Brown   PetscInt     m, row, col;
199c4762a1bSJed Brown   PetscReal    x, y, dx, dy;
200c4762a1bSJed Brown   PetscScalar *localptr;
201c4762a1bSJed Brown   PetscInt     i, mybase, myend, locsize;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   PetscFunctionBeginUser;
204c4762a1bSJed Brown   /* make the local  copies of parameters */
205c4762a1bSJed Brown   m  = data->m;
206c4762a1bSJed Brown   dx = data->dx;
207c4762a1bSJed Brown   dy = data->dy;
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* determine starting point of each processor */
2109566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(global, &mybase, &myend));
2119566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global, &locsize));
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* Initialize the array */
2149566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(global, &localptr));
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   for (i = 0; i < locsize; i++) {
217c4762a1bSJed Brown     row         = 1 + (mybase + i) - ((mybase + i) / m) * m;
218c4762a1bSJed Brown     col         = (mybase + i) / m + 1;
219c4762a1bSJed Brown     x           = dx * row;
220c4762a1bSJed Brown     y           = dy * col;
221c4762a1bSJed Brown     localptr[i] = f_ini(x, y);
222c4762a1bSJed Brown   }
223c4762a1bSJed Brown 
2249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(global, &localptr));
2253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
226c4762a1bSJed Brown }
227c4762a1bSJed Brown 
Monitor(TS ts,PetscInt step,PetscReal time,Vec global,PetscCtx ctx)228*2a8381b2SBarry Smith PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec global, PetscCtx ctx)
229d71ae5a4SJacob Faibussowitsch {
230c4762a1bSJed Brown   VecScatter         scatter;
231c4762a1bSJed Brown   IS                 from, to;
232c4762a1bSJed Brown   PetscInt           i, n, *idx, nsteps, maxsteps;
233c4762a1bSJed Brown   Vec                tmp_vec;
234c4762a1bSJed Brown   const PetscScalar *tmp;
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   PetscFunctionBeginUser;
2379566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &nsteps));
238c4762a1bSJed Brown   /* display output at selected time steps */
2399566063dSJacob Faibussowitsch   PetscCall(TSGetMaxSteps(ts, &maxsteps));
2403ba16761SJacob Faibussowitsch   if (nsteps % 10 != 0) PetscFunctionReturn(PETSC_SUCCESS);
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* Get the size of the vector */
2439566063dSJacob Faibussowitsch   PetscCall(VecGetSize(global, &n));
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /* Set the index sets */
2469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idx));
247c4762a1bSJed Brown   for (i = 0; i < n; i++) idx[i] = i;
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* Create local sequential vectors */
2509566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &tmp_vec));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /* Create scatter context */
2539566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &from));
2549566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &to));
2559566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global, from, tmp_vec, to, &scatter));
2569566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
2579566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_vec, &tmp));
26063a3b9bcSJacob 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])));
2619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_vec, &tmp));
262c4762a1bSJed Brown 
2639566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
2649566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
2659566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
2669566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
2679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_vec));
2683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
269c4762a1bSJed Brown }
270c4762a1bSJed Brown 
RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void * ptr)271d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat BB, void *ptr)
272d71ae5a4SJacob Faibussowitsch {
273c4762a1bSJed Brown   Data       *data = (Data *)ptr;
274c4762a1bSJed Brown   PetscScalar v[5];
275c4762a1bSJed Brown   PetscInt    idx[5], i, j, row;
276c4762a1bSJed Brown   PetscInt    m, n, mn;
277c4762a1bSJed Brown   PetscReal   dx, dy, a, epsilon, xc, xl, xr, yl, yr;
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   PetscFunctionBeginUser;
280c4762a1bSJed Brown   m       = data->m;
281c4762a1bSJed Brown   n       = data->n;
282c4762a1bSJed Brown   mn      = m * n;
283c4762a1bSJed Brown   dx      = data->dx;
284c4762a1bSJed Brown   dy      = data->dy;
285c4762a1bSJed Brown   a       = data->a;
286c4762a1bSJed Brown   epsilon = data->epsilon;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   xc = -2.0 * epsilon * (1.0 / (dx * dx) + 1.0 / (dy * dy));
289c4762a1bSJed Brown   xl = 0.5 * a / dx + epsilon / (dx * dx);
290c4762a1bSJed Brown   xr = -0.5 * a / dx + epsilon / (dx * dx);
291c4762a1bSJed Brown   yl = 0.5 * a / dy + epsilon / (dy * dy);
292c4762a1bSJed Brown   yr = -0.5 * a / dy + epsilon / (dy * dy);
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   row    = 0;
2959371c9d4SSatish Balay   v[0]   = xc;
2969371c9d4SSatish Balay   v[1]   = xr;
2979371c9d4SSatish Balay   v[2]   = yr;
2989371c9d4SSatish Balay   idx[0] = 0;
2999371c9d4SSatish Balay   idx[1] = 2;
3009371c9d4SSatish Balay   idx[2] = m;
3019566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 3, idx, v, INSERT_VALUES));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   row    = m - 1;
3049371c9d4SSatish Balay   v[0]   = 2.0 * xl;
3059371c9d4SSatish Balay   v[1]   = xc;
3069371c9d4SSatish Balay   v[2]   = yr;
3079371c9d4SSatish Balay   idx[0] = m - 2;
3089371c9d4SSatish Balay   idx[1] = m - 1;
3099371c9d4SSatish Balay   idx[2] = m - 1 + m;
3109566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 3, idx, v, INSERT_VALUES));
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   for (i = 1; i < m - 1; i++) {
313c4762a1bSJed Brown     row    = i;
3149371c9d4SSatish Balay     v[0]   = xl;
3159371c9d4SSatish Balay     v[1]   = xc;
3169371c9d4SSatish Balay     v[2]   = xr;
3179371c9d4SSatish Balay     v[3]   = yr;
3189371c9d4SSatish Balay     idx[0] = i - 1;
3199371c9d4SSatish Balay     idx[1] = i;
3209371c9d4SSatish Balay     idx[2] = i + 1;
3219371c9d4SSatish Balay     idx[3] = i + m;
3229566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
323c4762a1bSJed Brown   }
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   for (j = 1; j < n - 1; j++) {
326c4762a1bSJed Brown     row    = j * m;
3279371c9d4SSatish Balay     v[0]   = xc;
3289371c9d4SSatish Balay     v[1]   = xr;
3299371c9d4SSatish Balay     v[2]   = yl;
3309371c9d4SSatish Balay     v[3]   = yr;
3319371c9d4SSatish Balay     idx[0] = j * m;
3329371c9d4SSatish Balay     idx[1] = j * m;
3339371c9d4SSatish Balay     idx[2] = j * m - m;
3349371c9d4SSatish Balay     idx[3] = j * m + m;
3359566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
336c4762a1bSJed Brown 
337c4762a1bSJed Brown     row    = j * m + m - 1;
3389371c9d4SSatish Balay     v[0]   = xc;
3399371c9d4SSatish Balay     v[1]   = 2.0 * xl;
3409371c9d4SSatish Balay     v[2]   = yl;
3419371c9d4SSatish Balay     v[3]   = yr;
3429371c9d4SSatish Balay     idx[0] = j * m + m - 1;
3439371c9d4SSatish Balay     idx[1] = j * m + m - 1 - 1;
3449371c9d4SSatish Balay     idx[2] = j * m + m - 1 - m;
3459371c9d4SSatish Balay     idx[3] = j * m + m - 1 + m;
3469566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
347c4762a1bSJed Brown 
348c4762a1bSJed Brown     for (i = 1; i < m - 1; i++) {
349c4762a1bSJed Brown       row    = j * m + i;
3509371c9d4SSatish Balay       v[0]   = xc;
3519371c9d4SSatish Balay       v[1]   = xl;
3529371c9d4SSatish Balay       v[2]   = xr;
3539371c9d4SSatish Balay       v[3]   = yl;
3549371c9d4SSatish Balay       v[4]   = yr;
3559371c9d4SSatish Balay       idx[0] = j * m + i;
3569371c9d4SSatish Balay       idx[1] = j * m + i - 1;
3579371c9d4SSatish Balay       idx[2] = j * m + i + 1;
3589371c9d4SSatish Balay       idx[3] = j * m + i - m;
359c4762a1bSJed Brown       idx[4] = j * m + i + m;
3609566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &row, 5, idx, v, INSERT_VALUES));
361c4762a1bSJed Brown     }
362c4762a1bSJed Brown   }
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   row    = mn - m;
3659371c9d4SSatish Balay   v[0]   = xc;
3669371c9d4SSatish Balay   v[1]   = xr;
3679371c9d4SSatish Balay   v[2]   = 2.0 * yl;
3689371c9d4SSatish Balay   idx[0] = mn - m;
3699371c9d4SSatish Balay   idx[1] = mn - m + 1;
3709371c9d4SSatish Balay   idx[2] = mn - m - m;
3719566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 3, idx, v, INSERT_VALUES));
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   row    = mn - 1;
3749371c9d4SSatish Balay   v[0]   = xc;
3759371c9d4SSatish Balay   v[1]   = 2.0 * xl;
3769371c9d4SSatish Balay   v[2]   = 2.0 * yl;
3779371c9d4SSatish Balay   idx[0] = mn - 1;
3789371c9d4SSatish Balay   idx[1] = mn - 2;
3799371c9d4SSatish Balay   idx[2] = mn - 1 - m;
3809566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   for (i = 1; i < m - 1; i++) {
383c4762a1bSJed Brown     row    = mn - m + i;
3849371c9d4SSatish Balay     v[0]   = xl;
3859371c9d4SSatish Balay     v[1]   = xc;
3869371c9d4SSatish Balay     v[2]   = xr;
3879371c9d4SSatish Balay     v[3]   = 2.0 * yl;
3889371c9d4SSatish Balay     idx[0] = mn - m + i - 1;
3899371c9d4SSatish Balay     idx[1] = mn - m + i;
3909371c9d4SSatish Balay     idx[2] = mn - m + i + 1;
3919371c9d4SSatish Balay     idx[3] = mn - m + i - m;
3929566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
393c4762a1bSJed Brown   }
394c4762a1bSJed Brown 
3959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
398c4762a1bSJed Brown }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown /* globalout = -a*(u_x+u_y) + epsilon*(u_xx+u_yy) */
RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,PetscCtx ctx)401*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, PetscCtx ctx)
402d71ae5a4SJacob Faibussowitsch {
403c4762a1bSJed Brown   Data              *data = (Data *)ctx;
404c4762a1bSJed Brown   PetscInt           m, n, mn;
405c4762a1bSJed Brown   PetscReal          dx, dy;
406c4762a1bSJed Brown   PetscReal          xc, xl, xr, yl, yr;
407c4762a1bSJed Brown   PetscReal          a, epsilon;
408c4762a1bSJed Brown   PetscScalar       *outptr;
409c4762a1bSJed Brown   const PetscScalar *inptr;
410c4762a1bSJed Brown   PetscInt           i, j, len;
411c4762a1bSJed Brown   IS                 from, to;
412c4762a1bSJed Brown   PetscInt          *idx;
413c4762a1bSJed Brown   VecScatter         scatter;
414c4762a1bSJed Brown   Vec                tmp_in, tmp_out;
415c4762a1bSJed Brown 
416c4762a1bSJed Brown   PetscFunctionBeginUser;
417c4762a1bSJed Brown   m       = data->m;
418c4762a1bSJed Brown   n       = data->n;
419c4762a1bSJed Brown   mn      = m * n;
420c4762a1bSJed Brown   dx      = data->dx;
421c4762a1bSJed Brown   dy      = data->dy;
422c4762a1bSJed Brown   a       = data->a;
423c4762a1bSJed Brown   epsilon = data->epsilon;
424c4762a1bSJed Brown 
425c4762a1bSJed Brown   xc = -2.0 * epsilon * (1.0 / (dx * dx) + 1.0 / (dy * dy));
426c4762a1bSJed Brown   xl = 0.5 * a / dx + epsilon / (dx * dx);
427c4762a1bSJed Brown   xr = -0.5 * a / dx + epsilon / (dx * dx);
428c4762a1bSJed Brown   yl = 0.5 * a / dy + epsilon / (dy * dy);
429c4762a1bSJed Brown   yr = -0.5 * a / dy + epsilon / (dy * dy);
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   /* Get the length of parallel vector */
4329566063dSJacob Faibussowitsch   PetscCall(VecGetSize(globalin, &len));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /* Set the index sets */
4359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len, &idx));
436c4762a1bSJed Brown   for (i = 0; i < len; i++) idx[i] = i;
437c4762a1bSJed Brown 
438c4762a1bSJed Brown   /* Create local sequential vectors */
4399566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, len, &tmp_in));
4409566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tmp_in, &tmp_out));
441c4762a1bSJed Brown 
442c4762a1bSJed Brown   /* Create scatter context */
4439566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, len, idx, PETSC_COPY_VALUES, &from));
4449566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, len, idx, PETSC_COPY_VALUES, &to));
4459566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(globalin, from, tmp_in, to, &scatter));
4469566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
4479566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
4489566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   /*Extract income array - include ghost points */
4519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_in, &inptr));
452c4762a1bSJed Brown 
453c4762a1bSJed Brown   /* Extract outcome array*/
4549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(tmp_out, &outptr));
455c4762a1bSJed Brown 
456c4762a1bSJed Brown   outptr[0]     = xc * inptr[0] + xr * inptr[1] + yr * inptr[m];
457c4762a1bSJed Brown   outptr[m - 1] = 2.0 * xl * inptr[m - 2] + xc * inptr[m - 1] + yr * inptr[m - 1 + m];
458ad540459SPierre Jolivet   for (i = 1; i < m - 1; i++) outptr[i] = xc * inptr[i] + xl * inptr[i - 1] + xr * inptr[i + 1] + yr * inptr[i + m];
459c4762a1bSJed Brown 
460c4762a1bSJed Brown   for (j = 1; j < n - 1; j++) {
461c4762a1bSJed Brown     outptr[j * m]         = xc * inptr[j * m] + xr * inptr[j * m + 1] + yl * inptr[j * m - m] + yr * inptr[j * m + m];
462c4762a1bSJed 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];
463ad540459SPierre Jolivet     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];
464c4762a1bSJed Brown   }
465c4762a1bSJed Brown 
466c4762a1bSJed Brown   outptr[mn - m] = xc * inptr[mn - m] + xr * inptr[mn - m + 1] + 2.0 * yl * inptr[mn - m - m];
467c4762a1bSJed Brown   outptr[mn - 1] = 2.0 * xl * inptr[mn - 2] + xc * inptr[mn - 1] + 2.0 * yl * inptr[mn - 1 - m];
468ad540459SPierre Jolivet   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];
469c4762a1bSJed Brown 
4709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_in, &inptr));
4719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(tmp_out, &outptr));
472c4762a1bSJed Brown 
4739566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(tmp_out, from, globalout, to, &scatter));
4749566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
4759566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
476c4762a1bSJed Brown 
477aeebefc2SPierre Jolivet   /* Destroy idx and scatter */
4789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_in));
4799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_out));
4809566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
4819566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
4829566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
483c4762a1bSJed Brown 
4849566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
4853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
486c4762a1bSJed Brown }
487c4762a1bSJed Brown 
PostStep(TS ts)488d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStep(TS ts)
489d71ae5a4SJacob Faibussowitsch {
490c4762a1bSJed Brown   PetscReal t;
491c4762a1bSJed Brown 
492c4762a1bSJed Brown   PetscFunctionBeginUser;
4939566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
4949566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  PostStep, t: %g\n", (double)t));
4953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
496c4762a1bSJed Brown }
497c4762a1bSJed Brown 
498c4762a1bSJed Brown /*TEST
499c4762a1bSJed Brown 
500c4762a1bSJed Brown     test:
501c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
502c4762a1bSJed Brown       output_file: output/ex4.out
503c4762a1bSJed Brown 
504c4762a1bSJed Brown     test:
505c4762a1bSJed Brown       suffix: 2
506c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
507c4762a1bSJed Brown       nsize: 2
508c4762a1bSJed Brown       output_file: output/ex4.out
509c4762a1bSJed Brown 
510c4762a1bSJed Brown     test:
511c4762a1bSJed Brown       suffix: 3
512c4762a1bSJed Brown       args: -ts_fd -ts_type cn
513c4762a1bSJed Brown 
514c4762a1bSJed Brown     test:
515c4762a1bSJed Brown       suffix: 4
516c4762a1bSJed Brown       args: -ts_fd -ts_type cn
517c4762a1bSJed Brown       output_file: output/ex4_3.out
518c4762a1bSJed Brown       nsize: 2
519c4762a1bSJed Brown 
520c4762a1bSJed Brown     test:
521c4762a1bSJed Brown       suffix: 5
522c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
523c4762a1bSJed Brown       output_file: output/ex4.out
524c4762a1bSJed Brown 
525c4762a1bSJed Brown     test:
526c4762a1bSJed Brown       suffix: 6
527c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
528c4762a1bSJed Brown       output_file: output/ex4.out
529c4762a1bSJed Brown       nsize: 2
530c4762a1bSJed Brown 
531c4762a1bSJed Brown     test:
532c4762a1bSJed Brown       suffix: 7
533c4762a1bSJed Brown       requires: !single
534188af4bfSBarry Smith       args: -ts_fd -ts_type beuler -test_PostStep -ts_time_step .1
535c4762a1bSJed Brown 
536c4762a1bSJed Brown     test:
537c4762a1bSJed Brown       suffix: 8
538c4762a1bSJed Brown       requires: !single
539188af4bfSBarry Smith       args: -ts_type rk -ts_rk_type 5dp -ts_time_step .01 -ts_adapt_type none -ts_view
540c4762a1bSJed Brown 
541c4762a1bSJed Brown TEST*/
542