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 PetscBag bag; /* Holds problem parameters */ 38c4762a1bSJed Brown } AppCtx; 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* 41c4762a1bSJed Brown In 2D, plane Poiseuille flow has exact solution: 42c4762a1bSJed Brown 43c4762a1bSJed Brown u = \Delta/(2 \nu) y (1 - y) + u_0 44c4762a1bSJed Brown v = 0 45c4762a1bSJed Brown p = -\Delta x 46c4762a1bSJed Brown f = 0 47c4762a1bSJed Brown 48c4762a1bSJed Brown so that 49c4762a1bSJed Brown 50c4762a1bSJed Brown -\nu \Delta u + \nabla p + f = <\Delta, 0> + <-\Delta, 0> + <0, 0> = 0 51c4762a1bSJed Brown \nabla \cdot u = 0 + 0 = 0 52c4762a1bSJed Brown 53c4762a1bSJed Brown In 3D we use exact solution: 54c4762a1bSJed Brown 55c4762a1bSJed Brown u = \Delta/(4 \nu) (y (1 - y) + z (1 - z)) + u_0 56c4762a1bSJed Brown v = 0 57c4762a1bSJed Brown w = 0 58c4762a1bSJed Brown p = -\Delta x 59c4762a1bSJed Brown f = 0 60c4762a1bSJed Brown 61c4762a1bSJed Brown so that 62c4762a1bSJed Brown 63c4762a1bSJed Brown -\nu \Delta u + \nabla p + f = <Delta, 0, 0> + <-Delta, 0, 0> + <0, 0, 0> = 0 64c4762a1bSJed Brown \nabla \cdot u = 0 + 0 + 0 = 0 65c4762a1bSJed Brown 66c4762a1bSJed Brown Note that these functions use coordinates X in the global (rotated) frame 67c4762a1bSJed Brown */ 68*2a8381b2SBarry Smith PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, PetscCtx ctx) 69d71ae5a4SJacob Faibussowitsch { 70c4762a1bSJed Brown Parameter *param = (Parameter *)ctx; 71c4762a1bSJed Brown PetscReal Delta = param->Delta; 72c4762a1bSJed Brown PetscReal nu = param->nu; 73c4762a1bSJed Brown PetscReal u_0 = param->u_0; 74c4762a1bSJed Brown PetscReal fac = (PetscReal)(dim - 1); 75c4762a1bSJed Brown PetscInt d; 76c4762a1bSJed Brown 77c4762a1bSJed Brown u[0] = u_0; 78c4762a1bSJed Brown for (d = 1; d < dim; ++d) u[0] += Delta / (fac * 2.0 * nu) * X[d] * (1.0 - X[d]); 79c4762a1bSJed Brown for (d = 1; d < dim; ++d) u[d] = 0.0; 803ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83*2a8381b2SBarry Smith PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, PetscCtx ctx) 84d71ae5a4SJacob Faibussowitsch { 85c4762a1bSJed Brown Parameter *param = (Parameter *)ctx; 86c4762a1bSJed Brown PetscReal Delta = param->Delta; 87c4762a1bSJed Brown 88c4762a1bSJed Brown p[0] = -Delta * X[0]; 893ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92*2a8381b2SBarry Smith PetscErrorCode wall_velocity(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, PetscCtx ctx) 93d71ae5a4SJacob Faibussowitsch { 94c4762a1bSJed Brown Parameter *param = (Parameter *)ctx; 95c4762a1bSJed Brown PetscReal u_0 = param->u_0; 96c4762a1bSJed Brown PetscInt d; 97c4762a1bSJed Brown 98c4762a1bSJed Brown u[0] = u_0; 99c4762a1bSJed Brown for (d = 1; d < dim; ++d) u[d] = 0.0; 1003ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed 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} 104c4762a1bSJed Brown u[Ncomp] = {p} */ 105d71ae5a4SJacob Faibussowitsch 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[]) 106d71ae5a4SJacob Faibussowitsch { 107c4762a1bSJed Brown const PetscReal nu = PetscRealPart(constants[1]); 108c4762a1bSJed Brown const PetscInt Nc = dim; 109c4762a1bSJed Brown PetscInt c, d; 110c4762a1bSJed Brown 111c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 112c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 113c4762a1bSJed Brown /* f1[c*dim+d] = 0.5*nu*(u_x[c*dim+d] + u_x[d*dim+c]); */ 114c4762a1bSJed Brown f1[c * dim + d] = nu * u_x[c * dim + d]; 115c4762a1bSJed Brown } 116c4762a1bSJed Brown f1[c * dim + c] -= u[uOff[1]]; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c4762a1bSJed 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} */ 121d71ae5a4SJacob Faibussowitsch 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[]) 122d71ae5a4SJacob Faibussowitsch { 123c4762a1bSJed Brown PetscInt d; 124c4762a1bSJed Brown for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d]; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* Residual functions are in reference coordinates */ 128d71ae5a4SJacob Faibussowitsch static void f0_bd_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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 129d71ae5a4SJacob Faibussowitsch { 130c4762a1bSJed Brown const PetscReal Delta = PetscRealPart(constants[0]); 131c4762a1bSJed Brown PetscReal alpha = PetscRealPart(constants[3]); 132c4762a1bSJed Brown PetscReal X = PetscCosReal(alpha) * x[0] + PetscSinReal(alpha) * x[1]; 133c4762a1bSJed Brown PetscInt d; 134c4762a1bSJed Brown 135ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f0[d] = -Delta * X * n[d]; 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* < q, \nabla\cdot u > 139c4762a1bSJed Brown NcompI = 1, NcompJ = dim */ 140d71ae5a4SJacob Faibussowitsch 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[]) 141d71ae5a4SJacob Faibussowitsch { 142c4762a1bSJed Brown PetscInt d; 143c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */ 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* -< \nabla\cdot v, p > 147c4762a1bSJed Brown NcompI = dim, NcompJ = 1 */ 148d71ae5a4SJacob Faibussowitsch 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[]) 149d71ae5a4SJacob Faibussowitsch { 150c4762a1bSJed Brown PetscInt d; 151c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */ 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T > 155c4762a1bSJed Brown This just gives \nabla u, give the perdiagonal for the transpose */ 156d71ae5a4SJacob Faibussowitsch 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[]) 157d71ae5a4SJacob Faibussowitsch { 158c4762a1bSJed Brown const PetscReal nu = PetscRealPart(constants[1]); 159c4762a1bSJed Brown const PetscInt Nc = dim; 160c4762a1bSJed Brown PetscInt c, d; 161c4762a1bSJed Brown 162c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 163ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g3[((c * Nc + c) * dim + d) * dim + d] = nu; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupParameters(AppCtx *user) 168d71ae5a4SJacob Faibussowitsch { 169c4762a1bSJed Brown PetscBag bag; 170c4762a1bSJed Brown Parameter *p; 171c4762a1bSJed Brown 172c4762a1bSJed Brown PetscFunctionBeginUser; 173c4762a1bSJed Brown /* setup PETSc parameter bag */ 174*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, &p)); 1759566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters")); 176c4762a1bSJed Brown bag = user->bag; 1779566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->Delta, 1.0, "Delta", "Pressure drop per unit length")); 1789566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity")); 1799566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->u_0, 0.0, "u_0", "Tangential velocity at the wall")); 1809566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->alpha, 0.0, "alpha", "Angle of pipe wall to x-axis")); 1813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 185d71ae5a4SJacob Faibussowitsch { 186c4762a1bSJed Brown PetscFunctionBeginUser; 1879566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1889566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1899566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 190c4762a1bSJed Brown { 191c4762a1bSJed Brown Parameter *param; 192c4762a1bSJed Brown Vec coordinates; 193c4762a1bSJed Brown PetscScalar *coords; 194c4762a1bSJed Brown PetscReal alpha; 195c4762a1bSJed Brown PetscInt cdim, N, bs, i; 196c4762a1bSJed Brown 1979566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*dm, &cdim)); 1989566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(*dm, &coordinates)); 1999566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &N)); 2009566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coordinates, &bs)); 20163a3b9bcSJacob Faibussowitsch PetscCheck(bs == cdim, comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %" PetscInt_FMT " != embedding dimension %" PetscInt_FMT, bs, cdim); 2029566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 203*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, ¶m)); 204c4762a1bSJed Brown alpha = param->alpha; 205c4762a1bSJed Brown for (i = 0; i < N; i += cdim) { 206c4762a1bSJed Brown PetscScalar x = coords[i + 0]; 207c4762a1bSJed Brown PetscScalar y = coords[i + 1]; 208c4762a1bSJed Brown 209c4762a1bSJed Brown coords[i + 0] = PetscCosReal(alpha) * x - PetscSinReal(alpha) * y; 210c4762a1bSJed Brown coords[i + 1] = PetscSinReal(alpha) * x + PetscCosReal(alpha) * y; 211c4762a1bSJed Brown } 2129566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 2139566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(*dm, coordinates)); 214c4762a1bSJed Brown } 2159566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 2163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupProblem(DM dm, AppCtx *user) 220d71ae5a4SJacob Faibussowitsch { 22145480ffeSMatthew G. Knepley PetscDS ds; 22245480ffeSMatthew G. Knepley PetscWeakForm wf; 22345480ffeSMatthew G. Knepley DMLabel label; 224c4762a1bSJed Brown Parameter *ctx; 22545480ffeSMatthew G. Knepley PetscInt id, bd; 226c4762a1bSJed Brown 227c4762a1bSJed Brown PetscFunctionBeginUser; 228*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, &ctx)); 2299566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2309566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u)); 2319566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, NULL)); 2329566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 2339566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL)); 2349566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL)); 23545480ffeSMatthew G. Knepley 23645480ffeSMatthew G. Knepley id = 2; 2379566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2389566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd)); 2399566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2409566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL)); 241c4762a1bSJed Brown /* Setup constants */ 242c4762a1bSJed Brown { 243c4762a1bSJed Brown Parameter *param; 244c4762a1bSJed Brown PetscScalar constants[4]; 245c4762a1bSJed Brown 246*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, ¶m)); 247c4762a1bSJed Brown 248c4762a1bSJed Brown constants[0] = param->Delta; 249c4762a1bSJed Brown constants[1] = param->nu; 250c4762a1bSJed Brown constants[2] = param->u_0; 251c4762a1bSJed Brown constants[3] = param->alpha; 2529566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 4, constants)); 253c4762a1bSJed Brown } 254c4762a1bSJed Brown /* Setup Boundary Conditions */ 255c4762a1bSJed Brown id = 3; 25657d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)wall_velocity, NULL, ctx, NULL)); 257c4762a1bSJed Brown id = 1; 25857d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)wall_velocity, NULL, ctx, NULL)); 259c4762a1bSJed Brown /* Setup exact solution */ 2609566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, ctx)); 2619566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, linear_p, ctx)); 2623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown 265d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 266d71ae5a4SJacob Faibussowitsch { 267c4762a1bSJed Brown DM cdm = dm; 268c4762a1bSJed Brown PetscFE fe[2]; 269c4762a1bSJed Brown Parameter *param; 27030602db0SMatthew G. Knepley PetscBool simplex; 27130602db0SMatthew G. Knepley PetscInt dim; 272c4762a1bSJed Brown MPI_Comm comm; 273c4762a1bSJed Brown 274c4762a1bSJed Brown PetscFunctionBeginUser; 2759566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2769566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2779566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 2799566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity")); 2809566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 2819566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 2829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure")); 283c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2849566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 2859566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 2869566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2879566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 288*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, ¶m)); 289c4762a1bSJed Brown while (cdm) { 2909566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 2919566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBasisRotation(cdm, param->alpha, 0.0, 0.0)); 2929566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 293c4762a1bSJed Brown } 2949566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 2959566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 300d71ae5a4SJacob Faibussowitsch { 301c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 302c4762a1bSJed Brown DM dm; /* problem definition */ 303c4762a1bSJed Brown Vec u, r; /* solution and residual */ 304c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 305c4762a1bSJed Brown 306327415f7SBarry Smith PetscFunctionBeginUser; 3079566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3089566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag)); 3099566063dSJacob Faibussowitsch PetscCall(SetupParameters(&user)); 3109566063dSJacob Faibussowitsch PetscCall(PetscBagSetFromOptions(user.bag)); 3119566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 3129566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 3139566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 3149566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user)); 315c4762a1bSJed Brown /* Setup problem */ 3169566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 3179566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 318c4762a1bSJed Brown 3199566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3209566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 321c4762a1bSJed Brown 3226493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user)); 323c4762a1bSJed Brown 3249566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 325c4762a1bSJed Brown 326c4762a1bSJed Brown { 32730602db0SMatthew G. Knepley PetscDS ds; 3288434afd1SBarry Smith PetscSimplePointFn *exactFuncs[2]; 329c4762a1bSJed Brown void *ctxs[2]; 330c4762a1bSJed Brown 3319566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3329566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 3339566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 3349566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u)); 3359566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution")); 3369566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-exact_vec_view")); 337c4762a1bSJed Brown } 3389566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 3399566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 3409566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Solution")); 3419566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 3429566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 343c4762a1bSJed Brown 3449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 3469566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3479566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 3489566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 3499566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 350b122ec5aSJacob Faibussowitsch return 0; 351c4762a1bSJed Brown } 352c4762a1bSJed Brown 353c4762a1bSJed Brown /*TEST 354c4762a1bSJed Brown 355c4762a1bSJed Brown # Convergence 356c4762a1bSJed Brown test: 357c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv 358c4762a1bSJed Brown requires: !single 35930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 \ 360c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 361c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \ 362c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 363c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 364c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 365c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 366c4762a1bSJed Brown test: 367c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv_u0 368c4762a1bSJed Brown requires: !single 36930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 \ 370c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 371c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \ 372c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 373c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 374c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 375c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 376c4762a1bSJed Brown test: 377c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv_u0_alpha 378c4762a1bSJed Brown requires: !single 37930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 -alpha 0.3927 \ 380c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 381c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \ 382c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 383c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 384c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 385c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 386c4762a1bSJed Brown test: 387c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv_gmg_vanka 388c4762a1bSJed Brown requires: !single long_runtime 38930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_plex_box_faces 2,2 -dm_refine_hierarchy 1 \ 390c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 391c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 1 -snes_error_if_not_converged \ 392c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 393c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 394c4762a1bSJed Brown -fieldsplit_velocity_pc_type mg \ 395c4762a1bSJed Brown -fieldsplit_velocity_mg_levels_pc_type patch -fieldsplit_velocity_mg_levels_pc_patch_exclude_subspaces 1 \ 396c4762a1bSJed Brown -fieldsplit_velocity_mg_levels_pc_patch_construct_codim 0 -fieldsplit_velocity_mg_levels_pc_patch_construct_type vanka \ 397c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-5 -fieldsplit_pressure_pc_type jacobi 398c4762a1bSJed Brown test: 399c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv 400c4762a1bSJed Brown requires: triangle !single 401c4762a1bSJed Brown args: -dm_plex_separate_marker -dm_refine 1 \ 402c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 403c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 404c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 405c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 406c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 407c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 408c4762a1bSJed Brown test: 409c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv_u0_alpha 410c4762a1bSJed Brown requires: triangle !single 411c4762a1bSJed Brown args: -dm_plex_separate_marker -dm_refine 0 -u_0 0.125 -alpha 0.3927 \ 412c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 413c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 414c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 415c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 416c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 417c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 418c4762a1bSJed Brown test: 419c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv_gmg_vcycle 42096b316e5SStefano Zampini TODO: broken (requires subDMs hooks) 421c4762a1bSJed Brown requires: triangle !single 42230602db0SMatthew G. Knepley args: -dm_plex_separate_marker -dm_plex_box_faces 2,2 -dm_refine_hierarchy 1 \ 423c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 424c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 425c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 426c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 427c4762a1bSJed Brown -fieldsplit_velocity_pc_type mg \ 428c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 429c4762a1bSJed Brown TEST*/ 430