1c4762a1bSJed Brown static char help[] = "Time dependent Navier-Stokes problem in 2d and 3d with finite elements.\n\ 2c4762a1bSJed Brown We solve the Navier-Stokes in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4c4762a1bSJed Brown This example supports discretized auxiliary fields (Re) as well as\n\ 5c4762a1bSJed Brown multilevel nonlinear solvers.\n\ 6c4762a1bSJed Brown Contributed by: Julian Andrej <juan@tf.uni-kiel.de>\n\n\n"; 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscdmplex.h> 9c4762a1bSJed Brown #include <petscsnes.h> 10c4762a1bSJed Brown #include <petscts.h> 11c4762a1bSJed Brown #include <petscds.h> 12c4762a1bSJed Brown 13c4762a1bSJed Brown /* 14c4762a1bSJed Brown Navier-Stokes equation: 15c4762a1bSJed Brown 16c4762a1bSJed Brown du/dt + u . grad u - \Delta u - grad p = f 17c4762a1bSJed Brown div u = 0 18c4762a1bSJed Brown */ 19c4762a1bSJed Brown 20c4762a1bSJed Brown typedef struct { 21c4762a1bSJed Brown PetscInt mms; 22c4762a1bSJed Brown } AppCtx; 23c4762a1bSJed Brown 24c4762a1bSJed Brown #define REYN 400.0 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* MMS1 27c4762a1bSJed Brown 28c4762a1bSJed Brown u = t + x^2 + y^2; 29c4762a1bSJed Brown v = t + 2*x^2 - 2*x*y; 30c4762a1bSJed Brown p = x + y - 1; 31c4762a1bSJed Brown 32c4762a1bSJed Brown f_x = -2*t*(x + y) + 2*x*y^2 - 4*x^2*y - 2*x^3 + 4.0/Re - 1.0 33c4762a1bSJed Brown f_y = -2*t*x + 2*y^3 - 4*x*y^2 - 2*x^2*y + 4.0/Re - 1.0 34c4762a1bSJed Brown 35c4762a1bSJed Brown so that 36c4762a1bSJed Brown 37c4762a1bSJed Brown 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> 38c4762a1bSJed Brown + <-t (2x + 2y) + 2xy^2 - 4x^2y - 2x^3 + 4/Re - 1, -2xt + 2y^3 - 4xy^2 - 2x^2y + 4/Re - 1> = 0 39c4762a1bSJed Brown \nabla \cdot u = 2x - 2x = 0 40c4762a1bSJed Brown 41c4762a1bSJed Brown where 42c4762a1bSJed Brown 43c4762a1bSJed Brown <u, v> . <<u_x, v_x>, <u_y, v_y>> = <u u_x + v u_y, u v_x + v v_y> 44c4762a1bSJed Brown */ 459371c9d4SSatish Balay PetscErrorCode mms1_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 46c4762a1bSJed Brown u[0] = time + x[0] * x[0] + x[1] * x[1]; 47c4762a1bSJed Brown u[1] = time + 2.0 * x[0] * x[0] - 2.0 * x[0] * x[1]; 48c4762a1bSJed Brown return 0; 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 519371c9d4SSatish Balay PetscErrorCode mms1_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) { 52c4762a1bSJed Brown *p = x[0] + x[1] - 1.0; 53c4762a1bSJed Brown return 0; 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* MMS 2*/ 57c4762a1bSJed Brown 589371c9d4SSatish Balay static PetscErrorCode mms2_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 59c4762a1bSJed Brown u[0] = PetscSinReal(time + x[0]) * PetscSinReal(time + x[1]); 60c4762a1bSJed Brown u[1] = PetscCosReal(time + x[0]) * PetscCosReal(time + x[1]); 61c4762a1bSJed Brown return 0; 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 649371c9d4SSatish Balay static PetscErrorCode mms2_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) { 65c4762a1bSJed Brown *p = PetscSinReal(time + x[0] - x[1]); 66c4762a1bSJed Brown return 0; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 699371c9d4SSatish Balay 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[]) { 70c4762a1bSJed Brown const PetscReal Re = REYN; 71c4762a1bSJed Brown const PetscInt Ncomp = dim; 72c4762a1bSJed Brown PetscInt c, d; 73c4762a1bSJed Brown 74c4762a1bSJed Brown for (c = 0; c < Ncomp; ++c) { 75*ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d]; 76c4762a1bSJed Brown } 77c4762a1bSJed Brown f0[0] += u_t[0]; 78c4762a1bSJed Brown f0[1] += u_t[1]; 79c4762a1bSJed Brown 80c4762a1bSJed Brown 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; 81c4762a1bSJed Brown 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; 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 849371c9d4SSatish Balay 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[]) { 85c4762a1bSJed Brown const PetscReal Re = REYN; 86c4762a1bSJed Brown const PetscInt Ncomp = dim; 87c4762a1bSJed Brown PetscInt c, d; 88c4762a1bSJed Brown 89c4762a1bSJed Brown for (c = 0; c < Ncomp; ++c) { 90*ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d]; 91c4762a1bSJed Brown } 92c4762a1bSJed Brown f0[0] += u_t[0]; 93c4762a1bSJed Brown f0[1] += u_t[1]; 94c4762a1bSJed Brown 95c4762a1bSJed Brown 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; 96c4762a1bSJed Brown 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; 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 999371c9d4SSatish Balay 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[]) { 100c4762a1bSJed Brown const PetscReal Re = REYN; 101c4762a1bSJed Brown const PetscInt Ncomp = dim; 102c4762a1bSJed Brown PetscInt comp, d; 103c4762a1bSJed Brown 104c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 105*ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f1[comp * dim + d] = 1.0 / Re * u_x[comp * dim + d]; 106c4762a1bSJed Brown f1[comp * dim + comp] -= u[Ncomp]; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 1109371c9d4SSatish Balay 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[]) { 111c4762a1bSJed Brown PetscInt d; 112c4762a1bSJed Brown for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d]; 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 1159371c9d4SSatish Balay 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[]) { 116c4762a1bSJed Brown PetscInt d; 117c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = 0.0; 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* 121c4762a1bSJed Brown (psi_i, u_j grad_j u_i) ==> (\psi_i, \phi_j grad_j u_i) 122c4762a1bSJed Brown */ 1239371c9d4SSatish Balay 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[]) { 124c4762a1bSJed Brown PetscInt NcI = dim, NcJ = dim; 125c4762a1bSJed Brown PetscInt fc, gc; 126c4762a1bSJed Brown PetscInt d; 127c4762a1bSJed Brown 128*ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g0[d * dim + d] = u_tShift; 129c4762a1bSJed Brown 130c4762a1bSJed Brown for (fc = 0; fc < NcI; ++fc) { 131*ad540459SPierre Jolivet for (gc = 0; gc < NcJ; ++gc) g0[fc * NcJ + gc] += u_x[fc * NcJ + gc]; 132c4762a1bSJed Brown } 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* 136c4762a1bSJed Brown (psi_i, u_j grad_j u_i) ==> (\psi_i, \u_j grad_j \phi_i) 137c4762a1bSJed Brown */ 1389371c9d4SSatish Balay 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[]) { 139c4762a1bSJed Brown PetscInt NcI = dim; 140c4762a1bSJed Brown PetscInt NcJ = dim; 141c4762a1bSJed Brown PetscInt fc, gc, dg; 142c4762a1bSJed Brown for (fc = 0; fc < NcI; ++fc) { 143c4762a1bSJed Brown for (gc = 0; gc < NcJ; ++gc) { 144c4762a1bSJed Brown for (dg = 0; dg < dim; ++dg) { 145c4762a1bSJed Brown /* kronecker delta */ 146*ad540459SPierre Jolivet if (fc == gc) g1[(fc * NcJ + gc) * dim + dg] += u[dg]; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown } 149c4762a1bSJed Brown } 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* < q, \nabla\cdot u > 153c4762a1bSJed Brown NcompI = 1, NcompJ = dim */ 1549371c9d4SSatish Balay 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[]) { 155c4762a1bSJed Brown PetscInt d; 156c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */ 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* -< \nabla\cdot v, p > 160c4762a1bSJed Brown NcompI = dim, NcompJ = 1 */ 1619371c9d4SSatish Balay 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[]) { 162c4762a1bSJed Brown PetscInt d; 163c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */ 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T > 167c4762a1bSJed Brown This just gives \nabla u, give the perdiagonal for the transpose */ 1689371c9d4SSatish Balay 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[]) { 169c4762a1bSJed Brown const PetscReal Re = REYN; 170c4762a1bSJed Brown const PetscInt Ncomp = dim; 171c4762a1bSJed Brown PetscInt compI, d; 172c4762a1bSJed Brown 173c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 174*ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g3[((compI * Ncomp + compI) * dim + d) * dim + d] = 1.0 / Re; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 1789371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 179c4762a1bSJed Brown PetscFunctionBeginUser; 180c4762a1bSJed Brown options->mms = 1; 181c4762a1bSJed Brown 182d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Navier-Stokes Equation Options", "DMPLEX"); 1839566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mms", "The manufactured solution to use", "ex46.c", options->mms, &options->mms, NULL)); 184d0609cedSBarry Smith PetscOptionsEnd(); 185c4762a1bSJed Brown PetscFunctionReturn(0); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 1889371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) { 189c4762a1bSJed Brown PetscFunctionBeginUser; 1909566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1919566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1929566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1939566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 194c4762a1bSJed Brown PetscFunctionReturn(0); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 1979371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) { 19845480ffeSMatthew G. Knepley PetscDS ds; 19945480ffeSMatthew G. Knepley DMLabel label; 200c4762a1bSJed Brown const PetscInt id = 1; 20130602db0SMatthew G. Knepley PetscInt dim; 202c4762a1bSJed Brown 203c4762a1bSJed Brown PetscFunctionBeginUser; 2049566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2059566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2069566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 20730602db0SMatthew G. Knepley switch (dim) { 20830602db0SMatthew G. Knepley case 2: 209c4762a1bSJed Brown switch (ctx->mms) { 210c4762a1bSJed Brown case 1: 2119566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_mms1_u, f1_u)); 2129566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p)); 2139566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu)); 2149566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL)); 2159566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL)); 2169566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms1_u_2d, ctx)); 2179566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, mms1_p_2d, ctx)); 2189566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))mms1_u_2d, NULL, ctx, NULL)); 219c4762a1bSJed Brown break; 220c4762a1bSJed Brown case 2: 2219566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_mms2_u, f1_u)); 2229566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p)); 2239566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu)); 2249566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL)); 2259566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL)); 2269566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms2_u_2d, ctx)); 2279566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, mms2_p_2d, ctx)); 2289566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))mms2_u_2d, NULL, ctx, NULL)); 229c4762a1bSJed Brown break; 23063a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid MMS %" PetscInt_FMT, ctx->mms); 231c4762a1bSJed Brown } 232c4762a1bSJed Brown break; 23363a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown PetscFunctionReturn(0); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 2389371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx) { 239c4762a1bSJed Brown MPI_Comm comm; 24030602db0SMatthew G. Knepley DM cdm = dm; 24130602db0SMatthew G. Knepley PetscFE fe[2]; 24230602db0SMatthew G. Knepley PetscInt dim; 24330602db0SMatthew G. Knepley PetscBool simplex; 244c4762a1bSJed Brown 245c4762a1bSJed Brown PetscFunctionBeginUser; 2469566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2479566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2489566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2499566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 2509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity")); 2519566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 2529566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure")); 254c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2559566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 2569566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 2579566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2589566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, ctx)); 259c4762a1bSJed Brown while (cdm) { 260c4762a1bSJed Brown PetscObject pressure; 261c4762a1bSJed Brown MatNullSpace nsp; 262c4762a1bSJed Brown 2639566063dSJacob Faibussowitsch PetscCall(DMGetField(cdm, 1, NULL, &pressure)); 2649566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nsp)); 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nsp)); 2669566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp)); 267c4762a1bSJed Brown 2689566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 2699566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 270c4762a1bSJed Brown } 2719566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 2729566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 273c4762a1bSJed Brown PetscFunctionReturn(0); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 2769371c9d4SSatish Balay static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) { 27730602db0SMatthew G. Knepley PetscSimplePointFunc funcs[2]; 27830602db0SMatthew G. Knepley void *ctxs[2]; 279c4762a1bSJed Brown DM dm; 28030602db0SMatthew G. Knepley PetscDS ds; 281c4762a1bSJed Brown PetscReal ferrors[2]; 282c4762a1bSJed Brown 283c4762a1bSJed Brown PetscFunctionBeginUser; 2849566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 2859566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2869566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0])); 2879566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1])); 2889566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, crtime, funcs, ctxs, u, ferrors)); 2899566063dSJacob Faibussowitsch 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])); 290c4762a1bSJed Brown PetscFunctionReturn(0); 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 2939371c9d4SSatish Balay int main(int argc, char **argv) { 294c4762a1bSJed Brown AppCtx ctx; 295c4762a1bSJed Brown DM dm; 296c4762a1bSJed Brown TS ts; 297c4762a1bSJed Brown Vec u, r; 298c4762a1bSJed Brown 299327415f7SBarry Smith PetscFunctionBeginUser; 3009566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3019566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3029566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx)); 3039566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &ctx)); 3049566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &ctx)); 3059566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 306c4762a1bSJed Brown 3079566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 309c4762a1bSJed Brown 3109566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 3119566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MonitorError, &ctx, NULL)); 3129566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 3139566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx)); 3149566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx)); 3159566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx)); 3169566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 3179566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 3189566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 319c4762a1bSJed Brown 32030602db0SMatthew G. Knepley { 32130602db0SMatthew G. Knepley PetscSimplePointFunc funcs[2]; 32230602db0SMatthew G. Knepley void *ctxs[2]; 32330602db0SMatthew G. Knepley PetscDS ds; 32430602db0SMatthew G. Knepley 3259566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3269566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0])); 3279566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1])); 3289566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, u)); 32930602db0SMatthew G. Knepley } 3309566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 3319566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 332c4762a1bSJed Brown 3339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 3359566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 3369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3379566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 338b122ec5aSJacob Faibussowitsch return 0; 339c4762a1bSJed Brown } 340c4762a1bSJed Brown 341c4762a1bSJed Brown /*TEST 342c4762a1bSJed Brown 343c4762a1bSJed Brown # Full solves 344c4762a1bSJed Brown test: 345c4762a1bSJed Brown suffix: 2d_p2p1_r1 346c4762a1bSJed Brown requires: !single triangle 347c4762a1bSJed Brown filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 34830602db0SMatthew G. Knepley args: -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 34930602db0SMatthew G. Knepley -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \ 35030602db0SMatthew G. Knepley -snes_monitor_short -snes_converged_reason \ 35130602db0SMatthew G. Knepley -ksp_monitor_short -ksp_converged_reason \ 35230602db0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \ 35330602db0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu \ 35430602db0SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi 35530602db0SMatthew G. Knepley 356c4762a1bSJed Brown test: 357c4762a1bSJed Brown suffix: 2d_q2q1_r1 358c4762a1bSJed Brown requires: !single 359c4762a1bSJed Brown filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" -e "s~ 0\]~ 0.0\]~g" 36030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 36130602db0SMatthew G. Knepley -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \ 36230602db0SMatthew G. Knepley -snes_monitor_short -snes_converged_reason \ 36330602db0SMatthew G. Knepley -ksp_monitor_short -ksp_converged_reason \ 36430602db0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \ 36530602db0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu \ 36630602db0SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi 367c4762a1bSJed Brown 368c4762a1bSJed Brown TEST*/ 369