xref: /petsc/src/ts/tutorials/ex46.c (revision f97672e55eacc8688507b9471cd7ec2664d7f203)
1 static char help[] = "Time dependent Navier-Stokes problem in 2d and 3d with finite elements.\n\
2 We solve the Navier-Stokes in a rectangular\n\
3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4 This example supports discretized auxiliary fields (Re) as well as\n\
5 multilevel nonlinear solvers.\n\
6 Contributed by: Julian Andrej <juan@tf.uni-kiel.de>\n\n\n";
7 
8 #include <petscdmplex.h>
9 #include <petscsnes.h>
10 #include <petscts.h>
11 #include <petscds.h>
12 
13 /*
14   Navier-Stokes equation:
15 
16   du/dt + u . grad u - \Delta u - grad p = f
17   div u  = 0
18 */
19 
20 typedef struct {
21   PetscInt mms;
22 } AppCtx;
23 
24 #define REYN 400.0
25 
26 /* MMS1
27 
28   u = t + x^2 + y^2;
29   v = t + 2*x^2 - 2*x*y;
30   p = x + y - 1;
31 
32   f_x = -2*t*(x + y) + 2*x*y^2 - 4*x^2*y - 2*x^3 + 4.0/Re - 1.0
33   f_y = -2*t*x       + 2*y^3 - 4*x*y^2 - 2*x^2*y + 4.0/Re - 1.0
34 
35   so that
36 
37     u_t + u \cdot \nabla u - 1/Re \Delta u + \nabla p + f = <1, 1> + <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t 2x + 2x^2y + 4xy^2 - 2y^3> - 1/Re <4, 4> + <1, 1>
38                                                     + <-t (2x + 2y) + 2xy^2 - 4x^2y - 2x^3 + 4/Re - 1, -2xt + 2y^3 - 4xy^2 - 2x^2y + 4/Re - 1> = 0
39     \nabla \cdot u                                  = 2x - 2x = 0
40 
41   where
42 
43     <u, v> . <<u_x, v_x>, <u_y, v_y>> = <u u_x + v u_y, u v_x + v v_y>
44 */
45 PetscErrorCode mms1_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
46 {
47   u[0] = time + x[0]*x[0] + x[1]*x[1];
48   u[1] = time + 2.0*x[0]*x[0] - 2.0*x[0]*x[1];
49   return 0;
50 }
51 
52 PetscErrorCode mms1_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx)
53 {
54   *p = x[0] + x[1] - 1.0;
55   return 0;
56 }
57 
58 /* MMS 2*/
59 
60 static PetscErrorCode mms2_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
61 {
62   u[0] = PetscSinReal(time + x[0])*PetscSinReal(time + x[1]);
63   u[1] = PetscCosReal(time + x[0])*PetscCosReal(time + x[1]);
64   return 0;
65 }
66 
67 static PetscErrorCode mms2_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx)
68 {
69   *p = PetscSinReal(time + x[0] - x[1]);
70   return 0;
71 }
72 
73 static void f0_mms1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
74                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
75                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
76                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
77 {
78   const PetscReal Re    = REYN;
79   const PetscInt  Ncomp = dim;
80   PetscInt        c, d;
81 
82   for (c = 0; c < Ncomp; ++c) {
83     for (d = 0; d < dim; ++d) {
84       f0[c] += u[d] * u_x[c*dim+d];
85     }
86   }
87   f0[0] += u_t[0];
88   f0[1] += u_t[1];
89 
90   f0[0] += -2.0*t*(x[0] + x[1]) + 2.0*x[0]*x[1]*x[1] - 4.0*x[0]*x[0]*x[1] - 2.0*x[0]*x[0]*x[0] + 4.0/Re - 1.0;
91   f0[1] += -2.0*t*x[0]          + 2.0*x[1]*x[1]*x[1] - 4.0*x[0]*x[1]*x[1] - 2.0*x[0]*x[0]*x[1] + 4.0/Re - 1.0;
92 }
93 
94 static void f0_mms2_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
98 {
99   const PetscReal Re    = REYN;
100   const PetscInt  Ncomp = dim;
101   PetscInt        c, d;
102 
103   for (c = 0; c < Ncomp; ++c) {
104     for (d = 0; d < dim; ++d) {
105       f0[c] += u[d] * u_x[c*dim+d];
106     }
107   }
108   f0[0] += u_t[0];
109   f0[1] += u_t[1];
110 
111   f0[0] -= ( Re*((1.0L/2.0L)*PetscSinReal(2*t + 2*x[0]) + PetscSinReal(2*t + x[0] + x[1]) + PetscCosReal(t + x[0] - x[1])) + 2.0*PetscSinReal(t + x[0])*PetscSinReal(t + x[1]))/Re;
112   f0[1] -= (-Re*((1.0L/2.0L)*PetscSinReal(2*t + 2*x[1]) + PetscSinReal(2*t + x[0] + x[1]) + PetscCosReal(t + x[0] - x[1])) + 2.0*PetscCosReal(t + x[0])*PetscCosReal(t + x[1]))/Re;
113 }
114 
115 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
116                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
117                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
118                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
119 {
120   const PetscReal Re    = REYN;
121   const PetscInt  Ncomp = dim;
122   PetscInt        comp, d;
123 
124   for (comp = 0; comp < Ncomp; ++comp) {
125     for (d = 0; d < dim; ++d) {
126       f1[comp*dim+d] = 1.0/Re * u_x[comp*dim+d];
127     }
128     f1[comp*dim+comp] -= u[Ncomp];
129   }
130 }
131 
132 static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
133                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
134                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
135                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
136 {
137   PetscInt d;
138   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
139 }
140 
141 static void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
142                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
143                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
144                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
145 {
146   PetscInt d;
147   for (d = 0; d < dim; ++d) f1[d] = 0.0;
148 }
149 
150 /*
151   (psi_i, u_j grad_j u_i) ==> (\psi_i, \phi_j grad_j u_i)
152 */
153 static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
154                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
155                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
156                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
157 {
158   PetscInt NcI = dim, NcJ = dim;
159   PetscInt fc, gc;
160   PetscInt d;
161 
162   for (d = 0; d < dim; ++d) {
163     g0[d*dim+d] = u_tShift;
164   }
165 
166   for (fc = 0; fc < NcI; ++fc) {
167     for (gc = 0; gc < NcJ; ++gc) {
168       g0[fc*NcJ+gc] += u_x[fc*NcJ+gc];
169     }
170   }
171 }
172 
173 /*
174   (psi_i, u_j grad_j u_i) ==> (\psi_i, \u_j grad_j \phi_i)
175 */
176 static void g1_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
177                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
178                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
179                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
180 {
181   PetscInt NcI = dim;
182   PetscInt NcJ = dim;
183   PetscInt fc, gc, dg;
184   for (fc = 0; fc < NcI; ++fc) {
185     for (gc = 0; gc < NcJ; ++gc) {
186       for (dg = 0; dg < dim; ++dg) {
187         /* kronecker delta */
188         if (fc == gc) {
189           g1[(fc*NcJ+gc)*dim+dg] += u[dg];
190         }
191       }
192     }
193   }
194 }
195 
196 /* < q, \nabla\cdot u >
197    NcompI = 1, NcompJ = dim */
198 static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
199                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
200                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
201                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
202 {
203   PetscInt d;
204   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
205 }
206 
207 /* -< \nabla\cdot v, p >
208     NcompI = dim, NcompJ = 1 */
209 static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
210                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
211                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
212                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
213 {
214   PetscInt d;
215   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
216 }
217 
218 /* < \nabla v, \nabla u + {\nabla u}^T >
219    This just gives \nabla u, give the perdiagonal for the transpose */
220 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
221                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
222                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
223                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
224 {
225   const PetscReal Re    = REYN;
226   const PetscInt  Ncomp = dim;
227   PetscInt        compI, d;
228 
229   for (compI = 0; compI < Ncomp; ++compI) {
230     for (d = 0; d < dim; ++d) {
231       g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0/Re;
232     }
233   }
234 }
235 
236 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
237 {
238   PetscFunctionBeginUser;
239   options->mms = 1;
240 
241   PetscOptionsBegin(comm, "", "Navier-Stokes Equation Options", "DMPLEX");
242   PetscCall(PetscOptionsInt("-mms", "The manufactured solution to use", "ex46.c", options->mms, &options->mms, NULL));
243   PetscOptionsEnd();
244   PetscFunctionReturn(0);
245 }
246 
247 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
248 {
249   PetscFunctionBeginUser;
250   PetscCall(DMCreate(comm, dm));
251   PetscCall(DMSetType(*dm, DMPLEX));
252   PetscCall(DMSetFromOptions(*dm));
253   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
254   PetscFunctionReturn(0);
255 }
256 
257 static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
258 {
259   PetscDS        ds;
260   DMLabel        label;
261   const PetscInt id = 1;
262   PetscInt       dim;
263 
264   PetscFunctionBeginUser;
265   PetscCall(DMGetDimension(dm, &dim));
266   PetscCall(DMGetDS(dm, &ds));
267   PetscCall(DMGetLabel(dm, "marker", &label));
268   switch (dim) {
269   case 2:
270     switch (ctx->mms) {
271     case 1:
272       PetscCall(PetscDSSetResidual(ds, 0, f0_mms1_u, f1_u));
273       PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p));
274       PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL,  g3_uu));
275       PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
276       PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL));
277       PetscCall(PetscDSSetExactSolution(ds, 0, mms1_u_2d, ctx));
278       PetscCall(PetscDSSetExactSolution(ds, 1, mms1_p_2d, ctx));
279       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms1_u_2d, NULL, ctx, NULL));
280       break;
281     case 2:
282       PetscCall(PetscDSSetResidual(ds, 0, f0_mms2_u, f1_u));
283       PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p));
284       PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL,  g3_uu));
285       PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
286       PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL));
287       PetscCall(PetscDSSetExactSolution(ds, 0, mms2_u_2d, ctx));
288       PetscCall(PetscDSSetExactSolution(ds, 1, mms2_p_2d, ctx));
289       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms2_u_2d, NULL, ctx, NULL));
290       break;
291     default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid MMS %" PetscInt_FMT, ctx->mms);
292     }
293     break;
294   default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
295   }
296   PetscFunctionReturn(0);
297 }
298 
299 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
300 {
301   MPI_Comm        comm;
302   DM              cdm = dm;
303   PetscFE         fe[2];
304   PetscInt        dim;
305   PetscBool       simplex;
306 
307   PetscFunctionBeginUser;
308   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
309   PetscCall(DMGetDimension(dm, &dim));
310   PetscCall(DMPlexIsSimplex(dm, &simplex));
311   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
312   PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity"));
313   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
314   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
315   PetscCall(PetscObjectSetName((PetscObject) fe[1], "pressure"));
316   /* Set discretization and boundary conditions for each mesh */
317   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0]));
318   PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1]));
319   PetscCall(DMCreateDS(dm));
320   PetscCall(SetupProblem(dm, ctx));
321   while (cdm) {
322     PetscObject  pressure;
323     MatNullSpace nsp;
324 
325     PetscCall(DMGetField(cdm, 1, NULL, &pressure));
326     PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nsp));
327     PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject) nsp));
328     PetscCall(MatNullSpaceDestroy(&nsp));
329 
330     PetscCall(DMCopyDisc(dm, cdm));
331     PetscCall(DMGetCoarseDM(cdm, &cdm));
332   }
333   PetscCall(PetscFEDestroy(&fe[0]));
334   PetscCall(PetscFEDestroy(&fe[1]));
335   PetscFunctionReturn(0);
336 }
337 
338 static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
339 {
340   PetscSimplePointFunc funcs[2];
341   void                *ctxs[2];
342   DM                   dm;
343   PetscDS              ds;
344   PetscReal            ferrors[2];
345 
346   PetscFunctionBeginUser;
347   PetscCall(TSGetDM(ts, &dm));
348   PetscCall(DMGetDS(dm, &ds));
349   PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
350   PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
351   PetscCall(DMComputeL2FieldDiff(dm, crtime, funcs, ctxs, u, ferrors));
352   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g]\n", (int) step, (double) crtime, (double) ferrors[0], (double) ferrors[1]));
353   PetscFunctionReturn(0);
354 }
355 
356 int main(int argc, char **argv)
357 {
358   AppCtx         ctx;
359   DM             dm;
360   TS             ts;
361   Vec            u, r;
362 
363   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
364   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
365   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
366   PetscCall(DMSetApplicationContext(dm, &ctx));
367   PetscCall(SetupDiscretization(dm, &ctx));
368   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
369 
370   PetscCall(DMCreateGlobalVector(dm, &u));
371   PetscCall(VecDuplicate(u, &r));
372 
373   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
374   PetscCall(TSMonitorSet(ts, MonitorError, &ctx, NULL));
375   PetscCall(TSSetDM(ts, dm));
376   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
377   PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
378   PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
379   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
380   PetscCall(TSSetFromOptions(ts));
381   PetscCall(DMTSCheckFromOptions(ts, u));
382 
383   {
384     PetscSimplePointFunc funcs[2];
385     void                *ctxs[2];
386     PetscDS              ds;
387 
388     PetscCall(DMGetDS(dm, &ds));
389     PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
390     PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
391     PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, u));
392   }
393   PetscCall(TSSolve(ts, u));
394   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
395 
396   PetscCall(VecDestroy(&u));
397   PetscCall(VecDestroy(&r));
398   PetscCall(TSDestroy(&ts));
399   PetscCall(DMDestroy(&dm));
400   PetscCall(PetscFinalize());
401   return 0;
402 }
403 
404 /*TEST
405 
406   # Full solves
407   test:
408     suffix: 2d_p2p1_r1
409     requires: !single triangle
410     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
411     args: -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
412           -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \
413           -snes_monitor_short -snes_converged_reason \
414           -ksp_monitor_short -ksp_converged_reason \
415           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
416             -fieldsplit_velocity_pc_type lu \
417             -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi
418 
419   test:
420     suffix: 2d_q2q1_r1
421     requires: !single
422     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" -e "s~ 0\]~ 0.0\]~g"
423     args: -dm_plex_simplex 0 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
424           -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \
425           -snes_monitor_short -snes_converged_reason \
426           -ksp_monitor_short -ksp_converged_reason \
427           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
428             -fieldsplit_velocity_pc_type lu \
429             -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi
430 
431 TEST*/
432