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*9371c9d4SSatish Balay PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) { 69c4762a1bSJed Brown Parameter *param = (Parameter *)ctx; 70c4762a1bSJed Brown PetscReal Delta = param->Delta; 71c4762a1bSJed Brown PetscReal nu = param->nu; 72c4762a1bSJed Brown PetscReal u_0 = param->u_0; 73c4762a1bSJed Brown PetscReal fac = (PetscReal)(dim - 1); 74c4762a1bSJed Brown PetscInt d; 75c4762a1bSJed Brown 76c4762a1bSJed Brown u[0] = u_0; 77c4762a1bSJed Brown for (d = 1; d < dim; ++d) u[0] += Delta / (fac * 2.0 * nu) * X[d] * (1.0 - X[d]); 78c4762a1bSJed Brown for (d = 1; d < dim; ++d) u[d] = 0.0; 79c4762a1bSJed Brown return 0; 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82*9371c9d4SSatish Balay PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) { 83c4762a1bSJed Brown Parameter *param = (Parameter *)ctx; 84c4762a1bSJed Brown PetscReal Delta = param->Delta; 85c4762a1bSJed Brown 86c4762a1bSJed Brown p[0] = -Delta * X[0]; 87c4762a1bSJed Brown return 0; 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90*9371c9d4SSatish Balay PetscErrorCode wall_velocity(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) { 91c4762a1bSJed Brown Parameter *param = (Parameter *)ctx; 92c4762a1bSJed Brown PetscReal u_0 = param->u_0; 93c4762a1bSJed Brown PetscInt d; 94c4762a1bSJed Brown 95c4762a1bSJed Brown u[0] = u_0; 96c4762a1bSJed Brown for (d = 1; d < dim; ++d) u[d] = 0.0; 97c4762a1bSJed Brown return 0; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed 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} 101c4762a1bSJed Brown u[Ncomp] = {p} */ 102*9371c9d4SSatish Balay 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[]) { 103c4762a1bSJed Brown const PetscReal nu = PetscRealPart(constants[1]); 104c4762a1bSJed Brown const PetscInt Nc = dim; 105c4762a1bSJed Brown PetscInt c, d; 106c4762a1bSJed Brown 107c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 108c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 109c4762a1bSJed Brown /* f1[c*dim+d] = 0.5*nu*(u_x[c*dim+d] + u_x[d*dim+c]); */ 110c4762a1bSJed Brown f1[c * dim + d] = nu * u_x[c * dim + d]; 111c4762a1bSJed Brown } 112c4762a1bSJed Brown f1[c * dim + c] -= u[uOff[1]]; 113c4762a1bSJed Brown } 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 116c4762a1bSJed 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} */ 117*9371c9d4SSatish Balay 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[]) { 118c4762a1bSJed Brown PetscInt d; 119c4762a1bSJed Brown for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d]; 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* Residual functions are in reference coordinates */ 123*9371c9d4SSatish Balay 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[]) { 124c4762a1bSJed Brown const PetscReal Delta = PetscRealPart(constants[0]); 125c4762a1bSJed Brown PetscReal alpha = PetscRealPart(constants[3]); 126c4762a1bSJed Brown PetscReal X = PetscCosReal(alpha) * x[0] + PetscSinReal(alpha) * x[1]; 127c4762a1bSJed Brown PetscInt d; 128c4762a1bSJed Brown 129*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { f0[d] = -Delta * X * n[d]; } 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* < q, \nabla\cdot u > 133c4762a1bSJed Brown NcompI = 1, NcompJ = dim */ 134*9371c9d4SSatish Balay 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[]) { 135c4762a1bSJed Brown PetscInt d; 136c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */ 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* -< \nabla\cdot v, p > 140c4762a1bSJed Brown NcompI = dim, NcompJ = 1 */ 141*9371c9d4SSatish Balay 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[]) { 142c4762a1bSJed Brown PetscInt d; 143c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */ 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T > 147c4762a1bSJed Brown This just gives \nabla u, give the perdiagonal for the transpose */ 148*9371c9d4SSatish Balay 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[]) { 149c4762a1bSJed Brown const PetscReal nu = PetscRealPart(constants[1]); 150c4762a1bSJed Brown const PetscInt Nc = dim; 151c4762a1bSJed Brown PetscInt c, d; 152c4762a1bSJed Brown 153c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 154*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { g3[((c * Nc + c) * dim + d) * dim + d] = nu; } 155c4762a1bSJed Brown } 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158*9371c9d4SSatish Balay static PetscErrorCode SetupParameters(AppCtx *user) { 159c4762a1bSJed Brown PetscBag bag; 160c4762a1bSJed Brown Parameter *p; 161c4762a1bSJed Brown 162c4762a1bSJed Brown PetscFunctionBeginUser; 163c4762a1bSJed Brown /* setup PETSc parameter bag */ 1649566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)&p)); 1659566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters")); 166c4762a1bSJed Brown bag = user->bag; 1679566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->Delta, 1.0, "Delta", "Pressure drop per unit length")); 1689566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity")); 1699566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->u_0, 0.0, "u_0", "Tangential velocity at the wall")); 1709566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->alpha, 0.0, "alpha", "Angle of pipe wall to x-axis")); 171c4762a1bSJed Brown PetscFunctionReturn(0); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174*9371c9d4SSatish Balay PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 175c4762a1bSJed Brown PetscFunctionBeginUser; 1769566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1779566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1789566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 179c4762a1bSJed Brown { 180c4762a1bSJed Brown Parameter *param; 181c4762a1bSJed Brown Vec coordinates; 182c4762a1bSJed Brown PetscScalar *coords; 183c4762a1bSJed Brown PetscReal alpha; 184c4762a1bSJed Brown PetscInt cdim, N, bs, i; 185c4762a1bSJed Brown 1869566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*dm, &cdim)); 1879566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(*dm, &coordinates)); 1889566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &N)); 1899566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coordinates, &bs)); 19063a3b9bcSJacob Faibussowitsch PetscCheck(bs == cdim, comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %" PetscInt_FMT " != embedding dimension %" PetscInt_FMT, bs, cdim); 1919566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 1929566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 193c4762a1bSJed Brown alpha = param->alpha; 194c4762a1bSJed Brown for (i = 0; i < N; i += cdim) { 195c4762a1bSJed Brown PetscScalar x = coords[i + 0]; 196c4762a1bSJed Brown PetscScalar y = coords[i + 1]; 197c4762a1bSJed Brown 198c4762a1bSJed Brown coords[i + 0] = PetscCosReal(alpha) * x - PetscSinReal(alpha) * y; 199c4762a1bSJed Brown coords[i + 1] = PetscSinReal(alpha) * x + PetscCosReal(alpha) * y; 200c4762a1bSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 2029566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(*dm, coordinates)); 203c4762a1bSJed Brown } 2049566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 205c4762a1bSJed Brown PetscFunctionReturn(0); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208*9371c9d4SSatish Balay PetscErrorCode SetupProblem(DM dm, AppCtx *user) { 20945480ffeSMatthew G. Knepley PetscDS ds; 21045480ffeSMatthew G. Knepley PetscWeakForm wf; 21145480ffeSMatthew G. Knepley DMLabel label; 212c4762a1bSJed Brown Parameter *ctx; 21345480ffeSMatthew G. Knepley PetscInt id, bd; 214c4762a1bSJed Brown 215c4762a1bSJed Brown PetscFunctionBeginUser; 2169566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)&ctx)); 2179566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2189566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u)); 2199566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, NULL)); 2209566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 2219566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL)); 22345480ffeSMatthew G. Knepley 22445480ffeSMatthew G. Knepley id = 2; 2259566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2269566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd)); 2279566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2289566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL)); 229c4762a1bSJed Brown /* Setup constants */ 230c4762a1bSJed Brown { 231c4762a1bSJed Brown Parameter *param; 232c4762a1bSJed Brown PetscScalar constants[4]; 233c4762a1bSJed Brown 2349566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown constants[0] = param->Delta; 237c4762a1bSJed Brown constants[1] = param->nu; 238c4762a1bSJed Brown constants[2] = param->u_0; 239c4762a1bSJed Brown constants[3] = param->alpha; 2409566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 4, constants)); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown /* Setup Boundary Conditions */ 243c4762a1bSJed Brown id = 3; 2449566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall", label, 1, &id, 0, 0, NULL, (void (*)(void))wall_velocity, NULL, ctx, NULL)); 245c4762a1bSJed Brown id = 1; 2469566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall", label, 1, &id, 0, 0, NULL, (void (*)(void))wall_velocity, NULL, ctx, NULL)); 247c4762a1bSJed Brown /* Setup exact solution */ 2489566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, ctx)); 2499566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, linear_p, ctx)); 250c4762a1bSJed Brown PetscFunctionReturn(0); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253*9371c9d4SSatish Balay PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) { 254c4762a1bSJed Brown DM cdm = dm; 255c4762a1bSJed Brown PetscFE fe[2]; 256c4762a1bSJed Brown Parameter *param; 25730602db0SMatthew G. Knepley PetscBool simplex; 25830602db0SMatthew G. Knepley PetscInt dim; 259c4762a1bSJed Brown MPI_Comm comm; 260c4762a1bSJed Brown 261c4762a1bSJed Brown PetscFunctionBeginUser; 2629566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2639566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2659566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity")); 2679566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 2689566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure")); 270c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2719566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 2729566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 2739566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2749566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 2759566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 276c4762a1bSJed Brown while (cdm) { 2779566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 2789566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBasisRotation(cdm, param->alpha, 0.0, 0.0)); 2799566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 280c4762a1bSJed Brown } 2819566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 2829566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 283c4762a1bSJed Brown PetscFunctionReturn(0); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286*9371c9d4SSatish Balay int main(int argc, char **argv) { 287c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 288c4762a1bSJed Brown DM dm; /* problem definition */ 289c4762a1bSJed Brown Vec u, r; /* solution and residual */ 290c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 291c4762a1bSJed Brown 292327415f7SBarry Smith PetscFunctionBeginUser; 2939566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2949566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag)); 2959566063dSJacob Faibussowitsch PetscCall(SetupParameters(&user)); 2969566063dSJacob Faibussowitsch PetscCall(PetscBagSetFromOptions(user.bag)); 2979566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 2989566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2999566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 3009566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user)); 301c4762a1bSJed Brown /* Setup problem */ 3029566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 3039566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 304c4762a1bSJed Brown 3059566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 307c4762a1bSJed Brown 3089566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 309c4762a1bSJed Brown 3109566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 311c4762a1bSJed Brown 312c4762a1bSJed Brown { 31330602db0SMatthew G. Knepley PetscDS ds; 31430602db0SMatthew G. Knepley PetscSimplePointFunc exactFuncs[2]; 315c4762a1bSJed Brown void *ctxs[2]; 316c4762a1bSJed Brown 3179566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3189566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 3199566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 3209566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u)); 3219566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution")); 3229566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-exact_vec_view")); 323c4762a1bSJed Brown } 3249566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 3259566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Solution")); 3279566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 3289566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 329c4762a1bSJed Brown 3309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 3329566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3339566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 3349566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 3359566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 336b122ec5aSJacob Faibussowitsch return 0; 337c4762a1bSJed Brown } 338c4762a1bSJed Brown 339c4762a1bSJed Brown /*TEST 340c4762a1bSJed Brown 341c4762a1bSJed Brown # Convergence 342c4762a1bSJed Brown test: 343c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv 344c4762a1bSJed Brown requires: !single 34530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 \ 346c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 347c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \ 348c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 349c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 350c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 351c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 352c4762a1bSJed Brown test: 353c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv_u0 354c4762a1bSJed Brown requires: !single 35530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 \ 356c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 357c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \ 358c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 359c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 360c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 361c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 362c4762a1bSJed Brown test: 363c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv_u0_alpha 364c4762a1bSJed Brown requires: !single 36530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 -alpha 0.3927 \ 366c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 367c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \ 368c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 369c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 370c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 371c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 372c4762a1bSJed Brown test: 373c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv_gmg_vanka 374c4762a1bSJed Brown requires: !single long_runtime 37530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_plex_box_faces 2,2 -dm_refine_hierarchy 1 \ 376c4762a1bSJed Brown -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 377c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 1 -snes_error_if_not_converged \ 378c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 379c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 380c4762a1bSJed Brown -fieldsplit_velocity_pc_type mg \ 381c4762a1bSJed Brown -fieldsplit_velocity_mg_levels_pc_type patch -fieldsplit_velocity_mg_levels_pc_patch_exclude_subspaces 1 \ 382c4762a1bSJed Brown -fieldsplit_velocity_mg_levels_pc_patch_construct_codim 0 -fieldsplit_velocity_mg_levels_pc_patch_construct_type vanka \ 383c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-5 -fieldsplit_pressure_pc_type jacobi 384c4762a1bSJed Brown test: 385c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv 386c4762a1bSJed Brown requires: triangle !single 387c4762a1bSJed Brown args: -dm_plex_separate_marker -dm_refine 1 \ 388c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 389c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 390c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 391c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 392c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 393c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 394c4762a1bSJed Brown test: 395c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv_u0_alpha 396c4762a1bSJed Brown requires: triangle !single 397c4762a1bSJed Brown args: -dm_plex_separate_marker -dm_refine 0 -u_0 0.125 -alpha 0.3927 \ 398c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 399c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 400c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 401c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 402c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 403c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 404c4762a1bSJed Brown test: 405c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv_gmg_vcycle 406c4762a1bSJed Brown requires: triangle !single 40730602db0SMatthew G. Knepley args: -dm_plex_separate_marker -dm_plex_box_faces 2,2 -dm_refine_hierarchy 1 \ 408c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 409c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 410c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 411c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 412c4762a1bSJed Brown -fieldsplit_velocity_pc_type mg \ 413c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 414c4762a1bSJed Brown TEST*/ 415