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