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 3269b95281SMatthew G. Knepley typedef enum {SOL_LINEAR, SOL_QUADRATIC, SOL_QUARTIC, SOL_QUARTIC_NEUMANN, SOL_UNKNOWN, NUM_SOLTYPE} SolType; 3369b95281SMatthew G. Knepley const char *SolTypeNames[NUM_SOLTYPE+3] = {"linear", "quadratic", "quartic", "quartic_neumann", "unknown", "SolType", "SOL_", NULL}; 34c4762a1bSJed Brown 35c4762a1bSJed Brown typedef struct { 36c4762a1bSJed Brown SolType solType; /* The type of exact solution */ 37c4762a1bSJed Brown } AppCtx; 38c4762a1bSJed Brown 3969b95281SMatthew G. Knepley static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 4069b95281SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 4169b95281SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 4269b95281SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 4369b95281SMatthew G. Knepley { 4469b95281SMatthew G. Knepley PetscInt d; 4569b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0]+d*dim+d]; 4669b95281SMatthew G. Knepley } 4769b95281SMatthew G. Knepley 48c4762a1bSJed Brown /* 2D Dirichlet potential example 49c4762a1bSJed Brown 50c4762a1bSJed Brown u = x 51c4762a1bSJed Brown q = <1, 0> 52c4762a1bSJed Brown f = 0 53c4762a1bSJed Brown 54c4762a1bSJed Brown We will need a boundary integral of u over \Gamma. 55c4762a1bSJed Brown */ 56c4762a1bSJed Brown static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 57c4762a1bSJed Brown { 58c4762a1bSJed Brown u[0] = x[0]; 59c4762a1bSJed Brown return 0; 60c4762a1bSJed Brown } 61c4762a1bSJed Brown static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 62c4762a1bSJed Brown { 63c4762a1bSJed Brown PetscInt c; 64c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = c ? 0.0 : 1.0; 65c4762a1bSJed Brown return 0; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 6869b95281SMatthew G. Knepley static void f0_linear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 6969b95281SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 7069b95281SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 7169b95281SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 7269b95281SMatthew G. Knepley { 7369b95281SMatthew G. Knepley f0[0] = 0.0; 7469b95281SMatthew 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); 7569b95281SMatthew G. Knepley } 7669b95281SMatthew G. Knepley static void f0_bd_linear_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 7769b95281SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 7869b95281SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 7969b95281SMatthew G. Knepley PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 8069b95281SMatthew G. Knepley { 8169b95281SMatthew G. Knepley PetscScalar potential; 8269b95281SMatthew G. Knepley PetscInt d; 8369b95281SMatthew G. Knepley 8469b95281SMatthew G. Knepley linear_u(dim, t, x, dim, &potential, NULL); 8569b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[d] = -potential*n[d]; 8669b95281SMatthew G. Knepley } 8769b95281SMatthew G. Knepley 88c4762a1bSJed Brown /* 2D Dirichlet potential example 89c4762a1bSJed Brown 90c4762a1bSJed Brown u = x^2 + y^2 91c4762a1bSJed Brown q = <2x, 2y> 92c4762a1bSJed Brown f = 4 93c4762a1bSJed Brown 94c4762a1bSJed Brown We will need a boundary integral of u over \Gamma. 95c4762a1bSJed Brown */ 96c4762a1bSJed Brown static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 97c4762a1bSJed Brown { 98c4762a1bSJed Brown PetscInt d; 99c4762a1bSJed Brown 100c4762a1bSJed Brown u[0] = 0.0; 101c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] += x[d]*x[d]; 102c4762a1bSJed Brown return 0; 103c4762a1bSJed Brown } 104c4762a1bSJed Brown static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown PetscInt c; 107c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = 2.0*x[c]; 108c4762a1bSJed Brown return 0; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 11169b95281SMatthew G. Knepley static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 11269b95281SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 11369b95281SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 11469b95281SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 11569b95281SMatthew G. Knepley { 11669b95281SMatthew G. Knepley f0[0] = -4.0; 11769b95281SMatthew 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); 11869b95281SMatthew G. Knepley } 11969b95281SMatthew G. Knepley static void f0_bd_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 12069b95281SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 12169b95281SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 12269b95281SMatthew G. Knepley PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 12369b95281SMatthew G. Knepley { 12469b95281SMatthew G. Knepley PetscScalar potential; 12569b95281SMatthew G. Knepley PetscInt d; 12669b95281SMatthew G. Knepley 12769b95281SMatthew G. Knepley quadratic_u(dim, t, x, dim, &potential, NULL); 12869b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[d] = -potential*n[d]; 12969b95281SMatthew G. Knepley } 13069b95281SMatthew G. Knepley 131c4762a1bSJed Brown /* 2D Dirichlet potential example 132c4762a1bSJed Brown 133c4762a1bSJed Brown u = x (1-x) y (1-y) 134c4762a1bSJed Brown q = <(1-2x) y (1-y), x (1-x) (1-2y)> 135c4762a1bSJed Brown f = -y (1-y) - x (1-x) 136c4762a1bSJed Brown 137c4762a1bSJed Brown u|_\Gamma = 0 so that the boundary integral vanishes 138c4762a1bSJed Brown */ 139c4762a1bSJed Brown static PetscErrorCode quartic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 140c4762a1bSJed Brown { 141c4762a1bSJed Brown PetscInt d; 142c4762a1bSJed Brown 143c4762a1bSJed Brown u[0] = 1.0; 144c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] *= x[d]*(1.0 - x[d]); 145c4762a1bSJed Brown return 0; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown static PetscErrorCode quartic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 148c4762a1bSJed Brown { 149c4762a1bSJed Brown PetscInt c, d; 150c4762a1bSJed Brown 151c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 152c4762a1bSJed Brown u[c] = 1.0; 153c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 154c4762a1bSJed Brown if (c == d) u[c] *= 1 - 2.0*x[d]; 155c4762a1bSJed Brown else u[c] *= x[d]*(1.0 - x[d]); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown } 158c4762a1bSJed Brown return 0; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown static void f0_quartic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 162c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 163c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 164c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 165c4762a1bSJed Brown { 166c4762a1bSJed Brown PetscInt d; 167c4762a1bSJed Brown f0[0] = 0.0; 16869b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += 2.0*x[d]*(1.0 - x[d]); 16969b95281SMatthew 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); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 17269b95281SMatthew G. Knepley /* 2D Dirichlet potential example 17369b95281SMatthew G. Knepley 17469b95281SMatthew G. Knepley u = x (1-x) y (1-y) 17569b95281SMatthew G. Knepley q = <(1-2x) y (1-y), x (1-x) (1-2y)> 17669b95281SMatthew G. Knepley f = -y (1-y) - x (1-x) 17769b95281SMatthew G. Knepley 17869b95281SMatthew G. Knepley du/dn_\Gamma = {(1-2x) y (1-y), x (1-x) (1-2y)} . n produces an essential condition on q 17969b95281SMatthew G. Knepley */ 18069b95281SMatthew G. Knepley 181c4762a1bSJed Brown static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 182c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 183c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 184c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 185c4762a1bSJed Brown { 186c4762a1bSJed Brown PetscInt c; 18769b95281SMatthew G. Knepley for (c = 0; c < dim; ++c) f0[c] = u[uOff[0]+c]; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* <\nabla\cdot w, u> = <\nabla w, Iu> */ 191c4762a1bSJed Brown static void f1_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 192c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 193c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 194c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 195c4762a1bSJed Brown { 19669b95281SMatthew G. Knepley PetscInt c; 19769b95281SMatthew G. Knepley for (c = 0; c < dim; ++c) f1[c*dim+c] = u[uOff[1]]; 198c4762a1bSJed Brown } 199c4762a1bSJed Brown 20069b95281SMatthew G. Knepley /* <t, q> */ 201c4762a1bSJed Brown static void g0_qq(PetscInt dim, PetscInt Nf, PetscInt NfAux, 202c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 203c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 204c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 205c4762a1bSJed Brown { 206c4762a1bSJed Brown PetscInt c; 207c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c*dim+c] = 1.0; 208c4762a1bSJed Brown } 209c4762a1bSJed Brown 21069b95281SMatthew G. Knepley /* <\nabla\cdot t, u> = <\nabla t, Iu> */ 211c4762a1bSJed Brown static void g2_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 212c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 213c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 214c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 215c4762a1bSJed Brown { 216c4762a1bSJed Brown PetscInt d; 217c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d*dim+d] = 1.0; 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 22069b95281SMatthew G. Knepley /* <v, \nabla\cdot q> */ 221c4762a1bSJed Brown static void g1_uq(PetscInt dim, PetscInt Nf, PetscInt NfAux, 222c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 223c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 224c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 225c4762a1bSJed Brown { 226c4762a1bSJed Brown PetscInt d; 22769b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 231c4762a1bSJed Brown { 232c4762a1bSJed Brown PetscFunctionBeginUser; 233c4762a1bSJed Brown options->solType = SOL_LINEAR; 234d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 2359566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sol_type", "Type of exact solution", "ex24.c", SolTypeNames, (PetscEnum) options->solType, (PetscEnum *) &options->solType, NULL)); 236d0609cedSBarry Smith PetscOptionsEnd(); 237c4762a1bSJed Brown PetscFunctionReturn(0); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 241c4762a1bSJed Brown { 242c4762a1bSJed Brown PetscFunctionBeginUser; 2439566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2449566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Example Mesh")); 2469566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 2479566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2489566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 249c4762a1bSJed Brown PetscFunctionReturn(0); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 253c4762a1bSJed Brown { 25445480ffeSMatthew G. Knepley PetscDS ds; 25545480ffeSMatthew G. Knepley DMLabel label; 25645480ffeSMatthew G. Knepley PetscWeakForm wf; 257c4762a1bSJed Brown const PetscInt id = 1; 25845480ffeSMatthew G. Knepley PetscInt bd; 259c4762a1bSJed Brown 260c4762a1bSJed Brown PetscFunctionBeginUser; 2619566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2629566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2639566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_q, f1_q)); 2649566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qu, NULL)); 2669566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_uq, NULL, NULL)); 267c4762a1bSJed Brown switch (user->solType) 268c4762a1bSJed Brown { 269c4762a1bSJed Brown case SOL_LINEAR: 2709566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_linear_u, NULL)); 2719566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 2729566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2739566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_linear_q, 0, NULL)); 2749566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, linear_q, user)); 2759566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, linear_u, user)); 276c4762a1bSJed Brown break; 277c4762a1bSJed Brown case SOL_QUADRATIC: 2789566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_quadratic_u, NULL)); 2799566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 2809566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2819566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_quadratic_q, 0, NULL)); 2829566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_q, user)); 2839566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quadratic_u, user)); 284c4762a1bSJed Brown break; 285c4762a1bSJed Brown case SOL_QUARTIC: 2869566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL)); 2879566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user)); 2889566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user)); 289c4762a1bSJed Brown break; 29069b95281SMatthew G. Knepley case SOL_QUARTIC_NEUMANN: 2919566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL)); 2929566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user)); 2939566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user)); 2949566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "Flux condition", label, 1, &id, 0, 0, NULL, (void (*)(void)) quartic_q, NULL, user, NULL)); 29569b95281SMatthew G. Knepley break; 29698921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown PetscFunctionReturn(0); 299c4762a1bSJed Brown } 300c4762a1bSJed Brown 301c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 302c4762a1bSJed Brown { 303c4762a1bSJed Brown DM cdm = dm; 304c4762a1bSJed Brown PetscFE feq, feu; 30569b95281SMatthew G. Knepley DMPolytopeType ct; 30669b95281SMatthew G. Knepley PetscBool simplex; 30769b95281SMatthew G. Knepley PetscInt dim, cStart; 308c4762a1bSJed Brown 309c4762a1bSJed Brown PetscFunctionBeginUser; 3109566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3119566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 3129566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 31369b95281SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 314c4762a1bSJed Brown /* Create finite element */ 3159566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field_", -1, &feq)); 3169566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) feq, "field")); 3179566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "potential_", -1, &feu)); 3189566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) feu, "potential")); 3199566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(feq, feu)); 320c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 3219566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) feq)); 3229566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject) feu)); 3239566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 3249566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 325c4762a1bSJed Brown while (cdm) { 3269566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm,cdm)); 3279566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 328c4762a1bSJed Brown } 3299566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feq)); 3309566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feu)); 331c4762a1bSJed Brown PetscFunctionReturn(0); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown 334c4762a1bSJed Brown int main(int argc, char **argv) 335c4762a1bSJed Brown { 336c4762a1bSJed Brown DM dm; /* Problem specification */ 337c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 338c4762a1bSJed Brown Vec u; /* Solutions */ 339c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 340c4762a1bSJed Brown 341*327415f7SBarry Smith PetscFunctionBeginUser; 3429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3439566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 344c4762a1bSJed Brown /* Primal system */ 3459566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 3469566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 3479566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 3489566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, SetupPrimalProblem, &user)); 3499566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3509566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 3519566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 3529566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 3539566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 3549566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 3559566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 3569566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 3579566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 358c4762a1bSJed Brown /* Cleanup */ 3599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3609566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 3619566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 363b122ec5aSJacob Faibussowitsch return 0; 364c4762a1bSJed Brown } 365c4762a1bSJed Brown 366c4762a1bSJed Brown /*TEST 367c4762a1bSJed Brown 368c4762a1bSJed Brown test: 36969b95281SMatthew G. Knepley suffix: 2d_bdm1_p0 370c4762a1bSJed Brown requires: triangle 371c4762a1bSJed Brown args: -sol_type linear \ 37269b95281SMatthew G. Knepley -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 1 \ 373c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 374c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 375c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 376c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 377c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 378c4762a1bSJed Brown test: 37969b95281SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.0, 1.0] 38069b95281SMatthew G. Knepley # Using -sol_type quadratic -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.9, 1.0] 38169b95281SMatthew G. Knepley suffix: 2d_bdm1_p0_conv 382c4762a1bSJed Brown requires: triangle 383c4762a1bSJed Brown args: -sol_type quartic \ 384c4762a1bSJed Brown -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 385c4762a1bSJed Brown -snes_error_if_not_converged \ 386c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 387c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 388c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 389c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 390c4762a1bSJed Brown 391c4762a1bSJed Brown test: 39269b95281SMatthew 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 39369b95281SMatthew G. Knepley # VTK output: -potential_view vtk: -exact_vec_view vtk: 39469b95281SMatthew G. Knepley # VTU output: -potential_view vtk:multifield.vtu -exact_vec_view vtk:exact.vtu 39569b95281SMatthew G. Knepley suffix: 2d_q2_p0 396c4762a1bSJed Brown requires: triangle 39730602db0SMatthew G. Knepley args: -sol_type linear -dm_plex_simplex 0 \ 39869b95281SMatthew G. Knepley -field_petscspace_degree 2 -dm_refine 0 \ 399c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 400c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 401c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 402c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 403c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 404c4762a1bSJed Brown test: 40569b95281SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [0.0057, 1.0] 40669b95281SMatthew G. Knepley suffix: 2d_q2_p0_conv 407c4762a1bSJed Brown requires: triangle 40830602db0SMatthew G. Knepley args: -sol_type linear -dm_plex_simplex 0 \ 409c4762a1bSJed Brown -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 41069b95281SMatthew G. Knepley -snes_error_if_not_converged \ 411c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 412c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 413c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 414c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 41569b95281SMatthew G. Knepley test: 41669b95281SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [-0.014, 0.0066] 41769b95281SMatthew G. Knepley suffix: 2d_q2_p0_neumann_conv 41869b95281SMatthew G. Knepley requires: triangle 41930602db0SMatthew G. Knepley args: -sol_type quartic_neumann -dm_plex_simplex 0 \ 42069b95281SMatthew G. Knepley -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 42169b95281SMatthew G. Knepley -snes_error_if_not_converged \ 42269b95281SMatthew G. Knepley -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 42369b95281SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 42469b95281SMatthew G. Knepley -fieldsplit_field_pc_type lu \ 42569b95281SMatthew G. Knepley -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type svd 42669b95281SMatthew G. Knepley 427c4762a1bSJed Brown TEST*/ 428c4762a1bSJed Brown 429c4762a1bSJed Brown /* These tests will be active once tensor P^- is working 430c4762a1bSJed Brown 431c4762a1bSJed Brown test: 432c4762a1bSJed Brown suffix: 2d_bdmq1_p0_0 433c4762a1bSJed Brown requires: triangle 43430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -sol_type linear \ 435c4762a1bSJed 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 \ 436c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 437c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 438c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 439c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 440c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 441c4762a1bSJed Brown test: 442c4762a1bSJed Brown suffix: 2d_bdmq1_p0_2 443c4762a1bSJed Brown requires: triangle 44430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -sol_type quartic \ 445c4762a1bSJed 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 \ 446c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 447c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 448c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 449c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 450c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 451c4762a1bSJed Brown 452c4762a1bSJed Brown */ 453