1c4762a1bSJed Brown static char help[] = "Poisson Problem in mixed form with 2d and 3d with finite elements.\n\ 2c4762a1bSJed Brown We solve the Poisson problem in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4c4762a1bSJed Brown This example supports automatic convergence estimation\n\ 5c4762a1bSJed Brown and Hdiv elements.\n\n\n"; 6c4762a1bSJed Brown 769b95281SMatthew G. Knepley /* 869b95281SMatthew G. Knepley 969b95281SMatthew G. Knepley The mixed form of Poisson's equation, e.g. https://firedrakeproject.org/demos/poisson_mixed.py.html, is given 1069b95281SMatthew G. Knepley in the strong form by 1169b95281SMatthew G. Knepley \begin{align} 1269b95281SMatthew G. Knepley \vb{q} - \nabla u &= 0 \\ 1369b95281SMatthew G. Knepley \nabla \cdot \vb{q} &= f 1469b95281SMatthew G. Knepley \end{align} 1569b95281SMatthew G. Knepley where $u$ is the potential, as in the original problem, but we also discretize the gradient of potential $\vb{q}$, 1669b95281SMatthew G. Knepley or flux, directly. The weak form is then 1769b95281SMatthew G. Knepley \begin{align} 1869b95281SMatthew G. Knepley <t, \vb{q}> + <\nabla \cdot t, u> - <t_n, u>_\Gamma &= 0 \\ 1969b95281SMatthew G. Knepley <v, \nabla \cdot \vb{q}> &= <v, f> 2069b95281SMatthew G. Knepley \end{align} 2169b95281SMatthew G. Knepley For the original Poisson problem, the Dirichlet boundary forces an essential boundary condition on the potential space, 2269b95281SMatthew G. Knepley and the Neumann boundary gives a natural boundary condition in the weak form. In the mixed formulation, the Neumann 2369b95281SMatthew G. Knepley boundary gives an essential boundary condition on the flux space, $\vb{q} \cdot \vu{n} = h$, and the Dirichlet condition 2469b95281SMatthew G. Knepley becomes a natural condition in the weak form, <t_n, g>_\Gamma. 2569b95281SMatthew G. Knepley */ 2669b95281SMatthew G. Knepley 27c4762a1bSJed Brown #include <petscdmplex.h> 28c4762a1bSJed Brown #include <petscsnes.h> 29c4762a1bSJed Brown #include <petscds.h> 30c4762a1bSJed Brown #include <petscconvest.h> 31c4762a1bSJed Brown 329371c9d4SSatish Balay typedef enum { 339371c9d4SSatish Balay SOL_LINEAR, 349371c9d4SSatish Balay SOL_QUADRATIC, 359371c9d4SSatish Balay SOL_QUARTIC, 369371c9d4SSatish Balay SOL_QUARTIC_NEUMANN, 379371c9d4SSatish Balay SOL_UNKNOWN, 389371c9d4SSatish Balay NUM_SOLTYPE 399371c9d4SSatish Balay } SolType; 4069b95281SMatthew G. Knepley const char *SolTypeNames[NUM_SOLTYPE + 3] = {"linear", "quadratic", "quartic", "quartic_neumann", "unknown", "SolType", "SOL_", NULL}; 41c4762a1bSJed Brown 42c4762a1bSJed Brown typedef struct { 43c4762a1bSJed Brown SolType solType; /* The type of exact solution */ 44c4762a1bSJed Brown } AppCtx; 45c4762a1bSJed Brown 46d71ae5a4SJacob Faibussowitsch static void f0_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 f0[]) 47d71ae5a4SJacob Faibussowitsch { 4869b95281SMatthew G. Knepley PetscInt d; 4969b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0] + d * dim + d]; 5069b95281SMatthew G. Knepley } 5169b95281SMatthew G. Knepley 52c4762a1bSJed Brown /* 2D Dirichlet potential example 53c4762a1bSJed Brown 54c4762a1bSJed Brown u = x 55c4762a1bSJed Brown q = <1, 0> 56c4762a1bSJed Brown f = 0 57c4762a1bSJed Brown 58c4762a1bSJed Brown We will need a boundary integral of u over \Gamma. 59c4762a1bSJed Brown */ 60d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 61d71ae5a4SJacob Faibussowitsch { 62c4762a1bSJed Brown u[0] = x[0]; 633ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 64c4762a1bSJed Brown } 65d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 66d71ae5a4SJacob Faibussowitsch { 67c4762a1bSJed Brown PetscInt c; 68c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = c ? 0.0 : 1.0; 693ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 70c4762a1bSJed Brown } 71c4762a1bSJed Brown 72d71ae5a4SJacob Faibussowitsch static void f0_linear_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 f0[]) 73d71ae5a4SJacob Faibussowitsch { 7469b95281SMatthew G. Knepley f0[0] = 0.0; 7569b95281SMatthew G. Knepley f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0); 7669b95281SMatthew G. Knepley } 77d71ae5a4SJacob Faibussowitsch static void f0_bd_linear_q(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[]) 78d71ae5a4SJacob Faibussowitsch { 7969b95281SMatthew G. Knepley PetscScalar potential; 8069b95281SMatthew G. Knepley PetscInt d; 8169b95281SMatthew G. Knepley 823ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, linear_u(dim, t, x, dim, &potential, NULL)); 8369b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[d] = -potential * n[d]; 8469b95281SMatthew G. Knepley } 8569b95281SMatthew G. Knepley 86c4762a1bSJed Brown /* 2D Dirichlet potential example 87c4762a1bSJed Brown 88c4762a1bSJed Brown u = x^2 + y^2 89c4762a1bSJed Brown q = <2x, 2y> 90c4762a1bSJed Brown f = 4 91c4762a1bSJed Brown 92c4762a1bSJed Brown We will need a boundary integral of u over \Gamma. 93c4762a1bSJed Brown */ 94d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 95d71ae5a4SJacob Faibussowitsch { 96c4762a1bSJed Brown PetscInt d; 97c4762a1bSJed Brown 98c4762a1bSJed Brown u[0] = 0.0; 99c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] += x[d] * x[d]; 1003ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 101c4762a1bSJed Brown } 102d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 103d71ae5a4SJacob Faibussowitsch { 104c4762a1bSJed Brown PetscInt c; 105c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = 2.0 * x[c]; 1063ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109d71ae5a4SJacob Faibussowitsch static void f0_quadratic_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 f0[]) 110d71ae5a4SJacob Faibussowitsch { 11169b95281SMatthew G. Knepley f0[0] = -4.0; 11269b95281SMatthew G. Knepley f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0); 11369b95281SMatthew G. Knepley } 114d71ae5a4SJacob Faibussowitsch static void f0_bd_quadratic_q(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[]) 115d71ae5a4SJacob Faibussowitsch { 11669b95281SMatthew G. Knepley PetscScalar potential; 11769b95281SMatthew G. Knepley PetscInt d; 11869b95281SMatthew G. Knepley 1193ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, quadratic_u(dim, t, x, dim, &potential, NULL)); 12069b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[d] = -potential * n[d]; 12169b95281SMatthew G. Knepley } 12269b95281SMatthew G. Knepley 123c4762a1bSJed Brown /* 2D Dirichlet potential example 124c4762a1bSJed Brown 125c4762a1bSJed Brown u = x (1-x) y (1-y) 126c4762a1bSJed Brown q = <(1-2x) y (1-y), x (1-x) (1-2y)> 127c4762a1bSJed Brown f = -y (1-y) - x (1-x) 128c4762a1bSJed Brown 129c4762a1bSJed Brown u|_\Gamma = 0 so that the boundary integral vanishes 130c4762a1bSJed Brown */ 131d71ae5a4SJacob Faibussowitsch static PetscErrorCode quartic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 132d71ae5a4SJacob Faibussowitsch { 133c4762a1bSJed Brown PetscInt d; 134c4762a1bSJed Brown 135c4762a1bSJed Brown u[0] = 1.0; 136c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] *= x[d] * (1.0 - x[d]); 1373ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 138c4762a1bSJed Brown } 139d71ae5a4SJacob Faibussowitsch static PetscErrorCode quartic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 140d71ae5a4SJacob Faibussowitsch { 141c4762a1bSJed Brown PetscInt c, d; 142c4762a1bSJed Brown 143c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 144c4762a1bSJed Brown u[c] = 1.0; 145c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 146c4762a1bSJed Brown if (c == d) u[c] *= 1 - 2.0 * x[d]; 147c4762a1bSJed Brown else u[c] *= x[d] * (1.0 - x[d]); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown } 1503ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153d71ae5a4SJacob Faibussowitsch static void f0_quartic_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 f0[]) 154d71ae5a4SJacob Faibussowitsch { 155c4762a1bSJed Brown PetscInt d; 156c4762a1bSJed Brown f0[0] = 0.0; 15769b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += 2.0 * x[d] * (1.0 - x[d]); 15869b95281SMatthew G. Knepley f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 16169b95281SMatthew G. Knepley /* 2D Dirichlet potential example 16269b95281SMatthew G. Knepley 16369b95281SMatthew G. Knepley u = x (1-x) y (1-y) 16469b95281SMatthew G. Knepley q = <(1-2x) y (1-y), x (1-x) (1-2y)> 16569b95281SMatthew G. Knepley f = -y (1-y) - x (1-x) 16669b95281SMatthew G. Knepley 16769b95281SMatthew G. Knepley du/dn_\Gamma = {(1-2x) y (1-y), x (1-x) (1-2y)} . n produces an essential condition on q 16869b95281SMatthew G. Knepley */ 16969b95281SMatthew G. Knepley 170d71ae5a4SJacob Faibussowitsch static void f0_q(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[]) 171d71ae5a4SJacob Faibussowitsch { 172c4762a1bSJed Brown PetscInt c; 17369b95281SMatthew G. Knepley for (c = 0; c < dim; ++c) f0[c] = u[uOff[0] + c]; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* <\nabla\cdot w, u> = <\nabla w, Iu> */ 177d71ae5a4SJacob Faibussowitsch static void f1_q(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[]) 178d71ae5a4SJacob Faibussowitsch { 17969b95281SMatthew G. Knepley PetscInt c; 18069b95281SMatthew G. Knepley for (c = 0; c < dim; ++c) f1[c * dim + c] = u[uOff[1]]; 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 18369b95281SMatthew G. Knepley /* <t, q> */ 184d71ae5a4SJacob Faibussowitsch static void g0_qq(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 g0[]) 185d71ae5a4SJacob Faibussowitsch { 186c4762a1bSJed Brown PetscInt c; 187c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 19069b95281SMatthew G. Knepley /* <\nabla\cdot t, u> = <\nabla t, Iu> */ 191d71ae5a4SJacob Faibussowitsch static void g2_qu(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[]) 192d71ae5a4SJacob Faibussowitsch { 193c4762a1bSJed Brown PetscInt d; 194c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d * dim + d] = 1.0; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 19769b95281SMatthew G. Knepley /* <v, \nabla\cdot q> */ 198d71ae5a4SJacob Faibussowitsch static void g1_uq(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[]) 199d71ae5a4SJacob Faibussowitsch { 200c4762a1bSJed Brown PetscInt d; 20169b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 204d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 205d71ae5a4SJacob Faibussowitsch { 206c4762a1bSJed Brown PetscFunctionBeginUser; 207c4762a1bSJed Brown options->solType = SOL_LINEAR; 208d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 2099566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sol_type", "Type of exact solution", "ex24.c", SolTypeNames, (PetscEnum)options->solType, (PetscEnum *)&options->solType, NULL)); 210d0609cedSBarry Smith PetscOptionsEnd(); 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 214d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 215d71ae5a4SJacob Faibussowitsch { 216c4762a1bSJed Brown PetscFunctionBeginUser; 2179566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2189566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Example Mesh")); 2209566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 2219566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2229566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 227d71ae5a4SJacob Faibussowitsch { 22845480ffeSMatthew G. Knepley PetscDS ds; 22945480ffeSMatthew G. Knepley DMLabel label; 23045480ffeSMatthew G. Knepley PetscWeakForm wf; 231c4762a1bSJed Brown const PetscInt id = 1; 23245480ffeSMatthew G. Knepley PetscInt bd; 233c4762a1bSJed Brown 234c4762a1bSJed Brown PetscFunctionBeginUser; 2359566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2369566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2379566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_q, f1_q)); 2389566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL)); 2399566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qu, NULL)); 2409566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_uq, NULL, NULL)); 2419371c9d4SSatish Balay switch (user->solType) { 242c4762a1bSJed Brown case SOL_LINEAR: 2439566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_linear_u, NULL)); 2449566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 2459566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2469566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_linear_q, 0, NULL)); 2479566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, linear_q, user)); 2489566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, linear_u, user)); 249c4762a1bSJed Brown break; 250c4762a1bSJed Brown case SOL_QUADRATIC: 2519566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_quadratic_u, NULL)); 2529566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 2539566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2549566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_quadratic_q, 0, NULL)); 2559566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_q, user)); 2569566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quadratic_u, user)); 257c4762a1bSJed Brown break; 258c4762a1bSJed Brown case SOL_QUARTIC: 2599566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL)); 2609566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user)); 2619566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user)); 262c4762a1bSJed Brown break; 26369b95281SMatthew G. Knepley case SOL_QUARTIC_NEUMANN: 2649566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user)); 2669566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user)); 2679566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "Flux condition", label, 1, &id, 0, 0, NULL, (void (*)(void))quartic_q, NULL, user, NULL)); 26869b95281SMatthew G. Knepley break; 269d71ae5a4SJacob Faibussowitsch default: 270d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]); 271c4762a1bSJed Brown } 2723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 276d71ae5a4SJacob Faibussowitsch { 277c4762a1bSJed Brown DM cdm = dm; 278c4762a1bSJed Brown PetscFE feq, feu; 27969b95281SMatthew G. Knepley DMPolytopeType ct; 28069b95281SMatthew G. Knepley PetscBool simplex; 28169b95281SMatthew G. Knepley PetscInt dim, cStart; 282c4762a1bSJed Brown 283c4762a1bSJed Brown PetscFunctionBeginUser; 2849566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2859566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 2869566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 28769b95281SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 288c4762a1bSJed Brown /* Create finite element */ 2899566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field_", -1, &feq)); 2909566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feq, "field")); 2919566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "potential_", -1, &feu)); 2929566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feu, "potential")); 2939566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(feq, feu)); 294c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2959566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)feq)); 2969566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)feu)); 2979566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2989566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 299c4762a1bSJed Brown while (cdm) { 3009566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 3019566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 302c4762a1bSJed Brown } 3039566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feq)); 3049566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feu)); 3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 309d71ae5a4SJacob Faibussowitsch { 310c4762a1bSJed Brown DM dm; /* Problem specification */ 311c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 312c4762a1bSJed Brown Vec u; /* Solutions */ 313c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 314c4762a1bSJed Brown 315327415f7SBarry Smith PetscFunctionBeginUser; 3169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3179566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 318c4762a1bSJed Brown /* Primal system */ 3199566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 3209566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 3219566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 3229566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, SetupPrimalProblem, &user)); 3239566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3249566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 3259566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "potential")); 3266493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user)); 3279566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 3289566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 3299566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 3309566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 3319566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 332c4762a1bSJed Brown /* Cleanup */ 3339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3349566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 3359566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3369566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 337b122ec5aSJacob Faibussowitsch return 0; 338c4762a1bSJed Brown } 339c4762a1bSJed Brown 340c4762a1bSJed Brown /*TEST 341c4762a1bSJed Brown 342c4762a1bSJed Brown test: 3436c2660fcSDarsh Nathawani suffix:2d_rt0_tri 3446c2660fcSDarsh Nathawani requires: triangle 3456c2660fcSDarsh Nathawani args: -sol_type linear -dmsnes_check 0.001 \ 3466c2660fcSDarsh Nathawani -potential_petscspace_degree 0 \ 3476c2660fcSDarsh Nathawani -potential_petscdualspace_lagrange_continuity 0 \ 3486c2660fcSDarsh Nathawani -field_petscspace_type ptrimmed \ 3496c2660fcSDarsh Nathawani -field_petscspace_components 2 \ 3506c2660fcSDarsh Nathawani -field_petscspace_ptrimmed_form_degree -1 \ 3516c2660fcSDarsh Nathawani -field_petscdualspace_order 1 \ 3526c2660fcSDarsh Nathawani -field_petscdualspace_form_degree -1 \ 3536c2660fcSDarsh Nathawani -field_petscdualspace_lagrange_trimmed true \ 3546c2660fcSDarsh Nathawani -field_petscfe_default_quadrature_order 2 \ 3556c2660fcSDarsh Nathawani -snes_error_if_not_converged \ 3566c2660fcSDarsh Nathawani -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 3576c2660fcSDarsh Nathawani -pc_type fieldsplit -pc_fieldsplit_type schur \ 3586c2660fcSDarsh Nathawani -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 3596c2660fcSDarsh Nathawani -fieldsplit_field_pc_type lu \ 3606c2660fcSDarsh Nathawani -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 3616c2660fcSDarsh Nathawani 3626c2660fcSDarsh Nathawani test: 3636c2660fcSDarsh Nathawani suffix:2d_rt0_quad 3646c2660fcSDarsh Nathawani requires: triangle 3656c2660fcSDarsh Nathawani args: -dm_plex_simplex 0 -sol_type linear -dmsnes_check 0.001 \ 3666c2660fcSDarsh Nathawani -potential_petscspace_degree 0 \ 3676c2660fcSDarsh Nathawani -potential_petscdualspace_lagrange_continuity 0 \ 3686c2660fcSDarsh Nathawani -field_petscspace_degree 1 \ 3696c2660fcSDarsh Nathawani -field_petscspace_type sum \ 3706c2660fcSDarsh Nathawani -field_petscspace_variables 2 \ 3716c2660fcSDarsh Nathawani -field_petscspace_components 2 \ 3726c2660fcSDarsh Nathawani -field_petscspace_sum_spaces 2 \ 3736c2660fcSDarsh Nathawani -field_petscspace_sum_concatenate true \ 3746c2660fcSDarsh Nathawani -field_sumcomp_0_petscspace_variables 2 \ 3756c2660fcSDarsh Nathawani -field_sumcomp_0_petscspace_type tensor \ 3766c2660fcSDarsh Nathawani -field_sumcomp_0_petscspace_tensor_spaces 2 \ 3776c2660fcSDarsh Nathawani -field_sumcomp_0_petscspace_tensor_uniform false \ 3786c2660fcSDarsh Nathawani -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \ 3796c2660fcSDarsh Nathawani -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \ 3806c2660fcSDarsh Nathawani -field_sumcomp_1_petscspace_variables 2 \ 3816c2660fcSDarsh Nathawani -field_sumcomp_1_petscspace_type tensor \ 3826c2660fcSDarsh Nathawani -field_sumcomp_1_petscspace_tensor_spaces 2 \ 3836c2660fcSDarsh Nathawani -field_sumcomp_1_petscspace_tensor_uniform false \ 3846c2660fcSDarsh Nathawani -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \ 3856c2660fcSDarsh Nathawani -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \ 3866c2660fcSDarsh Nathawani -field_petscdualspace_form_degree -1 \ 3876c2660fcSDarsh Nathawani -field_petscdualspace_order 1 \ 3886c2660fcSDarsh Nathawani -field_petscdualspace_lagrange_trimmed true \ 3896c2660fcSDarsh Nathawani -field_petscfe_default_quadrature_order 2 \ 3906c2660fcSDarsh Nathawani -pc_type fieldsplit -pc_fieldsplit_type schur \ 3916c2660fcSDarsh Nathawani -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 3926c2660fcSDarsh Nathawani -fieldsplit_field_pc_type lu \ 3936c2660fcSDarsh Nathawani -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 3946c2660fcSDarsh Nathawani 3956c2660fcSDarsh Nathawani test: 39669b95281SMatthew G. Knepley suffix: 2d_bdm1_p0 397c4762a1bSJed Brown requires: triangle 398c4762a1bSJed Brown args: -sol_type linear \ 39969b95281SMatthew G. Knepley -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 1 \ 400c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 401c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 402c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 403c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 404c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 405c4762a1bSJed Brown test: 40607a75e27SStefano Zampini nsize: 4 40707a75e27SStefano Zampini suffix: 2d_bdm1_p0_bddc 408910b42cbSStefano Zampini requires: triangle !single 40907a75e27SStefano Zampini args: -sol_type linear \ 41007a75e27SStefano Zampini -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 1 \ 41107a75e27SStefano Zampini -dmsnes_check .001 -snes_error_if_not_converged \ 412*79ab67a3SMatthew G. Knepley -ksp_error_if_not_converged -ksp_type cg -ksp_norm_type natural -ksp_divtol 1e10 \ 41307a75e27SStefano Zampini -petscpartitioner_type simple -dm_mat_type is \ 41407a75e27SStefano Zampini -pc_type bddc -pc_bddc_use_local_mat_graph 0 \ 41507a75e27SStefano Zampini -pc_bddc_benign_trick -pc_bddc_nonetflux -pc_bddc_detect_disconnected -pc_bddc_use_qr_single \ 41607a75e27SStefano Zampini -pc_bddc_coarse_redundant_pc_type svd -pc_bddc_neumann_pc_type svd -pc_bddc_dirichlet_pc_type svd 41707a75e27SStefano Zampini 41807a75e27SStefano Zampini test: 41907a75e27SStefano Zampini nsize: 9 420910b42cbSStefano Zampini requires: !single 42107a75e27SStefano Zampini suffix: 2d_rt1_p0_bddc 42207a75e27SStefano Zampini args: -sol_type quadratic \ 42307a75e27SStefano Zampini -potential_petscspace_degree 0 \ 42407a75e27SStefano Zampini -potential_petscdualspace_lagrange_use_moments \ 42507a75e27SStefano Zampini -potential_petscdualspace_lagrange_moment_order 2 \ 42607a75e27SStefano Zampini -field_petscfe_default_quadrature_order 1 \ 42707a75e27SStefano Zampini -field_petscspace_degree 1 \ 42807a75e27SStefano Zampini -field_petscspace_type sum \ 42907a75e27SStefano Zampini -field_petscspace_variables 2 \ 43007a75e27SStefano Zampini -field_petscspace_components 2 \ 43107a75e27SStefano Zampini -field_petscspace_sum_spaces 2 \ 43207a75e27SStefano Zampini -field_petscspace_sum_concatenate true \ 43307a75e27SStefano Zampini -field_sumcomp_0_petscspace_variables 2 \ 43407a75e27SStefano Zampini -field_sumcomp_0_petscspace_type tensor \ 43507a75e27SStefano Zampini -field_sumcomp_0_petscspace_tensor_spaces 2 \ 43607a75e27SStefano Zampini -field_sumcomp_0_petscspace_tensor_uniform false \ 43707a75e27SStefano Zampini -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \ 43807a75e27SStefano Zampini -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \ 43907a75e27SStefano Zampini -field_sumcomp_1_petscspace_variables 2 \ 44007a75e27SStefano Zampini -field_sumcomp_1_petscspace_type tensor \ 44107a75e27SStefano Zampini -field_sumcomp_1_petscspace_tensor_spaces 2 \ 44207a75e27SStefano Zampini -field_sumcomp_1_petscspace_tensor_uniform false \ 44307a75e27SStefano Zampini -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \ 44407a75e27SStefano Zampini -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \ 44507a75e27SStefano Zampini -field_petscdualspace_form_degree -1 \ 44607a75e27SStefano Zampini -field_petscdualspace_order 1 \ 44707a75e27SStefano Zampini -field_petscdualspace_lagrange_trimmed true \ 44807a75e27SStefano Zampini -dm_plex_box_faces 3,3 dm_refine 1 -dm_plex_simplex 0 \ 44907a75e27SStefano Zampini -dmsnes_check .001 -snes_error_if_not_converged \ 45007a75e27SStefano Zampini -ksp_error_if_not_converged -ksp_type cg \ 45107a75e27SStefano Zampini -petscpartitioner_type simple -dm_mat_type is \ 45207a75e27SStefano Zampini -pc_type bddc -pc_bddc_use_local_mat_graph 0 \ 45307a75e27SStefano Zampini -pc_bddc_benign_trick -pc_bddc_nonetflux -pc_bddc_detect_disconnected -pc_bddc_use_qr_single \ 45407a75e27SStefano Zampini -pc_bddc_coarse_redundant_pc_type svd -pc_bddc_neumann_pc_type svd -pc_bddc_dirichlet_pc_type svd 45507a75e27SStefano Zampini 45607a75e27SStefano Zampini test: 45769b95281SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.0, 1.0] 45869b95281SMatthew G. Knepley # Using -sol_type quadratic -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.9, 1.0] 45969b95281SMatthew G. Knepley suffix: 2d_bdm1_p0_conv 460c4762a1bSJed Brown requires: triangle 461c4762a1bSJed Brown args: -sol_type quartic \ 462c4762a1bSJed Brown -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 463c4762a1bSJed Brown -snes_error_if_not_converged \ 464c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 465c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 466c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 467c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 468c4762a1bSJed Brown 469c4762a1bSJed Brown test: 47069b95281SMatthew G. Knepley # HDF5 output: -dm_view hdf5:${PETSC_DIR}/sol.h5 -potential_view hdf5:${PETSC_DIR}/sol.h5::append -exact_vec_view hdf5:${PETSC_DIR}/sol.h5::append 47169b95281SMatthew G. Knepley # VTK output: -potential_view vtk: -exact_vec_view vtk: 47269b95281SMatthew G. Knepley # VTU output: -potential_view vtk:multifield.vtu -exact_vec_view vtk:exact.vtu 47369b95281SMatthew G. Knepley suffix: 2d_q2_p0 474c4762a1bSJed Brown requires: triangle 47530602db0SMatthew G. Knepley args: -sol_type linear -dm_plex_simplex 0 \ 47669b95281SMatthew G. Knepley -field_petscspace_degree 2 -dm_refine 0 \ 477c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 478c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 479c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 480c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 481c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 482c4762a1bSJed Brown test: 48369b95281SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [0.0057, 1.0] 48469b95281SMatthew G. Knepley suffix: 2d_q2_p0_conv 485c4762a1bSJed Brown requires: triangle 48630602db0SMatthew G. Knepley args: -sol_type linear -dm_plex_simplex 0 \ 487c4762a1bSJed Brown -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 48869b95281SMatthew G. Knepley -snes_error_if_not_converged \ 489c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 490c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 491c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 492c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 49369b95281SMatthew G. Knepley test: 49469b95281SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [-0.014, 0.0066] 49569b95281SMatthew G. Knepley suffix: 2d_q2_p0_neumann_conv 49669b95281SMatthew G. Knepley requires: triangle 49730602db0SMatthew G. Knepley args: -sol_type quartic_neumann -dm_plex_simplex 0 \ 49869b95281SMatthew G. Knepley -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 49969b95281SMatthew G. Knepley -snes_error_if_not_converged \ 50069b95281SMatthew G. Knepley -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 50169b95281SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 50269b95281SMatthew G. Knepley -fieldsplit_field_pc_type lu \ 50369b95281SMatthew G. Knepley -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type svd 50469b95281SMatthew G. Knepley 505c4762a1bSJed Brown TEST*/ 506c4762a1bSJed Brown 507c4762a1bSJed Brown /* These tests will be active once tensor P^- is working 508c4762a1bSJed Brown 509c4762a1bSJed Brown test: 510c4762a1bSJed Brown suffix: 2d_bdmq1_p0_0 511c4762a1bSJed Brown requires: triangle 51230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -sol_type linear \ 513c4762a1bSJed Brown -field_petscspace_poly_type pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \ 514c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 515c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 516c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 517c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 518c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 519c4762a1bSJed Brown test: 520c4762a1bSJed Brown suffix: 2d_bdmq1_p0_2 521c4762a1bSJed Brown requires: triangle 52230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -sol_type quartic \ 523c4762a1bSJed Brown -field_petscspace_poly_type_no pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \ 524c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 525c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 526c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 527c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 528c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 529c4762a1bSJed Brown 530c4762a1bSJed Brown */ 531