xref: /petsc/src/ts/tests/ex4.c (revision 750b007cd8d816cecd9de99077bb0a703b4cf61a)
1 /*
2        The Problem:
3            Solve the convection-diffusion equation:
4 
5              u_t+a*(u_x+u_y)=epsilon*(u_xx+u_yy)
6              u=0   at x=0, y=0
7              u_x=0 at x=1
8              u_y=0 at y=1
9              u = exp(-20.0*(pow(x-0.5,2.0)+pow(y-0.5,2.0))) at t=0
10 
11        This program tests the routine of computing the Jacobian by the
12        finite difference method as well as PETSc.
13 
14 */
15 
16 static char help[] = "Solve the convection-diffusion equation. \n\n";
17 
18 #include <petscts.h>
19 
20 typedef struct {
21   PetscInt  m;       /* the number of mesh points in x-direction */
22   PetscInt  n;       /* the number of mesh points in y-direction */
23   PetscReal dx;      /* the grid space in x-direction */
24   PetscReal dy;      /* the grid space in y-direction */
25   PetscReal a;       /* the convection coefficient    */
26   PetscReal epsilon; /* the diffusion coefficient     */
27   PetscReal tfinal;
28 } Data;
29 
30 extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
31 extern PetscErrorCode Initial(Vec, void *);
32 extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
33 extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
34 extern PetscErrorCode PostStep(TS);
35 
36 int main(int argc, char **argv) {
37   PetscInt      time_steps = 100, iout, NOUT = 1;
38   Vec           global;
39   PetscReal     dt, ftime, ftime_original;
40   TS            ts;
41   PetscViewer   viewfile;
42   Mat           J = 0;
43   Vec           x;
44   Data          data;
45   PetscInt      mn;
46   PetscBool     flg;
47   MatColoring   mc;
48   ISColoring    iscoloring;
49   MatFDColoring matfdcoloring        = 0;
50   PetscBool     fd_jacobian_coloring = PETSC_FALSE;
51   SNES          snes;
52   KSP           ksp;
53   PC            pc;
54 
55   PetscFunctionBeginUser;
56   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
57 
58   /* set data */
59   data.m       = 9;
60   data.n       = 9;
61   data.a       = 1.0;
62   data.epsilon = 0.1;
63   data.dx      = 1.0 / (data.m + 1.0);
64   data.dy      = 1.0 / (data.n + 1.0);
65   mn           = (data.m) * (data.n);
66   PetscCall(PetscOptionsGetInt(NULL, NULL, "-time", &time_steps, NULL));
67 
68   /* set initial conditions */
69   PetscCall(VecCreate(PETSC_COMM_WORLD, &global));
70   PetscCall(VecSetSizes(global, PETSC_DECIDE, mn));
71   PetscCall(VecSetFromOptions(global));
72   PetscCall(Initial(global, &data));
73   PetscCall(VecDuplicate(global, &x));
74 
75   /* create timestep context */
76   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
77   PetscCall(TSMonitorSet(ts, Monitor, &data, NULL));
78   PetscCall(TSSetType(ts, TSEULER));
79   dt             = 0.1;
80   ftime_original = data.tfinal = 1.0;
81 
82   PetscCall(TSSetTimeStep(ts, dt));
83   PetscCall(TSSetMaxSteps(ts, time_steps));
84   PetscCall(TSSetMaxTime(ts, ftime_original));
85   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
86   PetscCall(TSSetSolution(ts, global));
87 
88   /* set user provided RHSFunction and RHSJacobian */
89   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &data));
90   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
91   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, mn, mn));
92   PetscCall(MatSetFromOptions(J));
93   PetscCall(MatSeqAIJSetPreallocation(J, 5, NULL));
94   PetscCall(MatMPIAIJSetPreallocation(J, 5, NULL, 5, NULL));
95 
96   PetscCall(PetscOptionsHasName(NULL, NULL, "-ts_fd", &flg));
97   if (!flg) {
98     PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, &data));
99   } else {
100     PetscCall(TSGetSNES(ts, &snes));
101     PetscCall(PetscOptionsHasName(NULL, NULL, "-fd_color", &fd_jacobian_coloring));
102     if (fd_jacobian_coloring) { /* Use finite differences with coloring */
103       /* Get data structure of J */
104       PetscBool pc_diagonal;
105       PetscCall(PetscOptionsHasName(NULL, NULL, "-pc_diagonal", &pc_diagonal));
106       if (pc_diagonal) { /* the preconditioner of J is a diagonal matrix */
107         PetscInt    rstart, rend, i;
108         PetscScalar zero = 0.0;
109         PetscCall(MatGetOwnershipRange(J, &rstart, &rend));
110         for (i = rstart; i < rend; i++) PetscCall(MatSetValues(J, 1, &i, 1, &i, &zero, INSERT_VALUES));
111         PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
112         PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
113       } else {
114         /* Fill the structure using the expensive SNESComputeJacobianDefault. Temporarily set up the TS so we can call this function */
115         PetscCall(TSSetType(ts, TSBEULER));
116         PetscCall(TSSetUp(ts));
117         PetscCall(SNESComputeJacobianDefault(snes, x, J, J, ts));
118       }
119 
120       /* create coloring context */
121       PetscCall(MatColoringCreate(J, &mc));
122       PetscCall(MatColoringSetType(mc, MATCOLORINGSL));
123       PetscCall(MatColoringSetFromOptions(mc));
124       PetscCall(MatColoringApply(mc, &iscoloring));
125       PetscCall(MatColoringDestroy(&mc));
126       PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring));
127       PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode(*)(void))SNESTSFormFunction, ts));
128       PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
129       PetscCall(MatFDColoringSetUp(J, iscoloring, matfdcoloring));
130       PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, matfdcoloring));
131       PetscCall(ISColoringDestroy(&iscoloring));
132     } else { /* Use finite differences (slow) */
133       PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefault, NULL));
134     }
135   }
136 
137   /* Pick up a Petsc preconditioner */
138   /* one can always set method or preconditioner during the run time */
139   PetscCall(TSGetSNES(ts, &snes));
140   PetscCall(SNESGetKSP(snes, &ksp));
141   PetscCall(KSPGetPC(ksp, &pc));
142   PetscCall(PCSetType(pc, PCJACOBI));
143   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
144 
145   PetscCall(TSSetFromOptions(ts));
146   PetscCall(TSSetUp(ts));
147 
148   /* Test TSSetPostStep() */
149   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_PostStep", &flg));
150   if (flg) PetscCall(TSSetPostStep(ts, PostStep));
151 
152   PetscCall(PetscOptionsGetInt(NULL, NULL, "-NOUT", &NOUT, NULL));
153   for (iout = 1; iout <= NOUT; iout++) {
154     PetscCall(TSSetMaxSteps(ts, time_steps));
155     PetscCall(TSSetMaxTime(ts, iout * ftime_original / NOUT));
156     PetscCall(TSSolve(ts, global));
157     PetscCall(TSGetSolveTime(ts, &ftime));
158     PetscCall(TSSetTime(ts, ftime));
159     PetscCall(TSSetTimeStep(ts, dt));
160   }
161   /* Interpolate solution at tfinal */
162   PetscCall(TSGetSolution(ts, &global));
163   PetscCall(TSInterpolate(ts, ftime_original, global));
164 
165   PetscCall(PetscOptionsHasName(NULL, NULL, "-matlab_view", &flg));
166   if (flg) { /* print solution into a MATLAB file */
167     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "out.m", &viewfile));
168     PetscCall(PetscViewerPushFormat(viewfile, PETSC_VIEWER_ASCII_MATLAB));
169     PetscCall(VecView(global, viewfile));
170     PetscCall(PetscViewerPopFormat(viewfile));
171     PetscCall(PetscViewerDestroy(&viewfile));
172   }
173 
174   /* free the memories */
175   PetscCall(TSDestroy(&ts));
176   PetscCall(VecDestroy(&global));
177   PetscCall(VecDestroy(&x));
178   PetscCall(MatDestroy(&J));
179   if (fd_jacobian_coloring) PetscCall(MatFDColoringDestroy(&matfdcoloring));
180   PetscCall(PetscFinalize());
181   return 0;
182 }
183 
184 /* -------------------------------------------------------------------*/
185 /* the initial function */
186 PetscReal f_ini(PetscReal x, PetscReal y) {
187   PetscReal f;
188 
189   f = PetscExpReal(-20.0 * (PetscPowRealInt(x - 0.5, 2) + PetscPowRealInt(y - 0.5, 2)));
190   return f;
191 }
192 
193 PetscErrorCode Initial(Vec global, void *ctx) {
194   Data        *data = (Data *)ctx;
195   PetscInt     m, row, col;
196   PetscReal    x, y, dx, dy;
197   PetscScalar *localptr;
198   PetscInt     i, mybase, myend, locsize;
199 
200   PetscFunctionBeginUser;
201   /* make the local  copies of parameters */
202   m  = data->m;
203   dx = data->dx;
204   dy = data->dy;
205 
206   /* determine starting point of each processor */
207   PetscCall(VecGetOwnershipRange(global, &mybase, &myend));
208   PetscCall(VecGetLocalSize(global, &locsize));
209 
210   /* Initialize the array */
211   PetscCall(VecGetArrayWrite(global, &localptr));
212 
213   for (i = 0; i < locsize; i++) {
214     row         = 1 + (mybase + i) - ((mybase + i) / m) * m;
215     col         = (mybase + i) / m + 1;
216     x           = dx * row;
217     y           = dy * col;
218     localptr[i] = f_ini(x, y);
219   }
220 
221   PetscCall(VecRestoreArrayWrite(global, &localptr));
222   PetscFunctionReturn(0);
223 }
224 
225 PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec global, void *ctx) {
226   VecScatter         scatter;
227   IS                 from, to;
228   PetscInt           i, n, *idx, nsteps, maxsteps;
229   Vec                tmp_vec;
230   const PetscScalar *tmp;
231 
232   PetscFunctionBeginUser;
233   PetscCall(TSGetStepNumber(ts, &nsteps));
234   /* display output at selected time steps */
235   PetscCall(TSGetMaxSteps(ts, &maxsteps));
236   if (nsteps % 10 != 0) PetscFunctionReturn(0);
237 
238   /* Get the size of the vector */
239   PetscCall(VecGetSize(global, &n));
240 
241   /* Set the index sets */
242   PetscCall(PetscMalloc1(n, &idx));
243   for (i = 0; i < n; i++) idx[i] = i;
244 
245   /* Create local sequential vectors */
246   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &tmp_vec));
247 
248   /* Create scatter context */
249   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &from));
250   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, idx, PETSC_COPY_VALUES, &to));
251   PetscCall(VecScatterCreate(global, from, tmp_vec, to, &scatter));
252   PetscCall(VecScatterBegin(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
253   PetscCall(VecScatterEnd(scatter, global, tmp_vec, INSERT_VALUES, SCATTER_FORWARD));
254 
255   PetscCall(VecGetArrayRead(tmp_vec, &tmp));
256   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])));
257   PetscCall(VecRestoreArrayRead(tmp_vec, &tmp));
258 
259   PetscCall(PetscFree(idx));
260   PetscCall(ISDestroy(&from));
261   PetscCall(ISDestroy(&to));
262   PetscCall(VecScatterDestroy(&scatter));
263   PetscCall(VecDestroy(&tmp_vec));
264   PetscFunctionReturn(0);
265 }
266 
267 PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat BB, void *ptr) {
268   Data       *data = (Data *)ptr;
269   PetscScalar v[5];
270   PetscInt    idx[5], i, j, row;
271   PetscInt    m, n, mn;
272   PetscReal   dx, dy, a, epsilon, xc, xl, xr, yl, yr;
273 
274   PetscFunctionBeginUser;
275   m       = data->m;
276   n       = data->n;
277   mn      = m * n;
278   dx      = data->dx;
279   dy      = data->dy;
280   a       = data->a;
281   epsilon = data->epsilon;
282 
283   xc = -2.0 * epsilon * (1.0 / (dx * dx) + 1.0 / (dy * dy));
284   xl = 0.5 * a / dx + epsilon / (dx * dx);
285   xr = -0.5 * a / dx + epsilon / (dx * dx);
286   yl = 0.5 * a / dy + epsilon / (dy * dy);
287   yr = -0.5 * a / dy + epsilon / (dy * dy);
288 
289   row    = 0;
290   v[0]   = xc;
291   v[1]   = xr;
292   v[2]   = yr;
293   idx[0] = 0;
294   idx[1] = 2;
295   idx[2] = m;
296   PetscCall(MatSetValues(A, 1, &row, 3, idx, v, INSERT_VALUES));
297 
298   row    = m - 1;
299   v[0]   = 2.0 * xl;
300   v[1]   = xc;
301   v[2]   = yr;
302   idx[0] = m - 2;
303   idx[1] = m - 1;
304   idx[2] = m - 1 + m;
305   PetscCall(MatSetValues(A, 1, &row, 3, idx, v, INSERT_VALUES));
306 
307   for (i = 1; i < m - 1; i++) {
308     row    = i;
309     v[0]   = xl;
310     v[1]   = xc;
311     v[2]   = xr;
312     v[3]   = yr;
313     idx[0] = i - 1;
314     idx[1] = i;
315     idx[2] = i + 1;
316     idx[3] = i + m;
317     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
318   }
319 
320   for (j = 1; j < n - 1; j++) {
321     row    = j * m;
322     v[0]   = xc;
323     v[1]   = xr;
324     v[2]   = yl;
325     v[3]   = yr;
326     idx[0] = j * m;
327     idx[1] = j * m;
328     idx[2] = j * m - m;
329     idx[3] = j * m + m;
330     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
331 
332     row    = j * m + m - 1;
333     v[0]   = xc;
334     v[1]   = 2.0 * xl;
335     v[2]   = yl;
336     v[3]   = yr;
337     idx[0] = j * m + m - 1;
338     idx[1] = j * m + m - 1 - 1;
339     idx[2] = j * m + m - 1 - m;
340     idx[3] = j * m + m - 1 + m;
341     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
342 
343     for (i = 1; i < m - 1; i++) {
344       row    = j * m + i;
345       v[0]   = xc;
346       v[1]   = xl;
347       v[2]   = xr;
348       v[3]   = yl;
349       v[4]   = yr;
350       idx[0] = j * m + i;
351       idx[1] = j * m + i - 1;
352       idx[2] = j * m + i + 1;
353       idx[3] = j * m + i - m;
354       idx[4] = j * m + i + m;
355       PetscCall(MatSetValues(A, 1, &row, 5, idx, v, INSERT_VALUES));
356     }
357   }
358 
359   row    = mn - m;
360   v[0]   = xc;
361   v[1]   = xr;
362   v[2]   = 2.0 * yl;
363   idx[0] = mn - m;
364   idx[1] = mn - m + 1;
365   idx[2] = mn - m - m;
366   PetscCall(MatSetValues(A, 1, &row, 3, idx, v, INSERT_VALUES));
367 
368   row    = mn - 1;
369   v[0]   = xc;
370   v[1]   = 2.0 * xl;
371   v[2]   = 2.0 * yl;
372   idx[0] = mn - 1;
373   idx[1] = mn - 2;
374   idx[2] = mn - 1 - m;
375   PetscCall(MatSetValues(A, 1, &i, 3, idx, v, INSERT_VALUES));
376 
377   for (i = 1; i < m - 1; i++) {
378     row    = mn - m + i;
379     v[0]   = xl;
380     v[1]   = xc;
381     v[2]   = xr;
382     v[3]   = 2.0 * yl;
383     idx[0] = mn - m + i - 1;
384     idx[1] = mn - m + i;
385     idx[2] = mn - m + i + 1;
386     idx[3] = mn - m + i - m;
387     PetscCall(MatSetValues(A, 1, &row, 4, idx, v, INSERT_VALUES));
388   }
389 
390   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
391   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
392 
393   PetscFunctionReturn(0);
394 }
395 
396 /* globalout = -a*(u_x+u_y) + epsilon*(u_xx+u_yy) */
397 PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx) {
398   Data              *data = (Data *)ctx;
399   PetscInt           m, n, mn;
400   PetscReal          dx, dy;
401   PetscReal          xc, xl, xr, yl, yr;
402   PetscReal          a, epsilon;
403   PetscScalar       *outptr;
404   const PetscScalar *inptr;
405   PetscInt           i, j, len;
406   IS                 from, to;
407   PetscInt          *idx;
408   VecScatter         scatter;
409   Vec                tmp_in, tmp_out;
410 
411   PetscFunctionBeginUser;
412   m       = data->m;
413   n       = data->n;
414   mn      = m * n;
415   dx      = data->dx;
416   dy      = data->dy;
417   a       = data->a;
418   epsilon = data->epsilon;
419 
420   xc = -2.0 * epsilon * (1.0 / (dx * dx) + 1.0 / (dy * dy));
421   xl = 0.5 * a / dx + epsilon / (dx * dx);
422   xr = -0.5 * a / dx + epsilon / (dx * dx);
423   yl = 0.5 * a / dy + epsilon / (dy * dy);
424   yr = -0.5 * a / dy + epsilon / (dy * dy);
425 
426   /* Get the length of parallel vector */
427   PetscCall(VecGetSize(globalin, &len));
428 
429   /* Set the index sets */
430   PetscCall(PetscMalloc1(len, &idx));
431   for (i = 0; i < len; i++) idx[i] = i;
432 
433   /* Create local sequential vectors */
434   PetscCall(VecCreateSeq(PETSC_COMM_SELF, len, &tmp_in));
435   PetscCall(VecDuplicate(tmp_in, &tmp_out));
436 
437   /* Create scatter context */
438   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, len, idx, PETSC_COPY_VALUES, &from));
439   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, len, idx, PETSC_COPY_VALUES, &to));
440   PetscCall(VecScatterCreate(globalin, from, tmp_in, to, &scatter));
441   PetscCall(VecScatterBegin(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
442   PetscCall(VecScatterEnd(scatter, globalin, tmp_in, INSERT_VALUES, SCATTER_FORWARD));
443   PetscCall(VecScatterDestroy(&scatter));
444 
445   /*Extract income array - include ghost points */
446   PetscCall(VecGetArrayRead(tmp_in, &inptr));
447 
448   /* Extract outcome array*/
449   PetscCall(VecGetArrayWrite(tmp_out, &outptr));
450 
451   outptr[0]     = xc * inptr[0] + xr * inptr[1] + yr * inptr[m];
452   outptr[m - 1] = 2.0 * xl * inptr[m - 2] + xc * inptr[m - 1] + yr * inptr[m - 1 + m];
453   for (i = 1; i < m - 1; i++) outptr[i] = xc * inptr[i] + xl * inptr[i - 1] + xr * inptr[i + 1] + yr * inptr[i + m];
454 
455   for (j = 1; j < n - 1; j++) {
456     outptr[j * m]         = xc * inptr[j * m] + xr * inptr[j * m + 1] + yl * inptr[j * m - m] + yr * inptr[j * m + m];
457     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];
458     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];
459   }
460 
461   outptr[mn - m] = xc * inptr[mn - m] + xr * inptr[mn - m + 1] + 2.0 * yl * inptr[mn - m - m];
462   outptr[mn - 1] = 2.0 * xl * inptr[mn - 2] + xc * inptr[mn - 1] + 2.0 * yl * inptr[mn - 1 - m];
463   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];
464 
465   PetscCall(VecRestoreArrayRead(tmp_in, &inptr));
466   PetscCall(VecRestoreArrayWrite(tmp_out, &outptr));
467 
468   PetscCall(VecScatterCreate(tmp_out, from, globalout, to, &scatter));
469   PetscCall(VecScatterBegin(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
470   PetscCall(VecScatterEnd(scatter, tmp_out, globalout, INSERT_VALUES, SCATTER_FORWARD));
471 
472   /* Destroy idx aand scatter */
473   PetscCall(VecDestroy(&tmp_in));
474   PetscCall(VecDestroy(&tmp_out));
475   PetscCall(ISDestroy(&from));
476   PetscCall(ISDestroy(&to));
477   PetscCall(VecScatterDestroy(&scatter));
478 
479   PetscCall(PetscFree(idx));
480   PetscFunctionReturn(0);
481 }
482 
483 PetscErrorCode PostStep(TS ts) {
484   PetscReal t;
485 
486   PetscFunctionBeginUser;
487   PetscCall(TSGetTime(ts, &t));
488   PetscCall(PetscPrintf(PETSC_COMM_SELF, "  PostStep, t: %g\n", (double)t));
489   PetscFunctionReturn(0);
490 }
491 
492 /*TEST
493 
494     test:
495       args: -ts_fd -ts_type beuler
496       output_file: output/ex4.out
497 
498     test:
499       suffix: 2
500       args: -ts_fd -ts_type beuler
501       nsize: 2
502       output_file: output/ex4.out
503 
504     test:
505       suffix: 3
506       args: -ts_fd -ts_type cn
507 
508     test:
509       suffix: 4
510       args: -ts_fd -ts_type cn
511       output_file: output/ex4_3.out
512       nsize: 2
513 
514     test:
515       suffix: 5
516       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
517       output_file: output/ex4.out
518 
519     test:
520       suffix: 6
521       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
522       output_file: output/ex4.out
523       nsize: 2
524 
525     test:
526       suffix: 7
527       requires: !single
528       args: -ts_fd -ts_type beuler -test_PostStep -ts_dt .1
529 
530     test:
531       suffix: 8
532       requires: !single
533       args: -ts_type rk -ts_rk_type 5dp -ts_dt .01 -ts_adapt_type none -ts_view
534 
535 TEST*/
536