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, PetscCtx 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, PetscCtx 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, PetscCtx 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, PetscCtx 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 d;
79: for (PetscInt 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 d;
95: for (PetscInt 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 d;
111: for (PetscInt 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;
136: for (PetscInt d = 0; d < dim; ++d) g0[d * dim + d] = u_tShift;
138: for (PetscInt fc = 0; fc < NcI; ++fc) {
139: for (PetscInt gc = 0; gc < NcJ; ++gc) g0[fc * NcJ + gc] += u_x[fc * NcJ + gc];
140: }
141: }
143: /*
144: (psi_i, u_j grad_j u_i) ==> (\psi_i, \u_j grad_j \phi_i)
145: */
146: 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[])
147: {
148: PetscInt NcI = dim;
149: PetscInt NcJ = dim;
150: for (PetscInt fc = 0; fc < NcI; ++fc) {
151: for (PetscInt gc = 0; gc < NcJ; ++gc) {
152: for (PetscInt dg = 0; dg < dim; ++dg) {
153: /* kronecker delta */
154: if (fc == gc) g1[(fc * NcJ + gc) * dim + dg] += u[dg];
155: }
156: }
157: }
158: }
160: /* < q, \nabla\cdot u >
161: NcompI = 1, NcompJ = dim */
162: 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[])
163: {
164: for (PetscInt d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
165: }
167: /* -< \nabla\cdot v, p >
168: NcompI = dim, NcompJ = 1 */
169: 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[])
170: {
171: for (PetscInt d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
172: }
174: /* < \nabla v, \nabla u + {\nabla u}^T >
175: This just gives \nabla u, give the perdiagonal for the transpose */
176: 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[])
177: {
178: const PetscReal Re = REYN;
179: const PetscInt Ncomp = dim;
181: for (PetscInt compI = 0; compI < Ncomp; ++compI)
182: for (PetscInt d = 0; d < dim; ++d) g3[((compI * Ncomp + compI) * dim + d) * dim + d] = 1.0 / Re;
183: }
185: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
186: {
187: PetscFunctionBeginUser;
188: options->mms = 1;
190: PetscOptionsBegin(comm, "", "Navier-Stokes Equation Options", "DMPLEX");
191: PetscCall(PetscOptionsInt("-mms", "The manufactured solution to use", "ex46.c", options->mms, &options->mms, NULL));
192: PetscOptionsEnd();
193: PetscFunctionReturn(PETSC_SUCCESS);
194: }
196: static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
197: {
198: PetscFunctionBeginUser;
199: PetscCall(DMCreate(comm, dm));
200: PetscCall(DMSetType(*dm, DMPLEX));
201: PetscCall(DMSetFromOptions(*dm));
202: PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
203: PetscFunctionReturn(PETSC_SUCCESS);
204: }
206: static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
207: {
208: PetscDS ds;
209: DMLabel label;
210: const PetscInt id = 1;
211: PetscInt dim;
213: PetscFunctionBeginUser;
214: PetscCall(DMGetDimension(dm, &dim));
215: PetscCall(DMGetDS(dm, &ds));
216: PetscCall(DMGetLabel(dm, "marker", &label));
217: switch (dim) {
218: case 2:
219: switch (ctx->mms) {
220: case 1:
221: PetscCall(PetscDSSetResidual(ds, 0, f0_mms1_u, f1_u));
222: PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p));
223: PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu));
224: PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
225: PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL));
226: PetscCall(PetscDSSetExactSolution(ds, 0, mms1_u_2d, ctx));
227: PetscCall(PetscDSSetExactSolution(ds, 1, mms1_p_2d, ctx));
228: PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)mms1_u_2d, NULL, ctx, NULL));
229: break;
230: case 2:
231: PetscCall(PetscDSSetResidual(ds, 0, f0_mms2_u, f1_u));
232: PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p));
233: PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu));
234: PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
235: PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL));
236: PetscCall(PetscDSSetExactSolution(ds, 0, mms2_u_2d, ctx));
237: PetscCall(PetscDSSetExactSolution(ds, 1, mms2_p_2d, ctx));
238: PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)mms2_u_2d, NULL, ctx, NULL));
239: break;
240: default:
241: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid MMS %" PetscInt_FMT, ctx->mms);
242: }
243: break;
244: default:
245: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
246: }
247: PetscFunctionReturn(PETSC_SUCCESS);
248: }
250: static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
251: {
252: MPI_Comm comm;
253: DM cdm = dm;
254: PetscFE fe[2];
255: PetscInt dim;
256: PetscBool simplex;
258: PetscFunctionBeginUser;
259: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
260: PetscCall(DMGetDimension(dm, &dim));
261: PetscCall(DMPlexIsSimplex(dm, &simplex));
262: PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
263: PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity"));
264: PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
265: PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
266: PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure"));
267: /* Set discretization and boundary conditions for each mesh */
268: PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
269: PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
270: PetscCall(DMCreateDS(dm));
271: PetscCall(SetupProblem(dm, ctx));
272: while (cdm) {
273: PetscObject pressure;
274: MatNullSpace nsp;
276: PetscCall(DMGetField(cdm, 1, NULL, &pressure));
277: PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nsp));
278: PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nsp));
279: PetscCall(MatNullSpaceDestroy(&nsp));
281: PetscCall(DMCopyDisc(dm, cdm));
282: PetscCall(DMGetCoarseDM(cdm, &cdm));
283: }
284: PetscCall(PetscFEDestroy(&fe[0]));
285: PetscCall(PetscFEDestroy(&fe[1]));
286: PetscFunctionReturn(PETSC_SUCCESS);
287: }
289: static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, PetscCtx ctx)
290: {
291: PetscSimplePointFn *funcs[2];
292: void *ctxs[2];
293: DM dm;
294: PetscDS ds;
295: PetscReal ferrors[2];
297: PetscFunctionBeginUser;
298: PetscCall(TSGetDM(ts, &dm));
299: PetscCall(DMGetDS(dm, &ds));
300: PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
301: PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
302: PetscCall(DMComputeL2FieldDiff(dm, crtime, funcs, ctxs, u, ferrors));
303: 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]));
304: PetscFunctionReturn(PETSC_SUCCESS);
305: }
307: int main(int argc, char **argv)
308: {
309: AppCtx ctx;
310: DM dm;
311: TS ts;
312: Vec u, r;
314: PetscFunctionBeginUser;
315: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
316: PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
317: PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
318: PetscCall(DMSetApplicationContext(dm, &ctx));
319: PetscCall(SetupDiscretization(dm, &ctx));
320: PetscCall(DMPlexCreateClosureIndex(dm, NULL));
322: PetscCall(DMCreateGlobalVector(dm, &u));
323: PetscCall(VecDuplicate(u, &r));
325: PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
326: PetscCall(TSMonitorSet(ts, MonitorError, &ctx, NULL));
327: PetscCall(TSSetDM(ts, dm));
328: PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
329: PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
330: PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
331: PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
332: PetscCall(TSSetFromOptions(ts));
333: PetscCall(DMTSCheckFromOptions(ts, u));
335: {
336: PetscSimplePointFn *funcs[2];
337: void *ctxs[2];
338: PetscDS ds;
340: PetscCall(DMGetDS(dm, &ds));
341: PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
342: PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
343: PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, u));
344: }
345: PetscCall(TSSolve(ts, u));
346: PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
348: PetscCall(VecDestroy(&u));
349: PetscCall(VecDestroy(&r));
350: PetscCall(TSDestroy(&ts));
351: PetscCall(DMDestroy(&dm));
352: PetscCall(PetscFinalize());
353: return 0;
354: }
356: /*TEST
358: # Full solves
359: test:
360: suffix: 2d_p2p1_r1
361: requires: !single triangle
362: filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
363: args: -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
364: -ts_type beuler -ts_max_steps 10 -ts_time_step 0.1 -ts_monitor -dmts_check \
365: -snes_monitor_short -snes_converged_reason \
366: -ksp_monitor_short -ksp_converged_reason \
367: -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
368: -fieldsplit_velocity_pc_type lu \
369: -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi
371: test:
372: suffix: 2d_q2q1_r1
373: requires: !single
374: filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" -e "s~ 0\]~ 0.0\]~g"
375: args: -dm_plex_simplex 0 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
376: -ts_type beuler -ts_max_steps 10 -ts_time_step 0.1 -ts_monitor -dmts_check \
377: -snes_monitor_short -snes_converged_reason \
378: -ksp_monitor_short -ksp_converged_reason \
379: -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
380: -fieldsplit_velocity_pc_type lu \
381: -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi
383: TEST*/