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 */ 45c4762a1bSJed Brown PetscErrorCode mms1_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 46c4762a1bSJed Brown { 47c4762a1bSJed Brown u[0] = time + x[0]*x[0] + x[1]*x[1]; 48c4762a1bSJed Brown u[1] = time + 2.0*x[0]*x[0] - 2.0*x[0]*x[1]; 49c4762a1bSJed Brown return 0; 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52c4762a1bSJed Brown PetscErrorCode mms1_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) 53c4762a1bSJed Brown { 54c4762a1bSJed Brown *p = x[0] + x[1] - 1.0; 55c4762a1bSJed Brown return 0; 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* MMS 2*/ 59c4762a1bSJed Brown 60c4762a1bSJed Brown static PetscErrorCode mms2_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 61c4762a1bSJed Brown { 62c4762a1bSJed Brown u[0] = PetscSinReal(time + x[0])*PetscSinReal(time + x[1]); 63c4762a1bSJed Brown u[1] = PetscCosReal(time + x[0])*PetscCosReal(time + x[1]); 64c4762a1bSJed Brown return 0; 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown static PetscErrorCode mms2_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown *p = PetscSinReal(time + x[0] - x[1]); 70c4762a1bSJed Brown return 0; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown static void f0_mms1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 74c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 75c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 76c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 77c4762a1bSJed Brown { 78c4762a1bSJed Brown const PetscReal Re = REYN; 79c4762a1bSJed Brown const PetscInt Ncomp = dim; 80c4762a1bSJed Brown PetscInt c, d; 81c4762a1bSJed Brown 82c4762a1bSJed Brown for (c = 0; c < Ncomp; ++c) { 83c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 84c4762a1bSJed Brown f0[c] += u[d] * u_x[c*dim+d]; 85c4762a1bSJed Brown } 86c4762a1bSJed Brown } 87c4762a1bSJed Brown f0[0] += u_t[0]; 88c4762a1bSJed Brown f0[1] += u_t[1]; 89c4762a1bSJed Brown 90c4762a1bSJed 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; 91c4762a1bSJed 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; 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 94c4762a1bSJed Brown static void f0_mms2_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 95c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 96c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 97c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 98c4762a1bSJed Brown { 99c4762a1bSJed Brown const PetscReal Re = REYN; 100c4762a1bSJed Brown const PetscInt Ncomp = dim; 101c4762a1bSJed Brown PetscInt c, d; 102c4762a1bSJed Brown 103c4762a1bSJed Brown for (c = 0; c < Ncomp; ++c) { 104c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 105c4762a1bSJed Brown f0[c] += u[d] * u_x[c*dim+d]; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown } 108c4762a1bSJed Brown f0[0] += u_t[0]; 109c4762a1bSJed Brown f0[1] += u_t[1]; 110c4762a1bSJed Brown 111c4762a1bSJed 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; 112c4762a1bSJed 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; 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 116c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 117c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 118c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 119c4762a1bSJed Brown { 120c4762a1bSJed Brown const PetscReal Re = REYN; 121c4762a1bSJed Brown const PetscInt Ncomp = dim; 122c4762a1bSJed Brown PetscInt comp, d; 123c4762a1bSJed Brown 124c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 125c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 126c4762a1bSJed Brown f1[comp*dim+d] = 1.0/Re * u_x[comp*dim+d]; 127c4762a1bSJed Brown } 128c4762a1bSJed Brown f1[comp*dim+comp] -= u[Ncomp]; 129c4762a1bSJed Brown } 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 133c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 134c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 135c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 136c4762a1bSJed Brown { 137c4762a1bSJed Brown PetscInt d; 138c4762a1bSJed Brown for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown static void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 142c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 143c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 144c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 145c4762a1bSJed Brown { 146c4762a1bSJed Brown PetscInt d; 147c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = 0.0; 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* 151c4762a1bSJed Brown (psi_i, u_j grad_j u_i) ==> (\psi_i, \phi_j grad_j u_i) 152c4762a1bSJed Brown */ 153c4762a1bSJed Brown static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 154c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 155c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 156c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 157c4762a1bSJed Brown { 158c4762a1bSJed Brown PetscInt NcI = dim, NcJ = dim; 159c4762a1bSJed Brown PetscInt fc, gc; 160c4762a1bSJed Brown PetscInt d; 161c4762a1bSJed Brown 162c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 163c4762a1bSJed Brown g0[d*dim+d] = u_tShift; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown for (fc = 0; fc < NcI; ++fc) { 167c4762a1bSJed Brown for (gc = 0; gc < NcJ; ++gc) { 168c4762a1bSJed Brown g0[fc*NcJ+gc] += u_x[fc*NcJ+gc]; 169c4762a1bSJed Brown } 170c4762a1bSJed Brown } 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* 174c4762a1bSJed Brown (psi_i, u_j grad_j u_i) ==> (\psi_i, \u_j grad_j \phi_i) 175c4762a1bSJed Brown */ 176c4762a1bSJed Brown static void g1_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 177c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 178c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 179c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 180c4762a1bSJed Brown { 181c4762a1bSJed Brown PetscInt NcI = dim; 182c4762a1bSJed Brown PetscInt NcJ = dim; 183c4762a1bSJed Brown PetscInt fc, gc, dg; 184c4762a1bSJed Brown for (fc = 0; fc < NcI; ++fc) { 185c4762a1bSJed Brown for (gc = 0; gc < NcJ; ++gc) { 186c4762a1bSJed Brown for (dg = 0; dg < dim; ++dg) { 187c4762a1bSJed Brown /* kronecker delta */ 188c4762a1bSJed Brown if (fc == gc) { 189c4762a1bSJed Brown g1[(fc*NcJ+gc)*dim+dg] += u[dg]; 190c4762a1bSJed Brown } 191c4762a1bSJed Brown } 192c4762a1bSJed Brown } 193c4762a1bSJed Brown } 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* < q, \nabla\cdot u > 197c4762a1bSJed Brown NcompI = 1, NcompJ = dim */ 198c4762a1bSJed Brown static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 199c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 200c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 201c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 202c4762a1bSJed Brown { 203c4762a1bSJed Brown PetscInt d; 204c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */ 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* -< \nabla\cdot v, p > 208c4762a1bSJed Brown NcompI = dim, NcompJ = 1 */ 209c4762a1bSJed Brown static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, 210c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 211c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 212c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 213c4762a1bSJed Brown { 214c4762a1bSJed Brown PetscInt d; 215c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */ 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T > 219c4762a1bSJed Brown This just gives \nabla u, give the perdiagonal for the transpose */ 220c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 221c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 222c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 223c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 224c4762a1bSJed Brown { 225c4762a1bSJed Brown const PetscReal Re = REYN; 226c4762a1bSJed Brown const PetscInt Ncomp = dim; 227c4762a1bSJed Brown PetscInt compI, d; 228c4762a1bSJed Brown 229c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 230c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 231c4762a1bSJed Brown g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0/Re; 232c4762a1bSJed Brown } 233c4762a1bSJed Brown } 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 237c4762a1bSJed Brown { 238c4762a1bSJed Brown PetscFunctionBeginUser; 239c4762a1bSJed Brown options->mms = 1; 240c4762a1bSJed Brown 241*d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Navier-Stokes Equation Options", "DMPLEX"); 2429566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mms", "The manufactured solution to use", "ex46.c", options->mms, &options->mms, NULL)); 243*d0609cedSBarry Smith PetscOptionsEnd(); 244c4762a1bSJed Brown PetscFunctionReturn(0); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) 248c4762a1bSJed Brown { 249c4762a1bSJed Brown PetscFunctionBeginUser; 2509566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2519566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2529566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2539566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 254c4762a1bSJed Brown PetscFunctionReturn(0); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) 258c4762a1bSJed Brown { 25945480ffeSMatthew G. Knepley PetscDS ds; 26045480ffeSMatthew G. Knepley DMLabel label; 261c4762a1bSJed Brown const PetscInt id = 1; 26230602db0SMatthew G. Knepley PetscInt dim; 263c4762a1bSJed Brown 264c4762a1bSJed Brown PetscFunctionBeginUser; 2659566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2669566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2679566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 26830602db0SMatthew G. Knepley switch (dim) { 26930602db0SMatthew G. Knepley case 2: 270c4762a1bSJed Brown switch (ctx->mms) { 271c4762a1bSJed Brown case 1: 2729566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_mms1_u, f1_u)); 2739566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p)); 2749566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu)); 2759566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL)); 2769566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms1_u_2d, ctx)); 2789566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, mms1_p_2d, ctx)); 2799566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms1_u_2d, NULL, ctx, NULL)); 280c4762a1bSJed Brown break; 281c4762a1bSJed Brown case 2: 2829566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_mms2_u, f1_u)); 2839566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p)); 2849566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu)); 2859566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL)); 2869566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL)); 2879566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms2_u_2d, ctx)); 2889566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, mms2_p_2d, ctx)); 2899566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms2_u_2d, NULL, ctx, NULL)); 290c4762a1bSJed Brown break; 29198921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid MMS %D", ctx->mms); 292c4762a1bSJed Brown } 293c4762a1bSJed Brown break; 29498921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %D", dim); 295c4762a1bSJed Brown } 296c4762a1bSJed Brown PetscFunctionReturn(0); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx) 300c4762a1bSJed Brown { 301c4762a1bSJed Brown MPI_Comm comm; 30230602db0SMatthew G. Knepley DM cdm = dm; 30330602db0SMatthew G. Knepley PetscFE fe[2]; 30430602db0SMatthew G. Knepley PetscInt dim; 30530602db0SMatthew G. Knepley PetscBool simplex; 306c4762a1bSJed Brown 307c4762a1bSJed Brown PetscFunctionBeginUser; 3089566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 3099566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3109566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 3119566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 3129566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity")); 3139566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 3149566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 3159566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[1], "pressure")); 316c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 3179566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 3189566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 3199566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 3209566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, ctx)); 321c4762a1bSJed Brown while (cdm) { 322c4762a1bSJed Brown PetscObject pressure; 323c4762a1bSJed Brown MatNullSpace nsp; 324c4762a1bSJed Brown 3259566063dSJacob Faibussowitsch PetscCall(DMGetField(cdm, 1, NULL, &pressure)); 3269566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nsp)); 3279566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject) nsp)); 3289566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp)); 329c4762a1bSJed Brown 3309566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 3319566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 332c4762a1bSJed Brown } 3339566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 3349566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 335c4762a1bSJed Brown PetscFunctionReturn(0); 336c4762a1bSJed Brown } 337c4762a1bSJed Brown 338c4762a1bSJed Brown static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 339c4762a1bSJed Brown { 34030602db0SMatthew G. Knepley PetscSimplePointFunc funcs[2]; 34130602db0SMatthew G. Knepley void *ctxs[2]; 342c4762a1bSJed Brown DM dm; 34330602db0SMatthew G. Knepley PetscDS ds; 344c4762a1bSJed Brown PetscReal ferrors[2]; 345c4762a1bSJed Brown 346c4762a1bSJed Brown PetscFunctionBeginUser; 3479566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 3489566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3499566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0])); 3509566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1])); 3519566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, crtime, funcs, ctxs, u, ferrors)); 3529566063dSJacob 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])); 353c4762a1bSJed Brown PetscFunctionReturn(0); 354c4762a1bSJed Brown } 355c4762a1bSJed Brown 356c4762a1bSJed Brown int main(int argc, char **argv) 357c4762a1bSJed Brown { 358c4762a1bSJed Brown AppCtx ctx; 359c4762a1bSJed Brown DM dm; 360c4762a1bSJed Brown TS ts; 361c4762a1bSJed Brown Vec u, r; 362c4762a1bSJed Brown 3639566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3649566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 3659566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx)); 3669566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &ctx)); 3679566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &ctx)); 3689566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 369c4762a1bSJed Brown 3709566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 372c4762a1bSJed Brown 3739566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 3749566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MonitorError, &ctx, NULL)); 3759566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm)); 3769566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx)); 3779566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx)); 3789566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx)); 3799566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 3809566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 3819566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u)); 382c4762a1bSJed Brown 38330602db0SMatthew G. Knepley { 38430602db0SMatthew G. Knepley PetscSimplePointFunc funcs[2]; 38530602db0SMatthew G. Knepley void *ctxs[2]; 38630602db0SMatthew G. Knepley PetscDS ds; 38730602db0SMatthew G. Knepley 3889566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3899566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0])); 3909566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1])); 3919566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, u)); 39230602db0SMatthew G. Knepley } 3939566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 3949566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 395c4762a1bSJed Brown 3969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 3989566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 3999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4009566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 401b122ec5aSJacob Faibussowitsch return 0; 402c4762a1bSJed Brown } 403c4762a1bSJed Brown 404c4762a1bSJed Brown /*TEST 405c4762a1bSJed Brown 406c4762a1bSJed Brown # Full solves 407c4762a1bSJed Brown test: 408c4762a1bSJed Brown suffix: 2d_p2p1_r1 409c4762a1bSJed Brown requires: !single triangle 410c4762a1bSJed Brown filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" 41130602db0SMatthew G. Knepley args: -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 41230602db0SMatthew G. Knepley -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \ 41330602db0SMatthew G. Knepley -snes_monitor_short -snes_converged_reason \ 41430602db0SMatthew G. Knepley -ksp_monitor_short -ksp_converged_reason \ 41530602db0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \ 41630602db0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu \ 41730602db0SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi 41830602db0SMatthew G. Knepley 419c4762a1bSJed Brown test: 420c4762a1bSJed Brown suffix: 2d_q2q1_r1 421c4762a1bSJed Brown requires: !single 422c4762a1bSJed Brown filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" -e "s~ 0\]~ 0.0\]~g" 42330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 42430602db0SMatthew G. Knepley -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \ 42530602db0SMatthew G. Knepley -snes_monitor_short -snes_converged_reason \ 42630602db0SMatthew G. Knepley -ksp_monitor_short -ksp_converged_reason \ 42730602db0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \ 42830602db0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu \ 42930602db0SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi 430c4762a1bSJed Brown 431c4762a1bSJed Brown TEST*/ 432