1c4762a1bSJed Brown static char help[] = "Poiseuille Flow in 2d and 3d channels with finite elements.\n\ 2c4762a1bSJed Brown We solve the Poiseuille flow problem in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown /*F 6c4762a1bSJed Brown A Poiseuille flow is a steady-state isoviscous Stokes flow in a pipe of constant cross-section. We discretize using the 7c4762a1bSJed Brown finite element method on an unstructured mesh. The weak form equations are 8c4762a1bSJed Brown \begin{align*} 9c4762a1bSJed Brown < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p > + < v, \Delta \hat n >_{\Gamma_o} = 0 10c4762a1bSJed Brown < q, \nabla\cdot u > = 0 11c4762a1bSJed Brown \end{align*} 12c4762a1bSJed Brown where $\nu$ is the kinematic viscosity, $\Delta$ is the pressure drop per unit length, assuming that pressure is 0 on 13c4762a1bSJed Brown the left edge, and $\Gamma_o$ is the outlet boundary at the right edge of the pipe. The normal velocity will be zero at 14c4762a1bSJed Brown the wall, but we will allow a fixed tangential velocity $u_0$. 15c4762a1bSJed Brown 16c4762a1bSJed Brown In order to test our global to local basis transformation, we will allow the pipe to be at an angle $\alpha$ to the 17c4762a1bSJed Brown coordinate axes. 18c4762a1bSJed Brown 19c4762a1bSJed Brown For visualization, use 20c4762a1bSJed Brown 21c4762a1bSJed Brown -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append 22c4762a1bSJed Brown F*/ 23c4762a1bSJed Brown 24c4762a1bSJed Brown #include <petscdmplex.h> 25c4762a1bSJed Brown #include <petscsnes.h> 26c4762a1bSJed Brown #include <petscds.h> 27c4762a1bSJed Brown #include <petscbag.h> 28c4762a1bSJed Brown 29c4762a1bSJed Brown typedef struct { 30c4762a1bSJed Brown PetscReal Delta; /* Pressure drop per unit length */ 31c4762a1bSJed Brown PetscReal nu; /* Kinematic viscosity */ 32c4762a1bSJed Brown PetscReal u_0; /* Tangential velocity at the wall */ 33c4762a1bSJed Brown PetscReal alpha; /* Angle of pipe wall to x-axis */ 34c4762a1bSJed Brown } Parameter; 35c4762a1bSJed Brown 36c4762a1bSJed Brown typedef struct { 37c4762a1bSJed Brown /* Domain and mesh definition */ 38c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 39c4762a1bSJed Brown PetscBool simplex; /* Use simplices or tensor product cells */ 40c4762a1bSJed Brown PetscInt cells[3]; /* The initial domain division */ 41c4762a1bSJed Brown /* Problem definition */ 42c4762a1bSJed Brown PetscBag bag; /* Holds problem parameters */ 43c4762a1bSJed Brown PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 44c4762a1bSJed Brown } AppCtx; 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* 47c4762a1bSJed Brown In 2D, plane Poiseuille flow has exact solution: 48c4762a1bSJed Brown 49c4762a1bSJed Brown u = \Delta/(2 \nu) y (1 - y) + u_0 50c4762a1bSJed Brown v = 0 51c4762a1bSJed Brown p = -\Delta x 52c4762a1bSJed Brown f = 0 53c4762a1bSJed Brown 54c4762a1bSJed Brown so that 55c4762a1bSJed Brown 56c4762a1bSJed Brown -\nu \Delta u + \nabla p + f = <\Delta, 0> + <-\Delta, 0> + <0, 0> = 0 57c4762a1bSJed Brown \nabla \cdot u = 0 + 0 = 0 58c4762a1bSJed Brown 59c4762a1bSJed Brown In 3D we use exact solution: 60c4762a1bSJed Brown 61c4762a1bSJed Brown u = \Delta/(4 \nu) (y (1 - y) + z (1 - z)) + u_0 62c4762a1bSJed Brown v = 0 63c4762a1bSJed Brown w = 0 64c4762a1bSJed Brown p = -\Delta x 65c4762a1bSJed Brown f = 0 66c4762a1bSJed Brown 67c4762a1bSJed Brown so that 68c4762a1bSJed Brown 69c4762a1bSJed Brown -\nu \Delta u + \nabla p + f = <Delta, 0, 0> + <-Delta, 0, 0> + <0, 0, 0> = 0 70c4762a1bSJed Brown \nabla \cdot u = 0 + 0 + 0 = 0 71c4762a1bSJed Brown 72c4762a1bSJed Brown Note that these functions use coordinates X in the global (rotated) frame 73c4762a1bSJed Brown */ 74c4762a1bSJed Brown PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 75c4762a1bSJed Brown { 76c4762a1bSJed Brown Parameter *param = (Parameter *) ctx; 77c4762a1bSJed Brown PetscReal Delta = param->Delta; 78c4762a1bSJed Brown PetscReal nu = param->nu; 79c4762a1bSJed Brown PetscReal u_0 = param->u_0; 80c4762a1bSJed Brown PetscReal fac = (PetscReal) (dim - 1); 81c4762a1bSJed Brown PetscInt d; 82c4762a1bSJed Brown 83c4762a1bSJed Brown u[0] = u_0; 84c4762a1bSJed Brown for (d = 1; d < dim; ++d) u[0] += Delta/(fac * 2.0*nu) * X[d] * (1.0 - X[d]); 85c4762a1bSJed Brown for (d = 1; d < dim; ++d) u[d] = 0.0; 86c4762a1bSJed Brown return 0; 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89c4762a1bSJed Brown PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 90c4762a1bSJed Brown { 91c4762a1bSJed Brown Parameter *param = (Parameter *) ctx; 92c4762a1bSJed Brown PetscReal Delta = param->Delta; 93c4762a1bSJed Brown 94c4762a1bSJed Brown p[0] = -Delta * X[0]; 95c4762a1bSJed Brown return 0; 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown PetscErrorCode wall_velocity(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 99c4762a1bSJed Brown { 100c4762a1bSJed Brown Parameter *param = (Parameter *) ctx; 101c4762a1bSJed Brown PetscReal u_0 = param->u_0; 102c4762a1bSJed Brown PetscInt d; 103c4762a1bSJed Brown 104c4762a1bSJed Brown u[0] = u_0; 105c4762a1bSJed Brown for (d = 1; d < dim; ++d) u[d] = 0.0; 106c4762a1bSJed Brown return 0; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z} 110c4762a1bSJed Brown u[Ncomp] = {p} */ 111c4762a1bSJed Brown void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 112c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 113c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 114c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 115c4762a1bSJed Brown { 116c4762a1bSJed Brown const PetscReal nu = PetscRealPart(constants[1]); 117c4762a1bSJed Brown const PetscInt Nc = dim; 118c4762a1bSJed Brown PetscInt c, d; 119c4762a1bSJed Brown 120c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 121c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 122c4762a1bSJed Brown /* f1[c*dim+d] = 0.5*nu*(u_x[c*dim+d] + u_x[d*dim+c]); */ 123c4762a1bSJed Brown f1[c*dim+d] = nu*u_x[c*dim+d]; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown f1[c*dim+c] -= u[uOff[1]]; 126c4762a1bSJed Brown } 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z} */ 130c4762a1bSJed Brown void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 131c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 132c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 133c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 134c4762a1bSJed Brown { 135c4762a1bSJed Brown PetscInt d; 136c4762a1bSJed Brown for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Residual functions are in reference coordinates */ 140c4762a1bSJed Brown static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 141c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 142c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 143c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 144c4762a1bSJed Brown { 145c4762a1bSJed Brown const PetscReal Delta = PetscRealPart(constants[0]); 146c4762a1bSJed Brown PetscReal alpha = PetscRealPart(constants[3]); 147c4762a1bSJed Brown PetscReal X = PetscCosReal(alpha)*x[0] + PetscSinReal(alpha)*x[1]; 148c4762a1bSJed Brown PetscInt d; 149c4762a1bSJed Brown 150c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 151c4762a1bSJed Brown f0[d] = -Delta * X * n[d]; 152c4762a1bSJed Brown } 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* < q, \nabla\cdot u > 156c4762a1bSJed Brown NcompI = 1, NcompJ = dim */ 157c4762a1bSJed Brown void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 158c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 159c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 160c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 161c4762a1bSJed Brown { 162c4762a1bSJed Brown PetscInt d; 163c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */ 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* -< \nabla\cdot v, p > 167c4762a1bSJed Brown NcompI = dim, NcompJ = 1 */ 168c4762a1bSJed Brown void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, 169c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 170c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 171c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 172c4762a1bSJed Brown { 173c4762a1bSJed Brown PetscInt d; 174c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */ 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T > 178c4762a1bSJed Brown This just gives \nabla u, give the perdiagonal for the transpose */ 179c4762a1bSJed Brown void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 180c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 181c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 182c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 183c4762a1bSJed Brown { 184c4762a1bSJed Brown const PetscReal nu = PetscRealPart(constants[1]); 185c4762a1bSJed Brown const PetscInt Nc = dim; 186c4762a1bSJed Brown PetscInt c, d; 187c4762a1bSJed Brown 188c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 189c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 190c4762a1bSJed Brown g3[((c*Nc+c)*dim+d)*dim+d] = nu; 191c4762a1bSJed Brown } 192c4762a1bSJed Brown } 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 196c4762a1bSJed Brown { 197c4762a1bSJed Brown PetscInt n = 3; 198c4762a1bSJed Brown PetscErrorCode ierr; 199c4762a1bSJed Brown 200c4762a1bSJed Brown PetscFunctionBeginUser; 201c4762a1bSJed Brown options->dim = 2; 202c4762a1bSJed Brown options->simplex = PETSC_TRUE; 203c4762a1bSJed Brown options->cells[0] = 3; 204c4762a1bSJed Brown options->cells[1] = 3; 205c4762a1bSJed Brown options->cells[2] = 3; 206c4762a1bSJed Brown 207c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Poiseuille Flow Options", "DMPLEX");CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex62.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex62.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = PetscOptionsIntArray("-cells", "The initial mesh division", "ex62.c", options->cells, &n, NULL);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = PetscOptionsEnd(); 212c4762a1bSJed Brown PetscFunctionReturn(0); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215c4762a1bSJed Brown static PetscErrorCode SetupParameters(AppCtx *user) 216c4762a1bSJed Brown { 217c4762a1bSJed Brown PetscBag bag; 218c4762a1bSJed Brown Parameter *p; 219c4762a1bSJed Brown PetscErrorCode ierr; 220c4762a1bSJed Brown 221c4762a1bSJed Brown PetscFunctionBeginUser; 222c4762a1bSJed Brown /* setup PETSc parameter bag */ 223c4762a1bSJed Brown ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = PetscBagSetName(user->bag, "par", "Poiseuille flow parameters");CHKERRQ(ierr); 225c4762a1bSJed Brown bag = user->bag; 226c4762a1bSJed Brown ierr = PetscBagRegisterReal(bag, &p->Delta, 1.0, "Delta", "Pressure drop per unit length");CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity");CHKERRQ(ierr); 228c4762a1bSJed Brown ierr = PetscBagRegisterReal(bag, &p->u_0, 0.0, "u_0", "Tangential velocity at the wall");CHKERRQ(ierr); 229c4762a1bSJed Brown ierr = PetscBagRegisterReal(bag, &p->alpha, 0.0, "alpha", "Angle of pipe wall to x-axis");CHKERRQ(ierr); 230c4762a1bSJed Brown PetscFunctionReturn(0); 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 234c4762a1bSJed Brown { 235c4762a1bSJed Brown PetscInt dim = user->dim; 236c4762a1bSJed Brown PetscErrorCode ierr; 237c4762a1bSJed Brown 238c4762a1bSJed Brown PetscFunctionBeginUser; 239c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, user->simplex, user->cells, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 240c4762a1bSJed Brown { 241c4762a1bSJed Brown Parameter *param; 242c4762a1bSJed Brown Vec coordinates; 243c4762a1bSJed Brown PetscScalar *coords; 244c4762a1bSJed Brown PetscReal alpha; 245c4762a1bSJed Brown PetscInt cdim, N, bs, i; 246c4762a1bSJed Brown 247c4762a1bSJed Brown ierr = DMGetCoordinateDim(*dm, &cdim);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = DMGetCoordinates(*dm, &coordinates);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = VecGetBlockSize(coordinates, &bs);CHKERRQ(ierr); 251c4762a1bSJed Brown if (bs != cdim) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %D != embedding dimension %D", bs, cdim); 252c4762a1bSJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 253c4762a1bSJed Brown ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 254c4762a1bSJed Brown alpha = param->alpha; 255c4762a1bSJed Brown for (i = 0; i < N; i += cdim) { 256c4762a1bSJed Brown PetscScalar x = coords[i+0]; 257c4762a1bSJed Brown PetscScalar y = coords[i+1]; 258c4762a1bSJed Brown 259c4762a1bSJed Brown coords[i+0] = PetscCosReal(alpha)*x - PetscSinReal(alpha)*y; 260c4762a1bSJed Brown coords[i+1] = PetscSinReal(alpha)*x + PetscCosReal(alpha)*y; 261c4762a1bSJed Brown } 262c4762a1bSJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 263c4762a1bSJed Brown ierr = DMSetCoordinates(*dm, coordinates);CHKERRQ(ierr); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown { 266c4762a1bSJed Brown DM pdm = NULL; 267c4762a1bSJed Brown PetscPartitioner part; 268c4762a1bSJed Brown 269c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr); 272c4762a1bSJed Brown if (pdm) { 273c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 274c4762a1bSJed Brown *dm = pdm; 275c4762a1bSJed Brown } 276c4762a1bSJed Brown } 277c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 278c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 279c4762a1bSJed Brown PetscFunctionReturn(0); 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 282c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, AppCtx *user) 283c4762a1bSJed Brown { 284*45480ffeSMatthew G. Knepley PetscDS ds; 285*45480ffeSMatthew G. Knepley PetscWeakForm wf; 286*45480ffeSMatthew G. Knepley DMLabel label; 287c4762a1bSJed Brown Parameter *ctx; 288*45480ffeSMatthew G. Knepley PetscInt id, bd; 289c4762a1bSJed Brown PetscErrorCode ierr; 290c4762a1bSJed Brown 291c4762a1bSJed Brown PetscFunctionBeginUser; 292*45480ffeSMatthew G. Knepley ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 293*45480ffeSMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 294*45480ffeSMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, NULL, f1_u);CHKERRQ(ierr); 295*45480ffeSMatthew G. Knepley ierr = PetscDSSetResidual(ds, 1, f0_p, NULL);CHKERRQ(ierr); 296*45480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 297*45480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL);CHKERRQ(ierr); 298*45480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL);CHKERRQ(ierr); 299*45480ffeSMatthew G. Knepley 300*45480ffeSMatthew G. Knepley id = 2; 301*45480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 302*45480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd);CHKERRQ(ierr); 303*45480ffeSMatthew G. Knepley ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 304*45480ffeSMatthew G. Knepley ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, f0_bd_u, 0, NULL);CHKERRQ(ierr); 305c4762a1bSJed Brown /* Setup constants */ 306c4762a1bSJed Brown { 307c4762a1bSJed Brown Parameter *param; 308c4762a1bSJed Brown PetscScalar constants[4]; 309c4762a1bSJed Brown 310c4762a1bSJed Brown ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 311c4762a1bSJed Brown 312c4762a1bSJed Brown constants[0] = param->Delta; 313c4762a1bSJed Brown constants[1] = param->nu; 314c4762a1bSJed Brown constants[2] = param->u_0; 315c4762a1bSJed Brown constants[3] = param->alpha; 316*45480ffeSMatthew G. Knepley ierr = PetscDSSetConstants(ds, 4, constants);CHKERRQ(ierr); 317c4762a1bSJed Brown } 318c4762a1bSJed Brown /* Setup Boundary Conditions */ 319c4762a1bSJed Brown id = 3; 320*45480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) wall_velocity, NULL, ctx, NULL);CHKERRQ(ierr); 321c4762a1bSJed Brown id = 1; 322*45480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) wall_velocity, NULL, ctx, NULL);CHKERRQ(ierr); 323c4762a1bSJed Brown /* Setup exact solution */ 324c4762a1bSJed Brown user->exactFuncs[0] = quadratic_u; 325c4762a1bSJed Brown user->exactFuncs[1] = linear_p; 326*45480ffeSMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, user->exactFuncs[0], ctx);CHKERRQ(ierr); 327*45480ffeSMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 1, user->exactFuncs[1], ctx);CHKERRQ(ierr); 328c4762a1bSJed Brown PetscFunctionReturn(0); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 332c4762a1bSJed Brown { 333c4762a1bSJed Brown DM cdm = dm; 334c4762a1bSJed Brown const PetscInt dim = user->dim; 335c4762a1bSJed Brown PetscFE fe[2]; 336c4762a1bSJed Brown Parameter *param; 337c4762a1bSJed Brown MPI_Comm comm; 338c4762a1bSJed Brown PetscErrorCode ierr; 339c4762a1bSJed Brown 340c4762a1bSJed Brown PetscFunctionBeginUser; 341c4762a1bSJed Brown /* Create finite element */ 342c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 343c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, dim, user->simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 344c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 345c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, user->simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 346c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 347c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 348c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 349c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 350c4762a1bSJed Brown ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 351c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 352c4762a1bSJed Brown ierr = SetupProblem(dm, user);CHKERRQ(ierr); 353c4762a1bSJed Brown ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 354c4762a1bSJed Brown while (cdm) { 355c4762a1bSJed Brown ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 356c4762a1bSJed Brown ierr = DMPlexCreateBasisRotation(cdm, param->alpha, 0.0, 0.0);CHKERRQ(ierr); 357c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 358c4762a1bSJed Brown } 359c4762a1bSJed Brown ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 360c4762a1bSJed Brown ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 361c4762a1bSJed Brown PetscFunctionReturn(0); 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 364c4762a1bSJed Brown int main(int argc, char **argv) 365c4762a1bSJed Brown { 366c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 367c4762a1bSJed Brown DM dm; /* problem definition */ 368c4762a1bSJed Brown Vec u, r; /* solution and residual */ 369c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 370c4762a1bSJed Brown PetscErrorCode ierr; 371c4762a1bSJed Brown 372c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 373c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr); 375c4762a1bSJed Brown ierr = SetupParameters(&user);CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 377c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 378c4762a1bSJed Brown ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 380c4762a1bSJed Brown /* Setup problem */ 381c4762a1bSJed Brown ierr = PetscMalloc(2 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &user.exactFuncs);CHKERRQ(ierr); 382c4762a1bSJed Brown ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 383c4762a1bSJed Brown ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 384c4762a1bSJed Brown 385c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 386c4762a1bSJed Brown ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 387c4762a1bSJed Brown 388c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr); 389c4762a1bSJed Brown 390c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 391c4762a1bSJed Brown 392c4762a1bSJed Brown { 393c4762a1bSJed Brown Parameter *param; 394c4762a1bSJed Brown void *ctxs[2]; 395c4762a1bSJed Brown 396c4762a1bSJed Brown ierr = PetscBagGetData(user.bag, (void **) ¶m);CHKERRQ(ierr); 397c4762a1bSJed Brown ctxs[0] = ctxs[1] = param; 398c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, user.exactFuncs, ctxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 399c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "Exact Solution");CHKERRQ(ierr); 400c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-exact_vec_view");CHKERRQ(ierr); 401c4762a1bSJed Brown } 402348a1646SMatthew G. Knepley ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 403c4762a1bSJed Brown ierr = VecSet(u, 0.0);CHKERRQ(ierr); 404c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr); 405c4762a1bSJed Brown ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 406c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 407c4762a1bSJed Brown 408c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 409c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 410c4762a1bSJed Brown ierr = PetscFree(user.exactFuncs);CHKERRQ(ierr); 411c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 412c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 413c4762a1bSJed Brown ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr); 414c4762a1bSJed Brown ierr = PetscFinalize(); 415c4762a1bSJed Brown return ierr; 416c4762a1bSJed Brown } 417c4762a1bSJed Brown 418c4762a1bSJed Brown /*TEST 419c4762a1bSJed Brown 420c4762a1bSJed Brown # Convergence 421c4762a1bSJed Brown test: 422c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv 423c4762a1bSJed Brown requires: !single 424c4762a1bSJed Brown args: -simplex 0 -dm_plex_separate_marker -dm_refine 1 \ 425c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 426c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \ 427c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 428c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 429c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 430c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 431c4762a1bSJed Brown test: 432c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv_u0 433c4762a1bSJed Brown requires: !single 434c4762a1bSJed Brown args: -simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 \ 435c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 436c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \ 437c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 438c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 439c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 440c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 441c4762a1bSJed Brown test: 442c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv_u0_alpha 443c4762a1bSJed Brown requires: !single 444c4762a1bSJed Brown args: -simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 -alpha 0.3927 \ 445c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 446c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \ 447c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 448c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 449c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 450c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 451c4762a1bSJed Brown test: 452c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv_gmg_vanka 453c4762a1bSJed Brown requires: !single long_runtime 454c4762a1bSJed Brown args: -simplex 0 -dm_plex_separate_marker -cells 2,2 -dm_refine_hierarchy 1 \ 455c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 456c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 1 -snes_error_if_not_converged \ 457c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 458c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 459c4762a1bSJed Brown -fieldsplit_velocity_pc_type mg \ 460c4762a1bSJed Brown -fieldsplit_velocity_mg_levels_pc_type patch -fieldsplit_velocity_mg_levels_pc_patch_exclude_subspaces 1 \ 461c4762a1bSJed Brown -fieldsplit_velocity_mg_levels_pc_patch_construct_codim 0 -fieldsplit_velocity_mg_levels_pc_patch_construct_type vanka \ 462c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-5 -fieldsplit_pressure_pc_type jacobi 463c4762a1bSJed Brown test: 464c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv 465c4762a1bSJed Brown requires: triangle !single 466c4762a1bSJed Brown args: -dm_plex_separate_marker -dm_refine 1 \ 467c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 468c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 469c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 470c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 471c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 472c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 473c4762a1bSJed Brown test: 474c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv_u0_alpha 475c4762a1bSJed Brown requires: triangle !single 476c4762a1bSJed Brown args: -dm_plex_separate_marker -dm_refine 0 -u_0 0.125 -alpha 0.3927 \ 477c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 478c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 479c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 480c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 481c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 482c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 483c4762a1bSJed Brown test: 484c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv_gmg_vcycle 485c4762a1bSJed Brown requires: triangle !single 486c4762a1bSJed Brown args: -dm_plex_separate_marker -cells 2,2 -dm_refine_hierarchy 1 \ 487c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 488c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 489c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 490c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 491c4762a1bSJed Brown -fieldsplit_velocity_pc_type mg \ 492c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 493c4762a1bSJed Brown TEST*/ 494