Actual source code: ex46.c

  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";

  8: #include <petscdmplex.h>
  9: #include <petscsnes.h>
 10: #include <petscts.h>
 11: #include <petscds.h>

 13: /*
 14:   Navier-Stokes equation:

 16:   du/dt + u . grad u - \Delta u - grad p = f
 17:   div u  = 0
 18: */

 20: typedef struct {
 21:   PetscInt mms;
 22: } AppCtx;

 24: #define REYN 400.0

 26: /* MMS1

 28:   u = t + x^2 + y^2;
 29:   v = t + 2*x^2 - 2*x*y;
 30:   p = x + y - 1;

 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

 35:   so that

 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

 41:   where

 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 PETSC_SUCCESS;
 50: }

 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 PETSC_SUCCESS;
 56: }

 58: /* MMS 2*/

 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 PETSC_SUCCESS;
 65: }

 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 PETSC_SUCCESS;
 71: }

 73: static void f0_mms1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
 74: {
 75:   const PetscReal Re    = REYN;
 76:   const PetscInt  Ncomp = dim;
 77:   PetscInt        c, d;

 79:   for (c = 0; c < Ncomp; ++c) {
 80:     for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d];
 81:   }
 82:   f0[0] += u_t[0];
 83:   f0[1] += u_t[1];

 85:   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;
 86:   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;
 87: }

 89: static void f0_mms2_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
 90: {
 91:   const PetscReal Re    = REYN;
 92:   const PetscInt  Ncomp = dim;
 93:   PetscInt        c, d;

 95:   for (c = 0; c < Ncomp; ++c) {
 96:     for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d];
 97:   }
 98:   f0[0] += u_t[0];
 99:   f0[1] += u_t[1];

101:   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;
102:   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;
103: }

105: static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
106: {
107:   const PetscReal Re    = REYN;
108:   const PetscInt  Ncomp = dim;
109:   PetscInt        comp, d;

111:   for (comp = 0; comp < Ncomp; ++comp) {
112:     for (d = 0; d < dim; ++d) f1[comp * dim + d] = 1.0 / Re * u_x[comp * dim + d];
113:     f1[comp * dim + comp] -= u[Ncomp];
114:   }
115: }

117: static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
118: {
119:   PetscInt d;
120:   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d];
121: }

123: static void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
124: {
125:   PetscInt d;
126:   for (d = 0; d < dim; ++d) f1[d] = 0.0;
127: }

129: /*
130:   (psi_i, u_j grad_j u_i) ==> (\psi_i, \phi_j grad_j u_i)
131: */
132: static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
133: {
134:   PetscInt NcI = dim, NcJ = dim;
135:   PetscInt fc, gc;
136:   PetscInt d;

138:   for (d = 0; d < dim; ++d) g0[d * dim + d] = u_tShift;

140:   for (fc = 0; fc < NcI; ++fc) {
141:     for (gc = 0; gc < NcJ; ++gc) g0[fc * NcJ + gc] += u_x[fc * NcJ + gc];
142:   }
143: }

145: /*
146:   (psi_i, u_j grad_j u_i) ==> (\psi_i, \u_j grad_j \phi_i)
147: */
148: static void g1_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
149: {
150:   PetscInt NcI = dim;
151:   PetscInt NcJ = dim;
152:   PetscInt fc, gc, dg;
153:   for (fc = 0; fc < NcI; ++fc) {
154:     for (gc = 0; gc < NcJ; ++gc) {
155:       for (dg = 0; dg < dim; ++dg) {
156:         /* kronecker delta */
157:         if (fc == gc) g1[(fc * NcJ + gc) * dim + dg] += u[dg];
158:       }
159:     }
160:   }
161: }

163: /* < q, \nabla\cdot u >
164:    NcompI = 1, NcompJ = dim */
165: static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
166: {
167:   PetscInt d;
168:   for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
169: }

171: /* -< \nabla\cdot v, p >
172:     NcompI = dim, NcompJ = 1 */
173: static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
174: {
175:   PetscInt d;
176:   for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
177: }

179: /* < \nabla v, \nabla u + {\nabla u}^T >
180:    This just gives \nabla u, give the perdiagonal for the transpose */
181: static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
182: {
183:   const PetscReal Re    = REYN;
184:   const PetscInt  Ncomp = dim;
185:   PetscInt        compI, d;

187:   for (compI = 0; compI < Ncomp; ++compI) {
188:     for (d = 0; d < dim; ++d) g3[((compI * Ncomp + compI) * dim + d) * dim + d] = 1.0 / Re;
189:   }
190: }

192: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
193: {
194:   PetscFunctionBeginUser;
195:   options->mms = 1;

197:   PetscOptionsBegin(comm, "", "Navier-Stokes Equation Options", "DMPLEX");
198:   PetscCall(PetscOptionsInt("-mms", "The manufactured solution to use", "ex46.c", options->mms, &options->mms, NULL));
199:   PetscOptionsEnd();
200:   PetscFunctionReturn(PETSC_SUCCESS);
201: }

203: static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
204: {
205:   PetscFunctionBeginUser;
206:   PetscCall(DMCreate(comm, dm));
207:   PetscCall(DMSetType(*dm, DMPLEX));
208:   PetscCall(DMSetFromOptions(*dm));
209:   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
210:   PetscFunctionReturn(PETSC_SUCCESS);
211: }

213: static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
214: {
215:   PetscDS        ds;
216:   DMLabel        label;
217:   const PetscInt id = 1;
218:   PetscInt       dim;

220:   PetscFunctionBeginUser;
221:   PetscCall(DMGetDimension(dm, &dim));
222:   PetscCall(DMGetDS(dm, &ds));
223:   PetscCall(DMGetLabel(dm, "marker", &label));
224:   switch (dim) {
225:   case 2:
226:     switch (ctx->mms) {
227:     case 1:
228:       PetscCall(PetscDSSetResidual(ds, 0, f0_mms1_u, f1_u));
229:       PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p));
230:       PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu));
231:       PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
232:       PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL));
233:       PetscCall(PetscDSSetExactSolution(ds, 0, mms1_u_2d, ctx));
234:       PetscCall(PetscDSSetExactSolution(ds, 1, mms1_p_2d, ctx));
235:       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))mms1_u_2d, NULL, ctx, NULL));
236:       break;
237:     case 2:
238:       PetscCall(PetscDSSetResidual(ds, 0, f0_mms2_u, f1_u));
239:       PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p));
240:       PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu));
241:       PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
242:       PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL));
243:       PetscCall(PetscDSSetExactSolution(ds, 0, mms2_u_2d, ctx));
244:       PetscCall(PetscDSSetExactSolution(ds, 1, mms2_p_2d, ctx));
245:       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))mms2_u_2d, NULL, ctx, NULL));
246:       break;
247:     default:
248:       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid MMS %" PetscInt_FMT, ctx->mms);
249:     }
250:     break;
251:   default:
252:     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
253:   }
254:   PetscFunctionReturn(PETSC_SUCCESS);
255: }

257: static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
258: {
259:   MPI_Comm  comm;
260:   DM        cdm = dm;
261:   PetscFE   fe[2];
262:   PetscInt  dim;
263:   PetscBool simplex;

265:   PetscFunctionBeginUser;
266:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
267:   PetscCall(DMGetDimension(dm, &dim));
268:   PetscCall(DMPlexIsSimplex(dm, &simplex));
269:   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
270:   PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity"));
271:   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
272:   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
273:   PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure"));
274:   /* Set discretization and boundary conditions for each mesh */
275:   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
276:   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
277:   PetscCall(DMCreateDS(dm));
278:   PetscCall(SetupProblem(dm, ctx));
279:   while (cdm) {
280:     PetscObject  pressure;
281:     MatNullSpace nsp;

283:     PetscCall(DMGetField(cdm, 1, NULL, &pressure));
284:     PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nsp));
285:     PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nsp));
286:     PetscCall(MatNullSpaceDestroy(&nsp));

288:     PetscCall(DMCopyDisc(dm, cdm));
289:     PetscCall(DMGetCoarseDM(cdm, &cdm));
290:   }
291:   PetscCall(PetscFEDestroy(&fe[0]));
292:   PetscCall(PetscFEDestroy(&fe[1]));
293:   PetscFunctionReturn(PETSC_SUCCESS);
294: }

296: static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
297: {
298:   PetscSimplePointFn *funcs[2];
299:   void               *ctxs[2];
300:   DM                  dm;
301:   PetscDS             ds;
302:   PetscReal           ferrors[2];

304:   PetscFunctionBeginUser;
305:   PetscCall(TSGetDM(ts, &dm));
306:   PetscCall(DMGetDS(dm, &ds));
307:   PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
308:   PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
309:   PetscCall(DMComputeL2FieldDiff(dm, crtime, funcs, ctxs, u, ferrors));
310:   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]));
311:   PetscFunctionReturn(PETSC_SUCCESS);
312: }

314: int main(int argc, char **argv)
315: {
316:   AppCtx ctx;
317:   DM     dm;
318:   TS     ts;
319:   Vec    u, r;

321:   PetscFunctionBeginUser;
322:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
323:   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
324:   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
325:   PetscCall(DMSetApplicationContext(dm, &ctx));
326:   PetscCall(SetupDiscretization(dm, &ctx));
327:   PetscCall(DMPlexCreateClosureIndex(dm, NULL));

329:   PetscCall(DMCreateGlobalVector(dm, &u));
330:   PetscCall(VecDuplicate(u, &r));

332:   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
333:   PetscCall(TSMonitorSet(ts, MonitorError, &ctx, NULL));
334:   PetscCall(TSSetDM(ts, dm));
335:   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
336:   PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
337:   PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
338:   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
339:   PetscCall(TSSetFromOptions(ts));
340:   PetscCall(DMTSCheckFromOptions(ts, u));

342:   {
343:     PetscSimplePointFn *funcs[2];
344:     void               *ctxs[2];
345:     PetscDS             ds;

347:     PetscCall(DMGetDS(dm, &ds));
348:     PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
349:     PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
350:     PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, u));
351:   }
352:   PetscCall(TSSolve(ts, u));
353:   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));

355:   PetscCall(VecDestroy(&u));
356:   PetscCall(VecDestroy(&r));
357:   PetscCall(TSDestroy(&ts));
358:   PetscCall(DMDestroy(&dm));
359:   PetscCall(PetscFinalize());
360:   return 0;
361: }

363: /*TEST

365:   # Full solves
366:   test:
367:     suffix: 2d_p2p1_r1
368:     requires: !single triangle
369:     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
370:     args: -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
371:           -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \
372:           -snes_monitor_short -snes_converged_reason \
373:           -ksp_monitor_short -ksp_converged_reason \
374:           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
375:             -fieldsplit_velocity_pc_type lu \
376:             -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi

378:   test:
379:     suffix: 2d_q2q1_r1
380:     requires: !single
381:     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" -e "s~ 0\]~ 0.0\]~g"
382:     args: -dm_plex_simplex 0 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
383:           -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \
384:           -snes_monitor_short -snes_converged_reason \
385:           -ksp_monitor_short -ksp_converged_reason \
386:           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
387:             -fieldsplit_velocity_pc_type lu \
388:             -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi

390: TEST*/